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; 22a4efd8eaSMatthew Knepley PetscTruth defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */ 2330ad9308SMatthew Knepley PetscInt bs; /* Block size for IS and Mat structures */ 2430ad9308SMatthew Knepley PetscInt nsplits; /* Number of field divisions defined */ 2579416396SBarry Smith Vec *x,*y,w1,w2; 2630ad9308SMatthew Knepley Mat *pmat; /* The diagonal block for each split */ 2730ad9308SMatthew Knepley Mat *Afield; /* The rows of the matrix associated with each split */ 28704ba839SBarry Smith PetscTruth issetup; 2930ad9308SMatthew Knepley /* Only used when Schur complement preconditioning is used */ 3030ad9308SMatthew Knepley Mat B; /* The (0,1) block */ 3130ad9308SMatthew Knepley Mat C; /* The (1,0) block */ 3230ad9308SMatthew Knepley Mat schur; /* The Schur complement S = D - C A^{-1} B */ 3330ad9308SMatthew Knepley KSP kspschur; /* The solver for S */ 34e69d4d44SBarry Smith PetscTruth schurpre; /* preconditioner for the Schur complement is built from pmat[1] == D */ 3597bbdb24SBarry Smith PC_FieldSplitLink head; 360971522cSBarry Smith } PC_FieldSplit; 370971522cSBarry Smith 3816913363SBarry Smith /* 3916913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 4016913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 4116913363SBarry Smith PC you could change this. 4216913363SBarry Smith */ 430971522cSBarry Smith #undef __FUNCT__ 440971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 450971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 460971522cSBarry Smith { 470971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 480971522cSBarry Smith PetscErrorCode ierr; 490971522cSBarry Smith PetscTruth iascii; 500971522cSBarry Smith PetscInt i,j; 515a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 520971522cSBarry Smith 530971522cSBarry Smith PetscFunctionBegin; 540971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 550971522cSBarry Smith if (iascii) { 5651f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 5769a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 580971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 590971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 601ab39975SBarry Smith if (ilink->fields) { 610971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 6279416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 635a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 6479416396SBarry Smith if (j > 0) { 6579416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 6679416396SBarry Smith } 675a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 680971522cSBarry Smith } 690971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 7079416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 711ab39975SBarry Smith } else { 721ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 731ab39975SBarry Smith } 745a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 755a9f2f41SSatish Balay ilink = ilink->next; 760971522cSBarry Smith } 770971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 780971522cSBarry Smith } else { 790971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 800971522cSBarry Smith } 810971522cSBarry Smith PetscFunctionReturn(0); 820971522cSBarry Smith } 830971522cSBarry Smith 840971522cSBarry Smith #undef __FUNCT__ 853b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 863b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 873b224e63SBarry Smith { 883b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 893b224e63SBarry Smith PetscErrorCode ierr; 903b224e63SBarry Smith PetscTruth iascii; 913b224e63SBarry Smith PetscInt i,j; 923b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 933b224e63SBarry Smith KSP ksp; 943b224e63SBarry Smith 953b224e63SBarry Smith PetscFunctionBegin; 963b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 973b224e63SBarry Smith if (iascii) { 983b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 993b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 1003b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1013b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 1023b224e63SBarry Smith if (ilink->fields) { 1033b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1043b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1053b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1063b224e63SBarry Smith if (j > 0) { 1073b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1083b224e63SBarry Smith } 1093b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1103b224e63SBarry Smith } 1113b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1123b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1133b224e63SBarry Smith } else { 1143b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1153b224e63SBarry Smith } 1163b224e63SBarry Smith ilink = ilink->next; 1173b224e63SBarry Smith } 1183b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1193b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1203b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1213b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 1223b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1233b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1243b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1253b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 1263b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1273b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1283b224e63SBarry Smith } else { 1293b224e63SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1303b224e63SBarry Smith } 1313b224e63SBarry Smith PetscFunctionReturn(0); 1323b224e63SBarry Smith } 1333b224e63SBarry Smith 1343b224e63SBarry Smith #undef __FUNCT__ 13569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 13669a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1370971522cSBarry Smith { 1380971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1390971522cSBarry Smith PetscErrorCode ierr; 1405a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 141d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 142d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 143d32f9abdSBarry Smith char optionname[128]; 1440971522cSBarry Smith 1450971522cSBarry Smith PetscFunctionBegin; 146d32f9abdSBarry Smith if (!ilink) { 147d32f9abdSBarry Smith 148521d7252SBarry Smith if (jac->bs <= 0) { 149704ba839SBarry Smith if (pc->pmat) { 150521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 151704ba839SBarry Smith } else { 152704ba839SBarry Smith jac->bs = 1; 153704ba839SBarry Smith } 154521d7252SBarry Smith } 155d32f9abdSBarry Smith 156d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 157d32f9abdSBarry Smith if (!flg) { 158d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 159d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 160d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 161d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 162d32f9abdSBarry Smith while (PETSC_TRUE) { 163d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 164d32f9abdSBarry Smith nfields = jac->bs; 165d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 166d32f9abdSBarry Smith if (!flg2) break; 167d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 168d32f9abdSBarry Smith flg = PETSC_FALSE; 169d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 170d32f9abdSBarry Smith } 171d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 172d32f9abdSBarry Smith } 173d32f9abdSBarry Smith 174d32f9abdSBarry Smith if (flg) { 175d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 17679416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 17779416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1785a9f2f41SSatish Balay while (ilink) { 1795a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1805a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 181521d7252SBarry Smith } 1825a9f2f41SSatish Balay ilink = ilink->next; 18379416396SBarry Smith } 18497bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 18579416396SBarry Smith for (i=0; i<jac->bs; i++) { 18679416396SBarry Smith if (!fields[i]) { 18779416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 18879416396SBarry Smith } else { 18979416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 19079416396SBarry Smith } 19179416396SBarry Smith } 19268dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 193521d7252SBarry Smith } 194edf189efSBarry Smith } else if (jac->nsplits == 1) { 195edf189efSBarry Smith if (ilink->is) { 196edf189efSBarry Smith IS is2; 197edf189efSBarry Smith PetscInt nmin,nmax; 198edf189efSBarry Smith 199edf189efSBarry Smith ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr); 200edf189efSBarry Smith ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr); 201edf189efSBarry Smith ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr); 202edf189efSBarry Smith ierr = ISDestroy(is2);CHKERRQ(ierr); 203edf189efSBarry Smith } else { 204edf189efSBarry Smith SETERRQ(PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()"); 205edf189efSBarry Smith } 206d32f9abdSBarry Smith } 20769a612a9SBarry Smith PetscFunctionReturn(0); 20869a612a9SBarry Smith } 20969a612a9SBarry Smith 21069a612a9SBarry Smith 21169a612a9SBarry Smith #undef __FUNCT__ 21269a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 21369a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 21469a612a9SBarry Smith { 21569a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 21669a612a9SBarry Smith PetscErrorCode ierr; 2175a9f2f41SSatish Balay PC_FieldSplitLink ilink; 218cf502942SBarry Smith PetscInt i,nsplit,ccsize; 21969a612a9SBarry Smith MatStructure flag = pc->flag; 22068dd23aaSBarry Smith PetscTruth sorted,getsub; 22169a612a9SBarry Smith 22269a612a9SBarry Smith PetscFunctionBegin; 22369a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 22497bbdb24SBarry Smith nsplit = jac->nsplits; 2255a9f2f41SSatish Balay ilink = jac->head; 22697bbdb24SBarry Smith 22797bbdb24SBarry Smith /* get the matrices for each split */ 228704ba839SBarry Smith if (!jac->issetup) { 2291b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 23097bbdb24SBarry Smith 231704ba839SBarry Smith jac->issetup = PETSC_TRUE; 232704ba839SBarry Smith 233704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 23451f519a2SBarry Smith bs = jac->bs; 23597bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 236cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2371b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2381b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2391b9fc7fcSBarry Smith if (jac->defaultsplit) { 240704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 241704ba839SBarry Smith } else if (!ilink->is) { 242ccb205f8SBarry Smith if (ilink->nfields > 1) { 2435a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2445a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2451b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2461b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2471b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 24897bbdb24SBarry Smith } 24997bbdb24SBarry Smith } 250704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 251ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 252ccb205f8SBarry Smith } else { 253704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 254ccb205f8SBarry Smith } 2553e197d65SBarry Smith } 2563e197d65SBarry Smith ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr); 257704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 2581b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 259704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 260704ba839SBarry Smith ilink = ilink->next; 2611b9fc7fcSBarry Smith } 2621b9fc7fcSBarry Smith } 2631b9fc7fcSBarry Smith 264704ba839SBarry Smith ilink = jac->head; 26597bbdb24SBarry Smith if (!jac->pmat) { 266cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 267cf502942SBarry Smith for (i=0; i<nsplit; i++) { 268704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 269704ba839SBarry Smith ilink = ilink->next; 270cf502942SBarry Smith } 27197bbdb24SBarry Smith } else { 272cf502942SBarry Smith for (i=0; i<nsplit; i++) { 273704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 274704ba839SBarry Smith ilink = ilink->next; 275cf502942SBarry Smith } 27697bbdb24SBarry Smith } 27797bbdb24SBarry Smith 27868dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 27968dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 2803b224e63SBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 28168dd23aaSBarry Smith ilink = jac->head; 28268dd23aaSBarry Smith if (!jac->Afield) { 28368dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 28468dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 28511755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 28668dd23aaSBarry Smith ilink = ilink->next; 28768dd23aaSBarry Smith } 28868dd23aaSBarry Smith } else { 28968dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 29011755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 29168dd23aaSBarry Smith ilink = ilink->next; 29268dd23aaSBarry Smith } 29368dd23aaSBarry Smith } 29468dd23aaSBarry Smith } 29568dd23aaSBarry Smith 2963b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 2973b224e63SBarry Smith IS ccis; 298e69d4d44SBarry Smith PetscInt N,nlocal,nis; 2993b224e63SBarry Smith if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 30068dd23aaSBarry Smith 3013b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 3023b224e63SBarry Smith if (jac->schur) { 3033b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3043b224e63SBarry Smith ilink = jac->head; 305edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 306e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 307e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 308e69d4d44SBarry Smith nlocal = nlocal - nis; 309e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 3103b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3113b224e63SBarry Smith ilink = ilink->next; 312edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 313e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 314e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 315e69d4d44SBarry Smith nlocal = nlocal - nis; 316e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 3173b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3183b224e63SBarry Smith ierr = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 319e69d4d44SBarry Smith if (jac->schurpre) { 320e69d4d44SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],pc->flag);CHKERRQ(ierr); 321e69d4d44SBarry Smith } else { 3223b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,pc->flag);CHKERRQ(ierr); 323e69d4d44SBarry Smith } 3243b224e63SBarry Smith 3253b224e63SBarry Smith } else { 3261cee3971SBarry Smith KSP ksp; 3273b224e63SBarry Smith 3283b224e63SBarry Smith /* extract the B and C matrices */ 3293b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3303b224e63SBarry Smith ilink = jac->head; 331edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 332e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 333e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 334e69d4d44SBarry Smith nlocal = nlocal - nis; 335e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3363b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3373b224e63SBarry Smith ilink = ilink->next; 338edf189efSBarry Smith ierr = ISComplement(ilink->cis,0,N,&ccis);CHKERRQ(ierr); 339e69d4d44SBarry Smith ierr = ISGetLocalSize(ilink->is,&nis);CHKERRQ(ierr); 340e69d4d44SBarry Smith ierr = MatGetLocalSize(pc->mat,PETSC_NULL,&nlocal);CHKERRQ(ierr); 341e69d4d44SBarry Smith nlocal = nlocal - nis; 342e69d4d44SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,nlocal,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3433b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3443b224e63SBarry Smith ierr = MatCreateSchurComplement(jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr); 3451cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 346e69d4d44SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr); 3473b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3483b224e63SBarry Smith 3493b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 3501cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 351e69d4d44SBarry Smith if (jac->schurpre) { 352e69d4d44SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->pmat[1],DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 353e69d4d44SBarry Smith } else { 3543b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 355e69d4d44SBarry Smith } 3563b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 357edf189efSBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr); 3583b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3593b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3603b224e63SBarry Smith 3613b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3623b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3633b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3643b224e63SBarry Smith ilink = jac->head; 3653b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3663b224e63SBarry Smith ilink = ilink->next; 3673b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3683b224e63SBarry Smith } 3693b224e63SBarry Smith } else { 37097bbdb24SBarry Smith /* set up the individual PCs */ 37197bbdb24SBarry Smith i = 0; 3725a9f2f41SSatish Balay ilink = jac->head; 3735a9f2f41SSatish Balay while (ilink) { 3745a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 3753b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3765a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 3775a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 37897bbdb24SBarry Smith i++; 3795a9f2f41SSatish Balay ilink = ilink->next; 3800971522cSBarry Smith } 38197bbdb24SBarry Smith 38297bbdb24SBarry Smith /* create work vectors for each split */ 3831b9fc7fcSBarry Smith if (!jac->x) { 38497bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 3855a9f2f41SSatish Balay ilink = jac->head; 38697bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 387906ed7ccSBarry Smith Vec *vl,*vr; 388906ed7ccSBarry Smith 389906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 390906ed7ccSBarry Smith ilink->x = *vr; 391906ed7ccSBarry Smith ilink->y = *vl; 392906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 393906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 3945a9f2f41SSatish Balay jac->x[i] = ilink->x; 3955a9f2f41SSatish Balay jac->y[i] = ilink->y; 3965a9f2f41SSatish Balay ilink = ilink->next; 39797bbdb24SBarry Smith } 3983b224e63SBarry Smith } 3993b224e63SBarry Smith } 4003b224e63SBarry Smith 4013b224e63SBarry Smith 402d0f46423SBarry Smith if (!jac->head->sctx) { 4033b224e63SBarry Smith Vec xtmp; 4043b224e63SBarry Smith 40579416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 4061b9fc7fcSBarry Smith 4075a9f2f41SSatish Balay ilink = jac->head; 4081b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 4091b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 410704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 4115a9f2f41SSatish Balay ilink = ilink->next; 4121b9fc7fcSBarry Smith } 4131b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 4141b9fc7fcSBarry Smith } 4150971522cSBarry Smith PetscFunctionReturn(0); 4160971522cSBarry Smith } 4170971522cSBarry Smith 4185a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 419ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 420ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 4215a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 422ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 423ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 42479416396SBarry Smith 4250971522cSBarry Smith #undef __FUNCT__ 4263b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 4273b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 4283b224e63SBarry Smith { 4293b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4303b224e63SBarry Smith PetscErrorCode ierr; 4313b224e63SBarry Smith KSP ksp; 4323b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 4333b224e63SBarry Smith 4343b224e63SBarry Smith PetscFunctionBegin; 4353b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4363b224e63SBarry Smith 4373b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4383b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4393b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4403b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4413b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4423b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4433b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4443b224e63SBarry Smith 4453b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4463b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4473b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4483b224e63SBarry Smith 4493b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4503b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4513b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4523b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4533b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4543b224e63SBarry Smith 4553b224e63SBarry Smith PetscFunctionReturn(0); 4563b224e63SBarry Smith } 4573b224e63SBarry Smith 4583b224e63SBarry Smith #undef __FUNCT__ 4590971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4600971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4610971522cSBarry Smith { 4620971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4630971522cSBarry Smith PetscErrorCode ierr; 4645a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4653e197d65SBarry Smith PetscInt bs,cnt; 4660971522cSBarry Smith 4670971522cSBarry Smith PetscFunctionBegin; 46851f519a2SBarry Smith CHKMEMQ; 469e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 47051f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 47151f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 47251f519a2SBarry Smith 47379416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4741b9fc7fcSBarry Smith if (jac->defaultsplit) { 4750971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 4765a9f2f41SSatish Balay while (ilink) { 4775a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4785a9f2f41SSatish Balay ilink = ilink->next; 4790971522cSBarry Smith } 4800971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 4811b9fc7fcSBarry Smith } else { 482efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4835a9f2f41SSatish Balay while (ilink) { 4845a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4855a9f2f41SSatish Balay ilink = ilink->next; 4861b9fc7fcSBarry Smith } 4871b9fc7fcSBarry Smith } 48816913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 48979416396SBarry Smith if (!jac->w1) { 49079416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 49179416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 49279416396SBarry Smith } 493efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4945a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4953e197d65SBarry Smith cnt = 1; 4965a9f2f41SSatish Balay while (ilink->next) { 4975a9f2f41SSatish Balay ilink = ilink->next; 4983e197d65SBarry Smith if (jac->Afield) { 4993e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 5003e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 5013e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 5023e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5033e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5043e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 5053e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5063e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5073e197d65SBarry Smith } else { 5083e197d65SBarry Smith /* compute the residual over the entire vector */ 5091ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 510efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 5115a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 51279416396SBarry Smith } 5133e197d65SBarry Smith } 51451f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 51511755939SBarry Smith cnt -= 2; 51651f519a2SBarry Smith while (ilink->previous) { 51751f519a2SBarry Smith ilink = ilink->previous; 51811755939SBarry Smith if (jac->Afield) { 51911755939SBarry Smith /* compute the residual only over the part of the vector needed */ 52011755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 52111755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 52211755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52311755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52411755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 52511755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 52611755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 52711755939SBarry Smith } else { 5281ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 52951f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 53051f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 53179416396SBarry Smith } 53251f519a2SBarry Smith } 53311755939SBarry Smith } 53416913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 53551f519a2SBarry Smith CHKMEMQ; 5360971522cSBarry Smith PetscFunctionReturn(0); 5370971522cSBarry Smith } 5380971522cSBarry Smith 539421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 540ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 541ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 542421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 543ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 544ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 545421e10b8SBarry Smith 546421e10b8SBarry Smith #undef __FUNCT__ 547421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 548421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 549421e10b8SBarry Smith { 550421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 551421e10b8SBarry Smith PetscErrorCode ierr; 552421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 553421e10b8SBarry Smith PetscInt bs; 554421e10b8SBarry Smith 555421e10b8SBarry Smith PetscFunctionBegin; 556421e10b8SBarry Smith CHKMEMQ; 557421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 558421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 559421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 560421e10b8SBarry Smith 561421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 562421e10b8SBarry Smith if (jac->defaultsplit) { 563421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 564421e10b8SBarry Smith while (ilink) { 565421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 566421e10b8SBarry Smith ilink = ilink->next; 567421e10b8SBarry Smith } 568421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 569421e10b8SBarry Smith } else { 570421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 571421e10b8SBarry Smith while (ilink) { 572421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 573421e10b8SBarry Smith ilink = ilink->next; 574421e10b8SBarry Smith } 575421e10b8SBarry Smith } 576421e10b8SBarry Smith } else { 577421e10b8SBarry Smith if (!jac->w1) { 578421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 579421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 580421e10b8SBarry Smith } 581421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 582421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 583421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 584421e10b8SBarry Smith while (ilink->next) { 585421e10b8SBarry Smith ilink = ilink->next; 5869989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 587421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 588421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 589421e10b8SBarry Smith } 590421e10b8SBarry Smith while (ilink->previous) { 591421e10b8SBarry Smith ilink = ilink->previous; 5929989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 593421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 594421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 595421e10b8SBarry Smith } 596421e10b8SBarry Smith } else { 597421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 598421e10b8SBarry Smith ilink = ilink->next; 599421e10b8SBarry Smith } 600421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 601421e10b8SBarry Smith while (ilink->previous) { 602421e10b8SBarry Smith ilink = ilink->previous; 6039989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 604421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 605421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 606421e10b8SBarry Smith } 607421e10b8SBarry Smith } 608421e10b8SBarry Smith } 609421e10b8SBarry Smith CHKMEMQ; 610421e10b8SBarry Smith PetscFunctionReturn(0); 611421e10b8SBarry Smith } 612421e10b8SBarry Smith 6130971522cSBarry Smith #undef __FUNCT__ 6140971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 6150971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 6160971522cSBarry Smith { 6170971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6180971522cSBarry Smith PetscErrorCode ierr; 6195a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 6200971522cSBarry Smith 6210971522cSBarry Smith PetscFunctionBegin; 6225a9f2f41SSatish Balay while (ilink) { 6235a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 6245a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 6255a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 6265a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 627704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 628704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 6295a9f2f41SSatish Balay next = ilink->next; 630704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 631704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 6325a9f2f41SSatish Balay ilink = next; 6330971522cSBarry Smith } 63405b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 63597bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 63668dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 63779416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 63879416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6393b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 640d0f46423SBarry Smith if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);} 6413b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6423b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6430971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6440971522cSBarry Smith PetscFunctionReturn(0); 6450971522cSBarry Smith } 6460971522cSBarry Smith 6470971522cSBarry Smith #undef __FUNCT__ 6480971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6490971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6500971522cSBarry Smith { 6511b9fc7fcSBarry Smith PetscErrorCode ierr; 65251f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 653e69d4d44SBarry Smith PetscTruth flg,set; 6541b9fc7fcSBarry Smith char optionname[128]; 6559dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6563b224e63SBarry Smith PCCompositeType ctype; 6571b9fc7fcSBarry Smith 6580971522cSBarry Smith PetscFunctionBegin; 6591b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 66051f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 66151f519a2SBarry Smith if (flg) { 66251f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 66351f519a2SBarry Smith } 664704ba839SBarry Smith 6653b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6663b224e63SBarry Smith if (flg) { 6673b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6683b224e63SBarry Smith } 669d32f9abdSBarry Smith 670c30613efSMatthew Knepley /* Only setup fields once */ 671c30613efSMatthew Knepley if ((jac->bs > 0) && (jac->nsplits == 0)) { 672d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 673d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 67451f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6751b9fc7fcSBarry Smith while (PETSC_TRUE) { 67613f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 67751f519a2SBarry Smith nfields = jac->bs; 6781b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6791b9fc7fcSBarry Smith if (!flg) break; 6801b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 6811b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 6821b9fc7fcSBarry Smith } 68351f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 684d32f9abdSBarry Smith } 685e69d4d44SBarry Smith ierr = PetscOptionsTruth("-pc_fieldsplit_schur_precondition","Build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",jac->schurpre,&flg,&set);CHKERRQ(ierr); 686e69d4d44SBarry Smith if (set) { 687e69d4d44SBarry Smith ierr = PCFieldSplitSchurPrecondition(pc,flg);CHKERRQ(ierr); 688e69d4d44SBarry Smith } 6891b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6900971522cSBarry Smith PetscFunctionReturn(0); 6910971522cSBarry Smith } 6920971522cSBarry Smith 6930971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 6940971522cSBarry Smith 6950971522cSBarry Smith EXTERN_C_BEGIN 6960971522cSBarry Smith #undef __FUNCT__ 6970971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 698dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 6990971522cSBarry Smith { 70097bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7010971522cSBarry Smith PetscErrorCode ierr; 7025a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 70369a612a9SBarry Smith char prefix[128]; 70451f519a2SBarry Smith PetscInt i; 7050971522cSBarry Smith 7060971522cSBarry Smith PetscFunctionBegin; 7070971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 70851f519a2SBarry Smith for (i=0; i<n; i++) { 70951f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 71051f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 71151f519a2SBarry Smith } 712704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 713704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 7145a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 7155a9f2f41SSatish Balay ilink->nfields = n; 7165a9f2f41SSatish Balay ilink->next = PETSC_NULL; 7177adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7181cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7195a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 72069a612a9SBarry Smith 7217adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 7227adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 72369a612a9SBarry Smith } else { 72413f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 72569a612a9SBarry Smith } 7265a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 7270971522cSBarry Smith 7280971522cSBarry Smith if (!next) { 7295a9f2f41SSatish Balay jac->head = ilink; 73051f519a2SBarry Smith ilink->previous = PETSC_NULL; 7310971522cSBarry Smith } else { 7320971522cSBarry Smith while (next->next) { 7330971522cSBarry Smith next = next->next; 7340971522cSBarry Smith } 7355a9f2f41SSatish Balay next->next = ilink; 73651f519a2SBarry Smith ilink->previous = next; 7370971522cSBarry Smith } 7380971522cSBarry Smith jac->nsplits++; 7390971522cSBarry Smith PetscFunctionReturn(0); 7400971522cSBarry Smith } 7410971522cSBarry Smith EXTERN_C_END 7420971522cSBarry Smith 743e69d4d44SBarry Smith EXTERN_C_BEGIN 744e69d4d44SBarry Smith #undef __FUNCT__ 745e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur" 746e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp) 747e69d4d44SBarry Smith { 748e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 749e69d4d44SBarry Smith PetscErrorCode ierr; 750e69d4d44SBarry Smith 751e69d4d44SBarry Smith PetscFunctionBegin; 752e69d4d44SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 753e69d4d44SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr); 754e69d4d44SBarry Smith (*subksp)[1] = jac->kspschur; 755e69d4d44SBarry Smith PetscFunctionReturn(0); 756e69d4d44SBarry Smith } 757e69d4d44SBarry Smith EXTERN_C_END 7580971522cSBarry Smith 7590971522cSBarry Smith EXTERN_C_BEGIN 7600971522cSBarry Smith #undef __FUNCT__ 76169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 762dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7630971522cSBarry Smith { 7640971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7650971522cSBarry Smith PetscErrorCode ierr; 7660971522cSBarry Smith PetscInt cnt = 0; 7675a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7680971522cSBarry Smith 7690971522cSBarry Smith PetscFunctionBegin; 77069a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7715a9f2f41SSatish Balay while (ilink) { 7725a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7735a9f2f41SSatish Balay ilink = ilink->next; 7740971522cSBarry Smith } 7750971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7760971522cSBarry Smith *n = jac->nsplits; 7770971522cSBarry Smith PetscFunctionReturn(0); 7780971522cSBarry Smith } 7790971522cSBarry Smith EXTERN_C_END 7800971522cSBarry Smith 781704ba839SBarry Smith EXTERN_C_BEGIN 782704ba839SBarry Smith #undef __FUNCT__ 783704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 784704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 785704ba839SBarry Smith { 786704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 787704ba839SBarry Smith PetscErrorCode ierr; 788704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 789704ba839SBarry Smith char prefix[128]; 790704ba839SBarry Smith 791704ba839SBarry Smith PetscFunctionBegin; 79216913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 7931ab39975SBarry Smith ilink->is = is; 7941ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 795704ba839SBarry Smith ilink->next = PETSC_NULL; 796704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 7971cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 798704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 799704ba839SBarry Smith 800704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 801704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 802704ba839SBarry Smith } else { 803704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 804704ba839SBarry Smith } 805704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 806704ba839SBarry Smith 807704ba839SBarry Smith if (!next) { 808704ba839SBarry Smith jac->head = ilink; 809704ba839SBarry Smith ilink->previous = PETSC_NULL; 810704ba839SBarry Smith } else { 811704ba839SBarry Smith while (next->next) { 812704ba839SBarry Smith next = next->next; 813704ba839SBarry Smith } 814704ba839SBarry Smith next->next = ilink; 815704ba839SBarry Smith ilink->previous = next; 816704ba839SBarry Smith } 817704ba839SBarry Smith jac->nsplits++; 818704ba839SBarry Smith 819704ba839SBarry Smith PetscFunctionReturn(0); 820704ba839SBarry Smith } 821704ba839SBarry Smith EXTERN_C_END 822704ba839SBarry Smith 8230971522cSBarry Smith #undef __FUNCT__ 8240971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 8250971522cSBarry Smith /*@ 8260971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 8270971522cSBarry Smith 8280971522cSBarry Smith Collective on PC 8290971522cSBarry Smith 8300971522cSBarry Smith Input Parameters: 8310971522cSBarry Smith + pc - the preconditioner context 8320971522cSBarry Smith . n - the number of fields in this split 8330971522cSBarry Smith . fields - the fields in this split 8340971522cSBarry Smith 8350971522cSBarry Smith Level: intermediate 8360971522cSBarry Smith 837d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 838d32f9abdSBarry Smith 839d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 840d32f9abdSBarry 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 841d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 842d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 843d32f9abdSBarry Smith 844d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 8450971522cSBarry Smith 8460971522cSBarry Smith @*/ 847dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 8480971522cSBarry Smith { 8490971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 8500971522cSBarry Smith 8510971522cSBarry Smith PetscFunctionBegin; 8520971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8530971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 8540971522cSBarry Smith if (f) { 8550971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8560971522cSBarry Smith } 8570971522cSBarry Smith PetscFunctionReturn(0); 8580971522cSBarry Smith } 8590971522cSBarry Smith 8600971522cSBarry Smith #undef __FUNCT__ 861704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 862704ba839SBarry Smith /*@ 863704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 864704ba839SBarry Smith 865704ba839SBarry Smith Collective on PC 866704ba839SBarry Smith 867704ba839SBarry Smith Input Parameters: 868704ba839SBarry Smith + pc - the preconditioner context 869704ba839SBarry Smith . is - the index set that defines the vector elements in this field 870704ba839SBarry Smith 871d32f9abdSBarry Smith 872d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 873d32f9abdSBarry Smith 874704ba839SBarry Smith Level: intermediate 875704ba839SBarry Smith 876704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 877704ba839SBarry Smith 878704ba839SBarry Smith @*/ 879704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 880704ba839SBarry Smith { 881704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 882704ba839SBarry Smith 883704ba839SBarry Smith PetscFunctionBegin; 884704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 885704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 886704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 887704ba839SBarry Smith if (f) { 888704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 889704ba839SBarry Smith } 890704ba839SBarry Smith PetscFunctionReturn(0); 891704ba839SBarry Smith } 892704ba839SBarry Smith 893704ba839SBarry Smith #undef __FUNCT__ 89451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 89551f519a2SBarry Smith /*@ 89651f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 89751f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 89851f519a2SBarry Smith 89951f519a2SBarry Smith Collective on PC 90051f519a2SBarry Smith 90151f519a2SBarry Smith Input Parameters: 90251f519a2SBarry Smith + pc - the preconditioner context 90351f519a2SBarry Smith - bs - the block size 90451f519a2SBarry Smith 90551f519a2SBarry Smith Level: intermediate 90651f519a2SBarry Smith 90751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 90851f519a2SBarry Smith 90951f519a2SBarry Smith @*/ 91051f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 91151f519a2SBarry Smith { 91251f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 91351f519a2SBarry Smith 91451f519a2SBarry Smith PetscFunctionBegin; 91551f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 91651f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 91751f519a2SBarry Smith if (f) { 91851f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 91951f519a2SBarry Smith } 92051f519a2SBarry Smith PetscFunctionReturn(0); 92151f519a2SBarry Smith } 92251f519a2SBarry Smith 92351f519a2SBarry Smith #undef __FUNCT__ 92469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 9250971522cSBarry Smith /*@C 92669a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 9270971522cSBarry Smith 92869a612a9SBarry Smith Collective on KSP 9290971522cSBarry Smith 9300971522cSBarry Smith Input Parameter: 9310971522cSBarry Smith . pc - the preconditioner context 9320971522cSBarry Smith 9330971522cSBarry Smith Output Parameters: 9340971522cSBarry Smith + n - the number of split 93569a612a9SBarry Smith - pc - the array of KSP contexts 9360971522cSBarry Smith 9370971522cSBarry Smith Note: 938d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 939d32f9abdSBarry Smith (not the KSP just the array that contains them). 9400971522cSBarry Smith 94169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 9420971522cSBarry Smith 9430971522cSBarry Smith Level: advanced 9440971522cSBarry Smith 9450971522cSBarry Smith .seealso: PCFIELDSPLIT 9460971522cSBarry Smith @*/ 947dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 9480971522cSBarry Smith { 94969a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 9500971522cSBarry Smith 9510971522cSBarry Smith PetscFunctionBegin; 9520971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9530971522cSBarry Smith PetscValidIntPointer(n,2); 95469a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 9550971522cSBarry Smith if (f) { 95669a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 9570971522cSBarry Smith } else { 95869a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9590971522cSBarry Smith } 9600971522cSBarry Smith PetscFunctionReturn(0); 9610971522cSBarry Smith } 9620971522cSBarry Smith 963e69d4d44SBarry Smith #undef __FUNCT__ 964e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition" 965e69d4d44SBarry Smith /*@ 966e69d4d44SBarry Smith PCFieldSplitSchurPrecondition - Indicates if the Schur complement is preconditioned by a preconditioner constructed by the 967e69d4d44SBarry Smith D matrix. Otherwise no preconditioner is used. 968e69d4d44SBarry Smith 969e69d4d44SBarry Smith Collective on PC 970e69d4d44SBarry Smith 971e69d4d44SBarry Smith Input Parameters: 972e69d4d44SBarry Smith + pc - the preconditioner context 973e69d4d44SBarry Smith - flg - build the preconditioner 974e69d4d44SBarry Smith 975e69d4d44SBarry Smith Options Database: 976e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 977e69d4d44SBarry Smith 978e69d4d44SBarry Smith Level: intermediate 979e69d4d44SBarry Smith 980e69d4d44SBarry 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 981e69d4d44SBarry Smith 982e69d4d44SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFIELDSPLIT 983e69d4d44SBarry Smith 984e69d4d44SBarry Smith @*/ 985e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PetscTruth flg) 986e69d4d44SBarry Smith { 987e69d4d44SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscTruth); 988e69d4d44SBarry Smith 989e69d4d44SBarry Smith PetscFunctionBegin; 990e69d4d44SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 991e69d4d44SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr); 992e69d4d44SBarry Smith if (f) { 993e69d4d44SBarry Smith ierr = (*f)(pc,flg);CHKERRQ(ierr); 994e69d4d44SBarry Smith } 995e69d4d44SBarry Smith PetscFunctionReturn(0); 996e69d4d44SBarry Smith } 997e69d4d44SBarry Smith 998e69d4d44SBarry Smith EXTERN_C_BEGIN 999e69d4d44SBarry Smith #undef __FUNCT__ 1000e69d4d44SBarry Smith #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit" 1001e69d4d44SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PetscTruth flg) 1002e69d4d44SBarry Smith { 1003e69d4d44SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1004e69d4d44SBarry Smith 1005e69d4d44SBarry Smith PetscFunctionBegin; 1006e69d4d44SBarry Smith jac->schurpre = flg; 1007e69d4d44SBarry Smith PetscFunctionReturn(0); 1008e69d4d44SBarry Smith } 1009e69d4d44SBarry Smith EXTERN_C_END 1010e69d4d44SBarry Smith 101130ad9308SMatthew Knepley #undef __FUNCT__ 101230ad9308SMatthew Knepley #define __FUNCT__ "PCFieldSplitGetSchurBlocks" 101330ad9308SMatthew Knepley /*@C 101430ad9308SMatthew Knepley PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement 101530ad9308SMatthew Knepley 101630ad9308SMatthew Knepley Collective on KSP 101730ad9308SMatthew Knepley 101830ad9308SMatthew Knepley Input Parameter: 101930ad9308SMatthew Knepley . pc - the preconditioner context 102030ad9308SMatthew Knepley 102130ad9308SMatthew Knepley Output Parameters: 102230ad9308SMatthew Knepley + A - the (0,0) block 102330ad9308SMatthew Knepley . B - the (0,1) block 102430ad9308SMatthew Knepley . C - the (1,0) block 102530ad9308SMatthew Knepley - D - the (1,1) block 102630ad9308SMatthew Knepley 102730ad9308SMatthew Knepley Level: advanced 102830ad9308SMatthew Knepley 102930ad9308SMatthew Knepley .seealso: PCFIELDSPLIT 103030ad9308SMatthew Knepley @*/ 103130ad9308SMatthew Knepley PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D) 103230ad9308SMatthew Knepley { 103330ad9308SMatthew Knepley PC_FieldSplit *jac = (PC_FieldSplit *) pc->data; 103430ad9308SMatthew Knepley 103530ad9308SMatthew Knepley PetscFunctionBegin; 103630ad9308SMatthew Knepley PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1037cf8ad460SMatthew Knepley if (jac->type != PC_COMPOSITE_SCHUR) {SETERRQ(PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");} 103830ad9308SMatthew Knepley if (A) *A = jac->pmat[0]; 103930ad9308SMatthew Knepley if (B) *B = jac->B; 104030ad9308SMatthew Knepley if (C) *C = jac->C; 104130ad9308SMatthew Knepley if (D) *D = jac->pmat[1]; 104230ad9308SMatthew Knepley PetscFunctionReturn(0); 104330ad9308SMatthew Knepley } 104430ad9308SMatthew Knepley 104579416396SBarry Smith EXTERN_C_BEGIN 104679416396SBarry Smith #undef __FUNCT__ 104779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 1048dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 104979416396SBarry Smith { 105079416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1051e69d4d44SBarry Smith PetscErrorCode ierr; 105279416396SBarry Smith 105379416396SBarry Smith PetscFunctionBegin; 105479416396SBarry Smith jac->type = type; 10553b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 10563b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 10573b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 1058e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr); 1059e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr); 1060e69d4d44SBarry Smith 10613b224e63SBarry Smith } else { 10623b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 10633b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 1064e69d4d44SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 10659e7d6b0aSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr); 10663b224e63SBarry Smith } 106779416396SBarry Smith PetscFunctionReturn(0); 106879416396SBarry Smith } 106979416396SBarry Smith EXTERN_C_END 107079416396SBarry Smith 107151f519a2SBarry Smith EXTERN_C_BEGIN 107251f519a2SBarry Smith #undef __FUNCT__ 107351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 107451f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 107551f519a2SBarry Smith { 107651f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 107751f519a2SBarry Smith 107851f519a2SBarry Smith PetscFunctionBegin; 107951f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 108051f519a2SBarry 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); 108151f519a2SBarry Smith jac->bs = bs; 108251f519a2SBarry Smith PetscFunctionReturn(0); 108351f519a2SBarry Smith } 108451f519a2SBarry Smith EXTERN_C_END 108551f519a2SBarry Smith 108679416396SBarry Smith #undef __FUNCT__ 108779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 1088bc08b0f1SBarry Smith /*@ 108979416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 109079416396SBarry Smith 109179416396SBarry Smith Collective on PC 109279416396SBarry Smith 109379416396SBarry Smith Input Parameter: 109479416396SBarry Smith . pc - the preconditioner context 109581540f2fSBarry Smith . type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR 109679416396SBarry Smith 109779416396SBarry Smith Options Database Key: 1098a4efd8eaSMatthew Knepley . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type 109979416396SBarry Smith 110079416396SBarry Smith Level: Developer 110179416396SBarry Smith 110279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 110379416396SBarry Smith 110479416396SBarry Smith .seealso: PCCompositeSetType() 110579416396SBarry Smith 110679416396SBarry Smith @*/ 1107dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 110879416396SBarry Smith { 110979416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 111079416396SBarry Smith 111179416396SBarry Smith PetscFunctionBegin; 111279416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 111379416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 111479416396SBarry Smith if (f) { 111579416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 111679416396SBarry Smith } 111779416396SBarry Smith PetscFunctionReturn(0); 111879416396SBarry Smith } 111979416396SBarry Smith 11200971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 11210971522cSBarry Smith /*MC 1122a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 11230971522cSBarry Smith fields or groups of fields 11240971522cSBarry Smith 11250971522cSBarry Smith 1126edf189efSBarry Smith To set options on the solvers for each block append -fieldsplit_ to all the PC 1127edf189efSBarry Smith options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1 11280971522cSBarry Smith 1129a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 113069a612a9SBarry Smith and set the options directly on the resulting KSP object 11310971522cSBarry Smith 11320971522cSBarry Smith Level: intermediate 11330971522cSBarry Smith 113479416396SBarry Smith Options Database Keys: 113581540f2fSBarry Smith + -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 113681540f2fSBarry Smith . -pc_fieldsplit_default - automatically add any fields to additional splits that have not 113781540f2fSBarry Smith been supplied explicitly by -pc_fieldsplit_%d_fields 113881540f2fSBarry Smith . -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 113981540f2fSBarry Smith . -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative> 1140e69d4d44SBarry Smith . -pc_fieldsplit_schur_precondition <true,false> default is true 114179416396SBarry Smith 1142edf189efSBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_ 11433b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 11443b224e63SBarry Smith 11453b224e63SBarry Smith 1146d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1147d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1148d32f9abdSBarry Smith 1149d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1150d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1151d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1152d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1153d32f9abdSBarry Smith 1154d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1155d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1156d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1157d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1158d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1159d32f9abdSBarry Smith 1160e69d4d44SBarry Smith For the Schur complement preconditioner if J = ( A B ) 1161e69d4d44SBarry Smith ( C D ) 1162e69d4d44SBarry Smith the preconditioner is 1163e69d4d44SBarry Smith (I -B inv(A)) ( inv(A) 0 ) (I 0 ) 1164e69d4d44SBarry Smith (0 I ) ( 0 inv(S) ) (-C inv(A) I ) 1165edf189efSBarry Smith where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of 1166e69d4d44SBarry Smith inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is 1167edf189efSBarry Smith 0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default 1168e69d4d44SBarry Smith D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this 1169e69d4d44SBarry Smith option. 1170e69d4d44SBarry Smith 1171edf189efSBarry Smith If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS 1172edf189efSBarry Smith is used automatically for a second block. 1173edf189efSBarry Smith 1174a541d17aSBarry Smith Concepts: physics based preconditioners, block preconditioners 11750971522cSBarry Smith 1176a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners 1177e69d4d44SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 11780971522cSBarry Smith M*/ 11790971522cSBarry Smith 11800971522cSBarry Smith EXTERN_C_BEGIN 11810971522cSBarry Smith #undef __FUNCT__ 11820971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1183dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 11840971522cSBarry Smith { 11850971522cSBarry Smith PetscErrorCode ierr; 11860971522cSBarry Smith PC_FieldSplit *jac; 11870971522cSBarry Smith 11880971522cSBarry Smith PetscFunctionBegin; 118938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 11900971522cSBarry Smith jac->bs = -1; 11910971522cSBarry Smith jac->nsplits = 0; 11923e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 1193e69d4d44SBarry Smith jac->schurpre = PETSC_TRUE; 119451f519a2SBarry Smith 11950971522cSBarry Smith pc->data = (void*)jac; 11960971522cSBarry Smith 11970971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1198421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 11990971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 12000971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 12010971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 12020971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 12030971522cSBarry Smith pc->ops->applyrichardson = 0; 12040971522cSBarry Smith 120569a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 120669a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 12070971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 12080971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1209704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1210704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 121179416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 121279416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 121351f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 121451f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 12150971522cSBarry Smith PetscFunctionReturn(0); 12160971522cSBarry Smith } 12170971522cSBarry Smith EXTERN_C_END 12180971522cSBarry Smith 12190971522cSBarry Smith 1220a541d17aSBarry Smith /*MC 1221a541d17aSBarry Smith Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an 1222a541d17aSBarry Smith overview of these methods. 1223a541d17aSBarry Smith 1224a541d17aSBarry Smith Consider the solution to ( A B ) (x_1) = (b_1) 1225a541d17aSBarry Smith ( C D ) (x_2) (b_2) 1226a541d17aSBarry Smith 1227a541d17aSBarry Smith Important special cases, the Stokes equation: C = B' and D = 0 (A B) (x_1) = (b_1) 1228a541d17aSBarry Smith B' 0) (x_2) (b_2) 1229a541d17aSBarry Smith 1230a541d17aSBarry Smith One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners 1231a541d17aSBarry Smith for this block system. 1232a541d17aSBarry Smith 1233a541d17aSBarry Smith Consider an additional matrix (Ap Bp) 1234a541d17aSBarry Smith (Cp Dp) where some or all of the entries may be the same as 1235a541d17aSBarry Smith in the original matrix (for example Ap == A). 1236a541d17aSBarry Smith 1237a541d17aSBarry Smith In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the 1238a541d17aSBarry Smith approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...) 1239a541d17aSBarry Smith 1240a541d17aSBarry Smith Block Jacobi: x_1 = A^ b_1 1241a541d17aSBarry Smith x_2 = D^ b_2 1242a541d17aSBarry Smith 1243a541d17aSBarry Smith Lower block Gauss-Seidel: x_1 = A^ b_1 1244a541d17aSBarry Smith x_2 = D^ (b_2 - C x_1) variant x_2 = D^ (b_2 - Cp x_1) 1245a541d17aSBarry Smith 1246a541d17aSBarry 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) 1247a541d17aSBarry Smith Interestingly this form is not actually a symmetric matrix, the symmetric version is 1248a541d17aSBarry Smith x_1 = A^(b_1 - B x_2) variant x_1 = A^(b_1 - Bp x_2) 1249a541d17aSBarry Smith 1250*0bc0a719SBarry Smith Level: intermediate 1251*0bc0a719SBarry Smith 1252a541d17aSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT 1253a541d17aSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition() 1254a541d17aSBarry Smith M*/ 1255