OGS
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int DisplacementDim>
class ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >

Definition at line 32 of file SmallDeformationLocalAssemblerFracture.h.

#include <SmallDeformationLocalAssemblerFracture.h>

Inheritance diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >:
[legend]

Public Types

using ShapeMatricesType
 
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using HMatricesType = HMatrixPolicyType<ShapeFunction, DisplacementDim>
 
using HMatrixType = typename HMatricesType::HMatrixType
 
using StiffnessMatrixType = typename HMatricesType::StiffnessMatrixType
 
using NodalForceVectorType = typename HMatricesType::NodalForceVectorType
 
using NodalDisplacementVectorType
 
using ForceVectorType = typename HMatricesType::ForceVectorType
 
using GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>
 

Public Member Functions

 SmallDeformationLocalAssemblerFracture (SmallDeformationLocalAssemblerFracture const &)=delete
 
 SmallDeformationLocalAssemblerFracture (SmallDeformationLocalAssemblerFracture &&)=delete
 
 SmallDeformationLocalAssemblerFracture (MeshLib::Element const &e, std::size_t const n_variables, std::size_t const local_matrix_size, std::vector< unsigned > const &dofIndex_to_localIndex, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SmallDeformationProcessData< DisplacementDim > &process_data)
 
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
 
void assembleWithJacobian (double const t, double const dt, Eigen::VectorXd const &local_u, Eigen::VectorXd &local_b, Eigen::MatrixXd &local_J) override
 
void preTimestepConcrete (std::vector< double > const &, double const, double const) override
 
void computeSecondaryVariableConcreteWithVector (const double t, Eigen::VectorXd const &local_u) override
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point.
 
