33 template <
typename Mat>
36 for (
unsigned r = 0;
r < mat.rows(); ++
r)
40 case MatOutType::OGS5:
47 case MatOutType::PYTHON:
56 for (
unsigned c = 0;
c < mat.cols(); ++
c)
60 case MatOutType::OGS5:
61 std::printf(
" %.16e", mat(
r,
c));
63 case MatOutType::PYTHON:
68 std::printf(
" %23.16g", mat(
r,
c));
75 case MatOutType::OGS5:
78 case MatOutType::PYTHON:
86 template <
typename Vec>
89 for (
unsigned r = 0;
r < vec.size(); ++
r)
93 case MatOutType::OGS5:
98 std::printf(
"| %.16e | ", vec[
r]);
100 case MatOutType::PYTHON:
105 std::printf(
"[ %23.16g ]", vec[
r]);
118 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
122 bool is_axially_symmetric,
123 unsigned const integration_order,
126 _integration_method(integration_order),
129 e, is_axially_symmetric, _integration_method)),
130 _d(asm_params, _integration_method.getNumberOfPoints(), GlobalDim)
134 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
136 double const ,
double const , std::vector<double>
const& local_x,
137 std::vector<double>
const& ,
138 std::vector<double>& local_M_data, std::vector<double>& local_K_data,
139 std::vector<double>& local_b_data)
141 auto const local_matrix_size = local_x.size();
143 assert(local_matrix_size == ShapeFunction::NPOINTS *
NODAL_DOF);
145 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
146 local_M_data, local_matrix_size, local_matrix_size);
147 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
148 local_K_data, local_matrix_size, local_matrix_size);
149 auto local_b = MathLib::createZeroedVector<NodalVectorType>(local_b_data,
152 unsigned const n_integration_points =
153 _integration_method.getNumberOfPoints();
155 _d.preEachAssemble();
157 for (
unsigned ip = 0; ip < n_integration_points; ip++)
159 auto const& sm = _shape_matrices[ip];
160 auto const& wp = _integration_method.getWeightedPoint(ip);
161 auto const weight = wp.getWeight();
163 _d.assembleIntegrationPoint(ip, local_x, sm, weight, local_M, local_K,
167 if (_d.getAssemblyParameters().output_element_matrices)
169 std::puts(
"### Element: ?");
171 std::puts(
"---Velocity of water");
172 for (
auto const& vs : _d.getData().velocity)
177 std::printf(
"%23.16e ", v);
182 std::printf(
"\n---Mass matrix: \n");
186 std::printf(
"---Laplacian + Advective + Content matrix: \n");
190 std::printf(
"---RHS: \n");
196 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
197 std::vector<double>
const&
201 std::vector<GlobalVector*>
const& ,
202 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
203 std::vector<double>& )
const
205 return _d.getData().solid_density;
208 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
209 std::vector<double>
const&
213 std::vector<GlobalVector*>
const& ,
214 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
215 std::vector<double>& cache)
const
217 auto const rho_SR = _d.getData().solid_density;
218 auto const rho_SR_dry = _d.getAssemblyParameters().rho_SR_dry;
221 cache.reserve(rho_SR.size());
223 transform(cbegin(rho_SR), cend(rho_SR), back_inserter(cache),
224 [&](
auto const rho) {
return rho / rho_SR_dry - 1.0; });
229 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
230 std::vector<double>
const&
233 const double , std::vector<GlobalVector*>
const& ,
234 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
235 std::vector<double>& cache)
const
237 auto const fac = _d.getData().reaction_adaptor->getReactionDampingFactor();
238 auto const num_integration_points = _d.getData().solid_density.size();
241 cache.resize(num_integration_points, fac);
246 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
247 std::vector<double>
const&
251 std::vector<GlobalVector*>
const& ,
252 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
253 std::vector<double>& )
const
255 return _d.getData().reaction_rate;
258 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
259 std::vector<double>
const&
263 std::vector<GlobalVector*>
const& x,
264 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
265 std::vector<double>& cache)
const
267 auto const n_integration_points = _integration_method.getNumberOfPoints();
269 constexpr
int process_id = 0;
272 assert(!indices.empty());
273 auto const local_x = x[process_id]->get(indices);
276 Eigen::Matrix<double, NODAL_DOF, Eigen::Dynamic, Eigen::RowMajor>>(
277 local_x,
NODAL_DOF, ShapeFunction_::NPOINTS);
281 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
282 cache, GlobalDim, n_integration_points);
284 for (
unsigned i = 0; i < n_integration_points; ++i)
291 auto const& k = _d.getAssemblyParameters().solid_perm_tensor;
292 cache_mat.col(i).noalias() =
293 k * (_shape_matrices[i].dNdx *
301 template <
typename ShapeFunction_,
typename IntegrationMethod_,
int GlobalDim>
304 std::vector<double>
const& local_x_prev_ts)
306 return _d.getReactionAdaptor().checkBounds(local_x, local_x_prev_ts);
std::vector< double > const & getIntPtReactionRate(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
bool checkBounds(std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts) override
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
ShapeFunction_ ShapeFunction
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
TESLocalAssembler(MeshLib::Element const &e, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, AssemblyParams const &asm_params)
std::vector< double > const & getIntPtLoading(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 & getIntPtSolidDensity(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 & getIntPtReactionDampingFactor(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
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 > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
double fluid_viscosity(const double p, const double T, const double x)
const unsigned COMPONENT_ID_PRESSURE
void ogs5OutMat(const Mat &mat)
void ogs5OutVec(const Vec &vec)
const MatOutType MATRIX_OUTPUT_FORMAT