43 int main(
int argc,
char *argv[])
48 const char *mesh_file =
"../data/star.mesh";
50 bool visualization = 1;
53 args.
AddOption(&mesh_file,
"-m",
"--mesh",
56 "Finite element order (polynomial degree).");
57 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
59 "Enable or disable GLVis visualization.");
71 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
80 (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
81 for (
int l = 0; l < ref_levels; l++)
100 block_offsets[1] = R_space->
GetVSize();
101 block_offsets[2] = W_space->
GetVSize();
104 std::cout <<
"***********************************************************\n";
105 std::cout <<
"dim(R) = " << block_offsets[1] - block_offsets[0] <<
"\n";
106 std::cout <<
"dim(W) = " << block_offsets[2] - block_offsets[1] <<
"\n";
107 std::cout <<
"dim(R+W) = " << block_offsets.
Last() <<
"\n";
108 std::cout <<
"***********************************************************\n";
176 for (
int i = 0; i < Md.Size(); i++)
184 #ifndef MFEM_USE_SUITESPARSE
218 <<
" iterations with a residual norm of " << solver.
GetFinalNorm() <<
".\n";
221 <<
" iterations. Residual norm is " << solver.
GetFinalNorm() <<
".\n";
222 std::cout <<
"MINRES solver took " << chrono.
RealTime() <<
"s. \n";
226 u.
MakeRef(R_space, x.GetBlock(0), 0);
227 p.
MakeRef(W_space, x.GetBlock(1), 0);
229 int order_quad = max(2, 2*order+1);
231 for (
int i=0; i < Geometry::NumGeom; ++i)
241 std::cout <<
"|| u_h - u_ex || / || u_ex || = " << err_u / norm_u <<
"\n";
242 std::cout <<
"|| p_h - p_ex || / || p_ex || = " << err_p / norm_p <<
"\n";
248 ofstream mesh_ofs(
"ex5.mesh");
249 mesh_ofs.precision(8);
250 mesh->
Print(mesh_ofs);
252 ofstream u_ofs(
"sol_u.gf");
256 ofstream p_ofs(
"sol_p.gf");
270 char vishost[] =
"localhost";
274 u_sock <<
"solution\n" << *mesh << u <<
"window_title 'Velocity'" << endl;
277 p_sock <<
"solution\n" << *mesh << p <<
"window_title 'Pressure'" << endl;
310 u(0) = - exp(xi)*sin(yi)*cos(zi);
311 u(1) = - exp(xi)*cos(yi)*cos(zi);
315 u(2) = exp(xi)*sin(yi)*sin(zi);
331 return exp(xi)*sin(yi)*cos(zi);
Class for domain integration L(v) := (f, v)
Class for integration rule.
Class for grid function - Vector with associated FE space.
int GetNumIterations() const
Data type for scaled Jacobi-type smoother of sparse matrix.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void fFun(const Vector &x, Vector &f)
void MakeRef(FiniteElementSpace *f, Vector &v, int v_offset)
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Subclass constant coefficient.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
virtual void Save()
Save the collection and a VisIt root file.
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
int main(int argc, char *argv[])
void ScaleRow(const int row, const double scale)
double GetFinalNorm() const
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void SetBlock(int i, int j, SparseMatrix *mat)
Set A(i,j) = mat.
virtual void Mult(const Vector &b, Vector &x) const
Operator application.
double ComputeLpNorm(double p, Coefficient &coeff, Mesh &mesh, const IntegrationRule *irs[])
IntegrationRules IntRules(0)
A global object with all integration rules (defined in intrules.cpp)
void SetPrintLevel(int print_lvl)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
double f_natural(const Vector &x)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
void Transpose(const Table &A, Table &At, int _ncols_A)
Transpose a Table.
double pFun_ex(const Vector &x)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void PrintUsage(std::ostream &out) const
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
virtual void Print(std::ostream &out=std::cout) const
Print the mesh to the given stream using the default MFEM mesh format.
void SetAbsTol(double atol)
void SetRelTol(double rtol)
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)
void PartialSum()
Partial Sum.
void uFun_ex(const Vector &x, Vector &u)
double gFun(const Vector &x)
void PrintOptions(std::ostream &out) const
T & Last()
Return the last element in the array.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
class for C-function coefficient
Vector & GetBlock(int i)
Get the i-th vector in the block.
Integrator for (Q u, v) for VectorFiniteElements.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Arbitrary order "L2-conforming" discontinuous finite elements.