v0.13.2
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Private Attributes | Friends | List of all members
MoFEM::KspCtx Struct Reference

Interface for linear (KSP) solver. More...

#include <src/petsc/KspCtx.hpp>

Collaboration diagram for MoFEM::KspCtx:
[legend]

Public Types

typedef MoFEM::PairNameFEMethodPtr PairNameFEMethodPtr
 
typedef MoFEM::FEMethodsSequence FEMethodsSequence
 
typedef MoFEM::BasicMethodsSequence BasicMethodsSequence
 

Public Member Functions

 KspCtx (MoFEM::Interface &m_field, const std::string &_problem_name)
 
virtual ~KspCtx ()=default
 
FEMethodsSequenceget_loops_to_do_Mat ()
 
FEMethodsSequenceget_loops_to_do_Rhs ()
 
BasicMethodsSequenceget_preProcess_to_do_Rhs ()
 
BasicMethodsSequenceget_postProcess_to_do_Rhs ()
 
BasicMethodsSequenceget_preProcess_to_do_Mat ()
 
BasicMethodsSequenceget_postProcess_to_do_Mat ()
 
MoFEMErrorCode clearLoops ()
 Clear loops. More...
 

Public Attributes

MoFEM::InterfacemField
 
moab::Interface & moab
 
std::string problemName
 Problem name. More...
 
MoFEMTypes bH
 If set to MF_EXIST check if element exist. More...
 
FEMethodsSequence loops_to_do_Mat
 
FEMethodsSequence loops_to_do_Rhs
 
BasicMethodsSequence preProcess_Mat
 
BasicMethodsSequence postProcess_Mat
 
BasicMethodsSequence preProcess_Rhs
 
BasicMethodsSequence postProcess_Rhs
 

Private Attributes

PetscLogEvent MOFEM_EVENT_KspRhs
 
PetscLogEvent MOFEM_EVENT_KspMat
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 

Friends

PetscErrorCode KspRhs (KSP ksp, Vec f, void *ctx)
 Run over elements in the lists. More...
 
PetscErrorCode KspMat (KSP ksp, Mat A, Mat B, void *ctx)
 Run over elenents in the list. More...
 

Detailed Description

Interface for linear (KSP) solver.

Definition at line 15 of file KspCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

Definition at line 25 of file KspCtx.hpp.

◆ FEMethodsSequence

Definition at line 24 of file KspCtx.hpp.

◆ PairNameFEMethodPtr

Definition at line 23 of file KspCtx.hpp.

Constructor & Destructor Documentation

◆ KspCtx()

MoFEM::KspCtx::KspCtx ( MoFEM::Interface m_field,
const std::string &  _problem_name 
)
inline

Definition at line 40 of file KspCtx.hpp.

41 : mField(m_field), moab(m_field.get_moab()), problemName(_problem_name),
42 bH(MF_EXIST) {
43 PetscLogEventRegister("LoopKSPRhs", 0, &MOFEM_EVENT_KspRhs);
44 PetscLogEventRegister("LoopKSPMat", 0, &MOFEM_EVENT_KspMat);
45 }
@ MF_EXIST
Definition: definitions.h:100
virtual moab::Interface & get_moab()=0
moab::Interface & moab
Definition: KspCtx.hpp:18
PetscLogEvent MOFEM_EVENT_KspRhs
Definition: KspCtx.hpp:101
PetscLogEvent MOFEM_EVENT_KspMat
Definition: KspCtx.hpp:102
MoFEM::Interface & mField
Definition: KspCtx.hpp:17
MoFEMTypes bH
If set to MF_EXIST check if element exist.
Definition: KspCtx.hpp:21
std::string problemName
Problem name.
Definition: KspCtx.hpp:20

◆ ~KspCtx()

virtual MoFEM::KspCtx::~KspCtx ( )
virtualdefault

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode MoFEM::KspCtx::clearLoops ( )

Clear loops.

Returns
MoFEMErrorCode

Definition at line 11 of file KspCtx.cpp.

