Loading [MathJax]/extensions/tex2jax.js
OGS
ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

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

Definition at line 23 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, 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::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () 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)
 
bool requiresNormalization () const override
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
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, GlobalVector &b, GlobalMatrix &Jac) override
 
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) 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
 
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id) override
 
- Private Member Functions inherited from ProcessLib::AssemblyMixin< SmallDeformationProcess< DisplacementDim > >
void initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
 
void updateActiveElements ()
 
void assemble (double const t, double const dt, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const process_id, GlobalMatrix &, GlobalMatrix &, GlobalVector &)
 
void assembleWithJacobian (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac)
 

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
 

Friends

class AssemblyMixin< SmallDeformationProcess< DisplacementDim > >
 

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
 
CellAverageData cell_average_data_
 
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 49 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 32 of file SmallDeformationProcess.cpp.

42 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
43 integration_order, std::move(process_variables),
44 std::move(secondary_variables)),
45 AssemblyMixin<SmallDeformationProcess<DisplacementDim>>{
47 process_data_(std::move(process_data))
48{
50 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
51
53 mesh, "MaterialForces", MeshLib::MeshItemType::Node, DisplacementDim);
54
55 process_data_.principal_stress_vector[0] =
57 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_1",
59
60 process_data_.principal_stress_vector[1] =
62 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_2",
64
65 process_data_.principal_stress_vector[2] =
67 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_3",
69
70 process_data_.principal_stress_values =
72 const_cast<MeshLib::Mesh&>(mesh), "principal_stress_values",
74
76 DisplacementDim>(SmallDeformationLocalAssemblerInterface<
77 DisplacementDim>::getReflectionDataForOutput(),
78 _integration_point_writer, integration_order,
80}
std::string const name
Definition Process.h:362
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > _jacobian_assembler
Definition Process.h:376
SmallDeformationProcessData< DisplacementDim > process_data_
std::vector< std::unique_ptr< LocalAssemblerInterface > > local_assemblers_
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
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::_jacobian_assembler.

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 138 of file SmallDeformationProcess.cpp.

142{
143 DBUG("Assemble SmallDeformationProcess.");
144
145 AssemblyMixin<SmallDeformationProcess<DisplacementDim>>::assemble(
146 t, dt, x, x_prev, process_id, M, K, b);
147}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
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
Definition Process.cpp:265

References DBUG().

◆ 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,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 150 of file SmallDeformationProcess.cpp.

155{
156 DBUG("AssembleWithJacobian SmallDeformationProcess.");
157
158 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
160 // Call global assembler for each local assembly item.
161 try
162 {
163 AssemblyMixin<SmallDeformationProcess<DisplacementDim>>::
164 assembleWithJacobian(t, dt, x, x_prev, process_id, b, Jac);
165 }
166 catch (NumLib::AssemblyException const&)
167 {
170 std::negate<double>());
171 throw;
172 }
174 *nodal_forces_, std::negate<double>());
175}
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
void assembleWithJacobian(const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac) final
Definition Process.cpp:286
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 DBUG().

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

235{
236 DBUG("Compute the secondary variables for SmallDeformationProcess.");
237 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
238 dof_tables.reserve(x.size());
239 std::generate_n(std::back_inserter(dof_tables), x.size(),
240 [&]() { return _local_to_global_index_map.get(); });
241
244 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id);
245
247}
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)
CellAverageData cell_average_data_
Definition Process.h:372
std::vector< std::size_t > const & getActiveElementIDs() const
Definition Process.h:167
void computeCellAverages(CellAverageData &cell_average_data, std::vector< std::unique_ptr< LAIntf > > const &local_assemblers)
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::computeCellAverages(), ProcessLib::LocalAssemblerInterface::computeSecondaryVariable(), DBUG(), and NumLib::SerialExecutor::executeSelectedMemberOnDereferenced().

◆ initializeAssemblyOnSubmeshes()

template<int DisplacementDim>
std::vector< std::vector< std::string > > ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & meshes)
overrideprivatevirtual

Initializes the assembly on submeshes

