47#ifndef MFEM_USE_SUPERLU
48#error This example requires that MFEM is built with MFEM_USE_SUPERLU=YES
54int main(
int argc,
char *argv[])
63 const char *mesh_file =
"../../data/star.mesh";
65 const char *device_config =
"cpu";
66 bool visualization =
true;
73 args.
AddOption(&mesh_file,
"-m",
"--mesh",
76 "Finite element order (polynomial degree) or -1 for"
77 " isoparametric space.");
78 args.
AddOption(&device_config,
"-d",
"--device",
79 "Device configuration string, see Device::Configure().");
80 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
82 "Enable or disable GLVis visualization.");
83 args.
AddOption(&slu_colperm,
"-cp",
"--colperm",
84 "SuperLU Column Permutation Method: 0-NATURAL, 1-MMD-ATA "
85 "2-MMD_AT_PLUS_A, 3-COLAMD, 4-METIS_AT_PLUS_A, 5-PARMETIS "
87 args.
AddOption(&slu_rowperm,
"-rp",
"--rowperm",
88 "SuperLU Row Permutation Method: 0-NOROWPERM, 1-LargeDiag");
89 args.
AddOption(&slu_iterref,
"-ir",
"--iterref",
90 "SuperLU Iterative Refinement: 0-NOREFINE, 1-Single, "
92 args.
AddOption(&slu_npdep,
"-npdep",
"--npdepth",
93 "Depth of 3D parition for SuperLU (>= 7.2.0)");
111 Device device(device_config);
112 if (myid == 0) { device.
Print(); }
117 Mesh mesh(mesh_file, 1, 1);
126 (int)floor(log(1000./mesh.
GetNE())/log(2.)/
dim);
127 for (
int l = 0; l < ref_levels; l++)
136 ParMesh pmesh(MPI_COMM_WORLD, mesh);
139 int par_ref_levels = 2;
140 for (
int l = 0; l < par_ref_levels; l++)
162 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
174 cout <<
"Number of finite element unknowns: " << size << endl;
217 a.FormLinearSystem(ess_tdof_list, x,
b, A, X, B);
225 if (slu_colperm == 0)
229 else if (slu_colperm == 1)
233 else if (slu_colperm == 2)
237 else if (slu_colperm == 3)
241 else if (slu_colperm == 4)
245 else if (slu_colperm == 5)
249 else if (slu_colperm == 6)
254 if (slu_rowperm == 0)
258 else if (slu_rowperm == 1)
260#ifdef MFEM_USE_SUPERLU5
267 if (slu_iterref == 0)
271 else if (slu_iterref == 1)
275 else if (slu_iterref == 2)
279 else if (slu_iterref == 3)
298 ostringstream mesh_name, sol_name;
299 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
300 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
302 ofstream mesh_ofs(mesh_name.str().c_str());
303 mesh_ofs.precision(8);
304 pmesh.
Print(mesh_ofs);
306 ofstream sol_ofs(sol_name.str().c_str());
307 sol_ofs.precision(8);
317 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
318 sol_sock.precision(8);
319 sol_sock <<
"solution\n" << pmesh << x << flush;
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
Arbitrary order H1-conforming (continuous) finite elements.
Wrapper for hypre's ParCSR matrix class.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
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.
void GetNodes(Vector &node_coord) const
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of 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().
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
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.
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
Class for parallel grid function.
void Save(std::ostream &out) const override
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
void SetColumnPermutation(superlu::ColPerm col_perm)
Specify how to permute the columns of the matrix.
void SetRowPermutation(superlu::RowPerm row_perm)
Specify how to permute the rows of the matrix.
void Mult(const Vector &x, Vector &y) const
Factor and solve the linear system .
void SetSymmetricPattern(bool sym)
Specify whether the matrix has a symmetric pattern to avoid extra work (default false)
void SetOperator(const Operator &op)
Set the operator/matrix.
void SetIterativeRefine(superlu::IterRefine iter_ref)
Specify how to handle iterative refinement.
void SetPrintStatistics(bool print_stat)
Specify whether to print the solver statistics (default true)
@ NATURAL
Natural ordering.
@ MMD_ATA
Minimum degree ordering on structure of .
@ METIS_AT_PLUS_A
Sequential ordering on structure of using the METIS package.
@ MMD_AT_PLUS_A
Minimum degree ordering on structure of .
@ PARMETIS
Sequential ordering on structure of using the PARMETIS package.
@ ZOLTAN
Use the Zoltan library from Sandia to define the column ordering.
@ COLAMD
Approximate minimum degree column ordering.
@ SLU_SINGLE
Iterative refinement accumulating residuals in a float.
@ SLU_DOUBLE
Iterative refinement accumulating residuals in a double.
@ SLU_EXTRA
Iterative refinement accumulating residuals in a higher precision variable.
@ NOREFINE
No interative refinement.
@ NOROWPERM
No row permutation.
@ LargeDiag_MC64
Duff/Koster algorithm to make the diagonals large compared to the off-diagonals. Use LargeDiag for Su...