85static const char *enum_str[] =
 
  109int main(
int argc, 
char *argv[])
 
  116   const char *mesh_file = 
"../../data/inline-quad.mesh";
 
  122   bool static_cond = 
false;
 
  125   bool visualization = 
true;
 
  127   bool paraview = 
false;
 
  130   args.
AddOption(&mesh_file, 
"-m", 
"--mesh",
 
  131                  "Mesh file to use.");
 
  133                  "Finite element order (polynomial degree).");
 
  134   args.
AddOption(&delta_order, 
"-do", 
"--delta-order",
 
  135                  "Order enrichment for DPG test space.");
 
  137                  "Epsilon coefficient");
 
  138   args.
AddOption(&ref, 
"-ref", 
"--num-refinements",
 
  139                  "Number of uniform refinements");
 
  140   args.
AddOption(&theta, 
"-theta", 
"--theta",
 
  141                  "Theta parameter for AMR");
 
  142   args.
AddOption(&iprob, 
"-prob", 
"--problem", 
"Problem case" 
  143                  " 0: lshape, 1: General");
 
  145                  "Vector Coefficient beta");
 
  146   args.
AddOption(&static_cond, 
"-sc", 
"--static-condensation", 
"-no-sc",
 
  147                  "--no-static-condensation", 
"Enable static condensation.");
 
  148   args.
AddOption(&visualization, 
"-vis", 
"--visualization", 
"-no-vis",
 
  149                  "--no-visualization",
 
  150                  "Enable or disable GLVis visualization.");
 
  151   args.
AddOption(¶view, 
"-paraview", 
"--paraview", 
"-no-paraview",
 
  153                  "Enable or disable ParaView visualization.");
 
  154   args.
AddOption(&visport, 
"-p", 
"--send-port", 
"Socket for GLVis.");
 
  165   if (iprob > 3) { iprob = 3; }
 
  168   if (
prob == prob_type::EJ || 
prob == prob_type::curved_streamlines ||
 
  169       prob == prob_type::bdr_layer)
 
  171      mesh_file = 
"../../data/inline-quad.mesh";
 
  174   Mesh mesh(mesh_file, 1, 1);
 
  176   MFEM_VERIFY(
dim > 1, 
"Dimension = 1 is not supported in this example");
 
  212   ParMesh pmesh(MPI_COMM_WORLD, mesh);
 
  246   int test_order = order+delta_order;
 
  268   trial_fes.
Append(sigma_fes);
 
  269   trial_fes.
Append(hatu_fes);
 
  270   trial_fes.
Append(hatf_fes);
 
  275   a->StoreMatrices(
true);
 
  279                         TrialSpace::u_space, TestSpace::v_space);
 
  283                         TrialSpace::sigma_space, TestSpace::v_space);
 
  287                         TrialSpace::u_space, TestSpace::tau_space);
 
  291                         TrialSpace::sigma_space, TestSpace::tau_space);
 
  295                         TrialSpace::hatu_space, TestSpace::tau_space);
 
  299                         TrialSpace::hatf_space, TestSpace::v_space);
 
  314                        TestSpace::v_space, TestSpace::v_space);
 
  317                        TestSpace::v_space, TestSpace::v_space);
 
  320                        TestSpace::v_space, TestSpace::v_space);
 
  323                        TestSpace::tau_space, TestSpace::tau_space);
 
  326                        TestSpace::tau_space, TestSpace::tau_space);
 
  329   if (
prob == prob_type::sinusoidal ||
 
  330       prob == prob_type::curved_streamlines)
 
  352      std::cout << 
"  Ref |" 
  356         std::cout << 
"  L2 Error  |" 
  359      std::cout << 
"  Residual  |" 
  361                << 
" CG it  |" << std::endl;
 
  362      std::cout << std::string((
exact_known) ? 72 : 50,
'-')
 
  366   if (static_cond) { 
a->EnableStaticCondensation(); }
 
  386   for (
int it = 0; it<=ref; it++)
 
  399         if (
prob == prob_type::EJ)
 
  417      int n = ess_tdof_list_uhat.
Size();
 
  418      int m = ess_tdof_list_fhat.
Size();
 
  420      for (
int j = 0; j < n; j++)
 
  422         ess_tdof_list[j] = ess_tdof_list_uhat[j]
 
  426      for (
int j = 0; j < m; j++)
 
  428         ess_tdof_list[j+n] = ess_tdof_list_fhat[j]
 
  452      a->FormLinearSystem(ess_tdof_list,x,Ah,X,B);
 
  465         M.SetDiagonalBlock(0,amg0);
 
  466         M.SetDiagonalBlock(1,amg1);
 
  472      M.SetDiagonalBlock(skip,amg2);
 
  485      M.SetDiagonalBlock(skip+1,prec);
 
  496      a->RecoverFEMSolution(X,x);
 
  497      Vector & residuals = 
a->ComputeResidual(x);
 
  502      real_t gresidual = residual * residual;
 
  505                    MPI_MAX, MPI_COMM_WORLD);
 
  507                    MPI_SUM, MPI_COMM_WORLD);
 
  509      gresidual = sqrt(gresidual);
 
  512      for (
int iel = 0; iel<pmesh.
GetNE(); iel++)
 
  514         if (residuals[iel] > theta * maxresidual)
 
  516            elements_to_refine.
