33template <
typename Mat>
36 for (
unsigned r = 0; r < mat.rows(); ++r)
56 for (
unsigned c = 0; c < mat.cols(); ++c)
61 std::printf(
" %.16e", mat(r, c));
68 std::printf(
" %23.16g", mat(r, c));
86template <
typename Vec>
89 for (
unsigned r = 0; r < vec.size(); ++r)
98 std::printf(
"| %.16e | ", vec[r]);
105 std::printf(
"[ %23.16g ]", vec[r]);
118template <
typename ShapeFunction_,
int GlobalDim>
123 bool const is_axially_symmetric,
126 _integration_method(integration_method),
129 e, is_axially_symmetric, _integration_method)),
130 _d(asm_params, _integration_method.getNumberOfPoints(), GlobalDim)
134template <
typename ShapeFunction_,
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();
144 assert(local_matrix_size == ShapeFunction::NPOINTS *
NODAL_DOF);
147 local_M_data, local_matrix_size, local_matrix_size);
149 local_K_data, local_matrix_size, local_matrix_size);
151 local_b_data, local_matrix_size);
153 unsigned const n_integration_points =
154 _integration_method.getNumberOfPoints();
156 _d.preEachAssemble();
158 for (
unsigned ip = 0; ip < n_integration_points; ip++)
160 auto const& sm = _shape_matrices[ip];
161 auto const& wp = _integration_method.getWeightedPoint(ip);
162 auto const weight = wp.getWeight();
164 _d.assembleIntegrationPoint(ip, local_x, sm, weight, local_M, local_K,
168 if (_d.getAssemblyParameters().output_element_matrices)
170 std::puts(
"### Element: ?");
172 std::puts(
"---Velocity of water");
173 for (
auto const& vs : _d.getData().velocity)
178 std::printf(
"%23.16e ", v);
183 std::printf(
"\n---Mass matrix: \n");
187 std::printf(
"---Laplacian + Advective + Content matrix: \n");
191 std::printf(
"---RHS: \n");
197template <
typename ShapeFunction_,
int GlobalDim>
198std::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;
208template <
typename ShapeFunction_,
int GlobalDim>
209std::vector<double>
const&
212 std::vector<GlobalVector*>
const& ,
213 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
214 std::vector<double>& cache)
const
216 auto const rho_SR = _d.getData().solid_density;
217 auto const rho_SR_dry = _d.getAssemblyParameters().rho_SR_dry;
220 cache.reserve(rho_SR.size());
222 transform(cbegin(rho_SR), cend(rho_SR), back_inserter(cache),
223 [&](
auto const rho) {
return rho / rho_SR_dry - 1.0; });
228template <
typename ShapeFunction_,
int GlobalDim>
229std::vector<double>
const&
231 const double , std::vector<GlobalVector*>
const& ,
232 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
233 std::vector<double>& cache)
const
235 auto const fac = _d.getData().reaction_adaptor->getReactionDampingFactor();
236 auto const num_integration_points = _d.getData().solid_density.size();
239 cache.resize(num_integration_points, fac);
244template <
typename ShapeFunction_,
int GlobalDim>
245std::vector<double>
const&
248 std::vector<GlobalVector*>
const& ,
249 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& ,
250 std::vector<double>& )
const
252 return _d.getData().reaction_rate;
255template <
typename ShapeFunction_,
int GlobalDim>
256std::vector<double>
const&
259 std::vector<GlobalVector*>
const& x,
260 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
261 std::vector<double>& cache)
const
263 auto const n_integration_points = _integration_method.getNumberOfPoints();
265 constexpr int process_id = 0;
268 assert(!indices.empty());
269 auto const local_x = x[process_id]->get(indices);
272 Eigen::Matrix<double, NODAL_DOF, Eigen::Dynamic, Eigen::RowMajor>>(
273 local_x,
NODAL_DOF, ShapeFunction_::NPOINTS);
277 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
278 cache, GlobalDim, n_integration_points);
280 for (
unsigned i = 0; i < n_integration_points; ++i)
287 auto const& k = _d.getAssemblyParameters().solid_perm_tensor;
288 cache_mat.col(i).noalias() =
290 (_shape_matrices[i].dNdx *
298template <
typename ShapeFunction_,
int GlobalDim>
300 std::vector<double>
const& local_x,
301 std::vector<double>
const& local_x_prev_ts)
303 return _d.getReactionAdaptor().checkBounds(local_x, local_x_prev_ts);
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
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
ShapeFunction_ ShapeFunction
TESLocalAssembler(MeshLib::Element const &e, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, AssemblyParams const &asm_params)
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
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 & getIntPtReactionRate(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
bool checkBounds(std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts) override
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
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 &, const ShapeMatrix &)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
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