11 {
13 loops_to_do_Mat.clear();
14 loops_to_do_Rhs.clear();
15 preProcess_Mat.clear();
16 postProcess_Mat.clear();
17 preProcess_Rhs.clear();
19}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
FEMethodsSequence loops_to_do_Rhs
Definition: KspCtx.hpp:29
BasicMethodsSequence preProcess_Rhs
Definition: KspCtx.hpp:35
BasicMethodsSequence postProcess_Mat
Definition: KspCtx.hpp:33
BasicMethodsSequence preProcess_Mat
Definition: KspCtx.hpp:31
FEMethodsSequence loops_to_do_Mat
Definition: KspCtx.hpp:27

◆ get_loops_to_do_Mat()

FEMethodsSequence & MoFEM::KspCtx::get_loops_to_do_Mat ( )
inline
Returns
return reference to vector with FEMethod to calculate matrix

Definition at line 51 of file KspCtx.hpp.

51{ return loops_to_do_Mat; }

◆ get_loops_to_do_Rhs()

FEMethodsSequence & MoFEM::KspCtx::get_loops_to_do_Rhs ( )
inline
Returns
return vector to vector with FEMethod to vector

Definition at line 56 of file KspCtx.hpp.

56{ return loops_to_do_Rhs; }

◆ get_postProcess_to_do_Mat()

BasicMethodsSequence & MoFEM::KspCtx::get_postProcess_to_do_Mat ( )
inline

The sequence of BasicMethod is executed after tangent matrix is calculated. It can be used to setup data structures, e.g. aggregate data from processors or to apply essential boundary conditions.

Returns
reference to BasicMethod for postprocessing

Definition at line 88 of file KspCtx.hpp.

88{ return postProcess_Mat; }

◆ get_postProcess_to_do_Rhs()

BasicMethodsSequence & MoFEM::KspCtx::get_postProcess_to_do_Rhs ( )
inline

The sequence of BasicMethod is executed after residual is calculated. It can be used to setup data structures, e.g. aggregate data from processors or to apply essential boundary conditions.

Returns
reference to BasicMethod for postprocessing

Definition at line 74 of file KspCtx.hpp.

74{ return postProcess_Rhs; }
BasicMethodsSequence postProcess_Rhs
Definition: KspCtx.hpp:37

◆ get_preProcess_to_do_Mat()

BasicMethodsSequence & MoFEM::KspCtx::get_preProcess_to_do_Mat ( )
inline
Returns
reference to BasicMethod for preprocessing

Definition at line 79 of file KspCtx.hpp.

79{ return preProcess_Mat; }

◆ get_preProcess_to_do_Rhs()

BasicMethodsSequence & MoFEM::KspCtx::get_preProcess_to_do_Rhs ( )
inline

The sequence of BasicMethod is executed before residual is calculated. It can be used to setup data structures, e.g. zero global variable which is integrated in domain, e.g. for calculation of strain energy.

Returns
reference to BasicMethod for preprocessing

Definition at line 65 of file KspCtx.hpp.

65{ return preProcess_Rhs; }

Friends And Related Function Documentation

◆ KspMat

PetscErrorCode KspMat ( KSP  ksp,
Mat  A,
Mat  B,
void *  ctx 
)
friend

Run over elenents in the list.

Parameters
kspKSP solver
Amatrix
BPreconditioner matrix
ctxdata context, i.e. KspCtx
Returns
error code

Definition at line 85 of file KspCtx.cpp.

