OGS
ProcessLib::Process Class Referenceabstract

Detailed Description

Definition at line 39 of file Process.h.

#include <Process.h>

Inheritance diagram for ProcessLib::Process:
[legend]
Collaboration diagram for ProcessLib::Process:
[legend]

Public Member Functions

 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
 
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable (const int) const
 
MeshLib::MeshgetMesh () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::string > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Public Attributes

std::string const name
 

Static Public Attributes

static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 

Protected Member Functions

NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
virtual void constructDofTable ()
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 

Protected Attributes

MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Private Member Functions

virtual void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order)=0
 Process specific initialization called by initialize().
 
virtual void initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
virtual void setInitialConditionsConcreteProcess (std::vector< GlobalVector * > &, double const, int const)
 
virtual void preAssembleConcreteProcess (const double, double const, GlobalVector const &)
 
virtual void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)=0
 
virtual void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)=0
 
virtual void preTimestepConcreteProcess (std::vector< GlobalVector * > const &, const double, const double, const int)
 
virtual void postTimestepConcreteProcess (std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, const double, int const)
 
virtual void postNonLinearSolverConcreteProcess (std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, double const, int const)
 
virtual void preIterationConcreteProcess (const unsigned, GlobalVector const &)
 
virtual void computeSecondaryVariableConcrete (double const, double const, std::vector< GlobalVector * > const &, GlobalVector const &, int const)
 
virtual NumLib::IterationResult postIterationConcreteProcess (GlobalVector const &)
 
virtual void preOutputConcreteProcess (const double, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)
 
void initializeExtrapolator ()
 
void computeSparsityPattern ()
 

Private Attributes

std::vector< SourceTermCollection_source_term_collections
 
ExtrapolatorData _extrapolator_data
 
std::vector< std::size_t > _ids_of_active_elements
 Union of active element ids per process variable.
 

Constructor & Destructor Documentation

◆ Process()

ProcessLib::Process::Process ( std::string name_,
MeshLib::Mesh & mesh,
std::unique_ptr< AbstractJacobianAssembler > && jacobian_assembler,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
unsigned const integration_order,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && process_variables,
SecondaryVariableCollection && secondary_variables,
const bool use_monolithic_scheme = true )

Definition at line 28 of file Process.cpp.

38 : name(std::move(name_)),
39 _mesh(mesh),
40 _secondary_variables(std::move(secondary_variables)),
41 _jacobian_assembler(std::move(jacobian_assembler)),
43 _use_monolithic_scheme(use_monolithic_scheme),
44 _integration_order(integration_order),
45 _process_variables(std::move(process_variables)),
47 [&](const std::size_t number_of_process_variables)
48 -> std::vector<BoundaryConditionCollection>
49 {
50 std::vector<BoundaryConditionCollection> pcs_BCs;
51 pcs_BCs.reserve(number_of_process_variables);
52 for (std::size_t i = 0; i < number_of_process_variables; i++)
53 {
54 pcs_BCs.emplace_back(BoundaryConditionCollection(parameters));
55 }
56 return pcs_BCs;
57 }(_process_variables.size())),
59 [&](const std::size_t number_of_processes)
60 -> std::vector<SourceTermCollection>
61 {
62 std::vector<SourceTermCollection> pcs_sts;
63 pcs_sts.reserve(number_of_processes);
64 for (std::size_t i = 0; i < number_of_processes; i++)
65 {
66 pcs_sts.emplace_back(SourceTermCollection(parameters));
67 }
68 return pcs_sts;
69 }(_process_variables.size()))
70{
71}
std::string const name
Definition Process.h:351
std::vector< BoundaryConditionCollection > _boundary_conditions
Definition Process.h:392
MeshLib::Mesh & _mesh
Definition Process.h:354
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
Definition Process.h:387
SecondaryVariableCollection _secondary_variables
Definition Process.h:359
VectorMatrixAssembler _global_assembler
Definition Process.h:364
unsigned const _integration_order
Definition Process.h:371
std::vector< SourceTermCollection > _source_term_collections
Definition Process.h:398
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:363
const bool _use_monolithic_scheme
Definition Process.h:366

Member Function Documentation

◆ assemble()

void ProcessLib::Process::assemble ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
final

Definition at line 249 of file Process.cpp.

254{
255 assert(x.size() == x_prev.size());
256 for (std::size_t i = 0; i < x.size(); i++)
257 {
260 }
261
262 assembleConcreteProcess(t, dt, x, x_prev, process_id, M, K, b);
263
264 // the last argument is for the jacobian, nullptr is for a unused jacobian
265 _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, K, b,
266 nullptr);
267
268 // the last argument is for the jacobian, nullptr is for a unused jacobian
269 _source_term_collections[process_id].integrate(t, *x[process_id], b,
270 nullptr);
271}
virtual void assembleConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)=0
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27

References _boundary_conditions, _source_term_collections, assembleConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ assembleConcreteProcess()

virtual void ProcessLib::Process::assembleConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
privatepure virtual

