v0.14.0
ConstrainMatrixCtx.cpp
Go to the documentation of this file.
1 /** \file ConstrainMatrixCtx.cpp
2  * \brief Implementation of projection matrix
3  *
4  *
5  */
6 
7 #include <ConstrainMatrixCtx.hpp>
8 namespace MoFEM {
9 const static bool debug = false;
10 
11 #define INIT_DATA_CONSTRAINMATRIXCTX \
12  C(PETSC_NULL), CT(PETSC_NULL), CCT(PETSC_NULL), CTC(PETSC_NULL), \
13  K(PETSC_NULL), Cx(PETSC_NULL), CCTm1_Cx(PETSC_NULL), \
14  CT_CCTm1_Cx(PETSC_NULL), CTCx(PETSC_NULL), Qx(PETSC_NULL), \
15  KQx(PETSC_NULL), initQorP(true), initQTKQ(true), createKSP(create_ksp), \
16  createScatter(true), cancelKSPMonitor(true), \
17  ownConstrainMatrix(own_contrain_matrix)
18 
20  string x_problem, string y_problem,
21  bool create_ksp,
22  bool own_contrain_matrix)
23  : mField(m_field), INIT_DATA_CONSTRAINMATRIXCTX, xProblem(x_problem),
24  yProblem(y_problem) {
25  PetscLogEventRegister("ProjectionInit", 0, &MOFEM_EVENT_projInit);
26  PetscLogEventRegister("ProjectionQ", 0, &MOFEM_EVENT_projQ);
27  PetscLogEventRegister("ProjectionP", 0, &MOFEM_EVENT_projP);
28  PetscLogEventRegister("ProjectionR", 0, &MOFEM_EVENT_projR);
29  PetscLogEventRegister("ProjectionRT", 0, &MOFEM_EVENT_projRT);
30  PetscLogEventRegister("ProjectionCTC_QTKQ", 0, &MOFEM_EVENT_projCTC_QTKQ);
31 }
32 
34  VecScatter scatter, bool create_ksp,
35  bool own_contrain_matrix)
36  : mField(m_field), INIT_DATA_CONSTRAINMATRIXCTX, sCatter(scatter) {
37  PetscLogEventRegister("ProjectionInit", 0, &MOFEM_EVENT_projInit);
38  PetscLogEventRegister("ProjectionQ", 0, &MOFEM_EVENT_projQ);
39  PetscLogEventRegister("ProjectionP", 0, &MOFEM_EVENT_projP);
40  PetscLogEventRegister("ProjectionR", 0, &MOFEM_EVENT_projR);
41  PetscLogEventRegister("ProjectionRT", 0, &MOFEM_EVENT_projRT);
42  PetscLogEventRegister("ProjectionCTC_QTKQ", 0, &MOFEM_EVENT_projCTC_QTKQ);
43 }
44 
47  if (initQorP) {
48  initQorP = false;
49 
50  PetscLogEventBegin(MOFEM_EVENT_projInit, 0, 0, 0, 0);
51  CHKERR MatTranspose(C, MAT_INITIAL_MATRIX, &CT);
52  // need to be calculated when C is changed
53  CHKERR MatTransposeMatMult(CT, CT, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &CCT);
54  if (createKSP) {
55  CHKERR KSPCreate(mField.get_comm(), &kSP);
56  // need to be recalculated when C is changed
57  CHKERR KSPSetOperators(kSP, CCT, CCT);
58  CHKERR KSPSetFromOptions(kSP);
59  CHKERR KSPSetInitialGuessKnoll(kSP, PETSC_TRUE);
60  CHKERR KSPGetTolerances(kSP, &rTol, &absTol, &dTol, &maxIts);
61  CHKERR KSPSetUp(kSP);
62  if (cancelKSPMonitor) {
63  CHKERR KSPMonitorCancel(kSP);
64  }
65  }
66 #if PETSC_VERSION_GE(3, 5, 3)
67  CHKERR MatCreateVecs(C, &X, PETSC_NULL);
68  CHKERR MatCreateVecs(C, PETSC_NULL, &Cx);
69  CHKERR MatCreateVecs(CCT, PETSC_NULL, &CCTm1_Cx);
70 #else
71  CHKERR MatGetVecs(C, &X, PETSC_NULL);
72  CHKERR MatGetVecs(C, PETSC_NULL, &Cx);
73  CHKERR MatGetVecs(CCT, PETSC_NULL, &CCTm1_Cx);
74 #endif
75  CHKERR VecDuplicate(X, &CT_CCTm1_Cx);
76  if (createScatter) {
77  CHKERR mField.getInterface<VecManager>()->vecScatterCreate(
78  x, xProblem, ROW, X, yProblem, COL, &sCatter);
79  }
80  PetscLogEventEnd(MOFEM_EVENT_projInit, 0, 0, 0, 0);
81  }
83 }
84 
87  if (initQorP)
89  CHKERR MatTranspose(C, MAT_REUSE_MATRIX, &CT);
90  CHKERR MatTransposeMatMult(CT, CT, MAT_REUSE_MATRIX, PETSC_DEFAULT, &CCT);
92 }
93 
96  if (initQorP)
98  CHKERR MatDestroy(&CT);
99  CHKERR MatDestroy(&CCT);
100  if (createKSP) {
101  CHKERR KSPDestroy(&kSP);
102  }
103  CHKERR VecDestroy(&X);
104  CHKERR VecDestroy(&Cx);
105  CHKERR VecDestroy(&CCTm1_Cx);
106  CHKERR VecDestroy(&CT_CCTm1_Cx);
107  if (createScatter) {
108  CHKERR VecScatterDestroy(&sCatter);
109  }
110  initQorP = true;
112 }
113 
116  if (initQTKQ) {
117  initQTKQ = false;
118  PetscLogEventBegin(MOFEM_EVENT_projInit, 0, 0, 0, 0);
119  // need to be recalculated when C is changed
120  CHKERR MatTransposeMatMult(C, C, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &CTC);
121  if (debug) {
122  // MatView(CCT,PETSC_VIEWER_DRAW_WORLD);
123  int m, n;
124  MatGetSize(CCT, &m, &n);
125  PetscPrintf(mField.get_comm(), "CCT size (%d,%d)\n", m, n);
126  // std::string wait;
127  // std::cin >> wait;
128  }
129 #if PETSC_VERSION_GE(3, 5, 3)
130  CHKERR MatCreateVecs(K, &Qx, PETSC_NULL);
131  CHKERR MatCreateVecs(K, PETSC_NULL, &KQx);
132  CHKERR MatCreateVecs(CTC, PETSC_NULL, &CTCx);
133 #else
134  CHKERR MatGetVecs(K, &Qx, PETSC_NULL);
135  CHKERR MatGetVecs(K, PETSC_NULL, &KQx);
136  CHKERR MatGetVecs(CTC, PETSC_NULL, &CTCx);
137 #endif
138  PetscLogEventEnd(MOFEM_EVENT_projInit, 0, 0, 0, 0);
139  }
141 }
142 
145  if (initQTKQ)
147  CHKERR MatTransposeMatMult(C, C, MAT_REUSE_MATRIX, PETSC_DEFAULT, &CTC);
149 }
150 
153  if (initQTKQ)
155  CHKERR MatDestroy(&CTC);
156  CHKERR VecDestroy(&Qx);
157  CHKERR VecDestroy(&KQx);
158  CHKERR VecDestroy(&CTCx);
159  initQTKQ = true;
161 }
162 
165  void *void_ctx;
166  CHKERR MatShellGetContext(Q, &void_ctx);
167  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
168  PetscLogEventBegin(ctx->MOFEM_EVENT_projQ, 0, 0, 0, 0);
169  CHKERR ctx->initializeQorP(x);
170  CHKERR VecCopy(x, f);
171  CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
172  SCATTER_FORWARD);
173  CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
174  if (debug) {
175  // CHKERR VecView(ctx->X,PETSC_VIEWER_STDOUT_WORLD);
176  CHKERR VecScatterBegin(ctx->sCatter, ctx->X, f, INSERT_VALUES,
177  SCATTER_REVERSE);
178  CHKERR VecScatterEnd(ctx->sCatter, ctx->X, f, INSERT_VALUES,
179  SCATTER_REVERSE);
180  PetscBool flg;
181  CHKERR VecEqual(x, f, &flg);
182  if (flg == PETSC_FALSE)
183  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "scatter is not working");
184  }
185  CHKERR MatMult(ctx->C, ctx->X, ctx->Cx);
186  CHKERR KSPSolve(ctx->kSP, ctx->Cx, ctx->CCTm1_Cx);
187  CHKERR MatMult(ctx->CT, ctx->CCTm1_Cx, ctx->CT_CCTm1_Cx);
188  CHKERR VecScale(ctx->CT_CCTm1_Cx, -1);
189  CHKERR VecScatterBegin(ctx->sCatter, ctx->CT_CCTm1_Cx, f, ADD_VALUES,
190  SCATTER_REVERSE);
191  CHKERR VecScatterEnd(ctx->sCatter, ctx->CT_CCTm1_Cx, f, ADD_VALUES,
192  SCATTER_REVERSE);
193  PetscLogEventEnd(ctx->MOFEM_EVENT_projQ, 0, 0, 0, 0);
195 }
196 
199  void *void_ctx;
200  CHKERR MatShellGetContext(P, &void_ctx);
201  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
202  PetscLogEventBegin(ctx->MOFEM_EVENT_projP, 0, 0, 0, 0);
203  CHKERR ctx->initializeQorP(x);
204  CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
205  SCATTER_FORWARD);
206  CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
207  CHKERR MatMult(ctx->C, ctx->X, ctx->Cx);
208  CHKERR KSPSolve(ctx->kSP, ctx->Cx, ctx->CCTm1_Cx);
209  CHKERR MatMult(ctx->CT, ctx->CCTm1_Cx, ctx->CT_CCTm1_Cx);
210  CHKERR VecZeroEntries(f);
211  CHKERR VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD);
212  CHKERR VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD);
213  CHKERR VecScatterBegin(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
214  SCATTER_REVERSE);
215  CHKERR VecScatterEnd(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
216  SCATTER_REVERSE);
217  PetscLogEventEnd(ctx->MOFEM_EVENT_projP, 0, 0, 0, 0);
219 }
220 
223  void *void_ctx;
224  CHKERR MatShellGetContext(R, &void_ctx);
225  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
226  PetscLogEventBegin(ctx->MOFEM_EVENT_projR, 0, 0, 0, 0);
227  if (ctx->initQorP)
228  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
229  "you have to call first initQorP or use Q matrix");
230  CHKERR KSPSolve(ctx->kSP, x, ctx->CCTm1_Cx);
231  CHKERR MatMult(ctx->CT, ctx->CCTm1_Cx, ctx->CT_CCTm1_Cx);
232  CHKERR VecZeroEntries(f);
233  CHKERR VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD);
234  CHKERR VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD);
235  CHKERR VecScatterBegin(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
236  SCATTER_REVERSE);
237  CHKERR VecScatterEnd(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
238  SCATTER_REVERSE);
239  PetscLogEventEnd(ctx->MOFEM_EVENT_projR, 0, 0, 0, 0);
241 }
242 
245  void *void_ctx;
246  CHKERR MatShellGetContext(RT, &void_ctx);
247  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
248  PetscLogEventBegin(ctx->MOFEM_EVENT_projRT, 0, 0, 0, 0);
249  if (ctx->initQorP)
250  SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
251  "you have to call first initQorP or use Q matrix");
252  CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
253  SCATTER_FORWARD);
254  CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
255  CHKERR MatMult(ctx->C, ctx->X, ctx->Cx);
256  CHKERR KSPSolve(ctx->kSP, ctx->Cx, f);
257  PetscLogEventEnd(ctx->MOFEM_EVENT_projRT, 0, 0, 0, 0);
259 }
260 
263  void *void_ctx;
264  CHKERR MatShellGetContext(CTC_QTKQ, &void_ctx);
265  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
266  PetscLogEventBegin(ctx->MOFEM_EVENT_projCTC_QTKQ, 0, 0, 0, 0);
267  Mat Q;
268  int M, N, m, n;
269  CHKERR MatGetSize(ctx->K, &M, &N);
270  CHKERR MatGetLocalSize(ctx->K, &m, &n);
271  CHKERR MatCreateShell(ctx->mField.get_comm(), m, n, M, N, ctx, &Q);
272  CHKERR MatShellSetOperation(Q, MATOP_MULT,
273  (void (*)(void))ProjectionMatrixMultOpQ);
274  CHKERR ctx->initializeQTKQ();
275  CHKERR MatMult(Q, x, ctx->Qx);
276  CHKERR MatMult(ctx->K, ctx->Qx, ctx->KQx);
277  CHKERR MatMult(Q, ctx->KQx, f);
278  CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
279  SCATTER_FORWARD);
280  CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
281  CHKERR MatMult(ctx->CTC, ctx->X, ctx->CTCx);
282  CHKERR VecScatterBegin(ctx->sCatter, ctx->CTCx, f, ADD_VALUES,
283  SCATTER_REVERSE);
284  CHKERR VecScatterEnd(ctx->sCatter, ctx->CTCx, f, ADD_VALUES, SCATTER_REVERSE);
285  CHKERR MatDestroy(&Q);
286  PetscLogEventEnd(ctx->MOFEM_EVENT_projCTC_QTKQ, 0, 0, 0, 0);
288 }
289 
292  void *void_ctx;
293  CHKERR MatShellGetContext(Q, &void_ctx);
294  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
295  CHKERR ctx->destroyQorP();
297 }
300  void *void_ctx;
301  CHKERR MatShellGetContext(QTKQ, &void_ctx);
302  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
303  CHKERR ctx->destroyQTKQ();
305 }
306 
307 } // namespace MoFEM
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::ConstrainMatrixCtx::yProblem
string yProblem
Definition: ConstrainMatrixCtx.hpp:36
MoFEM::ConstrainMatrixCtx::initializeQorP
MoFEMErrorCode initializeQorP(Vec x)
initialize vectors and matrices for Q and P shell matrices, scattering is set based on x_problem and ...
Definition: ConstrainMatrixCtx.cpp:45
MoFEM::ConstrainMatrixCtx::maxIts
PetscInt maxIts
Definition: ConstrainMatrixCtx.hpp:73
ConstrainMatrixCtx::initializeQorP
MoFEMErrorCode initializeQorP(Vec x)
initialize vectors and matrices for Q and P shell matrices, scattering is set based on x_problem and ...
MoFEM::debug
const static bool debug
Definition: ConstrainMatrixCtx.cpp:9
ConstrainMatrixCtx::KQx
Vec KQx
Definition: ConstrainMatrixCtx.hpp:23
ConstrainMatrixCtx::destroyQorP
MoFEMErrorCode destroyQorP()
destroy sub-matrices used for shell matrices P, Q, R, RT
MoFEM::ConstrainMatrixCtx::Cx
Vec Cx
Definition: ConstrainMatrixCtx.hpp:24
MoFEM::ConstrainMatrixMultOpRT
MoFEMErrorCode ConstrainMatrixMultOpRT(Mat RT, Vec x, Vec f)
Multiplication operator for RT = (CCT)^-TC.
Definition: ConstrainMatrixCtx.cpp:243
MoFEM::ConstrainMatrixCtx::recalculateCTC
MoFEMErrorCode recalculateCTC()
re-calculate CTC matrix has been changed since initialization
Definition: ConstrainMatrixCtx.cpp:143
MoFEM::ConstrainMatrixCtx::destroyQorP
MoFEMErrorCode destroyQorP()
destroy sub-matrices used for shell matrices P, Q, R, RT
Definition: ConstrainMatrixCtx.cpp:94
MoFEM::CoreInterface::get_comm
virtual MPI_Comm & get_comm() const =0
ConstrainMatrixCtx::MOFEM_EVENT_projCTC_QTKQ
PetscLogEvent MOFEM_EVENT_projCTC_QTKQ
Definition: ConstrainMatrixCtx.hpp:41
MoFEM::ConstrainMatrixDestroyOpPorQ
MoFEMErrorCode ConstrainMatrixDestroyOpPorQ(Mat Q)
Destroy shell matrix Q.
Definition: ConstrainMatrixCtx.cpp:290
MoFEM::ConstrainMatrixCtx::CCT
Mat CCT
Definition: ConstrainMatrixCtx.hpp:23
MoFEM::ConstrainMatrixCtx::MOFEM_EVENT_projP
PetscLogEvent MOFEM_EVENT_projP
Definition: ConstrainMatrixCtx.hpp:40
MoFEM::ConstrainMatrixMultOpP
MoFEMErrorCode ConstrainMatrixMultOpP(Mat P, Vec x, Vec f)
Multiplication operator for P = CT(CCT)^-1C.
Definition: ConstrainMatrixCtx.cpp:197
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
ConstrainMatrixCtx::X
Vec X
Definition: ConstrainMatrixCtx.hpp:23
ConstrainMatrixCtx::C
Mat C
Definition: ConstrainMatrixCtx.hpp:21
MoFEM::ConstrainMatrixCtx::MOFEM_EVENT_projCTC_QTKQ
PetscLogEvent MOFEM_EVENT_projCTC_QTKQ
Definition: ConstrainMatrixCtx.hpp:43
ConstrainMatrixCtx::K
Mat K
Definition: ConstrainMatrixCtx.hpp:21
MoFEM::ConstrainMatrixCtx::absTol
PetscReal absTol
Definition: ConstrainMatrixCtx.hpp:72
MOFEM_IMPOSSIBLE_CASE
@ MOFEM_IMPOSSIBLE_CASE
Definition: definitions.h:35
ConstrainMatrixCtx::kSP
KSP kSP
Definition: ConstrainMatrixCtx.hpp:20
ConstrainMatrixCtx::initQorP
bool initQorP
Definition: ConstrainMatrixCtx.hpp:24
MoFEM::ConstrainMatrixCtx::initQTKQ
bool initQTKQ
Definition: ConstrainMatrixCtx.hpp:26
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::ConstrainMatrixCtx::C
Mat C
Definition: ConstrainMatrixCtx.hpp:23
MoFEM::ConstrainMatrixCtx::CTC
Mat CTC
Definition: ConstrainMatrixCtx.hpp:23
ConstrainMatrixCtx::destroyQTKQ
MoFEMErrorCode destroyQTKQ()
destroy sub-matrices used for shell matrix QTKQ
MoFEM::ConstrainMatrixCtx::CTCx
Vec CTCx
Definition: ConstrainMatrixCtx.hpp:24
ROW
@ ROW
Definition: definitions.h:136
MoFEM::ConstrainMatrixCtx::destroyQTKQ
MoFEMErrorCode destroyQTKQ()
destroy sub-matrices used for shell matrix QTKQ
Definition: ConstrainMatrixCtx.cpp:151
MoFEM::ConstrainMatrixCtx::initializeQTKQ
MoFEMErrorCode initializeQTKQ()
initialize vectors and matrices for CTC+QTKQ shell matrices, scattering is set based on x_problem and...
Definition: ConstrainMatrixCtx.cpp:114
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
ConstrainMatrixCtx::Cx
Vec Cx
Definition: ConstrainMatrixCtx.hpp:22
MoFEM::ConstrainMatrixCtx::cancelKSPMonitor
bool cancelKSPMonitor
Definition: ConstrainMatrixCtx.hpp:29
ConstrainMatrixCtx::MOFEM_EVENT_projR
PetscLogEvent MOFEM_EVENT_projR
Definition: ConstrainMatrixCtx.hpp:39
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::ConstrainMatrixCtx::MOFEM_EVENT_projQ
PetscLogEvent MOFEM_EVENT_projQ
Definition: ConstrainMatrixCtx.hpp:39
ConstrainMatrixCtx::Qx
Vec Qx
Definition: ConstrainMatrixCtx.hpp:23
MoFEM::ProjectionMatrixMultOpQ
MoFEMErrorCode ProjectionMatrixMultOpQ(Mat Q, Vec x, Vec f)
Multiplication operator for Q = I-CTC(CCT)^-1C.
Definition: ConstrainMatrixCtx.cpp:163
R
@ R
Definition: free_surface.cpp:394
MoFEM::ConstrainMatrixMultOpCTC_QTKQ
MoFEMErrorCode ConstrainMatrixMultOpCTC_QTKQ(Mat CTC_QTKQ, Vec x, Vec f)
Multiplication operator for RT = (CCT)^-TC.
Definition: ConstrainMatrixCtx.cpp:261
MoFEM::ConstrainMatrixCtx::CT
Mat CT
Definition: ConstrainMatrixCtx.hpp:23
MoFEM::ConstrainMatrixCtx::sCatter
VecScatter sCatter
Definition: ConstrainMatrixCtx.hpp:35
EshelbianPlasticity::P
@ P
Definition: EshelbianContact.cpp:197
COL
@ COL
Definition: definitions.h:136
ConstrainMatrixCtx::CT_CCTm1_Cx
Vec CT_CCTm1_Cx
Definition: ConstrainMatrixCtx.hpp:22
MoFEM::ConstrainMatrixCtx::dTol
PetscReal dTol
Definition: ConstrainMatrixCtx.hpp:72
MoFEM::ConstrainMatrixCtx::recalculateCTandCCT
MoFEMErrorCode recalculateCTandCCT()
re-calculate CT and CCT if C matrix has been changed since initialization
Definition: ConstrainMatrixCtx.cpp:85
ConstrainMatrixCtx::sCatter
VecScatter sCatter
Definition: ConstrainMatrixCtx.hpp:33
MoFEM::ConstrainMatrixCtx::xProblem
string xProblem
Definition: ConstrainMatrixCtx.hpp:36
MoFEM::ConstrainMatrixCtx::rTol
PetscReal rTol
Definition: ConstrainMatrixCtx.hpp:70
ConstrainMatrixCtx::MOFEM_EVENT_projRT
PetscLogEvent MOFEM_EVENT_projRT
Definition: ConstrainMatrixCtx.hpp:40
ConstrainMatrixCtx::CTC
Mat CTC
Definition: ConstrainMatrixCtx.hpp:21
MoFEM::ConstrainMatrixCtx::K
Mat K
Definition: ConstrainMatrixCtx.hpp:23
MoFEM::ConstrainMatrixCtx::MOFEM_EVENT_projR
PetscLogEvent MOFEM_EVENT_projR
Definition: ConstrainMatrixCtx.hpp:41
MoFEM::ConstrainMatrixCtx::createKSP
bool createKSP
Definition: ConstrainMatrixCtx.hpp:27
ConstrainMatrixCtx::CCTm1_Cx
Vec CCTm1_Cx
Definition: ConstrainMatrixCtx.hpp:22
MoFEM::ConstrainMatrixCtx::CT_CCTm1_Cx
Vec CT_CCTm1_Cx
Definition: ConstrainMatrixCtx.hpp:24
MoFEM::VecManager
Vector manager is used to create vectors \mofem_vectors.
Definition: VecManager.hpp:23
MoFEM::ConstrainMatrixDestroyOpQTKQ
MoFEMErrorCode ConstrainMatrixDestroyOpQTKQ(Mat QTKQ)
Destroy shell matrix.
Definition: ConstrainMatrixCtx.cpp:298
MoFEM::ConstrainMatrixCtx::createScatter
bool createScatter
Definition: ConstrainMatrixCtx.hpp:28
MoFEM::ConstrainMatrixCtx::mField
MoFEM::Interface & mField
Definition: ConstrainMatrixCtx.hpp:20
convert.n
n
Definition: convert.py:82
ConstrainMatrixCtx::MOFEM_EVENT_projP
PetscLogEvent MOFEM_EVENT_projP
Definition: ConstrainMatrixCtx.hpp:38
MoFEM::ConstrainMatrixCtx::MOFEM_EVENT_projInit
PetscLogEvent MOFEM_EVENT_projInit
Definition: ConstrainMatrixCtx.hpp:38
N
const int N
Definition: speed_test.cpp:3
MoFEM::ConstrainMatrixCtx::X
Vec X
Definition: ConstrainMatrixCtx.hpp:25
ConstrainMatrixCtx::CTCx
Vec CTCx
Definition: ConstrainMatrixCtx.hpp:22
MoFEM::ConstrainMatrixCtx::initQorP
bool initQorP
Definition: ConstrainMatrixCtx.hpp:26
MoFEM::ConstrainMatrixCtx::Qx
Vec Qx
Definition: ConstrainMatrixCtx.hpp:25
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
ConstrainMatrixCtx::mField
MoFEM::Interface & mField
Definition: ConstrainMatrixCtx.hpp:18
MoFEM::ConstrainMatrixMultOpR
MoFEMErrorCode ConstrainMatrixMultOpR(Mat R, Vec x, Vec f)
Multiplication operator for R = CT(CCT)^-1.
Definition: ConstrainMatrixCtx.cpp:221
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
ConstrainMatrixCtx::MOFEM_EVENT_projQ
PetscLogEvent MOFEM_EVENT_projQ
Definition: ConstrainMatrixCtx.hpp:37
MoFEM::ConstrainMatrixCtx::MOFEM_EVENT_projRT
PetscLogEvent MOFEM_EVENT_projRT
Definition: ConstrainMatrixCtx.hpp:42
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::ConstrainMatrixCtx::KQx
Vec KQx
Definition: ConstrainMatrixCtx.hpp:25
ConstrainMatrixCtx::CT
Mat CT
Definition: ConstrainMatrixCtx.hpp:21
INIT_DATA_CONSTRAINMATRIXCTX
#define INIT_DATA_CONSTRAINMATRIXCTX
Definition: ConstrainMatrixCtx.cpp:11
ConstrainMatrixCtx
structure for projection matrices
Definition: ConstrainMatrixCtx.hpp:16
ConstrainMatrixCtx::initializeQTKQ
MoFEMErrorCode initializeQTKQ()
initialize vectors and matrices for CTC+QTKQ shell matrices, scattering is set based on x_problem and...
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::ConstrainMatrixCtx::CCTm1_Cx
Vec CCTm1_Cx
Definition: ConstrainMatrixCtx.hpp:24
MoFEM::ConstrainMatrixCtx::ConstrainMatrixCtx
ConstrainMatrixCtx(MoFEM::Interface &m_field, string x_problem, string y_problem, bool create_ksp=true, bool own_contrain_matrix=false)
Definition: ConstrainMatrixCtx.cpp:19
MoFEM::ConstrainMatrixCtx::kSP
KSP kSP
Definition: ConstrainMatrixCtx.hpp:22
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359