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++) { 31bc8c1f72SJose Roman if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 32111ade9eSPeter Brune if (nasm->xl){ ierr = VecDestroy(&nasm->xl[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" 68111ade9eSPeter Brune PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) { 69111ade9eSPeter Brune PetscErrorCode ierr; 70111ade9eSPeter Brune Vec bcs = (Vec)ctx; 71111ade9eSPeter Brune PetscFunctionBegin; 72111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 73111ade9eSPeter Brune PetscFunctionReturn(0); 74111ade9eSPeter Brune } 75111ade9eSPeter Brune 76111ade9eSPeter Brune #undef __FUNCT__ 77eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM" 78eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes) 79eaedb033SPeter Brune { 80eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 81eaedb033SPeter Brune PetscErrorCode ierr; 82*0a696f66SPeter Brune DM dm,ddm; 83111ade9eSPeter Brune DM *subdms; 84111ade9eSPeter Brune PetscInt i; 85eaedb033SPeter Brune const char *optionsprefix; 86111ade9eSPeter Brune Vec F; 87eaedb033SPeter Brune 88eaedb033SPeter Brune PetscFunctionBegin; 89eaedb033SPeter Brune 90eaedb033SPeter Brune if (!nasm->subsnes) { 91eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 92*0a696f66SPeter Brune if (dm) { 93eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 94111ade9eSPeter Brune ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 95*0a696f66SPeter Brune if (!subdms) { 96*0a696f66SPeter Brune ierr = DMCreateDomainDecompositionDM(dm,"default",&ddm);CHKERRQ(ierr); 97*0a696f66SPeter Brune if (!ddm) { 98*0a696f66SPeter Brune SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 99*0a696f66SPeter Brune } 100*0a696f66SPeter Brune ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr); 101*0a696f66SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 102*0a696f66SPeter Brune ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 103*0a696f66SPeter Brune } 104*0a696f66SPeter Brune if (!subdms) { 105*0a696f66SPeter Brune SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 106*0a696f66SPeter Brune } 107111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 108eaedb033SPeter Brune 109eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 110111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 111111ade9eSPeter Brune 112111ade9eSPeter Brune for (i=0;i<nasm->n;i++) { 113cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 114cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 115cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 116cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 117cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 118111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 119111ade9eSPeter Brune } 120111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 121eaedb033SPeter Brune } else { 122eaedb033SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 123eaedb033SPeter Brune } 124eaedb033SPeter Brune } else { 125eaedb033SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 126eaedb033SPeter Brune } 127111ade9eSPeter Brune /* allocate the global vectors */ 128111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 129111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 130111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 131111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 132111ade9eSPeter Brune 133111ade9eSPeter Brune for (i=0;i<nasm->n;i++) { 134111ade9eSPeter Brune DM subdm; 135111ade9eSPeter Brune ierr = SNESGetFunction(nasm->subsnes[i],&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 136111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr); 137111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr); 138111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr); 139111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 140111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 141111ade9eSPeter Brune ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,PETSC_NULL,nasm->xl[i]);CHKERRQ(ierr); 142111ade9eSPeter Brune } 143111ade9eSPeter Brune 144eaedb033SPeter Brune PetscFunctionReturn(0); 145eaedb033SPeter Brune } 146eaedb033SPeter Brune 147eaedb033SPeter Brune #undef __FUNCT__ 148eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 149eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 150eaedb033SPeter Brune { 151eaedb033SPeter Brune PetscErrorCode ierr; 152*0a696f66SPeter Brune DM dm,ddm; 153*0a696f66SPeter Brune char ddm_name[1024]; 154111ade9eSPeter Brune PCASMType asmtype; 155111ade9eSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 156111ade9eSPeter Brune PetscBool flg; 157111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 158eaedb033SPeter Brune PetscFunctionBegin; 159111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 160111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 161111ade9eSPeter Brune if (flg) {nasm->type = asmtype;} 162*0a696f66SPeter Brune ierr = PetscOptionsString("-snes_nasm_decomposition", "Name of the DM defining the composition", "SNESSetDM", ddm_name, ddm_name,1024,&flg);CHKERRQ(ierr); 163*0a696f66SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 164*0a696f66SPeter Brune if (flg) { 165*0a696f66SPeter Brune if (dm) { 166*0a696f66SPeter Brune ierr = DMCreateDomainDecompositionDM(dm, ddm_name, &ddm);CHKERRQ(ierr); 167*0a696f66SPeter Brune if (!ddm) { 168*0a696f66SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Unknown DM decomposition name %s", ddm_name); 169*0a696f66SPeter Brune } 170*0a696f66SPeter Brune ierr = PetscInfo(snes,"Using domain decomposition DM defined using options database\n");CHKERRQ(ierr); 171*0a696f66SPeter Brune ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr); 172*0a696f66SPeter Brune } else { 173*0a696f66SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "No DM to decompose"); 174*0a696f66SPeter Brune } 175*0a696f66SPeter Brune } 176eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 177eaedb033SPeter Brune PetscFunctionReturn(0); 178eaedb033SPeter Brune } 179eaedb033SPeter Brune 180eaedb033SPeter Brune #undef __FUNCT__ 181eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 182eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 183eaedb033SPeter Brune { 184eaedb033SPeter Brune PetscFunctionBegin; 185*0a696f66SPeter Brune 186*0a696f66SPeter Brune 187eaedb033SPeter Brune PetscFunctionReturn(0); 188eaedb033SPeter Brune } 189eaedb033SPeter Brune 190eaedb033SPeter Brune #undef __FUNCT__ 191eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 192111ade9eSPeter Brune PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 193eaedb033SPeter Brune PetscErrorCode ierr; 194111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 195eaedb033SPeter Brune PetscFunctionBegin; 196eaedb033SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 197111ade9eSPeter Brune ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 198eaedb033SPeter Brune PetscFunctionReturn(0); 199eaedb033SPeter Brune } 200eaedb033SPeter Brune 201eaedb033SPeter Brune EXTERN_C_BEGIN 202eaedb033SPeter Brune #undef __FUNCT__ 203eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 204111ade9eSPeter Brune PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 205eaedb033SPeter Brune PetscInt i; 206eaedb033SPeter Brune PetscErrorCode ierr; 207eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 208eaedb033SPeter Brune PetscFunctionBegin; 209eaedb033SPeter Brune if (snes->setupcalled)SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 210eaedb033SPeter Brune 211111ade9eSPeter Brune /* tear down the previously set things */ 212111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 213111ade9eSPeter Brune 214eaedb033SPeter Brune nasm->n = n; 215111ade9eSPeter Brune if (oscatter) { 216111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 217eaedb033SPeter Brune } 218111ade9eSPeter Brune if (iscatter) { 219111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 220eaedb033SPeter Brune } 221111ade9eSPeter Brune if (gscatter) { 222111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 223111ade9eSPeter Brune } 224111ade9eSPeter Brune if (oscatter) { 225111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 226eaedb033SPeter Brune for (i=0; i<n; i++) { 227111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 228eaedb033SPeter Brune } 229111ade9eSPeter Brune } 230111ade9eSPeter Brune if (iscatter) { 231111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 232eaedb033SPeter Brune for (i=0; i<n; i++) { 233111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 234eaedb033SPeter Brune } 235eaedb033SPeter Brune } 236111ade9eSPeter Brune if (gscatter) { 237111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 238eaedb033SPeter Brune for (i=0; i<n; i++) { 239111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 240eaedb033SPeter Brune } 241eaedb033SPeter Brune } 242111ade9eSPeter Brune 243eaedb033SPeter Brune if (subsnes) { 244eaedb033SPeter Brune ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 245eaedb033SPeter Brune for (i=0; i<n; i++) { 246eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 247eaedb033SPeter Brune } 248eaedb033SPeter Brune } 249eaedb033SPeter Brune PetscFunctionReturn(0); 250eaedb033SPeter Brune } 251eaedb033SPeter Brune EXTERN_C_END 252eaedb033SPeter Brune 253eaedb033SPeter Brune #undef __FUNCT__ 254eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 255111ade9eSPeter Brune PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) { 256eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 257258e1594SPeter Brune SNES subsnes; 258eaedb033SPeter Brune PetscInt i; 259eaedb033SPeter Brune PetscErrorCode ierr; 260111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 261111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 262111ade9eSPeter Brune DM dm,subdm; 263eaedb033SPeter Brune PetscFunctionBegin; 264eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 265111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 266eaedb033SPeter Brune for(i=0;i<nasm->n;i++) { 267258e1594SPeter Brune subsnes = nasm->subsnes[i]; 268258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 269111ade9eSPeter Brune iscat = nasm->iscatter[i]; 270111ade9eSPeter Brune oscat = nasm->oscatter[i]; 271111ade9eSPeter Brune gscat = nasm->gscatter[i]; 272111ade9eSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 273111ade9eSPeter Brune 274111ade9eSPeter Brune /* scatter the solution to the local solution */ 275eaedb033SPeter Brune Xl = nasm->x[i]; 276111ade9eSPeter Brune Xlloc = nasm->xl[i]; 277111ade9eSPeter Brune Yl = nasm->y[i]; 278111ade9eSPeter Brune 279111ade9eSPeter Brune ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 280111ade9eSPeter Brune ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 281111ade9eSPeter Brune 282111ade9eSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 283111ade9eSPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 284111ade9eSPeter Brune 285111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 286111ade9eSPeter Brune 287eaedb033SPeter Brune if (B) { 288111ade9eSPeter Brune /* scatter the RHS to the local RHS */ 289eaedb033SPeter Brune Bl = nasm->b[i]; 290111ade9eSPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 291111ade9eSPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 292eaedb033SPeter Brune } else { 293eaedb033SPeter Brune Bl = PETSC_NULL; 294eaedb033SPeter Brune } 295258e1594SPeter Brune ierr = SNESSolve(subsnes,Bl,Yl);CHKERRQ(ierr); 296111ade9eSPeter Brune 297111ade9eSPeter Brune ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr); 298111ade9eSPeter Brune 299111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 300111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 301111ade9eSPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 302111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 303111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 304111ade9eSPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 305111ade9eSPeter Brune } else { 306111ade9eSPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Only basic and interpolate types are supported for SNESNASM"); 307111ade9eSPeter Brune } 308eaedb033SPeter Brune } 309eaedb033SPeter Brune 310cd939e56SPeter Brune ierr = VecAssemblyBegin(Y);CHKERRQ(ierr); 311cd939e56SPeter Brune ierr = VecAssemblyEnd(Y);CHKERRQ(ierr); 312cd939e56SPeter Brune 313111ade9eSPeter Brune ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr); 314111ade9eSPeter Brune 315eaedb033SPeter Brune PetscFunctionReturn(0); 316eaedb033SPeter Brune } 317eaedb033SPeter Brune 318eaedb033SPeter Brune #undef __FUNCT__ 319eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 320eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 321eaedb033SPeter Brune { 322eaedb033SPeter Brune Vec F; 323eaedb033SPeter Brune Vec X; 324eaedb033SPeter Brune Vec B; 325111ade9eSPeter Brune Vec Y; 326eaedb033SPeter Brune PetscInt i; 327eaedb033SPeter Brune PetscReal fnorm; 328eaedb033SPeter Brune PetscErrorCode ierr; 329eaedb033SPeter Brune SNESNormType normtype; 330eaedb033SPeter Brune 331eaedb033SPeter Brune PetscFunctionBegin; 332eaedb033SPeter Brune 333eaedb033SPeter Brune X = snes->vec_sol; 334111ade9eSPeter Brune Y = snes->vec_sol_update; 335eaedb033SPeter Brune F = snes->vec_func; 336eaedb033SPeter Brune B = snes->vec_rhs; 337eaedb033SPeter Brune 338eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 339eaedb033SPeter Brune snes->iter = 0; 340eaedb033SPeter Brune snes->norm = 0.; 341eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 342eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 343eaedb033SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 344eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 345eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 346eaedb033SPeter Brune if (!snes->vec_func_init_set) { 347eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 348eaedb033SPeter Brune if (snes->domainerror) { 349eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 350eaedb033SPeter Brune PetscFunctionReturn(0); 351eaedb033SPeter Brune } 352eaedb033SPeter Brune } else { 353eaedb033SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 354eaedb033SPeter Brune } 355eaedb033SPeter Brune 356eaedb033SPeter Brune /* convergence test */ 357eaedb033SPeter Brune if (!snes->norm_init_set) { 358eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 359eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 360eaedb033SPeter Brune } else { 361eaedb033SPeter Brune fnorm = snes->norm_init; 362eaedb033SPeter Brune snes->norm_init_set = PETSC_FALSE; 363eaedb033SPeter Brune } 364eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 365eaedb033SPeter Brune snes->iter = 0; 366eaedb033SPeter Brune snes->norm = fnorm; 367eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 368eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 369eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 370eaedb033SPeter Brune 371eaedb033SPeter Brune /* set parameter for default relative tolerance convergence test */ 372eaedb033SPeter Brune snes->ttol = fnorm*snes->rtol; 373eaedb033SPeter Brune 374eaedb033SPeter Brune /* test convergence */ 375eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 376eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 377eaedb033SPeter Brune } else { 378eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 379eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 380eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 381eaedb033SPeter Brune } 382eaedb033SPeter Brune 383eaedb033SPeter Brune /* Call general purpose update function */ 384eaedb033SPeter Brune if (snes->ops->update) { 385eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 386eaedb033SPeter Brune } 387eaedb033SPeter Brune 388eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 389111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 390eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 391eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 392eaedb033SPeter Brune if (snes->domainerror) { 393eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 394eaedb033SPeter Brune PetscFunctionReturn(0); 395eaedb033SPeter Brune } 396eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 397eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 398eaedb033SPeter Brune } 399eaedb033SPeter Brune /* Monitor convergence */ 400eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 401eaedb033SPeter Brune snes->iter = i+1; 402eaedb033SPeter Brune snes->norm = fnorm; 403eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 404eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 405eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 406eaedb033SPeter Brune /* Test for convergence */ 407eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 408eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 409eaedb033SPeter Brune /* Call general purpose update function */ 410eaedb033SPeter Brune if (snes->ops->update) { 411eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 412eaedb033SPeter Brune } 413eaedb033SPeter Brune } 414eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 415eaedb033SPeter Brune if (i == snes->max_its) { 416eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 417eaedb033SPeter Brune if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 418eaedb033SPeter Brune } 419eaedb033SPeter Brune } else { 420eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 421eaedb033SPeter Brune } 422eaedb033SPeter Brune PetscFunctionReturn(0); 423eaedb033SPeter Brune } 424eaedb033SPeter Brune 425eaedb033SPeter Brune /*MC 426eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 427eaedb033SPeter Brune 428eaedb033SPeter Brune Level: advanced 429eaedb033SPeter Brune 430eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 431eaedb033SPeter Brune M*/ 432eaedb033SPeter Brune 433eaedb033SPeter Brune EXTERN_C_BEGIN 434eaedb033SPeter Brune #undef __FUNCT__ 435eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 436eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes) 437eaedb033SPeter Brune { 438eaedb033SPeter Brune SNES_NASM *nasm; 439eaedb033SPeter Brune PetscErrorCode ierr; 440eaedb033SPeter Brune 441eaedb033SPeter Brune PetscFunctionBegin; 442eaedb033SPeter Brune 443eaedb033SPeter Brune 444eaedb033SPeter Brune ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 445eaedb033SPeter Brune snes->data = (void*)nasm; 446eaedb033SPeter Brune 447eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 448eaedb033SPeter Brune nasm->subsnes = 0; 449eaedb033SPeter Brune nasm->x = 0; 450111ade9eSPeter Brune nasm->xl = 0; 451111ade9eSPeter Brune nasm->y = 0; 452eaedb033SPeter Brune nasm->b = 0; 453111ade9eSPeter Brune nasm->oscatter = 0; 454111ade9eSPeter Brune nasm->iscatter = 0; 455111ade9eSPeter Brune nasm->gscatter = 0; 456111ade9eSPeter Brune 457111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 458eaedb033SPeter Brune 459eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 460eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 461eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 462eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 463eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 464eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 465eaedb033SPeter Brune 466eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 467eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 468eaedb033SPeter Brune 469eaedb033SPeter Brune if (!snes->tolerancesset) { 470eaedb033SPeter Brune snes->max_its = 10000; 471eaedb033SPeter Brune snes->max_funcs = 10000; 472eaedb033SPeter Brune } 473eaedb033SPeter Brune 474eaedb033SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM", 475eaedb033SPeter Brune SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 476eaedb033SPeter Brune 477eaedb033SPeter Brune PetscFunctionReturn(0); 478eaedb033SPeter Brune } 479eaedb033SPeter Brune EXTERN_C_END 480