OGS 6.1.0-1699-ge946d4c5f
SmallDeformationProcess.cpp
Go to the documentation of this file.
1 
11 
12 #include <cassert>
13 #include <nlohmann/json.hpp>
14 
15 #include "BaseLib/Functional.h"
16 #include "ProcessLib/Process.h"
18 
19 #include "SmallDeformationFEM.h"
20 
21 namespace ProcessLib
22 {
23 namespace SmallDeformation
24 {
25 template <int DisplacementDim>
27  MeshLib::Mesh& mesh,
28  std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
29  std::vector<std::unique_ptr<ParameterBase>> const& parameters,
30  unsigned const integration_order,
31  std::vector<std::vector<std::reference_wrapper<ProcessVariable>>>&&
32  process_variables,
34  SecondaryVariableCollection&& secondary_variables,
35  NumLib::NamedFunctionCaller&& named_function_caller)
36  : Process(mesh, std::move(jacobian_assembler), parameters,
37  integration_order, std::move(process_variables),
38  std::move(secondary_variables), std::move(named_function_caller)),
39  _process_data(std::move(process_data))
40 {
41  _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
42  mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
43 
44  _material_forces = MeshLib::getOrCreateMeshProperty<double>(
45  mesh, "MaterialForces", MeshLib::MeshItemType::Node, DisplacementDim);
46 
47  _integration_point_writer.emplace_back(
48  std::make_unique<SigmaIntegrationPointWriter>(
49  static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
50  2 /*integration order*/, [this]() {
51  // Result containing integration point data for each local
52  // assembler.
53  std::vector<std::vector<double>> result;
54  result.resize(_local_assemblers.size());
55 
56  for (std::size_t i = 0; i < _local_assemblers.size(); ++i)
57  {
58  auto const& local_asm = *_local_assemblers[i];
59 
60  result[i] = local_asm.getSigma();
61  }
62 
63  return result;
64  }));
65 }
66 
67 template <int DisplacementDim>
69 {
70  return false;
71 }
72 
73 template <int DisplacementDim>
75  NumLib::LocalToGlobalIndexMap const& dof_table,
76  MeshLib::Mesh const& mesh,
77  unsigned const integration_order)
78 {
79  using nlohmann::json;
80 
82  DisplacementDim, SmallDeformationLocalAssembler>(
83  mesh.getElements(), dof_table, _local_assemblers,
84  mesh.isAxiallySymmetric(), integration_order, _process_data);
85 
86  // TODO move the two data members somewhere else.
87  // for extrapolation of secondary variables
88  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
91  std::make_unique<NumLib::LocalToGlobalIndexMap>(
92  std::move(all_mesh_subsets_single_component),
93  // by location order is needed for output
95 
97  "free_energy_density",
98  makeExtrapolator(1, getExtrapolator(), _local_assemblers,
100 
102  "sigma",
104  DisplacementDim>::RowsAtCompileTime,
105  getExtrapolator(), _local_assemblers,
107 
109  "epsilon",
111  DisplacementDim>::RowsAtCompileTime,
112  getExtrapolator(), _local_assemblers,
114 
115  //
116  // enable output of internal variables defined by material models
117  //
118 
119  // Collect the internal variables for all solid materials.
120  std::vector<typename MaterialLib::Solids::MechanicsBase<
121  DisplacementDim>::InternalVariable>
122  internal_variables;
123  for (auto const& material_id__solid_material :
124  _process_data.solid_materials)
125  {
126  auto const variables =
127  material_id__solid_material.second->getInternalVariables();
128  copy(begin(variables), end(variables),
129  back_inserter(internal_variables));
130  }
131 
132  // Register the internal variables.
133  for (auto const& internal_variable : internal_variables)
134  {
135  auto const& name = internal_variable.name;
136  auto const& fct = internal_variable.getter;
137  auto const num_components = internal_variable.num_components;
138  DBUG("Registering internal variable %s.", name.c_str());
139 
140  auto getIntPtValues = BaseLib::easyBind(
141  [fct, num_components](
142  LocalAssemblerInterface const& loc_asm,
143  const double /*t*/,
144  GlobalVector const& /*current_solution*/,
145  NumLib::LocalToGlobalIndexMap const& /*dof_table*/,
146  std::vector<double>& cache) -> std::vector<double> const& {
147  const unsigned num_int_pts =
149 
150  cache.clear();
151  auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
152  double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>(
153  cache, num_components, num_int_pts);
154 
155  // TODO avoid the heap allocation (one per finite element)
156  std::vector<double> cache_column(num_int_pts);
157 
158  for (unsigned i = 0; i < num_int_pts; ++i)
159  {
160  auto const& state = loc_asm.getMaterialStateVariablesAt(i);
161 
162  auto const& int_pt_values = fct(state, cache_column);
163  assert(int_pt_values.size() == num_components);
164  auto const int_pt_values_vec =
165  MathLib::toVector(int_pt_values);
166 
167  cache_mat.col(i).noalias() = int_pt_values_vec;
168  }
169 
170  return cache;
171  });
172 
174  name,
175  makeExtrapolator(num_components, getExtrapolator(),
176  _local_assemblers, std::move(getIntPtValues)));
177  }
178 
179  // Set initial conditions for integration point data.
180  for (auto const& ip_writer : _integration_point_writer)
181  {
182  // Find the mesh property with integration point writer's name.
183  auto const& name = ip_writer->name();
184  if (!mesh.getProperties().existsPropertyVector<double>(name))
185  {
186  continue;
187  }
188  auto const& mesh_property =
189  *mesh.getProperties().template getPropertyVector<double>(name);
190 
191  // The mesh property must be defined on integration points.
192  if (mesh_property.getMeshItemType() !=
194  {
195  continue;
196  }
197 
198  auto const ip_meta_data = getIntegrationPointMetaData(mesh, name);
199 
200  // Check the number of components.
201  if (ip_meta_data.n_components != mesh_property.getNumberOfComponents())
202  {
203  OGS_FATAL(
204  "Different number of components in meta data (%d) than in "
205  "the integration point field data for '%s': %d.",
206  ip_meta_data.n_components, name.c_str(),
207  mesh_property.getNumberOfComponents());
208  }
209 
210  // Now we have a properly named vtk's field data array and the
211  // corresponding meta data.
212  std::size_t position = 0;
213  for (auto& local_asm : _local_assemblers)
214  {
215  std::size_t const integration_points_read =
216  local_asm->setIPDataInitialConditions(
217  name, &mesh_property[position],
218  ip_meta_data.integration_order);
219  if (integration_points_read == 0)
220  {
221  OGS_FATAL(
222  "No integration points read in the integration point "
223  "initial conditions set function.");
224  }
225  position += integration_points_read * ip_meta_data.n_components;
226  }
227  }
228 }
229 
230 template <int DisplacementDim>
232  const double t, GlobalVector const& x, GlobalMatrix& M, GlobalMatrix& K,
233  GlobalVector& b)
234 {
235  DBUG("Assemble SmallDeformationProcess.");
236 
237  const int process_id = 0;
238  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
239 
240  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
241  dof_table = {std::ref(*_local_to_global_index_map)};
242  // Call global assembler for each local assembly item.
245  pv.getActiveElementIDs(), dof_table, t, x, M, K, b,
247 }
248 
249 template <int DisplacementDim>
251  assembleWithJacobianConcreteProcess(const double t, GlobalVector const& x,
252  GlobalVector const& xdot,
253  const double dxdot_dx,
254  const double dx_dx, GlobalMatrix& M,
255  GlobalMatrix& K, GlobalVector& b,
256  GlobalMatrix& Jac)
257 {
258  DBUG("AssembleWithJacobian SmallDeformationProcess.");
259 
260  const int process_id = 0;
261  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
262 
263  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
264  dof_table = {std::ref(*_local_to_global_index_map)};
265  // Call global assembler for each local assembly item.
268  _local_assemblers, pv.getActiveElementIDs(), dof_table, t, x,
269  xdot, dxdot_dx, dx_dx, M, K, b, Jac, _coupled_solutions);
270 
272  *_nodal_forces, std::negate<double>());
273 }
274 
275 template <int DisplacementDim>
277  GlobalVector const& x, double const t, double const dt,
278  const int process_id)
279 {
280  DBUG("PreTimestep SmallDeformationProcess.");
281 
282  _process_data.dt = dt;
283  _process_data.t = t;
284 
285  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
286 
290  x, t, dt);
291 }
292 
293 template <int DisplacementDim>
295  GlobalVector const& x, const double /*t*/, const double /*delta_t*/,
296  int const process_id)
297 {
298  DBUG("PostTimestep SmallDeformationProcess.");
299 
300  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
301 
305 
306  std::unique_ptr<GlobalVector> material_forces;
308  material_forces, _local_assemblers, *_local_to_global_index_map, x);
309 
310  material_forces->copyValues(*_material_forces);
311 }
312 
313 template class SmallDeformationProcess<2>;
314 template class SmallDeformationProcess<3>;
315 
316 } // namespace SmallDeformation
317 } // namespace ProcessLib
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Mesh const &mesh, std::string const &name)
SmallDeformationProcess(MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterBase >> const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable >>> &&process_variables, SmallDeformationProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, NumLib::NamedFunctionCaller &&named_function_caller)
virtual std::vector< InternalVariable > getInternalVariables() const
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:262
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:123
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
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, GlobalVector const &x, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, CoupledSolutionsForStaggeredScheme const *const cpl_xs)
MeshLib::Properties & getProperties()
Definition: Mesh.h:131
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)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:150
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, GlobalVector const &x, GlobalVector const &xdot, const double dxdot_dx, const double dx_dx, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac, CoupledSolutionsForStaggeredScheme const *const cpl_xs)
bool isAxiallySymmetric() const
Definition: Mesh.h:134
void preTimestepConcreteProcess(GlobalVector const &x, double const t, double const dt, const int process_id) override
virtual std::vector< double > const & getIntPtFreeEnergyDensity(const double, GlobalVector const &, NumLib::LocalToGlobalIndexMap const &, std::vector< double > &cache) const =0
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:31
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:264
Builds expression trees of named functions dynamically at runtime.
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 std::vector< double > const & getIntPtEpsilon(const double, GlobalVector const &, NumLib::LocalToGlobalIndexMap const &, std::vector< double > &cache) const =0
void assembleConcreteProcess(const double t, GlobalVector const &x, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
SecondaryVariableCollection _secondary_variables
Definition: Process.h:266
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
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)
std::enable_if< std::is_same< MethodClass, typename std::remove_cv< typename std::remove_pointer< typename std::decay< Object >::type >::type >::type >::value, std::function< ReturnType(Args...)> >::type easyBind(ReturnType(MethodClass::*method)(Args...), Object &&obj)
Definition: Functional.h:180
virtual MaterialLib::Solids::MechanicsBase< DisplacementDim >::MaterialStateVariables const & getMaterialStateVariablesAt(unsigned) const =0
void postTimestepConcreteProcess(GlobalVector const &x, const double t, const double delta_t, int const process_id) override
SmallDeformationProcessData< DisplacementDim > _process_data
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition: Mesh.h:108
bool existsPropertyVector(std::string const &name) const
Definition: Properties.h:79
virtual void postTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x)
virtual std::vector< double > const & getIntPtSigma(const double, GlobalVector const &, NumLib::LocalToGlobalIndexMap const &, std::vector< double > &cache) const =0
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
CoupledSolutionsForStaggeredScheme * _coupled_solutions
Definition: Process.h:279
void initializeConcreteProcess(NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
Process specific initialization called by initialize().
Handles configuration of several secondary variables from the project file.
Ordering data by spatial location.
#define OGS_FATAL(fmt,...)
Definition: Error.h:71
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
void copy(MatrixOrVector const &x, MatrixOrVector &y)
Copies x to y.
Definition: LinAlg.h:36
void assembleWithJacobianConcreteProcess(const double t, GlobalVector const &x, GlobalVector const &xdot, const double dxdot_dx, const double dx_dx, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
std::vector< std::unique_ptr< IntegrationPointWriter > > _integration_point_writer
Definition: Process.h:290
VectorMatrixAssembler _global_assembler
Definition: Process.h:273
Eigen::Matrix< double, KelvinVectorDimensions< DisplacementDim >::value, 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:49
std::vector< std::size_t > & getActiveElementIDs() const