Actual source code: schurm.c

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

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



 17: PetscErrorCode MatCreateVecs_SchurComplement(Mat N,Vec *right,Vec *left)
 18: {
 19:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 20:   PetscErrorCode      ierr;

 23:   if (Na->D) {
 24:     MatCreateVecs(Na->D,right,left);
 25:     return(0);
 26:   }
 27:   if (right) {
 28:     MatCreateVecs(Na->B,right,NULL);
 29:   }
 30:   if (left) {
 31:     MatCreateVecs(Na->C,NULL,left);
 32:   }
 33:   return(0);
 34: }

 38: PetscErrorCode MatView_SchurComplement(Mat N,PetscViewer viewer)
 39: {
 40:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 41:   PetscErrorCode      ierr;

 44:   PetscViewerASCIIPrintf(viewer,"Schur complement A11 - A10 inv(A00) A01\n");
 45:   if (Na->D) {
 46:     PetscViewerASCIIPrintf(viewer,"A11\n");
 47:     PetscViewerASCIIPushTab(viewer);
 48:     MatView(Na->D,viewer);
 49:     PetscViewerASCIIPopTab(viewer);
 50:   } else {
 51:     PetscViewerASCIIPrintf(viewer,"A11 = 0\n");
 52:   }
 53:   PetscViewerASCIIPrintf(viewer,"A10\n");
 54:   PetscViewerASCIIPushTab(viewer);
 55:   MatView(Na->C,viewer);
 56:   PetscViewerASCIIPopTab(viewer);
 57:   PetscViewerASCIIPrintf(viewer,"KSP of A00\n");
 58:   PetscViewerASCIIPushTab(viewer);
 59:   KSPView(Na->ksp,viewer);
 60:   PetscViewerASCIIPopTab(viewer);
 61:   PetscViewerASCIIPrintf(viewer,"A01\n");
 62:   PetscViewerASCIIPushTab(viewer);
 63:   MatView(Na->B,viewer);
 64:   PetscViewerASCIIPopTab(viewer);
 65:   return(0);
 66: }

 68: /*
 69:            A11^T - A01^T ksptrans(A00,Ap00) A10^T
 70: */
 73: PetscErrorCode MatMultTranspose_SchurComplement(Mat N,Vec x,Vec y)
 74: {
 75:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 76:   PetscErrorCode      ierr;

 79:   if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
 80:   if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
 81:   MatMultTranspose(Na->C,x,Na->work1);
 82:   KSPSolveTranspose(Na->ksp,Na->work1,Na->work2);
 83:   MatMultTranspose(Na->B,Na->work2,y);
 84:   VecScale(y,-1.0);
 85:   if (Na->D) {
 86:     MatMultTransposeAdd(Na->D,x,y,y);
 87:   }
 88:   return(0);
 89: }

 91: /*
 92:            A11 - A10 ksp(A00,Ap00) A01
 93: */
 96: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
 97: {
 98:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 99:   PetscErrorCode      ierr;

102:   if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
103:   if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
104:   MatMult(Na->B,x,Na->work1);
105:   KSPSolve(Na->ksp,Na->work1,Na->work2);
106:   MatMult(Na->C,Na->work2,y);
107:   VecScale(y,-1.0);
108:   if (Na->D) {
109:     MatMultAdd(Na->D,x,y,y);
110:   }
111:   return(0);
112: }

114: /*
115:            A11 - A10 ksp(A00,Ap00) A01
116: */
119: PetscErrorCode MatMultAdd_SchurComplement(Mat N,Vec x,Vec y,Vec z)
120: {
121:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
122:   PetscErrorCode      ierr;

125:   if (!Na->work1) {MatCreateVecs(Na->A,&Na->work1,NULL);}
126:   if (!Na->work2) {MatCreateVecs(Na->A,&Na->work2,NULL);}
127:   MatMult(Na->B,x,Na->work1);
128:   KSPSolve(Na->ksp,Na->work1,Na->work2);
129:   if (y == z) {
130:     VecScale(Na->work2,-1.0);
131:     MatMultAdd(Na->C,Na->work2,z,z);
132:   } else {
133:     MatMult(Na->C,Na->work2,z);
134:     VecAYPX(z,-1.0,y);
135:   }
136:   if (Na->D) {
137:     MatMultAdd(Na->D,x,z,z);
138:   }
139:   return(0);
140: }

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