85 {
86 // PetscValidHeaderSpecific(ksp,KSP_CLASSID,1);
88 KspCtx *ksp_ctx = static_cast<KspCtx *>(ctx);
89 PetscLogEventBegin(ksp_ctx->MOFEM_EVENT_KspMat, 0, 0, 0, 0);
90
91 ksp_ctx->matAssembleSwitch = boost::movelib::make_unique<bool>(true);
92
93 auto set = [&](auto &fe) {
94 fe.ksp = ksp;
95 fe.ksp_A = A;
96 fe.ksp_B = B;
97 fe.ksp_ctx = KspMethod::CTX_OPERATORS;
99 };
100
101 auto unset = [&](auto &fe) {
102 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
103 fe.data_ctx = PetscData::CtxSetNone;
104 };
105
106 auto ent_data_cache = boost::make_shared<std::vector<EntityCacheDofs>>();
107 auto ent_row_cache =
108 boost::make_shared<std::vector<EntityCacheNumeredDofs>>();
109 auto ent_col_cache =
110 boost::make_shared<std::vector<EntityCacheNumeredDofs>>();
111
112 // pre-procsess
113 for (auto &bit : ksp_ctx->preProcess_Mat) {
114 bit->matAssembleSwitch = boost::move(ksp_ctx->matAssembleSwitch);
115 set(*bit);
116 CHKERR ksp_ctx->mField.problem_basic_method_preProcess(ksp_ctx->problemName,
117 *bit);
118 unset(*bit);
119 ksp_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
120 }
121
122 auto cache_ptr = boost::make_shared<CacheTuple>();
123 CHKERR ksp_ctx->mField.cache_problem_entities(ksp_ctx->problemName,
124 cache_ptr);
125
126 // operators
127 for (auto &lit : ksp_ctx->loops_to_do_Mat) {
128 lit.second->matAssembleSwitch = boost::move(ksp_ctx->matAssembleSwitch);
129 set(*lit.second);
130 CHKERR ksp_ctx->mField.loop_finite_elements(ksp_ctx->problemName, lit.first,
131 *(lit.second), nullptr,
132 ksp_ctx->bH, cache_ptr);
133 unset(*lit.second);
134 ksp_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
135 }
136
137 // post-process
138 for (auto &bit : ksp_ctx->postProcess_Mat) {
139 bit->matAssembleSwitch = boost::move(ksp_ctx->matAssembleSwitch);
140 set(*bit);
141 CHKERR ksp_ctx->mField.problem_basic_method_postProcess(
142 ksp_ctx->problemName, *bit);
143 unset(*bit);
144 ksp_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
145 }
146
147 if (ksp_ctx->matAssembleSwitch) {
148 CHKERR MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
149 CHKERR MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
150 }
151 PetscLogEventEnd(ksp_ctx->MOFEM_EVENT_KspMat, 0, 0, 0, 0);
153}
#define CHKERR
Inline error check.
Definition: definitions.h:535
auto bit
set bit
constexpr AssemblyType A
KspCtx(MoFEM::Interface &m_field, const std::string &_problem_name)
Definition: KspCtx.hpp:40
static constexpr Switches CtxSetA
Definition: LoopMethods.hpp:37
static constexpr Switches CtxSetNone
Definition: LoopMethods.hpp:35
static constexpr Switches CtxSetB
Definition: LoopMethods.hpp:38

◆ KspRhs

PetscErrorCode KspRhs ( KSP  ksp,
Vec  f,
void *  ctx 
)
friend

Run over elements in the lists.

Parameters
kspKSP solver
fthe right hand side vector
ctxdata context, i.e. KspCtx
Returns
error code

Definition at line 21 of file KspCtx.cpp.