Implemented in ProcessLib::ComponentTransport::ComponentTransportProcess, ProcessLib::HeatTransportBHE::HeatTransportBHEProcess, ProcessLib::HT::HTProcess, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >, ProcessLib::LiquidFlow::LiquidFlowProcess, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >, ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess, ProcessLib::RichardsFlow::RichardsFlowProcess, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion, ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >, ProcessLib::TES::TESProcess, ProcessLib::TH2M::TH2MProcess< DisplacementDim >, ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >, ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess, ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >, ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess, ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess, ProcessLib::HeatConduction::HeatConductionProcess, ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >, ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >, and ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >.

Referenced by assemble().

◆ assembleWithJacobian()

void ProcessLib::Process::assembleWithJacobian ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
GlobalMatrix & Jac )
final

Definition at line 273 of file Process.cpp.

279{
280 assert(x.size() == x_prev.size());
281 for (std::size_t i = 0; i < x.size(); i++)
282 {
285 }
286
287 assembleWithJacobianConcreteProcess(t, dt, x, x_prev, process_id, M, K, b,
288 Jac);
289
290 // TODO: apply BCs to Jacobian.
291 _boundary_conditions[process_id].applyNaturalBC(t, x, process_id, K, b,
292 &Jac);
293
294 // the last argument is for the jacobian, nullptr is for a unused jacobian
295 _source_term_collections[process_id].integrate(t, *x[process_id], b, &Jac);
296}
virtual void assembleWithJacobianConcreteProcess(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)=0

References _boundary_conditions, _source_term_collections, assembleWithJacobianConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ assembleWithJacobianConcreteProcess()

virtual void ProcessLib::Process::assembleWithJacobianConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
GlobalMatrix & Jac )
privatepure virtual

Implemented in ProcessLib::ComponentTransport::ComponentTransportProcess, ProcessLib::HeatTransportBHE::HeatTransportBHEProcess, ProcessLib::HT::HTProcess, ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >, ProcessLib::LiquidFlow::LiquidFlowProcess, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >, ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess, ProcessLib::RichardsFlow::RichardsFlowProcess, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >, ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion, ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >, ProcessLib::TES::TESProcess, ProcessLib::TH2M::TH2MProcess< DisplacementDim >, ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >, ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess, ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >, ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess, ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess, ProcessLib::HeatConduction::HeatConductionProcess, and ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >.

Referenced by assembleWithJacobian().

◆ computeSecondaryVariable()

void ProcessLib::Process::computeSecondaryVariable ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
GlobalVector const & x_prev,
int const process_id )

compute secondary variables for the coupled equations or for output.

Definition at line 470 of file Process.cpp.

475{
476 for (auto const* solution : x)
479
480 computeSecondaryVariableConcrete(t, dt, x, x_prev, process_id);
481}
virtual void computeSecondaryVariableConcrete(double const, double const, std::vector< GlobalVector * > const &, GlobalVector const &, int const)
Definition Process.h:282

References computeSecondaryVariableConcrete(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ computeSecondaryVariableConcrete()

◆ computeSparsityPattern()

void ProcessLib::Process::computeSparsityPattern ( )
private

Computes and stores global matrix' sparsity pattern from given DOF-table.

Definition at line 418 of file Process.cpp.

419{
422}
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:357
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:379
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References _local_to_global_index_map, _mesh, _sparsity_pattern, and NumLib::computeSparsityPattern().

Referenced by initialize().

◆ constructDofTable()

void ProcessLib::Process::constructDofTable ( )
protectedvirtual

This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map. Other cases can be considered by overloading this member function in the derived class.

Reimplemented in ProcessLib::HeatTransportBHE::HeatTransportBHEProcess, ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >, ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >, ProcessLib::TH2M::TH2MProcess< DisplacementDim >, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >, and ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >.

Definition at line 298 of file Process.cpp.

299{
301 {
303
304 return;
305 }
306
307 // For staggered scheme:
308 const int specified_process_id = 0;
310}
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:345
void constructMonolithicProcessDofTable()
Definition Process.cpp:312

References _use_monolithic_scheme, constructDofTableOfSpecifiedProcessStaggeredScheme(), and constructMonolithicProcessDofTable().

Referenced by initialize().

◆ constructDofTableOfSpecifiedProcessStaggeredScheme()

void ProcessLib::Process::constructDofTableOfSpecifiedProcessStaggeredScheme ( const int specified_process_id)
protected

Construct the DOF table for a specified process in the staggered scheme, which is stored in the member of this class, _local_to_global_index_map.

Definition at line 345 of file Process.cpp.

347{
348 // Create single component dof in every of the mesh nodes.
350 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
351
352 // Vector of mesh subsets.
353 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
354
355 // Vector of the number of variable components
356 std::vector<int> vec_var_n_components;
357 // Collect the mesh subsets in a vector for each variables' components.
358 std::generate_n(std::back_inserter(all_mesh_subsets),
359 _process_variables[specified_process_id][0]
360 .get()
361 .getNumberOfGlobalComponents(),
362 [&]() { return *_mesh_subset_all_nodes; });
363
364 // Create a vector of the number of variable components.
365 vec_var_n_components.push_back(_process_variables[specified_process_id][0]
366 .get()
367 .getNumberOfGlobalComponents());
369 std::make_unique<NumLib::LocalToGlobalIndexMap>(
370 std::move(all_mesh_subsets), vec_var_n_components,
372
374}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:355
@ BY_LOCATION
Ordering data by spatial location.
auto & get(Tuples &... ts)
Definition Get.h:67

References _local_to_global_index_map, _mesh, _mesh_subset_all_nodes, _process_variables, NumLib::BY_LOCATION, and MeshLib::Mesh::getNodes().

Referenced by constructDofTable().

◆ constructMonolithicProcessDofTable()

void ProcessLib::Process::constructMonolithicProcessDofTable ( )
protected

Construct the DOF table for the monolithic scheme, which is stored in the member of this class, _local_to_global_index_map.

Definition at line 312 of file Process.cpp.

313{
314 // Create single component dof in every of the mesh nodes.
316 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
317
318 // Vector of mesh subsets.
319 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
320
321 // Collect the mesh subsets in a vector for the components of each
322 // variables.
323 for (ProcessVariable const& pv : _process_variables[0])
324 {
325 std::generate_n(std::back_inserter(all_mesh_subsets),
326 pv.getNumberOfGlobalComponents(),
327 [&]() { return *_mesh_subset_all_nodes; });
328 }
329
330 // Create a vector of the number of variable components
331 std::vector<int> vec_var_n_components;
332 transform(cbegin(_process_variables[0]), cend(_process_variables[0]),
333 back_inserter(vec_var_n_components),
334 [](ProcessVariable const& pv)
335 { return pv.getNumberOfGlobalComponents(); });
336
338 std::make_unique<NumLib::LocalToGlobalIndexMap>(
339 std::move(all_mesh_subsets), vec_var_n_components,
341
343}

References _local_to_global_index_map, _mesh, _mesh_subset_all_nodes, _process_variables, NumLib::BY_LOCATION, MeshLib::Mesh::getNodes(), and ProcessLib::ProcessVariable::getNumberOfGlobalComponents().

Referenced by constructDofTable().

◆ extrapolateIntegrationPointValuesToNodes()

virtual void ProcessLib::Process::extrapolateIntegrationPointValuesToNodes ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< GlobalVector * > &  )
inlinevirtual

Definition at line 106 of file Process.h.

110 {
111 }

◆ getActiveElementIDs()

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

Definition at line 160 of file Process.h.

161 {
163 }
std::vector< std::size_t > _ids_of_active_elements
Union of active element ids per process variable.
Definition Process.h:403

References _ids_of_active_elements.

Referenced by updateDeactivatedSubdomains().

◆ getDOFTable()

◆ getDOFTableForExtrapolatorData()

std::tuple< NumLib::LocalToGlobalIndexMap *, bool > ProcessLib::Process::getDOFTableForExtrapolatorData ( ) const
protectedvirtual

Get the address of a LocalToGlobalIndexMap, and the status of its memory. If the LocalToGlobalIndexMap is created as new in this function, the function also returns a true boolean value to let Extrapolator manage the memory by the address returned by this function.

Returns
Address of a LocalToGlobalIndexMap and its memory status.

Reimplemented in ProcessLib::HT::HTProcess, ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >, ProcessLib::TH2M::TH2MProcess< DisplacementDim >, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >, and ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >.

Definition at line 377 of file Process.cpp.

378{
379 if (_local_to_global_index_map->getNumberOfGlobalComponents() == 1)
380 {
381 // For single-variable-single-component processes reuse the existing DOF
382 // table.
383 const bool manage_storage = false;
384 return std::make_tuple(_local_to_global_index_map.get(),
385 manage_storage);
386 }
387
388 // Otherwise construct a new DOF table.
389 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component;
390 all_mesh_subsets_single_component.emplace_back(*_mesh_subset_all_nodes);
391
392 const bool manage_storage = true;
393
394 return std::make_tuple(new NumLib::LocalToGlobalIndexMap(
395 std::move(all_mesh_subsets_single_component),
396 // by location order is needed for output
398 manage_storage);
399}

References _local_to_global_index_map, _mesh_subset_all_nodes, and NumLib::BY_LOCATION.

Referenced by initializeExtrapolator().

◆ getExtrapolator()

NumLib::Extrapolator & ProcessLib::Process::getExtrapolator ( ) const
inlineprotected

◆ getFlux()

virtual Eigen::Vector3d ProcessLib::Process::getFlux ( std::size_t ,
MathLib::Point3d const & ,
double const ,
std::vector< GlobalVector * > const &  ) const
inlinevirtual

◆ getIndicesOfResiduumWithoutInitialCompensation()

std::vector< GlobalIndexType > ProcessLib::Process::getIndicesOfResiduumWithoutInitialCompensation ( ) const
overrideprotected
Returns
The global indices for the entries of the global residuum vector that do not need initial non-equilibrium compensation.

Definition at line 507 of file Process.cpp.

508{
509 std::vector<GlobalIndexType> indices;
510
511 for (std::size_t process_id = 0; process_id < _process_variables.size();
512 process_id++)
513 {
514 auto const& dof_table_of_process = getDOFTable(process_id);
515
516 auto const& per_process_variables = _process_variables[process_id];
517 for (std::size_t variable_id = 0;
518 variable_id < per_process_variables.size();
519 variable_id++)
520 {
521 auto const& pv = per_process_variables[variable_id];
522 DBUG("Set the initial condition of variable {:s} of process {:d}.",
523 pv.get().getName().data(), process_id);
524
525 if ((pv.get().compensateNonEquilibriumInitialResiduum()))
526 {
527 continue;
528 }
529
530 auto const num_comp = pv.get().getNumberOfGlobalComponents();
531
532 for (int component_id = 0; component_id < num_comp; ++component_id)
533 {
534 auto const& mesh_subset = dof_table_of_process.getMeshSubset(
535 variable_id, component_id);
536 auto const mesh_id = mesh_subset.getMeshID();
537 for (auto const* node : mesh_subset.getNodes())
538 {
539 MeshLib::Location const l(
540 mesh_id, MeshLib::MeshItemType::Node, node->getID());
541
542 auto global_index =
543 std::abs(dof_table_of_process.getGlobalIndex(
544 l, variable_id, component_id));
545
546 indices.push_back(global_index);
547 }
548 }
549 }
550 }
551
552 return indices;
553}
std::vector< std::size_t > getNodes(GeoLib::Point const &pnt, std::vector< MeshLib::Node * > const &nodes, MeshLib::PropertyVector< int > const &mat_ids, std::pair< int, int > const &mat_limits, std::pair< double, double > const &elevation_limits, MeshLib::Mesh const &mesh)
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
virtual NumLib::LocalToGlobalIndexMap const & getDOFTable(const int) const
Definition Process.h:147

References _process_variables, DBUG(), getDOFTable(), and MeshLib::Node.

◆ getIntegrationPointWriters()

std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & ProcessLib::Process::getIntegrationPointWriters ( ) const
inline

Definition at line 171 of file Process.h.

172 {
174 }
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:377

References _integration_point_writer.

Referenced by anonymous_namespace{ProcessOutputData.cpp}::getIntegrationPointWriters().

◆ getKnownSolutions()

std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * ProcessLib::Process::getKnownSolutions ( double const t,
GlobalVector const & x,
int const process_id ) const
inlinefinal

Definition at line 141 of file Process.h.

143 {
144 return _boundary_conditions[process_id].getKnownSolutions(t, x);
145 }

References _boundary_conditions.

◆ getMatrixSpecifications()

◆ getMesh()

◆ getProcessVariables()

std::vector< std::reference_wrapper< ProcessVariable > > const & ProcessLib::Process::getProcessVariables ( const int process_id) const
inline

Definition at line 155 of file Process.h.

156 {
157 return _process_variables[process_id];
158 }

References _process_variables.

Referenced by ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >::HydroMechanicsProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleConcreteProcess(), ProcessLib::TES::TESProcess::assembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess(), ProcessLib::AssemblyMixin< Process >::assembleGeneric(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleWithJacobianConcreteProcess(), ProcessLib::TES::TESProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::computeSecondaryVariableConcrete(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::computeSecondaryVariableConcrete(), ProcessLib::LiquidFlow::LiquidFlowProcess::computeSecondaryVariableConcrete(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::ComponentTransportProcess::computeSecondaryVariableConcrete(), ProcessLib::createConstraintDirichletBoundaryCondition(), anonymous_namespace{ProcessOutputData.cpp}::getProcessVariablesOfAllProcesses(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::initializeConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::postTimestepConcreteProcess(), ProcessLib::HT::HTProcess::postTimestepConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::postTimestepConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::postTimestepConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::postTimestepConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::preTimestepConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation(), ProcessLib::AssemblyMixin< Process >::updateActiveElements(), and updateDeactivatedSubdomains().

◆ getSecondaryVariables()

SecondaryVariableCollection const & ProcessLib::Process::getSecondaryVariables ( ) const
inline

Definition at line 165 of file Process.h.

166 {
168 }

References _secondary_variables.

Referenced by ProcessLib::createProcessOutputData().

◆ getSingleComponentDOFTable()

NumLib::LocalToGlobalIndexMap const & ProcessLib::Process::getSingleComponentDOFTable ( ) const
inlineprotected

◆ initialize()

void ProcessLib::Process::initialize ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)

Definition at line 102 of file Process.cpp.

104{
105 DBUG("Initialize process.");
106
107 DBUG("Construct dof mappings.");
109
110 DBUG("Compute sparsity pattern");
112
113 DBUG("Initialize the extrapolator");
115
118
119 DBUG("Initialize boundary conditions.");
121}
virtual void initializeBoundaryConditions(std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:88
virtual void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order)=0
Process specific initialization called by initialize().
void initializeExtrapolator()
Definition Process.cpp:401
virtual void constructDofTable()
Definition Process.cpp:298
void computeSparsityPattern()
Definition Process.cpp:418

References _integration_order, _local_to_global_index_map, _mesh, computeSparsityPattern(), constructDofTable(), DBUG(), initializeBoundaryConditions(), initializeConcreteProcess(), and initializeExtrapolator().

◆ initializeBoundaryConditions()

void ProcessLib::Process::initializeBoundaryConditions ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
privatevirtual

Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().

Reimplemented in ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >, ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >, ProcessLib::TH2M::TH2MProcess< DisplacementDim >, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >, and ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >.

Definition at line 88 of file Process.cpp.

90{
91 // The number of processes is identical to the size of _process_variables,
92 // the vector contains variables for different processes. See the
93 // documentation of _process_variables.
94 const std::size_t number_of_processes = _process_variables.size();
95 for (std::size_t pcs_id = 0; pcs_id < number_of_processes; pcs_id++)
96 {
98 *_local_to_global_index_map, pcs_id, media);
99 }
100}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:73

References _local_to_global_index_map, _process_variables, and initializeProcessBoundaryConditionsAndSourceTerms().

Referenced by initialize().

◆ initializeConcreteProcess()

virtual void ProcessLib::Process::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
privatepure virtual

Process specific initialization called by initialize().

Implemented in ProcessLib::ComponentTransport::ComponentTransportProcess, ProcessLib::HeatConduction::HeatConductionProcess, ProcessLib::HeatTransportBHE::HeatTransportBHEProcess, ProcessLib::HT::HTProcess, ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >, ProcessLib::LargeDeformation::LargeDeformationProcess< DisplacementDim >, ProcessLib::LIE::HydroMechanics::HydroMechanicsProcess< GlobalDim >, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >, ProcessLib::LiquidFlow::LiquidFlowProcess, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >, ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess, ProcessLib::RichardsFlow::RichardsFlowProcess, ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >, ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >, ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion, ProcessLib::StokesFlow::StokesFlowProcess< GlobalDim >, ProcessLib::TES::TESProcess, ProcessLib::TH2M::TH2MProcess< DisplacementDim >, ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess, ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >, ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >, ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess, ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim, ConstitutiveTraits >, ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess, and ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess.

Referenced by initialize().

◆ initializeExtrapolator()

void ProcessLib::Process::initializeExtrapolator ( )
private

Definition at line 401 of file Process.cpp.

402{
403 NumLib::LocalToGlobalIndexMap* dof_table_single_component;
404 bool manage_storage;
405
406 std::tie(dof_table_single_component, manage_storage) =
408
409 std::unique_ptr<NumLib::Extrapolator> extrapolator(
411 *dof_table_single_component));
412
413 // TODO: Later on the DOF table can change during the simulation!
414 _extrapolator_data = ExtrapolatorData(
415 std::move(extrapolator), dof_table_single_component, manage_storage);
416}
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData() const
Definition Process.cpp:377

References _extrapolator_data, and getDOFTableForExtrapolatorData().

Referenced by initialize().

◆ initializeProcessBoundaryConditionsAndSourceTerms()

void ProcessLib::Process::initializeProcessBoundaryConditionsAndSourceTerms ( const NumLib::LocalToGlobalIndexMap & dof_table,
const int process_id,
std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media )
protected

Initialize the boundary conditions for a single process or coupled processes modelled by the monolithic scheme. It is called by initializeBoundaryConditions().

Definition at line 73 of file Process.cpp.

76{
77 auto const& per_process_variables = _process_variables[process_id];
78 auto& per_process_BCs = _boundary_conditions[process_id];
79
80 per_process_BCs.addBCsForProcessVariables(per_process_variables, dof_table,
81 _integration_order, *this, media);
82
83 auto& per_process_sts = _source_term_collections[process_id];
84 per_process_sts.addSourceTermsForProcessVariables(
85 per_process_variables, dof_table, _integration_order);
86}

References _boundary_conditions, _integration_order, _process_variables, and _source_term_collections.

Referenced by initializeBoundaryConditions().

◆ isMonolithicSchemeUsed()

virtual bool ProcessLib::Process::isMonolithicSchemeUsed ( ) const
inlinevirtual

◆ postIteration()

NumLib::IterationResult ProcessLib::Process::postIteration ( GlobalVector const & x)
final

Definition at line 489 of file Process.cpp.

490{
493}
virtual NumLib::IterationResult postIterationConcreteProcess(GlobalVector const &)
Definition Process.h:291

References postIterationConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ postIterationConcreteProcess()

virtual NumLib::IterationResult ProcessLib::Process::postIterationConcreteProcess ( GlobalVector const & )
inlineprivatevirtual

◆ postNonLinearSolver()

void ProcessLib::Process::postNonLinearSolver ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
double const dt,
int const process_id )

Calculates secondary variables, e.g. stress and strain for deformation analysis, only after nonlinear solver being successfully conducted.

Definition at line 453 of file Process.cpp.

457{
458 for (auto* const solution : x)
459 {
461 }
462 for (auto* const solution : x_prev)
463 {
465 }
466
467 postNonLinearSolverConcreteProcess(x, x_prev, t, dt, process_id);
468}
virtual void postNonLinearSolverConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, double const, int const)
Definition Process.h:270

References postNonLinearSolverConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ postNonLinearSolverConcreteProcess()

virtual void ProcessLib::Process::postNonLinearSolverConcreteProcess ( std::vector< GlobalVector * > const & ,
std::vector< GlobalVector * > const & ,
const double ,
double const ,
int const  )
inlineprivatevirtual

◆ postTimestep()

void ProcessLib::Process::postTimestep ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
const double delta_t,
int const process_id )

Postprocessing after a complete timestep.

Definition at line 434 of file Process.cpp.

438{
439 for (auto* const solution : x)
440 {
442 }
443 for (auto* const solution : x_prev)
444 {
446 }
447
448 postTimestepConcreteProcess(x, x_prev, t, delta_t, process_id);
449
450 _boundary_conditions[process_id].postTimestep(t, x, process_id);
451}
virtual void postTimestepConcreteProcess(std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, const double, const double, int const)
Definition Process.h:261

References _boundary_conditions, postTimestepConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ postTimestepConcreteProcess()

virtual void ProcessLib::Process::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & ,
std::vector< GlobalVector * > const & ,
const double ,
const double ,
int const  )
inlineprivatevirtual

