https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
MultiSpeciesDiffusionCG Class Reference

Creates all the objects needed to solve diffusion equations for multiple species with a continuous Galerkin finite element discretization. More...

#include <MultiSpeciesDiffusionCG.h>

Inheritance diagram for MultiSpeciesDiffusionCG:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 MultiSpeciesDiffusionCG (const InputParameters &parameters)
 
void addComponent (const ActionComponent &component) override
 
virtual InputParameters getAdditionalRMParams () const
 
virtual void act () override final
 
virtual void actOnAdditionalTasks ()
 
void addBlocks (const std::vector< SubdomainName > &blocks)
 
void addBlocksById (const std::vector< SubdomainID > &block_ids)
 
const std::vector< SubdomainName > & blocks () const
 
bool checkBlockRestrictionIdentical (const std::string &object_name, const std::vector< SubdomainName > &blocks, const bool error_if_not_identical=true) const
 
const T * getCoupledPhysics (const PhysicsName &phys_name, const bool allow_fail=false) const
 
const std::vector< T *> getCoupledPhysics (const bool allow_fail=false) const
 
unsigned int dimension () const
 
const ActionComponentgetActionComponent (const ComponentName &comp_name) const
 
void checkComponentType (const ActionComponent &component) const
 
const std::vector< VariableName > & solverVariableNames () const
 
const std::vector< VariableName > & auxVariableNames () const
 
void timedAct ()
 
MooseObjectName uniqueActionName () const
 
const std::string & specificTaskName () const
 
const std::set< std::string > & getAllTasks () const
 
void appendTask (const std::string &task)
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
PerfGraphperfGraph ()
 
void assertParamDefined (const std::string &libmesh_dbg_var(param)) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 

Public Attributes

const ConsoleStream _console
 

Static Public Attributes

static constexpr auto SYSTEM
 
static constexpr auto NAME
 

Protected Member Functions

virtual void addSolverVariables () override
 
virtual void addFEKernels () override
 
virtual void addFEBCs () override
 
void assertParamDefined (const std::string &param) const
 
bool isTransient () const
 
FactorygetFactory ()
 
FactorygetFactory () const
 
virtual FEProblemBasegetProblem ()
 
virtual const FEProblemBasegetProblem () const
 
void prepareCopyVariablesFromMesh () const
 
void copyVariablesFromMesh (const std::vector< VariableName > &variables_to_copy, bool are_nonlinear=true)
 
std::string prefix () const
 
void saveSolverVariableName (const VariableName &var_name)
 
void saveAuxVariableName (const VariableName &var_name)
 
bool variableExists (const VariableName &var_name, bool error_if_aux) const
 
bool solverVariableExists (const VariableName &var_name) const
 
const SolverSystemName & getSolverSystem (unsigned int variable_index) const
 
const SolverSystemName & getSolverSystem (const VariableName &variable_name) const
 
void addRequiredPhysicsTask (const std::string &task)
 
void assignBlocks (InputParameters &params, const std::vector< SubdomainName > &blocks) const
 
bool allMeshBlocks (const std::vector< SubdomainName > &blocks) const
 
bool allMeshBlocks (const std::set< SubdomainName > &blocks) const
 
std::set< SubdomainIDgetSubdomainIDs (const std::set< SubdomainName > &blocks) const
 
std::vector< std::string > getSubdomainNamesAndIDs (const std::set< SubdomainID > &blocks) const
 
void addPetscPairsToPetscOptions (const std::vector< std::pair< MooseEnumItem, std::string >> &petsc_pair_options)
 
bool isVariableFV (const VariableName &var_name) const
 
bool isVariableScalar (const VariableName &var_name) const
 
bool shouldCreateVariable (const VariableName &var_name, const std::vector< SubdomainName > &blocks, const bool error_if_aux)
 
bool shouldCreateIC (const VariableName &var_name, const std::vector< SubdomainName > &blocks, const bool ic_is_default_ic, const bool error_if_already_defined) const
 
bool shouldCreateTimeDerivative (const VariableName &var_name, const std::vector< SubdomainName > &blocks, const bool error_if_already_defined) const
 
