49 int main(
int argc,
char *argv[])
53 MPI_Init(&argc, &argv);
54 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
55 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
58 const char *mesh_file =
"../data/beam-tri.mesh";
62 bool visualization = 1;
67 args.
AddOption(&mesh_file,
"-m",
"--mesh",
70 "Finite element order (polynomial degree).");
72 "Number of desired eigenmodes.");
74 "Random seed used to initialize LOBPCG.");
75 args.
AddOption(&amg_elast,
"-elast",
"--amg-for-elasticity",
"-sys",
77 "Use the special AMG elasticity solver (GM/LN approaches), "
78 "or standard AMG for systems (unknown approach).");
79 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
81 "Enable or disable GLVis visualization.");
82 args.
AddOption(&adios2,
"-adios2",
"--adios2-streams",
"-no-adios2",
83 "--no-adios2-streams",
84 "Save data using adios2 streams.");
103 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
109 cerr <<
"\nInput mesh should have at least two materials!"
110 <<
" (See schematic in ex12p.cpp)\n"
129 (int)floor(log(1000./mesh->
GetNE())/log(2.)/
dim);
130 for (
int l = 0; l < ref_levels; l++)
142 int par_ref_levels = 1;
143 for (
int l = 0; l < par_ref_levels; l++)
157 const bool use_nodal_fespace = pmesh->
NURBSext && !amg_elast;
158 if (use_nodal_fespace)
171 cout <<
"Number of unknowns: " << size << endl
172 <<
"Assembling: " << flush;
186 lambda(0) = lambda(1)*50;
201 cout <<
"matrix ... " << flush;
215 cout <<
"done." << endl;
264 if (!use_nodal_fespace)
272 ostringstream mesh_name, mode_name;
273 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
275 ofstream mesh_ofs(mesh_name.str().c_str());
276 mesh_ofs.precision(8);
277 pmesh->
Print(mesh_ofs);
279 for (
int i=0; i<nev; i++)
284 mode_name <<
"mode_" << setfill(
'0') << setw(2) << i <<
"."
285 << setfill(
'0') << setw(6) << myid;
287 ofstream mode_ofs(mode_name.str().c_str());
288 mode_ofs.precision(8);
296 #ifdef MFEM_USE_ADIOS2
299 std::string postfix(mesh_file);
300 postfix.erase(0, std::string(
"../data/").size() );
301 postfix +=
"_o" + std::to_string(order);
305 pmesh->
Print(adios2output);
306 for (
int i=0; i<nev; i++)
310 x.
Save(adios2output,
"mode_" + std::to_string(i));
323 for (
int i=0; i<nev; i++)
327 cout <<
"Eigenmode " << i+1 <<
'/' << nev
328 <<
", Lambda = " << eigenvalues[i] << endl;
334 mode_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n"
335 <<
"solution\n" << *pmesh << x << flush
336 <<
"window_title 'Eigenmode " << i+1 <<
'/' << nev
337 <<
", Lambda = " << eigenvalues[i] <<
"'" << endl;
342 cout <<
"press (q)uit or (c)ontinue --> " << flush;
345 MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
HYPRE_BigInt GlobalTrueVSize() const
int GetNE() const
Returns number of elements.
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 GetEigenvalues(Array< double > &eigenvalues) const
Collect the converged eigenvalues.
int close()
Close the socketstream.
The BoomerAMG solver in hypre.
void SetMaxIter(int max_iter)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
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.
virtual void Print(std::ostream &out=mfem::out) const
void SetElasticityOptions(ParFiniteElementSpace *fespace)
void SetNodalFESpace(FiniteElementSpace *nfes)
void PrintUsage(std::ostream &out) const
Print the usage message.
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 SetRandomSeed(int s)
NURBSExtension * NURBSext
Optional NURBS mesh extension.
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
void DegreeElevate(int rel_degree, int degree=16)
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 GetNodes(Vector &node_coord) const
void SetPrecondUsageMode(int pcg_mode)
void SetSystemsOptions(int dim, bool order_bynodes=false)
Arbitrary order H1-conforming (continuous) finite elements.
Class for parallel grid function.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Wrapper for hypre's ParCSR matrix class.
void SetNumModes(int num_eigs)
Class for parallel meshes.
void Solve()
Solve the eigenproblem.
Array< int > attributes
A list of all unique element attributes used by the Mesh.
bool Good() const
Return true if the command line options were parsed successfully.