Actual source code: schurm.c

petsc-3.3-p7 2013-05-11
  2: #include <petsc-private/matimpl.h>          /*I "petscmat.h" I*/
  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;


 14: PetscErrorCode MatView_SchurComplement(Mat N,PetscViewer viewer)
 15: {
 16:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 17:   PetscErrorCode       ierr;

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


 45: /*
 46:            A11 - A10 ksp(A00,Ap00)  A01 
 47: */
 50: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
 51: {
 52:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 53:   PetscErrorCode       ierr;

 56:   if (!Na->work1) {MatGetVecs(Na->A,&Na->work1,PETSC_NULL);}
 57:   if (!Na->work2) {MatGetVecs(Na->A,&Na->work2,PETSC_NULL);}
 58:   MatMult(Na->B,x,Na->work1);
 59:   KSPSolve(Na->ksp,Na->work1,Na->work2);
 60:   MatMult(Na->C,Na->work2,y);
 61:   VecScale(y,-1.0);
 62:   if (Na->D) {
 63:     MatMultAdd(Na->D,x,y,y);
 64:   }
 65:   return(0);
 66: }

 70: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N)
 71: {
 72:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 73:   PetscErrorCode       ierr;

 76:   KSPSetFromOptions(Na->ksp);
 77:   return(0);
 78: }
 79: 
 82: PetscErrorCode MatDestroy_SchurComplement(Mat N)
 83: {
 84:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 85:   PetscErrorCode       ierr;

 88:   MatDestroy(&Na->A);
 89:   MatDestroy(&Na->Ap);
 90:   MatDestroy(&Na->B);
 91:   MatDestroy(&Na->C);
 92:   MatDestroy(&Na->D);
 93:   VecDestroy(&Na->work1);
 94:   VecDestroy(&Na->work2);
 95:   KSPDestroy(&Na->ksp);
 96:   PetscFree(N->data);
 97:   return(0);
 98: }

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

105:    Collective on Mat

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

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

113:    Level: intermediate

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

119:           All four matrices must have the same MPI communicator

121:           A00 and  A11 must be square matrices

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

125: @*/
126: PetscErrorCode  MatCreateSchurComplement(Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11,Mat *N)
127: {
128:   PetscErrorCode       ierr;
129:   PetscInt             m,n;
130:   Mat_SchurComplement  *Na;

140:   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);
141:   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);
142:   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);
143:   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);
144:   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);
145:   if (A11) {
148:     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);
149:     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);
150:   }

152:   MatGetLocalSize(A01,PETSC_NULL,&n);
153:   MatGetLocalSize(A10,&m,PETSC_NULL);
154:   MatCreate(((PetscObject)A00)->comm,N);
155:   MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);
156:   PetscObjectChangeTypeName((PetscObject)*N,MATSCHURCOMPLEMENT);
157: 
158:   PetscNewLog(*N,Mat_SchurComplement,&Na);
159:   (*N)->data = (void*) Na;
160:   PetscObjectReference((PetscObject)A00);
161:   PetscObjectReference((PetscObject)Ap00);
162:   PetscObjectReference((PetscObject)A01);
163:   PetscObjectReference((PetscObject)A10);
164:   Na->A     = A00;
165:   Na->Ap    = Ap00;
166:   Na->B     = A01;
167:   Na->C     = A10;
168:   Na->D     = A11;
169:   if (A11) {
170:     PetscObjectReference((PetscObject)A11);
171:   }

173:   (*N)->ops->destroy        = MatDestroy_SchurComplement;
174:   (*N)->ops->view           = MatView_SchurComplement;
175:   (*N)->ops->mult           = MatMult_SchurComplement;
176:   (*N)->ops->setfromoptions = MatSetFromOptions_SchurComplement;
177:   (*N)->assembled           = PETSC_TRUE;
178:   (*N)->preallocated        = PETSC_TRUE;

180:   PetscLayoutSetUp((*N)->rmap);
181:   PetscLayoutSetUp((*N)->cmap);

183:   KSPCreate(((PetscObject)A00)->comm,&Na->ksp);
184:   KSPSetOperators(Na->ksp,A00,Ap00,SAME_NONZERO_PATTERN);
185:   return(0);
186: }

190: /*@
191:       MatSchurComplementGetKSP - Creates gets the KSP object that is used in the Schur complement matrix

193:    Not Collective

195:    Input Parameter:
196: .   A - matrix created with MatCreateSchurComplement()

198:    Output Parameter:
199: .   ksp - the linear solver object

201:    Options Database:
202: -     -fieldsplit_0_XXX sets KSP and PC options for the A block solver inside the Schur complement

204:    Level: intermediate

206:    Notes: 
207: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()

209: @*/
210: PetscErrorCode  MatSchurComplementGetKSP(Mat A,KSP *ksp)
211: {
212:   Mat_SchurComplement  *Na;

217:   Na = (Mat_SchurComplement*)A->data;
218:   *ksp = Na->ksp;
219:   return(0);
220: }

224: /*@
225:       MatSchurComplementUpdate - Updates the Schur complement matrix object with new submatrices

227:    Collective on Mat

229:    Input Parameters:
230: +   N - the matrix obtained with MatCreateSchurComplement()
231: .   A,B,C,D  - the four parts of the original matrix (D is optional)
232: -   str - either SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER

234:  
235:    Level: intermediate

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

239:           A and  D must be square matrices

241:           All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement()
242:           though they need not be the same matrices

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

246: @*/
247: PetscErrorCode  MatSchurComplementUpdate(Mat N,Mat A,Mat Ap,Mat B,Mat C,Mat D,MatStructure str)
248: {
249:   PetscErrorCode       ierr;
250:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;

259:   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);
260:   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);
261:   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);
262:   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);
263:   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);
264:   if (D) {
267:     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);
268:     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);
269:   }

271:   PetscObjectReference((PetscObject)A);
272:   PetscObjectReference((PetscObject)Ap);
273:   PetscObjectReference((PetscObject)B);
274:   PetscObjectReference((PetscObject)C);
275:   if (D) {
276:     PetscObjectReference((PetscObject)D);
277:   }

279:   MatDestroy(&Na->A);
280:   MatDestroy(&Na->Ap);
281:   MatDestroy(&Na->B);
282:   MatDestroy(&Na->C);
283:   MatDestroy(&Na->D);

285:   Na->A  = A;
286:   Na->Ap = Ap;
287:   Na->B  = B;
288:   Na->C  = C;
289:   Na->D  = D;

291:   KSPSetOperators(Na->ksp,A,Ap,str);
292:   return(0);
293: }


298: /*@C
299:   MatSchurComplementGetSubmatrices - Get the individual submatrices in the Schur complement

301:   Collective on Mat

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

307:   Note:
308:   D is optional, and thus can be PETSC_NULL

310:   Level: intermediate

312: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdate()
313: @*/
314: PetscErrorCode  MatSchurComplementGetSubmatrices(Mat N,Mat *A,Mat *Ap,Mat *B,Mat *C,Mat *D)
315: {
316:   Mat_SchurComplement *Na = (Mat_SchurComplement *) N->data;
317:   PetscErrorCode       ierr;
318:   PetscBool            flg;

322:   PetscObjectTypeCompare((PetscObject)N,MATSCHURCOMPLEMENT,&flg);
323:   if (flg) {
324:     if (A)  *A  = Na->A;
325:     if (Ap) *Ap = Na->Ap;
326:     if (B)  *B  = Na->B;
327:     if (C)  *C  = Na->C;
328:     if (D)  *D  = Na->D;
329:   } else {
330:     if (A)  *A  = 0;
331:     if (Ap) *Ap = 0;
332:     if (B)  *B  = 0;
333:     if (C)  *C  = 0;
334:     if (D)  *D  = 0;
335:   }
336:   return(0);
337: }

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

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

358:   if (mreuse != MAT_IGNORE_MATRIX) {
359:     /* Use MatSchurComplement */
360:     if (mreuse == MAT_REUSE_MATRIX) {
361:       MatSchurComplementGetSubmatrices(*newmat,&A,&Ap,&B,&C,&D);
362:       if (!A || !Ap || !B || !C) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset");
363:       if (A != Ap) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator");
364:       MatDestroy(&Ap); /* get rid of extra reference */
365:     }
366:     MatGetSubMatrix(mat,isrow0,iscol0,mreuse,&A);
367:     MatGetSubMatrix(mat,isrow0,iscol1,mreuse,&B);
368:     MatGetSubMatrix(mat,isrow1,iscol0,mreuse,&C);
369:     MatGetSubMatrix(mat,isrow1,iscol1,mreuse,&D);
370:     switch (mreuse) {
371:     case MAT_INITIAL_MATRIX:
372:       MatCreateSchurComplement(A,A,B,C,D,newmat);
373:       break;
374:     case MAT_REUSE_MATRIX:
375:       MatSchurComplementUpdate(*newmat,A,A,B,C,D,DIFFERENT_NONZERO_PATTERN);
376:       break;
377:     default:
378:       SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Unrecognized value of mreuse");
379:     }
380:   }
381:   if (preuse != MAT_IGNORE_MATRIX) {
382:     /* Use the diagonal part of A to form D - C inv(diag(A)) B */
383:     Mat Ad,AdB,S;
384:     Vec diag;
385:     PetscInt i,m,n,mstart,mend;
386:     PetscScalar *x;

388:     /* We could compose these with newpmat so that the matrices can be reused. */
389:     if (!A) {MatGetSubMatrix(mat,isrow0,iscol0,MAT_INITIAL_MATRIX,&A);}
390:     if (!B) {MatGetSubMatrix(mat,isrow0,iscol1,MAT_INITIAL_MATRIX,&B);}
391:     if (!C) {MatGetSubMatrix(mat,isrow1,iscol0,MAT_INITIAL_MATRIX,&C);}
392:     if (!D) {MatGetSubMatrix(mat,isrow1,iscol1,MAT_INITIAL_MATRIX,&D);}

394:     MatGetVecs(A,&diag,PETSC_NULL);
395:     MatGetDiagonal(A,diag);
396:     VecReciprocal(diag);
397:     MatGetLocalSize(A,&m,&n);
398:     /* 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
399:      * for column-oriented formats the columns of C can be scaled.  Would skip creating a silly diagonal matrix. */
400:     MatCreate(((PetscObject)A)->comm,&Ad);
401:     MatSetSizes(Ad,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
402:     MatSetOptionsPrefix(Ad,((PetscObject)mat)->prefix);
403:     MatAppendOptionsPrefix(Ad,"diag_");
404:     MatSetFromOptions(Ad);
405:     MatSeqAIJSetPreallocation(Ad,1,PETSC_NULL);
406:     MatMPIAIJSetPreallocation(Ad,1,PETSC_NULL,0,PETSC_NULL);
407:     MatGetOwnershipRange(Ad,&mstart,&mend);
408:     VecGetArray(diag,&x);
409:     for (i=mstart; i<mend; i++) {
410:       MatSetValue(Ad,i,i,x[i-mstart],INSERT_VALUES);
411:     }
412:     VecRestoreArray(diag,&x);
413:     MatAssemblyBegin(Ad,MAT_FINAL_ASSEMBLY);
414:     MatAssemblyEnd(Ad,MAT_FINAL_ASSEMBLY);
415:     VecDestroy(&diag);

417:     MatMatMult(Ad,B,MAT_INITIAL_MATRIX,1,&AdB);
418:     S = (preuse == MAT_REUSE_MATRIX) ? *newpmat : (Mat)0;
419:     MatMatMult(C,AdB,preuse,PETSC_DEFAULT,&S);
420:     MatAYPX(S,-1,D,DIFFERENT_NONZERO_PATTERN);
421:     *newpmat = S;
422:     MatDestroy(&Ad);
423:     MatDestroy(&AdB);
424:   }
425:   MatDestroy(&A);
426:   MatDestroy(&B);
427:   MatDestroy(&C);
428:   MatDestroy(&D);
429:   return(0);
430: }

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

437:     Collective on Mat

439:     Input Parameters:
440: +   mat - Matrix in which the complement is to be taken
441: .   isrow0 - rows to eliminate
442: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
443: .   isrow1 - rows in which the Schur complement is formed
444: .   iscol1 - columns in which the Schur complement is formed
445: .   mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newmat
446: -   preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newpmat

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

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

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

462:     Level: advanced

464:     Concepts: matrices^submatrices

466: .seealso: MatGetSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement()
467: @*/
468: PetscErrorCode  MatGetSchurComplement(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat)
469: {
470:   PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*);

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

483:   PetscObjectQueryFunction((PetscObject)mat,"MatGetSchurComplement_C",(void(**)(void))&f);
484:   if (f) {
485:     (*f)(mat,isrow0,iscol0,isrow1,iscol1,mreuse,newmat,preuse,newpmat);
486:   } else {
487:     MatGetSchurComplement_Basic(mat,isrow0,iscol0,isrow1,iscol1,mreuse,newmat,preuse,newpmat);
488:   }
489:   return(0);
490: }