OGS
ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface Class Referenceabstract

Detailed Description

Definition at line 56 of file ComponentTransportFEM.h.

#include <ComponentTransportFEM.h>

Inheritance diagram for ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface:
[legend]
Collaboration diagram for ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface:
[legend]

Public Member Functions

 ComponentTransportLocalAssemblerInterface ()=default
virtual void setChemicalSystemID (std::size_t const)=0
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)
virtual void postSpeciationCalculation (std::size_t const ele_id, double const t, double const dt)=0
virtual void computeReactionRelatedSecondaryVariable (std::size_t const ele_id)=0
virtual std::vector< double > const & getIntPtLiquidDensity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
virtual 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 =0
virtual std::vector< double > const & getIntPtMolarFlux (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache, int const component_id) const =0
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
virtual ~LocalAssemblerInterface ()=default
virtual void setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, 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 assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
virtual void assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_x_prev, 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_x_prev, 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_x_prev, int const process_id, 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_prev, 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, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
void postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
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.
virtual int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix (const unsigned integration_point) const =0
 Provides the shape matrix at the given integration point.
virtual ~ExtrapolatableElement ()=default

Private Member Functions

virtual void initializeChemicalSystemConcrete (Eigen::VectorXd const &, double const)=0
virtual void setChemicalSystemConcrete (Eigen::VectorXd const &, double const, double const)=0
virtual 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)=0

Constructor & Destructor Documentation

◆ ComponentTransportLocalAssemblerInterface()

ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::ComponentTransportLocalAssemblerInterface ( )
default

Member Function Documentation

◆ assembleReactionEquation()

void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::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 )
inline

Definition at line 111 of file ComponentTransportFEM.h.

116 {
117 std::vector<double> local_x_vec;
118
119 auto const n_processes = x.size();
120 for (std::size_t pcs_id = 0; pcs_id < n_processes; ++pcs_id)
121 {
122 auto const indices =
123 NumLib::getIndices(mesh_item_id, *dof_tables[pcs_id]);
124 assert(!indices.empty());
125 auto const local_solution = x[pcs_id]->get(indices);
126 local_x_vec.insert(std::end(local_x_vec),
127 std::begin(local_solution),
128 std::end(local_solution));
129 }
130 auto const local_x = MathLib::toVector(local_x_vec);
131
132 auto const indices =
133 NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
134 auto const num_r_c = indices.size();
135
136 std::vector<double> local_M_data;
137 local_M_data.reserve(num_r_c * num_r_c);
138 std::vector<double> local_K_data;
139 local_K_data.reserve(num_r_c * num_r_c);
140 std::vector<double> local_b_data;
141 local_b_data.reserve(num_r_c);
142
143 assembleReactionEquationConcrete(t, dt, local_x, local_M_data,
144 local_K_data, local_b_data,
145 process_id);
146
147 auto const r_c_indices =
149 if (!local_M_data.empty())
150 {
151 auto const local_M =
152 MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
153 M.add(r_c_indices, local_M);
154 }
155 if (!local_K_data.empty())
156 {
157 auto const local_K =
158 MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
159 K.add(r_c_indices, local_K);
160 }
161 if (!local_b_data.empty())
162 {
163 b.add(indices, local_b_data);
164 }
165 }
int add(IndexType row, IndexType col, double val)
Definition EigenMatrix.h:80
void add(IndexType rowId, double v)
add entry
Definition EigenVector.h:70
MathLib::RowColumnIndices< GlobalIndexType > RowColumnIndices
virtual 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)=0
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &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 MathLib::EigenMatrix::add(), MathLib::EigenVector::add(), assembleReactionEquationConcrete(), NumLib::getIndices(), MathLib::toMatrix(), and MathLib::toVector().

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation().

◆ assembleReactionEquationConcrete()

virtual void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::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 )
privatepure virtual

◆ computeReactionRelatedSecondaryVariable()

virtual void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::computeReactionRelatedSecondaryVariable ( std::size_t const ele_id)
pure virtual

◆ getIntPtDarcyVelocity()

virtual std::vector< double > const & ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtLiquidDensity()

virtual std::vector< double > const & ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtLiquidDensity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache ) const
pure virtual

◆ getIntPtMolarFlux()

virtual std::vector< double > const & ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getIntPtMolarFlux ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & cache,
int const component_id ) const
pure virtual

◆ initializeChemicalSystem()

void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::initializeChemicalSystem ( std::size_t const mesh_item_id,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_tables,
std::vector< GlobalVector * > const & x,
double const t )
inline

Definition at line 65 of file ComponentTransportFEM.h.

69 {
70 std::vector<double> local_x_vec;
71
72 auto const n_processes = x.size();
73 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
74 {
75 auto const indices =
76 NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
77 assert(!indices.empty());
78 auto const local_solution = x[process_id]->get(indices);
79 local_x_vec.insert(std::end(local_x_vec),
80 std::begin(local_solution),
81 std::end(local_solution));
82 }
83 auto const local_x = MathLib::toVector(local_x_vec);
84
86 }
virtual void initializeChemicalSystemConcrete(Eigen::VectorXd const &, double const)=0

References NumLib::getIndices(), initializeChemicalSystemConcrete(), and MathLib::toVector().

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::setInitialConditionsConcreteProcess().

◆ initializeChemicalSystemConcrete()

virtual void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::initializeChemicalSystemConcrete ( Eigen::VectorXd const & ,
double const  )
privatepure virtual

◆ postSpeciationCalculation()

virtual void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::postSpeciationCalculation ( std::size_t const ele_id,
double const t,
double const dt )
pure virtual

◆ setChemicalSystem()

void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::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 )
inline

Definition at line 88 of file ComponentTransportFEM.h.

92 {
93 std::vector<double> local_x_vec;
94
95 auto const n_processes = x.size();
96 for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
97 {
98 auto const indices =
99 NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
100 assert(!indices.empty());
101 auto const local_solution = x[process_id]->get(indices);
102 local_x_vec.insert(std::end(local_x_vec),
103 std::begin(local_solution),
104 std::end(local_solution));
105 }
106 auto const local_x = MathLib::toVector(local_x_vec);
107
108 setChemicalSystemConcrete(local_x, t, dt);
109 }
virtual void setChemicalSystemConcrete(Eigen::VectorXd const &, double const, double const)=0

References NumLib::getIndices(), setChemicalSystemConcrete(), and MathLib::toVector().

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::solveReactionEquation().

◆ setChemicalSystemConcrete()

virtual void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystemConcrete ( Eigen::VectorXd const & ,
double const ,
double const  )
privatepure virtual

◆ setChemicalSystemID()

virtual void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setChemicalSystemID ( std::size_t const )
pure virtual

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