Actual source code: schurm.c

  1: #include <../src/ksp/ksp/utils/schurm/schurm.h>

  3: const char *const MatSchurComplementAinvTypes[] = {"DIAG","LUMP","BLOCKDIAG","MatSchurComplementAinvType","MAT_SCHUR_COMPLEMENT_AINV_",NULL};

  5: PetscErrorCode MatCreateVecs_SchurComplement(Mat N,Vec *right,Vec *left)
  6: {
  7:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

  9:   if (Na->D) {
 10:     MatCreateVecs(Na->D,right,left);
 11:     return 0;
 12:   }
 13:   if (right) {
 14:     MatCreateVecs(Na->B,right,NULL);
 15:   }
 16:   if (left) {
 17:     MatCreateVecs(Na->C,NULL,left);
 18:   }
 19:   return 0;
 20: }

 22: PetscErrorCode MatView_SchurComplement(Mat N,PetscViewer viewer)
 23: {
 24:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

 26:   PetscViewerASCIIPrintf(viewer,"Schur complement A11 - A10 inv(A00) A01\n");
 27:   if (Na->D) {
 28:     PetscViewerASCIIPrintf(viewer,"A11\n");
 29:     PetscViewerASCIIPushTab(viewer);
 30:     MatView(Na->D,viewer);
 31:     PetscViewerASCIIPopTab(viewer);
 32:   } else {
 33:     PetscViewerASCIIPrintf(viewer,"A11 = 0\n");
 34:   }
 35:   PetscViewerASCIIPrintf(viewer,"A10\n");
 36:   PetscViewerASCIIPushTab(viewer);
 37:   MatView(Na->C,viewer);
 38:   PetscViewerASCIIPopTab(viewer);
 39:   PetscViewerASCIIPrintf(viewer,"KSP of A00\n");
 40:   PetscViewerASCIIPushTab(viewer);
 41:   KSPView(Na->ksp,viewer);
 42:   PetscViewerASCIIPopTab(viewer);
 43:   PetscViewerASCIIPrintf(viewer,"A01\n");
 44:   PetscViewerASCIIPushTab(viewer);
 45:   MatView(Na->B,viewer);
 46:   PetscViewerASCIIPopTab(viewer);
 47:   return 0;
 48: }

 50: /*
 51:            A11^T - A01^T ksptrans(A00,Ap00) A10^T
 52: */
 53: PetscErrorCode MatMultTranspose_SchurComplement(Mat N,Vec x,Vec y)
 54: {
 55:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

 57:   if (!Na->work1) MatCreateVecs(Na->A,&Na->work1,NULL);
 58:   if (!Na->work2) MatCreateVecs(Na->A,&Na->work2,NULL);
 59:   MatMultTranspose(Na->C,x,Na->work1);
 60:   KSPSolveTranspose(Na->ksp,Na->work1,Na->work2);
 61:   MatMultTranspose(Na->B,Na->work2,y);
 62:   VecScale(y,-1.0);
 63:   if (Na->D) {
 64:     MatMultTransposeAdd(Na->D,x,y,y);
 65:   }
 66:   return 0;
 67: }

 69: /*
 70:            A11 - A10 ksp(A00,Ap00) A01
 71: */
 72: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
 73: {
 74:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

 76:   if (!Na->work1) MatCreateVecs(Na->A,&Na->work1,NULL);
 77:   if (!Na->work2) MatCreateVecs(Na->A,&Na->work2,NULL);
 78:   MatMult(Na->B,x,Na->work1);
 79:   KSPSolve(Na->ksp,Na->work1,Na->work2);
 80:   MatMult(Na->C,Na->work2,y);
 81:   VecScale(y,-1.0);
 82:   if (Na->D) {
 83:     MatMultAdd(Na->D,x,y,y);
 84:   }
 85:   return 0;
 86: }

 88: /*
 89:            A11 - A10 ksp(A00,Ap00) A01
 90: */
 91: PetscErrorCode MatMultAdd_SchurComplement(Mat N,Vec x,Vec y,Vec z)
 92: {
 93:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

 95:   if (!Na->work1) MatCreateVecs(Na->A,&Na->work1,NULL);
 96:   if (!Na->work2) MatCreateVecs(Na->A,&Na->work2,NULL);
 97:   MatMult(Na->B,x,Na->work1);
 98:   KSPSolve(Na->ksp,Na->work1,Na->work2);
 99:   if (y == z) {
100:     VecScale(Na->work2,-1.0);
101:     MatMultAdd(Na->C,Na->work2,z,z);
102:   } else {
103:     MatMult(Na->C,Na->work2,z);
104:     VecAYPX(z,-1.0,y);
105:   }
106:   if (Na->D) {
107:     MatMultAdd(Na->D,x,z,z);
108:   }
109:   return 0;
110: }

