Actual source code: schurm.c

petsc-master 2017-07-26
Report Typos and Errors
  1:  #include <petsc/private/matimpl.h>
  2:  #include <petscksp.h>
  3: const char *const MatSchurComplementAinvTypes[] = {"DIAG","LUMP","MatSchurComplementAinvType","MAT_SCHUR_COMPLEMENT_AINV_",0};

  5: typedef struct {
  6:   Mat                        A,Ap,B,C,D;
  7:   KSP                        ksp;
  8:   Vec                        work1,work2;
  9:   MatSchurComplementAinvType ainvtype;
 10: } Mat_SchurComplement;

 12: PetscErrorCode MatCreateVecs_SchurComplement(Mat N,Vec *right,Vec *left)
 13: {
 14:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 15:   PetscErrorCode      ierr;

 18:   if (Na->D) {
 19:     MatCreateVecs(Na->D,right,left);
 20:     return(0);
 21:   }
 22:   if (right) {
 23:     MatCreateVecs(Na->B,right,NULL);
 24:   }
 25:   if (left) {
 26:     MatCreateVecs(Na->C,NULL,left);
 27:   }
 28:   return(0);
 29: }

 31: PetscErrorCode MatView_SchurComplement(Mat N,PetscViewer viewer)
 32: {
 33:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 34:   PetscErrorCode      ierr;

 37:   PetscViewerASCIIPrintf(viewer,"Schur complement A11 - A10 inv(A00) A01\n");
 38:   if (Na->D) {
 39:     PetscViewerASCIIPrintf(viewer,"A11\n");
 40:     PetscViewerASCIIPushTab(viewer);
 41:     MatView(Na->D,viewer);
 42:     PetscViewerASCIIPopTab(viewer);
 43:   } else {
 44:     PetscViewerASCIIPrintf(viewer,"A11 = 0\n");
 45:   }
 46:   PetscViewerASCIIPrintf(viewer,"A10\n");
 47:   PetscViewerASCIIPushTab(viewer);
 48:   MatView(Na->C,viewer);
 49:   PetscViewerASCIIPopTab(viewer);
 50:   PetscViewerASCIIPrintf(viewer,"KSP of A00\n");
 51:   PetscViewerASCIIPushTab(viewer);
 52:   KSPView(Na->ksp,viewer);
 53:   PetscViewerASCIIPopTab(viewer);
 54:   PetscViewerASCIIPrintf(viewer,"A01\n");
 55:   PetscViewerASCIIPushTab(viewer);
 56:   MatView(Na->B,viewer);
 57:   PetscViewerASCIIPopTab(viewer);
 58:   return(0);
 59: }

 61: /*
 62:            A11^T - A01^T ksptrans(A00,Ap00) A10^T
 63: */
 64: PetscErrorCode MatMultTranspose_SchurComplement(Mat N,Vec x,Vec y)
 65: {
 66:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 67:   PetscErrorCode      ierr;

 70:   if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
 71:   if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
 72:   MatMultTranspose(Na->C,x,Na->work1);
 73:   KSPSolveTranspose(Na->ksp,Na->work1,Na->work2);
 74:   MatMultTranspose(Na->B,Na->work2,y);
 75:   VecScale(y,-1.0);
 76:   if (Na->D) {
 77:     MatMultTransposeAdd(Na->D,x,y,y);
 78:   }
 79:   return(0);
 80: }

 82: /*
 83:            A11 - A10 ksp(A00,Ap00) A01
 84: */
 85: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
 86: {
 87:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 88:   PetscErrorCode      ierr;

 91:   if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
 92:   if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
 93:   MatMult(Na->B,x,Na->work1);
 94:   KSPSolve(Na->ksp,Na->work1,Na->work2);
 95:   MatMult(Na->C,Na->work2,y);
 96:   VecScale(y,-1.0);
 97:   if (Na->D) {
 98:     MatMultAdd(Na->D,x,y,y);
 99:   }
100:   return(0);
101: }

103: /*
104:            A11 - A10 ksp(A00,Ap00) A01
105: */
106: PetscErrorCode MatMultAdd_SchurComplement(Mat N,Vec x,Vec y,Vec z)
107: {
108:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
109:   PetscErrorCode      ierr;

112:   if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
113:   if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
114:   MatMult(Na->B,x,Na->work1);
115:   KSPSolve(Na->ksp,Na->work1,Na->work2);
116:   if (y == z) {
117:     VecScale(Na->work2,-1.0);
118:     MatMultAdd(Na->C,Na->work2,z,z);
119:   } else {
120:     MatMult(Na->C,Na->work2,z);
121:     VecAYPX(z,-1.0,y);
122:   }
123:   if (Na->D) {
124:     MatMultAdd(Na->D,x,z,z);
125:   }
126:   return(0);
127: }

129: PetscErrorCode MatSetFromOptions_SchurComplement(PetscOptionItems *PetscOptionsObject,Mat N)
130: {
131:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
132:   PetscErrorCode      ierr;

135:   PetscOptionsHead(PetscOptionsObject,"MatSchurComplementOptions");
136:   Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
137:   PetscOptionsEnum("-mat_schur_complement_ainv_type","Type of approximation for inv(A00) used when assembling Sp = A11 - A10 inv(A00) A01","MatSchurComplementSetAinvType",MatSchurComplementAinvTypes,(PetscEnum)Na->ainvtype,(PetscEnum*)&Na->ainvtype,NULL);
138:   PetscOptionsTail();
139:   KSPSetFromOptions(Na->ksp);
140:   return(0);
141: }

143: PetscErrorCode MatDestroy_SchurComplement(Mat N)
144: {
145:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
146:   PetscErrorCode      ierr;

149:   MatDestroy(&Na->A);
150:   MatDestroy(&Na->Ap);
151:   MatDestroy(&Na->B);
152:   MatDestroy(&Na->C);
153:   MatDestroy(&Na->D);
154:   VecDestroy(&Na->work1);
155:   VecDestroy(&Na->work2);
156:   KSPDestroy(&Na->ksp);
157:   PetscFree(N->data);
158:   return(0);
159: }

161: /*@
162:       MatCreateSchurComplement - Creates a new matrix object that behaves like the Schur complement of a matrix

164:    Collective on Mat

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

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

173:    Level: intermediate

175:    Notes: The Schur complement is NOT actually formed! Rather, this
176:           object performs the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
177:           for Schur complement S and a KSP solver to approximate the action of A^{-1}.

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

181:           A00 and  A11 must be square matrices.

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

185: @*/
186: PetscErrorCode  MatCreateSchurComplement(Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11,Mat *S)
187: {

191:   KSPInitializePackage();
192:   MatCreate(((PetscObject)A00)->comm,S);
193:   MatSetType(*S,MATSCHURCOMPLEMENT);
194:   MatSchurComplementSetSubMatrices(*S,A00,Ap00,A01,A10,A11);
195:   return(0);
196: }

198: /*@
199:       MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement

201:    Collective on Mat

203:    Input Parameter:
204: +   S                - matrix obtained with MatCreateSchurComplement (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
205: .   A00,A01,A10,A11  - the four parts of A = [A00 A01; A10 A11] (A11 is optional)
206: -   Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.

208:    Level: intermediate

210:    Notes: The Schur complement is NOT actually formed! Rather, this
211:           object performs the matrix-vector product by using formula S = A11 - A10 A^{-1} A01
212:           for Schur complement S and a KSP solver to approximate the action of A^{-1}.

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

216:           A00 and  A11 must be square matrices.

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

220: @*/
221: PetscErrorCode  MatSchurComplementSetSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
222: {
223:   PetscErrorCode      ierr;
224:   PetscInt            m,n;
225:   Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;

228:   if (S->assembled) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementUpdateSubMatrices() for already used matrix");
236:   if (A00->rmap->n != A00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local columns %D",A00->rmap->n,A00->cmap->n);
237:   if (A00->rmap->n != Ap00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local rows of Ap00 %D",A00->rmap->n,Ap00->rmap->n);
238:   if (Ap00->rmap->n != Ap00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of Ap00 %D do not equal local columns %D",Ap00->rmap->n,Ap00->cmap->n);
239:   if (A00->cmap->n != A01->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A00 %D do not equal local rows of A01 %D",A00->cmap->n,A01->rmap->n);
240:   if (A10->cmap->n != A00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A10 %D do not equal local rows of A00 %D",A10->cmap->n,A00->rmap->n);
241:   if (A11) {
244:     if (A10->rmap->n != A11->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A10 %D do not equal local rows A11 %D",A10->rmap->n,A11->rmap->n);
245:   }

247:   MatGetLocalSize(A01,NULL,&n);
248:   MatGetLocalSize(A10,&m,NULL);
249:   MatSetSizes(S,m,n,PETSC_DECIDE,PETSC_DECIDE);
250:   PetscObjectReference((PetscObject)A00);
251:   PetscObjectReference((PetscObject)Ap00);
252:   PetscObjectReference((PetscObject)A01);
253:   PetscObjectReference((PetscObject)A10);
254:   Na->A  = A00;
255:   Na->Ap = Ap00;
256:   Na->B  = A01;
257:   Na->C  = A10;
258:   Na->D  = A11;
259:   if (A11) {
260:     PetscObjectReference((PetscObject)A11);
261:   }
262:   S->assembled    = PETSC_TRUE;
263:   S->preallocated = PETSC_TRUE;

265:   PetscLayoutSetUp((S)->rmap);
266:   PetscLayoutSetUp((S)->cmap);
267:   KSPSetOperators(Na->ksp,A00,Ap00);
268:   return(0);
269: }

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

274:   Not Collective

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

279:   Output Parameter:
280: . ksp - the linear solver object

282:   Options Database:
283: . -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.

285:   Level: intermediate

287: .seealso: MatSchurComplementSetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate()
288: @*/
289: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
290: {
291:   Mat_SchurComplement *Na;

296:   Na   = (Mat_SchurComplement*) S->data;
297:   *ksp = Na->ksp;
298:   return(0);
299: }

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

304:   Not Collective

306:   Input Parameters:
307: + S   - matrix created with MatCreateSchurComplement()
308: - ksp - the linear solver object

310:   Level: developer

312:   Developer Notes:
313:     This is used in PCFieldSplit to reuse the 0-split KSP to implement ksp(A00,Ap00) in S.

315: .seealso: MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate(), MATSCHURCOMPLEMENT
316: @*/
317: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
318: {
319:   Mat_SchurComplement *Na;
320:   PetscErrorCode      ierr;

325:   Na      = (Mat_SchurComplement*) S->data;
326:   PetscObjectReference((PetscObject)ksp);
327:   KSPDestroy(&Na->ksp);
328:   Na->ksp = ksp;
329:   KSPSetOperators(Na->ksp, Na->A, Na->Ap);
330:   return(0);
331: }

333: /*@
334:       MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices

336:    Collective on Mat

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

343:    Level: intermediate

345:    Notes: All four matrices must have the same MPI communicator

347:           A00 and  A11 must be square matrices

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

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

354: @*/
355: PetscErrorCode  MatSchurComplementUpdateSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
356: {
357:   PetscErrorCode      ierr;
358:   Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;

361:   if (!S->assembled) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementSetSubMatrices() for a new matrix");
368:   if (A00->rmap->n != A00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local columns %D",A00->rmap->n,A00->cmap->n);
369:   if (A00->rmap->n != Ap00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A00 %D do not equal local rows of Ap00 %D",A00->rmap->n,Ap00->rmap->n);
370:   if (Ap00->rmap->n != Ap00->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of Ap00 %D do not equal local columns %D",Ap00->rmap->n,Ap00->cmap->n);
371:   if (A00->cmap->n != A01->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A00 %D do not equal local rows of A01 %D",A00->cmap->n,A01->rmap->n);
372:   if (A10->cmap->n != A00->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A10 %D do not equal local rows of A00 %D",A10->cmap->n,A00->rmap->n);
373:   if (A11) {
376:     if (A10->rmap->n != A11->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A10 %D do not equal local rows A11 %D",A10->rmap->n,A11->rmap->n);
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: }


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

407:   Collective on Mat

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

412:   Output Paramters:
413: + A00,A01,A10,A11  - the four parts of the original matrix A = [A00 A01; A10 A11] (A11 is optional)
414: - Ap00             - preconditioning matrix for use in ksp(A00,Ap00) to approximate the action of A^{-1}.

416:   Note: A11 is optional, and thus can be NULL.  The submatrices are not increfed before they are returned and should not be modified or destroyed.

418:   Level: intermediate

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

430:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&flg);
431:   if (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:   } else {
438:     if (A00) *A00 = 0;
439:     if (Ap00) *Ap00 = 0;
440:     if (A01) *A01 = 0;
441:     if (A10) *A10 = 0;
442:     if (A11) *A11 = 0;
443:   }
444:   return(0);
445: }

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

449: /*@
450:   MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly

452:   Collective on Mat

454:   Input Parameter:
455: . M - the matrix obtained with MatCreateSchurComplement()

457:   Output Parameter:
458: . S - the Schur complement matrix

460:   Note: This can be expensive, so it is mainly for testing

462:   Level: advanced

464: .seealso: MatCreateSchurComplement(), MatSchurComplementUpdate()
465: @*/
466: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat M, Mat *S)
467: {
468:   Mat            B, C, D;
469:   KSP            ksp;
470:   PC             pc;
471:   PetscBool      isLU, isILU;
472:   PetscReal      fill = 2.0;

476:   MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);
477:   MatSchurComplementGetKSP(M, &ksp);
478:   KSPGetPC(ksp, &pc);
479:   PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);
480:   PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);
481:   if (isLU || isILU) {
482:     Mat       fact, Bd, AinvB, AinvBd;
483:     PetscReal eps = 1.0e-10;

485:     /* This can be sped up for banded LU */
486:     KSPSetUp(ksp);
487:     PCFactorGetMatrix(pc, &fact);
488:     MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd);
489:     MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
490:     MatMatSolve(fact, Bd, AinvBd);
491:     MatDestroy(&Bd);
492:     MatChop(AinvBd, eps);
493:     MatConvert(AinvBd, MATAIJ, MAT_INITIAL_MATRIX, &AinvB);
494:     MatDestroy(&AinvBd);
495:     MatMatMult(C, AinvB, MAT_INITIAL_MATRIX, fill, S);
496:     MatDestroy(&AinvB);
497:   } else {
498:     Mat Ainvd, Ainv;

500:     PCComputeExplicitOperator(pc, &Ainvd);
501:     MatConvert(Ainvd, MATAIJ, MAT_INITIAL_MATRIX, &Ainv);
502:     MatDestroy(&Ainvd);
503: #if 0
504:     /* Symmetric version */
505:     MatPtAP(Ainv, B, MAT_INITIAL_MATRIX, fill, S);
506: #else
507:     /* Nonsymmetric version */
508:     MatMatMatMult(C, Ainv, B, MAT_INITIAL_MATRIX, fill, S);
509: #endif
510:     MatDestroy(&Ainv);
511:   }
512:   if (D) {
513:     MatAXPY(*S, -1.0, D, DIFFERENT_NONZERO_PATTERN);
514:    }
515:   MatScale(*S,-1.0);
516:   return(0);
517: }

519: /* Developer Notes: This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
520: PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatSchurComplementAinvType ainvtype, MatReuse preuse,Mat *newpmat)
521: {
523:   Mat            A=0,Ap=0,B=0,C=0,D=0;
524:   MatReuse       reuse;

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

537:   if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

539:   reuse = MAT_INITIAL_MATRIX;
540:   if (mreuse == MAT_REUSE_MATRIX) {
541:     MatSchurComplementGetSubMatrices(*newmat,&A,&Ap,&B,&C,&D);
542:     if (!A || !Ap || !B || !C) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset");
543:     if (A != Ap) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator");
544:     MatDestroy(&Ap); /* get rid of extra reference */
545:     reuse = MAT_REUSE_MATRIX;
546:   }
547:   MatCreateSubMatrix(mat,isrow0,iscol0,reuse,&A);
548:   MatCreateSubMatrix(mat,isrow0,iscol1,reuse,&B);
549:   MatCreateSubMatrix(mat,isrow1,iscol0,reuse,&C);
550:   MatCreateSubMatrix(mat,isrow1,iscol1,reuse,&D);
551:   switch (mreuse) {
552:   case MAT_INITIAL_MATRIX:
553:     MatCreateSchurComplement(A,A,B,C,D,newmat);
554:     break;
555:   case MAT_REUSE_MATRIX:
556:     MatSchurComplementUpdateSubMatrices(*newmat,A,A,B,C,D);
557:     break;
558:   default:
559:     if (mreuse != MAT_IGNORE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Unrecognized value of mreuse");
560:   }
561:   if (preuse != MAT_IGNORE_MATRIX) {
562:     MatCreateSchurComplementPmat(A,B,C,D,ainvtype,preuse,newpmat);
563:   }
564:   MatDestroy(&A);
565:   MatDestroy(&B);
566:   MatDestroy(&C);
567:   MatDestroy(&D);
568:   return(0);
569: }

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

574:     Collective on Mat

576:     Input Parameters:
577: +   A      - matrix in which the complement is to be taken
578: .   isrow0 - rows to eliminate
579: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
580: .   isrow1 - rows in which the Schur complement is formed
581: .   iscol1 - columns in which the Schur complement is formed
582: .   mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in S
583: .   plump  - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
584:                        MAT_SCHUR_COMPLEMENT_AINV_DIAG or MAT_SCHUR_COMPLEMENT_AINV_LUMP
585: -   preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in Sp

587:     Output Parameters:
588: +   S      - exact Schur complement, often of type MATSCHURCOMPLEMENT which is difficult to use for preconditioning
589: -   Sp     - approximate Schur complement suitable for preconditioning

591:     Note:
592:     Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
593:     application-specific information.  The default for assembled matrices is to use the inverse of the diagonal of
594:     the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
595:     before forming inv(diag(A00)).

597:     Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
598:     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
599:     "MatGetSchurComplement_C" to their function.  If their function needs to fall back to the default implementation, it
600:     should call MatGetSchurComplement_Basic().

602:     Level: advanced

604:     Concepts: matrices^submatrices

606: .seealso: MatCreateSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement(), MatSchurComplementAinvType
607: @*/
608: PetscErrorCode  MatGetSchurComplement(Mat A,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *S,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Sp)
609: {
610:   PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*) = NULL;

621:   if (A->factortype) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
622:   f = NULL;
623:   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. */
624:     PetscObjectQueryFunction((PetscObject)*S,"MatGetSchurComplement_C",&f);
625:   }
626:   if (f) {
627:       (*f)(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,preuse,Sp);
628:   } else {
629:     MatGetSchurComplement_Basic(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,ainvtype,preuse,Sp);
630:   }
631:   return(0);
632: }

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

637:     Not collective.

639:     Input Parameters:
640: +   S        - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
641: -   ainvtype - type of approximation used to form A00inv from A00 when assembling Sp = A11 - A10 A00inv A01:
642:                       MAT_SCHUR_COMPLEMENT_AINV_DIAG or MAT_SCHUR_COMPLEMENT_AINV_LUMP

644:     Options database:
645:     -mat_schur_complement_ainv_type diag | lump

647:     Note:
648:     Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
649:     application-specific information.  The default for assembled matrices is to use the inverse of the diagonal of
650:     the (0,0) block A00 in place of A00^{-1}. This rarely produces a scalable algorithm. Optionally, A00 can be lumped
651:     before forming inv(diag(A00)).

653:     Level: advanced

655:     Concepts: matrices^submatrices

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

668:   PetscObjectGetType((PetscObject)S,&t);
669:   PetscStrcmp(t,MATSCHURCOMPLEMENT,&isschur);
670:   if (!isschur) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected Mat of type MATSCHURCOMPLEMENT, got %s instead",t);
671:   schur = (Mat_SchurComplement*)S->data;
672:   if (ainvtype != MAT_SCHUR_COMPLEMENT_AINV_DIAG && ainvtype != MAT_SCHUR_COMPLEMENT_AINV_LUMP) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatSchurComplementAinvType: %D",ainvtype);
673:   schur->ainvtype = ainvtype;
674:   return(0);
675: }

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

680:     Not collective.

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

685:     Output Parameter:
686: .   ainvtype - type of approximation used to form A00inv from A00 when assembling Sp = A11 - A10 A00inv A01:
687:                       MAT_SCHUR_COMPLEMENT_AINV_DIAG or MAT_SCHUR_COMPLEMENT_AINV_LUMP

689:     Note:
690:     Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
691:     application-specific information.  The default for assembled matrices is to use the inverse of the diagonal of
692:     the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
693:     before forming inv(diag(A00)).

695:     Level: advanced

697:     Concepts: matrices^submatrices

699: .seealso: MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementSetAinvType()
700: @*/
701: PetscErrorCode  MatSchurComplementGetAinvType(Mat S,MatSchurComplementAinvType *ainvtype)
702: {
703:   PetscErrorCode      ierr;
704:   const char*         t;
705:   PetscBool           isschur;
706:   Mat_SchurComplement *schur;

710:   PetscObjectGetType((PetscObject)S,&t);
711:   PetscStrcmp(t,MATSCHURCOMPLEMENT,&isschur);
712:   if (!isschur) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected Mat of type MATSCHURCOMPLEMENT, got %s instead",t);
713:   schur = (Mat_SchurComplement*)S->data;
714:   if (ainvtype) *ainvtype = schur->ainvtype;
715:   return(0);
716: }

718: /*@
719:     MatCreateSchurComplementPmat - create a preconditioning matrix for the Schur complement by assembling Sp = A11 - A10 inv(diag(A00)) A01

721:     Collective on Mat

723:     Input Parameters:
724: +   A00,A01,A10,A11      - the four parts of the original matrix A = [A00 A01; A10 A11] (A01,A10, and A11 are optional, implying zero matrices)
725: .   ainvtype             - type of approximation for inv(A00) used when forming Sp = A11 - A10 inv(A00) A01
726: -   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

728:     Output Parameter:
729: -   Spmat                - approximate Schur complement suitable for preconditioning S = A11 - A10 inv(diag(A00)) A01

731:     Note:
732:     Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
733:     application-specific information.  The default for assembled matrices is to use the inverse of the diagonal of
734:     the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
735:     before forming inv(diag(A00)).

737:     Level: advanced

739:     Concepts: matrices^submatrices

741: .seealso: MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementAinvType
742: @*/
743: PetscErrorCode  MatCreateSchurComplementPmat(Mat A00,Mat A01,Mat A10,Mat A11,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Spmat)
744: {

747:   PetscInt       N00;

750:   /* Use an appropriate approximate inverse of A00 to form A11 - A10 inv(diag(A00)) A01; a NULL A01, A10 or A11 indicates a zero matrix. */
751:   /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
752:   if ((!A01 || !A10) & !A11) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot assemble Spmat: A01, A10 and A11 are all NULL.");

754:   if (preuse == MAT_IGNORE_MATRIX) return(0);

756:   /* A zero size A00 or empty A01 or A10 imply S = A11. */
757:   MatGetSize(A00,&N00,NULL);
758:   if (!A01 || !A10 || !N00) {
759:     if (preuse == MAT_INITIAL_MATRIX) {
760:       MatDuplicate(A11,MAT_COPY_VALUES,Spmat);
761:     } else { /* MAT_REUSE_MATRIX */
762:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
763:       MatCopy(A11,*Spmat,DIFFERENT_NONZERO_PATTERN);
764:     }

766:   } else {
767:     Mat         AdB;
768:     Vec         diag;

770:     MatCreateVecs(A00,&diag,NULL);
771:     if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
772:       MatGetRowSum(A00,diag);
773:     } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
774:       MatGetDiagonal(A00,diag);
775:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatSchurComplementAinvType: %D", ainvtype);
776:     VecReciprocal(diag);
777:     MatDuplicate(A01,MAT_COPY_VALUES,&AdB);
778:     MatDiagonalScale(AdB,diag,NULL);
779:     VecDestroy(&diag);
780:     /* Cannot really reuse Spmat in MatMatMult() because of MatAYPX() -->
781:          MatAXPY() --> MatHeaderReplace() --> MatDestroy_XXX_MatMatMult()  */
782:     MatDestroy(Spmat);
783:     MatMatMult(A10,AdB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,Spmat);
784:     if (!A11) {
785:       MatScale(*Spmat,-1.0);
786:     } else {
787:       /* TODO: when can we pass SAME_NONZERO_PATTERN? */
788:       MatAYPX(*Spmat,-1,A11,DIFFERENT_NONZERO_PATTERN);
789:     }
790:     MatDestroy(&AdB);
791:   }
792:   return(0);
793: }

795: PetscErrorCode  MatSchurComplementGetPmat_Basic(Mat S,MatReuse preuse,Mat *Spmat)
796: {
797:   Mat A,B,C,D;
798:   Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;
799:   PetscErrorCode      ierr;

802:   if (preuse == MAT_IGNORE_MATRIX) return(0);

804:   MatSchurComplementGetSubMatrices(S,&A,NULL,&B,&C,&D);
805:   if (!A) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Schur complement component matrices unset");
806:   MatCreateSchurComplementPmat(A,B,C,D,schur->ainvtype,preuse,Spmat);
807:   return(0);
808: }

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

813:     Collective on Mat

815:     Input Parameters:
816: +   S      - matrix obtained with MatCreateSchurComplement() (or equivalent) and implementing the action of A11 - A10 ksp(A00,Ap00) A01
817: -   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

819:     Output Parameter:
820: -   Sp     - approximate Schur complement suitable for preconditioning S = A11 - A10 inv(diag(A00)) A01

822:     Note:
823:     Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
824:     application-specific information.  The default for assembled matrices is to use the inverse of the diagonal of
825:     the (0,0) block A00 in place of A00^{-1}. This rarely produce a scalable algorithm. Optionally, A00 can be lumped
826:     before forming inv(diag(A00)).

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

833:     Level: advanced

835:     Concepts: matrices^submatrices

837: .seealso: MatCreateSubMatrix(), PCFIELDSPLIT, MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementSetAinvType()
838: @*/
839: PetscErrorCode  MatSchurComplementGetPmat(Mat S,MatReuse preuse,Mat *Sp)
840: {
841:   PetscErrorCode ierr,(*f)(Mat,MatReuse,Mat*);

847:   if (S->factortype) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

849:   PetscObjectQueryFunction((PetscObject)S,"MatSchurComplementGetPmat_C",&f);
850:   if (f) {
851:     (*f)(S,preuse,Sp);
852:   } else {
853:     MatSchurComplementGetPmat_Basic(S,preuse,Sp);
854:   }
855:   return(0);
856: }

858: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
859: {
860:   PetscErrorCode      ierr;
861:   Mat_SchurComplement *Na;

864:   PetscNewLog(N,&Na);
865:   N->data = (void*) Na;

867:   N->ops->destroy        = MatDestroy_SchurComplement;
868:   N->ops->getvecs        = MatCreateVecs_SchurComplement;
869:   N->ops->view           = MatView_SchurComplement;
870:   N->ops->mult           = MatMult_SchurComplement;
871:   N->ops->multtranspose  = MatMultTranspose_SchurComplement;
872:   N->ops->multadd        = MatMultAdd_SchurComplement;
873:   N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
874:   N->assembled           = PETSC_FALSE;
875:   N->preallocated        = PETSC_FALSE;

877:   KSPCreate(PetscObjectComm((PetscObject)N),&Na->ksp);
878:   PetscObjectChangeTypeName((PetscObject)N,MATSCHURCOMPLEMENT);
879:   return(0);
880: }

882: static PetscBool KSPMatRegisterAllCalled;

884: /*@C
885:   KSPMatRegisterAll - Registers all matrix implementations in the KSP package.

887:   Not Collective

889:   Level: advanced

891: .keywords: Mat, KSP, register, all

893: .seealso: MatRegisterAll(),  KSPInitializePackage()
894: @*/
895: PetscErrorCode KSPMatRegisterAll(void)
896: {

900:   if (KSPMatRegisterAllCalled) return(0);
901:   KSPMatRegisterAllCalled = PETSC_TRUE;
902:   MatRegister(MATSCHURCOMPLEMENT,MatCreate_SchurComplement);
903:   return(0);
904: }