void reportPotentiallyMissedParameters (const std::vector< std::string > &param_names, const std::string &object_type) const
 
bool addRelationshipManagers (Moose::RelationshipManagerType when_type, const InputParameters &moose_object_pars)
 
void associateWithParameter (const std::string &param_name, InputParameters &params) const
 
void associateWithParameter (const InputParameters &from_params, const std::string &param_name, InputParameters &params) const
 
const T & getMeshProperty (const std::string &data_name, const std::string &prefix)
 
const T & getMeshProperty (const std::string &data_name)
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name) const
 
bool hasMeshProperty (const std::string &data_name) const
 
std::string meshPropertyName (const std::string &data_name) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 
std::string timedSectionName (const std::string &section_name) const
 
void checkParamsBothSetOrNotSet (const std::string &param1, const std::string &param2) const
 
void checkSecondParamSetOnlyIfFirstOneTrue (const std::string &param1, const std::string &param2) const
 
void checkSecondParamSetOnlyIfFirstOneSet (const std::string &param1, const std::string &param2) const
 
void checkSecondParamNotSetIfFirstOneSet (const std::string &param1, const std::string &param2) const
 
void checkVectorParamsSameLength (const std::string &param1, const std::string &param2) const
 
void checkVectorParamAndMultiMooseEnumLength (const std::string &param1, const std::string &param2) const
 
void checkTwoDVectorParamsSameLength (const std::string &param1, const std::string &param2) const
 
void checkVectorParamsNoOverlap (const std::vector< std::string > &param_vecs) const
 
void checkTwoDVectorParamsNoRespectiveOverlap (const std::vector< std::string > &param_vecs) const
 
void checkTwoDVectorParamInnerSameLengthAsOneDVector (const std::string &param1, const std::string &param2) const
 
void checkTwoDVectorParamMultiMooseEnumSameLength (const std::string &param1, const std::string &param2, const bool error_for_param2) const
 
void checkVectorParamNotEmpty (const std::string &param1) const
 
void checkVectorParamsSameLengthIfSet (const std::string &param1, const std::string &param2, const bool ignore_empty_default_param2=false) const
 
void checkVectorParamLengthSameAsCombinedOthers (const std::string &param1, const std::string &param2, const std::string &param3) const
 
void checkBlockwiseConsistency (const std::string &block_param_name, const std::vector< std::string > &parameter_names) const
 
bool parameterConsistent (const InputParameters &other_param, const std::string &param_name) const
 
void warnInconsistent (const InputParameters &parameters, const std::string &param_name) const
 
void errorDependentParameter (const std::string &param1, const std::string &value_not_set, const std::vector< std::string > &dependent_params) const
 
void errorInconsistentDependentParameter (const std::string &param1, const std::string &value_set, const std::vector< std::string > &dependent_params) const
 

Static Protected Member Functions

static std::string meshPropertyName (const std::string &data_name, const std::string &prefix)
 

Protected Attributes

const std::vector< VariableName > & _species_names
 Name of the diffused variables. More...
 
const unsigned int _num_species
 Number of species. More...
 
const std::vector< std::vector< BoundaryName > > & _neumann_boundaries
 Boundaries on which a Neumann boundary condition is applied. Outer indexing is variables. More...
 
const std::vector< std::vector< BoundaryName > > & _dirichlet_boundaries
 Boundaries on which a Dirichlet boundary condition is applied. Outer indexing is variables. More...
 
const bool _use_ad
 Whether to use automatic differentiation or not. More...
 
std::vector< SolverSystemName > _system_names
 
std::vector< unsigned int_system_numbers
 
const bool _verbose
 
const MooseEnum_preconditioning
 
std::vector< SubdomainName > _blocks
 
std::string _registered_identifier
 
std::string _specific_task_name
 
std::set< std::string > _all_tasks
 
ActionWarehouse_awh
 
const std::string & _current_task
 
std::shared_ptr< MooseMesh > & _mesh
 
std::shared_ptr< MooseMesh > & _displaced_mesh
 
std::shared_ptr< FEProblemBase > & _problem
 
PerfID _act_timer
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
MooseApp_pg_moose_app
 
