https://mooseframework.inl.gov
Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Moose::MFEM::TimeDependentEquationSystem Class Reference

Class to store weak form components for time dependent PDEs. More...

#include <EquationSystem.h>

Inheritance diagram for Moose::MFEM::TimeDependentEquationSystem:
[legend]

Public Member Functions

 TimeDependentEquationSystem ()
 
void AddTrialVariableNameIfMissing (const std::string &trial_var_name) override
 Add trial variable to EquationSystem. More...
 
virtual void SetTimeStep (double dt)
 
virtual void UpdateEquationSystem ()
 
virtual void AddKernel (std::shared_ptr< MFEMKernel > kernel) override
 Add kernels. More...
 
virtual void BuildBilinearForms () override
 
virtual void FormLegacySystem (mfem::OperatorHandle &op, mfem::BlockVector &truedXdt, mfem::BlockVector &trueRHS) override
 
virtual void FormSystem (mfem::OperatorHandle &op, mfem::BlockVector &truedXdt, mfem::BlockVector &trueRHS) override
 
const std::vector< std::string > & TrialVarTimeDerivativeNames () const
 
virtual void AddTestVariableNameIfMissing (const std::string &test_var_name)
 Add test variable to EquationSystem. More...
 
virtual void AddIntegratedBC (std::shared_ptr< MFEMIntegratedBC > kernel)
 
virtual void AddEssentialBC (std::shared_ptr< MFEMEssentialBC > bc)
 
virtual void ApplyEssentialBCs ()
 
virtual void Init (Moose::MFEM::GridFunctions &gridfunctions, const Moose::MFEM::FESpaces &fespaces, mfem::AssemblyLevel assembly_level)
 Build forms. More...
 
virtual void BuildLinearForms ()
 
virtual void BuildMixedBilinearForms ()
 
virtual void BuildEquationSystem ()
 
virtual void FormLinearSystem (mfem::OperatorHandle &op, mfem::BlockVector &trueX, mfem::BlockVector &trueRHS)
 Form linear system, with essential boundary conditions accounted for. More...
 
virtual void BuildJacobian (mfem::BlockVector &trueX, mfem::BlockVector &trueRHS)
 Build linear system, with essential boundary conditions accounted for. More...
 
void Mult (const mfem::Vector &u, mfem::Vector &residual) const override
 Compute residual y = Mu. More...
 
mfem::Operator & GetGradient (const mfem::Vector &u) const override
 Compute J = M + grad_H(u) More...
 
virtual void RecoverFEMSolution (mfem::BlockVector &trueX, Moose::MFEM::GridFunctions &gridfunctions)
 Update variable from solution vector after solve. More...
 
const std::vector< std::string > & TrialVarNames () const
 
const std::vector< std::string > & TestVarNames () const
 

Public Attributes

std::vector< mfem::Array< int > > _ess_tdof_lists
 

Protected Member Functions

void DeleteAllBlocks ()
 Deletes the HypreParMatrix associated with any pointer stored in _h_blocks, and then proceeds to delete all dynamically allocated memory for _h_blocks itself, resetting all dimensions to zero. More...
 
bool VectorContainsName (const std::vector< std::string > &the_vector, const std::string &name) const
 
template<class FormType >
void ApplyDomainBLFIntegrators (const std::string &trial_var_name, const std::string &test_var_name, std::shared_ptr< FormType > form, Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel >>>> &kernels_map)
 Template method for applying BilinearFormIntegrators on domains from kernels to a BilinearForm, or MixedBilinearForm. More...
 
void ApplyDomainLFIntegrators (const std::string &test_var_name, std::shared_ptr< mfem::ParLinearForm > form, Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel >>>> &kernels_map)
 
template<class FormType >
void ApplyBoundaryBLFIntegrators (const std::string &trial_var_name, const std::string &test_var_name, std::shared_ptr< FormType > form, Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC >>>> &integrated_bc_map)
 
void ApplyBoundaryLFIntegrators (const std::string &test_var_name, std::shared_ptr< mfem::ParLinearForm > form, Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC >>>> &integrated_bc_map)
 

Protected Attributes

mfem::ConstantCoefficient _dt_coef
 Coefficient for timestep scaling. More...
 
std::vector< std::string > _trial_var_time_derivative_names
 
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel > > > > _td_kernels_map
 
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _td_blfs
 Container to store contributions to weak form of the form (F du/dt, v) More...
 
std::vector< std::string > _trial_var_names
 Names of all variables corresponding to gridfunctions. More...
 
Moose::MFEM::GridFunctions _trial_variables
 Pointers to trial variables. More...
 
std::vector< std::string > _test_var_names
 Names of all test variables corresponding to linear forms in this equation system. More...
 
std::vector< mfem::ParFiniteElementSpace * > _test_pfespaces
 Pointers to finite element spaces associated with test variables. More...
 
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _blfs
 
Moose::MFEM::NamedFieldsMap< mfem::ParLinearForm > _lfs
 
Moose::MFEM::NamedFieldsMap< mfem::ParNonlinearForm > _nlfs
 
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< mfem::ParMixedBilinearForm > > _mblfs
 
std::vector< std::unique_ptr< mfem::ParGridFunction > > _xs
 Gridfunctions for setting Dirichlet BCs. More...
 
std::vector< std::unique_ptr< mfem::ParGridFunction > > _dxdts
 
mfem::Array2D< const mfem::HypreParMatrix * > _h_blocks
 
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel > > > > _kernels_map
 Arrays to store kernels to act on each component of weak form. More...
 
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC > > > > _integrated_bc_map
 Arrays to store integrated BCs to act on each component of weak form. More...
 
Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMEssentialBC > > > _essential_bc_map
 Arrays to store essential BCs to act on each component of weak form. More...
 
mfem::OperatorHandle _jacobian
 
mfem::AssemblyLevel _assembly_level
 

Detailed Description

Class to store weak form components for time dependent PDEs.

