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 */ 17b20c023fSPeter Brune 18b20c023fSPeter Brune /* logging events */ 19b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 20b20c023fSPeter Brune PetscLogEvent eventsubsolve; 21eaedb033SPeter Brune } SNES_NASM; 22eaedb033SPeter Brune 23b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 24b20c023fSPeter Brune 25eaedb033SPeter Brune #undef __FUNCT__ 26eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM" 27eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes) 28eaedb033SPeter Brune { 29eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 30eaedb033SPeter Brune PetscErrorCode ierr; 31eaedb033SPeter Brune PetscInt i; 326e111a19SKarl Rupp 33eaedb033SPeter Brune PetscFunctionBegin; 34eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 35111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 36f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 37111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 38bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 39eaedb033SPeter Brune 40bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 41111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 42111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 43111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 44eaedb033SPeter Brune } 45111ade9eSPeter Brune 46111ade9eSPeter Brune if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 47111ade9eSPeter Brune if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 48111ade9eSPeter Brune if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 49111ade9eSPeter Brune if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 50111ade9eSPeter Brune 51111ade9eSPeter Brune if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 52111ade9eSPeter Brune if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 53111ade9eSPeter Brune if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 54111ade9eSPeter Brune if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 55b20c023fSPeter Brune 56b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 57b20c023fSPeter Brune nasm->eventsubsolve = 0; 58eaedb033SPeter Brune PetscFunctionReturn(0); 59eaedb033SPeter Brune } 60eaedb033SPeter Brune 61eaedb033SPeter Brune #undef __FUNCT__ 62eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM" 63eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes) 64eaedb033SPeter Brune { 65eaedb033SPeter Brune PetscErrorCode ierr; 666e111a19SKarl Rupp 67eaedb033SPeter Brune PetscFunctionBegin; 68eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 6922d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 70eaedb033SPeter Brune PetscFunctionReturn(0); 71eaedb033SPeter Brune } 72eaedb033SPeter Brune 73eaedb033SPeter Brune #undef __FUNCT__ 74111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 750adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 760adebc6cSBarry Smith { 77111ade9eSPeter Brune PetscErrorCode ierr; 78111ade9eSPeter Brune Vec bcs = (Vec)ctx; 796e111a19SKarl Rupp 80111ade9eSPeter Brune PetscFunctionBegin; 81111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 82111ade9eSPeter Brune PetscFunctionReturn(0); 83111ade9eSPeter Brune } 84111ade9eSPeter Brune 85111ade9eSPeter Brune #undef __FUNCT__ 86eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM" 87eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes) 88eaedb033SPeter Brune { 89eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 90eaedb033SPeter Brune PetscErrorCode ierr; 9176857b2aSPeter Brune DM dm,subdm; 92111ade9eSPeter Brune DM *subdms; 93111ade9eSPeter Brune PetscInt i; 94eaedb033SPeter Brune const char *optionsprefix; 95111ade9eSPeter Brune Vec F; 96ed3c11a9SPeter Brune PetscMPIInt size; 97ed3c11a9SPeter Brune KSP ksp; 98ed3c11a9SPeter Brune PC pc; 99eaedb033SPeter Brune 100eaedb033SPeter Brune PetscFunctionBegin; 101eaedb033SPeter Brune if (!nasm->subsnes) { 102eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1030a696f66SPeter Brune if (dm) { 104eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1050298fd71SBarry Smith ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 106ce94432eSBarry Smith if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 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 for (i=0; i<nasm->n; i++) { 112cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 113cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 114cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 115cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 116ed3c11a9SPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 117ed3c11a9SPeter Brune if (size == 1) { 118ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 119ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 120ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 121ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 122ed3c11a9SPeter Brune } 123cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 124111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 125111ade9eSPeter Brune } 126111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 127ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 128ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 129111ade9eSPeter Brune /* allocate the global vectors */ 130e245e0aaSPeter Brune if (!nasm->x) { 131111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 132e245e0aaSPeter Brune ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr); 133e245e0aaSPeter Brune } 134e245e0aaSPeter Brune if (!nasm->xl) { 135111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 136e245e0aaSPeter Brune ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr); 137e245e0aaSPeter Brune } 138e245e0aaSPeter Brune if (!nasm->y) { 139111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 140e245e0aaSPeter Brune ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr); 141e245e0aaSPeter Brune } 142e245e0aaSPeter Brune if (!nasm->b) { 143111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 144e245e0aaSPeter Brune ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr); 145e245e0aaSPeter Brune } 146111ade9eSPeter Brune 147111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1480298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 14976857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 15076857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 15176857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 15276857b2aSPeter Brune if (!nasm->xl[i]) { 153111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 154111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 15576857b2aSPeter Brune } 1560298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 157111ade9eSPeter Brune } 158d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);} 159eaedb033SPeter Brune PetscFunctionReturn(0); 160eaedb033SPeter Brune } 161eaedb033SPeter Brune 162eaedb033SPeter Brune #undef __FUNCT__ 163eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 164eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 165eaedb033SPeter Brune { 166eaedb033SPeter Brune PetscErrorCode ierr; 167111ade9eSPeter Brune PCASMType asmtype; 168b20c023fSPeter Brune PetscBool flg,monflg; 169111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1706e111a19SKarl Rupp 171eaedb033SPeter Brune PetscFunctionBegin; 172111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 173111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 1741aa26658SKarl Rupp if (flg) nasm->type = asmtype; 175b20c023fSPeter Brune flg = PETSC_FALSE; 176b20c023fSPeter Brune monflg = PETSC_TRUE; 177b20c023fSPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 178d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 179b20c023fSPeter Brune if (flg) { 180b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 181b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 182b20c023fSPeter Brune } 183eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 184eaedb033SPeter Brune PetscFunctionReturn(0); 185eaedb033SPeter Brune } 186eaedb033SPeter Brune 187eaedb033SPeter Brune #undef __FUNCT__ 188eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 189eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 190eaedb033SPeter Brune { 191b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 192b20c023fSPeter Brune PetscErrorCode ierr; 193b20c023fSPeter Brune PetscMPIInt rank; 194b20c023fSPeter Brune PetscInt i,N; 195b20c023fSPeter Brune PetscBool iascii,isstring; 196b20c023fSPeter Brune PetscViewer sviewer; 197ce94432eSBarry Smith MPI_Comm comm; 198b20c023fSPeter Brune 199b20c023fSPeter Brune PetscFunctionBegin; 200ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 201b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 202b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 203b20c023fSPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 204b20c023fSPeter Brune ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr); 205b20c023fSPeter Brune if (iascii) { 206b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr); 207b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Nonlinear Additive Schwarz: restriction/interpolation type - %s\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 208b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Nonlinear Additive Schwarz: subSNES iterations: %D subKSP iterations: %D\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 209b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 210b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 211b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 212b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local SNES objects:\n");CHKERRQ(ierr); 213b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 214b20c023fSPeter Brune if (!rank) { 215b20c023fSPeter Brune for (i=0; i<nasm->n; i++) { 216b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 217b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 218b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 219b20c023fSPeter Brune } 220b20c023fSPeter Brune } 221b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 222b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 223b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 224b20c023fSPeter Brune } else if (isstring) { 225b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 226b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 227b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 228b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 229b20c023fSPeter Brune } 230eaedb033SPeter Brune PetscFunctionReturn(0); 231eaedb033SPeter Brune } 232eaedb033SPeter Brune 233eaedb033SPeter Brune #undef __FUNCT__ 234eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 23576857b2aSPeter Brune /*@ 23676857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 23776857b2aSPeter Brune 23876857b2aSPeter Brune Not Collective 23976857b2aSPeter Brune 24076857b2aSPeter Brune Input Parameters: 24176857b2aSPeter Brune + SNES - the SNES context 24276857b2aSPeter Brune . n - the number of local subdomains 24376857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 24476857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 24576857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 24676857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 24776857b2aSPeter Brune 24876857b2aSPeter Brune Level: intermediate 24976857b2aSPeter Brune 25076857b2aSPeter Brune .keywords: SNES, NASM 25176857b2aSPeter Brune 25276857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 25376857b2aSPeter Brune @*/ 254a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 255a6dfd86eSKarl Rupp { 256eaedb033SPeter Brune PetscErrorCode ierr; 257111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 2586e111a19SKarl Rupp 259eaedb033SPeter Brune PetscFunctionBegin; 260eaedb033SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 261111ade9eSPeter Brune ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 262eaedb033SPeter Brune PetscFunctionReturn(0); 263eaedb033SPeter Brune } 264eaedb033SPeter Brune 265eaedb033SPeter Brune EXTERN_C_BEGIN 266eaedb033SPeter Brune #undef __FUNCT__ 267eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 268a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 269a6dfd86eSKarl Rupp { 270eaedb033SPeter Brune PetscInt i; 271eaedb033SPeter Brune PetscErrorCode ierr; 272eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 2736e111a19SKarl Rupp 274eaedb033SPeter Brune PetscFunctionBegin; 275ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 276eaedb033SPeter Brune 277111ade9eSPeter Brune /* tear down the previously set things */ 278111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 279111ade9eSPeter Brune 280eaedb033SPeter Brune nasm->n = n; 281111ade9eSPeter Brune if (oscatter) { 282111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 283eaedb033SPeter Brune } 284111ade9eSPeter Brune if (iscatter) { 285111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 286eaedb033SPeter Brune } 287111ade9eSPeter Brune if (gscatter) { 288111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 289111ade9eSPeter Brune } 290111ade9eSPeter Brune if (oscatter) { 291111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 292eaedb033SPeter Brune for (i=0; i<n; i++) { 293111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 294eaedb033SPeter Brune } 295111ade9eSPeter Brune } 296111ade9eSPeter Brune if (iscatter) { 297111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 298eaedb033SPeter Brune for (i=0; i<n; i++) { 299111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 300eaedb033SPeter Brune } 301eaedb033SPeter Brune } 302111ade9eSPeter Brune if (gscatter) { 303111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 304eaedb033SPeter Brune for (i=0; i<n; i++) { 305111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 306eaedb033SPeter Brune } 307eaedb033SPeter Brune } 308111ade9eSPeter Brune 309eaedb033SPeter Brune if (subsnes) { 310eaedb033SPeter Brune ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 311eaedb033SPeter Brune for (i=0; i<n; i++) { 312eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 313eaedb033SPeter Brune } 314eaedb033SPeter Brune } 315eaedb033SPeter Brune PetscFunctionReturn(0); 316eaedb033SPeter Brune } 317eaedb033SPeter Brune EXTERN_C_END 318eaedb033SPeter Brune 319eaedb033SPeter Brune #undef __FUNCT__ 32076857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains" 32176857b2aSPeter Brune /*@ 32276857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 32376857b2aSPeter Brune 32476857b2aSPeter Brune Not Collective 32576857b2aSPeter Brune 32676857b2aSPeter Brune Input Parameters: 32776857b2aSPeter Brune . SNES - the SNES context 32876857b2aSPeter Brune 32976857b2aSPeter Brune Output Parameters: 33076857b2aSPeter Brune + n - the number of local subdomains 33176857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 33276857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 33376857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 33476857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 33576857b2aSPeter Brune 33676857b2aSPeter Brune Level: intermediate 33776857b2aSPeter Brune 33876857b2aSPeter Brune .keywords: SNES, NASM 33976857b2aSPeter Brune 34076857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 34176857b2aSPeter Brune @*/ 34276857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 34376857b2aSPeter Brune { 34476857b2aSPeter Brune PetscErrorCode ierr; 34576857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 34676857b2aSPeter Brune 34776857b2aSPeter Brune PetscFunctionBegin; 34876857b2aSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 34976857b2aSPeter Brune ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 35076857b2aSPeter Brune PetscFunctionReturn(0); 35176857b2aSPeter Brune } 35276857b2aSPeter Brune 35376857b2aSPeter Brune EXTERN_C_BEGIN 35476857b2aSPeter Brune #undef __FUNCT__ 35576857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 35676857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 35776857b2aSPeter Brune { 35876857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 35976857b2aSPeter Brune 36076857b2aSPeter Brune PetscFunctionBegin; 36176857b2aSPeter Brune if (n) *n = nasm->n; 36276857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 36376857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 36476857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 36576857b2aSPeter Brune if (subsnes) *subsnes = nasm->subsnes; 36676857b2aSPeter Brune PetscFunctionReturn(0); 36776857b2aSPeter Brune } 36876857b2aSPeter Brune EXTERN_C_END 36976857b2aSPeter Brune 37076857b2aSPeter Brune #undef __FUNCT__ 37176857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs" 37276857b2aSPeter Brune /*@ 37376857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 37476857b2aSPeter Brune 37576857b2aSPeter Brune Not Collective 37676857b2aSPeter Brune 37776857b2aSPeter Brune Input Parameters: 37876857b2aSPeter Brune . SNES - the SNES context 37976857b2aSPeter Brune 38076857b2aSPeter Brune Output Parameters: 38176857b2aSPeter Brune + n - the number of local subdomains 38276857b2aSPeter Brune . x - The subdomain solution vector 38376857b2aSPeter Brune . y - The subdomain step vector 38476857b2aSPeter Brune . b - The subdomain RHS vector 38576857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 38676857b2aSPeter Brune 38776857b2aSPeter Brune Level: developer 38876857b2aSPeter Brune 38976857b2aSPeter Brune .keywords: SNES, NASM 39076857b2aSPeter Brune 39176857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 39276857b2aSPeter Brune @*/ 39376857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 39476857b2aSPeter Brune { 39576857b2aSPeter Brune PetscErrorCode ierr; 39676857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 39776857b2aSPeter Brune 39876857b2aSPeter Brune PetscFunctionBegin; 39976857b2aSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",(void (**)(void))&f);CHKERRQ(ierr); 40076857b2aSPeter Brune ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr); 40176857b2aSPeter Brune PetscFunctionReturn(0); 40276857b2aSPeter Brune } 40376857b2aSPeter Brune 40476857b2aSPeter Brune EXTERN_C_BEGIN 40576857b2aSPeter Brune #undef __FUNCT__ 40676857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 40776857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 40876857b2aSPeter Brune { 40976857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 41076857b2aSPeter Brune 41176857b2aSPeter Brune PetscFunctionBegin; 41276857b2aSPeter Brune if (n) *n = nasm->n; 41376857b2aSPeter Brune if (x) *x = nasm->x; 41476857b2aSPeter Brune if (y) *y = nasm->y; 41576857b2aSPeter Brune if (b) *b = nasm->b; 41676857b2aSPeter Brune if (xl) *xl = nasm->xl; 41776857b2aSPeter Brune PetscFunctionReturn(0); 41876857b2aSPeter Brune } 41976857b2aSPeter Brune EXTERN_C_END 42076857b2aSPeter Brune 421d728fb7dSPeter Brune #undef __FUNCT__ 422d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 423d728fb7dSPeter Brune /*@ 424d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 425d728fb7dSPeter Brune 426d728fb7dSPeter Brune Collective on SNES 427d728fb7dSPeter Brune 428d728fb7dSPeter Brune Input Parameters: 429d728fb7dSPeter Brune + SNES - the SNES context 430d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 431d728fb7dSPeter Brune 432d728fb7dSPeter Brune Level: developer 433d728fb7dSPeter Brune 434d728fb7dSPeter Brune Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 435d728fb7dSPeter Brune is needed at each linear iteration. 436d728fb7dSPeter Brune 437d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN 438d728fb7dSPeter Brune 439d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 440d728fb7dSPeter Brune @*/ 441d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 442d728fb7dSPeter Brune { 443d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 444d728fb7dSPeter Brune PetscErrorCode ierr; 445d728fb7dSPeter Brune 446d728fb7dSPeter Brune PetscFunctionBegin; 447d728fb7dSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",(void (**)(void))&f);CHKERRQ(ierr); 448d728fb7dSPeter Brune ierr = (f)(snes,flg);CHKERRQ(ierr); 449d728fb7dSPeter Brune PetscFunctionReturn(0); 450d728fb7dSPeter Brune } 451d728fb7dSPeter Brune 452d728fb7dSPeter Brune EXTERN_C_BEGIN 453d728fb7dSPeter Brune #undef __FUNCT__ 454d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 455d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 456d728fb7dSPeter Brune { 457d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 458d728fb7dSPeter Brune 459d728fb7dSPeter Brune PetscFunctionBegin; 460d728fb7dSPeter Brune nasm->finaljacobian = flg; 461d728fb7dSPeter Brune PetscFunctionReturn(0); 462d728fb7dSPeter Brune } 463d728fb7dSPeter Brune EXTERN_C_END 464d728fb7dSPeter Brune 46576857b2aSPeter Brune 46676857b2aSPeter Brune #undef __FUNCT__ 467eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 4680adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 4690adebc6cSBarry Smith { 470eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 471258e1594SPeter Brune SNES subsnes; 472eaedb033SPeter Brune PetscInt i; 473eaedb033SPeter Brune PetscErrorCode ierr; 474111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 475111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 476111ade9eSPeter Brune DM dm,subdm; 4770adebc6cSBarry Smith 478eaedb033SPeter Brune PetscFunctionBegin; 479eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 480111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 481b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 482eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 48370c78f05SPeter Brune /* scatter the solution to the local solution */ 48470c78f05SPeter Brune Xlloc = nasm->xl[i]; 48570c78f05SPeter Brune gscat = nasm->gscatter[i]; 48670c78f05SPeter Brune oscat = nasm->oscatter[i]; 48770c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 48870c78f05SPeter Brune if (B) { 48970c78f05SPeter Brune /* scatter the RHS to the local RHS */ 49070c78f05SPeter Brune Bl = nasm->b[i]; 49170c78f05SPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 49270c78f05SPeter Brune } 49370c78f05SPeter Brune } 494b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 495b20c023fSPeter Brune 496b20c023fSPeter Brune 497b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 49870c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 49970c78f05SPeter Brune Xl = nasm->x[i]; 50070c78f05SPeter Brune Xlloc = nasm->xl[i]; 50170c78f05SPeter Brune Yl = nasm->y[i]; 502258e1594SPeter Brune subsnes = nasm->subsnes[i]; 503258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 504111ade9eSPeter Brune iscat = nasm->iscatter[i]; 505111ade9eSPeter Brune oscat = nasm->oscatter[i]; 506111ade9eSPeter Brune gscat = nasm->gscatter[i]; 507ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 50824b7f281SPeter Brune if (B) { 50924b7f281SPeter Brune Bl = nasm->b[i]; 510ed3c11a9SPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 511ed3c11a9SPeter Brune } else Bl = NULL; 512ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 51370c78f05SPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 51470c78f05SPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 515111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 516ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 517ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 518111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 519111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 520111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 521111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 522ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 523eaedb033SPeter Brune } 524ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 525ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 52670c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 52770c78f05SPeter Brune Yl = nasm->y[i]; 52870c78f05SPeter Brune iscat = nasm->iscatter[i]; 52970c78f05SPeter Brune oscat = nasm->oscatter[i]; 53070c78f05SPeter Brune if (nasm->type == PC_ASM_BASIC) { 53170c78f05SPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 53270c78f05SPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 53370c78f05SPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 534ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 53570c78f05SPeter Brune } 536b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 537111ade9eSPeter Brune ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr); 538eaedb033SPeter Brune PetscFunctionReturn(0); 539eaedb033SPeter Brune } 540eaedb033SPeter Brune 541eaedb033SPeter Brune #undef __FUNCT__ 542d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 543d728fb7dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec X) 544d728fb7dSPeter Brune { 545d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 546d728fb7dSPeter Brune SNES subsnes; 547ca44f815SPeter Brune PetscInt i,lag = 1; 548d728fb7dSPeter Brune PetscErrorCode ierr; 549*e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 550d728fb7dSPeter Brune VecScatter oscat,gscat; 551d728fb7dSPeter Brune DM dm,subdm; 552d728fb7dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 553d728fb7dSPeter Brune PetscFunctionBegin; 554*e59f0a30SPeter Brune F = snes->vec_func; 555*e59f0a30SPeter Brune if (snes->normtype == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 556d728fb7dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 557d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 558d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 559d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 560d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 561*e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 562d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 563d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 564d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 565*e59f0a30SPeter Brune ierr = VecScatterBegin(oscat,F,Fl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 566d728fb7dSPeter Brune } 567d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 568d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 569*e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 570d728fb7dSPeter Brune Xl = nasm->x[i]; 571d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 572d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 573d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 574d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 575ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 576*e59f0a30SPeter Brune ierr = VecScatterEnd(oscat,F,Fl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 577ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 578d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 579d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 580d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 581ca44f815SPeter Brune 582ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 583ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 5842f543b25SPeter Brune ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr); 585ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 586d728fb7dSPeter Brune ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr); 587d728fb7dSPeter Brune } 588d728fb7dSPeter Brune PetscFunctionReturn(0); 589d728fb7dSPeter Brune } 590d728fb7dSPeter Brune 591d728fb7dSPeter Brune #undef __FUNCT__ 592eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 593eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 594eaedb033SPeter Brune { 595eaedb033SPeter Brune Vec F; 596eaedb033SPeter Brune Vec X; 597eaedb033SPeter Brune Vec B; 598111ade9eSPeter Brune Vec Y; 599eaedb033SPeter Brune PetscInt i; 600ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 601eaedb033SPeter Brune PetscErrorCode ierr; 602eaedb033SPeter Brune SNESNormType normtype; 603d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 604eaedb033SPeter Brune 605eaedb033SPeter Brune PetscFunctionBegin; 606eaedb033SPeter Brune X = snes->vec_sol; 607111ade9eSPeter Brune Y = snes->vec_sol_update; 608eaedb033SPeter Brune F = snes->vec_func; 609eaedb033SPeter Brune B = snes->vec_rhs; 610eaedb033SPeter Brune 611eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 612eaedb033SPeter Brune snes->iter = 0; 613eaedb033SPeter Brune snes->norm = 0.; 614eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 615eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 616eaedb033SPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 617eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 618eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 619eaedb033SPeter Brune if (!snes->vec_func_init_set) { 620eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 621eaedb033SPeter Brune if (snes->domainerror) { 622eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 623eaedb033SPeter Brune PetscFunctionReturn(0); 624eaedb033SPeter Brune } 6251aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 626eaedb033SPeter Brune 627eaedb033SPeter Brune /* convergence test */ 628eaedb033SPeter Brune if (!snes->norm_init_set) { 629eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 630eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 631eaedb033SPeter Brune } else { 632eaedb033SPeter Brune fnorm = snes->norm_init; 633eaedb033SPeter Brune snes->norm_init_set = PETSC_FALSE; 634eaedb033SPeter Brune } 635eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 636eaedb033SPeter Brune snes->iter = 0; 637eaedb033SPeter Brune snes->norm = fnorm; 638eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 639a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 640eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 641eaedb033SPeter Brune 642eaedb033SPeter Brune /* set parameter for default relative tolerance convergence test */ 643eaedb033SPeter Brune snes->ttol = fnorm*snes->rtol; 644eaedb033SPeter Brune 645eaedb033SPeter Brune /* test convergence */ 646eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 647eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 648eaedb033SPeter Brune } else { 649eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 650a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 651eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 652eaedb033SPeter Brune } 653eaedb033SPeter Brune 654eaedb033SPeter Brune /* Call general purpose update function */ 655eaedb033SPeter Brune if (snes->ops->update) { 656eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 657eaedb033SPeter Brune } 658eaedb033SPeter Brune 659eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 660111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 661eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 662eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 663eaedb033SPeter Brune if (snes->domainerror) { 664eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 665d728fb7dSPeter Brune break; 666eaedb033SPeter Brune } 667eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 668eaedb033SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 669eaedb033SPeter Brune } 670eaedb033SPeter Brune /* Monitor convergence */ 671eaedb033SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 672eaedb033SPeter Brune snes->iter = i+1; 673eaedb033SPeter Brune snes->norm = fnorm; 674eaedb033SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 675a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 676eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 677eaedb033SPeter Brune /* Test for convergence */ 678*e59f0a30SPeter Brune if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 679d728fb7dSPeter Brune if (snes->reason) break; 680eaedb033SPeter Brune /* Call general purpose update function */ 681*e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 682eaedb033SPeter Brune } 683d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 684eaedb033SPeter Brune if (normtype == SNES_NORM_FUNCTION) { 685eaedb033SPeter Brune if (i == snes->max_its) { 686eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 687eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 688eaedb033SPeter Brune } 6891aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 690eaedb033SPeter Brune PetscFunctionReturn(0); 691eaedb033SPeter Brune } 692eaedb033SPeter Brune 693eaedb033SPeter Brune /*MC 694eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 695eaedb033SPeter Brune 69670c78603SPeter Brune Options Database: 69770c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 69870c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 69970c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 70070c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 70170c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 70270c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 70370c78603SPeter Brune 704eaedb033SPeter Brune Level: advanced 705eaedb033SPeter Brune 706eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 707eaedb033SPeter Brune M*/ 708eaedb033SPeter Brune 709eaedb033SPeter Brune EXTERN_C_BEGIN 710eaedb033SPeter Brune #undef __FUNCT__ 711eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 712eaedb033SPeter Brune PetscErrorCode SNESCreate_NASM(SNES snes) 713eaedb033SPeter Brune { 714eaedb033SPeter Brune SNES_NASM *nasm; 715eaedb033SPeter Brune PetscErrorCode ierr; 716eaedb033SPeter Brune 717eaedb033SPeter Brune PetscFunctionBegin; 718eaedb033SPeter Brune ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 719eaedb033SPeter Brune snes->data = (void*)nasm; 720eaedb033SPeter Brune 721eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 722eaedb033SPeter Brune nasm->subsnes = 0; 723eaedb033SPeter Brune nasm->x = 0; 724111ade9eSPeter Brune nasm->xl = 0; 725111ade9eSPeter Brune nasm->y = 0; 726eaedb033SPeter Brune nasm->b = 0; 727111ade9eSPeter Brune nasm->oscatter = 0; 728111ade9eSPeter Brune nasm->iscatter = 0; 729111ade9eSPeter Brune nasm->gscatter = 0; 730111ade9eSPeter Brune 731111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 732d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 733eaedb033SPeter Brune 734eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 735eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 736eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 737eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 738eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 739eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 740eaedb033SPeter Brune 741eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 742eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 743eaedb033SPeter Brune 7440298fd71SBarry Smith nasm->eventrestrictinterp = 0; 7450298fd71SBarry Smith nasm->eventsubsolve = 0; 746b20c023fSPeter Brune 747eaedb033SPeter Brune if (!snes->tolerancesset) { 748eaedb033SPeter Brune snes->max_its = 10000; 749eaedb033SPeter Brune snes->max_funcs = 10000; 750eaedb033SPeter Brune } 751eaedb033SPeter Brune 752eaedb033SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM", 753eaedb033SPeter Brune SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 75476857b2aSPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMGetSubdomains_C","SNESNASMGetSubdomains_NASM", 75576857b2aSPeter Brune SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 75676857b2aSPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMGetSubdomainVecs_C","SNESNASMGetSubdomainVecs_NASM", 75776857b2aSPeter Brune SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 758d728fb7dSPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C","SNESNASMSetComputeFinalJacobian_NASM", 759d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 760eaedb033SPeter Brune PetscFunctionReturn(0); 761eaedb033SPeter Brune } 762eaedb033SPeter Brune EXTERN_C_END 763