OGS
ProcessLib::Deformation Namespace Reference

Functions

template<int DisplacementDim, int NPOINTS, typename DNDX_Type >
double divergence (const Eigen::Ref< Eigen::Matrix< double, NPOINTS *DisplacementDim, 1 > const > &u, DNDX_Type const &dNdx)
 Divergence of displacement, the volumetric strain. More...
 
template<int DisplacementDim, int NPOINTS, typename N_Type , typename DNDX_Type , typename GMatrixType >
void computeGMatrix (DNDX_Type const &dNdx, GMatrixType &g_matrix, const bool is_axially_symmetric, N_Type const &N, const double radius)
 Fills a G-matrix based on given shape function dN/dx values. More...
 
template<typename LocalAssemblerInterface , typename AddSecondaryVariableCallback , int DisplacementDim>
void solidMaterialInternalToSecondaryVariables (std::map< int, std::unique_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim >>> const &solid_materials, AddSecondaryVariableCallback const &add_secondary_variable)
 
template<typename LocalAssemblerInterface , typename IntegrationPointWriter , int DisplacementDim>
void solidMaterialInternalVariablesToIntegrationPointWriter (std::map< int, std::unique_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim >>> const &solid_materials, std::vector< std::unique_ptr< LocalAssemblerInterface >> const &local_assemblers, std::vector< std::unique_ptr< IntegrationPointWriter >> &integration_point_writer, int const integration_order)
 

Function Documentation

◆ computeGMatrix()

template<int DisplacementDim, int NPOINTS, typename N_Type , typename DNDX_Type , typename GMatrixType >
void ProcessLib::Deformation::computeGMatrix ( DNDX_Type const &  dNdx,
GMatrixType &  g_matrix,
const bool  is_axially_symmetric,
N_Type const &  N,
const double  radius 
)

Fills a G-matrix based on given shape function dN/dx values.

Definition at line 25 of file GMatrix.h.

30 {
31  static_assert(0 < DisplacementDim && DisplacementDim <= 3,
32  "LinearGMatrix::computeGMatrix: DisplacementDim must be in "
33  "range [1,3].");
34 
35  g_matrix.setZero();
36 
37  switch (DisplacementDim)
38  {
39  case 3:
40  // The gradient coordinates are organized in the following order:
41  // (1,1), (1,2), (1,3)
42  // (2,1), (2,2), (2,3)
43  // (3,1), (3,2), (3,3)
44  for (int d = 0; d < DisplacementDim; ++d)
45  {
46  for (int i = 0; i < NPOINTS; ++i)
47  {
48  g_matrix(d + 0 * DisplacementDim, i + 0 * NPOINTS) =
49  dNdx(d, i);
50  g_matrix(d + 1 * DisplacementDim, i + 1 * NPOINTS) =
51  dNdx(d, i);
52  g_matrix(d + 2 * DisplacementDim, i + 2 * NPOINTS) =
53  dNdx(d, i);
54  }
55  }
56  break;
57  case 2:
58  // The gradient coordinates are organized in the following order:
59  // (1,1), (1,2)
60  // (2,1), (2,2)
61  // (3,3)
62  for (int d = 0; d < DisplacementDim; ++d)
63  {
64  for (int i = 0; i < NPOINTS; ++i)
65  {
66  g_matrix(d, i) = dNdx(d, i);
67  g_matrix(d + DisplacementDim, i + NPOINTS) = dNdx(d, i);
68  }
69  }
70  if (is_axially_symmetric)
71  {
72  for (int i = 0; i < NPOINTS; ++i)
73  {
74  g_matrix(4, i) = N[i] / radius;
75  }
76  }
77  break;
78  default:
79  break;
80  }
81 }

◆ divergence()

template<int DisplacementDim, int NPOINTS, typename DNDX_Type >
double ProcessLib::Deformation::divergence ( const Eigen::Ref< Eigen::Matrix< double, NPOINTS *DisplacementDim, 1 > const > &  u,
DNDX_Type const &  dNdx 
)

Divergence of displacement, the volumetric strain.

Definition at line 19 of file Divergence.h.

23 {
24  double divergence = 0;
25  for (int i = 0; i < DisplacementDim; ++i)
26  {
27  divergence += dNdx.template block<1, NPOINTS>(i, 0) *
28  u.template segment<NPOINTS>(i * NPOINTS);
29  }
30  return divergence;
31 }
double divergence(const Eigen::Ref< Eigen::Matrix< double, NPOINTS *DisplacementDim, 1 > const > &u, DNDX_Type const &dNdx)
Divergence of displacement, the volumetric strain.
Definition: Divergence.h:19

Referenced by ProcessLib::SmallDeformationNonlocal::SmallDeformationNonlocalLocalAssembler< ShapeFunction, IntegrationMethod, DisplacementDim >::computeCrackIntegral().

◆ solidMaterialInternalToSecondaryVariables()

template<typename LocalAssemblerInterface , typename AddSecondaryVariableCallback , int DisplacementDim>
void ProcessLib::Deformation::solidMaterialInternalToSecondaryVariables ( std::map< int, std::unique_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim >>> const &  solid_materials,
AddSecondaryVariableCallback const &  add_secondary_variable 
)

