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 54 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 &&)
std::string const & getName () const
MeshLib::Mesh const & getMesh () const
 Returns a mesh on which the process variable is defined.
std::vector< DeactivatedSubdomain > 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.
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::reference_wrapper< ProcessVariable > > const &all_process_variables_for_this_process, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
std::vector< std::unique_ptr< SourceTermBase > > createSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int variable_id, unsigned const integration_order, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, std::vector< std::reference_wrapper< ProcessVariable > > const &all_process_variables_for_this_process, const MeshLib::Mesh &bulk_mesh)
ParameterLib::Parameter< double > const & getInitialCondition () const
unsigned getShapeFunctionOrder () const
bool compensateNonEquilibriumInitialResiduum () const
 ~ProcessVariable ()

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< DeactivatedSubdomain_deactivated_subdomains
std::vector< std::size_t > _ids_of_active_elements
MeshLib::PropertyVector< unsigned char > * _is_active = nullptr
ParameterLib::Parameter< double > const & _initial_condition
std::vector< BoundaryConditionConfig_bc_configs
std::vector< SourceTermConfig_source_term_configs
const bool _compensate_non_equilibrium_initial_residuum

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
Input File Parameter
prj__process_variables__process_variable__source_terms__source_term__type

Definition at line 105 of file ProcessVariable.cpp.

112 :
113 _name(config.getConfigParameter<std::string>("name")),
114 _mesh(mesh),
116 _n_components(config.getConfigParameter<int>("components")),
118 _shapefunction_order(config.getConfigParameter<unsigned>("order")),
120 createDeactivatedSubdomains(config, mesh, parameters, curves)),
123 config.getConfigParameter<std::string>("initial_condition"),
124 parameters, _n_components, &mesh)),
127 config.getConfigParameter<bool>(
128 "compensate_non_equilibrium_initial_residuum", false))
129{
130 DBUG("Constructing process variable {:s}", _name);
131
133 {
134 OGS_FATAL("The given shape function order {:d} is not supported",
136 }
137
138 // Boundary conditions
139 if (auto bcs_config =
141 config.getConfigSubtreeOptional("boundary_conditions"))
142 {
143 for (
144 auto bc_config :
146 bcs_config->getConfigSubtreeList("boundary_condition"))
147 {
148 auto const bc_meshes = findMeshInConfig(bc_config, meshes);
149 auto component_id =
151 bc_config.getConfigParameterOptional<int>("component");
152
153 if (!component_id && _n_components == 1)
154 {
155 // default value for single component vars.
156 component_id = 0;
157 }
158
159 _bc_configs.emplace_back(
160 std::move(bc_config), bc_meshes, component_id,
162 }
163 }
164 else
165 {
166 INFO("No boundary conditions for process variable '{:s}' found.",
167 _name);
168 }
169
170 // Source terms
172 if (auto sts_config = config.getConfigSubtreeOptional("source_terms"))
173 {
174 for (
175 auto st_config :
177 sts_config->getConfigSubtreeList("source_term"))
178 {
179 auto st_meshes = findMeshInConfig(st_config, meshes);
180 auto component_id =
182 st_config.getConfigParameterOptional<int>("component");
183 auto const type =
185 st_config.peekConfigParameter<std::string>("type");
186
187 if (!component_id)
188 {
189 if (_n_components == 1)
190 {
191 // default value for single component vars.
192 component_id = 0;
193 }
194 else if (type == "Anchor")
195 {
196 // dummy value
197 component_id = 0;
198 }
199 else if (type == "EmbeddedAnchor")
200 {
201 // dummy value
202 component_id = 0;
203 }
204 else
205 {
206 OGS_FATAL(
207 "Specifying the component id (<component>) for a "
208 "source term for a non-scalar process variable is "
209 "mandatory.");
210 }
211 }
212
213 _source_term_configs.emplace_back(std::move(st_config),
214 st_meshes.front(), *component_id);
215 }
216 }
217 else
218 {
219 INFO("No source terms for process variable '{:s}' found.", _name);
220 }
221
222 if (!_deactivated_subdomains.empty())
223 {
224 _is_active = getOrCreateMeshProperty<unsigned char>(
225 _mesh, _name + "_active", MeshLib::MeshItemType::Cell, 1);
226 std::fill(std::begin(*_is_active), std::end(*_is_active), 1u);
227 }
228}
#define OGS_FATAL(...)
Definition Error.h:19
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
ParameterLib::Parameter< double > const & _initial_condition
std::vector< BoundaryConditionConfig > _bc_configs
std::vector< DeactivatedSubdomain > _deactivated_subdomains
const bool _compensate_non_equilibrium_initial_residuum
MeshLib::PropertyVector< unsigned char > * _is_active
std::vector< SourceTermConfig > _source_term_configs
OGS_NO_DANGLING Parameter< ParameterDataType > & findParameter(std::string const &parameter_name, std::vector< std::unique_ptr< ParameterBase > > const &parameters, int const num_components, MeshLib::Mesh const *const mesh=nullptr)
std::vector< DeactivatedSubdomain > 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)
std::vector< std::reference_wrapper< const MeshLib::Mesh > > findMeshInConfig(BaseLib::ConfigTree const &config, std::vector< std::unique_ptr< MeshLib::Mesh > > const &meshes)

References _bc_configs, _compensate_non_equilibrium_initial_residuum, _deactivated_subdomains, _initial_condition, _is_active, _mesh, _n_components, _name, _shapefunction_order, _source_term_configs, MeshLib::Cell, ProcessLib::createDeactivatedSubdomains(), DBUG(), BaseLib::ConfigTree::getConfigSubtreeOptional(), INFO(), and OGS_FATAL.

Referenced by ProcessVariable().

◆ ProcessVariable() [2/2]

ProcessLib::ProcessVariable::ProcessVariable ( ProcessVariable && )
default

References ProcessVariable().

◆ ~ProcessVariable()

ProcessLib::ProcessVariable::~ProcessVariable ( )
default

Member Function Documentation

◆ compensateNonEquilibriumInitialResiduum()

bool ProcessLib::ProcessVariable::compensateNonEquilibriumInitialResiduum ( ) const
inline

◆ 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,
std::vector< std::reference_wrapper< ProcessVariable > > const & all_process_variables_for_this_process,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )

