Actual source code: schurm.c
petsc-dev 2014-02-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: KSPInitializePackage();
180: MatCreate(((PetscObject)A00)->comm,N);
181: MatSetType(*N,MATSCHURCOMPLEMENT);
182: MatSchurComplementSet(*N,A00,Ap00,A01,A10,A11);
183: return(0);
184: }
188: /*@
189: MatSchurComplementSet - Sets the matrices that define the Schur complement
191: Collective on Mat
193: Input Parameter:
194: + N - matrix obtained with MatCreate() and MatSetType(MATSCHURCOMPLEMENT);
195: - A00,A01,A10,A11 - the four parts of the original matrix (A00 is optional)
197: Level: intermediate
199: Notes: The Schur complement is NOT actually formed! Rather this
200: object performs the matrix-vector product by using the the formula for
201: the Schur complement and a KSP solver to approximate the action of inv(A)
203: All four matrices must have the same MPI communicator
205: A00 and A11 must be square matrices
207: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatSchurComplementUpdate(), MatCreateTranspose(), MatGetSchurComplement()
209: @*/
210: PetscErrorCode MatSchurComplementSet(Mat N,Mat A00,Mat Ap00,Mat A01,Mat A10,Mat A11)
211: {
212: PetscErrorCode ierr;
213: PetscInt m,n;
214: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
217: if (N->assembled) SETERRQ(PetscObjectComm((PetscObject)N),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementUpdate() for already used matrix");
225: 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);
226: 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);
227: 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);
228: 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);
229: 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);
230: if (A11) {
233: 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);
234: 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);
235: }
237: MatGetLocalSize(A01,NULL,&n);
238: MatGetLocalSize(A10,&m,NULL);
239: MatSetSizes(N,m,n,PETSC_DECIDE,PETSC_DECIDE);
240: PetscObjectReference((PetscObject)A00);
241: PetscObjectReference((PetscObject)Ap00);
242: PetscObjectReference((PetscObject)A01);
243: PetscObjectReference((PetscObject)A10);
244: Na->A = A00;
245: Na->Ap = Ap00;
246: Na->B = A01;
247: Na->C = A10;
248: Na->D = A11;
249: if (A11) {
250: PetscObjectReference((PetscObject)A11);
251: }
252: N->assembled = PETSC_TRUE;
253: N->preallocated = PETSC_TRUE;
255: PetscLayoutSetUp((N)->rmap);
256: PetscLayoutSetUp((N)->cmap);
257: KSPSetOperators(Na->ksp,A00,Ap00,SAME_NONZERO_PATTERN);
258: return(0);
259: }
263: /*@
264: MatSchurComplementGetKSP - Gets the KSP object that is used to invert A in the Schur complement matrix S = C A^{-1} B
266: Not Collective
268: Input Parameter:
269: . A - matrix created with MatCreateSchurComplement()
271: Output Parameter:
272: . ksp - the linear solver object
274: Options Database:
275: . -fieldsplit_0_XXX sets KSP and PC options for the A block solver inside the Schur complement
277: Level: intermediate
279: .seealso: MatSchurComplementSetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate()
280: @*/
281: PetscErrorCode MatSchurComplementGetKSP(Mat A, KSP *ksp)
282: {
283: Mat_SchurComplement *Na;
288: Na = (Mat_SchurComplement*) A->data;
289: *ksp = Na->ksp;
290: return(0);
291: }
295: /*@
296: MatSchurComplementSetKSP - Sets the KSP object that is used to invert A in the Schur complement matrix S = C A^{-1} B
298: Not Collective
300: Input Parameters:
301: + A - matrix created with MatCreateSchurComplement()
302: - ksp - the linear solver object
304: Level: developer
306: Developer Notes:
307: There is likely no use case for this function.
309: .seealso: MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatCreateNormal(), MatMult(), MatCreate(), MATSCHURCOMPLEMENT
310: @*/
311: PetscErrorCode MatSchurComplementSetKSP(Mat A, KSP ksp)
312: {
313: Mat_SchurComplement *Na;
314: PetscErrorCode ierr;
319: Na = (Mat_SchurComplement*) A->data;
320: PetscObjectReference((PetscObject)ksp);
321: KSPDestroy(&Na->ksp);
322: Na->ksp = ksp;
323: KSPSetOperators(Na->ksp, Na->A, Na->Ap, SAME_NONZERO_PATTERN);
324: return(0);
325: }
329: /*@
330: MatSchurComplementUpdate - Updates the Schur complement matrix object with new submatrices
332: Collective on Mat
334: Input Parameters:
335: + N - the matrix obtained with MatCreateSchurComplement()
336: . A,B,C,D - the four parts of the original matrix (D is optional)
337: - str - either SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER
340: Level: intermediate
342: Notes: All four matrices must have the same MPI communicator
344: A and D must be square matrices
346: All of the matrices provided must have the same sizes as was used with MatCreateSchurComplement() or MatSchurComplementSet()
347: though they need not be the same matrices
349: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement()
351: @*/
352: PetscErrorCode MatSchurComplementUpdate(Mat N,Mat A,Mat Ap,Mat B,Mat C,Mat D,MatStructure str)
353: {
354: PetscErrorCode ierr;
355: Mat_SchurComplement *Na = (Mat_SchurComplement*)N->data;
358: if (!N->assembled) SETERRQ(PetscObjectComm((PetscObject)N),PETSC_ERR_ARG_WRONGSTATE,"Use MatSchurComplementSet() for new matrix");
365: 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);
366: 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);
367: 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);
368: 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);
369: 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);
370: if (D) {
373: 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);
374: 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);
375: }
377: PetscObjectReference((PetscObject)A);
378: PetscObjectReference((PetscObject)Ap);
379: PetscObjectReference((PetscObject)B);
380: PetscObjectReference((PetscObject)C);
381: if (D) {
382: PetscObjectReference((PetscObject)D);
383: }
385: MatDestroy(&Na->A);
386: MatDestroy(&Na->Ap);
387: MatDestroy(&Na->B);
388: MatDestroy(&Na->C);
389: MatDestroy(&Na->D);
391: Na->A = A;
392: Na->Ap = Ap;
393: Na->B = B;
394: Na->C = C;
395: Na->D = D;
397: KSPSetOperators(Na->ksp,A,Ap,str);
398: return(0);
399: }
404: /*@C
405: MatSchurComplementGetSubmatrices - Get the individual submatrices in the Schur complement
407: Collective on Mat
409: Input Parameters:
410: + N - the matrix obtained with MatCreateSchurComplement()
411: - A,B,C,D - the four parts of the original matrix (D is optional)
413: Note:
414: D is optional, and thus can be NULL
416: Level: intermediate
418: .seealso: MatCreateNormal(), MatMult(), MatCreate(), MatSchurComplementGetKSP(), MatCreateSchurComplement(), MatSchurComplementUpdate()
419: @*/
420: PetscErrorCode MatSchurComplementGetSubmatrices(Mat N,Mat *A,Mat *Ap,Mat *B,Mat *C,Mat *D)
421: {
422: Mat_SchurComplement *Na = (Mat_SchurComplement*) N->data;
423: PetscErrorCode ierr;
424: PetscBool flg;
428: PetscObjectTypeCompare((PetscObject)N,MATSCHURCOMPLEMENT,&flg);
429: if (flg) {
430: if (A) *A = Na->A;
431: if (Ap) *Ap = Na->Ap;
432: if (B) *B = Na->B;
433: if (C) *C = Na->C;
434: if (D) *D = Na->D;
435: } else {
436: if (A) *A = 0;
437: if (Ap) *Ap = 0;
438: if (B) *B = 0;
439: if (C) *C = 0;
440: if (D) *D = 0;
441: }
442: return(0);
443: }
447: /* Developer Notes: This should be implemented with a MatCreate_SchurComplement() as that is the standard design for new Mat classes. */
448: PetscErrorCode MatGetSchurComplement_Basic(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat)
449: {
451: Mat A=0,Ap=0,B=0,C=0,D=0;
462: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
464: if (mreuse != MAT_IGNORE_MATRIX) {
465: /* Use MatSchurComplement */
466: if (mreuse == MAT_REUSE_MATRIX) {
467: MatSchurComplementGetSubmatrices(*newmat,&A,&Ap,&B,&C,&D);
468: if (!A || !Ap || !B || !C) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Attempting to reuse matrix but Schur complement matrices unset");
469: if (A != Ap) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Preconditioning matrix does not match operator");
470: MatDestroy(&Ap); /* get rid of extra reference */
471: }
472: MatGetSubMatrix(mat,isrow0,iscol0,mreuse,&A);
473: MatGetSubMatrix(mat,isrow0,iscol1,mreuse,&B);
474: MatGetSubMatrix(mat,isrow1,iscol0,mreuse,&C);
475: MatGetSubMatrix(mat,isrow1,iscol1,mreuse,&D);
476: switch (mreuse) {
477: case MAT_INITIAL_MATRIX:
478: MatCreateSchurComplement(A,A,B,C,D,newmat);
479: break;
480: case MAT_REUSE_MATRIX:
481: MatSchurComplementUpdate(*newmat,A,A,B,C,D,DIFFERENT_NONZERO_PATTERN);
482: break;
483: default:
484: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Unrecognized value of mreuse");
485: }
486: }
487: if (preuse != MAT_IGNORE_MATRIX) {
488: /* Use the diagonal part of A to form D - C inv(diag(A)) B */
489: Mat Ad,AdB,S;
490: Vec diag;
491: PetscInt i,m,n,mstart,mend;
492: PetscScalar *x;
494: /* We could compose these with newpmat so that the matrices can be reused. */
495: if (!A) {MatGetSubMatrix(mat,isrow0,iscol0,MAT_INITIAL_MATRIX,&A);}
496: if (!B) {MatGetSubMatrix(mat,isrow0,iscol1,MAT_INITIAL_MATRIX,&B);}
497: if (!C) {MatGetSubMatrix(mat,isrow1,iscol0,MAT_INITIAL_MATRIX,&C);}
498: if (!D) {MatGetSubMatrix(mat,isrow1,iscol1,MAT_INITIAL_MATRIX,&D);}
500: MatGetVecs(A,&diag,NULL);
501: MatGetDiagonal(A,diag);
502: VecReciprocal(diag);
503: MatGetLocalSize(A,&m,&n);
504: /* 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
505: * for column-oriented formats the columns of C can be scaled. Would skip creating a silly diagonal matrix. */
506: MatCreate(PetscObjectComm((PetscObject)A),&Ad);
507: MatSetSizes(Ad,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
508: MatSetOptionsPrefix(Ad,((PetscObject)mat)->prefix);
509: MatAppendOptionsPrefix(Ad,"diag_");
510: MatSetFromOptions(Ad);
511: MatSeqAIJSetPreallocation(Ad,1,NULL);
512: MatMPIAIJSetPreallocation(Ad,1,NULL,0,NULL);
513: MatGetOwnershipRange(Ad,&mstart,&mend);
514: VecGetArray(diag,&x);
515: for (i=mstart; i<mend; i++) {
516: MatSetValue(Ad,i,i,x[i-mstart],INSERT_VALUES);
517: }
518: VecRestoreArray(diag,&x);
519: MatAssemblyBegin(Ad,MAT_FINAL_ASSEMBLY);
520: MatAssemblyEnd(Ad,MAT_FINAL_ASSEMBLY);
521: VecDestroy(&diag);
523: MatMatMult(Ad,B,MAT_INITIAL_MATRIX,1,&AdB);
524: S = (preuse == MAT_REUSE_MATRIX) ? *newpmat : (Mat)0;
525: MatMatMult(C,AdB,preuse,PETSC_DEFAULT,&S);
526: MatAYPX(S,-1,D,DIFFERENT_NONZERO_PATTERN);
527: *newpmat = S;
528: MatDestroy(&Ad);
529: MatDestroy(&AdB);
530: }
531: MatDestroy(&A);
532: MatDestroy(&B);
533: MatDestroy(&C);
534: MatDestroy(&D);
535: return(0);
536: }
540: /*@
541: MatGetSchurComplement - Obtain the Schur complement from eliminating part of the matrix in another part.
543: Collective on Mat
545: Input Parameters:
546: + mat - Matrix in which the complement is to be taken
547: . isrow0 - rows to eliminate
548: . iscol0 - columns to eliminate, (isrow0,iscol0) should be square and nonsingular
549: . isrow1 - rows in which the Schur complement is formed
550: . iscol1 - columns in which the Schur complement is formed
551: . mreuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newmat
552: - preuse - MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX, use MAT_IGNORE_MATRIX to put nothing in newpmat
554: Output Parameters:
555: + newmat - exact Schur complement, often of type MATSCHURCOMPLEMENT which is difficult to use for preconditioning
556: - newpmat - approximate Schur complement suitable for preconditioning
558: Note:
559: Since the real Schur complement is usually dense, providing a good approximation to newpmat usually requires
560: application-specific information. The default for assembled matrices is to use the diagonal of the (0,0) block
561: which will rarely produce a scalable algorithm.
563: Sometimes users would like to provide problem-specific data in the Schur complement, usually only for special row
564: and column index sets. In that case, the user should call PetscObjectComposeFunction() to set
565: "MatNestGetSubMat_C" to their function. If their function needs to fall back to the default implementation, it
566: should call MatGetSchurComplement_Basic().
568: Level: advanced
570: Concepts: matrices^submatrices
572: .seealso: MatGetSubMatrix(), PCFIELDSPLIT, MatCreateSchurComplement()
573: @*/
574: PetscErrorCode MatGetSchurComplement(Mat mat,IS isrow0,IS iscol0,IS isrow1,IS iscol1,MatReuse mreuse,Mat *newmat,MatReuse preuse,Mat *newpmat)
575: {
576: PetscErrorCode ierr,(*f)(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatReuse,Mat*);
587: if (mat->factortype) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
589: PetscObjectQueryFunction((PetscObject)mat,"MatGetSchurComplement_C",&f);
590: if (f) {
591: (*f)(mat,isrow0,iscol0,isrow1,iscol1,mreuse,newmat,preuse,newpmat);
592: } else {
593: MatGetSchurComplement_Basic(mat,isrow0,iscol0,isrow1,iscol1,mreuse,newmat,preuse,newpmat);
594: }
595: return(0);
596: }
600: PETSC_EXTERN PetscErrorCode MatCreate_SchurComplement(Mat N)
601: {
602: PetscErrorCode ierr;
603: Mat_SchurComplement *Na;
606: PetscNewLog(N,&Na);
607: N->data = (void*) Na;
609: N->ops->destroy = MatDestroy_SchurComplement;
610: N->ops->getvecs = MatGetVecs_SchurComplement;
611: N->ops->view = MatView_SchurComplement;
612: N->ops->mult = MatMult_SchurComplement;
613: N->ops->multadd = MatMultAdd_SchurComplement;
614: N->ops->setfromoptions = MatSetFromOptions_SchurComplement;
615: N->assembled = PETSC_FALSE;
616: N->preallocated = PETSC_FALSE;
618: KSPCreate(PetscObjectComm((PetscObject)N),&Na->ksp);
619: PetscObjectChangeTypeName((PetscObject)N,MATSCHURCOMPLEMENT);
620: return(0);
621: }
623: static PetscBool KSPMatRegisterAllCalled;
627: /*@C
628: KSPMatRegisterAll - Registers all matrix implementations in the KSP package.
630: Not Collective
632: Level: advanced
634: .keywords: Mat, KSP, register, all
636: .seealso: MatRegisterAll(), MatRegisterDestroy(), KSPInitializePackage()
637: @*/
638: PetscErrorCode KSPMatRegisterAll()
639: {
643: if (KSPMatRegisterAllCalled) return(0);
644: KSPMatRegisterAllCalled = PETSC_TRUE;
645: MatRegister(MATSCHURCOMPLEMENT,MatCreate_SchurComplement);
646: return(0);
647: }