150:   PetscOptionsHead("MatSchurComplementOptions");
151:   Na->ainvtype = MAT_SCHUR_COMPLEMENT_AINV_DIAG;
152:   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);
153:   PetscOptionsTail();
154:   KSPSetFromOptions(Na->ksp);
155:   return(0);
156: }

160: PetscErrorCode MatDestroy_SchurComplement(Mat N)
161: {
162:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
163:   PetscErrorCode      ierr;

166:   MatDestroy(&Na->A);
167:   MatDestroy(&Na->Ap);
168:   MatDestroy(&Na->B);
169:   MatDestroy(&Na->C);
170:   MatDestroy(&Na->D);
171:   VecDestroy(&Na->work1);
172:   VecDestroy(&Na->work2);
173:   KSPDestroy(&Na->ksp);
174:   PetscFree(N->data);
175:   return(0);
176: }

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

183:    Collective on Mat

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

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

192:    Level: intermediate

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

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

200:           A00 and  A11 must be square matrices.

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

204: @*/
205: PetscErrorCode  MatCreateSchurComplement(Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11,Mat *S)
206: {

210:   KSPInitializePackage();
211:   MatCreate(((PetscObject)A00)->comm,S);
212:   MatSetType(*S,MATSCHURCOMPLEMENT);
213:   MatSchurComplementSetSubMatrices(*S,A00,Ap00,A01,A10,A11);
214:   return(0);
215: }

219: /*@
220:       MatSchurComplementSetSubMatrices - Sets the matrices that define the Schur complement

222:    Collective on Mat

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

229:    Level: intermediate

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

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

237:           A00 and  A11 must be square matrices.

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

241: @*/
242: PetscErrorCode  MatSchurComplementSetSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
243: {
244:   PetscErrorCode      ierr;
245:   PetscInt            m,n;
246:   Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;

249:   if (S->assembled) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementUpdateSubMatrices() for already used matrix");
257:   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);
258:   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);
259:   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);
260:   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);
261:   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);
262:   if (A11) {
265:     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);
266:   }

268:   MatGetLocalSize(A01,NULL,&n);
269:   MatGetLocalSize(A10,&m,NULL);
270:   MatSetSizes(S,m,n,PETSC_DECIDE,PETSC_DECIDE);
271:   PetscObjectReference((PetscObject)A00);
272:   PetscObjectReference((PetscObject)Ap00);
273:   PetscObjectReference((PetscObject)A01);
274:   PetscObjectReference((PetscObject)A10);
275:   Na->A  = A00;
276:   Na->Ap = Ap00;
277:   Na->B  = A01;
278:   Na->C  = A10;
279:   Na->D  = A11;
280:   if (A11) {
281:     PetscObjectReference((PetscObject)A11);
282:   }
283:   S->assembled    = PETSC_TRUE;
284:   S->preallocated = PETSC_TRUE;

286:   PetscLayoutSetUp((S)->rmap);
287:   PetscLayoutSetUp((S)->cmap);
288:   KSPSetOperators(Na->ksp,A00,Ap00);
289:   return(0);
290: }

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

297:   Not Collective

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

302:   Output Parameter:
303: . ksp - the linear solver object

305:   Options Database:
306: . -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.

308:   Level: intermediate

310: .seealso: MatSchurComplementSetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate()
311: @*/
312: PetscErrorCode MatSchurComplementGetKSP(Mat S, KSP *ksp)
313: {
314:   Mat_SchurComplement *Na;

319:   Na   = (Mat_SchurComplement*) S->data;
320:   *ksp = Na->ksp;
321:   return(0);
322: }

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

329:   Not Collective

331:   Input Parameters:
332: + S   - matrix created with MatCreateSchurComplement()
333: - ksp - the linear solver object

335:   Level: developer

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

