43 #include "utilities/MPIDataTypes.hpp"
46 using namespace parelag;
52 int main(
int argc,
char *argv[])
60 Timer total_timer = TimeManager::AddTimer(
"Program Execution -- Total");
61 Timer init_timer = TimeManager::AddTimer(
"Initial Setup");
64 cout <<
"-- This is an example of using a geometric-like multilevel "
65 "hierarchy, constructed by ParELAG,\n"
66 "-- to solve respective finite element H(curl) and H(div) forms: \n"
67 "(alpha curl u, curl v) + (beta u, v);\n"
68 "(alpha div u, div v) + (beta u, v).\n\n";
71 const char *xml_file_c = NULL;
76 args.
AddOption(&xml_file_c,
"-f",
"--xml-file",
77 "XML parameter list (an XML file with detailed parameters).");
78 args.
AddOption(&hcurl,
"-curl",
"--hcurl",
"-div",
"--hdiv",
79 "Whether the H(curl) or H(div) form is being solved.");
80 args.
AddOption(&visualize,
"-v",
"--visualize",
"-nv",
"--no-visualize",
81 "Use GLVis to visualize the final solution and the "
83 args.
AddOption(&tolSVD,
"-s",
"--svd-tol",
84 "SVD tolerance. It is used for filtering out local linear "
85 "dependencies in the basis construction and extension "
86 "process in ParELAG. Namely, right singular vectors with "
87 "singular values smaller than this tolerance are removed.");
101 xml_file_c =
"MultilevelHcurlSolver_cube_example_parameters.xml";
105 xml_file_c =
"MultilevelHdivSolver_cube_example_parameters.xml";
109 cout <<
"No XML parameter list provided! Using default "
110 << xml_file_c <<
"." << endl;
117 string xml_file(xml_file_c);
120 unique_ptr<ParameterList> master_list;
121 ifstream xml_in(xml_file);
126 cerr <<
"ERROR: Cannot read from input file: " << xml_file <<
".\n";
130 SimpleXMLParameterListReader reader;
131 master_list = reader.GetParameterList(xml_in);
135 ParameterList& prob_list = master_list->Sublist(
"Problem parameters",
true);
138 const string meshfile = prob_list.Get(
"Mesh file",
"");
143 int ser_ref_levels = prob_list.Get(
"Serial refinement levels", -1);
147 const int par_ref_levels = prob_list.Get(
"Parallel refinement levels", 2);
152 const int amge_levels = prob_list.Get(
"AMGe levels", 0);
155 const int feorder = prob_list.Get(
"Finite element order", 0);
159 const int upscalingOrder = prob_list.Get(
"Upscaling order", 0);
165 vector<int> par_ess_attr = prob_list.Get(
"Essential attributes",
171 vector<double> alpha_vals = prob_list.Get(
"alpha values",
172 vector<double> {1.0});
177 vector<double> beta_vals = prob_list.Get(
"beta values",
178 vector<double> {1.0});
181 auto list_of_solvers = prob_list.Get<list<string>>(
"List of linear solvers");
183 ostringstream mesh_msg;
186 mesh_msg <<
'\n' << string(50,
'*') <<
'\n'
187 <<
"* Mesh: " << meshfile <<
"\n*\n"
188 <<
"* FE order: " << feorder <<
'\n'
189 <<
"* Upscaling order: " << upscalingOrder <<
"\n*\n";
193 shared_ptr<ParMesh> pmesh;
197 cout <<
"\nReading and refining serial mesh...\n";
198 cout <<
"Times to refine mesh in serial: " << ser_ref_levels <<
".\n";
201 ifstream imesh(meshfile);
206 cerr <<
"ERROR: Cannot open mesh file: " << meshfile <<
".\n";
211 auto mesh = make_unique<Mesh>(imesh,
true,
true);
214 for (
int l = 0; l < ser_ref_levels; ++l)
218 cout <<
"Refining mesh in serial: " << l + 1 <<
"...\n";
220 mesh->UniformRefinement();
223 if (ser_ref_levels < 0)
226 for (; mesh->GetNE() < 6 * num_ranks; ++ser_ref_levels)
230 cout <<
"Refining mesh in serial: " << ser_ref_levels + 1
233 mesh->UniformRefinement();
239 cout <<
"Times refined mesh in serial: " << ser_ref_levels <<
".\n";
240 cout <<
"Building and refining parallel mesh...\n";
241 cout <<
"Times to refine mesh in parallel: " << par_ref_levels
243 mesh_msg <<
"* Serial refinements: " << ser_ref_levels <<
'\n'
244 <<
"* Coarse mesh size: " << mesh->GetNE() <<
"\n*\n";
247 pmesh = make_shared<ParMesh>(MPI_COMM_WORLD, *mesh);
251 MFEM_VERIFY(par_ess_attr.size() <= 1 ||
252 par_ess_attr.size() == (unsigned) pmesh->bdr_attributes.Max(),
253 "Incorrect size of the essential attributes vector in parameters"
255 vector<Array<int>> ess_attr(1);
256 ess_attr[0].SetSize(pmesh->bdr_attributes.Max());
257 if (par_ess_attr.size() == 0)
261 else if (par_ess_attr.size() == 1)
263 ess_attr[0] = par_ess_attr[0];
267 for (
unsigned i = 0; i < par_ess_attr.size(); ++i)
269 ess_attr[0][i] = par_ess_attr[i];
274 MFEM_VERIFY(alpha_vals.size() <= 1 ||
275 alpha_vals.size() == (unsigned) pmesh->attributes.Max(),
276 "Incorrect size of the 'alpha' local values vector in parameters"
278 MFEM_VERIFY(alpha_vals.size() <= 1 ||
279 alpha_vals.size() == (unsigned) pmesh->attributes.Max(),
280 "Incorrect size of the 'alpha' local values vector in parameters"
285 if (alpha_vals.size() == 0)
289 else if (alpha_vals.size() == 1)
291 alpha = alpha_vals[0];
295 for (
unsigned i = 0; i < alpha_vals.size(); ++i)
297 alpha(i+1) = alpha_vals[i];
301 if (beta_vals.size() == 0)
305 else if (beta_vals.size() == 1)
311 for (
unsigned i = 0; i < beta_vals.size(); ++i)
313 beta(i+1) = beta_vals[i];
318 const int nDimensions = pmesh->Dimension();
325 MFEM_VERIFY(nDimensions == 3,
"Only 3D problems are currently supported.");
327 const int nLevels = amge_levels <= 0 ? par_ref_levels + 1 : amge_levels;
328 MFEM_VERIFY(nLevels <= par_ref_levels + 1,
329 "Number of AMGe levels too high relative to parallel"
331 vector<int> level_nElements(nLevels);
332 for (
int l = 0; l < par_ref_levels; ++l)
336 cout <<
"Refining mesh in parallel: " << l + 1
337 << (par_ref_levels - l > nLevels ?
" (not in hierarchy)"
341 if (par_ref_levels - l < nLevels)
343 level_nElements[par_ref_levels - l] = pmesh->GetNE();
345 pmesh->UniformRefinement();
347 level_nElements[0] = pmesh->GetNE();
351 cout <<
"Times refined mesh in parallel: " << par_ref_levels <<
".\n";
355 size_t local_num_elmts = pmesh->GetNE(), global_num_elmts;
356 MPI_Reduce(&local_num_elmts, &global_num_elmts, 1, GetMPIType<size_t>(0),
357 MPI_SUM, 0, MPI_COMM_WORLD);
360 mesh_msg <<
"* Parallel refinements: " << par_ref_levels <<
'\n'
361 <<
"* Fine mesh size: " << global_num_elmts <<
'\n'
362 <<
"* Total levels: " << nLevels <<
'\n'
363 << string(50,
'*') <<
"\n\n";
369 cout << mesh_msg.str();
374 Timer agg_timer = TimeManager::AddTimer(
"Mesh Agglomeration -- Total");
375 Timer agg0_timer = TimeManager::AddTimer(
"Mesh Agglomeration -- Level 0");
378 cout <<
"Agglomerating topology for " << nLevels - 1
379 <<
" coarse levels...\n";
386 MFEMRefinedMeshPartitioner partitioner(nDimensions);
387 vector<shared_ptr<AgglomeratedTopology>> topology(nLevels);
391 cout <<
"Agglomerating level: 0...\n";
394 topology[0] = make_shared<AgglomeratedTopology>(pmesh, nDimensions);
398 cout <<
"Level 0 global number of mesh entities: "
400 GetNumberGlobalTrueEntities((AgglomeratedTopology::Entity)0);
401 for (
int j = 1; j <= nDimensions; ++j)
402 cout <<
", " << topology[0]->
403 GetNumberGlobalTrueEntities((AgglomeratedTopology::Entity)j);
409 for (
int l = 0; l < nLevels - 1; ++l)
411 Timer aggl_timer = TimeManager::AddTimer(std::string(
"Mesh "
412 "Agglomeration -- Level ").
413 append(std::to_string(l+1)));
414 Array<int> partitioning(topology[l]->GetNumberLocalEntities(AT_elem));
415 partitioner.Partition(topology[l]->GetNumberLocalEntities(AT_elem),
416 level_nElements[l + 1], partitioning);
420 cout <<
"Agglomerating level: " << l + 1 <<
"...\n";
423 topology[l + 1] = topology[l]->CoarsenLocalPartitioning(partitioning,
427 cout <<
"Level " << l + 1 <<
" global number of mesh entities: "
429 GetNumberGlobalTrueEntities((AgglomeratedTopology::Entity)0);
430 for (
int j = 1; j <= nDimensions; ++j)
432 cout <<
", " << topology[l + 1]->
433 GetNumberGlobalTrueEntities((AgglomeratedTopology::Entity)j);
441 if (visualize && nDimensions <= 3)
443 for (
int l = 1; l < nLevels; ++l)
445 ShowTopologyAgglomeratedElements(topology[l].
get(), pmesh.get());
451 Timer derham_timer = TimeManager::AddTimer(
"DeRhamSequence Construction -- "
453 Timer derham0_timer = TimeManager::AddTimer(
"DeRhamSequence Construction -- "
457 cout <<
"Building the fine-level de Rham sequence...\n";
460 vector<shared_ptr<DeRhamSequence>> sequence(topology.size());
462 const int jform = DeRhamSequence::GetForm(nDimensions,
463 hcurl ? DeRhamSequence::HCURL :
464 DeRhamSequence::HDIV);
465 if (nDimensions == 3)
467 sequence[0] = make_shared<DeRhamSequence3D_FE>(topology[0], pmesh.get(),
468 feorder,
true,
false);
472 MFEM_VERIFY(nDimensions == 2,
"Only 2D or 3D problems are supported "
473 "by the utilized ParELAG.");
476 MFEM_ABORT(
"No H(curl) 2D interpretation of form 1 is implemented.");
478 sequence[0] = make_shared<DeRhamSequence2D_Hdiv_FE>(topology[0],
479 pmesh.get(), feorder,
493 sequence[0]->SetjformStart(0);
495 DeRhamSequenceFE *DRSequence_FE = sequence[0]->FemSequence();
496 MFEM_VERIFY(DRSequence_FE,
497 "Failed to obtain the fine-level de Rham sequence.");
501 cout <<
"Level 0 global number of dofs: "
502 << DRSequence_FE->GetDofHandler(0)->GetDofTrueDof().
504 for (
int j = 1; j <= nDimensions; ++j)
506 cout <<
", " << DRSequence_FE->GetDofHandler(j)->GetDofTrueDof().
514 cout <<
"Setting coefficients and computing fine-level local "
518 DRSequence_FE->ReplaceMassIntegrator(AT_elem, jform,
519 make_unique<VectorFEMassIntegrator>(beta),
false);
520 if (hcurl && nDimensions == 3)
522 DRSequence_FE->ReplaceMassIntegrator(AT_elem, jform + 1,
523 make_unique<VectorFEMassIntegrator>(
alpha),
true);
527 DRSequence_FE->ReplaceMassIntegrator(AT_elem, jform + 1,
528 make_unique<MassIntegrator>(
alpha),
true);
533 cout <<
"Interpolating and setting polynomial targets...\n";
536 DRSequence_FE->SetUpscalingTargets(nDimensions, upscalingOrder);
537 derham0_timer.Stop();
541 cout <<
"Building the coarse-level de Rham sequences...\n";
544 for (
int l = 0; l < nLevels - 1; ++l)
546 Timer derhaml_timer = TimeManager::AddTimer(std::string(
"DeRhamSequence "
547 "Construction -- Level ").
548 append(std::to_string(l+1)));
551 cout <<
"Building the level " << l + 1 <<
" de Rham sequences...\n";
554 sequence[l]->SetSVDTol(tolSVD);
555 sequence[l + 1] = sequence[l]->Coarsen();
559 auto DRSequence = sequence[l + 1];
560 cout <<
"Level " << l + 1 <<
" global number of dofs: "
561 << DRSequence->GetDofHandler(0)->GetDofTrueDof().
563 for (
int j = 1; j <= nDimensions; ++j)
565 cout <<
", " << DRSequence->GetDofHandler(j)->GetDofTrueDof().
573 Timer assemble_timer = TimeManager::AddTimer(
"Fine Matrix Assembly");
576 cout <<
"Assembling the fine-level system...\n";
586 auto rhsform = make_unique<LinearForm>(fespace);
589 unique_ptr<Vector> rhs = move(rhsform);
593 auto solgf = make_unique<GridFunction>(fespace);
594 solgf->ProjectCoefficient(solcoeff);
595 unique_ptr<Vector> sol = move(solgf);
598 const SharingMap& hcurlhdiv_dofTrueDof =
599 sequence[0]->GetDofHandler(jform)->GetDofTrueDof();
602 Vector B(hcurlhdiv_dofTrueDof.GetTrueLocalSize());
605 shared_ptr<HypreParMatrix> A;
624 auto M1 = sequence[0]->ComputeMassOperator(jform),
625 M2 = sequence[0]->ComputeMassOperator(jform + 1);
626 auto D1 = sequence[0]->GetDerivativeOperator(jform);
633 auto spA = ToUnique(
Add(*M1, *ToUnique(
RAP(*D1, *M2, *D1))));
638 sequence[0]->GetDofHandler(jform)->MarkDofsOnSelectedBndr(ess_attr[0],
641 for (
int i = 0; i < spA->Height(); ++i)
645 spA->EliminateRowCol(i, sol->
Elem(i), *rhs);
649 A = Assemble(hcurlhdiv_dofTrueDof, *spA, hcurlhdiv_dofTrueDof);
650 hcurlhdiv_dofTrueDof.Assemble(*rhs, B);
654 cout <<
"A size: " << A->GetGlobalNumRows() <<
'x'
655 << A->GetGlobalNumCols() <<
'\n' <<
" A NNZ: " << A->NNZ() <<
'\n';
657 MFEM_VERIFY(B.Size() == A->Height(),
658 "Matrix and vector size are incompatible.");
659 assemble_timer.Stop();
662 Timer solvers_timer = TimeManager::AddTimer(
"Solvers -- Total");
665 cout <<
"\nRunning fine-level solvers...\n\n";
669 auto lib = SolverLibrary::CreateLibrary(
670 master_list->Sublist(
"Preconditioner Library"));
673 for (
const auto& solver_name : list_of_solvers)
675 Timer solver_timer = TimeManager::AddTimer(std::string(
"Solver \"").
677 append(
"\" -- Total"));
679 auto solver_factory = lib->GetSolverFactory(solver_name);
680 auto solver_state = solver_factory->GetDefaultState();
681 solver_state->SetDeRhamSequence(sequence[0]);
682 solver_state->SetBoundaryLabels(ess_attr);
683 solver_state->SetForms({jform});
686 Timer build_timer = TimeManager::AddTimer(std::string(
"Solver \"").
688 append(
"\" -- Build"));
691 cout <<
"Building solver \"" << solver_name <<
"\"...\n";
693 unique_ptr<Solver> solver = solver_factory->BuildSolver(A, *solver_state);
697 Timer pre_timer = TimeManager::AddTimer(std::string(
"Solver \"").
699 append(
"\" -- Pre-solve"));
702 cout <<
"Solving system with \"" << solver_name <<
"\"...\n";
707 Vector X(A->Width()), x(sequence[0]->GetNumberOfDofs(jform));
716 double local_norm = tmp.
Norml2();
717 local_norm *= local_norm;
719 MPI_Reduce(&local_norm, &global_norm, 1, GetMPIType(local_norm),
720 MPI_SUM, 0, MPI_COMM_WORLD);
723 cout <<
"Initial residual l2 norm: " << sqrt(global_norm) <<
'\n';
729 Timer solve_timer = TimeManager::AddTimer(std::string(
"Solver \"").
731 append(
"\" -- Solve"));
735 Timer post_timer = TimeManager::AddTimer(std::string(
"Solver \"").
737 append(
"\" -- Post-solve"));
744 double local_norm = tmp.
Norml2();
745 local_norm *= local_norm;
747 MPI_Reduce(&local_norm, &global_norm, 1, GetMPIType(local_norm),
748 MPI_SUM, 0, MPI_COMM_WORLD);
751 cout <<
"Final residual l2 norm: " << sqrt(global_norm) <<
'\n';
757 cout <<
"Solver \"" << solver_name <<
"\" finished.\n";
763 hcurlhdiv_dofTrueDof.Distribute(X, x);
764 MultiVector tmp(x.GetData(), 1, x.Size());
765 sequence[0]->show(jform, tmp);
769 solvers_timer.Stop();
772 TimeManager::Print();
776 cout <<
"\nFinished.\n";
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(), hypre will be finalized automatically at program exit.
double & Elem(int i)
Access Vector entries. Index i = 0 .. size-1.
double Norml2() const
Returns the l2 norm of the vector.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
static int WorldSize()
Return the size of MPI_COMM_WORLD.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
void rhsfunc(const Vector &, Vector &)
double f(const Vector &xvec)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
static void Init()
Singleton creation with Mpi::Init();.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
void PrintUsage(std::ostream &out) const
Print the usage message.
A general vector function coefficient.
double p(const Vector &x, double t)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
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...
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
void PrintOptions(std::ostream &out) const
Print the options.
void bdrfunc(const Vector &, Vector &)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
bool Good() const
Return true if the command line options were parsed successfully.