OGS
ProcessLib::AssemblyMixin< Process > Class Template Reference

Detailed Description

template<typename Process>
class ProcessLib::AssemblyMixin< Process >

A mixin providing assembly functionality to a specific Process.

The process must be derived from this class (CRTP).

Definition at line 99 of file AssemblyMixin.h.

#include <AssemblyMixin.h>

Inheritance diagram for ProcessLib::AssemblyMixin< Process >:
[legend]
Collaboration diagram for ProcessLib::AssemblyMixin< Process >:
[legend]

Public Member Functions

void initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names)
 
void updateActiveElements ()
 
void assemble (double const t, double const dt, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const process_id, GlobalMatrix &, GlobalMatrix &, GlobalVector &)
 
void assembleWithJacobian (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector &b, GlobalMatrix &Jac)
 

Private Member Functions

Processderived ()
 
Process const & derived () const
 
 AssemblyMixinBase (AbstractJacobianAssembler &jacobian_assembler)
 
- Private Member Functions inherited from ProcessLib::AssemblyMixinBase
 AssemblyMixinBase (AbstractJacobianAssembler &jacobian_assembler)
 
void initializeAssemblyOnSubmeshes (MeshLib::Mesh &bulk_mesh, std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const &pvs)
 
void updateActiveElements (ProcessLib::Process const &process)
 

Private Attributes

friend Process
 
- Private Attributes inherited from ProcessLib::AssemblyMixinBase
std::vector< SubmeshAssemblyDatasubmesh_assembly_data_
 
std::vector< std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > > residuum_vectors_bulk_
 
std::size_t b_submesh_id_ = 0
 ID of the b vector on submeshes, cf. NumLib::VectorProvider.
 
Assembly::ParallelVectorMatrixAssembler pvma_
 

Additional Inherited Members

- Static Private Member Functions inherited from ProcessLib::AssemblyMixinBase
static void copyResiduumVectorsToBulkMesh (GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > residuum_vectors)
 
static void copyResiduumVectorsToSubmesh (int const process_id, GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, SubmeshAssemblyData const &sad)
 

Member Function Documentation

◆ assemble()

template<typename Process >
void ProcessLib::AssemblyMixin< Process >::assemble ( double const t,
double const dt,
std::vector< GlobalVector * > const & ,
std::vector< GlobalVector * > const & ,
int const process_id,
GlobalMatrix & ,
GlobalMatrix & ,
GlobalVector &  )
inline

Implementation similar to assembleWithJacobian calling Assembly::ParallelVectorMatrixAssembler::assemble function. Residuum must be correctly computed.

Definition at line 143 of file AssemblyMixin.h.

148 {
149 DBUG("AssemblyMixin assemble(t={}, dt={}, process_id={}).", t, dt,
150 process_id);
151
155 OGS_FATAL("AssemblyMixin for Picard scheme is not yet implemented.");
156 }
#define OGS_FATAL(...)
Definition Error.h:26
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30

References DBUG(), and OGS_FATAL.

◆ assembleWithJacobian()

template<typename Process >
void ProcessLib::AssemblyMixin< Process >::assembleWithJacobian ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalVector & b,
GlobalMatrix & Jac )
inline

Definition at line 158 of file AssemblyMixin.h.

163 {
164 DBUG("AssemblyMixin assembleWithJacobian(t={}, dt={}, process_id={}).",
165 t, dt, process_id);
166
167 std::vector<NumLib::LocalToGlobalIndexMap const*> const dof_tables =
168 derived().getDOFTables(x.size());
169
170 auto const& loc_asms = derived().local_assemblers_;
171
172 if (!submesh_assembly_data_.empty())
173 {
175 b, b_submesh_id_);
176
177 for (auto const& sad : submesh_assembly_data_)
178 {
179 b_submesh.setZero();
180
181 pvma_.assembleWithJacobian(loc_asms, sad.active_element_ids,
182 dof_tables, t, dt, x, x_prev,
183 process_id, b_submesh, Jac);
184
185 MathLib::LinAlg::axpy(b, 1.0, b_submesh);
186
188 process_id, b_submesh, *(dof_tables[process_id]), sad);
189 }
190
192 }
193 else
194 {
196 loc_asms, derived().getActiveElementIDs(), dof_tables, t, dt, x,
197 x_prev, process_id, b, Jac);
198 }
199
201 b, *(dof_tables[process_id]), residuum_vectors_bulk_[process_id]);
202 }
virtual GlobalVector & getVector(std::size_t &id)=0
Get an uninitialized vector with the given id.
virtual void releaseVector(GlobalVector const &x)=0
Assembly::ParallelVectorMatrixAssembler pvma_
static void copyResiduumVectorsToBulkMesh(GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > residuum_vectors)
std::size_t b_submesh_id_
ID of the b vector on submeshes, cf. NumLib::VectorProvider.
static void copyResiduumVectorsToSubmesh(int const process_id, GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, SubmeshAssemblyData const &sad)
std::vector< SubmeshAssemblyData > submesh_assembly_data_
std::vector< std::vector< std::reference_wrapper< MeshLib::PropertyVector< double > > > > residuum_vectors_bulk_
void assembleWithJacobian(BaseLib::PolymorphicRandomAccessContainerView< LocalAssemblerInterface > const &local_assemblers, std::vector< std::size_t > const &active_elements, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &xs, std::vector< GlobalVector * > const &x_prevs, int const process_id, GlobalVector &b, GlobalMatrix &Jac)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:385
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:57
static NUMLIB_EXPORT VectorProvider & provider

References ProcessLib::Assembly::ParallelVectorMatrixAssembler::assembleWithJacobian(), MathLib::LinAlg::axpy(), ProcessLib::AssemblyMixinBase::b_submesh_id_, ProcessLib::AssemblyMixinBase::copyResiduumVectorsToBulkMesh(), ProcessLib::AssemblyMixinBase::copyResiduumVectorsToSubmesh(), DBUG(), ProcessLib::AssemblyMixin< Process >::derived(), ProcessLib::Process::getDOFTables(), NumLib::VectorProvider::getVector(), NumLib::GlobalVectorProvider::provider, ProcessLib::AssemblyMixinBase::pvma_, NumLib::VectorProvider::releaseVector(), ProcessLib::AssemblyMixinBase::residuum_vectors_bulk_, and ProcessLib::AssemblyMixinBase::submesh_assembly_data_.

◆ AssemblyMixinBase()

template<typename Process >
ProcessLib::AssemblyMixinBase::AssemblyMixinBase ( AbstractJacobianAssembler & jacobian_assembler)
inlineexplicitprivate

Definition at line 51 of file AssemblyMixin.h.

52 : pvma_{jacobian_assembler} {};

◆ derived() [1/2]

◆ derived() [2/2]

template<typename Process >
Process const & ProcessLib::AssemblyMixin< Process >::derived ( ) const
inlineprivate

Definition at line 206 of file AssemblyMixin.h.

207 {
208 return static_cast<Process const&>(*this);
209 }

◆ initializeAssemblyOnSubmeshes()

template<typename Process >
void ProcessLib::AssemblyMixin< Process >::initializeAssemblyOnSubmeshes ( std::vector< std::reference_wrapper< MeshLib::Mesh > > const & submeshes,
std::vector< std::vector< std::string > > const & residuum_names )
inline

Specifies that the assembly of the process should take place on the passed submeshes.

Note
This method must be called at most once.
Attention
The passed submeshes must be non-overlapping and must cover the entire simulation domain. The caller is responsible to meet this requirement.
The order of the passed residuum_names matters and must match the order of the residuum vectors in the process. Otherwise residuum output on submeshes will be wrong.
Note
As of January 2023 assembly on submeshes has neither been tested with domain decomposition, nor with domain deactivation nor staggered coupling. Use at your own risk in these cases.
Todo
Implement and test with domain decomposition, domain deactivation and staggered coupling.

Definition at line 128 of file AssemblyMixin.h.

131 {
133 derived().getMesh(), submeshes, residuum_names,
134 derived().getProcessVariables());
135 }
void initializeAssemblyOnSubmeshes(MeshLib::Mesh &bulk_mesh, std::vector< std::reference_wrapper< MeshLib::Mesh > > const &submeshes, std::vector< std::vector< std::string > > const &residuum_names, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > const &pvs)
OGSMesh getMesh(std::string const &name)

References ProcessLib::AssemblyMixin< Process >::derived(), getMesh(), and ProcessLib::AssemblyMixinBase::initializeAssemblyOnSubmeshes().

◆ updateActiveElements()

template<typename Process >
void ProcessLib::AssemblyMixin< Process >::updateActiveElements ( )
inline

Member Data Documentation

◆ Process

template<typename Process >
friend ProcessLib::AssemblyMixin< Process >::Process
private

Definition at line 104 of file AssemblyMixin.h.


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