92int main(
int argc,
char *argv[])
95 const char *mesh_file =
"../../data/inline-quad.mesh";
99 bool visualization =
true;
101 bool static_cond =
false;
104 args.
AddOption(&mesh_file,
"-m",
"--mesh",
105 "Mesh file to use.");
107 "Finite element order (polynomial degree).");
108 args.
AddOption(&delta_order,
"-do",
"--delta_order",
109 "Order enrichment for DPG test space.");
110 args.
AddOption(&ref,
"-ref",
"--num_refinements",
111 "Number of uniform refinements");
112 args.
AddOption(&iprob,
"-prob",
"--problem",
"Problem case"
113 " 0: manufactured, 1: general");
114 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
115 "--no-static-condensation",
"Enable static condensation.");
116 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
117 "--no-visualization",
118 "Enable or disable GLVis visualization.");
127 if (iprob > 1) { iprob = 1; }
130 Mesh mesh(mesh_file, 1, 1);
132 MFEM_VERIFY(
dim > 1,
"Dimension = 1 is not supported in this example");
164 int test_order = order+delta_order;
172 trial_fes.
Append(sigma_fes);
173 trial_fes.
Append(hatu_fes);
174 trial_fes.
Append(hatsigma_fes);
193 TrialSpace::u_space,TestSpace::tau_space);
197 negone)), TrialSpace::sigma_space, TestSpace::tau_space);
201 TrialSpace::sigma_space,TestSpace::v_space);
205 TrialSpace::hatu_space,TestSpace::tau_space);
209 TrialSpace::hatsigma_space, TestSpace::v_space);
214 TestSpace::tau_space, TestSpace::tau_space);
217 TestSpace::tau_space, TestSpace::tau_space);
220 TestSpace::v_space, TestSpace::v_space);
223 TestSpace::v_space, TestSpace::v_space);
226 if (
prob == prob_type::manufactured)
242 if (
prob == prob_type::manufactured)
244 std::cout <<
"\n Ref |"
248 <<
" PCG it |" << endl;
249 std::cout << std::string(50,
'-')
255 if (static_cond) {
a->EnableStaticCondensation(); }
256 for (
int it = 0; it<=ref; it++)
270 for (
int i = 0; i < ess_tdof_list.
Size(); i++)
280 offsets[4] = hatsigma_fes->
GetVSize();
285 if (
prob == prob_type::manufactured)
293 a->FormLinearSystem(ess_tdof_list,x,Ah,X,B);
312 a->RecoverFEMSolution(X,x);
318 if (
prob == prob_type::manufactured)
323 real_t L2Error = sqrt(u_err*u_err + sigma_err*sigma_err);
324 real_t rate_err = (it) ?
dim*log(err0/L2Error)/log((
real_t)dof0/l2dofs) : 0.0;
328 std::ios oldState(
nullptr);
329 oldState.copyfmt(std::cout);
330 std::cout << std::right << std::setw(5) << it <<
" | "
331 << std::setw(10) << dof0 <<
" | "
332 << std::setprecision(3)
333 << std::setw(10) << std::scientific << err0 <<
" | "
334 << std::setprecision(2)
335 << std::setw(6) << std::fixed << rate_err <<
" | "
338 std::cout.copyfmt(oldState);
343 const char * keys = (it == 0 &&
dim == 2) ?
"jRcm\n" :
nullptr;
347 "Numerical u", 0,0, 500, 500, keys);
349 "Numerical flux", 500,0,500, 500, keys);
352 if (it == ref) {
break; }
355 for (
int i =0; i<trial_fes.
Size(); i++)
357 trial_fes[i]->Update(
false);
388 for (
int i = 0; i<du.
Size(); i++)
390 du[i] = M_PI * cos(
alpha);
398 return - M_PI*M_PI *
u * X.
Size();
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
int Append(const T &el)
Append element 'el' to array, resize if necessary.
A class to handle Block diagonal preconditioners in a matrix-free implementation.
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
Array< int > & RowOffsets()
Return the row offsets for block starts.
int NumRowBlocks() const
Return the number of row blocks.
A class to handle Vectors in a block fashion.
Vector & GetBlock(int i)
Get the i-th vector in the block.
Conjugate gradient method.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
A coefficient that is constant across space and time.
for Raviart-Thomas elements
Class for domain integration .
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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 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 ...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
A general function coefficient.
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
virtual void MakeRef(FiniteElementSpace *f, real_t *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
void ProjectBdrCoefficient(Coefficient &coeff, const Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Arbitrary order H1-conforming (continuous) finite elements.
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Arbitrary order "L2-conforming" discontinuous 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.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
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.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
A general vector function coefficient.
int Size() const
Returns the size of the vector.
real_t Sum() const
Return the sum of the vector entries.
void SetSize(int s)
Resize the vector to size s.
real_t sigma(const Vector &x)
void exact_gradu(const Vector &X, Vector &gradu)
real_t exact_hatu(const Vector &X)
real_t exact_laplacian_u(const Vector &X)
real_t exact_u(const Vector &X)
void exact_sigma(const Vector &X, Vector &sigma)
void exact_hatsigma(const Vector &X, Vector &hatsigma)
real_t f_exact(const Vector &X)
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
real_t u(const Vector &xvec)
std::function< real_t(const Vector &)> f(real_t mass_coeff)