const std::string _prefix
 
const Parallel::Communicator & _communicator
 

Detailed Description

Creates all the objects needed to solve diffusion equations for multiple species with a continuous Galerkin finite element discretization.

Definition at line 18 of file MultiSpeciesDiffusionCG.h.

Constructor & Destructor Documentation

◆ MultiSpeciesDiffusionCG()

MultiSpeciesDiffusionCG::MultiSpeciesDiffusionCG ( const InputParameters parameters)

Definition at line 30 of file MultiSpeciesDiffusionCG.C.

32 {
33 }
MultiSpeciesDiffusionPhysicsBase(const InputParameters &parameters)
const InputParameters & parameters() const

Member Function Documentation

◆ addComponent()

void MultiSpeciesDiffusionPhysicsBase::addComponent ( const ActionComponent component)
overridevirtualinherited

Reimplemented from PhysicsBase.

Definition at line 158 of file MultiSpeciesDiffusionPhysicsBase.C.

159 {
160  for (const auto & block : component.blocks())
161  _blocks.push_back(block);
162 }
static const std::string component
Definition: NS.h:153
std::vector< SubdomainName > _blocks

◆ addFEBCs()

void MultiSpeciesDiffusionCG::addFEBCs ( )
overrideprotectedvirtual

Reimplemented from PhysicsBase.

Definition at line 132 of file MultiSpeciesDiffusionCG.C.

