OGS
ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >

Definition at line 24 of file SmallDeformationProcess.h.

#include <SmallDeformationProcess.h>

Inheritance diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >:
[legend]

Public Member Functions

 SmallDeformationProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::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, SmallDeformationProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables)
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id) override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from 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)
 
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. More...
 
void postTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep. More...
 
void postNonLinearSolver (GlobalVector const &x, GlobalVector const &xdot, 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_dot, int const process_id)
 compute secondary variables for the coupled equations or for output. More...
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize ()
 
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 setCoupledSolutionsForStaggeredScheme (CoupledSolutionsForStaggeredScheme *const coupled_solutions)
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
bool isMonolithicSchemeUsed () const
 
virtual void setCoupledTermForTheStaggeredSchemeToLocalAssemblers (int 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 &xdot, 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 &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
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
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< 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)
 

Private Types

using LocalAssemblerInterface = SmallDeformationLocalAssemblerInterface
 

Private Member Functions

void constructDofTable () override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize(). More...
 
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
 
void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
 

Private Attributes

SmallDeformationProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
std::vector< MeshLib::Element * > _vec_matrix_elements
 
std::vector< int > _vec_fracture_mat_IDs
 
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_elements
 
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_matrix_elements
 
std::vector< std::vector< MeshLib::Node * > > _vec_fracture_nodes
 
std::vector< MeshLib::Node * > _vec_junction_nodes
 
std::vector< std::vector< MeshLib::Element * > > _vec_junction_fracture_matrix_elements
 
std::vector< std::unique_ptr< MeshLib::MeshSubset const > > _mesh_subset_fracture_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_junction_nodes
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_matrix_nodes
 

Additional Inherited Members

- Public Types inherited from ProcessLib::Process
using NonlinearSolver = NumLib::NonlinearSolverBase
 
using TimeDiscretization = NumLib::TimeDiscretization
 
- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Protected Member Functions inherited from ProcessLib::Process
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
 
virtual void constructDofTable ()
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_prosess_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
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
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
CoupledSolutionsForStaggeredScheme_coupled_solutions
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Member Typedef Documentation

◆ LocalAssemblerInterface

template<int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::LocalAssemblerInterface = SmallDeformationLocalAssemblerInterface
private

Definition at line 55 of file SmallDeformationProcess.h.

Constructor & Destructor Documentation

◆ SmallDeformationProcess()

template<int DisplacementDim>
ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::SmallDeformationProcess ( std::string  name,
MeshLib::Mesh mesh,
std::unique_ptr< ProcessLib::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,
SmallDeformationProcessData< DisplacementDim > &&  process_data,
SecondaryVariableCollection &&  secondary_variables 
)

Definition at line 31 of file SmallDeformationProcess.cpp.