Definition at line 238 of file ProcessVariable.cpp.

247{
248 std::vector<std::unique_ptr<BoundaryCondition>> bcs;
249 bcs.reserve(_bc_configs.size());
250
251 for (auto const& config : _bc_configs)
252 {
253 auto bc_vec = createBoundaryCondition(
254 config, dof_table, _mesh, variable_id, integration_order,
255 _shapefunction_order, parameters, process,
256 all_process_variables_for_this_process, media);
257 for (auto& bc : bc_vec)
258 {
259 // bc can be a nullptr if in parallel execution the bc isn't
260 // associated to the particular subdomain
261 if (bc != nullptr)
262 {
263 bcs.push_back(std::move(bc));
264 }
265 }
266 }
267
269 parameters, bcs);
270
271 // Clear configs to trigger ConfigTree destruction and validation of
272 // unprocessed parameters before the simulation starts.
273 _bc_configs.clear();
274
275 return bcs;
276}
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::vector< 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, std::vector< std::reference_wrapper< ProcessVariable > > const &all_process_variables_for_this_process, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)

References _bc_configs, _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 278 of file ProcessVariable.cpp.

282{
283 for (auto const& deactivated_subdomain : _deactivated_subdomains)
284 {
285 auto const& deactivated_subdomain_mesh =
286 deactivated_subdomain.deactivated_subdomain_mesh;
287 auto const* parameter = &ParameterLib::findParameter<double>(
289 bool const set_outer_nodes_dirichlet_values =
290 deactivated_subdomain.boundary_value_parameter != nullptr;
291 if (set_outer_nodes_dirichlet_values)
292 {
293 parameter = deactivated_subdomain.boundary_value_parameter;
294 }
295
296 for (int component_id = 0;
297 component_id <
298 dof_table.getNumberOfVariableComponents(variable_id);
299 component_id++)
300 {
301 auto bc = std::make_unique<DeactivatedSubdomainDirichlet>(
302 *_is_active, deactivated_subdomain.time_interval, *parameter,
303 set_outer_nodes_dirichlet_values, deactivated_subdomain_mesh,
304 dof_table, variable_id, component_id);
305
306#ifdef USE_PETSC
307 // TODO: make it work under PETSc too.
308 if (bc == nullptr)
309 {
310 continue;
311 }
312#endif // USE_PETSC
313 bcs.push_back(std::move(bc));
314 }
315 }
316}
int getNumberOfVariableComponents(int variable_id) const
static PROCESSLIB_EXPORT const std::string zero_parameter_name

References _deactivated_subdomains, _is_active, ParameterLib::findParameter(), NumLib::LocalToGlobalIndexMap::getNumberOfVariableComponents(), and ProcessLib::DeactivatedSubdomain::zero_parameter_name.

Referenced by createBoundaryConditions().

◆ createSourceTerms()

std::vector< std::unique_ptr< SourceTermBase > > 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,
std::vector< std::reference_wrapper< ProcessVariable > > const & all_process_variables_for_this_process,
const MeshLib::Mesh & bulk_mesh )

