OGS
ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
class ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >

Definition at line 213 of file ComponentTransportFEM.h.

#include <ComponentTransportFEM.h>

Inheritance diagram for ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >:
[legend]

Public Member Functions

 LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, ComponentTransportProcessData const &process_data, std::vector< std::reference_wrapper< ProcessVariable >> const &transport_process_variables)
 
void setChemicalSystemID (std::size_t const) override
 
void initializeChemicalSystemConcrete (Eigen::VectorXd const &local_x, double const t) override
 
void setChemicalSystemConcrete (Eigen::VectorXd const &local_x, double const t, double dt) override
 
void postSpeciationCalculation (std::size_t const ele_id, double const t, double const dt) override
 
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
 
void assembleBlockMatrices (int const component_id, double const t, double const dt, Eigen::Ref< const NodalVectorType > const &C_nodal_values, Eigen::Ref< const NodalVectorType > const &p_nodal_values, Eigen::Ref< LocalBlockMatrixType > KCC, Eigen::Ref< LocalBlockMatrixType > MCC, Eigen::Ref< LocalBlockMatrixType > MCp, Eigen::Ref< LocalBlockMatrixType > MpC, Eigen::Ref< LocalBlockMatrixType > Kpp, Eigen::Ref< LocalBlockMatrixType > Mpp, Eigen::Ref< LocalSegmentVectorType > Bp)
 
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) override
 
void assembleHydraulicEquation (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
void assembleComponentTransportEquation (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &, int const transport_process_id)
 
void assembleReactionEquationConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, int const transport_process_id) 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
 
std::vector< double > const & calculateIntPtDarcyVelocity (const double t, Eigen::Ref< const NodalVectorType > const &p_nodal_values, Eigen::Ref< const NodalVectorType > const &C_nodal_values, std::vector< double > &cache) const
 
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const override
 Provides the shape matrix at the given integration point. More...
 
Eigen::Vector3d getFlux (MathLib::Point3d const &pnt_local_coords, double const t, std::vector< double > const &local_x) const override
 
std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &local_x) override
 