Definition at line 22 of file SolidMaterialInternalToSecondaryVariables.h.

26 {
27  // Collect the internal variables for all solid materials.
28  std::vector<typename MaterialLib::Solids::MechanicsBase<
29  DisplacementDim>::InternalVariable>
30  internal_variables;
31  for (auto const& material_id__solid_material : solid_materials)
32  {
33  auto const variables =
34  material_id__solid_material.second->getInternalVariables();
35  copy(begin(variables), end(variables),
36  back_inserter(internal_variables));
37  }
38 
39  // Register the internal variables.
40  for (auto const& internal_variable : internal_variables)
41  {
42  auto const& name = internal_variable.name;
43  auto const& fct = internal_variable.getter;
44  auto const num_components = internal_variable.num_components;
45  DBUG("Registering internal variable {:s}.", name);
46 
47  auto getIntPtValues =
48  [fct, num_components](
49  LocalAssemblerInterface const& loc_asm,
50  const double /*t*/,
51  std::vector<GlobalVector*> const& /*x*/,
52  std::vector<
53  NumLib::LocalToGlobalIndexMap const*> const& /*dof_table*/,
54  std::vector<double>& cache) -> std::vector<double> const& {
55  const unsigned num_int_pts = loc_asm.getNumberOfIntegrationPoints();
56 
57  cache.clear();
58  auto cache_mat = MathLib::createZeroedMatrix<Eigen::Matrix<
59  double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>>(
60  cache, num_components, num_int_pts);
61 
62  // TODO avoid the heap allocation (one per finite element)
63  std::vector<double> cache_column(num_int_pts);
64 
65  for (unsigned i = 0; i < num_int_pts; ++i)
66  {
67  auto const& state = loc_asm.getMaterialStateVariablesAt(i);
68 
69  auto const& int_pt_values = fct(state, cache_column);
70  assert(int_pt_values.size() ==
71  static_cast<std::size_t>(num_components));
72  auto const int_pt_values_vec = MathLib::toVector(int_pt_values);
73 
74  cache_mat.col(i).noalias() = int_pt_values_vec;
75  }
76 
77  return cache;
78  };
79 
80  add_secondary_variable(name, num_components, std::move(getIntPtValues));
81  }
82 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32
virtual std::vector< InternalVariable > getInternalVariables() const

References MathLib::LinAlg::copy(), MathLib::createZeroedMatrix(), DBUG(), MaterialLib::Solids::MechanicsBase< DisplacementDim >::getInternalVariables(), MaterialPropertyLib::name, and MathLib::toVector().

Referenced by ProcessLib::HydroMechanics::HydroMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::SmallDeformation::SmallDeformationProcess< DisplacementDim >::initializeConcreteProcess(), ProcessLib::ThermoMechanics::ThermoMechanicsProcess< DisplacementDim >::initializeConcreteProcess(), and ProcessLib::ThermoRichardsMechanics::ThermoRichardsMechanicsProcess< DisplacementDim >::initializeConcreteProcess().

◆ solidMaterialInternalVariablesToIntegrationPointWriter()

template<typename LocalAssemblerInterface , typename IntegrationPointWriter , int DisplacementDim>
void ProcessLib::Deformation::solidMaterialInternalVariablesToIntegrationPointWriter ( std::map< int, std::unique_ptr< MaterialLib::Solids::MechanicsBase< DisplacementDim >>> const &  solid_materials,
std::vector< std::unique_ptr< LocalAssemblerInterface >> const &  local_assemblers,
std::vector< std::unique_ptr< IntegrationPointWriter >> &  integration_point_writer,
int const  integration_order 
)

Definition at line 86 of file SolidMaterialInternalToSecondaryVariables.h.

94 {
95  // Collect the internal variables for all solid materials.
96  std::vector<typename MaterialLib::Solids::MechanicsBase<
97  DisplacementDim>::InternalVariable>
98  internal_variables;
99  for (auto const& solid_material : solid_materials)
100  {
101  auto const ivs = solid_material.second->getInternalVariables();
102  copy(begin(ivs), end(ivs), back_inserter(internal_variables));
103  }
104 
105  // Create integration point writers for each of the internal variables.
106  for (auto const& iv : internal_variables)
107  {
108  DBUG("Creating integration point writer for internal variable {:s}.",
109  iv.name);
110 
111  integration_point_writer.emplace_back(
112  std::make_unique<IntegrationPointWriter>(
113  "material_state_variable_" + iv.name + "_ip", iv.num_components,
114  integration_order, local_assemblers,
115  &LocalAssemblerInterface::getMaterialStateVariableInternalState,
116  iv.reference, iv.num_components));
117  }
118 }

References MathLib::LinAlg::copy(), DBUG(), and MaterialLib::Solids::MechanicsBase< DisplacementDim >::getInternalVariables().

Referenced by ProcessLib::RichardsMechanics::RichardsMechanicsProcess< DisplacementDim >::initializeConcreteProcess().