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; 2297bbdb24SBarry Smith PetscTruth defaultsplit; 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; 3197bbdb24SBarry Smith PC_FieldSplitLink head; 320971522cSBarry Smith } PC_FieldSplit; 330971522cSBarry Smith 3416913363SBarry Smith /* 3516913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 3616913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 3716913363SBarry Smith PC you could change this. 3816913363SBarry Smith */ 390971522cSBarry Smith #undef __FUNCT__ 400971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 410971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 420971522cSBarry Smith { 430971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 440971522cSBarry Smith PetscErrorCode ierr; 450971522cSBarry Smith PetscTruth iascii; 460971522cSBarry Smith PetscInt i,j; 475a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 480971522cSBarry Smith 490971522cSBarry Smith PetscFunctionBegin; 500971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 510971522cSBarry Smith if (iascii) { 5251f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 5369a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 540971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 550971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 561ab39975SBarry Smith if (ilink->fields) { 570971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 5879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 595a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 6079416396SBarry Smith if (j > 0) { 6179416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 6279416396SBarry Smith } 635a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 640971522cSBarry Smith } 650971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6679416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 671ab39975SBarry Smith } else { 681ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 691ab39975SBarry Smith } 705a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 715a9f2f41SSatish Balay ilink = ilink->next; 720971522cSBarry Smith } 730971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 740971522cSBarry Smith } else { 750971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 760971522cSBarry Smith } 770971522cSBarry Smith PetscFunctionReturn(0); 780971522cSBarry Smith } 790971522cSBarry Smith 800971522cSBarry Smith #undef __FUNCT__ 813b224e63SBarry Smith #define __FUNCT__ "PCView_FieldSplit_Schur" 823b224e63SBarry Smith static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer) 833b224e63SBarry Smith { 843b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 853b224e63SBarry Smith PetscErrorCode ierr; 863b224e63SBarry Smith PetscTruth iascii; 873b224e63SBarry Smith PetscInt i,j; 883b224e63SBarry Smith PC_FieldSplitLink ilink = jac->head; 893b224e63SBarry Smith KSP ksp; 903b224e63SBarry Smith 913b224e63SBarry Smith PetscFunctionBegin; 923b224e63SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 933b224e63SBarry Smith if (iascii) { 943b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr); 953b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Split info:\n");CHKERRQ(ierr); 963b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 973b224e63SBarry Smith for (i=0; i<jac->nsplits; i++) { 983b224e63SBarry Smith if (ilink->fields) { 993b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 1003b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 1013b224e63SBarry Smith for (j=0; j<ilink->nfields; j++) { 1023b224e63SBarry Smith if (j > 0) { 1033b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 1043b224e63SBarry Smith } 1053b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 1063b224e63SBarry Smith } 1073b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 1083b224e63SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 1093b224e63SBarry Smith } else { 1103b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 1113b224e63SBarry Smith } 1123b224e63SBarry Smith ilink = ilink->next; 1133b224e63SBarry Smith } 1143b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr); 1153b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 1163b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1173b224e63SBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 1183b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1193b224e63SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr); 1203b224e63SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1213b224e63SBarry Smith ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr); 1223b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1233b224e63SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1243b224e63SBarry Smith } else { 1253b224e63SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 1263b224e63SBarry Smith } 1273b224e63SBarry Smith PetscFunctionReturn(0); 1283b224e63SBarry Smith } 1293b224e63SBarry Smith 1303b224e63SBarry Smith #undef __FUNCT__ 13169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 13269a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 1330971522cSBarry Smith { 1340971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1350971522cSBarry Smith PetscErrorCode ierr; 1365a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 137d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 138d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 139d32f9abdSBarry Smith char optionname[128]; 1400971522cSBarry Smith 1410971522cSBarry Smith PetscFunctionBegin; 142d32f9abdSBarry Smith if (!ilink) { 143d32f9abdSBarry Smith 144521d7252SBarry Smith if (jac->bs <= 0) { 145704ba839SBarry Smith if (pc->pmat) { 146521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 147704ba839SBarry Smith } else { 148704ba839SBarry Smith jac->bs = 1; 149704ba839SBarry Smith } 150521d7252SBarry Smith } 151d32f9abdSBarry Smith 152d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 153d32f9abdSBarry Smith if (!flg) { 154d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 155d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 156d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 157d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 158d32f9abdSBarry Smith while (PETSC_TRUE) { 159d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 160d32f9abdSBarry Smith nfields = jac->bs; 161d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 162d32f9abdSBarry Smith if (!flg2) break; 163d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 164d32f9abdSBarry Smith flg = PETSC_FALSE; 165d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 166d32f9abdSBarry Smith } 167d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 168d32f9abdSBarry Smith } 169d32f9abdSBarry Smith 170d32f9abdSBarry Smith if (flg) { 171d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 17279416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 17379416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1745a9f2f41SSatish Balay while (ilink) { 1755a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1765a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 177521d7252SBarry Smith } 1785a9f2f41SSatish Balay ilink = ilink->next; 17979416396SBarry Smith } 18097bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 18179416396SBarry Smith for (i=0; i<jac->bs; i++) { 18279416396SBarry Smith if (!fields[i]) { 18379416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 18479416396SBarry Smith } else { 18579416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 18679416396SBarry Smith } 18779416396SBarry Smith } 18868dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 189521d7252SBarry Smith } 190d32f9abdSBarry Smith } 19169a612a9SBarry Smith PetscFunctionReturn(0); 19269a612a9SBarry Smith } 19369a612a9SBarry Smith 19469a612a9SBarry Smith 19569a612a9SBarry Smith #undef __FUNCT__ 19669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 19769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 19869a612a9SBarry Smith { 19969a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 20069a612a9SBarry Smith PetscErrorCode ierr; 2015a9f2f41SSatish Balay PC_FieldSplitLink ilink; 202cf502942SBarry Smith PetscInt i,nsplit,ccsize; 20369a612a9SBarry Smith MatStructure flag = pc->flag; 20468dd23aaSBarry Smith PetscTruth sorted,getsub; 20569a612a9SBarry Smith 20669a612a9SBarry Smith PetscFunctionBegin; 20769a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 20897bbdb24SBarry Smith nsplit = jac->nsplits; 2095a9f2f41SSatish Balay ilink = jac->head; 21097bbdb24SBarry Smith 21197bbdb24SBarry Smith /* get the matrices for each split */ 212704ba839SBarry Smith if (!jac->issetup) { 2131b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 21497bbdb24SBarry Smith 215704ba839SBarry Smith jac->issetup = PETSC_TRUE; 216704ba839SBarry Smith 217704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 21851f519a2SBarry Smith bs = jac->bs; 21997bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 220cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 2211b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 2221b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 2231b9fc7fcSBarry Smith if (jac->defaultsplit) { 224704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 225704ba839SBarry Smith } else if (!ilink->is) { 226ccb205f8SBarry Smith if (ilink->nfields > 1) { 2275a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 2285a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 2291b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 2301b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 2311b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 23297bbdb24SBarry Smith } 23397bbdb24SBarry Smith } 234704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 235ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 236ccb205f8SBarry Smith } else { 237704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 238ccb205f8SBarry Smith } 2393e197d65SBarry Smith } 2403e197d65SBarry Smith ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr); 241704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 2421b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 243704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 244704ba839SBarry Smith ilink = ilink->next; 2451b9fc7fcSBarry Smith } 2461b9fc7fcSBarry Smith } 2471b9fc7fcSBarry Smith 248704ba839SBarry Smith ilink = jac->head; 24997bbdb24SBarry Smith if (!jac->pmat) { 250cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 251cf502942SBarry Smith for (i=0; i<nsplit; i++) { 252704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 253704ba839SBarry Smith ilink = ilink->next; 254cf502942SBarry Smith } 25597bbdb24SBarry Smith } else { 256cf502942SBarry Smith for (i=0; i<nsplit; i++) { 257704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 258704ba839SBarry Smith ilink = ilink->next; 259cf502942SBarry Smith } 26097bbdb24SBarry Smith } 26197bbdb24SBarry Smith 26268dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 26368dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 2643b224e63SBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE && jac->type != PC_COMPOSITE_SCHUR) { 26568dd23aaSBarry Smith ilink = jac->head; 26668dd23aaSBarry Smith if (!jac->Afield) { 26768dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 26868dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 26911755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 27068dd23aaSBarry Smith ilink = ilink->next; 27168dd23aaSBarry Smith } 27268dd23aaSBarry Smith } else { 27368dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 27411755939SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 27568dd23aaSBarry Smith ilink = ilink->next; 27668dd23aaSBarry Smith } 27768dd23aaSBarry Smith } 27868dd23aaSBarry Smith } 27968dd23aaSBarry Smith 2803b224e63SBarry Smith if (jac->type == PC_COMPOSITE_SCHUR) { 2813b224e63SBarry Smith IS ccis; 2823b224e63SBarry Smith PetscInt N; 2833b224e63SBarry Smith if (nsplit != 2) SETERRQ(PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields"); 28468dd23aaSBarry Smith 2853b224e63SBarry Smith /* need to handle case when one is resetting up the preconditioner */ 2863b224e63SBarry Smith if (jac->schur) { 2873b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 2883b224e63SBarry Smith ilink = jac->head; 2893b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 2903b224e63SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr); 2913b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 2923b224e63SBarry Smith ilink = ilink->next; 2933b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 2943b224e63SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr); 2953b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 2963b224e63SBarry Smith ierr = MatSchurComplementUpdate(jac->schur,jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr); 2973b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,pc->flag);CHKERRQ(ierr); 2983b224e63SBarry Smith 2993b224e63SBarry Smith } else { 300*1cee3971SBarry Smith KSP ksp; 3013b224e63SBarry Smith 3023b224e63SBarry Smith /* extract the B and C matrices */ 3033b224e63SBarry Smith ierr = MatGetSize(pc->mat,PETSC_NULL,&N);CHKERRQ(ierr); 3043b224e63SBarry Smith ilink = jac->head; 3053b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 3063b224e63SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr); 3073b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3083b224e63SBarry Smith ilink = ilink->next; 3093b224e63SBarry Smith ierr = ISComplement(ilink->cis,N,&ccis);CHKERRQ(ierr); 3103b224e63SBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,ccis,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr); 3113b224e63SBarry Smith ierr = ISDestroy(ccis);CHKERRQ(ierr); 3123b224e63SBarry Smith ierr = MatCreateSchurComplement(jac->pmat[0],jac->B,jac->C,jac->pmat[1],&jac->schur);CHKERRQ(ierr); 313*1cee3971SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 314*1cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3153b224e63SBarry Smith ierr = MatSetFromOptions(jac->schur);CHKERRQ(ierr); 3163b224e63SBarry Smith 3173b224e63SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr); 318*1cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr); 3193b224e63SBarry Smith ierr = KSPSetOperators(jac->kspschur,jac->schur,jac->schur,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 3203b224e63SBarry Smith ierr = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr); 3213b224e63SBarry Smith ierr = KSPAppendOptionsPrefix(jac->kspschur,"schur_");CHKERRQ(ierr); 3223b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3233b224e63SBarry Smith ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr); 3243b224e63SBarry Smith 3253b224e63SBarry Smith ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr); 3263b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr); 3273b224e63SBarry Smith ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr); 3283b224e63SBarry Smith ilink = jac->head; 3293b224e63SBarry Smith ilink->x = jac->x[0]; ilink->y = jac->y[0]; 3303b224e63SBarry Smith ilink = ilink->next; 3313b224e63SBarry Smith ilink->x = jac->x[1]; ilink->y = jac->y[1]; 3323b224e63SBarry Smith } 3333b224e63SBarry Smith } else { 33497bbdb24SBarry Smith /* set up the individual PCs */ 33597bbdb24SBarry Smith i = 0; 3365a9f2f41SSatish Balay ilink = jac->head; 3375a9f2f41SSatish Balay while (ilink) { 3385a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 3393b224e63SBarry Smith /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */ 3405a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 3415a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 34297bbdb24SBarry Smith i++; 3435a9f2f41SSatish Balay ilink = ilink->next; 3440971522cSBarry Smith } 34597bbdb24SBarry Smith 34697bbdb24SBarry Smith /* create work vectors for each split */ 3471b9fc7fcSBarry Smith if (!jac->x) { 34879416396SBarry Smith Vec xtmp; 34997bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 3505a9f2f41SSatish Balay ilink = jac->head; 35197bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 352906ed7ccSBarry Smith Vec *vl,*vr; 353906ed7ccSBarry Smith 354906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 355906ed7ccSBarry Smith ilink->x = *vr; 356906ed7ccSBarry Smith ilink->y = *vl; 357906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 358906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 3595a9f2f41SSatish Balay jac->x[i] = ilink->x; 3605a9f2f41SSatish Balay jac->y[i] = ilink->y; 3615a9f2f41SSatish Balay ilink = ilink->next; 36297bbdb24SBarry Smith } 3633b224e63SBarry Smith } 3643b224e63SBarry Smith } 3653b224e63SBarry Smith 3663b224e63SBarry Smith 3673b224e63SBarry Smith if (1) { 3683b224e63SBarry Smith Vec xtmp; 3693b224e63SBarry Smith 37079416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 3711b9fc7fcSBarry Smith 3725a9f2f41SSatish Balay ilink = jac->head; 3731b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 3741b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 375704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 3765a9f2f41SSatish Balay ilink = ilink->next; 3771b9fc7fcSBarry Smith } 3781b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 3791b9fc7fcSBarry Smith } 3800971522cSBarry Smith PetscFunctionReturn(0); 3810971522cSBarry Smith } 3820971522cSBarry Smith 3835a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 384ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 385ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 3865a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 387ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 388ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 38979416396SBarry Smith 3900971522cSBarry Smith #undef __FUNCT__ 3913b224e63SBarry Smith #define __FUNCT__ "PCApply_FieldSplit_Schur" 3923b224e63SBarry Smith static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y) 3933b224e63SBarry Smith { 3943b224e63SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3953b224e63SBarry Smith PetscErrorCode ierr; 3963b224e63SBarry Smith KSP ksp; 3973b224e63SBarry Smith PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next; 3983b224e63SBarry Smith 3993b224e63SBarry Smith PetscFunctionBegin; 4003b224e63SBarry Smith ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr); 4013b224e63SBarry Smith 4023b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4033b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4043b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4053b224e63SBarry Smith ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr); 4063b224e63SBarry Smith ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr); 4073b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4083b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4093b224e63SBarry Smith 4103b224e63SBarry Smith ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr); 4113b224e63SBarry Smith ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4123b224e63SBarry Smith ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4133b224e63SBarry Smith 4143b224e63SBarry Smith ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr); 4153b224e63SBarry Smith ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr); 4163b224e63SBarry Smith ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr); 4173b224e63SBarry Smith ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4183b224e63SBarry Smith ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4193b224e63SBarry Smith 4203b224e63SBarry Smith PetscFunctionReturn(0); 4213b224e63SBarry Smith } 4223b224e63SBarry Smith 4233b224e63SBarry Smith #undef __FUNCT__ 4240971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 4250971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 4260971522cSBarry Smith { 4270971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4280971522cSBarry Smith PetscErrorCode ierr; 4295a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4303e197d65SBarry Smith PetscInt bs,cnt; 4310971522cSBarry Smith 4320971522cSBarry Smith PetscFunctionBegin; 43351f519a2SBarry Smith CHKMEMQ; 434e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 43551f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 43651f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 43751f519a2SBarry Smith 43879416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 4391b9fc7fcSBarry Smith if (jac->defaultsplit) { 4400971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 4415a9f2f41SSatish Balay while (ilink) { 4425a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4435a9f2f41SSatish Balay ilink = ilink->next; 4440971522cSBarry Smith } 4450971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 4461b9fc7fcSBarry Smith } else { 447efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4485a9f2f41SSatish Balay while (ilink) { 4495a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4505a9f2f41SSatish Balay ilink = ilink->next; 4511b9fc7fcSBarry Smith } 4521b9fc7fcSBarry Smith } 45316913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 45479416396SBarry Smith if (!jac->w1) { 45579416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 45679416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 45779416396SBarry Smith } 458efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 4595a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 4603e197d65SBarry Smith cnt = 1; 4615a9f2f41SSatish Balay while (ilink->next) { 4625a9f2f41SSatish Balay ilink = ilink->next; 4633e197d65SBarry Smith if (jac->Afield) { 4643e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 4653e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 4663e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 4673e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4683e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4693e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 4703e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4713e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4723e197d65SBarry Smith } else { 4733e197d65SBarry Smith /* compute the residual over the entire vector */ 4741ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 475efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 4765a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 47779416396SBarry Smith } 4783e197d65SBarry Smith } 47951f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 48011755939SBarry Smith cnt -= 2; 48151f519a2SBarry Smith while (ilink->previous) { 48251f519a2SBarry Smith ilink = ilink->previous; 48311755939SBarry Smith if (jac->Afield) { 48411755939SBarry Smith /* compute the residual only over the part of the vector needed */ 48511755939SBarry Smith ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr); 48611755939SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 48711755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48811755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48911755939SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 49011755939SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 49111755939SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 49211755939SBarry Smith } else { 4931ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 49451f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 49551f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 49679416396SBarry Smith } 49751f519a2SBarry Smith } 49811755939SBarry Smith } 49916913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 50051f519a2SBarry Smith CHKMEMQ; 5010971522cSBarry Smith PetscFunctionReturn(0); 5020971522cSBarry Smith } 5030971522cSBarry Smith 504421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 505ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 506ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 507421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 508ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 509ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 510421e10b8SBarry Smith 511421e10b8SBarry Smith #undef __FUNCT__ 512421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 513421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 514421e10b8SBarry Smith { 515421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 516421e10b8SBarry Smith PetscErrorCode ierr; 517421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 518421e10b8SBarry Smith PetscInt bs; 519421e10b8SBarry Smith 520421e10b8SBarry Smith PetscFunctionBegin; 521421e10b8SBarry Smith CHKMEMQ; 522421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 523421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 524421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 525421e10b8SBarry Smith 526421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 527421e10b8SBarry Smith if (jac->defaultsplit) { 528421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 529421e10b8SBarry Smith while (ilink) { 530421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 531421e10b8SBarry Smith ilink = ilink->next; 532421e10b8SBarry Smith } 533421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 534421e10b8SBarry Smith } else { 535421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 536421e10b8SBarry Smith while (ilink) { 537421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 538421e10b8SBarry Smith ilink = ilink->next; 539421e10b8SBarry Smith } 540421e10b8SBarry Smith } 541421e10b8SBarry Smith } else { 542421e10b8SBarry Smith if (!jac->w1) { 543421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 544421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 545421e10b8SBarry Smith } 546421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 547421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 548421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 549421e10b8SBarry Smith while (ilink->next) { 550421e10b8SBarry Smith ilink = ilink->next; 5519989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 552421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 553421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 554421e10b8SBarry Smith } 555421e10b8SBarry Smith while (ilink->previous) { 556421e10b8SBarry Smith ilink = ilink->previous; 5579989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 558421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 559421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 560421e10b8SBarry Smith } 561421e10b8SBarry Smith } else { 562421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 563421e10b8SBarry Smith ilink = ilink->next; 564421e10b8SBarry Smith } 565421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 566421e10b8SBarry Smith while (ilink->previous) { 567421e10b8SBarry Smith ilink = ilink->previous; 5689989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 569421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 570421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 571421e10b8SBarry Smith } 572421e10b8SBarry Smith } 573421e10b8SBarry Smith } 574421e10b8SBarry Smith CHKMEMQ; 575421e10b8SBarry Smith PetscFunctionReturn(0); 576421e10b8SBarry Smith } 577421e10b8SBarry Smith 5780971522cSBarry Smith #undef __FUNCT__ 5790971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 5800971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 5810971522cSBarry Smith { 5820971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5830971522cSBarry Smith PetscErrorCode ierr; 5845a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 5850971522cSBarry Smith 5860971522cSBarry Smith PetscFunctionBegin; 5875a9f2f41SSatish Balay while (ilink) { 5885a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 5895a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 5905a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 5915a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 592704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 593704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 5945a9f2f41SSatish Balay next = ilink->next; 595704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 596704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 5975a9f2f41SSatish Balay ilink = next; 5980971522cSBarry Smith } 59905b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 60097bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 60168dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 60279416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 60379416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 6043b224e63SBarry Smith if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);} 6053b224e63SBarry Smith if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);} 6063b224e63SBarry Smith if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);} 6070971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 6080971522cSBarry Smith PetscFunctionReturn(0); 6090971522cSBarry Smith } 6100971522cSBarry Smith 6110971522cSBarry Smith #undef __FUNCT__ 6120971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 6130971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 6140971522cSBarry Smith { 6151b9fc7fcSBarry Smith PetscErrorCode ierr; 61651f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 6171b9fc7fcSBarry Smith PetscTruth flg; 6181b9fc7fcSBarry Smith char optionname[128]; 6199dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6203b224e63SBarry Smith PCCompositeType ctype; 6211b9fc7fcSBarry Smith 6220971522cSBarry Smith PetscFunctionBegin; 6231b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 62451f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 62551f519a2SBarry Smith if (flg) { 62651f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 62751f519a2SBarry Smith } 628704ba839SBarry Smith 6293b224e63SBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr); 6303b224e63SBarry Smith if (flg) { 6313b224e63SBarry Smith ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr); 6323b224e63SBarry Smith } 633d32f9abdSBarry Smith 634d32f9abdSBarry Smith if (jac->bs > 0) { 635d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 636d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 63751f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 6381b9fc7fcSBarry Smith while (PETSC_TRUE) { 63913f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 64051f519a2SBarry Smith nfields = jac->bs; 6411b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 6421b9fc7fcSBarry Smith if (!flg) break; 6431b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 6441b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 6451b9fc7fcSBarry Smith } 64651f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 647d32f9abdSBarry Smith } 6481b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6490971522cSBarry Smith PetscFunctionReturn(0); 6500971522cSBarry Smith } 6510971522cSBarry Smith 6520971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 6530971522cSBarry Smith 6540971522cSBarry Smith EXTERN_C_BEGIN 6550971522cSBarry Smith #undef __FUNCT__ 6560971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 657dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 6580971522cSBarry Smith { 65997bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 6600971522cSBarry Smith PetscErrorCode ierr; 6615a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 66269a612a9SBarry Smith char prefix[128]; 66351f519a2SBarry Smith PetscInt i; 6640971522cSBarry Smith 6650971522cSBarry Smith PetscFunctionBegin; 6660971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 66751f519a2SBarry Smith for (i=0; i<n; i++) { 66851f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 66951f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 67051f519a2SBarry Smith } 671704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 672704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 6735a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 6745a9f2f41SSatish Balay ilink->nfields = n; 6755a9f2f41SSatish Balay ilink->next = PETSC_NULL; 6767adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 677*1cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 6785a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 67969a612a9SBarry Smith 6807adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 6817adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 68269a612a9SBarry Smith } else { 68313f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 68469a612a9SBarry Smith } 6855a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 6860971522cSBarry Smith 6870971522cSBarry Smith if (!next) { 6885a9f2f41SSatish Balay jac->head = ilink; 68951f519a2SBarry Smith ilink->previous = PETSC_NULL; 6900971522cSBarry Smith } else { 6910971522cSBarry Smith while (next->next) { 6920971522cSBarry Smith next = next->next; 6930971522cSBarry Smith } 6945a9f2f41SSatish Balay next->next = ilink; 69551f519a2SBarry Smith ilink->previous = next; 6960971522cSBarry Smith } 6970971522cSBarry Smith jac->nsplits++; 6980971522cSBarry Smith PetscFunctionReturn(0); 6990971522cSBarry Smith } 7000971522cSBarry Smith EXTERN_C_END 7010971522cSBarry Smith 7020971522cSBarry Smith 7030971522cSBarry Smith EXTERN_C_BEGIN 7040971522cSBarry Smith #undef __FUNCT__ 70569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 706dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 7070971522cSBarry Smith { 7080971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 7090971522cSBarry Smith PetscErrorCode ierr; 7100971522cSBarry Smith PetscInt cnt = 0; 7115a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7120971522cSBarry Smith 7130971522cSBarry Smith PetscFunctionBegin; 71469a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 7155a9f2f41SSatish Balay while (ilink) { 7165a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 7175a9f2f41SSatish Balay ilink = ilink->next; 7180971522cSBarry Smith } 7190971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 7200971522cSBarry Smith *n = jac->nsplits; 7210971522cSBarry Smith PetscFunctionReturn(0); 7220971522cSBarry Smith } 7230971522cSBarry Smith EXTERN_C_END 7240971522cSBarry Smith 725704ba839SBarry Smith EXTERN_C_BEGIN 726704ba839SBarry Smith #undef __FUNCT__ 727704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 728704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 729704ba839SBarry Smith { 730704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 731704ba839SBarry Smith PetscErrorCode ierr; 732704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 733704ba839SBarry Smith char prefix[128]; 734704ba839SBarry Smith 735704ba839SBarry Smith PetscFunctionBegin; 73616913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 7371ab39975SBarry Smith ilink->is = is; 7381ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 739704ba839SBarry Smith ilink->next = PETSC_NULL; 740704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 741*1cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 742704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 743704ba839SBarry Smith 744704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 745704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 746704ba839SBarry Smith } else { 747704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 748704ba839SBarry Smith } 749704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 750704ba839SBarry Smith 751704ba839SBarry Smith if (!next) { 752704ba839SBarry Smith jac->head = ilink; 753704ba839SBarry Smith ilink->previous = PETSC_NULL; 754704ba839SBarry Smith } else { 755704ba839SBarry Smith while (next->next) { 756704ba839SBarry Smith next = next->next; 757704ba839SBarry Smith } 758704ba839SBarry Smith next->next = ilink; 759704ba839SBarry Smith ilink->previous = next; 760704ba839SBarry Smith } 761704ba839SBarry Smith jac->nsplits++; 762704ba839SBarry Smith 763704ba839SBarry Smith PetscFunctionReturn(0); 764704ba839SBarry Smith } 765704ba839SBarry Smith EXTERN_C_END 766704ba839SBarry Smith 7670971522cSBarry Smith #undef __FUNCT__ 7680971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 7690971522cSBarry Smith /*@ 7700971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 7710971522cSBarry Smith 7720971522cSBarry Smith Collective on PC 7730971522cSBarry Smith 7740971522cSBarry Smith Input Parameters: 7750971522cSBarry Smith + pc - the preconditioner context 7760971522cSBarry Smith . n - the number of fields in this split 7770971522cSBarry Smith . fields - the fields in this split 7780971522cSBarry Smith 7790971522cSBarry Smith Level: intermediate 7800971522cSBarry Smith 781d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 782d32f9abdSBarry Smith 783d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 784d32f9abdSBarry 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 785d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 786d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 787d32f9abdSBarry Smith 788d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 7890971522cSBarry Smith 7900971522cSBarry Smith @*/ 791dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 7920971522cSBarry Smith { 7930971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 7940971522cSBarry Smith 7950971522cSBarry Smith PetscFunctionBegin; 7960971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7970971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 7980971522cSBarry Smith if (f) { 7990971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 8000971522cSBarry Smith } 8010971522cSBarry Smith PetscFunctionReturn(0); 8020971522cSBarry Smith } 8030971522cSBarry Smith 8040971522cSBarry Smith #undef __FUNCT__ 805704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 806704ba839SBarry Smith /*@ 807704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 808704ba839SBarry Smith 809704ba839SBarry Smith Collective on PC 810704ba839SBarry Smith 811704ba839SBarry Smith Input Parameters: 812704ba839SBarry Smith + pc - the preconditioner context 813704ba839SBarry Smith . is - the index set that defines the vector elements in this field 814704ba839SBarry Smith 815d32f9abdSBarry Smith 816d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 817d32f9abdSBarry Smith 818704ba839SBarry Smith Level: intermediate 819704ba839SBarry Smith 820704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 821704ba839SBarry Smith 822704ba839SBarry Smith @*/ 823704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 824704ba839SBarry Smith { 825704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 826704ba839SBarry Smith 827704ba839SBarry Smith PetscFunctionBegin; 828704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 829704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 830704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 831704ba839SBarry Smith if (f) { 832704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 833704ba839SBarry Smith } 834704ba839SBarry Smith PetscFunctionReturn(0); 835704ba839SBarry Smith } 836704ba839SBarry Smith 837704ba839SBarry Smith #undef __FUNCT__ 83851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 83951f519a2SBarry Smith /*@ 84051f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 84151f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 84251f519a2SBarry Smith 84351f519a2SBarry Smith Collective on PC 84451f519a2SBarry Smith 84551f519a2SBarry Smith Input Parameters: 84651f519a2SBarry Smith + pc - the preconditioner context 84751f519a2SBarry Smith - bs - the block size 84851f519a2SBarry Smith 84951f519a2SBarry Smith Level: intermediate 85051f519a2SBarry Smith 85151f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 85251f519a2SBarry Smith 85351f519a2SBarry Smith @*/ 85451f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 85551f519a2SBarry Smith { 85651f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 85751f519a2SBarry Smith 85851f519a2SBarry Smith PetscFunctionBegin; 85951f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 86051f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 86151f519a2SBarry Smith if (f) { 86251f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 86351f519a2SBarry Smith } 86451f519a2SBarry Smith PetscFunctionReturn(0); 86551f519a2SBarry Smith } 86651f519a2SBarry Smith 86751f519a2SBarry Smith #undef __FUNCT__ 86869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 8690971522cSBarry Smith /*@C 87069a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 8710971522cSBarry Smith 87269a612a9SBarry Smith Collective on KSP 8730971522cSBarry Smith 8740971522cSBarry Smith Input Parameter: 8750971522cSBarry Smith . pc - the preconditioner context 8760971522cSBarry Smith 8770971522cSBarry Smith Output Parameters: 8780971522cSBarry Smith + n - the number of split 87969a612a9SBarry Smith - pc - the array of KSP contexts 8800971522cSBarry Smith 8810971522cSBarry Smith Note: 882d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 883d32f9abdSBarry Smith (not the KSP just the array that contains them). 8840971522cSBarry Smith 88569a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 8860971522cSBarry Smith 8870971522cSBarry Smith Level: advanced 8880971522cSBarry Smith 8890971522cSBarry Smith .seealso: PCFIELDSPLIT 8900971522cSBarry Smith @*/ 891dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 8920971522cSBarry Smith { 89369a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 8940971522cSBarry Smith 8950971522cSBarry Smith PetscFunctionBegin; 8960971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8970971522cSBarry Smith PetscValidIntPointer(n,2); 89869a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 8990971522cSBarry Smith if (f) { 90069a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 9010971522cSBarry Smith } else { 90269a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 9030971522cSBarry Smith } 9040971522cSBarry Smith PetscFunctionReturn(0); 9050971522cSBarry Smith } 9060971522cSBarry Smith 90779416396SBarry Smith EXTERN_C_BEGIN 90879416396SBarry Smith #undef __FUNCT__ 90979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 910dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 91179416396SBarry Smith { 91279416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 91379416396SBarry Smith 91479416396SBarry Smith PetscFunctionBegin; 91579416396SBarry Smith jac->type = type; 9163b224e63SBarry Smith if (type == PC_COMPOSITE_SCHUR) { 9173b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit_Schur; 9183b224e63SBarry Smith pc->ops->view = PCView_FieldSplit_Schur; 9193b224e63SBarry Smith } else { 9203b224e63SBarry Smith pc->ops->apply = PCApply_FieldSplit; 9213b224e63SBarry Smith pc->ops->view = PCView_FieldSplit; 9223b224e63SBarry Smith } 92379416396SBarry Smith PetscFunctionReturn(0); 92479416396SBarry Smith } 92579416396SBarry Smith EXTERN_C_END 92679416396SBarry Smith 92751f519a2SBarry Smith EXTERN_C_BEGIN 92851f519a2SBarry Smith #undef __FUNCT__ 92951f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 93051f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 93151f519a2SBarry Smith { 93251f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 93351f519a2SBarry Smith 93451f519a2SBarry Smith PetscFunctionBegin; 93551f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 93651f519a2SBarry 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); 93751f519a2SBarry Smith jac->bs = bs; 93851f519a2SBarry Smith PetscFunctionReturn(0); 93951f519a2SBarry Smith } 94051f519a2SBarry Smith EXTERN_C_END 94151f519a2SBarry Smith 94279416396SBarry Smith #undef __FUNCT__ 94379416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 944bc08b0f1SBarry Smith /*@ 94579416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 94679416396SBarry Smith 94779416396SBarry Smith Collective on PC 94879416396SBarry Smith 94979416396SBarry Smith Input Parameter: 95079416396SBarry Smith . pc - the preconditioner context 951ce9499c7SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 95279416396SBarry Smith 95379416396SBarry Smith Options Database Key: 954ce9499c7SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type 95579416396SBarry Smith 95679416396SBarry Smith Level: Developer 95779416396SBarry Smith 95879416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 95979416396SBarry Smith 96079416396SBarry Smith .seealso: PCCompositeSetType() 96179416396SBarry Smith 96279416396SBarry Smith @*/ 963dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 96479416396SBarry Smith { 96579416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 96679416396SBarry Smith 96779416396SBarry Smith PetscFunctionBegin; 96879416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 96979416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 97079416396SBarry Smith if (f) { 97179416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 97279416396SBarry Smith } 97379416396SBarry Smith PetscFunctionReturn(0); 97479416396SBarry Smith } 97579416396SBarry Smith 9760971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 9770971522cSBarry Smith /*MC 978a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 9790971522cSBarry Smith fields or groups of fields 9800971522cSBarry Smith 9810971522cSBarry Smith 9820971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 983d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 9840971522cSBarry Smith 985a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 98669a612a9SBarry Smith and set the options directly on the resulting KSP object 9870971522cSBarry Smith 9880971522cSBarry Smith Level: intermediate 9890971522cSBarry Smith 99079416396SBarry Smith Options Database Keys: 9912e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 9922e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 9932e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 99451f519a2SBarry Smith . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 9952e70eadcSBarry Smith - -pc_splitfield_type <additive,multiplicative> 99679416396SBarry Smith 9973b224e63SBarry Smith - Options prefix for inner solvers when using Schur complement preconditioner are -schurAblock_ and -schur, 9983b224e63SBarry Smith for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields 9993b224e63SBarry Smith 10003b224e63SBarry Smith 1001d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 1002d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 1003d32f9abdSBarry Smith 1004d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 1005d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 1006d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 1007d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 1008d32f9abdSBarry Smith 1009d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 1010d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 1011d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 1012d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 1013d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 1014d32f9abdSBarry Smith 10150971522cSBarry Smith Concepts: physics based preconditioners 10160971522cSBarry Smith 10170971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1018d32f9abdSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS() 10190971522cSBarry Smith M*/ 10200971522cSBarry Smith 10210971522cSBarry Smith EXTERN_C_BEGIN 10220971522cSBarry Smith #undef __FUNCT__ 10230971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 1024dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 10250971522cSBarry Smith { 10260971522cSBarry Smith PetscErrorCode ierr; 10270971522cSBarry Smith PC_FieldSplit *jac; 10280971522cSBarry Smith 10290971522cSBarry Smith PetscFunctionBegin; 103038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 10310971522cSBarry Smith jac->bs = -1; 10320971522cSBarry Smith jac->nsplits = 0; 10333e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 103451f519a2SBarry Smith 10350971522cSBarry Smith pc->data = (void*)jac; 10360971522cSBarry Smith 10370971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 1038421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 10390971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 10400971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 10410971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 10420971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 10430971522cSBarry Smith pc->ops->applyrichardson = 0; 10440971522cSBarry Smith 104569a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 104669a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 10470971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 10480971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 1049704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 1050704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 105179416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 105279416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 105351f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 105451f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 10550971522cSBarry Smith PetscFunctionReturn(0); 10560971522cSBarry Smith } 10570971522cSBarry Smith EXTERN_C_END 10580971522cSBarry Smith 10590971522cSBarry Smith 1060