v0.15.0
Loading...
Searching...
No Matches
ConstrainMatrixCtx.cpp
Go to the documentation of this file.
1/** \file ConstrainMatrixCtx.cpp
2 * \brief Implementation of projection matrix
3 *
4 *
5 */
6
8namespace MoFEM {
9const static bool debug = false;
10
11#define INIT_DATA_CONSTRAINMATRIXCTX \
12 C(PETSC_NULLPTR), CT(PETSC_NULLPTR), CCT(PETSC_NULLPTR), CTC(PETSC_NULLPTR), \
13 K(PETSC_NULLPTR), Cx(PETSC_NULLPTR), CCTm1_Cx(PETSC_NULLPTR), \
14 CT_CCTm1_Cx(PETSC_NULLPTR), CTCx(PETSC_NULLPTR), Qx(PETSC_NULLPTR), \
15 KQx(PETSC_NULLPTR), 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_NULLPTR);
68 CHKERR MatCreateVecs(C, PETSC_NULLPTR, &Cx);
69 CHKERR MatCreateVecs(CCT, PETSC_NULLPTR, &CCTm1_Cx);
70#else
71 CHKERR MatGetVecs(C, &X, PETSC_NULLPTR);
72 CHKERR MatGetVecs(C, PETSC_NULLPTR, &Cx);
73 CHKERR MatGetVecs(CCT, PETSC_NULLPTR, &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_NULLPTR);
131 CHKERR MatCreateVecs(K, PETSC_NULLPTR, &KQx);
132 CHKERR MatCreateVecs(CTC, PETSC_NULLPTR, &CTCx);
133#else
134 CHKERR MatGetVecs(K, &Qx, PETSC_NULLPTR);
135 CHKERR MatGetVecs(K, PETSC_NULLPTR, &KQx);
136 CHKERR MatGetVecs(CTC, PETSC_NULLPTR, &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
#define INIT_DATA_CONSTRAINMATRIXCTX
@ COL
@ ROW
#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 ...
@ MOFEM_IMPOSSIBLE_CASE
Definition definitions.h:35
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
@ R
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.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
static const bool debug
FTensor::Index< 'm', 3 > m
const int N
Definition speed_test.cpp:3
structure for projection matrices
MoFEMErrorCode destroyQTKQ()
destroy sub-matrices used for shell matrix QTKQ
ConstrainMatrixCtx(MoFEM::Interface &m_field, string x_problem, string y_problem, bool create_ksp=true, bool own_contrain_matrix=false)
MoFEMErrorCode recalculateCTC()
re-calculate CTC matrix has been changed since initialization
MoFEMErrorCode destroyQorP()
destroy sub-matrices used for shell matrices P, Q, R, RT
MoFEMErrorCode initializeQTKQ()
initialize vectors and matrices for CTC+QTKQ shell matrices, scattering is set based on x_problem and...
MoFEMErrorCode initializeQorP(Vec x)
initialize vectors and matrices for Q and P shell matrices, scattering is set based on x_problem and ...
MoFEMErrorCode recalculateCTandCCT()
re-calculate CT and CCT if C matrix has been changed since initialization
virtual MPI_Comm & get_comm() const =0
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.