Parameters
meshesthe submeshes on whom the assembly shall proceed.
Attention
meshes must be a must be a non-overlapping cover of the entire simulation domain (bulk mesh)!
Returns
The names of the residuum vectors that will be assembled for each process: outer vector of size 1 for monolithic schemes and greater for staggered schemes.

Reimplemented from ProcessLib::SubmeshAssemblySupport.

Definition at line 219 of file SmallDeformationProcess.cpp.

221{
222 INFO("SmallDeformation process initializeSubmeshOutput().");
223 std::vector<std::vector<std::string>> residuum_names{{"NodalForces"}};
224
225 AssemblyMixin<SmallDeformationProcess<DisplacementDim>>::
226 initializeAssemblyOnSubmeshes(meshes, residuum_names);
227
228 return residuum_names;
229}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
std::vector< std::vector< std::string > > initializeAssemblyOnSubmeshes(std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes) override

References INFO().

◆ 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 89 of file SmallDeformationProcess.cpp.

93{
95 DisplacementDim, SmallDeformationLocalAssembler>(
96 mesh.getElements(), dof_table, local_assemblers_,
97 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
99
100 auto add_secondary_variable = [&](std::string const& name,
101 int const num_components,
102 auto get_ip_values_function)
103 {
105 name,
106 makeExtrapolator(num_components, getExtrapolator(),
108 std::move(get_ip_values_function)));
109 };
110
112 SmallDeformationLocalAssemblerInterface<
113 DisplacementDim>::getReflectionDataForOutput(),
115
116 //
117 // enable output of internal variables defined by material models
118 //
120 LocalAssemblerInterface>(process_data_.solid_materials,
121 add_secondary_variable);
122
125 process_data_.solid_materials, local_assemblers_,
126 _integration_point_writer, integration_order);
127
130
131 // Initialize local assemblers after all variables have been set.
135}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
SmallDeformationLocalAssemblerInterface< DisplacementDim > LocalAssemblerInterface
void solidMaterialInternalVariablesToIntegrationPointWriter(std::map< int, std::shared_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::shared_ptr< SolidMaterial > > const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
void addReflectedSecondaryVariables(ReflData const &reflection_data, SecondaryVariableCollection &secondary_variables, NumLib::Extrapolator &extrapolator, std::vector< std::unique_ptr< LocAsmIF > > const &local_assemblers)
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::Reflection::addReflectedSecondaryVariables(), 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 83 of file SmallDeformationProcess.cpp.

84{
85 return false;
86}

◆ 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 194 of file SmallDeformationProcess.cpp.

198{
199 DBUG("PostTimestep SmallDeformationProcess.");
200 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
201 dof_tables.reserve(x.size());
202 std::generate_n(std::back_inserter(dof_tables), x.size(),
203 [&]() { return _local_to_global_index_map.get(); });
204
207 getActiveElementIDs(), dof_tables, x, x_prev, t, dt, process_id);
208
209 std::unique_ptr<GlobalVector> material_forces;
212 *x[process_id]);
213
214 material_forces->copyValues(*material_forces_);
215}
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::LocalAssemblerInterface::postTimestep(), and ProcessLib::SmallDeformation::writeMaterialForces().

◆ preTimestepConcreteProcess()

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

181{
182 DBUG("PreTimestep SmallDeformationProcess.");
183
186 getActiveElementIDs(), *_local_to_global_index_map, *x[process_id], t,
187 dt);
188
189 AssemblyMixin<
191}
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)
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)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Friends And Related Symbol Documentation

◆ AssemblyMixin< SmallDeformationProcess< DisplacementDim > >

template<int DisplacementDim>
friend class AssemblyMixin< SmallDeformationProcess< DisplacementDim > >
friend

Definition at line 96 of file SmallDeformationProcess.h.

Member Data Documentation

◆ local_assemblers_

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

Definition at line 89 of file SmallDeformationProcess.h.

◆ material_forces_

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

Definition at line 92 of file SmallDeformationProcess.h.

◆ nodal_forces_

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

Definition at line 91 of file SmallDeformationProcess.h.

◆ process_data_

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

Definition at line 87 of file SmallDeformationProcess.h.


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