112: PetscErrorCode MatSetFromOptions_SchurComplement(PetscOptionItems *PetscOptionsObject,Mat N)
113: {
114:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

116:   PetscOptionsHead(PetscOptionsObject,"MatSchurComplementOptions");
117:   Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
118:   PetscOptionsEnum("-mat_schur_complement_ainv_type","Type of approximation for DIAGFORM(A00) used when assembling Sp = A11 - A10 inv(DIAGFORM(A00)) A01","MatSchurComplementSetAinvType",MatSchurComplementAinvTypes,(PetscEnum)Na->ainvtype,(PetscEnum*)&Na->ainvtype,NULL);
119:   PetscOptionsTail();
120:   KSPSetFromOptions(Na->ksp);
121:   return 0;
122: }

124: PetscErrorCode MatDestroy_SchurComplement(Mat N)
125: {
126:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

128:   MatDestroy(&Na->A);
129:   MatDestroy(&Na->Ap);
130:   MatDestroy(&Na->B);
131:   MatDestroy(&Na->C);
132:   MatDestroy(&Na->D);
133:   VecDestroy(&Na->work1);
134:   VecDestroy(&Na->work2);
135:   KSPDestroy(&Na->ksp);
136:   PetscFree(N->data);
137:   return 0;
138: }

140: /*@
141:       MatCreateSchurComplement - Creates a new Mat that behaves like the Schur complement of a matrix

143:    Collective on A00

145:    Input Parameters:
146: +   A00,A01,A10,A11  - the four parts of the original matrix A = [A00 A01; A10 A11] (A11 is optional)
147: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}

149:    Output Parameter:
150: .   S - the matrix that the Schur complement S = A11 - A10 ksp(A00,Ap00) A01

152:    Level: intermediate

154:    Notes:
155:     The Schur complement is NOT explicitly formed! Rather, this function returns a virtual Schur complement
156:     that can compute the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
157:     for Schur complement S and a KSP solver to approximate the action of A^{-1}.

159:     All four matrices must have the same MPI communicator.

161:     A00 and  A11 must be square matrices.

163:     MatGetSchurComplement() takes as arguments the index sets for the submatrices and returns both the virtual Schur complement (what this returns) plus
164:     a sparse approximation to the Schur complement (useful for building a preconditioner for the Schur complement) which can be obtained from this
165:     matrix with MatSchurComplementGetPmat()

167:     Developer Notes:
168:     The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
169:     remove redundancy and be clearer and simpler.

171: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdateSubMatrices(), MatCreateTranspose(), MatGetSchurComplement(),
172:           MatSchurComplementGetPmat(), MatSchurComplementSetSubMatrices()

174: @*/
175: PetscErrorCode  MatCreateSchurComplement(Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11,Mat *S)
176: {
177:   KSPInitializePackage();
178:   MatCreate(PetscObjectComm((PetscObject)A00),S);
179:   MatSetType(*S,MATSCHURCOMPLEMENT);
180:   MatSchurComplementSetSubMatrices(*S,A00,Ap00,A01,A10,A11);
181:   return 0;
182: }

184: /*@
185:       MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement

187:    Collective on S

189:    Input Parameters:
190: +   S                - matrix obtained with MatSetType(S,MATSCHURCOMPLEMENT)
191: .   A00,A01,A10,A11  - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
192: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.

194:    Level: intermediate

196:    Notes:
197:      The Schur complement is NOT explicitly formed! Rather, this
198:      object performs the matrix-vector product of the Schur complement by using formula S = A11 - A10 ksp(A00,Ap00) A01

200:      All four matrices must have the same MPI communicator.

202:      A00 and  A11 must be square matrices.

204:      This is to be used in the context of code such as
205: $        MatSetType(S,MATSCHURCOMPLEMENT);
206: $        MatSchurComplementSetSubMatrices(S,...);

208:     while MatSchurComplementUpdateSubMatrices() should only be called after MatCreateSchurComplement() or MatSchurComplementSetSubMatrices()

210: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdateSubMatrices(), MatCreateTranspose(), MatCreateSchurComplement(), MatGetSchurComplement()

212: @*/
213: PetscErrorCode  MatSchurComplementSetSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
214: {
215:   Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;
216:   PetscBool           isschur;

218:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&isschur);
219:   if (!isschur) return 0;
233:   if (A11) {
237:   }

239:   MatSetSizes(S,A10->rmap->n,A01->cmap->n,A10->rmap->N,A01->cmap->N);
240:   PetscObjectReference((PetscObject)A00);
241:   PetscObjectReference((PetscObject)Ap00);
242:   PetscObjectReference((PetscObject)A01);
243:   PetscObjectReference((PetscObject)A10);
244:   Na->A  = A00;
245:   Na->Ap = Ap00;
246:   Na->B  = A01;
247:   Na->C  = A10;
248:   Na->D  = A11;
249:   if (A11) {
250:     PetscObjectReference((PetscObject)A11);
251:   }
252:   MatSetUp(S);
253:   KSPSetOperators(Na->ksp,A00,Ap00);
254:   S->assembled = PETSC_TRUE;
255:   return 0;
256: }

258: /*@
259:   MatSchurComplementGetKSP - Gets the KSP object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01

261:   Not Collective

263:   Input Parameter:
264: . S - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01

266:   Output Parameter:
267: . ksp - the linear solver object

269:   Options Database:
270: . -fieldsplit_<splitname_0>_XXX sets KSP and PC options for the 0-split solver inside the Schur complement used in PCFieldSplit; default <splitname_0> is 0.

272:   Level: intermediate

274: .seealso: MatSchurComplementSetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate()
275: @*/
276: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
277: {
278:   Mat_SchurComplement *Na;
279:   PetscBool           isschur;

282:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&isschur);
285:   Na   = (Mat_SchurComplement*) S->data;
286:   *ksp = Na->ksp;
287:   return 0;
288: }

290: /*@
291:   MatSchurComplementSetKSP - Sets the KSP object that is used to solve with A00 in the Schur complement matrix S = A11 - A10 ksp(A00,Ap00) A01

293:   Not Collective

295:   Input Parameters:
296: + S   - matrix created with MatCreateSchurComplement()
297: - ksp - the linear solver object

299:   Level: developer

301:   Developer Notes:
302:     This is used in PCFieldSplit to reuse the 0-split KSP to implement ksp(A00,Ap00) in S.
303:     The KSP operators are overwritten with A00 and Ap00 currently set in S.

305: .seealso: MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate(), MATSCHURCOMPLEMENT
306: @*/
307: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
308: {
309:   Mat_SchurComplement *Na;
310:   PetscBool           isschur;

313:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&isschur);
314:   if (!isschur) return 0;
316:   Na      = (Mat_SchurComplement*) S->data;
317:   PetscObjectReference((PetscObject)ksp);
318:   KSPDestroy(&Na->ksp);
319:   Na->ksp = ksp;
320:   KSPSetOperators(Na->ksp, Na->A, Na->Ap);
321:   return 0;
322: }

324: /*@
325:       MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices

327:    Collective on S

329:    Input Parameters:
330: +   S                - matrix obtained with MatCreateSchurComplement() (or MatSchurSetSubMatrices()) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
331: .   A00,A01,A10,A11  - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
332: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.

334:    Level: intermediate

336:    Notes:
337:      All four matrices must have the same MPI communicator

339:      A00 and  A11 must be square matrices

341:      All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement() or MatSchurComplementSetSubMatrices()
342:      though they need not be the same matrices.

344:      This can only be called after MatCreateSchurComplement() or MatSchurComplementSetSubMatrices(), it cannot be called immediately after MatSetType(S,MATSCHURCOMPLEMENT);

346:    Developer Notes:
347:      This code is almost identical to MatSchurComplementSetSubMatrices(). The API should be refactored.

349: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()

351: @*/
352: PetscErrorCode  MatSchurComplementUpdateSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
353: {
354:   Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;
355:   PetscBool           isschur;

358:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&isschur);
359:   if (!isschur) return 0;
373:   if (A11) {
377:   }

379:   PetscObjectReference((PetscObject)A00);
380:   PetscObjectReference((PetscObject)Ap00);
381:   PetscObjectReference((PetscObject)A01);
382:   PetscObjectReference((PetscObject)A10);
383:   if (A11) {
384:     PetscObjectReference((PetscObject)A11);
385:   }

387:   MatDestroy(&Na->A);
388:   MatDestroy(&Na->Ap);
389:   MatDestroy(&Na->B);
390:   MatDestroy(&Na->C);
391:   MatDestroy(&Na->D);

393:   Na->A  = A00;
394:   Na->Ap = Ap00;
395:   Na->B  = A01;
396:   Na->C  = A10;
397:   Na->D  = A11;

399:   KSPSetOperators(Na->ksp,A00,Ap00);
400:   return 0;
401: }

403: /*@C
404:   MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement

406:   Collective on S

408:   Input Parameter:
409: . S    - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01

411:   Output Parameters:
412: + A00  - the upper-left block of the original matrix A = [A00 A01; A10 A11]
413: . Ap00 - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}
414: . A01  - the upper-right block of the original matrix A = [A00 A01; A10 A11]
415: . A10  - the lower-left block of the original matrix A = [A00 A01; A10 A11]
416: - A11  - (optional) the lower-right block of the original matrix A = [A00 A01; A10 A11]

418:   Note: A11 is optional, and thus can be NULL.  The reference counts of the submatrices are not increased before they are returned and the matrices should not be modified or destroyed.

420:   Level: intermediate

422: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdateSubMatrices()
423: @*/
424: PetscErrorCode  MatSchurComplementGetSubMatrices(Mat S,Mat *A00,Mat *Ap00,Mat *A01,Mat *A10,Mat *A11)
425: {
426:   Mat_SchurComplement *Na = (Mat_SchurComplement*) S->data;
427:   PetscBool           flg;

430:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&flg);
432:   if (A00) *A00 = Na->A;
433:   if (Ap00) *Ap00 = Na->Ap;
434:   if (A01) *A01 = Na->B;
435:   if (A10) *A10 = Na->C;
436:   if (A11) *A11 = Na->D;
437:   return 0;
438: }

440: #include <petsc/private/kspimpl.h>

442: /*@
443:   MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly

445:   Collective on M

447:   Input Parameter:
448: . M - the matrix obtained with MatCreateSchurComplement()

450:   Output Parameter:
451: . S - the Schur complement matrix

453:   Notes:
454:     This can be expensive, so it is mainly for testing

456:     Use MatSchurComplementGetPmat() to get a sparse approximation for the Schur complement suitable for use in building a preconditioner

458:   Level: advanced

460: .seealso: MatCreateSchurComplement(), MatSchurComplementUpdate(), MatSchurComplementGetPmat()
461: @*/
462: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat A, Mat *S)
463: {
464:   Mat            B, C, D, Bd, AinvBd;
465:   KSP            ksp;
466:   PetscInt       n,N,m,M;

468:   MatSchurComplementGetSubMatrices(A, NULL, NULL, &B, &C, &D);
469:   MatSchurComplementGetKSP(A, &ksp);
470:   KSPSetUp(ksp);
471:   MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd);
472:   MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
473:   KSPMatSolve(ksp, Bd, AinvBd);
474:   MatDestroy(&Bd);
475:   MatChop(AinvBd, PETSC_SMALL);
476:   if (D) {
477:     MatGetLocalSize(D, &m, &n);
478:     MatGetSize(D, &M, &N);
479:     MatCreateDense(PetscObjectComm((PetscObject)A), m, n, M, N, NULL, S);
480:   }
481:   MatMatMult(C, AinvBd, D ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX, PETSC_DEFAULT, S);
482:   MatDestroy(&AinvBd);
483:   if (D) {
484:     MatAXPY(*S, -1.0, D, DIFFERENT_NONZERO_PATTERN);
485:   }
486:   MatConvert(*S, MATAIJ, MAT_INPLACE_MATRIX, S);
487:   MatScale(*S, -1.0);
488:   return 0;
489: }

491: /* Developer Notes:
492:     This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
493: PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *S,MatSchurComplementAinvType ainvtype, MatReuse preuse,Mat *Sp)
494: {
495:   Mat            A=NULL,Ap=NULL,B=NULL,C=NULL,D=NULL;
496:   MatReuse       reuse;

507:   if (mreuse == MAT_IGNORE_MATRIX && preuse == MAT_IGNORE_MATRIX) return 0;


513:   reuse = MAT_INITIAL_MATRIX;
514:   if (mreuse == MAT_REUSE_MATRIX) {
515:     MatSchurComplementGetSubMatrices(*S,&A,&Ap,&B,&C,&D);
518:     MatDestroy(&Ap); /* get rid of extra reference */
519:     reuse = MAT_REUSE_MATRIX;
520:   }
521:   MatCreateSubMatrix(mat,isrow0,iscol0,reuse,&A);
522:   MatCreateSubMatrix(mat,isrow0,iscol1,reuse,&B);
523:   MatCreateSubMatrix(mat,isrow1,iscol0,reuse,&C);
524:   MatCreateSubMatrix(mat,isrow1,iscol1,reuse,&D);
525:   switch (mreuse) {
526:   case MAT_INITIAL_MATRIX:
527:     MatCreateSchurComplement(A,A,B,C,D,S);
528:     break;
529:   case MAT_REUSE_MATRIX:
530:     MatSchurComplementUpdateSubMatrices(*S,A,A,B,C,D);
531:     break;
532:   default:
534:   }
535:   if (preuse != MAT_IGNORE_MATRIX) {
536:     MatCreateSchurComplementPmat(A,B,C,D,ainvtype,preuse,Sp);
537:   }
538:   MatDestroy(&A);
539:   MatDestroy(&B);
540:   MatDestroy(&C);
541:   MatDestroy(&D);
542:   return 0;
543: }

545: /*@
546:     MatGetSchurComplement - Obtain the Schur complement from eliminating part of the matrix in another part.

548:     Collective on A

550:     Input Parameters:
551: +   A      - matrix in which the complement is to be taken
552: .   isrow0 - rows to eliminate
553: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
554: .   isrow1 - rows in which the Schur complement is formed
555: .   iscol1 - columns in which the Schur complement is formed
556: .   mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in S
557: .   ainvtype - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
558:                        MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, or MAT_SCHUR_COMPLEMENT_AINV_LUMP
559: -   preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in Sp

561:     Output Parameters:
562: +   S      - exact Schur complement, often of type MATSCHURCOMPLEMENT which is difficult to use for preconditioning
563: -   Sp     - approximate Schur complement from which a preconditioner can be built A11 - A10 inv(DIAGFORM(A00)) A01

565:     Note:
566:     Since the real Schur complement is usually dense, providing a good approximation to Sp usually requires
567:     application-specific information.

569:     Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
570:     and column index sets.  In that case, the user should call PetscObjectComposeFunction() on the *S matrix and pass mreuse of MAT_REUSE_MATRIX to set
571:     "MatGetSchurComplement_C" to their function.  If their function needs to fall back to the default implementation, it
572:     should call MatGetSchurComplement_Basic().

574:     MatCreateSchurComplement() takes as arguments the four submatrices and returns the virtual Schur complement (what this function returns in S).

576:     MatSchurComplementGetPmat() takes the virtual Schur complement and returns an explicit approximate Schur complement (what this returns in Sp).

578:     In other words calling MatCreateSchurComplement() followed by MatSchurComplementGetPmat() produces the same output as this function but with slightly different
579:     inputs. The actually submatrices of the original block matrix instead of index sets to the submatrices.

581:     Developer Notes:
582:     The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
583:     remove redundancy and be clearer and simpler.

585:     Level: advanced

587: .seealso: MatCreateSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement(), MatSchurComplementAinvType
588: @*/
589: PetscErrorCode  MatGetSchurComplement(Mat A,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *S,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Sp)
590: {
591:   PetscErrorCode (*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*) = NULL;

605:   if (mreuse == MAT_REUSE_MATRIX) { /* This is the only situation, in which we can demand that the user pass a non-NULL pointer to non-garbage in S. */
606:     PetscObjectQueryFunction((PetscObject)*S,"MatGetSchurComplement_C",&f);
607:   }
608:   if (f) (*f)(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,preuse,Sp);
609:   else MatGetSchurComplement_Basic(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,ainvtype,preuse,Sp);
610:   return 0;
611: }

613: /*@
614:     MatSchurComplementSetAinvType - set the type of approximation used for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()

616:     Not collective.

618:     Input Parameters:
619: +   S        - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
620: -   ainvtype - type of approximation to be used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
621:                       MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, or MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG

623:     Options database:
624:     -mat_schur_complement_ainv_type diag | lump | blockdiag

626:     Level: advanced

628: .seealso: MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementGetAinvType()
629: @*/
630: PetscErrorCode  MatSchurComplementSetAinvType(Mat S,MatSchurComplementAinvType ainvtype)
631: {
632:   PetscBool           isschur;
633:   Mat_SchurComplement *schur;

636:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&isschur);
637:   if (!isschur) return 0;
639:   schur = (Mat_SchurComplement*)S->data;
641:   schur->ainvtype = ainvtype;
642:   return 0;
643: }

645: /*@
646:     MatSchurComplementGetAinvType - get the type of approximation for the inverse of the (0,0) block used in forming Sp in MatSchurComplementGetPmat()

648:     Not collective.

650:     Input Parameter:
651: .   S      - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01

653:     Output Parameter:
654: .   ainvtype - type of approximation used to form approximate Schur complement Sp = A11 - A10 inv(DIAGFORM(A00)) A01:
655:                       MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, or MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG

657:     Level: advanced

659: .seealso: MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementSetAinvType()
660: @*/
661: PetscErrorCode  MatSchurComplementGetAinvType(Mat S,MatSchurComplementAinvType *ainvtype)
662: {
663:   PetscBool           isschur;
664:   Mat_SchurComplement *schur;

667:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&isschur);
669:   schur = (Mat_SchurComplement*)S->data;
670:   if (ainvtype) *ainvtype = schur->ainvtype;
671:   return 0;
672: }

674: /*@
675:     MatCreateSchurComplementPmat - create a preconditioning matrix for the Schur complement by explicitly assembling the sparse matrix
676:         Sp = A11 - A10 inv(DIAGFORM(A00)) A01

678:     Collective on A00

680:     Input Parameters:
681: +   A00      - the upper-left part of the original matrix A = [A00 A01; A10 A11]
682: .   A01      - (optional) the upper-right part of the original matrix A = [A00 A01; A10 A11]
683: .   A10      - (optional) the lower-left part of the original matrix A = [A00 A01; A10 A11]
684: .   A11      - (optional) the lower-right part of the original matrix A = [A00 A01; A10 A11]
685: .   ainvtype - type of approximation for DIAGFORM(A00) used when forming Sp = A11 - A10 inv(DIAGFORM(A00)) A01. See MatSchurComplementAinvType.
686: -   preuse   - MAT_INITIAL_MATRIX for a new Sp, or MAT_REUSE_MATRIX to reuse an existing Sp, or MAT_IGNORE_MATRIX to put nothing in Sp

688:     Output Parameter:
689: -   Sp    - approximate Schur complement suitable for preconditioning the true Schur complement S = A11 - A10 inv(A00) A01

691:     Level: advanced

693: .seealso: MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementAinvType
694: @*/
695: PetscErrorCode  MatCreateSchurComplementPmat(Mat A00,Mat A01,Mat A10,Mat A11,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Sp)
696: {
697:   PetscInt       N00;

699:   /* Use an appropriate approximate inverse of A00 to form A11 - A10 inv(DIAGFORM(A00)) A01; a NULL A01, A10 or A11 indicates a zero matrix. */
700:   /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
703:   if (preuse == MAT_IGNORE_MATRIX) return 0;

705:   /* A zero size A00 or empty A01 or A10 imply S = A11. */
706:   MatGetSize(A00,&N00,NULL);
707:   if (!A01 || !A10 || !N00) {
708:     if (preuse == MAT_INITIAL_MATRIX) {
709:       MatDuplicate(A11,MAT_COPY_VALUES,Sp);
710:     } else { /* MAT_REUSE_MATRIX */
711:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
712:       MatCopy(A11,*Sp,DIFFERENT_NONZERO_PATTERN);
713:     }
714:   } else {
715:     Mat AdB;
716:     Vec diag;

718:     if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP || ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
719:       MatDuplicate(A01,MAT_COPY_VALUES,&AdB);
720:       MatCreateVecs(A00,&diag,NULL);
721:       if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
722:         MatGetRowSum(A00,diag);
723:       } else {
724:         MatGetDiagonal(A00,diag);
725:       }
726:       VecReciprocal(diag);
727:       MatDiagonalScale(AdB,diag,NULL);
728:       VecDestroy(&diag);
729:     } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG) {
730:       Mat      A00_inv;
731:       MatType  type;
732:       MPI_Comm comm;

734:       PetscObjectGetComm((PetscObject)A00,&comm);
735:       MatGetType(A00,&type);
736:       MatCreate(comm,&A00_inv);
737:       MatSetType(A00_inv,type);
738:       MatInvertBlockDiagonalMat(A00,A00_inv);
739:       MatMatMult(A00_inv,A01,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AdB);
740:       MatDestroy(&A00_inv);
741:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatSchurComplementAinvType: %d", ainvtype);
742:     /* Cannot really reuse Sp in MatMatMult() because of MatAYPX() -->
743:          MatAXPY() --> MatHeaderReplace() --> MatDestroy_XXX_MatMatMult()  */
744:     MatDestroy(Sp);
745:     MatMatMult(A10,AdB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,Sp);
746:     if (!A11) {
747:       MatScale(*Sp,-1.0);
748:     } else {
749:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
750:       MatAYPX(*Sp,-1,A11,DIFFERENT_NONZERO_PATTERN);
751:     }
752:     MatDestroy(&AdB);
753:   }
754:   return 0;
755: }

