1eaedb033SPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2111ade9eSPeter Brune #include <petscdm.h> 3eaedb033SPeter Brune 4eaedb033SPeter Brune typedef struct { 5eaedb033SPeter Brune PetscInt n; /* local subdomains */ 6eaedb033SPeter Brune SNES *subsnes; /* nonlinear solvers for each subdomain */ 7eaedb033SPeter Brune 8eaedb033SPeter Brune Vec *x; /* solution vectors */ 9111ade9eSPeter Brune Vec *xl; /* solution local vectors */ 10111ade9eSPeter Brune Vec *y; /* step vectors */ 11eaedb033SPeter Brune Vec *b; /* rhs vectors */ 12eaedb033SPeter Brune 13111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 14111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 15111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 16111ade9eSPeter Brune 17111ade9eSPeter Brune PCASMType type; /* ASM type */ 18111ade9eSPeter Brune 19111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 20eaedb033SPeter Brune } SNES_NASM; 21eaedb033SPeter Brune 22eaedb033SPeter Brune #undef __FUNCT__ 23eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM" 24eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes) 25eaedb033SPeter Brune { 26eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 27eaedb033SPeter Brune PetscErrorCode ierr; 28eaedb033SPeter Brune PetscInt i; 29eaedb033SPeter Brune PetscFunctionBegin; 30eaedb033SPeter Brune for (i=0;i<nasm->n;i++) { 31111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 32*f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 33111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 34bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 35eaedb033SPeter Brune 36bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 37111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 38111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 39111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 40eaedb033SPeter Brune } 41111ade9eSPeter Brune 42111ade9eSPeter Brune if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 43111ade9eSPeter Brune if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 44111ade9eSPeter Brune if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 45111ade9eSPeter Brune if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 46111ade9eSPeter Brune 47111ade9eSPeter Brune if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 48111ade9eSPeter Brune if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 49111ade9eSPeter Brune if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 50111ade9eSPeter Brune if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 51111ade9eSPeter Brune 52eaedb033SPeter Brune PetscFunctionReturn(0); 53eaedb033SPeter Brune } 54eaedb033SPeter Brune 55eaedb033SPeter Brune #undef __FUNCT__ 56eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM" 57eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes) 58eaedb033SPeter Brune { 59eaedb033SPeter Brune PetscErrorCode ierr; 60eaedb033SPeter Brune PetscFunctionBegin; 61eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 6222d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 63eaedb033SPeter Brune PetscFunctionReturn(0); 64eaedb033SPeter Brune } 65eaedb033SPeter Brune 66eaedb033SPeter Brune #undef __FUNCT__ 67111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 680adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 690adebc6cSBarry Smith { 70111ade9eSPeter Brune PetscErrorCode ierr; 71111ade9eSPeter Brune Vec bcs = (Vec)ctx; 72111ade9eSPeter Brune PetscFunctionBegin; 73111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 74111ade9eSPeter Brune PetscFunctionReturn(0); 75111ade9eSPeter Brune } 76111ade9eSPeter Brune 77111ade9eSPeter Brune #undef __FUNCT__ 78eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM" 79eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes) 80eaedb033SPeter Brune { 81eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 82eaedb033SPeter Brune PetscErrorCode ierr; 830a696f66SPeter Brune DM dm,ddm; 84111ade9eSPeter Brune DM *subdms; 85111ade9eSPeter Brune PetscInt i; 86eaedb033SPeter Brune const char *optionsprefix; 87111ade9eSPeter Brune Vec F; 88eaedb033SPeter Brune 89eaedb033SPeter Brune PetscFunctionBegin; 90eaedb033SPeter Brune 91eaedb033SPeter Brune if (!nasm->subsnes) { 92eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 930a696f66SPeter Brune if (dm) { 94eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 95111ade9eSPeter Brune ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 960a696f66SPeter Brune if (!subdms) { 970a696f66SPeter Brune ierr = DMCreateDomainDecompositionDM(dm,"default",&ddm);CHKERRQ(ierr); 980a696f66SPeter Brune if (!ddm) { 990a696f66SPeter Brune SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 1000a696f66SPeter Brune } 1010a696f66SPeter Brune ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr); 1020a696f66SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1030a696f66SPeter Brune ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 1040a696f66SPeter Brune } 1050a696f66SPeter Brune if (!subdms) { 1060a696f66SPeter Brune SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 1070a696f66SPeter Brune } 108111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 109eaedb033SPeter Brune 110eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 111111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 112111ade9eSPeter Brune 113111ade9eSPeter Brune for (i=0;i<nasm->n;i++) { 114cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 115cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 116cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 117cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 118cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 119111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 120111ade9eSPeter Brune } 121111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 122eaedb033SPeter Brune } else { 123eaedb033SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 124eaedb033SPeter Brune } 125eaedb033SPeter Brune } else { 126eaedb033SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 127eaedb033SPeter Brune } 128111ade9eSPeter Brune /* allocate the global vectors */ 129111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 130111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 131111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 132111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 133111ade9eSPeter Brune 134111ade9eSPeter Brune for (i=0;i<nasm->n;i++) { 135111ade9eSPeter Brune DM subdm; 136111ade9eSPeter Brune ierr = SNESGetFunction(nasm->subsnes[i],&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 137111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr); 138111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr); 139111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr); 140111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 141111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 142111ade9eSPeter Brune ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,PETSC_NULL,nasm->xl[i]);CHKERRQ(ierr); 143111ade9eSPeter Brune } 144111ade9eSPeter Brune 145eaedb033SPeter Brune PetscFunctionReturn(0); 146eaedb033SPeter Brune } 147eaedb033SPeter Brune 148eaedb033SPeter Brune #undef __FUNCT__ 149eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 150eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 151eaedb033SPeter Brune { 152eaedb033SPeter Brune PetscErrorCode ierr; 1530a696f66SPeter Brune DM dm,ddm; 1540a696f66SPeter Brune char ddm_name[1024]; 155111ade9eSPeter Brune PCASMType asmtype; 156111ade9eSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 157111ade9eSPeter Brune PetscBool flg; 158111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 159eaedb033SPeter Brune PetscFunctionBegin; 160111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 161111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 162111ade9eSPeter Brune if (flg) {nasm->type = asmtype;} 1630a696f66SPeter Brune ierr = PetscOptionsString("-snes_nasm_decomposition", "Name of the DM defining the composition", "SNESSetDM", ddm_name, ddm_name,1024,&flg);CHKERRQ(ierr); 1640a696f66SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1650a696f66SPeter Brune if (flg) { 1660a696f66SPeter Brune if (dm) { 1670a696f66SPeter Brune ierr = DMCreateDomainDecompositionDM(dm, ddm_name, &ddm);CHKERRQ(ierr); 1680a696f66SPeter Brune if (!ddm) { 1690a696f66SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Unknown DM decomposition name %s", ddm_name); 1700a696f66SPeter Brune } 1710a696f66SPeter Brune ierr = PetscInfo(snes,"Using domain decomposition DM defined using options database\n");CHKERRQ(ierr); 1720a696f66SPeter Brune ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr); 1730a696f66SPeter Brune } else { 1740a696f66SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "No DM to decompose"); 1750a696f66SPeter Brune } 1760a696f66SPeter Brune } 177eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 178eaedb033SPeter Brune PetscFunctionReturn(0); 179eaedb033SPeter Brune } 180eaedb033SPeter Brune 181eaedb033SPeter Brune #undef __FUNCT__ 182eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 183eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 184eaedb033SPeter Brune { 185eaedb033SPeter Brune PetscFunctionBegin; 1860a696f66SPeter Brune 1870a696f66SPeter Brune 188eaedb033SPeter Brune PetscFunctionReturn(0); 189eaedb033SPeter Brune } 190eaedb033SPeter Brune 191eaedb033SPeter Brune #undef __FUNCT__ 192eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 193111ade9eSPeter Brune PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 194eaedb033SPeter Brune PetscErrorCode ierr; 195111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 196eaedb033SPeter Brune PetscFunctionBegin; 197eaedb033SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 198111ade9eSPeter Brune ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 199eaedb033SPeter Brune PetscFunctionReturn(0); 200eaedb033SPeter Brune } 201eaedb033SPeter Brune 202eaedb033SPeter Brune EXTERN_C_BEGIN 203eaedb033SPeter Brune #undef __FUNCT__ 204eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 205111ade9eSPeter Brune PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 206eaedb033SPeter Brune PetscInt i; 207eaedb033SPeter Brune PetscErrorCode ierr; 208eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 209eaedb033SPeter Brune PetscFunctionBegin; 210eaedb033SPeter Brune if (snes->setupcalled)SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 211eaedb033SPeter Brune 212111ade9eSPeter Brune /* tear down the previously set things */ 213111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 214111ade9eSPeter Brune 215eaedb033SPeter Brune nasm->n = n; 216111ade9eSPeter Brune if (oscatter) { 217111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 218eaedb033SPeter Brune } 219111ade9eSPeter Brune if (iscatter) { 220111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 221eaedb033SPeter Brune } 222111ade9eSPeter Brune if (gscatter) { 223111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 224111ade9eSPeter Brune } 225111ade9eSPeter Brune if (oscatter) { 226111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 227eaedb033SPeter Brune for (i=0; i<n; i++) { 228111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 229eaedb033SPeter Brune } 230111ade9eSPeter Brune } 231111ade9eSPeter Brune if (iscatter) { 232111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 233eaedb033SPeter Brune for (i=0; i<n; i++) { 234111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 235eaedb033SPeter Brune } 236eaedb033SPeter Brune } 237111ade9eSPeter Brune if (gscatter) { 238111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 239eaedb033SPeter Brune for (i=0; i<n; i++) { 240111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 241eaedb033SPeter Brune } 242eaedb033SPeter Brune } 243111ade9eSPeter Brune 244eaedb033SPeter Brune if (subsnes) { 245eaedb033SPeter Brune ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 246eaedb033SPeter Brune for (i=0; i<n; i++) { 247eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 248eaedb033SPeter Brune } 249eaedb033SPeter Brune } 250eaedb033SPeter Brune PetscFunctionReturn(0); 251eaedb033SPeter Brune } 252eaedb033SPeter Brune EXTERN_C_END 253eaedb033SPeter Brune 254eaedb033SPeter Brune #undef __FUNCT__ 255eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 2560adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 2570adebc6cSBarry Smith { 258eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 259258e1594SPeter Brune SNES subsnes; 260eaedb033SPeter Brune PetscInt i; 261eaedb033SPeter Brune PetscErrorCode ierr; 262111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 263111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 264111ade9eSPeter Brune DM dm,subdm; 2650adebc6cSBarry Smith 266eaedb033SPeter Brune PetscFunctionBegin; 267eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 268111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 269eaedb033SPeter Brune for(i=0;i<nasm->n;i++) { 270258e1594SPeter Brune subsnes = nasm->subsnes[i]; 271258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 272111ade9eSPeter Brune iscat = nasm->iscatter[i]; 273111ade9eSPeter Brune oscat = nasm->oscatter[i]; 274111ade9eSPeter Brune gscat = nasm->gscatter[i]; 275111ade9eSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 276111ade9eSPeter Brune 277111ade9eSPeter Brune /* scatter the solution to the local solution */ 278eaedb033SPeter Brune Xl = nasm->x[i]; 279111ade9eSPeter Brune Xlloc = nasm->xl[i]; 280111ade9eSPeter Brune Yl = nasm->y[i]; 281111ade9eSPeter Brune 282111ade9eSPeter Brune ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 283111ade9eSPeter Brune ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 284111ade9eSPeter Brune 285111ade9eSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 286111ade9eSPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 287111ade9eSPeter Brune 288111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 289111ade9eSPeter Brune 290eaedb033SPeter Brune if (B) { 291111ade9eSPeter Brune /* scatter the RHS to the local RHS */ 292eaedb033SPeter Brune Bl = nasm->b[i]; 293111ade9eSPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 294111ade9eSPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 295eaedb033SPeter Brune } else { 296eaedb033SPeter Brune Bl = PETSC_NULL; 297eaedb033SPeter Brune } 298258e1594SPeter Brune ierr = SNESSolve(subsnes,Bl,Yl);CHKERRQ(ierr); 299111ade9eSPeter Brune 300111ade9eSPeter Brune ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr); 301111ade9eSPeter Brune 302111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 303111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 304111ade9eSPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 305111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 306111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 307111ade9eSPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 308111ade9eSPeter Brune } else { 309111ade9eSPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Only basic and interpolate types are supported for SNESNASM"); 310111ade9eSPeter Brune } 311eaedb033SPeter Brune } 312eaedb033SPeter Brune 313cd939e56SPeter Brune ierr = VecAssemblyBegin(Y);CHKERRQ(ierr); 314cd939e56SPeter Brune ierr = VecAssemblyEnd(Y);CHKERRQ(ierr); 315cd939e56SPeter Brune 316111ade9eSPeter Brune ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr); 317111ade9eSPeter Brune 318eaedb033SPeter Brune PetscFunctionReturn(0); 319eaedb033SPeter Brune } 320eaedb033SPeter Brune 321eaedb033SPeter Brune #undef __FUNCT__ 322eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 323eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 324eaedb033SPeter Brune { 325eaedb033SPeter Brune Vec F; 326eaedb033SPeter Brune Vec X; 327eaedb033SPeter Brune Vec B; 328111ade9eSPeter Brune Vec Y; 329eaedb033SPeter Brune PetscInt i; 330eaedb033SPeter Brune PetscReal fnorm; 331eaedb033SPeter Brune PetscErrorCode ierr; 332eaedb033SPeter Brune SNESNormType normtype; 333eaedb033SPeter Brune 334eaedb033SPeter Brune PetscFunctionBegin; 335eaedb033SPeter Brune 336eaedb033SPeter Brune X = snes->vec_sol; 337111ade9eSPeter Brune Y = snes->vec_sol_update; 338eaedb033SPeter Brune F = snes->vec_func; 339eaedb033SPeter Brune B = snes->vec_rhs; 340eaedb033SPeter Brune 341eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 342eaedb033SPeter Brune snes->iter = 0; 343eaedb033SPeter Brune snes->norm = 0.; 344eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 345eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 346eaedb033SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 347eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 348eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 349eaedb033SPeter Brune if (!snes->vec_func_init_set) { 350eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 351eaedb033SPeter Brune if (snes->domainerror) { 352eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 353eaedb033SPeter Brune PetscFunctionReturn(0); 354eaedb033SPeter Brune } 355eaedb033SPeter Brune } else { 356eaedb033SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 357eaedb033SPeter Brune } 358eaedb033SPeter Brune 359eaedb033SPeter Brune /* convergence test */ 360eaedb033SPeter Brune if (!snes->norm_init_set) { 361eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 362eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 363eaedb033SPeter Brune } else { 364eaedb033SPeter Brune fnorm = snes->norm_init; 365eaedb033SPeter Brune snes->norm_init_set = PETSC_FALSE; 366eaedb033SPeter Brune } 367eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 368eaedb033SPeter Brune snes->iter = 0; 369eaedb033SPeter Brune snes->norm = fnorm; 370eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 371eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 372eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 373eaedb033SPeter Brune 374eaedb033SPeter Brune /* set parameter for default relative tolerance convergence test */ 375eaedb033SPeter Brune snes->ttol = fnorm*snes->rtol; 376eaedb033SPeter Brune 377eaedb033SPeter Brune /* test convergence */ 378eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 379eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 380eaedb033SPeter Brune } else { 381eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 382eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 383eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 384eaedb033SPeter Brune } 385eaedb033SPeter Brune 386eaedb033SPeter Brune /* Call general purpose update function */ 387eaedb033SPeter Brune if (snes->ops->update) { 388eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 389eaedb033SPeter Brune } 390eaedb033SPeter Brune 391eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 392111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 393eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 394eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 395eaedb033SPeter Brune if (snes->domainerror) { 396eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 397eaedb033SPeter Brune PetscFunctionReturn(0); 398eaedb033SPeter Brune } 399eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 400eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 401eaedb033SPeter Brune } 402eaedb033SPeter Brune /* Monitor convergence */ 403eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 404eaedb033SPeter Brune snes->iter = i+1; 405eaedb033SPeter Brune snes->norm = fnorm; 406eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 407eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 408eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 409eaedb033SPeter Brune /* Test for convergence */ 410eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 411eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 412eaedb033SPeter Brune /* Call general purpose update function */ 413eaedb033SPeter Brune if (snes->ops->update) { 414eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 415eaedb033SPeter Brune } 416eaedb033SPeter Brune } 417eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 418eaedb033SPeter Brune if (i == snes->max_its) { 419eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 420eaedb033SPeter Brune if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 421eaedb033SPeter Brune } 422eaedb033SPeter Brune } else { 423eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 424eaedb033SPeter Brune } 425eaedb033SPeter Brune PetscFunctionReturn(0); 426eaedb033SPeter Brune } 427eaedb033SPeter Brune 428eaedb033SPeter Brune /*MC 429eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 430eaedb033SPeter Brune 431eaedb033SPeter Brune Level: advanced 432eaedb033SPeter Brune 433eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 434eaedb033SPeter Brune M*/ 435eaedb033SPeter Brune 436eaedb033SPeter Brune EXTERN_C_BEGIN 437eaedb033SPeter Brune #undef __FUNCT__ 438eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 439eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes) 440eaedb033SPeter Brune { 441eaedb033SPeter Brune SNES_NASM *nasm; 442eaedb033SPeter Brune PetscErrorCode ierr; 443eaedb033SPeter Brune 444eaedb033SPeter Brune PetscFunctionBegin; 445eaedb033SPeter Brune 446eaedb033SPeter Brune 447eaedb033SPeter Brune ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 448eaedb033SPeter Brune snes->data = (void*)nasm; 449eaedb033SPeter Brune 450eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 451eaedb033SPeter Brune nasm->subsnes = 0; 452eaedb033SPeter Brune nasm->x = 0; 453111ade9eSPeter Brune nasm->xl = 0; 454111ade9eSPeter Brune nasm->y = 0; 455eaedb033SPeter Brune nasm->b = 0; 456111ade9eSPeter Brune nasm->oscatter = 0; 457111ade9eSPeter Brune nasm->iscatter = 0; 458111ade9eSPeter Brune nasm->gscatter = 0; 459111ade9eSPeter Brune 460111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 461eaedb033SPeter Brune 462eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 463eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 464eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 465eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 466eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 467eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 468eaedb033SPeter Brune 469eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 470eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 471eaedb033SPeter Brune 472eaedb033SPeter Brune if (!snes->tolerancesset) { 473eaedb033SPeter Brune snes->max_its = 10000; 474eaedb033SPeter Brune snes->max_funcs = 10000; 475eaedb033SPeter Brune } 476eaedb033SPeter Brune 477eaedb033SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM", 478eaedb033SPeter Brune SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 479eaedb033SPeter Brune 480eaedb033SPeter Brune PetscFunctionReturn(0); 481eaedb033SPeter Brune } 482eaedb033SPeter Brune EXTERN_C_END 483