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. 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, const double dxdot_dx, const double dx_dx, 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 * getIntegrationPointWriter (MeshLib::Mesh const &mesh) 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< DisplacementDim >
 

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(). 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, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
 
void computeSecondaryVariableConcrete (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, const int 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
 

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
 
- 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::SmallDeformation::SmallDeformationProcess< DisplacementDim >::LocalAssemblerInterface = SmallDeformationLocalAssemblerInterface<DisplacementDim>
private

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

37  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
38  integration_order, std::move(process_variables),
39  std::move(secondary_variables)),
40  _process_data(std::move(process_data))
41 {
42  _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
43  mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
44 
45  _material_forces = MeshLib::getOrCreateMeshProperty<double>(
46  mesh, "MaterialForces", MeshLib::MeshItemType::Node, DisplacementDim);
47 
48  _process_data.principal_stress_vector[0] =
49  MeshLib::getOrCreateMeshProperty<double>(
50  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_1",
52 
53  _process_data.principal_stress_vector[1] =
54  MeshLib::getOrCreateMeshProperty<double>(
55  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_2",
57 
58  _process_data.principal_stress_vector[2] =
59  MeshLib::getOrCreateMeshProperty<double>(
60  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_vector_3",
62 
63  _process_data.principal_stress_values =
64  MeshLib::getOrCreateMeshProperty<double>(
65  const_cast<MeshLib::Mesh&>(mesh), "principal_stress_values",
67 
68  // TODO (naumov) remove ip suffix. Probably needs modification of the mesh
69  // properties, s.t. there is no "overlapping" with cell/point data.
70  // See getOrCreateMeshProperty.
71  _integration_point_writer.emplace_back(
72  std::make_unique<IntegrationPointWriter>(
73  "sigma_ip",
74  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
75  integration_order, _local_assemblers,
77 }
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::string const name
Definition: Process.h:323
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:350
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
SmallDeformationProcessData< DisplacementDim > _process_data
std::vector< std::unique_ptr< LocalAssemblerInterface > > _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, MeshLib::Cell, MeshLib::Mesh::getDimension(), ProcessLib::SmallDeformation::SmallDeformationLocalAssemblerInterface< DisplacementDim >::getSigma(), and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

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

192 {
193  DBUG("Assemble SmallDeformationProcess.");
194 
195  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
196 
197  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
198  dof_table = {std::ref(*_local_to_global_index_map)};
199  // Call global assembler for each local assembly item.
202  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
203  b);
204 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
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 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 &  xdot,
const double  dxdot_dx,
const double  dx_dx,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 207 of file SmallDeformationProcess.cpp.

213 {
214  DBUG("AssembleWithJacobian SmallDeformationProcess.");
215 
216  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
217 
218  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
219  dof_table = {std::ref(*_local_to_global_index_map)};
220  // Call global assembler for each local assembly item.
223  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot,
224  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
225 
227  *_nodal_forces, std::negate<double>());
228 }
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, const double dxdot_dx, const double dx_dx, 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 mapFunction)
Definition: DOFTableUtil.h:59

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

◆ computeSecondaryVariableConcrete()

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

Reimplemented from ProcessLib::Process.

Definition at line 255 of file SmallDeformationProcess.cpp.

258 {
259  DBUG("Compute the secondary variables for SmallDeformationProcess.");
260  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
261  dof_tables.reserve(x.size());
262  std::generate_n(std::back_inserter(dof_tables), x.size(),
263  [&]() { return _local_to_global_index_map.get(); });
264 
265  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
268  pv.getActiveElementIDs(), dof_tables, t, dt, x, x_dot, process_id);
269 }
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(), 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 86 of file SmallDeformationProcess.cpp.

90 {
91  using nlohmann::json;
92 
94  DisplacementDim, SmallDeformationLocalAssembler>(
95  mesh.getElements(), dof_table, _local_assemblers,
96  mesh.isAxiallySymmetric(), integration_order, _process_data);
97 
98  auto add_secondary_variable = [&](std::string const& name,
99  int const num_components,
100  auto get_ip_values_function)
101  {
103  name,
104  makeExtrapolator(num_components, getExtrapolator(),
106  std::move(get_ip_values_function)));
107  };
108 
109  add_secondary_variable("free_energy_density",
110  1,
112 
113  add_secondary_variable("sigma",
115  DisplacementDim>::RowsAtCompileTime,
117 
118  add_secondary_variable("epsilon",
120  DisplacementDim>::RowsAtCompileTime,
122 
123  //
124  // enable output of internal variables defined by material models
125  //
127  LocalAssemblerInterface>(_process_data.solid_materials,
128  add_secondary_variable);
129 
130  // Set initial conditions for integration point data.
131  for (auto const& ip_writer : _integration_point_writer)
132  {
133  // Find the mesh property with integration point writer's name.
134  auto const& name = ip_writer->name();
135  if (!mesh.getProperties().existsPropertyVector<double>(name))
136  {
137  continue;
138  }
139  auto const& mesh_property =
140  *mesh.getProperties().template getPropertyVector<double>(name);
141 
142  // The mesh property must be defined on integration points.
143  if (mesh_property.getMeshItemType() !=
145  {
146  continue;
147  }
148 
149  auto const ip_meta_data = getIntegrationPointMetaData(mesh, name);
150 
151  // Check the number of components.
152  if (ip_meta_data.n_components !=
153  mesh_property.getNumberOfGlobalComponents())
154  {
155  OGS_FATAL(
156  "Different number of components in meta data ({:d}) than in "
157  "the integration point field data for '{:s}': {:d}.",
158  ip_meta_data.n_components, name,
159  mesh_property.getNumberOfGlobalComponents());
160  }
161 
162  // Now we have a properly named vtk's field data array and the
163  // corresponding meta data.
164  std::size_t position = 0;
165  for (auto& local_asm : _local_assemblers)
166  {
167  std::size_t const integration_points_read =
168  local_asm->setIPDataInitialConditions(
169  name, &mesh_property[position],
170  ip_meta_data.integration_order);
171  if (integration_points_read == 0)
172  {
173  OGS_FATAL(
174  "No integration points read in the integration point "
175  "initial conditions set function.");
176  }
177  position += integration_points_read * ip_meta_data.n_components;
178  }
179  }
180 
181  // Initialize local assemblers after all variables have been set.
185 }
#define OGS_FATAL(...)
Definition: Error.h:26
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
SmallDeformationLocalAssemblerInterface< DisplacementDim > LocalAssemblerInterface
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
void solidMaterialInternalToSecondaryVariables(std::map< int, std::unique_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim >>> 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, ExtraCtorArgs &&... extra_ctor_args)
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Mesh const &mesh, std::string const &name)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual std::vector< double > const & getIntPtFreeEnergyDensity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0

References ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Properties::existsPropertyVector(), MeshLib::Mesh::getElements(), ProcessLib::getIntegrationPointMetaData(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::IntegrationPoint, MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), MaterialPropertyLib::name, OGS_FATAL, and ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables().

◆ isLinear()

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

Definition at line 80 of file SmallDeformationProcess.cpp.

81 {
82  return false;
83 }

◆ postTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 231 of file SmallDeformationProcess.cpp.

234 {
235  DBUG("PostTimestep SmallDeformationProcess.");
236  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
237  dof_tables.reserve(x.size());
238  std::generate_n(std::back_inserter(dof_tables), x.size(),
239  [&]() { return _local_to_global_index_map.get(); });
240 
241  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
244  pv.getActiveElementIDs(), dof_tables, x, t, dt);
245 
246  std::unique_ptr<GlobalVector> material_forces;
249  *x[process_id]);
250 
251  material_forces->copyValues(*_material_forces);
252 }
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
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

◆ _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: