87static const char *enum_str[] =
 
  103int main(
int argc, 
char *argv[])
 
  111   const char *mesh_file = 
"../../data/inline-quad.mesh";
 
  117   bool static_cond = 
false;
 
  119   bool visualization = 
true;
 
  121   bool paraview = 
false;
 
  124   args.
AddOption(&mesh_file, 
"-m", 
"--mesh",
 
  125                  "Mesh file to use.");
 
  127                  "Finite element order (polynomial degree).");
 
  128   args.
AddOption(&delta_order, 
"-do", 
"--delta_order",
 
  129                  "Order enrichment for DPG test space.");
 
  130   args.
AddOption(&sref, 
"-sref", 
"--num-serial-refinements",
 
  131                  "Number of initial serial uniform refinements");
 
  132   args.
AddOption(&pref, 
"-pref", 
"--num-parallel-refinements",
 
  133                  "Number of AMR refinements");
 
  134   args.
AddOption(&theta, 
"-theta", 
"--theta-factor",
 
  135                  "Refinement factor (0 indicates uniform refinements) ");
 
  136   args.
AddOption(&iprob, 
"-prob", 
"--problem", 
"Problem case" 
  137                  " 0: manufactured, 1: L-shape");
 
  138   args.
AddOption(&static_cond, 
"-sc", 
"--static-condensation", 
"-no-sc",
 
  139                  "--no-static-condensation", 
"Enable static condensation.");
 
  140   args.
AddOption(&visualization, 
"-vis", 
"--visualization", 
"-no-vis",
 
  141                  "--no-visualization",
 
  142                  "Enable or disable GLVis visualization.");
 
  143   args.
AddOption(¶view, 
"-paraview", 
"--paraview", 
"-no-paraview",
 
  145                  "Enable or disable ParaView visualization.");
 
  146   args.
AddOption(&visport, 
"-p", 
"--send-port", 
"Socket for GLVis.");
 
  157   if (iprob > 1) { iprob = 1; }
 
  162      mesh_file = 
"../../data/l-shape.mesh";
 
  170   Mesh mesh(mesh_file, 1, 1);
 
  172   MFEM_VERIFY(
dim > 1, 
"Dimension = 1 is not supported in this example");
 
  180      int size = 
nodes->Size()/2;
 
  181      for (
int i = 0; i<size; i++)
 
  184         (*nodes)[2*i] =  2*(*nodes)[2*i+1]-1;
 
  185         (*nodes)[2*i+1] = -2*x+1;
 
  190   for (
int i = 0; i<sref; i++)
 
  197   ParMesh pmesh(MPI_COMM_WORLD, mesh);
 
  232   int test_order = order+delta_order;
 
  240   trial_fes.
Append(sigma_fes);
 
  241   trial_fes.
Append(hatu_fes);
 
  242   trial_fes.
Append(hatsigma_fes);
 
  257   a->StoreMatrices(
true); 
 
  261                         TrialSpace::u_space,TestSpace::tau_space);
 
  265                                                    negone)), TrialSpace::sigma_space, TestSpace::tau_space);
 
  269                         TrialSpace::sigma_space,TestSpace::v_space);
 
  273                         TrialSpace::hatu_space,TestSpace::tau_space);
 
  277                         TrialSpace::hatsigma_space, TestSpace::v_space);
 
  282                        TestSpace::tau_space, TestSpace::tau_space);
 
  285                        TestSpace::tau_space, TestSpace::tau_space);
 
  288                        TestSpace::v_space, TestSpace::v_space);
 
  291                        TestSpace::v_space, TestSpace::v_space);
 
  294   if (
prob == prob_type::manufactured)
 
  308      std::cout << 
"\n  Ref |" 
  314                << 
" PCG it |" << endl;
 
  315      std::cout << std::string(72,
'-') << endl;
 
  343   if (static_cond) { 
a->EnableStaticCondensation(); }
 
  344   for (
int it = 0; it<=pref; it++)
 
  358      for (
int i = 0; i < ess_tdof_list.
Size(); i++)
 
  368      offsets[4] = hatsigma_fes->
GetVSize();
 
  390         M.SetDiagonalBlock(0,amg0);
 
  391         M.SetDiagonalBlock(1,amg1);
 
  397      M.SetDiagonalBlock(skip,amg2);
 
  409      M.SetDiagonalBlock(skip+1,prec);
 
  419      a->RecoverFEMSolution(X,x);
 
  421      Vector & residuals = 
