90int main(
int argc,
char *argv[])
93 cout <<
"This example is not supported in single precision.\n\n";
94 return MFEM_SKIP_RETURN_VALUE;
98 const char *mesh_file =
"../data/star.mesh";
102 bool visualization =
true;
103 bool verification =
false;
106 args.
AddOption(&mesh_file,
"-m",
"--mesh",
107 "Mesh file to use.");
109 "Finite element order (polynomial degree) or -1 for"
110 " isoparametric space.");
111 args.
AddOption(&num_refs,
"-r",
"--refs",
112 "Number of uniform refinements");
114 "Fractional exponent");
115 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
116 "--no-visualization",
117 "Enable or disable GLVis visualization.");
118 args.
AddOption(&verification,
"-ver",
"--verification",
"-no-ver",
120 "Use sinusoidal function (f) for manufactured "
131 int progress_steps = 1;
135 const int power_of_laplace = (int)floor(
alpha);
136 real_t exponent_to_approximate =
alpha - power_of_laplace;
137 bool integer_order =
false;
139 if (abs(exponent_to_approximate) > 1e-12)
141 mfem::out <<
"Approximating the fractional exponent "
142 << exponent_to_approximate
149 alpha = exponent_to_approximate + power_of_laplace;
153 integer_order =
true;
154 mfem::out <<
"Treating integer order PDE." << endl;
158 Mesh mesh(mesh_file, 1, 1);
162 for (
int i = 0; i < num_refs; i++)
170 cout <<
"Number of degrees of freedom: "
186 for (
int i=0; i<x.Size(); i++)
188 val *= sin(M_PI*x(i));
190 return pow(x.Size()*pow(M_PI,2),
alpha) * val;
226 if (power_of_laplace > 0)
247 mfem::out <<
"\nComputing (-Δ) ^ -" << power_of_laplace
248 <<
" ( f ) " << endl;
249 for (
int i = 0; i < power_of_laplace; i++)
252 PCG(*Op, M, B, X, 3, 300, 1e-12, 0.0);
255 if (i == power_of_laplace - 1)
259 if (integer_order && verification)
270 oss_f.str(
""); oss_f.clear();
271 oss_f <<
"Step " << progress_steps++ <<
": Solution of PDE -Δ ^ "
274 fout <<
"solution\n" << mesh << g
275 <<
"window_title '" << oss_f.str() <<
"'" << flush;
305 ostringstream oss_x, oss_u;
315 for (
int i = 0; i < coeffs.
Size(); i++)
317 mfem::out <<
"\nSolving PDE -Δ u + " << -poles[i]
318 <<
" u = " << coeffs[i] <<
" g " << endl;
334 a.FormLinearSystem(ess_tdof_list, x,
b, A, X, B);
339 PCG(*A, M, B, X, 3, 300, 1e-12, 0.0);
342 a.RecoverFEMSolution(X,
b, x);
351 oss_x.str(
""); oss_x.clear();
352 oss_x <<
"Step " << progress_steps
353 <<
": Solution of PDE -Δ u + " << -poles[i]
354 <<
" u = " << coeffs[i] <<
" g";
355 xout <<
"solution\n" << mesh << x
356 <<
"window_title '" << oss_x.str() <<
"'" << flush;
358 oss_u.str(
""); oss_u.clear();
359 oss_u <<
"Step " << progress_steps + 1
360 <<
": Solution of fractional PDE (-Δ)^" <<
alpha
362 uout <<
"solution\n" << mesh <<
u
363 <<
"window_title '" << oss_u.str() <<
"'"
374 auto solution = [] (
const Vector &x)
377 for (
int i=0; i<x.
Size(); i++)
379 val *= sin(M_PI*x(i));
384 real_t l2_error =
u.ComputeL2Error(sol);
386 string manufactured_solution,expected_mesh;
390 manufactured_solution =
"sin(π x)";
391 expected_mesh =
"inline_segment.mesh";
394 manufactured_solution =
"sin(π x) sin(π y)";
395 expected_mesh =
"inline_quad.mesh";
398 manufactured_solution =
"sin(π x) sin(π y) sin(π z)";
399 expected_mesh =
"inline_hex.mesh";
404 <<
"\n\nSolution Verification in "<<
dim <<
"D \n\n"
405 <<
"Manufactured solution : " << manufactured_solution <<
"\n"
406 <<
"Expected mesh : " << expected_mesh <<
"\n"
407 <<
"Your mesh : " << mesh_file <<
"\n"
408 <<
"L2 error : " << l2_error <<
"\n\n"
409 << string(80,
'=') << endl;
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int Size() const
Return the logical size of the array.
A coefficient that is constant across space and time.
Class for domain integration .
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
virtual const SparseMatrix * GetRestrictionMatrix() const
The returned SparseMatrix is owned by the FiniteElementSpace.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
A general function coefficient.
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Arbitrary order H1-conforming (continuous) finite elements.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
int Dimension() const
Dimension of the reference space used within the elements.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Pointer to an Operator of a specified type.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
bool Good() const
Return true if the command line options were parsed successfully.
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
int Size() const
Returns the size of the vector.
void SetSubVectorComplement(const Array< int > &dofs, const real_t val)
Set all vector entries NOT in the dofs Array to the given val.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
void ComputePartialFractionApproximation(real_t &alpha, Array< real_t > &coeffs, Array< real_t > &poles, real_t lmax=1000., real_t tol=1e-10, int npoints=1000, int max_order=100)
real_t u(const Vector &xvec)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
std::function< real_t(const Vector &)> f(real_t mass_coeff)