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 94 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 &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
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, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)
void preOutput (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)

Private Member Functions

Processderived ()
Process const & derived () const
std::exception_ptr assembleOnBulkMeshOrOnSubmeshCommon (Assembly::CommonAssemblyData const &assembly_data, AssembledMatrixCache &mat_cache, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset)
 Common code for Picard assembly on the bulk mesh or on submeshes.
std::exception_ptr assembleOnBulkMesh (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset)
std::exception_ptr assembleOnSubmeshes (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset, bool const copy_residua_to_mesh)
std::exception_ptr assembleWithJacobianOnBulkMeshOrOnSubmeshCommon (Assembly::CommonAssemblyData const &assembly_data, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset)
std::exception_ptr assembleWithJacobianOnSubmeshes (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset, bool const)
void preOutputPicard (double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
void preOutputNewton (double const, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)
 AssemblyMixinBase (AbstractJacobianAssembler &jacobian_assembler, bool const is_linear, bool const use_monolithic_scheme)
Private Member Functions inherited from ProcessLib::AssemblyMixinBase
 AssemblyMixinBase (AbstractJacobianAssembler &jacobian_assembler, bool const is_linear, bool const use_monolithic_scheme)
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)
bool isLinear () const

Private Attributes

friend Process
Private Attributes inherited from ProcessLib::AssemblyMixinBase
std::vector< ProcessLib::Assembly::SubmeshAssemblyDatasubmesh_assembly_data_
 SubmeshAssemblyData for each submesh.
std::vector< AssembledMatrixCachesubmesh_matrix_cache_
 AssembledMatrixCache for each submesh.
std::optional< ProcessLib::Assembly::BulkMeshAssemblyDatabulk_mesh_assembly_data_
 Empty if submesh assembly is used.
AssembledMatrixCache bulk_mesh_matrix_cache_
Assembly::ParallelVectorMatrixAssembler pvma_
bool is_linear_
std::optional< NumLib::NonlinearSolverTaglast_assembly_was_

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, ProcessLib::Assembly::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 & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
std::vector< std::size_t > const *const sorted_element_subset = nullptr,
bool const copy_residua_to_mesh = false )
inline

Definition at line 140 of file AssemblyMixin.h.

146 {
147 DBUG("AssemblyMixin assemble(t={:g}, dt={:g}, process_id={}).", t, dt,
148 process_id);
149
151
153 derived().getDOFTables(x.size());
154
162
164
166 {
167 // Note: computeResiduum() uses bwd/fwd Euler scheme.
169 dt, *x[process_id], *x_prev[process_id], M, K, b);
170
173 bulk_mesh_assembly_data_->residuum_vectors[process_id]);
174 }
175 }
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:22
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::optional< NumLib::NonlinearSolverTag > last_assembly_was_
std::vector< ProcessLib::Assembly::SubmeshAssemblyData > submesh_assembly_data_
SubmeshAssemblyData for each submesh.
std::optional< ProcessLib::Assembly::BulkMeshAssemblyData > bulk_mesh_assembly_data_
Empty if submesh assembly is used.
std::exception_ptr assembleOnSubmeshes(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset, bool const copy_residua_to_mesh)
std::exception_ptr assembleOnBulkMesh(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:376
void allRanksThrowOrNone(std::exception_ptr const &exception, auto &&warning_callback)
Definition MPI.h:236
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)

References BaseLib::MPI::allRanksThrowOrNone(), assembleOnBulkMesh(), assembleOnSubmeshes(), ProcessLib::AssemblyMixinBase::bulk_mesh_assembly_data_, ProcessLib::computeResiduum(), ProcessLib::AssemblyMixinBase::copyResiduumVectorsToBulkMesh(), DBUG(), derived(), ProcessLib::Process::getDOFTables(), ProcessLib::AssemblyMixinBase::last_assembly_was_, NumLib::Picard, and ProcessLib::AssemblyMixinBase::submesh_assembly_data_.

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::assembleConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleConcreteProcess(), and preOutputPicard().

