template<typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
class ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >
Definition at line 79 of file RichardsFlowFEM.h.
|
| LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, RichardsFlowProcessData const &process_data) |
|
void | assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override |
|
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
| Provides the shape matrix at the given integration point. More...
|
|
std::vector< double > const & | getIntPtSaturation (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const 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 |
|
virtual | ~LocalAssemblerInterface ()=default |
|
void | setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, 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_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) |
|
virtual void | assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) |
|
virtual void | assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, 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_dot, 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, double const t, double const dt) |
|
void | postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id) |
|
virtual std::vector< double > | interpolateNodalValuesToIntegrationPoints (std::vector< double > const &) |
|
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. More...
|
|
template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assemble |
( |
double const |
t, |
|
|
double const |
dt, |
|
|
std::vector< double > const & |
local_x, |
|
|
std::vector< double > const & |
, |
|
|
std::vector< double > & |
local_M_data, |
|
|
std::vector< double > & |
local_K_data, |
|
|
std::vector< double > & |
local_b_data |
|
) |
| |
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 134 of file RichardsFlowFEM.h.
141 auto const local_matrix_size = local_x.size();
144 assert(local_matrix_size == ShapeFunction::NPOINTS *
NUM_NODAL_DOF);
146 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
147 local_M_data, local_matrix_size, local_matrix_size);
148 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
149 local_K_data, local_matrix_size, local_matrix_size);
150 auto local_b = MathLib::createZeroedVector<NodalVectorType>(
151 local_b_data, local_matrix_size);
153 unsigned const n_integration_points =
160 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
166 .template value<double>(vars, pos, t, dt);
168 for (
unsigned ip = 0; ip < n_integration_points; ip++)
171 double p_int_pt = 0.0;
174 vars[
static_cast<int>(
176 vars[
static_cast<int>(
180 MaterialPropertyLib::formEigenTensor<GlobalDim>(
182 .value(vars, pos, t, dt));
186 .template value<double>(vars, pos, t, dt);
191 .template value<double>(vars, pos, t, dt);
193 vars[
static_cast<int>(
196 double const dSw_dpc =
199 .template dValue<double>(
203 auto const drhow_dp =
206 .template dValue<double>(
211 .template value<double>(vars, pos, t, dt);
212 double const mass_mat_coeff =
215 local_M.noalias() += mass_mat_coeff *
_ip_data[ip].mass_operator;
221 .template value<double>(vars, pos, t, dt);
225 .template value<double>(vars, pos, t, dt);
228 _ip_data[ip].integration_weight * (k_rel / mu);
235 .template value<double>(vars, pos, t, dt);
237 assert(body_force.size() == GlobalDim);
241 local_b.noalias() += (k_rel / mu) * rho_w * gravity_operator;
246 for (
int idx_ml = 0; idx_ml < local_M.cols(); idx_ml++)
248 double const mass_lump_val = local_M.col(idx_ml).sum();
249 local_M.col(idx_ml).setZero();
250 local_M(idx_ml, idx_ml) = mass_lump_val;
virtual std::size_t getID() const final
Returns the ID of the element.
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
typename ShapeMatricesType::NodalVectorType NodalVectorType
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
std::unique_ptr< MaterialPropertyLib::MaterialSpatialDistributionMap > media_map
bool const has_mass_lumping
Eigen::VectorXd const specific_body_force
References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_saturation, MaterialPropertyLib::capillary_pressure, MaterialPropertyLib::density, MeshLib::Element::getID(), ProcessLib::RichardsFlow::RichardsFlowProcessData::has_gravity, ProcessLib::RichardsFlow::RichardsFlowProcessData::has_mass_lumping, MaterialPropertyLib::liquid_saturation, ProcessLib::RichardsFlow::RichardsFlowProcessData::media_map, ProcessLib::RichardsFlow::NUM_NODAL_DOF, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), ProcessLib::RichardsFlow::RichardsFlowProcessData::specific_body_force, MaterialPropertyLib::storage, MaterialPropertyLib::temperature, and MaterialPropertyLib::viscosity.
template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Implements ProcessLib::RichardsFlow::RichardsFlowLocalAssemblerInterface.
Definition at line 274 of file RichardsFlowFEM.h.
282 double const dt = std::numeric_limits<double>::quiet_NaN();
284 constexpr
int process_id = 0;
287 assert(!indices.empty());
288 auto const local_x = x[process_id]->get(indices);
295 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
302 .template value<double>(vars, pos, t, dt);
304 auto const p_nodal_values = Eigen::Map<const NodalVectorType>(
305 &local_x[0], ShapeFunction::NPOINTS);
307 unsigned const n_integration_points =
312 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
313 cache, GlobalDim, n_integration_points);
315 for (
unsigned ip = 0; ip < n_integration_points; ++ip)
317 double p_int_pt = 0.0;
319 vars[
static_cast<int>(
321 vars[
static_cast<int>(
327 .template value<double>(vars, pos, t, dt);
328 vars[
static_cast<int>(
332 MaterialPropertyLib::formEigenTensor<GlobalDim>(
334 .value(vars, pos, t, dt));
340 .template value<double>(vars, pos, t, dt);
343 .template value<double>(vars, pos, t, dt);
345 cache_vec.col(ip).noalias() =
346 -K_mat_coeff *
_ip_data[ip].dNdx * p_nodal_values;
352 .template value<double>(vars, pos, t, dt);
354 assert(body_force.size() == GlobalDim);
357 cache_vec.col(ip).noalias() += K_mat_coeff * rho_w * body_force;
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
References ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::RichardsFlow::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, MaterialPropertyLib::capillary_pressure, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MeshLib::Element::getID(), NumLib::getIndices(), ProcessLib::RichardsFlow::RichardsFlowProcessData::has_gravity, MaterialPropertyLib::liquid_saturation, ProcessLib::RichardsFlow::RichardsFlowProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, MaterialPropertyLib::relative_permeability, MaterialPropertyLib::saturation, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), ProcessLib::RichardsFlow::RichardsFlowProcessData::specific_body_force, MaterialPropertyLib::temperature, and MaterialPropertyLib::viscosity.