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 101 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 145 of file AssemblyMixin.h.

150 {
151 DBUG("AssemblyMixin assemble(t={}, dt={}, process_id={}).", t, dt,
152 process_id);
153
157 OGS_FATAL("AssemblyMixin for Picard scheme is not yet implemented.");
158 }
#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 160 of file AssemblyMixin.h.

165 {
166 DBUG("AssemblyMixin assembleWithJacobian(t={}, dt={}, process_id={}).",
167 t, dt, process_id);
168
170 derived().getDOFTables(x.size());
171
172 auto const& loc_asms = derived().local_assemblers_;
173
175 if (!submesh_assembly_data_.empty())
176 {
179
180 for (auto const& sad : submesh_assembly_data_)
181 {
182 b_submesh.setZero();
183
184 try
185 {
186 pvma_.assembleWithJacobian(loc_asms, sad.active_element_ids,
187 dof_tables, t, dt, x, x_prev,
189 }
190 catch (NumLib::AssemblyException const&)
191 {
193 }
194
196
199 }
200
202 }
203 else
204 {
205 try
206 {
207 pvma_.assembleWithJacobian(
209 dt, x, x_prev, process_id, b, Jac);
210 }
211 catch (NumLib::AssemblyException const&)
212 {
214 }
215 }
216
219
220 if (BaseLib::MPI::anyOf(exception != nullptr))
221 {
222 if (exception) // Only the rank with the exception rethrows...
223 {
225 }
226 // ... but all ranks quit.
227 return;
228 }
229
232 }
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_
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:383
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:198
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:57

References BaseLib::MPI::anyOf(), MathLib::LinAlg::axpy(), ProcessLib::AssemblyMixinBase::b_submesh_id_, ProcessLib::AssemblyMixinBase::copyResiduumVectorsToBulkMesh(), ProcessLib::AssemblyMixinBase::copyResiduumVectorsToSubmesh(), DBUG(), derived(), MathLib::LinAlg::finalizeAssembly(), ProcessLib::Process::getDOFTables(), NumLib::GlobalVectorProvider::provider, ProcessLib::AssemblyMixinBase::pvma_, 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 53 of file AssemblyMixin.h.

◆ derived() [1/2]

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

Definition at line 235 of file AssemblyMixin.h.

235{ return static_cast<Process&>(*this); }

References Process.

Referenced by assembleWithJacobian(), initializeAssemblyOnSubmeshes(), and updateActiveElements().

◆ derived() [2/2]

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

Definition at line 236 of file AssemblyMixin.h.

237 {
238 return static_cast<Process const&>(*this);
239 }

References Process.

◆ 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 130 of file AssemblyMixin.h.

133 {
137 }
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)

References derived(), and ProcessLib::AssemblyMixinBase::initializeAssemblyOnSubmeshes().

◆ updateActiveElements()

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

Definition at line 139 of file AssemblyMixin.h.

140 {
142 }
void updateActiveElements(ProcessLib::Process const &process)

References derived(), and ProcessLib::AssemblyMixinBase::updateActiveElements().

Member Data Documentation

◆ Process

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

Definition at line 106 of file AssemblyMixin.h.

Referenced by derived(), and derived().


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