Definition at line 283 of file EquationSystem.h.

Constructor & Destructor Documentation

◆ TimeDependentEquationSystem()

Moose::MFEM::TimeDependentEquationSystem::TimeDependentEquationSystem ( )

Definition at line 403 of file EquationSystem.C.

403 : _dt_coef(1.0) {}
mfem::ConstantCoefficient _dt_coef
Coefficient for timestep scaling.

Member Function Documentation

◆ AddEssentialBC()

void Moose::MFEM::EquationSystem::AddEssentialBC ( std::shared_ptr< MFEMEssentialBC bc)
virtualinherited

Definition at line 104 of file EquationSystem.C.

105 {
106  AddTestVariableNameIfMissing(bc->getTestVariableName());
107  auto test_var_name = bc->getTestVariableName();
108  if (!_essential_bc_map.Has(test_var_name))
109  {
110  auto bcs = std::make_shared<std::vector<std::shared_ptr<MFEMEssentialBC>>>();
111  _essential_bc_map.Register(test_var_name, std::move(bcs));
112  }
113  _essential_bc_map.GetRef(test_var_name).push_back(std::move(bc));
114 }
virtual void AddTestVariableNameIfMissing(const std::string &test_var_name)
Add test variable to EquationSystem.
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMEssentialBC > > > _essential_bc_map
Arrays to store essential BCs to act on each component of weak form.
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ AddIntegratedBC()

void Moose::MFEM::EquationSystem::AddIntegratedBC ( std::shared_ptr< MFEMIntegratedBC kernel)
virtualinherited

Definition at line 81 of file EquationSystem.C.

82 {
83  AddTestVariableNameIfMissing(bc->getTestVariableName());
84  AddTrialVariableNameIfMissing(bc->getTrialVariableName());
85  auto trial_var_name = bc->getTrialVariableName();
86  auto test_var_name = bc->getTestVariableName();
87  if (!_integrated_bc_map.Has(test_var_name))
88  {
89  auto integrated_bc_field_map = std::make_shared<
91  _integrated_bc_map.Register(test_var_name, std::move(integrated_bc_field_map));
92  }
93  // Register new integrated bc map if not present for the test/trial variable
94  // pair
95  if (!_integrated_bc_map.Get(test_var_name)->Has(trial_var_name))
96  {
97  auto bcs = std::make_shared<std::vector<std::shared_ptr<MFEMIntegratedBC>>>();
98  _integrated_bc_map.Get(test_var_name)->Register(trial_var_name, std::move(bcs));
99  }
100  _integrated_bc_map.GetRef(test_var_name).Get(trial_var_name)->push_back(std::move(bc));
101 }
virtual void AddTestVariableNameIfMissing(const std::string &test_var_name)
Add test variable to EquationSystem.
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
Lightweight adaptor over an std::map from strings to pointer to T.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
virtual void AddTrialVariableNameIfMissing(const std::string &trial_var_name)
Add trial variable to EquationSystem.
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC > > > > _integrated_bc_map
Arrays to store integrated BCs to act on each component of weak form.
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ AddKernel()

void Moose::MFEM::TimeDependentEquationSystem::AddKernel ( std::shared_ptr< MFEMKernel kernel)
overridevirtual

Add kernels.

Reimplemented from Moose::MFEM::EquationSystem.

Definition at line 435 of file EquationSystem.C.

436 {
437  if (kernel->getTrialVariableName() == GetTimeDerivativeName(kernel->getTestVariableName()))
438  {
439  auto trial_var_name = kernel->getTrialVariableName();
440  auto test_var_name = kernel->getTestVariableName();
441  AddTestVariableNameIfMissing(test_var_name);
442  AddTrialVariableNameIfMissing(test_var_name);
443  if (!_td_kernels_map.Has(test_var_name))
444  {
445  auto kernel_field_map =
446  std::make_shared<Moose::MFEM::NamedFieldsMap<std::vector<std::shared_ptr<MFEMKernel>>>>();
447  _td_kernels_map.Register(test_var_name, std::move(kernel_field_map));
448  }
449  // Register new kernels map if not present for the test variable
450  if (!_td_kernels_map.Get(test_var_name)->Has(test_var_name))
451  {
452  auto kernels = std::make_shared<std::vector<std::shared_ptr<MFEMKernel>>>();
453  _td_kernels_map.Get(test_var_name)->Register(test_var_name, std::move(kernels));
454  }
455  _td_kernels_map.GetRef(test_var_name).Get(test_var_name)->push_back(std::move(kernel));
456  }
457  else
458  {
460  }
461 }
virtual void AddTestVariableNameIfMissing(const std::string &test_var_name)
Add test variable to EquationSystem.
virtual void AddKernel(std::shared_ptr< MFEMKernel > kernel)
Add kernels.
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel > > > > _td_kernels_map
void AddTrialVariableNameIfMissing(const std::string &trial_var_name) override
Add trial variable to EquationSystem.
std::string GetTimeDerivativeName(std::string name)
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ AddTestVariableNameIfMissing()

void Moose::MFEM::EquationSystem::AddTestVariableNameIfMissing ( const std::string &  test_var_name)
virtualinherited

Add test variable to EquationSystem.

Definition at line 49 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::AddEssentialBC(), Moose::MFEM::EquationSystem::AddIntegratedBC(), Moose::MFEM::EquationSystem::AddKernel(), and AddKernel().

50 {
51  if (!VectorContainsName(_test_var_names, test_var_name))
52  {
53  _test_var_names.push_back(test_var_name);
54  }
55 }
bool VectorContainsName(const std::vector< std::string > &the_vector, const std::string &name) const
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.

◆ AddTrialVariableNameIfMissing()

void Moose::MFEM::TimeDependentEquationSystem::AddTrialVariableNameIfMissing ( const std::string &  trial_var_name)
overridevirtual

Add trial variable to EquationSystem.