Append(iel);
 
  534         L2Error = sqrt(u_err*u_err + sigma_err*sigma_err);
 
  535         rate_err = (it) ? 
dim*log(err0/L2Error)/log((
real_t)dof0/dofs) : 0.0;
 
  538      real_t rate_res = (it) ? 
dim*log(res0/gresidual)/log((
real_t)dof0/dofs) : 0.0;
 
  545         std::ios oldState(
nullptr);
 
  546         oldState.copyfmt(std::cout);
 
  547         std::cout << std::right << std::setw(5) << it << 
" | " 
  548                   << std::setw(10) <<  dof0 << 
" | ";
 
  551            std::cout << std::setprecision(3) << std::setw(10)
 
  552                      << std::scientific <<  err0 << 
" | " 
  553                      << std::setprecision(2)
 
  554                      << std::setw(6) << std::fixed << rate_err << 
" | " ;
 
  556         std::cout << std::setprecision(3)
 
  557                   << std::setw(10) << std::scientific <<  res0 << 
" | " 
  558                   << std::setprecision(2)
 
  559                   << std::setw(6) << std::fixed << rate_res << 
" | " 
  560                   << std::setw(6) << std::fixed << num_iter << 
" | " 
  562         std::cout.copyfmt(oldState);
 
  567         const char * keys = (it == 0 && 
dim == 2) ? 
"cgRjmlk\n" : 
nullptr;
 
  570                        "Numerical u", 0,0, 500, 500, keys);
 
  572                        "Numerical flux", 501,0,500, 500, keys);
 
  588      for (
int i =0; i<trial_fes.
Size(); i++)
 
  590         trial_fes[i]->Update(
false);
 
 
  627   if (X.
Size() == 3) { z = X[2]; }
 
  641         real_t denom = exp(-r2) - exp(-r1);
 
  643         real_t g1 = exp(r2*(x-1.));
 
  644         real_t g2 = exp(r1*(x-1.));
 
  646         return g * cos(M_PI * y)/denom;
 
  656         MFEM_ABORT(
"Wrong code path");
 
 
  667   if (X.
Size() == 3) { z = X[2]; }
 
  675         for (
int i = 0; i<du.
Size(); i++)
 
  677            du[i] = M_PI * cos(
alpha);
 
  686         real_t denom = exp(-r2) - exp(-r1);
 
  688         real_t g1 = exp(r2*(x-1.));
 
  690         real_t g2 = exp(r1*(x-1.));
 
  695         real_t u_x = g_x * cos(M_PI * y)/denom;
 
  696         real_t u_y = -M_PI * g * sin(M_PI*y)/denom;
 
  711         MFEM_ABORT(
"Wrong code path");
 
 
  721   if (X.
Size() == 3) { z = X[2]; }
 
  728         return - M_PI*M_PI * 
u * X.
Size();
 
  736         real_t denom = exp(-r2) - exp(-r1);
 
  738         real_t g1 = exp(r2*(x-1.));
 
  741         real_t g2 = exp(r1*(x-1.));
 
  745         real_t g_xx = g1_xx - g2_xx;
 
  747         real_t u = g * cos(M_PI * y)/denom;
 
  748         real_t u_xx = g_xx * cos(M_PI * y)/denom;
 
  749         real_t u_yy = -M_PI * M_PI * 
u;
 
  761         MFEM_ABORT(
"Wrong code path");
 
 
  787   for (
int i = 0; i<hatf.
Size(); i++)
 
  789      hatf[i] = beta_val[i] * 
u - 
sigma[i];
 
 
  804   for (
int i = 0; i<du.
Size(); i++)
 
  806      s += beta_val[i] * du[i];
 
 
  813   if (
prob == prob_type::bdr_layer)
 
 
  840   if (
prob == prob_type::curved_streamlines)
 
  844      beta_val(0) = exp(x)*sin(y);
 
  845      beta_val(1) = exp(x)*cos(y);
 
 
  858   for (
int i = 0; i < pmesh->
GetNE(); i++)
 
 
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.
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
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 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.
real_t GetElementVolume(int i)
void EnsureNCMesh(bool simplices_nonconforming=false)
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 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.
Matrix coefficient defined as the outer product of two vector coefficients.
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.
void GetElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be al...
ParMesh * GetParMesh() const
void Update(bool want_transform=true) override
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)
ParFiniteElementSpace * ParFESpace() const
void MakeRef(FiniteElementSpace *f, real_t *v) override
Make the ParGridFunction reference external data on a new FiniteElementSpace.
void Update() override
Transform by the Space UpdateMatrix (e.g., on Mesh change).
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
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...
Vector coefficient defined as a product of scalar and vector coefficients.
A general vector function coefficient.
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
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.
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)
real_t exact_hatu(const Vector &X)
real_t exact_laplacian_u(const Vector &X)
void exact_gradu(const Vector &X, Vector &du)
real_t bdr_data(const Vector &X)
void exact_hatf(const Vector &X, Vector &hatf)
real_t exact_u(const Vector &X)
void beta_function(const Vector &X, Vector &beta_val)
void exact_sigma(const Vector &X, Vector &sigma)
real_t f_exact(const Vector &X)
void setup_test_norm_coeffs(ParGridFunction &c1_gf, ParGridFunction &c2_gf)
Helper struct to convert a C++ type to an MPI type.