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; 2997bbdb24SBarry Smith PC_FieldSplitLink head; 300971522cSBarry Smith } PC_FieldSplit; 310971522cSBarry Smith 3216913363SBarry Smith /* 3316913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 3416913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 3516913363SBarry Smith PC you could change this. 3616913363SBarry Smith */ 370971522cSBarry Smith #undef __FUNCT__ 380971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 390971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 400971522cSBarry Smith { 410971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 420971522cSBarry Smith PetscErrorCode ierr; 430971522cSBarry Smith PetscTruth iascii; 440971522cSBarry Smith PetscInt i,j; 455a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 460971522cSBarry Smith 470971522cSBarry Smith PetscFunctionBegin; 480971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 490971522cSBarry Smith if (iascii) { 5051f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 5169a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 520971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 530971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 541ab39975SBarry Smith if (ilink->fields) { 550971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 5679416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 575a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 5879416396SBarry Smith if (j > 0) { 5979416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 6079416396SBarry Smith } 615a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 620971522cSBarry Smith } 630971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6479416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 651ab39975SBarry Smith } else { 661ab39975SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr); 671ab39975SBarry Smith } 685a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 695a9f2f41SSatish Balay ilink = ilink->next; 700971522cSBarry Smith } 710971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 720971522cSBarry Smith } else { 730971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 740971522cSBarry Smith } 750971522cSBarry Smith PetscFunctionReturn(0); 760971522cSBarry Smith } 770971522cSBarry Smith 780971522cSBarry Smith #undef __FUNCT__ 7969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 8069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 810971522cSBarry Smith { 820971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 830971522cSBarry Smith PetscErrorCode ierr; 845a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 85d32f9abdSBarry Smith PetscInt i = 0,*ifields,nfields; 86d32f9abdSBarry Smith PetscTruth flg = PETSC_FALSE,*fields,flg2; 87d32f9abdSBarry Smith char optionname[128]; 880971522cSBarry Smith 890971522cSBarry Smith PetscFunctionBegin; 90d32f9abdSBarry Smith if (!ilink) { 91d32f9abdSBarry Smith 92521d7252SBarry Smith if (jac->bs <= 0) { 93704ba839SBarry Smith if (pc->pmat) { 94521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 95704ba839SBarry Smith } else { 96704ba839SBarry Smith jac->bs = 1; 97704ba839SBarry Smith } 98521d7252SBarry Smith } 99d32f9abdSBarry Smith 100d32f9abdSBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 101d32f9abdSBarry Smith if (!flg) { 102d32f9abdSBarry Smith /* Allow user to set fields from command line, if bs was known at the time of PCSetFromOptions_FieldSplit() 103d32f9abdSBarry Smith then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */ 104d32f9abdSBarry Smith flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */ 105d32f9abdSBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr); 106d32f9abdSBarry Smith while (PETSC_TRUE) { 107d32f9abdSBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 108d32f9abdSBarry Smith nfields = jac->bs; 109d32f9abdSBarry Smith ierr = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr); 110d32f9abdSBarry Smith if (!flg2) break; 111d32f9abdSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 112d32f9abdSBarry Smith flg = PETSC_FALSE; 113d32f9abdSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr); 114d32f9abdSBarry Smith } 115d32f9abdSBarry Smith ierr = PetscFree(ifields);CHKERRQ(ierr); 116d32f9abdSBarry Smith } 117d32f9abdSBarry Smith 118d32f9abdSBarry Smith if (flg) { 119d32f9abdSBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 12079416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 12179416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 1225a9f2f41SSatish Balay while (ilink) { 1235a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 1245a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 125521d7252SBarry Smith } 1265a9f2f41SSatish Balay ilink = ilink->next; 12779416396SBarry Smith } 12897bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 12979416396SBarry Smith for (i=0; i<jac->bs; i++) { 13079416396SBarry Smith if (!fields[i]) { 13179416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 13279416396SBarry Smith } else { 13379416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 13479416396SBarry Smith } 13579416396SBarry Smith } 13668dd23aaSBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 137521d7252SBarry Smith } 138d32f9abdSBarry Smith } 13969a612a9SBarry Smith PetscFunctionReturn(0); 14069a612a9SBarry Smith } 14169a612a9SBarry Smith 14269a612a9SBarry Smith 14369a612a9SBarry Smith #undef __FUNCT__ 14469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 14569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 14669a612a9SBarry Smith { 14769a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 14869a612a9SBarry Smith PetscErrorCode ierr; 1495a9f2f41SSatish Balay PC_FieldSplitLink ilink; 150cf502942SBarry Smith PetscInt i,nsplit,ccsize; 15169a612a9SBarry Smith MatStructure flag = pc->flag; 15268dd23aaSBarry Smith PetscTruth sorted,getsub; 15369a612a9SBarry Smith 15469a612a9SBarry Smith PetscFunctionBegin; 15569a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 15697bbdb24SBarry Smith nsplit = jac->nsplits; 1575a9f2f41SSatish Balay ilink = jac->head; 15897bbdb24SBarry Smith 15997bbdb24SBarry Smith /* get the matrices for each split */ 160704ba839SBarry Smith if (!jac->issetup) { 1611b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 16297bbdb24SBarry Smith 163704ba839SBarry Smith jac->issetup = PETSC_TRUE; 164704ba839SBarry Smith 165704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 16651f519a2SBarry Smith bs = jac->bs; 16797bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 168cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 1691b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 1701b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1711b9fc7fcSBarry Smith if (jac->defaultsplit) { 172704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 173704ba839SBarry Smith } else if (!ilink->is) { 174ccb205f8SBarry Smith if (ilink->nfields > 1) { 1755a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 1765a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 1771b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 1781b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 1791b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 18097bbdb24SBarry Smith } 18197bbdb24SBarry Smith } 182704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 183ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 184ccb205f8SBarry Smith } else { 185704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 186ccb205f8SBarry Smith } 187*3e197d65SBarry Smith } 188*3e197d65SBarry Smith ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr); 189*3e197d65SBarry Smith printf("csize %d\n",ilink->csize);CHKERRQ(ierr); 190704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 1911b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 192704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 193*3e197d65SBarry Smith ierr = ISView(ilink->is,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 194*3e197d65SBarry Smith ierr = ISView(ilink->cis,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 195704ba839SBarry Smith ilink = ilink->next; 1961b9fc7fcSBarry Smith } 1971b9fc7fcSBarry Smith } 1981b9fc7fcSBarry Smith 199704ba839SBarry Smith ilink = jac->head; 20097bbdb24SBarry Smith if (!jac->pmat) { 201cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 202cf502942SBarry Smith for (i=0; i<nsplit; i++) { 203704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 204704ba839SBarry Smith ilink = ilink->next; 205cf502942SBarry Smith } 20697bbdb24SBarry Smith } else { 207cf502942SBarry Smith for (i=0; i<nsplit; i++) { 208704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 209704ba839SBarry Smith ilink = ilink->next; 210cf502942SBarry Smith } 21197bbdb24SBarry Smith } 21297bbdb24SBarry Smith 21368dd23aaSBarry Smith /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */ 21468dd23aaSBarry Smith ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr); 21568dd23aaSBarry Smith if (getsub && jac->type != PC_COMPOSITE_ADDITIVE) { 21668dd23aaSBarry Smith ilink = jac->head; 21768dd23aaSBarry Smith if (!jac->Afield) { 21868dd23aaSBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr); 21968dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 22068dd23aaSBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,ilink->csize,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 22168dd23aaSBarry Smith ilink = ilink->next; 22268dd23aaSBarry Smith } 22368dd23aaSBarry Smith } else { 22468dd23aaSBarry Smith for (i=0; i<nsplit; i++) { 22568dd23aaSBarry Smith ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,ilink->csize,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr); 22668dd23aaSBarry Smith ilink = ilink->next; 22768dd23aaSBarry Smith } 22868dd23aaSBarry Smith } 22968dd23aaSBarry Smith } 23068dd23aaSBarry Smith 23168dd23aaSBarry Smith 23297bbdb24SBarry Smith /* set up the individual PCs */ 23397bbdb24SBarry Smith i = 0; 2345a9f2f41SSatish Balay ilink = jac->head; 2355a9f2f41SSatish Balay while (ilink) { 2365a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 2375a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 2385a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 23997bbdb24SBarry Smith i++; 2405a9f2f41SSatish Balay ilink = ilink->next; 2410971522cSBarry Smith } 24297bbdb24SBarry Smith 24397bbdb24SBarry Smith /* create work vectors for each split */ 2441b9fc7fcSBarry Smith if (!jac->x) { 24579416396SBarry Smith Vec xtmp; 24697bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 2475a9f2f41SSatish Balay ilink = jac->head; 24897bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 249906ed7ccSBarry Smith Vec *vl,*vr; 250906ed7ccSBarry Smith 251906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 252906ed7ccSBarry Smith ilink->x = *vr; 253906ed7ccSBarry Smith ilink->y = *vl; 254906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 255906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 2565a9f2f41SSatish Balay jac->x[i] = ilink->x; 2575a9f2f41SSatish Balay jac->y[i] = ilink->y; 2585a9f2f41SSatish Balay ilink = ilink->next; 25997bbdb24SBarry Smith } 26079416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 2611b9fc7fcSBarry Smith 2625a9f2f41SSatish Balay ilink = jac->head; 2631b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 2641b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 265704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 2665a9f2f41SSatish Balay ilink = ilink->next; 2671b9fc7fcSBarry Smith } 2681b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 2691b9fc7fcSBarry Smith } 2700971522cSBarry Smith PetscFunctionReturn(0); 2710971522cSBarry Smith } 2720971522cSBarry Smith 2735a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 274ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 275ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 2765a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 277ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 278ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 27979416396SBarry Smith 2800971522cSBarry Smith #undef __FUNCT__ 2810971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 2820971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 2830971522cSBarry Smith { 2840971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2850971522cSBarry Smith PetscErrorCode ierr; 2865a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 287*3e197d65SBarry Smith PetscInt bs,cnt; 2880971522cSBarry Smith 2890971522cSBarry Smith PetscFunctionBegin; 29051f519a2SBarry Smith CHKMEMQ; 291e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 29251f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 29351f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 29451f519a2SBarry Smith 29579416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 2961b9fc7fcSBarry Smith if (jac->defaultsplit) { 2970971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 2985a9f2f41SSatish Balay while (ilink) { 2995a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 3005a9f2f41SSatish Balay ilink = ilink->next; 3010971522cSBarry Smith } 3020971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 3031b9fc7fcSBarry Smith } else { 304efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 3055a9f2f41SSatish Balay while (ilink) { 3065a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 3075a9f2f41SSatish Balay ilink = ilink->next; 3081b9fc7fcSBarry Smith } 3091b9fc7fcSBarry Smith } 31016913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 31179416396SBarry Smith if (!jac->w1) { 31279416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 31379416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 31479416396SBarry Smith } 315efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 3165a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 317*3e197d65SBarry Smith cnt = 1; 3185a9f2f41SSatish Balay while (ilink->next) { 3195a9f2f41SSatish Balay ilink = ilink->next; 320*3e197d65SBarry Smith if (jac->Afield) { 321*3e197d65SBarry Smith /* compute the residual only over the part of the vector needed */ 322*3e197d65SBarry Smith ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr); 323*3e197d65SBarry Smith ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr); 324*3e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 325*3e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 326*3e197d65SBarry Smith ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 327*3e197d65SBarry Smith ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 328*3e197d65SBarry Smith ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 329*3e197d65SBarry Smith } else { 330*3e197d65SBarry Smith /* compute the residual over the entire vector */ 3311ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 332efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 3335a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 33479416396SBarry Smith } 335*3e197d65SBarry Smith } 33651f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 33751f519a2SBarry Smith while (ilink->previous) { 33851f519a2SBarry Smith ilink = ilink->previous; 3391ab39975SBarry Smith ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr); 34051f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 34151f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 34279416396SBarry Smith } 34351f519a2SBarry Smith } 34416913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 34551f519a2SBarry Smith CHKMEMQ; 3460971522cSBarry Smith PetscFunctionReturn(0); 3470971522cSBarry Smith } 3480971522cSBarry Smith 349421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 350ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 351ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 352421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 353ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 354ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 355421e10b8SBarry Smith 356421e10b8SBarry Smith #undef __FUNCT__ 357421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 358421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 359421e10b8SBarry Smith { 360421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 361421e10b8SBarry Smith PetscErrorCode ierr; 362421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 363421e10b8SBarry Smith PetscInt bs; 364421e10b8SBarry Smith 365421e10b8SBarry Smith PetscFunctionBegin; 366421e10b8SBarry Smith CHKMEMQ; 367421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 368421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 369421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 370421e10b8SBarry Smith 371421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 372421e10b8SBarry Smith if (jac->defaultsplit) { 373421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 374421e10b8SBarry Smith while (ilink) { 375421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 376421e10b8SBarry Smith ilink = ilink->next; 377421e10b8SBarry Smith } 378421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 379421e10b8SBarry Smith } else { 380421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 381421e10b8SBarry Smith while (ilink) { 382421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 383421e10b8SBarry Smith ilink = ilink->next; 384421e10b8SBarry Smith } 385421e10b8SBarry Smith } 386421e10b8SBarry Smith } else { 387421e10b8SBarry Smith if (!jac->w1) { 388421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 389421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 390421e10b8SBarry Smith } 391421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 392421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 393421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 394421e10b8SBarry Smith while (ilink->next) { 395421e10b8SBarry Smith ilink = ilink->next; 3969989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 397421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 398421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 399421e10b8SBarry Smith } 400421e10b8SBarry Smith while (ilink->previous) { 401421e10b8SBarry Smith ilink = ilink->previous; 4029989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 403421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 404421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 405421e10b8SBarry Smith } 406421e10b8SBarry Smith } else { 407421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 408421e10b8SBarry Smith ilink = ilink->next; 409421e10b8SBarry Smith } 410421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 411421e10b8SBarry Smith while (ilink->previous) { 412421e10b8SBarry Smith ilink = ilink->previous; 4139989ab13SBarry Smith ierr = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr); 414421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 415421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 416421e10b8SBarry Smith } 417421e10b8SBarry Smith } 418421e10b8SBarry Smith } 419421e10b8SBarry Smith CHKMEMQ; 420421e10b8SBarry Smith PetscFunctionReturn(0); 421421e10b8SBarry Smith } 422421e10b8SBarry Smith 4230971522cSBarry Smith #undef __FUNCT__ 4240971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 4250971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 4260971522cSBarry Smith { 4270971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4280971522cSBarry Smith PetscErrorCode ierr; 4295a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 4300971522cSBarry Smith 4310971522cSBarry Smith PetscFunctionBegin; 4325a9f2f41SSatish Balay while (ilink) { 4335a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 4345a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 4355a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 4365a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 437704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 438704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 4395a9f2f41SSatish Balay next = ilink->next; 440704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 441704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 4425a9f2f41SSatish Balay ilink = next; 4430971522cSBarry Smith } 44405b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 44597bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 44668dd23aaSBarry Smith if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);} 44779416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 44879416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 4490971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 4500971522cSBarry Smith PetscFunctionReturn(0); 4510971522cSBarry Smith } 4520971522cSBarry Smith 4530971522cSBarry Smith #undef __FUNCT__ 4540971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 4550971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 4560971522cSBarry Smith { 4571b9fc7fcSBarry Smith PetscErrorCode ierr; 45851f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 4591b9fc7fcSBarry Smith PetscTruth flg; 4601b9fc7fcSBarry Smith char optionname[128]; 4619dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4621b9fc7fcSBarry Smith 4630971522cSBarry Smith PetscFunctionBegin; 4641b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 46551f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 46651f519a2SBarry Smith if (flg) { 46751f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 46851f519a2SBarry Smith } 469704ba839SBarry Smith 4709dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 471d32f9abdSBarry Smith 472d32f9abdSBarry Smith if (jac->bs > 0) { 473d32f9abdSBarry Smith /* only allow user to set fields from command line if bs is already known. 474d32f9abdSBarry Smith otherwise user can set them in PCFieldSplitSetDefaults() */ 47551f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 4761b9fc7fcSBarry Smith while (PETSC_TRUE) { 47713f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 47851f519a2SBarry Smith nfields = jac->bs; 4791b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 4801b9fc7fcSBarry Smith if (!flg) break; 4811b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 4821b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 4831b9fc7fcSBarry Smith } 48451f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 485d32f9abdSBarry Smith } 4861b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4870971522cSBarry Smith PetscFunctionReturn(0); 4880971522cSBarry Smith } 4890971522cSBarry Smith 4900971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 4910971522cSBarry Smith 4920971522cSBarry Smith EXTERN_C_BEGIN 4930971522cSBarry Smith #undef __FUNCT__ 4940971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 495dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 4960971522cSBarry Smith { 49797bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4980971522cSBarry Smith PetscErrorCode ierr; 4995a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 50069a612a9SBarry Smith char prefix[128]; 50151f519a2SBarry Smith PetscInt i; 5020971522cSBarry Smith 5030971522cSBarry Smith PetscFunctionBegin; 5040971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 50551f519a2SBarry Smith for (i=0; i<n; i++) { 50651f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 50751f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 50851f519a2SBarry Smith } 509704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 510704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 5115a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 5125a9f2f41SSatish Balay ilink->nfields = n; 5135a9f2f41SSatish Balay ilink->next = PETSC_NULL; 5147adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 5155a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 51669a612a9SBarry Smith 5177adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 5187adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 51969a612a9SBarry Smith } else { 52013f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 52169a612a9SBarry Smith } 5225a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 5230971522cSBarry Smith 5240971522cSBarry Smith if (!next) { 5255a9f2f41SSatish Balay jac->head = ilink; 52651f519a2SBarry Smith ilink->previous = PETSC_NULL; 5270971522cSBarry Smith } else { 5280971522cSBarry Smith while (next->next) { 5290971522cSBarry Smith next = next->next; 5300971522cSBarry Smith } 5315a9f2f41SSatish Balay next->next = ilink; 53251f519a2SBarry Smith ilink->previous = next; 5330971522cSBarry Smith } 5340971522cSBarry Smith jac->nsplits++; 5350971522cSBarry Smith PetscFunctionReturn(0); 5360971522cSBarry Smith } 5370971522cSBarry Smith EXTERN_C_END 5380971522cSBarry Smith 5390971522cSBarry Smith 5400971522cSBarry Smith EXTERN_C_BEGIN 5410971522cSBarry Smith #undef __FUNCT__ 54269a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 543dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 5440971522cSBarry Smith { 5450971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 5460971522cSBarry Smith PetscErrorCode ierr; 5470971522cSBarry Smith PetscInt cnt = 0; 5485a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 5490971522cSBarry Smith 5500971522cSBarry Smith PetscFunctionBegin; 55169a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 5525a9f2f41SSatish Balay while (ilink) { 5535a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 5545a9f2f41SSatish Balay ilink = ilink->next; 5550971522cSBarry Smith } 5560971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 5570971522cSBarry Smith *n = jac->nsplits; 5580971522cSBarry Smith PetscFunctionReturn(0); 5590971522cSBarry Smith } 5600971522cSBarry Smith EXTERN_C_END 5610971522cSBarry Smith 562704ba839SBarry Smith EXTERN_C_BEGIN 563704ba839SBarry Smith #undef __FUNCT__ 564704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 565704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 566704ba839SBarry Smith { 567704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 568704ba839SBarry Smith PetscErrorCode ierr; 569704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 570704ba839SBarry Smith char prefix[128]; 571704ba839SBarry Smith 572704ba839SBarry Smith PetscFunctionBegin; 57316913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 5741ab39975SBarry Smith ilink->is = is; 5751ab39975SBarry Smith ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 576704ba839SBarry Smith ilink->next = PETSC_NULL; 577704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 578704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 579704ba839SBarry Smith 580704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 581704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 582704ba839SBarry Smith } else { 583704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 584704ba839SBarry Smith } 585704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 586704ba839SBarry Smith 587704ba839SBarry Smith if (!next) { 588704ba839SBarry Smith jac->head = ilink; 589704ba839SBarry Smith ilink->previous = PETSC_NULL; 590704ba839SBarry Smith } else { 591704ba839SBarry Smith while (next->next) { 592704ba839SBarry Smith next = next->next; 593704ba839SBarry Smith } 594704ba839SBarry Smith next->next = ilink; 595704ba839SBarry Smith ilink->previous = next; 596704ba839SBarry Smith } 597704ba839SBarry Smith jac->nsplits++; 598704ba839SBarry Smith 599704ba839SBarry Smith PetscFunctionReturn(0); 600704ba839SBarry Smith } 601704ba839SBarry Smith EXTERN_C_END 602704ba839SBarry Smith 6030971522cSBarry Smith #undef __FUNCT__ 6040971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 6050971522cSBarry Smith /*@ 6060971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 6070971522cSBarry Smith 6080971522cSBarry Smith Collective on PC 6090971522cSBarry Smith 6100971522cSBarry Smith Input Parameters: 6110971522cSBarry Smith + pc - the preconditioner context 6120971522cSBarry Smith . n - the number of fields in this split 6130971522cSBarry Smith . fields - the fields in this split 6140971522cSBarry Smith 6150971522cSBarry Smith Level: intermediate 6160971522cSBarry Smith 617d32f9abdSBarry Smith Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field. 618d32f9abdSBarry Smith 619d32f9abdSBarry Smith The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block 620d32f9abdSBarry 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 621d32f9abdSBarry Smith 0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x.... 622d32f9abdSBarry Smith where the numbered entries indicate what is in the field. 623d32f9abdSBarry Smith 624d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS() 6250971522cSBarry Smith 6260971522cSBarry Smith @*/ 627dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 6280971522cSBarry Smith { 6290971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 6300971522cSBarry Smith 6310971522cSBarry Smith PetscFunctionBegin; 6320971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6330971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 6340971522cSBarry Smith if (f) { 6350971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 6360971522cSBarry Smith } 6370971522cSBarry Smith PetscFunctionReturn(0); 6380971522cSBarry Smith } 6390971522cSBarry Smith 6400971522cSBarry Smith #undef __FUNCT__ 641704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 642704ba839SBarry Smith /*@ 643704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 644704ba839SBarry Smith 645704ba839SBarry Smith Collective on PC 646704ba839SBarry Smith 647704ba839SBarry Smith Input Parameters: 648704ba839SBarry Smith + pc - the preconditioner context 649704ba839SBarry Smith . is - the index set that defines the vector elements in this field 650704ba839SBarry Smith 651d32f9abdSBarry Smith 652d32f9abdSBarry Smith Notes: Use PCFieldSplitSetFields(), for fields defined by strided types. 653d32f9abdSBarry Smith 654704ba839SBarry Smith Level: intermediate 655704ba839SBarry Smith 656704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 657704ba839SBarry Smith 658704ba839SBarry Smith @*/ 659704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 660704ba839SBarry Smith { 661704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 662704ba839SBarry Smith 663704ba839SBarry Smith PetscFunctionBegin; 664704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 665704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 666704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 667704ba839SBarry Smith if (f) { 668704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 669704ba839SBarry Smith } 670704ba839SBarry Smith PetscFunctionReturn(0); 671704ba839SBarry Smith } 672704ba839SBarry Smith 673704ba839SBarry Smith #undef __FUNCT__ 67451f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 67551f519a2SBarry Smith /*@ 67651f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 67751f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 67851f519a2SBarry Smith 67951f519a2SBarry Smith Collective on PC 68051f519a2SBarry Smith 68151f519a2SBarry Smith Input Parameters: 68251f519a2SBarry Smith + pc - the preconditioner context 68351f519a2SBarry Smith - bs - the block size 68451f519a2SBarry Smith 68551f519a2SBarry Smith Level: intermediate 68651f519a2SBarry Smith 68751f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 68851f519a2SBarry Smith 68951f519a2SBarry Smith @*/ 69051f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 69151f519a2SBarry Smith { 69251f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 69351f519a2SBarry Smith 69451f519a2SBarry Smith PetscFunctionBegin; 69551f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 69651f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 69751f519a2SBarry Smith if (f) { 69851f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 69951f519a2SBarry Smith } 70051f519a2SBarry Smith PetscFunctionReturn(0); 70151f519a2SBarry Smith } 70251f519a2SBarry Smith 70351f519a2SBarry Smith #undef __FUNCT__ 70469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 7050971522cSBarry Smith /*@C 70669a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 7070971522cSBarry Smith 70869a612a9SBarry Smith Collective on KSP 7090971522cSBarry Smith 7100971522cSBarry Smith Input Parameter: 7110971522cSBarry Smith . pc - the preconditioner context 7120971522cSBarry Smith 7130971522cSBarry Smith Output Parameters: 7140971522cSBarry Smith + n - the number of split 71569a612a9SBarry Smith - pc - the array of KSP contexts 7160971522cSBarry Smith 7170971522cSBarry Smith Note: 718d32f9abdSBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user 719d32f9abdSBarry Smith (not the KSP just the array that contains them). 7200971522cSBarry Smith 72169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 7220971522cSBarry Smith 7230971522cSBarry Smith Level: advanced 7240971522cSBarry Smith 7250971522cSBarry Smith .seealso: PCFIELDSPLIT 7260971522cSBarry Smith @*/ 727dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 7280971522cSBarry Smith { 72969a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 7300971522cSBarry Smith 7310971522cSBarry Smith PetscFunctionBegin; 7320971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7330971522cSBarry Smith PetscValidIntPointer(n,2); 73469a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 7350971522cSBarry Smith if (f) { 73669a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 7370971522cSBarry Smith } else { 73869a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 7390971522cSBarry Smith } 7400971522cSBarry Smith PetscFunctionReturn(0); 7410971522cSBarry Smith } 7420971522cSBarry Smith 74379416396SBarry Smith EXTERN_C_BEGIN 74479416396SBarry Smith #undef __FUNCT__ 74579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 746dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 74779416396SBarry Smith { 74879416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 74979416396SBarry Smith 75079416396SBarry Smith PetscFunctionBegin; 75179416396SBarry Smith jac->type = type; 75279416396SBarry Smith PetscFunctionReturn(0); 75379416396SBarry Smith } 75479416396SBarry Smith EXTERN_C_END 75579416396SBarry Smith 75651f519a2SBarry Smith EXTERN_C_BEGIN 75751f519a2SBarry Smith #undef __FUNCT__ 75851f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 75951f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 76051f519a2SBarry Smith { 76151f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 76251f519a2SBarry Smith 76351f519a2SBarry Smith PetscFunctionBegin; 76451f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 76551f519a2SBarry 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); 76651f519a2SBarry Smith jac->bs = bs; 76751f519a2SBarry Smith PetscFunctionReturn(0); 76851f519a2SBarry Smith } 76951f519a2SBarry Smith EXTERN_C_END 77051f519a2SBarry Smith 77179416396SBarry Smith #undef __FUNCT__ 77279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 773bc08b0f1SBarry Smith /*@ 77479416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 77579416396SBarry Smith 77679416396SBarry Smith Collective on PC 77779416396SBarry Smith 77879416396SBarry Smith Input Parameter: 77979416396SBarry Smith . pc - the preconditioner context 780ce9499c7SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 78179416396SBarry Smith 78279416396SBarry Smith Options Database Key: 783ce9499c7SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type 78479416396SBarry Smith 78579416396SBarry Smith Level: Developer 78679416396SBarry Smith 78779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 78879416396SBarry Smith 78979416396SBarry Smith .seealso: PCCompositeSetType() 79079416396SBarry Smith 79179416396SBarry Smith @*/ 792dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 79379416396SBarry Smith { 79479416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 79579416396SBarry Smith 79679416396SBarry Smith PetscFunctionBegin; 79779416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 79879416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 79979416396SBarry Smith if (f) { 80079416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 80179416396SBarry Smith } 80279416396SBarry Smith PetscFunctionReturn(0); 80379416396SBarry Smith } 80479416396SBarry Smith 8050971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 8060971522cSBarry Smith /*MC 807a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 8080971522cSBarry Smith fields or groups of fields 8090971522cSBarry Smith 8100971522cSBarry Smith 8110971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 812d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 8130971522cSBarry Smith 814a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 81569a612a9SBarry Smith and set the options directly on the resulting KSP object 8160971522cSBarry Smith 8170971522cSBarry Smith Level: intermediate 8180971522cSBarry Smith 81979416396SBarry Smith Options Database Keys: 8202e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 8212e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 8222e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 82351f519a2SBarry Smith . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 8242e70eadcSBarry Smith - -pc_splitfield_type <additive,multiplicative> 82579416396SBarry Smith 826d32f9abdSBarry Smith Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS() 827d32f9abdSBarry Smith to define a field by an arbitrary collection of entries. 828d32f9abdSBarry Smith 829d32f9abdSBarry Smith If no fields are set the default is used. The fields are defined by entries strided by bs, 830d32f9abdSBarry Smith beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(), 831d32f9abdSBarry Smith if this is not called the block size defaults to the blocksize of the second matrix passed 832d32f9abdSBarry Smith to KSPSetOperators()/PCSetOperators(). 833d32f9abdSBarry Smith 834d32f9abdSBarry Smith Currently for the multiplicative version, the updated residual needed for the next field 835d32f9abdSBarry Smith solve is computed via a matrix vector product over the entire array. An optimization would be 836d32f9abdSBarry Smith to update the residual only for the part of the right hand side associated with the next field 837d32f9abdSBarry Smith solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the 838d32f9abdSBarry Smith part of the matrix needed to just update part of the residual). 839d32f9abdSBarry Smith 8400971522cSBarry Smith Concepts: physics based preconditioners 8410971522cSBarry Smith 8420971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 843d32f9abdSBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS() 8440971522cSBarry Smith M*/ 8450971522cSBarry Smith 8460971522cSBarry Smith EXTERN_C_BEGIN 8470971522cSBarry Smith #undef __FUNCT__ 8480971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 849dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 8500971522cSBarry Smith { 8510971522cSBarry Smith PetscErrorCode ierr; 8520971522cSBarry Smith PC_FieldSplit *jac; 8530971522cSBarry Smith 8540971522cSBarry Smith PetscFunctionBegin; 85538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 8560971522cSBarry Smith jac->bs = -1; 8570971522cSBarry Smith jac->nsplits = 0; 858*3e197d65SBarry Smith jac->type = PC_COMPOSITE_MULTIPLICATIVE; 85951f519a2SBarry Smith 8600971522cSBarry Smith pc->data = (void*)jac; 8610971522cSBarry Smith 8620971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 863421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 8640971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 8650971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 8660971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 8670971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 8680971522cSBarry Smith pc->ops->applyrichardson = 0; 8690971522cSBarry Smith 87069a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 87169a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 8720971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 8730971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 874704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 875704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 87679416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 87779416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 87851f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 87951f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 8800971522cSBarry Smith PetscFunctionReturn(0); 8810971522cSBarry Smith } 8820971522cSBarry Smith EXTERN_C_END 8830971522cSBarry Smith 8840971522cSBarry Smith 885