OGS
SteadyStateDiffusionFEM.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <vector>
14
30
31namespace ProcessLib
32{
33namespace SteadyStateDiffusion
34{
35const unsigned NUM_NODAL_DOF = 1;
36
40{
41public:
42 virtual std::vector<double> const& getIntPtDarcyVelocity(
43 const double t,
44 std::vector<GlobalVector*> const& x,
45 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
46 std::vector<double>& cache) const = 0;
47};
48
49template <typename ShapeFunction, int GlobalDim>
51{
54
56 ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim>;
57
61
62public:
66 MeshLib::Element const& element,
67 std::size_t const /*local_matrix_size*/,
68 NumLib::GenericIntegrationMethod const& integration_method,
69 bool const is_axially_symmetric,
70 SteadyStateDiffusionData const& process_data)
71 : _element(element),
72 _process_data(process_data),
73 _integration_method(integration_method),
75 NumLib::initShapeMatrices<ShapeFunction, ShapeMatricesType,
76 GlobalDim>(
77 element, is_axially_symmetric, _integration_method))
78 {
79 }
80
81 void assemble(double const t, double const dt,
82 std::vector<double> const& local_x,
83 std::vector<double> const& /*local_x_prev*/,
84 std::vector<double>& /*local_M_data*/,
85 std::vector<double>& local_K_data,
86 std::vector<double>& /*local_b_data*/) override
87 {
88 auto const local_matrix_size = local_x.size();
89 // This assertion is valid only if all nodal d.o.f. use the same shape
90 // matrices.
91 assert(local_matrix_size == ShapeFunction::NPOINTS * NUM_NODAL_DOF);
92
93 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
94 local_K_data, local_matrix_size, local_matrix_size);
95
96 unsigned const n_integration_points =
98
101
102 auto const& medium =
105 vars.temperature =
106 medium
107 .property(
109 .template value<double>(vars, pos, t, dt);
110
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 pos.setIntegrationPoint(ip);
114 auto const& sm = _shape_matrices[ip];
115 auto const& wp = _integration_method.getWeightedPoint(ip);
116
117 double p_int_pt = 0.0;
118 NumLib::shapeFunctionInterpolate(local_x, sm.N, p_int_pt);
119 vars.liquid_phase_pressure = p_int_pt;
120 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
122 .value(vars, pos, t, dt));
123
124 local_K.noalias() += sm.dNdx.transpose() * k * sm.dNdx * sm.detJ *
125 sm.integralMeasure * wp.getWeight();
126 }
127 }
128
131 Eigen::Vector3d getFlux(MathLib::Point3d const& p_local_coords,
132 double const t,
133 std::vector<double> const& local_x) const override
134 {
135 // TODO (tf) Temporary value not used by current material models. Need
136 // extension of getFlux interface
137 double const dt = std::numeric_limits<double>::quiet_NaN();
138
139 // Eval shape matrices at given point
140 // Note: Axial symmetry is set to false here, because we only need dNdx
141 // here, which is not affected by axial symmetry.
142 auto const shape_matrices =
144 GlobalDim>(
145 _element, false /*is_axially_symmetric*/,
146 std::array{p_local_coords})[0];
147
148 // fetch hydraulic conductivity
151 auto const& medium =
153
155 vars.temperature =
156 medium
157 .property(
159 .template value<double>(vars, pos, t, dt);
160 double pressure = 0.0;
161 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, pressure);
162 vars.liquid_phase_pressure = pressure;
163
164 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
166 .value(vars, pos, t, dt));
167
168 Eigen::Vector3d flux(0.0, 0.0, 0.0);
169 flux.head<GlobalDim>() =
170 -k * shape_matrices.dNdx *
171 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
172
173 return flux;
174 }
175
176 Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
177 const unsigned integration_point) const override
178 {
179 auto const& N = _shape_matrices[integration_point].N;
180
181 // assumes N is stored contiguously in memory
182 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
183 }
184
185 std::vector<double> const& getIntPtDarcyVelocity(
186 const double t,
187 std::vector<GlobalVector*> const& x,
188 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
189 std::vector<double>& cache) const override
190 {
191 // TODO (tf) Temporary value not used by current material models. Need
192 // extension of secondary variable interface.
193 double const dt = std::numeric_limits<double>::quiet_NaN();
194
195 auto const n_integration_points =
197
198 int const process_id = 0; // monolithic scheme
199 auto const indices =
200 NumLib::getIndices(_element.getID(), *dof_table[process_id]);
201 assert(!indices.empty());
202 auto const local_x = x[process_id]->get(indices);
203 auto const local_x_vec =
204 MathLib::toVector<Eigen::Matrix<double, ShapeFunction::NPOINTS, 1>>(
205 local_x, ShapeFunction::NPOINTS);
206
207 cache.clear();
208 auto cache_mat = MathLib::createZeroedMatrix<
209 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
210 cache, GlobalDim, n_integration_points);
211
214
215 auto const& medium =
217
219 vars.temperature =
220 medium
221 .property(
223 .template value<double>(vars, pos, t, dt);
224 double pressure = 0.0;
225 for (unsigned i = 0; i < n_integration_points; ++i)
226 {
227 pos.setIntegrationPoint(i);
229 pressure);
230 vars.liquid_phase_pressure = pressure;
231
232 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
234 .value(vars, pos, t, dt));
235 // dimensions: (d x 1) = (d x n) * (n x 1)
236 cache_mat.col(i).noalias() =
237 -k * _shape_matrices[i].dNdx * local_x_vec;
238 }
239
240 return cache;
241 }
242
243private:
246
248 std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>>
250};
251
252} // namespace SteadyStateDiffusion
253} // namespace ProcessLib
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)
void setIntegrationPoint(unsigned integration_point)
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &, std::vector< double > &local_K_data, std::vector< double > &) override
LocalAssemblerData(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SteadyStateDiffusionData const &process_data)
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
Eigen::Vector3d getFlux(MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const override
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
typename LocalAssemblerTraits::LocalVector NodalVectorType
NumLib::GenericIntegrationMethod const & _integration_method
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
VectorType< GlobalDim > GlobalDimVectorType
MaterialPropertyLib::MaterialSpatialDistributionMap media_map
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix