75 #include "../common/mfem-common.hpp"
84 int main(
int argc,
char *argv[])
86 #ifdef HYPRE_USING_GPU
87 cout <<
"\nAs of mfem-4.3 and hypre-2.22.0 (July 2021) this miniapp\n"
88 <<
"is NOT supported with the GPU version of hypre.\n\n";
98 const char *mesh_file =
"../../data/inline-quad.mesh";
100 bool visualization =
true;
101 int ser_ref_levels = 0;
102 int dirichlet_level_set_type = -1;
103 int neumann_level_set_type = -1;
104 bool dirichlet_combo =
false;
107 bool include_cut_cell =
false;
110 args.
AddOption(&mesh_file,
"-m",
"--mesh",
111 "Mesh file to use.");
113 "Finite element order (polynomial degree) or -1 for"
114 " isoparametric space.");
115 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
116 "--no-visualization",
117 "Enable or disable GLVis visualization.");
118 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
119 "Number of times to refine the mesh uniformly in serial.");
120 args.
AddOption(&dirichlet_level_set_type,
"-lst",
"--level-set-type",
122 args.
AddOption(&neumann_level_set_type,
"-nlst",
"--neumann-level-set-type",
123 "neumann-level-set-type.");
124 args.
AddOption(&dirichlet_combo,
"-dc",
"--dcombo",
125 "no-dc",
"--no-dcombo",
126 "Combination of two Dirichlet level sets.");
127 args.
AddOption(&ho_terms,
"-ho",
"--high-order",
128 "Additional high-order terms to include");
129 args.
AddOption(&alpha,
"-alpha",
"--alpha",
130 "Nitsche penalty parameter (~1 for 2D, ~10 for 3D).");
131 args.
AddOption(&include_cut_cell,
"-cut",
"--cut",
"-no-cut-cell",
133 "Include or not include elements cut by true boundary.");
146 if (dirichlet_level_set_type <= 0 && neumann_level_set_type <= 0)
148 dirichlet_level_set_type = 1;
150 MFEM_VERIFY((neumann_level_set_type > 0 && ho_terms < 1) ==
false,
151 "Shifted Neumann BC requires extra terms, i.e., -ho >= 1.");
156 if (myid == 0) { device.
Print(); }
159 Mesh mesh(mesh_file, 1, 1);
164 std::cout <<
"Number of elements: " << mesh.
GetNE() << std::endl;
166 MFEM_VERIFY(mesh.
Conforming(),
"AMR capability is not implemented yet!");
169 ParMesh pmesh(MPI_COMM_WORLD, mesh);
174 if (order < 1) { order = 1; }
184 pmesh.SetNodalFESpace(&pfespace_mesh);
186 pmesh.SetNodalGridFunction(&x_mesh);
187 vxyz = *pmesh.GetNodes();
188 int nodes_cnt = vxyz.
Size()/
dim;
189 if (dirichlet_level_set_type == 3)
191 for (
int i = 0; i < nodes_cnt; i++)
194 vxyz(i+nodes_cnt) = (1.+1.e-4)*vxyz(i+nodes_cnt)-1.e-4;
197 pmesh.SetNodes(vxyz);
202 cout <<
"Number of finite element unknowns: " << gtsize << endl;
222 if (dirichlet_level_set_type > 0)
227 filter.
Filter(*dirichlet_dist_coef, level_set_gf);
241 MFEM_VERIFY(dirichlet_level_set_type == 5,
242 "The combo level set example has been only set for"
243 " dirichlet_level_set_type == 5.");
252 if (neumann_level_set_type > 0)
267 for (
int i = 0; i < elem_marker_gf.
Size(); i++)
269 elem_marker_gf(i) = (double)elem_marker[i];
275 "Element Flags", 0, 0,
s,
s,
"Rjmpc");
287 for (
int i = 0; i < sb_dofs.
Size(); i++)
289 face_dofs(sb_dofs[i]) = 1.0;
295 "Shifted Face Dofs", 0,
s,
s,
s,
"Rjmplo");
313 if (dirichlet_level_set_type == 1 || dirichlet_level_set_type == 2 ||
314 dirichlet_level_set_type == 3)
326 filter.
Filter(combo_dist_coef, filt_gf);
335 "Input Level Set", 0, 2*
s,
s,
s,
"Rjmm");
352 "Distance Vector",
s,
s,
s,
s,
"Rjmmpcvv", 1);
357 const int max_elem_attr = pmesh.attributes.Max();
360 bool inactive_elements =
false;
361 for (
int i = 0; i < pmesh.GetNE(); i++)
363 if (!include_cut_cell &&
364 (elem_marker[i] == ShiftedFaceMarker::SBElementType::OUTSIDE ||
365 elem_marker[i] >= ShiftedFaceMarker::SBElementType::CUT))
367 pmesh.SetAttribute(i, max_elem_attr+1);
368 inactive_elements =
true;
370 if (include_cut_cell &&
371 elem_marker[i] == ShiftedFaceMarker::SBElementType::OUTSIDE)
373 pmesh.SetAttribute(i, max_elem_attr+1);
374 inactive_elements =
true;
377 bool inactive_elements_global;
378 MPI_Allreduce(&inactive_elements, &inactive_elements_global, 1, MPI_C_BOOL,
379 MPI_LOR, MPI_COMM_WORLD);
380 if (inactive_elements_global) { ess_elem.
Append(0); }
381 pmesh.SetAttributes();
387 if (dirichlet_level_set_type == 1 || dirichlet_level_set_type == 4 ||
388 dirichlet_level_set_type == 5 || dirichlet_level_set_type == 6 ||
389 dirichlet_level_set_type == 8 ||
390 neumann_level_set_type == 1 || neumann_level_set_type == 7)
394 else if (dirichlet_level_set_type == 2 || neumann_level_set_type == 2)
398 else if (dirichlet_level_set_type == 3)
402 else { MFEM_ABORT(
"RHS function not set for level set type.\n"); }
409 if (dirichlet_level_set_type == 1 || dirichlet_level_set_type >= 4)
414 else if (dirichlet_level_set_type == 2)
419 else if (dirichlet_level_set_type == 3)
434 if (neumann_level_set_type == 1)
439 else if (neumann_level_set_type == 2)
445 else if (neumann_level_set_type == 7)
450 else if (neumann_level_set_type > 0)
452 MFEM_ABORT(
" Normal vector coefficient not implemented for level set.");
459 int ls_cut_marker = ShiftedFaceMarker::SBElementType::CUT;
462 Array<int> bf_dirichlet_marker(0), bf_neumann_marker(0);
464 if (dirichlet_level_set_type > 0)
479 bf_dirichlet_marker.Append(ls_cut_marker);
491 bf_dirichlet_marker.Append(ls_cut_marker);
497 if (neumann_level_set_type > 0)
499 MFEM_VERIFY(!include_cut_cell,
"include_cut_cell option must be set to"
500 " false for Neumann boundary conditions.");
502 &pmesh, *nbcCoef, *dist_vec,
503 *normalbcCoef, elem_marker,
504 ho_terms, include_cut_cell,
506 bf_neumann_marker.
Append(ls_cut_marker);
518 if (dirichlet_level_set_type > 0)
530 ho_terms), ess_shift_bdr);
534 if (neumann_level_set_type > 0)
574 ofstream mesh_ofs(
"diffusion.mesh");
575 mesh_ofs.precision(8);
576 pmesh.PrintAsOne(mesh_ofs);
577 ofstream sol_ofs(
"diffusion.gf");
578 sol_ofs.precision(8);
598 "Solution",
s, 0,
s,
s,
"Rj");
602 if (dirichlet_level_set_type == 2 || dirichlet_level_set_type == 3 ||
603 (dirichlet_level_set_type == -1 && neumann_level_set_type == 2))
608 for (
int i = 0; i < nodes_cnt; i++)
611 pxyz(1) = vxyz(i+nodes_cnt);
612 double exact_val = 0.;
613 if (dirichlet_level_set_type == 2 || neumann_level_set_type == 2)
617 else if (dirichlet_level_set_type == 3)
621 error(i) = std::fabs(x(i) - exact_val);
630 "Error", 2*
s, 0,
s,
s,
"Rj");
636 std::cout <<
"Global L2 error: " << global_error << endl;
641 if (myid == 0) { std::cout << setprecision(10) << norm << std::endl; }
652 delete neumann_dist_coef;
653 delete dirichlet_dist_coef;
654 delete dirichlet_dist_coef_2;
int Size() const
Return the logical size of the array.
Class for domain integration L(v) := (f, v)
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
double rhs_fun_circle(const Vector &x)
f for the Poisson problem (-nabla^2 u = f).
void ExchangeFaceNbrData()
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
Helper class for ParaView visualization data.
Pointer to an Operator of a specified type.
HYPRE_BigInt GlobalTrueVSize() const
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
virtual void ProjectDiscCoefficient(VectorCoefficient &coeff)
Project a discontinuous vector coefficient as a grid function on a continuous finite element space...
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
double dirichlet_velocity_xy_exponent(const Vector &x)
The BoomerAMG solver in hypre.
double homogeneous(const Vector &x)
Boundary conditions - Dirichlet.
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
void ExchangeFaceNbrData()
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
int Append(const T &el)
Append element 'el' to array, resize if necessary.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void ListShiftedFaceDofs(const Array< int > &elem_marker, Array< int > &sface_dof_list) const
void ListEssentialTDofs(const Array< int > &elem_marker, const Array< int > &sface_dof_list, Array< int > &ess_tdof_list, Array< int > &ess_shift_bdr) const
void SetMaxIter(int max_it)
Distance vector to the zero level-set.
Level set coefficient: +1 inside the true domain, -1 outside.
static void Init()
Singleton creation with Mpi::Init();.
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetTime(double t)
Set physical time (for time-dependent simulations)
double rhs_fun_xy_sinusoidal(const Vector &x)
virtual double ComputeL1Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
void SetRelTol(double rtol)
Combination of level sets: +1 inside the true domain, -1 outside.
void Add_Level_Set_Coefficient(Dist_Level_Set_Coefficient &dls_)
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
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...
virtual void ComputeVectorDistance(Coefficient &zero_level_set, ParGridFunction &distance)
double rhs_fun_xy_exponent(const Vector &x)
void SaveAsOne(const char *fname, int precision=16) const
double dirichlet_velocity_xy_sinusoidal(const Vector &x)
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
void normal_vector_1(const Vector &x, Vector &p)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void PrintOptions(std::ostream &out) const
Print the options.
void normal_vector_2(const Vector &x, Vector &p)
Normal vector for level_set_type = 7. Circle centered at [0.5 , 0.6].
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
void SetLevelsOfDetail(int levels_of_detail_)
double traction_xy_exponent(const Vector &x)
Neumann condition for exponent based solution.
void Clear()
Clear the contents of the Mesh.
A general function coefficient.
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)
Vector coefficient defined by a vector GridFunction.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
virtual void Save() override
Class for parallel grid function.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Class for parallel meshes.
void MarkElements(const ParGridFunction &ls_func, Array< int > &elem_marker)
void Filter(ParGridFunction &func, ParGridFunction &ffield)
Arbitrary order "L2-conforming" discontinuous finite elements.
double AvgElementSize(ParMesh &pmesh)
bool Good() const
Return true if the command line options were parsed successfully.