◆ preAssemble()

void ProcessLib::Process::preAssemble ( const double t,
double const dt,
GlobalVector const & x )
final

Definition at line 243 of file Process.cpp.

245{
247}
virtual void preAssembleConcreteProcess(const double, double const, GlobalVector const &)
Definition Process.h:236

References preAssembleConcreteProcess().

◆ preAssembleConcreteProcess()

virtual void ProcessLib::Process::preAssembleConcreteProcess ( const double ,
double const ,
GlobalVector const &  )
inlineprivatevirtual

Reimplemented in ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalProcess< DisplacementDim >.

Definition at line 236 of file Process.h.

239 {
240 }

Referenced by preAssemble().

◆ preIteration()

void ProcessLib::Process::preIteration ( const unsigned iter,
GlobalVector const & x )
final

Definition at line 483 of file Process.cpp.

484{
487}
virtual void preIterationConcreteProcess(const unsigned, GlobalVector const &)
Definition Process.h:277

References preIterationConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ preIterationConcreteProcess()

virtual void ProcessLib::Process::preIterationConcreteProcess ( const unsigned ,
GlobalVector const &  )
inlineprivatevirtual

Reimplemented in ProcessLib::TES::TESProcess.

Definition at line 277 of file Process.h.

279 {
280 }

