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 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_NASM(SNES snes) 34d71ae5a4SJacob Faibussowitsch { 35eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 36eaedb033SPeter Brune PetscInt i; 376e111a19SKarl Rupp 38eaedb033SPeter Brune PetscFunctionBegin; 39eaedb033SPeter Brune for (i = 0; i < nasm->n; i++) { 409566063dSJacob Faibussowitsch if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i])); 419566063dSJacob Faibussowitsch if (nasm->x) PetscCall(VecDestroy(&nasm->x[i])); 429566063dSJacob Faibussowitsch if (nasm->y) PetscCall(VecDestroy(&nasm->y[i])); 439566063dSJacob Faibussowitsch if (nasm->b) PetscCall(VecDestroy(&nasm->b[i])); 44eaedb033SPeter Brune 459566063dSJacob Faibussowitsch if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i])); 469566063dSJacob Faibussowitsch if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i])); 479566063dSJacob Faibussowitsch if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i])); 489566063dSJacob Faibussowitsch if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i])); 499566063dSJacob Faibussowitsch if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i])); 50eaedb033SPeter Brune } 51111ade9eSPeter Brune 529566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->x)); 539566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->xl)); 549566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->y)); 559566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->b)); 56111ade9eSPeter Brune 579566063dSJacob Faibussowitsch if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit)); 58602bec5dSPeter Brune 599566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->subsnes)); 609566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->oscatter)); 619566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->oscatter_copy)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->iscatter)); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->gscatter)); 64b20c023fSPeter Brune 6548a46eb9SPierre Jolivet if (nasm->weight_set) PetscCall(VecDestroy(&nasm->weight)); 66f10b3e88SPatrick Farrell 67b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 68b20c023fSPeter Brune nasm->eventsubsolve = 0; 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70eaedb033SPeter Brune } 71eaedb033SPeter Brune 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NASM(SNES snes) 73d71ae5a4SJacob Faibussowitsch { 74eaedb033SPeter Brune PetscFunctionBegin; 759566063dSJacob Faibussowitsch PetscCall(SNESReset_NASM(snes)); 762e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", NULL)); 772e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", NULL)); 782e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", NULL)); 792e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", NULL)); 802e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", NULL)); 812e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", NULL)); 822e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", NULL)); 832e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", NULL)); 849566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86eaedb033SPeter Brune } 87eaedb033SPeter Brune 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx) 89d71ae5a4SJacob Faibussowitsch { 90111ade9eSPeter Brune Vec bcs = (Vec)ctx; 916e111a19SKarl Rupp 92111ade9eSPeter Brune PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(VecCopy(bcs, l)); 943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95111ade9eSPeter Brune } 96111ade9eSPeter Brune 97d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NASM(SNES snes) 98d71ae5a4SJacob Faibussowitsch { 99eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 10076857b2aSPeter Brune DM dm, subdm; 101111ade9eSPeter Brune DM *subdms; 102111ade9eSPeter Brune PetscInt i; 103eaedb033SPeter Brune const char *optionsprefix; 104111ade9eSPeter Brune Vec F; 105ed3c11a9SPeter Brune PetscMPIInt size; 106ed3c11a9SPeter Brune KSP ksp; 107ed3c11a9SPeter Brune PC pc; 108eaedb033SPeter Brune 109eaedb033SPeter Brune PetscFunctionBegin; 110eaedb033SPeter Brune if (!nasm->subsnes) { 1119566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1120a696f66SPeter Brune if (dm) { 113eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1149566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms)); 11528b400f6SJacob Faibussowitsch PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 1169566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter)); 1179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy)); 11848a46eb9SPierre Jolivet for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i])); 119eaedb033SPeter Brune 1209566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 1219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes)); 122111ade9eSPeter Brune for (i = 0; i < nasm->n; i++) { 1239566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_SELF, &nasm->subsnes[i])); 1249566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1)); 1259566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix)); 1269566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_")); 1279566063dSJacob Faibussowitsch PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i])); 1289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]), &size)); 129ed3c11a9SPeter Brune if (size == 1) { 1309566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(nasm->subsnes[i], &ksp)); 1319566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1329566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY)); 1339566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCLU)); 134ed3c11a9SPeter Brune } 1359566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(nasm->subsnes[i])); 1369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdms[i])); 137111ade9eSPeter Brune } 1389566063dSJacob Faibussowitsch PetscCall(PetscFree(subdms)); 139ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!"); 140ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!"); 141111ade9eSPeter Brune /* allocate the global vectors */ 14248a46eb9SPierre Jolivet if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x)); 14348a46eb9SPierre Jolivet if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl)); 14448a46eb9SPierre Jolivet if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y)); 14548a46eb9SPierre Jolivet if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b)); 146111ade9eSPeter Brune 147111ade9eSPeter Brune for (i = 0; i < nasm->n; i++) { 1489566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL)); 1499566063dSJacob Faibussowitsch if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i])); 1509566063dSJacob Faibussowitsch if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i])); 1519566063dSJacob Faibussowitsch if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i])); 15276857b2aSPeter Brune if (!nasm->xl[i]) { 1539566063dSJacob Faibussowitsch PetscCall(SNESGetDM(nasm->subsnes[i], &subdm)); 1549566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i])); 1559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i])); 156111ade9eSPeter Brune } 15761ba4676SBarry Smith } 158602bec5dSPeter Brune if (nasm->finaljacobian) { 1599566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 16048a46eb9SPierre Jolivet if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit)); 16148a46eb9SPierre Jolivet for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i])); 162602bec5dSPeter Brune } 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 164eaedb033SPeter Brune } 165eaedb033SPeter Brune 166d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems *PetscOptionsObject) 167d71ae5a4SJacob Faibussowitsch { 168111ade9eSPeter Brune PCASMType asmtype; 16983dc3634SPierre Jolivet PetscBool flg, monflg; 170111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 1716e111a19SKarl Rupp 172eaedb033SPeter Brune PetscFunctionBegin; 173d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options"); 1749566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg)); 1759566063dSJacob Faibussowitsch if (flg) PetscCall(SNESNASMSetType(snes, asmtype)); 176b20c023fSPeter Brune flg = PETSC_FALSE; 177b20c023fSPeter Brune monflg = PETSC_TRUE; 1789566063dSJacob Faibussowitsch PetscCall(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)); 1799566063dSJacob Faibussowitsch if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping)); 1809566063dSJacob Faibussowitsch PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail")); 1819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL)); 1829566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL)); 1839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg)); 184b20c023fSPeter Brune if (flg) { 1859566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve)); 1869566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp)); 187b20c023fSPeter Brune } 188d0609cedSBarry Smith PetscOptionsHeadEnd(); 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 190eaedb033SPeter Brune } 191eaedb033SPeter Brune 192d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 193d71ae5a4SJacob Faibussowitsch { 194b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 195a4f17876SPeter Brune PetscMPIInt rank, size; 196dd2fa690SBarry Smith PetscInt i, N, bsz; 197b20c023fSPeter Brune PetscBool iascii, isstring; 198b20c023fSPeter Brune PetscViewer sviewer; 199ce94432eSBarry Smith MPI_Comm comm; 20083dc3634SPierre Jolivet PetscViewerFormat format; 20183dc3634SPierre Jolivet const char *prefix; 202b20c023fSPeter Brune 203b20c023fSPeter Brune PetscFunctionBegin; 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 2059566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2069566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 2079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2091c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm)); 210b20c023fSPeter Brune if (iascii) { 21163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %" PetscInt_FMT "\n", N)); 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 21383dc3634SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 214a4f17876SPeter Brune if (nasm->subsnes) { 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block on rank 0:\n")); 2169566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 2179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 2189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2199566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 220dd400576SPatrick Sanan if (rank == 0) { 2219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2229566063dSJacob Faibussowitsch PetscCall(SNESView(nasm->subsnes[0], sviewer)); 2239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 224a4f17876SPeter Brune } 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 227a4f17876SPeter Brune } 228a4f17876SPeter Brune } else { 229a4f17876SPeter Brune /* print the solver on each block */ 2309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 23163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, nasm->n)); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following SNES objects:\n")); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 238a4f17876SPeter Brune for (i = 0; i < nasm->n; i++) { 2399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(nasm->x[i], &bsz)); 24063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz)); 2419566063dSJacob Faibussowitsch PetscCall(SNESView(nasm->subsnes[i], sviewer)); 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n")); 243b20c023fSPeter Brune } 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 247a4f17876SPeter Brune } 248b20c023fSPeter Brune } else if (isstring) { 24963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type])); 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2519566063dSJacob Faibussowitsch if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer)); 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 253b20c023fSPeter Brune } 2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 255eaedb033SPeter Brune } 256eaedb033SPeter Brune 257e0331734SPeter Brune /*@ 258420bcc1bSBarry Smith SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM` 259e0331734SPeter Brune 260c3339decSBarry Smith Logically Collective 261e0331734SPeter Brune 262e0331734SPeter Brune Input Parameters: 263f6dfbefdSBarry Smith + snes - the `SNES` context 264f6dfbefdSBarry Smith - type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT` 265e0331734SPeter Brune 266420bcc1bSBarry Smith Options Database Key: 267420bcc1bSBarry Smith . -snes_nasm_type <basic,restrict> - type of subdomain update used 268420bcc1bSBarry Smith 269e0331734SPeter Brune Level: intermediate 270e0331734SPeter Brune 271420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 272e0331734SPeter Brune @*/ 273d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type) 274d71ae5a4SJacob Faibussowitsch { 275e0331734SPeter Brune PetscErrorCode (*f)(SNES, PCASMType); 276e0331734SPeter Brune 277e0331734SPeter Brune PetscFunctionBegin; 2789566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f)); 2799566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, type)); 2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 281e0331734SPeter Brune } 282e0331734SPeter Brune 283d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type) 284d71ae5a4SJacob Faibussowitsch { 285e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 286e0331734SPeter Brune 287e0331734SPeter Brune PetscFunctionBegin; 2880b121fc5SBarry Smith PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types"); 289e0331734SPeter Brune nasm->type = type; 2903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 291e0331734SPeter Brune } 292e0331734SPeter Brune 293e0331734SPeter Brune /*@ 294420bcc1bSBarry Smith SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM` 295e0331734SPeter Brune 296c3339decSBarry Smith Logically Collective 297e0331734SPeter Brune 2982fe279fdSBarry Smith Input Parameter: 299f6dfbefdSBarry Smith . snes - the `SNES` context 300e0331734SPeter Brune 3012fe279fdSBarry Smith Output Parameter: 302e0331734SPeter Brune . type - the type of update 303e0331734SPeter Brune 304e0331734SPeter Brune Level: intermediate 305e0331734SPeter Brune 306420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 307e0331734SPeter Brune @*/ 308d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type) 309d71ae5a4SJacob Faibussowitsch { 310e0331734SPeter Brune PetscFunctionBegin; 311cac4c232SBarry Smith PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type)); 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313e0331734SPeter Brune } 314e0331734SPeter Brune 315d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type) 316d71ae5a4SJacob Faibussowitsch { 317e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 318e0331734SPeter Brune 319e0331734SPeter Brune PetscFunctionBegin; 320e0331734SPeter Brune *type = nasm->type; 3213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 322e0331734SPeter Brune } 323e0331734SPeter Brune 32476857b2aSPeter Brune /*@ 325f6dfbefdSBarry Smith SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver 32676857b2aSPeter Brune 327f6dfbefdSBarry Smith Logically Collective 32876857b2aSPeter Brune 32976857b2aSPeter Brune Input Parameters: 330420bcc1bSBarry Smith + snes - the `SNES` context 33176857b2aSPeter Brune . n - the number of local subdomains 33276857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 33376857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 33476857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 33576857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 33676857b2aSPeter Brune 33776857b2aSPeter Brune Level: intermediate 33876857b2aSPeter Brune 339420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 34076857b2aSPeter Brune @*/ 341d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 342d71ae5a4SJacob Faibussowitsch { 343111ade9eSPeter Brune PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *); 3446e111a19SKarl Rupp 345eaedb033SPeter Brune PetscFunctionBegin; 3469566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f)); 3479566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 349eaedb033SPeter Brune } 350eaedb033SPeter Brune 351d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 352d71ae5a4SJacob Faibussowitsch { 353eaedb033SPeter Brune PetscInt i; 354eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 3556e111a19SKarl Rupp 356eaedb033SPeter Brune PetscFunctionBegin; 35728b400f6SJacob Faibussowitsch PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 358eaedb033SPeter Brune 359111ade9eSPeter Brune /* tear down the previously set things */ 3609566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 361111ade9eSPeter Brune 362eaedb033SPeter Brune nasm->n = n; 363111ade9eSPeter Brune if (oscatter) { 3649566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i])); 365eaedb033SPeter Brune } 366111ade9eSPeter Brune if (iscatter) { 3679566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i])); 368eaedb033SPeter Brune } 369111ade9eSPeter Brune if (gscatter) { 3709566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i])); 371111ade9eSPeter Brune } 372111ade9eSPeter Brune if (oscatter) { 3739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->oscatter)); 3749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->oscatter_copy)); 375eaedb033SPeter Brune for (i = 0; i < n; i++) { 376111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 3779566063dSJacob Faibussowitsch PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i])); 378eaedb033SPeter Brune } 379111ade9eSPeter Brune } 380111ade9eSPeter Brune if (iscatter) { 3819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->iscatter)); 382ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i]; 383eaedb033SPeter Brune } 384111ade9eSPeter Brune if (gscatter) { 3859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->gscatter)); 386ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i]; 387eaedb033SPeter Brune } 388111ade9eSPeter Brune 389eaedb033SPeter Brune if (subsnes) { 3909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->subsnes)); 391ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i]; 392eaedb033SPeter Brune } 3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 394eaedb033SPeter Brune } 395eaedb033SPeter Brune 39676857b2aSPeter Brune /*@ 397f6dfbefdSBarry Smith SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver 39876857b2aSPeter Brune 399f6dfbefdSBarry Smith Not Collective but some of the objects returned will be parallel 40076857b2aSPeter Brune 401f899ff85SJose E. Roman Input Parameter: 402f6dfbefdSBarry Smith . snes - the `SNES` context 40376857b2aSPeter Brune 40476857b2aSPeter Brune Output Parameters: 40576857b2aSPeter Brune + n - the number of local subdomains 40676857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 40776857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 40876857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 40976857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 41076857b2aSPeter Brune 41176857b2aSPeter Brune Level: intermediate 41276857b2aSPeter Brune 413420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()` 41476857b2aSPeter Brune @*/ 415d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 416d71ae5a4SJacob Faibussowitsch { 41776857b2aSPeter Brune PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **); 41876857b2aSPeter Brune 41976857b2aSPeter Brune PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f)); 4219566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 4223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42376857b2aSPeter Brune } 42476857b2aSPeter Brune 425d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 426d71ae5a4SJacob Faibussowitsch { 42776857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 42876857b2aSPeter Brune 42976857b2aSPeter Brune PetscFunctionBegin; 43076857b2aSPeter Brune if (n) *n = nasm->n; 43176857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 43276857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 43376857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 43483dc3634SPierre Jolivet if (subsnes) *subsnes = nasm->subsnes; 4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43676857b2aSPeter Brune } 43776857b2aSPeter Brune 43876857b2aSPeter Brune /*@ 439f6dfbefdSBarry Smith SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver 44076857b2aSPeter Brune 44176857b2aSPeter Brune Not Collective 44276857b2aSPeter Brune 443f899ff85SJose E. Roman Input Parameter: 444f6dfbefdSBarry Smith . snes - the `SNES` context 44576857b2aSPeter Brune 44676857b2aSPeter Brune Output Parameters: 44776857b2aSPeter Brune + n - the number of local subdomains 44876857b2aSPeter Brune . x - The subdomain solution vector 44976857b2aSPeter Brune . y - The subdomain step vector 45076857b2aSPeter Brune . b - The subdomain RHS vector 45176857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 45276857b2aSPeter Brune 45376857b2aSPeter Brune Level: developer 45476857b2aSPeter Brune 455420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 45676857b2aSPeter Brune @*/ 457d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 458d71ae5a4SJacob Faibussowitsch { 45976857b2aSPeter Brune PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **); 46076857b2aSPeter Brune 46176857b2aSPeter Brune PetscFunctionBegin; 4629566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f)); 4639566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, x, y, b, xl)); 4643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46576857b2aSPeter Brune } 46676857b2aSPeter Brune 467d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 468d71ae5a4SJacob Faibussowitsch { 46976857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 47076857b2aSPeter Brune 47176857b2aSPeter Brune PetscFunctionBegin; 47276857b2aSPeter Brune if (n) *n = nasm->n; 47376857b2aSPeter Brune if (x) *x = nasm->x; 47476857b2aSPeter Brune if (y) *y = nasm->y; 47576857b2aSPeter Brune if (b) *b = nasm->b; 47676857b2aSPeter Brune if (xl) *xl = nasm->xl; 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47876857b2aSPeter Brune } 47976857b2aSPeter Brune 480d728fb7dSPeter Brune /*@ 481f6dfbefdSBarry Smith SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the 482f6dfbefdSBarry Smith nonlinear additive Schwarz solver 483d728fb7dSPeter Brune 484c3339decSBarry Smith Collective 485d728fb7dSPeter Brune 486d728fb7dSPeter Brune Input Parameters: 487f6dfbefdSBarry Smith + snes - the SNES context 488f6dfbefdSBarry Smith - flg - `PETSC_TRUE` to compute the Jacobians 489d728fb7dSPeter Brune 490d728fb7dSPeter Brune Level: developer 491d728fb7dSPeter Brune 49295452b02SPatrick Sanan Notes: 493f6dfbefdSBarry Smith This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian 494d728fb7dSPeter Brune is needed at each linear iteration. 495d728fb7dSPeter Brune 496420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 497d728fb7dSPeter Brune @*/ 498d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg) 499d71ae5a4SJacob Faibussowitsch { 500d728fb7dSPeter Brune PetscErrorCode (*f)(SNES, PetscBool); 501d728fb7dSPeter Brune 502d728fb7dSPeter Brune PetscFunctionBegin; 5039566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f)); 5049566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg)); 5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 506d728fb7dSPeter Brune } 507d728fb7dSPeter Brune 508d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg) 509d71ae5a4SJacob Faibussowitsch { 510d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 511d728fb7dSPeter Brune 512d728fb7dSPeter Brune PetscFunctionBegin; 513d728fb7dSPeter Brune nasm->finaljacobian = flg; 5143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 515d728fb7dSPeter Brune } 51676857b2aSPeter Brune 517610116beSPeter Brune /*@ 518f6dfbefdSBarry Smith SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 519610116beSPeter Brune 52020f4b53cSBarry Smith Logically Collective 521610116beSPeter Brune 522610116beSPeter Brune Input Parameters: 523f6dfbefdSBarry Smith + snes - the `SNES` context 524610116beSPeter Brune - dmp - damping 525610116beSPeter Brune 526420bcc1bSBarry Smith Options Database Key: 527420bcc1bSBarry Smith . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains) 528420bcc1bSBarry Smith 529610116beSPeter Brune Level: intermediate 530610116beSPeter Brune 531420bcc1bSBarry Smith Note: 53295452b02SPatrick Sanan The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5335dfa0f3bSBarry Smith 534420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()` 535610116beSPeter Brune @*/ 536d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp) 537d71ae5a4SJacob Faibussowitsch { 538610116beSPeter Brune PetscErrorCode (*f)(SNES, PetscReal); 539610116beSPeter Brune 540610116beSPeter Brune PetscFunctionBegin; 5419566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f)); 5429566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, dmp)); 5433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 544610116beSPeter Brune } 545610116beSPeter Brune 546d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp) 547d71ae5a4SJacob Faibussowitsch { 548610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 549610116beSPeter Brune 550610116beSPeter Brune PetscFunctionBegin; 551610116beSPeter Brune nasm->damping = dmp; 5523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 553610116beSPeter Brune } 554610116beSPeter Brune 555610116beSPeter Brune /*@ 556f6dfbefdSBarry Smith SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 557610116beSPeter Brune 558610116beSPeter Brune Not Collective 559610116beSPeter Brune 560f6dfbefdSBarry Smith Input Parameter: 561420bcc1bSBarry Smith . snes - the `SNES` context 562f6dfbefdSBarry Smith 563f6dfbefdSBarry Smith Output Parameter: 564f6dfbefdSBarry Smith . dmp - damping 565610116beSPeter Brune 566610116beSPeter Brune Level: intermediate 567610116beSPeter Brune 568420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()` 569610116beSPeter Brune @*/ 570d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp) 571d71ae5a4SJacob Faibussowitsch { 572610116beSPeter Brune PetscFunctionBegin; 573cac4c232SBarry Smith PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp)); 5743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 575610116beSPeter Brune } 576610116beSPeter Brune 577d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp) 578d71ae5a4SJacob Faibussowitsch { 579610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 580610116beSPeter Brune 581610116beSPeter Brune PetscFunctionBegin; 582610116beSPeter Brune *dmp = nasm->damping; 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 584610116beSPeter Brune } 585610116beSPeter Brune 58614eb1c5cSMatthew G. Knepley /* 58714eb1c5cSMatthew G. Knepley Input Parameters: 58814eb1c5cSMatthew G. Knepley + snes - The solver 58914eb1c5cSMatthew G. Knepley . B - The RHS vector 59014eb1c5cSMatthew G. Knepley - X - The initial guess 59114eb1c5cSMatthew G. Knepley 5922fe279fdSBarry Smith Output Parameter: 59314eb1c5cSMatthew G. Knepley . Y - The solution update 59414eb1c5cSMatthew G. Knepley 59514eb1c5cSMatthew G. Knepley TODO: All scatters should be packed into one 59614eb1c5cSMatthew G. Knepley */ 59766976f2fSJacob Faibussowitsch static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X) 598d71ae5a4SJacob Faibussowitsch { 599eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 600258e1594SPeter Brune SNES subsnes; 601eaedb033SPeter Brune PetscInt i; 602610116beSPeter Brune PetscReal dmp; 603f10b3e88SPatrick Farrell Vec Xl, Bl, Yl, Xlloc; 604f10b3e88SPatrick Farrell VecScatter iscat, oscat, gscat, oscat_copy; 605111ade9eSPeter Brune DM dm, subdm; 606e0331734SPeter Brune PCASMType type; 6070adebc6cSBarry Smith 608eaedb033SPeter Brune PetscFunctionBegin; 6099566063dSJacob Faibussowitsch PetscCall(SNESNASMGetType(snes, &type)); 6109566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 6119566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0)); 6129566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 613eaedb033SPeter Brune for (i = 0; i < nasm->n; i++) { 614f10b3e88SPatrick Farrell /* scatter the solution to the global solution and the local solution */ 615f10b3e88SPatrick Farrell Xl = nasm->x[i]; 61670c78f05SPeter Brune Xlloc = nasm->xl[i]; 61770c78f05SPeter Brune oscat = nasm->oscatter[i]; 618f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 619f10b3e88SPatrick Farrell gscat = nasm->gscatter[i]; 6209566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 6219566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 62270c78f05SPeter Brune if (B) { 62370c78f05SPeter Brune /* scatter the RHS to the local RHS */ 62470c78f05SPeter Brune Bl = nasm->b[i]; 6259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 62670c78f05SPeter Brune } 62770c78f05SPeter Brune } 6289566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 629b20c023fSPeter Brune 6309566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0)); 63170c78f05SPeter Brune for (i = 0; i < nasm->n; i++) { 63270c78f05SPeter Brune Xl = nasm->x[i]; 63370c78f05SPeter Brune Xlloc = nasm->xl[i]; 63470c78f05SPeter Brune Yl = nasm->y[i]; 635258e1594SPeter Brune subsnes = nasm->subsnes[i]; 6369566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes, &subdm)); 637111ade9eSPeter Brune iscat = nasm->iscatter[i]; 638111ade9eSPeter Brune oscat = nasm->oscatter[i]; 639f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 640111ade9eSPeter Brune gscat = nasm->gscatter[i]; 6419566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 6429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 64324b7f281SPeter Brune if (B) { 64424b7f281SPeter Brune Bl = nasm->b[i]; 6459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 646ed3c11a9SPeter Brune } else Bl = NULL; 647f10b3e88SPatrick Farrell 6489566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 6499566063dSJacob Faibussowitsch PetscCall(VecCopy(Xl, Yl)); 6509566063dSJacob Faibussowitsch PetscCall(SNESSolve(subsnes, Bl, Xl)); 6519566063dSJacob Faibussowitsch PetscCall(VecAYPX(Yl, -1.0, Xl)); 6529566063dSJacob Faibussowitsch PetscCall(VecScale(Yl, nasm->damping)); 653e0331734SPeter Brune if (type == PC_ASM_BASIC) { 6549566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 6559566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 656e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 6579566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 6589566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 659ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM"); 660eaedb033SPeter Brune } 6619566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0)); 6629566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 6631baa6e33SBarry Smith if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight)); 6649566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 6659566063dSJacob Faibussowitsch PetscCall(SNESNASMGetDamping(snes, &dmp)); 6669566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, dmp, Y)); 6673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 668eaedb033SPeter Brune } 669eaedb033SPeter Brune 670d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 671d71ae5a4SJacob Faibussowitsch { 672602bec5dSPeter Brune Vec X = Xfinal; 673d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 674d728fb7dSPeter Brune SNES subsnes; 675ca44f815SPeter Brune PetscInt i, lag = 1; 676e59f0a30SPeter Brune Vec Xlloc, Xl, Fl, F; 677d728fb7dSPeter Brune VecScatter oscat, gscat; 678d728fb7dSPeter Brune DM dm, subdm; 679d1e9a80fSBarry Smith 680d728fb7dSPeter Brune PetscFunctionBegin; 681602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 682e59f0a30SPeter Brune F = snes->vec_func; 6839566063dSJacob Faibussowitsch if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F)); 6849566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 6859566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 6869566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 687602bec5dSPeter Brune if (nasm->fjtype != 1) { 688d728fb7dSPeter Brune for (i = 0; i < nasm->n; i++) { 689d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 690d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 6919566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 692602bec5dSPeter Brune } 693d728fb7dSPeter Brune } 6949566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 695d728fb7dSPeter Brune for (i = 0; i < nasm->n; i++) { 696e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 697d728fb7dSPeter Brune Xl = nasm->x[i]; 698d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 699d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 700d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 701d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 7029566063dSJacob Faibussowitsch if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 7039566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes, &subdm)); 7049566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 705602bec5dSPeter Brune if (nasm->fjtype != 1) { 7069566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl)); 7079566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl)); 708602bec5dSPeter Brune } 709ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 710ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 7119566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(subsnes, Xl, Fl)); 7129566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre)); 713ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 714d728fb7dSPeter Brune } 7153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 716d728fb7dSPeter Brune } 717d728fb7dSPeter Brune 718d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NASM(SNES snes) 719d71ae5a4SJacob Faibussowitsch { 720eaedb033SPeter Brune Vec F; 721eaedb033SPeter Brune Vec X; 722eaedb033SPeter Brune Vec B; 723111ade9eSPeter Brune Vec Y; 724eaedb033SPeter Brune PetscInt i; 725ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 726365a6726SPeter Brune SNESNormSchedule normschedule; 727d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 728eaedb033SPeter Brune 729eaedb033SPeter Brune PetscFunctionBegin; 730c579b300SPatrick Farrell 7310b121fc5SBarry Smith PetscCheck(!snes->xl & !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 732c579b300SPatrick Farrell 7339566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 734eaedb033SPeter Brune X = snes->vec_sol; 735111ade9eSPeter Brune Y = snes->vec_sol_update; 736eaedb033SPeter Brune F = snes->vec_func; 737eaedb033SPeter Brune B = snes->vec_rhs; 738eaedb033SPeter Brune 7399566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 740eaedb033SPeter Brune snes->iter = 0; 741eaedb033SPeter Brune snes->norm = 0.; 7429566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 743eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 7449566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 745365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 746eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 747eaedb033SPeter Brune if (!snes->vec_func_init_set) { 7489566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 7491aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 750eaedb033SPeter Brune 7519566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 752422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 7539566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 754eaedb033SPeter Brune snes->iter = 0; 755eaedb033SPeter Brune snes->norm = fnorm; 7569566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7579566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 758eaedb033SPeter Brune 759eaedb033SPeter Brune /* test convergence */ 7602d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 7612d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, snes->norm)); 7623ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 763eaedb033SPeter Brune } else { 7649566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7659566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7662d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 767eaedb033SPeter Brune } 768eaedb033SPeter Brune 769eaedb033SPeter Brune /* Call general purpose update function */ 770dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 771602bec5dSPeter Brune /* copy the initial solution over for later */ 7729566063dSJacob Faibussowitsch if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit)); 773eaedb033SPeter Brune 774eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 7759566063dSJacob Faibussowitsch PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X)); 776365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 7779566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 7789566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 779422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 780eaedb033SPeter Brune } 781eaedb033SPeter Brune /* Monitor convergence */ 7829566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 783eaedb033SPeter Brune snes->iter = i + 1; 784eaedb033SPeter Brune snes->norm = fnorm; 7859566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7869566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 787eaedb033SPeter Brune /* Test for convergence */ 7882d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 7892d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 790d728fb7dSPeter Brune if (snes->reason) break; 791eaedb033SPeter Brune /* Call general purpose update function */ 792dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 793eaedb033SPeter Brune } 79407b62357SFande Kong if (nasm->finaljacobian) { 7959566063dSJacob Faibussowitsch PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X)); 79607b62357SFande Kong SNESCheckJacobianDomainerror(snes); 79707b62357SFande Kong } 7983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 799eaedb033SPeter Brune } 800eaedb033SPeter Brune 801eaedb033SPeter Brune /*MC 802*1d27aa22SBarry Smith SNESNASM - Nonlinear Additive Schwarz solver {cite}`ck02`, {cite}`bruneknepleysmithtu15` 803eaedb033SPeter Brune 804f6dfbefdSBarry Smith Options Database Keys: 80570c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 80670c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 807420bcc1bSBarry Smith . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 808*1d27aa22SBarry Smith . -snes_nasm_finaljacobian - compute the local and global Jacobians of the final iterate 809*1d27aa22SBarry Smith . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the Jacobian is calculated at 81070c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 81170c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 81270c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 81370c78603SPeter Brune 814eaedb033SPeter Brune Level: advanced 815eaedb033SPeter Brune 816f6dfbefdSBarry Smith Note: 817f6dfbefdSBarry Smith This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for 818f6dfbefdSBarry Smith the `SNESASPIN` solver 819f6dfbefdSBarry Smith 820f6dfbefdSBarry Smith Developer Note: 821f6dfbefdSBarry Smith This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to 822f6dfbefdSBarry Smith false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if 823f6dfbefdSBarry Smith `SNESNASM` is used as a nonlinear preconditioner for `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the 824f6dfbefdSBarry Smith Jacobian (needed by `SNESASPIN`) 825f6dfbefdSBarry Smith and this utilizes the inner `KSP` object for storing the matrices, but the `KSP` is never used for solving a linear system. When `SNESNASM` is 826f6dfbefdSBarry Smith used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES` 827f6dfbefdSBarry Smith object (in this case `SNESNASM`) to inherit the outer Jacobian matrices. 828951fe5abSBarry Smith 829420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`, 830420bcc1bSBarry Smith `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`, 831420bcc1bSBarry Smith `SNESNASMGetSNES()`, `SNESNASMGetNumber()` 832eaedb033SPeter Brune M*/ 833eaedb033SPeter Brune 834d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 835d71ae5a4SJacob Faibussowitsch { 836eaedb033SPeter Brune SNES_NASM *nasm; 837eaedb033SPeter Brune 838eaedb033SPeter Brune PetscFunctionBegin; 8394dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&nasm)); 840eaedb033SPeter Brune snes->data = (void *)nasm; 841eaedb033SPeter Brune 842eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 8439e5d0892SLisandro Dalcin nasm->subsnes = NULL; 8449e5d0892SLisandro Dalcin nasm->x = NULL; 8459e5d0892SLisandro Dalcin nasm->xl = NULL; 8469e5d0892SLisandro Dalcin nasm->y = NULL; 8479e5d0892SLisandro Dalcin nasm->b = NULL; 8489e5d0892SLisandro Dalcin nasm->oscatter = NULL; 8499e5d0892SLisandro Dalcin nasm->oscatter_copy = NULL; 8509e5d0892SLisandro Dalcin nasm->iscatter = NULL; 8519e5d0892SLisandro Dalcin nasm->gscatter = NULL; 852610116beSPeter Brune nasm->damping = 1.; 853111ade9eSPeter Brune 854111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 855d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 856f10b3e88SPatrick Farrell nasm->weight_set = PETSC_FALSE; 857eaedb033SPeter Brune 858eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 859eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 860eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 861eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 862eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 863eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 864eaedb033SPeter Brune 865eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 866efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 867eaedb033SPeter Brune 8684fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8694fc747eaSLawrence Mitchell 870602bec5dSPeter Brune nasm->fjtype = 0; 871602bec5dSPeter Brune nasm->xinit = NULL; 8720298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8730298fd71SBarry Smith nasm->eventsubsolve = 0; 874b20c023fSPeter Brune 875eaedb033SPeter Brune if (!snes->tolerancesset) { 876eaedb033SPeter Brune snes->max_its = 10000; 877eaedb033SPeter Brune snes->max_funcs = 10000; 878eaedb033SPeter Brune } 879eaedb033SPeter Brune 8809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM)); 8819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM)); 8829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM)); 8839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM)); 8849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM)); 8859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM)); 8869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM)); 8879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM)); 8883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 889eaedb033SPeter Brune } 89099e0435eSBarry Smith 891448b6425SPatrick Farrell /*@ 892448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 893448b6425SPatrick Farrell 89420f4b53cSBarry Smith Not Collective 895448b6425SPatrick Farrell 896448b6425SPatrick Farrell Input Parameters: 89720f4b53cSBarry Smith + snes - the `SNES` context 898448b6425SPatrick Farrell - i - the number of the subsnes to get 899448b6425SPatrick Farrell 9002fe279fdSBarry Smith Output Parameter: 901448b6425SPatrick Farrell . subsnes - the subsolver context 902448b6425SPatrick Farrell 903448b6425SPatrick Farrell Level: intermediate 904448b6425SPatrick Farrell 905420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()` 906448b6425SPatrick Farrell @*/ 907d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes) 908d71ae5a4SJacob Faibussowitsch { 909448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 910448b6425SPatrick Farrell 911448b6425SPatrick Farrell PetscFunctionBegin; 9120b121fc5SBarry Smith PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver"); 913448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 9143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 915448b6425SPatrick Farrell } 916448b6425SPatrick Farrell 917448b6425SPatrick Farrell /*@ 918448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 919448b6425SPatrick Farrell 92020f4b53cSBarry Smith Not Collective 921448b6425SPatrick Farrell 9222fe279fdSBarry Smith Input Parameter: 92320f4b53cSBarry Smith . snes - the `SNES` context 924448b6425SPatrick Farrell 9252fe279fdSBarry Smith Output Parameter: 926448b6425SPatrick Farrell . n - the number of subsolvers 927448b6425SPatrick Farrell 928448b6425SPatrick Farrell Level: intermediate 929448b6425SPatrick Farrell 930420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()` 931448b6425SPatrick Farrell @*/ 932d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n) 933d71ae5a4SJacob Faibussowitsch { 934448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 935448b6425SPatrick Farrell 936448b6425SPatrick Farrell PetscFunctionBegin; 937448b6425SPatrick Farrell *n = nasm->n; 9383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 939448b6425SPatrick Farrell } 940448b6425SPatrick Farrell 941f10b3e88SPatrick Farrell /*@ 942f10b3e88SPatrick Farrell SNESNASMSetWeight - Sets weight to use when adding overlapping updates 943f10b3e88SPatrick Farrell 944f10b3e88SPatrick Farrell Collective 945f10b3e88SPatrick Farrell 946f10b3e88SPatrick Farrell Input Parameters: 94720f4b53cSBarry Smith + snes - the `SNES` context 948f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 949f10b3e88SPatrick Farrell 950f10b3e88SPatrick Farrell Level: intermediate 951f10b3e88SPatrick Farrell 952420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM` 953f10b3e88SPatrick Farrell @*/ 954d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight) 955d71ae5a4SJacob Faibussowitsch { 956f10b3e88SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 957f10b3e88SPatrick Farrell 958f10b3e88SPatrick Farrell PetscFunctionBegin; 959f10b3e88SPatrick Farrell 9609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nasm->weight)); 961f10b3e88SPatrick Farrell nasm->weight_set = PETSC_TRUE; 962f10b3e88SPatrick Farrell nasm->weight = weight; 9639566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)nasm->weight)); 964f10b3e88SPatrick Farrell 9653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 966f10b3e88SPatrick Farrell } 967