Actual source code: schurm.c

  1: #define PETSCMAT_DLL

 3:  #include private/matimpl.h
 4:  #include petscksp.h

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

 12: /*
 13:            D - C inv(A) B 
 14: */
 17: PetscErrorCode MatMult_SchurComplement(Mat N,Vec x,Vec y)
 18: {
 19:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 20:   PetscErrorCode       ierr;

 23:   if (!Na->work1) {MatGetVecs(Na->A,&Na->work1,PETSC_NULL);}
 24:   if (!Na->work2) {MatGetVecs(Na->A,&Na->work2,PETSC_NULL);}
 25:   MatMult(Na->B,x,Na->work1);
 26:   KSPSolve(Na->ksp,Na->work1,Na->work2);
 27:   MatMult(Na->C,Na->work2,y);
 28:   VecScale(y,-1.0);
 29:   if (Na->D) {
 30:     MatMultAdd(Na->D,x,y,y);
 31:   }
 32:   return(0);
 33: }

 37: PetscErrorCode MatSetFromOptions_SchurComplement(Mat N)
 38: {
 39:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 40:   PetscErrorCode       ierr;

 43:   KSPSetFromOptions(Na->ksp);
 44:   return(0);
 45: }
 46: 
 49: PetscErrorCode MatDestroy_SchurComplement(Mat N)
 50: {
 51:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;
 52:   PetscErrorCode       ierr;

 55:   if (Na->A)  {MatDestroy(Na->A);}
 56:   if (Na->Ap) {MatDestroy(Na->Ap);}
 57:   if (Na->B)  {MatDestroy(Na->B);}
 58:   if (Na->C)  {MatDestroy(Na->C);}
 59:   if (Na->D)  {MatDestroy(Na->D);}
 60:   if (Na->work1) {VecDestroy(Na->work1);}
 61:   if (Na->work2) {VecDestroy(Na->work2);}
 62:   KSPDestroy(Na->ksp);
 63:   PetscFree(Na);
 64:   return(0);
 65: }

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

 72:    Collective on Mat

 74:    Input Parameter:
 75: .   A,B,C,D  - the four parts of the original matrix (D is optional)

 77:    Output Parameter:
 78: .   N - the matrix that the Schur complement D - C inv(A) B

 80:    Level: intermediate

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

 86:           All four matrices must have the same MPI communicator

 88:           A and  D must be square matrices

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

 92: @*/
 93: PetscErrorCode  MatCreateSchurComplement(Mat A,Mat Ap,Mat B,Mat C,Mat D,Mat *N)
 94: {
 95:   PetscErrorCode       ierr;
 96:   PetscInt             m,n;
 97:   Mat_SchurComplement  *Na;

107:   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local columns %D",A->rmap->n,A->cmap->n);
108:   if (A->rmap->n != Ap->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local rows of Ap %D",A->rmap->n,Ap->rmap->n);
109:   if (Ap->rmap->n != Ap->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of Ap %D do not equal local columns %D",Ap->rmap->n,Ap->cmap->n);
110:   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local columns of A %D do not equal local rows of B %D",A->cmap->n,B->rmap->n);
111:   if (C->cmap->n != A->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local columns of C %D do not equal local rows of A %D",C->cmap->n,A->rmap->n);
112:   if (D) {
115:     if (D->rmap->n != D->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of D %D do not equal local columns %D",D->rmap->n,D->cmap->n);
116:     if (C->rmap->n != D->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of C %D do not equal local rows D %D",C->rmap->n,D->rmap->n);
117:   }

119:   MatGetLocalSize(B,PETSC_NULL,&n);
120:   MatGetLocalSize(C,&m,PETSC_NULL);
121:   MatCreate(((PetscObject)A)->comm,N);
122:   MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);
123:   PetscObjectChangeTypeName((PetscObject)*N,MATSCHURCOMPLEMENT);
124: 
125:   PetscNewLog(*N,Mat_SchurComplement,&Na);
126:   (*N)->data = (void*) Na;
127:   PetscObjectReference((PetscObject)A);
128:   PetscObjectReference((PetscObject)Ap);
129:   PetscObjectReference((PetscObject)B);
130:   PetscObjectReference((PetscObject)C);
131:   Na->A     = A;
132:   Na->Ap    = Ap;
133:   Na->B     = B;
134:   Na->C     = C;
135:   Na->D     = D;
136:   if (D) {
137:     PetscObjectReference((PetscObject)D);
138:   }

140:   (*N)->ops->destroy        = MatDestroy_SchurComplement;
141:   (*N)->ops->mult           = MatMult_SchurComplement;
142:   (*N)->ops->setfromoptions = MatSetFromOptions_SchurComplement;
143:   (*N)->assembled           = PETSC_TRUE;

145:   /* treats the new matrix as having block size of 1 which is most likely the case */
146:   PetscLayoutSetBlockSize((*N)->rmap,1);
147:   PetscLayoutSetBlockSize((*N)->cmap,1);
148:   PetscLayoutSetUp((*N)->rmap);
149:   PetscLayoutSetUp((*N)->cmap);

151:   KSPCreate(((PetscObject)A)->comm,&Na->ksp);
152:   KSPSetOptionsPrefix(Na->ksp,((PetscObject)A)->prefix);
153:   KSPAppendOptionsPrefix(Na->ksp,"fieldsplit_0_");
154:   KSPSetOperators(Na->ksp,A,Ap,SAME_NONZERO_PATTERN);
155:   return(0);
156: }

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

163:    Not Collective

165:    Input Parameter:
166: .   A - matrix created with MatCreateSchurComplement()

168:    Output Parameter:
169: .   ksp - the linear solver object

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

174:    Level: intermediate

176:    Notes: 
177: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()

179: @*/
180: PetscErrorCode  MatSchurComplementGetKSP(Mat A,KSP *ksp)
181: {
182:   Mat_SchurComplement  *Na;

187:   Na = (Mat_SchurComplement*)A->data;
188:   *ksp = Na->ksp;
189:   return(0);
190: }

194: /*@
195:       MatSchurComplementUpdate - Updates the Schur complement matrix object with new submatrices

197:    Collective on Mat

199:    Input Parameters:
200: +   N - the matrix obtained with MatCreateSchurComplement()
201: .   A,B,C,D  - the four parts of the original matrix (D is optional)
202: -   str - either SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER

204:  
205:    Level: intermediate

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

209:           A and  D must be square matrices

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

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

216: @*/
217: PetscErrorCode  MatSchurComplementUpdate(Mat N,Mat A,Mat Ap,Mat B,Mat C,Mat D,MatStructure str)
218: {
219:   PetscErrorCode       ierr;
220:   Mat_SchurComplement  *Na = (Mat_SchurComplement*)N->data;

229:   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local columns %D",A->rmap->n,A->cmap->n);
230:   if (A->rmap->n != Ap->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of A %D do not equal local rows of Ap %D",A->rmap->n,Ap->rmap->n);
231:   if (Ap->rmap->n != Ap->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of Ap %D do not equal local columns %D",Ap->rmap->n,Ap->cmap->n);
232:   if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local columns of A %D do not equal local rows of B %D",A->cmap->n,B->rmap->n);
233:   if (C->cmap->n != A->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local columns of C %D do not equal local rows of A %D",C->cmap->n,A->rmap->n);
234:   if (D) {
237:     if (D->rmap->n != D->cmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of D %D do not equal local columns %D",D->rmap->n,D->cmap->n);
238:     if (C->rmap->n != D->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local rows of C %D do not equal local rows D %D",C->rmap->n,D->rmap->n);
239:   }

241:   MatDestroy(Na->A);
242:   MatDestroy(Na->Ap);
243:   MatDestroy(Na->B);
244:   MatDestroy(Na->C);
245:   if (Na->D) {
246:     MatDestroy(Na->D);
247:   }
248:   PetscObjectReference((PetscObject)A);
249:   PetscObjectReference((PetscObject)Ap);
250:   PetscObjectReference((PetscObject)B);
251:   PetscObjectReference((PetscObject)C);
252:   Na->A     = A;
253:   Na->Ap    = Ap;
254:   Na->B     = B;
255:   Na->C     = C;
256:   Na->D     = D;
257:   if (D) {
258:     PetscObjectReference((PetscObject)D);
259:   }

261:   KSPSetOperators(Na->ksp,A,Ap,str);
262:   return(0);
263: }


268: /*@C
269:   MatSchurComplementGetSubmatrices - Get the individual submatrices in the Schur complement

271:   Collective on Mat

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

277:   Note:
278:   D is optional, and thus can be PETSC_NULL

280:   Level: intermediate

282: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdate()
283: @*/
284: PetscErrorCode  MatSchurComplementGetSubmatrices(Mat N,Mat *A,Mat *Ap,Mat *B,Mat *C,Mat *D)
285: {
286:   Mat_SchurComplement *Na = (Mat_SchurComplement *) N->data;
287:   PetscErrorCode       ierr;
288:   PetscTruth           flg;

292:   PetscTypeCompare((PetscObject)N,MATSCHURCOMPLEMENT,&flg);
293:   if (flg) {
294:     if (A)  *A  = Na->A;
295:     if (Ap) *Ap = Na->Ap;
296:     if (B)  *B  = Na->B;
297:     if (C)  *C  = Na->C;
298:     if (D)  *D  = Na->D;
299:   } else {
300:     if (A)  *A  = 0;
301:     if (Ap) *Ap = 0;
302:     if (B)  *B  = 0;
303:     if (C)  *C  = 0;
304:     if (D)  *D  = 0;
305:   }
306:   return(0);
307: }


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

315:     Collective on Mat

317:     Input Parameters:
318: +   mat - Matrix in which the complement is to be taken
319: .   isrow0 - rows to eliminate
320: .   iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
321: .   isrow1 - rows in which the Schur complement is formed
322: .   iscol1 - columns in which the Schur complement is formed
323: .   mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newmat
324: -   preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newpmat

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

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

335:     Level: advanced

337:     Concepts: matrices^submatrices

339: .seealso: MatGetSubMatrix(), PCFIELDSPLIT
340: @*/
341: PetscErrorCode  MatGetSchurComplement(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat)
342: {
343:   PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*);

354:   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");

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

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

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

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