v0.14.0
Functions
Constrain Projection Matrix
Collaboration diagram for Constrain Projection Matrix:

Functions

MoFEMErrorCode MoFEM::ProjectionMatrixMultOpQ (Mat Q, Vec x, Vec f)
 Multiplication operator for Q = I-CTC(CCT)^-1C. More...
 
MoFEMErrorCode MoFEM::ConstrainMatrixMultOpP (Mat P, Vec x, Vec f)
 Multiplication operator for P = CT(CCT)^-1C. More...
 
MoFEMErrorCode MoFEM::ConstrainMatrixMultOpR (Mat R, Vec x, Vec f)
 Multiplication operator for R = CT(CCT)^-1. More...
 
MoFEMErrorCode MoFEM::ConstrainMatrixMultOpRT (Mat RT, Vec x, Vec f)
 Multiplication operator for RT = (CCT)^-TC. More...
 
MoFEMErrorCode MoFEM::ConstrainMatrixMultOpCTC_QTKQ (Mat CTC_QTKQ, Vec x, Vec f)
 Multiplication operator for RT = (CCT)^-TC. More...
 
MoFEMErrorCode MoFEM::ConstrainMatrixDestroyOpPorQ (Mat Q)
 Destroy shell matrix Q. More...
 
MoFEMErrorCode MoFEM::ConstrainMatrixDestroyOpQTKQ (Mat QTKQ)
 Destroy shell matrix. More...
 
MoFEMErrorCode ProjectionMatrixMultOpQ (Mat Q, Vec x, Vec f)
 Multiplication operator for Q = I-CTC(CCT)^-1C. More...
 
MoFEMErrorCode ConstrainMatrixMultOpP (Mat P, Vec x, Vec f)
 Multiplication operator for P = CT(CCT)^-1C. More...
 
MoFEMErrorCode ConstrainMatrixMultOpR (Mat R, Vec x, Vec f)
 Multiplication operator for R = CT(CCT)^-1. More...
 
MoFEMErrorCode ConstrainMatrixMultOpRT (Mat RT, Vec x, Vec f)
 Multiplication operator for RT = (CCT)^-TC. More...
 
MoFEMErrorCode ConstrainMatrixMultOpCTC_QTKQ (Mat CTC_QTKQ, Vec x, Vec f)
 Multiplication operator for RT = (CCT)^-TC. More...
 
MoFEMErrorCode ConstrainMatrixDestroyOpPorQ (Mat Q)
 Destroy shell matrix Q. More...
 
MoFEMErrorCode ConstrainMatrixDestroyOpQTKQ (Mat QTKQ)
 Destroy shell matrix. More...
 

Detailed Description

Function Documentation

◆ ConstrainMatrixDestroyOpPorQ() [1/2]

MoFEMErrorCode ConstrainMatrixDestroyOpPorQ ( Mat  Q)

Destroy shell matrix Q.

Mat Q; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&Q);
MatShellSetOperation(Q,MATOP_MULT,(void(*)(void))ProjectionMatrixMultOpQ);
MatShellSetOperation(Q,MATOP_DESTROY,(void(*)(void))ConstrainMatrixDestroyOpPorQ);

Definition at line 290 of file ConstrainMatrixCtx.cpp.

290  {
292  void *void_ctx;
293  CHKERR MatShellGetContext(Q, &void_ctx);
294  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
295  CHKERR ctx->destroyQorP();
297 }

◆ ConstrainMatrixDestroyOpPorQ() [2/2]

MoFEMErrorCode MoFEM::ConstrainMatrixDestroyOpPorQ ( Mat  Q)

Destroy shell matrix Q.

Mat Q; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&Q);
MatShellSetOperation(Q,MATOP_MULT,(void(*)(void))ProjectionMatrixMultOpQ);
MatShellSetOperation(Q,MATOP_DESTROY,(void(*)(void))ConstrainMatrixDestroyOpPorQ);

Definition at line 290 of file ConstrainMatrixCtx.cpp.

290  {
292  void *void_ctx;
293  CHKERR MatShellGetContext(Q, &void_ctx);
294  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
295  CHKERR ctx->destroyQorP();
297 }

◆ ConstrainMatrixDestroyOpQTKQ() [1/2]

MoFEMErrorCode ConstrainMatrixDestroyOpQTKQ ( Mat  QTKQ)

Destroy shell matrix.

Mat CTC_QTKQ; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&Q);
MatShellSetOperation(Q,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpCTC_QTKQ);
MatShellSetOperation(Q,MATOP_DESTROY,(void(*)(void))mat_destroy_QTKQ);

Definition at line 298 of file ConstrainMatrixCtx.cpp.

298  {
300  void *void_ctx;
301  CHKERR MatShellGetContext(QTKQ, &void_ctx);
302  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
303  CHKERR ctx->destroyQTKQ();
305 }

◆ ConstrainMatrixDestroyOpQTKQ() [2/2]

MoFEMErrorCode MoFEM::ConstrainMatrixDestroyOpQTKQ ( Mat  QTKQ)

Destroy shell matrix.

Mat CTC_QTKQ; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&Q);
MatShellSetOperation(Q,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpCTC_QTKQ);
MatShellSetOperation(Q,MATOP_DESTROY,(void(*)(void))mat_destroy_QTKQ);

Definition at line 298 of file ConstrainMatrixCtx.cpp.

298  {
300  void *void_ctx;
301  CHKERR MatShellGetContext(QTKQ, &void_ctx);
302  ConstrainMatrixCtx *ctx = (ConstrainMatrixCtx *)void_ctx;
303  CHKERR ctx->destroyQTKQ();
305 }

◆ ConstrainMatrixMultOpCTC_QTKQ() [1/2]

MoFEMErrorCode ConstrainMatrixMultOpCTC_QTKQ ( Mat  CTC_QTKQ,
Vec  x,
Vec  f 
)

Multiplication operator for RT = (CCT)^-TC.

Mat CTC_QTKQ; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&CTC_QTKQ);
MatShellSetOperation(CTC_QTKQ,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpCTC_QTKQ);
MatShellSetOperation(CTC_QTKQ,MATOP_DESTROY,(void(*)(void))ConstrainMatrixDestroyOpQTKQ);

Definition at line 261 of file ConstrainMatrixCtx.cpp.

261  {
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 }

◆ ConstrainMatrixMultOpCTC_QTKQ() [2/2]

MoFEMErrorCode MoFEM::ConstrainMatrixMultOpCTC_QTKQ ( Mat  CTC_QTKQ,
Vec  x,
Vec  f 
)

Multiplication operator for RT = (CCT)^-TC.

Mat CTC_QTKQ; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&CTC_QTKQ);
MatShellSetOperation(CTC_QTKQ,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpCTC_QTKQ);
MatShellSetOperation(CTC_QTKQ,MATOP_DESTROY,(void(*)(void))ConstrainMatrixDestroyOpQTKQ);

Definition at line 261 of file ConstrainMatrixCtx.cpp.

261  {
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 }

◆ ConstrainMatrixMultOpP() [1/2]

MoFEMErrorCode ConstrainMatrixMultOpP ( Mat  P,
Vec  x,
Vec  f 
)

Multiplication operator for P = CT(CCT)^-1C.

Mat P; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&P);
MatShellSetOperation(P,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpP);

Definition at line 197 of file ConstrainMatrixCtx.cpp.

197  {
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 }

◆ ConstrainMatrixMultOpP() [2/2]

MoFEMErrorCode MoFEM::ConstrainMatrixMultOpP ( Mat  P,
Vec  x,
Vec  f 
)

Multiplication operator for P = CT(CCT)^-1C.

Mat P; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&P);
MatShellSetOperation(P,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpP);

Definition at line 197 of file ConstrainMatrixCtx.cpp.

197  {
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 }

◆ ConstrainMatrixMultOpR() [1/2]

MoFEMErrorCode ConstrainMatrixMultOpR ( Mat  R,
Vec  x,
Vec  f 
)

Multiplication operator for R = CT(CCT)^-1.

Mat R; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&R);
MatShellSetOperation(R,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpR);

Definition at line 221 of file ConstrainMatrixCtx.cpp.

221  {
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 }

◆ ConstrainMatrixMultOpR() [2/2]

MoFEMErrorCode MoFEM::ConstrainMatrixMultOpR ( Mat  R,
Vec  x,
Vec  f 
)

Multiplication operator for R = CT(CCT)^-1.

Mat R; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&R);
MatShellSetOperation(R,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpR);

Definition at line 221 of file ConstrainMatrixCtx.cpp.

221  {
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 }

◆ ConstrainMatrixMultOpRT() [1/2]

MoFEMErrorCode ConstrainMatrixMultOpRT ( Mat  RT,
Vec  x,
Vec  f 
)

Multiplication operator for RT = (CCT)^-TC.