Referenced by preIteration().

◆ preOutput()

void ProcessLib::Process::preOutput ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id )

Definition at line 495 of file Process.cpp.

499{
500 for (auto const* solution : x)
502
503 preOutputConcreteProcess(t, dt, x, x_prev, process_id);
504}
virtual void preOutputConcreteProcess(const double, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)
Definition Process.h:297

References preOutputConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ preOutputConcreteProcess()

virtual void ProcessLib::Process::preOutputConcreteProcess ( const double ,
double const ,
std::vector< GlobalVector * > const & ,
std::vector< GlobalVector * > const & ,
int const  )
inlineprivatevirtual

◆ preTimestep()

void ProcessLib::Process::preTimestep ( std::vector< GlobalVector * > const & x,
const double t,
const double delta_t,
const int process_id )

Preprocessing before starting assembly for new timestep.

Definition at line 424 of file Process.cpp.

426{
427 for (auto* const solution : x)
429 preTimestepConcreteProcess(x, t, delta_t, process_id);
430
431 _boundary_conditions[process_id].preTimestep(t, x, process_id);
432}
virtual void preTimestepConcreteProcess(std::vector< GlobalVector * > const &, const double, const double, const int)
Definition Process.h:253

References _boundary_conditions, preTimestepConcreteProcess(), and MathLib::LinAlg::setLocalAccessibleVector().