133 {
134  if (isParamSetByUser("neumann_boundaries"))
135  {
136  const auto & boundary_fluxes =
137  getParam<std::vector<std::vector<MooseFunctorName>>>("boundary_fluxes");
138 
139  for (const auto s : index_range(_species_names))
140  {
141  const auto & var_name = _species_names[s];
142 
143  for (const auto i : index_range(_neumann_boundaries[s]))
144  {
145  const auto & bc_flux = boundary_fluxes[s][i];
146  // Select the boundary type based on the user parameters and what we know to be most
147  // efficient We could actually just use the very last option for everything but the
148  // performance is better if one uses the specialized objects
149  std::string bc_type = "";
150  if (MooseUtils::parsesToReal(bc_flux))
151  bc_type = _use_ad ? "ADNeumannBC" : "NeumannBC";
152  else if (getProblem().hasVariable(bc_flux))
153  bc_type = "CoupledVarNeumannBC"; // not AD, but still perfect Jacobian
154  else if (getProblem().hasFunction(bc_flux))
155  bc_type = _use_ad ? "ADFunctionNeumannBC" : "FunctionNeumannBC";
156  else if (getProblem().hasPostprocessorValueByName(bc_flux))
157  bc_type = "PostprocessorNeumannBC";
158  else // this is AD, but we can mix AD and non-AD
159  bc_type = "FunctorNeumannBC";
160 
161  // Get the parameters for the object type chosen and set the common parameters
162  InputParameters params = getFactory().getValidParams(bc_type);
163  params.set<NonlinearVariableName>("variable") = var_name;
164  params.set<std::vector<BoundaryName>>("boundary") = {_neumann_boundaries[s][i]};
165 
166  // Set the flux parameter for the specific type of NeumannBC used
167  if (MooseUtils::parsesToReal(bc_flux))
168  params.set<Real>("value") = MooseUtils::convert<Real>(bc_flux);
169  else if (getProblem().hasVariable(bc_flux))
170  params.set<std::vector<VariableName>>("v") = {bc_flux};
171  else if (getProblem().hasFunction(bc_flux))
172  params.set<FunctionName>("function") = bc_flux;
173  else if (getProblem().hasPostprocessorValueByName(bc_flux))
174  params.set<PostprocessorName>("postprocessor") = bc_flux;
175  else
176  params.set<MooseFunctorName>("functor") = bc_flux;
177 
179  bc_type, prefix() + var_name + "_neumann_bc_" + _neumann_boundaries[s][i], params);
180  }
181  }
182  }
183  if (isParamSetByUser("dirichlet_boundaries"))
184  {
185  const auto & boundary_values =
186  getParam<std::vector<std::vector<MooseFunctorName>>>("boundary_values");
187  for (const auto s : index_range(_species_names))
188  {
189  const auto & var_name = _species_names[s];
190  for (const auto i : index_range(_dirichlet_boundaries[s]))
191  {
192  const auto & bc_value = boundary_values[s][i];
193  // Select the boundary type based on the user parameters and what we know to be most
194  // efficient
195  std::string bc_type = "";
196  if (MooseUtils::parsesToReal(bc_value))
197  bc_type = _use_ad ? "ADDirichletBC" : "DirichletBC";
198  else if (getProblem().hasVariable(bc_value))
199  bc_type = _use_ad ? "ADMatchedValueBC" : "MatchedValueBC";
200  else if (getProblem().hasFunction(bc_value))
201  bc_type = _use_ad ? "ADFunctionDirichletBC" : "FunctionDirichletBC";
202  else if (getProblem().hasPostprocessorValueByName(bc_value))
203  bc_type = "PostprocessorDirichletBC";
204  else // this is AD, but we can mix AD and non-AD
205  bc_type = "FunctorDirichletBC";
206 
207  InputParameters params = getFactory().getValidParams(bc_type);
208  params.set<NonlinearVariableName>("variable") = var_name;
209  params.set<std::vector<BoundaryName>>("boundary") = {_dirichlet_boundaries[s][i]};
210 
211  // Set the flux parameter for the specific type of DirichletBC used
212  if (MooseUtils::parsesToReal(bc_value))
213  params.set<Real>("value") = MooseUtils::convert<Real>(bc_value);
214  else if (getProblem().hasVariable(bc_value))
215  params.set<std::vector<VariableName>>("v") = {bc_value};
216  else if (getProblem().hasFunction(bc_value))
217  params.set<FunctionName>("function") = bc_value;
218  else if (getProblem().hasPostprocessorValueByName(bc_value))
219  params.set<PostprocessorName>("postprocessor") = bc_value;
220  else
221  params.set<MooseFunctorName>("functor") = bc_value;
222 
224  bc_type, prefix() + var_name + "_dirichlet_bc_" + _dirichlet_boundaries[s][i], params);
225  }
226  }
227  }
228 }
std::string prefix() const
virtual bool hasVariable(const std::string &var_name) const override
bool parsesToReal(const std::string &input)
Factory & getFactory()
const std::vector< std::vector< BoundaryName > > & _neumann_boundaries
Boundaries on which a Neumann boundary condition is applied. Outer indexing is variables.
T & set(const std::string &name, bool quiet_mode=false)
InputParameters getValidParams(const std::string &name) const
virtual void addBoundaryCondition(const std::string &bc_name, const std::string &name, InputParameters &parameters)
virtual FEProblemBase & getProblem()
bool hasPostprocessorValueByName(const PostprocessorName &name) const
const std::vector< VariableName > & _species_names
Name of the diffused variables.
bool isParamSetByUser(const std::string &nm) const
const std::vector< std::vector< BoundaryName > > & _dirichlet_boundaries
Boundaries on which a Dirichlet boundary condition is applied. Outer indexing is variables.
const bool _use_ad
Whether to use automatic differentiation or not.
virtual bool hasFunction(const std::string &name, const THREAD_ID tid=0)
auto index_range(const T &sizable)

◆ addFEKernels()

void MultiSpeciesDiffusionCG::addFEKernels ( )
overrideprotectedvirtual

Reimplemented from PhysicsBase.

Definition at line 36 of file MultiSpeciesDiffusionCG.C.