Mat RT; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&RT);
MatShellSetOperation(RT,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpRT);

Definition at line 243 of file ConstrainMatrixCtx.cpp.

243  {
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 }

◆ ConstrainMatrixMultOpRT() [2/2]

MoFEMErrorCode MoFEM::ConstrainMatrixMultOpRT ( Mat  RT,
Vec  x,
Vec  f 
)

Multiplication operator for RT = (CCT)^-TC.

Mat RT; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&RT);
MatShellSetOperation(RT,MATOP_MULT,(void(*)(void))ConstrainMatrixMultOpRT);

Definition at line 243 of file ConstrainMatrixCtx.cpp.

243  {
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 }

◆ ProjectionMatrixMultOpQ() [1/2]

MoFEMErrorCode ProjectionMatrixMultOpQ ( Mat  Q,
Vec  x,
Vec  f 
)

Multiplication operator for Q = I-CTC(CCT)^-1C.

Mat Q; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&Q);
MatShellSetOperation(Q,MATOP_MULT,(void(*)(void))ProjectionMatrixMultOpQ);
MatShellSetOperation(Q,MATOP_DESTROY,(void(*)(void))ConstrainMatrixDestroyOpPorQ);

Definition at line 163 of file ConstrainMatrixCtx.cpp.

163  {
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 }

◆ ProjectionMatrixMultOpQ() [2/2]

MoFEMErrorCode MoFEM::ProjectionMatrixMultOpQ ( Mat  Q,
Vec  x,
Vec  f 
)

Multiplication operator for Q = I-CTC(CCT)^-1C.

Mat Q; //for problem
projection_matrix_ctx(m_field,problem_name,contrains_problem_name);
CHKERR MatCreateShell(PETSC_COMM_WORLD,m,m,M,M,&projection_matrix_ctx,&Q);
MatShellSetOperation(Q,MATOP_MULT,(void(*)(void))ProjectionMatrixMultOpQ);
MatShellSetOperation(Q,MATOP_DESTROY,(void(*)(void))ConstrainMatrixDestroyOpPorQ);

Definition at line 163 of file ConstrainMatrixCtx.cpp.

163  {
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 }
ConstrainMatrixCtx::initializeQorP
MoFEMErrorCode initializeQorP(Vec x)
initialize vectors and matrices for Q and P shell matrices, scattering is set based on x_problem and ...
ConstrainMatrixCtx::KQx
Vec KQx
Definition: ConstrainMatrixCtx.hpp:23
ConstrainMatrixCtx::destroyQorP
MoFEMErrorCode destroyQorP()
destroy sub-matrices used for shell matrices P, Q, R, RT
MoFEM::ConstrainMatrixMultOpRT
MoFEMErrorCode ConstrainMatrixMultOpRT(Mat RT, Vec x, Vec f)
Multiplication operator for RT = (CCT)^-TC.
Definition: ConstrainMatrixCtx.cpp:243
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::ConstrainMatrixMultOpP
MoFEMErrorCode ConstrainMatrixMultOpP(Mat P, Vec x, Vec f)
Multiplication operator for P = CT(CCT)^-1C.
Definition: ConstrainMatrixCtx.cpp:197
ConstrainMatrixCtx::X
Vec X
Definition: ConstrainMatrixCtx.hpp:23
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
ConstrainMatrixCtx::destroyQTKQ
MoFEMErrorCode destroyQTKQ()
destroy sub-matrices used for shell matrix QTKQ
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::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
EshelbianPlasticity::P
@ P
Definition: EshelbianContact.cpp:197
ConstrainMatrixCtx::CT_CCTm1_Cx
Vec CT_CCTm1_Cx
Definition: ConstrainMatrixCtx.hpp:22
ConstrainMatrixCtx::sCatter
VecScatter sCatter
Definition: ConstrainMatrixCtx.hpp:33
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::ConstrainMatrixDestroyOpQTKQ
MoFEMErrorCode ConstrainMatrixDestroyOpQTKQ(Mat QTKQ)
Destroy shell matrix.
Definition: ConstrainMatrixCtx.cpp:298
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::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
debug
const static bool debug
Definition: ConstrainMatrixCtx.cpp:13
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
ConstrainMatrixCtx::MOFEM_EVENT_projQ
PetscLogEvent MOFEM_EVENT_projQ
Definition: ConstrainMatrixCtx.hpp:37
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
ConstrainMatrixCtx::CT
Mat CT
Definition: ConstrainMatrixCtx.hpp:21
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
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359