12 #ifndef SPDE_SOLVERS_HPP 13 #define SPDE_SOLVERS_HPP 16 #include <unordered_map> 71 double alpha,
double beta,
double gamma,
double &glb_err);
105 double l1 = 0.1,
double l2 = 0.1,
double l3 = 0.1,
106 double e1 = 0.0,
double e2 = 0.0,
double e3 = 0.0);
127 double l2,
double l3,
133 double e1,
double e2,
double e3,
144 double beta,
int exponent = 1);
156 void ActivateRepeatedSolve() { repeated_solve_ =
true; }
159 void DeactivateRepeatedSolve() { repeated_solve_ =
false; }
163 void UpdateRHS(ParLinearForm &
b)
const;
166 void ComputeRationalCoefficients(
double exponent);
171 HypreParMatrix stiffness_;
172 HypreParMatrix mass_bc_;
176 const SparseMatrix *restriction_matrix_ =
nullptr;
177 const Operator *prolongation_matrix_ =
nullptr;
184 Array<int> ess_tdof_list_;
185 ParFiniteElementSpace *fespace_ptr_;
189 Array<int> dbc_marker_;
190 Array<int> rbc_marker_;
193 Array<double> coeffs_;
194 Array<double> poles_;
199 int integer_order_of_exponent_ = 0;
210 int print_level_ = 1;
213 bool repeated_solve_ =
false;
214 bool integer_order_ =
false;
215 bool apply_lift_ =
false;
217 WhiteGaussianNoiseDomainLFIntegrator *integ =
nullptr;
218 ParLinearForm * b_wn =
nullptr;
224 #endif // SPDE_SOLVERS_HPP void ComputeBoundaryError(const ParGridFunction &solution)
void SetRobinCoefficient(double coefficient)
Set the robin coefficient for the boundary.
void SetPrintLevel(int print_level)
Set the print level.
void Solve(ParLinearForm &b, ParGridFunction &x)
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)
void UpdateIntegrationCoefficients(int i, double &alpha, double &beta, double &gamma)
Abstract parallel finite element space.
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.
double IntegrateBC(const ParGridFunction &x, const Array< int > &bdr, double alpha, double beta, double gamma, double &glb_err)
void GenerateRandomField(ParGridFunction &x)
void SetupRandomFieldGenerator(int seed=0)
void AddInhomogeneousDirichletBoundaryCondition(int boundary, double coefficient)
Add a inhomogeneous Dirichlet boundary condition to the boundary.
void AddHomogeneousBoundaryCondition(int boundary, BoundaryType type)
Add a homogeneous boundary condition to the boundary.
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...
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.
Class for parallel grid function.
void VerifyDefinedBoundaries(const Mesh &mesh) const
Wrapper for hypre's ParCSR matrix class.