◆ preTimestepConcreteProcess()

◆ setInitialConditions()

void ProcessLib::Process::setInitialConditions ( std::vector< GlobalVector * > & process_solutions,
std::vector< GlobalVector * > const & process_solutions_prev,
double const t,
int const process_id )

Definition at line 123 of file Process.cpp.

128{
129 auto& x = *process_solutions[process_id];
130 auto& x_prev = *process_solutions_prev[process_id];
131
132 // getDOFTableOfProcess can be overloaded by the specific process.
133 auto const& dof_table_of_process = getDOFTable(process_id);
134
135 auto const& per_process_variables = _process_variables[process_id];
136 for (std::size_t variable_id = 0;
137 variable_id < per_process_variables.size();
138 variable_id++)
139 {
142
143 auto const& pv = per_process_variables[variable_id];
144 DBUG("Set the initial condition of variable {:s} of process {:d}.",
145 pv.get().getName().data(), process_id);
146
147 auto const& ic = pv.get().getInitialCondition();
148
149 auto const num_comp = pv.get().getNumberOfGlobalComponents();
150
151 for (int component_id = 0; component_id < num_comp; ++component_id)
152 {
153 auto const& mesh_subset =
154 dof_table_of_process.getMeshSubset(variable_id, component_id);
155 auto const mesh_id = mesh_subset.getMeshID();
156 for (auto const* node : mesh_subset.getNodes())
157 {
159 node->getID());
160
161 pos.setNodeID(node->getID());
162 pos.setCoordinates(*node);
163 auto const& ic_value = ic(t, pos);
164
165 auto global_index =
166 std::abs(dof_table_of_process.getGlobalIndex(l, variable_id,
167 component_id));
168#ifdef USE_PETSC
169 // The global indices of the ghost entries of the global
170 // matrix or the global vectors need to be set as negative
171 // values for equation assembly, however the global indices
172 // start from zero. Therefore, any ghost entry with zero
173 // index is assigned an negative value of the vector size
174 // or the matrix dimension. To assign the initial value for
175 // the ghost entries, the negative indices of the ghost
176 // entries are restored to zero.
177 if (global_index == x.size())
178 global_index = 0;
179#endif
180 x.set(global_index, ic_value[component_id]);
181 }
182 }
183 }
184
186 MathLib::LinAlg::copy(x, x_prev); // pushState
187
190
191 setInitialConditionsConcreteProcess(process_solutions, t, process_id);
192}
void setNodeID(std::size_t node_id)
void setCoordinates(MathLib::Point3d const &coordinates)
virtual void setInitialConditionsConcreteProcess(std::vector< GlobalVector * > &, double const, int const)
Definition Process.h:229
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:170
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37

