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

Detailed Description

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

Definition at line 22 of file SmallDeformationProcess.h.

#include <SmallDeformationProcess.h>

Inheritance diagram for ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::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)
 
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.
 
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
 

Private Types

using LocalAssemblerInterface
 

Private Member Functions

void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
void assembleConcreteProcess (double const 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) override
 
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) override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) override
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
 

Private Attributes

SmallDeformationProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 
MeshLib::PropertyVector< double > * _material_forces = nullptr
 
Assembly::GlobalMatrixOutput _global_output
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
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 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
 
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
 

Member Typedef Documentation

◆ LocalAssemblerInterface

template<int DisplacementDim>
using ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::LocalAssemblerInterface
private
Initial value:
SmallDeformationLocalAssemblerInterface<DisplacementDim>

Definition at line 44 of file SmallDeformationProcess.h.

Constructor & Destructor Documentation

◆ SmallDeformationProcess()

template<int DisplacementDim>
ProcessLib::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 30 of file SmallDeformationProcess.cpp.

40 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
41 integration_order, std::move(process_variables),
42 std::move(secondary_variables)),
43 _process_data(std::move(process_data))
44{
45 _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
46 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
47
48 _material_forces = MeshLib::getOrCreateMeshProperty<double>(
49 mesh, "MaterialForces", MeshLib::MeshItemType::Node, DisplacementDim);
50
51 _process_data.principal_stress_vector[0] =
52 MeshLib::getOrCreateMeshProperty<double>(
53 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_1",
55
56 _process_data.principal_stress_vector[1] =
57 MeshLib::getOrCreateMeshProperty<double>(
58 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_2",
60
61 _process_data.principal_stress_vector[2] =
62 MeshLib::getOrCreateMeshProperty<double>(
63 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_3",
65
66 _process_data.principal_stress_values =
67 MeshLib::getOrCreateMeshProperty<double>(
68 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_values",
70
72 DisplacementDim>(SmallDeformationLocalAssemblerInterface<
73 DisplacementDim>::getReflectionDataForOutput(),
74 _integration_point_writer, integration_order,
76}
std::string const name
Definition Process.h:354
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:380
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:44
SmallDeformationProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
void addReflectedIntegrationPointWriters(ReflData const &reflection_data, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writers, unsigned const integration_order, std::vector< std::unique_ptr< LocAsmIF > > const &local_assemblers)

References ProcessLib::Process::_integration_point_writer, ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_local_assemblers, ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_material_forces, ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_nodal_forces, ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_process_data, ProcessLib::Reflection::addReflectedIntegrationPointWriters(), MeshLib::Cell, and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::assembleConcreteProcess ( double const 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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 134 of file SmallDeformationProcess.cpp.

138{
139 DBUG("Assemble SmallDeformationProcess.");
140
141 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
142
143 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
145 // Call global assembler for each local assembly item.
148 pv.getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, M, K,
149 b);
150
151 _global_output(t, process_id, M, K, b);
152}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
VectorMatrixAssembler _global_assembler
Definition Process.h:367
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:360
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const 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)
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(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::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 )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 155 of file SmallDeformationProcess.cpp.

160{
161 DBUG("AssembleWithJacobian SmallDeformationProcess.");
162
163 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
164
165 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
167 // Call global assembler for each local assembly item.
170 _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x,
171 x_prev, process_id, M, K, b, Jac);
172
174 *_nodal_forces, std::negate<double>());
175
176 _global_output(t, process_id, M, K, b, &Jac);
177}
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, 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)
void transformVariableFromGlobalVector(GlobalVector const &input_vector, int const variable_id, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, MeshLib::PropertyVector< double > &output_vector, Functor map_function)

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

◆ computeSecondaryVariableConcrete()

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

Reimplemented from ProcessLib::Process.

Definition at line 205 of file SmallDeformationProcess.cpp.

208{
209 DBUG("Compute the secondary variables for SmallDeformationProcess.");
210 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
211 dof_tables.reserve(x.size());
212 std::generate_n(std::back_inserter(dof_tables), x.size(),
213 [&]() { return _local_to_global_index_map.get(); });
214
215 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
218 pv.getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id);
219}
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_prev, 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(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::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 85 of file SmallDeformationProcess.cpp.

89{
91 DisplacementDim, SmallDeformationLocalAssembler>(
92 mesh.getElements(), dof_table, _local_assemblers,
93 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
95
96 auto add_secondary_variable = [&](std::string const& name,
97 int const num_components,
98 auto get_ip_values_function)
99 {
101 name,
102 makeExtrapolator(num_components, getExtrapolator(),
104 std::move(get_ip_values_function)));
105 };
106
107 ProcessLib::Reflection::addReflectedSecondaryVariables<DisplacementDim>(
108 SmallDeformationLocalAssemblerInterface<
109 DisplacementDim>::getReflectionDataForOutput(),
111
112 //
113 // enable output of internal variables defined by material models
114 //
116 LocalAssemblerInterface>(_process_data.solid_materials,
117 add_secondary_variable);
118
121 _process_data.solid_materials, _local_assemblers,
122 _integration_point_writer, integration_order);
123
126
127 // Initialize local assemblers after all variables have been set.
131}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:362
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
SmallDeformationLocalAssemblerInterface< DisplacementDim > LocalAssemblerInterface
void solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::unique_ptr< SolidMaterial > > const &solid_materials, std::vector< std::unique_ptr< LocalAssemblerInterface > > const &local_assemblers, std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > &integration_point_writer, int const integration_order)
void solidMaterialInternalToSecondaryVariables(std::map< int, std::unique_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)

References ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getElements(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), ProcessLib::setIPDataInitialConditions(), ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables(), and ProcessLib::Deformation::solidMaterialInternalVariablesToIntegrationPointWriter().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 79 of file SmallDeformationProcess.cpp.

80{
81 return false;
82}

◆ postTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 180 of file SmallDeformationProcess.cpp.

184{
185 DBUG("PostTimestep SmallDeformationProcess.");
186 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
187 dof_tables.reserve(x.size());
188 std::generate_n(std::back_inserter(dof_tables), x.size(),
189 [&]() { return _local_to_global_index_map.get(); });
190
191 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
194 pv.getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
195
196 std::unique_ptr<GlobalVector> material_forces;
199 *x[process_id]);
200
201 material_forces->copyValues(*_material_forces);
202}
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void writeMaterialForces(std::unique_ptr< GlobalVector > &material_forces, std::vector< std::unique_ptr< LocalAssemblerInterface > > const &local_assemblers, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, GlobalVector const &x)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), ProcessLib::LocalAssemblerInterface::postTimestep(), and ProcessLib::SmallDeformation::writeMaterialForces().

Member Data Documentation

◆ _global_output

template<int DisplacementDim>
Assembly::GlobalMatrixOutput ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_global_output
private

Definition at line 82 of file SmallDeformationProcess.h.

◆ _local_assemblers

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

◆ _material_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_material_forces = nullptr
private

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

template<int DisplacementDim>
SmallDeformationProcessData<DisplacementDim> ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::_process_data
private

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