58 int main(
int argc,
char *argv[])
61 Mpi::Init(argc, argv);
62 int num_procs = Mpi::WorldSize();
63 int myid = Mpi::WorldRank();
67 const char *mesh_file =
"../data/star.mesh";
68 int ser_ref_levels = 2;
69 int par_ref_levels = 1;
73 bool slu_solver =
false;
74 bool sp_solver =
false;
75 bool cpardiso_solver =
false;
76 bool visualization = 1;
79 args.
AddOption(&mesh_file,
"-m",
"--mesh",
81 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
82 "Number of times to refine the mesh uniformly in serial.");
83 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
84 "Number of times to refine the mesh uniformly in parallel.");
86 "Finite element order (polynomial degree) or -1 for"
87 " isoparametric space.");
89 "Number of desired eigenmodes.");
91 "Random seed used to initialize LOBPCG.");
92 #ifdef MFEM_USE_SUPERLU
93 args.
AddOption(&slu_solver,
"-slu",
"--superlu",
"-no-slu",
94 "--no-superlu",
"Use the SuperLU Solver.");
96 #ifdef MFEM_USE_STRUMPACK
97 args.
AddOption(&sp_solver,
"-sp",
"--strumpack",
"-no-sp",
98 "--no-strumpack",
"Use the STRUMPACK Solver.");
100 #ifdef MFEM_USE_MKL_CPARDISO
101 args.
AddOption(&cpardiso_solver,
"-cpardiso",
"--cpardiso",
"-no-cpardiso",
102 "--no-cpardiso",
"Use the MKL CPardiso Solver.");
104 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
105 "--no-visualization",
106 "Enable or disable GLVis visualization.");
108 if (slu_solver && sp_solver)
111 cout <<
"WARNING: Both SuperLU and STRUMPACK have been selected,"
112 <<
" please choose either one." << endl
113 <<
" Defaulting to SuperLU." << endl;
137 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
143 for (
int lev = 0; lev < ser_ref_levels; lev++)
154 for (
int lev = 0; lev < par_ref_levels; lev++)
179 cout <<
"Number of unknowns: " << size << endl;
220 #if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
222 #ifdef MFEM_USE_SUPERLU
228 #ifdef MFEM_USE_STRUMPACK
243 if (!slu_solver && !sp_solver && !cpardiso_solver)
251 #ifdef MFEM_USE_SUPERLU
262 #ifdef MFEM_USE_STRUMPACK
276 #ifdef MFEM_USE_MKL_CPARDISO
280 cpardiso->
SetMatrixType(CPardisoSolver::MatType::REAL_STRUCTURE_SYMMETRIC);
281 cpardiso->SetPrintLevel(1);
282 cpardiso->SetOperator(*A);
310 ostringstream mesh_name, mode_name;
311 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
313 ofstream mesh_ofs(mesh_name.str().c_str());
314 mesh_ofs.precision(8);
315 pmesh->
Print(mesh_ofs);
317 for (
int i=0; i<nev; i++)
322 mode_name <<
"mode_" << setfill(
'0') << setw(2) << i <<
"."
323 << setfill(
'0') << setw(6) << myid;
325 ofstream mode_ofs(mode_name.str().c_str());
326 mode_ofs.precision(8);
338 mode_sock.precision(8);
340 for (
int i=0; i<nev; i++)
344 cout <<
"Eigenmode " << i+1 <<
'/' << nev
345 <<
", Lambda = " << eigenvalues[i] << endl;
351 mode_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
352 <<
"solution\n" << *pmesh << x << flush
353 <<
"window_title 'Eigenmode " << i+1 <<
'/' << nev
354 <<
", Lambda = " << eigenvalues[i] <<
"'" << endl;
359 cout <<
"press (q)uit or (c)ontinue --> " << flush;
362 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
377 #if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
int Size() const
Return the logical size of the array.
MPI_Comm GetComm() const
MPI communicator.
A coefficient that is constant across space and time.
HYPRE_BigInt GlobalTrueVSize() const
void SetFromCommandLine()
virtual void Save(std::ostream &out) const
void SetPreconditioner(Solver &precond)
void SetMassMatrix(Operator &M)
Abstract parallel finite element space.
void SetPrintLevel(int logging)
void SetSymmetricPattern(bool sym)
void GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
int close()
Close the socketstream.
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
The BoomerAMG solver in hypre.
void SetPrintFactorStatistics(bool print_stat)
void SetMaxIter(int max_iter)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void SetColumnPermutation(superlu::ColPerm col_perm)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetPrintLevel(int print_level)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetPrintSolveStatistics(bool print_stat)
void SetPrintStatistics(bool print_stat)
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
void SetRandomSeed(int s)
void SetKrylovSolver(strumpack::KrylovSolver method)
void SetOperator(Operator &A)
void PrintOptions(std::ostream &out) const
Print the options.
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
void GetNodes(Vector &node_coord) const
void SetPrecondUsageMode(int pcg_mode)
Arbitrary order H1-conforming (continuous) finite elements.
void SetMatrixType(MatType mat_type)
Set the matrix type.
void Print(std::ostream &out=mfem::out) const override
Class for parallel grid function.
Wrapper for hypre's ParCSR matrix class.
void SetNumModes(int num_eigs)
Class for parallel meshes.
void Solve()
Solve the eigenproblem.
MKL Parallel Direct Sparse Solver for Clusters.
bool Good() const
Return true if the command line options were parsed successfully.