References _process_variables, MathLib::LinAlg::copy(), DBUG(), MathLib::LinAlg::finalizeAssembly(), getDOFTable(), MeshLib::Node, ParameterLib::SpatialPosition::setCoordinates(), setInitialConditionsConcreteProcess(), MathLib::LinAlg::setLocalAccessibleVector(), and ParameterLib::SpatialPosition::setNodeID().

◆ setInitialConditionsConcreteProcess()

◆ solveReactionEquation()

virtual void ProcessLib::Process::solveReactionEquation ( std::vector< GlobalVector * > & ,
std::vector< GlobalVector * > const & ,
double const ,
double const ,
NumLib::EquationSystem & ,
int const  )
inlinevirtual

Reimplemented in ProcessLib::ComponentTransport::ComponentTransportProcess.

Definition at line 187 of file Process.h.

192 {
193 }

◆ updateDeactivatedSubdomains()

void ProcessLib::Process::updateDeactivatedSubdomains ( double const time,
const int process_id )

Definition at line 202 of file Process.cpp.

204{
205 auto const& variables_per_process = getProcessVariables(process_id);
206 for (auto const& variable : variables_per_process)
207 {
208 variable.get().updateDeactivatedSubdomains(time);
209 }
210
212
213 auto active_elements_ids = ranges::views::transform(
214 [](auto const& variable)
215 { return variable.get().getActiveElementIDs(); });
216
217 // Early return if there's any process variable with all elements active.
218 if (ranges::any_of(variables_per_process | active_elements_ids,
219 [](auto const& vector) { return vector.empty(); }))
220 {
221 return;
222 }
223
224 // Some process variable has deactivated elements. Create union of active
225 // ids.
226
227 _ids_of_active_elements = // there is at least one process variable.
228 variables_per_process[0].get().getActiveElementIDs();
229
230 for (auto const& pv_active_element_ids :
231 variables_per_process | ranges::views::drop(1) | active_elements_ids)
232 {
233 std::vector<std::size_t> new_active_elements;
234 new_active_elements.reserve(_ids_of_active_elements.size() +
235 pv_active_element_ids.size());
236 ranges::set_union(_ids_of_active_elements, pv_active_element_ids,
237 std::back_inserter(new_active_elements));
238
239 _ids_of_active_elements = std::move(new_active_elements);
240 }
241}
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:160
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155

References _ids_of_active_elements, getActiveElementIDs(), and getProcessVariables().

Member Data Documentation

◆ _boundary_conditions

std::vector<BoundaryConditionCollection> ProcessLib::Process::_boundary_conditions
protected

Vector for boundary conditions. For the monolithic scheme or a single process, the size of the vector is one. For the staggered scheme, the size of vector is the number of the coupled processes.

Definition at line 392 of file Process.h.