340: .seealso: MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate(), MATSCHURCOMPLEMENT
341: @*/
342: PetscErrorCode MatSchurComplementSetKSP(Mat S, KSP ksp)
343: {
344:   Mat_SchurComplement *Na;
345:   PetscErrorCode      ierr;

350:   Na      = (Mat_SchurComplement*) S->data;
351:   PetscObjectReference((PetscObject)ksp);
352:   KSPDestroy(&Na->ksp);
353:   Na->ksp = ksp;
354:   KSPSetOperators(Na->ksp, Na->A, Na->Ap);
355:   return(0);
356: }

360: /*@
361:       MatSchurComplementUpdateSubMatrices - Updates the Schur complement matrix object with new submatrices

363:    Collective on Mat

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

370:    Level: intermediate

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

374:           A00 and  A11 must be square matrices

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

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

381: @*/
382: PetscErrorCode  MatSchurComplementUpdateSubMatrices(Mat S,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
383: {
384:   PetscErrorCode      ierr;
385:   Mat_SchurComplement *Na = (Mat_SchurComplement*)S->data;

388:   if (!S->assembled) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementSetSubMatrices() for a new matrix");
395:   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);
396:   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);
397:   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);
398:   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);
399:   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);
400:   if (A11) {
403:     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);
404:   }

406:   PetscObjectReference((PetscObject)A00);
407:   PetscObjectReference((PetscObject)Ap00);
408:   PetscObjectReference((PetscObject)A01);
409:   PetscObjectReference((PetscObject)A10);
410:   if (A11) {
411:     PetscObjectReference((PetscObject)A11);
412:   }

414:   MatDestroy(&Na->A);
415:   MatDestroy(&Na->Ap);
416:   MatDestroy(&Na->B);
417:   MatDestroy(&Na->C);
418:   MatDestroy(&Na->D);

420:   Na->A  = A00;
421:   Na->Ap = Ap00;
422:   Na->B  = A01;
423:   Na->C  = A10;
424:   Na->D  = A11;

426:   KSPSetOperators(Na->ksp,A00,Ap00);
427:   return(0);
428: }


433: /*@C
434:   MatSchurComplementGetSubMatrices - Get the individual submatrices in the Schur complement

436:   Collective on Mat

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

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

445:   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.

447:   Level: intermediate

449: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdateSubMatrices()
450: @*/
451: PetscErrorCode  MatSchurComplementGetSubMatrices(Mat S,Mat *A00,Mat *Ap00,Mat *A01,Mat *A10,Mat *A11)
452: {
453:   Mat_SchurComplement *Na = (Mat_SchurComplement*) S->data;
454:   PetscErrorCode      ierr;
455:   PetscBool           flg;

459:   PetscObjectTypeCompare((PetscObject)S,MATSCHURCOMPLEMENT,&flg);
460:   if (flg) {
461:     if (A00) *A00 = Na->A;
462:     if (Ap00) *Ap00 = Na->Ap;
463:     if (A01) *A01 = Na->B;
464:     if (A10) *A10 = Na->C;
465:     if (A11) *A11 = Na->D;
466:   } else {
467:     if (A00) *A00 = 0;
468:     if (Ap00) *Ap00 = 0;
469:     if (A01) *A01 = 0;
470:     if (A10) *A10 = 0;
471:     if (A11) *A11 = 0;
472:   }
473:   return(0);
474: }

476: #include <petsc-private/kspimpl.h>

480: /*@
481:   MatSchurComplementComputeExplicitOperator - Compute the Schur complement matrix explicitly

483:   Collective on Mat

485:   Input Parameter:
486: . M - the matrix obtained with MatCreateSchurComplement()

488:   Output Parameter:
489: . S - the Schur complement matrix

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

493:   Level: advanced

495: .seealso: MatCreateSchurComplement(), MatSchurComplementUpdate()
496: @*/
497: PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat M, Mat *S)
498: {
499:   Mat            B, C, D;
500:   KSP            ksp;
501:   PC             pc;
502:   PetscBool      isLU, isILU;
503:   PetscReal      fill = 2.0;

507:   MatSchurComplementGetSubMatrices(M, NULL, NULL, &B, &C, &D);
508:   MatSchurComplementGetKSP(M, &ksp);
509:   KSPGetPC(ksp, &pc);
510:   PetscObjectTypeCompare((PetscObject) pc, PCLU, &isLU);
511:   PetscObjectTypeCompare((PetscObject) pc, PCILU, &isILU);
512:   if (isLU || isILU) {
513:     Mat       fact, Bd, AinvB, AinvBd;
514:     PetscReal eps = 1.0e-10;

516:     /* This can be sped up for banded LU */
517:     KSPSetUp(ksp);
518:     PCFactorGetMatrix(pc, &fact);
519:     MatConvert(B, MATDENSE, MAT_INITIAL_MATRIX, &Bd);
520:     MatDuplicate(Bd, MAT_DO_NOT_COPY_VALUES, &AinvBd);
521:     MatMatSolve(fact, Bd, AinvBd);
522:     MatDestroy(&Bd);
523:     MatChop(AinvBd, eps);
524:     MatConvert(AinvBd, MATAIJ, MAT_INITIAL_MATRIX, &AinvB);
525:     MatDestroy(&AinvBd);
526:     MatMatMult(C, AinvB, MAT_INITIAL_MATRIX, fill, S);
527:     MatDestroy(&AinvB);
528:   } else {
529:     Mat Ainvd, Ainv;

531:     PCComputeExplicitOperator(pc, &Ainvd);
532:     MatConvert(Ainvd, MATAIJ, MAT_INITIAL_MATRIX, &Ainv);
533:     MatDestroy(&Ainvd);
534: #if 0
535:     /* Symmetric version */
536:     MatPtAP(Ainv, B, MAT_INITIAL_MATRIX, fill, S);
537: #else
538:     /* Nonsymmetric version */
539:     MatMatMatMult(C, Ainv, B, MAT_INITIAL_MATRIX, fill, S);
540: #endif
541:     MatDestroy(&Ainv);
542:   }

544:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO);
545:   MatView(*S, PETSC_VIEWER_STDOUT_WORLD);
546:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);

548:   if (D) {
549:     MatInfo info;

551:     MatGetInfo(D, MAT_GLOBAL_SUM, &info);
552:     if (info.nz_used) SETERRQ(PetscObjectComm((PetscObject) M), PETSC_ERR_SUP, "Not yet implemented");
553:   }
554:   return(0);
555: }

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

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

576:   if (mreuse != MAT_IGNORE_MATRIX) {
577:     /* Use MatSchurComplement */
578:     if (mreuse == MAT_REUSE_MATRIX) {
579:       MatSchurComplementGetSubMatrices(*newmat,&A,&Ap,&B,&C,&D);
580:       if (!A || !Ap || !B || !C) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset");
581:       if (A != Ap) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator");
582:       MatDestroy(&Ap); /* get rid of extra reference */
583:     }
584:     MatGetSubMatrix(mat,isrow0,iscol0,mreuse,&A);
585:     MatGetSubMatrix(mat,isrow0,iscol1,mreuse,&B);
586:     MatGetSubMatrix(mat,isrow1,iscol0,mreuse,&C);
587:     MatGetSubMatrix(mat,isrow1,iscol1,mreuse,&D);
588:     switch (mreuse) {
589:     case MAT_INITIAL_MATRIX:
590:       MatCreateSchurComplement(A,A,B,C,D,newmat);
591:       break;
592:     case MAT_REUSE_MATRIX:
593:       MatSchurComplementUpdateSubMatrices(*newmat,A,A,B,C,D);
594:       break;
595:     default:
596:       SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Unrecognized value of mreuse");
597:     }
598:   }
599:   if (preuse != MAT_IGNORE_MATRIX) {
600:     MatCreateSchurComplementPmat(A,B,C,D,ainvtype,preuse,newpmat);
601:   }
602:   MatDestroy(&A);
603:   MatDestroy(&B);
604:   MatDestroy(&C);
605:   MatDestroy(&D);
606:   return(0);
607: }

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

614:     Collective on Mat

