1eaedb033SPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2*111ade9eSPeter 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 */ 9*111ade9eSPeter Brune Vec *xl; /* solution local vectors */ 10*111ade9eSPeter Brune Vec *y; /* step vectors */ 11eaedb033SPeter Brune Vec *b; /* rhs vectors */ 12eaedb033SPeter Brune 13*111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 14*111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 15*111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 16*111ade9eSPeter Brune 17*111ade9eSPeter Brune PCASMType type; /* ASM type */ 18*111ade9eSPeter Brune 19*111ade9eSPeter 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); } 32*111ade9eSPeter Brune if (nasm->xl){ ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 33*111ade9eSPeter 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); } 37*111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 38*111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 39*111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 40eaedb033SPeter Brune } 41*111ade9eSPeter Brune 42*111ade9eSPeter Brune if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 43*111ade9eSPeter Brune if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 44*111ade9eSPeter Brune if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 45*111ade9eSPeter Brune if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 46*111ade9eSPeter Brune 47*111ade9eSPeter Brune if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 48*111ade9eSPeter Brune if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 49*111ade9eSPeter Brune if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 50*111ade9eSPeter Brune if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 51*111ade9eSPeter 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); 62eaedb033SPeter Brune ierr = PetscFree(snes->data); 63eaedb033SPeter Brune PetscFunctionReturn(0); 64eaedb033SPeter Brune } 65eaedb033SPeter Brune 66eaedb033SPeter Brune #undef __FUNCT__ 67*111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 68*111ade9eSPeter Brune PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) { 69*111ade9eSPeter Brune PetscErrorCode ierr; 70*111ade9eSPeter Brune Vec bcs = (Vec)ctx; 71*111ade9eSPeter Brune PetscFunctionBegin; 72*111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 73*111ade9eSPeter Brune PetscFunctionReturn(0); 74*111ade9eSPeter Brune } 75*111ade9eSPeter Brune 76*111ade9eSPeter 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; 82eaedb033SPeter Brune DM dm; 83*111ade9eSPeter Brune DM *subdms; 84*111ade9eSPeter Brune PetscInt i; 85eaedb033SPeter Brune const char *optionsprefix; 86*111ade9eSPeter Brune Vec F; 87eaedb033SPeter Brune 88eaedb033SPeter Brune PetscFunctionBegin; 89eaedb033SPeter Brune 90eaedb033SPeter Brune if (!nasm->subsnes) { 91eaedb033SPeter Brune if (snes->dm) { 92eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 93eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 94*111ade9eSPeter Brune /* create the subdomain dms on this proc */ 95*111ade9eSPeter Brune ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 96*111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 97eaedb033SPeter Brune 98eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 99*111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 100*111ade9eSPeter Brune 101*111ade9eSPeter Brune for (i=0;i<nasm->n;i++) { 102*111ade9eSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[0]);CHKERRQ(ierr); 103eaedb033SPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[0],optionsprefix);CHKERRQ(ierr); 104eaedb033SPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[0],"sub_");CHKERRQ(ierr); 105*111ade9eSPeter Brune ierr = SNESSetDM(nasm->subsnes[0],subdms[i]);CHKERRQ(ierr); 106eaedb033SPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[0]);CHKERRQ(ierr); 107*111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 108*111ade9eSPeter Brune } 109*111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 110eaedb033SPeter Brune } else { 111eaedb033SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 112eaedb033SPeter Brune } 113eaedb033SPeter Brune } else { 114eaedb033SPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 115eaedb033SPeter Brune } 116*111ade9eSPeter Brune /* allocate the global vectors */ 117*111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 118*111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 119*111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 120*111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 121*111ade9eSPeter Brune 122*111ade9eSPeter Brune for (i=0;i<nasm->n;i++) { 123*111ade9eSPeter Brune DM subdm; 124*111ade9eSPeter Brune ierr = SNESGetFunction(nasm->subsnes[i],&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 125*111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr); 126*111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr); 127*111ade9eSPeter Brune ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr); 128*111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 129*111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 130*111ade9eSPeter Brune ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,PETSC_NULL,nasm->xl[i]);CHKERRQ(ierr); 131*111ade9eSPeter Brune } 132*111ade9eSPeter Brune 133eaedb033SPeter Brune PetscFunctionReturn(0); 134eaedb033SPeter Brune } 135eaedb033SPeter Brune 136eaedb033SPeter Brune #undef __FUNCT__ 137eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 138eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 139eaedb033SPeter Brune { 140eaedb033SPeter Brune PetscErrorCode ierr; 141*111ade9eSPeter Brune PCASMType asmtype; 142*111ade9eSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 143*111ade9eSPeter Brune PetscBool flg; 144*111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 145eaedb033SPeter Brune PetscFunctionBegin; 146*111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 147*111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 148*111ade9eSPeter Brune if (flg) {nasm->type = asmtype;} 149eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 150eaedb033SPeter Brune PetscFunctionReturn(0); 151eaedb033SPeter Brune } 152eaedb033SPeter Brune 153eaedb033SPeter Brune #undef __FUNCT__ 154eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 155eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 156eaedb033SPeter Brune { 157eaedb033SPeter Brune PetscFunctionBegin; 158eaedb033SPeter Brune PetscFunctionReturn(0); 159eaedb033SPeter Brune } 160eaedb033SPeter Brune 161eaedb033SPeter Brune #undef __FUNCT__ 162eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 163*111ade9eSPeter Brune PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 164eaedb033SPeter Brune PetscErrorCode ierr; 165*111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 166eaedb033SPeter Brune PetscFunctionBegin; 167eaedb033SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 168*111ade9eSPeter Brune ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 169eaedb033SPeter Brune PetscFunctionReturn(0); 170eaedb033SPeter Brune } 171eaedb033SPeter Brune 172eaedb033SPeter Brune EXTERN_C_BEGIN 173eaedb033SPeter Brune #undef __FUNCT__ 174eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 175*111ade9eSPeter Brune PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 176eaedb033SPeter Brune PetscInt i; 177eaedb033SPeter Brune PetscErrorCode ierr; 178eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 179eaedb033SPeter Brune PetscFunctionBegin; 180eaedb033SPeter Brune if (snes->setupcalled)SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 181eaedb033SPeter Brune 182*111ade9eSPeter Brune /* tear down the previously set things */ 183*111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 184*111ade9eSPeter Brune 185eaedb033SPeter Brune nasm->n = n; 186*111ade9eSPeter Brune if (oscatter) { 187*111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 188eaedb033SPeter Brune } 189*111ade9eSPeter Brune if (iscatter) { 190*111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 191eaedb033SPeter Brune } 192*111ade9eSPeter Brune if (gscatter) { 193*111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 194*111ade9eSPeter Brune } 195*111ade9eSPeter Brune if (oscatter) { 196*111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 197eaedb033SPeter Brune for (i=0; i<n; i++) { 198*111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 199eaedb033SPeter Brune } 200*111ade9eSPeter Brune } 201*111ade9eSPeter Brune if (iscatter) { 202*111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 203eaedb033SPeter Brune for (i=0; i<n; i++) { 204*111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 205eaedb033SPeter Brune } 206eaedb033SPeter Brune } 207*111ade9eSPeter Brune if (gscatter) { 208*111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 209eaedb033SPeter Brune for (i=0; i<n; i++) { 210*111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 211eaedb033SPeter Brune } 212eaedb033SPeter Brune } 213*111ade9eSPeter Brune 214eaedb033SPeter Brune if (subsnes) { 215eaedb033SPeter Brune ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 216eaedb033SPeter Brune for (i=0; i<n; i++) { 217eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 218eaedb033SPeter Brune } 219eaedb033SPeter Brune } 220eaedb033SPeter Brune PetscFunctionReturn(0); 221eaedb033SPeter Brune } 222eaedb033SPeter Brune EXTERN_C_END 223eaedb033SPeter Brune 224eaedb033SPeter Brune #undef __FUNCT__ 225eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 226*111ade9eSPeter Brune PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) { 227eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 228eaedb033SPeter Brune PetscInt i; 229eaedb033SPeter Brune PetscErrorCode ierr; 230*111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 231*111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 232*111ade9eSPeter Brune DM dm,subdm; 233eaedb033SPeter Brune PetscFunctionBegin; 234eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 235*111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 236eaedb033SPeter Brune for(i=0;i<nasm->n;i++) { 237*111ade9eSPeter Brune ierr = SNESGetDM(snes,&subdm);CHKERRQ(ierr); 238*111ade9eSPeter Brune iscat = nasm->iscatter[i]; 239*111ade9eSPeter Brune oscat = nasm->oscatter[i]; 240*111ade9eSPeter Brune gscat = nasm->gscatter[i]; 241*111ade9eSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 242*111ade9eSPeter Brune 243*111ade9eSPeter Brune /* scatter the solution to the local solution */ 244eaedb033SPeter Brune Xl = nasm->x[i]; 245*111ade9eSPeter Brune Xlloc = nasm->xl[i]; 246*111ade9eSPeter Brune Yl = nasm->y[i]; 247*111ade9eSPeter Brune 248*111ade9eSPeter Brune ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 249*111ade9eSPeter Brune ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 250*111ade9eSPeter Brune 251*111ade9eSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 252*111ade9eSPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 253*111ade9eSPeter Brune 254*111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 255*111ade9eSPeter Brune 256eaedb033SPeter Brune if (B) { 257*111ade9eSPeter Brune /* scatter the RHS to the local RHS */ 258eaedb033SPeter Brune Bl = nasm->b[i]; 259*111ade9eSPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 260*111ade9eSPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 261eaedb033SPeter Brune } else { 262eaedb033SPeter Brune Bl = PETSC_NULL; 263eaedb033SPeter Brune } 264*111ade9eSPeter Brune ierr = SNESSolve(nasm->subsnes[i],Bl,Yl);CHKERRQ(ierr); 265*111ade9eSPeter Brune 266*111ade9eSPeter Brune ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr); 267*111ade9eSPeter Brune 268*111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 269*111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 270*111ade9eSPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 271*111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 272*111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 273*111ade9eSPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 274*111ade9eSPeter Brune } else { 275*111ade9eSPeter Brune SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Only basic and interpolate types are supported for SNESNASM"); 276*111ade9eSPeter Brune } 277eaedb033SPeter Brune } 278eaedb033SPeter Brune 279*111ade9eSPeter Brune ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr); 280*111ade9eSPeter Brune 281eaedb033SPeter Brune PetscFunctionReturn(0); 282eaedb033SPeter Brune } 283eaedb033SPeter Brune 284eaedb033SPeter Brune #undef __FUNCT__ 285eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 286eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 287eaedb033SPeter Brune { 288eaedb033SPeter Brune Vec F; 289eaedb033SPeter Brune Vec X; 290eaedb033SPeter Brune Vec B; 291*111ade9eSPeter Brune Vec Y; 292eaedb033SPeter Brune PetscInt i; 293eaedb033SPeter Brune PetscReal fnorm; 294eaedb033SPeter Brune PetscErrorCode ierr; 295eaedb033SPeter Brune SNESNormType normtype; 296eaedb033SPeter Brune 297eaedb033SPeter Brune PetscFunctionBegin; 298eaedb033SPeter Brune 299eaedb033SPeter Brune X = snes->vec_sol; 300*111ade9eSPeter Brune Y = snes->vec_sol_update; 301eaedb033SPeter Brune F = snes->vec_func; 302eaedb033SPeter Brune B = snes->vec_rhs; 303eaedb033SPeter Brune 304eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 305eaedb033SPeter Brune snes->iter = 0; 306eaedb033SPeter Brune snes->norm = 0.; 307eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 308eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 309eaedb033SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 310eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 311eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 312eaedb033SPeter Brune if (!snes->vec_func_init_set) { 313eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 314eaedb033SPeter Brune if (snes->domainerror) { 315eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 316eaedb033SPeter Brune PetscFunctionReturn(0); 317eaedb033SPeter Brune } 318eaedb033SPeter Brune } else { 319eaedb033SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 320eaedb033SPeter Brune } 321eaedb033SPeter Brune 322eaedb033SPeter Brune /* convergence test */ 323eaedb033SPeter Brune if (!snes->norm_init_set) { 324eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 325eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 326eaedb033SPeter Brune } else { 327eaedb033SPeter Brune fnorm = snes->norm_init; 328eaedb033SPeter Brune snes->norm_init_set = PETSC_FALSE; 329eaedb033SPeter Brune } 330eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 331eaedb033SPeter Brune snes->iter = 0; 332eaedb033SPeter Brune snes->norm = fnorm; 333eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 334eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 335eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 336eaedb033SPeter Brune 337eaedb033SPeter Brune /* set parameter for default relative tolerance convergence test */ 338eaedb033SPeter Brune snes->ttol = fnorm*snes->rtol; 339eaedb033SPeter Brune 340eaedb033SPeter Brune /* test convergence */ 341eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 342eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 343eaedb033SPeter Brune } else { 344eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 345eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 346eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 347eaedb033SPeter Brune } 348eaedb033SPeter Brune 349eaedb033SPeter Brune /* Call general purpose update function */ 350eaedb033SPeter Brune if (snes->ops->update) { 351eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 352eaedb033SPeter Brune } 353eaedb033SPeter Brune 354eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 355*111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 356eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 357eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 358eaedb033SPeter Brune if (snes->domainerror) { 359eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 360eaedb033SPeter Brune PetscFunctionReturn(0); 361eaedb033SPeter Brune } 362eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 363eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 364eaedb033SPeter Brune } 365eaedb033SPeter Brune /* Monitor convergence */ 366eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 367eaedb033SPeter Brune snes->iter = i+1; 368eaedb033SPeter Brune snes->norm = fnorm; 369eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 370eaedb033SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 371eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 372eaedb033SPeter Brune /* Test for convergence */ 373eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 374eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 375eaedb033SPeter Brune /* Call general purpose update function */ 376eaedb033SPeter Brune if (snes->ops->update) { 377eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 378eaedb033SPeter Brune } 379eaedb033SPeter Brune } 380eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 381eaedb033SPeter Brune if (i == snes->max_its) { 382eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 383eaedb033SPeter Brune if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 384eaedb033SPeter Brune } 385eaedb033SPeter Brune } else { 386eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 387eaedb033SPeter Brune } 388eaedb033SPeter Brune PetscFunctionReturn(0); 389eaedb033SPeter Brune } 390eaedb033SPeter Brune 391eaedb033SPeter Brune /*MC 392eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 393eaedb033SPeter Brune 394eaedb033SPeter Brune Level: advanced 395eaedb033SPeter Brune 396eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 397eaedb033SPeter Brune M*/ 398eaedb033SPeter Brune 399eaedb033SPeter Brune EXTERN_C_BEGIN 400eaedb033SPeter Brune #undef __FUNCT__ 401eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 402eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes) 403eaedb033SPeter Brune { 404eaedb033SPeter Brune SNES_NASM *nasm; 405eaedb033SPeter Brune PetscErrorCode ierr; 406eaedb033SPeter Brune 407eaedb033SPeter Brune PetscFunctionBegin; 408eaedb033SPeter Brune 409eaedb033SPeter Brune 410eaedb033SPeter Brune ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 411eaedb033SPeter Brune snes->data = (void*)nasm; 412eaedb033SPeter Brune 413eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 414eaedb033SPeter Brune nasm->subsnes = 0; 415eaedb033SPeter Brune nasm->x = 0; 416*111ade9eSPeter Brune nasm->xl = 0; 417*111ade9eSPeter Brune nasm->y = 0; 418eaedb033SPeter Brune nasm->b = 0; 419*111ade9eSPeter Brune nasm->oscatter = 0; 420*111ade9eSPeter Brune nasm->iscatter = 0; 421*111ade9eSPeter Brune nasm->gscatter = 0; 422*111ade9eSPeter Brune 423*111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 424eaedb033SPeter Brune 425eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 426eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 427eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 428eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 429eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 430eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 431eaedb033SPeter Brune 432eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 433eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 434eaedb033SPeter Brune 435eaedb033SPeter Brune if (!snes->tolerancesset) { 436eaedb033SPeter Brune snes->max_its = 10000; 437eaedb033SPeter Brune snes->max_funcs = 10000; 438eaedb033SPeter Brune } 439eaedb033SPeter Brune 440eaedb033SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM", 441eaedb033SPeter Brune SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 442eaedb033SPeter Brune 443eaedb033SPeter Brune PetscFunctionReturn(0); 444eaedb033SPeter Brune } 445eaedb033SPeter Brune EXTERN_C_END 446