37 {
38  for (const auto s : index_range(_species_names))
39  {
40  const auto & var_name = _species_names[s];
41  // Diffusion term
42  if (isParamValid("diffusivity_matprops") || isParamValid("diffusivity_functors"))
43  {
44  // Select the kernel type based on the user parameters
45  std::string kernel_type;
46  if (isParamValid("diffusivity_matprops"))
47  kernel_type = _use_ad ? "ADMatDiffusion" : "MatDiffusion";
48  else if (isParamValid("diffusivity_functors"))
49  {
50  const auto & d = getParam<std::vector<MooseFunctorName>>("diffusivity_functors")[s];
51  if (getProblem().hasFunction(d))
52  kernel_type = "FunctionDiffusion";
53  else
54  paramError(
55  "diffusivity_functors", "No diffusion kernel implemented for the source type of", d);
56  }
57 
58  InputParameters params = getFactory().getValidParams(kernel_type);
59  params.set<NonlinearVariableName>("variable") = var_name;
60  assignBlocks(params, _blocks);
61 
62  // Transfer the diffusivity parameter from the Physics to the kernel
63  if (isParamValid("diffusivity_matprops"))
64  params.set<MaterialPropertyName>("diffusivity") =
65  getParam<std::vector<MaterialPropertyName>>("diffusivity_matprops")[s];
66  else if (isParamValid("diffusivity_functors"))
67  params.set<FunctionName>("function") =
68  getParam<std::vector<MooseFunctorName>>("diffusivity_functors")[s];
69 
70  getProblem().addKernel(kernel_type, prefix() + var_name + "_diffusion", params);
71  }
72 
73  // Source term
74  if (isParamValid("source_functors"))
75  {
76  // Select the kernel type based on the user parameters
77  std::string kernel_type;
78  const auto & sources = getParam<std::vector<MooseFunctorName>>("source_functors");
79  const auto & source = sources[s];
80  if (MooseUtils::parsesToReal(source) || getProblem().hasFunction(source) ||
82  kernel_type = _use_ad ? "ADBodyForce" : "BodyForce";
83  else if (getProblem().hasVariable(source))
84  kernel_type = _use_ad ? "ADCoupledForce" : "CoupledForce";
85  else
86  paramError("source_functors",
87  "No kernel defined for a source term in CG for the type of '",
88  source,
89  "'");
90 
91  InputParameters params = getFactory().getValidParams(kernel_type);
92  params.set<NonlinearVariableName>("variable") = var_name;
93  assignBlocks(params, _blocks);
94 
95  // Transfer the source and coefficient parameter from the Physics to the kernel
96  const auto coefs = getParam<std::vector<Real>>("source_coefs");
97  const auto coef = coefs[s];
98  if (MooseUtils::parsesToReal(source))
99  params.set<Real>("value") = MooseUtils::convert<Real>(source) * coef;
100  else if (getProblem().hasFunction(source))
101  {
102  params.set<Real>("value") = coef;
103  params.set<FunctionName>("function") = source;
104  }
105  else if (getProblem().hasPostprocessorValueByName(source))
106  {
107  params.set<Real>("value") = coef;
108  params.set<PostprocessorName>("postprocessor") = source;
109  }
110  else
111  {
112  params.set<Real>("coef") = coef;
113  params.set<std::vector<VariableName>>("v") = {source};
114  }
115 
116  getProblem().addKernel(kernel_type, prefix() + var_name + "_source", params);
117  }
118 
119  // Time derivative term
120  if (shouldCreateTimeDerivative(var_name, _blocks, false))
121  {
122  const std::string kernel_type = _use_ad ? "ADTimeDerivative" : "TimeDerivative";
123  InputParameters params = getFactory().getValidParams(kernel_type);
124  params.set<NonlinearVariableName>("variable") = var_name;
125  assignBlocks(params, _blocks);
126  getProblem().addKernel(kernel_type, prefix() + var_name + "_time", params);
127  }
128  }
129 }
std::string prefix() const
virtual bool hasVariable(const std::string &var_name) const override
bool parsesToReal(const std::string &input)
void assignBlocks(InputParameters &params, const std::vector< SubdomainName > &blocks) const
Factory & getFactory()
T & set(const std::string &name, bool quiet_mode=false)
InputParameters getValidParams(const std::string &name) const
bool shouldCreateTimeDerivative(const VariableName &var_name, const std::vector< SubdomainName > &blocks, const bool error_if_already_defined) const
virtual void addKernel(const std::string &kernel_name, const std::string &name, InputParameters &parameters)
std::vector< SubdomainName > _blocks
bool isParamValid(const std::string &name) const
virtual FEProblemBase & getProblem()
bool hasPostprocessorValueByName(const PostprocessorName &name) const
const T & getParam(const std::string &name) const
void paramError(const std::string &param, Args... args) const
const std::vector< VariableName > & _species_names
Name of the diffused variables.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const bool _use_ad
Whether to use automatic differentiation or not.
virtual bool hasFunction(const std::string &name, const THREAD_ID tid=0)
auto index_range(const T &sizable)

◆ addSolverVariables()

void MultiSpeciesDiffusionCG::addSolverVariables ( )
overrideprotectedvirtual

Reimplemented from PhysicsBase.

Definition at line 231 of file MultiSpeciesDiffusionCG.C.

232 {
233  for (const auto & var_name : _species_names)
234  {
235  // If the variable was added outside the Physics
236  if (variableExists(var_name, /*error_if_aux*/ true))
237  {
238  if (isParamValid("variable_order"))
239  paramError("variable_order",
240  "Cannot specify the variable order if variable " + var_name +
241  " is defined outside the Physics block");
242  else
243  return;
244  }
245 
246  const std::string variable_type = "MooseVariable";
247  InputParameters params = getFactory().getValidParams(variable_type);
248  params.set<MooseEnum>("order") = getParam<MooseEnum>("variable_order");
249  assignBlocks(params, _blocks);
250  params.set<SolverSystemName>("solver_sys") = getSolverSystem(var_name);
251 
252  getProblem().addVariable(variable_type, var_name, params);
253  }
254 }
void assignBlocks(InputParameters &params, const std::vector< SubdomainName > &blocks) const
Factory & getFactory()
T & set(const std::string &name, bool quiet_mode=false)
InputParameters getValidParams(const std::string &name) const
std::vector< SubdomainName > _blocks
bool isParamValid(const std::string &name) const
virtual FEProblemBase & getProblem()
const SolverSystemName & getSolverSystem(unsigned int variable_index) const
void paramError(const std::string &param, Args... args) const
virtual void addVariable(const std::string &var_type, const std::string &var_name, InputParameters &params)
const std::vector< VariableName > & _species_names
Name of the diffused variables.
bool variableExists(const VariableName &var_name, bool error_if_aux) const

◆ validParams()

InputParameters MultiSpeciesDiffusionCG::validParams ( )
static

Definition at line 20 of file MultiSpeciesDiffusionCG.C.

21 {
23  params.addClassDescription("Discretizes diffusion equations for several species with the "
24  "continuous Galerkin finite element method");
25  params.transferParam<MooseEnum>(MooseVariableBase::validParams(), "order", "variable_order");
26 
27  return params;
28 }
static InputParameters validParams()
void transferParam(const InputParameters &source_param, const std::string &name, const std::string &new_name="", const std::string &new_description="")
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _dirichlet_boundaries

const std::vector<std::vector<BoundaryName> >& MultiSpeciesDiffusionPhysicsBase::_dirichlet_boundaries
protectedinherited

Boundaries on which a Dirichlet boundary condition is applied. Outer indexing is variables.

Definition at line 43 of file MultiSpeciesDiffusionPhysicsBase.h.

Referenced by addFEBCs().

◆ _neumann_boundaries

const std::vector<std::vector<BoundaryName> >& MultiSpeciesDiffusionPhysicsBase::_neumann_boundaries
protectedinherited

Boundaries on which a Neumann boundary condition is applied. Outer indexing is variables.

Definition at line 41 of file MultiSpeciesDiffusionPhysicsBase.h.

Referenced by addFEBCs().

◆ _num_species

const unsigned int MultiSpeciesDiffusionPhysicsBase::_num_species
protectedinherited

Number of species.

Definition at line 39 of file MultiSpeciesDiffusionPhysicsBase.h.

◆ _species_names

const std::vector<VariableName>& MultiSpeciesDiffusionPhysicsBase::_species_names
protectedinherited

◆ _use_ad

const bool MultiSpeciesDiffusionPhysicsBase::_use_ad
protectedinherited

Whether to use automatic differentiation or not.

Definition at line 46 of file MultiSpeciesDiffusionPhysicsBase.h.

Referenced by addFEBCs(), addFEKernels(), and MultiSpeciesDiffusionPhysicsBase::addPostprocessors().


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