OGS
ThermoMechanicsProcess.cpp
Go to the documentation of this file.
1 
11 #include "ThermoMechanicsProcess.h"
12 
13 #include <cassert>
14 
20 #include "ThermoMechanicsFEM.h"
21 
22 namespace ProcessLib
23 {
24 namespace ThermoMechanics
25 {
26 template <int DisplacementDim>
28  std::string name, MeshLib::Mesh& mesh,
29  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
30  std::vector<std::unique_ptr<ParameterLib::ParameterBase>> const& parameters,
31  unsigned const integration_order,
32  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
33  process_variables,
35  SecondaryVariableCollection&& secondary_variables,
36  bool const use_monolithic_scheme)
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 }
69 
70 template <int DisplacementDim>
72 {
73  return false;
74 }
75 
76 template <int DisplacementDim>
79  const int process_id) const
80 {
81  // For the monolithic scheme or the M process (deformation) in the staggered
82  // scheme.
83  if (_use_monolithic_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.
92  auto const& l = *_local_to_global_index_map_single_component;
93  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
94  &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
95 }
96 
97 // TODO [WW]: remove if (_use_monolithic_scheme) during the refactoring of the
98 // coupling part.
99 template <int DisplacementDim>
101 {
102  // Note: the heat conduction process and the mechanical process use the same
103  // order of shape functions.
104 
105  if (_use_monolithic_scheme)
106  {
107  constructMonolithicProcessDofTable();
108  return;
109  }
110  constructDofTableOfSpecifiedProcessStaggeredScheme(
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{
116  *_mesh_subset_all_nodes};
117  _local_to_global_index_map_single_component.reset(
119  std::move(all_mesh_subsets_single_component),
120  // by location order is needed for output
122 
123  if (!_use_monolithic_scheme)
124  {
125  _sparsity_pattern_with_single_component =
127  *_local_to_global_index_map_single_component, _mesh);
128  }
129 }
130 
131 template <int DisplacementDim>
133  NumLib::LocalToGlobalIndexMap const& dof_table,
134  MeshLib::Mesh const& mesh,
135  unsigned const integration_order)
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  {
146  _secondary_variables.addSecondaryVariable(
147  name,
148  makeExtrapolator(num_components, getExtrapolator(),
149  _local_assemblers,
150  std::move(get_ip_values_function)));
151  };
152 
153  add_secondary_variable("sigma",
155  DisplacementDim>::RowsAtCompileTime,
156  &LocalAssemblerInterface::getIntPtSigma);
157 
158  add_secondary_variable("epsilon",
160  DisplacementDim>::RowsAtCompileTime,
161  &LocalAssemblerInterface::getIntPtEpsilon);
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.
223  &LocalAssemblerInterface::initialize, _local_assemblers,
224  *_local_to_global_index_map);
225 }
226 
227 template <int DisplacementDim>
229 {
230  if (_use_monolithic_scheme)
231  {
232  const int process_id_of_thermomechanics = 0;
233  initializeProcessBoundaryConditionsAndSourceTerms(
234  *_local_to_global_index_map, process_id_of_thermomechanics);
235  return;
236  }
237 
238  // Staggered scheme:
239  // for the equations of heat conduction
240  initializeProcessBoundaryConditionsAndSourceTerms(
241  *_local_to_global_index_map_single_component,
242  _process_data.heat_conduction_process_id);
243 
244  // for the equations of deformation.
245  initializeProcessBoundaryConditionsAndSourceTerms(
246  *_local_to_global_index_map, _process_data.mechanics_process_id);
247 }
248 
249 template <int DisplacementDim>
251  const double t, double const dt, std::vector<GlobalVector*> const& x,
252  std::vector<GlobalVector*> const& xdot, int const process_id,
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.
263  _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers,
264  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
265  b);
266 }
267 
268 template <int DisplacementDim>
271  const double t, double const dt, std::vector<GlobalVector*> const& x,
272  std::vector<GlobalVector*> const& xdot, const double dxdot_dx,
273  const double dx_dx, int const process_id, GlobalMatrix& M,
275 {
276  DBUG("AssembleJacobian ThermoMechanicsProcess.");
277 
278  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
279  dof_tables;
280  // For the monolithic scheme
281  if (_use_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(
309  *_local_to_global_index_map_single_component);
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(
316  *_local_to_global_index_map_single_component);
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  {
330  if (_use_monolithic_scheme)
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  };
343  if (_use_monolithic_scheme ||
344  process_id == _process_data.heat_conduction_process_id)
345  {
346  copyRhs(0, *_heat_flux);
347  }
348  if (_use_monolithic_scheme ||
349  process_id == _process_data.mechanics_process_id)
350  {
351  copyRhs(1, *_nodal_forces);
352  }
353 }
354 
355 template <int DisplacementDim>
357  std::vector<GlobalVector*> const& x, double const t, double const dt,
358  const int process_id)
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  {
369  &LocalAssemblerInterface::preTimestep, _local_assemblers,
370  pv.getActiveElementIDs(), *_local_to_global_index_map,
371  *x[process_id], t, dt);
372  return;
373  }
374 }
375 
376 template <int DisplacementDim>
378  std::vector<GlobalVector*> const& x, double const t, double const dt,
379  int const process_id)
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];
397  &LocalAssemblerInterface::postTimestep, _local_assemblers,
398  pv.getActiveElementIDs(), dof_tables, x, t, dt);
399 }
400 
401 template <int DisplacementDim>
404 {
405  if (_process_data.mechanics_process_id == process_id)
406  {
407  return *_local_to_global_index_map;
408  }
409 
410  // For the equation of pressure
411  return *_local_to_global_index_map_single_component;
412 }
413 
414 template class ThermoMechanicsProcess<2>;
415 template class ThermoMechanicsProcess<3>;
416 
417 } // namespace ThermoMechanics
418 } // namespace ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
Global vector based on Eigen vector.
Definition: EigenVector.h:26
bool isAxiallySymmetric() const
Definition: Mesh.h:126
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition: Mesh.h:71
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:98
Properties & getProperties()
Definition: Mesh.h:123
bool existsPropertyVector(std::string const &name) const
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)
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)
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:350
Handles configuration of several secondary variables from the project file.
Local assembler of ThermoMechanics process.
void postTimestepConcreteProcess(std::vector< GlobalVector * > const &x, const double t, const double dt, int const process_id) override
void preTimestepConcreteProcess(std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
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)
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
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const 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 initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
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 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)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
@ 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.
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
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 executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > getEpsilonMechanical() const =0