36#error This examples requires that MFEM is build with MFEM_USE_SLEPC=YES
42int main(
int argc,
char *argv[])
51 const char *mesh_file =
"../../data/star.mesh";
52 int ser_ref_levels = 2;
53 int par_ref_levels = 1;
57 bool slu_solver =
false;
58 bool sp_solver =
false;
59 bool visualization = 1;
60 bool use_slepc =
true;
61 const char *slepcrc_file =
"";
62 const char *device_config =
"cpu";
65 args.
AddOption(&mesh_file,
"-m",
"--mesh",
67 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
68 "Number of times to refine the mesh uniformly in serial.");
69 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
70 "Number of times to refine the mesh uniformly in parallel.");
72 "Finite element order (polynomial degree) or -1 for"
73 " isoparametric space.");
75 "Number of desired eigenmodes.");
77 "Random seed used to initialize LOBPCG.");
78#ifdef MFEM_USE_SUPERLU
79 args.
AddOption(&slu_solver,
"-slu",
"--superlu",
"-no-slu",
80 "--no-superlu",
"Use the SuperLU Solver.");
82#ifdef MFEM_USE_STRUMPACK
83 args.
AddOption(&sp_solver,
"-sp",
"--strumpack",
"-no-sp",
84 "--no-strumpack",
"Use the STRUMPACK Solver.");
86 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
88 "Enable or disable GLVis visualization.");
89 args.
AddOption(&use_slepc,
"-useslepc",
"--useslepc",
"-no-slepc",
90 "--no-slepc",
"Use or not SLEPc to solve the eigenvalue problem");
91 args.
AddOption(&slepcrc_file,
"-slepcopts",
"--slepcopts",
92 "SlepcOptions file to use.");
93 args.
AddOption(&device_config,
"-d",
"--device",
94 "Device configuration string, see Device::Configure().");
96 if (slu_solver && sp_solver)
99 cout <<
"WARNING: Both SuperLU and STRUMPACK have been selected,"
100 <<
" please choose either one." << endl
101 <<
" Defaulting to SuperLU." << endl;
124 Device device(device_config);
125 if (myid == 0) { device.
Print(); }
133 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
139 for (
int lev = 0; lev < ser_ref_levels; lev++)
150 for (
int lev = 0; lev < par_ref_levels; lev++)
175 cout <<
"Number of unknowns: " << size << endl;
203 a->EliminateEssentialBCDiag(ess_bdr, 1.0);
219 a->ParallelAssemble(Ah);
220 if (!use_slepc) { Ah.Get(A); }
222 Ah.SetOperatorOwner(
false);
225 if (!use_slepc) {Mh.
Get(M); }
229#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
231#ifdef MFEM_USE_SUPERLU
237#ifdef MFEM_USE_STRUMPACK
254 if (!slu_solver && !sp_solver)
262#ifdef MFEM_USE_SUPERLU
273#ifdef MFEM_USE_STRUMPACK
281 strumpack->
SetMatching(strumpack::MatchingJob::NONE);
330 for (
int i=0; i<nev; i++)
341 ostringstream mesh_name, mode_name;
342 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
344 ofstream mesh_ofs(mesh_name.str().c_str());
345 mesh_ofs.precision(8);
346 pmesh->
Print(mesh_ofs);
348 for (
int i=0; i<nev; i++)
361 mode_name <<
"mode_" << setfill(
'0') << setw(2) << i <<
"."
362 << setfill(
'0') << setw(6) << myid;
364 ofstream mode_ofs(mode_name.str().c_str());
365 mode_ofs.precision(8);
377 mode_sock.precision(8);
379 for (
int i=0; i<nev; i++)
383 cout <<
"Eigenmode " << i+1 <<
'/' << nev
384 <<
", Lambda = " << eigenvalues[i] << endl;
398 mode_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
399 <<
"solution\n" << *pmesh << x << flush
400 <<
"window_title 'Eigenmode " << i+1 <<
'/' << nev
401 <<
", Lambda = " << eigenvalues[i] <<
"'" << endl;
406 cout <<
"press (q)uit or (c)ontinue --> " << flush;
409 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
427#if defined(MFEM_USE_SUPERLU) || defined(MFEM_USE_STRUMPACK)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
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.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Arbitrary order H1-conforming (continuous) finite elements.
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_level)
void SetMassMatrix(Operator &M)
void SetPrintLevel(int logging)
void SetPreconditioner(Solver &precond)
void GetEigenvalues(Array< real_t > &eigenvalues) const
Collect the converged eigenvalues.
void SetOperator(Operator &A)
void SetNumModes(int num_eigs)
void Solve()
Solve the eigenproblem.
void SetPrecondUsageMode(int pcg_mode)
void SetRandomSeed(int s)
void SetMaxIter(int max_iter)
const HypreParVector & GetEigenvector(unsigned int i) const
Extract a single eigenvector.
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.
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.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Type
Enumeration defining IDs for some classes derived from Operator.
@ Hypre_ParCSR
ID for class HypreParMatrix.
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
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.
HYPRE_BigInt GlobalTrueVSize() const
int GetTrueVSize() const override
Return the number of local vector true dofs.
Class for parallel grid function.
void Save(std::ostream &out) const override
void Distribute(const Vector *tv)
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Wrapper for PETSc's matrix class.
void SetOperator(const Operator &op)
Set the operator/matrix.
void SetMatching(strumpack::MatchingJob job)
Configure static pivoting for stability.
void SetCompression(strumpack::CompressionType type)
Select compression for sparse data types.
void SetPrintFactorStatistics(bool print_stat)
Set up verbose printing during the factor step.
void SetKrylovSolver(strumpack::KrylovSolver method)
Set the Krylov solver method to use.
void SetFromCommandLine()
Set options that were captured from the command line.
void SetPrintSolveStatistics(bool print_stat)
Set up verbose printing during the solve step.
void SetReorderingStrategy(strumpack::ReorderingStrategy method)
Set matrix reordering strategy.
@ TARGET_REAL
The eigenvalues with the real component closest to the target value.
void SetOperators(const PetscParMatrix &op, const PetscParMatrix &opB)
Set operators for generalized eigenvalue problem.
void SetWhichEigenpairs(Which which)
Set the which eigenvalues the solver will target and the order they will be indexed in.
void SetTarget(real_t target)
Set the target value for the eigenpairs you want when using SlepcEigenSolver::TARGET_MAGNITUDE or Sle...
void GetEigenvalue(unsigned int i, real_t &lr) const
Get the ith eigenvalue after the system has been solved.
void SetNumModes(int num_eigs)
Set the number of eigenmodes to compute.
void SetSpectralTransformation(SpectralTransformation transformation)
Set the spectral transformation strategy for acceletating convergenvce. Both SlepcEigenSolver::SHIFT ...
@ SHIFT_INVERT
Utilize the shift and invert strategy.
void GetEigenvector(unsigned int i, Vector &vr) const
Get the ith eigenvector after the system has been solved.
void Solve()
Solve the eigenvalue problem for the specified number of eigenvalues.
void SetColumnPermutation(superlu::ColPerm col_perm)
Specify how to permute the columns of the matrix.
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 SetPrintStatistics(bool print_stat)
Specify whether to print the solver statistics (default true)
int close()
Close the socketstream.
@ PARMETIS
Sequential ordering on structure of using the PARMETIS package.
void MFEMInitializeSlepc()