757: PetscErrorCode  MatSchurComplementGetPmat_Basic(Mat S,MatReuse preuse,Mat *Sp)
758: {
759:   Mat A,B,C,D;
760:   Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;

762:   if (preuse == MAT_IGNORE_MATRIX) return 0;
763:   MatSchurComplementGetSubMatrices(S,&A,NULL,&B,&C,&D);
765:   MatCreateSchurComplementPmat(A,B,C,D,schur->ainvtype,preuse,Sp);
766:   return 0;
767: }

769: /*@
770:     MatSchurComplementGetPmat - Obtain a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(DIAGFORM(A00)) A01

772:     Collective on S

774:     Input Parameters:
775: +   S      - matrix obtained with MatCreateSchurComplement() (or equivalent) that implements the action of A11 - A10 ksp(A00,Ap00) A01
776: -   preuse - MAT_INITIAL_MATRIX for a new Sp, or MAT_REUSE_MATRIX to reuse an existing Sp, or MAT_IGNORE_MATRIX to put nothing in Sp

778:     Output Parameter:
779: -   Sp     - approximate Schur complement suitable for preconditioning the exact Schur complement S = A11 - A10 inv(A00) A01

781:     Note:
782:     The approximation of Sp depends on the the argument passed to to MatSchurComplementSetAinvType()
783:     MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, or MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG or
784:     -mat_schur_complement_ainv_type <diag,lump,blockdiag>

786:     Sometimes users would like to provide problem-specific data in the Schur complement, usually only
787:     for special row and column index sets.  In that case, the user should call PetscObjectComposeFunction() to set
788:     "MatSchurComplementGetPmat_C" to their function.  If their function needs to fall back to the default implementation,
789:     it should call MatSchurComplementGetPmat_Basic().

791:     Developer Notes:
792:     The API that includes MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementGetPmat() should be refactored to
793:     remove redundancy and be clearer and simpler.

795:     This routine should be called MatSchurComplementCreatePmat()

797:     Level: advanced

799: .seealso: MatCreateSubMatrix(), PCFIELDSPLIT, MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementSetAinvType()
800: @*/
801: PetscErrorCode  MatSchurComplementGetPmat(Mat S,MatReuse preuse,Mat *Sp)
802: {
803:   PetscErrorCode (*f)(Mat,MatReuse,Mat*);


812:   PetscObjectQueryFunction((PetscObject)S,"MatSchurComplementGetPmat_C",&f);
813:   if (f) (*f)(S,preuse,Sp);
814:   else MatSchurComplementGetPmat_Basic(S,preuse,Sp);
815:   return 0;
816: }

818: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
819: {
820:   Mat_SchurComplement *Na;

822:   PetscNewLog(N,&Na);
823:   N->data = (void*) Na;

825:   N->ops->destroy        = MatDestroy_SchurComplement;
826:   N->ops->getvecs        = MatCreateVecs_SchurComplement;
827:   N->ops->view           = MatView_SchurComplement;
828:   N->ops->mult           = MatMult_SchurComplement;
829:   N->ops->multtranspose  = MatMultTranspose_SchurComplement;
830:   N->ops->multadd        = MatMultAdd_SchurComplement;
831:   N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
832:   N->assembled           = PETSC_FALSE;
833:   N->preallocated        = PETSC_FALSE;

835:   KSPCreate(PetscObjectComm((PetscObject)N),&Na->ksp);
836:   PetscObjectChangeTypeName((PetscObject)N,MATSCHURCOMPLEMENT);
837:   return 0;
838: }