OGS
ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface Class Referenceabstract

Detailed Description

Definition at line 62 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
 
void setStaggeredCoupledSolutions (std::size_t const, CoupledSolutionsForStaggeredScheme *const coupling_term)
 
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 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 & getInterpolatedLocalSolution (const double, std::vector< GlobalVector * > const &int_pt_x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &cache) const =0
 
- 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 assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, 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_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
 
virtual void assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data)
 
virtual void computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id)
 
virtual void preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
 
virtual void postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)
 
void postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
 
virtual std::vector< double > interpolateNodalValuesToIntegrationPoints (std::vector< double > const &)
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const
 
virtual Eigen::Vector3d getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const
 Fits to staggered scheme. More...
 
- 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. More...
 
virtual ~ExtrapolatableElement ()=default
 

Protected Attributes

CoupledSolutionsForStaggeredScheme_coupled_solutions {nullptr}
 

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 124 of file ComponentTransportFEM.h.

129  {
130  std::vector<double> local_x_vec;
131 
132  auto const n_processes = x.size();
133  for (std::size_t pcs_id = 0; pcs_id < n_processes; ++pcs_id)
134  {
135  auto const indices =
136  NumLib::getIndices(mesh_item_id, *dof_tables[pcs_id]);
137  assert(!indices.empty());
138  auto const local_solution = x[pcs_id]->get(indices);
139  local_x_vec.insert(std::end(local_x_vec),
140  std::begin(local_solution),
141  std::end(local_solution));
142  }
143  auto const local_x = MathLib::toVector(local_x_vec);
144 
145  auto const indices =
146  NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
147  auto const num_r_c = indices.size();
148 
149  std::vector<double> local_M_data;
150  local_M_data.reserve(num_r_c * num_r_c);
151  std::vector<double> local_K_data;
152  local_K_data.reserve(num_r_c * num_r_c);
153  std::vector<double> local_b_data;
154  local_b_data.reserve(num_r_c);
155 
156  assembleReactionEquationConcrete(t, dt, local_x, local_M_data,
157  local_K_data, local_b_data,
158  process_id);
159 
160  auto const r_c_indices =
162  if (!local_M_data.empty())
163  {
164  auto const local_M =
165  MathLib::toMatrix(local_M_data, num_r_c, num_r_c);
166  M.add(r_c_indices, local_M);
167  }
168  if (!local_K_data.empty())
169  {
170  auto const local_K =
171  MathLib::toMatrix(local_K_data, num_r_c, num_r_c);
172  K.add(r_c_indices, local_K);
173  }
174  if (!local_b_data.empty())
175  {
176  b.add(indices, local_b_data);
177  }
178  }
int add(IndexType row, IndexType col, double val)
Definition: EigenMatrix.h:87
void add(IndexType rowId, double v)
add entry
Definition: EigenVector.h:80
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)
Definition: EigenMapTools.h:72
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

◆ getInterpolatedLocalSolution()

virtual std::vector<double> const& ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::getInterpolatedLocalSolution ( const double  ,
std::vector< GlobalVector * > const &  int_pt_x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const &  ,
std::vector< double > &  cache 
) const
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

◆ 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 78 of file ComponentTransportFEM.h.

82  {
83  std::vector<double> local_x_vec;
84 
85  auto const n_processes = x.size();
86  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
87  {
88  auto const indices =
89  NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
90  assert(!indices.empty());
91  auto const local_solution = x[process_id]->get(indices);
92  local_x_vec.insert(std::end(local_x_vec),
93  std::begin(local_solution),
94  std::end(local_solution));
95  }
96  auto const local_x = MathLib::toVector(local_x_vec);
97 
99  }
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 101 of file ComponentTransportFEM.h.

105  {
106  std::vector<double> local_x_vec;
107 
108  auto const n_processes = x.size();
109  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
110  {
111  auto const indices =
112  NumLib::getIndices(mesh_item_id, *dof_tables[process_id]);
113  assert(!indices.empty());
114  auto const local_solution = x[process_id]->get(indices);
115  local_x_vec.insert(std::end(local_x_vec),
116  std::begin(local_solution),
117  std::end(local_solution));
118  }
119  auto const local_x = MathLib::toVector(local_x_vec);
120 
121  setChemicalSystemConcrete(local_x, t, dt);
122  }
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

◆ setStaggeredCoupledSolutions()

void ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::setStaggeredCoupledSolutions ( std::size_t const  ,
CoupledSolutionsForStaggeredScheme *const  coupling_term 
)
inline

Definition at line 69 of file ComponentTransportFEM.h.

72  {
73  _coupled_solutions = coupling_term;
74  }

References _coupled_solutions.

Member Data Documentation

◆ _coupled_solutions

CoupledSolutionsForStaggeredScheme* ProcessLib::ComponentTransport::ComponentTransportLocalAssemblerInterface::_coupled_solutions {nullptr}
protected

Definition at line 209 of file ComponentTransportFEM.h.

Referenced by setStaggeredCoupledSolutions().


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