15 #include "../../examples/ex33.hpp" 28 return (fespace_ptr->
GetMyRank() == 0 && print_level > 0);
33 os <<
"\n<Boundary Info>\n" 34 <<
" Boundary Conditions:\n";
37 os <<
" Boundary " << it.first <<
": ";
58 bool first_print_statement =
true;
62 os <<
" Inhomogeneous Dirichlet defined on ";
65 if (!first_print_statement)
71 first_print_statement =
false;
73 os << it.first <<
"(=" << it.second <<
")";
77 os <<
"<Boundary Info>\n\n";
89 if (boundary.
Find(it.first) == -1)
91 MFEM_ABORT(
" Boundary " 93 <<
" is not defined on the mesh but in Boundary class." 101 std::vector<int> boundary_attributes_keys;
102 for (
int i = 0; i < boundary.
Size(); i++)
106 boundary_attributes_keys.push_back(boundary[i]);
109 if (!boundary_attributes_keys.empty())
112 for (
const auto &it : boundary_attributes_keys)
116 mfem::out <<
") are defined on the mesh but not in the";
117 mfem::out <<
" boundary attributes (Use Neumann).";
125 MFEM_ABORT(
" Periodic boundaries must be defined on the mesh" 126 <<
", not in Boundaries. Exiting...")
140 mfem::out <<
"<Boundary::ComputeBoundaryError>" 163 mfem::out <<
"->Boundary " << i + 1 <<
"\n";
166 mfem::out <<
" Gamma : " << gamma <<
"\n";
167 mfem::out <<
" Average : " << avg <<
"\n";
168 mfem::out <<
" Error : " << error <<
"\n\n";
174 mfem::out <<
"<Boundary::ComputeBoundaryError>" << std::endl;
243 double alpha,
double beta,
double gamma,
double &glb_err)
246 double &nrm = loc_vals[0];
247 double &avg = loc_vals[1];
248 double &error = loc_vals[2];
254 const bool a_is_zero =
alpha == 0.0;
255 const bool b_is_zero =
beta == 0.0;
258 MFEM_ASSERT(fes.
GetVDim() == 1,
"");
265 for (
int i = 0; i < mesh.GetNBE(); i++)
267 if (bdr[mesh.GetBdrAttribute(i) - 1] == 0)
280 const int int_order = 2 * fe.
GetOrder() + 3;
313 val +=
beta * (shape * loc_dofs);
317 nrm += ip.
weight * face_weight;
320 avg += val * ip.
weight * face_weight;
324 error += (val * val) * ip.
weight * face_weight;
329 MPI_Allreduce(loc_vals, glb_vals, 3, MPI_DOUBLE, MPI_SUM, fes.
GetComm());
331 double glb_nrm = glb_vals[0];
332 double glb_avg = glb_vals[1];
333 glb_err = glb_vals[2];
336 if (std::abs(glb_nrm) > 0.0)
344 glb_err = (glb_err >= 0.0) ? sqrt(glb_err) : -sqrt(-glb_err);
352 double l3,
double e1,
double e2,
double e3)
355 fespace_ptr_(fespace),
367 mfem::out <<
"<SPDESolver> Initialize Solver .." << std::endl;
376 if (bdr_attributes.Size() > 0)
379 nbc = bdr_attributes.
Max() - bdr_attributes.Min() + 1;
394 dbc_marker_[it.first - 1] = 1;
397 rbc_marker_[it.first - 1] = 1;
413 alpha_ = (nu_ +
dim / 2.0) / 2.0;
414 integer_order_of_exponent_ =
static_cast<int>(std::floor(alpha_));
415 double exponent_to_approximate = alpha_ - integer_order_of_exponent_;
418 ComputeRationalCoefficients(exponent_to_approximate);
423 auto diffusion_tensor =
446 if (prolongation_matrix_)
452 mfem::err <<
"<SPDESolver> prolongation matrix is not defined" << std::endl;
454 if (restriction_matrix_)
460 mfem::err <<
"<SPDESolver> restriction matrix is not defined" << std::endl;
467 <<
" [s]" << std::endl;
486 if (integer_order_of_exponent_ > 0)
490 mfem::out <<
"<SPDESolver> Solving PDE (A)^" << integer_order_of_exponent_
491 <<
" u = f" << std::endl;
493 ActivateRepeatedSolve();
494 Solve(
b, helper_gf, 1.0, 1.0, integer_order_of_exponent_);
507 DeactivateRepeatedSolve();
518 for (
int i = 0; i < coeffs_.
Size(); i++)
522 mfem::out <<
"\n<SPDESolver> Solving PDE -Δ u + " << -poles_[i]
523 <<
" u = " << coeffs_[i] <<
" g " << std::endl;
526 Solve(
b, helper_gf, 1.0 - poles_[i], coeffs_[i]);
541 <<
" [s]" << std::endl;
551 b_wn->AddDomainIntegrator(integ);
558 MFEM_ABORT(
"Need to call SPDESolver::SetupRandomFieldGenerator(...) first");
564 (*b_wn) *= normalization;
571 double l2,
double l3,
589 const double gamma1 = tgamma(nu + static_cast<double>(
dim) / 2.0);
590 const double gamma2 = tgamma(nu);
591 return sqrt(pow(2 * M_PI,
dim / 2.0) * det * gamma1 /
592 (gamma2 * pow(nu,
dim / 2.0)));
596 double l3,
double e1,
597 double e2,
double e3,
603 const double c1 = cos(e1);
604 const double s1 = sin(e1);
605 const double c2 = cos(e2);
606 const double s2 = sin(e2);
607 const double c3 = cos(e3);
608 const double s3 = sin(e3);
612 R(0, 0) = c1 * c3 - c2 * s1 * s3;
613 R(0, 1) = -c1 * s3 - c2 * c3 * s1;
615 R(1, 0) = c3 * s1 + c1 * c2 * s3;
616 R(1, 1) = c1 * c2 * c3 - s1 * s3;
624 l(0) = std::pow(l1, 2);
625 l(1) = std::pow(l2, 2);
626 l(2) = std::pow(l3, 2);
627 l *= (1 / (2.0 * nu));
637 const double c1 = cos(e1);
638 const double s1 = sin(e1);
645 l(0) = std::pow(l1, 2);
646 l(1) = std::pow(l2, 2);
647 l *= (1 / (2.0 * nu));
655 res(0, 0) = std::pow(l1, 2) / (2.0 * nu);
661 double beta,
int exponent)
666 if (prolongation_matrix_)
683 restriction_matrix_->
Mult(x, X_);
687 Add(1.0, stiffness_,
alpha, mass_bc_);
688 HypreParMatrix *Ae = Op->EliminateRowsCols(ess_tdof_list_);
689 Op->EliminateBC(*Ae, ess_tdof_list_, X_, B_);
691 for (
int i = 0; i < exponent; i++)
694 SolveLinearSystem(Op);
701 GridFunctionCoefficient gfc(&x);
702 ParLinearForm previous_solution(fespace_ptr_);
703 previous_solution.AddDomainIntegrator(
new DomainLFIntegrator(gfc));
704 previous_solution.Assemble();
706 Op->EliminateBC(*Ae, ess_tdof_list_, X_, B_);
713 void SPDESolver::LiftSolution(ParGridFunction &x)
722 mfem::out <<
"\n<SPDESolver> Applying inhomogeneous DBC" << std::endl;
726 ParGridFunction helper_gf(fespace_ptr_);
734 marker[bc.first - 1] = 1;
735 ConstantCoefficient cc(bc.second);
736 helper_gf.ProjectBdrCoefficient(cc, marker);
740 ParLinearForm
b(fespace_ptr_);
741 ConstantCoefficient zero(0.0);
742 b.AddDomainIntegrator(
new DomainLFIntegrator(zero));
746 Solve(
b, helper_gf, 1.0, 1.0);
755 void SPDESolver::UpdateRHS(ParLinearForm &
b)
const 757 if (!repeated_solve_)
762 if (restriction_matrix_)
774 void SPDESolver::SolveLinearSystem(
const HypreParMatrix *Op)
776 HypreBoomerAMG prec(*Op);
777 prec.SetPrintLevel(-1);
778 CGSolver cg(fespace_ptr_->
GetComm());
782 cg.SetPreconditioner(prec);
784 cg.SetPrintLevel(std::max(0, print_level_ - 1));
788 void SPDESolver::ComputeRationalCoefficients(
double exponent)
790 if (abs(exponent) > 1e-12)
794 mfem::out <<
"<SPDESolver> Approximating the fractional exponent " 795 << exponent << std::endl;
801 alpha_ = exponent + integer_order_of_exponent_;
805 integer_order_ =
true;
808 mfem::out <<
"<SPDESolver> Treating integer order PDE." << std::endl;
Abstract class for all finite elements.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
A matrix coefficient that is constant in space and time.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void ComputeBoundaryError(const ParGridFunction &solution)
void SetRobinCoefficient(double coefficient)
Set the robin coefficient for the boundary.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
int GetNPoints() const
Returns the number of the points in the integration rule.
Class for an integration rule - an Array of IntegrationPoint.
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void Solve(ParLinearForm &b, ParGridFunction &x)
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
A coefficient that is constant across space and time.
int Dimension() const
Dimension of the reference space used within the elements.
void SetSize(int s)
Resize the vector to size s.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Data type dense matrix using column-major storage.
static DenseMatrix ConstructMatrixCoefficient(double l1, double l2, double l3, double e1, double e2, double e3, double nu, int dim)
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
void UpdateIntegrationCoefficients(int i, double &alpha, double &beta, double &gamma)
Abstract parallel finite element space.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
void CalcOrtho(const DenseMatrix &J, Vector &n)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
void Stop()
Stop the stopwatch.
void PrintInfo(std::ostream &os=mfem::out) const
Print the information specifying the boundary conditions.
std::map< int, double > dirichlet_coefficients
Coefficient for inhomogeneous Dirichlet boundary conditions.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
double IntegrateBC(const ParGridFunction &x, const Array< int > &bdr, double alpha, double beta, double gamma, double &glb_err)
void GenerateRandomField(ParGridFunction &x)
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
int Find(const T &el) const
Return the first index where 'el' is found; return -1 if not found.
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
ParMesh * GetParMesh() const
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
void MultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
void SetupRandomFieldGenerator(int seed=0)
void AddInhomogeneousDirichletBoundaryCondition(int boundary, double coefficient)
Add a inhomogeneous Dirichlet boundary condition to the boundary.
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
bool PrintOutput(const ParFiniteElementSpace *fespace_ptr, int print_level)
void Start()
Start the stopwatch. The elapsed time is not cleared.
void AddHomogeneousBoundaryCondition(int boundary, BoundaryType type)
Add a homogeneous boundary condition to the boundary.
int GetVDim() const
Returns vector dimension.
virtual DofTransformation * GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void Transpose()
(*this) = (*this)^t
virtual const FiniteElement * GetFE(int i) const
SPDESolver(double nu, const Boundary &bc, ParFiniteElementSpace *fespace, double l1=0.1, double l2=0.1, double l3=0.1, double e1=0.0, double e2=0.0, double e3=0.0)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
ParFiniteElementSpace * ParFESpace() const
void ComputePartialFractionApproximation(double &alpha, Array< double > &coeffs, Array< double > &poles, double lmax=1000., double tol=1e-10, int npoints=1000, int max_order=100)
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int SpaceDimension() const
Dimension of the physical space containing the mesh.
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
Class for integration point with weight.
static double ConstructNormalizationCoefficient(double nu, double l1, double l2, double l3, int dim)
std::map< int, BoundaryType > boundary_attributes
Map to assign homogeneous boundary conditions to defined boundary types.
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
int Size() const
Return the logical size of the array.
Class for parallel grid function.
void VerifyDefinedBoundaries(const Mesh &mesh) const
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Class for parallel meshes.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...