std::vector< double > const & getIntPtSigma (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtEpsilon (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtFractureStress (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
std::vector< double > const & getIntPtFractureAperture (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
 
- Public Member Functions inherited from ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface
 SmallDeformationLocalAssemblerInterface ()
 
 SmallDeformationLocalAssemblerInterface (std::size_t n_local_size, std::vector< unsigned > dofIndex_to_localIndex)
 
void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x_, std::vector< double > const &, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) override
 
void computeSecondaryVariableConcrete (double const t, double const, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
 
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, int const process_id)
 
virtual void initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
 
virtual void preAssemble (double const, double const, std::vector< double > const &)
 
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, int const process_id, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
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_prev, int const process_id)
 
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 postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const
 Fits to staggered scheme.
 
virtual std::optional< VectorSegmentgetVectorDeformationSegment () const
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Types

using IntegrationPointDataType
 

Private Attributes

SmallDeformationProcessData< DisplacementDim > & _process_data
 
std::vector< FractureProperty * > _fracture_props
 
std::vector< JunctionProperty * > _junction_props
 
std::unordered_map< int, int > _fracID_to_local
 
FractureProperty const * _fracture_property = nullptr
 
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
 
NumLib::GenericIntegrationMethod const & _integration_method
 
std::vector< ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > _shape_matrices
 
MeshLib::Element const & _element
 
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
 

Additional Inherited Members

- Protected Member Functions inherited from ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface

Member Typedef Documentation

◆ ForceVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::ForceVectorType = typename HMatricesType::ForceVectorType

Definition at line 49 of file SmallDeformationLocalAssemblerFracture.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::GlobalDimVectorType = Eigen::Matrix<double, DisplacementDim, 1>

Definition at line 50 of file SmallDeformationLocalAssemblerFracture.h.

◆ HMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::HMatricesType = HMatrixPolicyType<ShapeFunction, DisplacementDim>

Definition at line 41 of file SmallDeformationLocalAssemblerFracture.h.

◆ HMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::HMatrixType = typename HMatricesType::HMatrixType

Definition at line 43 of file SmallDeformationLocalAssemblerFracture.h.

◆ IntegrationPointDataType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::IntegrationPointDataType
private
Initial value:
IntegrationPointDataFracture<HMatricesType, DisplacementDim>

Definition at line 147 of file SmallDeformationLocalAssemblerFracture.h.

◆ NodalDisplacementVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::NodalDisplacementVectorType
Initial value:
VectorType< _number_of_dof > NodalForceVectorType

Definition at line 46 of file SmallDeformationLocalAssemblerFracture.h.

◆ NodalForceVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::NodalForceVectorType = typename HMatricesType::NodalForceVectorType

Definition at line 45 of file SmallDeformationLocalAssemblerFracture.h.

◆ NodalMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 38 of file SmallDeformationLocalAssemblerFracture.h.

◆ NodalVectorType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 39 of file SmallDeformationLocalAssemblerFracture.h.

◆ ShapeMatrices

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 40 of file SmallDeformationLocalAssemblerFracture.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::ShapeMatricesType

◆ StiffnessMatrixType

template<typename ShapeFunction , int DisplacementDim>
using ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::StiffnessMatrixType = typename HMatricesType::StiffnessMatrixType

Definition at line 44 of file SmallDeformationLocalAssemblerFracture.h.

Constructor & Destructor Documentation

◆ SmallDeformationLocalAssemblerFracture() [1/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerFracture ( SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim > const & )
delete

◆ SmallDeformationLocalAssemblerFracture() [2/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerFracture ( SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim > && )
delete

◆ SmallDeformationLocalAssemblerFracture() [3/3]

template<typename ShapeFunction , int DisplacementDim>
ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::SmallDeformationLocalAssemblerFracture ( MeshLib::Element const & e,
std::size_t const n_variables,
std::size_t const local_matrix_size,
std::vector< unsigned > const & dofIndex_to_localIndex,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
SmallDeformationProcessData< DisplacementDim > & process_data )

Definition at line 31 of file SmallDeformationLocalAssemblerFracture-impl.h.

41 n_variables * ShapeFunction::NPOINTS * DisplacementDim,
42 dofIndex_to_localIndex),
43 _process_data(process_data),
44 _integration_method(integration_method),
47 DisplacementDim>(e, is_axially_symmetric,
49 _element(e)
50{
51 assert(_element.getDimension() == DisplacementDim - 1);
52
53 unsigned const n_integration_points =
55
56 _ip_data.reserve(n_integration_points);
57 _secondary_data.N.resize(n_integration_points);
58
59 auto mat_id = (*_process_data.mesh_prop_materialIDs)[e.getID()];
60 auto frac_id = _process_data.map_materialID_to_fractureID[mat_id];
61 _fracture_property = &_process_data.fracture_properties[frac_id];
62 for (auto fid : process_data.vec_ele_connected_fractureIDs[e.getID()])
63 {
64 _fracID_to_local.insert({fid, _fracture_props.size()});
65 _fracture_props.push_back(&_process_data.fracture_properties[fid]);
66 }
67
68 _junction_props = process_data.vec_ele_connected_junctionIDs[e.getID()] |
69 ranges::views::transform(
70 [&](auto const jid)
71 { return &_process_data.junction_properties[jid]; }) |
72 ranges::to<std::vector>;
73
75 x_position.setElementID(_element.getID());
76 for (unsigned ip = 0; ip < n_integration_points; ip++)
77 {
78 _ip_data.emplace_back(*_process_data.fracture_model);
79 auto const& sm = _shape_matrices[ip];
80 auto& ip_data = _ip_data[ip];
81 ip_data.integration_weight =
83 sm.integralMeasure * sm.detJ;
84 ip_data.h_matrices.setZero(DisplacementDim,
85 ShapeFunction::NPOINTS * DisplacementDim);
86
87 computeHMatrix<DisplacementDim, ShapeFunction::NPOINTS,
89 HMatrixType>(sm.N, ip_data.h_matrices);
90
91 // Initialize current time step values
92 ip_data.w.setZero(DisplacementDim);
93 ip_data.sigma.setZero(DisplacementDim);
94
95 // Previous time step values are not initialized and are set later.
96 ip_data.sigma_prev.resize(DisplacementDim);
97 ip_data.w_prev.resize(DisplacementDim);
98
99 ip_data.C.resize(DisplacementDim, DisplacementDim);
100
101 ip_data.aperture0 = _fracture_property->aperture0(0, x_position)[0];
102 ip_data.aperture_prev = ip_data.aperture0;
103
104 _secondary_data.N[ip] = sm.N;
105 }
106}
constexpr double getWeight() const
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
std::vector< IntegrationPointDataType, Eigen::aligned_allocator< IntegrationPointDataType > > _ip_data
std::vector< ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > _shape_matrices
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
void computeHMatrix(N_Type const &N, HMatrixType &H)
Fills a H-matrix based on given shape function.
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
ParameterLib::Parameter< double > const & aperture0
Initial aperture.
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_element, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_fracID_to_local, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_fracture_property, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_fracture_props, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_integration_method, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_ip_data, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_junction_props, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_process_data, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_secondary_data, ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_shape_matrices, ProcessLib::LIE::FractureProperty::aperture0, ProcessLib::computeHMatrix(), MeshLib::Element::getDimension(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), ProcessLib::LIE::SmallDeformation::SecondaryData< ShapeMatrixType >::N, ParameterLib::SpatialPosition::setElementID(), ProcessLib::LIE::SmallDeformation::SmallDeformationProcessData< DisplacementDim >::vec_ele_connected_fractureIDs, and ProcessLib::LIE::SmallDeformation::SmallDeformationProcessData< DisplacementDim >::vec_ele_connected_junctionIDs.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::assemble ( double const ,
double const ,
std::vector< double > const & ,
std::vector< double > const & ,
std::vector< double > & ,
std::vector< double > & ,
std::vector< double > &  )
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 66 of file SmallDeformationLocalAssemblerFracture.h.

72 {
74 "SmallDeformationLocalAssembler: assembly without jacobian is not "
75 "implemented.");
76 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ assembleWithJacobian()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::assembleWithJacobian ( double const t,
double const dt,
Eigen::VectorXd const & local_u,
Eigen::VectorXd & local_b,
Eigen::MatrixXd & local_J )
overridevirtual

Reimplemented from ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 109 of file SmallDeformationLocalAssemblerFracture-impl.h.

113{
114 auto const N_DOF_PER_VAR = ShapeFunction::NPOINTS * DisplacementDim;
115 auto const n_fractures = _fracture_props.size();
116 auto const n_junctions = _junction_props.size();
117 auto const n_enrich_var = n_fractures + n_junctions;
118
119 //--------------------------------------------------------------------------------------
120 // prepare sub vectors, matrices for regular displacement (u) and
121 // displacement jumps (g)
122 //
123 // example with two fractures with one intersection:
124 // b = |b(g1)|
125 // |b(g2)|
126 //
127 // J = |J(g1,g1) J(g1,g2)|
128 // |J(g2,g1) J(g2,g2)|
129 //--------------------------------------------------------------------------------------
130
131 using BlockVectorType =
132 typename Eigen::VectorXd::FixedSegmentReturnType<N_DOF_PER_VAR>::Type;
133 using BlockMatrixType =
134 Eigen::Block<Eigen::MatrixXd, N_DOF_PER_VAR, N_DOF_PER_VAR>;
135
136 std::vector<BlockVectorType> vec_local_b_g;
137 for (unsigned i = 0; i < n_enrich_var; i++)
138 {
139 vec_local_b_g.push_back(
140 local_b.segment<N_DOF_PER_VAR>(N_DOF_PER_VAR * i));
141 }
142 std::vector<std::vector<BlockMatrixType>> vec_local_J_gg(n_enrich_var);
143 for (unsigned i = 0; i < n_enrich_var; i++)
144 {
145 for (unsigned j = 0; j < n_enrich_var; j++)
146 {
147 auto sub_gg = local_J.block<N_DOF_PER_VAR, N_DOF_PER_VAR>(
148 N_DOF_PER_VAR * i, N_DOF_PER_VAR * j);
149 vec_local_J_gg[i].push_back(sub_gg);
150 }
151 }
152
153 std::vector<Eigen::VectorXd> vec_nodal_g;
154 for (unsigned i = 0; i < n_enrich_var; i++)
155 {
156 auto sub = const_cast<Eigen::VectorXd&>(local_u).segment<N_DOF_PER_VAR>(
157 N_DOF_PER_VAR * i);
158 vec_nodal_g.push_back(sub);
159 }
160
161 //------------------------------------------------
162 // integration
163 //------------------------------------------------
164 // the index of a normal (normal to a fracture plane) component
165 // in a displacement vector
166 int const index_normal = DisplacementDim - 1;
167 auto const& R = _fracture_property->R;
168
169 unsigned const n_integration_points =
171
173 x_position.setElementID(_element.getID());
174
175 for (unsigned ip = 0; ip < n_integration_points; ip++)
176 {
177 auto& ip_data = _ip_data[ip];
178 auto const& integration_weight = ip_data.integration_weight;
179 auto const& H = ip_data.h_matrices;
180 auto& mat = ip_data.fracture_material;
181 auto& sigma = ip_data.sigma;
182 auto const& sigma_prev = ip_data.sigma_prev;
183 auto& w = ip_data.w;
184 auto const& w_prev = ip_data.w_prev;
185 auto& C = ip_data.C;
186 auto& state = *ip_data.material_state_variables;
187 auto const& N = _secondary_data.N[ip];
188
189 auto const ip_physical_coords(computePhysicalCoordinates(_element, N));
190 std::vector<double> const levelsets(duGlobalEnrichments(
192 _fracID_to_local, ip_physical_coords));
193
194 // du = du^hat + sum_i(enrich^br_i(x) * [u]_i) + sum_i(enrich^junc_i(x)
195 // * [u]_i)
196 Eigen::VectorXd nodal_gap(N_DOF_PER_VAR);
197 nodal_gap.setZero();
198 for (unsigned i = 0; i < n_enrich_var; i++)
199 {
200 nodal_gap += levelsets[i] * vec_nodal_g[i];
201 }
202
203 // displacement jumps
204 w.noalias() = R * H * nodal_gap;
205
206 // total aperture
207 ip_data.aperture = ip_data.aperture0 + w[index_normal];
208
209 // local C, local stress
210 mat.computeConstitutiveRelation(
211 t, x_position, ip_data.aperture0,
212 Eigen::Matrix<double, DisplacementDim, 1>::Zero(), // TODO (naumov)
213 // Replace with
214 // initial
215 // stress values
216 w_prev, w, sigma_prev, sigma, C, state);
217
218 // r_[u] += H^T*Stress
219 for (unsigned i = 0; i < n_enrich_var; i++)
220 {
221 vec_local_b_g[i].noalias() -= levelsets[i] * H.transpose() *
222 R.transpose() * sigma *
223 integration_weight;
224 }
225
226 // J_[u][u] += H^T*C*H
227 for (unsigned i = 0; i < n_enrich_var; i++)
228 {
229 for (unsigned j = 0; j < n_enrich_var; j++)
230 {
231 // J_[u][u] += (levelset * B)^T * C * (levelset * B)
232 vec_local_J_gg[i][j].noalias() +=
233 (levelsets[i] * H.transpose() * R.transpose()) * C *
234 (levelsets[j] * R * H) * integration_weight;
235 }
236 }
237 }
238}
std::vector< double > duGlobalEnrichments(std::size_t this_frac_id, std::vector< FractureProperty * > const &frac_props, std::vector< JunctionProperty * > const &junction_props, std::unordered_map< int, int > const &fracID_to_local, Eigen::Vector3d const &x)
Eigen::Vector3d computePhysicalCoordinates(MeshLib::Element const &e, Eigen::MatrixBase< Derived > const &shape)
Definition Utils.h:24
Eigen::MatrixXd R
Rotation matrix from global to local coordinates.

References ProcessLib::LIE::computePhysicalCoordinates(), ProcessLib::LIE::duGlobalEnrichments(), and ParameterLib::SpatialPosition::setElementID().

◆ computeSecondaryVariableConcreteWithVector()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::computeSecondaryVariableConcreteWithVector ( const double t,
Eigen::VectorXd const & local_u )
overridevirtual

Implements ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 241 of file SmallDeformationLocalAssemblerFracture-impl.h.

244{
245 auto const N_DOF_PER_VAR = ShapeFunction::NPOINTS * DisplacementDim;
246 auto const n_fractures = _fracture_props.size();
247 auto const n_junctions = _junction_props.size();
248 auto const n_enrich_var = n_fractures + n_junctions;
249
250 auto const& R = _fracture_property->R;
251
252 // the index of a normal (normal to a fracture plane) component
253 // in a displacement vector
254 int const index_normal = DisplacementDim - 1;
255
256 unsigned const n_integration_points =
258
260 auto const e_id = _element.getID();
261 x_position.setElementID(e_id);
262
263 std::vector<Eigen::VectorXd> vec_nodal_g;
264 for (unsigned i = 0; i < n_enrich_var; i++)
265 {
266 auto sub = const_cast<Eigen::VectorXd&>(local_u).segment<N_DOF_PER_VAR>(
267 N_DOF_PER_VAR * i);
268 vec_nodal_g.push_back(sub);
269 }
270
271 for (unsigned ip = 0; ip < n_integration_points; ip++)
272 {
273 auto& ip_data = _ip_data[ip];
274 auto const& H = ip_data.h_matrices;
275 auto& mat = ip_data.fracture_material;
276 auto& sigma = ip_data.sigma;
277 auto const& sigma_prev = ip_data.sigma_prev;
278 auto& w = ip_data.w;
279 auto const& w_prev = ip_data.w_prev;
280 auto& C = ip_data.C;
281 auto& state = *ip_data.material_state_variables;
282 auto& b_m = ip_data.aperture;
283 auto const& N = _secondary_data.N[ip];
284
285 auto const ip_physical_coords(computePhysicalCoordinates(_element, N));
286 std::vector<double> const levelsets(duGlobalEnrichments(
288 _fracID_to_local, ip_physical_coords));
289
290 // du = du^hat + sum_i(enrich^br_i(x) * [u]_i) + sum_i(enrich^junc_i(x)
291 // * [u]_i)
292 Eigen::VectorXd nodal_gap(N_DOF_PER_VAR);
293 nodal_gap.setZero();
294 for (unsigned i = 0; i < n_enrich_var; i++)
295 {
296 nodal_gap += levelsets[i] * vec_nodal_g[i];
297 }
298
299 // displacement jumps in local coordinates
300 w.noalias() = R * H * nodal_gap;
301
302 // aperture
303 b_m = ip_data.aperture0 + w[index_normal];
304 if (b_m < 0.0)
305 {
306 OGS_FATAL(
307 "Element {:d}, gp {:d}: Fracture aperture is {:g}, but it must "
308 "be non-negative.",
309 _element.getID(), ip, b_m);
310 }
311
312 // local C, local stress
313 mat.computeConstitutiveRelation(
314 t, x_position, ip_data.aperture0,
315 Eigen::Matrix<double, DisplacementDim, 1>::Zero(), // TODO (naumov)
316 // Replace with
317 // initial
318 // stress values
319 w_prev, w, sigma_prev, sigma, C, state);
320 }
321
322 double ele_b = 0;
323 ForceVectorType ele_sigma = ForceVectorType::Zero(DisplacementDim);
324 ForceVectorType ele_w = ForceVectorType::Zero(DisplacementDim);
325
326 for (unsigned ip = 0; ip < n_integration_points; ip++)
327 {
328 auto const& ip_data = _ip_data[ip];
329 ele_b += ip_data.aperture;
330 ele_w += ip_data.w;
331 ele_sigma += ip_data.sigma;
332 }
333 ele_b /= n_integration_points;
334 ele_w /= n_integration_points;
335 ele_sigma /= n_integration_points;
336 (*_process_data.mesh_prop_b)[_element.getID()] = ele_b;
337
338 Eigen::Map<GlobalDimVectorType>(
339 &(*_process_data.element_fracture_stresses)[e_id * DisplacementDim]) =
340 ele_sigma;
341
342 Eigen::Map<GlobalDimVectorType>(
343 &(*_process_data.element_local_jumps)[e_id * DisplacementDim]) = ele_w;
344}

References ProcessLib::LIE::computePhysicalCoordinates(), ProcessLib::LIE::duGlobalEnrichments(), OGS_FATAL, and ParameterLib::SpatialPosition::setElementID().

◆ getIntPtEpsilon()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::getIntPtEpsilon ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 118 of file SmallDeformationLocalAssemblerFracture.h.

123 {
124 cache.resize(0);
125 return cache;
126 }

◆ getIntPtFractureAperture()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::getIntPtFractureAperture ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 371 of file SmallDeformationLocalAssemblerFracture-impl.h.

377{
380}
std::vector< double > const & getIntegrationPointScalarData(IntegrationPointDataVector const &ip_data_vector, MemberType IpData::*const member, std::vector< double > &cache)

References ProcessLib::getIntegrationPointScalarData().

◆ getIntPtFractureStress()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::getIntPtFractureStress ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
overridevirtual

Implements ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 348 of file SmallDeformationLocalAssemblerFracture-impl.h.

354{
355 unsigned const n_integration_points = _ip_data.size();
356 cache.clear();
357 auto cache_matrix = MathLib::createZeroedMatrix<Eigen::Matrix<
358 double, DisplacementDim, Eigen::Dynamic, Eigen::RowMajor>>(
359 cache, DisplacementDim, n_integration_points);
360
361 for (unsigned ip = 0; ip < n_integration_points; ip++)
362 {
363 cache_matrix.col(ip).noalias() = _ip_data[ip].sigma;
364 }
365
366 return cache;
367}
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References MathLib::createZeroedMatrix().

◆ getIntPtSigma()

template<typename ShapeFunction , int DisplacementDim>
std::vector< double > const & ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::getIntPtSigma ( const double ,
std::vector< GlobalVector * > const & ,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & ,
std::vector< double > & cache ) const
inlineoverridevirtual

Implements ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerInterface.

Definition at line 108 of file SmallDeformationLocalAssemblerFracture.h.

113 {
114 cache.resize(0);
115 return cache;
116 }

◆ getShapeMatrix()

template<typename ShapeFunction , int DisplacementDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 99 of file SmallDeformationLocalAssemblerFracture.h.

101 {
102 auto const& N = _secondary_data.N[integration_point];
103
104 // assumes N is stored contiguously in memory
105 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
106 }

References ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_secondary_data, and ProcessLib::LIE::SmallDeformation::SecondaryData< ShapeMatrixType >::N.

◆ preTimestepConcrete()

template<typename ShapeFunction , int DisplacementDim>
void ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::preTimestepConcrete ( std::vector< double > const & ,
double const ,
double const  )
inlineoverridevirtual

Member Data Documentation

◆ _element

◆ _fracID_to_local

template<typename ShapeFunction , int DisplacementDim>
std::unordered_map<int, int> ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_fracID_to_local
private

◆ _fracture_property

◆ _fracture_props

◆ _integration_method

◆ _ip_data

◆ _junction_props

◆ _process_data

◆ _secondary_data

◆ _shape_matrices

template<typename ShapeFunction , int DisplacementDim>
std::vector<ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices> > ProcessLib::LIE::SmallDeformation::SmallDeformationLocalAssemblerFracture< ShapeFunction, DisplacementDim >::_shape_matrices
private

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