21 {
22 // PetscValidHeaderSpecific(ksp,KSP_CLASSID,1);
24 KspCtx *ksp_ctx = static_cast<KspCtx *>(ctx);
25 PetscLogEventBegin(ksp_ctx->MOFEM_EVENT_KspRhs, 0, 0, 0, 0);
26
27 ksp_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
28
29 auto set = [&](auto &fe) {
30 fe.ksp = ksp;
31 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
32 fe.data_ctx = PetscData::CtxSetF;
33 fe.ksp_f = f;
34 };
35
36 auto unset = [&](auto &fe) {
37 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
38 fe.data_ctx = PetscData::CtxSetNone;
39 };
40
41 // pre-process
42 for (auto &bit : ksp_ctx->preProcess_Rhs) {
43 bit->vecAssembleSwitch = boost::move(ksp_ctx->vecAssembleSwitch);
44 set(*bit);
45 CHKERR ksp_ctx->mField.problem_basic_method_preProcess(ksp_ctx->problemName,
46 *bit);
47 unset(*bit);
48 ksp_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
49 }
50
51 auto cache_ptr = boost::make_shared<CacheTuple>();
52 CHKERR ksp_ctx->mField.cache_problem_entities(ksp_ctx->problemName,
53 cache_ptr);
54
55 // operators
56 for (auto &lit : ksp_ctx->loops_to_do_Rhs) {
57 lit.second->vecAssembleSwitch = boost::move(ksp_ctx->vecAssembleSwitch);
58 set(*lit.second);
59 CHKERR ksp_ctx->mField.loop_finite_elements(ksp_ctx->problemName, lit.first,
60 *(lit.second), nullptr,
61 ksp_ctx->bH, cache_ptr);
62 unset(*lit.second);
63 ksp_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
64 }
65
66 // post-process
67 for (auto &bit : ksp_ctx->postProcess_Rhs) {
68 bit->vecAssembleSwitch = boost::move(ksp_ctx->vecAssembleSwitch);
69 set(*bit);
70 CHKERR ksp_ctx->mField.problem_basic_method_postProcess(
71 ksp_ctx->problemName, *bit);
72 unset(*bit);
73 ksp_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
74 }
75
76 if (*ksp_ctx->vecAssembleSwitch) {
77 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
78 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
79 CHKERR VecAssemblyBegin(f);
80 CHKERR VecAssemblyEnd(f);
81 }
82 PetscLogEventEnd(ksp_ctx->MOFEM_EVENT_KspRhs, 0, 0, 0, 0);
84}
static constexpr Switches CtxSetF
Definition: LoopMethods.hpp:36

Member Data Documentation

◆ bH

MoFEMTypes MoFEM::KspCtx::bH

If set to MF_EXIST check if element exist.

Definition at line 21 of file KspCtx.hpp.

◆ loops_to_do_Mat

FEMethodsSequence MoFEM::KspCtx::loops_to_do_Mat

Sequence of finite elements instances assembling tangent matrix

Definition at line 27 of file KspCtx.hpp.

◆ loops_to_do_Rhs

FEMethodsSequence MoFEM::KspCtx::loops_to_do_Rhs

Sequence of finite elements instances assembling residual vector

Definition at line 29 of file KspCtx.hpp.

◆ matAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::KspCtx::matAssembleSwitch
private

Definition at line 105 of file KspCtx.hpp.

◆ mField

MoFEM::Interface& MoFEM::KspCtx::mField

Definition at line 17 of file KspCtx.hpp.

◆ moab

moab::Interface& MoFEM::KspCtx::moab

Definition at line 18 of file KspCtx.hpp.

◆ MOFEM_EVENT_KspMat

PetscLogEvent MoFEM::KspCtx::MOFEM_EVENT_KspMat
private

Definition at line 102 of file KspCtx.hpp.

◆ MOFEM_EVENT_KspRhs

PetscLogEvent MoFEM::KspCtx::MOFEM_EVENT_KspRhs
private

Definition at line 101 of file KspCtx.hpp.

◆ postProcess_Mat

BasicMethodsSequence MoFEM::KspCtx::postProcess_Mat

Sequence of methods run after tangent matrix is assembled

Definition at line 33 of file KspCtx.hpp.

◆ postProcess_Rhs

BasicMethodsSequence MoFEM::KspCtx::postProcess_Rhs

Sequence of methods run after residual is assembled

Definition at line 37 of file KspCtx.hpp.

◆ preProcess_Mat

BasicMethodsSequence MoFEM::KspCtx::preProcess_Mat

Sequence of methods run before tangent matrix is assembled

Definition at line 31 of file KspCtx.hpp.

◆ preProcess_Rhs

BasicMethodsSequence MoFEM::KspCtx::preProcess_Rhs

Sequence of methods run before residual is assembled

Definition at line 35 of file KspCtx.hpp.

◆ problemName

std::string MoFEM::KspCtx::problemName

Problem name.

Definition at line 20 of file KspCtx.hpp.

◆ vecAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::KspCtx::vecAssembleSwitch
private

Definition at line 104 of file KspCtx.hpp.


The documentation for this struct was generated from the following files: