OGS
ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >

Definition at line 23 of file ThermoMechanicsProcess.h.

#include <ThermoMechanicsProcess.h>

Inheritance diagram for ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >:
[legend]

Public Member Functions

 ThermoMechanicsProcess (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, ThermoMechanicsProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const 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, 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
 
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 = ThermoMechanicsLocalAssemblerInterface< DisplacementDim >
 

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 initializeBoundaryConditions () override
 
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 preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
 
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override
 

Private Attributes

ThermoMechanicsProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
GlobalSparsityPattern _sparsity_pattern_with_single_component
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 
MeshLib::PropertyVector< double > * _heat_flux = 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)
 
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::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::LocalAssemblerInterface = ThermoMechanicsLocalAssemblerInterface<DisplacementDim>
private

Definition at line 58 of file ThermoMechanicsProcess.h.

Constructor & Destructor Documentation

◆ ThermoMechanicsProcess()

template<int DisplacementDim>
ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::ThermoMechanicsProcess ( 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,
ThermoMechanicsProcessData< DisplacementDim > &&  process_data,
SecondaryVariableCollection &&  secondary_variables,
bool const  use_monolithic_scheme 
)

Definition at line 27 of file ThermoMechanicsProcess.cpp.

37  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
38  integration_order, std::move(process_variables),
39  std::move(secondary_variables), use_monolithic_scheme),
40  _process_data(std::move(process_data))
41 {
42  _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
43  mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
44 
45  _heat_flux = MeshLib::getOrCreateMeshProperty<double>(
46  mesh, "HeatFlux", MeshLib::MeshItemType::Node, 1);
47 
48  _integration_point_writer.emplace_back(
49  std::make_unique<IntegrationPointWriter>(
50  "sigma_ip",
51  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
52  integration_order, _local_assemblers,
54 
55  _integration_point_writer.emplace_back(
56  std::make_unique<IntegrationPointWriter>(
57  "epsilon_ip",
58  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
59  integration_order, _local_assemblers,
61 
62  _integration_point_writer.emplace_back(
63  std::make_unique<IntegrationPointWriter>(
64  "epsilon_m_ip",
65  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
66  integration_order, _local_assemblers,
68 }
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
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
ThermoMechanicsProcessData< DisplacementDim > _process_data
virtual std::vector< double > getEpsilonMechanical() const =0

References ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_heat_flux, ProcessLib::Process::_integration_point_writer, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_local_assemblers, ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_nodal_forces, MeshLib::Mesh::getDimension(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >::getEpsilon(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >::getEpsilonMechanical(), ProcessLib::ThermoMechanics::ThermoMechanicsLocalAssemblerInterface< DisplacementDim >::getSigma(), and MeshLib::Node.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< 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 250 of file ThermoMechanicsProcess.cpp.

254 {
255  DBUG("Assemble ThermoMechanicsProcess.");
256 
257  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
258  dof_table = {std::ref(*_local_to_global_index_map)};
259  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
260 
261  // Call global assembler for each local assembly item.
264  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
265  b);
266 }
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::ThermoMechanics::ThermoMechanicsProcess< 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 269 of file ThermoMechanicsProcess.cpp.

275 {
276  DBUG("AssembleJacobian ThermoMechanicsProcess.");
277 
278  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
279  dof_tables;
280  // For the monolithic scheme
282  {
283  DBUG(
284  "Assemble the Jacobian of ThermoMechanics for the monolithic"
285  " scheme.");
286  dof_tables.emplace_back(*_local_to_global_index_map);
287  }
288  else
289  {
290  // For the staggered scheme
291  if (process_id == _process_data.heat_conduction_process_id)
292  {
293  DBUG(
294  "Assemble the Jacobian equations of heat conduction process in "
295  "ThermoMechanics for the staggered scheme.");
296  }
297  else
298  {
299  DBUG(
300  "Assemble the Jacobian equations of mechanical process in "
301  "ThermoMechanics for the staggered scheme.");
302  }
303 
304  // For the flexible appearance order of processes in the coupling.
305  if (_process_data.heat_conduction_process_id ==
306  0) // First: the heat conduction process
307  {
308  dof_tables.emplace_back(
310  dof_tables.emplace_back(*_local_to_global_index_map);
311  }
312  else // vice versa
313  {
314  dof_tables.emplace_back(*_local_to_global_index_map);
315  dof_tables.emplace_back(
317  }
318  }
319 
320  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
321 
324  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot,
325  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
326 
327  // TODO (naumov): Refactor the copy rhs part. This is copy from HM.
328  auto copyRhs = [&](int const variable_id, auto& output_vector)
329  {
331  {
332  transformVariableFromGlobalVector(b, variable_id, dof_tables[0],
333  output_vector,
334  std::negate<double>());
335  }
336  else
337  {
338  transformVariableFromGlobalVector(b, 0, dof_tables[process_id],
339  output_vector,
340  std::negate<double>());
341  }
342  };
344  process_id == _process_data.heat_conduction_process_id)
345  {
346  copyRhs(0, *_heat_flux);
347  }
349  process_id == _process_data.mechanics_process_id)
350  {
351  copyRhs(1, *_nodal_forces);
352  }
353 }
const bool _use_monolithic_scheme
Definition: Process.h:335
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
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().

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< 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 100 of file ThermoMechanicsProcess.cpp.

101 {
102  // Note: the heat conduction process and the mechanical process use the same
103  // order of shape functions.
104 
106  {
108  return;
109  }
111  _process_data.mechanics_process_id);
112 
113  // TODO move the two data members somewhere else.
114  // for extrapolation of secondary variables of stress or strain
115  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
119  std::move(all_mesh_subsets_single_component),
120  // by location order is needed for output
122 
124  {
128  }
129 }
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_prosess_id)
Definition: Process.cpp:303
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:327
MeshLib::Mesh & _mesh
Definition: Process.h:326
void constructMonolithicProcessDofTable()
Definition: Process.cpp:270
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References NumLib::BY_LOCATION, and NumLib::computeSparsityPattern().

◆ getDOFTable()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::getDOFTable ( const int  process_id) const
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 403 of file ThermoMechanicsProcess.cpp.

404 {
405  if (_process_data.mechanics_process_id == process_id)
406  {
408  }
409 
410  // For the equation of pressure
412 }

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::getMatrixSpecifications ( const int  process_id) const
override

Get the size and the sparse pattern of the global matrix in order to create the global matrices and vectors for the system equations of this process.

Parameters
process_idProcess ID. If the monolithic scheme is applied, process_id = 0.
Returns
Matrix specifications including size and sparse pattern.

Definition at line 78 of file ThermoMechanicsProcess.cpp.

80 {
81  // For the monolithic scheme or the M process (deformation) in the staggered
82  // scheme.
84  process_id == _process_data.mechanics_process_id)
85  {
86  auto const& l = *_local_to_global_index_map;
87  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
88  &l.getGhostIndices(), &this->_sparsity_pattern};
89  }
90 
91  // For staggered scheme and T process.
93  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
94  &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
95 }
GlobalSparsityPattern _sparsity_pattern
Definition: Process.h:352

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeBoundaryConditions
overrideprivatevirtual

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

Reimplemented from ProcessLib::Process.

Definition at line 228 of file ThermoMechanicsProcess.cpp.

229 {
231  {
232  const int process_id_of_thermomechanics = 0;
234  *_local_to_global_index_map, process_id_of_thermomechanics);
235  return;
236  }
237 
238  // Staggered scheme:
239  // for the equations of heat conduction
242  _process_data.heat_conduction_process_id);
243 
244  // for the equations of deformation.
246  *_local_to_global_index_map, _process_data.mechanics_process_id);
247 }
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:67

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< 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 132 of file ThermoMechanicsProcess.cpp.

136 {
138  DisplacementDim, ThermoMechanicsLocalAssembler>(
139  mesh.getElements(), dof_table, _local_assemblers,
140  mesh.isAxiallySymmetric(), integration_order, _process_data);
141 
142  auto add_secondary_variable = [&](std::string const& name,
143  int const num_components,
144  auto get_ip_values_function)
145  {
147  name,
148  makeExtrapolator(num_components, getExtrapolator(),
150  std::move(get_ip_values_function)));
151  };
152 
153  add_secondary_variable("sigma",
155  DisplacementDim>::RowsAtCompileTime,
157 
158  add_secondary_variable("epsilon",
160  DisplacementDim>::RowsAtCompileTime,
162 
163  //
164  // enable output of internal variables defined by material models
165  //
167  LocalAssemblerInterface>(_process_data.solid_materials,
168  add_secondary_variable);
169 
170  // Set initial conditions for integration point data.
171  for (auto const& ip_writer : _integration_point_writer)
172  {
173  // Find the mesh property with integration point writer's name.
174  auto const& name = ip_writer->name();
175  if (!mesh.getProperties().existsPropertyVector<double>(name))
176  {
177  continue;
178  }
179  auto const& mesh_property =
180  *mesh.getProperties().template getPropertyVector<double>(name);
181 
182  // The mesh property must be defined on integration points.
183  if (mesh_property.getMeshItemType() !=
185  {
186  continue;
187  }
188 
189  auto const ip_meta_data = getIntegrationPointMetaData(mesh, name);
190 
191  // Check the number of components.
192  if (ip_meta_data.n_components !=
193  mesh_property.getNumberOfGlobalComponents())
194  {
195  OGS_FATAL(
196  "Different number of components in meta data ({:d}) than in "
197  "the integration point field data for '{:s}': {:d}.",
198  ip_meta_data.n_components, name,
199  mesh_property.getNumberOfGlobalComponents());
200  }
201 
202  // Now we have a properly named vtk's field data array and the
203  // corresponding meta data.
204  std::size_t position = 0;
205  for (auto& local_asm : _local_assemblers)
206  {
207  std::size_t const integration_points_read =
208  local_asm->setIPDataInitialConditions(
209  name, &mesh_property[position],
210  ip_meta_data.integration_order);
211  if (integration_points_read == 0)
212  {
213  OGS_FATAL(
214  "No integration points read in the integration point "
215  "initial conditions set function.");
216  }
217  position += integration_points_read * ip_meta_data.n_components;
218  }
219  }
220 
221  // Initialize local assemblers after all variables have been set.
225 }
#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)
ThermoMechanicsLocalAssemblerInterface< 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 & 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 & getIntPtEpsilon(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::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::isLinear
override

Definition at line 71 of file ThermoMechanicsProcess.cpp.

72 {
73  return false;
74 }

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< 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 377 of file ThermoMechanicsProcess.cpp.

380 {
381  if (process_id != 0)
382  {
383  return;
384  }
385 
386  DBUG("PostTimestep ThermoMechanicsProcess.");
387  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
388  auto const n_processes = x.size();
389  dof_tables.reserve(n_processes);
390  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
391  {
392  dof_tables.push_back(&getDOFTable(process_id));
393  }
394 
395  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
398  pv.getActiveElementIDs(), dof_tables, x, t, dt);
399 }
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)
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

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

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanics::ThermoMechanicsProcess< 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 356 of file ThermoMechanicsProcess.cpp.

359 {
360  DBUG("PreTimestep ThermoMechanicsProcess.");
361 
362  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
363 
364  assert(process_id < 2);
365 
366  if (process_id == _process_data.mechanics_process_id)
367  {
371  *x[process_id], t, dt);
372  return;
373  }
374 }
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(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Member Data Documentation

◆ _heat_flux

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_heat_flux = nullptr
private

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_local_assemblers
private

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 99 of file ThermoMechanicsProcess.h.

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

template<int DisplacementDim>
ThermoMechanicsProcessData<DisplacementDim> ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_process_data
private

Definition at line 94 of file ThermoMechanicsProcess.h.

◆ _sparsity_pattern_with_single_component

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::_sparsity_pattern_with_single_component
private

Sparsity pattern for the heat conduction equation, and it is initialized only if the staggered scheme is used.

Definition at line 103 of file ThermoMechanicsProcess.h.


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