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; 21779f8d099SPeter Brune PetscInt i,j,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); 253a4f17876SPeter Brune for (j=0; j<size; j++) { 254b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 255a4f17876SPeter Brune if (rank == j) { 256a4f17876SPeter Brune for (i=0; i<nasm->n; i++) { 257a4f17876SPeter Brune ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr); 258a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 259b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 260a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 261b20c023fSPeter Brune } 262b20c023fSPeter Brune } 263a4f17876SPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 264a4f17876SPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 265a4f17876SPeter Brune } 266b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 267a4f17876SPeter Brune } 268b20c023fSPeter Brune } else if (isstring) { 269b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 270b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 271b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 272b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 273b20c023fSPeter Brune } 274eaedb033SPeter Brune PetscFunctionReturn(0); 275eaedb033SPeter Brune } 276eaedb033SPeter Brune 277eaedb033SPeter Brune #undef __FUNCT__ 278eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 27976857b2aSPeter Brune /*@ 28076857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 28176857b2aSPeter Brune 28276857b2aSPeter Brune Not Collective 28376857b2aSPeter Brune 28476857b2aSPeter Brune Input Parameters: 28576857b2aSPeter Brune + SNES - the SNES context 28676857b2aSPeter Brune . n - the number of local subdomains 28776857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 28876857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 28976857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 29076857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 29176857b2aSPeter Brune 29276857b2aSPeter Brune Level: intermediate 29376857b2aSPeter Brune 29476857b2aSPeter Brune .keywords: SNES, NASM 29576857b2aSPeter Brune 29676857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 29776857b2aSPeter Brune @*/ 298a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 299a6dfd86eSKarl Rupp { 300eaedb033SPeter Brune PetscErrorCode ierr; 301111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 3026e111a19SKarl Rupp 303eaedb033SPeter Brune PetscFunctionBegin; 3040005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 3054b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 306eaedb033SPeter Brune PetscFunctionReturn(0); 307eaedb033SPeter Brune } 308eaedb033SPeter Brune 309eaedb033SPeter Brune #undef __FUNCT__ 310eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 311a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 312a6dfd86eSKarl Rupp { 313eaedb033SPeter Brune PetscInt i; 314eaedb033SPeter Brune PetscErrorCode ierr; 315eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3166e111a19SKarl Rupp 317eaedb033SPeter Brune PetscFunctionBegin; 318ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 319eaedb033SPeter Brune 320111ade9eSPeter Brune /* tear down the previously set things */ 321111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 322111ade9eSPeter Brune 323eaedb033SPeter Brune nasm->n = n; 324111ade9eSPeter Brune if (oscatter) { 325111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 326eaedb033SPeter Brune } 327111ade9eSPeter Brune if (iscatter) { 328111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 329eaedb033SPeter Brune } 330111ade9eSPeter Brune if (gscatter) { 331111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 332111ade9eSPeter Brune } 333111ade9eSPeter Brune if (oscatter) { 334785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 335eaedb033SPeter Brune for (i=0; i<n; i++) { 336111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 337eaedb033SPeter Brune } 338111ade9eSPeter Brune } 339111ade9eSPeter Brune if (iscatter) { 340785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 341eaedb033SPeter Brune for (i=0; i<n; i++) { 342111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 343eaedb033SPeter Brune } 344eaedb033SPeter Brune } 345111ade9eSPeter Brune if (gscatter) { 346785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 347eaedb033SPeter Brune for (i=0; i<n; i++) { 348111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 349eaedb033SPeter Brune } 350eaedb033SPeter Brune } 351111ade9eSPeter Brune 352eaedb033SPeter Brune if (subsnes) { 353785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 354eaedb033SPeter Brune for (i=0; i<n; i++) { 355eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 356eaedb033SPeter Brune } 357a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 358eaedb033SPeter Brune } 359eaedb033SPeter Brune PetscFunctionReturn(0); 360eaedb033SPeter Brune } 361eaedb033SPeter Brune 362eaedb033SPeter Brune #undef __FUNCT__ 36376857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains" 36476857b2aSPeter Brune /*@ 36576857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 36676857b2aSPeter Brune 36776857b2aSPeter Brune Not Collective 36876857b2aSPeter Brune 36976857b2aSPeter Brune Input Parameters: 37076857b2aSPeter Brune . SNES - the SNES context 37176857b2aSPeter Brune 37276857b2aSPeter Brune Output Parameters: 37376857b2aSPeter Brune + n - the number of local subdomains 37476857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 37576857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 37676857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 37776857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 37876857b2aSPeter Brune 37976857b2aSPeter Brune Level: intermediate 38076857b2aSPeter Brune 38176857b2aSPeter Brune .keywords: SNES, NASM 38276857b2aSPeter Brune 38376857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 38476857b2aSPeter Brune @*/ 38576857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 38676857b2aSPeter Brune { 38776857b2aSPeter Brune PetscErrorCode ierr; 38876857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 38976857b2aSPeter Brune 39076857b2aSPeter Brune PetscFunctionBegin; 3910005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 3924b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 39376857b2aSPeter Brune PetscFunctionReturn(0); 39476857b2aSPeter Brune } 39576857b2aSPeter Brune 39676857b2aSPeter Brune #undef __FUNCT__ 39776857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 39876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 39976857b2aSPeter Brune { 40076857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 40176857b2aSPeter Brune 40276857b2aSPeter Brune PetscFunctionBegin; 40376857b2aSPeter Brune if (n) *n = nasm->n; 40476857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 40576857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 40676857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 407a4f17876SPeter Brune if (subsnes) { 408a4f17876SPeter Brune *subsnes = nasm->subsnes; 409a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 410a4f17876SPeter Brune } 41176857b2aSPeter Brune PetscFunctionReturn(0); 41276857b2aSPeter Brune } 41376857b2aSPeter Brune 41476857b2aSPeter Brune #undef __FUNCT__ 41576857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs" 41676857b2aSPeter Brune /*@ 41776857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 41876857b2aSPeter Brune 41976857b2aSPeter Brune Not Collective 42076857b2aSPeter Brune 42176857b2aSPeter Brune Input Parameters: 42276857b2aSPeter Brune . SNES - the SNES context 42376857b2aSPeter Brune 42476857b2aSPeter Brune Output Parameters: 42576857b2aSPeter Brune + n - the number of local subdomains 42676857b2aSPeter Brune . x - The subdomain solution vector 42776857b2aSPeter Brune . y - The subdomain step vector 42876857b2aSPeter Brune . b - The subdomain RHS vector 42976857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 43076857b2aSPeter Brune 43176857b2aSPeter Brune Level: developer 43276857b2aSPeter Brune 43376857b2aSPeter Brune .keywords: SNES, NASM 43476857b2aSPeter Brune 43576857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 43676857b2aSPeter Brune @*/ 43776857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 43876857b2aSPeter Brune { 43976857b2aSPeter Brune PetscErrorCode ierr; 44076857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 44176857b2aSPeter Brune 44276857b2aSPeter Brune PetscFunctionBegin; 4430005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 4444b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 44576857b2aSPeter Brune PetscFunctionReturn(0); 44676857b2aSPeter Brune } 44776857b2aSPeter Brune 44876857b2aSPeter Brune #undef __FUNCT__ 44976857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 45076857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 45176857b2aSPeter Brune { 45276857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 45376857b2aSPeter Brune 45476857b2aSPeter Brune PetscFunctionBegin; 45576857b2aSPeter Brune if (n) *n = nasm->n; 45676857b2aSPeter Brune if (x) *x = nasm->x; 45776857b2aSPeter Brune if (y) *y = nasm->y; 45876857b2aSPeter Brune if (b) *b = nasm->b; 45976857b2aSPeter Brune if (xl) *xl = nasm->xl; 46076857b2aSPeter Brune PetscFunctionReturn(0); 46176857b2aSPeter Brune } 46276857b2aSPeter Brune 463d728fb7dSPeter Brune #undef __FUNCT__ 464d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 465d728fb7dSPeter Brune /*@ 466d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 467d728fb7dSPeter Brune 468d728fb7dSPeter Brune Collective on SNES 469d728fb7dSPeter Brune 470d728fb7dSPeter Brune Input Parameters: 471d728fb7dSPeter Brune + SNES - the SNES context 472d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 473d728fb7dSPeter Brune 474d728fb7dSPeter Brune Level: developer 475d728fb7dSPeter Brune 476d728fb7dSPeter Brune Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 477d728fb7dSPeter Brune is needed at each linear iteration. 478d728fb7dSPeter Brune 479d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN 480d728fb7dSPeter Brune 481d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 482d728fb7dSPeter Brune @*/ 483d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 484d728fb7dSPeter Brune { 485d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 486d728fb7dSPeter Brune PetscErrorCode ierr; 487d728fb7dSPeter Brune 488d728fb7dSPeter Brune PetscFunctionBegin; 4890005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 4904b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 491d728fb7dSPeter Brune PetscFunctionReturn(0); 492d728fb7dSPeter Brune } 493d728fb7dSPeter Brune 494d728fb7dSPeter Brune #undef __FUNCT__ 495d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 496d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 497d728fb7dSPeter Brune { 498d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 499d728fb7dSPeter Brune 500d728fb7dSPeter Brune PetscFunctionBegin; 501d728fb7dSPeter Brune nasm->finaljacobian = flg; 502602bec5dSPeter Brune if (flg) snes->usesksp = PETSC_TRUE; 503d728fb7dSPeter Brune PetscFunctionReturn(0); 504d728fb7dSPeter Brune } 50576857b2aSPeter Brune 50676857b2aSPeter Brune #undef __FUNCT__ 507610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping" 508610116beSPeter Brune /*@ 509610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 510610116beSPeter Brune 511610116beSPeter Brune Logically collective on SNES 512610116beSPeter Brune 513610116beSPeter Brune Input Parameters: 514610116beSPeter Brune + SNES - the SNES context 515610116beSPeter Brune - dmp - damping 516610116beSPeter Brune 517610116beSPeter Brune Level: intermediate 518610116beSPeter Brune 519610116beSPeter Brune .keywords: SNES, NASM, damping 520610116beSPeter Brune 521610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 522610116beSPeter Brune @*/ 523610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 524610116beSPeter Brune { 525610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 526610116beSPeter Brune PetscErrorCode ierr; 527610116beSPeter Brune 528610116beSPeter Brune PetscFunctionBegin; 529610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 530e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 531610116beSPeter Brune PetscFunctionReturn(0); 532610116beSPeter Brune } 533610116beSPeter Brune 534610116beSPeter Brune #undef __FUNCT__ 535610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping_NASM" 536610116beSPeter Brune PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 537610116beSPeter Brune { 538610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 539610116beSPeter Brune 540610116beSPeter Brune PetscFunctionBegin; 541610116beSPeter Brune nasm->damping = dmp; 542610116beSPeter Brune PetscFunctionReturn(0); 543610116beSPeter Brune } 544610116beSPeter Brune 545610116beSPeter Brune #undef __FUNCT__ 546610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping" 547610116beSPeter Brune /*@ 548610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 549610116beSPeter Brune 550610116beSPeter Brune Not Collective 551610116beSPeter Brune 552610116beSPeter Brune Input Parameters: 553610116beSPeter Brune + SNES - the SNES context 554610116beSPeter Brune - dmp - damping 555610116beSPeter Brune 556610116beSPeter Brune Level: intermediate 557610116beSPeter Brune 558610116beSPeter Brune .keywords: SNES, NASM, damping 559610116beSPeter Brune 560610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping() 561610116beSPeter Brune @*/ 562610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 563610116beSPeter Brune { 564610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal*); 565610116beSPeter Brune PetscErrorCode ierr; 566610116beSPeter Brune 567610116beSPeter Brune PetscFunctionBegin; 568610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 569e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 570610116beSPeter Brune PetscFunctionReturn(0); 571610116beSPeter Brune } 572610116beSPeter Brune 573610116beSPeter Brune #undef __FUNCT__ 574610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping_NASM" 575610116beSPeter Brune PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 576610116beSPeter Brune { 577610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 578610116beSPeter Brune 579610116beSPeter Brune PetscFunctionBegin; 580610116beSPeter Brune *dmp = nasm->damping; 581610116beSPeter Brune PetscFunctionReturn(0); 582610116beSPeter Brune } 583610116beSPeter Brune 584610116beSPeter Brune 585610116beSPeter Brune #undef __FUNCT__ 586eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 5870adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 5880adebc6cSBarry Smith { 589eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 590258e1594SPeter Brune SNES subsnes; 591eaedb033SPeter Brune PetscInt i; 592610116beSPeter Brune PetscReal dmp; 593eaedb033SPeter Brune PetscErrorCode ierr; 594111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 595111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 596111ade9eSPeter Brune DM dm,subdm; 5970adebc6cSBarry Smith 598eaedb033SPeter Brune PetscFunctionBegin; 599eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 600610116beSPeter Brune ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 601111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 602b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 603eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 60470c78f05SPeter Brune /* scatter the solution to the local solution */ 60570c78f05SPeter Brune Xlloc = nasm->xl[i]; 60670c78f05SPeter Brune gscat = nasm->gscatter[i]; 60770c78f05SPeter Brune oscat = nasm->oscatter[i]; 60870c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 60970c78f05SPeter Brune if (B) { 61070c78f05SPeter Brune /* scatter the RHS to the local RHS */ 61170c78f05SPeter Brune Bl = nasm->b[i]; 61270c78f05SPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 61370c78f05SPeter Brune } 61470c78f05SPeter Brune } 615b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 616b20c023fSPeter Brune 617b20c023fSPeter Brune 618b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 61970c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 62070c78f05SPeter Brune Xl = nasm->x[i]; 62170c78f05SPeter Brune Xlloc = nasm->xl[i]; 62270c78f05SPeter Brune Yl = nasm->y[i]; 623258e1594SPeter Brune subsnes = nasm->subsnes[i]; 624258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 625111ade9eSPeter Brune iscat = nasm->iscatter[i]; 626111ade9eSPeter Brune oscat = nasm->oscatter[i]; 627111ade9eSPeter Brune gscat = nasm->gscatter[i]; 628ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 62924b7f281SPeter Brune if (B) { 63024b7f281SPeter Brune Bl = nasm->b[i]; 631ed3c11a9SPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 632ed3c11a9SPeter Brune } else Bl = NULL; 633ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 63470c78f05SPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 63570c78f05SPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 636111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 637ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 638ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 639111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 640111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 641111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 642111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 643ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 644eaedb033SPeter Brune } 645ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 646ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 64770c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 64870c78f05SPeter Brune Yl = nasm->y[i]; 64970c78f05SPeter Brune iscat = nasm->iscatter[i]; 65070c78f05SPeter Brune oscat = nasm->oscatter[i]; 65170c78f05SPeter Brune if (nasm->type == PC_ASM_BASIC) { 65270c78f05SPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 65370c78f05SPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 65470c78f05SPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 655ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 65670c78f05SPeter Brune } 657b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 658610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 659eaedb033SPeter Brune PetscFunctionReturn(0); 660eaedb033SPeter Brune } 661eaedb033SPeter Brune 662eaedb033SPeter Brune #undef __FUNCT__ 663d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 664602bec5dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 665d728fb7dSPeter Brune { 666602bec5dSPeter Brune Vec X = Xfinal; 667d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 668d728fb7dSPeter Brune SNES subsnes; 669ca44f815SPeter Brune PetscInt i,lag = 1; 670d728fb7dSPeter Brune PetscErrorCode ierr; 671e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 672d728fb7dSPeter Brune VecScatter oscat,gscat; 673d728fb7dSPeter Brune DM dm,subdm; 674d728fb7dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 675d728fb7dSPeter Brune PetscFunctionBegin; 676602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 677e59f0a30SPeter Brune F = snes->vec_func; 678365a6726SPeter Brune if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 679d728fb7dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 680d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 681d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 682602bec5dSPeter Brune if (nasm->fjtype != 1) { 683d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 684d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 685d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 686d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 687d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 688602bec5dSPeter Brune } 689d728fb7dSPeter Brune } 690d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 691d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 692e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 693d728fb7dSPeter Brune Xl = nasm->x[i]; 694d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 695d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 696d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 697d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 698602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 699ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 700d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 701602bec5dSPeter Brune if (nasm->fjtype != 1) { 702d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 703d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 704602bec5dSPeter Brune } 705ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 706ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 707602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 7082f543b25SPeter Brune ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr); 709ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 710d728fb7dSPeter Brune ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr); 711d728fb7dSPeter Brune } 712d728fb7dSPeter Brune PetscFunctionReturn(0); 713d728fb7dSPeter Brune } 714d728fb7dSPeter Brune 715d728fb7dSPeter Brune #undef __FUNCT__ 716eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 717eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 718eaedb033SPeter Brune { 719eaedb033SPeter Brune Vec F; 720eaedb033SPeter Brune Vec X; 721eaedb033SPeter Brune Vec B; 722111ade9eSPeter Brune Vec Y; 723eaedb033SPeter Brune PetscInt i; 724ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 725eaedb033SPeter Brune PetscErrorCode ierr; 726365a6726SPeter Brune SNESNormSchedule normschedule; 727d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 728eaedb033SPeter Brune 729eaedb033SPeter Brune PetscFunctionBegin; 730eaedb033SPeter Brune X = snes->vec_sol; 731111ade9eSPeter Brune Y = snes->vec_sol_update; 732eaedb033SPeter Brune F = snes->vec_func; 733eaedb033SPeter Brune B = snes->vec_rhs; 734eaedb033SPeter Brune 735e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 736eaedb033SPeter Brune snes->iter = 0; 737eaedb033SPeter Brune snes->norm = 0.; 738e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 739eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 740365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 741365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 742eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 743eaedb033SPeter Brune if (!snes->vec_func_init_set) { 744eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 745eaedb033SPeter Brune if (snes->domainerror) { 746eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 747eaedb033SPeter Brune PetscFunctionReturn(0); 748eaedb033SPeter Brune } 7491aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 750eaedb033SPeter Brune 751eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 752189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 753189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 754189a9710SBarry Smith PetscFunctionReturn(0); 755189a9710SBarry Smith } 756e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 757eaedb033SPeter Brune snes->iter = 0; 758eaedb033SPeter Brune snes->norm = fnorm; 759e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 760a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 761eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 762eaedb033SPeter Brune 763eaedb033SPeter Brune /* test convergence */ 764eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 765eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 766eaedb033SPeter Brune } else { 767e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 768a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 769eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 770eaedb033SPeter Brune } 771eaedb033SPeter Brune 772eaedb033SPeter Brune /* Call general purpose update function */ 773eaedb033SPeter Brune if (snes->ops->update) { 774eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 775eaedb033SPeter Brune } 776602bec5dSPeter Brune /* copy the initial solution over for later */ 777602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 778eaedb033SPeter Brune 779eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 780111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 781365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 782eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 783eaedb033SPeter Brune if (snes->domainerror) { 784eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 785d728fb7dSPeter Brune break; 786eaedb033SPeter Brune } 787eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 788189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 789189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 790189a9710SBarry Smith break; 791189a9710SBarry Smith } 792eaedb033SPeter Brune } 793eaedb033SPeter Brune /* Monitor convergence */ 794e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 795eaedb033SPeter Brune snes->iter = i+1; 796eaedb033SPeter Brune snes->norm = fnorm; 797e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 798a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 799eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 800eaedb033SPeter Brune /* Test for convergence */ 801365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 802d728fb7dSPeter Brune if (snes->reason) break; 803eaedb033SPeter Brune /* Call general purpose update function */ 804e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 805eaedb033SPeter Brune } 806d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 807365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 808eaedb033SPeter Brune if (i == snes->max_its) { 809eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 810eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 811eaedb033SPeter Brune } 8121aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 813eaedb033SPeter Brune PetscFunctionReturn(0); 814eaedb033SPeter Brune } 815eaedb033SPeter Brune 816eaedb033SPeter Brune /*MC 817eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 818eaedb033SPeter Brune 81970c78603SPeter Brune Options Database: 82070c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 82170c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 82270c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 823602bec5dSPeter Brune . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at 82470c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 82570c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 82670c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 82770c78603SPeter Brune 828eaedb033SPeter Brune Level: advanced 829eaedb033SPeter Brune 830eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 831eaedb033SPeter Brune M*/ 832eaedb033SPeter Brune 833eaedb033SPeter Brune #undef __FUNCT__ 834eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 8358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 836eaedb033SPeter Brune { 837eaedb033SPeter Brune SNES_NASM *nasm; 838eaedb033SPeter Brune PetscErrorCode ierr; 839eaedb033SPeter Brune 840eaedb033SPeter Brune PetscFunctionBegin; 841*b00a9115SJed Brown ierr = PetscNewLog(snes,&nasm);CHKERRQ(ierr); 842eaedb033SPeter Brune snes->data = (void*)nasm; 843eaedb033SPeter Brune 844eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 845eaedb033SPeter Brune nasm->subsnes = 0; 846eaedb033SPeter Brune nasm->x = 0; 847111ade9eSPeter Brune nasm->xl = 0; 848111ade9eSPeter Brune nasm->y = 0; 849eaedb033SPeter Brune nasm->b = 0; 850111ade9eSPeter Brune nasm->oscatter = 0; 851111ade9eSPeter Brune nasm->iscatter = 0; 852111ade9eSPeter Brune nasm->gscatter = 0; 853610116beSPeter Brune nasm->damping = 1.; 854111ade9eSPeter Brune 855111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 856d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 857a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 858eaedb033SPeter Brune 859eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 860eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 861eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 862eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 863eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 864eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 865eaedb033SPeter Brune 866eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 867eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 868eaedb033SPeter Brune 869602bec5dSPeter Brune nasm->fjtype = 0; 870602bec5dSPeter Brune nasm->xinit = NULL; 8710298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8720298fd71SBarry Smith nasm->eventsubsolve = 0; 873b20c023fSPeter Brune 874eaedb033SPeter Brune if (!snes->tolerancesset) { 875eaedb033SPeter Brune snes->max_its = 10000; 876eaedb033SPeter Brune snes->max_funcs = 10000; 877eaedb033SPeter Brune } 878eaedb033SPeter Brune 879bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 880bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 88190bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 88290bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 883bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 884bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 885eaedb033SPeter Brune PetscFunctionReturn(0); 886eaedb033SPeter Brune } 88799e0435eSBarry Smith 888