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 */ 18a4f17876SPeter Brune PetscBool same_local_solves; /* flag to determine if the solvers have been individually modified */ 19b20c023fSPeter Brune 20b20c023fSPeter Brune /* logging events */ 21b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 22b20c023fSPeter Brune PetscLogEvent eventsubsolve; 23602bec5dSPeter Brune 24602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 25602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 26eaedb033SPeter Brune } SNES_NASM; 27eaedb033SPeter Brune 28b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 29602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 30b20c023fSPeter Brune 31eaedb033SPeter Brune #undef __FUNCT__ 32eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM" 33eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes) 34eaedb033SPeter Brune { 35eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 36eaedb033SPeter Brune PetscErrorCode ierr; 37eaedb033SPeter Brune PetscInt i; 386e111a19SKarl Rupp 39eaedb033SPeter Brune PetscFunctionBegin; 40eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 41111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 42f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 43111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 44bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 45eaedb033SPeter Brune 46bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 47111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 48111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 49111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 50eaedb033SPeter Brune } 51111ade9eSPeter Brune 52111ade9eSPeter Brune if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 53111ade9eSPeter Brune if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 54111ade9eSPeter Brune if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 55111ade9eSPeter Brune if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 56111ade9eSPeter Brune 57602bec5dSPeter Brune if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 58602bec5dSPeter Brune 59111ade9eSPeter Brune if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 60111ade9eSPeter Brune if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 61111ade9eSPeter Brune if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 62111ade9eSPeter Brune if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 63b20c023fSPeter Brune 64b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 65b20c023fSPeter Brune nasm->eventsubsolve = 0; 66eaedb033SPeter Brune PetscFunctionReturn(0); 67eaedb033SPeter Brune } 68eaedb033SPeter Brune 69eaedb033SPeter Brune #undef __FUNCT__ 70eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM" 71eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes) 72eaedb033SPeter Brune { 73eaedb033SPeter Brune PetscErrorCode ierr; 746e111a19SKarl Rupp 75eaedb033SPeter Brune PetscFunctionBegin; 76eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 7722d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 78eaedb033SPeter Brune PetscFunctionReturn(0); 79eaedb033SPeter Brune } 80eaedb033SPeter Brune 81eaedb033SPeter Brune #undef __FUNCT__ 82111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 830adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 840adebc6cSBarry Smith { 85111ade9eSPeter Brune PetscErrorCode ierr; 86111ade9eSPeter Brune Vec bcs = (Vec)ctx; 876e111a19SKarl Rupp 88111ade9eSPeter Brune PetscFunctionBegin; 89111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 90111ade9eSPeter Brune PetscFunctionReturn(0); 91111ade9eSPeter Brune } 92111ade9eSPeter Brune 93111ade9eSPeter Brune #undef __FUNCT__ 94eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM" 95eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes) 96eaedb033SPeter Brune { 97eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 98eaedb033SPeter Brune PetscErrorCode ierr; 9976857b2aSPeter Brune DM dm,subdm; 100111ade9eSPeter Brune DM *subdms; 101111ade9eSPeter Brune PetscInt i; 102eaedb033SPeter Brune const char *optionsprefix; 103111ade9eSPeter Brune Vec F; 104ed3c11a9SPeter Brune PetscMPIInt size; 105ed3c11a9SPeter Brune KSP ksp; 106ed3c11a9SPeter Brune PC pc; 107eaedb033SPeter Brune 108eaedb033SPeter Brune PetscFunctionBegin; 109eaedb033SPeter Brune if (!nasm->subsnes) { 110eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1110a696f66SPeter Brune if (dm) { 112eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1130298fd71SBarry Smith ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 114ce94432eSBarry Smith if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 115111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 116eaedb033SPeter Brune 117eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 118785e854fSJed Brown ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr); 119111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 120cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 121cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 122cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 123cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 124ed3c11a9SPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 125ed3c11a9SPeter Brune if (size == 1) { 126ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 127ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 128ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 129ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 130ed3c11a9SPeter Brune } 131cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 132111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 133111ade9eSPeter Brune } 134111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 135ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 136ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 137111ade9eSPeter Brune /* allocate the global vectors */ 138e245e0aaSPeter Brune if (!nasm->x) { 1391795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr); 140e245e0aaSPeter Brune } 141e245e0aaSPeter Brune if (!nasm->xl) { 1421795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr); 143e245e0aaSPeter Brune } 144e245e0aaSPeter Brune if (!nasm->y) { 1451795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr); 146e245e0aaSPeter Brune } 147e245e0aaSPeter Brune if (!nasm->b) { 1481795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr); 149e245e0aaSPeter Brune } 150111ade9eSPeter Brune 151111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1520298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 15376857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 15476857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 15576857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 15676857b2aSPeter Brune if (!nasm->xl[i]) { 157111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 158111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 15976857b2aSPeter Brune } 1600298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 161111ade9eSPeter Brune } 162602bec5dSPeter Brune if (nasm->finaljacobian) { 163602bec5dSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 164602bec5dSPeter Brune if (nasm->fjtype == 2) { 165602bec5dSPeter Brune ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 166602bec5dSPeter Brune } 167602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 168602bec5dSPeter Brune ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 169602bec5dSPeter Brune } 170602bec5dSPeter Brune } 171eaedb033SPeter Brune PetscFunctionReturn(0); 172eaedb033SPeter Brune } 173eaedb033SPeter Brune 174eaedb033SPeter Brune #undef __FUNCT__ 175eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 176eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 177eaedb033SPeter Brune { 178eaedb033SPeter Brune PetscErrorCode ierr; 179111ade9eSPeter Brune PCASMType asmtype; 180a4f17876SPeter Brune PetscBool flg,monflg,subviewflg; 181111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1826e111a19SKarl Rupp 183eaedb033SPeter Brune PetscFunctionBegin; 184111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 185111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 1861aa26658SKarl Rupp if (flg) nasm->type = asmtype; 187b20c023fSPeter Brune flg = PETSC_FALSE; 188b20c023fSPeter Brune monflg = PETSC_TRUE; 189610116beSPeter Brune ierr = PetscOptionsReal("-snes_nasm_damping","Log times for subSNES solves and restriction","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr); 190610116beSPeter Brune if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 191a4f17876SPeter Brune subviewflg = PETSC_FALSE; 192a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr); 193a4f17876SPeter Brune if (flg) { 194a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 195a4f17876SPeter Brune if (!subviewflg) { 196a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 197a4f17876SPeter Brune } 198a4f17876SPeter Brune } 199d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 200602bec5dSPeter Brune ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 201a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 202b20c023fSPeter Brune if (flg) { 203b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 204b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 205b20c023fSPeter Brune } 206eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 207eaedb033SPeter Brune PetscFunctionReturn(0); 208eaedb033SPeter Brune } 209eaedb033SPeter Brune 210eaedb033SPeter Brune #undef __FUNCT__ 211eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 212eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 213eaedb033SPeter Brune { 214b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 215b20c023fSPeter Brune PetscErrorCode ierr; 216a4f17876SPeter Brune PetscMPIInt rank,size; 217dd2fa690SBarry Smith PetscInt i,N,bsz; 218b20c023fSPeter Brune PetscBool iascii,isstring; 219b20c023fSPeter Brune PetscViewer sviewer; 220ce94432eSBarry Smith MPI_Comm comm; 221b20c023fSPeter Brune 222b20c023fSPeter Brune PetscFunctionBegin; 223ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 224b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 225b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 226b20c023fSPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 227a4f17876SPeter Brune ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 22879f8d099SPeter Brune ierr = MPI_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 229b20c023fSPeter Brune if (iascii) { 230b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr); 231a4f17876SPeter Brune if (nasm->same_local_solves) { 232a4f17876SPeter Brune if (nasm->subsnes) { 233a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve is the same for all blocks:\n");CHKERRQ(ierr); 234a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 235a4f17876SPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 236a4f17876SPeter Brune if (!rank) { 237a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 238a4f17876SPeter Brune ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr); 239a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 240a4f17876SPeter Brune } 241a4f17876SPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 242a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 243a4f17876SPeter Brune } 244a4f17876SPeter Brune } else { 245a4f17876SPeter Brune /* print the solver on each block */ 246b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 247b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 248b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 249a4f17876SPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 250a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr); 251b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 252a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 253b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 254a4f17876SPeter Brune for (i=0; i<nasm->n; i++) { 255a4f17876SPeter Brune ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr); 256a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 257b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 258a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 259b20c023fSPeter Brune } 260a4f17876SPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 261a4f17876SPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 262b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 263a4f17876SPeter Brune } 264b20c023fSPeter Brune } else if (isstring) { 265b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 266b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 267b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 268b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 269b20c023fSPeter Brune } 270eaedb033SPeter Brune PetscFunctionReturn(0); 271eaedb033SPeter Brune } 272eaedb033SPeter Brune 273eaedb033SPeter Brune #undef __FUNCT__ 274eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 27576857b2aSPeter Brune /*@ 27676857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 27776857b2aSPeter Brune 27876857b2aSPeter Brune Not Collective 27976857b2aSPeter Brune 28076857b2aSPeter Brune Input Parameters: 28176857b2aSPeter Brune + SNES - the SNES context 28276857b2aSPeter Brune . n - the number of local subdomains 28376857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 28476857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 28576857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 28676857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 28776857b2aSPeter Brune 28876857b2aSPeter Brune Level: intermediate 28976857b2aSPeter Brune 29076857b2aSPeter Brune .keywords: SNES, NASM 29176857b2aSPeter Brune 29276857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 29376857b2aSPeter Brune @*/ 294a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 295a6dfd86eSKarl Rupp { 296eaedb033SPeter Brune PetscErrorCode ierr; 297111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 2986e111a19SKarl Rupp 299eaedb033SPeter Brune PetscFunctionBegin; 3000005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 3014b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 302eaedb033SPeter Brune PetscFunctionReturn(0); 303eaedb033SPeter Brune } 304eaedb033SPeter Brune 305eaedb033SPeter Brune #undef __FUNCT__ 306eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 307a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 308a6dfd86eSKarl Rupp { 309eaedb033SPeter Brune PetscInt i; 310eaedb033SPeter Brune PetscErrorCode ierr; 311eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3126e111a19SKarl Rupp 313eaedb033SPeter Brune PetscFunctionBegin; 314ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 315eaedb033SPeter Brune 316111ade9eSPeter Brune /* tear down the previously set things */ 317111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 318111ade9eSPeter Brune 319eaedb033SPeter Brune nasm->n = n; 320111ade9eSPeter Brune if (oscatter) { 321111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 322eaedb033SPeter Brune } 323111ade9eSPeter Brune if (iscatter) { 324111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 325eaedb033SPeter Brune } 326111ade9eSPeter Brune if (gscatter) { 327111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 328111ade9eSPeter Brune } 329111ade9eSPeter Brune if (oscatter) { 330785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 331eaedb033SPeter Brune for (i=0; i<n; i++) { 332111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 333eaedb033SPeter Brune } 334111ade9eSPeter Brune } 335111ade9eSPeter Brune if (iscatter) { 336785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 337eaedb033SPeter Brune for (i=0; i<n; i++) { 338111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 339eaedb033SPeter Brune } 340eaedb033SPeter Brune } 341111ade9eSPeter Brune if (gscatter) { 342785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 343eaedb033SPeter Brune for (i=0; i<n; i++) { 344111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 345eaedb033SPeter Brune } 346eaedb033SPeter Brune } 347111ade9eSPeter Brune 348eaedb033SPeter Brune if (subsnes) { 349785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 350eaedb033SPeter Brune for (i=0; i<n; i++) { 351eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 352eaedb033SPeter Brune } 353a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 354eaedb033SPeter Brune } 355eaedb033SPeter Brune PetscFunctionReturn(0); 356eaedb033SPeter Brune } 357eaedb033SPeter Brune 358eaedb033SPeter Brune #undef __FUNCT__ 35976857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains" 36076857b2aSPeter Brune /*@ 36176857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 36276857b2aSPeter Brune 36376857b2aSPeter Brune Not Collective 36476857b2aSPeter Brune 36576857b2aSPeter Brune Input Parameters: 36676857b2aSPeter Brune . SNES - the SNES context 36776857b2aSPeter Brune 36876857b2aSPeter Brune Output Parameters: 36976857b2aSPeter Brune + n - the number of local subdomains 37076857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 37176857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 37276857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 37376857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 37476857b2aSPeter Brune 37576857b2aSPeter Brune Level: intermediate 37676857b2aSPeter Brune 37776857b2aSPeter Brune .keywords: SNES, NASM 37876857b2aSPeter Brune 37976857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 38076857b2aSPeter Brune @*/ 38176857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 38276857b2aSPeter Brune { 38376857b2aSPeter Brune PetscErrorCode ierr; 38476857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 38576857b2aSPeter Brune 38676857b2aSPeter Brune PetscFunctionBegin; 3870005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 3884b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 38976857b2aSPeter Brune PetscFunctionReturn(0); 39076857b2aSPeter Brune } 39176857b2aSPeter Brune 39276857b2aSPeter Brune #undef __FUNCT__ 39376857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 39476857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 39576857b2aSPeter Brune { 39676857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 39776857b2aSPeter Brune 39876857b2aSPeter Brune PetscFunctionBegin; 39976857b2aSPeter Brune if (n) *n = nasm->n; 40076857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 40176857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 40276857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 403a4f17876SPeter Brune if (subsnes) { 404a4f17876SPeter Brune *subsnes = nasm->subsnes; 405a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 406a4f17876SPeter Brune } 40776857b2aSPeter Brune PetscFunctionReturn(0); 40876857b2aSPeter Brune } 40976857b2aSPeter Brune 41076857b2aSPeter Brune #undef __FUNCT__ 41176857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs" 41276857b2aSPeter Brune /*@ 41376857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 41476857b2aSPeter Brune 41576857b2aSPeter Brune Not Collective 41676857b2aSPeter Brune 41776857b2aSPeter Brune Input Parameters: 41876857b2aSPeter Brune . SNES - the SNES context 41976857b2aSPeter Brune 42076857b2aSPeter Brune Output Parameters: 42176857b2aSPeter Brune + n - the number of local subdomains 42276857b2aSPeter Brune . x - The subdomain solution vector 42376857b2aSPeter Brune . y - The subdomain step vector 42476857b2aSPeter Brune . b - The subdomain RHS vector 42576857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 42676857b2aSPeter Brune 42776857b2aSPeter Brune Level: developer 42876857b2aSPeter Brune 42976857b2aSPeter Brune .keywords: SNES, NASM 43076857b2aSPeter Brune 43176857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 43276857b2aSPeter Brune @*/ 43376857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 43476857b2aSPeter Brune { 43576857b2aSPeter Brune PetscErrorCode ierr; 43676857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 43776857b2aSPeter Brune 43876857b2aSPeter Brune PetscFunctionBegin; 4390005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 4404b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 44176857b2aSPeter Brune PetscFunctionReturn(0); 44276857b2aSPeter Brune } 44376857b2aSPeter Brune 44476857b2aSPeter Brune #undef __FUNCT__ 44576857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 44676857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 44776857b2aSPeter Brune { 44876857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 44976857b2aSPeter Brune 45076857b2aSPeter Brune PetscFunctionBegin; 45176857b2aSPeter Brune if (n) *n = nasm->n; 45276857b2aSPeter Brune if (x) *x = nasm->x; 45376857b2aSPeter Brune if (y) *y = nasm->y; 45476857b2aSPeter Brune if (b) *b = nasm->b; 45576857b2aSPeter Brune if (xl) *xl = nasm->xl; 45676857b2aSPeter Brune PetscFunctionReturn(0); 45776857b2aSPeter Brune } 45876857b2aSPeter Brune 459d728fb7dSPeter Brune #undef __FUNCT__ 460d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 461d728fb7dSPeter Brune /*@ 462d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 463d728fb7dSPeter Brune 464d728fb7dSPeter Brune Collective on SNES 465d728fb7dSPeter Brune 466d728fb7dSPeter Brune Input Parameters: 467d728fb7dSPeter Brune + SNES - the SNES context 468d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 469d728fb7dSPeter Brune 470d728fb7dSPeter Brune Level: developer 471d728fb7dSPeter Brune 472d728fb7dSPeter Brune Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 473d728fb7dSPeter Brune is needed at each linear iteration. 474d728fb7dSPeter Brune 475d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN 476d728fb7dSPeter Brune 477d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 478d728fb7dSPeter Brune @*/ 479d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 480d728fb7dSPeter Brune { 481d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 482d728fb7dSPeter Brune PetscErrorCode ierr; 483d728fb7dSPeter Brune 484d728fb7dSPeter Brune PetscFunctionBegin; 4850005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 4864b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 487d728fb7dSPeter Brune PetscFunctionReturn(0); 488d728fb7dSPeter Brune } 489d728fb7dSPeter Brune 490d728fb7dSPeter Brune #undef __FUNCT__ 491d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 492d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 493d728fb7dSPeter Brune { 494d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 495d728fb7dSPeter Brune 496d728fb7dSPeter Brune PetscFunctionBegin; 497d728fb7dSPeter Brune nasm->finaljacobian = flg; 498602bec5dSPeter Brune if (flg) snes->usesksp = PETSC_TRUE; 499d728fb7dSPeter Brune PetscFunctionReturn(0); 500d728fb7dSPeter Brune } 50176857b2aSPeter Brune 50276857b2aSPeter Brune #undef __FUNCT__ 503610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping" 504610116beSPeter Brune /*@ 505610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 506610116beSPeter Brune 507610116beSPeter Brune Logically collective on SNES 508610116beSPeter Brune 509610116beSPeter Brune Input Parameters: 510610116beSPeter Brune + SNES - the SNES context 511610116beSPeter Brune - dmp - damping 512610116beSPeter Brune 513610116beSPeter Brune Level: intermediate 514610116beSPeter Brune 515610116beSPeter Brune .keywords: SNES, NASM, damping 516610116beSPeter Brune 517610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 518610116beSPeter Brune @*/ 519610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 520610116beSPeter Brune { 521610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 522610116beSPeter Brune PetscErrorCode ierr; 523610116beSPeter Brune 524610116beSPeter Brune PetscFunctionBegin; 525610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 526e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 527610116beSPeter Brune PetscFunctionReturn(0); 528610116beSPeter Brune } 529610116beSPeter Brune 530610116beSPeter Brune #undef __FUNCT__ 531610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping_NASM" 532610116beSPeter Brune PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 533610116beSPeter Brune { 534610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 535610116beSPeter Brune 536610116beSPeter Brune PetscFunctionBegin; 537610116beSPeter Brune nasm->damping = dmp; 538610116beSPeter Brune PetscFunctionReturn(0); 539610116beSPeter Brune } 540610116beSPeter Brune 541610116beSPeter Brune #undef __FUNCT__ 542610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping" 543610116beSPeter Brune /*@ 544610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 545610116beSPeter Brune 546610116beSPeter Brune Not Collective 547610116beSPeter Brune 548610116beSPeter Brune Input Parameters: 549610116beSPeter Brune + SNES - the SNES context 550610116beSPeter Brune - dmp - damping 551610116beSPeter Brune 552610116beSPeter Brune Level: intermediate 553610116beSPeter Brune 554610116beSPeter Brune .keywords: SNES, NASM, damping 555610116beSPeter Brune 556610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping() 557610116beSPeter Brune @*/ 558610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 559610116beSPeter Brune { 560610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal*); 561610116beSPeter Brune PetscErrorCode ierr; 562610116beSPeter Brune 563610116beSPeter Brune PetscFunctionBegin; 564610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 565e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 566610116beSPeter Brune PetscFunctionReturn(0); 567610116beSPeter Brune } 568610116beSPeter Brune 569610116beSPeter Brune #undef __FUNCT__ 570610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping_NASM" 571610116beSPeter Brune PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 572610116beSPeter Brune { 573610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 574610116beSPeter Brune 575610116beSPeter Brune PetscFunctionBegin; 576610116beSPeter Brune *dmp = nasm->damping; 577610116beSPeter Brune PetscFunctionReturn(0); 578610116beSPeter Brune } 579610116beSPeter Brune 580610116beSPeter Brune 581610116beSPeter Brune #undef __FUNCT__ 582eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 5830adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 5840adebc6cSBarry Smith { 585eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 586258e1594SPeter Brune SNES subsnes; 587eaedb033SPeter Brune PetscInt i; 588610116beSPeter Brune PetscReal dmp; 589eaedb033SPeter Brune PetscErrorCode ierr; 590111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 591111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 592111ade9eSPeter Brune DM dm,subdm; 5930adebc6cSBarry Smith 594eaedb033SPeter Brune PetscFunctionBegin; 595eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 596610116beSPeter Brune ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 597111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 598b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 599eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 60070c78f05SPeter Brune /* scatter the solution to the local solution */ 60170c78f05SPeter Brune Xlloc = nasm->xl[i]; 60270c78f05SPeter Brune gscat = nasm->gscatter[i]; 60370c78f05SPeter Brune oscat = nasm->oscatter[i]; 60470c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 60570c78f05SPeter Brune if (B) { 60670c78f05SPeter Brune /* scatter the RHS to the local RHS */ 60770c78f05SPeter Brune Bl = nasm->b[i]; 60870c78f05SPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 60970c78f05SPeter Brune } 61070c78f05SPeter Brune } 611b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 612b20c023fSPeter Brune 613b20c023fSPeter Brune 614b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 61570c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 61670c78f05SPeter Brune Xl = nasm->x[i]; 61770c78f05SPeter Brune Xlloc = nasm->xl[i]; 61870c78f05SPeter Brune Yl = nasm->y[i]; 619258e1594SPeter Brune subsnes = nasm->subsnes[i]; 620258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 621111ade9eSPeter Brune iscat = nasm->iscatter[i]; 622111ade9eSPeter Brune oscat = nasm->oscatter[i]; 623111ade9eSPeter Brune gscat = nasm->gscatter[i]; 624ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 62524b7f281SPeter Brune if (B) { 62624b7f281SPeter Brune Bl = nasm->b[i]; 627ed3c11a9SPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 628ed3c11a9SPeter Brune } else Bl = NULL; 629ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 63070c78f05SPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 63170c78f05SPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 632111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 633ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 634ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 635111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 636111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 637111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 638111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 639ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 640eaedb033SPeter Brune } 641ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 642ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 64370c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 64470c78f05SPeter Brune Yl = nasm->y[i]; 64570c78f05SPeter Brune iscat = nasm->iscatter[i]; 64670c78f05SPeter Brune oscat = nasm->oscatter[i]; 64770c78f05SPeter Brune if (nasm->type == PC_ASM_BASIC) { 64870c78f05SPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 64970c78f05SPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 65070c78f05SPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 651ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 65270c78f05SPeter Brune } 653b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 654610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 655eaedb033SPeter Brune PetscFunctionReturn(0); 656eaedb033SPeter Brune } 657eaedb033SPeter Brune 658eaedb033SPeter Brune #undef __FUNCT__ 659d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 660602bec5dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 661d728fb7dSPeter Brune { 662602bec5dSPeter Brune Vec X = Xfinal; 663d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 664d728fb7dSPeter Brune SNES subsnes; 665ca44f815SPeter Brune PetscInt i,lag = 1; 666d728fb7dSPeter Brune PetscErrorCode ierr; 667e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 668d728fb7dSPeter Brune VecScatter oscat,gscat; 669d728fb7dSPeter Brune DM dm,subdm; 670d728fb7dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 671d728fb7dSPeter Brune PetscFunctionBegin; 672602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 673e59f0a30SPeter Brune F = snes->vec_func; 674365a6726SPeter Brune if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 675d728fb7dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 676d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 677d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 678602bec5dSPeter Brune if (nasm->fjtype != 1) { 679d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 680d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 681d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 682d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 683d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 684602bec5dSPeter Brune } 685d728fb7dSPeter Brune } 686d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 687d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 688e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 689d728fb7dSPeter Brune Xl = nasm->x[i]; 690d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 691d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 692d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 693d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 694602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 695ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 696d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 697602bec5dSPeter Brune if (nasm->fjtype != 1) { 698d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 699d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 700602bec5dSPeter Brune } 701ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 702ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 703602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 7042f543b25SPeter Brune ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr); 705ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 706d728fb7dSPeter Brune ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr); 707d728fb7dSPeter Brune } 708d728fb7dSPeter Brune PetscFunctionReturn(0); 709d728fb7dSPeter Brune } 710d728fb7dSPeter Brune 711d728fb7dSPeter Brune #undef __FUNCT__ 712eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 713eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 714eaedb033SPeter Brune { 715eaedb033SPeter Brune Vec F; 716eaedb033SPeter Brune Vec X; 717eaedb033SPeter Brune Vec B; 718111ade9eSPeter Brune Vec Y; 719eaedb033SPeter Brune PetscInt i; 720ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 721eaedb033SPeter Brune PetscErrorCode ierr; 722365a6726SPeter Brune SNESNormSchedule normschedule; 723d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 724eaedb033SPeter Brune 725eaedb033SPeter Brune PetscFunctionBegin; 726*fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 727eaedb033SPeter Brune X = snes->vec_sol; 728111ade9eSPeter Brune Y = snes->vec_sol_update; 729eaedb033SPeter Brune F = snes->vec_func; 730eaedb033SPeter Brune B = snes->vec_rhs; 731eaedb033SPeter Brune 732e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 733eaedb033SPeter Brune snes->iter = 0; 734eaedb033SPeter Brune snes->norm = 0.; 735e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 736eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 737365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 738365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 739eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 740eaedb033SPeter Brune if (!snes->vec_func_init_set) { 741eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 742eaedb033SPeter Brune if (snes->domainerror) { 743eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 744eaedb033SPeter Brune PetscFunctionReturn(0); 745eaedb033SPeter Brune } 7461aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 747eaedb033SPeter Brune 748eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 749189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 750189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 751189a9710SBarry Smith PetscFunctionReturn(0); 752189a9710SBarry Smith } 753e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 754eaedb033SPeter Brune snes->iter = 0; 755eaedb033SPeter Brune snes->norm = fnorm; 756e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 757a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 758eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 759eaedb033SPeter Brune 760eaedb033SPeter Brune /* test convergence */ 761eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 762eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 763eaedb033SPeter Brune } else { 764e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 765a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 766eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 767eaedb033SPeter Brune } 768eaedb033SPeter Brune 769eaedb033SPeter Brune /* Call general purpose update function */ 770eaedb033SPeter Brune if (snes->ops->update) { 771eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 772eaedb033SPeter Brune } 773602bec5dSPeter Brune /* copy the initial solution over for later */ 774602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 775eaedb033SPeter Brune 776eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 777111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 778365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 779eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 780eaedb033SPeter Brune if (snes->domainerror) { 781eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 782d728fb7dSPeter Brune break; 783eaedb033SPeter Brune } 784eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 785189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 786189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 787189a9710SBarry Smith break; 788189a9710SBarry Smith } 789eaedb033SPeter Brune } 790eaedb033SPeter Brune /* Monitor convergence */ 791e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 792eaedb033SPeter Brune snes->iter = i+1; 793eaedb033SPeter Brune snes->norm = fnorm; 794e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 795a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 796eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 797eaedb033SPeter Brune /* Test for convergence */ 798365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 799d728fb7dSPeter Brune if (snes->reason) break; 800eaedb033SPeter Brune /* Call general purpose update function */ 801e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 802eaedb033SPeter Brune } 803d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 804365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 805eaedb033SPeter Brune if (i == snes->max_its) { 806eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 807eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 808eaedb033SPeter Brune } 8091aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 810eaedb033SPeter Brune PetscFunctionReturn(0); 811eaedb033SPeter Brune } 812eaedb033SPeter Brune 813eaedb033SPeter Brune /*MC 814eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 815eaedb033SPeter Brune 81670c78603SPeter Brune Options Database: 81770c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 81870c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 81970c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 820602bec5dSPeter Brune . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at 82170c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 82270c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 82370c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 82470c78603SPeter Brune 825eaedb033SPeter Brune Level: advanced 826eaedb033SPeter Brune 827eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 828eaedb033SPeter Brune M*/ 829eaedb033SPeter Brune 830eaedb033SPeter Brune #undef __FUNCT__ 831eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 8328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 833eaedb033SPeter Brune { 834eaedb033SPeter Brune SNES_NASM *nasm; 835eaedb033SPeter Brune PetscErrorCode ierr; 836eaedb033SPeter Brune 837eaedb033SPeter Brune PetscFunctionBegin; 838b00a9115SJed Brown ierr = PetscNewLog(snes,&nasm);CHKERRQ(ierr); 839eaedb033SPeter Brune snes->data = (void*)nasm; 840eaedb033SPeter Brune 841eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 842eaedb033SPeter Brune nasm->subsnes = 0; 843eaedb033SPeter Brune nasm->x = 0; 844111ade9eSPeter Brune nasm->xl = 0; 845111ade9eSPeter Brune nasm->y = 0; 846eaedb033SPeter Brune nasm->b = 0; 847111ade9eSPeter Brune nasm->oscatter = 0; 848111ade9eSPeter Brune nasm->iscatter = 0; 849111ade9eSPeter Brune nasm->gscatter = 0; 850610116beSPeter Brune nasm->damping = 1.; 851111ade9eSPeter Brune 852111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 853d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 854a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 855eaedb033SPeter Brune 856eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 857eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 858eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 859eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 860eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 861eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 862eaedb033SPeter Brune 863eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 864eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 865eaedb033SPeter Brune 866602bec5dSPeter Brune nasm->fjtype = 0; 867602bec5dSPeter Brune nasm->xinit = NULL; 8680298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8690298fd71SBarry Smith nasm->eventsubsolve = 0; 870b20c023fSPeter Brune 871eaedb033SPeter Brune if (!snes->tolerancesset) { 872eaedb033SPeter Brune snes->max_its = 10000; 873eaedb033SPeter Brune snes->max_funcs = 10000; 874eaedb033SPeter Brune } 875eaedb033SPeter Brune 876bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 877bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 87890bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 87990bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 880bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 881bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 882eaedb033SPeter Brune PetscFunctionReturn(0); 883eaedb033SPeter Brune } 88499e0435eSBarry Smith 885