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 /*@ 258f6dfbefdSBarry Smith SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver 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 266e0331734SPeter Brune Level: intermediate 267e0331734SPeter Brune 268f6dfbefdSBarry Smith .seealso: `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 269e0331734SPeter Brune @*/ 270d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type) 271d71ae5a4SJacob Faibussowitsch { 272e0331734SPeter Brune PetscErrorCode (*f)(SNES, PCASMType); 273e0331734SPeter Brune 274e0331734SPeter Brune PetscFunctionBegin; 2759566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f)); 2769566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, type)); 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278e0331734SPeter Brune } 279e0331734SPeter Brune 280d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type) 281d71ae5a4SJacob Faibussowitsch { 282e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 283e0331734SPeter Brune 284e0331734SPeter Brune PetscFunctionBegin; 2850b121fc5SBarry Smith PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types"); 286e0331734SPeter Brune nasm->type = type; 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 288e0331734SPeter Brune } 289e0331734SPeter Brune 290e0331734SPeter Brune /*@ 291f6dfbefdSBarry Smith SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver 292e0331734SPeter Brune 293c3339decSBarry Smith Logically Collective 294e0331734SPeter Brune 295e0331734SPeter Brune Input Parameters: 296f6dfbefdSBarry Smith . snes - the `SNES` context 297e0331734SPeter Brune 298e0331734SPeter Brune Output Parameters: 299e0331734SPeter Brune . type - the type of update 300e0331734SPeter Brune 301e0331734SPeter Brune Level: intermediate 302e0331734SPeter Brune 303f6dfbefdSBarry Smith .seealso: `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 304e0331734SPeter Brune @*/ 305d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type) 306d71ae5a4SJacob Faibussowitsch { 307e0331734SPeter Brune PetscFunctionBegin; 308cac4c232SBarry Smith PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310e0331734SPeter Brune } 311e0331734SPeter Brune 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type) 313d71ae5a4SJacob Faibussowitsch { 314e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 315e0331734SPeter Brune 316e0331734SPeter Brune PetscFunctionBegin; 317e0331734SPeter Brune *type = nasm->type; 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319e0331734SPeter Brune } 320e0331734SPeter Brune 32176857b2aSPeter Brune /*@ 322f6dfbefdSBarry Smith SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver 32376857b2aSPeter Brune 324f6dfbefdSBarry Smith Logically Collective 32576857b2aSPeter Brune 32676857b2aSPeter Brune Input Parameters: 327f6dfbefdSBarry Smith + snes - the SNES context 32876857b2aSPeter Brune . n - the number of local subdomains 32976857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 33076857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 33176857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 33276857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 33376857b2aSPeter Brune 33476857b2aSPeter Brune Level: intermediate 33576857b2aSPeter Brune 336db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()` 33776857b2aSPeter Brune @*/ 338d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 339d71ae5a4SJacob Faibussowitsch { 340111ade9eSPeter Brune PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *); 3416e111a19SKarl Rupp 342eaedb033SPeter Brune PetscFunctionBegin; 3439566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f)); 3449566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346eaedb033SPeter Brune } 347eaedb033SPeter Brune 348d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 349d71ae5a4SJacob Faibussowitsch { 350eaedb033SPeter Brune PetscInt i; 351eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 3526e111a19SKarl Rupp 353eaedb033SPeter Brune PetscFunctionBegin; 35428b400f6SJacob Faibussowitsch PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 355eaedb033SPeter Brune 356111ade9eSPeter Brune /* tear down the previously set things */ 3579566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 358111ade9eSPeter Brune 359eaedb033SPeter Brune nasm->n = n; 360111ade9eSPeter Brune if (oscatter) { 3619566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i])); 362eaedb033SPeter Brune } 363111ade9eSPeter Brune if (iscatter) { 3649566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i])); 365eaedb033SPeter Brune } 366111ade9eSPeter Brune if (gscatter) { 3679566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i])); 368111ade9eSPeter Brune } 369111ade9eSPeter Brune if (oscatter) { 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->oscatter)); 3719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->oscatter_copy)); 372eaedb033SPeter Brune for (i = 0; i < n; i++) { 373111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 3749566063dSJacob Faibussowitsch PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i])); 375eaedb033SPeter Brune } 376111ade9eSPeter Brune } 377111ade9eSPeter Brune if (iscatter) { 3789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->iscatter)); 379ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i]; 380eaedb033SPeter Brune } 381111ade9eSPeter Brune if (gscatter) { 3829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->gscatter)); 383ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i]; 384eaedb033SPeter Brune } 385111ade9eSPeter Brune 386eaedb033SPeter Brune if (subsnes) { 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->subsnes)); 388ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i]; 389eaedb033SPeter Brune } 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391eaedb033SPeter Brune } 392eaedb033SPeter Brune 39376857b2aSPeter Brune /*@ 394f6dfbefdSBarry Smith SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver 39576857b2aSPeter Brune 396f6dfbefdSBarry Smith Not Collective but some of the objects returned will be parallel 39776857b2aSPeter Brune 398f899ff85SJose E. Roman Input Parameter: 399f6dfbefdSBarry Smith . snes - the `SNES` context 40076857b2aSPeter Brune 40176857b2aSPeter Brune Output Parameters: 40276857b2aSPeter Brune + n - the number of local subdomains 40376857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 40476857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 40576857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 40676857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 40776857b2aSPeter Brune 40876857b2aSPeter Brune Level: intermediate 40976857b2aSPeter Brune 410db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetSubdomains()` 41176857b2aSPeter Brune @*/ 412d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 413d71ae5a4SJacob Faibussowitsch { 41476857b2aSPeter Brune PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **); 41576857b2aSPeter Brune 41676857b2aSPeter Brune PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f)); 4189566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42076857b2aSPeter Brune } 42176857b2aSPeter Brune 422d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 423d71ae5a4SJacob Faibussowitsch { 42476857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 42576857b2aSPeter Brune 42676857b2aSPeter Brune PetscFunctionBegin; 42776857b2aSPeter Brune if (n) *n = nasm->n; 42876857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 42976857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 43076857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 43183dc3634SPierre Jolivet if (subsnes) *subsnes = nasm->subsnes; 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43376857b2aSPeter Brune } 43476857b2aSPeter Brune 43576857b2aSPeter Brune /*@ 436f6dfbefdSBarry Smith SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver 43776857b2aSPeter Brune 43876857b2aSPeter Brune Not Collective 43976857b2aSPeter Brune 440f899ff85SJose E. Roman Input Parameter: 441f6dfbefdSBarry Smith . snes - the `SNES` context 44276857b2aSPeter Brune 44376857b2aSPeter Brune Output Parameters: 44476857b2aSPeter Brune + n - the number of local subdomains 44576857b2aSPeter Brune . x - The subdomain solution vector 44676857b2aSPeter Brune . y - The subdomain step vector 44776857b2aSPeter Brune . b - The subdomain RHS vector 44876857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 44976857b2aSPeter Brune 45076857b2aSPeter Brune Level: developer 45176857b2aSPeter Brune 452db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()` 45376857b2aSPeter Brune @*/ 454d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 455d71ae5a4SJacob Faibussowitsch { 45676857b2aSPeter Brune PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **); 45776857b2aSPeter Brune 45876857b2aSPeter Brune PetscFunctionBegin; 4599566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f)); 4609566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, x, y, b, xl)); 4613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46276857b2aSPeter Brune } 46376857b2aSPeter Brune 464d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 465d71ae5a4SJacob Faibussowitsch { 46676857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 46776857b2aSPeter Brune 46876857b2aSPeter Brune PetscFunctionBegin; 46976857b2aSPeter Brune if (n) *n = nasm->n; 47076857b2aSPeter Brune if (x) *x = nasm->x; 47176857b2aSPeter Brune if (y) *y = nasm->y; 47276857b2aSPeter Brune if (b) *b = nasm->b; 47376857b2aSPeter Brune if (xl) *xl = nasm->xl; 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47576857b2aSPeter Brune } 47676857b2aSPeter Brune 477d728fb7dSPeter Brune /*@ 478f6dfbefdSBarry Smith SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the 479f6dfbefdSBarry Smith nonlinear additive Schwarz solver 480d728fb7dSPeter Brune 481c3339decSBarry Smith Collective 482d728fb7dSPeter Brune 483d728fb7dSPeter Brune Input Parameters: 484f6dfbefdSBarry Smith + snes - the SNES context 485f6dfbefdSBarry Smith - flg - `PETSC_TRUE` to compute the Jacobians 486d728fb7dSPeter Brune 487d728fb7dSPeter Brune Level: developer 488d728fb7dSPeter Brune 48995452b02SPatrick Sanan Notes: 490f6dfbefdSBarry Smith This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian 491d728fb7dSPeter Brune is needed at each linear iteration. 492d728fb7dSPeter Brune 493db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSubdomains()` 494d728fb7dSPeter Brune @*/ 495d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg) 496d71ae5a4SJacob Faibussowitsch { 497d728fb7dSPeter Brune PetscErrorCode (*f)(SNES, PetscBool); 498d728fb7dSPeter Brune 499d728fb7dSPeter Brune PetscFunctionBegin; 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f)); 5019566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg)); 5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 503d728fb7dSPeter Brune } 504d728fb7dSPeter Brune 505d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg) 506d71ae5a4SJacob Faibussowitsch { 507d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 508d728fb7dSPeter Brune 509d728fb7dSPeter Brune PetscFunctionBegin; 510d728fb7dSPeter Brune nasm->finaljacobian = flg; 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 512d728fb7dSPeter Brune } 51376857b2aSPeter Brune 514610116beSPeter Brune /*@ 515f6dfbefdSBarry Smith SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 516610116beSPeter Brune 517*20f4b53cSBarry Smith Logically Collective 518610116beSPeter Brune 519610116beSPeter Brune Input Parameters: 520f6dfbefdSBarry Smith + snes - the `SNES` context 521610116beSPeter Brune - dmp - damping 522610116beSPeter Brune 523610116beSPeter Brune Level: intermediate 524610116beSPeter Brune 52595452b02SPatrick Sanan Notes: 52695452b02SPatrick Sanan The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5275dfa0f3bSBarry Smith 528db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetDamping()` 529610116beSPeter Brune @*/ 530d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp) 531d71ae5a4SJacob Faibussowitsch { 532610116beSPeter Brune PetscErrorCode (*f)(SNES, PetscReal); 533610116beSPeter Brune 534610116beSPeter Brune PetscFunctionBegin; 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f)); 5369566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, dmp)); 5373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 538610116beSPeter Brune } 539610116beSPeter Brune 540d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp) 541d71ae5a4SJacob Faibussowitsch { 542610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 543610116beSPeter Brune 544610116beSPeter Brune PetscFunctionBegin; 545610116beSPeter Brune nasm->damping = dmp; 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 547610116beSPeter Brune } 548610116beSPeter Brune 549610116beSPeter Brune /*@ 550f6dfbefdSBarry Smith SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 551610116beSPeter Brune 552610116beSPeter Brune Not Collective 553610116beSPeter Brune 554f6dfbefdSBarry Smith Input Parameter: 555f6dfbefdSBarry Smith . snes - the SNES context 556f6dfbefdSBarry Smith 557f6dfbefdSBarry Smith Output Parameter: 558f6dfbefdSBarry Smith . dmp - damping 559610116beSPeter Brune 560610116beSPeter Brune Level: intermediate 561610116beSPeter Brune 562db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMSetDamping()` 563610116beSPeter Brune @*/ 564d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp) 565d71ae5a4SJacob Faibussowitsch { 566610116beSPeter Brune PetscFunctionBegin; 567cac4c232SBarry Smith PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp)); 5683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 569610116beSPeter Brune } 570610116beSPeter Brune 571d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp) 572d71ae5a4SJacob Faibussowitsch { 573610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 574610116beSPeter Brune 575610116beSPeter Brune PetscFunctionBegin; 576610116beSPeter Brune *dmp = nasm->damping; 5773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 578610116beSPeter Brune } 579610116beSPeter Brune 58014eb1c5cSMatthew G. Knepley /* 58114eb1c5cSMatthew G. Knepley Input Parameters: 58214eb1c5cSMatthew G. Knepley + snes - The solver 58314eb1c5cSMatthew G. Knepley . B - The RHS vector 58414eb1c5cSMatthew G. Knepley - X - The initial guess 58514eb1c5cSMatthew G. Knepley 58614eb1c5cSMatthew G. Knepley Output Parameters: 58714eb1c5cSMatthew G. Knepley . Y - The solution update 58814eb1c5cSMatthew G. Knepley 58914eb1c5cSMatthew G. Knepley TODO: All scatters should be packed into one 59014eb1c5cSMatthew G. Knepley */ 591d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X) 592d71ae5a4SJacob Faibussowitsch { 593eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 594258e1594SPeter Brune SNES subsnes; 595eaedb033SPeter Brune PetscInt i; 596610116beSPeter Brune PetscReal dmp; 597f10b3e88SPatrick Farrell Vec Xl, Bl, Yl, Xlloc; 598f10b3e88SPatrick Farrell VecScatter iscat, oscat, gscat, oscat_copy; 599111ade9eSPeter Brune DM dm, subdm; 600e0331734SPeter Brune PCASMType type; 6010adebc6cSBarry Smith 602eaedb033SPeter Brune PetscFunctionBegin; 6039566063dSJacob Faibussowitsch PetscCall(SNESNASMGetType(snes, &type)); 6049566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 6059566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0)); 6069566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 607eaedb033SPeter Brune for (i = 0; i < nasm->n; i++) { 608f10b3e88SPatrick Farrell /* scatter the solution to the global solution and the local solution */ 609f10b3e88SPatrick Farrell Xl = nasm->x[i]; 61070c78f05SPeter Brune Xlloc = nasm->xl[i]; 61170c78f05SPeter Brune oscat = nasm->oscatter[i]; 612f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 613f10b3e88SPatrick Farrell gscat = nasm->gscatter[i]; 6149566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 6159566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 61670c78f05SPeter Brune if (B) { 61770c78f05SPeter Brune /* scatter the RHS to the local RHS */ 61870c78f05SPeter Brune Bl = nasm->b[i]; 6199566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 62070c78f05SPeter Brune } 62170c78f05SPeter Brune } 6229566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 623b20c023fSPeter Brune 6249566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0)); 62570c78f05SPeter Brune for (i = 0; i < nasm->n; i++) { 62670c78f05SPeter Brune Xl = nasm->x[i]; 62770c78f05SPeter Brune Xlloc = nasm->xl[i]; 62870c78f05SPeter Brune Yl = nasm->y[i]; 629258e1594SPeter Brune subsnes = nasm->subsnes[i]; 6309566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes, &subdm)); 631111ade9eSPeter Brune iscat = nasm->iscatter[i]; 632111ade9eSPeter Brune oscat = nasm->oscatter[i]; 633f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 634111ade9eSPeter Brune gscat = nasm->gscatter[i]; 6359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 6369566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 63724b7f281SPeter Brune if (B) { 63824b7f281SPeter Brune Bl = nasm->b[i]; 6399566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 640ed3c11a9SPeter Brune } else Bl = NULL; 641f10b3e88SPatrick Farrell 6429566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 6439566063dSJacob Faibussowitsch PetscCall(VecCopy(Xl, Yl)); 6449566063dSJacob Faibussowitsch PetscCall(SNESSolve(subsnes, Bl, Xl)); 6459566063dSJacob Faibussowitsch PetscCall(VecAYPX(Yl, -1.0, Xl)); 6469566063dSJacob Faibussowitsch PetscCall(VecScale(Yl, nasm->damping)); 647e0331734SPeter Brune if (type == PC_ASM_BASIC) { 6489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 6499566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 650e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 6519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 6529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 653ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM"); 654eaedb033SPeter Brune } 6559566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0)); 6569566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 6571baa6e33SBarry Smith if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight)); 6589566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 6599566063dSJacob Faibussowitsch PetscCall(SNESNASMGetDamping(snes, &dmp)); 6609566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, dmp, Y)); 6613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 662eaedb033SPeter Brune } 663eaedb033SPeter Brune 664d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 665d71ae5a4SJacob Faibussowitsch { 666602bec5dSPeter Brune Vec X = Xfinal; 667d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 668d728fb7dSPeter Brune SNES subsnes; 669ca44f815SPeter Brune PetscInt i, lag = 1; 670e59f0a30SPeter Brune Vec Xlloc, Xl, Fl, F; 671d728fb7dSPeter Brune VecScatter oscat, gscat; 672d728fb7dSPeter Brune DM dm, subdm; 673d1e9a80fSBarry Smith 674d728fb7dSPeter Brune PetscFunctionBegin; 675602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 676e59f0a30SPeter Brune F = snes->vec_func; 6779566063dSJacob Faibussowitsch if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F)); 6789566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 6799566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 6809566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 681602bec5dSPeter Brune if (nasm->fjtype != 1) { 682d728fb7dSPeter Brune for (i = 0; i < nasm->n; i++) { 683d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 684d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 6859566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 686602bec5dSPeter Brune } 687d728fb7dSPeter Brune } 6889566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 689d728fb7dSPeter Brune for (i = 0; i < nasm->n; i++) { 690e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 691d728fb7dSPeter Brune Xl = nasm->x[i]; 692d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 693d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 694d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 695d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 6969566063dSJacob Faibussowitsch if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 6979566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes, &subdm)); 6989566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 699602bec5dSPeter Brune if (nasm->fjtype != 1) { 7009566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl)); 7019566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl)); 702602bec5dSPeter Brune } 703ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 704ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 7059566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(subsnes, Xl, Fl)); 7069566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre)); 707ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 708d728fb7dSPeter Brune } 7093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 710d728fb7dSPeter Brune } 711d728fb7dSPeter Brune 712d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NASM(SNES snes) 713d71ae5a4SJacob Faibussowitsch { 714eaedb033SPeter Brune Vec F; 715eaedb033SPeter Brune Vec X; 716eaedb033SPeter Brune Vec B; 717111ade9eSPeter Brune Vec Y; 718eaedb033SPeter Brune PetscInt i; 719ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 720365a6726SPeter Brune SNESNormSchedule normschedule; 721d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 722eaedb033SPeter Brune 723eaedb033SPeter Brune PetscFunctionBegin; 724c579b300SPatrick Farrell 7250b121fc5SBarry 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); 726c579b300SPatrick Farrell 7279566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 728eaedb033SPeter Brune X = snes->vec_sol; 729111ade9eSPeter Brune Y = snes->vec_sol_update; 730eaedb033SPeter Brune F = snes->vec_func; 731eaedb033SPeter Brune B = snes->vec_rhs; 732eaedb033SPeter Brune 7339566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 734eaedb033SPeter Brune snes->iter = 0; 735eaedb033SPeter Brune snes->norm = 0.; 7369566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 737eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 7389566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 739365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 740eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 741eaedb033SPeter Brune if (!snes->vec_func_init_set) { 7429566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 7431aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 744eaedb033SPeter Brune 7459566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 746422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 7479566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 748eaedb033SPeter Brune snes->iter = 0; 749eaedb033SPeter Brune snes->norm = fnorm; 7509566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7519566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7529566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, snes->norm)); 753eaedb033SPeter Brune 754eaedb033SPeter Brune /* test convergence */ 755dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 7563ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 757eaedb033SPeter Brune } else { 7589566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7599566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7609566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, 0, snes->norm)); 761eaedb033SPeter Brune } 762eaedb033SPeter Brune 763eaedb033SPeter Brune /* Call general purpose update function */ 764dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 765602bec5dSPeter Brune /* copy the initial solution over for later */ 7669566063dSJacob Faibussowitsch if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit)); 767eaedb033SPeter Brune 768eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 7699566063dSJacob Faibussowitsch PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X)); 770365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 7719566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 7729566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 773422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 774eaedb033SPeter Brune } 775eaedb033SPeter Brune /* Monitor convergence */ 7769566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 777eaedb033SPeter Brune snes->iter = i + 1; 778eaedb033SPeter Brune snes->norm = fnorm; 7799566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7809566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7819566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 782eaedb033SPeter Brune /* Test for convergence */ 783dbbe0bcdSBarry Smith if (normschedule == SNES_NORM_ALWAYS) PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 784d728fb7dSPeter Brune if (snes->reason) break; 785eaedb033SPeter Brune /* Call general purpose update function */ 786dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 787eaedb033SPeter Brune } 78807b62357SFande Kong if (nasm->finaljacobian) { 7899566063dSJacob Faibussowitsch PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X)); 79007b62357SFande Kong SNESCheckJacobianDomainerror(snes); 79107b62357SFande Kong } 792365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 793eaedb033SPeter Brune if (i == snes->max_its) { 79463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", snes->max_its)); 795eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 796eaedb033SPeter Brune } 797f6dfbefdSBarry Smith } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a nonlinear preconditioner */ 7983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 799eaedb033SPeter Brune } 800eaedb033SPeter Brune 801eaedb033SPeter Brune /*MC 802f6dfbefdSBarry Smith SNESNASM - Nonlinear Additive Schwarz solver 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 8075dfa0f3bSBarry Smith . -snes_asm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 80870c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 809150d49b7SHong Zhang . -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 8294f02bc6aSBarry Smith References: 830606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8314f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8324f02bc6aSBarry Smith 833db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`, `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()` 834eaedb033SPeter Brune M*/ 835eaedb033SPeter Brune 836d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 837d71ae5a4SJacob Faibussowitsch { 838eaedb033SPeter Brune SNES_NASM *nasm; 839eaedb033SPeter Brune 840eaedb033SPeter Brune PetscFunctionBegin; 8414dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&nasm)); 842eaedb033SPeter Brune snes->data = (void *)nasm; 843eaedb033SPeter Brune 844eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 8459e5d0892SLisandro Dalcin nasm->subsnes = NULL; 8469e5d0892SLisandro Dalcin nasm->x = NULL; 8479e5d0892SLisandro Dalcin nasm->xl = NULL; 8489e5d0892SLisandro Dalcin nasm->y = NULL; 8499e5d0892SLisandro Dalcin nasm->b = NULL; 8509e5d0892SLisandro Dalcin nasm->oscatter = NULL; 8519e5d0892SLisandro Dalcin nasm->oscatter_copy = NULL; 8529e5d0892SLisandro Dalcin nasm->iscatter = NULL; 8539e5d0892SLisandro Dalcin nasm->gscatter = NULL; 854610116beSPeter Brune nasm->damping = 1.; 855111ade9eSPeter Brune 856111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 857d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 858f10b3e88SPatrick Farrell nasm->weight_set = PETSC_FALSE; 859eaedb033SPeter Brune 860eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 861eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 862eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 863eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 864eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 865eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 866eaedb033SPeter Brune 867eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 868efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 869eaedb033SPeter Brune 8704fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8714fc747eaSLawrence Mitchell 872602bec5dSPeter Brune nasm->fjtype = 0; 873602bec5dSPeter Brune nasm->xinit = NULL; 8740298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8750298fd71SBarry Smith nasm->eventsubsolve = 0; 876b20c023fSPeter Brune 877eaedb033SPeter Brune if (!snes->tolerancesset) { 878eaedb033SPeter Brune snes->max_its = 10000; 879eaedb033SPeter Brune snes->max_funcs = 10000; 880eaedb033SPeter Brune } 881eaedb033SPeter Brune 8829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM)); 8839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM)); 8849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM)); 8859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM)); 8869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM)); 8879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM)); 8889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM)); 8899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM)); 8903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 891eaedb033SPeter Brune } 89299e0435eSBarry Smith 893448b6425SPatrick Farrell /*@ 894448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 895448b6425SPatrick Farrell 896*20f4b53cSBarry Smith Not Collective 897448b6425SPatrick Farrell 898448b6425SPatrick Farrell Input Parameters: 899*20f4b53cSBarry Smith + snes - the `SNES` context 900448b6425SPatrick Farrell - i - the number of the subsnes to get 901448b6425SPatrick Farrell 902448b6425SPatrick Farrell Output Parameters: 903448b6425SPatrick Farrell . subsnes - the subsolver context 904448b6425SPatrick Farrell 905448b6425SPatrick Farrell Level: intermediate 906448b6425SPatrick Farrell 907db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetNumber()` 908448b6425SPatrick Farrell @*/ 909d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes) 910d71ae5a4SJacob Faibussowitsch { 911448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 912448b6425SPatrick Farrell 913448b6425SPatrick Farrell PetscFunctionBegin; 9140b121fc5SBarry Smith PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver"); 915448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 9163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 917448b6425SPatrick Farrell } 918448b6425SPatrick Farrell 919448b6425SPatrick Farrell /*@ 920448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 921448b6425SPatrick Farrell 922*20f4b53cSBarry Smith Not Collective 923448b6425SPatrick Farrell 924448b6425SPatrick Farrell Input Parameters: 925*20f4b53cSBarry Smith . snes - the `SNES` context 926448b6425SPatrick Farrell 927448b6425SPatrick Farrell Output Parameters: 928448b6425SPatrick Farrell . n - the number of subsolvers 929448b6425SPatrick Farrell 930448b6425SPatrick Farrell Level: intermediate 931448b6425SPatrick Farrell 932db781477SPatrick Sanan .seealso: `SNESNASM`, `SNESNASMGetSNES()` 933448b6425SPatrick Farrell @*/ 934d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n) 935d71ae5a4SJacob Faibussowitsch { 936448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 937448b6425SPatrick Farrell 938448b6425SPatrick Farrell PetscFunctionBegin; 939448b6425SPatrick Farrell *n = nasm->n; 9403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 941448b6425SPatrick Farrell } 942448b6425SPatrick Farrell 943f10b3e88SPatrick Farrell /*@ 944f10b3e88SPatrick Farrell SNESNASMSetWeight - Sets weight to use when adding overlapping updates 945f10b3e88SPatrick Farrell 946f10b3e88SPatrick Farrell Collective 947f10b3e88SPatrick Farrell 948f10b3e88SPatrick Farrell Input Parameters: 949*20f4b53cSBarry Smith + snes - the `SNES` context 950f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 951f10b3e88SPatrick Farrell 952f10b3e88SPatrick Farrell Level: intermediate 953f10b3e88SPatrick Farrell 954db781477SPatrick Sanan .seealso: `SNESNASM` 955f10b3e88SPatrick Farrell @*/ 956d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight) 957d71ae5a4SJacob Faibussowitsch { 958f10b3e88SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 959f10b3e88SPatrick Farrell 960f10b3e88SPatrick Farrell PetscFunctionBegin; 961f10b3e88SPatrick Farrell 9629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nasm->weight)); 963f10b3e88SPatrick Farrell nasm->weight_set = PETSC_TRUE; 964f10b3e88SPatrick Farrell nasm->weight = weight; 9659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)nasm->weight)); 966f10b3e88SPatrick Farrell 9673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 968f10b3e88SPatrick Farrell } 969