34#error This example requires that MFEM is built with MFEM_USE_AMGX=YES
40int main(
int argc,
char *argv[])
43 const char *mesh_file =
"../../data/star.mesh";
45 bool static_cond =
false;
47 const char *device_config =
"cpu";
48 bool visualization =
true;
50 bool amgx_solver =
true;
51 const char* amgx_json_file =
"";
54 args.
AddOption(&mesh_file,
"-m",
"--mesh",
57 "Finite element order (polynomial degree) or -1 for"
58 " isoparametric space.");
59 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
60 "--no-static-condensation",
"Enable static condensation.");
61 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
62 "--no-partial-assembly",
"Enable Partial Assembly.");
63 args.
AddOption(&amgx_lib,
"-amgx",
"--amgx-lib",
"-no-amgx",
64 "--no-amgx-lib",
"Use AmgX in example.");
65 args.
AddOption(&amgx_json_file,
"--amgx-file",
"--amgx-file",
66 "AMGX solver config file (overrides --amgx-solver, --amgx-verbose)");
67 args.
AddOption(&amgx_solver,
"--amgx-solver",
"--amgx-solver",
68 "--amgx-preconditioner",
"--amgx-preconditioner",
69 "Configure AMGX as solver or preconditioner.");
70 args.
AddOption(&device_config,
"-d",
"--device",
71 "Device configuration string, see Device::Configure().");
72 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
74 "Enable or disable GLVis visualization.");
85 Device device(device_config);
91 Mesh mesh(mesh_file, 1, 1);
100 (int)floor(log(50000./mesh.
GetNE())/log(2.)/
dim);
101 for (
int l = 0; l < ref_levels; l++)
121 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
129 cout <<
"Number of finite element unknowns: "
162 if (pa) {
a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
169 if (static_cond) {
a.EnableStaticCondensation(); }
174 a.FormLinearSystem(ess_tdof_list, x,
b, A, X, B);
176 cout <<
"Size of linear system: " << A->
Height() << endl;
185 PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
189 CG(*A, B, X, 1, 400, 1e-12, 0.0);
192 else if (amgx_lib && strcmp(amgx_json_file,
"") == 0)
194 bool amgx_verbose =
false;
197 PCG(*A, amgx, B, X, 1, 200, 1e-12, 0.0);
199 else if (amgx_lib && strcmp(amgx_json_file,
"") != 0)
221#ifndef MFEM_USE_SUITESPARSE
224 PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
228 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
230 umf_solver.
Mult(B, X);
235 a.RecoverFEMSolution(X,
b, x);
239 ofstream mesh_ofs(
"refined.mesh");
240 mesh_ofs.precision(8);
241 mesh.
Print(mesh_ofs);
242 ofstream sol_ofs(
"sol.gf");
243 sol_ofs.precision(8);
252 sol_sock.precision(8);
253 sol_sock <<
"solution\n" << mesh << x << flush;
40int main(
int argc,
char *argv[]) {
…}
MFEM wrapper for Nvidia's multigrid library, AmgX (github.com/NVIDIA/AMGX)
@ EXTERNAL
Configure will be read from a specified file.
void InitSerial()
Initialize the AmgX library for serial execution once the solver configuration has been established t...
void SetConvergenceCheck(bool setConvergenceCheck_=true)
Add a check for convergence after applying Mult.
void ReadParameters(const std::string config, CONFIG_SRC source)
Read in the AmgX parameters either through a file or directly through a properly formated string....
void Mult(const Vector &b, Vector &x) const override
Utilize the AmgX library to solve the linear system where the "matrix" is the AMG approximation to th...
void SetOperator(const Operator &op) override
Sets the Operator that is going to be solved via AmgX. Supports operators based on either an MFEM Spa...
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.
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 ...
Data type for Gauss-Seidel smoother of sparse matrix.
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.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
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.
Direct sparse solver using UMFPACK.
real_t Control[UMFPACK_CONTROL]
void SetOperator(const Operator &op) override
Factorize the given Operator op which must be a SparseMatrix.
void Mult(const Vector &b, Vector &x) const override
Direct solution of the linear system using UMFPACK.
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.