23 yProblem(y_problem) {
24 PetscLogEventRegister("ProjectionInit", 0, &MOFEM_EVENT_projInit);
25 PetscLogEventRegister("ProjectionQ", 0, &MOFEM_EVENT_projQ);
26 PetscLogEventRegister("ProjectionP", 0, &MOFEM_EVENT_projP);
27 PetscLogEventRegister("ProjectionR", 0, &MOFEM_EVENT_projR);
28 PetscLogEventRegister("ProjectionRT", 0, &MOFEM_EVENT_projRT);
29 PetscLogEventRegister("ProjectionCTC_QTKQ", 0, &MOFEM_EVENT_projCTC_QTKQ);
30}
31
33 VecScatter scatter, bool create_ksp,
34 bool own_contrain_matrix)
36 PetscLogEventRegister("ProjectionInit", 0, &MOFEM_EVENT_projInit);
37 PetscLogEventRegister("ProjectionQ", 0, &MOFEM_EVENT_projQ);
38 PetscLogEventRegister("ProjectionP", 0, &MOFEM_EVENT_projP);
39 PetscLogEventRegister("ProjectionR", 0, &MOFEM_EVENT_projR);
40 PetscLogEventRegister("ProjectionRT", 0, &MOFEM_EVENT_projRT);
41 PetscLogEventRegister("ProjectionCTC_QTKQ", 0, &MOFEM_EVENT_projCTC_QTKQ);
42}
43
46 if (initQorP) {
47 initQorP = false;
48
49 PetscLogEventBegin(MOFEM_EVENT_projInit, 0, 0, 0, 0);
50 CHKERR MatTranspose(C, MAT_INITIAL_MATRIX, &CT);
51
52 CHKERR MatTransposeMatMult(CT, CT, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &CCT);
53 if (createKSP) {
54 CHKERR KSPCreate(mField.get_comm(), &kSP);
55
56 CHKERR KSPSetOperators(kSP, CCT, CCT);
57 CHKERR KSPSetFromOptions(kSP);
58 CHKERR KSPSetInitialGuessKnoll(kSP, PETSC_TRUE);
59 CHKERR KSPGetTolerances(kSP, &rTol, &absTol, &dTol, &maxIts);
61 if (cancelKSPMonitor) {
62 CHKERR KSPMonitorCancel(kSP);
63 }
64 }
65#if PETSC_VERSION_GE(3, 5, 3)
66 CHKERR MatCreateVecs(C, &X, PETSC_NULLPTR);
67 CHKERR MatCreateVecs(C, PETSC_NULLPTR, &Cx);
68 CHKERR MatCreateVecs(CCT, PETSC_NULLPTR, &CCTm1_Cx);
69#else
70 CHKERR MatGetVecs(C, &X, PETSC_NULLPTR);
71 CHKERR MatGetVecs(C, PETSC_NULLPTR, &Cx);
72 CHKERR MatGetVecs(CCT, PETSC_NULLPTR, &CCTm1_Cx);
73#endif
74 CHKERR VecDuplicate(X, &CT_CCTm1_Cx);
75 if (createScatter) {
76 CHKERR mField.getInterface<VecManager>()->vecScatterCreate(
77 x, xProblem,
ROW, X, yProblem,
COL, &sCatter);
78 }
79 PetscLogEventEnd(MOFEM_EVENT_projInit, 0, 0, 0, 0);
80 }
82}
83
86 if (initQorP)
88 CHKERR MatTranspose(C, MAT_REUSE_MATRIX, &CT);
89 CHKERR MatTransposeMatMult(CT, CT, MAT_REUSE_MATRIX, PETSC_DEFAULT, &CCT);
91}
92
95 if (initQorP)
99 if (createKSP) {
101 }
104 CHKERR VecDestroy(&CCTm1_Cx);
105 CHKERR VecDestroy(&CT_CCTm1_Cx);
106 if (createScatter) {
107 CHKERR VecScatterDestroy(&sCatter);
108 }
109 initQorP = true;
111}
112
115 if (initQTKQ) {
116 initQTKQ = false;
117 PetscLogEventBegin(MOFEM_EVENT_projInit, 0, 0, 0, 0);
118
119 CHKERR MatTransposeMatMult(C, C, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &CTC);
121
123 MatGetSize(CCT, &
m, &
n);
124 PetscPrintf(mField.get_comm(),
"CCT size (%d,%d)\n",
m,
n);
125
126
127 }
128#if PETSC_VERSION_GE(3, 5, 3)
129 CHKERR MatCreateVecs(K, &Qx, PETSC_NULLPTR);
130 CHKERR MatCreateVecs(K, PETSC_NULLPTR, &KQx);
131 CHKERR MatCreateVecs(CTC, PETSC_NULLPTR, &CTCx);
132#else
133 CHKERR MatGetVecs(K, &Qx, PETSC_NULLPTR);
134 CHKERR MatGetVecs(K, PETSC_NULLPTR, &KQx);
135 CHKERR MatGetVecs(CTC, PETSC_NULLPTR, &CTCx);
136#endif
137 PetscLogEventEnd(MOFEM_EVENT_projInit, 0, 0, 0, 0);
138 }
140}
141
144 if (initQTKQ)
146 CHKERR MatTransposeMatMult(C, C, MAT_REUSE_MATRIX, PETSC_DEFAULT, &CTC);
148}
149
152 if (initQTKQ)
158 initQTKQ = true;
160}
161
164 void *void_ctx;
165 CHKERR MatShellGetContext(Q, &void_ctx);
166 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
167 PetscLogEventBegin(ctx->MOFEM_EVENT_projQ, 0, 0, 0, 0);
168 CHKERR ctx->initializeQorP(x);
170 CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
171 SCATTER_FORWARD);
172 CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
174
175 CHKERR VecScatterBegin(ctx->sCatter, ctx->X, f, INSERT_VALUES,
176 SCATTER_REVERSE);
177 CHKERR VecScatterEnd(ctx->sCatter, ctx->X, f, INSERT_VALUES,
178 SCATTER_REVERSE);
179 PetscBool flg;
180 CHKERR VecEqual(x, f, &flg);
181 if (flg == PETSC_FALSE)
183 }
184 CHKERR MatMult(ctx->C, ctx->X, ctx->Cx);
185 CHKERR KSPSolve(ctx->kSP, ctx->Cx, ctx->CCTm1_Cx);
186 CHKERR MatMult(ctx->CT, ctx->CCTm1_Cx, ctx->CT_CCTm1_Cx);
187 CHKERR VecScale(ctx->CT_CCTm1_Cx, -1);
188 CHKERR VecScatterBegin(ctx->sCatter, ctx->CT_CCTm1_Cx, f, ADD_VALUES,
189 SCATTER_REVERSE);
190 CHKERR VecScatterEnd(ctx->sCatter, ctx->CT_CCTm1_Cx, f, ADD_VALUES,
191 SCATTER_REVERSE);
192 PetscLogEventEnd(ctx->MOFEM_EVENT_projQ, 0, 0, 0, 0);
194}
195
198 void *void_ctx;
199 CHKERR MatShellGetContext(P, &void_ctx);
200 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
201 PetscLogEventBegin(ctx->MOFEM_EVENT_projP, 0, 0, 0, 0);
202 CHKERR ctx->initializeQorP(x);
203 CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
204 SCATTER_FORWARD);
205 CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
206 CHKERR MatMult(ctx->C, ctx->X, ctx->Cx);
207 CHKERR KSPSolve(ctx->kSP, ctx->Cx, ctx->CCTm1_Cx);
208 CHKERR MatMult(ctx->CT, ctx->CCTm1_Cx, ctx->CT_CCTm1_Cx);
210 CHKERR VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD);
211 CHKERR VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD);
212 CHKERR VecScatterBegin(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
213 SCATTER_REVERSE);
214 CHKERR VecScatterEnd(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
215 SCATTER_REVERSE);
216 PetscLogEventEnd(ctx->MOFEM_EVENT_projP, 0, 0, 0, 0);
218}
219
222 void *void_ctx;
223 CHKERR MatShellGetContext(
R, &void_ctx);
224 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
225 PetscLogEventBegin(ctx->MOFEM_EVENT_projR, 0, 0, 0, 0);
226 if (ctx->initQorP)
228 "you have to call first initQorP or use Q matrix");
229 CHKERR KSPSolve(ctx->kSP, x, ctx->CCTm1_Cx);
230 CHKERR MatMult(ctx->CT, ctx->CCTm1_Cx, ctx->CT_CCTm1_Cx);
232 CHKERR VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD);
233 CHKERR VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD);
234 CHKERR VecScatterBegin(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
235 SCATTER_REVERSE);
236 CHKERR VecScatterEnd(ctx->sCatter, ctx->CT_CCTm1_Cx, f, INSERT_VALUES,
237 SCATTER_REVERSE);
238 PetscLogEventEnd(ctx->MOFEM_EVENT_projR, 0, 0, 0, 0);
240}
241
244 void *void_ctx;
245 CHKERR MatShellGetContext(RT, &void_ctx);
246 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
247 PetscLogEventBegin(ctx->MOFEM_EVENT_projRT, 0, 0, 0, 0);
248 if (ctx->initQorP)
250 "you have to call first initQorP or use Q matrix");
251 CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
252 SCATTER_FORWARD);
253 CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
254 CHKERR MatMult(ctx->C, ctx->X, ctx->Cx);
255 CHKERR KSPSolve(ctx->kSP, ctx->Cx, f);
256 PetscLogEventEnd(ctx->MOFEM_EVENT_projRT, 0, 0, 0, 0);
258}
259
262 void *void_ctx;
263 CHKERR MatShellGetContext(CTC_QTKQ, &void_ctx);
264 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
265 PetscLogEventBegin(ctx->MOFEM_EVENT_projCTC_QTKQ, 0, 0, 0, 0);
266 Mat Q;
268 CHKERR MatGetSize(ctx->K, &M, &
N);
269 CHKERR MatGetLocalSize(ctx->K, &
m, &
n);
270 CHKERR MatCreateShell(ctx->mField.get_comm(),
m,
n, M,
N, ctx, &Q);
271 CHKERR MatShellSetOperation(Q, MATOP_MULT,
272 (void (*)(void))ProjectionMatrixMultOpQ);
273 CHKERR ctx->initializeQTKQ();
274 CHKERR MatMult(Q, x, ctx->Qx);
275 CHKERR MatMult(ctx->K, ctx->Qx, ctx->KQx);
276 CHKERR MatMult(Q, ctx->KQx, f);
277 CHKERR VecScatterBegin(ctx->sCatter, x, ctx->X, INSERT_VALUES,
278 SCATTER_FORWARD);
279 CHKERR VecScatterEnd(ctx->sCatter, x, ctx->X, INSERT_VALUES, SCATTER_FORWARD);
280 CHKERR MatMult(ctx->CTC, ctx->X, ctx->CTCx);
281 CHKERR VecScatterBegin(ctx->sCatter, ctx->CTCx, f, ADD_VALUES,
282 SCATTER_REVERSE);
283 CHKERR VecScatterEnd(ctx->sCatter, ctx->CTCx, f, ADD_VALUES, SCATTER_REVERSE);
285 PetscLogEventEnd(ctx->MOFEM_EVENT_projCTC_QTKQ, 0, 0, 0, 0);
287}
288
291 void *void_ctx;
292 CHKERR MatShellGetContext(Q, &void_ctx);
293 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
294 CHKERR ctx->destroyQorP();
296}
299 void *void_ctx;
300 CHKERR MatShellGetContext(QTKQ, &void_ctx);
301 ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
302 CHKERR ctx->destroyQTKQ();
304}
305
306}
#define INIT_DATA_CONSTRAINMATRIXCTX
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MoFEMErrorCode ConstrainMatrixMultOpCTC_QTKQ(Mat CTC_QTKQ, Vec x, Vec f)
Multiplication operator for RT = (CCT)^-TC.
MoFEMErrorCode ConstrainMatrixMultOpP(Mat P, Vec x, Vec f)
Multiplication operator for P = CT(CCT)^-1C.
MoFEMErrorCode ConstrainMatrixDestroyOpPorQ(Mat Q)
Destroy shell matrix Q.
MoFEMErrorCode ProjectionMatrixMultOpQ(Mat Q, Vec x, Vec f)
Multiplication operator for Q = I-CTC(CCT)^-1C.
MoFEMErrorCode ConstrainMatrixMultOpR(Mat R, Vec x, Vec f)
Multiplication operator for R = CT(CCT)^-1.
MoFEMErrorCode ConstrainMatrixDestroyOpQTKQ(Mat QTKQ)
Destroy shell matrix.
MoFEMErrorCode ConstrainMatrixMultOpRT(Mat RT, Vec x, Vec f)
Multiplication operator for RT = (CCT)^-TC.
const double n
refractive index of diffusive medium
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
FTensor::Index< 'M', 3 > M
FTensor::Index< 'm', 3 > m
Deprecated interface functions.