Actual source code: schurm.c

petsc-3.4.2 2013-07-02
  2: #include <petsc-private/matimpl.h>
  3: #include <petscksp.h>                 /*I "petscksp.h" I*/

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

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

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

 34: PetscErrorCode MatView_SchurComplement(Mat N,PetscViewer viewer)
 35: {
 36:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 37:   PetscErrorCode      ierr;

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


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

 76:   if (!Na->work1) {MatGetVecs(Na->A,&Na->work1,NULL);}
 77:   if (!Na->work2) {MatGetVecs(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: */
 93: PetscErrorCode MatMultAdd_SchurComplement(Mat N,Vec x,Vec y,Vec z)
 94: {
 95:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
 96:   PetscErrorCode      ierr;

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

118: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N)
119: {
120:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
121:   PetscErrorCode      ierr;

124:   KSPSetFromOptions(Na->ksp);
125:   return(0);
126: }

130: PetscErrorCode MatDestroy_SchurComplement(Mat N)
131: {
132:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
133:   PetscErrorCode      ierr;

136:   MatDestroy(&Na->A);
137:   MatDestroy(&Na->Ap);
138:   MatDestroy(&Na->B);
139:   MatDestroy(&Na->C);
140:   MatDestroy(&Na->D);
141:   VecDestroy(&Na->work1);
142:   VecDestroy(&Na->work2);
143:   KSPDestroy(&Na->ksp);
144:   PetscFree(N->data);
145:   return(0);
146: }

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

153:    Collective on Mat

155:    Input Parameter:
156: .   A00,A01,A10,A11  - the four parts of the original matrix (A00 is optional)

158:    Output Parameter:
159: .   N - the matrix that the Schur complement A11 - A10 ksp(A00) A01

161:    Level: intermediate

163:    Notes: The Schur complement is NOT actually formed! Rather this
164:           object performs the matrix-vector product by using the the formula for
165:           the Schur complement and a KSP solver to approximate the action of inv(A)

167:           All four matrices must have the same MPI communicator

169:           A00 and  A11 must be square matrices

171: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdate(), MatCreateTranspose(), MatGetSchurComplement()

173: @*/
174: PetscErrorCode  MatCreateSchurComplement(Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11,Mat *N)
175: {

179:   MatCreate(((PetscObject)A00)->comm,N);
180:   MatSetType(*N,MATSCHURCOMPLEMENT);
181:   MatSchurComplementSet(*N,A00,Ap00,A01,A10,A11);
182:   return(0);
183: }

187: /*@
188:       MatSchurComplementSet - Sets the matrices that define the Schur complement

190:    Collective on Mat

192:    Input Parameter:
193: +   N - matrix obtained with MatCreate() and MatSetType(MATSCHURCOMPLEMENT);
194: -   A00,A01,A10,A11  - the four parts of the original matrix (A00 is optional)

196:    Level: intermediate

198:    Notes: The Schur complement is NOT actually formed! Rather this
199:           object performs the matrix-vector product by using the the formula for
200:           the Schur complement and a KSP solver to approximate the action of inv(A)

202:           All four matrices must have the same MPI communicator

204:           A00 and  A11 must be square matrices

206: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdate(), MatCreateTranspose(), MatGetSchurComplement()

208: @*/
209: PetscErrorCode  MatSchurComplementSet(Mat N,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
210: {
211:   PetscErrorCode      ierr;
212:   PetscInt            m,n;
213:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

216:   if (N->assembled) SETERRQ(PetscObjectComm((PetscObject)N),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementUpdate() for already used matrix");
224:   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);
225:   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);
226:   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);
227:   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);
228:   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);
229:   if (A11) {
232:     if (A11->rmap->n != A11->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A11 %D do not equal local columns %D",A11->rmap->n,A11->cmap->n);
233:     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);
234:   }

236:   MatGetLocalSize(A01,NULL,&n);
237:   MatGetLocalSize(A10,&m,NULL);
238:   MatSetSizes(N,m,n,PETSC_DECIDE,PETSC_DECIDE);
239:   PetscObjectReference((PetscObject)A00);
240:   PetscObjectReference((PetscObject)Ap00);
241:   PetscObjectReference((PetscObject)A01);
242:   PetscObjectReference((PetscObject)A10);
243:   Na->A  = A00;
244:   Na->Ap = Ap00;
245:   Na->B  = A01;
246:   Na->C  = A10;
247:   Na->D  = A11;
248:   if (A11) {
249:     PetscObjectReference((PetscObject)A11);
250:   }
251:   N->assembled    = PETSC_TRUE;
252:   N->preallocated = PETSC_TRUE;

254:   PetscLayoutSetUp((N)->rmap);
255:   PetscLayoutSetUp((N)->cmap);
256:   KSPSetOperators(Na->ksp,A00,Ap00,SAME_NONZERO_PATTERN);
257:   return(0);
258: }

262: /*@
263:   MatSchurComplementGetKSP - Gets the KSP object that is used to invert A in the Schur complement matrix S = C A^{-1} B

265:   Not Collective

267:   Input Parameter:
268: . A - matrix created with MatCreateSchurComplement()

270:   Output Parameter:
271: . ksp - the linear solver object

273:   Options Database:
274: . -fieldsplit_0_XXX sets KSP and PC options for the A block solver inside the Schur complement

276:   Level: intermediate

278: .seealso: MatSchurComplementSetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate()
279: @*/
280: PetscErrorCode MatSchurComplementGetKSP(Mat A, KSP *ksp)
281: {
282:   Mat_SchurComplement *Na;

287:   Na   = (Mat_SchurComplement*) A->data;
288:   *ksp = Na->ksp;
289:   return(0);
290: }

294: /*@
295:   MatSchurComplementSetKSP - Sets the KSP object that is used to invert A in the Schur complement matrix S = C A^{-1} B

297:   Not Collective

299:   Input Parameters:
300: + A - matrix created with MatCreateSchurComplement()
301: - ksp - the linear solver object

303:   Level: developer

305:   Developer Notes:
306:     There is likely no use case for this function.

308: .seealso: MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate(), MATSCHURCOMPLEMENT
309: @*/
310: PetscErrorCode MatSchurComplementSetKSP(Mat A, KSP ksp)
311: {
312:   Mat_SchurComplement *Na;
313:   PetscErrorCode      ierr;

318:   Na      = (Mat_SchurComplement*) A->data;
319:   PetscObjectReference((PetscObject)ksp);
320:   KSPDestroy(&Na->ksp);
321:   Na->ksp = ksp;
322:   KSPSetOperators(Na->ksp, Na->A, Na->Ap, SAME_NONZERO_PATTERN);
323:   return(0);
324: }

328: /*@
329:       MatSchurComplementUpdate - Updates the Schur complement matrix object with new submatrices

331:    Collective on Mat

333:    Input Parameters:
334: +   N - the matrix obtained with MatCreateSchurComplement()
335: .   A,B,C,D  - the four parts of the original matrix (D is optional)
336: -   str - either SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER


339:    Level: intermediate

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

343:           A and  D must be square matrices

345:           All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement() or MatSchurComplementSet()
346:           though they need not be the same matrices

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

350: @*/
351: PetscErrorCode  MatSchurComplementUpdate(Mat N,Mat A,Mat Ap,Mat B,Mat C,Mat D,MatStructure str)
352: {
353:   PetscErrorCode      ierr;
354:   Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;

357:   if (!N->assembled) SETERRQ(PetscObjectComm((PetscObject)N),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementSet() for new matrix");
364:   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local columns %D",A->rmap->n,A->cmap->n);
365:   if (A->rmap->n != Ap->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local rows of Ap %D",A->rmap->n,Ap->rmap->n);
366:   if (Ap->rmap->n != Ap->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of Ap %D do not equal local columns %D",Ap->rmap->n,Ap->cmap->n);
367:   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of A %D do not equal local rows of B %D",A->cmap->n,B->rmap->n);
368:   if (C->cmap->n != A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local columns of C %D do not equal local rows of A %D",C->cmap->n,A->rmap->n);
369:   if (D) {
372:     if (D->rmap->n != D->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of D %D do not equal local columns %D",D->rmap->n,D->cmap->n);
373:     if (C->rmap->n != D->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local rows of C %D do not equal local rows D %D",C->rmap->n,D->rmap->n);
374:   }

376:   PetscObjectReference((PetscObject)A);
377:   PetscObjectReference((PetscObject)Ap);
378:   PetscObjectReference((PetscObject)B);
379:   PetscObjectReference((PetscObject)C);
380:   if (D) {
381:     PetscObjectReference((PetscObject)D);
382:   }

384:   MatDestroy(&Na->A);
385:   MatDestroy(&Na->Ap);
386:   MatDestroy(&Na->B);
387:   MatDestroy(&Na->C);
388:   MatDestroy(&Na->D);

390:   Na->A  = A;
391:   Na->Ap = Ap;
392:   Na->B  = B;
393:   Na->C  = C;
394:   Na->D  = D;

396:   KSPSetOperators(Na->ksp,A,Ap,str);
397:   return(0);
398: }


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

406:   Collective on Mat

408:   Input Parameters:
409: + N - the matrix obtained with MatCreateSchurComplement()
410: - A,B,C,D  - the four parts of the original matrix (D is optional)

412:   Note:
413:   D is optional, and thus can be NULL

415:   Level: intermediate

417: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdate()
418: @*/
419: PetscErrorCode  MatSchurComplementGetSubmatrices(Mat N,Mat *A,Mat *Ap,Mat *B,Mat *C,Mat *D)
420: {
421:   Mat_SchurComplement *Na = (Mat_SchurComplement*) N->data;
422:   PetscErrorCode      ierr;
423:   PetscBool           flg;

427:   PetscObjectTypeCompare((PetscObject)N,MATSCHURCOMPLEMENT,&flg);
428:   if (flg) {
429:     if (A) *A = Na->A;
430:     if (Ap) *Ap = Na->Ap;
431:     if (B) *B = Na->B;
432:     if (C) *C = Na->C;
433:     if (D) *D = Na->D;
434:   } else {
435:     if (A) *A = 0;
436:     if (Ap) *Ap = 0;
437:     if (B) *B = 0;
438:     if (C) *C = 0;
439:     if (D) *D = 0;
440:   }
441:   return(0);
442: }

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

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

463:   if (mreuse != MAT_IGNORE_MATRIX) {
464:     /* Use MatSchurComplement */
465:     if (mreuse == MAT_REUSE_MATRIX) {
466:       MatSchurComplementGetSubmatrices(*newmat,&A,&Ap,&B,&C,&D);
467:       if (!A || !Ap || !B || !C) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset");
468:       if (A != Ap) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator");
469:       MatDestroy(&Ap); /* get rid of extra reference */
470:     }
471:     MatGetSubMatrix(mat,isrow0,iscol0,mreuse,&A);
472:     MatGetSubMatrix(mat,isrow0,iscol1,mreuse,&B);
473:     MatGetSubMatrix(mat,isrow1,iscol0,mreuse,&C);
474:     MatGetSubMatrix(mat,isrow1,iscol1,mreuse,&D);
475:     switch (mreuse) {
476:     case MAT_INITIAL_MATRIX:
477:       MatCreateSchurComplement(A,A,B,C,D,newmat);
478:       break;
479:     case MAT_REUSE_MATRIX:
480:       MatSchurComplementUpdate(*newmat,A,A,B,C,D,DIFFERENT_NONZERO_PATTERN);
481:       break;
482:     default:
483:       SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Unrecognized value of mreuse");
484:     }
485:   }
486:   if (preuse != MAT_IGNORE_MATRIX) {
487:     /* Use the diagonal part of A to form D - C inv(diag(A)) B */
488:     Mat         Ad,AdB,S;
489:     Vec         diag;
490:     PetscInt    i,m,n,mstart,mend;
491:     PetscScalar *x;

493:     /* We could compose these with newpmat so that the matrices can be reused. */
494:     if (!A) {MatGetSubMatrix(mat,isrow0,iscol0,MAT_INITIAL_MATRIX,&A);}
495:     if (!B) {MatGetSubMatrix(mat,isrow0,iscol1,MAT_INITIAL_MATRIX,&B);}
496:     if (!C) {MatGetSubMatrix(mat,isrow1,iscol0,MAT_INITIAL_MATRIX,&C);}
497:     if (!D) {MatGetSubMatrix(mat,isrow1,iscol1,MAT_INITIAL_MATRIX,&D);}

499:     MatGetVecs(A,&diag,NULL);
500:     MatGetDiagonal(A,diag);
501:     VecReciprocal(diag);
502:     MatGetLocalSize(A,&m,&n);
503:     /* We need to compute S = D - C inv(diag(A)) B.  For row-oriented formats, it is easy to scale the rows of B and
504:      * for column-oriented formats the columns of C can be scaled.  Would skip creating a silly diagonal matrix. */
505:     MatCreate(PetscObjectComm((PetscObject)A),&Ad);
506:     MatSetSizes(Ad,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
507:     MatSetOptionsPrefix(Ad,((PetscObject)mat)->prefix);
508:     MatAppendOptionsPrefix(Ad,"diag_");
509:     MatSetFromOptions(Ad);
510:     MatSeqAIJSetPreallocation(Ad,1,NULL);
511:     MatMPIAIJSetPreallocation(Ad,1,NULL,0,NULL);
512:     MatGetOwnershipRange(Ad,&mstart,&mend);
513:     VecGetArray(diag,&x);
514:     for (i=mstart; i<mend; i++) {
515:       MatSetValue(Ad,i,i,x[i-mstart],INSERT_VALUES);
516:     }
517:     VecRestoreArray(diag,&x);
518:     MatAssemblyBegin(Ad,MAT_FINAL_ASSEMBLY);
519:     MatAssemblyEnd(Ad,MAT_FINAL_ASSEMBLY);
520:     VecDestroy(&diag);

522:     MatMatMult(Ad,B,MAT_INITIAL_MATRIX,1,&AdB);
523:     S        = (preuse == MAT_REUSE_MATRIX) ? *newpmat : (Mat)0;
524:     MatMatMult(C,AdB,preuse,PETSC_DEFAULT,&S);
525:     MatAYPX(S,-1,D,DIFFERENT_NONZERO_PATTERN);
526:     *newpmat = S;
527:     MatDestroy(&Ad);
528:     MatDestroy(&AdB);
529:   }
530:   MatDestroy(&A);
531:   MatDestroy(&B);
532:   MatDestroy(&C);
533:   MatDestroy(&D);
534:   return(0);
535: }

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

542:     Collective on Mat

544:     Input Parameters:
545: +   mat - Matrix in which the complement is to be taken
546: .   isrow0 - rows to eliminate
547: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
548: .   isrow1 - rows in which the Schur complement is formed
549: .   iscol1 - columns in which the Schur complement is formed
550: .   mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newmat
551: -   preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newpmat

553:     Output Parameters:
554: +   newmat - exact Schur complement, often of type MATSCHURCOMPLEMENT which is difficult to use for preconditioning
555: -   newpmat - approximate Schur complement suitable for preconditioning

557:     Note:
558:     Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
559:     application-specific information.  The default for assembled matrices is to use the diagonal of the (0,0) block
560:     which will rarely produce a scalable algorithm.

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

567:     Level: advanced

569:     Concepts: matrices^submatrices

571: .seealso: MatGetSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement()
572: @*/
573: PetscErrorCode  MatGetSchurComplement(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat)
574: {
575:   PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*);

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

588:   PetscObjectQueryFunction((PetscObject)mat,"MatGetSchurComplement_C",&f);
589:   if (f) {
590:     (*f)(mat,isrow0,iscol0,isrow1,iscol1,mreuse,newmat,preuse,newpmat);
591:   } else {
592:     MatGetSchurComplement_Basic(mat,isrow0,iscol0,isrow1,iscol1,mreuse,newmat,preuse,newpmat);
593:   }
594:   return(0);
595: }

599: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
600: {
601:   PetscErrorCode      ierr;
602:   Mat_SchurComplement *Na;

605:   PetscNewLog(N,Mat_SchurComplement,&Na);
606:   N->data = (void*) Na;

608:   N->ops->destroy        = MatDestroy_SchurComplement;
609:   N->ops->getvecs        = MatGetVecs_SchurComplement;
610:   N->ops->view           = MatView_SchurComplement;
611:   N->ops->mult           = MatMult_SchurComplement;
612:   N->ops->multadd        = MatMultAdd_SchurComplement;
613:   N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
614:   N->assembled           = PETSC_FALSE;
615:   N->preallocated        = PETSC_FALSE;

617:   KSPCreate(PetscObjectComm((PetscObject)N),&Na->ksp);
618:   PetscObjectChangeTypeName((PetscObject)N,MATSCHURCOMPLEMENT);
619:   return(0);
620: }

622: static PetscBool KSPMatRegisterAllCalled;

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

629:   Not Collective

631:   Level: advanced

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

635: .seealso: MatRegisterAll(), MatRegisterDestroy(), KSPInitializePackage()
636: @*/
637: PetscErrorCode KSPMatRegisterAll()
638: {

642:   if (KSPMatRegisterAllCalled) return(0);
643:   KSPMatRegisterAllCalled = PETSC_TRUE;
644:   MatRegister(MATSCHURCOMPLEMENT,MatCreate_SchurComplement);
645:   return(0);
646: }