a->ComputeResidual(x);
 
  426      real_t globalresidual = residual * residual;
 
  429                    MPI_MAX, MPI_COMM_WORLD);
 
  430      MPI_Allreduce(MPI_IN_PLACE, &globalresidual, 1,
 
  433      globalresidual = sqrt(globalresidual);
 
  443      real_t L2Error = sqrt(u_err*u_err + sigma_err*sigma_err);
 
  444      real_t rate_err = (it) ? 
dim*log(err0/L2Error)/log((
real_t)dof0/dofs) : 0.0;
 
  445      real_t rate_res = (it) ? 
dim*log(res0/globalresidual)/log((
 
  448      res0 = globalresidual;
 
  453         std::ios oldState(
nullptr);
 
  454         oldState.copyfmt(std::cout);
 
  455         std::cout << std::right << std::setw(5) << it << 
" | " 
  456                   << std::setw(10) <<  dof0 << 
" | " 
  457                   << std::setprecision(3)
 
  458                   << std::setw(10) << std::scientific <<  err0 << 
" | " 
  459                   << std::setprecision(2)
 
  460                   << std::setw(6) << std::fixed << rate_err << 
" | " 
  461                   << std::setprecision(3)
 
  462                   << std::setw(10) << std::scientific <<  res0 << 
" | " 
  463                   << std::setprecision(2)
 
  464                   << std::setw(6) << std::fixed << rate_res << 
" | " 
  467         std::cout.copyfmt(oldState);
 
  472         const char * keys = (it == 0 && 
dim == 2) ? 
"jRcm\n" : 
nullptr;
 
  476                        "Numerical u", 0,0,500,500,keys);
 
  478                        "Numerical flux", 500,0,500,500,keys);
 
  488      if (it == pref) { 
break; }
 
  491      for (
int iel = 0; iel<pmesh.
GetNE(); iel++)
 
  493         if (residuals[iel] >= theta * maxresidual)
 
  495            elements_to_refine.
Append(iel);
 
  501      for (
int i =0; i<trial_fes.
Size(); i++)
 
  503         trial_fes[i]->Update(
false);
 
 
  536         real_t r = sqrt(x*x + y*y);
 
  539         if (phi < 0) { phi += 2*M_PI; }
 
 
  561         real_t r = sqrt(x*x + y*y);
 
  564         if (phi < 0) { phi += 2*M_PI; }
 
  568         real_t phi_x = - y / (r*r);
 
  579         for (
int i = 0; i<du.
Size(); i++)
 
  581            du[i] = M_PI * cos(
alpha);
 
 
  592      case prob_type::manufactured:
 
  596         return - M_PI*M_PI * 
u * X.
Size();
 
  600         MFEM_ABORT(
"Should be unreachable");
 
 
  626               "f_exact should not be called for l-shape benchmark problem, i.e., f = 0")
 
 
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.
A class to handle Block systems in a matrix-free implementation.
Array< int > & RowOffsets()
Return the row offsets for block starts.
Operator & GetBlock(int i, int j)
Return a reference to block i,j.
A class to handle Vectors in a block fashion.
Vector & GetBlock(int i)
Get the i-th vector in the block.
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
A coefficient that is constant across space and time.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetTime(real_t t)
Set physical time (for time-dependent simulations)
void SetPrefixPath(const std::string &prefix)
Set the path where the DataCollection will be saved.
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...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
A general function coefficient.
Class for grid function - Vector with associated FE space.
Arbitrary order H1-conforming (continuous) finite elements.
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
The Auxiliary-space Divergence Solver in hypre.
The Auxiliary-space Maxwell Solver in hypre.
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_level)
Wrapper for hypre's ParCSR matrix class.
Abstract class for hypre's solvers and preconditioners.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
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.
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
void EnsureNodes()
Make sure that the mesh has valid nodes, i.e. its geometry is described by a vector finite element gr...
void Clear()
Clear the contents of the Mesh.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
void EnsureNCMesh(bool simplices_nonconforming=false)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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 FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
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.
Abstract parallel finite element space.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
int GetTrueVSize() const override
Return the number of local vector true dofs.
Class for parallel grid function.
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Returns ||u_ex - u_h||_L2 in parallel for H1 or L2 elements.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Class for parallel meshes.
Helper class for ParaView visualization data.
void SetHighOrderOutput(bool high_order_output_)
void SetLevelsOfDetail(int levels_of_detail_)
void SetDataFormat(VTKFormat fmt)
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.
real_t Norml2() const
Returns the l2 norm of the vector.
real_t Max() const
Returns the maximal element of the vector.
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 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)
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)
Helper struct to convert a C++ type to an MPI type.
std::array< int, NCMesh::MaxFaceNodes > nodes