1af0996ceSBarry Smith #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 */ 11f10b3e88SPatrick Farrell Vec weight; /* weighting for adding updates on overlaps, in global space */ 12111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 13f10b3e88SPatrick Farrell VecScatter *oscatter_copy; /* copy of the above */ 14111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 15111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 16111ade9eSPeter Brune PCASMType type; /* ASM type */ 17111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 18d728fb7dSPeter Brune PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 19610116beSPeter Brune PetscReal damping; /* damping parameter for updates from the blocks */ 20a4f17876SPeter Brune PetscBool same_local_solves; /* flag to determine if the solvers have been individually modified */ 21f10b3e88SPatrick Farrell PetscBool weight_set; /* use a weight in the overlap updates */ 22b20c023fSPeter Brune 23b20c023fSPeter Brune /* logging events */ 24b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 25b20c023fSPeter Brune PetscLogEvent eventsubsolve; 26602bec5dSPeter Brune 27602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 28602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 29eaedb033SPeter Brune } SNES_NASM; 30eaedb033SPeter Brune 31b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 32602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 33b20c023fSPeter Brune 3440244768SBarry Smith static PetscErrorCode SNESReset_NASM(SNES snes) 35eaedb033SPeter Brune { 36eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 37eaedb033SPeter Brune PetscErrorCode ierr; 38eaedb033SPeter Brune PetscInt i; 396e111a19SKarl Rupp 40eaedb033SPeter Brune PetscFunctionBegin; 41eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 42111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 43f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 44111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 45bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 46eaedb033SPeter Brune 47bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 48111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 49f10b3e88SPatrick Farrell if (nasm->oscatter_copy) { ierr = VecScatterDestroy(&nasm->oscatter_copy[i]);CHKERRQ(ierr); } 50111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 51111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 52eaedb033SPeter Brune } 53111ade9eSPeter Brune 540677ede6SBarry Smith ierr = PetscFree(nasm->x);CHKERRQ(ierr); 550677ede6SBarry Smith ierr = PetscFree(nasm->xl);CHKERRQ(ierr); 560677ede6SBarry Smith ierr = PetscFree(nasm->y);CHKERRQ(ierr); 570677ede6SBarry Smith ierr = PetscFree(nasm->b);CHKERRQ(ierr); 58111ade9eSPeter Brune 59602bec5dSPeter Brune if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 60602bec5dSPeter Brune 610677ede6SBarry Smith ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr); 620677ede6SBarry Smith ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr); 63f10b3e88SPatrick Farrell ierr = PetscFree(nasm->oscatter_copy);CHKERRQ(ierr); 640677ede6SBarry Smith ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr); 650677ede6SBarry Smith ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr); 66b20c023fSPeter Brune 67f10b3e88SPatrick Farrell if (nasm->weight_set) { 68f10b3e88SPatrick Farrell ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr); 69f10b3e88SPatrick Farrell } 70f10b3e88SPatrick Farrell 71b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 72b20c023fSPeter Brune nasm->eventsubsolve = 0; 73eaedb033SPeter Brune PetscFunctionReturn(0); 74eaedb033SPeter Brune } 75eaedb033SPeter Brune 7640244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes) 77eaedb033SPeter Brune { 78eaedb033SPeter Brune PetscErrorCode ierr; 796e111a19SKarl Rupp 80eaedb033SPeter Brune PetscFunctionBegin; 81eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 8222d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 83eaedb033SPeter Brune PetscFunctionReturn(0); 84eaedb033SPeter Brune } 85eaedb033SPeter Brune 8640244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 870adebc6cSBarry Smith { 88111ade9eSPeter Brune PetscErrorCode ierr; 89111ade9eSPeter Brune Vec bcs = (Vec)ctx; 906e111a19SKarl Rupp 91111ade9eSPeter Brune PetscFunctionBegin; 92111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 93111ade9eSPeter Brune PetscFunctionReturn(0); 94111ade9eSPeter Brune } 95111ade9eSPeter Brune 9640244768SBarry Smith static PetscErrorCode SNESSetUp_NASM(SNES snes) 97eaedb033SPeter Brune { 98eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 99eaedb033SPeter Brune PetscErrorCode ierr; 10076857b2aSPeter Brune DM dm,subdm; 101111ade9eSPeter Brune DM *subdms; 102111ade9eSPeter Brune PetscInt i; 103eaedb033SPeter Brune const char *optionsprefix; 104111ade9eSPeter Brune Vec F; 105ed3c11a9SPeter Brune PetscMPIInt size; 106ed3c11a9SPeter Brune KSP ksp; 107ed3c11a9SPeter Brune PC pc; 108eaedb033SPeter Brune 109eaedb033SPeter Brune PetscFunctionBegin; 110eaedb033SPeter Brune if (!nasm->subsnes) { 111eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1120a696f66SPeter Brune if (dm) { 113eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1140298fd71SBarry Smith ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 115ce94432eSBarry Smith if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 116111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 117f10b3e88SPatrick Farrell ierr = PetscMalloc1(nasm->n, &nasm->oscatter_copy);CHKERRQ(ierr); 118f10b3e88SPatrick Farrell for (i=0; i<nasm->n; i++) { 119f10b3e88SPatrick Farrell ierr = VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]);CHKERRQ(ierr); 120f10b3e88SPatrick Farrell } 121eaedb033SPeter Brune 122eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 123785e854fSJed Brown ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr); 124111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 125cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 1267a2f0ea0SPatrick Farrell ierr = PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1);CHKERRQ(ierr); 127cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 128cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 129cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 130ed3c11a9SPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 131ed3c11a9SPeter Brune if (size == 1) { 132ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 133ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 134ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 135ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 136ed3c11a9SPeter Brune } 137cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 138111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 139111ade9eSPeter Brune } 140111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 141ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 142ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 143111ade9eSPeter Brune /* allocate the global vectors */ 144e245e0aaSPeter Brune if (!nasm->x) { 1451795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr); 146e245e0aaSPeter Brune } 147e245e0aaSPeter Brune if (!nasm->xl) { 1481795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr); 149e245e0aaSPeter Brune } 150e245e0aaSPeter Brune if (!nasm->y) { 1511795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr); 152e245e0aaSPeter Brune } 153e245e0aaSPeter Brune if (!nasm->b) { 1541795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr); 155e245e0aaSPeter Brune } 156111ade9eSPeter Brune 157111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1580298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 15976857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 16076857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 16176857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 16276857b2aSPeter Brune if (!nasm->xl[i]) { 163111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 164111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 1650298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 166111ade9eSPeter Brune } 16761ba4676SBarry Smith } 168602bec5dSPeter Brune if (nasm->finaljacobian) { 169602bec5dSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 170602bec5dSPeter Brune if (nasm->fjtype == 2) { 171602bec5dSPeter Brune ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 172602bec5dSPeter Brune } 173602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 174602bec5dSPeter Brune ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 175602bec5dSPeter Brune } 176602bec5dSPeter Brune } 177eaedb033SPeter Brune PetscFunctionReturn(0); 178eaedb033SPeter Brune } 179eaedb033SPeter Brune 18040244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes) 181eaedb033SPeter Brune { 182eaedb033SPeter Brune PetscErrorCode ierr; 183111ade9eSPeter Brune PCASMType asmtype; 184a4f17876SPeter Brune PetscBool flg,monflg,subviewflg; 185111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1866e111a19SKarl Rupp 187eaedb033SPeter Brune PetscFunctionBegin; 188e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwartz options");CHKERRQ(ierr); 189111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 190e0331734SPeter Brune if (flg) {ierr = SNESNASMSetType(snes,asmtype);CHKERRQ(ierr);} 191b20c023fSPeter Brune flg = PETSC_FALSE; 192b20c023fSPeter Brune monflg = PETSC_TRUE; 1935dfa0f3bSBarry Smith ierr = PetscOptionsReal("-snes_nasm_damping","The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr); 194610116beSPeter Brune if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 195a4f17876SPeter Brune subviewflg = PETSC_FALSE; 196a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr); 197a4f17876SPeter Brune if (flg) { 198a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 199a4f17876SPeter Brune if (!subviewflg) { 200a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 201a4f17876SPeter Brune } 202a4f17876SPeter Brune } 203d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 204602bec5dSPeter Brune ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 205a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 206b20c023fSPeter Brune if (flg) { 207b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 208b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 209b20c023fSPeter Brune } 210eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 211eaedb033SPeter Brune PetscFunctionReturn(0); 212eaedb033SPeter Brune } 213eaedb033SPeter Brune 21440244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 215eaedb033SPeter Brune { 216b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 217b20c023fSPeter Brune PetscErrorCode ierr; 218a4f17876SPeter Brune PetscMPIInt rank,size; 219dd2fa690SBarry Smith PetscInt i,N,bsz; 220b20c023fSPeter Brune PetscBool iascii,isstring; 221b20c023fSPeter Brune PetscViewer sviewer; 222ce94432eSBarry Smith MPI_Comm comm; 223b20c023fSPeter Brune 224b20c023fSPeter Brune PetscFunctionBegin; 225ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 226b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 227b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 228b20c023fSPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 229a4f17876SPeter Brune ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 230b2566f29SBarry Smith ierr = MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 231b20c023fSPeter Brune if (iascii) { 232efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %D\n",N);CHKERRQ(ierr); 233a4f17876SPeter Brune if (nasm->same_local_solves) { 234a4f17876SPeter Brune if (nasm->subsnes) { 235a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve is the same for all blocks:\n");CHKERRQ(ierr); 236a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2373f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 238a4f17876SPeter Brune if (!rank) { 239a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 240a4f17876SPeter Brune ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr); 241a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 242a4f17876SPeter Brune } 2433f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 244a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 245a4f17876SPeter Brune } 246a4f17876SPeter Brune } else { 247a4f17876SPeter Brune /* print the solver on each block */ 2481575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 249b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 250b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2511575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 252a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr); 253b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 254a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2553f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 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 } 2623f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 263a4f17876SPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 264b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 265a4f17876SPeter Brune } 266b20c023fSPeter Brune } else if (isstring) { 267b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 2683f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 269b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 2703f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 271b20c023fSPeter Brune } 272eaedb033SPeter Brune PetscFunctionReturn(0); 273eaedb033SPeter Brune } 274eaedb033SPeter Brune 275e0331734SPeter Brune /*@ 276e0331734SPeter Brune SNESNASMSetType - Set the type of subdomain update used 277e0331734SPeter Brune 278e0331734SPeter Brune Logically Collective on SNES 279e0331734SPeter Brune 280e0331734SPeter Brune Input Parameters: 281e0331734SPeter Brune + SNES - the SNES context 282e0331734SPeter Brune - type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT 283e0331734SPeter Brune 284e0331734SPeter Brune Level: intermediate 285e0331734SPeter Brune 286e0331734SPeter Brune .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType() 287e0331734SPeter Brune @*/ 288e0331734SPeter Brune PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type) 289e0331734SPeter Brune { 290e0331734SPeter Brune PetscErrorCode ierr; 291e0331734SPeter Brune PetscErrorCode (*f)(SNES,PCASMType); 292e0331734SPeter Brune 293e0331734SPeter Brune PetscFunctionBegin; 294e0331734SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f);CHKERRQ(ierr); 295e0331734SPeter Brune if (f) {ierr = (f)(snes,type);CHKERRQ(ierr);} 296e0331734SPeter Brune PetscFunctionReturn(0); 297e0331734SPeter Brune } 298e0331734SPeter Brune 29940244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type) 300e0331734SPeter Brune { 301e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 302e0331734SPeter Brune 303e0331734SPeter Brune PetscFunctionBegin; 304e0331734SPeter Brune if (type != PC_ASM_BASIC && type != PC_ASM_RESTRICT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESNASM only supports basic and restrict types"); 305e0331734SPeter Brune nasm->type = type; 306e0331734SPeter Brune PetscFunctionReturn(0); 307e0331734SPeter Brune } 308e0331734SPeter Brune 309e0331734SPeter Brune /*@ 310e0331734SPeter Brune SNESNASMGetType - Get the type of subdomain update used 311e0331734SPeter Brune 312e0331734SPeter Brune Logically Collective on SNES 313e0331734SPeter Brune 314e0331734SPeter Brune Input Parameters: 315e0331734SPeter Brune . SNES - the SNES context 316e0331734SPeter Brune 317e0331734SPeter Brune Output Parameters: 318e0331734SPeter Brune . type - the type of update 319e0331734SPeter Brune 320e0331734SPeter Brune Level: intermediate 321e0331734SPeter Brune 322e0331734SPeter Brune .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType() 323e0331734SPeter Brune @*/ 324e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type) 325e0331734SPeter Brune { 326e0331734SPeter Brune PetscErrorCode ierr; 327e0331734SPeter Brune 328e0331734SPeter Brune PetscFunctionBegin; 3292a808120SBarry Smith ierr = PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));CHKERRQ(ierr); 330e0331734SPeter Brune PetscFunctionReturn(0); 331e0331734SPeter Brune } 332e0331734SPeter Brune 33340244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type) 334e0331734SPeter Brune { 335e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 336e0331734SPeter Brune 337e0331734SPeter Brune PetscFunctionBegin; 338e0331734SPeter Brune *type = nasm->type; 339e0331734SPeter Brune PetscFunctionReturn(0); 340e0331734SPeter Brune } 341e0331734SPeter Brune 34276857b2aSPeter Brune /*@ 34376857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 34476857b2aSPeter Brune 34576857b2aSPeter Brune Not Collective 34676857b2aSPeter Brune 34776857b2aSPeter Brune Input Parameters: 34876857b2aSPeter Brune + SNES - the SNES context 34976857b2aSPeter Brune . n - the number of local subdomains 35076857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 35176857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 35276857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 35376857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 35476857b2aSPeter Brune 35576857b2aSPeter Brune Level: intermediate 35676857b2aSPeter Brune 35776857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 35876857b2aSPeter Brune @*/ 359a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 360a6dfd86eSKarl Rupp { 361eaedb033SPeter Brune PetscErrorCode ierr; 362111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 3636e111a19SKarl Rupp 364eaedb033SPeter Brune PetscFunctionBegin; 3650005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 3664b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 367eaedb033SPeter Brune PetscFunctionReturn(0); 368eaedb033SPeter Brune } 369eaedb033SPeter Brune 37040244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 371a6dfd86eSKarl Rupp { 372eaedb033SPeter Brune PetscInt i; 373eaedb033SPeter Brune PetscErrorCode ierr; 374eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3756e111a19SKarl Rupp 376eaedb033SPeter Brune PetscFunctionBegin; 377ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 378eaedb033SPeter Brune 379111ade9eSPeter Brune /* tear down the previously set things */ 380111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 381111ade9eSPeter Brune 382eaedb033SPeter Brune nasm->n = n; 383111ade9eSPeter Brune if (oscatter) { 384111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 385eaedb033SPeter Brune } 386111ade9eSPeter Brune if (iscatter) { 387111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 388eaedb033SPeter Brune } 389111ade9eSPeter Brune if (gscatter) { 390111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 391111ade9eSPeter Brune } 392111ade9eSPeter Brune if (oscatter) { 393785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 394f10b3e88SPatrick Farrell ierr = PetscMalloc1(n,&nasm->oscatter_copy);CHKERRQ(ierr); 395eaedb033SPeter Brune for (i=0; i<n; i++) { 396111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 397f10b3e88SPatrick Farrell ierr = VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]);CHKERRQ(ierr); 398eaedb033SPeter Brune } 399111ade9eSPeter Brune } 400111ade9eSPeter Brune if (iscatter) { 401785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 402eaedb033SPeter Brune for (i=0; i<n; i++) { 403111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 404eaedb033SPeter Brune } 405eaedb033SPeter Brune } 406111ade9eSPeter Brune if (gscatter) { 407785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 408eaedb033SPeter Brune for (i=0; i<n; i++) { 409111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 410eaedb033SPeter Brune } 411eaedb033SPeter Brune } 412111ade9eSPeter Brune 413eaedb033SPeter Brune if (subsnes) { 414785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 415eaedb033SPeter Brune for (i=0; i<n; i++) { 416eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 417eaedb033SPeter Brune } 418a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 419eaedb033SPeter Brune } 420eaedb033SPeter Brune PetscFunctionReturn(0); 421eaedb033SPeter Brune } 422eaedb033SPeter Brune 42376857b2aSPeter Brune /*@ 42476857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 42576857b2aSPeter Brune 42676857b2aSPeter Brune Not Collective 42776857b2aSPeter Brune 42876857b2aSPeter Brune Input Parameters: 42976857b2aSPeter Brune . SNES - the SNES context 43076857b2aSPeter Brune 43176857b2aSPeter Brune Output Parameters: 43276857b2aSPeter Brune + n - the number of local subdomains 43376857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 43476857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 43576857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 43676857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 43776857b2aSPeter Brune 43876857b2aSPeter Brune Level: intermediate 43976857b2aSPeter Brune 44076857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 44176857b2aSPeter Brune @*/ 44276857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 44376857b2aSPeter Brune { 44476857b2aSPeter Brune PetscErrorCode ierr; 44576857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 44676857b2aSPeter Brune 44776857b2aSPeter Brune PetscFunctionBegin; 4480005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 4494b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 45076857b2aSPeter Brune PetscFunctionReturn(0); 45176857b2aSPeter Brune } 45276857b2aSPeter Brune 45340244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 45476857b2aSPeter Brune { 45576857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 45676857b2aSPeter Brune 45776857b2aSPeter Brune PetscFunctionBegin; 45876857b2aSPeter Brune if (n) *n = nasm->n; 45976857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 46076857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 46176857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 462a4f17876SPeter Brune if (subsnes) { 463a4f17876SPeter Brune *subsnes = nasm->subsnes; 464a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 465a4f17876SPeter Brune } 46676857b2aSPeter Brune PetscFunctionReturn(0); 46776857b2aSPeter Brune } 46876857b2aSPeter Brune 46976857b2aSPeter Brune /*@ 47076857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 47176857b2aSPeter Brune 47276857b2aSPeter Brune Not Collective 47376857b2aSPeter Brune 47476857b2aSPeter Brune Input Parameters: 47576857b2aSPeter Brune . SNES - the SNES context 47676857b2aSPeter Brune 47776857b2aSPeter Brune Output Parameters: 47876857b2aSPeter Brune + n - the number of local subdomains 47976857b2aSPeter Brune . x - The subdomain solution vector 48076857b2aSPeter Brune . y - The subdomain step vector 48176857b2aSPeter Brune . b - The subdomain RHS vector 48276857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 48376857b2aSPeter Brune 48476857b2aSPeter Brune Level: developer 48576857b2aSPeter Brune 48676857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 48776857b2aSPeter Brune @*/ 48876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 48976857b2aSPeter Brune { 49076857b2aSPeter Brune PetscErrorCode ierr; 49176857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 49276857b2aSPeter Brune 49376857b2aSPeter Brune PetscFunctionBegin; 4940005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 4954b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 49676857b2aSPeter Brune PetscFunctionReturn(0); 49776857b2aSPeter Brune } 49876857b2aSPeter Brune 49940244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 50076857b2aSPeter Brune { 50176857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 50276857b2aSPeter Brune 50376857b2aSPeter Brune PetscFunctionBegin; 50476857b2aSPeter Brune if (n) *n = nasm->n; 50576857b2aSPeter Brune if (x) *x = nasm->x; 50676857b2aSPeter Brune if (y) *y = nasm->y; 50776857b2aSPeter Brune if (b) *b = nasm->b; 50876857b2aSPeter Brune if (xl) *xl = nasm->xl; 50976857b2aSPeter Brune PetscFunctionReturn(0); 51076857b2aSPeter Brune } 51176857b2aSPeter Brune 512d728fb7dSPeter Brune /*@ 513d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 514d728fb7dSPeter Brune 515d728fb7dSPeter Brune Collective on SNES 516d728fb7dSPeter Brune 517d728fb7dSPeter Brune Input Parameters: 518d728fb7dSPeter Brune + SNES - the SNES context 519d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 520d728fb7dSPeter Brune 521d728fb7dSPeter Brune Level: developer 522d728fb7dSPeter Brune 52395452b02SPatrick Sanan Notes: 52495452b02SPatrick Sanan This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 525d728fb7dSPeter Brune is needed at each linear iteration. 526d728fb7dSPeter Brune 527d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 528d728fb7dSPeter Brune @*/ 529d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 530d728fb7dSPeter Brune { 531d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 532d728fb7dSPeter Brune PetscErrorCode ierr; 533d728fb7dSPeter Brune 534d728fb7dSPeter Brune PetscFunctionBegin; 5350005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 5364b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 537d728fb7dSPeter Brune PetscFunctionReturn(0); 538d728fb7dSPeter Brune } 539d728fb7dSPeter Brune 54040244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 541d728fb7dSPeter Brune { 542d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 543d728fb7dSPeter Brune 544d728fb7dSPeter Brune PetscFunctionBegin; 545d728fb7dSPeter Brune nasm->finaljacobian = flg; 546d728fb7dSPeter Brune PetscFunctionReturn(0); 547d728fb7dSPeter Brune } 54876857b2aSPeter Brune 549610116beSPeter Brune /*@ 550610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 551610116beSPeter Brune 552610116beSPeter Brune Logically collective on SNES 553610116beSPeter Brune 554610116beSPeter Brune Input Parameters: 555610116beSPeter Brune + SNES - the SNES context 556610116beSPeter Brune - dmp - damping 557610116beSPeter Brune 558610116beSPeter Brune Level: intermediate 559610116beSPeter Brune 56095452b02SPatrick Sanan Notes: 56195452b02SPatrick Sanan The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5625dfa0f3bSBarry Smith 563610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 564610116beSPeter Brune @*/ 565610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 566610116beSPeter Brune { 567610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 568610116beSPeter Brune PetscErrorCode ierr; 569610116beSPeter Brune 570610116beSPeter Brune PetscFunctionBegin; 571610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 572e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 573610116beSPeter Brune PetscFunctionReturn(0); 574610116beSPeter Brune } 575610116beSPeter Brune 57640244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 577610116beSPeter Brune { 578610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 579610116beSPeter Brune 580610116beSPeter Brune PetscFunctionBegin; 581610116beSPeter Brune nasm->damping = dmp; 582610116beSPeter Brune PetscFunctionReturn(0); 583610116beSPeter Brune } 584610116beSPeter Brune 585610116beSPeter Brune /*@ 586610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 587610116beSPeter Brune 588610116beSPeter Brune Not Collective 589610116beSPeter Brune 590610116beSPeter Brune Input Parameters: 591610116beSPeter Brune + SNES - the SNES context 592610116beSPeter Brune - dmp - damping 593610116beSPeter Brune 594610116beSPeter Brune Level: intermediate 595610116beSPeter Brune 596610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping() 597610116beSPeter Brune @*/ 598610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 599610116beSPeter Brune { 600610116beSPeter Brune PetscErrorCode ierr; 601610116beSPeter Brune 602610116beSPeter Brune PetscFunctionBegin; 6034a2f8832SBarry Smith ierr = PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));CHKERRQ(ierr); 604610116beSPeter Brune PetscFunctionReturn(0); 605610116beSPeter Brune } 606610116beSPeter Brune 60740244768SBarry Smith static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 608610116beSPeter Brune { 609610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 610610116beSPeter Brune 611610116beSPeter Brune PetscFunctionBegin; 612610116beSPeter Brune *dmp = nasm->damping; 613610116beSPeter Brune PetscFunctionReturn(0); 614610116beSPeter Brune } 615610116beSPeter Brune 616610116beSPeter Brune 61714eb1c5cSMatthew G. Knepley /* 61814eb1c5cSMatthew G. Knepley Input Parameters: 61914eb1c5cSMatthew G. Knepley + snes - The solver 62014eb1c5cSMatthew G. Knepley . B - The RHS vector 62114eb1c5cSMatthew G. Knepley - X - The initial guess 62214eb1c5cSMatthew G. Knepley 62314eb1c5cSMatthew G. Knepley Output Parameters: 62414eb1c5cSMatthew G. Knepley . Y - The solution update 62514eb1c5cSMatthew G. Knepley 62614eb1c5cSMatthew G. Knepley TODO: All scatters should be packed into one 62714eb1c5cSMatthew G. Knepley */ 6280adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 6290adebc6cSBarry Smith { 630eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 631258e1594SPeter Brune SNES subsnes; 632eaedb033SPeter Brune PetscInt i; 633610116beSPeter Brune PetscReal dmp; 634eaedb033SPeter Brune PetscErrorCode ierr; 635f10b3e88SPatrick Farrell Vec Xl,Bl,Yl,Xlloc; 636f10b3e88SPatrick Farrell VecScatter iscat,oscat,gscat,oscat_copy; 637111ade9eSPeter Brune DM dm,subdm; 638e0331734SPeter Brune PCASMType type; 6390adebc6cSBarry Smith 640eaedb033SPeter Brune PetscFunctionBegin; 641e0331734SPeter Brune ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr); 642eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 643111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 644b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 645eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 646f10b3e88SPatrick Farrell /* scatter the solution to the global solution and the local solution */ 647f10b3e88SPatrick Farrell Xl = nasm->x[i]; 64870c78f05SPeter Brune Xlloc = nasm->xl[i]; 64970c78f05SPeter Brune oscat = nasm->oscatter[i]; 650f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 651f10b3e88SPatrick Farrell gscat = nasm->gscatter[i]; 652f10b3e88SPatrick Farrell ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 65370c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 65470c78f05SPeter Brune if (B) { 65570c78f05SPeter Brune /* scatter the RHS to the local RHS */ 65670c78f05SPeter Brune Bl = nasm->b[i]; 657f10b3e88SPatrick Farrell ierr = VecScatterBegin(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 65870c78f05SPeter Brune } 65970c78f05SPeter Brune } 660b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 661b20c023fSPeter Brune 662b20c023fSPeter Brune 663b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 66470c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 66570c78f05SPeter Brune Xl = nasm->x[i]; 66670c78f05SPeter Brune Xlloc = nasm->xl[i]; 66770c78f05SPeter Brune Yl = nasm->y[i]; 668258e1594SPeter Brune subsnes = nasm->subsnes[i]; 669258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 670111ade9eSPeter Brune iscat = nasm->iscatter[i]; 671111ade9eSPeter Brune oscat = nasm->oscatter[i]; 672f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 673111ade9eSPeter Brune gscat = nasm->gscatter[i]; 674f10b3e88SPatrick Farrell ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 675ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 67624b7f281SPeter Brune if (B) { 67724b7f281SPeter Brune Bl = nasm->b[i]; 678f10b3e88SPatrick Farrell ierr = VecScatterEnd(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 679ed3c11a9SPeter Brune } else Bl = NULL; 680f10b3e88SPatrick Farrell 681ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 682111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 683ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 684ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 685f10b3e88SPatrick Farrell ierr = VecScale(Yl, nasm->damping);CHKERRQ(ierr); 686e0331734SPeter Brune if (type == PC_ASM_BASIC) { 687111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 688*1ed9ada7SJunchao Zhang ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 689e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 690111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 691*1ed9ada7SJunchao Zhang ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 692ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 693eaedb033SPeter Brune } 694ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 695ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 696f10b3e88SPatrick Farrell if (nasm->weight_set) { 697f10b3e88SPatrick Farrell ierr = VecPointwiseMult(Y,Y,nasm->weight);CHKERRQ(ierr); 698f10b3e88SPatrick Farrell } 699b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 7005dfa0f3bSBarry Smith ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 701610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 702eaedb033SPeter Brune PetscFunctionReturn(0); 703eaedb033SPeter Brune } 704eaedb033SPeter Brune 70540244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 706d728fb7dSPeter Brune { 707602bec5dSPeter Brune Vec X = Xfinal; 708d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 709d728fb7dSPeter Brune SNES subsnes; 710ca44f815SPeter Brune PetscInt i,lag = 1; 711d728fb7dSPeter Brune PetscErrorCode ierr; 712e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 713d728fb7dSPeter Brune VecScatter oscat,gscat; 714d728fb7dSPeter Brune DM dm,subdm; 715d1e9a80fSBarry Smith 716d728fb7dSPeter Brune PetscFunctionBegin; 717602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 718e59f0a30SPeter Brune F = snes->vec_func; 719365a6726SPeter Brune if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 720d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 721d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 722d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 723602bec5dSPeter Brune if (nasm->fjtype != 1) { 724d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 725d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 726d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 727d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 728602bec5dSPeter Brune } 729d728fb7dSPeter Brune } 730d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 731d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 732e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 733d728fb7dSPeter Brune Xl = nasm->x[i]; 734d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 735d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 736d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 737d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 738602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 739ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 740d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 741602bec5dSPeter Brune if (nasm->fjtype != 1) { 742d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 743d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 744602bec5dSPeter Brune } 745ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 746ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 747602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 748d1e9a80fSBarry Smith ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr); 749ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 750d728fb7dSPeter Brune } 751d728fb7dSPeter Brune PetscFunctionReturn(0); 752d728fb7dSPeter Brune } 753d728fb7dSPeter Brune 75440244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes) 755eaedb033SPeter Brune { 756eaedb033SPeter Brune Vec F; 757eaedb033SPeter Brune Vec X; 758eaedb033SPeter Brune Vec B; 759111ade9eSPeter Brune Vec Y; 760eaedb033SPeter Brune PetscInt i; 761ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 762eaedb033SPeter Brune PetscErrorCode ierr; 763365a6726SPeter Brune SNESNormSchedule normschedule; 764d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 765eaedb033SPeter Brune 766eaedb033SPeter Brune PetscFunctionBegin; 767c579b300SPatrick Farrell 7686c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 769c579b300SPatrick Farrell 770fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 771eaedb033SPeter Brune X = snes->vec_sol; 772111ade9eSPeter Brune Y = snes->vec_sol_update; 773eaedb033SPeter Brune F = snes->vec_func; 774eaedb033SPeter Brune B = snes->vec_rhs; 775eaedb033SPeter Brune 776e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 777eaedb033SPeter Brune snes->iter = 0; 778eaedb033SPeter Brune snes->norm = 0.; 779e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 780eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 781365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 782365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 783eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 784eaedb033SPeter Brune if (!snes->vec_func_init_set) { 785eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 7861aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 787eaedb033SPeter Brune 788eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 789422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 790e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 791eaedb033SPeter Brune snes->iter = 0; 792eaedb033SPeter Brune snes->norm = fnorm; 793e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 794a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 795eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 796eaedb033SPeter Brune 797eaedb033SPeter Brune /* test convergence */ 798eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 799eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 800eaedb033SPeter Brune } else { 801e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 802a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 803eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 804eaedb033SPeter Brune } 805eaedb033SPeter Brune 806eaedb033SPeter Brune /* Call general purpose update function */ 807eaedb033SPeter Brune if (snes->ops->update) { 808eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 809eaedb033SPeter Brune } 810602bec5dSPeter Brune /* copy the initial solution over for later */ 811602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 812eaedb033SPeter Brune 813eaedb033SPeter Brune for (i=0; i < snes->max_its; i++) { 814111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 815365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 816eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 817eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 818422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 819eaedb033SPeter Brune } 820eaedb033SPeter Brune /* Monitor convergence */ 821e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 822eaedb033SPeter Brune snes->iter = i+1; 823eaedb033SPeter Brune snes->norm = fnorm; 824e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 825a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 826eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 827eaedb033SPeter Brune /* Test for convergence */ 828365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 829d728fb7dSPeter Brune if (snes->reason) break; 830eaedb033SPeter Brune /* Call general purpose update function */ 831e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 832eaedb033SPeter Brune } 83307b62357SFande Kong if (nasm->finaljacobian) { 83407b62357SFande Kong ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr); 83507b62357SFande Kong SNESCheckJacobianDomainerror(snes); 83607b62357SFande Kong } 837365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 838eaedb033SPeter Brune if (i == snes->max_its) { 839eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 840eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 841eaedb033SPeter Brune } 8421aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 843eaedb033SPeter Brune PetscFunctionReturn(0); 844eaedb033SPeter Brune } 845eaedb033SPeter Brune 846eaedb033SPeter Brune /*MC 847eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 848eaedb033SPeter Brune 84970c78603SPeter Brune Options Database: 85070c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 85170c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 8525dfa0f3bSBarry Smith . -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 85370c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 854150d49b7SHong Zhang . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at 85570c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 85670c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 85770c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 85870c78603SPeter Brune 859eaedb033SPeter Brune Level: advanced 860eaedb033SPeter Brune 861951fe5abSBarry Smith Developer Note: This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to 862951fe5abSBarry Smith false and SNESView() and -snes_view do not display a KSP object. However the flag nasm->finaljacobian is set (for example if 863951fe5abSBarry Smith NASM is used as a nonlinear preconditioner for KSPASPIN) then SNESSetUpMatrices() is called to generate the Jacobian (needed by KSPASPIN) 864951fe5abSBarry Smith and this utilizes the KSP for storing the matrices, but the KSP is never used for solving a linear system. Note that when SNESNASM is 865951fe5abSBarry Smith used by SNESASPIN they share the same Jacobian matrices because SNESSetUp() (called on the outer SNES KSPASPIN) causes the inner SNES 866951fe5abSBarry Smith object (in this case SNESNASM) to inherit the outer Jacobian matrices. 867951fe5abSBarry Smith 8684f02bc6aSBarry Smith References: 86996a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8704f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8714f02bc6aSBarry Smith 8725dfa0f3bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping() 873eaedb033SPeter Brune M*/ 874eaedb033SPeter Brune 8758cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 876eaedb033SPeter Brune { 877eaedb033SPeter Brune SNES_NASM *nasm; 878eaedb033SPeter Brune PetscErrorCode ierr; 879eaedb033SPeter Brune 880eaedb033SPeter Brune PetscFunctionBegin; 881b00a9115SJed Brown ierr = PetscNewLog(snes,&nasm);CHKERRQ(ierr); 882eaedb033SPeter Brune snes->data = (void*)nasm; 883eaedb033SPeter Brune 884eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 885eaedb033SPeter Brune nasm->subsnes = 0; 886eaedb033SPeter Brune nasm->x = 0; 887111ade9eSPeter Brune nasm->xl = 0; 888111ade9eSPeter Brune nasm->y = 0; 889eaedb033SPeter Brune nasm->b = 0; 890111ade9eSPeter Brune nasm->oscatter = 0; 891f10b3e88SPatrick Farrell nasm->oscatter_copy = 0; 892111ade9eSPeter Brune nasm->iscatter = 0; 893111ade9eSPeter Brune nasm->gscatter = 0; 894610116beSPeter Brune nasm->damping = 1.; 895111ade9eSPeter Brune 896111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 897d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 898a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 899f10b3e88SPatrick Farrell nasm->weight_set = PETSC_FALSE; 900eaedb033SPeter Brune 901eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 902eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 903eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 904eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 905eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 906eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 907eaedb033SPeter Brune 908eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 909efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 910eaedb033SPeter Brune 9114fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 9124fc747eaSLawrence Mitchell 913602bec5dSPeter Brune nasm->fjtype = 0; 914602bec5dSPeter Brune nasm->xinit = NULL; 9150298fd71SBarry Smith nasm->eventrestrictinterp = 0; 9160298fd71SBarry Smith nasm->eventsubsolve = 0; 917b20c023fSPeter Brune 918eaedb033SPeter Brune if (!snes->tolerancesset) { 919eaedb033SPeter Brune snes->max_its = 10000; 920eaedb033SPeter Brune snes->max_funcs = 10000; 921eaedb033SPeter Brune } 922eaedb033SPeter Brune 923e0331734SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr); 924e0331734SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr); 925bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 926bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 92790bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 92890bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 929bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 930bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 931eaedb033SPeter Brune PetscFunctionReturn(0); 932eaedb033SPeter Brune } 93399e0435eSBarry Smith 934448b6425SPatrick Farrell /*@ 935448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 936448b6425SPatrick Farrell 937448b6425SPatrick Farrell Not collective 938448b6425SPatrick Farrell 939448b6425SPatrick Farrell Input Parameters: 940448b6425SPatrick Farrell + snes - the SNES context 941448b6425SPatrick Farrell - i - the number of the subsnes to get 942448b6425SPatrick Farrell 943448b6425SPatrick Farrell Output Parameters: 944448b6425SPatrick Farrell . subsnes - the subsolver context 945448b6425SPatrick Farrell 946448b6425SPatrick Farrell Level: intermediate 947448b6425SPatrick Farrell 948448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetNumber() 949448b6425SPatrick Farrell @*/ 950448b6425SPatrick Farrell PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes) 951448b6425SPatrick Farrell { 952448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 953448b6425SPatrick Farrell 954448b6425SPatrick Farrell PetscFunctionBegin; 955448b6425SPatrick Farrell if (i < 0 || i >= nasm->n) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver"); 956448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 957448b6425SPatrick Farrell PetscFunctionReturn(0); 958448b6425SPatrick Farrell } 959448b6425SPatrick Farrell 960448b6425SPatrick Farrell /*@ 961448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 962448b6425SPatrick Farrell 963448b6425SPatrick Farrell Not collective 964448b6425SPatrick Farrell 965448b6425SPatrick Farrell Input Parameters: 966448b6425SPatrick Farrell . snes - the SNES context 967448b6425SPatrick Farrell 968448b6425SPatrick Farrell Output Parameters: 969448b6425SPatrick Farrell . n - the number of subsolvers 970448b6425SPatrick Farrell 971448b6425SPatrick Farrell Level: intermediate 972448b6425SPatrick Farrell 973448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetSNES() 974448b6425SPatrick Farrell @*/ 975448b6425SPatrick Farrell PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n) 976448b6425SPatrick Farrell { 977448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 978448b6425SPatrick Farrell 979448b6425SPatrick Farrell PetscFunctionBegin; 980448b6425SPatrick Farrell *n = nasm->n; 981448b6425SPatrick Farrell PetscFunctionReturn(0); 982448b6425SPatrick Farrell } 983448b6425SPatrick Farrell 984f10b3e88SPatrick Farrell /*@ 985f10b3e88SPatrick Farrell SNESNASMSetWeight - Sets weight to use when adding overlapping updates 986f10b3e88SPatrick Farrell 987f10b3e88SPatrick Farrell Collective 988f10b3e88SPatrick Farrell 989f10b3e88SPatrick Farrell Input Parameters: 990f10b3e88SPatrick Farrell + snes - the SNES context 991f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 992f10b3e88SPatrick Farrell 993f10b3e88SPatrick Farrell Level: intermediate 994f10b3e88SPatrick Farrell 995f10b3e88SPatrick Farrell .seealso: SNESNASM 996f10b3e88SPatrick Farrell @*/ 997f10b3e88SPatrick Farrell PetscErrorCode SNESNASMSetWeight(SNES snes,Vec weight) 998f10b3e88SPatrick Farrell { 999f10b3e88SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 1000f10b3e88SPatrick Farrell PetscErrorCode ierr; 1001f10b3e88SPatrick Farrell 1002f10b3e88SPatrick Farrell PetscFunctionBegin; 1003f10b3e88SPatrick Farrell 1004f10b3e88SPatrick Farrell ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr); 1005f10b3e88SPatrick Farrell nasm->weight_set = PETSC_TRUE; 1006f10b3e88SPatrick Farrell nasm->weight = weight; 1007f10b3e88SPatrick Farrell ierr = PetscObjectReference((PetscObject)nasm->weight);CHKERRQ(ierr); 1008f10b3e88SPatrick Farrell 1009f10b3e88SPatrick Farrell PetscFunctionReturn(0); 1010f10b3e88SPatrick Farrell } 1011