Reimplemented from Moose::MFEM::EquationSystem.

Definition at line 406 of file EquationSystem.C.

Referenced by AddKernel().

407 {
408  // The TimeDependentEquationSystem operator expects to act on a vector of variable time
409  // derivatives
410  std::string var_time_derivative_name = GetTimeDerivativeName(var_name);
411 
412  if (!VectorContainsName(_trial_var_names, var_time_derivative_name))
413  {
414  _trial_var_names.push_back(var_time_derivative_name);
415  _trial_var_time_derivative_names.push_back(var_time_derivative_name);
416  }
417 }
bool VectorContainsName(const std::vector< std::string > &the_vector, const std::string &name) const
std::string GetTimeDerivativeName(std::string name)
std::vector< std::string > _trial_var_names
Names of all variables corresponding to gridfunctions.
std::vector< std::string > _trial_var_time_derivative_names

◆ ApplyBoundaryBLFIntegrators()

template<class FormType >
void Moose::MFEM::EquationSystem::ApplyBoundaryBLFIntegrators ( const std::string &  trial_var_name,
const std::string &  test_var_name,
std::shared_ptr< FormType >  form,
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC >>>> &  integrated_bc_map 
)
protectedinherited

Definition at line 230 of file EquationSystem.h.

237 {
238  if (integrated_bc_map.Has(test_var_name) &&
239  integrated_bc_map.Get(test_var_name)->Has(trial_var_name))
240  {
241  auto bcs = integrated_bc_map.GetRef(test_var_name).GetRef(trial_var_name);
242  for (auto & bc : bcs)
243  {
244  mfem::BilinearFormIntegrator * integ = bc->createBFIntegrator();
245  if (integ != nullptr)
246  {
247  bc->isBoundaryRestricted()
248  ? form->AddBoundaryIntegrator(std::move(integ), bc->getBoundaryMarkers())
249  : form->AddBoundaryIntegrator(std::move(integ));
250  }
251  }
252  }
253 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ ApplyBoundaryLFIntegrators()

void Moose::MFEM::EquationSystem::ApplyBoundaryLFIntegrators ( const std::string &  test_var_name,
std::shared_ptr< mfem::ParLinearForm >  form,
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC >>>> &  integrated_bc_map 
)
inlineprotectedinherited

Definition at line 256 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystem::BuildLinearForms().

262 {
263  if (integrated_bc_map.Has(test_var_name) &&
264  integrated_bc_map.Get(test_var_name)->Has(test_var_name))
265  {
266  auto bcs = integrated_bc_map.GetRef(test_var_name).GetRef(test_var_name);
267  for (auto & bc : bcs)
268  {
269  mfem::LinearFormIntegrator * integ = bc->createLFIntegrator();
270  if (integ != nullptr)
271  {
272  bc->isBoundaryRestricted()
273  ? form->AddBoundaryIntegrator(std::move(integ), bc->getBoundaryMarkers())
274  : form->AddBoundaryIntegrator(std::move(integ));
275  }
276  }
277  }
278 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ ApplyDomainBLFIntegrators()

template<class FormType >
void Moose::MFEM::EquationSystem::ApplyDomainBLFIntegrators ( const std::string &  trial_var_name,
const std::string &  test_var_name,
std::shared_ptr< FormType >  form,
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel >>>> &  kernels_map 
)
protectedinherited

Template method for applying BilinearFormIntegrators on domains from kernels to a BilinearForm, or MixedBilinearForm.

Definition at line 182 of file EquationSystem.h.

188 {
189  if (kernels_map.Has(test_var_name) && kernels_map.Get(test_var_name)->Has(trial_var_name))
190  {
191  auto kernels = kernels_map.GetRef(test_var_name).GetRef(trial_var_name);
192  for (auto & kernel : kernels)
193  {
194  mfem::BilinearFormIntegrator * integ = kernel->createBFIntegrator();
195  if (integ != nullptr)
196  {
197  kernel->isSubdomainRestricted()
198  ? form->AddDomainIntegrator(std::move(integ), kernel->getSubdomainMarkers())
199  : form->AddDomainIntegrator(std::move(integ));
200  }
201  }
202  }
203 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ ApplyDomainLFIntegrators()

void Moose::MFEM::EquationSystem::ApplyDomainLFIntegrators ( const std::string &  test_var_name,
std::shared_ptr< mfem::ParLinearForm >  form,
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel >>>> &  kernels_map 
)
inlineprotectedinherited

Definition at line 206 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystem::BuildLinearForms().

211 {
212  if (kernels_map.Has(test_var_name) && kernels_map.Get(test_var_name)->Has(test_var_name))
213  {
214  auto kernels = kernels_map.GetRef(test_var_name).GetRef(test_var_name);
215  for (auto & kernel : kernels)
216  {
217  mfem::LinearFormIntegrator * integ = kernel->createLFIntegrator();
218  if (integ != nullptr)
219  {
220  kernel->isSubdomainRestricted()
221  ? form->AddDomainIntegrator(std::move(integ), kernel->getSubdomainMarkers())
222  : form->AddDomainIntegrator(std::move(integ));
223  }
224  }
225  }
226 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
T & GetRef(const std::string &field_name) const
Returns a reference to a field.

◆ ApplyEssentialBCs()

void Moose::MFEM::EquationSystem::ApplyEssentialBCs ( )
virtualinherited

Definition at line 117 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::BuildLinearForms().

118 {
119  _ess_tdof_lists.resize(_test_var_names.size());
120  for (const auto i : index_range(_test_var_names))
121  {
122  auto test_var_name = _test_var_names.at(i);
123  if (!_essential_bc_map.Has(test_var_name))
124  continue;
125 
126  // Set default value of gridfunction used in essential BC. Values
127  // overwritten in applyEssentialBCs
128  mfem::ParGridFunction & trial_gf(*(_xs.at(i)));
129  auto * const pmesh = _test_pfespaces.at(i)->GetParMesh();
130  mooseAssert(pmesh, "parallel mesh is null");
131 
132  auto bcs = _essential_bc_map.GetRef(test_var_name);
133  mfem::Array<int> global_ess_markers(pmesh->bdr_attributes.Max());
134  global_ess_markers = 0;
135  for (auto & bc : bcs)
136  {
137  bc->ApplyBC(trial_gf);
138 
139  mfem::Array<int> ess_bdrs(bc->getBoundaryMarkers());
140  for (auto it = 0; it != pmesh->bdr_attributes.Max(); ++it)
141  {
142  global_ess_markers[it] = std::max(global_ess_markers[it], ess_bdrs[it]);
143  }
144  }
145  trial_gf.FESpace()->GetEssentialTrueDofs(global_ess_markers, _ess_tdof_lists.at(i));
146  }
147 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
std::vector< std::unique_ptr< mfem::ParGridFunction > > _xs
Gridfunctions for setting Dirichlet BCs.
std::vector< mfem::Array< int > > _ess_tdof_lists
auto max(const L &left, const R &right)
Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMEssentialBC > > > _essential_bc_map
Arrays to store essential BCs to act on each component of weak form.
std::vector< mfem::ParFiniteElementSpace * > _test_pfespaces
Pointers to finite element spaces associated with test variables.
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
T & GetRef(const std::string &field_name) const
Returns a reference to a field.
auto index_range(const T &sizable)

◆ BuildBilinearForms()

void Moose::MFEM::TimeDependentEquationSystem::BuildBilinearForms ( )
overridevirtual

Reimplemented from Moose::MFEM::EquationSystem.

Definition at line 464 of file EquationSystem.C.

Referenced by UpdateEquationSystem().

465 {
467 
468  // Build and assemble bilinear forms acting on time derivatives
469  for (const auto i : index_range(_test_var_names))
470  {
471  auto test_var_name = _test_var_names.at(i);
472 
473  _td_blfs.Register(test_var_name,
474  std::make_shared<mfem::ParBilinearForm>(_test_pfespaces.at(i)));
475 
476  // Apply kernels to td_blf
477  auto td_blf = _td_blfs.GetShared(test_var_name);
478  td_blf->SetAssemblyLevel(_assembly_level);
479  ApplyBoundaryBLFIntegrators<mfem::ParBilinearForm>(
480  test_var_name, test_var_name, td_blf, _integrated_bc_map);
481  ApplyDomainBLFIntegrators<mfem::ParBilinearForm>(
482  test_var_name, test_var_name, td_blf, _td_kernels_map);
483 
484  // Recover and scale integrators from blf. This is to apply the dt*du/dt contributions from the
485  // operator on the trial variable in the implicit integration scheme
486  auto blf = _blfs.Get(test_var_name);
487  auto integs = blf->GetDBFI();
488  auto b_integs = blf->GetBBFI();
489  auto markers = blf->GetBBFI_Marker();
490 
491  mfem::SumIntegrator * sum = new mfem::SumIntegrator(false);
492  ScaleIntegrator * scaled_sum = new ScaleIntegrator(sum, _dt_coef.constant, true);
493 
494  for (int i = 0; i < integs->Size(); ++i)
495  {
496  sum->AddIntegrator(*integs[i]);
497  }
498 
499  for (int i = 0; i < b_integs->Size(); ++i)
500  {
501  td_blf->AddBoundaryIntegrator(new ScaleIntegrator(*b_integs[i], _dt_coef.constant, false),
502  *(*markers[i]));
503  }
504 
505  // scaled_sum is owned by td_blf
506  td_blf->AddDomainIntegrator(scaled_sum);
507 
508  // Assemble form
509  td_blf->Assemble();
510  }
511 }
mfem::ConstantCoefficient _dt_coef
Coefficient for timestep scaling.
virtual void BuildBilinearForms()
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel > > > > _td_kernels_map
mfem::AssemblyLevel _assembly_level
std::vector< mfem::ParFiniteElementSpace * > _test_pfespaces
Pointers to finite element spaces associated with test variables.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
std::shared_ptr< T > GetShared(const std::string &field_name) const
Returns a shared pointer to the field. This is guaranteed to return a non-null shared pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC > > > > _integrated_bc_map
Arrays to store integrated BCs to act on each component of weak form.
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _td_blfs
Container to store contributions to weak form of the form (F du/dt, v)
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _blfs
auto index_range(const T &sizable)

◆ BuildEquationSystem()

void Moose::MFEM::EquationSystem::BuildEquationSystem ( )
virtualinherited

◆ BuildJacobian()

void Moose::MFEM::EquationSystem::BuildJacobian ( mfem::BlockVector &  trueX,
mfem::BlockVector &  trueRHS 
)
virtualinherited

Build linear system, with essential boundary conditions accounted for.

Definition at line 254 of file EquationSystem.C.

Referenced by Moose::MFEM::TimeDomainEquationSystemProblemOperator::BuildEquationSystemOperator(), and Moose::MFEM::EquationSystemProblemOperator::Solve().

255 {
256  height = trueX.Size();
257  width = trueRHS.Size();
258  FormLinearSystem(_jacobian, trueX, trueRHS);
259 }
mfem::OperatorHandle _jacobian
virtual void FormLinearSystem(mfem::OperatorHandle &op, mfem::BlockVector &trueX, mfem::BlockVector &trueRHS)
Form linear system, with essential boundary conditions accounted for.

◆ BuildLinearForms()

void Moose::MFEM::EquationSystem::BuildLinearForms ( )
virtualinherited

Definition at line 314 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::BuildEquationSystem(), and UpdateEquationSystem().

315 {
316  // Register linear forms
317  for (const auto i : index_range(_test_var_names))
318  {
319  auto test_var_name = _test_var_names.at(i);
320  _lfs.Register(test_var_name, std::make_shared<mfem::ParLinearForm>(_test_pfespaces.at(i)));
321  _lfs.GetRef(test_var_name) = 0.0;
322  }
323  // Apply boundary conditions
325 
326  for (auto & test_var_name : _test_var_names)
327  {
328  // Apply kernels
329  auto lf = _lfs.GetShared(test_var_name);
330  ApplyDomainLFIntegrators(test_var_name, lf, _kernels_map);
332  lf->Assemble();
333  }
334 }
void ApplyBoundaryLFIntegrators(const std::string &test_var_name, std::shared_ptr< mfem::ParLinearForm > form, Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC >>>> &integrated_bc_map)
std::vector< mfem::ParFiniteElementSpace * > _test_pfespaces
Pointers to finite element spaces associated with test variables.
std::shared_ptr< T > GetShared(const std::string &field_name) const
Returns a shared pointer to the field. This is guaranteed to return a non-null shared pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMIntegratedBC > > > > _integrated_bc_map
Arrays to store integrated BCs to act on each component of weak form.
void ApplyDomainLFIntegrators(const std::string &test_var_name, std::shared_ptr< mfem::ParLinearForm > form, Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel >>>> &kernels_map)
virtual void ApplyEssentialBCs()
Moose::MFEM::NamedFieldsMap< mfem::ParLinearForm > _lfs
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel > > > > _kernels_map
Arrays to store kernels to act on each component of weak form.
T & GetRef(const std::string &field_name) const
Returns a reference to a field.
auto index_range(const T &sizable)

◆ BuildMixedBilinearForms()

void Moose::MFEM::EquationSystem::BuildMixedBilinearForms ( )
virtualinherited

Definition at line 358 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::BuildEquationSystem(), and UpdateEquationSystem().

359 {
360  // Register mixed bilinear forms. Note that not all combinations may
361  // have a kernel
362 
363  // Create mblf for each test/trial pair
364  for (const auto i : index_range(_test_var_names))
365  {
366  auto test_var_name = _test_var_names.at(i);
367  auto test_mblfs = std::make_shared<Moose::MFEM::NamedFieldsMap<mfem::ParMixedBilinearForm>>();
368  for (const auto j : index_range(_test_var_names))
369  {
370  auto trial_var_name = _trial_var_names.at(j);
371  auto mblf = std::make_shared<mfem::ParMixedBilinearForm>(_test_pfespaces.at(j),
372  _test_pfespaces.at(i));
373  // Register MixedBilinearForm if kernels exist for it, and assemble
374  // kernels
375  if (_kernels_map.Has(test_var_name) && _kernels_map.Get(test_var_name)->Has(trial_var_name) &&
376  test_var_name != trial_var_name)
377  {
378  mblf->SetAssemblyLevel(_assembly_level);
379  // Apply all mixed kernels with this test/trial pair
380  ApplyDomainBLFIntegrators<mfem::ParMixedBilinearForm>(
381  trial_var_name, test_var_name, mblf, _kernels_map);
382  // Assemble mixed bilinear forms
383  mblf->Assemble();
384  // Register mixed bilinear forms associated with a single trial variable
385  // for the current test variable
386  test_mblfs->Register(trial_var_name, mblf);
387  }
388  }
389  // Register all mixed bilinear form sets associated with a single test
390  // variable
391  _mblfs.Register(test_var_name, test_mblfs);
392  }
393 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< mfem::ParMixedBilinearForm > > _mblfs
mfem::AssemblyLevel _assembly_level
std::vector< mfem::ParFiniteElementSpace * > _test_pfespaces
Pointers to finite element spaces associated with test variables.
std::vector< std::string > _trial_var_names
Names of all variables corresponding to gridfunctions.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap< std::vector< std::shared_ptr< MFEMKernel > > > > _kernels_map
Arrays to store kernels to act on each component of weak form.
auto index_range(const T &sizable)

◆ DeleteAllBlocks()

void Moose::MFEM::EquationSystem::DeleteAllBlocks ( )
protectedinherited

Deletes the HypreParMatrix associated with any pointer stored in _h_blocks, and then proceeds to delete all dynamically allocated memory for _h_blocks itself, resetting all dimensions to zero.

Definition at line 21 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::FormLegacySystem(), FormLegacySystem(), and Moose::MFEM::EquationSystem::~EquationSystem().

22 {
23  for (const auto i : make_range(_h_blocks.NumRows()))
24  for (const auto j : make_range(_h_blocks.NumCols()))
25  delete _h_blocks(i, j);
26  _h_blocks.DeleteAll();
27 }
mfem::Array2D< const mfem::HypreParMatrix * > _h_blocks
IntRange< T > make_range(T beg, T end)

◆ FormLegacySystem()

void Moose::MFEM::TimeDependentEquationSystem::FormLegacySystem ( mfem::OperatorHandle &  op,
mfem::BlockVector &  truedXdt,
mfem::BlockVector &  trueRHS 
)
overridevirtual

Reimplemented from Moose::MFEM::EquationSystem.

Definition at line 514 of file EquationSystem.C.

517 {
518 
519  // Allocate block operator
520  DeleteAllBlocks();
521  _h_blocks.SetSize(_test_var_names.size(), _test_var_names.size());
522  // Form diagonal blocks.
523  for (const auto i : index_range(_test_var_names))
524  {
525  auto & test_var_name = _test_var_names.at(i);
526  auto td_blf = _td_blfs.Get(test_var_name);
527  auto blf = _blfs.Get(test_var_name);
528  auto lf = _lfs.Get(test_var_name);
529  // if implicit, add contribution to linear form from terms involving state
530  // variable at previous timestep: {
531  blf->AddMult(*_trial_variables.Get(test_var_name), *lf, -1.0);
532  // }
533  mfem::Vector aux_x, aux_rhs;
534  // Update solution values on Dirichlet values to be in terms of du/dt instead of u
535  mfem::Vector bc_x = *(_xs.at(i).get());
536  bc_x -= *_trial_variables.Get(test_var_name);
537  bc_x /= _dt_coef.constant;
538 
539  // Form linear system for operator acting on vector of du/dt
540  mfem::HypreParMatrix * aux_a = new mfem::HypreParMatrix;
541  td_blf->FormLinearSystem(_ess_tdof_lists.at(i), bc_x, *lf, *aux_a, aux_x, aux_rhs);
542  _h_blocks(i, i) = aux_a;
543  truedXdt.GetBlock(i) = aux_x;
544  trueRHS.GetBlock(i) = aux_rhs;
545  }
546 
547  truedXdt.SyncFromBlocks();
548  trueRHS.SyncFromBlocks();
549 
550  // Create monolithic matrix
551  op.Reset(mfem::HypreParMatrixFromBlocks(_h_blocks));
552 }
mfem::ConstantCoefficient _dt_coef
Coefficient for timestep scaling.
std::vector< std::unique_ptr< mfem::ParGridFunction > > _xs
Gridfunctions for setting Dirichlet BCs.
std::vector< mfem::Array< int > > _ess_tdof_lists
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
mfem::Array2D< const mfem::HypreParMatrix * > _h_blocks
void DeleteAllBlocks()
Deletes the HypreParMatrix associated with any pointer stored in _h_blocks, and then proceeds to dele...
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _td_blfs
Container to store contributions to weak form of the form (F du/dt, v)
Moose::MFEM::GridFunctions _trial_variables
Pointers to trial variables.
Moose::MFEM::NamedFieldsMap< mfem::ParLinearForm > _lfs
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _blfs
auto index_range(const T &sizable)

◆ FormLinearSystem()

void Moose::MFEM::EquationSystem::FormLinearSystem ( mfem::OperatorHandle &  op,
mfem::BlockVector &  trueX,
mfem::BlockVector &  trueRHS 
)
virtualinherited

Form linear system, with essential boundary conditions accounted for.

Definition at line 150 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::BuildJacobian().

153 {
154 
155  switch (_assembly_level)
156  {
157  case mfem::AssemblyLevel::LEGACY:
158  FormLegacySystem(op, trueX, trueRHS);
159  break;
160  default:
161  mooseAssert(_test_var_names.size() == 1,
162  "Non-legacy assembly is only supported for single-variable systems");
163  mooseAssert(
164  _test_var_names.size() == _trial_var_names.size(),
165  "Non-legacy assembly is only supported for single test and trial variable systems");
166  FormSystem(op, trueX, trueRHS);
167  }
168 }
mfem::AssemblyLevel _assembly_level
std::vector< std::string > _trial_var_names
Names of all variables corresponding to gridfunctions.
virtual void FormLegacySystem(mfem::OperatorHandle &op, mfem::BlockVector &trueX, mfem::BlockVector &trueRHS)
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
virtual void FormSystem(mfem::OperatorHandle &op, mfem::BlockVector &trueX, mfem::BlockVector &trueRHS)

◆ FormSystem()

void Moose::MFEM::TimeDependentEquationSystem::FormSystem ( mfem::OperatorHandle &  op,
mfem::BlockVector &  truedXdt,
mfem::BlockVector &  trueRHS 
)
overridevirtual

Reimplemented from Moose::MFEM::EquationSystem.

Definition at line 555 of file EquationSystem.C.

558 {
559  auto & test_var_name = _test_var_names.at(0);
560  auto td_blf = _td_blfs.Get(test_var_name);
561  auto blf = _blfs.Get(test_var_name);
562  auto lf = _lfs.Get(test_var_name);
563  // if implicit, add contribution to linear form from terms involving state
564  // variable at previous timestep: {
565 
566  // The AddMult method in mfem::BilinearForm is not defined for non-legacy assembly
567  mfem::Vector lf_prev(lf->Size());
568  blf->Mult(*_trial_variables.Get(test_var_name), lf_prev);
569  *lf -= lf_prev;
570  // }
571  mfem::Vector aux_x, aux_rhs;
572  // Update solution values on Dirichlet values to be in terms of du/dt instead of u
573  mfem::Vector bc_x = *(_xs.at(0).get());
574  bc_x -= *_trial_variables.Get(test_var_name);
575  bc_x /= _dt_coef.constant;
576 
577  // Form linear system for operator acting on vector of du/dt
578  mfem::OperatorPtr aux_a;
579  td_blf->FormLinearSystem(_ess_tdof_lists.at(0), bc_x, *lf, aux_a, aux_x, aux_rhs);
580 
581  truedXdt.GetBlock(0) = aux_x;
582  trueRHS.GetBlock(0) = aux_rhs;
583  truedXdt.SyncFromBlocks();
584  trueRHS.SyncFromBlocks();
585 
586  // Create monolithic matrix
587  op.Reset(aux_a.Ptr());
588  aux_a.SetOperatorOwner(false);
589 }
mfem::ConstantCoefficient _dt_coef
Coefficient for timestep scaling.
std::vector< std::unique_ptr< mfem::ParGridFunction > > _xs
Gridfunctions for setting Dirichlet BCs.
std::vector< mfem::Array< int > > _ess_tdof_lists
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _td_blfs
Container to store contributions to weak form of the form (F du/dt, v)
Moose::MFEM::GridFunctions _trial_variables
Pointers to trial variables.
Moose::MFEM::NamedFieldsMap< mfem::ParLinearForm > _lfs
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _blfs

◆ GetGradient()

mfem::Operator & Moose::MFEM::EquationSystem::GetGradient ( const mfem::Vector &  u) const
overrideinherited

Compute J = M + grad_H(u)

Definition at line 270 of file EquationSystem.C.

271 {
272  return *_jacobian;
273 }
mfem::OperatorHandle _jacobian

◆ Init()

void Moose::MFEM::EquationSystem::Init ( Moose::MFEM::GridFunctions gridfunctions,
const Moose::MFEM::FESpaces fespaces,
mfem::AssemblyLevel  assembly_level 
)
virtualinherited

Build forms.

Definition at line 288 of file EquationSystem.C.

291 {
292  _assembly_level = assembly_level;
293 
294  for (auto & test_var_name : _test_var_names)
295  {
296  if (!gridfunctions.Has(test_var_name))
297  {
298  MFEM_ABORT("Test variable " << test_var_name
299  << " requested by equation system during initialisation was "
300  "not found in gridfunctions");
301  }
302  // Store pointers to variable FESpaces
303  _test_pfespaces.push_back(gridfunctions.Get(test_var_name)->ParFESpace());
304  // Create auxiliary gridfunctions for applying Dirichlet conditions
305  _xs.emplace_back(
306  std::make_unique<mfem::ParGridFunction>(gridfunctions.Get(test_var_name)->ParFESpace()));
307  _dxdts.emplace_back(
308  std::make_unique<mfem::ParGridFunction>(gridfunctions.Get(test_var_name)->ParFESpace()));
309  _trial_variables.Register(test_var_name, gridfunctions.GetShared(test_var_name));
310  }
311 }
bool Has(const std::string &field_name) const
Predicate to check if a field is registered with name field_name.
std::vector< std::unique_ptr< mfem::ParGridFunction > > _xs
Gridfunctions for setting Dirichlet BCs.
std::vector< std::unique_ptr< mfem::ParGridFunction > > _dxdts
mfem::AssemblyLevel _assembly_level
std::vector< mfem::ParFiniteElementSpace * > _test_pfespaces
Pointers to finite element spaces associated with test variables.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
std::shared_ptr< T > GetShared(const std::string &field_name) const
Returns a shared pointer to the field. This is guaranteed to return a non-null shared pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
void Register(const std::string &field_name, FieldArgs &&... args)
Construct new field with name field_name and register.
Moose::MFEM::GridFunctions _trial_variables
Pointers to trial variables.

◆ Mult()

void Moose::MFEM::EquationSystem::Mult ( const mfem::Vector &  u,
mfem::Vector &  residual 
) const
overrideinherited

Compute residual y = Mu.

Definition at line 262 of file EquationSystem.C.

263 {
264  _jacobian->Mult(x, residual);
265  x.HostRead();
266  residual.HostRead();
267 }
mfem::OperatorHandle _jacobian

◆ RecoverFEMSolution()

void Moose::MFEM::EquationSystem::RecoverFEMSolution ( mfem::BlockVector &  trueX,
Moose::MFEM::GridFunctions gridfunctions 
)
virtualinherited

Update variable from solution vector after solve.

Definition at line 276 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystemProblemOperator::Solve().

278 {
279  for (const auto i : index_range(_trial_var_names))
280  {
281  auto & trial_var_name = _trial_var_names.at(i);
282  trueX.GetBlock(i).SyncAliasMemory(trueX);
283  gridfunctions.Get(trial_var_name)->Distribute(&(trueX.GetBlock(i)));
284  }
285 }
std::vector< std::string > _trial_var_names
Names of all variables corresponding to gridfunctions.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
auto index_range(const T &sizable)

◆ SetTimeStep()

void Moose::MFEM::TimeDependentEquationSystem::SetTimeStep ( double  dt)
virtual

Definition at line 420 of file EquationSystem.C.

Referenced by Moose::MFEM::TimeDomainEquationSystemProblemOperator::BuildEquationSystemOperator().

421 {
422  if (fabs(dt - _dt_coef.constant) > 1.0e-12 * dt)
423  {
424  _dt_coef.constant = dt;
425  for (auto test_var_name : _test_var_names)
426  {
427  auto blf = _blfs.Get(test_var_name);
428  blf->Update();
429  blf->Assemble();
430  }
431  }
432 }
mfem::ConstantCoefficient _dt_coef
Coefficient for timestep scaling.
T * Get(const std::string &field_name) const
Returns a non-owning pointer to the field. This is guaranteed to return a non-null pointer...
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.
Moose::MFEM::NamedFieldsMap< mfem::ParBilinearForm > _blfs

◆ TestVarNames()

const std::vector<std::string>& Moose::MFEM::EquationSystem::TestVarNames ( ) const
inlineinherited

Definition at line 87 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystemProblemOperator::SetGridFunctions(), and Moose::MFEM::TimeDomainEquationSystemProblemOperator::SetGridFunctions().

87 { return _test_var_names; }
std::vector< std::string > _test_var_names
Names of all test variables corresponding to linear forms in this equation system.

◆ TrialVarNames()

const std::vector<std::string>& Moose::MFEM::EquationSystem::TrialVarNames ( ) const
inlineinherited

Definition at line 86 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystemProblemOperator::SetGridFunctions(), and Moose::MFEM::TimeDomainEquationSystemProblemOperator::SetGridFunctions().

86 { return _trial_var_names; }
std::vector< std::string > _trial_var_names
Names of all variables corresponding to gridfunctions.

◆ TrialVarTimeDerivativeNames()

const std::vector< std::string > & Moose::MFEM::TimeDependentEquationSystem::TrialVarTimeDerivativeNames ( ) const
inline

Definition at line 316 of file EquationSystem.h.

317 {
319 }
std::vector< std::string > _trial_var_time_derivative_names

◆ UpdateEquationSystem()

void Moose::MFEM::TimeDependentEquationSystem::UpdateEquationSystem ( )
virtual

◆ VectorContainsName()

bool Moose::MFEM::EquationSystem::VectorContainsName ( const std::vector< std::string > &  the_vector,
const std::string &  name 
) const
protectedinherited

Definition at line 30 of file EquationSystem.C.

Referenced by Moose::MFEM::EquationSystem::AddTestVariableNameIfMissing(), Moose::MFEM::EquationSystem::AddTrialVariableNameIfMissing(), and AddTrialVariableNameIfMissing().

32 {
33 
34  auto iter = std::find(the_vector.begin(), the_vector.end(), name);
35 
36  return (iter != the_vector.end());
37 }
std::string name(const ElemQuality q)

Member Data Documentation

◆ _assembly_level

mfem::AssemblyLevel Moose::MFEM::EquationSystem::_assembly_level
protectedinherited

◆ _blfs

Moose::MFEM::NamedFieldsMap<mfem::ParBilinearForm> Moose::MFEM::EquationSystem::_blfs
protectedinherited

◆ _dt_coef

mfem::ConstantCoefficient Moose::MFEM::TimeDependentEquationSystem::_dt_coef
protected

Coefficient for timestep scaling.

Definition at line 306 of file EquationSystem.h.

Referenced by BuildBilinearForms(), FormLegacySystem(), FormSystem(), and SetTimeStep().

◆ _dxdts

std::vector<std::unique_ptr<mfem::ParGridFunction> > Moose::MFEM::EquationSystem::_dxdts
protectedinherited

Definition at line 158 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystem::Init().

◆ _ess_tdof_lists

std::vector<mfem::Array<int> > Moose::MFEM::EquationSystem::_ess_tdof_lists
inherited

◆ _essential_bc_map

Moose::MFEM::NamedFieldsMap<std::vector<std::shared_ptr<MFEMEssentialBC> > > Moose::MFEM::EquationSystem::_essential_bc_map
protectedinherited

Arrays to store essential BCs to act on each component of weak form.

Named according to test variable.

Definition at line 173 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystem::AddEssentialBC(), and Moose::MFEM::EquationSystem::ApplyEssentialBCs().

◆ _h_blocks

mfem::Array2D<const mfem::HypreParMatrix *> Moose::MFEM::EquationSystem::_h_blocks
protectedinherited

◆ _integrated_bc_map

Moose::MFEM::NamedFieldsMap< Moose::MFEM::NamedFieldsMap<std::vector<std::shared_ptr<MFEMIntegratedBC> > > > Moose::MFEM::EquationSystem::_integrated_bc_map
protectedinherited

Arrays to store integrated BCs to act on each component of weak form.

Named according to test and trial variables.

Definition at line 170 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystem::AddIntegratedBC(), Moose::MFEM::EquationSystem::BuildBilinearForms(), BuildBilinearForms(), and Moose::MFEM::EquationSystem::BuildLinearForms().

◆ _jacobian

mfem::OperatorHandle Moose::MFEM::EquationSystem::_jacobian
mutableprotectedinherited

◆ _kernels_map

Moose::MFEM::NamedFieldsMap<Moose::MFEM::NamedFieldsMap<std::vector<std::shared_ptr<MFEMKernel> > > > Moose::MFEM::EquationSystem::_kernels_map
protectedinherited

Arrays to store kernels to act on each component of weak form.

Named according to test and trial variables.

Definition at line 165 of file EquationSystem.h.

Referenced by Moose::MFEM::EquationSystem::AddKernel(), Moose::MFEM::EquationSystem::BuildBilinearForms(), Moose::MFEM::EquationSystem::BuildLinearForms(), and Moose::MFEM::EquationSystem::BuildMixedBilinearForms().

◆ _lfs

Moose::MFEM::NamedFieldsMap<mfem::ParLinearForm> Moose::MFEM::EquationSystem::_lfs
protectedinherited

◆ _mblfs

Moose::MFEM::NamedFieldsMap<Moose::MFEM::NamedFieldsMap<mfem::ParMixedBilinearForm> > Moose::MFEM::EquationSystem::_mblfs
protectedinherited

◆ _nlfs

Moose::MFEM::NamedFieldsMap<mfem::ParNonlinearForm> Moose::MFEM::EquationSystem::_nlfs
protectedinherited

Definition at line 118 of file EquationSystem.h.

◆ _td_blfs

Moose::MFEM::NamedFieldsMap<mfem::ParBilinearForm> Moose::MFEM::TimeDependentEquationSystem::_td_blfs
protected

Container to store contributions to weak form of the form (F du/dt, v)

Definition at line 312 of file EquationSystem.h.

Referenced by BuildBilinearForms(), FormLegacySystem(), and FormSystem().

◆ _td_kernels_map

Moose::MFEM::NamedFieldsMap<Moose::MFEM::NamedFieldsMap<std::vector<std::shared_ptr<MFEMKernel> > > > Moose::MFEM::TimeDependentEquationSystem::_td_kernels_map
protected

Definition at line 310 of file EquationSystem.h.

Referenced by AddKernel(), and BuildBilinearForms().

◆ _test_pfespaces

std::vector<mfem::ParFiniteElementSpace *> Moose::MFEM::EquationSystem::_test_pfespaces
protectedinherited

◆ _test_var_names

std::vector<std::string> Moose::MFEM::EquationSystem::_test_var_names
protectedinherited

◆ _trial_var_names

std::vector<std::string> Moose::MFEM::EquationSystem::_trial_var_names
protectedinherited

◆ _trial_var_time_derivative_names

std::vector<std::string> Moose::MFEM::TimeDependentEquationSystem::_trial_var_time_derivative_names
protected

Definition at line 307 of file EquationSystem.h.

Referenced by AddTrialVariableNameIfMissing(), and TrialVarTimeDerivativeNames().

◆ _trial_variables

Moose::MFEM::GridFunctions Moose::MFEM::EquationSystem::_trial_variables
protectedinherited

Pointers to trial variables.

Definition at line 109 of file EquationSystem.h.

Referenced by FormLegacySystem(), FormSystem(), and Moose::MFEM::EquationSystem::Init().

◆ _xs

std::vector<std::unique_ptr<mfem::ParGridFunction> > Moose::MFEM::EquationSystem::_xs
protectedinherited

The documentation for this class was generated from the following files: