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 */ 11111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 12111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 13111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 14111ade9eSPeter Brune PCASMType type; /* ASM type */ 15111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 16d728fb7dSPeter Brune PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 17610116beSPeter Brune PetscReal damping; /* damping parameter for updates from the blocks */ 18a4f17876SPeter Brune PetscBool same_local_solves; /* flag to determine if the solvers have been individually modified */ 19b20c023fSPeter Brune 20b20c023fSPeter Brune /* logging events */ 21b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 22b20c023fSPeter Brune PetscLogEvent eventsubsolve; 23602bec5dSPeter Brune 24602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 25602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 26eaedb033SPeter Brune } SNES_NASM; 27eaedb033SPeter Brune 28b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 29602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 30b20c023fSPeter Brune 3140244768SBarry Smith static PetscErrorCode SNESReset_NASM(SNES snes) 32eaedb033SPeter Brune { 33eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 34eaedb033SPeter Brune PetscErrorCode ierr; 35eaedb033SPeter Brune PetscInt i; 366e111a19SKarl Rupp 37eaedb033SPeter Brune PetscFunctionBegin; 38eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 39111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 40f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 41111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 42bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 43eaedb033SPeter Brune 44bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 45111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 46111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 47111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 48eaedb033SPeter Brune } 49111ade9eSPeter Brune 500677ede6SBarry Smith ierr = PetscFree(nasm->x);CHKERRQ(ierr); 510677ede6SBarry Smith ierr = PetscFree(nasm->xl);CHKERRQ(ierr); 520677ede6SBarry Smith ierr = PetscFree(nasm->y);CHKERRQ(ierr); 530677ede6SBarry Smith ierr = PetscFree(nasm->b);CHKERRQ(ierr); 54111ade9eSPeter Brune 55602bec5dSPeter Brune if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 56602bec5dSPeter Brune 570677ede6SBarry Smith ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr); 580677ede6SBarry Smith ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr); 590677ede6SBarry Smith ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr); 600677ede6SBarry Smith ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr); 61b20c023fSPeter Brune 62b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 63b20c023fSPeter Brune nasm->eventsubsolve = 0; 64eaedb033SPeter Brune PetscFunctionReturn(0); 65eaedb033SPeter Brune } 66eaedb033SPeter Brune 6740244768SBarry Smith static PetscErrorCode SNESDestroy_NASM(SNES snes) 68eaedb033SPeter Brune { 69eaedb033SPeter Brune PetscErrorCode ierr; 706e111a19SKarl Rupp 71eaedb033SPeter Brune PetscFunctionBegin; 72eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 7322d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 74eaedb033SPeter Brune PetscFunctionReturn(0); 75eaedb033SPeter Brune } 76eaedb033SPeter Brune 7740244768SBarry Smith static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 780adebc6cSBarry Smith { 79111ade9eSPeter Brune PetscErrorCode ierr; 80111ade9eSPeter Brune Vec bcs = (Vec)ctx; 816e111a19SKarl Rupp 82111ade9eSPeter Brune PetscFunctionBegin; 83111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 84111ade9eSPeter Brune PetscFunctionReturn(0); 85111ade9eSPeter Brune } 86111ade9eSPeter Brune 8740244768SBarry Smith static PetscErrorCode SNESSetUp_NASM(SNES snes) 88eaedb033SPeter Brune { 89eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 90eaedb033SPeter Brune PetscErrorCode ierr; 9176857b2aSPeter Brune DM dm,subdm; 92111ade9eSPeter Brune DM *subdms; 93111ade9eSPeter Brune PetscInt i; 94eaedb033SPeter Brune const char *optionsprefix; 95111ade9eSPeter Brune Vec F; 96ed3c11a9SPeter Brune PetscMPIInt size; 97ed3c11a9SPeter Brune KSP ksp; 98ed3c11a9SPeter Brune PC pc; 99eaedb033SPeter Brune 100eaedb033SPeter Brune PetscFunctionBegin; 101eaedb033SPeter Brune if (!nasm->subsnes) { 102eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1030a696f66SPeter Brune if (dm) { 104eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1050298fd71SBarry Smith ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 106ce94432eSBarry Smith if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 107111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 108eaedb033SPeter Brune 109eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 110785e854fSJed Brown ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr); 111111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 112cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 1137a2f0ea0SPatrick Farrell ierr = PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1);CHKERRQ(ierr); 114cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 115cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 116cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 117ed3c11a9SPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 118ed3c11a9SPeter Brune if (size == 1) { 119ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 120ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 121ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 122ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 123ed3c11a9SPeter Brune } 124cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 125111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 126111ade9eSPeter Brune } 127111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 128ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 129ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 130111ade9eSPeter Brune /* allocate the global vectors */ 131e245e0aaSPeter Brune if (!nasm->x) { 1321795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr); 133e245e0aaSPeter Brune } 134e245e0aaSPeter Brune if (!nasm->xl) { 1351795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr); 136e245e0aaSPeter Brune } 137e245e0aaSPeter Brune if (!nasm->y) { 1381795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr); 139e245e0aaSPeter Brune } 140e245e0aaSPeter Brune if (!nasm->b) { 1411795a4d1SJed Brown ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr); 142e245e0aaSPeter Brune } 143111ade9eSPeter Brune 144111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1450298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 14676857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 14776857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 14876857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 14976857b2aSPeter Brune if (!nasm->xl[i]) { 150111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 151111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 1520298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 153111ade9eSPeter Brune } 15461ba4676SBarry Smith } 155602bec5dSPeter Brune if (nasm->finaljacobian) { 156602bec5dSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 157602bec5dSPeter Brune if (nasm->fjtype == 2) { 158602bec5dSPeter Brune ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 159602bec5dSPeter Brune } 160602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 161602bec5dSPeter Brune ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 162602bec5dSPeter Brune } 163602bec5dSPeter Brune } 164eaedb033SPeter Brune PetscFunctionReturn(0); 165eaedb033SPeter Brune } 166eaedb033SPeter Brune 16740244768SBarry Smith static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes) 168eaedb033SPeter Brune { 169eaedb033SPeter Brune PetscErrorCode ierr; 170111ade9eSPeter Brune PCASMType asmtype; 171a4f17876SPeter Brune PetscBool flg,monflg,subviewflg; 172111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1736e111a19SKarl Rupp 174eaedb033SPeter Brune PetscFunctionBegin; 175e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwartz options");CHKERRQ(ierr); 176111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 177e0331734SPeter Brune if (flg) {ierr = SNESNASMSetType(snes,asmtype);CHKERRQ(ierr);} 178b20c023fSPeter Brune flg = PETSC_FALSE; 179b20c023fSPeter Brune monflg = PETSC_TRUE; 1805dfa0f3bSBarry 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); 181610116beSPeter Brune if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 182a4f17876SPeter Brune subviewflg = PETSC_FALSE; 183a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr); 184a4f17876SPeter Brune if (flg) { 185a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 186a4f17876SPeter Brune if (!subviewflg) { 187a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 188a4f17876SPeter Brune } 189a4f17876SPeter Brune } 190d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 191602bec5dSPeter Brune ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 192a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 193b20c023fSPeter Brune if (flg) { 194b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 195b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 196b20c023fSPeter Brune } 197eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 198eaedb033SPeter Brune PetscFunctionReturn(0); 199eaedb033SPeter Brune } 200eaedb033SPeter Brune 20140244768SBarry Smith static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 202eaedb033SPeter Brune { 203b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 204b20c023fSPeter Brune PetscErrorCode ierr; 205a4f17876SPeter Brune PetscMPIInt rank,size; 206dd2fa690SBarry Smith PetscInt i,N,bsz; 207b20c023fSPeter Brune PetscBool iascii,isstring; 208b20c023fSPeter Brune PetscViewer sviewer; 209ce94432eSBarry Smith MPI_Comm comm; 210b20c023fSPeter Brune 211b20c023fSPeter Brune PetscFunctionBegin; 212ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 213b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 214b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 215b20c023fSPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 216a4f17876SPeter Brune ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 217b2566f29SBarry Smith ierr = MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 218b20c023fSPeter Brune if (iascii) { 219*efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %D\n",N);CHKERRQ(ierr); 220a4f17876SPeter Brune if (nasm->same_local_solves) { 221a4f17876SPeter Brune if (nasm->subsnes) { 222a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve is the same for all blocks:\n");CHKERRQ(ierr); 223a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2243f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 225a4f17876SPeter Brune if (!rank) { 226a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 227a4f17876SPeter Brune ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr); 228a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 229a4f17876SPeter Brune } 2303f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 231a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 232a4f17876SPeter Brune } 233a4f17876SPeter Brune } else { 234a4f17876SPeter Brune /* print the solver on each block */ 2351575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 236b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 237b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2381575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 239a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr); 240b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 241a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2423f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 243a4f17876SPeter Brune for (i=0; i<nasm->n; i++) { 244a4f17876SPeter Brune ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr); 245a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 246b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 247a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 248b20c023fSPeter Brune } 2493f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 250a4f17876SPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 251b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 252a4f17876SPeter Brune } 253b20c023fSPeter Brune } else if (isstring) { 254b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 2553f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 256b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 2573f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 258b20c023fSPeter Brune } 259eaedb033SPeter Brune PetscFunctionReturn(0); 260eaedb033SPeter Brune } 261eaedb033SPeter Brune 262e0331734SPeter Brune /*@ 263e0331734SPeter Brune SNESNASMSetType - Set the type of subdomain update used 264e0331734SPeter Brune 265e0331734SPeter Brune Logically Collective on SNES 266e0331734SPeter Brune 267e0331734SPeter Brune Input Parameters: 268e0331734SPeter Brune + SNES - the SNES context 269e0331734SPeter Brune - type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT 270e0331734SPeter Brune 271e0331734SPeter Brune Level: intermediate 272e0331734SPeter Brune 273e0331734SPeter Brune .keywords: SNES, NASM 274e0331734SPeter Brune 275e0331734SPeter Brune .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType() 276e0331734SPeter Brune @*/ 277e0331734SPeter Brune PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type) 278e0331734SPeter Brune { 279e0331734SPeter Brune PetscErrorCode ierr; 280e0331734SPeter Brune PetscErrorCode (*f)(SNES,PCASMType); 281e0331734SPeter Brune 282e0331734SPeter Brune PetscFunctionBegin; 283e0331734SPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f);CHKERRQ(ierr); 284e0331734SPeter Brune if (f) {ierr = (f)(snes,type);CHKERRQ(ierr);} 285e0331734SPeter Brune PetscFunctionReturn(0); 286e0331734SPeter Brune } 287e0331734SPeter Brune 28840244768SBarry Smith static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type) 289e0331734SPeter Brune { 290e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 291e0331734SPeter Brune 292e0331734SPeter Brune PetscFunctionBegin; 293e0331734SPeter 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"); 294e0331734SPeter Brune nasm->type = type; 295e0331734SPeter Brune PetscFunctionReturn(0); 296e0331734SPeter Brune } 297e0331734SPeter Brune 298e0331734SPeter Brune /*@ 299e0331734SPeter Brune SNESNASMGetType - Get the type of subdomain update used 300e0331734SPeter Brune 301e0331734SPeter Brune Logically Collective on SNES 302e0331734SPeter Brune 303e0331734SPeter Brune Input Parameters: 304e0331734SPeter Brune . SNES - the SNES context 305e0331734SPeter Brune 306e0331734SPeter Brune Output Parameters: 307e0331734SPeter Brune . type - the type of update 308e0331734SPeter Brune 309e0331734SPeter Brune Level: intermediate 310e0331734SPeter Brune 311e0331734SPeter Brune .keywords: SNES, NASM 312e0331734SPeter Brune 313e0331734SPeter Brune .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType() 314e0331734SPeter Brune @*/ 315e0331734SPeter Brune PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type) 316e0331734SPeter Brune { 317e0331734SPeter Brune PetscErrorCode ierr; 318e0331734SPeter Brune 319e0331734SPeter Brune PetscFunctionBegin; 3202a808120SBarry Smith ierr = PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));CHKERRQ(ierr); 321e0331734SPeter Brune PetscFunctionReturn(0); 322e0331734SPeter Brune } 323e0331734SPeter Brune 32440244768SBarry Smith static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type) 325e0331734SPeter Brune { 326e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 327e0331734SPeter Brune 328e0331734SPeter Brune PetscFunctionBegin; 329e0331734SPeter Brune *type = nasm->type; 330e0331734SPeter Brune PetscFunctionReturn(0); 331e0331734SPeter Brune } 332e0331734SPeter Brune 33376857b2aSPeter Brune /*@ 33476857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 33576857b2aSPeter Brune 33676857b2aSPeter Brune Not Collective 33776857b2aSPeter Brune 33876857b2aSPeter Brune Input Parameters: 33976857b2aSPeter Brune + SNES - the SNES context 34076857b2aSPeter Brune . n - the number of local subdomains 34176857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 34276857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 34376857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 34476857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 34576857b2aSPeter Brune 34676857b2aSPeter Brune Level: intermediate 34776857b2aSPeter Brune 34876857b2aSPeter Brune .keywords: SNES, NASM 34976857b2aSPeter Brune 35076857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 35176857b2aSPeter Brune @*/ 352a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 353a6dfd86eSKarl Rupp { 354eaedb033SPeter Brune PetscErrorCode ierr; 355111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 3566e111a19SKarl Rupp 357eaedb033SPeter Brune PetscFunctionBegin; 3580005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 3594b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 360eaedb033SPeter Brune PetscFunctionReturn(0); 361eaedb033SPeter Brune } 362eaedb033SPeter Brune 36340244768SBarry Smith static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 364a6dfd86eSKarl Rupp { 365eaedb033SPeter Brune PetscInt i; 366eaedb033SPeter Brune PetscErrorCode ierr; 367eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3686e111a19SKarl Rupp 369eaedb033SPeter Brune PetscFunctionBegin; 370ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 371eaedb033SPeter Brune 372111ade9eSPeter Brune /* tear down the previously set things */ 373111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 374111ade9eSPeter Brune 375eaedb033SPeter Brune nasm->n = n; 376111ade9eSPeter Brune if (oscatter) { 377111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 378eaedb033SPeter Brune } 379111ade9eSPeter Brune if (iscatter) { 380111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 381eaedb033SPeter Brune } 382111ade9eSPeter Brune if (gscatter) { 383111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 384111ade9eSPeter Brune } 385111ade9eSPeter Brune if (oscatter) { 386785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 387eaedb033SPeter Brune for (i=0; i<n; i++) { 388111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 389eaedb033SPeter Brune } 390111ade9eSPeter Brune } 391111ade9eSPeter Brune if (iscatter) { 392785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 393eaedb033SPeter Brune for (i=0; i<n; i++) { 394111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 395eaedb033SPeter Brune } 396eaedb033SPeter Brune } 397111ade9eSPeter Brune if (gscatter) { 398785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 399eaedb033SPeter Brune for (i=0; i<n; i++) { 400111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 401eaedb033SPeter Brune } 402eaedb033SPeter Brune } 403111ade9eSPeter Brune 404eaedb033SPeter Brune if (subsnes) { 405785e854fSJed Brown ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 406eaedb033SPeter Brune for (i=0; i<n; i++) { 407eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 408eaedb033SPeter Brune } 409a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 410eaedb033SPeter Brune } 411eaedb033SPeter Brune PetscFunctionReturn(0); 412eaedb033SPeter Brune } 413eaedb033SPeter Brune 41476857b2aSPeter Brune /*@ 41576857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 41676857b2aSPeter Brune 41776857b2aSPeter Brune Not Collective 41876857b2aSPeter Brune 41976857b2aSPeter Brune Input Parameters: 42076857b2aSPeter Brune . SNES - the SNES context 42176857b2aSPeter Brune 42276857b2aSPeter Brune Output Parameters: 42376857b2aSPeter Brune + n - the number of local subdomains 42476857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 42576857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 42676857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 42776857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 42876857b2aSPeter Brune 42976857b2aSPeter Brune Level: intermediate 43076857b2aSPeter Brune 43176857b2aSPeter Brune .keywords: SNES, NASM 43276857b2aSPeter Brune 43376857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 43476857b2aSPeter Brune @*/ 43576857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 43676857b2aSPeter Brune { 43776857b2aSPeter Brune PetscErrorCode ierr; 43876857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 43976857b2aSPeter Brune 44076857b2aSPeter Brune PetscFunctionBegin; 4410005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 4424b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 44376857b2aSPeter Brune PetscFunctionReturn(0); 44476857b2aSPeter Brune } 44576857b2aSPeter Brune 44640244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 44776857b2aSPeter Brune { 44876857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 44976857b2aSPeter Brune 45076857b2aSPeter Brune PetscFunctionBegin; 45176857b2aSPeter Brune if (n) *n = nasm->n; 45276857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 45376857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 45476857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 455a4f17876SPeter Brune if (subsnes) { 456a4f17876SPeter Brune *subsnes = nasm->subsnes; 457a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 458a4f17876SPeter Brune } 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 .keywords: SNES, NASM 48076857b2aSPeter Brune 48176857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 48276857b2aSPeter Brune @*/ 48376857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 48476857b2aSPeter Brune { 48576857b2aSPeter Brune PetscErrorCode ierr; 48676857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 48776857b2aSPeter Brune 48876857b2aSPeter Brune PetscFunctionBegin; 4890005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 4904b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 49176857b2aSPeter Brune PetscFunctionReturn(0); 49276857b2aSPeter Brune } 49376857b2aSPeter Brune 49440244768SBarry Smith static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 49576857b2aSPeter Brune { 49676857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 49776857b2aSPeter Brune 49876857b2aSPeter Brune PetscFunctionBegin; 49976857b2aSPeter Brune if (n) *n = nasm->n; 50076857b2aSPeter Brune if (x) *x = nasm->x; 50176857b2aSPeter Brune if (y) *y = nasm->y; 50276857b2aSPeter Brune if (b) *b = nasm->b; 50376857b2aSPeter Brune if (xl) *xl = nasm->xl; 50476857b2aSPeter Brune PetscFunctionReturn(0); 50576857b2aSPeter Brune } 50676857b2aSPeter Brune 507d728fb7dSPeter Brune /*@ 508d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 509d728fb7dSPeter Brune 510d728fb7dSPeter Brune Collective on SNES 511d728fb7dSPeter Brune 512d728fb7dSPeter Brune Input Parameters: 513d728fb7dSPeter Brune + SNES - the SNES context 514d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 515d728fb7dSPeter Brune 516d728fb7dSPeter Brune Level: developer 517d728fb7dSPeter Brune 518d728fb7dSPeter Brune Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 519d728fb7dSPeter Brune is needed at each linear iteration. 520d728fb7dSPeter Brune 521d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN 522d728fb7dSPeter Brune 523d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 524d728fb7dSPeter Brune @*/ 525d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 526d728fb7dSPeter Brune { 527d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 528d728fb7dSPeter Brune PetscErrorCode ierr; 529d728fb7dSPeter Brune 530d728fb7dSPeter Brune PetscFunctionBegin; 5310005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 5324b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 533d728fb7dSPeter Brune PetscFunctionReturn(0); 534d728fb7dSPeter Brune } 535d728fb7dSPeter Brune 53640244768SBarry Smith static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 537d728fb7dSPeter Brune { 538d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 539d728fb7dSPeter Brune 540d728fb7dSPeter Brune PetscFunctionBegin; 541d728fb7dSPeter Brune nasm->finaljacobian = flg; 542602bec5dSPeter Brune if (flg) snes->usesksp = PETSC_TRUE; 543d728fb7dSPeter Brune PetscFunctionReturn(0); 544d728fb7dSPeter Brune } 54576857b2aSPeter Brune 546610116beSPeter Brune /*@ 547610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 548610116beSPeter Brune 549610116beSPeter Brune Logically collective on SNES 550610116beSPeter Brune 551610116beSPeter Brune Input Parameters: 552610116beSPeter Brune + SNES - the SNES context 553610116beSPeter Brune - dmp - damping 554610116beSPeter Brune 555610116beSPeter Brune Level: intermediate 556610116beSPeter Brune 5575dfa0f3bSBarry Smith Notes: The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5585dfa0f3bSBarry Smith 559610116beSPeter Brune .keywords: SNES, NASM, damping 560610116beSPeter Brune 561610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 562610116beSPeter Brune @*/ 563610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 564610116beSPeter Brune { 565610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 566610116beSPeter Brune PetscErrorCode ierr; 567610116beSPeter Brune 568610116beSPeter Brune PetscFunctionBegin; 569610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 570e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 571610116beSPeter Brune PetscFunctionReturn(0); 572610116beSPeter Brune } 573610116beSPeter Brune 57440244768SBarry Smith static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 575610116beSPeter Brune { 576610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 577610116beSPeter Brune 578610116beSPeter Brune PetscFunctionBegin; 579610116beSPeter Brune nasm->damping = dmp; 580610116beSPeter Brune PetscFunctionReturn(0); 581610116beSPeter Brune } 582610116beSPeter Brune 583610116beSPeter Brune /*@ 584610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 585610116beSPeter Brune 586610116beSPeter Brune Not Collective 587610116beSPeter Brune 588610116beSPeter Brune Input Parameters: 589610116beSPeter Brune + SNES - the SNES context 590610116beSPeter Brune - dmp - damping 591610116beSPeter Brune 592610116beSPeter Brune Level: intermediate 593610116beSPeter Brune 594610116beSPeter Brune .keywords: SNES, NASM, damping 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; 635111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 636111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 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++) { 64670c78f05SPeter Brune /* scatter the solution to the local solution */ 64770c78f05SPeter Brune Xlloc = nasm->xl[i]; 64870c78f05SPeter Brune gscat = nasm->gscatter[i]; 64970c78f05SPeter Brune oscat = nasm->oscatter[i]; 65070c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 65170c78f05SPeter Brune if (B) { 65270c78f05SPeter Brune /* scatter the RHS to the local RHS */ 65370c78f05SPeter Brune Bl = nasm->b[i]; 65470c78f05SPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 65570c78f05SPeter Brune } 65670c78f05SPeter Brune } 657b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 658b20c023fSPeter Brune 659b20c023fSPeter Brune 660b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 66170c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 66270c78f05SPeter Brune Xl = nasm->x[i]; 66370c78f05SPeter Brune Xlloc = nasm->xl[i]; 66470c78f05SPeter Brune Yl = nasm->y[i]; 665258e1594SPeter Brune subsnes = nasm->subsnes[i]; 666258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 667111ade9eSPeter Brune iscat = nasm->iscatter[i]; 668111ade9eSPeter Brune oscat = nasm->oscatter[i]; 669111ade9eSPeter Brune gscat = nasm->gscatter[i]; 670ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 67124b7f281SPeter Brune if (B) { 67224b7f281SPeter Brune Bl = nasm->b[i]; 673ed3c11a9SPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 674ed3c11a9SPeter Brune } else Bl = NULL; 675ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 67614eb1c5cSMatthew G. Knepley /* Could scatter directly from X */ 67770c78f05SPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 67870c78f05SPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 679111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 680ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 681ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 682e0331734SPeter Brune if (type == PC_ASM_BASIC) { 683111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 684e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 685111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 686ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 687eaedb033SPeter Brune } 688ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 689ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 69070c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 69170c78f05SPeter Brune Yl = nasm->y[i]; 69270c78f05SPeter Brune iscat = nasm->iscatter[i]; 69370c78f05SPeter Brune oscat = nasm->oscatter[i]; 694e0331734SPeter Brune if (type == PC_ASM_BASIC) { 69570c78f05SPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 696e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 69770c78f05SPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 698ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 69970c78f05SPeter Brune } 700b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 7015dfa0f3bSBarry Smith ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 702610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 703eaedb033SPeter Brune PetscFunctionReturn(0); 704eaedb033SPeter Brune } 705eaedb033SPeter Brune 70640244768SBarry Smith static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 707d728fb7dSPeter Brune { 708602bec5dSPeter Brune Vec X = Xfinal; 709d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 710d728fb7dSPeter Brune SNES subsnes; 711ca44f815SPeter Brune PetscInt i,lag = 1; 712d728fb7dSPeter Brune PetscErrorCode ierr; 713e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 714d728fb7dSPeter Brune VecScatter oscat,gscat; 715d728fb7dSPeter Brune DM dm,subdm; 716d1e9a80fSBarry Smith 717d728fb7dSPeter Brune PetscFunctionBegin; 718602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 719e59f0a30SPeter Brune F = snes->vec_func; 720365a6726SPeter Brune if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 721d1e9a80fSBarry Smith ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 722d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 723d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 724602bec5dSPeter Brune if (nasm->fjtype != 1) { 725d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 726d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 727d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 728d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 729602bec5dSPeter Brune } 730d728fb7dSPeter Brune } 731d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 732d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 733e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 734d728fb7dSPeter Brune Xl = nasm->x[i]; 735d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 736d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 737d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 738d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 739602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 740ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 741d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 742602bec5dSPeter Brune if (nasm->fjtype != 1) { 743d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 744d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 745602bec5dSPeter Brune } 746ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 747ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 748602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 749d1e9a80fSBarry Smith ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr); 750ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 751d728fb7dSPeter Brune } 752d728fb7dSPeter Brune PetscFunctionReturn(0); 753d728fb7dSPeter Brune } 754d728fb7dSPeter Brune 75540244768SBarry Smith static PetscErrorCode SNESSolve_NASM(SNES snes) 756eaedb033SPeter Brune { 757eaedb033SPeter Brune Vec F; 758eaedb033SPeter Brune Vec X; 759eaedb033SPeter Brune Vec B; 760111ade9eSPeter Brune Vec Y; 761eaedb033SPeter Brune PetscInt i; 762ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 763eaedb033SPeter Brune PetscErrorCode ierr; 764365a6726SPeter Brune SNESNormSchedule normschedule; 765d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 766eaedb033SPeter Brune 767eaedb033SPeter Brune PetscFunctionBegin; 768c579b300SPatrick Farrell 7696c4ed002SBarry 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); 770c579b300SPatrick Farrell 771fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 772eaedb033SPeter Brune X = snes->vec_sol; 773111ade9eSPeter Brune Y = snes->vec_sol_update; 774eaedb033SPeter Brune F = snes->vec_func; 775eaedb033SPeter Brune B = snes->vec_rhs; 776eaedb033SPeter Brune 777e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 778eaedb033SPeter Brune snes->iter = 0; 779eaedb033SPeter Brune snes->norm = 0.; 780e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 781eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 782365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 783365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 784eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 785eaedb033SPeter Brune if (!snes->vec_func_init_set) { 786eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 7871aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 788eaedb033SPeter Brune 789eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 790422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 791e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 792eaedb033SPeter Brune snes->iter = 0; 793eaedb033SPeter Brune snes->norm = fnorm; 794e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 795a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 796eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 797eaedb033SPeter Brune 798eaedb033SPeter Brune /* test convergence */ 799eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 800eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 801eaedb033SPeter Brune } else { 802e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 803a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 804eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 805eaedb033SPeter Brune } 806eaedb033SPeter Brune 807eaedb033SPeter Brune /* Call general purpose update function */ 808eaedb033SPeter Brune if (snes->ops->update) { 809eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 810eaedb033SPeter Brune } 811602bec5dSPeter Brune /* copy the initial solution over for later */ 812602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 813eaedb033SPeter Brune 814eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 815111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 816365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 817eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 818eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 819422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 820eaedb033SPeter Brune } 821eaedb033SPeter Brune /* Monitor convergence */ 822e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 823eaedb033SPeter Brune snes->iter = i+1; 824eaedb033SPeter Brune snes->norm = fnorm; 825e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 826a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 827eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 828eaedb033SPeter Brune /* Test for convergence */ 829365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 830d728fb7dSPeter Brune if (snes->reason) break; 831eaedb033SPeter Brune /* Call general purpose update function */ 832e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 833eaedb033SPeter Brune } 834d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 835365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 836eaedb033SPeter Brune if (i == snes->max_its) { 837eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 838eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 839eaedb033SPeter Brune } 8401aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 841eaedb033SPeter Brune PetscFunctionReturn(0); 842eaedb033SPeter Brune } 843eaedb033SPeter Brune 844eaedb033SPeter Brune /*MC 845eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 846eaedb033SPeter Brune 84770c78603SPeter Brune Options Database: 84870c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 84970c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 8505dfa0f3bSBarry Smith . -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 85170c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 852150d49b7SHong Zhang . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at 85370c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 85470c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 85570c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 85670c78603SPeter Brune 857eaedb033SPeter Brune Level: advanced 858eaedb033SPeter Brune 8594f02bc6aSBarry Smith References: 86096a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8614f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8624f02bc6aSBarry Smith 8635dfa0f3bSBarry Smith .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types), SNESNASMSetType(), SNESNASMGetType(), SNESNASMSetSubdomains(), SNESNASMGetSubdomains(), SNESNASMGetSubdomainVecs(), SNESNASMSetComputeFinalJacobian(), SNESNASMSetDamping(), SNESNASMGetDamping() 864eaedb033SPeter Brune M*/ 865eaedb033SPeter Brune 8668cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 867eaedb033SPeter Brune { 868eaedb033SPeter Brune SNES_NASM *nasm; 869eaedb033SPeter Brune PetscErrorCode ierr; 870eaedb033SPeter Brune 871eaedb033SPeter Brune PetscFunctionBegin; 872b00a9115SJed Brown ierr = PetscNewLog(snes,&nasm);CHKERRQ(ierr); 873eaedb033SPeter Brune snes->data = (void*)nasm; 874eaedb033SPeter Brune 875eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 876eaedb033SPeter Brune nasm->subsnes = 0; 877eaedb033SPeter Brune nasm->x = 0; 878111ade9eSPeter Brune nasm->xl = 0; 879111ade9eSPeter Brune nasm->y = 0; 880eaedb033SPeter Brune nasm->b = 0; 881111ade9eSPeter Brune nasm->oscatter = 0; 882111ade9eSPeter Brune nasm->iscatter = 0; 883111ade9eSPeter Brune nasm->gscatter = 0; 884610116beSPeter Brune nasm->damping = 1.; 885111ade9eSPeter Brune 886111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 887d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 888a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 889eaedb033SPeter Brune 890eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 891eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 892eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 893eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 894eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 895eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 896eaedb033SPeter Brune 897eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 898*efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 899eaedb033SPeter Brune 9004fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 9014fc747eaSLawrence Mitchell 902602bec5dSPeter Brune nasm->fjtype = 0; 903602bec5dSPeter Brune nasm->xinit = NULL; 9040298fd71SBarry Smith nasm->eventrestrictinterp = 0; 9050298fd71SBarry Smith nasm->eventsubsolve = 0; 906b20c023fSPeter Brune 907eaedb033SPeter Brune if (!snes->tolerancesset) { 908eaedb033SPeter Brune snes->max_its = 10000; 909eaedb033SPeter Brune snes->max_funcs = 10000; 910eaedb033SPeter Brune } 911eaedb033SPeter Brune 912e0331734SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr); 913e0331734SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr); 914bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 915bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 91690bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 91790bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 918bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 919bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 920eaedb033SPeter Brune PetscFunctionReturn(0); 921eaedb033SPeter Brune } 92299e0435eSBarry Smith 923448b6425SPatrick Farrell /*@ 924448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 925448b6425SPatrick Farrell 926448b6425SPatrick Farrell Not collective 927448b6425SPatrick Farrell 928448b6425SPatrick Farrell Input Parameters: 929448b6425SPatrick Farrell + snes - the SNES context 930448b6425SPatrick Farrell - i - the number of the subsnes to get 931448b6425SPatrick Farrell 932448b6425SPatrick Farrell Output Parameters: 933448b6425SPatrick Farrell . subsnes - the subsolver context 934448b6425SPatrick Farrell 935448b6425SPatrick Farrell Level: intermediate 936448b6425SPatrick Farrell 937448b6425SPatrick Farrell .keywords: SNES, NASM 938448b6425SPatrick Farrell 939448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetNumber() 940448b6425SPatrick Farrell @*/ 941448b6425SPatrick Farrell PetscErrorCode SNESNASMGetSNES(SNES snes,PetscInt i,SNES *subsnes) 942448b6425SPatrick Farrell { 943448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 944448b6425SPatrick Farrell 945448b6425SPatrick Farrell PetscFunctionBegin; 946448b6425SPatrick Farrell if (i < 0 || i >= nasm->n) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"No such subsolver"); 947448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 948448b6425SPatrick Farrell PetscFunctionReturn(0); 949448b6425SPatrick Farrell } 950448b6425SPatrick Farrell 951448b6425SPatrick Farrell /*@ 952448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 953448b6425SPatrick Farrell 954448b6425SPatrick Farrell Not collective 955448b6425SPatrick Farrell 956448b6425SPatrick Farrell Input Parameters: 957448b6425SPatrick Farrell . snes - the SNES context 958448b6425SPatrick Farrell 959448b6425SPatrick Farrell Output Parameters: 960448b6425SPatrick Farrell . n - the number of subsolvers 961448b6425SPatrick Farrell 962448b6425SPatrick Farrell Level: intermediate 963448b6425SPatrick Farrell 964448b6425SPatrick Farrell .keywords: SNES, NASM 965448b6425SPatrick Farrell 966448b6425SPatrick Farrell .seealso: SNESNASM, SNESNASMGetSNES() 967448b6425SPatrick Farrell @*/ 968448b6425SPatrick Farrell PetscErrorCode SNESNASMGetNumber(SNES snes,PetscInt *n) 969448b6425SPatrick Farrell { 970448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM*)snes->data; 971448b6425SPatrick Farrell 972448b6425SPatrick Farrell PetscFunctionBegin; 973448b6425SPatrick Farrell *n = nasm->n; 974448b6425SPatrick Farrell PetscFunctionReturn(0); 975448b6425SPatrick Farrell } 976448b6425SPatrick Farrell 977