Referenced by assemble(), assembleWithJacobian(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::createBHEBoundaryConditionTopBottom(), getKnownSolutions(), initializeProcessBoundaryConditionsAndSourceTerms(), postTimestep(), and preTimestep().

◆ _extrapolator_data

ExtrapolatorData ProcessLib::Process::_extrapolator_data
private

Definition at line 400 of file Process.h.

Referenced by getExtrapolator(), getSingleComponentDOFTable(), and initializeExtrapolator().

◆ _global_assembler

VectorMatrixAssembler ProcessLib::Process::_global_assembler
protected

Definition at line 364 of file Process.h.

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::assembleConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleConcreteProcess(), ProcessLib::TES::TESProcess::assembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleWithJacobianConcreteProcess(), ProcessLib::TES::TESProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleWithJacobianConcreteProcess(), and ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess().

◆ _ids_of_active_elements

std::vector<std::size_t> ProcessLib::Process::_ids_of_active_elements
private

Union of active element ids per process variable.

Definition at line 403 of file Process.h.

Referenced by getActiveElementIDs(), and updateDeactivatedSubdomains().

◆ _integration_order

unsigned const ProcessLib::Process::_integration_order
protected

◆ _integration_point_writer

◆ _jacobian_assembler

◆ _local_to_global_index_map

std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::Process::_local_to_global_index_map
protected

Definition at line 357 of file Process.h.

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::assembleConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleConcreteProcess(), ProcessLib::TES::TESProcess::assembleConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess(), ProcessLib::AssemblyMixin< Process >::assembleGeneric(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), ProcessLib::LiquidFlow::LiquidFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsComponentTransport::RichardsComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::RichardsFlow::RichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::assembleWithJacobianConcreteProcess(), ProcessLib::TES::TESProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPP::TwoPhaseFlowWithPPProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TwoPhaseFlowWithPrho::TwoPhaseFlowWithPrhoProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess(), ProcessLib::TES::TESProcess::computeEquilibriumLoading(), ProcessLib::TES::TESProcess::computeRelativeHumidity(), computeSparsityPattern(), ProcessLib::TES::TESProcess::computeVapourPartialPressure(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::constructDofTable(), constructDofTableOfSpecifiedProcessStaggeredScheme(), constructMonolithicProcessDofTable(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::createBHEBoundaryConditionTopBottom(), getDOFTable(), getDOFTableForExtrapolatorData(), ProcessLib::HT::HTProcess::getDOFTableForExtrapolatorData(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::getDOFTables(), ProcessLib::ComponentTransport::ComponentTransportProcess::getFlux(), ProcessLib::LiquidFlow::LiquidFlowProcess::getFlux(), ProcessLib::HT::HTProcess::getFlux(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusion::getFlux(), getMatrixSpecifications(), initialize(), initializeBoundaryConditions(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeConcreteProcess(), ProcessLib::TES::TESProcess::postIterationConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::preOutputConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::preOutputConcreteProcess(), ProcessLib::ThermalTwoPhaseFlowWithPP::ThermalTwoPhaseFlowWithPPProcess::preTimestepConcreteProcess(), and ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::setInitialConditionsConcreteProcess().

◆ _mesh

◆ _mesh_subset_all_nodes

◆ _process_variables

◆ _secondary_variables

◆ _source_term_collections

std::vector<SourceTermCollection> ProcessLib::Process::_source_term_collections
private

Vector for nodal source term collections. For the monolithic scheme or a single process, the size of the vector is one. For the staggered scheme, the size of vector is the number of the coupled processes.

Definition at line 398 of file Process.h.

Referenced by assemble(), assembleWithJacobian(), and initializeProcessBoundaryConditionsAndSourceTerms().

◆ _sparsity_pattern

GlobalSparsityPattern ProcessLib::Process::_sparsity_pattern
protected

Definition at line 379 of file Process.h.

Referenced by computeSparsityPattern(), and getMatrixSpecifications().

◆ _use_monolithic_scheme

const bool ProcessLib::Process::_use_monolithic_scheme
protected

Definition at line 366 of file Process.h.

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::ComponentTransportProcess(), ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::HydroMechanicsProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleConcreteProcess(), ProcessLib::HT::HTProcess::assembleConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HT::HTProcess::assembleWithJacobianConcreteProcess(), constructDofTable(), ProcessLib::HT::HTProcess::getDOFTableForExtrapolatorData(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::hasMechanicalProcess(), ProcessLib::TH2M::TH2MProcess< DisplacementDim >::hasMechanicalProcess(), ProcessLib::ThermoHydroMechanics::ThermoHydroMechanicsProcess< DisplacementDim >::hasMechanicalProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::initializeConcreteProcess(), ProcessLib::HT::HTProcess::initializeConcreteProcess(), isMonolithicSchemeUsed(), ProcessLib::HT::HTProcess::postTimestepConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::postTimestepConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::postTimestepConcreteProcess(), ProcessLib::ComponentTransport::ComponentTransportProcess::preOutputConcreteProcess(), and ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::setInitialConditionsConcreteProcess().

◆ constant_one_parameter_name

const std::string ProcessLib::Process::constant_one_parameter_name = "constant_one"
static

◆ name

std::string const ProcessLib::Process::name

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