56#ifndef MFEM_USE_GINKGO 
   57#error This example requires that MFEM is built with MFEM_USE_GINKGO=YES 
   63int main(
int argc, 
char *argv[])
 
   66   const char *mesh_file = 
"../../data/star.mesh";
 
   68   bool static_cond = 
false;
 
   70   const char *device_config = 
"cpu";
 
   71   bool visualization = 
true;
 
   72   int solver_config = 0;
 
   76   args.
AddOption(&mesh_file, 
"-m", 
"--mesh",
 
   79                  "Finite element order (polynomial degree) or -1 for" 
   80                  " isoparametric space.");
 
   81   args.
AddOption(&static_cond, 
"-sc", 
"--static-condensation", 
"-no-sc",
 
   82                  "--no-static-condensation", 
"Enable static condensation.");
 
   83   args.
AddOption(&pa, 
"-pa", 
"--partial-assembly", 
"-no-pa",
 
   84                  "--no-partial-assembly", 
"Enable Partial Assembly.");
 
   85   args.
AddOption(&device_config, 
"-d", 
"--device",
 
   86                  "Device configuration string, see Device::Configure().");
 
   87   args.
AddOption(&visualization, 
"-vis", 
"--visualization", 
"-no-vis",
 
   89                  "Enable or disable GLVis visualization.");
 
   90   args.
AddOption(&solver_config, 
"-s", 
"--solver-config",
 
   91                  "Solver and preconditioner combination: \n\t" 
   92                  "   0 - Ginkgo solver and Ginkgo preconditioner, \n\t" 
   93                  "   1 - Ginkgo solver and MFEM preconditioner, \n\t" 
   94                  "   2 - MFEM solver and Ginkgo preconditioner, \n\t" 
   95                  "   3 - MFEM solver and MFEM preconditioner.");
 
   96   args.
AddOption(&print_lvl, 
"-pl", 
"--print-level",
 
   97                  "Print level for iterative solver (1 prints every iteration).");
 
  108   Device device(device_config);
 
  114   Mesh *mesh = 
new Mesh(mesh_file, 1, 1);
 
  123         (int)floor(log(50000./mesh->
GetNE())/log(2.)/
dim);
 
  124      for (
int l = 0; l < ref_levels; l++)
 
  141      cout << 
"Using isoparametric FEs: " << fec->
Name() << endl;
 
  148   cout << 
"Number of finite element unknowns: " 
  181   if (pa) { 
a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
 
  188   if (static_cond) { 
a->EnableStaticCondensation(); }
 
  193   a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
 
  195   cout << 
"Size of linear system: " << A->
Height() << endl;
 
  200      switch (solver_config)
 
  205            cout << 
"Using Ginkgo solver + preconditioner...\n";
 
  214            ginkgo_solver.
Mult(B, X);
 
  221            cout << 
"Using Ginkgo solver + MFEM preconditioner...\n";
 
  232            ginkgo_solver.
Mult(B, X);
 
  239            cout << 
"Using MFEM solver + Ginkgo preconditioner...\n";
 
  243            PCG(*A, M, B, X, print_lvl, 400, 1e-12, 0.0);
 
  250            cout << 
"Using MFEM solver + MFEM preconditioner...\n";
 
  253            PCG(*A, M, B, X, print_lvl, 400, 1e-12, 0.0);
 
  266         switch (solver_config)
 
  271               cout << 
"Using Ginkgo solver + preconditioner...\n";
 
  272               MFEM_ABORT(
"Cannot use Ginkgo preconditioner in partial assembly mode.\n" 
  273                          "            Try -s 1 to test Ginkgo solver with an MFEM preconditioner.");
 
  280               cout << 
"Using Ginkgo solver + MFEM preconditioner...\n";
 
  290               ginkgo_solver.
Mult(B, X);
 
  297               cout << 
"Using MFEM solver + Ginkgo preconditioner...\n";
 
  298               MFEM_ABORT(
"Cannot use Ginkgo preconditioner in partial assembly mode.\n" 
  299                          "            Try -s 1 to test Ginkgo solver with an MFEM preconditioner.");
 
  306               cout << 
"Using MFEM solver + MFEM preconditioner...\n";
 
  307               PCG(*A, M, B, X, print_lvl, 400, 1e-12, 0.0);
 
  314         cout << 
"Using MFEM solver + no preconditioner...\n";
 
  315         CG(*A, B, X, print_lvl, 400, 1e-12, 0.0);
 
  320   a->RecoverFEMSolution(X, *
b, x);
 
  324   ofstream mesh_ofs(
"refined.mesh");
 
  325   mesh_ofs.precision(8);
 
  326   mesh->
Print(mesh_ofs);
 
  327   ofstream sol_ofs(
"sol.gf");
 
  328   sol_ofs.precision(8);
 
  337      sol_sock.precision(8);
 
  338      sol_sock << 
"solution\n" << *mesh << x << flush;
 
  345   if (order > 0) { 
delete fec; }
 
 
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
 
int Size() const
Return the logical size of the array.
 
A coefficient that is constant across space and time.
 
Data type for scaled Jacobi-type smoother of sparse matrix.
 
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
 
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
 
Class for domain integration .
 
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
 
virtual const char * Name() const
 
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 ...
 
void SetRelTol(real_t rtol)
 
void SetAbsTol(real_t atol)
 
void SetMaxIter(int max_it)
 
void SetOperator(const Operator &op) override
 
void Mult(const Vector &x, Vector &y) const override
 
void SetPrintLevel(int print_lvl)
 
void SetOperator(const Operator &op) override
 
Class for grid function - Vector with associated FE space.
 
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
 
Arbitrary order H1-conforming (continuous) finite elements.
 
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
 
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
 
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 UniformRefinement(int i, const DSTable &, int *, int *, int *)
 
Pointer to an Operator of a specified type.
 
Operator * Ptr() const
Access the underlying Operator pointer.
 
Jacobi smoothing for a given bilinear form (no matrix necessary).
 
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
 
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.
 
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
 
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
 
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.