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
163 [[unlikely]] if (BaseLib::MPI::anyOf(exception != nullptr))
164 {
165 if (exception) // Only the rank with the exception rethrows...
166 {
168 }
169 // TODO should other ranks rather throw a new exception?
170 // ... but all ranks quit.
171 return;
172 }
173
175 {
176 // Note: computeResiduum() uses bwd/fwd Euler scheme.
178 dt, *x[process_id], *x_prev[process_id], M, K, b);
179
182 bulk_mesh_assembly_data_->residuum_vectors[process_id]);
183 }
184 }
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
GlobalVector computeResiduum(double const dt, GlobalVector const &x, GlobalVector const &x_prev, GlobalMatrix const &M, GlobalMatrix const &K, GlobalVector const &b)

References BaseLib::MPI::anyOf(), 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 327 of file AssemblyMixin.h.

333 {
334 if (!bulk_mesh_matrix_cache_.hasMKb())
335 {
339 }
340
341 DBUG("Reusing saved global M, K, b.");
342
344 time_restore.start();
345
346 auto const [M_, K_, b_] = bulk_mesh_matrix_cache_.MKb();
350
351 INFO("[time] Restoring global M, K, b took {:g} s",
352 time_restore.elapsed());
353
354 return nullptr;
355 }
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 281 of file AssemblyMixin.h.

289 {
291
292 try
293 {
294 auto const& loc_asms = derived().local_assemblers_;
295
296 pvma_.assemble(
297 loc_asms,
298 assembly_data.activeElementIDsSorted(sorted_element_subset)
299 .get(),
300 dof_tables, t, dt, x, x_prev, process_id, M, K, b);
301 }
302 catch (NumLib::AssemblyException const&)
303 {
305 }
306
310
311 if (is_linear_)
312 {
313 DBUG("Saving global M, K, b for later reuse.");
314
316 time_save.start();
317
318 mat_cache.storeMKb(M, K, b);
319
320 INFO("[time] Saving global M, K, b took {:g} s",
321 time_save.elapsed());
322 }
323
324 return exception;
325 }
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 357 of file AssemblyMixin.h.

364 {
366
367 // Temporary matrices for assembly (reused across submeshes).
368 // They will remain nullptr if no assembly is necessary (linear case &
369 // cache populated).
373
374 auto getMKb = [&M_tmp, &K_tmp, &b_tmp,
376 {
377 if (cache.hasMKb())
378 {
379 return cache.MKb();
380 }
381
382 if (!(M_tmp && K_tmp && b_tmp))
383 {
385 mat_spec);
387 mat_spec);
389 mat_spec);
390 }
391
392 return std::tie(*M_tmp, *K_tmp, *b_tmp);
393 };
394
396 for (auto const& [sad, mat_cache] :
399 {
401
402 if (!mat_cache.hasMKb())
403 {
404 M_submesh.setZero();
405 K_submesh.setZero();
406 b_submesh.setZero();
407
411 }
412
413 // TODO: This operates on the entire global vector b once
414 // for each submesh. If the number of submeshes is high,
415 // that might not perform well.
419
421 {
422 // Note: computeResiduum() uses bwd/fwd Euler scheme.
426
429 sad);
430 }
431 }
432
433 return exception;
434 }
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 189 of file AssemblyMixin.h.

195 {
196 DBUG(
197 "AssemblyMixin assembleWithJacobian(t={:g}, dt={:g}, "
198 "process_id={}).",
199 t, dt, process_id);
200
202
203 // Note: We do not cache the Jacobian or b vector. It does not make
204 // sense, because even in the linear case the b vector changes if
205 // the solution changes (unless steady state is reached).
206 // Todo: If we separated Jacobian and residuum assembly, we could
207 // cache at least the Jacobian and only update the residuum.
208 if (is_linear_)
209 {
210 DBUG(
211 "You specified that this process is linear. The assembly for "
212 "the Newton-Raphson method will be run anyways.");
213 }
214
216 derived().getDOFTables(x.size());
217
226
227 [[unlikely]] if (BaseLib::MPI::anyOf(exception != nullptr))
228 {
229 if (exception) // Only the rank with the exception rethrows...
230 {
232 }
233 // ... but all ranks quit.
234 return;
235 }
236
238
239 // TODO better for submesh case, too?
240 // TODO check copy_residua_to_mesh, too. But that needs a refactoring of
241 // postTimestep(), computeSecondaryVariable() etc. hooks.
242 if (/*copy_residua_to_mesh &&*/ bulk_mesh_assembly_data_)
243 {
246 bulk_mesh_assembly_data_->residuum_vectors[process_id]);
247 }
248 }
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::anyOf(), 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 439 of file AssemblyMixin.h.

446 {
447 auto const& loc_asms = derived().local_assemblers_;
449
450 try
451 {
452 pvma_.assembleWithJacobian(
453 loc_asms,
454 assembly_data.activeElementIDsSorted(sorted_element_subset)
455 .get(),
457 }
458 catch (NumLib::AssemblyException const&)
459 {
461 }
462
464
465 return exception;
466 }

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

475 {
477
479 auto b_submesh =
481
482 for (auto const& sad : submesh_assembly_data_)
483 {
484 b_submesh->setZero();
485
489
491
492 // TODO enabling this check needs a refactoring of postTimestep(),
493 // computeSecondaryVariable() etc. hooks.
494 if (/*copy_residua_to_mesh*/ true)
495 {
498 }
499 }
500
501 return exception;
502 }

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

276 {
277 return static_cast<Process const&>(*this);
278 }

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

255 {
257 {
258 WARN("AssemblyMixin. Skipping preOutput() step.");
259 return;
260 }
261
262 switch (*last_assembly_was_)
263 {
266 return;
269 return;
270 }
271 }
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 527 of file AssemblyMixin.h.

532 {
533 // TODO enabling this needs a refactoring of postTimestep(),
534 // computeSecondaryVariable() etc. hooks.
535#if 0
536 auto const matrix_specification =
538
543
544 Jac->setZero();
545 b->setZero();
546
547 assembleWithJacobian(t, dt, x, x_prev, process_id, *b, *Jac, nullptr,
548 true);
549#endif
550 }
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 504 of file AssemblyMixin.h.

509 {
510 auto const matrix_specification =
512
519
520 M->setZero();
521 K->setZero();
522 b->setZero();
523
524 assemble(t, dt, x, x_prev, process_id, *M, *K, *b, nullptr, true);
525 }
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: