1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 30971522cSBarry Smith /* 40971522cSBarry Smith 50971522cSBarry Smith */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 70971522cSBarry Smith 80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 90971522cSBarry Smith struct _PC_FieldSplitLink { 1069a612a9SBarry Smith KSP ksp; 110971522cSBarry Smith Vec x,y; 120971522cSBarry Smith PetscInt nfields; 130971522cSBarry Smith PetscInt *fields; 141b9fc7fcSBarry Smith VecScatter sctx; 15704ba839SBarry Smith IS is,cis; 16704ba839SBarry Smith PetscInt csize; 1751f519a2SBarry Smith PC_FieldSplitLink next,previous; 180971522cSBarry Smith }; 190971522cSBarry Smith 200971522cSBarry Smith typedef struct { 2168dd23aaSBarry Smith PCCompositeType type; 22*a4efd8eaSMatthew Knepley PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 230971522cSBarry Smith PetscInt bs; 24704ba839SBarry Smith PetscInt nsplits; 2579416396SBarry Smith Vec *x,*y,w1,w2; 2697bbdb24SBarry Smith Mat *pmat; 2768dd23aaSBarry Smith Mat *Afield; /* the rows of the matrix associated with each field */ 28704ba839SBarry Smith PetscTruth issetup; 293b224e63SBarry Smith Mat B,C,schur; /* only used when Schur complement preconditioning is used */ 303b224e63SBarry Smith KSP kspschur; 31e69d4d44SBarry Smith PetscTruth schurpre; /* preconditioner for the Schur complement is built from pmat[1] == D */ 3297bbdb24SBarry Smith PC_FieldSplitLink head; 330971522cSBarry Smith } PC_FieldSplit; 340971522cSBarry Smith 3516913363SBarry Smith /* 3616913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 3716913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 3816913363SBarry Smith PC you could change this. 3916913363SBarry Smith */ 400971522cSBarry Smith #undef __FUNCT__ 410971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 420971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 430971522cSBarry Smith { 440971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 450971522cSBarry Smith PetscErrorCode ierr; 460971522cSBarry Smith PetscTruth iascii; 470971522cSBarry Smith PetscInt i,j; 485a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 490971522cSBarry Smith 500971522cSBarry Smith PetscFunctionBegin; 510971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 520971522cSBarry Smith if (iascii) { 5351f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 5469a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 550971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 560971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 571ab39975SBarry Smith if (ilink->fields) { 580971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 5979416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 605a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 6179416396SBarry Smith if (j > 0) { 6279416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 6379416396SBarry Smith } 645a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 650971522cSBarry Smith } 660971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6779416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 681ab39975SBarry Smith } else { 691ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 701ab39975SBarry Smith } 715a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 725a9f2f41SSatish Balay ilink = ilink->next; 730971522cSBarry Smith } 740971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 750971522cSBarry Smith } else { 760971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 770971522cSBarry Smith } 780971522cSBarry Smith PetscFunctionReturn(0); 790971522cSBarry Smith } 800971522cSBarry Smith 810971522cSBarry Smith #undef __FUNCT__ 823b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 833b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 843b224e63SBarry Smith { 853b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 863b224e63SBarry Smith PetscErrorCode ierr; 873b224e63SBarry Smith PetscTruth iascii; 883b224e63SBarry Smith PetscInt i,j; 893b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 903b224e63SBarry Smith KSP ksp; 913b224e63SBarry Smith 923b224e63SBarry Smith PetscFunctionBegin; 933b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 943b224e63SBarry Smith if (iascii) { 953b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 963b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 973b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 983b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 993b224e63SBarry Smith if (ilink->fields) { 1003b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1013b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1023b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1033b224e63SBarry Smith if (j > 0) { 1043b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1053b224e63SBarry Smith } 1063b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1073b224e63SBarry Smith } 1083b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1093b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1103b224e63SBarry Smith } else { 1113b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1123b224e63SBarry Smith } 1133b224e63SBarry Smith ilink = ilink->next; 1143b224e63SBarry Smith } 1153b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1163b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1173b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1183b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 1193b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1203b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1213b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1223b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 1233b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1243b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1253b224e63SBarry Smith } else { 1263b224e63SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1273b224e63SBarry Smith } 1283b224e63SBarry Smith PetscFunctionReturn(0); 1293b224e63SBarry Smith } 1303b224e63SBarry Smith 1313b224e63SBarry Smith #undef __FUNCT__ 13269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 13369a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1340971522cSBarry Smith { 1350971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1360971522cSBarry Smith PetscErrorCode ierr; 1375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 138d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 139d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 140d32f9abdSBarry Smith char optionname[128]; 1410971522cSBarry Smith 1420971522cSBarry Smith PetscFunctionBegin; 143d32f9abdSBarry Smith if (!ilink) { 144d32f9abdSBarry Smith 145521d7252SBarry Smith if (jac->bs <= 0) { 146704ba839SBarry Smith if (pc->pmat) { 147521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 148704ba839SBarry Smith } else { 149704ba839SBarry Smith jac->bs = 1; 150704ba839SBarry Smith } 151521d7252SBarry Smith } 152d32f9abdSBarry Smith 153d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 154d32f9abdSBarry Smith if (!flg) { 155d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 156d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 157d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 158d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 159d32f9abdSBarry Smith while (PETSC_TRUE) { 160d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 161d32f9abdSBarry Smith nfields = jac->bs; 162d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 163d32f9abdSBarry Smith if (!flg2) break; 164d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 165d32f9abdSBarry Smith flg = PETSC_FALSE; 166d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 167d32f9abdSBarry Smith } 168d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 169d32f9abdSBarry Smith } 170d32f9abdSBarry Smith 171d32f9abdSBarry Smith if (flg) { 172d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 17379416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 17479416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1755a9f2f41SSatish Balay while (ilink) { 1765a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1775a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 178521d7252SBarry Smith } 1795a9f2f41SSatish Balay ilink = ilink->next; 18079416396SBarry Smith } 18197bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 18279416396SBarry Smith for (i=0; i<jac->bs; i++) { 18379416396SBarry Smith if (!fields[i]) { 18479416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 18579416396SBarry Smith } else { 18679416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 18779416396SBarry Smith } 18879416396SBarry Smith } 18968dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 190521d7252SBarry Smith } 191edf189efSBarry Smith } else if (jac->nsplits == 1) { 192edf189efSBarry Smith if (ilink->is) { 193edf189efSBarry Smith IS is2; 194edf189efSBarry Smith PetscInt nmin,nmax; 195edf189efSBarry Smith 196edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 197edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 198edf189efSBarry Smith ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr); 199edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 200edf189efSBarry Smith } else { 201edf189efSBarry Smith SETERRQ(PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 202edf189efSBarry Smith } 203d32f9abdSBarry Smith } 20469a612a9SBarry Smith PetscFunctionReturn(0); 20569a612a9SBarry Smith } 20669a612a9SBarry Smith 20769a612a9SBarry Smith 20869a612a9SBarry Smith #undef __FUNCT__ 20969a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 21069a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 21169a612a9SBarry Smith { 21269a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 21369a612a9SBarry Smith PetscErrorCode ierr; 2145a9f2f41SSatish Balay PC_FieldSplitLink ilink; 215cf502942SBarry Smith PetscInt i,nsplit,ccsize; 21669a612a9SBarry Smith MatStructure flag = pc->flag; 21768dd23aaSBarry Smith PetscTruth sorted,getsub; 21869a612a9SBarry Smith 21969a612a9SBarry Smith PetscFunctionBegin; 22069a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 22197bbdb24SBarry Smith nsplit = jac->nsplits; 2225a9f2f41SSatish Balay ilink = jac->head; 22397bbdb24SBarry Smith 22497bbdb24SBarry Smith /* get the matrices for each split */ 225704ba839SBarry Smith if (!jac->issetup) { 2261b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 22797bbdb24SBarry Smith 228704ba839SBarry Smith jac->issetup = PETSC_TRUE; 229704ba839SBarry Smith 230704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 23151f519a2SBarry Smith bs = jac->bs; 23297bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 233cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2341b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2351b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2361b9fc7fcSBarry Smith if (jac->defaultsplit) { 237704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 238704ba839SBarry Smith } else if (!ilink->is) { 239ccb205f8SBarry Smith if (ilink->nfields > 1) { 2405a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2415a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2421b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2431b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2441b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 24597bbdb24SBarry Smith } 24697bbdb24SBarry Smith } 247704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 248ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 249ccb205f8SBarry Smith } else { 250704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 251ccb205f8SBarry Smith } 2523e197d65SBarry Smith } 2533e197d65SBarry Smith ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr); 254704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 2551b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 256704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 257704ba839SBarry Smith ilink = ilink->next; 2581b9fc7fcSBarry Smith } 2591b9fc7fcSBarry Smith } 2601b9fc7fcSBarry Smith 261704ba839SBarry Smith ilink = jac->head; 26297bbdb24SBarry Smith if (!jac->pmat) { 263cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 264cf502942SBarry Smith for (i=0; i<nsplit; i++) { 265704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 266704ba839SBarry Smith ilink = ilink->next; 267cf502942SBarry Smith } 26897bbdb24SBarry Smith } else { 269cf502942SBarry Smith for (i=0; i<nsplit; i++) { 270704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 271704ba839SBarry Smith ilink = ilink->next; 272cf502942SBarry Smith } 27397bbdb24SBarry Smith } 27497bbdb24SBarry Smith 27568dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 27668dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 2773b224e63SBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 27868dd23aaSBarry Smith ilink = jac->head; 27968dd23aaSBarry Smith if (!jac->Afield) { 28068dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 28168dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 28211755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 28368dd23aaSBarry Smith ilink = ilink->next; 28468dd23aaSBarry Smith } 28568dd23aaSBarry Smith } else { 28668dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 28711755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 28868dd23aaSBarry Smith ilink = ilink->next; 28968dd23aaSBarry Smith } 29068dd23aaSBarry Smith } 29168dd23aaSBarry Smith } 29268dd23aaSBarry Smith 2933b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 2943b224e63SBarry Smith IS ccis; 295e69d4d44SBarry Smith PetscInt N,nlocal,nis; 2963b224e63SBarry Smith if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 29768dd23aaSBarry Smith 2983b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 2993b224e63SBarry Smith if (jac->schur) { 3003b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3013b224e63SBarry Smith ilink = jac->head; 302edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 303e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 304e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 305e69d4d44SBarry Smith nlocal = nlocal - nis; 306e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3073b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3083b224e63SBarry Smith ilink = ilink->next; 309edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 310e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 311e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 312e69d4d44SBarry Smith nlocal = nlocal - nis; 313e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3143b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3153b224e63SBarry Smith ierr = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 316e69d4d44SBarry Smith if (jac->schurpre) { 317e69d4d44SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],pc->flag);CHKERRQ(ierr); 318e69d4d44SBarry Smith } else { 3193b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,pc->flag);CHKERRQ(ierr); 320e69d4d44SBarry Smith } 3213b224e63SBarry Smith 3223b224e63SBarry Smith } else { 3231cee3971SBarry Smith KSP ksp; 3243b224e63SBarry Smith 3253b224e63SBarry Smith /* extract the B and C matrices */ 3263b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3273b224e63SBarry Smith ilink = jac->head; 328edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 329e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 330e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 331e69d4d44SBarry Smith nlocal = nlocal - nis; 332e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3333b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3343b224e63SBarry Smith ilink = ilink->next; 335edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 336e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 337e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 338e69d4d44SBarry Smith nlocal = nlocal - nis; 339e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3403b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3413b224e63SBarry Smith ierr = MatCreateSchurComplement(jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr); 3421cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 343e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3443b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3453b224e63SBarry Smith 3463b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3471cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 348e69d4d44SBarry Smith if (jac->schurpre) { 349e69d4d44SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 350e69d4d44SBarry Smith } else { 3513b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 352e69d4d44SBarry Smith } 3533b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 354edf189efSBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 3553b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3563b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3573b224e63SBarry Smith 3583b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3593b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3603b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3613b224e63SBarry Smith ilink = jac->head; 3623b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3633b224e63SBarry Smith ilink = ilink->next; 3643b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3653b224e63SBarry Smith } 3663b224e63SBarry Smith } else { 36797bbdb24SBarry Smith /* set up the individual PCs */ 36897bbdb24SBarry Smith i = 0; 3695a9f2f41SSatish Balay ilink = jac->head; 3705a9f2f41SSatish Balay while (ilink) { 3715a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 3723b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3735a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 3745a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 37597bbdb24SBarry Smith i++; 3765a9f2f41SSatish Balay ilink = ilink->next; 3770971522cSBarry Smith } 37897bbdb24SBarry Smith 37997bbdb24SBarry Smith /* create work vectors for each split */ 3801b9fc7fcSBarry Smith if (!jac->x) { 38197bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 3825a9f2f41SSatish Balay ilink = jac->head; 38397bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 384906ed7ccSBarry Smith Vec *vl,*vr; 385906ed7ccSBarry Smith 386906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 387906ed7ccSBarry Smith ilink->x = *vr; 388906ed7ccSBarry Smith ilink->y = *vl; 389906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 390906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 3915a9f2f41SSatish Balay jac->x[i] = ilink->x; 3925a9f2f41SSatish Balay jac->y[i] = ilink->y; 3935a9f2f41SSatish Balay ilink = ilink->next; 39497bbdb24SBarry Smith } 3953b224e63SBarry Smith } 3963b224e63SBarry Smith } 3973b224e63SBarry Smith 3983b224e63SBarry Smith 399d0f46423SBarry Smith if (!jac->head->sctx) { 4003b224e63SBarry Smith Vec xtmp; 4013b224e63SBarry Smith 40279416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4031b9fc7fcSBarry Smith 4045a9f2f41SSatish Balay ilink = jac->head; 4051b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4061b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 407704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4085a9f2f41SSatish Balay ilink = ilink->next; 4091b9fc7fcSBarry Smith } 4101b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4111b9fc7fcSBarry Smith } 4120971522cSBarry Smith PetscFunctionReturn(0); 4130971522cSBarry Smith } 4140971522cSBarry Smith 4155a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 416ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 417ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4185a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 419ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 420ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 42179416396SBarry Smith 4220971522cSBarry Smith #undef __FUNCT__ 4233b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4243b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4253b224e63SBarry Smith { 4263b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4273b224e63SBarry Smith PetscErrorCode ierr; 4283b224e63SBarry Smith KSP ksp; 4293b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4303b224e63SBarry Smith 4313b224e63SBarry Smith PetscFunctionBegin; 4323b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4333b224e63SBarry Smith 4343b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4353b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4363b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4373b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4383b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4393b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4403b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4413b224e63SBarry Smith 4423b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4433b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4443b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4453b224e63SBarry Smith 4463b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4473b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4483b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4493b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4503b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4513b224e63SBarry Smith 4523b224e63SBarry Smith PetscFunctionReturn(0); 4533b224e63SBarry Smith } 4543b224e63SBarry Smith 4553b224e63SBarry Smith #undef __FUNCT__ 4560971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4570971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4580971522cSBarry Smith { 4590971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4600971522cSBarry Smith PetscErrorCode ierr; 4615a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4623e197d65SBarry Smith PetscInt bs,cnt; 4630971522cSBarry Smith 4640971522cSBarry Smith PetscFunctionBegin; 46551f519a2SBarry Smith CHKMEMQ; 466e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 46751f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 46851f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 46951f519a2SBarry Smith 47079416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4711b9fc7fcSBarry Smith if (jac->defaultsplit) { 4720971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 4735a9f2f41SSatish Balay while (ilink) { 4745a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4755a9f2f41SSatish Balay ilink = ilink->next; 4760971522cSBarry Smith } 4770971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 4781b9fc7fcSBarry Smith } else { 479efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4805a9f2f41SSatish Balay while (ilink) { 4815a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4825a9f2f41SSatish Balay ilink = ilink->next; 4831b9fc7fcSBarry Smith } 4841b9fc7fcSBarry Smith } 48516913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 48679416396SBarry Smith if (!jac->w1) { 48779416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 48879416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 48979416396SBarry Smith } 490efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4915a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4923e197d65SBarry Smith cnt = 1; 4935a9f2f41SSatish Balay while (ilink->next) { 4945a9f2f41SSatish Balay ilink = ilink->next; 4953e197d65SBarry Smith if (jac->Afield) { 4963e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 4973e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 4983e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 4993e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5003e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5013e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5023e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5033e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5043e197d65SBarry Smith } else { 5053e197d65SBarry Smith /* compute the residual over the entire vector */ 5061ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 507efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 5085a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 50979416396SBarry Smith } 5103e197d65SBarry Smith } 51151f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 51211755939SBarry Smith cnt -= 2; 51351f519a2SBarry Smith while (ilink->previous) { 51451f519a2SBarry Smith ilink = ilink->previous; 51511755939SBarry Smith if (jac->Afield) { 51611755939SBarry Smith /* compute the residual only over the part of the vector needed */ 51711755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 51811755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 51911755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52011755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52111755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 52211755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 52311755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 52411755939SBarry Smith } else { 5251ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 52651f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 52751f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 52879416396SBarry Smith } 52951f519a2SBarry Smith } 53011755939SBarry Smith } 53116913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 53251f519a2SBarry Smith CHKMEMQ; 5330971522cSBarry Smith PetscFunctionReturn(0); 5340971522cSBarry Smith } 5350971522cSBarry Smith 536421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 537ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 538ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 539421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 540ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 541ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 542421e10b8SBarry Smith 543421e10b8SBarry Smith #undef __FUNCT__ 544421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 545421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 546421e10b8SBarry Smith { 547421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 548421e10b8SBarry Smith PetscErrorCode ierr; 549421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 550421e10b8SBarry Smith PetscInt bs; 551421e10b8SBarry Smith 552421e10b8SBarry Smith PetscFunctionBegin; 553421e10b8SBarry Smith CHKMEMQ; 554421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 555421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 556421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 557421e10b8SBarry Smith 558421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 559421e10b8SBarry Smith if (jac->defaultsplit) { 560421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 561421e10b8SBarry Smith while (ilink) { 562421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 563421e10b8SBarry Smith ilink = ilink->next; 564421e10b8SBarry Smith } 565421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 566421e10b8SBarry Smith } else { 567421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 568421e10b8SBarry Smith while (ilink) { 569421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 570421e10b8SBarry Smith ilink = ilink->next; 571421e10b8SBarry Smith } 572421e10b8SBarry Smith } 573421e10b8SBarry Smith } else { 574421e10b8SBarry Smith if (!jac->w1) { 575421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 576421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 577421e10b8SBarry Smith } 578421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 579421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 580421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 581421e10b8SBarry Smith while (ilink->next) { 582421e10b8SBarry Smith ilink = ilink->next; 5839989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 584421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 585421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 586421e10b8SBarry Smith } 587421e10b8SBarry Smith while (ilink->previous) { 588421e10b8SBarry Smith ilink = ilink->previous; 5899989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 590421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 591421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 592421e10b8SBarry Smith } 593421e10b8SBarry Smith } else { 594421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 595421e10b8SBarry Smith ilink = ilink->next; 596421e10b8SBarry Smith } 597421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 598421e10b8SBarry Smith while (ilink->previous) { 599421e10b8SBarry Smith ilink = ilink->previous; 6009989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 601421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 602421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 603421e10b8SBarry Smith } 604421e10b8SBarry Smith } 605421e10b8SBarry Smith } 606421e10b8SBarry Smith CHKMEMQ; 607421e10b8SBarry Smith PetscFunctionReturn(0); 608421e10b8SBarry Smith } 609421e10b8SBarry Smith 6100971522cSBarry Smith #undef __FUNCT__ 6110971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6120971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6130971522cSBarry Smith { 6140971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6150971522cSBarry Smith PetscErrorCode ierr; 6165a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6170971522cSBarry Smith 6180971522cSBarry Smith PetscFunctionBegin; 6195a9f2f41SSatish Balay while (ilink) { 6205a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6215a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6225a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6235a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 624704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 625704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 6265a9f2f41SSatish Balay next = ilink->next; 627704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 628704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6295a9f2f41SSatish Balay ilink = next; 6300971522cSBarry Smith } 63105b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 63297bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 63368dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 63479416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 63579416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6363b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 637d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 6383b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6393b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6400971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6410971522cSBarry Smith PetscFunctionReturn(0); 6420971522cSBarry Smith } 6430971522cSBarry Smith 6440971522cSBarry Smith #undef __FUNCT__ 6450971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6460971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6470971522cSBarry Smith { 6481b9fc7fcSBarry Smith PetscErrorCode ierr; 64951f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 650e69d4d44SBarry Smith PetscTruth flg,set; 6511b9fc7fcSBarry Smith char optionname[128]; 6529dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6533b224e63SBarry Smith PCCompositeType ctype; 6541b9fc7fcSBarry Smith 6550971522cSBarry Smith PetscFunctionBegin; 6561b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 65751f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 65851f519a2SBarry Smith if (flg) { 65951f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 66051f519a2SBarry Smith } 661704ba839SBarry Smith 6623b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6633b224e63SBarry Smith if (flg) { 6643b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6653b224e63SBarry Smith } 666d32f9abdSBarry Smith 667d32f9abdSBarry Smith if (jac->bs > 0) { 668d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 669d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 67051f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6711b9fc7fcSBarry Smith while (PETSC_TRUE) { 67213f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 67351f519a2SBarry Smith nfields = jac->bs; 6741b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6751b9fc7fcSBarry Smith if (!flg) break; 6761b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 6771b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 6781b9fc7fcSBarry Smith } 67951f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 680d32f9abdSBarry Smith } 681e69d4d44SBarry Smith ierr = PetscOptionsTruth("-pc_fieldsplit_schur_precondition","Build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",jac->schurpre,&flg,&set);CHKERRQ(ierr); 682e69d4d44SBarry Smith if (set) { 683e69d4d44SBarry Smith ierr = PCFieldSplitSchurPrecondition(pc,flg);CHKERRQ(ierr); 684e69d4d44SBarry Smith } 6851b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6860971522cSBarry Smith PetscFunctionReturn(0); 6870971522cSBarry Smith } 6880971522cSBarry Smith 6890971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 6900971522cSBarry Smith 6910971522cSBarry Smith EXTERN_C_BEGIN 6920971522cSBarry Smith #undef __FUNCT__ 6930971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 694dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 6950971522cSBarry Smith { 69697bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6970971522cSBarry Smith PetscErrorCode ierr; 6985a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 69969a612a9SBarry Smith char prefix[128]; 70051f519a2SBarry Smith PetscInt i; 7010971522cSBarry Smith 7020971522cSBarry Smith PetscFunctionBegin; 7030971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 70451f519a2SBarry Smith for (i=0; i<n; i++) { 70551f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 70651f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 70751f519a2SBarry Smith } 708704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 709704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 7105a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 7115a9f2f41SSatish Balay ilink->nfields = n; 7125a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7137adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7141cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7155a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 71669a612a9SBarry Smith 7177adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 7187adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 71969a612a9SBarry Smith } else { 72013f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 72169a612a9SBarry Smith } 7225a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7230971522cSBarry Smith 7240971522cSBarry Smith if (!next) { 7255a9f2f41SSatish Balay jac->head = ilink; 72651f519a2SBarry Smith ilink->previous = PETSC_NULL; 7270971522cSBarry Smith } else { 7280971522cSBarry Smith while (next->next) { 7290971522cSBarry Smith next = next->next; 7300971522cSBarry Smith } 7315a9f2f41SSatish Balay next->next = ilink; 73251f519a2SBarry Smith ilink->previous = next; 7330971522cSBarry Smith } 7340971522cSBarry Smith jac->nsplits++; 7350971522cSBarry Smith PetscFunctionReturn(0); 7360971522cSBarry Smith } 7370971522cSBarry Smith EXTERN_C_END 7380971522cSBarry Smith 739e69d4d44SBarry Smith EXTERN_C_BEGIN 740e69d4d44SBarry Smith #undef __FUNCT__ 741e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 742e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 743e69d4d44SBarry Smith { 744e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 745e69d4d44SBarry Smith PetscErrorCode ierr; 746e69d4d44SBarry Smith 747e69d4d44SBarry Smith PetscFunctionBegin; 748e69d4d44SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 749e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 750e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 751e69d4d44SBarry Smith PetscFunctionReturn(0); 752e69d4d44SBarry Smith } 753e69d4d44SBarry Smith EXTERN_C_END 7540971522cSBarry Smith 7550971522cSBarry Smith EXTERN_C_BEGIN 7560971522cSBarry Smith #undef __FUNCT__ 75769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 758dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7590971522cSBarry Smith { 7600971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7610971522cSBarry Smith PetscErrorCode ierr; 7620971522cSBarry Smith PetscInt cnt = 0; 7635a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7640971522cSBarry Smith 7650971522cSBarry Smith PetscFunctionBegin; 76669a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7675a9f2f41SSatish Balay while (ilink) { 7685a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7695a9f2f41SSatish Balay ilink = ilink->next; 7700971522cSBarry Smith } 7710971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7720971522cSBarry Smith *n = jac->nsplits; 7730971522cSBarry Smith PetscFunctionReturn(0); 7740971522cSBarry Smith } 7750971522cSBarry Smith EXTERN_C_END 7760971522cSBarry Smith 777704ba839SBarry Smith EXTERN_C_BEGIN 778704ba839SBarry Smith #undef __FUNCT__ 779704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 780704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 781704ba839SBarry Smith { 782704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 783704ba839SBarry Smith PetscErrorCode ierr; 784704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 785704ba839SBarry Smith char prefix[128]; 786704ba839SBarry Smith 787704ba839SBarry Smith PetscFunctionBegin; 78816913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 7891ab39975SBarry Smith ilink->is = is; 7901ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 791704ba839SBarry Smith ilink->next = PETSC_NULL; 792704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7931cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 794704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 795704ba839SBarry Smith 796704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 797704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 798704ba839SBarry Smith } else { 799704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 800704ba839SBarry Smith } 801704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 802704ba839SBarry Smith 803704ba839SBarry Smith if (!next) { 804704ba839SBarry Smith jac->head = ilink; 805704ba839SBarry Smith ilink->previous = PETSC_NULL; 806704ba839SBarry Smith } else { 807704ba839SBarry Smith while (next->next) { 808704ba839SBarry Smith next = next->next; 809704ba839SBarry Smith } 810704ba839SBarry Smith next->next = ilink; 811704ba839SBarry Smith ilink->previous = next; 812704ba839SBarry Smith } 813704ba839SBarry Smith jac->nsplits++; 814704ba839SBarry Smith 815704ba839SBarry Smith PetscFunctionReturn(0); 816704ba839SBarry Smith } 817704ba839SBarry Smith EXTERN_C_END 818704ba839SBarry Smith 8190971522cSBarry Smith #undef __FUNCT__ 8200971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8210971522cSBarry Smith /*@ 8220971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8230971522cSBarry Smith 8240971522cSBarry Smith Collective on PC 8250971522cSBarry Smith 8260971522cSBarry Smith Input Parameters: 8270971522cSBarry Smith + pc - the preconditioner context 8280971522cSBarry Smith . n - the number of fields in this split 8290971522cSBarry Smith . fields - the fields in this split 8300971522cSBarry Smith 8310971522cSBarry Smith Level: intermediate 8320971522cSBarry Smith 833d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 834d32f9abdSBarry Smith 835d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 836d32f9abdSBarry Smith size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean 837d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 838d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 839d32f9abdSBarry Smith 840d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8410971522cSBarry Smith 8420971522cSBarry Smith @*/ 843dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 8440971522cSBarry Smith { 8450971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 8460971522cSBarry Smith 8470971522cSBarry Smith PetscFunctionBegin; 8480971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8490971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8500971522cSBarry Smith if (f) { 8510971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8520971522cSBarry Smith } 8530971522cSBarry Smith PetscFunctionReturn(0); 8540971522cSBarry Smith } 8550971522cSBarry Smith 8560971522cSBarry Smith #undef __FUNCT__ 857704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 858704ba839SBarry Smith /*@ 859704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 860704ba839SBarry Smith 861704ba839SBarry Smith Collective on PC 862704ba839SBarry Smith 863704ba839SBarry Smith Input Parameters: 864704ba839SBarry Smith + pc - the preconditioner context 865704ba839SBarry Smith . is - the index set that defines the vector elements in this field 866704ba839SBarry Smith 867d32f9abdSBarry Smith 868d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 869d32f9abdSBarry Smith 870704ba839SBarry Smith Level: intermediate 871704ba839SBarry Smith 872704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 873704ba839SBarry Smith 874704ba839SBarry Smith @*/ 875704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 876704ba839SBarry Smith { 877704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 878704ba839SBarry Smith 879704ba839SBarry Smith PetscFunctionBegin; 880704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 881704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 882704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 883704ba839SBarry Smith if (f) { 884704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 885704ba839SBarry Smith } 886704ba839SBarry Smith PetscFunctionReturn(0); 887704ba839SBarry Smith } 888704ba839SBarry Smith 889704ba839SBarry Smith #undef __FUNCT__ 89051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 89151f519a2SBarry Smith /*@ 89251f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 89351f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 89451f519a2SBarry Smith 89551f519a2SBarry Smith Collective on PC 89651f519a2SBarry Smith 89751f519a2SBarry Smith Input Parameters: 89851f519a2SBarry Smith + pc - the preconditioner context 89951f519a2SBarry Smith - bs - the block size 90051f519a2SBarry Smith 90151f519a2SBarry Smith Level: intermediate 90251f519a2SBarry Smith 90351f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 90451f519a2SBarry Smith 90551f519a2SBarry Smith @*/ 90651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 90751f519a2SBarry Smith { 90851f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 90951f519a2SBarry Smith 91051f519a2SBarry Smith PetscFunctionBegin; 91151f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 91251f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 91351f519a2SBarry Smith if (f) { 91451f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 91551f519a2SBarry Smith } 91651f519a2SBarry Smith PetscFunctionReturn(0); 91751f519a2SBarry Smith } 91851f519a2SBarry Smith 91951f519a2SBarry Smith #undef __FUNCT__ 92069a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9210971522cSBarry Smith /*@C 92269a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9230971522cSBarry Smith 92469a612a9SBarry Smith Collective on KSP 9250971522cSBarry Smith 9260971522cSBarry Smith Input Parameter: 9270971522cSBarry Smith . pc - the preconditioner context 9280971522cSBarry Smith 9290971522cSBarry Smith Output Parameters: 9300971522cSBarry Smith + n - the number of split 93169a612a9SBarry Smith - pc - the array of KSP contexts 9320971522cSBarry Smith 9330971522cSBarry Smith Note: 934d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 935d32f9abdSBarry Smith (not the KSP just the array that contains them). 9360971522cSBarry Smith 93769a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9380971522cSBarry Smith 9390971522cSBarry Smith Level: advanced 9400971522cSBarry Smith 9410971522cSBarry Smith .seealso: PCFIELDSPLIT 9420971522cSBarry Smith @*/ 943dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9440971522cSBarry Smith { 94569a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9460971522cSBarry Smith 9470971522cSBarry Smith PetscFunctionBegin; 9480971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9490971522cSBarry Smith PetscValidIntPointer(n,2); 95069a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9510971522cSBarry Smith if (f) { 95269a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 9530971522cSBarry Smith } else { 95469a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9550971522cSBarry Smith } 9560971522cSBarry Smith PetscFunctionReturn(0); 9570971522cSBarry Smith } 9580971522cSBarry Smith 959e69d4d44SBarry Smith #undef __FUNCT__ 960e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 961e69d4d44SBarry Smith /*@ 962e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 963e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 964e69d4d44SBarry Smith 965e69d4d44SBarry Smith Collective on PC 966e69d4d44SBarry Smith 967e69d4d44SBarry Smith Input Parameters: 968e69d4d44SBarry Smith + pc - the preconditioner context 969e69d4d44SBarry Smith - flg - build the preconditioner 970e69d4d44SBarry Smith 971e69d4d44SBarry Smith Options Database: 972e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 973e69d4d44SBarry Smith 974e69d4d44SBarry Smith Level: intermediate 975e69d4d44SBarry Smith 976e69d4d44SBarry Smith Notes: What should we do if the user wants to provide a different matrix (like a mass matrix) to use to build the preconditioner 977e69d4d44SBarry Smith 978e69d4d44SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFIELDSPLIT 979e69d4d44SBarry Smith 980e69d4d44SBarry Smith @*/ 981e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PetscTruth flg) 982e69d4d44SBarry Smith { 983e69d4d44SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 984e69d4d44SBarry Smith 985e69d4d44SBarry Smith PetscFunctionBegin; 986e69d4d44SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 987e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 988e69d4d44SBarry Smith if (f) { 989e69d4d44SBarry Smith ierr = (*f)(pc,flg);CHKERRQ(ierr); 990e69d4d44SBarry Smith } 991e69d4d44SBarry Smith PetscFunctionReturn(0); 992e69d4d44SBarry Smith } 993e69d4d44SBarry Smith 994e69d4d44SBarry Smith EXTERN_C_BEGIN 995e69d4d44SBarry Smith #undef __FUNCT__ 996e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 997e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PetscTruth flg) 998e69d4d44SBarry Smith { 999e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1000e69d4d44SBarry Smith 1001e69d4d44SBarry Smith PetscFunctionBegin; 1002e69d4d44SBarry Smith jac->schurpre = flg; 1003e69d4d44SBarry Smith PetscFunctionReturn(0); 1004e69d4d44SBarry Smith } 1005e69d4d44SBarry Smith EXTERN_C_END 1006e69d4d44SBarry Smith 100779416396SBarry Smith EXTERN_C_BEGIN 100879416396SBarry Smith #undef __FUNCT__ 100979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1010dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 101179416396SBarry Smith { 101279416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1013e69d4d44SBarry Smith PetscErrorCode ierr; 101479416396SBarry Smith 101579416396SBarry Smith PetscFunctionBegin; 101679416396SBarry Smith jac->type = type; 10173b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10183b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10193b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1020e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1021e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1022e69d4d44SBarry Smith 10233b224e63SBarry Smith } else { 10243b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10253b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1026e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 10279e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 10283b224e63SBarry Smith } 102979416396SBarry Smith PetscFunctionReturn(0); 103079416396SBarry Smith } 103179416396SBarry Smith EXTERN_C_END 103279416396SBarry Smith 103351f519a2SBarry Smith EXTERN_C_BEGIN 103451f519a2SBarry Smith #undef __FUNCT__ 103551f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 103651f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 103751f519a2SBarry Smith { 103851f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 103951f519a2SBarry Smith 104051f519a2SBarry Smith PetscFunctionBegin; 104151f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 104251f519a2SBarry Smith if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs); 104351f519a2SBarry Smith jac->bs = bs; 104451f519a2SBarry Smith PetscFunctionReturn(0); 104551f519a2SBarry Smith } 104651f519a2SBarry Smith EXTERN_C_END 104751f519a2SBarry Smith 104879416396SBarry Smith #undef __FUNCT__ 104979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1050bc08b0f1SBarry Smith /*@ 105179416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 105279416396SBarry Smith 105379416396SBarry Smith Collective on PC 105479416396SBarry Smith 105579416396SBarry Smith Input Parameter: 105679416396SBarry Smith . pc - the preconditioner context 1057*a4efd8eaSMatthew Knepley . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 105879416396SBarry Smith 105979416396SBarry Smith Options Database Key: 1060*a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 106179416396SBarry Smith 106279416396SBarry Smith Level: Developer 106379416396SBarry Smith 106479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 106579416396SBarry Smith 106679416396SBarry Smith .seealso: PCCompositeSetType() 106779416396SBarry Smith 106879416396SBarry Smith @*/ 1069dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 107079416396SBarry Smith { 107179416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 107279416396SBarry Smith 107379416396SBarry Smith PetscFunctionBegin; 107479416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 107579416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 107679416396SBarry Smith if (f) { 107779416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 107879416396SBarry Smith } 107979416396SBarry Smith PetscFunctionReturn(0); 108079416396SBarry Smith } 108179416396SBarry Smith 10820971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 10830971522cSBarry Smith /*MC 1084a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 10850971522cSBarry Smith fields or groups of fields 10860971522cSBarry Smith 10870971522cSBarry Smith 1088edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1089edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 10900971522cSBarry Smith 1091a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 109269a612a9SBarry Smith and set the options directly on the resulting KSP object 10930971522cSBarry Smith 10940971522cSBarry Smith Level: intermediate 10950971522cSBarry Smith 109679416396SBarry Smith Options Database Keys: 10972e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 10982e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 10992e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 110051f519a2SBarry Smith . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 1101e69d4d44SBarry Smith . -pc_splitfield_type <additive,multiplicative,schur,symmetric_multiplicative> 1102e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 110379416396SBarry Smith 1104edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 11053b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 11063b224e63SBarry Smith 11073b224e63SBarry Smith 1108d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1109d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1110d32f9abdSBarry Smith 1111d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1112d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1113d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1114d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1115d32f9abdSBarry Smith 1116d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1117d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1118d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1119d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1120d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1121d32f9abdSBarry Smith 1122e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1123e69d4d44SBarry Smith ( C D ) 1124e69d4d44SBarry Smith the preconditioner is 1125e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1126e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1127edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1128e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1129edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1130e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1131e69d4d44SBarry Smith option. 1132e69d4d44SBarry Smith 1133edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1134edf189efSBarry Smith is used automatically for a second block. 1135edf189efSBarry Smith 1136a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 11370971522cSBarry Smith 1138a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1139e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 11400971522cSBarry Smith M*/ 11410971522cSBarry Smith 11420971522cSBarry Smith EXTERN_C_BEGIN 11430971522cSBarry Smith #undef __FUNCT__ 11440971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1145dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 11460971522cSBarry Smith { 11470971522cSBarry Smith PetscErrorCode ierr; 11480971522cSBarry Smith PC_FieldSplit *jac; 11490971522cSBarry Smith 11500971522cSBarry Smith PetscFunctionBegin; 115138f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 11520971522cSBarry Smith jac->bs = -1; 11530971522cSBarry Smith jac->nsplits = 0; 11543e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1155e69d4d44SBarry Smith jac->schurpre = PETSC_TRUE; 115651f519a2SBarry Smith 11570971522cSBarry Smith pc->data = (void*)jac; 11580971522cSBarry Smith 11590971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1160421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 11610971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 11620971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 11630971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 11640971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 11650971522cSBarry Smith pc->ops->applyrichardson = 0; 11660971522cSBarry Smith 116769a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 116869a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 11690971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 11700971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1171704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1172704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 117379416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 117479416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 117551f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 117651f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 11770971522cSBarry Smith PetscFunctionReturn(0); 11780971522cSBarry Smith } 11790971522cSBarry Smith EXTERN_C_END 11800971522cSBarry Smith 11810971522cSBarry Smith 1182a541d17aSBarry Smith /*MC 1183a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1184a541d17aSBarry Smith overview of these methods. 1185a541d17aSBarry Smith 1186a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1187a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1188a541d17aSBarry Smith 1189a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1190a541d17aSBarry Smith B' 0) (x_2) (b_2) 1191a541d17aSBarry Smith 1192a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1193a541d17aSBarry Smith for this block system. 1194a541d17aSBarry Smith 1195a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1196a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1197a541d17aSBarry Smith in the original matrix (for example Ap == A). 1198a541d17aSBarry Smith 1199a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1200a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1201a541d17aSBarry Smith 1202a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1203a541d17aSBarry Smith x_2 = D^ b_2 1204a541d17aSBarry Smith 1205a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1206a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1207a541d17aSBarry Smith 1208a541d17aSBarry Smith Symmetric Gauss-Seidel: x_1 = x_1 + A^(b_1 - A x_1 - B x_2) variant x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2) 1209a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1210a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1211a541d17aSBarry Smith 1212a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1213a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1214a541d17aSBarry Smith M*/ 1215