template<typename ShapeFunction, typename IntegrationMethod>
class ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction, IntegrationMethod >
Definition at line 29 of file HeatTransportBHELocalAssemblerSoil.h.
|
| HeatTransportBHELocalAssemblerSoil (HeatTransportBHELocalAssemblerSoil const &)=delete |
|
| HeatTransportBHELocalAssemblerSoil (HeatTransportBHELocalAssemblerSoil &&)=delete |
|
| HeatTransportBHELocalAssemblerSoil (MeshLib::Element const &e, bool is_axially_symmetric, unsigned const integration_order, HeatTransportBHEProcessData &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 |
|
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
| Provides the shape matrix at the given integration point. More...
|
|
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 >
void ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction, IntegrationMethod >::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 > & |
|
|
) |
| |
|
overridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 71 of file HeatTransportBHELocalAssemblerSoil-impl.h.
79 assert(local_x.size() == ShapeFunction::NPOINTS);
82 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
83 local_M_data, ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
84 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
85 local_K_data, ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
91 auto const& solid_phase = medium.phase(
"Solid");
92 auto const& liquid_phase = medium.phase(
"AqueousLiquid");
96 unsigned const n_integration_points =
99 for (
unsigned ip = 0; ip < n_integration_points; ip++)
103 auto const&
N = ip_data.N;
104 auto const& dNdx = ip_data.dNdx;
105 auto const& w = ip_data.integration_weight;
107 double T_int_pt = 0.0;
114 auto const density_s =
116 .template value<double>(vars, pos, t, dt);
118 auto const heat_capacity_s =
122 .template value<double>(vars, pos, t, dt);
124 auto const density_f =
126 .template value<double>(vars, pos, t, dt);
128 auto const heat_capacity_f =
132 .template value<double>(vars, pos, t, dt);
136 .template value<double>(vars, pos, t, dt);
138 auto const velocity =
141 .template value<Eigen::Vector3d>(vars, pos, t, dt);
149 .value(vars, pos, t, dt));
153 double const velocity_magnitude = velocity.norm();
155 if (velocity_magnitude >= std::numeric_limits<double>::epsilon())
157 auto const thermal_dispersivity_longitudinal =
161 .template value<double>();
162 auto const thermal_dispersivity_transversal =
166 .template value<double>();
168 auto const thermal_dispersivity =
169 density_f * heat_capacity_f *
170 (thermal_dispersivity_transversal * velocity_magnitude *
171 Eigen::Matrix3d::Identity() +
172 (thermal_dispersivity_longitudinal -
173 thermal_dispersivity_transversal) /
174 velocity_magnitude * velocity * velocity.transpose());
175 thermal_conductivity_dispersivity += thermal_dispersivity;
180 (dNdx.transpose() * thermal_conductivity_dispersivity * dNdx +
181 N.transpose() * velocity.transpose() * dNdx * density_f *
186 local_M.noalias() +=
N.transpose() *
N * w *
187 (density_s * heat_capacity_s * (1 -
porosity) +
188 density_f * heat_capacity_f * porosity);
void setIntegrationPoint(unsigned integration_point)
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
@ thermal_transversal_dispersivity
@ thermal_longitudinal_dispersivity
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
std::unique_ptr< MaterialPropertyLib::MaterialSpatialDistributionMap > media_map
References MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor< 3 >(), MaterialPropertyLib::phase_velocity, MaterialPropertyLib::porosity, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::temperature, MaterialPropertyLib::thermal_conductivity, MaterialPropertyLib::thermal_longitudinal_dispersivity, and MaterialPropertyLib::thermal_transversal_dispersivity.