41 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
42 integration_order, std::move(process_variables),
43 std::move(secondary_variables)),
44 _process_data(std::move(process_data))
45{
46 std::vector<std::pair<std::size_t, std::vector<int>>>
47 vec_branch_nodeID_matIDs;
48 std::vector<std::pair<std::size_t, std::vector<int>>>
49 vec_junction_nodeID_matIDs;
53 _vec_fracture_nodes, vec_branch_nodeID_matIDs,
54 vec_junction_nodeID_matIDs);
55
56 if (_vec_fracture_mat_IDs.size() !=
57 _process_data.fracture_properties.size())
58 {
60 "The number of the given fracture properties ({:d}) are not "
61 "consistent with the number of fracture groups in a mesh ({:d}).",
62 _process_data.fracture_properties.size(),
64 }
65
66 // create a map from a material ID to a fracture ID
67 auto max_frac_mat_id = std::max_element(_vec_fracture_mat_IDs.begin(),
69 _process_data._map_materialID_to_fractureID.resize(*max_frac_mat_id + 1);
70 for (unsigned i = 0; i < _vec_fracture_mat_IDs.size(); i++)
71 {
72 _process_data._map_materialID_to_fractureID[_vec_fracture_mat_IDs[i]] =
73 i;
74 }
75
76 // create a table of connected fracture IDs for each element
77 _process_data._vec_ele_connected_fractureIDs.resize(
78 mesh.getNumberOfElements());
79 for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
80 {
81 for (auto e : _vec_fracture_matrix_elements[i])
82 {
83 _process_data._vec_ele_connected_fractureIDs[e->getID()].push_back(
84 i);
85 }
86 }
87
88 // set fracture property
89 for (auto& fracture_prop : _process_data.fracture_properties)
90 {
91 // based on the 1st element assuming a fracture forms a straight line
93 DisplacementDim,
94 *_vec_fracture_elements[fracture_prop.fracture_id][0],
95 fracture_prop);
96 }
97
98 // set branches
99 for (auto& vec_branch_nodeID_matID : vec_branch_nodeID_matIDs)
100 {
101 auto master_matId = vec_branch_nodeID_matID.second[0];
102 auto slave_matId = vec_branch_nodeID_matID.second[1];
103 auto& master_frac =
104 _process_data.fracture_properties
105 [_process_data._map_materialID_to_fractureID[master_matId]];
106 auto& slave_frac =
107 _process_data.fracture_properties
108 [_process_data._map_materialID_to_fractureID[slave_matId]];
109
110 master_frac.branches_master.push_back(
111 createBranchProperty(*mesh.getNode(vec_branch_nodeID_matID.first),
112 master_frac, slave_frac));
113
114 slave_frac.branches_slave.push_back(
115 createBranchProperty(*mesh.getNode(vec_branch_nodeID_matID.first),
116 master_frac, slave_frac));
117 }
118
119 // set junctions
120 transform(cbegin(vec_junction_nodeID_matIDs),
121 cend(vec_junction_nodeID_matIDs),
122 back_inserter(_vec_junction_nodes),
123 [&](auto& vec_junction_nodeID_matID)
124 {
125 return const_cast<MeshLib::Node*>(
126 _mesh.getNode(vec_junction_nodeID_matID.first));
127 });
128
129 for (std::size_t i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
130 {
131 auto const& material_ids = vec_junction_nodeID_matIDs[i].second;
132 assert(material_ids.size() == 2);
133 std::array<int, 2> fracture_ids{
134 {_process_data._map_materialID_to_fractureID[material_ids[0]],
135 _process_data._map_materialID_to_fractureID[material_ids[1]]}};
136
137 _process_data.junction_properties.emplace_back(
138 i, *mesh.getNode(vec_junction_nodeID_matIDs[i].first),
139 fracture_ids);
140 }
141
142 // create a table of connected junction IDs for each element
143 _process_data._vec_ele_connected_junctionIDs.resize(
144 mesh.getNumberOfElements());
145 for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
146 {
147 auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
148 for (auto e : mesh.getElementsConnectedToNode(*node))
149 {
150 _process_data._vec_ele_connected_junctionIDs[e->getID()].push_back(
151 i);
152 }
153 }
154
155 // create a table of junction node and connected elements
157 vec_junction_nodeID_matIDs.size());
158 for (unsigned i = 0; i < vec_junction_nodeID_matIDs.size(); i++)
159 {
160 auto node = mesh.getNode(vec_junction_nodeID_matIDs[i].first);
161 for (auto e : mesh.getElementsConnectedToNode(*node))
162 {
164 const_cast<MeshLib::Element*>(e));
165 }
166 }
167
168 //
169 // If Neumann BCs for the displacement_jump variable are required they need
170 // special treatment because of the levelset function. The implementation
171 // exists in the version 6.1.0 (e54815cc07ee89c81f953a4955b1c788595dd725)
172 // and was removed due to lack of applications.
173 //
174
175 MeshLib::PropertyVector<int> const* material_ids(
176 mesh.getProperties().getPropertyVector<int>("MaterialIDs"));
177 _process_data._mesh_prop_materialIDs = material_ids;
178}
#define OGS_FATAL(...)
Definition: Error.h:26
const Node * getNode(std::size_t idx) const
Get the node with the given index.
Definition: Mesh.h:79
Properties & getProperties()
Definition: Mesh.h:128
std::vector< Element const * > const & getElementsConnectedToNode(std::size_t node_id) const
Definition: Mesh.cpp:238
std::size_t getNumberOfElements() const
Get the number of elements.
Definition: Mesh.h:91
PropertyVector< T > const * getPropertyVector(std::string const &name) const
std::vector< std::vector< MeshLib::Element * > > _vec_junction_fracture_matrix_elements
std::vector< std::vector< MeshLib::Node * > > _vec_fracture_nodes
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_elements
std::vector< std::vector< MeshLib::Element * > > _vec_fracture_matrix_elements
SmallDeformationProcessData< DisplacementDim > _process_data
std::string const name
Definition: Process.h:325
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: Process.cpp:22
MeshLib::Mesh & _mesh
Definition: Process.h:328
void setFractureProperty(int const dim, MeshLib::Element const &e, FractureProperty &frac_prop)
void getFractureMatrixDataInMesh(MeshLib::Mesh const &mesh, std::vector< MeshLib::Element * > &vec_matrix_elements, std::vector< int > &vec_fracture_mat_IDs, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_elements, std::vector< std::vector< MeshLib::Element * > > &vec_fracture_matrix_elements, std::vector< std::vector< MeshLib::Node * > > &vec_fracture_nodes, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_branch_nodeID_matIDs, std::vector< std::pair< std::size_t, std::vector< int > > > &vec_junction_nodeID_matIDs)
Definition: MeshUtils.cpp:213
BranchProperty createBranchProperty(MeshLib::Node const &branchNode, FractureProperty const &master_frac, FractureProperty const &slave_frac)

References ProcessLib::Process::_mesh, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_process_data, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_elements, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_mat_IDs, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_matrix_elements, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_nodes, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_junction_fracture_matrix_elements, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_junction_nodes, ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_matrix_elements, ProcessLib::LIE::createBranchProperty(), MeshLib::Mesh::getElementsConnectedToNode(), ProcessLib::LIE::getFractureMatrixDataInMesh(), MeshLib::Mesh::getNode(), MeshLib::Mesh::getNumberOfElements(), MeshLib::Mesh::getProperties(), MeshLib::Properties::getPropertyVector(), OGS_FATAL, and ProcessLib::LIE::setFractureProperty().

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleConcreteProcess ( const double  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
std::vector< GlobalVector * > const &  xdot,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 553 of file SmallDeformationProcess.cpp.

557{
558 DBUG("Assemble SmallDeformationProcess.");
559
560 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
561
562 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
563 dof_table = {std::ref(*_local_to_global_index_map)};
564 // Call global assembler for each local assembly item.
567 pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
568 b);
569}
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:29
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:146
VectorMatrixAssembler _global_assembler
Definition: Process.h:335
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:331
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
static const double t
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and MathLib::t.

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( const double  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
std::vector< GlobalVector * > const &  xdot,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 571 of file SmallDeformationProcess.cpp.

578{
579 DBUG("AssembleWithJacobian SmallDeformationProcess.");
580
581 // Call global assembler for each local assembly item.
582 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
583
584 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
585 dof_table = {std::ref(*_local_to_global_index_map)};
588 _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
589 process_id, M, K, b, Jac);
590}
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and MathLib::t.

◆ computeSecondaryVariableConcrete()

template<int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::computeSecondaryVariableConcrete ( double const  t,
double const  dt,
std::vector< GlobalVector * > const &  x,
GlobalVector const &  x_dot,
int const  process_id 
)
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 529 of file SmallDeformationProcess.cpp.