Definition at line 374 of file ProcessVariable.cpp.

381{
382 std::vector<std::unique_ptr<SourceTermBase>> source_terms;
383
384 transform(cbegin(_source_term_configs), cend(_source_term_configs),
385 back_inserter(source_terms),
386 [&](auto const& config)
387 {
388 return createSourceTerm(
389 config, dof_table, config.mesh, variable_id,
390 integration_order, _shapefunction_order, parameters,
391 all_process_variables_for_this_process, bulk_mesh);
392 });
393
394 // Clear configs to trigger ConfigTree destruction and validation of
395 // unprocessed parameters before the simulation starts.
396 _source_term_configs.clear();
397
398 return source_terms;
399}
std::unique_ptr< SourceTermBase > 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, std::vector< std::reference_wrapper< ProcessVariable > > const &all_process_variables_for_this_process, const MeshLib::Mesh &bulk_mesh)

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 80 of file ProcessVariable.h.

81 {
83 }
std::vector< std::size_t > _ids_of_active_elements

References _ids_of_active_elements.

◆ getDeactivatedSubdomains()

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

Definition at line 73 of file ProcessVariable.h.

74 {
76 }

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 232 of file ProcessVariable.cpp.

233{
234 return _mesh;
235}

References _mesh.

◆ getName()

◆ getNumberOfGlobalComponents()

◆ getShapeFunctionOrder()

◆ updateDeactivatedSubdomains()

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

Definition at line 318 of file ProcessVariable.cpp.

319{
320 if (_deactivated_subdomains.empty())
321 {
322 return;
323 }
324
326
327 // If none of the deactivated subdomains is active at current time, then the
328 // _ids_of_active_elements remain empty.
329 if (std::none_of(begin(_deactivated_subdomains),
330 end(_deactivated_subdomains), [&](auto const& ds)
331 { return ds.isInTimeSupportInterval(time); }))
332 {
333 // Also mark all of the elements as active.
334 assert(_is_active != nullptr); // guaranteed by constructor
335 std::fill(std::begin(*_is_active), std::end(*_is_active), 1u);
336
337 return;
338 }
339
340 auto is_active_in_subdomain = [&](std::size_t const element_id,
341 DeactivatedSubdomain const& ds) -> bool
342 {
343 return (!ds.isInTimeSupportInterval(time)) ||
344 !ds.isDeactivated(*_mesh.getElement(element_id), time);
345 };
346
347 auto is_active_in_all_subdomains = [&](std::size_t const element_id) -> bool
348 {
349 return std::all_of(begin(_deactivated_subdomains),
351 [&](auto const& ds)
352 { return is_active_in_subdomain(element_id, ds); });
353 };
354
355 auto const number_of_elements = _mesh.getNumberOfElements();
356 for (std::size_t element_id = 0; element_id < number_of_elements;
357 element_id++)
358 {
359 if (is_active_in_all_subdomains(element_id))
360 {
361 _ids_of_active_elements.push_back(element_id);
362 }
363 }
364
365 // all elements are deactivated
366 std::fill(std::begin(*_is_active), std::end(*_is_active), 0u);
367
368 for (auto const id : _ids_of_active_elements)
369 {
370 (*_is_active)[id] = 1u;
371 }
372}

References _deactivated_subdomains, _ids_of_active_elements, _is_active, and _mesh.

Member Data Documentation

◆ _bc_configs

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

Definition at line 156 of file ProcessVariable.h.

Referenced by ProcessVariable(), and createBoundaryConditions().

◆ _compensate_non_equilibrium_initial_residuum

const bool ProcessLib::ProcessVariable::_compensate_non_equilibrium_initial_residuum
private

Definition at line 159 of file ProcessVariable.h.

Referenced by ProcessVariable(), and compensateNonEquilibriumInitialResiduum().

◆ _deactivated_subdomains

std::vector<DeactivatedSubdomain> 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 145 of file ProcessVariable.h.

Referenced by getActiveElementIDs(), and updateDeactivatedSubdomains().

◆ _initial_condition

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

Definition at line 154 of file ProcessVariable.h.

Referenced by ProcessVariable(), and getInitialCondition().

◆ _is_active

MeshLib::PropertyVector<unsigned char>* ProcessLib::ProcessVariable::_is_active = nullptr
private

◆ _mesh

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

◆ _n_components

const int ProcessLib::ProcessVariable::_n_components
private

Definition at line 125 of file ProcessVariable.h.

Referenced by ProcessVariable(), and getNumberOfGlobalComponents().

◆ _name

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

Definition at line 123 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 139 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 157 of file ProcessVariable.h.

Referenced by ProcessVariable(), and createSourceTerms().


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