616:     Input Parameters:
617: +   A      - matrix in which the complement is to be taken
618: .   isrow0 - rows to eliminate
619: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
620: .   isrow1 - rows in which the Schur complement is formed
621: .   iscol1 - columns in which the Schur complement is formed
622: .   mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newmat
623: .   plump  - the type of approximation used for the inverse of the (0,0) block used in forming Sp:
624:                        MAT_SCHUR_COMPLEMENT_AINV_DIAG | MAT_SCHUR_COMPLEMENT_AINV_LUMP
625: -   preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newpmat

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

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

637:     Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
638:     and column index sets.  In that case, the user should call PetscObjectComposeFunction() to set
639:     "MatNestGetSubMat_C" to their function.  If their function needs to fall back to the default implementation, it
640:     should call MatGetSchurComplement_Basic().

642:     Level: advanced

644:     Concepts: matrices^submatrices

646: .seealso: MatGetSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement(), MatSchurComplementAinvType
647: @*/
648: PetscErrorCode  MatGetSchurComplement(Mat A,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *S,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Sp)
649: {
650:   PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*);

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

663:   PetscObjectQueryFunction((PetscObject)S,"MatGetSchurComplement_C",&f);
664:   if (f) {
665:     (*f)(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,preuse,Sp);
666:   } else {
667:     MatGetSchurComplement_Basic(A,isrow0,iscol0,isrow1,iscol1,mreuse,S,ainvtype,preuse,Sp);
668:   }
669:   return(0);
670: }

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

677:     Not collective.

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

684:     Options database:
685:     -mat_schur_complement_ainv_type diag | lump

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

693:     Level: advanced

695:     Concepts: matrices^submatrices

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

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

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

722:     Not collective.

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

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

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: MatSchurComplementAinvType, MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementSetAinvType()
742: @*/
743: PetscErrorCode  MatSchurComplementGetAinvType(Mat S,MatSchurComplementAinvType *ainvtype)
744: {
745:   PetscErrorCode      ierr;
746:   const char*         t;
747:   PetscBool           isschur;
748:   Mat_SchurComplement *schur;

752:   PetscObjectGetType((PetscObject)S,&t);
753:   PetscStrcmp(t,MATSCHURCOMPLEMENT,&isschur);
754:   if (!isschur) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Expected Mat of type MATSCHURCOMPLEMENT, got %s instead",t);
755:   schur = (Mat_SchurComplement*)S->data;
756:   if (ainvtype) *ainvtype = schur->ainvtype;
757:   return(0);
758: }

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

765:     Collective on Mat

767:     Input Parameters:
768: +   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)
769: .   ainvtype             - type of approximation for inv(A00) used when forming Sp = A11 - A10 inv(A00) A01
770: -   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

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

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

781:     Level: advanced

783:     Concepts: matrices^submatrices

785: .seealso: MatCreateSchurComplement(), MatGetSchurComplement(), MatSchurComplementGetPmat(), MatSchurComplementAinvType
786: @*/
787: PetscErrorCode  MatCreateSchurComplementPmat(Mat A00,Mat A01,Mat A10,Mat A11,MatSchurComplementAinvType ainvtype,MatReuse preuse,Mat *Spmat)
788: {

791:   PetscInt       N00;

794:   /* 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. */
795:   /* TODO: Perhaps should create an appropriately-sized zero matrix of the same type as A00? */
796:   if ((!A01 || !A10) & !A11) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Cannot assemble Spmat: A01, A10 and A11 are all NULL.");

798:   /* A zero size A00 or empty A01 or A10 imply S = A11. */
799:   MatGetSize(A00,&N00,NULL);
800:   if (!A01 || !A10 || !N00) {
801:     if (!preuse) {
802:       MatDuplicate(A11,MAT_COPY_VALUES,Spmat);
803:     } else {
804:       MatCopy(A11,*Spmat,DIFFERENT_NONZERO_PATTERN);
805:     }

807:   } else {
808:     Mat         AdB,Sp;
809:     Vec         diag;

811:     MatCreateVecs(A00,&diag,NULL);
812:     if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_LUMP) {
813:       MatGetRowSum(A00,diag);
814:     } else if (ainvtype == MAT_SCHUR_COMPLEMENT_AINV_DIAG) {
815:       MatGetDiagonal(A00,diag);
816:     } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatSchurComplementAinvType: %D", ainvtype);
817:     VecReciprocal(diag);
818:     MatDuplicate(A01,MAT_COPY_VALUES,&AdB);
819:     MatDiagonalScale(AdB,diag,NULL);
820:     VecDestroy(&diag);
821:     Sp       = (preuse == MAT_REUSE_MATRIX) ? *Spmat : (Mat)0;
822:     MatMatMult(A10,AdB,preuse,PETSC_DEFAULT,&Sp);
823:     if (!A11) {
824:       MatScale(Sp,-1.0);
825:     } else {
826:       MatAYPX(Sp,-1,A11,DIFFERENT_NONZERO_PATTERN);
827:     }
828:     *Spmat    = Sp;
829:     MatDestroy(&AdB);
830:   }
831:   return(0);
832: }

836: PetscErrorCode  MatSchurComplementGetPmat_Basic(Mat S,MatReuse preuse,Mat *Spmat)
837: {
838:   Mat A,B,C,D;
839:   Mat_SchurComplement *schur = (Mat_SchurComplement *)S->data;
840:   PetscErrorCode      ierr;

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

845:   MatSchurComplementGetSubMatrices(S,&A,NULL,&B,&C,&D);
846:   if (!A) SETERRQ(PetscObjectComm((PetscObject)S),PETSC_ERR_ARG_WRONGSTATE,"Schur complement component matrices unset");
847:   MatCreateSchurComplementPmat(A,B,C,D,schur->ainvtype,preuse,Spmat);
848:   return(0);
849: }

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

856:     Collective on Mat

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

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

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

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

876:     Level: advanced

878:     Concepts: matrices^submatrices

880: .seealso: MatGetSubMatrix(), PCFIELDSPLIT, MatGetSchurComplement(), MatCreateSchurComplement(), MatSchurComplementSetAinvType()
881: @*/
882: PetscErrorCode  MatSchurComplementGetPmat(Mat S,MatReuse preuse,Mat *Sp)
883: {
884:   PetscErrorCode ierr,(*f)(Mat,MatReuse,Mat*);

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

892:   PetscObjectQueryFunction((PetscObject)S,"MatSchurComplementGetPmat_C",&f);
893:   if (f) {
894:     (*f)(S,preuse,Sp);
895:   } else {
896:     MatSchurComplementGetPmat_Basic(S,preuse,Sp);
897:   }
898:   return(0);
899: }

903: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
904: {
905:   PetscErrorCode      ierr;
906:   Mat_SchurComplement *Na;

909:   PetscNewLog(N,&Na);
910:   N->data = (void*) Na;

912:   N->ops->destroy        = MatDestroy_SchurComplement;
913:   N->ops->getvecs        = MatCreateVecs_SchurComplement;
914:   N->ops->view           = MatView_SchurComplement;
915:   N->ops->mult           = MatMult_SchurComplement;
916:   N->ops->multtranspose  = MatMultTranspose_SchurComplement;
917:   N->ops->multadd        = MatMultAdd_SchurComplement;
918:   N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
919:   N->assembled           = PETSC_FALSE;
920:   N->preallocated        = PETSC_FALSE;

922:   KSPCreate(PetscObjectComm((PetscObject)N),&Na->ksp);
923:   PetscObjectChangeTypeName((PetscObject)N,MATSCHURCOMPLEMENT);
924:   return(0);
925: }

927: static PetscBool KSPMatRegisterAllCalled;

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

934:   Not Collective

936:   Level: advanced

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

940: .seealso: MatRegisterAll(), MatRegisterDestroy(), KSPInitializePackage()
941: @*/
942: PetscErrorCode KSPMatRegisterAll()
943: {

947:   if (KSPMatRegisterAllCalled) return(0);
948:   KSPMatRegisterAllCalled = PETSC_TRUE;
949:   MatRegister(MATSCHURCOMPLEMENT,MatCreate_SchurComplement);
950:   return(0);
951: }