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 */ 20f10b3e88SPatrick Farrell PetscBool weight_set; /* use a weight in the overlap updates */ 21b20c023fSPeter Brune 22b20c023fSPeter Brune /* logging events */ 23b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 24b20c023fSPeter Brune PetscLogEvent eventsubsolve; 25602bec5dSPeter Brune 26602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 27602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 28eaedb033SPeter Brune } SNES_NASM; 29eaedb033SPeter Brune 309e5d0892SLisandro Dalcin const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",NULL}; 31602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 32b20c023fSPeter Brune 3340244768SBarry Smith static 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); } 48f10b3e88SPatrick Farrell if (nasm->oscatter_copy) { ierr = VecScatterDestroy(&nasm->oscatter_copy[i]);CHKERRQ(ierr); } 49111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 50111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 51eaedb033SPeter Brune } 52111ade9eSPeter Brune 530677ede6SBarry Smith ierr = PetscFree(nasm->x);CHKERRQ(ierr); 540677ede6SBarry Smith ierr = PetscFree(nasm->xl);CHKERRQ(ierr); 550677ede6SBarry Smith ierr = PetscFree(nasm->y);CHKERRQ(ierr); 560677ede6SBarry Smith ierr = PetscFree(nasm->b);CHKERRQ(ierr); 57111ade9eSPeter Brune 58602bec5dSPeter Brune if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 59602bec5dSPeter Brune 600677ede6SBarry Smith ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr); 610677ede6SBarry Smith ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr); 62f10b3e88SPatrick Farrell ierr = PetscFree(nasm->oscatter_copy);CHKERRQ(ierr); 630677ede6SBarry Smith ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr); 640677ede6SBarry Smith ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr); 65b20c023fSPeter Brune 66f10b3e88SPatrick Farrell if (nasm->weight_set) { 67f10b3e88SPatrick Farrell ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr); 68f10b3e88SPatrick Farrell } 69f10b3e88SPatrick Farrell 70b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 71b20c023fSPeter Brune nasm->eventsubsolve = 0; 72eaedb033SPeter Brune PetscFunctionReturn(0); 73eaedb033SPeter Brune } 74eaedb033SPeter Brune 7540244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes) 76eaedb033SPeter Brune { 77eaedb033SPeter Brune PetscErrorCode ierr; 786e111a19SKarl Rupp 79eaedb033SPeter Brune PetscFunctionBegin; 80eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 8122d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 82eaedb033SPeter Brune PetscFunctionReturn(0); 83eaedb033SPeter Brune } 84eaedb033SPeter Brune 8540244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 860adebc6cSBarry Smith { 87111ade9eSPeter Brune PetscErrorCode ierr; 88111ade9eSPeter Brune Vec bcs = (Vec)ctx; 896e111a19SKarl Rupp 90111ade9eSPeter Brune PetscFunctionBegin; 91111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 92111ade9eSPeter Brune PetscFunctionReturn(0); 93111ade9eSPeter Brune } 94111ade9eSPeter Brune 9540244768SBarry Smith static 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); 116f10b3e88SPatrick Farrell ierr = PetscMalloc1(nasm->n, &nasm->oscatter_copy);CHKERRQ(ierr); 117f10b3e88SPatrick Farrell for (i=0; i<nasm->n; i++) { 118f10b3e88SPatrick Farrell ierr = VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i]);CHKERRQ(ierr); 119f10b3e88SPatrick Farrell } 120eaedb033SPeter Brune 121eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 122785e854fSJed Brown ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr); 123111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 124cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 1257a2f0ea0SPatrick Farrell ierr = PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1);CHKERRQ(ierr); 126cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 127cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 128cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 129ffc4695bSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRMPI(ierr); 130ed3c11a9SPeter Brune if (size == 1) { 131ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 132ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 133ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 134ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 135ed3c11a9SPeter Brune } 136cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 137111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 138111ade9eSPeter Brune } 139111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 140ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 141ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 142111ade9eSPeter Brune /* allocate the global vectors */ 143e245e0aaSPeter Brune if (!nasm->x) { 1441795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr); 145e245e0aaSPeter Brune } 146e245e0aaSPeter Brune if (!nasm->xl) { 1471795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr); 148e245e0aaSPeter Brune } 149e245e0aaSPeter Brune if (!nasm->y) { 1501795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr); 151e245e0aaSPeter Brune } 152e245e0aaSPeter Brune if (!nasm->b) { 1531795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr); 154e245e0aaSPeter Brune } 155111ade9eSPeter Brune 156111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1570298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 15876857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 15976857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 16076857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 16176857b2aSPeter Brune if (!nasm->xl[i]) { 162111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 163111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 1640298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 165111ade9eSPeter Brune } 16661ba4676SBarry Smith } 167602bec5dSPeter Brune if (nasm->finaljacobian) { 168602bec5dSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 169602bec5dSPeter Brune if (nasm->fjtype == 2) { 170602bec5dSPeter Brune ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 171602bec5dSPeter Brune } 172602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 173602bec5dSPeter Brune ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 174602bec5dSPeter Brune } 175602bec5dSPeter Brune } 176eaedb033SPeter Brune PetscFunctionReturn(0); 177eaedb033SPeter Brune } 178eaedb033SPeter Brune 17940244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes) 180eaedb033SPeter Brune { 181eaedb033SPeter Brune PetscErrorCode ierr; 182111ade9eSPeter Brune PCASMType asmtype; 183*83dc3634SPierre Jolivet PetscBool flg,monflg; 184111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1856e111a19SKarl Rupp 186eaedb033SPeter Brune PetscFunctionBegin; 187f3f228e0SPierre Jolivet ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwarz options");CHKERRQ(ierr); 188111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 189e0331734SPeter Brune if (flg) {ierr = SNESNASMSetType(snes,asmtype);CHKERRQ(ierr);} 190b20c023fSPeter Brune flg = PETSC_FALSE; 191b20c023fSPeter Brune monflg = PETSC_TRUE; 1925dfa0f3bSBarry 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); 193610116beSPeter Brune if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 194*83dc3634SPierre Jolivet ierr = PetscOptionsDeprecated("-snes_nasm_sub_view",NULL,"3.15","Use -snes_view ::ascii_info_detail");CHKERRQ(ierr); 195d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 196602bec5dSPeter Brune ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 197a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 198b20c023fSPeter Brune if (flg) { 199b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 200b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 201b20c023fSPeter Brune } 202eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 203eaedb033SPeter Brune PetscFunctionReturn(0); 204eaedb033SPeter Brune } 205eaedb033SPeter Brune 20640244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 207eaedb033SPeter Brune { 208b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 209b20c023fSPeter Brune PetscErrorCode ierr; 210a4f17876SPeter Brune PetscMPIInt rank,size; 211dd2fa690SBarry Smith PetscInt i,N,bsz; 212b20c023fSPeter Brune PetscBool iascii,isstring; 213b20c023fSPeter Brune PetscViewer sviewer; 214ce94432eSBarry Smith MPI_Comm comm; 215*83dc3634SPierre Jolivet PetscViewerFormat format; 216*83dc3634SPierre Jolivet const char *prefix; 217b20c023fSPeter Brune 218b20c023fSPeter Brune PetscFunctionBegin; 219ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 220b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 221b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 222ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 223ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 224b2566f29SBarry Smith ierr = MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 225b20c023fSPeter Brune if (iascii) { 226efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %D\n",N);CHKERRQ(ierr); 227*83dc3634SPierre Jolivet ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 228*83dc3634SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 229a4f17876SPeter Brune if (nasm->subsnes) { 230*83dc3634SPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," Local solver information for first block on rank 0:\n");CHKERRQ(ierr); 231*83dc3634SPierre Jolivet ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 232*83dc3634SPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:"");CHKERRQ(ierr); 233a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2343f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 235a4f17876SPeter Brune if (!rank) { 236a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 237a4f17876SPeter Brune ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr); 238a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 239a4f17876SPeter Brune } 2403f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 241a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 242a4f17876SPeter Brune } 243a4f17876SPeter Brune } else { 244a4f17876SPeter Brune /* print the solver on each block */ 2451575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 246b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 247b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2481575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 249*83dc3634SPierre Jolivet ierr = PetscViewerASCIIPrintf(viewer," Local solver information for each block is in the following SNES objects:\n");CHKERRQ(ierr); 250b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 251a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2523f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 253a4f17876SPeter Brune for (i=0; i<nasm->n; i++) { 254a4f17876SPeter Brune ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr); 255a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 256b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 257a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 258b20c023fSPeter Brune } 2593f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 260a4f17876SPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 261b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 262a4f17876SPeter Brune } 263b20c023fSPeter Brune } else if (isstring) { 264b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 2653f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 266b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 2673f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 268b20c023fSPeter Brune } 269eaedb033SPeter Brune PetscFunctionReturn(0); 270eaedb033SPeter Brune } 271eaedb033SPeter Brune 272e0331734SPeter Brune /*@ 273e0331734SPeter Brune SNESNASMSetType - Set the type of subdomain update used 274e0331734SPeter Brune 275e0331734SPeter Brune Logically Collective on SNES 276e0331734SPeter Brune 277e0331734SPeter Brune Input Parameters: 278e0331734SPeter Brune + SNES - the SNES context 279e0331734SPeter Brune - type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT 280e0331734SPeter Brune 281e0331734SPeter Brune Level: intermediate 282e0331734SPeter Brune 283e0331734SPeter Brune .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType() 284e0331734SPeter Brune @*/ 285e0331734SPeter Brune PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type) 286e0331734SPeter Brune { 287e0331734SPeter Brune PetscErrorCode ierr; 288e0331734SPeter Brune PetscErrorCode (*f)(SNES,PCASMType); 289e0331734SPeter Brune 290e0331734SPeter Brune PetscFunctionBegin; 291e0331734SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f);CHKERRQ(ierr); 292e0331734SPeter Brune if (f) {ierr = (f)(snes,type);CHKERRQ(ierr);} 293e0331734SPeter Brune PetscFunctionReturn(0); 294e0331734SPeter Brune } 295e0331734SPeter Brune 29640244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type) 297e0331734SPeter Brune { 298e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 299e0331734SPeter Brune 300e0331734SPeter Brune PetscFunctionBegin; 301e0331734SPeter 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"); 302e0331734SPeter Brune nasm->type = type; 303e0331734SPeter Brune PetscFunctionReturn(0); 304e0331734SPeter Brune } 305e0331734SPeter Brune 306e0331734SPeter Brune /*@ 307e0331734SPeter Brune SNESNASMGetType - Get the type of subdomain update used 308e0331734SPeter Brune 309e0331734SPeter Brune Logically Collective on SNES 310e0331734SPeter Brune 311e0331734SPeter Brune Input Parameters: 312e0331734SPeter Brune . SNES - the SNES context 313e0331734SPeter Brune 314e0331734SPeter Brune Output Parameters: 315e0331734SPeter Brune . type - the type of update 316e0331734SPeter Brune 317e0331734SPeter Brune Level: intermediate 318e0331734SPeter Brune 319e0331734SPeter Brune .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType() 320e0331734SPeter Brune @*/ 321e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type) 322e0331734SPeter Brune { 323e0331734SPeter Brune PetscErrorCode ierr; 324e0331734SPeter Brune 325e0331734SPeter Brune PetscFunctionBegin; 3262a808120SBarry Smith ierr = PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));CHKERRQ(ierr); 327e0331734SPeter Brune PetscFunctionReturn(0); 328e0331734SPeter Brune } 329e0331734SPeter Brune 33040244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type) 331e0331734SPeter Brune { 332e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 333e0331734SPeter Brune 334e0331734SPeter Brune PetscFunctionBegin; 335e0331734SPeter Brune *type = nasm->type; 336e0331734SPeter Brune PetscFunctionReturn(0); 337e0331734SPeter Brune } 338e0331734SPeter Brune 33976857b2aSPeter Brune /*@ 34076857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 34176857b2aSPeter Brune 34276857b2aSPeter Brune Not Collective 34376857b2aSPeter Brune 34476857b2aSPeter Brune Input Parameters: 34576857b2aSPeter Brune + SNES - the SNES context 34676857b2aSPeter Brune . n - the number of local subdomains 34776857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 34876857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 34976857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 35076857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 35176857b2aSPeter Brune 35276857b2aSPeter Brune Level: intermediate 35376857b2aSPeter Brune 35476857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 35576857b2aSPeter Brune @*/ 356a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 357a6dfd86eSKarl Rupp { 358eaedb033SPeter Brune PetscErrorCode ierr; 359111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 3606e111a19SKarl Rupp 361eaedb033SPeter Brune PetscFunctionBegin; 3620005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 3634b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 364eaedb033SPeter Brune PetscFunctionReturn(0); 365eaedb033SPeter Brune } 366eaedb033SPeter Brune 36740244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 368a6dfd86eSKarl Rupp { 369eaedb033SPeter Brune PetscInt i; 370eaedb033SPeter Brune PetscErrorCode ierr; 371eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3726e111a19SKarl Rupp 373eaedb033SPeter Brune PetscFunctionBegin; 374ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 375eaedb033SPeter Brune 376111ade9eSPeter Brune /* tear down the previously set things */ 377111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 378111ade9eSPeter Brune 379eaedb033SPeter Brune nasm->n = n; 380111ade9eSPeter Brune if (oscatter) { 381111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 382eaedb033SPeter Brune } 383111ade9eSPeter Brune if (iscatter) { 384111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 385eaedb033SPeter Brune } 386111ade9eSPeter Brune if (gscatter) { 387111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 388111ade9eSPeter Brune } 389111ade9eSPeter Brune if (oscatter) { 390785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 391f10b3e88SPatrick Farrell ierr = PetscMalloc1(n,&nasm->oscatter_copy);CHKERRQ(ierr); 392eaedb033SPeter Brune for (i=0; i<n; i++) { 393111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 394f10b3e88SPatrick Farrell ierr = VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i]);CHKERRQ(ierr); 395eaedb033SPeter Brune } 396111ade9eSPeter Brune } 397111ade9eSPeter Brune if (iscatter) { 398785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 399eaedb033SPeter Brune for (i=0; i<n; i++) { 400111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 401eaedb033SPeter Brune } 402eaedb033SPeter Brune } 403111ade9eSPeter Brune if (gscatter) { 404785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 405eaedb033SPeter Brune for (i=0; i<n; i++) { 406111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 407eaedb033SPeter Brune } 408eaedb033SPeter Brune } 409111ade9eSPeter Brune 410eaedb033SPeter Brune if (subsnes) { 411785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 412eaedb033SPeter Brune for (i=0; i<n; i++) { 413eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 414eaedb033SPeter Brune } 415eaedb033SPeter Brune } 416eaedb033SPeter Brune PetscFunctionReturn(0); 417eaedb033SPeter Brune } 418eaedb033SPeter Brune 41976857b2aSPeter Brune /*@ 42076857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 42176857b2aSPeter Brune 42276857b2aSPeter Brune Not Collective 42376857b2aSPeter Brune 42476857b2aSPeter Brune Input Parameters: 42576857b2aSPeter Brune . SNES - the SNES context 42676857b2aSPeter Brune 42776857b2aSPeter Brune Output Parameters: 42876857b2aSPeter Brune + n - the number of local subdomains 42976857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 43076857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 43176857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 43276857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 43376857b2aSPeter Brune 43476857b2aSPeter Brune Level: intermediate 43576857b2aSPeter Brune 43676857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 43776857b2aSPeter Brune @*/ 43876857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 43976857b2aSPeter Brune { 44076857b2aSPeter Brune PetscErrorCode ierr; 44176857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 44276857b2aSPeter Brune 44376857b2aSPeter Brune PetscFunctionBegin; 4440005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 4454b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 44676857b2aSPeter Brune PetscFunctionReturn(0); 44776857b2aSPeter Brune } 44876857b2aSPeter Brune 44940244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 45076857b2aSPeter Brune { 45176857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 45276857b2aSPeter Brune 45376857b2aSPeter Brune PetscFunctionBegin; 45476857b2aSPeter Brune if (n) *n = nasm->n; 45576857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 45676857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 45776857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 458*83dc3634SPierre Jolivet if (subsnes) *subsnes = nasm->subsnes; 45976857b2aSPeter Brune PetscFunctionReturn(0); 46076857b2aSPeter Brune } 46176857b2aSPeter Brune 46276857b2aSPeter Brune /*@ 46376857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 46476857b2aSPeter Brune 46576857b2aSPeter Brune Not Collective 46676857b2aSPeter Brune 46776857b2aSPeter Brune Input Parameters: 46876857b2aSPeter Brune . SNES - the SNES context 46976857b2aSPeter Brune 47076857b2aSPeter Brune Output Parameters: 47176857b2aSPeter Brune + n - the number of local subdomains 47276857b2aSPeter Brune . x - The subdomain solution vector 47376857b2aSPeter Brune . y - The subdomain step vector 47476857b2aSPeter Brune . b - The subdomain RHS vector 47576857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 47676857b2aSPeter Brune 47776857b2aSPeter Brune Level: developer 47876857b2aSPeter Brune 47976857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 48076857b2aSPeter Brune @*/ 48176857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 48276857b2aSPeter Brune { 48376857b2aSPeter Brune PetscErrorCode ierr; 48476857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 48576857b2aSPeter Brune 48676857b2aSPeter Brune PetscFunctionBegin; 4870005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 4884b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 48976857b2aSPeter Brune PetscFunctionReturn(0); 49076857b2aSPeter Brune } 49176857b2aSPeter Brune 49240244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 49376857b2aSPeter Brune { 49476857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 49576857b2aSPeter Brune 49676857b2aSPeter Brune PetscFunctionBegin; 49776857b2aSPeter Brune if (n) *n = nasm->n; 49876857b2aSPeter Brune if (x) *x = nasm->x; 49976857b2aSPeter Brune if (y) *y = nasm->y; 50076857b2aSPeter Brune if (b) *b = nasm->b; 50176857b2aSPeter Brune if (xl) *xl = nasm->xl; 50276857b2aSPeter Brune PetscFunctionReturn(0); 50376857b2aSPeter Brune } 50476857b2aSPeter Brune 505d728fb7dSPeter Brune /*@ 5068bf45196SRichard Tran Mills SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence 507d728fb7dSPeter Brune 508d728fb7dSPeter Brune Collective on SNES 509d728fb7dSPeter Brune 510d728fb7dSPeter Brune Input Parameters: 511d728fb7dSPeter Brune + SNES - the SNES context 5128bf45196SRichard Tran Mills - flg - indication of whether to compute the Jacobians or not 513d728fb7dSPeter Brune 514d728fb7dSPeter Brune Level: developer 515d728fb7dSPeter Brune 51695452b02SPatrick Sanan Notes: 5178bf45196SRichard Tran Mills This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global Jacobian 518d728fb7dSPeter Brune is needed at each linear iteration. 519d728fb7dSPeter Brune 520d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 521d728fb7dSPeter Brune @*/ 522d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 523d728fb7dSPeter Brune { 524d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 525d728fb7dSPeter Brune PetscErrorCode ierr; 526d728fb7dSPeter Brune 527d728fb7dSPeter Brune PetscFunctionBegin; 5280005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 5294b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 530d728fb7dSPeter Brune PetscFunctionReturn(0); 531d728fb7dSPeter Brune } 532d728fb7dSPeter Brune 53340244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 534d728fb7dSPeter Brune { 535d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 536d728fb7dSPeter Brune 537d728fb7dSPeter Brune PetscFunctionBegin; 538d728fb7dSPeter Brune nasm->finaljacobian = flg; 539d728fb7dSPeter Brune PetscFunctionReturn(0); 540d728fb7dSPeter Brune } 54176857b2aSPeter Brune 542610116beSPeter Brune /*@ 543610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 544610116beSPeter Brune 545610116beSPeter Brune Logically collective on SNES 546610116beSPeter Brune 547610116beSPeter Brune Input Parameters: 548610116beSPeter Brune + SNES - the SNES context 549610116beSPeter Brune - dmp - damping 550610116beSPeter Brune 551610116beSPeter Brune Level: intermediate 552610116beSPeter Brune 55395452b02SPatrick Sanan Notes: 55495452b02SPatrick Sanan The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5555dfa0f3bSBarry Smith 556610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 557610116beSPeter Brune @*/ 558610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 559610116beSPeter Brune { 560610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 561610116beSPeter Brune PetscErrorCode ierr; 562610116beSPeter Brune 563610116beSPeter Brune PetscFunctionBegin; 564610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 565e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 566610116beSPeter Brune PetscFunctionReturn(0); 567610116beSPeter Brune } 568610116beSPeter Brune 56940244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 570610116beSPeter Brune { 571610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 572610116beSPeter Brune 573610116beSPeter Brune PetscFunctionBegin; 574610116beSPeter Brune nasm->damping = dmp; 575610116beSPeter Brune PetscFunctionReturn(0); 576610116beSPeter Brune } 577610116beSPeter Brune 578610116beSPeter Brune /*@ 579610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 580610116beSPeter Brune 581610116beSPeter Brune Not Collective 582610116beSPeter Brune 583610116beSPeter Brune Input Parameters: 584610116beSPeter Brune + SNES - the SNES context 585610116beSPeter Brune - dmp - damping 586610116beSPeter Brune 587610116beSPeter Brune Level: intermediate 588610116beSPeter Brune 589610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping() 590610116beSPeter Brune @*/ 591610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 592610116beSPeter Brune { 593610116beSPeter Brune PetscErrorCode ierr; 594610116beSPeter Brune 595610116beSPeter Brune PetscFunctionBegin; 5964a2f8832SBarry Smith ierr = PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));CHKERRQ(ierr); 597610116beSPeter Brune PetscFunctionReturn(0); 598610116beSPeter Brune } 599610116beSPeter Brune 60040244768SBarry Smith static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 601610116beSPeter Brune { 602610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 603610116beSPeter Brune 604610116beSPeter Brune PetscFunctionBegin; 605610116beSPeter Brune *dmp = nasm->damping; 606610116beSPeter Brune PetscFunctionReturn(0); 607610116beSPeter Brune } 608610116beSPeter Brune 609610116beSPeter Brune 61014eb1c5cSMatthew G. Knepley /* 61114eb1c5cSMatthew G. Knepley Input Parameters: 61214eb1c5cSMatthew G. Knepley + snes - The solver 61314eb1c5cSMatthew G. Knepley . B - The RHS vector 61414eb1c5cSMatthew G. Knepley - X - The initial guess 61514eb1c5cSMatthew G. Knepley 61614eb1c5cSMatthew G. Knepley Output Parameters: 61714eb1c5cSMatthew G. Knepley . Y - The solution update 61814eb1c5cSMatthew G. Knepley 61914eb1c5cSMatthew G. Knepley TODO: All scatters should be packed into one 62014eb1c5cSMatthew G. Knepley */ 6210adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 6220adebc6cSBarry Smith { 623eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 624258e1594SPeter Brune SNES subsnes; 625eaedb033SPeter Brune PetscInt i; 626610116beSPeter Brune PetscReal dmp; 627eaedb033SPeter Brune PetscErrorCode ierr; 628f10b3e88SPatrick Farrell Vec Xl,Bl,Yl,Xlloc; 629f10b3e88SPatrick Farrell VecScatter iscat,oscat,gscat,oscat_copy; 630111ade9eSPeter Brune DM dm,subdm; 631e0331734SPeter Brune PCASMType type; 6320adebc6cSBarry Smith 633eaedb033SPeter Brune PetscFunctionBegin; 634e0331734SPeter Brune ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr); 635eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 636111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 637b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 638eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 639f10b3e88SPatrick Farrell /* scatter the solution to the global solution and the local solution */ 640f10b3e88SPatrick Farrell Xl = nasm->x[i]; 64170c78f05SPeter Brune Xlloc = nasm->xl[i]; 64270c78f05SPeter Brune oscat = nasm->oscatter[i]; 643f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 644f10b3e88SPatrick Farrell gscat = nasm->gscatter[i]; 645f10b3e88SPatrick Farrell ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 64670c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 64770c78f05SPeter Brune if (B) { 64870c78f05SPeter Brune /* scatter the RHS to the local RHS */ 64970c78f05SPeter Brune Bl = nasm->b[i]; 650f10b3e88SPatrick Farrell ierr = VecScatterBegin(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 65170c78f05SPeter Brune } 65270c78f05SPeter Brune } 653b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 654b20c023fSPeter Brune 655b20c023fSPeter Brune 656b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 65770c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 65870c78f05SPeter Brune Xl = nasm->x[i]; 65970c78f05SPeter Brune Xlloc = nasm->xl[i]; 66070c78f05SPeter Brune Yl = nasm->y[i]; 661258e1594SPeter Brune subsnes = nasm->subsnes[i]; 662258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 663111ade9eSPeter Brune iscat = nasm->iscatter[i]; 664111ade9eSPeter Brune oscat = nasm->oscatter[i]; 665f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 666111ade9eSPeter Brune gscat = nasm->gscatter[i]; 667f10b3e88SPatrick Farrell ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 668ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 66924b7f281SPeter Brune if (B) { 67024b7f281SPeter Brune Bl = nasm->b[i]; 671f10b3e88SPatrick Farrell ierr = VecScatterEnd(oscat_copy,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 672ed3c11a9SPeter Brune } else Bl = NULL; 673f10b3e88SPatrick Farrell 674ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 675111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 676ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 677ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 678f10b3e88SPatrick Farrell ierr = VecScale(Yl, nasm->damping);CHKERRQ(ierr); 679e0331734SPeter Brune if (type == PC_ASM_BASIC) { 680111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6811ed9ada7SJunchao Zhang ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 682e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 683111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6841ed9ada7SJunchao Zhang ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 685ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 686eaedb033SPeter Brune } 687ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 688ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 689f10b3e88SPatrick Farrell if (nasm->weight_set) { 690f10b3e88SPatrick Farrell ierr = VecPointwiseMult(Y,Y,nasm->weight);CHKERRQ(ierr); 691f10b3e88SPatrick Farrell } 692b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 6935dfa0f3bSBarry Smith ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 694610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 695eaedb033SPeter Brune PetscFunctionReturn(0); 696eaedb033SPeter Brune } 697eaedb033SPeter Brune 69840244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 699d728fb7dSPeter Brune { 700602bec5dSPeter Brune Vec X = Xfinal; 701d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 702d728fb7dSPeter Brune SNES subsnes; 703ca44f815SPeter Brune PetscInt i,lag = 1; 704d728fb7dSPeter Brune PetscErrorCode ierr; 705e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 706d728fb7dSPeter Brune VecScatter oscat,gscat; 707d728fb7dSPeter Brune DM dm,subdm; 708d1e9a80fSBarry Smith 709d728fb7dSPeter Brune PetscFunctionBegin; 710602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 711e59f0a30SPeter Brune F = snes->vec_func; 712365a6726SPeter Brune if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 713d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 714d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 715d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 716602bec5dSPeter Brune if (nasm->fjtype != 1) { 717d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 718d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 719d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 720d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 721602bec5dSPeter Brune } 722d728fb7dSPeter Brune } 723d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 724d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 725e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 726d728fb7dSPeter Brune Xl = nasm->x[i]; 727d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 728d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 729d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 730d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 731602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 732ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 733d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 734602bec5dSPeter Brune if (nasm->fjtype != 1) { 735d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 736d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 737602bec5dSPeter Brune } 738ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 739ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 740602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 741d1e9a80fSBarry Smith ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr); 742ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 743d728fb7dSPeter Brune } 744d728fb7dSPeter Brune PetscFunctionReturn(0); 745d728fb7dSPeter Brune } 746d728fb7dSPeter Brune 74740244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes) 748eaedb033SPeter Brune { 749eaedb033SPeter Brune Vec F; 750eaedb033SPeter Brune Vec X; 751eaedb033SPeter Brune Vec B; 752111ade9eSPeter Brune Vec Y; 753eaedb033SPeter Brune PetscInt i; 754ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 755eaedb033SPeter Brune PetscErrorCode ierr; 756365a6726SPeter Brune SNESNormSchedule normschedule; 757d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 758eaedb033SPeter Brune 759eaedb033SPeter Brune PetscFunctionBegin; 760c579b300SPatrick Farrell 7616c4ed002SBarry 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); 762c579b300SPatrick Farrell 763fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 764eaedb033SPeter Brune X = snes->vec_sol; 765111ade9eSPeter Brune Y = snes->vec_sol_update; 766eaedb033SPeter Brune F = snes->vec_func; 767eaedb033SPeter Brune B = snes->vec_rhs; 768eaedb033SPeter Brune 769e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 770eaedb033SPeter Brune snes->iter = 0; 771eaedb033SPeter Brune snes->norm = 0.; 772e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 773eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 774365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 775365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 776eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 777eaedb033SPeter Brune if (!snes->vec_func_init_set) { 778eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 7791aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 780eaedb033SPeter Brune 781eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 782422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 783e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 784eaedb033SPeter Brune snes->iter = 0; 785eaedb033SPeter Brune snes->norm = fnorm; 786e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 787a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 788eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 789eaedb033SPeter Brune 790eaedb033SPeter Brune /* test convergence */ 791eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 792eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 793eaedb033SPeter Brune } else { 794e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 795a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 796eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 797eaedb033SPeter Brune } 798eaedb033SPeter Brune 799eaedb033SPeter Brune /* Call general purpose update function */ 800eaedb033SPeter Brune if (snes->ops->update) { 801eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 802eaedb033SPeter Brune } 803602bec5dSPeter Brune /* copy the initial solution over for later */ 804602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 805eaedb033SPeter Brune 806eaedb033SPeter Brune for (i=0; i < snes->max_its; i++) { 807111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 808365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 809eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 810eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 811422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 812eaedb033SPeter Brune } 813eaedb033SPeter Brune /* Monitor convergence */ 814e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 815eaedb033SPeter Brune snes->iter = i+1; 816eaedb033SPeter Brune snes->norm = fnorm; 817e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 818a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 819eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 820eaedb033SPeter Brune /* Test for convergence */ 821365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 822d728fb7dSPeter Brune if (snes->reason) break; 823eaedb033SPeter Brune /* Call general purpose update function */ 824e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 825eaedb033SPeter Brune } 82607b62357SFande Kong if (nasm->finaljacobian) { 82707b62357SFande Kong ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr); 82807b62357SFande Kong SNESCheckJacobianDomainerror(snes); 82907b62357SFande Kong } 830365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 831eaedb033SPeter Brune if (i == snes->max_its) { 832eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 833eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 834eaedb033SPeter Brune } 8351aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 836eaedb033SPeter Brune PetscFunctionReturn(0); 837eaedb033SPeter Brune } 838eaedb033SPeter Brune 839eaedb033SPeter Brune /*MC 840f3f228e0SPierre Jolivet SNESNASM - Nonlinear Additive Schwarz 841eaedb033SPeter Brune 84270c78603SPeter Brune Options Database: 84370c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 84470c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 8455dfa0f3bSBarry Smith . -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 84670c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 847150d49b7SHong Zhang . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at 84870c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 84970c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 85070c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 85170c78603SPeter Brune 852eaedb033SPeter Brune Level: advanced 853eaedb033SPeter Brune 854951fe5abSBarry 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 8557f851d37SRichard Tran Mills false and SNESView() and -snes_view do not display a KSP object. However, if the flag nasm->finaljacobian is set (for example, if 856951fe5abSBarry Smith NASM is used as a nonlinear preconditioner for KSPASPIN) then SNESSetUpMatrices() is called to generate the Jacobian (needed by KSPASPIN) 857951fe5abSBarry 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 858951fe5abSBarry Smith used by SNESASPIN they share the same Jacobian matrices because SNESSetUp() (called on the outer SNES KSPASPIN) causes the inner SNES 859951fe5abSBarry Smith object (in this case SNESNASM) to inherit the outer Jacobian matrices. 860951fe5abSBarry Smith 8614f02bc6aSBarry Smith References: 86296a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8634f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8644f02bc6aSBarry Smith 8655dfa0f3bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping() 866eaedb033SPeter Brune M*/ 867eaedb033SPeter Brune 8688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 869eaedb033SPeter Brune { 870eaedb033SPeter Brune SNES_NASM *nasm; 871eaedb033SPeter Brune PetscErrorCode ierr; 872eaedb033SPeter Brune 873eaedb033SPeter Brune PetscFunctionBegin; 874b00a9115SJed Brown ierr = PetscNewLog(snes,&nasm);CHKERRQ(ierr); 875eaedb033SPeter Brune snes->data = (void*)nasm; 876eaedb033SPeter Brune 877eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 8789e5d0892SLisandro Dalcin nasm->subsnes = NULL; 8799e5d0892SLisandro Dalcin nasm->x = NULL; 8809e5d0892SLisandro Dalcin nasm->xl = NULL; 8819e5d0892SLisandro Dalcin nasm->y = NULL; 8829e5d0892SLisandro Dalcin nasm->b = NULL; 8839e5d0892SLisandro Dalcin nasm->oscatter = NULL; 8849e5d0892SLisandro Dalcin nasm->oscatter_copy = NULL; 8859e5d0892SLisandro Dalcin nasm->iscatter = NULL; 8869e5d0892SLisandro Dalcin nasm->gscatter = NULL; 887610116beSPeter Brune nasm->damping = 1.; 888111ade9eSPeter Brune 889111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 890d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 891f10b3e88SPatrick Farrell nasm->weight_set = PETSC_FALSE; 892eaedb033SPeter Brune 893eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 894eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 895eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 896eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 897eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 898eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 899eaedb033SPeter Brune 900eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 901efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 902eaedb033SPeter Brune 9034fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 9044fc747eaSLawrence Mitchell 905602bec5dSPeter Brune nasm->fjtype = 0; 906602bec5dSPeter Brune nasm->xinit = NULL; 9070298fd71SBarry Smith nasm->eventrestrictinterp = 0; 9080298fd71SBarry Smith nasm->eventsubsolve = 0; 909b20c023fSPeter Brune 910eaedb033SPeter Brune if (!snes->tolerancesset) { 911eaedb033SPeter Brune snes->max_its = 10000; 912eaedb033SPeter Brune snes->max_funcs = 10000; 913eaedb033SPeter Brune } 914eaedb033SPeter Brune 915e0331734SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr); 916e0331734SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr); 917bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 918bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 91990bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 92090bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 921bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 922bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 923eaedb033SPeter Brune PetscFunctionReturn(0); 924eaedb033SPeter Brune } 92599e0435eSBarry Smith 926448b6425SPatrick Farrell /*@ 927448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 928448b6425SPatrick Farrell 929448b6425SPatrick Farrell Not collective 930448b6425SPatrick Farrell 931448b6425SPatrick Farrell Input Parameters: 932448b6425SPatrick Farrell + snes - the SNES context 933448b6425SPatrick Farrell - i - the number of the subsnes to get 934448b6425SPatrick Farrell 935448b6425SPatrick Farrell Output Parameters: 936448b6425SPatrick Farrell . subsnes - the subsolver context 937448b6425SPatrick Farrell 938448b6425SPatrick Farrell Level: intermediate 939448b6425SPatrick Farrell 940448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetNumber() 941448b6425SPatrick Farrell @*/ 942448b6425SPatrick Farrell PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes) 943448b6425SPatrick Farrell { 944448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 945448b6425SPatrick Farrell 946448b6425SPatrick Farrell PetscFunctionBegin; 947448b6425SPatrick Farrell if (i < 0 || i >= nasm->n) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver"); 948448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 949448b6425SPatrick Farrell PetscFunctionReturn(0); 950448b6425SPatrick Farrell } 951448b6425SPatrick Farrell 952448b6425SPatrick Farrell /*@ 953448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 954448b6425SPatrick Farrell 955448b6425SPatrick Farrell Not collective 956448b6425SPatrick Farrell 957448b6425SPatrick Farrell Input Parameters: 958448b6425SPatrick Farrell . snes - the SNES context 959448b6425SPatrick Farrell 960448b6425SPatrick Farrell Output Parameters: 961448b6425SPatrick Farrell . n - the number of subsolvers 962448b6425SPatrick Farrell 963448b6425SPatrick Farrell Level: intermediate 964448b6425SPatrick Farrell 965448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetSNES() 966448b6425SPatrick Farrell @*/ 967448b6425SPatrick Farrell PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n) 968448b6425SPatrick Farrell { 969448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 970448b6425SPatrick Farrell 971448b6425SPatrick Farrell PetscFunctionBegin; 972448b6425SPatrick Farrell *n = nasm->n; 973448b6425SPatrick Farrell PetscFunctionReturn(0); 974448b6425SPatrick Farrell } 975448b6425SPatrick Farrell 976f10b3e88SPatrick Farrell /*@ 977f10b3e88SPatrick Farrell SNESNASMSetWeight - Sets weight to use when adding overlapping updates 978f10b3e88SPatrick Farrell 979f10b3e88SPatrick Farrell Collective 980f10b3e88SPatrick Farrell 981f10b3e88SPatrick Farrell Input Parameters: 982f10b3e88SPatrick Farrell + snes - the SNES context 983f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 984f10b3e88SPatrick Farrell 985f10b3e88SPatrick Farrell Level: intermediate 986f10b3e88SPatrick Farrell 987f10b3e88SPatrick Farrell .seealso: SNESNASM 988f10b3e88SPatrick Farrell @*/ 989f10b3e88SPatrick Farrell PetscErrorCode SNESNASMSetWeight(SNES snes,Vec weight) 990f10b3e88SPatrick Farrell { 991f10b3e88SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 992f10b3e88SPatrick Farrell PetscErrorCode ierr; 993f10b3e88SPatrick Farrell 994f10b3e88SPatrick Farrell PetscFunctionBegin; 995f10b3e88SPatrick Farrell 996f10b3e88SPatrick Farrell ierr = VecDestroy(&nasm->weight);CHKERRQ(ierr); 997f10b3e88SPatrick Farrell nasm->weight_set = PETSC_TRUE; 998f10b3e88SPatrick Farrell nasm->weight = weight; 999f10b3e88SPatrick Farrell ierr = PetscObjectReference((PetscObject)nasm->weight);CHKERRQ(ierr); 1000f10b3e88SPatrick Farrell 1001f10b3e88SPatrick Farrell PetscFunctionReturn(0); 1002f10b3e88SPatrick Farrell } 1003