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;
126 bool paraview =
false;
129 args.
AddOption(&mesh_file,
"-m",
"--mesh",
130 "Mesh file to use.");
132 "Finite element order (polynomial degree).");
133 args.
AddOption(&delta_order,
"-do",
"--delta-order",
134 "Order enrichment for DPG test space.");
136 "Epsilon coefficient");
137 args.
AddOption(&ref,
"-ref",
"--num-refinements",
138 "Number of uniform refinements");
139 args.
AddOption(&theta,
"-theta",
"--theta",
140 "Theta parameter for AMR");
141 args.
AddOption(&iprob,
"-prob",
"--problem",
"Problem case"
142 " 0: lshape, 1: General");
144 "Vector Coefficient beta");
145 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
146 "--no-static-condensation",
"Enable static condensation.");
147 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
148 "--no-visualization",
149 "Enable or disable GLVis visualization.");
150 args.
AddOption(¶view,
"-paraview",
"--paraview",
"-no-paraview",
152 "Enable or disable ParaView visualization.");
163 if (iprob > 3) { iprob = 3; }
166 if (
prob == prob_type::EJ ||
prob == prob_type::curved_streamlines ||
167 prob == prob_type::bdr_layer)
169 mesh_file =
"../../data/inline-quad.mesh";
172 Mesh mesh(mesh_file, 1, 1);
174 MFEM_VERIFY(
dim > 1,
"Dimension = 1 is not supported in this example");
210 ParMesh pmesh(MPI_COMM_WORLD, mesh);
244 int test_order = order+delta_order;
266 trial_fes.
Append(sigma_fes);
267 trial_fes.
Append(hatu_fes);
268 trial_fes.
Append(hatf_fes);
273 a->StoreMatrices(
true);
277 TrialSpace::u_space, TestSpace::v_space);
281 TrialSpace::sigma_space, TestSpace::v_space);
285 TrialSpace::u_space, TestSpace::tau_space);
289 TrialSpace::sigma_space, TestSpace::tau_space);
293 TrialSpace::hatu_space, TestSpace::tau_space);
297 TrialSpace::hatf_space, TestSpace::v_space);
312 TestSpace::v_space, TestSpace::v_space);
315 TestSpace::v_space, TestSpace::v_space);
318 TestSpace::v_space, TestSpace::v_space);
321 TestSpace::tau_space, TestSpace::tau_space);
324 TestSpace::tau_space, TestSpace::tau_space);
327 if (
prob == prob_type::sinusoidal ||
328 prob == prob_type::curved_streamlines)
350 std::cout <<
" Ref |"
357 std::cout <<
" Residual |"
359 <<
" CG it |" << endl;
360 std::cout << std::string((
exact_known) ? 72 : 50,
'-')
364 if (static_cond) {
a->EnableStaticCondensation(); }
384 for (
int it = 0; it<=ref; it++)
397 if (
prob == prob_type::EJ)
415 int n = ess_tdof_list_uhat.
Size();
416 int m = ess_tdof_list_fhat.
Size();
418 for (
int j = 0; j < n; j++)
420 ess_tdof_list[j] = ess_tdof_list_uhat[j]
424 for (
int j = 0; j < m; j++)
426 ess_tdof_list[j+n] = ess_tdof_list_fhat[j]
450 a->FormLinearSystem(ess_tdof_list,x,Ah,X,B);
463 M.SetDiagonalBlock(0,amg0);
464 M.SetDiagonalBlock(1,amg1);
470 M.SetDiagonalBlock(skip,amg2);
483 M.SetDiagonalBlock(skip+1,prec);
494 a->RecoverFEMSolution(X,x);
495 Vector & residuals =
a->ComputeResidual(x);
500 real_t gresidual = residual * residual;
503 MPI_MAX, MPI_COMM_WORLD);
505 MPI_SUM, MPI_COMM_WORLD);
507 gresidual = sqrt(gresidual);
510 for (
int iel = 0; iel<pmesh.
GetNE(); iel++)
512 if (residuals[iel] > theta * maxresidual)
514 elements_to_refine.
Append(iel);
532 L2Error = sqrt(u_err*u_err + sigma_err*sigma_err);
533 rate_err = (it) ?
dim*log(err0/L2Error)/log((
real_t)dof0/dofs) : 0.0;
536 real_t rate_res = (it) ?
dim*log(res0/gresidual)/log((
real_t)dof0/dofs) : 0.0;
543 std::ios oldState(
nullptr);
544 oldState.copyfmt(std::cout);
545 std::cout << std::right << std::setw(5) << it <<
" | "
546 << std::setw(10) << dof0 <<
" | ";
549 std::cout << std::setprecision(3) << std::setw(10)
550 << std::scientific << err0 <<
" | "
551 << std::setprecision(2)
552 << std::setw(6) << std::fixed << rate_err <<
" | " ;
554 std::cout << std::setprecision(3)
555 << std::setw(10) << std::scientific << res0 <<
" | "
556 << std::setprecision(2)
557 << std::setw(6) << std::fixed << rate_res <<
" | "
558 << std::setw(6) << std::fixed << num_iter <<
" | "
560 std::cout.copyfmt(oldState);
565 const char * keys = (it == 0 &&
dim == 2) ?
"cgRjmlk\n" :
nullptr;
569 "Numerical u", 0,0, 500, 500, keys);
571 "Numerical flux", 501,0,500, 500, keys);
587 for (
int i =0; i<trial_fes.
Size(); i++)
589 trial_fes[i]->Update(
false);
626 if (X.
Size() == 3) { z = X[2]; }
640 real_t denom = exp(-r2) - exp(-r1);
642 real_t g1 = exp(r2*(x-1.));
643 real_t g2 = exp(r1*(x-1.));
645 return g * cos(M_PI * y)/denom;
655 MFEM_ABORT(
"Wrong code path");
666 if (X.
Size() == 3) { z = X[2]; }
674 for (
int i = 0; i<du.
Size(); i++)
676 du[i] = M_PI * cos(
alpha);
685 real_t denom = exp(-r2) - exp(-r1);
687 real_t g1 = exp(r2*(x-1.));
689 real_t g2 = exp(r1*(x-1.));
694 real_t u_x = g_x * cos(M_PI * y)/denom;
695 real_t u_y = -M_PI * g * sin(M_PI*y)/denom;
710 MFEM_ABORT(
"Wrong code path");
720 if (X.
Size() == 3) { z = X[2]; }
727 return - M_PI*M_PI *
u * X.
Size();
735 real_t denom = exp(-r2) - exp(-r1);
737 real_t g1 = exp(r2*(x-1.));
740 real_t g2 = exp(r1*(x-1.));
744 real_t g_xx = g1_xx - g2_xx;
746 real_t u = g * cos(M_PI * y)/denom;
747 real_t u_xx = g_xx * cos(M_PI * y)/denom;
748 real_t u_yy = -M_PI * M_PI *
u;
760 MFEM_ABORT(
"Wrong code path");
786 for (
int i = 0; i<hatf.
Size(); i++)
788 hatf[i] = beta_val[i] *
u -
sigma[i];
803 for (
int i = 0; i<du.
Size(); i++)
805 s += beta_val[i] * du[i];
812 if (
prob == prob_type::bdr_layer)
839 if (
prob == prob_type::curved_streamlines)
843 beta_val(0) = exp(x)*sin(y);
844 beta_val(1) = exp(x)*cos(y);
857 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.
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.
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
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_)
virtual void Save() override
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)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
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.