532{
533 DBUG("Compute the secondary variables for SmallDeformationProcess.");
534 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
535 dof_tables.reserve(x.size());
536 std::generate_n(std::back_inserter(dof_tables), x.size(),
537 [&]() { return _local_to_global_index_map.get(); });
538
539 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
542 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
543 x_dot, process_id);
544}
virtual void computeSecondaryVariable(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and MathLib::t.

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::constructDofTable
overrideprivatevirtual

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 from ProcessLib::Process.

Definition at line 181 of file SmallDeformationProcess.cpp.

182{
183 //------------------------------------------------------------
184 // prepare mesh subsets to define DoFs
185 //------------------------------------------------------------
186 // for extrapolation
188 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
189 // regular u
191 std::make_unique<MeshLib::MeshSubset>(_mesh, _mesh.getNodes());
192 // u jump
193 for (unsigned i = 0; i < _vec_fracture_nodes.size(); i++)
194 {
196 std::make_unique<MeshLib::MeshSubset const>(
198 }
199 // enrichment for junctions
201 std::make_unique<MeshLib::MeshSubset>(_mesh, _vec_junction_nodes);
202
203 // Collect the mesh subsets in a vector.
204 std::vector<MeshLib::MeshSubset> all_mesh_subsets;
205 std::generate_n(std::back_inserter(all_mesh_subsets), DisplacementDim,
206 [&]() { return *_mesh_subset_matrix_nodes; });
207 for (auto& ms : _mesh_subset_fracture_nodes)
208 {
209 std::generate_n(std::back_inserter(all_mesh_subsets),
210 DisplacementDim,
211 [&]() { return *ms; });
212 }
213 std::generate_n(std::back_inserter(all_mesh_subsets),
214 DisplacementDim,
215 [&]() { return *_mesh_subset_junction_nodes; });
216
217 std::vector<int> const vec_n_components(
218 1 + _vec_fracture_mat_IDs.size() + _vec_junction_nodes.size(),
219 DisplacementDim);
220
221 std::vector<std::vector<MeshLib::Element*> const*> vec_var_elements;
222 vec_var_elements.push_back(&_vec_matrix_elements);
223 for (unsigned i = 0; i < _vec_fracture_matrix_elements.size(); i++)
224 {
225 vec_var_elements.push_back(&_vec_fracture_matrix_elements[i]);
226 }
227 for (unsigned i = 0; i < _vec_junction_fracture_matrix_elements.size(); i++)
228 {
229 vec_var_elements.push_back(&_vec_junction_fracture_matrix_elements[i]);
230 }
231
233 std::make_unique<NumLib::LocalToGlobalIndexMap>(
234 std::move(all_mesh_subsets),
235 vec_n_components,
236 vec_var_elements,
238}
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition: Mesh.h:100
std::vector< std::unique_ptr< MeshLib::MeshSubset const > > _mesh_subset_fracture_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_junction_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_matrix_nodes
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:329
@ BY_COMPONENT
Ordering data by component type.

References NumLib::BY_COMPONENT.

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const &  dof_table,
MeshLib::Mesh const &  mesh,
unsigned const  integration_order 
)
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 241 of file SmallDeformationProcess.cpp.

245{
247 DisplacementDim, SmallDeformationLocalAssemblerMatrix,
248 SmallDeformationLocalAssemblerMatrixNearFracture,
249 SmallDeformationLocalAssemblerFracture>(
250 mesh.getElements(), dof_table, _local_assemblers,
251 mesh.isAxiallySymmetric(), integration_order, _process_data);
252
253 // TODO move the two data members somewhere else.
254 // for extrapolation of secondary variables
255 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
258 std::make_unique<NumLib::LocalToGlobalIndexMap>(
259 std::move(all_mesh_subsets_single_component),
260 // by location order is needed for output
262
264 "sigma_xx",
268
270 "sigma_yy",
274
276 "sigma_zz",
280
282 "sigma_xy",
286
287 if (DisplacementDim == 3)
288 {
290 "sigma_xz",
294
296 "sigma_yz",
300 }
301
303 "epsilon_xx",
307
309 "epsilon_yy",
313
315 "epsilon_zz",
319
321 "epsilon_xy",
325
326 if (DisplacementDim == 3)
327 {
329 "epsilon_xz",
333
335 "epsilon_yz",
339 }
340
341 auto mesh_prop_sigma_xx = MeshLib::getOrCreateMeshProperty<double>(
342 const_cast<MeshLib::Mesh&>(mesh), "stress_xx",
344 mesh_prop_sigma_xx->resize(mesh.getNumberOfElements());
345 _process_data._mesh_prop_stress_xx = mesh_prop_sigma_xx;
346
347 auto mesh_prop_sigma_yy = MeshLib::getOrCreateMeshProperty<double>(
348 const_cast<MeshLib::Mesh&>(mesh), "stress_yy",
350 mesh_prop_sigma_yy->resize(mesh.getNumberOfElements());
351 _process_data._mesh_prop_stress_yy = mesh_prop_sigma_yy;
352
353 auto mesh_prop_sigma_zz = MeshLib::getOrCreateMeshProperty<double>(
354 const_cast<MeshLib::Mesh&>(mesh), "stress_zz",
356 mesh_prop_sigma_zz->resize(mesh.getNumberOfElements());
357 _process_data._mesh_prop_stress_zz = mesh_prop_sigma_zz;
358
359 auto mesh_prop_sigma_xy = MeshLib::getOrCreateMeshProperty<double>(
360 const_cast<MeshLib::Mesh&>(mesh), "stress_xy",
362 mesh_prop_sigma_xy->resize(mesh.getNumberOfElements());
363 _process_data._mesh_prop_stress_xy = mesh_prop_sigma_xy;
364
365 if (DisplacementDim == 3)
366 {
367 auto mesh_prop_sigma_xz = MeshLib::getOrCreateMeshProperty<double>(
368 const_cast<MeshLib::Mesh&>(mesh), "stress_xz",
370 mesh_prop_sigma_xz->resize(mesh.getNumberOfElements());
371 _process_data._mesh_prop_stress_xz = mesh_prop_sigma_xz;
372
373 auto mesh_prop_sigma_yz = MeshLib::getOrCreateMeshProperty<double>(
374 const_cast<MeshLib::Mesh&>(mesh), "stress_yz",
376 mesh_prop_sigma_yz->resize(mesh.getNumberOfElements());
377 _process_data._mesh_prop_stress_yz = mesh_prop_sigma_yz;
378 }
379
380 auto mesh_prop_epsilon_xx = MeshLib::getOrCreateMeshProperty<double>(
381 const_cast<MeshLib::Mesh&>(mesh), "strain_xx",
383 mesh_prop_epsilon_xx->resize(mesh.getNumberOfElements());
384 _process_data._mesh_prop_strain_xx = mesh_prop_epsilon_xx;
385
386 auto mesh_prop_epsilon_yy = MeshLib::getOrCreateMeshProperty<double>(
387 const_cast<MeshLib::Mesh&>(mesh), "strain_yy",
389 mesh_prop_epsilon_yy->resize(mesh.getNumberOfElements());
390 _process_data._mesh_prop_strain_yy = mesh_prop_epsilon_yy;
391
392 auto mesh_prop_epsilon_zz = MeshLib::getOrCreateMeshProperty<double>(
393 const_cast<MeshLib::Mesh&>(mesh), "strain_zz",
395 mesh_prop_epsilon_zz->resize(mesh.getNumberOfElements());
396 _process_data._mesh_prop_strain_zz = mesh_prop_epsilon_zz;
397
398 auto mesh_prop_epsilon_xy = MeshLib::getOrCreateMeshProperty<double>(
399 const_cast<MeshLib::Mesh&>(mesh), "strain_xy",
401 mesh_prop_epsilon_xy->resize(mesh.getNumberOfElements());
402 _process_data._mesh_prop_strain_xy = mesh_prop_epsilon_xy;
403
404 if (DisplacementDim == 3)
405 {
406 auto mesh_prop_epsilon_xz = MeshLib::getOrCreateMeshProperty<double>(
407 const_cast<MeshLib::Mesh&>(mesh), "strain_xz",
409 mesh_prop_epsilon_xz->resize(mesh.getNumberOfElements());
410 _process_data._mesh_prop_strain_xz = mesh_prop_epsilon_xz;
411
412 auto mesh_prop_epsilon_yz = MeshLib::getOrCreateMeshProperty<double>(
413 const_cast<MeshLib::Mesh&>(mesh), "strain_yz",
415 mesh_prop_epsilon_yz->resize(mesh.getNumberOfElements());
416 _process_data._mesh_prop_strain_yz = mesh_prop_epsilon_yz;
417 }
418
419 for (MeshLib::Element const* e : _mesh.getElements())
420 {
421 if (e->getDimension() < DisplacementDim)
422 {
423 continue;
424 }
425
426 Eigen::Vector3d const pt(getCenterOfGravity(*e).asEigenVector3d());
427 std::vector<FractureProperty*> e_fracture_props;
428 std::unordered_map<int, int> e_fracID_to_local;
429 unsigned tmpi = 0;
430 for (auto fid :
431 _process_data._vec_ele_connected_fractureIDs[e->getID()])
432 {
433 e_fracture_props.push_back(&_process_data.fracture_properties[fid]);
434 e_fracID_to_local.insert({fid, tmpi++});
435 }
436 std::vector<JunctionProperty*> e_junction_props;
437 std::unordered_map<int, int> e_juncID_to_local;
438 tmpi = 0;
439 for (auto fid :
440 _process_data._vec_ele_connected_junctionIDs[e->getID()])
441 {
442 e_junction_props.push_back(&_process_data.junction_properties[fid]);
443 e_juncID_to_local.insert({fid, tmpi++});
444 }
445 std::vector<double> const levelsets(uGlobalEnrichments(
446 e_fracture_props, e_junction_props, e_fracID_to_local, pt));
447
448 for (unsigned i = 0; i < e_fracture_props.size(); i++)
449 {
450 auto mesh_prop_levelset = MeshLib::getOrCreateMeshProperty<double>(
451 const_cast<MeshLib::Mesh&>(mesh),
452 "levelset" +
453 std::to_string(e_fracture_props[i]->fracture_id + 1),
455 mesh_prop_levelset->resize(mesh.getNumberOfElements());
456 (*mesh_prop_levelset)[e->getID()] = levelsets[i];
457 }
458 for (unsigned i = 0; i < e_junction_props.size(); i++)
459 {
460 auto mesh_prop_levelset = MeshLib::getOrCreateMeshProperty<double>(
461 const_cast<MeshLib::Mesh&>(mesh),
462 "levelset" +
463 std::to_string(e_junction_props[i]->junction_id + 1 +
464 _process_data.fracture_properties.size()),
466 mesh_prop_levelset->resize(mesh.getNumberOfElements());
467 (*mesh_prop_levelset)[e->getID()] =
468 levelsets[i + e_fracture_props.size()];
469 }
470 }
471
472 auto mesh_prop_w_n = MeshLib::getOrCreateMeshProperty<double>(
473 const_cast<MeshLib::Mesh&>(mesh), "w_n", MeshLib::MeshItemType::Cell,
474 1);
475 mesh_prop_w_n->resize(mesh.getNumberOfElements());
476 _process_data._mesh_prop_w_n = mesh_prop_w_n;
477
478 auto mesh_prop_w_s = MeshLib::getOrCreateMeshProperty<double>(
479 const_cast<MeshLib::Mesh&>(mesh), "w_s", MeshLib::MeshItemType::Cell,
480 1);
481 mesh_prop_w_s->resize(mesh.getNumberOfElements());
482 _process_data._mesh_prop_w_s = mesh_prop_w_s;
483
484 auto mesh_prop_b = MeshLib::getOrCreateMeshProperty<double>(
485 const_cast<MeshLib::Mesh&>(mesh), "aperture",
487 mesh_prop_b->resize(mesh.getNumberOfElements());
488 auto const& mesh_prop_matid = *_process_data._mesh_prop_materialIDs;
489 for (auto const& fracture_prop : _process_data.fracture_properties)
490 {
491 for (MeshLib::Element const* e : _mesh.getElements())
492 {
493 if (e->getDimension() == DisplacementDim)
494 {
495 continue;
496 }
497 if (mesh_prop_matid[e->getID()] != fracture_prop.mat_id)
498 {
499 continue;
500 }
501 // Mean value for the element. This allows usage of node based
502 // properties for aperture.
503 (*mesh_prop_b)[e->getID()] =
504 fracture_prop.aperture0
505 .getNodalValuesOnElement(*e, /*time independent*/ 0)
506 .mean();
507 }
508 }
509 _process_data._mesh_prop_b = mesh_prop_b;
510
511 auto mesh_prop_fracture_stress_shear =
512 MeshLib::getOrCreateMeshProperty<double>(
513 const_cast<MeshLib::Mesh&>(mesh), "f_stress_s",
515 mesh_prop_fracture_stress_shear->resize(mesh.getNumberOfElements());
516 _process_data._mesh_prop_fracture_stress_shear =
517 mesh_prop_fracture_stress_shear;
518
519 auto mesh_prop_fracture_stress_normal =
520 MeshLib::getOrCreateMeshProperty<double>(
521 const_cast<MeshLib::Mesh&>(mesh), "f_stress_n",
523 mesh_prop_fracture_stress_normal->resize(mesh.getNumberOfElements());
524 _process_data._mesh_prop_fracture_stress_normal =
525 mesh_prop_fracture_stress_normal;
526}
Eigen::Vector3d const & asEigenVector3d() const
Definition: Point3d.h:67
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:103
virtual std::vector< double > const & getIntPtEpsilonXY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaYZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonXX(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaXY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaXZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonYZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaXX(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonYY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonZZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtEpsilonXZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaZZ(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigmaYY(const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
SecondaryVariableCollection _secondary_variables
Definition: Process.h:333
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:182
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
MathLib::Point3d getCenterOfGravity(Element const &element)
Calculates the center of gravity for the mesh element.
Definition: Element.cpp:124
@ BY_LOCATION
Ordering data by spatial location.
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ExtraCtorArgs &&... extra_ctor_args)
std::vector< double > uGlobalEnrichments(std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)

References MathLib::Point3d::asEigenVector3d(), NumLib::BY_LOCATION, MeshLib::Cell, ProcessLib::LIE::SmallDeformation::createLocalAssemblers(), MeshLib::getCenterOfGravity(), MeshLib::Mesh::getElements(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtEpsilonXX(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtEpsilonXY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtEpsilonXZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtEpsilonYY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtEpsilonYZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtEpsilonZZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtSigmaXX(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtSigmaXY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtSigmaXZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtSigmaYY(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtSigmaYZ(), ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface::getIntPtSigmaZZ(), MeshLib::Mesh::getNumberOfElements(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::LIE::uGlobalEnrichments().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::isLinear
override

Definition at line 547 of file SmallDeformationProcess.cpp.

548{
549 return false;
550}

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::preTimestepConcreteProcess ( std::vector< GlobalVector * > const &  x,
double const  t,
double const  dt,
const int  process_id 
)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 592 of file SmallDeformationProcess.cpp.

595{
596 DBUG("PreTimestep SmallDeformationProcess.");
597
598 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
602 *_local_to_global_index_map, *x[process_id], t, dt);
603}
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), ProcessLib::LocalAssemblerInterface::preTimestep(), and MathLib::t.

Member Data Documentation

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_local_assemblers
private

Definition at line 83 of file SmallDeformationProcess.h.

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 86 of file SmallDeformationProcess.h.

◆ _mesh_subset_fracture_nodes

template<int DisplacementDim>
std::vector<std::unique_ptr<MeshLib::MeshSubset const> > ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_mesh_subset_fracture_nodes
private

Definition at line 98 of file SmallDeformationProcess.h.

◆ _mesh_subset_junction_nodes

template<int DisplacementDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_mesh_subset_junction_nodes
private

Definition at line 99 of file SmallDeformationProcess.h.

◆ _mesh_subset_matrix_nodes

template<int DisplacementDim>
std::unique_ptr<MeshLib::MeshSubset const> ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_mesh_subset_matrix_nodes
private

Definition at line 100 of file SmallDeformationProcess.h.

◆ _process_data

◆ _vec_fracture_elements

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Element*> > ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_elements
private

◆ _vec_fracture_mat_IDs

template<int DisplacementDim>
std::vector<int> ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_mat_IDs
private

◆ _vec_fracture_matrix_elements

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Element*> > ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_matrix_elements
private

◆ _vec_fracture_nodes

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Node*> > ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_fracture_nodes
private

◆ _vec_junction_fracture_matrix_elements

template<int DisplacementDim>
std::vector<std::vector<MeshLib::Element*> > ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_junction_fracture_matrix_elements
private

◆ _vec_junction_nodes

template<int DisplacementDim>
std::vector<MeshLib::Node*> ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_junction_nodes
private

◆ _vec_matrix_elements

template<int DisplacementDim>
std::vector<MeshLib::Element*> ProcessLib::LIE::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_vec_matrix_elements
private

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