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 Vec *x; /* solution vectors */ 8111ade9eSPeter Brune Vec *xl; /* solution local vectors */ 9111ade9eSPeter Brune Vec *y; /* step vectors */ 10eaedb033SPeter Brune Vec *b; /* rhs vectors */ 11111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 12111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 13111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 14111ade9eSPeter Brune PCASMType type; /* ASM type */ 15111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 16d728fb7dSPeter Brune PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 17610116beSPeter Brune PetscReal damping; /* damping parameter for updates from the blocks */ 18b20c023fSPeter Brune 19b20c023fSPeter Brune /* logging events */ 20b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 21b20c023fSPeter Brune PetscLogEvent eventsubsolve; 22602bec5dSPeter Brune 23602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 24602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 25eaedb033SPeter Brune } SNES_NASM; 26eaedb033SPeter Brune 27b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 28602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 29b20c023fSPeter Brune 30eaedb033SPeter Brune #undef __FUNCT__ 31eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM" 32eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes) 33eaedb033SPeter Brune { 34eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 35eaedb033SPeter Brune PetscErrorCode ierr; 36eaedb033SPeter Brune PetscInt i; 376e111a19SKarl Rupp 38eaedb033SPeter Brune PetscFunctionBegin; 39eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 40111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 41f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 42111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 43bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 44eaedb033SPeter Brune 45bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 46111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 47111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 48111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 49eaedb033SPeter Brune } 50111ade9eSPeter Brune 51111ade9eSPeter Brune if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 52111ade9eSPeter Brune if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 53111ade9eSPeter Brune if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 54111ade9eSPeter Brune if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 55111ade9eSPeter Brune 56602bec5dSPeter Brune if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 57602bec5dSPeter Brune 58111ade9eSPeter Brune if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 59111ade9eSPeter Brune if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 60111ade9eSPeter Brune if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 61111ade9eSPeter Brune if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 62b20c023fSPeter Brune 63b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 64b20c023fSPeter Brune nasm->eventsubsolve = 0; 65eaedb033SPeter Brune PetscFunctionReturn(0); 66eaedb033SPeter Brune } 67eaedb033SPeter Brune 68eaedb033SPeter Brune #undef __FUNCT__ 69eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM" 70eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes) 71eaedb033SPeter Brune { 72eaedb033SPeter Brune PetscErrorCode ierr; 736e111a19SKarl Rupp 74eaedb033SPeter Brune PetscFunctionBegin; 75eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 7622d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 77eaedb033SPeter Brune PetscFunctionReturn(0); 78eaedb033SPeter Brune } 79eaedb033SPeter Brune 80eaedb033SPeter Brune #undef __FUNCT__ 81111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 820adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 830adebc6cSBarry Smith { 84111ade9eSPeter Brune PetscErrorCode ierr; 85111ade9eSPeter Brune Vec bcs = (Vec)ctx; 866e111a19SKarl Rupp 87111ade9eSPeter Brune PetscFunctionBegin; 88111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 89111ade9eSPeter Brune PetscFunctionReturn(0); 90111ade9eSPeter Brune } 91111ade9eSPeter Brune 92111ade9eSPeter Brune #undef __FUNCT__ 93eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM" 94eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes) 95eaedb033SPeter Brune { 96eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 97eaedb033SPeter Brune PetscErrorCode ierr; 9876857b2aSPeter Brune DM dm,subdm; 99111ade9eSPeter Brune DM *subdms; 100111ade9eSPeter Brune PetscInt i; 101eaedb033SPeter Brune const char *optionsprefix; 102111ade9eSPeter Brune Vec F; 103ed3c11a9SPeter Brune PetscMPIInt size; 104ed3c11a9SPeter Brune KSP ksp; 105ed3c11a9SPeter Brune PC pc; 106eaedb033SPeter Brune 107eaedb033SPeter Brune PetscFunctionBegin; 108eaedb033SPeter Brune if (!nasm->subsnes) { 109eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1100a696f66SPeter Brune if (dm) { 111eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1120298fd71SBarry Smith ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 113ce94432eSBarry Smith if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 114111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 115eaedb033SPeter Brune 116eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 117111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 118111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 119cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 120cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 121cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 122cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 123ed3c11a9SPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 124ed3c11a9SPeter Brune if (size == 1) { 125ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 126ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 127ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 128ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 129ed3c11a9SPeter Brune } 130cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 131111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 132111ade9eSPeter Brune } 133111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 134ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 135ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 136111ade9eSPeter Brune /* allocate the global vectors */ 137e245e0aaSPeter Brune if (!nasm->x) { 138111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 139e245e0aaSPeter Brune ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr); 140e245e0aaSPeter Brune } 141e245e0aaSPeter Brune if (!nasm->xl) { 142111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 143e245e0aaSPeter Brune ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr); 144e245e0aaSPeter Brune } 145e245e0aaSPeter Brune if (!nasm->y) { 146111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 147e245e0aaSPeter Brune ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr); 148e245e0aaSPeter Brune } 149e245e0aaSPeter Brune if (!nasm->b) { 150111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 151e245e0aaSPeter Brune ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr); 152e245e0aaSPeter Brune } 153111ade9eSPeter Brune 154111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1550298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 15676857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 15776857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 15876857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 15976857b2aSPeter Brune if (!nasm->xl[i]) { 160111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 161111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 16276857b2aSPeter Brune } 1630298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 164111ade9eSPeter Brune } 165602bec5dSPeter Brune if (nasm->finaljacobian) { 166602bec5dSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 167602bec5dSPeter Brune if (nasm->fjtype == 2) { 168602bec5dSPeter Brune ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 169602bec5dSPeter Brune } 170602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 171602bec5dSPeter Brune ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 172602bec5dSPeter Brune } 173602bec5dSPeter Brune } 174eaedb033SPeter Brune PetscFunctionReturn(0); 175eaedb033SPeter Brune } 176eaedb033SPeter Brune 177eaedb033SPeter Brune #undef __FUNCT__ 178eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 179eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 180eaedb033SPeter Brune { 181eaedb033SPeter Brune PetscErrorCode ierr; 182111ade9eSPeter Brune PCASMType asmtype; 183b20c023fSPeter Brune PetscBool flg,monflg; 184111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1856e111a19SKarl Rupp 186eaedb033SPeter Brune PetscFunctionBegin; 187111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 188111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 1891aa26658SKarl Rupp if (flg) nasm->type = asmtype; 190b20c023fSPeter Brune flg = PETSC_FALSE; 191b20c023fSPeter Brune monflg = PETSC_TRUE; 192610116beSPeter Brune ierr = PetscOptionsReal("-snes_nasm_damping","Log times for subSNES solves and restriction","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr); 193610116beSPeter Brune if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 194b20c023fSPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 195d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 196602bec5dSPeter Brune ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 197b20c023fSPeter Brune if (flg) { 198b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 199b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 200b20c023fSPeter Brune } 201eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 202eaedb033SPeter Brune PetscFunctionReturn(0); 203eaedb033SPeter Brune } 204eaedb033SPeter Brune 205eaedb033SPeter Brune #undef __FUNCT__ 206eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 207eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 208eaedb033SPeter Brune { 209b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 210b20c023fSPeter Brune PetscErrorCode ierr; 211b20c023fSPeter Brune PetscMPIInt rank; 212b20c023fSPeter Brune PetscInt i,N; 213b20c023fSPeter Brune PetscBool iascii,isstring; 214b20c023fSPeter Brune PetscViewer sviewer; 215ce94432eSBarry Smith MPI_Comm comm; 216b20c023fSPeter Brune 217b20c023fSPeter Brune PetscFunctionBegin; 218ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 219b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 220b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 221b20c023fSPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 222b20c023fSPeter Brune ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr); 223b20c023fSPeter Brune if (iascii) { 224b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr); 225b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Nonlinear Additive Schwarz: restriction/interpolation type - %s\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 226b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Nonlinear Additive Schwarz: subSNES iterations: %D subKSP iterations: %D\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 227b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 228b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 229b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 230b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local SNES objects:\n");CHKERRQ(ierr); 231b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 232b20c023fSPeter Brune if (!rank) { 233b20c023fSPeter Brune for (i=0; i<nasm->n; i++) { 234b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 235b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 236b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 237b20c023fSPeter Brune } 238b20c023fSPeter Brune } 239b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 240b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 241b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 242b20c023fSPeter Brune } else if (isstring) { 243b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 244b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 245b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 246b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 247b20c023fSPeter Brune } 248eaedb033SPeter Brune PetscFunctionReturn(0); 249eaedb033SPeter Brune } 250eaedb033SPeter Brune 251eaedb033SPeter Brune #undef __FUNCT__ 252eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 25376857b2aSPeter Brune /*@ 25476857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 25576857b2aSPeter Brune 25676857b2aSPeter Brune Not Collective 25776857b2aSPeter Brune 25876857b2aSPeter Brune Input Parameters: 25976857b2aSPeter Brune + SNES - the SNES context 26076857b2aSPeter Brune . n - the number of local subdomains 26176857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 26276857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 26376857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 26476857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 26576857b2aSPeter Brune 26676857b2aSPeter Brune Level: intermediate 26776857b2aSPeter Brune 26876857b2aSPeter Brune .keywords: SNES, NASM 26976857b2aSPeter Brune 27076857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 27176857b2aSPeter Brune @*/ 272a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 273a6dfd86eSKarl Rupp { 274eaedb033SPeter Brune PetscErrorCode ierr; 275111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 2766e111a19SKarl Rupp 277eaedb033SPeter Brune PetscFunctionBegin; 278eaedb033SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 279*4b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 280eaedb033SPeter Brune PetscFunctionReturn(0); 281eaedb033SPeter Brune } 282eaedb033SPeter Brune 283eaedb033SPeter Brune #undef __FUNCT__ 284eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 285a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 286a6dfd86eSKarl Rupp { 287eaedb033SPeter Brune PetscInt i; 288eaedb033SPeter Brune PetscErrorCode ierr; 289eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 2906e111a19SKarl Rupp 291eaedb033SPeter Brune PetscFunctionBegin; 292ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 293eaedb033SPeter Brune 294111ade9eSPeter Brune /* tear down the previously set things */ 295111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 296111ade9eSPeter Brune 297eaedb033SPeter Brune nasm->n = n; 298111ade9eSPeter Brune if (oscatter) { 299111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 300eaedb033SPeter Brune } 301111ade9eSPeter Brune if (iscatter) { 302111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 303eaedb033SPeter Brune } 304111ade9eSPeter Brune if (gscatter) { 305111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 306111ade9eSPeter Brune } 307111ade9eSPeter Brune if (oscatter) { 308111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 309eaedb033SPeter Brune for (i=0; i<n; i++) { 310111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 311eaedb033SPeter Brune } 312111ade9eSPeter Brune } 313111ade9eSPeter Brune if (iscatter) { 314111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 315eaedb033SPeter Brune for (i=0; i<n; i++) { 316111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 317eaedb033SPeter Brune } 318eaedb033SPeter Brune } 319111ade9eSPeter Brune if (gscatter) { 320111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 321eaedb033SPeter Brune for (i=0; i<n; i++) { 322111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 323eaedb033SPeter Brune } 324eaedb033SPeter Brune } 325111ade9eSPeter Brune 326eaedb033SPeter Brune if (subsnes) { 327eaedb033SPeter Brune ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 328eaedb033SPeter Brune for (i=0; i<n; i++) { 329eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 330eaedb033SPeter Brune } 331eaedb033SPeter Brune } 332eaedb033SPeter Brune PetscFunctionReturn(0); 333eaedb033SPeter Brune } 334eaedb033SPeter Brune 335eaedb033SPeter Brune #undef __FUNCT__ 33676857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains" 33776857b2aSPeter Brune /*@ 33876857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 33976857b2aSPeter Brune 34076857b2aSPeter Brune Not Collective 34176857b2aSPeter Brune 34276857b2aSPeter Brune Input Parameters: 34376857b2aSPeter Brune . SNES - the SNES context 34476857b2aSPeter Brune 34576857b2aSPeter Brune Output Parameters: 34676857b2aSPeter Brune + n - the number of local subdomains 34776857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 34876857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 34976857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 35076857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 35176857b2aSPeter Brune 35276857b2aSPeter Brune Level: intermediate 35376857b2aSPeter Brune 35476857b2aSPeter Brune .keywords: SNES, NASM 35576857b2aSPeter Brune 35676857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 35776857b2aSPeter Brune @*/ 35876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 35976857b2aSPeter Brune { 36076857b2aSPeter Brune PetscErrorCode ierr; 36176857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 36276857b2aSPeter Brune 36376857b2aSPeter Brune PetscFunctionBegin; 36476857b2aSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 365*4b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 36676857b2aSPeter Brune PetscFunctionReturn(0); 36776857b2aSPeter Brune } 36876857b2aSPeter Brune 36976857b2aSPeter Brune #undef __FUNCT__ 37076857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 37176857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 37276857b2aSPeter Brune { 37376857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 37476857b2aSPeter Brune 37576857b2aSPeter Brune PetscFunctionBegin; 37676857b2aSPeter Brune if (n) *n = nasm->n; 37776857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 37876857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 37976857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 38076857b2aSPeter Brune if (subsnes) *subsnes = nasm->subsnes; 38176857b2aSPeter Brune PetscFunctionReturn(0); 38276857b2aSPeter Brune } 38376857b2aSPeter Brune 38476857b2aSPeter Brune #undef __FUNCT__ 38576857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs" 38676857b2aSPeter Brune /*@ 38776857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 38876857b2aSPeter Brune 38976857b2aSPeter Brune Not Collective 39076857b2aSPeter Brune 39176857b2aSPeter Brune Input Parameters: 39276857b2aSPeter Brune . SNES - the SNES context 39376857b2aSPeter Brune 39476857b2aSPeter Brune Output Parameters: 39576857b2aSPeter Brune + n - the number of local subdomains 39676857b2aSPeter Brune . x - The subdomain solution vector 39776857b2aSPeter Brune . y - The subdomain step vector 39876857b2aSPeter Brune . b - The subdomain RHS vector 39976857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 40076857b2aSPeter Brune 40176857b2aSPeter Brune Level: developer 40276857b2aSPeter Brune 40376857b2aSPeter Brune .keywords: SNES, NASM 40476857b2aSPeter Brune 40576857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 40676857b2aSPeter Brune @*/ 40776857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 40876857b2aSPeter Brune { 40976857b2aSPeter Brune PetscErrorCode ierr; 41076857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 41176857b2aSPeter Brune 41276857b2aSPeter Brune PetscFunctionBegin; 41376857b2aSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",(void (**)(void))&f);CHKERRQ(ierr); 414*4b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 41576857b2aSPeter Brune PetscFunctionReturn(0); 41676857b2aSPeter Brune } 41776857b2aSPeter Brune 41876857b2aSPeter Brune #undef __FUNCT__ 41976857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 42076857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 42176857b2aSPeter Brune { 42276857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 42376857b2aSPeter Brune 42476857b2aSPeter Brune PetscFunctionBegin; 42576857b2aSPeter Brune if (n) *n = nasm->n; 42676857b2aSPeter Brune if (x) *x = nasm->x; 42776857b2aSPeter Brune if (y) *y = nasm->y; 42876857b2aSPeter Brune if (b) *b = nasm->b; 42976857b2aSPeter Brune if (xl) *xl = nasm->xl; 43076857b2aSPeter Brune PetscFunctionReturn(0); 43176857b2aSPeter Brune } 43276857b2aSPeter Brune 433d728fb7dSPeter Brune #undef __FUNCT__ 434d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 435d728fb7dSPeter Brune /*@ 436d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 437d728fb7dSPeter Brune 438d728fb7dSPeter Brune Collective on SNES 439d728fb7dSPeter Brune 440d728fb7dSPeter Brune Input Parameters: 441d728fb7dSPeter Brune + SNES - the SNES context 442d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 443d728fb7dSPeter Brune 444d728fb7dSPeter Brune Level: developer 445d728fb7dSPeter Brune 446d728fb7dSPeter Brune Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 447d728fb7dSPeter Brune is needed at each linear iteration. 448d728fb7dSPeter Brune 449d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN 450d728fb7dSPeter Brune 451d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 452d728fb7dSPeter Brune @*/ 453d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 454d728fb7dSPeter Brune { 455d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 456d728fb7dSPeter Brune PetscErrorCode ierr; 457d728fb7dSPeter Brune 458d728fb7dSPeter Brune PetscFunctionBegin; 459d728fb7dSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",(void (**)(void))&f);CHKERRQ(ierr); 460*4b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 461d728fb7dSPeter Brune PetscFunctionReturn(0); 462d728fb7dSPeter Brune } 463d728fb7dSPeter Brune 464d728fb7dSPeter Brune #undef __FUNCT__ 465d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 466d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 467d728fb7dSPeter Brune { 468d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 469d728fb7dSPeter Brune 470d728fb7dSPeter Brune PetscFunctionBegin; 471d728fb7dSPeter Brune nasm->finaljacobian = flg; 472602bec5dSPeter Brune if (flg) snes->usesksp = PETSC_TRUE; 473d728fb7dSPeter Brune PetscFunctionReturn(0); 474d728fb7dSPeter Brune } 47576857b2aSPeter Brune 47676857b2aSPeter Brune #undef __FUNCT__ 477610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping" 478610116beSPeter Brune /*@ 479610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 480610116beSPeter Brune 481610116beSPeter Brune Logically collective on SNES 482610116beSPeter Brune 483610116beSPeter Brune Input Parameters: 484610116beSPeter Brune + SNES - the SNES context 485610116beSPeter Brune - dmp - damping 486610116beSPeter Brune 487610116beSPeter Brune Level: intermediate 488610116beSPeter Brune 489610116beSPeter Brune .keywords: SNES, NASM, damping 490610116beSPeter Brune 491610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 492610116beSPeter Brune @*/ 493610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 494610116beSPeter Brune { 495610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 496610116beSPeter Brune PetscErrorCode ierr; 497610116beSPeter Brune 498610116beSPeter Brune PetscFunctionBegin; 499610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 500e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 501610116beSPeter Brune PetscFunctionReturn(0); 502610116beSPeter Brune } 503610116beSPeter Brune 504610116beSPeter Brune #undef __FUNCT__ 505610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping_NASM" 506610116beSPeter Brune PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 507610116beSPeter Brune { 508610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 509610116beSPeter Brune 510610116beSPeter Brune PetscFunctionBegin; 511610116beSPeter Brune nasm->damping = dmp; 512610116beSPeter Brune PetscFunctionReturn(0); 513610116beSPeter Brune } 514610116beSPeter Brune 515610116beSPeter Brune #undef __FUNCT__ 516610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping" 517610116beSPeter Brune /*@ 518610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 519610116beSPeter Brune 520610116beSPeter Brune Not Collective 521610116beSPeter Brune 522610116beSPeter Brune Input Parameters: 523610116beSPeter Brune + SNES - the SNES context 524610116beSPeter Brune - dmp - damping 525610116beSPeter Brune 526610116beSPeter Brune Level: intermediate 527610116beSPeter Brune 528610116beSPeter Brune .keywords: SNES, NASM, damping 529610116beSPeter Brune 530610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping() 531610116beSPeter Brune @*/ 532610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 533610116beSPeter Brune { 534610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal*); 535610116beSPeter Brune PetscErrorCode ierr; 536610116beSPeter Brune 537610116beSPeter Brune PetscFunctionBegin; 538610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 539e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 540610116beSPeter Brune PetscFunctionReturn(0); 541610116beSPeter Brune } 542610116beSPeter Brune 543610116beSPeter Brune #undef __FUNCT__ 544610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping_NASM" 545610116beSPeter Brune PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 546610116beSPeter Brune { 547610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 548610116beSPeter Brune 549610116beSPeter Brune PetscFunctionBegin; 550610116beSPeter Brune *dmp = nasm->damping; 551610116beSPeter Brune PetscFunctionReturn(0); 552610116beSPeter Brune } 553610116beSPeter Brune 554610116beSPeter Brune 555610116beSPeter Brune #undef __FUNCT__ 556eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 5570adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 5580adebc6cSBarry Smith { 559eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 560258e1594SPeter Brune SNES subsnes; 561eaedb033SPeter Brune PetscInt i; 562610116beSPeter Brune PetscReal dmp; 563eaedb033SPeter Brune PetscErrorCode ierr; 564111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 565111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 566111ade9eSPeter Brune DM dm,subdm; 5670adebc6cSBarry Smith 568eaedb033SPeter Brune PetscFunctionBegin; 569eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 570610116beSPeter Brune ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 571111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 572b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 573eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 57470c78f05SPeter Brune /* scatter the solution to the local solution */ 57570c78f05SPeter Brune Xlloc = nasm->xl[i]; 57670c78f05SPeter Brune gscat = nasm->gscatter[i]; 57770c78f05SPeter Brune oscat = nasm->oscatter[i]; 57870c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 57970c78f05SPeter Brune if (B) { 58070c78f05SPeter Brune /* scatter the RHS to the local RHS */ 58170c78f05SPeter Brune Bl = nasm->b[i]; 58270c78f05SPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 58370c78f05SPeter Brune } 58470c78f05SPeter Brune } 585b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 586b20c023fSPeter Brune 587b20c023fSPeter Brune 588b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 58970c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 59070c78f05SPeter Brune Xl = nasm->x[i]; 59170c78f05SPeter Brune Xlloc = nasm->xl[i]; 59270c78f05SPeter Brune Yl = nasm->y[i]; 593258e1594SPeter Brune subsnes = nasm->subsnes[i]; 594258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 595111ade9eSPeter Brune iscat = nasm->iscatter[i]; 596111ade9eSPeter Brune oscat = nasm->oscatter[i]; 597111ade9eSPeter Brune gscat = nasm->gscatter[i]; 598ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 59924b7f281SPeter Brune if (B) { 60024b7f281SPeter Brune Bl = nasm->b[i]; 601ed3c11a9SPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 602ed3c11a9SPeter Brune } else Bl = NULL; 603ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 60470c78f05SPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 60570c78f05SPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 606111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 607ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 608ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 609111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 610111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 611111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 612111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 613ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 614eaedb033SPeter Brune } 615ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 616ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 61770c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 61870c78f05SPeter Brune Yl = nasm->y[i]; 61970c78f05SPeter Brune iscat = nasm->iscatter[i]; 62070c78f05SPeter Brune oscat = nasm->oscatter[i]; 62170c78f05SPeter Brune if (nasm->type == PC_ASM_BASIC) { 62270c78f05SPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 62370c78f05SPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 62470c78f05SPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 625ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 62670c78f05SPeter Brune } 627b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 628610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 629eaedb033SPeter Brune PetscFunctionReturn(0); 630eaedb033SPeter Brune } 631eaedb033SPeter Brune 632eaedb033SPeter Brune #undef __FUNCT__ 633d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 634602bec5dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 635d728fb7dSPeter Brune { 636602bec5dSPeter Brune Vec X = Xfinal; 637d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 638d728fb7dSPeter Brune SNES subsnes; 639ca44f815SPeter Brune PetscInt i,lag = 1; 640d728fb7dSPeter Brune PetscErrorCode ierr; 641e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 642d728fb7dSPeter Brune VecScatter oscat,gscat; 643d728fb7dSPeter Brune DM dm,subdm; 644d728fb7dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 645d728fb7dSPeter Brune PetscFunctionBegin; 646602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 647e59f0a30SPeter Brune F = snes->vec_func; 648e59f0a30SPeter Brune if (snes->normtype == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 649d728fb7dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 650d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 651d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 652602bec5dSPeter Brune if (nasm->fjtype != 1) { 653d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 654d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 655d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 656d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 657d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 658602bec5dSPeter Brune } 659d728fb7dSPeter Brune } 660d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 661d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 662e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 663d728fb7dSPeter Brune Xl = nasm->x[i]; 664d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 665d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 666d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 667d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 668602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 669ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 670d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 671602bec5dSPeter Brune if (nasm->fjtype != 1) { 672d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 673d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 674602bec5dSPeter Brune } 675ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 676ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 677602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 6782f543b25SPeter Brune ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr); 679ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 680d728fb7dSPeter Brune ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr); 681d728fb7dSPeter Brune } 682d728fb7dSPeter Brune PetscFunctionReturn(0); 683d728fb7dSPeter Brune } 684d728fb7dSPeter Brune 685d728fb7dSPeter Brune #undef __FUNCT__ 686eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 687eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 688eaedb033SPeter Brune { 689eaedb033SPeter Brune Vec F; 690eaedb033SPeter Brune Vec X; 691eaedb033SPeter Brune Vec B; 692111ade9eSPeter Brune Vec Y; 693eaedb033SPeter Brune PetscInt i; 694ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 695eaedb033SPeter Brune PetscErrorCode ierr; 696eaedb033SPeter Brune SNESNormType normtype; 697d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 698eaedb033SPeter Brune 699eaedb033SPeter Brune PetscFunctionBegin; 700eaedb033SPeter Brune X = snes->vec_sol; 701111ade9eSPeter Brune Y = snes->vec_sol_update; 702eaedb033SPeter Brune F = snes->vec_func; 703eaedb033SPeter Brune B = snes->vec_rhs; 704eaedb033SPeter Brune 705eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 706eaedb033SPeter Brune snes->iter = 0; 707eaedb033SPeter Brune snes->norm = 0.; 708eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 709eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 710eaedb033SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 711eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 712eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 713eaedb033SPeter Brune if (!snes->vec_func_init_set) { 714eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 715eaedb033SPeter Brune if (snes->domainerror) { 716eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 717eaedb033SPeter Brune PetscFunctionReturn(0); 718eaedb033SPeter Brune } 7191aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 720eaedb033SPeter Brune 721eaedb033SPeter Brune /* convergence test */ 722eaedb033SPeter Brune if (!snes->norm_init_set) { 723eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 724189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 725189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 726189a9710SBarry Smith PetscFunctionReturn(0); 727189a9710SBarry Smith } 728eaedb033SPeter Brune } else { 729eaedb033SPeter Brune fnorm = snes->norm_init; 730eaedb033SPeter Brune snes->norm_init_set = PETSC_FALSE; 731eaedb033SPeter Brune } 732eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 733eaedb033SPeter Brune snes->iter = 0; 734eaedb033SPeter Brune snes->norm = fnorm; 735eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 736a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 737eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 738eaedb033SPeter Brune 739eaedb033SPeter Brune /* set parameter for default relative tolerance convergence test */ 740eaedb033SPeter Brune snes->ttol = fnorm*snes->rtol; 741eaedb033SPeter Brune 742eaedb033SPeter Brune /* test convergence */ 743eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 744eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 745eaedb033SPeter Brune } else { 746eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 747a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 748eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 749eaedb033SPeter Brune } 750eaedb033SPeter Brune 751eaedb033SPeter Brune /* Call general purpose update function */ 752eaedb033SPeter Brune if (snes->ops->update) { 753eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 754eaedb033SPeter Brune } 755602bec5dSPeter Brune /* copy the initial solution over for later */ 756602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 757eaedb033SPeter Brune 758eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 759111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 760eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 761eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 762eaedb033SPeter Brune if (snes->domainerror) { 763eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 764d728fb7dSPeter Brune break; 765eaedb033SPeter Brune } 766eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 767189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 768189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 769189a9710SBarry Smith break; 770189a9710SBarry Smith } 771eaedb033SPeter Brune } 772eaedb033SPeter Brune /* Monitor convergence */ 773eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 774eaedb033SPeter Brune snes->iter = i+1; 775eaedb033SPeter Brune snes->norm = fnorm; 776eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 777a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 778eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 779eaedb033SPeter Brune /* Test for convergence */ 780e59f0a30SPeter Brune if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 781d728fb7dSPeter Brune if (snes->reason) break; 782eaedb033SPeter Brune /* Call general purpose update function */ 783e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 784eaedb033SPeter Brune } 785d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 786eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 787eaedb033SPeter Brune if (i == snes->max_its) { 788eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 789eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 790eaedb033SPeter Brune } 7911aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 792eaedb033SPeter Brune PetscFunctionReturn(0); 793eaedb033SPeter Brune } 794eaedb033SPeter Brune 795eaedb033SPeter Brune /*MC 796eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 797eaedb033SPeter Brune 79870c78603SPeter Brune Options Database: 79970c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 80070c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 80170c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 802602bec5dSPeter Brune . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at 80370c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 80470c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 80570c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 80670c78603SPeter Brune 807eaedb033SPeter Brune Level: advanced 808eaedb033SPeter Brune 809eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 810eaedb033SPeter Brune M*/ 811eaedb033SPeter Brune 812eaedb033SPeter Brune #undef __FUNCT__ 813eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 8148cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 815eaedb033SPeter Brune { 816eaedb033SPeter Brune SNES_NASM *nasm; 817eaedb033SPeter Brune PetscErrorCode ierr; 818eaedb033SPeter Brune 819eaedb033SPeter Brune PetscFunctionBegin; 820eaedb033SPeter Brune ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 821eaedb033SPeter Brune snes->data = (void*)nasm; 822eaedb033SPeter Brune 823eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 824eaedb033SPeter Brune nasm->subsnes = 0; 825eaedb033SPeter Brune nasm->x = 0; 826111ade9eSPeter Brune nasm->xl = 0; 827111ade9eSPeter Brune nasm->y = 0; 828eaedb033SPeter Brune nasm->b = 0; 829111ade9eSPeter Brune nasm->oscatter = 0; 830111ade9eSPeter Brune nasm->iscatter = 0; 831111ade9eSPeter Brune nasm->gscatter = 0; 832610116beSPeter Brune nasm->damping = 1.; 833111ade9eSPeter Brune 834111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 835d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 836eaedb033SPeter Brune 837eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 838eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 839eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 840eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 841eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 842eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 843eaedb033SPeter Brune 844eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 845eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 846eaedb033SPeter Brune 847602bec5dSPeter Brune nasm->fjtype = 0; 848602bec5dSPeter Brune nasm->xinit = NULL; 8490298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8500298fd71SBarry Smith nasm->eventsubsolve = 0; 851b20c023fSPeter Brune 852eaedb033SPeter Brune if (!snes->tolerancesset) { 853eaedb033SPeter Brune snes->max_its = 10000; 854eaedb033SPeter Brune snes->max_funcs = 10000; 855eaedb033SPeter Brune } 856eaedb033SPeter Brune 85700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 85800de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C","SNESNASMGetSubdomains_NASM",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 859610116beSPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C","SNESNASMSetDamping_NASM",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 860610116beSPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C","SNESNASMGetDamping_NASM",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 86100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C","SNESNASMGetSubdomainVecs_NASM",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 86200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C","SNESNASMSetComputeFinalJacobian_NASM",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 863eaedb033SPeter Brune PetscFunctionReturn(0); 864eaedb033SPeter Brune } 86599e0435eSBarry Smith 866