std::vector< double > const & getInterpolatedLocalSolution (const double, std::vector< GlobalVector * > const &int_pt_x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const override
 
void computeSecondaryVariableConcrete (double const t, double const, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override
 
void postTimestepConcrete (Eigen::VectorXd const &, double const, double const) override
 
- Public Member Functions inherited from ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface
 ComponentTransportLocalAssemblerInterface ()=default
 
void setStaggeredCoupledSolutions (std::size_t const, CoupledSolutionsForStaggeredScheme *const coupling_term)
 
void initializeChemicalSystem (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t)
 
void setChemicalSystem (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 assembleReactionEquation (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, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, int const process_id)
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
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 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 Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
 Fits to staggered scheme. More...
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default
 

Private Types

using ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, GlobalDim >
 
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
 
using LocalBlockMatrixType = typename ShapeMatricesType::template MatrixType< pressure_size, pressure_size >
 
using LocalSegmentVectorType = typename ShapeMatricesType::template VectorType< pressure_size >
 
using LocalMatrixType = Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor >
 
using LocalVectorType = Eigen::Matrix< double, Eigen::Dynamic, 1 >
 
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
 
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
 
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
 
using GlobalDimNodalMatrixType = typename ShapeMatricesType::GlobalDimNodalMatrixType
 
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
 

Private Attributes

MeshLib::Element const & _element
 
ComponentTransportProcessData const & _process_data
 
IntegrationMethod const _integration_method
 
std::vector< std::reference_wrapper< ProcessVariable > > const _transport_process_variables
 
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType > > > _ip_data
 

Static Private Attributes

static const int pressure_index = 0
 
static const int first_concentration_index = ShapeFunction::NPOINTS
 
static const int pressure_size = ShapeFunction::NPOINTS
 
static const int concentration_size
 

Additional Inherited Members

- Protected Attributes inherited from ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface
CoupledSolutionsForStaggeredScheme_coupled_solutions {nullptr}
 

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 243 of file ComponentTransportFEM.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::GlobalDimNodalMatrixType = typename ShapeMatricesType::GlobalDimNodalMatrixType
private

Definition at line 241 of file ComponentTransportFEM.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 240 of file ComponentTransportFEM.h.

◆ LocalBlockMatrixType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalBlockMatrixType = typename ShapeMatricesType::template MatrixType<pressure_size, pressure_size>
private

Definition at line 227 of file ComponentTransportFEM.h.

◆ LocalMatrixType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalMatrixType = Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>
private

Definition at line 233 of file ComponentTransportFEM.h.

◆ LocalSegmentVectorType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalSegmentVectorType = typename ShapeMatricesType::template VectorType<pressure_size>
private

Definition at line 230 of file ComponentTransportFEM.h.

◆ LocalVectorType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalVectorType = Eigen::Matrix<double, Eigen::Dynamic, 1>
private

Definition at line 235 of file ComponentTransportFEM.h.

◆ NodalRowVectorType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 238 of file ComponentTransportFEM.h.

◆ NodalVectorType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::NodalVectorType = typename ShapeMatricesType::NodalVectorType
private

Definition at line 237 of file ComponentTransportFEM.h.

◆ ShapeMatrices

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 225 of file ComponentTransportFEM.h.

◆ ShapeMatricesType

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
using ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 224 of file ComponentTransportFEM.h.

Constructor & Destructor Documentation

◆ LocalAssemblerData()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalAssemblerData ( MeshLib::Element const &  element,
std::size_t const  local_matrix_size,
bool  is_axially_symmetric,
unsigned const  integration_order,
ComponentTransportProcessData const &  process_data,
std::vector< std::reference_wrapper< ProcessVariable >> const &  transport_process_variables 
)
inline

Definition at line 246 of file ComponentTransportFEM.h.

254  : _element(element),
255  _process_data(process_data),
256  _integration_method(integration_order),
257  _transport_process_variables(transport_process_variables)
258  {
259  (void)local_matrix_size;
260 
261  unsigned const n_integration_points =
262  _integration_method.getNumberOfPoints();
263  _ip_data.reserve(n_integration_points);
264 
266  pos.setElementID(_element.getID());
267 
268  auto const shape_matrices =
270  GlobalDim>(element, is_axially_symmetric,
272  auto const& medium =
273  *_process_data.media_map->getMedium(_element.getID());
274  for (unsigned ip = 0; ip < n_integration_points; ip++)
275  {
276  _ip_data.emplace_back(
277  shape_matrices[ip].N, shape_matrices[ip].dNdx,
278  _integration_method.getWeightedPoint(ip).getWeight() *
279  shape_matrices[ip].integralMeasure *
280  shape_matrices[ip].detJ);
281 
282  pos.setIntegrationPoint(ip);
283 
284  _ip_data[ip].porosity =
286  .template initialValue<double>(
287  pos, std::numeric_limits<double>::quiet_NaN() /*t*/);
288 
289  _ip_data[ip].pushBackState();
290  }
291  }
virtual std::size_t getID() const final
Returns the ID of the element.
Definition: Element.h:82
void setElementID(std::size_t element_id)
void setIntegrationPoint(unsigned integration_point)
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType > > > _ip_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
ComponentTransportProcessData const & _process_data
std::vector< std::reference_wrapper< ProcessVariable > > const _transport_process_variables
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)
std::unique_ptr< MaterialPropertyLib::MaterialSpatialDistributionMap > media_map

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, MeshLib::Element::getID(), NumLib::initShapeMatrices(), ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, MaterialPropertyLib::porosity, ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::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 443 of file ComponentTransportFEM.h.

449  {
450  auto const local_matrix_size = local_x.size();
451  // Nodal DOFs include pressure
452  int const num_nodal_dof = 1 + _transport_process_variables.size();
453  // This assertion is valid only if all nodal d.o.f. use the same shape
454  // matrices.
455  assert(local_matrix_size == ShapeFunction::NPOINTS * num_nodal_dof);
456 
457  auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>(
458  local_M_data, local_matrix_size, local_matrix_size);
459  auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>(
460  local_K_data, local_matrix_size, local_matrix_size);
461  auto local_b = MathLib::createZeroedVector<LocalVectorType>(
462  local_b_data, local_matrix_size);
463 
464  // Get block matrices
465  auto Kpp = local_K.template block<pressure_size, pressure_size>(
467  auto Mpp = local_M.template block<pressure_size, pressure_size>(
469  auto Bp = local_b.template segment<pressure_size>(pressure_index);
470 
471  auto local_p = Eigen::Map<const NodalVectorType>(
472  &local_x[pressure_index], pressure_size);
473 
474  auto const number_of_components = num_nodal_dof - 1;
475  for (int component_id = 0; component_id < number_of_components;
476  ++component_id)
477  {
478  /* Partitioned assembler matrix
479  * | pp | pc1 | pc2 | pc3 |
480  * |-----|-----|-----|-----|
481  * | c1p | c1c1| 0 | 0 |
482  * |-----|-----|-----|-----|
483  * | c2p | 0 | c2c2| 0 |
484  * |-----|-----|-----|-----|
485  * | c3p | 0 | 0 | c3c3|
486  */
487  auto concentration_index =
488  pressure_size + component_id * concentration_size;
489 
490  auto KCC =
491  local_K.template block<concentration_size, concentration_size>(
492  concentration_index, concentration_index);
493  auto MCC =
494  local_M.template block<concentration_size, concentration_size>(
495  concentration_index, concentration_index);
496  auto MCp =
497  local_M.template block<concentration_size, pressure_size>(
498  concentration_index, pressure_index);
499  auto MpC =
500  local_M.template block<pressure_size, concentration_size>(
501  pressure_index, concentration_index);
502 
503  auto local_C = Eigen::Map<const NodalVectorType>(
504  &local_x[concentration_index], concentration_size);
505 
506  assembleBlockMatrices(component_id, t, dt, local_C, local_p, KCC,
507  MCC, MCp, MpC, Kpp, Mpp, Bp);
508  }
509  }
void assembleBlockMatrices(int const component_id, double const t, double const dt, Eigen::Ref< const NodalVectorType > const &C_nodal_values, Eigen::Ref< const NodalVectorType > const &p_nodal_values, Eigen::Ref< LocalBlockMatrixType > KCC, Eigen::Ref< LocalBlockMatrixType > MCC, Eigen::Ref< LocalBlockMatrixType > MCp, Eigen::Ref< LocalBlockMatrixType > MpC, Eigen::Ref< LocalBlockMatrixType > Kpp, Eigen::Ref< LocalBlockMatrixType > Mpp, Eigen::Ref< LocalSegmentVectorType > Bp)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_transport_process_variables, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleBlockMatrices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_index, and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_size.

◆ assembleBlockMatrices()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleBlockMatrices ( int const  component_id,
double const  t,
double const  dt,
Eigen::Ref< const NodalVectorType > const &  C_nodal_values,
Eigen::Ref< const NodalVectorType > const &  p_nodal_values,
Eigen::Ref< LocalBlockMatrixType KCC,
Eigen::Ref< LocalBlockMatrixType MCC,
Eigen::Ref< LocalBlockMatrixType MCp,
Eigen::Ref< LocalBlockMatrixType MpC,
Eigen::Ref< LocalBlockMatrixType Kpp,
Eigen::Ref< LocalBlockMatrixType Mpp,
Eigen::Ref< LocalSegmentVectorType Bp 
)
inline

Definition at line 511 of file ComponentTransportFEM.h.

522  {
523  unsigned const n_integration_points =
524  _integration_method.getNumberOfPoints();
525 
527  pos.setElementID(_element.getID());
528 
529  auto const& b = _process_data.specific_body_force;
530 
532 
533  GlobalDimMatrixType const& I(
534  GlobalDimMatrixType::Identity(GlobalDim, GlobalDim));
535 
536  // Get material properties
537  auto const& medium =
538  *_process_data.media_map->getMedium(_element.getID());
539  // Select the only valid for component transport liquid phase.
540  auto const& phase = medium.phase("AqueousLiquid");
541 
542  // Assume that the component name is the same as the process variable
543  // name. Components are shifted by one because the first one is always
544  // pressure.
545  auto const& component = phase.component(
546  _transport_process_variables[component_id].get().getName());
547 
548  for (unsigned ip(0); ip < n_integration_points; ++ip)
549  {
550  pos.setIntegrationPoint(ip);
551 
552  auto& ip_data = _ip_data[ip];
553  auto const& N = ip_data.N;
554  auto const& dNdx = ip_data.dNdx;
555  auto const& w = ip_data.integration_weight;
556  auto& porosity = ip_data.porosity;
557 
558  double C_int_pt = 0.0;
559  double p_int_pt = 0.0;
560 
561  NumLib::shapeFunctionInterpolate(C_nodal_values, N, C_int_pt);
562  NumLib::shapeFunctionInterpolate(p_nodal_values, N, p_int_pt);
563 
564  vars[static_cast<int>(
566  vars[static_cast<int>(
568 
569  // update according to a particular porosity model
571  .template value<double>(vars, pos, t, dt);
572  vars[static_cast<int>(MaterialPropertyLib::Variable::porosity)] =
573  porosity;
574 
575  auto const& retardation_factor =
577  .template value<double>(vars, pos, t, dt);
578 
579  auto const& solute_dispersivity_transverse = medium.template value<
580  double>(
582 
583  auto const& solute_dispersivity_longitudinal =
584  medium.template value<double>(
587 
588  // Use the fluid density model to compute the density
589  // TODO (renchao): concentration of which component as the argument
590  // for calculation of fluid density
591  auto const density =
593  .template value<double>(vars, pos, t, dt);
594 
595  auto const decay_rate =
597  .template value<double>(vars, pos, t, dt);
598 
599  auto const& pore_diffusion_coefficient =
600  MaterialPropertyLib::formEigenTensor<GlobalDim>(
602  .value(vars, pos, t, dt));
603 
604  auto const& K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
606  vars, pos, t, dt));
607 
608  // Use the viscosity model to compute the viscosity
609  auto const mu = phase[MaterialPropertyLib::PropertyType::viscosity]
610  .template value<double>(vars, pos, t, dt);
611 
612  GlobalDimMatrixType const K_over_mu = K / mu;
613  GlobalDimVectorType const velocity =
615  ? GlobalDimVectorType(-K_over_mu *
616  (dNdx * p_nodal_values - density * b))
617  : GlobalDimVectorType(-K_over_mu * dNdx * p_nodal_values);
618 
619  const double drho_dp =
621  .template dValue<double>(
623  pos, t, dt);
624 
625  const double drho_dC =
627  .template dValue<double>(
629  t, dt);
630 
631  double const velocity_magnitude = velocity.norm();
632  GlobalDimMatrixType const hydrodynamic_dispersion =
633  velocity_magnitude != 0.0
634  ? GlobalDimMatrixType(porosity *
635  pore_diffusion_coefficient +
636  solute_dispersivity_transverse *
637  velocity_magnitude * I +
638  (solute_dispersivity_longitudinal -
639  solute_dispersivity_transverse) /
640  velocity_magnitude * velocity *
641  velocity.transpose())
643  pore_diffusion_coefficient +
644  solute_dispersivity_transverse *
645  velocity_magnitude * I);
646  const double R_times_phi(retardation_factor * porosity);
647  GlobalDimVectorType const mass_density_flow = velocity * density;
648  auto const N_t_N = (N.transpose() * N).eval();
650  {
651  MCp.noalias() += N_t_N * (C_int_pt * R_times_phi * drho_dp * w);
652  MCC.noalias() += N_t_N * (C_int_pt * R_times_phi * drho_dC * w);
653  KCC.noalias() -= dNdx.transpose() * mass_density_flow * N * w;
654  }
655  else
656  {
657  KCC.noalias() +=
658  N.transpose() * mass_density_flow.transpose() * dNdx * w;
659  }
660  MCC.noalias() += N_t_N * (R_times_phi * density * w);
661  KCC.noalias() += dNdx.transpose() * hydrodynamic_dispersion * dNdx *
662  (density * w) +
663  N_t_N * (decay_rate * R_times_phi * density * w);
664 
665  MpC.noalias() += N_t_N * (porosity * drho_dC * w);
666 
667  // Calculate Mpp, Kpp, and bp in the first loop over components
668  if (component_id == 0)
669  {
670  Mpp.noalias() += N_t_N * (porosity * drho_dp * w);
671  Kpp.noalias() +=
672  dNdx.transpose() * K_over_mu * dNdx * (density * w);
673 
675  {
676  Bp.noalias() += dNdx.transpose() * K_over_mu * b *
677  (density * density * w);
678  }
679  }
680  }
681  }
std::string getName(std::string const &line)
Returns the name/title from the "Zone"-description.
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
std::array< VariableType, static_cast< int >(Variable::number_of_variables)> VariableArray
Definition: VariableType.h:108
@ longitudinal_dispersivity
used to compute the hydrodynamic dispersion tensor.
Definition: PropertyType.h:58
@ transversal_dispersivity
used to compute the hydrodynamic dispersion tensor.
Definition: PropertyType.h:100
@ retardation_factor
specify retardation factor used in component transport process.
Definition: PropertyType.h:81
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_transport_process_variables, MaterialPropertyLib::concentration, MaterialPropertyLib::decay_rate, MaterialPropertyLib::density, MeshLib::Element::getID(), getName(), ProcessLib::ComponentTransport::ComponentTransportProcessData::has_gravity, MaterialPropertyLib::longitudinal_dispersivity, ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, ProcessLib::ComponentTransport::ComponentTransportProcessData::non_advective_form, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::pore_diffusion, MaterialPropertyLib::porosity, MaterialPropertyLib::retardation_factor, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), ProcessLib::ComponentTransport::ComponentTransportProcessData::specific_body_force, MaterialPropertyLib::transversal_dispersivity, and MaterialPropertyLib::viscosity.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assemble().

◆ assembleComponentTransportEquation()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleComponentTransportEquation ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &  local_xdot,
std::vector< double > &  local_M_data,
std::vector< double > &  local_K_data,
std::vector< double > &  ,
int const  transport_process_id 
)
inline

Definition at line 830 of file ComponentTransportFEM.h.

835  {
836  auto const local_p =
837  local_x.template segment<pressure_size>(pressure_index);
838  auto const local_C = local_x.template segment<concentration_size>(
840  (transport_process_id - 1) * concentration_size);
841  auto const local_pdot =
842  local_xdot.segment<pressure_size>(pressure_index);
843 
844  NodalVectorType local_T;
846  {
847  local_T =
849  }
850 
851  auto local_M = MathLib::createZeroedMatrix<LocalBlockMatrixType>(
852  local_M_data, concentration_size, concentration_size);
853  auto local_K = MathLib::createZeroedMatrix<LocalBlockMatrixType>(
854  local_K_data, concentration_size, concentration_size);
855 
856  unsigned const n_integration_points =
857  _integration_method.getNumberOfPoints();
858 
860  pos.setElementID(_element.getID());
861 
862  auto const& b = _process_data.specific_body_force;
863 
866 
867  GlobalDimMatrixType const& I(
868  GlobalDimMatrixType::Identity(GlobalDim, GlobalDim));
869 
870  auto const& medium =
871  *_process_data.media_map->getMedium(_element.getID());
872  auto const& phase = medium.phase("AqueousLiquid");
873  // Hydraulic process id is 0 and thus transport process id starts
874  // from 1.
875  auto const component_id = transport_process_id - 1;
876  auto const& component = phase.component(
877  _transport_process_variables[component_id].get().getName());
878 
879  for (unsigned ip(0); ip < n_integration_points; ++ip)
880  {
881  pos.setIntegrationPoint(ip);
882 
883  auto& ip_data = _ip_data[ip];
884  auto const& N = ip_data.N;
885  auto const& dNdx = ip_data.dNdx;
886  auto const& w = ip_data.integration_weight;
887  auto& porosity = ip_data.porosity;
888  auto const& porosity_prev = ip_data.porosity_prev;
889 
890  double C_int_pt = 0.0;
891  double p_int_pt = 0.0;
892 
893  NumLib::shapeFunctionInterpolate(local_C, N, C_int_pt);
894  NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
895 
896  vars[static_cast<int>(
898  vars[static_cast<int>(
900 
902  {
903  vars[static_cast<int>(
905  N.dot(local_T);
906  }
907 
908  // porosity
909  {
910  vars_prev[static_cast<int>(
911  MaterialPropertyLib::Variable::porosity)] = porosity_prev;
912 
913  porosity =
915  ? porosity_prev
917  .template value<double>(vars, vars_prev, pos, t,
918  dt);
919 
920  vars[static_cast<int>(
922  }
923 
924  auto const& retardation_factor =
926  .template value<double>(vars, pos, t, dt);
927 
928  auto const& solute_dispersivity_transverse = medium.template value<
929  double>(
931  auto const& solute_dispersivity_longitudinal =
932  medium.template value<double>(
935 
936  // Use the fluid density model to compute the density
937  auto const density =
939  .template value<double>(vars, pos, t, dt);
940  auto const decay_rate =
942  .template value<double>(vars, pos, t, dt);
943 
944  auto const& pore_diffusion_coefficient =
945  MaterialPropertyLib::formEigenTensor<GlobalDim>(
947  .value(vars, pos, t, dt));
948 
949  auto const& K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
951  vars, pos, t, dt));
952  // Use the viscosity model to compute the viscosity
953  auto const mu = phase[MaterialPropertyLib::PropertyType::viscosity]
954  .template value<double>(vars, pos, t, dt);
955 
956  GlobalDimMatrixType const K_over_mu = K / mu;
957  GlobalDimVectorType const velocity =
959  ? GlobalDimVectorType(-K_over_mu *
960  (dNdx * local_p - density * b))
961  : GlobalDimVectorType(-K_over_mu * dNdx * local_p);
962 
963  double const velocity_magnitude = velocity.norm();
964  GlobalDimMatrixType const hydrodynamic_dispersion =
965  velocity_magnitude != 0.0
966  ? GlobalDimMatrixType(porosity *
967  pore_diffusion_coefficient +
968  solute_dispersivity_transverse *
969  velocity_magnitude * I +
970  (solute_dispersivity_longitudinal -
971  solute_dispersivity_transverse) /
972  velocity_magnitude * velocity *
973  velocity.transpose())
975  pore_diffusion_coefficient +
976  solute_dispersivity_transverse *
977  velocity_magnitude * I);
978 
979  double const R_times_phi = retardation_factor * porosity;
980  auto const N_t_N = (N.transpose() * N).eval();
981 
983  {
984  const double drho_dC =
986  .template dValue<double>(
988  pos, t, dt);
989  local_M.noalias() +=
990  N_t_N * (R_times_phi * C_int_pt * drho_dC * w);
991  }
992 
993  local_M.noalias() += N_t_N * (R_times_phi * density * w);
994 
995  // coupling term
996 
998  {
999  double dot_p_int_pt = 0.0;
1000 
1001  NumLib::shapeFunctionInterpolate(local_pdot, N, dot_p_int_pt);
1002  const double drho_dp =
1004  .template dValue<double>(
1006  pos, t, dt);
1007 
1008  local_K.noalias() +=
1009  N_t_N * ((R_times_phi * drho_dp * dot_p_int_pt) * w) -
1010  dNdx.transpose() * velocity * N * (density * w);
1011  }
1012  else
1013  {
1014  local_K.noalias() +=
1015  N.transpose() * velocity.transpose() * dNdx * (density * w);
1016  }
1017  local_K.noalias() +=
1018  dNdx.transpose() * hydrodynamic_dispersion * dNdx *
1019  (density * w) +
1020  N_t_N * (decay_rate * R_times_phi * density * w);
1021  }
1022  }
typename ShapeMatricesType::NodalVectorType NodalVectorType
virtual Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > getNodalValuesOnElement(MeshLib::Element const &element, double const t) const
Returns a matrix of values for all nodes of the given element.
Definition: Parameter.h:163

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_transport_process_variables, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemically_induced_porosity_change, MaterialPropertyLib::concentration, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, MaterialPropertyLib::decay_rate, MaterialPropertyLib::density, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), getName(), ParameterLib::Parameter< T >::getNodalValuesOnElement(), ProcessLib::ComponentTransport::ComponentTransportProcessData::has_gravity, MaterialPropertyLib::longitudinal_dispersivity, ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, ProcessLib::ComponentTransport::ComponentTransportProcessData::non_advective_form, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::pore_diffusion, MaterialPropertyLib::porosity, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_index, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_size, MaterialPropertyLib::retardation_factor, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), ProcessLib::ComponentTransport::ComponentTransportProcessData::specific_body_force, MaterialPropertyLib::temperature, ProcessLib::ComponentTransport::ComponentTransportProcessData::temperature, MaterialPropertyLib::transversal_dispersivity, and MaterialPropertyLib::viscosity.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleForStaggeredScheme().

◆ assembleForStaggeredScheme()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::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 
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 683 of file ComponentTransportFEM.h.

690  {
691  if (process_id == _process_data.hydraulic_process_id)
692  {
693  assembleHydraulicEquation(t, dt, local_x, local_xdot, local_M_data,
694  local_K_data, local_b_data);
695  }
696  else
697  {
698  // Go for assembling in an order of transport process id.
699  assembleComponentTransportEquation(t, dt, local_x, local_xdot,
700  local_M_data, local_K_data,
701  local_b_data, process_id);
702  }
703  }
void assembleHydraulicEquation(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void assembleComponentTransportEquation(double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &, int const transport_process_id)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleComponentTransportEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleHydraulicEquation(), and ProcessLib::ComponentTransport::ComponentTransportProcessData::hydraulic_process_id.

◆ assembleHydraulicEquation()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleHydraulicEquation ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &  local_xdot,
std::vector< double > &  local_M_data,
std::vector< double > &  local_K_data,
std::vector< double > &  local_b_data 
)
inline

Definition at line 705 of file ComponentTransportFEM.h.

712  {
713  auto const local_p =
714  local_x.template segment<pressure_size>(pressure_index);
715  auto const local_C = local_x.template segment<concentration_size>(
717  auto const local_Cdot =
718  local_xdot.segment<concentration_size>(first_concentration_index);
719 
720  auto local_M = MathLib::createZeroedMatrix<LocalBlockMatrixType>(
721  local_M_data, pressure_size, pressure_size);
722  auto local_K = MathLib::createZeroedMatrix<LocalBlockMatrixType>(
723  local_K_data, pressure_size, pressure_size);
724  auto local_b = MathLib::createZeroedVector<LocalSegmentVectorType>(
725  local_b_data, pressure_size);
726 
727  unsigned const n_integration_points =
728  _integration_method.getNumberOfPoints();
729 
731  pos.setElementID(_element.getID());
732 
733  auto const& b = _process_data.specific_body_force;
734 
735  auto const& medium =
736  *_process_data.media_map->getMedium(_element.getID());
737  auto const& phase = medium.phase("AqueousLiquid");
738 
741 
742  for (unsigned ip(0); ip < n_integration_points; ++ip)
743  {
744  pos.setIntegrationPoint(ip);
745 
746  auto& ip_data = _ip_data[ip];
747  auto const& N = ip_data.N;
748  auto const& dNdx = ip_data.dNdx;
749  auto const& w = ip_data.integration_weight;
750  auto& porosity = ip_data.porosity;
751  auto const& porosity_prev = ip_data.porosity_prev;
752 
753  double C_int_pt = 0.0;
754  double p_int_pt = 0.0;
755 
756  NumLib::shapeFunctionInterpolate(local_C, N, C_int_pt);
757  NumLib::shapeFunctionInterpolate(local_p, N, p_int_pt);
758 
759  vars[static_cast<int>(
761  vars[static_cast<int>(
763 
764  // porosity
765  {
766  vars_prev[static_cast<int>(
767  MaterialPropertyLib::Variable::porosity)] = porosity_prev;
768 
769  porosity =
771  ? porosity_prev
773  .template value<double>(vars, vars_prev, pos, t,
774  dt);
775 
776  vars[static_cast<int>(
778  }
779 
780  // Use the fluid density model to compute the density
781  // TODO: Concentration of which component as one of arguments for
782  // calculation of fluid density
783  auto const density =
785  .template value<double>(vars, pos, t, dt);
786 
787  auto const& K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
789  vars, pos, t, dt));
790 
791  // Use the viscosity model to compute the viscosity
792  auto const mu = phase[MaterialPropertyLib::PropertyType::viscosity]
793  .template value<double>(vars, pos, t, dt);
794 
795  GlobalDimMatrixType const K_over_mu = K / mu;
796 
797  const double drho_dp =
799  .template dValue<double>(
801  pos, t, dt);
802  const double drho_dC =
804  .template dValue<double>(
806  t, dt);
807 
808  // matrix assembly
809  local_M.noalias() += w * N.transpose() * porosity * drho_dp * N;
810  local_K.noalias() +=
811  w * dNdx.transpose() * density * K_over_mu * dNdx;
812 
814  {
815  local_b.noalias() +=
816  w * density * density * dNdx.transpose() * K_over_mu * b;
817  }
818 
819  // coupling term
820  {
821  double dot_C_int_pt = 0.0;
822  NumLib::shapeFunctionInterpolate(local_Cdot, N, dot_C_int_pt);
823 
824  local_b.noalias() -=
825  w * N.transpose() * porosity * drho_dC * dot_C_int_pt;
826  }
827  }
828  }

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemically_induced_porosity_change, MaterialPropertyLib::concentration, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, MaterialPropertyLib::density, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), ProcessLib::ComponentTransport::ComponentTransportProcessData::has_gravity, ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::porosity, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_index, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_size, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::shapeFunctionInterpolate(), ProcessLib::ComponentTransport::ComponentTransportProcessData::specific_body_force, and MaterialPropertyLib::viscosity.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleForStaggeredScheme().

◆ assembleReactionEquationConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleReactionEquationConcrete ( double const  t,
double const  dt,
Eigen::VectorXd const &  local_x,
std::vector< double > &  local_M_data,
std::vector< double > &  local_K_data,
std::vector< double > &  local_b_data,
int const  transport_process_id 
)
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 1024 of file ComponentTransportFEM.h.

1029  {
1030  auto const local_C = local_x.template segment<concentration_size>(
1032  (transport_process_id - 1) * concentration_size);
1033 
1034  auto local_M = MathLib::createZeroedMatrix<LocalBlockMatrixType>(
1035  local_M_data, concentration_size, concentration_size);
1036  auto local_K = MathLib::createZeroedMatrix<LocalBlockMatrixType>(
1037  local_K_data, concentration_size, concentration_size);
1038  auto local_b = MathLib::createZeroedVector<LocalVectorType>(
1039  local_b_data, concentration_size);
1040 
1041  unsigned const n_integration_points =
1042  _integration_method.getNumberOfPoints();
1043 
1045  pos.setElementID(_element.getID());
1046 
1049 
1050  auto const& medium =
1051  *_process_data.media_map->getMedium(_element.getID());
1052  auto const component_id = transport_process_id - 1;
1053  for (unsigned ip(0); ip < n_integration_points; ++ip)
1054  {
1055  pos.setIntegrationPoint(ip);
1056 
1057  auto& ip_data = _ip_data[ip];
1058  auto const& N = ip_data.N;
1059  auto const w = ip_data.integration_weight;
1060  auto& porosity = ip_data.porosity;
1061  auto const& porosity_prev = ip_data.porosity_prev;
1062  auto const chemical_system_id = ip_data.chemical_system_id;
1063 
1064  double C_int_pt = 0.0;
1065  NumLib::shapeFunctionInterpolate(local_C, N, C_int_pt);
1066 
1067  vars[static_cast<int>(
1069 
1070  auto const porosity_dot = (porosity - porosity_prev) / dt;
1071 
1072  // porosity
1073  {
1074  vars_prev[static_cast<int>(
1075  MaterialPropertyLib::Variable::porosity)] = porosity_prev;
1076 
1077  porosity =
1079  ? porosity_prev
1081  .template value<double>(vars, vars_prev, pos, t,
1082  dt);
1083  }
1084 
1085  local_M.noalias() += w * N.transpose() * porosity * N;
1086 
1087  local_K.noalias() += w * N.transpose() * porosity_dot * N;
1088 
1089  auto const C_post_int_pt =
1091  component_id, chemical_system_id);
1092 
1093  local_b.noalias() +=
1094  w * N.transpose() * porosity * (C_post_int_pt - C_int_pt) / dt;
1095  }
1096  }
virtual double getConcentration(int const, GlobalIndexType const) const
ChemistryLib::ChemicalSolverInterface *const chemical_solver_interface

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemically_induced_porosity_change, MaterialPropertyLib::concentration, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, ChemistryLib::ChemicalSolverInterface::getConcentration(), MeshLib::Element::getID(), ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, MaterialPropertyLib::porosity, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), and NumLib::shapeFunctionInterpolate().

◆ calculateIntPtDarcyVelocity()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector<double> const& ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::calculateIntPtDarcyVelocity ( const double  t,
Eigen::Ref< const NodalVectorType > const &  p_nodal_values,
Eigen::Ref< const NodalVectorType > const &  C_nodal_values,
std::vector< double > &  cache 
) const
inline

Definition at line 1140 of file ComponentTransportFEM.h.

1145  {
1146  auto const n_integration_points =
1147  _integration_method.getNumberOfPoints();
1148 
1149  cache.clear();
1150  auto cache_mat = MathLib::createZeroedMatrix<
1151  Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
1152  cache, GlobalDim, n_integration_points);
1153 
1155  pos.setElementID(_element.getID());
1156 
1158 
1159  auto const& medium =
1160  *_process_data.media_map->getMedium(_element.getID());
1161  auto const& phase = medium.phase("AqueousLiquid");
1162 
1163  for (unsigned ip = 0; ip < n_integration_points; ++ip)
1164  {
1165  auto const& ip_data = _ip_data[ip];
1166  auto const& N = ip_data.N;
1167  auto const& dNdx = ip_data.dNdx;
1168  auto const& porosity = ip_data.porosity;
1169 
1170  pos.setIntegrationPoint(ip);
1171 
1172  double C_int_pt = 0.0;
1173  double p_int_pt = 0.0;
1174 
1175  NumLib::shapeFunctionInterpolate(C_nodal_values, N, C_int_pt);
1176  NumLib::shapeFunctionInterpolate(p_nodal_values, N, p_int_pt);
1177 
1178  vars[static_cast<int>(
1180  vars[static_cast<int>(
1182  vars[static_cast<int>(MaterialPropertyLib::Variable::porosity)] =
1183  porosity;
1184 
1185  // TODO (naumov) Temporary value not used by current material
1186  // models. Need extension of secondary variables interface.
1187  double const dt = std::numeric_limits<double>::quiet_NaN();
1188  auto const& K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
1190  vars, pos, t, dt));
1191  auto const mu = phase[MaterialPropertyLib::PropertyType::viscosity]
1192  .template value<double>(vars, pos, t, dt);
1193  GlobalDimMatrixType const K_over_mu = K / mu;
1194 
1195  cache_mat.col(ip).noalias() = -K_over_mu * dNdx * p_nodal_values;
1197  {
1198  auto const rho_w =
1200  .template value<double>(vars, pos, t, dt);
1201  auto const b = _process_data.specific_body_force;
1202  // here it is assumed that the vector b is directed 'downwards'
1203  cache_mat.col(ip).noalias() += K_over_mu * rho_w * b;
1204  }
1205  }
1206 
1207  return cache;
1208  }
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:32

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, MaterialPropertyLib::concentration, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MeshLib::Element::getID(), ProcessLib::ComponentTransport::ComponentTransportProcessData::has_gravity, ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::porosity, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), ProcessLib::ComponentTransport::ComponentTransportProcessData::specific_body_force, and MaterialPropertyLib::viscosity.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete(), and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getIntPtDarcyVelocity().

◆ computeSecondaryVariableConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete ( double const  t,
double const  ,
Eigen::VectorXd const &  local_x,
Eigen::VectorXd const &   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1319 of file ComponentTransportFEM.h.

1324  {
1325  auto const local_p =
1326  local_x.template segment<pressure_size>(pressure_index);
1327  auto const local_C = local_x.template segment<concentration_size>(
1329 
1330  std::vector<double> ele_velocity;
1331  calculateIntPtDarcyVelocity(t, local_p, local_C, ele_velocity);
1332 
1333  auto const n_integration_points =
1334  _integration_method.getNumberOfPoints();
1335  auto const ele_velocity_mat =
1336  MathLib::toMatrix(ele_velocity, GlobalDim, n_integration_points);
1337 
1338  auto const ele_id = _element.getID();
1339  Eigen::Map<LocalVectorType>(
1340  &(*_process_data.mesh_prop_velocity)[ele_id * GlobalDim],
1341  GlobalDim) =
1342  ele_velocity_mat.rowwise().sum() / n_integration_points;
1343 
1345  {
1347  {
1348  auto const& medium =
1349  *_process_data.media_map->getMedium(ele_id);
1350 
1352  pos.setElementID(ele_id);
1353 
1354  for (auto& ip_data : _ip_data)
1355  {
1356  ip_data.porosity = ip_data.porosity_prev;
1357 
1359  ->updatePorosityPostReaction(ip_data.chemical_system_id,
1360  medium, ip_data.porosity);
1361  }
1362  }
1363 
1364  std::vector<GlobalIndexType> chemical_system_indices;
1365  chemical_system_indices.reserve(n_integration_points);
1366  std::transform(_ip_data.begin(), _ip_data.end(),
1367  std::back_inserter(chemical_system_indices),
1368  [](auto const& ip_data)
1369  { return ip_data.chemical_system_id; });
1370 
1372  ele_id, chemical_system_indices);
1373  }
1374  }
virtual void updatePorosityPostReaction(GlobalIndexType const &, MaterialPropertyLib::Medium const &, double &)
virtual void computeSecondaryVariable(std::size_t const, std::vector< GlobalIndexType > const &)
std::vector< double > const & calculateIntPtDarcyVelocity(const double t, Eigen::Ref< const NodalVectorType > const &p_nodal_values, Eigen::Ref< const NodalVectorType > const &C_nodal_values, std::vector< double > &cache) const
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
Definition: EigenMapTools.h:72

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::calculateIntPtDarcyVelocity(), ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemically_induced_porosity_change, ChemistryLib::ChemicalSolverInterface::computeSecondaryVariable(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, ProcessLib::ComponentTransport::ComponentTransportProcessData::mesh_prop_velocity, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_index, ParameterLib::SpatialPosition::setElementID(), MathLib::toMatrix(), and ChemistryLib::ChemicalSolverInterface::updatePorosityPostReaction().

◆ getFlux()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Eigen::Vector3d ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getFlux ( MathLib::Point3d const &  ,
double const  ,
std::vector< double > const &   
) const
inlineoverridevirtual

Computes the flux in the point p_local_coords that is given in local coordinates using the values from local_x. Fits to monolithic scheme.

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1219 of file ComponentTransportFEM.h.

1222  {
1223  auto const local_p = Eigen::Map<const NodalVectorType>(
1224  &local_x[pressure_index], pressure_size);
1225  auto const local_C = Eigen::Map<const NodalVectorType>(
1227 
1228  // Eval shape matrices at given point
1229  // Note: Axial symmetry is set to false here, because we only need dNdx
1230  // here, which is not affected by axial symmetry.
1231  auto const shape_matrices =
1233  GlobalDim>(
1234  _element, false /*is_axially_symmetric*/,
1235  std::array{pnt_local_coords})[0];
1236 
1238  pos.setElementID(_element.getID());
1239 
1241 
1242  auto const& medium =
1243  *_process_data.media_map->getMedium(_element.getID());
1244  auto const& phase = medium.phase("AqueousLiquid");
1245 
1246  // local_x contains the local concentration and pressure values
1247  double c_int_pt;
1248  NumLib::shapeFunctionInterpolate(local_C, shape_matrices.N, c_int_pt);
1249  vars[static_cast<int>(MaterialPropertyLib::Variable::concentration)]
1250  .emplace<double>(c_int_pt);
1251 
1252  double p_int_pt;
1253  NumLib::shapeFunctionInterpolate(local_p, shape_matrices.N, p_int_pt);
1254  vars[static_cast<int>(MaterialPropertyLib::Variable::phase_pressure)]
1255  .emplace<double>(p_int_pt);
1256 
1257  // TODO (naumov) Temporary value not used by current material models.
1258  // Need extension of secondary variables interface.
1259  double const dt = std::numeric_limits<double>::quiet_NaN();
1260  auto const K = MaterialPropertyLib::formEigenTensor<GlobalDim>(
1262  vars, pos, t, dt));
1263 
1264  auto const mu = phase[MaterialPropertyLib::PropertyType::viscosity]
1265  .template value<double>(vars, pos, t, dt);
1266  GlobalDimMatrixType const K_over_mu = K / mu;
1267 
1268  GlobalDimVectorType q = -K_over_mu * shape_matrices.dNdx * local_p;
1269  auto const rho_w = phase[MaterialPropertyLib::PropertyType::density]
1270  .template value<double>(vars, pos, t, dt);
1272  {
1273  auto const b = _process_data.specific_body_force;
1274  q += K_over_mu * rho_w * b;
1275  }
1276  Eigen::Vector3d flux(0.0, 0.0, 0.0);
1277  flux.head<GlobalDim>() = rho_w * q;
1278  return flux;
1279  }
static const double q
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::concentration, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, MaterialPropertyLib::density, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), ProcessLib::ComponentTransport::ComponentTransportProcessData::has_gravity, ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, MaterialPropertyLib::permeability, MaterialPropertyLib::phase_pressure, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_index, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_size, MathLib::q, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), ProcessLib::ComponentTransport::ComponentTransportProcessData::specific_body_force, and MaterialPropertyLib::viscosity.

◆ getInterpolatedLocalSolution()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector<double> const& ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getInterpolatedLocalSolution ( const double  ,
std::vector< GlobalVector * > const &  int_pt_x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  ,
std::vector< double > &  cache 
) const
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 1296 of file ComponentTransportFEM.h.

1301  {
1303  assert(int_pt_x.size() == 1);
1304 
1305  cache.clear();
1306 
1307  unsigned const n_integration_points =
1308  _integration_method.getNumberOfPoints();
1309  for (unsigned ip(0); ip < n_integration_points; ++ip)
1310  {
1311  auto const& chemical_system_id = _ip_data[ip].chemical_system_id;
1312  auto const c_int_pt = int_pt_x[0]->get(chemical_system_id);
1313  cache.push_back(c_int_pt);
1314  }
1315 
1316  return cache;
1317  }

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, and ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface.

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector<double> const& ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getIntPtDarcyVelocity ( const double  t,
std::vector< GlobalVector * > const &  x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  dof_table,
std::vector< double > &  cache 
) const
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 1098 of file ComponentTransportFEM.h.

1103  {
1104  assert(x.size() == dof_table.size());
1105 
1106  auto const n_processes = x.size();
1107  std::vector<std::vector<double>> local_x;
1108  local_x.reserve(n_processes);
1109 
1110  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
1111  {
1112  auto const indices =
1113  NumLib::getIndices(_element.getID(), *dof_table[process_id]);
1114  assert(!indices.empty());
1115  local_x.push_back(x[process_id]->get(indices));
1116  }
1117 
1118  // only one process, must be monolithic.
1119  if (n_processes == 1)
1120  {
1121  auto const local_p = Eigen::Map<const NodalVectorType>(
1122  &local_x[0][pressure_index], pressure_size);
1123  auto const local_C = Eigen::Map<const NodalVectorType>(
1125  return calculateIntPtDarcyVelocity(t, local_p, local_C, cache);
1126  }
1127 
1128  // multiple processes, must be staggered.
1129  {
1130  constexpr int pressure_process_id = 0;
1131  constexpr int concentration_process_id = 1;
1132  auto const local_p = Eigen::Map<const NodalVectorType>(
1133  &local_x[pressure_process_id][0], pressure_size);
1134  auto const local_C = Eigen::Map<const NodalVectorType>(
1135  &local_x[concentration_process_id][0], concentration_size);
1136  return calculateIntPtDarcyVelocity(t, local_p, local_C, cache);
1137  }
1138  }
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::calculateIntPtDarcyVelocity(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), NumLib::getIndices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_index, and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::pressure_size.

◆ getShapeMatrix()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Eigen::Map<const Eigen::RowVectorXd> ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getShapeMatrix ( const unsigned  integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 1210 of file ComponentTransportFEM.h.

1212  {
1213  auto const& N = _ip_data[integration_point].N;
1214 
1215  // assumes N is stored contiguously in memory
1216  return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
1217  }

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data.

◆ initializeChemicalSystemConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::initializeChemicalSystemConcrete ( Eigen::VectorXd const &  local_x,
double const  t 
)
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 313 of file ComponentTransportFEM.h.

315  {
317 
318  auto const& medium =
319  *_process_data.media_map->getMedium(_element.getID());
320 
322  pos.setElementID(_element.getID());
323 
324  unsigned const n_integration_points =
325  _integration_method.getNumberOfPoints();
326  for (unsigned ip = 0; ip < n_integration_points; ip++)
327  {
328  auto& ip_data = _ip_data[ip];
329  auto const& N = ip_data.N;
330  auto const& chemical_system_id = ip_data.chemical_system_id;
331 
332  auto const n_component = _transport_process_variables.size();
333  std::vector<double> C_int_pt(n_component);
334  for (unsigned component_id = 0; component_id < n_component;
335  ++component_id)
336  {
337  auto const concentration_index =
339  component_id * concentration_size;
340  auto const local_C =
341  local_x.template segment<concentration_size>(
342  concentration_index);
343 
345  C_int_pt[component_id]);
346  }
347 
349  ->initializeChemicalSystemConcrete(C_int_pt, chemical_system_id,
350  medium, pos, t);
351  }
352  }
virtual void initializeChemicalSystemConcrete(std::vector< double > const &, GlobalIndexType const &, MaterialPropertyLib::Medium const &, ParameterLib::SpatialPosition const &, double const)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_transport_process_variables, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), ChemistryLib::ChemicalSolverInterface::initializeChemicalSystemConcrete(), ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, ParameterLib::SpatialPosition::setElementID(), and NumLib::shapeFunctionInterpolate().

◆ interpolateNodalValuesToIntegrationPoints()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector<double> ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::interpolateNodalValuesToIntegrationPoints ( std::vector< double > const &  local_x)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1281 of file ComponentTransportFEM.h.

1283  {
1284  unsigned const n_integration_points =
1285  _integration_method.getNumberOfPoints();
1286 
1287  std::vector<double> interpolated_values(n_integration_points);
1288  for (unsigned ip(0); ip < n_integration_points; ++ip)
1289  {
1291  interpolated_values[ip]);
1292  }
1293  return interpolated_values;
1294  }

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, and NumLib::shapeFunctionInterpolate().

◆ postSpeciationCalculation()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::postSpeciationCalculation ( std::size_t const  ele_id,
double const  t,
double const  dt 
)
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 416 of file ComponentTransportFEM.h.

418  {
420  {
421  return;
422  }
423 
424  auto const& medium = *_process_data.media_map->getMedium(ele_id);
425 
427  pos.setElementID(ele_id);
428 
429  for (auto& ip_data : _ip_data)
430  {
431  ip_data.porosity = ip_data.porosity_prev;
432 
434  ->updateVolumeFractionPostReaction(ip_data.chemical_system_id,
435  medium, pos,
436  ip_data.porosity, t, dt);
437 
439  ip_data.chemical_system_id, medium, ip_data.porosity);
440  }
441  }
virtual void updateVolumeFractionPostReaction(GlobalIndexType const &, MaterialPropertyLib::Medium const &, ParameterLib::SpatialPosition const &, double const, double const, double const)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemically_induced_porosity_change, ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, ParameterLib::SpatialPosition::setElementID(), ChemistryLib::ChemicalSolverInterface::updatePorosityPostReaction(), and ChemistryLib::ChemicalSolverInterface::updateVolumeFractionPostReaction().

◆ postTimestepConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::postTimestepConcrete ( Eigen::VectorXd const &  ,
double const  ,
double const   
)
inlineoverridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 1376 of file ComponentTransportFEM.h.

1378  {
1379  unsigned const n_integration_points =
1380  _integration_method.getNumberOfPoints();
1381 
1382  for (unsigned ip = 0; ip < n_integration_points; ip++)
1383  {
1384  _ip_data[ip].pushBackState();
1385  }
1386  }

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data.

◆ setChemicalSystemConcrete()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemConcrete ( Eigen::VectorXd const &  local_x,
double const  t,
double  dt 
)
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 354 of file ComponentTransportFEM.h.

356  {
358 
359  auto const& medium =
360  _process_data.media_map->getMedium(_element.getID());
361 
364 
366  pos.setElementID(_element.getID());
367 
368  unsigned const n_integration_points =
369  _integration_method.getNumberOfPoints();
370  for (unsigned ip = 0; ip < n_integration_points; ip++)
371  {
372  auto& ip_data = _ip_data[ip];
373  auto const& N = ip_data.N;
374  auto& porosity = ip_data.porosity;
375  auto const& porosity_prev = ip_data.porosity_prev;
376  auto const& chemical_system_id = ip_data.chemical_system_id;
377 
378  auto const n_component = _transport_process_variables.size();
379  std::vector<double> C_int_pt(n_component);
380  for (unsigned component_id = 0; component_id < n_component;
381  ++component_id)
382  {
383  auto const concentration_index =
385  component_id * concentration_size;
386  auto const local_C =
387  local_x.template segment<concentration_size>(
388  concentration_index);
389 
391  C_int_pt[component_id]);
392  }
393 
394  {
395  vars_prev[static_cast<int>(
396  MaterialPropertyLib::Variable::porosity)] = porosity_prev;
397 
398  porosity =
400  ? porosity_prev
401  : medium
402  ->property(
404  .template value<double>(vars, vars_prev, pos, t,
405  dt);
406 
407  vars[static_cast<int>(
409  }
410 
412  C_int_pt, chemical_system_id, medium, vars, pos, t, dt);
413  }
414  }
virtual void setChemicalSystemConcrete(std::vector< double > const &, GlobalIndexType const &, MaterialPropertyLib::Medium const *, MaterialPropertyLib::VariableArray const &, ParameterLib::SpatialPosition const &, double const, double const)

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_transport_process_variables, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemically_induced_porosity_change, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::concentration_size, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::first_concentration_index, MeshLib::Element::getID(), ProcessLib::ComponentTransport::ComponentTransportProcessData::media_map, MaterialPropertyLib::porosity, ChemistryLib::ChemicalSolverInterface::setChemicalSystemConcrete(), ParameterLib::SpatialPosition::setElementID(), and NumLib::shapeFunctionInterpolate().

◆ setChemicalSystemID()

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
void ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemID ( std::size_t const  )
inlineoverridevirtual

Implements ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface.

Definition at line 293 of file ComponentTransportFEM.h.

294  {
296  // chemical system index map
297  auto& chemical_system_index_map =
299 
300  unsigned const n_integration_points =
301  _integration_method.getNumberOfPoints();
302  for (unsigned ip = 0; ip < n_integration_points; ip++)
303  {
304  _ip_data[ip].chemical_system_id =
305  chemical_system_index_map.empty()
306  ? 0
307  : chemical_system_index_map.back() + 1;
308  chemical_system_index_map.push_back(
309  _ip_data[ip].chemical_system_id);
310  }
311  }
std::vector< GlobalIndexType > chemical_system_index_map

References ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data, ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::ComponentTransport::ComponentTransportProcessData::chemical_solver_interface, and ChemistryLib::ChemicalSolverInterface::chemical_system_index_map.

Member Data Documentation

◆ _element

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
MeshLib::Element const& ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element
private

◆ _integration_method

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
IntegrationMethod const ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method
private

Definition at line 1392 of file ComponentTransportFEM.h.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalAssemblerData(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleBlockMatrices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleComponentTransportEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleHydraulicEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleReactionEquationConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::calculateIntPtDarcyVelocity(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getInterpolatedLocalSolution(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::initializeChemicalSystemConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::interpolateNodalValuesToIntegrationPoints(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::postTimestepConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemConcrete(), and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemID().

◆ _ip_data

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
std::vector< IntegrationPointData<NodalRowVectorType, GlobalDimNodalMatrixType>, Eigen::aligned_allocator< IntegrationPointData<NodalRowVectorType, GlobalDimNodalMatrixType> > > ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_ip_data
private

Definition at line 1400 of file ComponentTransportFEM.h.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalAssemblerData(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleBlockMatrices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleComponentTransportEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleHydraulicEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleReactionEquationConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::calculateIntPtDarcyVelocity(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getInterpolatedLocalSolution(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getShapeMatrix(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::initializeChemicalSystemConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::interpolateNodalValuesToIntegrationPoints(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::postSpeciationCalculation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::postTimestepConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemConcrete(), and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemID().

◆ _process_data

template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
ComponentTransportProcessData const& ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data
private

Definition at line 1390 of file ComponentTransportFEM.h.

Referenced by ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::LocalAssemblerData(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleBlockMatrices(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleComponentTransportEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleForStaggeredScheme(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleHydraulicEquation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleReactionEquationConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::calculateIntPtDarcyVelocity(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getFlux(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getInterpolatedLocalSolution(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::initializeChemicalSystemConcrete(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::postSpeciationCalculation(), ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemConcrete(), and ProcessLib::ComponentTransport::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::setChemicalSystemID().

◆ _transport_process_variables

◆ concentration_size

◆ first_concentration_index

◆ pressure_index

◆ pressure_size


The documentation for this class was generated from the following file: