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