OGS
ProcessLib::ProcessVariable Class Reference

Detailed Description

A named process variable. Its properties includes the mesh, and the initial and boundary conditions as well as the source terms.

Definition at line 45 of file ProcessVariable.h.

#include <ProcessVariable.h>

Collaboration diagram for ProcessLib::ProcessVariable:
[legend]

Public Member Functions

 ProcessVariable (BaseLib::ConfigTree const &config, MeshLib::Mesh &mesh, std::vector< std::unique_ptr< MeshLib::Mesh >> const &meshes, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation >> const &curves)
 
 ProcessVariable (ProcessVariable &&other)
 
std::string const & getName () const
 
MeshLib::Mesh const & getMesh () const
 Returns a mesh on which the process variable is defined. More...
 
std::vector< std::unique_ptr< DeactivatedSubdomain const > > const & getDeactivatedSubdomains () const
 
void updateDeactivatedSubdomains (double const time)
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
int getNumberOfGlobalComponents () const
 Returns the number of components of the process variable. More...
 
std::vector< std::unique_ptr< BoundaryCondition > > createBoundaryConditions (const NumLib::LocalToGlobalIndexMap &dof_table, const int variable_id, unsigned const integration_order, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters, Process const &process)
 
std::vector< std::unique_ptr< SourceTerm > > createSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int variable_id, unsigned const integration_order, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters)
 
ParameterLib::Parameter< double > const & getInitialCondition () const
 
unsigned getShapeFunctionOrder () const
 

Private Member Functions

void createBoundaryConditionsForDeactivatedSubDomains (const NumLib::LocalToGlobalIndexMap &dof_table, const int variable_id, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters, std::vector< std::unique_ptr< BoundaryCondition >> &bcs)
 

Private Attributes

std::string const _name
 
MeshLib::Mesh_mesh
 
const int _n_components
 
unsigned _shapefunction_order
 
std::vector< std::unique_ptr< DeactivatedSubdomain const > > _deactivated_subdomains
 
std::vector< std::size_t > _ids_of_active_elements
 
ParameterLib::Parameter< double > const & _initial_condition
 
std::vector< BoundaryConditionConfig_bc_configs
 
std::vector< SourceTermConfig_source_term_configs
 

Constructor & Destructor Documentation

◆ ProcessVariable() [1/2]

ProcessLib::ProcessVariable::ProcessVariable ( BaseLib::ConfigTree const &  config,
MeshLib::Mesh mesh,
std::vector< std::unique_ptr< MeshLib::Mesh >> const &  meshes,
std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &  parameters,
std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation >> const &  curves 
)
Input File Parameter:
prj__process_variables__process_variable__boundary_conditions
Input File Parameter:
prj__process_variables__process_variable__boundary_conditions__boundary_condition
Input File Parameter:
prj__process_variables__process_variable__boundary_conditions__boundary_condition__component
Input File Parameter:
prj__process_variables__process_variable__source_terms
Input File Parameter:
prj__process_variables__process_variable__source_terms__source_term
Input File Parameter:
prj__process_variables__process_variable__source_terms__source_term__component

Definition at line 93 of file ProcessVariable.cpp.

101  :
102  _name(config.getConfigParameter<std::string>("name")),
103  _mesh(mesh),
105  _n_components(config.getConfigParameter<int>("components")),
107  _shapefunction_order(config.getConfigParameter<unsigned>("order")),
109  createDeactivatedSubdomains(config, mesh, parameters, curves)),
110  _initial_condition(ParameterLib::findParameter<double>(
112  config.getConfigParameter<std::string>("initial_condition"),
113  parameters, _n_components, &mesh))
114 {
115  DBUG("Constructing process variable {:s}", _name);
116 
118  {
119  OGS_FATAL("The given shape function order {:d} is not supported",
121  }
122 
123  // Boundary conditions
124  if (auto bcs_config =
126  config.getConfigSubtreeOptional("boundary_conditions"))
127  {
128  for (
129  auto bc_config :
131  bcs_config->getConfigSubtreeList("boundary_condition"))
132  {
133  auto const& bc_mesh = findMeshInConfig(bc_config, meshes);
134  auto component_id =
136  bc_config.getConfigParameterOptional<int>("component");
137 
138  if (!component_id && _n_components == 1)
139  {
140  // default value for single component vars.
141  component_id = 0;
142  }
143 
144  _bc_configs.emplace_back(std::move(bc_config), bc_mesh,
145  component_id);
146  }
147  }
148  else
149  {
150  INFO("No boundary conditions for process variable '{:s}' found.",
151  _name);
152  }
153 
154  // Source terms
156  if (auto sts_config = config.getConfigSubtreeOptional("source_terms"))
157  {
158  for (
159  auto st_config :
161  sts_config->getConfigSubtreeList("source_term"))
162  {
163  MeshLib::Mesh const& st_mesh = findMeshInConfig(st_config, meshes);
164  auto component_id =
166  st_config.getConfigParameterOptional<int>("component");
167 
168  if (!component_id && _n_components == 1)
169  {
170  // default value for single component vars.
171  component_id = 0;
172  }
173 
174  _source_term_configs.emplace_back(std::move(st_config), st_mesh,
175  component_id);
176  }
177  }
178  else
179  {
180  INFO("No source terms for process variable '{:s}' found.", _name);
181  }
182 }
#define OGS_FATAL(...)
Definition: Error.h:26
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
ParameterLib::Parameter< double > const & _initial_condition
std::vector< BoundaryConditionConfig > _bc_configs
std::vector< SourceTermConfig > _source_term_configs
std::vector< std::unique_ptr< DeactivatedSubdomain const > > _deactivated_subdomains
std::vector< std::unique_ptr< DeactivatedSubdomain const > > createDeactivatedSubdomains(BaseLib::ConfigTree const &config, MeshLib::Mesh const &mesh, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters, std::map< std::string, std::unique_ptr< MathLib::PiecewiseLinearInterpolation >> const &curves)
MeshLib::Mesh const & findMeshInConfig(BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh >> const &meshes)

References _bc_configs, _n_components, _name, _shapefunction_order, _source_term_configs, DBUG(), anonymous_namespace{ProcessVariable.cpp}::findMeshInConfig(), BaseLib::ConfigTree::getConfigSubtreeOptional(), INFO(), and OGS_FATAL.

◆ ProcessVariable() [2/2]

ProcessLib::ProcessVariable::ProcessVariable ( ProcessVariable &&  other)

Definition at line 184 of file ProcessVariable.cpp.

185  : _name(std::move(other._name)),
186  _mesh(other._mesh),
187  _n_components(other._n_components),
188  _shapefunction_order(other._shapefunction_order),
189  _deactivated_subdomains(std::move(other._deactivated_subdomains)),
190  _initial_condition(std::move(other._initial_condition)),
191  _bc_configs(std::move(other._bc_configs)),
192  _source_term_configs(std::move(other._source_term_configs))
193 {
194 }

Member Function Documentation

◆ createBoundaryConditions()

std::vector< std::unique_ptr< BoundaryCondition > > ProcessLib::ProcessVariable::createBoundaryConditions ( const NumLib::LocalToGlobalIndexMap dof_table,
const int  variable_id,
unsigned const  integration_order,
std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &  parameters,
Process const &  process 
)

Definition at line 207 of file ProcessVariable.cpp.

213 {
214  std::vector<std::unique_ptr<BoundaryCondition>> bcs;
215  bcs.reserve(_bc_configs.size());
216 
217  for (auto& config : _bc_configs)
218  {
219  auto bc = createBoundaryCondition(
220  config, dof_table, _mesh, variable_id, integration_order,
221  _shapefunction_order, parameters, process);
222 #ifdef USE_PETSC
223  if (bc == nullptr)
224  {
225  continue;
226  }
227 #endif // USE_PETSC
228  bcs.push_back(std::move(bc));
229  }
230 
231  if (_deactivated_subdomains.empty())
232  {
233  return bcs;
234  }
235 
237  parameters, bcs);
238  return bcs;
239 }
void createBoundaryConditionsForDeactivatedSubDomains(const NumLib::LocalToGlobalIndexMap &dof_table, const int variable_id, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters, std::vector< std::unique_ptr< BoundaryCondition >> &bcs)
std::unique_ptr< BoundaryCondition > createBoundaryCondition(const BoundaryConditionConfig &config, const NumLib::LocalToGlobalIndexMap &dof_table, const MeshLib::Mesh &bulk_mesh, const int variable_id, const unsigned integration_order, const unsigned shapefunction_order, const std::vector< std::unique_ptr< ParameterLib::ParameterBase >> &parameters, const Process &process)

References _bc_configs, _deactivated_subdomains, _mesh, _shapefunction_order, ProcessLib::createBoundaryCondition(), and createBoundaryConditionsForDeactivatedSubDomains().

Referenced by ProcessLib::BoundaryConditionCollection::addBCsForProcessVariables().

◆ createBoundaryConditionsForDeactivatedSubDomains()

void ProcessLib::ProcessVariable::createBoundaryConditionsForDeactivatedSubDomains ( const NumLib::LocalToGlobalIndexMap dof_table,
const int  variable_id,
std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &  parameters,
std::vector< std::unique_ptr< BoundaryCondition >> &  bcs 
)
private

Definition at line 241 of file ProcessVariable.cpp.

245 {
246  for (auto const& deactivated_subdomain : _deactivated_subdomains)
247  {
248  auto const& deactivated_subdomain_meshes =
249  deactivated_subdomain->deactivated_subdomain_meshes;
250  for (auto const& deactivated_subdomain_mesh :
251  deactivated_subdomain_meshes)
252  {
253  auto const* parameter = &ParameterLib::findParameter<double>(
255  bool const set_outer_nodes_dirichlet_values =
256  deactivated_subdomain->boundary_value_parameter != nullptr;
257  if (set_outer_nodes_dirichlet_values)
258  {
259  parameter = deactivated_subdomain->boundary_value_parameter;
260  }
261 
262  for (int component_id = 0;
263  component_id <
264  dof_table.getNumberOfVariableComponents(variable_id);
265  component_id++)
266  {
267  auto bc = std::make_unique<DeactivatedSubdomainDirichlet>(
269  deactivated_subdomain->time_interval, *parameter,
270  set_outer_nodes_dirichlet_values,
271  *deactivated_subdomain_mesh, dof_table, variable_id,
272  component_id);
273 
274 #ifdef USE_PETSC
275  // TODO: make it work under PETSc too.
276  if (bc == nullptr)
277  {
278  continue;
279  }
280 #endif // USE_PETSC
281  bcs.push_back(std::move(bc));
282  }
283  }
284  }
285 }
int getNumberOfVariableComponents(int variable_id) const
std::vector< std::size_t > _ids_of_active_elements
static const std::string zero_parameter_name

References _deactivated_subdomains, _ids_of_active_elements, NumLib::LocalToGlobalIndexMap::getNumberOfVariableComponents(), and ProcessLib::DeactivatedSubdomain::zero_parameter_name.

Referenced by createBoundaryConditions().

◆ createSourceTerms()

std::vector< std::unique_ptr< SourceTerm > > ProcessLib::ProcessVariable::createSourceTerms ( const NumLib::LocalToGlobalIndexMap dof_table,
const int  variable_id,
unsigned const  integration_order,
std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &  parameters 
)

Definition at line 335 of file ProcessVariable.cpp.

340 {
341  std::vector<std::unique_ptr<SourceTerm>> source_terms;
342 
343  transform(cbegin(_source_term_configs), cend(_source_term_configs),
344  back_inserter(source_terms),
345  [&](auto const& config)
346  {
347  return createSourceTerm(config, dof_table, config.mesh,
348  variable_id, integration_order,
349  _shapefunction_order, parameters);
350  });
351 
352  return source_terms;
353 }
std::unique_ptr< SourceTerm > createSourceTerm(const SourceTermConfig &config, const NumLib::LocalToGlobalIndexMap &dof_table_bulk, const MeshLib::Mesh &source_term_mesh, const int variable_id, const unsigned integration_order, const unsigned shapefunction_order, std::vector< std::unique_ptr< ParameterLib::ParameterBase >> const &parameters)

References _shapefunction_order, _source_term_configs, and ProcessLib::createSourceTerm().

Referenced by ProcessLib::SourceTermCollection::addSourceTermsForProcessVariables().

◆ getActiveElementIDs()

std::vector<std::size_t> const& ProcessLib::ProcessVariable::getActiveElementIDs ( ) const
inline

Definition at line 72 of file ProcessVariable.h.

73  {
75  }

References _ids_of_active_elements.

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::assembleConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleConcreteProcess(), ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::assembleConcreteProcess(), ProcessLib::TES::TESProcess::assembleConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::assembleConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::assembleWithJacobianConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleWithJacobianConcreteProcess(), ProcessLib::TES::TESProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::assembleWithJacobianConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::HeatConduction::HeatConductionProcess::computeSecondaryVariableConcrete(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::computeSecondaryVariableConcrete(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::LiquidFlow::LiquidFlowProcess::computeSecondaryVariableConcrete(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim >::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::ComponentTransportProcess::computeSecondaryVariableConcrete(), ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::ComponentTransportProcess::initializeConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::postIterationConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::postNonLinearSolverConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::postNonLinearSolverConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::postNonLinearSolverConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::postNonLinearSolverConcreteProcess(), ProcessLib::HT::HTProcess::postTimestepConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::postTimestepConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::postTimestepConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >::postTimestepConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::postTimestepConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::postTimestepConcreteProcess(), ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::postTimestepConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::postTimestepConcreteProcess(), ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >::preAssembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::preTimestepConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::preTimestepConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::preTimestepConcreteProcess(), ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::preTimestepConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::setCoupledTermForTheStaggeredSchemeToLocalAssemblers(), ProcessLib::HT::HTProcess::setCoupledTermForTheStaggeredSchemeToLocalAssemblers(), ProcessLib::ComponentTransport::ComponentTransportProcess::setInitialConditionsConcreteProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::setInitialConditionsConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::setInitialConditionsConcreteProcess(), and ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation().

◆ getDeactivatedSubdomains()

std::vector<std::unique_ptr<DeactivatedSubdomain const> > const& ProcessLib::ProcessVariable::getDeactivatedSubdomains ( ) const
inline

Definition at line 65 of file ProcessVariable.h.

66  {
68  }

References _deactivated_subdomains.

◆ getInitialCondition()

ParameterLib::Parameter<double> const& ProcessLib::ProcessVariable::getInitialCondition ( ) const
inline

◆ getMesh()

MeshLib::Mesh const & ProcessLib::ProcessVariable::getMesh ( ) const

Returns a mesh on which the process variable is defined.

Definition at line 201 of file ProcessVariable.cpp.

202 {
203  return _mesh;
204 }

References _mesh.

◆ getName()

◆ getNumberOfGlobalComponents()

◆ getShapeFunctionOrder()

◆ updateDeactivatedSubdomains()

void ProcessLib::ProcessVariable::updateDeactivatedSubdomains ( double const  time)

Definition at line 287 of file ProcessVariable.cpp.

288 {
289  _ids_of_active_elements.clear();
290 
291  // If none of the deactivated subdomains is active at current time, then the
292  // _ids_of_active_elements remain empty.
293  if (std::none_of(
295  [&](auto const& ds) { return ds->isInTimeSupportInterval(time); }))
296  {
297  return;
298  }
299 
300  auto const* const material_ids = MeshLib::materialIDs(_mesh);
301 
302  auto is_active_in_subdomain = [&](std::size_t const i,
303  DeactivatedSubdomain const& ds) -> bool
304  {
305  if (!ds.isInTimeSupportInterval(time))
306  {
307  return true;
308  }
309 
310  auto const& deactivated_materialIDs = ds.materialIDs;
311 
312  auto const& element_center = getCenterOfGravity(*_mesh.getElement(i));
313  if (std::binary_search(deactivated_materialIDs.begin(),
314  deactivated_materialIDs.end(),
315  (*material_ids)[i]) &&
316  ds.isDeactivated(element_center, time))
317  {
318  return false;
319  }
320  return true;
321  };
322 
323  auto const number_of_elements = _mesh.getNumberOfElements();
324  for (std::size_t i = 0; i < number_of_elements; i++)
325  {
326  if (std::all_of(
328  [&](auto const& ds) { return is_active_in_subdomain(i, *ds); }))
329  {
331  }
332  }
333 }
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
const Element * getElement(std::size_t idx) const
Get the element with the given index.
Definition: Mesh.h:77
std::size_t getNumberOfElements() const
Get the number of elements.
Definition: Mesh.h:86
PropertyVector< int > const * materialIDs(Mesh const &mesh)
Definition: Mesh.cpp:258
MeshLib::Node getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition: Element.cpp:126

References _deactivated_subdomains, _ids_of_active_elements, _mesh, MeshLib::getCenterOfGravity(), MeshLib::Mesh::getElement(), MeshLib::Element::getID(), MeshLib::Mesh::getNumberOfElements(), and MeshLib::materialIDs().

Member Data Documentation

◆ _bc_configs

std::vector<BoundaryConditionConfig> ProcessLib::ProcessVariable::_bc_configs
private

Definition at line 133 of file ProcessVariable.h.

Referenced by ProcessVariable(), and createBoundaryConditions().

◆ _deactivated_subdomains

std::vector<std::unique_ptr<DeactivatedSubdomain const> > ProcessLib::ProcessVariable::_deactivated_subdomains
private

◆ _ids_of_active_elements

std::vector<std::size_t> ProcessLib::ProcessVariable::_ids_of_active_elements
mutableprivate

IDs of the active elements. It is initialized only if there are deactivated subdomains.

Definition at line 123 of file ProcessVariable.h.

Referenced by createBoundaryConditionsForDeactivatedSubDomains(), getActiveElementIDs(), and updateDeactivatedSubdomains().

◆ _initial_condition

ParameterLib::Parameter<double> const& ProcessLib::ProcessVariable::_initial_condition
private

Definition at line 131 of file ProcessVariable.h.

Referenced by getInitialCondition().

◆ _mesh

MeshLib::Mesh& ProcessLib::ProcessVariable::_mesh
private

◆ _n_components

const int ProcessLib::ProcessVariable::_n_components
private

Definition at line 102 of file ProcessVariable.h.

Referenced by ProcessVariable(), and getNumberOfGlobalComponents().

◆ _name

std::string const ProcessLib::ProcessVariable::_name
private

Definition at line 100 of file ProcessVariable.h.

Referenced by ProcessVariable(), and getName().

◆ _shapefunction_order

unsigned ProcessLib::ProcessVariable::_shapefunction_order
private

The polynomial order of the process variable's shape functions.

Requires an appropriate mesh.

The order of the shape functions can not be higher than the maximum available order for the underlying geometric elements. For example the second order shape functions for a hexahedron are only possible if the geometric element is at least a 20-node hexahedron element (MeshLib::TemplateElement<MeshLib::HexRule20>), whereas linear shape functions are also available on the 8-node hexahedron (MeshLib::TemplateElement<MeshLib::HexRule8>).

See also
MeshLib::CellRule MeshLib::FaceRule MeshLib::EdgeRule.

Definition at line 116 of file ProcessVariable.h.

Referenced by ProcessVariable(), createBoundaryConditions(), createSourceTerms(), and getShapeFunctionOrder().

◆ _source_term_configs

std::vector<SourceTermConfig> ProcessLib::ProcessVariable::_source_term_configs
private

Definition at line 134 of file ProcessVariable.h.

Referenced by ProcessVariable(), and createSourceTerms().


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