Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
3 : //*
4 : //* All rights reserved, see COPYRIGHT for full restrictions
5 : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 : //*
7 : //* Licensed under LGPL 2.1, please see LICENSE for details
8 : //* https://www.gnu.org/licenses/lgpl-2.1.html
9 :
10 : #ifdef MOOSE_MFEM_ENABLED
11 :
12 : #include "MFEMProblem.h"
13 : #include "MFEMVariable.h"
14 : #include "MFEMIndicator.h"
15 : #include "MFEMSubMesh.h"
16 : #include "MFEMFunctorMaterial.h"
17 : #include "MFEMExecutedObject.h"
18 : #include "MFEMVectorUtils.h"
19 : #include "libmesh/string_to_enum.h"
20 :
21 : #include <vector>
22 : #include <algorithm>
23 : #include <map>
24 : #include <set>
25 : #include <deque>
26 : #include <sstream>
27 :
28 : registerMooseObject("MooseApp", MFEMProblem);
29 :
30 : InputParameters
31 7668 : MFEMProblem::validParams()
32 : {
33 7668 : InputParameters params = ExternalProblem::validParams();
34 15336 : params.addClassDescription("Problem type for building and solving the finite element problem "
35 : "using the MFEM finite element library.");
36 30672 : MooseEnum numeric_types("real complex", "real");
37 23004 : params.addParam<MooseEnum>("numeric_type", numeric_types, "Number type used for the problem");
38 :
39 15336 : return params;
40 7668 : }
41 :
42 1736 : MFEMProblem::MFEMProblem(const InputParameters & params)
43 3472 : : ExternalProblem(params), _num_type{static_cast<int>(getParam<MooseEnum>("numeric_type"))}
44 : {
45 : // Initialise Hypre for all MFEM problems.
46 1736 : mfem::Hypre::Init();
47 : // Disable multithreading for all MFEM problems (including any libMesh or MFEM subapps).
48 1736 : libMesh::libMeshPrivateData::_n_threads = 1;
49 1736 : setMesh();
50 1736 : }
51 :
52 : void
53 1400 : MFEMProblem::initialSetup()
54 : {
55 1400 : ExternalProblem::initialSetup();
56 :
57 : // MFEM indicators create their estimators during addIndicator(); markers still need an explicit
58 : // setup pass because they are no longer initialized through the libMesh/MOOSE user-object path.
59 1400 : std::vector<MFEMRefinementMarker *> markers;
60 1400 : theWarehouse().query().condition<AttribSystem>("Marker").queryInto(markers);
61 1426 : for (auto marker : markers)
62 26 : marker->initialSetup();
63 1400 : }
64 :
65 : void
66 13846 : MFEMProblem::execute(const ExecFlagType & exec_type)
67 : {
68 13846 : setCurrentExecuteOnFlag(exec_type);
69 13846 : executeMFEMObjects(exec_type);
70 :
71 13846 : ExternalProblem::execute(exec_type);
72 13846 : }
73 :
74 : void
75 1736 : MFEMProblem::setMesh()
76 : {
77 1736 : auto pmesh = mesh().getMFEMParMeshPtr();
78 1736 : getProblemData().pmesh = pmesh;
79 1736 : getProblemData().comm = pmesh->GetComm();
80 1736 : getProblemData().num_procs = pmesh->GetNRanks();
81 1736 : getProblemData().myid = pmesh->GetMyRank();
82 1736 : }
83 :
84 : void
85 1070 : MFEMProblem::addMFEMPreconditioner(const std::string & user_object_name,
86 : const std::string & name,
87 : InputParameters & parameters)
88 : {
89 2140 : addObject<MFEMSolverBase>(user_object_name, name, parameters);
90 1070 : }
91 :
92 : void
93 26 : MFEMProblem::addIndicator(const std::string & user_object_name,
94 : const std::string & name,
95 : InputParameters & parameters)
96 : {
97 52 : auto estimator = addObject<MFEMIndicator>(user_object_name, name, parameters).front();
98 :
99 : // construct the estimator itself
100 26 : estimator->createEstimator();
101 26 : }
102 :
103 : void
104 26 : MFEMProblem::addMarker(const std::string & user_object_name,
105 : const std::string & name,
106 : InputParameters & parameters)
107 : {
108 26 : getProblemData().refiner =
109 78 : addObject<MFEMRefinementMarker>(user_object_name, name, parameters).front();
110 26 : }
111 :
112 : void
113 892 : MFEMProblem::addMFEMSolver(const std::string & user_object_name,
114 : const std::string & name,
115 : InputParameters & parameters)
116 : {
117 892 : getProblemData().jacobian_solver =
118 2676 : addObject<MFEMSolverBase>(user_object_name, name, parameters).front();
119 892 : }
120 :
121 : void
122 1361 : MFEMProblem::addMFEMNonlinearSolver(unsigned int nl_max_its,
123 : mfem::real_t nl_abs_tol,
124 : mfem::real_t nl_rel_tol,
125 : unsigned int print_level)
126 : {
127 : // TODO: allow users to specify other mfem::IterativeSolvers
128 1361 : auto nl_solver = std::make_shared<mfem::NewtonSolver>(getComm());
129 :
130 : // Defaults to one iteration, without further nonlinear iterations
131 1361 : nl_solver->SetRelTol(nl_rel_tol);
132 1361 : nl_solver->SetAbsTol(nl_abs_tol);
133 1361 : nl_solver->SetMaxIter(nl_max_its);
134 1361 : nl_solver->SetPrintLevel(print_level);
135 1361 : getProblemData().nonlinear_solver = nl_solver;
136 1361 : }
137 :
138 : void
139 1216 : MFEMProblem::addBoundaryCondition(const std::string & bc_name,
140 : const std::string & name,
141 : InputParameters & parameters)
142 : {
143 2432 : auto bc = addObject<MFEMBoundaryCondition>(bc_name, name, parameters).front();
144 1216 : const auto & mfem_bc = *bc;
145 :
146 1216 : if (dynamic_cast<const MFEMIntegratedBC *>(&mfem_bc))
147 : {
148 56 : auto integrated_bc = std::dynamic_pointer_cast<MFEMIntegratedBC>(bc);
149 : auto eqsys =
150 56 : std::dynamic_pointer_cast<Moose::MFEM::EquationSystem>(getProblemData().eqn_system);
151 56 : if (eqsys)
152 56 : eqsys->AddIntegratedBC(std::move(integrated_bc));
153 : else
154 0 : mooseError("Cannot add integrated BC with name '" + name +
155 : "' because there is no corresponding equation system.");
156 56 : }
157 1160 : else if (dynamic_cast<const MFEMComplexIntegratedBC *>(&mfem_bc))
158 : {
159 12 : auto integrated_bc = std::dynamic_pointer_cast<MFEMComplexIntegratedBC>(bc);
160 : auto eqsys =
161 12 : std::dynamic_pointer_cast<Moose::MFEM::ComplexEquationSystem>(getProblemData().eqn_system);
162 12 : if (eqsys)
163 12 : eqsys->AddComplexIntegratedBC(std::move(integrated_bc));
164 : else
165 0 : mooseError("Cannot add complex integrated BC with name '" + name +
166 : "' because there is no corresponding equation system.");
167 12 : }
168 1148 : else if (dynamic_cast<const MFEMComplexEssentialBC *>(&mfem_bc))
169 : {
170 62 : auto essential_bc = std::dynamic_pointer_cast<MFEMComplexEssentialBC>(bc);
171 : auto eqsys =
172 62 : std::dynamic_pointer_cast<Moose::MFEM::ComplexEquationSystem>(getProblemData().eqn_system);
173 62 : if (eqsys)
174 62 : eqsys->AddComplexEssentialBCs(std::move(essential_bc));
175 : else
176 0 : mooseError("Cannot add boundary condition with name '" + name +
177 : "' because there is no corresponding equation system.");
178 62 : }
179 1086 : else if (dynamic_cast<const MFEMEssentialBC *>(&mfem_bc))
180 : {
181 1086 : auto essential_bc = std::dynamic_pointer_cast<MFEMEssentialBC>(bc);
182 : auto eqsys =
183 1086 : std::dynamic_pointer_cast<Moose::MFEM::EquationSystem>(getProblemData().eqn_system);
184 1086 : if (eqsys)
185 1086 : eqsys->AddEssentialBC(std::move(essential_bc));
186 : else
187 0 : mooseError("Cannot add boundary condition with name '" + name +
188 : "' because there is no corresponding equation system.");
189 1086 : }
190 : else
191 : {
192 0 : mooseError("Unsupported bc of type '", bc_name, "' and name '", name, "' detected.");
193 : }
194 1216 : }
195 :
196 : void
197 0 : MFEMProblem::addMaterial(const std::string &, const std::string &, InputParameters &)
198 : {
199 0 : mooseError(
200 : "MFEM materials must be added through the 'FunctorMaterials' block and not 'Materials'");
201 : }
202 :
203 : void
204 277 : MFEMProblem::addFunctorMaterial(const std::string & material_name,
205 : const std::string & name,
206 : InputParameters & parameters)
207 : {
208 570 : addObject<MFEMFunctorMaterial>(material_name, name, parameters);
209 261 : }
210 :
211 : void
212 2378 : MFEMProblem::addFESpace(const std::string & type,
213 : const std::string & name,
214 : InputParameters & parameters)
215 : {
216 4756 : auto & mfem_fespace = *addObject<MFEMFESpace>(type, name, parameters).front();
217 :
218 : // Register fespace and associated fe collection.
219 2378 : getProblemData().fecs.Register(name, mfem_fespace.getFEC());
220 2378 : getProblemData().fespaces.Register(name, mfem_fespace.getFESpace());
221 2378 : }
222 :
223 : void
224 1650 : MFEMProblem::addVariable(const std::string & var_type,
225 : const std::string & var_name,
226 : InputParameters & parameters)
227 : {
228 1650 : addGridFunction(var_type, var_name, parameters);
229 : // MOOSE variables store DoFs for the trial variable and its time derivatives up to second order;
230 : // MFEM GridFunctions store data for only one set of DoFs each, so we must add additional
231 : // GridFunctions for time derivatives.
232 1650 : if (isTransient())
233 : {
234 : const auto time_derivative_var_name =
235 151 : getMFEMObject<MFEMVariable>("MooseVariableBase", var_name).getTimeDerivativeName();
236 151 : getProblemData().time_derivative_map.addTimeDerivativeAssociation(var_name,
237 : time_derivative_var_name);
238 151 : addGridFunction(var_type, time_derivative_var_name, parameters);
239 151 : }
240 1650 : }
241 :
242 : void
243 3414 : MFEMProblem::addGridFunction(const std::string & var_type,
244 : const std::string & var_name,
245 : InputParameters & parameters)
246 : {
247 :
248 3414 : if (var_type == "MFEMVariable" || var_type == "MFEMComplexVariable")
249 : {
250 : // Add MFEM variable directly.
251 3378 : if (var_type == "MFEMComplexVariable")
252 186 : addObject<MFEMComplexVariable>(var_type, var_name, parameters);
253 : else
254 9948 : addObject<MFEMVariable>(var_type, var_name, parameters);
255 : }
256 : else
257 : {
258 : // Add MOOSE variable.
259 36 : ExternalProblem::addVariable(var_type, var_name, parameters);
260 :
261 : // Add MFEM variable indirectly ("gridfunction").
262 36 : InputParameters mfem_variable_params = addMFEMFESpaceFromMOOSEVariable(parameters);
263 144 : addObject<MFEMVariable>("MFEMVariable", var_name, mfem_variable_params);
264 36 : }
265 :
266 : // Register gridfunction.
267 3414 : if (var_type == "MFEMComplexVariable")
268 : {
269 : MFEMComplexVariable & mfem_variable =
270 62 : getMFEMObject<MFEMComplexVariable>("MooseVariableBase", var_name);
271 62 : getProblemData().cmplx_gridfunctions.Register(var_name, mfem_variable.getComplexGridFunction());
272 62 : mfem_variable.declareCoefficients();
273 : }
274 : else // must be real, but may have been set up indirectly from a MOOSE variable
275 : {
276 3352 : MFEMVariable & mfem_variable = getMFEMObject<MFEMVariable>("MooseVariableBase", var_name);
277 3352 : getProblemData().gridfunctions.Register(var_name, mfem_variable.getGridFunction());
278 3352 : mfem_variable.declareCoefficients();
279 : }
280 3414 : }
281 :
282 : void
283 1457 : MFEMProblem::addAuxVariable(const std::string & var_type,
284 : const std::string & var_name,
285 : InputParameters & parameters)
286 : {
287 : // We handle MFEM AuxVariables just like MFEM Variables, except
288 : // we do not add additional GridFunctions for time derivatives.
289 1457 : addGridFunction(var_type, var_name, parameters);
290 1457 : }
291 :
292 : void
293 430 : MFEMProblem::addAuxKernel(const std::string & kernel_name,
294 : const std::string & name,
295 : InputParameters & parameters)
296 : {
297 860 : addObject<MFEMExecutedObject>(kernel_name, name, parameters);
298 430 : }
299 :
300 : void
301 1697 : MFEMProblem::addKernel(const std::string & kernel_name,
302 : const std::string & name,
303 : InputParameters & parameters)
304 : {
305 3394 : auto kernel = addObject<MFEMKernel>(kernel_name, name, parameters).front();
306 1697 : const auto & kernel_object = *kernel;
307 :
308 1697 : if (dynamic_cast<const MFEMComplexKernel *>(&kernel_object))
309 : {
310 61 : auto complex_kernel = std::dynamic_pointer_cast<MFEMComplexKernel>(kernel);
311 : auto eqsys =
312 61 : std::dynamic_pointer_cast<Moose::MFEM::ComplexEquationSystem>(getProblemData().eqn_system);
313 61 : if (eqsys)
314 61 : eqsys->AddComplexKernel(std::move(complex_kernel));
315 : else
316 0 : mooseError("Cannot add complex kernel with name '" + name +
317 : "' because there is no corresponding equation system.");
318 61 : }
319 : else
320 : {
321 : auto eqsys =
322 1636 : std::dynamic_pointer_cast<Moose::MFEM::EquationSystem>(getProblemData().eqn_system);
323 1636 : if (eqsys)
324 1636 : eqsys->AddKernel(std::move(kernel));
325 : else
326 0 : mooseError("Cannot add kernel with name '" + name +
327 : "' because there is no corresponding equation system.");
328 1636 : }
329 1697 : }
330 :
331 : void
332 61 : MFEMProblem::addRealComponentToKernel(const std::string & kernel_name,
333 : const std::string & name,
334 : InputParameters & parameters)
335 : {
336 : auto parent_ptr = std::dynamic_pointer_cast<MFEMComplexKernel>(
337 61 : getMFEMObject<MFEMComplexKernel>("Kernel", name).getSharedPtr());
338 244 : parameters.set<VariableName>("variable") = parent_ptr->getParam<VariableName>("variable");
339 122 : auto kernel_ptr = addObject<MFEMKernel>(kernel_name, name + "_real", parameters).front();
340 61 : parent_ptr->setRealKernel(kernel_ptr);
341 61 : }
342 :
343 : void
344 55 : MFEMProblem::addImagComponentToKernel(const std::string & kernel_name,
345 : const std::string & name,
346 : InputParameters & parameters)
347 : {
348 : auto parent_ptr = std::dynamic_pointer_cast<MFEMComplexKernel>(
349 55 : getMFEMObject<MFEMComplexKernel>("Kernel", name).getSharedPtr());
350 220 : parameters.set<VariableName>("variable") = parent_ptr->getParam<VariableName>("variable");
351 110 : auto kernel_ptr = addObject<MFEMKernel>(kernel_name, name + "_imag", parameters).front();
352 55 : parent_ptr->setImagKernel(kernel_ptr);
353 55 : }
354 :
355 : void
356 0 : MFEMProblem::addRealComponentToBC(const std::string & kernel_name,
357 : const std::string & name,
358 : InputParameters & parameters)
359 : {
360 : auto parent_ptr = std::dynamic_pointer_cast<MFEMComplexIntegratedBC>(
361 0 : getMFEMObject<MFEMComplexIntegratedBC>("BoundaryCondition", name).getSharedPtr());
362 0 : parameters.set<VariableName>("variable") = parent_ptr->getParam<VariableName>("variable");
363 0 : parameters.set<std::vector<BoundaryName>>("boundary") =
364 0 : parent_ptr->getParam<std::vector<BoundaryName>>("boundary");
365 : auto bc_ptr = std::dynamic_pointer_cast<MFEMIntegratedBC>(
366 0 : addObject<MFEMBoundaryCondition>(kernel_name, name + "_real", parameters).front());
367 0 : parent_ptr->setRealBC(bc_ptr);
368 0 : }
369 :
370 : void
371 0 : MFEMProblem::addImagComponentToBC(const std::string & kernel_name,
372 : const std::string & name,
373 : InputParameters & parameters)
374 : {
375 : auto parent_ptr = std::dynamic_pointer_cast<MFEMComplexIntegratedBC>(
376 0 : getMFEMObject<MFEMComplexIntegratedBC>("BoundaryCondition", name).getSharedPtr());
377 0 : parameters.set<VariableName>("variable") = parent_ptr->getParam<VariableName>("variable");
378 0 : parameters.set<std::vector<BoundaryName>>("boundary") =
379 0 : parent_ptr->getParam<std::vector<BoundaryName>>("boundary");
380 : auto bc_ptr = std::dynamic_pointer_cast<MFEMIntegratedBC>(
381 0 : addObject<MFEMBoundaryCondition>(kernel_name, name + "_imag", parameters).front());
382 0 : parent_ptr->setImagBC(bc_ptr);
383 0 : }
384 :
385 : int
386 314 : vectorFunctionDim(const std::string & type, const InputParameters & parameters)
387 : {
388 628 : if (parameters.isParamSetByUser("expression_z"))
389 284 : return 3;
390 90 : if (parameters.isParamSetByUser("expression_y") || type == "LevelSetOlssonVortex")
391 28 : return 2;
392 4 : if (parameters.isParamSetByUser("expression_x"))
393 2 : return 1;
394 :
395 0 : return 3;
396 : }
397 :
398 : const std::vector<std::string> SCALAR_FUNCS = {"Axisymmetric2D3DSolutionFunction",
399 : "BicubicSplineFunction",
400 : "CoarsenedPiecewiseLinear",
401 : "CompositeFunction",
402 : "ConstantFunction",
403 : "ImageFunction",
404 : "ParsedFunction",
405 : "ParsedGradFunction",
406 : "PeriodicFunction",
407 : "PiecewiseBilinear",
408 : "PiecewiseConstant",
409 : "PiecewiseConstantFromCSV",
410 : "PiecewiseLinear",
411 : "PiecewiseLinearFromVectorPostprocessor",
412 : "PiecewiseMultiInterpolation",
413 : "PiecewiseMulticonstant",
414 : "SolutionFunction",
415 : "SplineFunction",
416 : "FunctionSeries",
417 : "LevelSetOlssonBubble",
418 : "LevelSetOlssonPlane",
419 : "NearestReporterCoordinatesFunction",
420 : "ParameterMeshFunction",
421 : "ParsedOptimizationFunction",
422 : "FourierNoise",
423 : "MovingPlanarFront",
424 : "MultiControlDrumFunction",
425 : "Grad2ParsedFunction",
426 : "GradParsedFunction",
427 : "ScaledAbsDifferenceDRLRewardFunction",
428 : "CircularAreaHydraulicDiameterFunction",
429 : "CosineHumpFunction",
430 : "CosineTransitionFunction",
431 : "CubicTransitionFunction",
432 : "GeneralizedCircumference",
433 : "PiecewiseFunction",
434 : "TimeRampFunction"},
435 : VECTOR_FUNCS = {"ParsedVectorFunction", "LevelSetOlssonVortex"};
436 :
437 : void
438 1319 : MFEMProblem::addFunction(const std::string & type,
439 : const std::string & name,
440 : InputParameters & parameters)
441 : {
442 1319 : ExternalProblem::addFunction(type, name, parameters);
443 1319 : auto & func = getFunction(name);
444 : // FIXME: Do we want to have optimised versions for when functions
445 : // are only of space or only of time.
446 1319 : if (std::find(SCALAR_FUNCS.begin(), SCALAR_FUNCS.end(), type) != SCALAR_FUNCS.end())
447 : {
448 935 : getCoefficients().declareScalar<mfem::FunctionCoefficient>(
449 : name,
450 935 : [&func](const mfem::Vector & p, mfem::real_t t) -> mfem::real_t
451 1851354 : { return func.value(t, Moose::MFEM::libMeshPointFromMFEMVector(p)); });
452 : }
453 384 : else if (std::find(VECTOR_FUNCS.begin(), VECTOR_FUNCS.end(), type) != VECTOR_FUNCS.end())
454 : {
455 314 : int dim = vectorFunctionDim(type, parameters);
456 314 : getCoefficients().declareVector<mfem::VectorFunctionCoefficient>(
457 : name,
458 : dim,
459 314 : [&func, dim](const mfem::Vector & p, mfem::real_t t, mfem::Vector & u)
460 : {
461 : libMesh::RealVectorValue vector_value =
462 1033366 : func.vectorValue(t, Moose::MFEM::libMeshPointFromMFEMVector(p));
463 3815218 : for (int i = 0; i < dim; i++)
464 : {
465 2781852 : u[i] = vector_value(i);
466 : }
467 1033366 : });
468 : }
469 70 : else if ("MFEMParsedFunction" != type)
470 : {
471 2 : mooseWarning("Could not identify whether function ",
472 : type,
473 : " is scalar or vector; no MFEM coefficient object created.");
474 : }
475 1317 : }
476 :
477 : void
478 575 : MFEMProblem::addPostprocessor(const std::string & type,
479 : const std::string & name,
480 : InputParameters & parameters)
481 : {
482 575 : if (parameters.getSystemAttributeName() == "MFEMExecutedObject")
483 : {
484 1102 : checkUserObjectNameCollision(name, "Postprocessor");
485 1102 : addObject<MFEMExecutedObject>(type, name, parameters);
486 551 : const PostprocessorValue & val = getPostprocessorValueByName(name);
487 551 : getCoefficients().declareScalar<mfem::FunctionCoefficient>(
488 553 : name, [&val](const mfem::Vector &) -> mfem::real_t { return val; });
489 : }
490 : else
491 24 : ExternalProblem::addPostprocessor(type, name, parameters);
492 575 : }
493 :
494 : void
495 306 : MFEMProblem::addVectorPostprocessor(const std::string & type,
496 : const std::string & name,
497 : InputParameters & parameters)
498 : {
499 306 : if (parameters.getSystemAttributeName() == "MFEMExecutedObject")
500 : {
501 612 : checkUserObjectNameCollision(name, "VectorPostprocessor");
502 918 : addObject<MFEMExecutedObject>(type, name, parameters);
503 : }
504 : else
505 0 : ExternalProblem::addVectorPostprocessor(type, name, parameters);
506 306 : }
507 :
508 : InputParameters
509 36 : MFEMProblem::addMFEMFESpaceFromMOOSEVariable(InputParameters & parameters)
510 : {
511 :
512 72 : InputParameters fespace_params = _factory.getValidParams("MFEMGenericFESpace");
513 72 : InputParameters mfem_variable_params = _factory.getValidParams("MFEMVariable");
514 :
515 : auto moose_fe_type =
516 72 : FEType(Utility::string_to_enum<Order>(parameters.get<MooseEnum>("order")),
517 108 : Utility::string_to_enum<FEFamily>(parameters.get<MooseEnum>("family")));
518 :
519 36 : std::string mfem_family;
520 36 : int mfem_vdim = 1;
521 :
522 36 : switch (moose_fe_type.family)
523 : {
524 12 : case FEFamily::LAGRANGE:
525 12 : mfem_family = "H1";
526 12 : mfem_vdim = 1;
527 12 : break;
528 12 : case FEFamily::LAGRANGE_VEC:
529 12 : mfem_family = "H1";
530 12 : mfem_vdim = 3;
531 12 : break;
532 6 : case FEFamily::MONOMIAL:
533 6 : mfem_family = "L2";
534 6 : mfem_vdim = 1;
535 6 : break;
536 6 : case FEFamily::MONOMIAL_VEC:
537 6 : mfem_family = "L2";
538 6 : mfem_vdim = 3;
539 6 : break;
540 0 : default:
541 0 : mooseError("Unable to set MFEM FESpace for MOOSE variable");
542 : break;
543 : }
544 :
545 : // Create fespace name. If this already exists, we will reuse this for
546 : // the mfem variable ("gridfunction").
547 72 : const std::string fespace_name = mfem_family + "_" +
548 144 : std::to_string(mesh().getMFEMParMesh().Dimension()) + "D_P" +
549 144 : std::to_string(moose_fe_type.order.get_order());
550 :
551 : // Set all fespace parameters.
552 36 : fespace_params.set<std::string>("fec_name") = fespace_name;
553 108 : fespace_params.set<int>("vdim") = mfem_vdim;
554 :
555 72 : if (!hasMFEMObject("MFEMFESpace", fespace_name)) // Create the fespace (implicit).
556 : {
557 24 : addFESpace("MFEMGenericFESpace", fespace_name, fespace_params);
558 : }
559 :
560 108 : mfem_variable_params.set<MFEMFESpaceName>("fespace") = fespace_name;
561 :
562 72 : return mfem_variable_params;
563 36 : }
564 :
565 : void
566 2371 : MFEMProblem::displaceMesh()
567 : {
568 : // Displace mesh
569 2371 : if (mesh().shouldDisplace())
570 : {
571 10 : mesh().displace(static_cast<mfem::GridFunction const &>(*getMeshDisplacementGridFunction()));
572 : // TODO: update FESpaces GridFunctions etc for transient solves
573 : }
574 2371 : }
575 :
576 : std::optional<std::reference_wrapper<mfem::ParGridFunction const>>
577 1410 : MFEMProblem::getMeshDisplacementGridFunction()
578 : {
579 : // If C++23 transform were available this would be easier
580 1410 : auto const displacement_variable = mesh().getMeshDisplacementVariable();
581 1410 : if (displacement_variable)
582 : {
583 20 : return *_problem_data.gridfunctions.Get(displacement_variable.value());
584 : }
585 : else
586 : {
587 1390 : return std::nullopt;
588 : }
589 : }
590 :
591 : void
592 0 : MFEMProblem::rebalanceMesh(mfem::ParMesh & pmesh)
593 : {
594 0 : if (pmesh.Nonconforming())
595 : {
596 0 : pmesh.Rebalance();
597 0 : updateFESpaces();
598 0 : updateGridFunctions();
599 : }
600 0 : }
601 :
602 : void
603 13 : MFEMProblem::updateFESpaces()
604 : {
605 39 : for (const auto & fe_space_pair : _problem_data.fespaces)
606 26 : fe_space_pair.second->Update();
607 13 : }
608 :
609 : void
610 26 : MFEMProblem::updateGridFunctions()
611 : {
612 78 : for (const auto & gridfunction_pair : _problem_data.gridfunctions)
613 52 : gridfunction_pair.second->Update();
614 26 : }
615 :
616 : std::vector<VariableName>
617 0 : MFEMProblem::getAuxVariableNames()
618 : {
619 0 : return systemBaseAuxiliary().getVariableNames();
620 : }
621 :
622 : MFEMMesh &
623 33214 : MFEMProblem::mesh()
624 : {
625 : mooseAssert(ExternalProblem::mesh().type() == "MFEMMesh",
626 : "Please choose the MFEMMesh mesh type for an MFEMProblem\n");
627 33214 : return static_cast<MFEMMesh &>(_mesh);
628 : }
629 :
630 : const MFEMMesh &
631 2136 : MFEMProblem::mesh() const
632 : {
633 2136 : return const_cast<MFEMProblem *>(this)->mesh();
634 : }
635 :
636 : void
637 167 : MFEMProblem::addSubMesh(const std::string & var_type,
638 : const std::string & var_name,
639 : InputParameters & parameters)
640 : {
641 334 : auto & mfem_submesh = *addObject<MFEMSubMesh>(var_type, var_name, parameters).front();
642 : // Register submesh.
643 167 : getProblemData().submeshes.Register(var_name, mfem_submesh.getSubMesh());
644 167 : }
645 :
646 : void
647 635 : MFEMProblem::addTransfer(const std::string & transfer_name,
648 : const std::string & name,
649 : InputParameters & parameters)
650 : {
651 635 : if (parameters.getBase() == "MFEMSubMeshTransfer")
652 624 : addObject<MFEMExecutedObject>(transfer_name, name, parameters);
653 : else
654 427 : ExternalProblem::addTransfer(transfer_name, name, parameters);
655 635 : }
656 :
657 : void
658 1079 : MFEMProblem::addInitialCondition(const std::string & ic_name,
659 : const std::string & name,
660 : InputParameters & parameters)
661 : {
662 2158 : addObject<MFEMExecutedObject>(ic_name, name, parameters);
663 1079 : }
664 :
665 : void
666 13852 : MFEMProblem::executeMFEMObjects(const ExecFlagType & exec_type)
667 : {
668 13852 : std::vector<MFEMExecutedObject *> objects;
669 13852 : theWarehouse()
670 13852 : .query()
671 13852 : .condition<AttribSystem>("MFEMExecutedObject")
672 13852 : .condition<AttribExecOns>(exec_type)
673 27704 : .condition<AttribThread>(0)
674 13852 : .queryInto(objects);
675 :
676 13852 : std::map<std::string, const MFEMExecutedObject *> suppliers;
677 16910 : for (auto * const object : objects)
678 6118 : for (const auto & item : object->getSuppliedItems())
679 : {
680 3060 : const auto [it, inserted] = suppliers.emplace(item, object);
681 3060 : if (!inserted && it->second != object)
682 2 : mooseError("MFEM executed-object dependency ambiguity on ",
683 : exec_type,
684 : ": both '",
685 2 : it->second->name(),
686 : "' and '",
687 2 : object->name(),
688 : "' supply '",
689 : item,
690 : "'.");
691 : }
692 :
693 16906 : for (auto * const object : objects)
694 : {
695 3056 : object->initialize();
696 3056 : object->execute();
697 3056 : object->finalize();
698 :
699 3056 : if (auto * const pp = dynamic_cast<const Postprocessor *>(object))
700 : {
701 591 : _reporter_data.finalize(pp->PPName());
702 591 : setPostprocessorValueByName(pp->PPName(), pp->getValue());
703 : }
704 :
705 3056 : if (auto * const vpp = dynamic_cast<VectorPostprocessor *>(object))
706 624 : _reporter_data.finalize(vpp->PPName());
707 : }
708 13854 : }
709 :
710 : std::string
711 1384 : MFEMProblem::solverTypeString(const unsigned int libmesh_dbg_var(solver_sys_num))
712 : {
713 : mooseAssert(solver_sys_num == 0, "No support for multi-system with MFEM right now");
714 1384 : return MooseUtils::prettyCppType(getProblemData().jacobian_solver.get());
715 : }
716 :
717 : bool
718 36 : MFEMProblem::hasMFEMObject(const std::string & system, const std::string & name) const
719 : {
720 36 : std::vector<MooseObject *> objs;
721 36 : theWarehouse()
722 36 : .query()
723 36 : .condition<AttribSystem>(system)
724 72 : .condition<AttribThread>(0)
725 36 : .condition<AttribName>(name)
726 36 : .queryInto(objs);
727 72 : return !objs.empty();
728 36 : }
729 :
730 : #endif
|