◆ assembleOnBulkMesh()

template<typename Process>
std::exception_ptr ProcessLib::AssemblyMixin< Process >::assembleOnBulkMesh ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_tables,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
std::vector< std::size_t > const *const sorted_element_subset )
inlineprivate

Definition at line 310 of file AssemblyMixin.h.

316 {
317 if (!bulk_mesh_matrix_cache_.hasMKb())
318 {
322 }
323
324 DBUG("Reusing saved global M, K, b.");
325
327 time_restore.start();
328
329 auto const [M_, K_, b_] = bulk_mesh_matrix_cache_.MKb();
333
334 INFO("[time] Restoring global M, K, b took {:g} s",
335 time_restore.elapsed());
336
337 return nullptr;
338 }
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
AssembledMatrixCache bulk_mesh_matrix_cache_
std::exception_ptr assembleOnBulkMeshOrOnSubmeshCommon(Assembly::CommonAssemblyData const &assembly_data, AssembledMatrixCache &mat_cache, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset)
Common code for Picard assembly on the bulk mesh or on submeshes.
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30

References assembleOnBulkMeshOrOnSubmeshCommon(), ProcessLib::AssemblyMixinBase::bulk_mesh_assembly_data_, ProcessLib::AssemblyMixinBase::bulk_mesh_matrix_cache_, MathLib::LinAlg::copy(), DBUG(), BaseLib::RunTime::elapsed(), INFO(), and BaseLib::RunTime::start().

Referenced by assemble().

◆ assembleOnBulkMeshOrOnSubmeshCommon()

template<typename Process>
std::exception_ptr ProcessLib::AssemblyMixin< Process >::assembleOnBulkMeshOrOnSubmeshCommon ( Assembly::CommonAssemblyData const & assembly_data,
AssembledMatrixCache & mat_cache,
double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_tables,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
std::vector< std::size_t > const *const sorted_element_subset )
inlinenodiscardprivate

Common code for Picard assembly on the bulk mesh or on submeshes.

Definition at line 264 of file AssemblyMixin.h.

272 {
274
275 try
276 {
277 auto const& loc_asms = derived().local_assemblers_;
278
279 pvma_.assemble(
280 loc_asms,
281 assembly_data.activeElementIDsSorted(sorted_element_subset)
282 .get(),
283 dof_tables, t, dt, x, x_prev, process_id, M, K, b);
284 }
285 catch (NumLib::AssemblyException const&)
286 {
288 }
289
293
294 if (is_linear_)
295 {
296 DBUG("Saving global M, K, b for later reuse.");
297
299 time_save.start();
300
301 mat_cache.storeMKb(M, K, b);
302
303 INFO("[time] Saving global M, K, b took {:g} s",
304 time_save.elapsed());
305 }
306
307 return exception;
308 }
Assembly::ParallelVectorMatrixAssembler pvma_
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:191

References ProcessLib::Assembly::CommonAssemblyData::activeElementIDsSorted(), DBUG(), derived(), BaseLib::RunTime::elapsed(), MathLib::LinAlg::finalizeAssembly(), INFO(), ProcessLib::AssemblyMixinBase::is_linear_, ProcessLib::AssemblyMixinBase::pvma_, BaseLib::RunTime::start(), and ProcessLib::AssembledMatrixCache::storeMKb().

Referenced by assembleOnBulkMesh(), and assembleOnSubmeshes().

◆ assembleOnSubmeshes()

template<typename Process>
std::exception_ptr ProcessLib::AssemblyMixin< Process >::assembleOnSubmeshes ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_tables,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
std::vector< std::size_t > const *const sorted_element_subset,
bool const copy_residua_to_mesh )
inlinenodiscardprivate

Definition at line 340 of file AssemblyMixin.h.

347 {
349
350 // Temporary matrices for assembly (reused across submeshes).
351 // They will remain nullptr if no assembly is necessary (linear case &
352 // cache populated).
356
357 auto getMKb = [&M_tmp, &K_tmp, &b_tmp,
359 {
360 if (cache.hasMKb())
361 {
362 return cache.MKb();
363 }
364
365 if (!(M_tmp && K_tmp && b_tmp))
366 {
368 mat_spec);
370 mat_spec);
372 mat_spec);
373 }
374
375 return std::tie(*M_tmp, *K_tmp, *b_tmp);
376 };
377
379 for (auto const& [sad, mat_cache] :
382 {
384
385 if (!mat_cache.hasMKb())
386 {
387 M_submesh.setZero();
388 K_submesh.setZero();
389 b_submesh.setZero();
390
394 }
395
396 // TODO: This operates on the entire global vector b once
397 // for each submesh. If the number of submeshes is high,
398 // that might not perform well.
402
404 {
405 // Note: computeResiduum() uses bwd/fwd Euler scheme.
409
412 sad);
413 }
414 }
415
416 return exception;
417 }
std::vector< AssembledMatrixCache > submesh_matrix_cache_
AssembledMatrixCache for each submesh.
static void copyResiduumVectorsToSubmesh(int const process_id, GlobalVector const &rhs, NumLib::LocalToGlobalIndexMap const &local_to_global_index_map, ProcessLib::Assembly::SubmeshAssemblyData const &sad)
MathLib::MatrixSpecifications getMatrixSpecifications(const int process_id) const override
Definition Process.cpp:203
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50

References assembleOnBulkMeshOrOnSubmeshCommon(), MathLib::LinAlg::axpy(), ProcessLib::computeResiduum(), ProcessLib::AssemblyMixinBase::copyResiduumVectorsToSubmesh(), derived(), ProcessLib::Process::getMatrixSpecifications(), ProcessLib::AssemblyMixinBase::submesh_assembly_data_, and ProcessLib::AssemblyMixinBase::submesh_matrix_cache_.

Referenced by assemble().

◆ 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,
std::vector< std::size_t > const *const sorted_element_subset = nullptr,
bool const copy_residua_to_mesh = false )
inline

Definition at line 180 of file AssemblyMixin.h.

186 {
187 DBUG(
188 "AssemblyMixin assembleWithJacobian(t={:g}, dt={:g}, "
189 "process_id={}).",
190 t, dt, process_id);
191
193
194 // Note: We do not cache the Jacobian or b vector. It does not make
195 // sense, because even in the linear case the b vector changes if
196 // the solution changes (unless steady state is reached).
197 // Todo: If we separated Jacobian and residuum assembly, we could
198 // cache at least the Jacobian and only update the residuum.
199 if (is_linear_)
200 {
201 DBUG(
202 "You specified that this process is linear. The assembly for "
203 "the Newton-Raphson method will be run anyways.");
204 }
205
207 derived().getDOFTables(x.size());
208
217
219
221
222 // TODO better for submesh case, too?
223 // TODO check copy_residua_to_mesh, too. But that needs a refactoring of
224 // postTimestep(), computeSecondaryVariable() etc. hooks.
225 if (/*copy_residua_to_mesh &&*/ bulk_mesh_assembly_data_)
226 {
229 bulk_mesh_assembly_data_->residuum_vectors[process_id]);
230 }
231 }
std::exception_ptr assembleWithJacobianOnBulkMeshOrOnSubmeshCommon(Assembly::CommonAssemblyData const &assembly_data, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset)
std::exception_ptr assembleWithJacobianOnSubmeshes(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, int const process_id, GlobalVector &b, GlobalMatrix &Jac, std::vector< std::size_t > const *const sorted_element_subset, bool const)

References BaseLib::MPI::allRanksThrowOrNone(), assembleWithJacobianOnBulkMeshOrOnSubmeshCommon(), assembleWithJacobianOnSubmeshes(), ProcessLib::AssemblyMixinBase::bulk_mesh_assembly_data_, ProcessLib::AssemblyMixinBase::copyResiduumVectorsToBulkMesh(), DBUG(), derived(), MathLib::LinAlg::finalizeAssembly(), ProcessLib::Process::getDOFTables(), ProcessLib::AssemblyMixinBase::is_linear_, ProcessLib::AssemblyMixinBase::last_assembly_was_, NumLib::Newton, and ProcessLib::AssemblyMixinBase::submesh_assembly_data_.

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatConduction::HeatConductionProcess::assembleWithJacobianConcreteProcess(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::assembleWithJacobianConcreteProcess(), ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::assembleWithJacobianConcreteProcess(), and preOutputNewton().

◆ assembleWithJacobianOnBulkMeshOrOnSubmeshCommon()

template<typename Process>
std::exception_ptr ProcessLib::AssemblyMixin< Process >::assembleWithJacobianOnBulkMeshOrOnSubmeshCommon ( Assembly::CommonAssemblyData const & assembly_data,
double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_tables,
int const process_id,
GlobalVector & b,
GlobalMatrix & Jac,
std::vector< std::size_t > const *const sorted_element_subset )
inlinenodiscardprivate

Common code for the Newton-Raphson assembly on the bulk mesh or on submeshes.

Definition at line 422 of file AssemblyMixin.h.

429 {
430 auto const& loc_asms = derived().local_assemblers_;
432
433 try
434 {
435 pvma_.assembleWithJacobian(
436 loc_asms,
437 assembly_data.activeElementIDsSorted(sorted_element_subset)
438 .get(),
440 }
441 catch (NumLib::AssemblyException const&)
442 {
444 }
445
447
448 return exception;
449 }

References ProcessLib::Assembly::CommonAssemblyData::activeElementIDsSorted(), derived(), MathLib::LinAlg::finalizeAssembly(), and ProcessLib::AssemblyMixinBase::pvma_.

Referenced by assembleWithJacobian(), and assembleWithJacobianOnSubmeshes().

◆ assembleWithJacobianOnSubmeshes()

template<typename Process>
std::exception_ptr ProcessLib::AssemblyMixin< Process >::assembleWithJacobianOnSubmeshes ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_tables,
int const process_id,
GlobalVector & b,
GlobalMatrix & Jac,
std::vector< std::size_t > const *const sorted_element_subset,
bool const  )
inlinenodiscardprivate

Definition at line 451 of file AssemblyMixin.h.

458 {
460
462 auto b_submesh =
464
465 for (auto const& sad : submesh_assembly_data_)
466 {
467 b_submesh->setZero();
468
472
474
475 // TODO enabling this check needs a refactoring of postTimestep(),
476 // computeSecondaryVariable() etc. hooks.
477 if (/*copy_residua_to_mesh*/ true)
478 {
481 }
482 }
483
484 return exception;
485 }

References assembleWithJacobianOnBulkMeshOrOnSubmeshCommon(), MathLib::LinAlg::axpy(), ProcessLib::AssemblyMixinBase::copyResiduumVectorsToSubmesh(), derived(), ProcessLib::Process::getMatrixSpecifications(), and ProcessLib::AssemblyMixinBase::submesh_assembly_data_.

Referenced by assembleWithJacobian().

◆ AssemblyMixinBase()

template<typename Process>
ProcessLib::AssemblyMixinBase::AssemblyMixinBase ( AbstractJacobianAssembler & jacobian_assembler,
bool const is_linear,
bool const use_monolithic_scheme )
explicitprivate

Definition at line 32 of file AssemblyMixin.cpp.

102{
104 {
105 OGS_FATAL(
106 "You requested to assemble only once in combination with staggered "
107 "coupling. This use case is not yet implemented.");
108 }
109
110 if (is_linear_)
111 {
112 WARN(
113 "You specified that the process simulated by OGS is linear. With "
114 "this optimization, the system is assembled only once, and the "
115 "non-linear solver performs a single iteration per time step. No "
116 "non-linearities will be resolved, and OGS will not detect whether "
117 "any exist. It is your responsibility to ensure that the assembled "
118 "systems are indeed linear; there is no safety net.\n"
119 "\n"
120 "Be aware that non-linearities may arise in subtle ways. For "
121 "example, in an HT problem with a steady-state velocity field, if "
122 "the initial pressure field does not match the steady-state "
123 "condition (e.g., initialized to zero), then the steady-state "
124 "velocity field only emerges during the first iteration—making the "
125 "system effectively non-linear under this optimization.");
126 }
127}
#define OGS_FATAL(...)
Definition Error.h:19
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34

References OGS_FATAL.

◆ derived() [1/2]

◆ derived() [2/2]

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

Definition at line 258 of file AssemblyMixin.h.

259 {
260 return static_cast<Process const&>(*this);
261 }

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

126 {
130 }
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().

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::initializeAssemblyOnSubmeshes(), ProcessLib::HeatConduction::HeatConductionProcess::initializeAssemblyOnSubmeshes(), ProcessLib::HeatTransportBHE::HeatTransportBHEProcess::initializeAssemblyOnSubmeshes(), and ProcessLib::ThermoRichardsFlow::ThermoRichardsFlowProcess::initializeAssemblyOnSubmeshes().

◆ preOutput()

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

Definition at line 233 of file AssemblyMixin.h.

238 {
240 {
241 WARN("AssemblyMixin. Skipping preOutput() step.");
242 return;
243 }
244
245 switch (*last_assembly_was_)
246 {
249 return;
252 return;
253 }
254 }
void preOutputPicard(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
void preOutputNewton(double const, double const, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > const &, int const)

References ProcessLib::AssemblyMixinBase::last_assembly_was_, NumLib::Newton, NumLib::Picard, preOutputNewton(), preOutputPicard(), and WARN().

Referenced by ProcessLib::ComponentTransport::ComponentTransportProcess::preOutputConcreteProcess(), and ProcessLib::HeatConduction::HeatConductionProcess::preOutputConcreteProcess().

◆ preOutputNewton()

template<typename Process>
void ProcessLib::AssemblyMixin< Process >::preOutputNewton ( double const ,
double const ,
std::vector< GlobalVector * > const & ,
std::vector< GlobalVector * > const & ,
int const  )
inlineprivate

Definition at line 510 of file AssemblyMixin.h.

515 {
516 // TODO enabling this needs a refactoring of postTimestep(),
517 // computeSecondaryVariable() etc. hooks.
518#if 0
519 auto const matrix_specification =
521
526
527 Jac->setZero();
528 b->setZero();
529
530 assembleWithJacobian(t, dt, x, x_prev, process_id, *b, *Jac, nullptr,
531 true);
532#endif
533 }
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, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)

References assembleWithJacobian(), derived(), and ProcessLib::Process::getMatrixSpecifications().

Referenced by preOutput().

◆ preOutputPicard()

template<typename Process>
void ProcessLib::AssemblyMixin< Process >::preOutputPicard ( double const t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id )
inlineprivate

Definition at line 487 of file AssemblyMixin.h.

492 {
493 auto const matrix_specification =
495
502
503 M->setZero();
504 K->setZero();
505 b->setZero();
506
507 assemble(t, dt, x, x_prev, process_id, *M, *K, *b, nullptr, true);
508 }
void assemble(double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, std::vector< std::size_t > const *const sorted_element_subset=nullptr, bool const copy_residua_to_mesh=false)

References assemble(), derived(), and ProcessLib::Process::getMatrixSpecifications().

Referenced by preOutput().

◆ updateActiveElements()

Member Data Documentation

◆ Process

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

Definition at line 99 of file AssemblyMixin.h.

Referenced by derived(), and derived().


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