1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 324cdb843SMatthew G. Knepley #include <petscds.h> 4af0996ceSBarry Smith #include <petsc/private/petscimpl.h> 5af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h> 6552f7358SJed Brown 7d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 8d2b2dc1eSMatthew G. Knepley #include <petscdmceed.h> 9d2b2dc1eSMatthew G. Knepley #include <petscdmplexceed.h> 10d2b2dc1eSMatthew G. Knepley #endif 11d2b2dc1eSMatthew G. Knepley 12d71ae5a4SJacob Faibussowitsch static void pressure_Private(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar p[]) 13d71ae5a4SJacob Faibussowitsch { 14649ef022SMatthew Knepley p[0] = u[uOff[1]]; 15649ef022SMatthew Knepley } 16649ef022SMatthew Knepley 17649ef022SMatthew Knepley /* 1820f4b53cSBarry Smith SNESCorrectDiscretePressure_Private - Add a vector in the nullspace to make the continuum integral of the pressure field equal to zero. 1920f4b53cSBarry Smith This is normally used only to evaluate convergence rates for the pressure accurately. 20649ef022SMatthew Knepley 21c3339decSBarry Smith Collective 22649ef022SMatthew Knepley 23649ef022SMatthew Knepley Input Parameters: 24420bcc1bSBarry Smith + snes - The `SNES` 25649ef022SMatthew Knepley . pfield - The field number for pressure 26649ef022SMatthew Knepley . nullspace - The pressure nullspace 27649ef022SMatthew Knepley . u - The solution vector 28649ef022SMatthew Knepley - ctx - An optional user context 29649ef022SMatthew Knepley 30649ef022SMatthew Knepley Output Parameter: 31649ef022SMatthew Knepley . u - The solution with a continuum pressure integral of zero 32649ef022SMatthew Knepley 332fe279fdSBarry Smith Level: developer 342fe279fdSBarry Smith 35420bcc1bSBarry Smith Note: 36649ef022SMatthew Knepley If int(u) = a and int(n) = b, then int(u - a/b n) = a - a/b b = 0. We assume that the nullspace is a single vector given explicitly. 37649ef022SMatthew Knepley 38420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESConvergedCorrectPressure()` 39649ef022SMatthew Knepley */ 40d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESCorrectDiscretePressure_Private(SNES snes, PetscInt pfield, MatNullSpace nullspace, Vec u, void *ctx) 41d71ae5a4SJacob Faibussowitsch { 42649ef022SMatthew Knepley DM dm; 43649ef022SMatthew Knepley PetscDS ds; 44649ef022SMatthew Knepley const Vec *nullvecs; 45649ef022SMatthew Knepley PetscScalar pintd, *intc, *intn; 46649ef022SMatthew Knepley MPI_Comm comm; 47649ef022SMatthew Knepley PetscInt Nf, Nv; 48649ef022SMatthew Knepley 49649ef022SMatthew Knepley PetscFunctionBegin; 509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 519566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 5228b400f6SJacob Faibussowitsch PetscCheck(dm, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a SNES DM"); 5328b400f6SJacob Faibussowitsch PetscCheck(nullspace, comm, PETSC_ERR_ARG_WRONG, "Cannot compute test without a Jacobian nullspace"); 549566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 559566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, pfield, pressure_Private)); 569566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(nullspace, NULL, &Nv, &nullvecs)); 5763a3b9bcSJacob Faibussowitsch PetscCheck(Nv == 1, comm, PETSC_ERR_ARG_OUTOFRANGE, "Can only handle a single null vector for pressure, not %" PetscInt_FMT, Nv); 589566063dSJacob Faibussowitsch PetscCall(VecDot(nullvecs[0], u, &pintd)); 5908401ef6SPierre Jolivet PetscCheck(PetscAbsScalar(pintd) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Discrete integral of pressure: %g", (double)PetscRealPart(pintd)); 609566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 619566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nf, &intc, Nf, &intn)); 629566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, nullvecs[0], intn, ctx)); 639566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx)); 649566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -intc[pfield] / intn[pfield], nullvecs[0])); 65649ef022SMatthew Knepley #if defined(PETSC_USE_DEBUG) 669566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, intc, ctx)); 6708401ef6SPierre Jolivet PetscCheck(PetscAbsScalar(intc[pfield]) <= PETSC_SMALL, comm, PETSC_ERR_ARG_WRONG, "Continuum integral of pressure after correction: %g", (double)PetscRealPart(intc[pfield])); 68649ef022SMatthew Knepley #endif 699566063dSJacob Faibussowitsch PetscCall(PetscFree2(intc, intn)); 703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71649ef022SMatthew Knepley } 72649ef022SMatthew Knepley 73649ef022SMatthew Knepley /*@C 74ceaaa498SBarry Smith SNESConvergedCorrectPressure - The regular `SNES` convergence test that, up on convergence, adds a vector in the nullspace 75ceaaa498SBarry Smith to make the continuum integral of the pressure field equal to zero. 76649ef022SMatthew Knepley 77c3339decSBarry Smith Logically Collective 78649ef022SMatthew Knepley 79649ef022SMatthew Knepley Input Parameters: 80f6dfbefdSBarry Smith + snes - the `SNES` context 81649ef022SMatthew Knepley . it - the iteration (0 indicates before any Newton steps) 82649ef022SMatthew Knepley . xnorm - 2-norm of current iterate 83e4094ef1SJacob Faibussowitsch . gnorm - 2-norm of current step 84e4094ef1SJacob Faibussowitsch . f - 2-norm of function at current iterate 85649ef022SMatthew Knepley - ctx - Optional user context 86649ef022SMatthew Knepley 87649ef022SMatthew Knepley Output Parameter: 88f6dfbefdSBarry Smith . reason - `SNES_CONVERGED_ITERATING`, `SNES_CONVERGED_ITS`, or `SNES_DIVERGED_FNORM_NAN` 89649ef022SMatthew Knepley 9020f4b53cSBarry Smith Options Database Key: 9120f4b53cSBarry Smith . -snes_convergence_test correct_pressure - see `SNESSetFromOptions()` 9220f4b53cSBarry Smith 9320f4b53cSBarry Smith Level: advanced 9420f4b53cSBarry Smith 95649ef022SMatthew Knepley Notes: 96f6dfbefdSBarry Smith In order to use this convergence test, you must set up several PETSc structures. First fields must be added to the `DM`, and a `PetscDS` 97f6dfbefdSBarry Smith must be created with discretizations of those fields. We currently assume that the pressure field has index 1. 98f6dfbefdSBarry Smith The pressure field must have a nullspace, likely created using the `DMSetNullSpaceConstructor()` interface. 99f6dfbefdSBarry Smith Last we must be able to integrate the pressure over the domain, so the `DM` attached to the SNES `must` be a `DMPLEX` at this time. 100f362779dSJed Brown 101ceaaa498SBarry Smith Developer Note: 102ceaaa498SBarry Smith This is a total misuse of the `SNES` convergence test handling system. It should be removed. Perhaps a `SNESSetPostSolve()` could 103ceaaa498SBarry Smith be constructed to handle this process. 104ceaaa498SBarry Smith 105420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DM`, `SNESConvergedDefault()`, `SNESSetConvergenceTest()`, `DMSetNullSpaceConstructor()` 106649ef022SMatthew Knepley @*/ 107d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESConvergedCorrectPressure(SNES snes, PetscInt it, PetscReal xnorm, PetscReal gnorm, PetscReal f, SNESConvergedReason *reason, void *ctx) 108d71ae5a4SJacob Faibussowitsch { 109649ef022SMatthew Knepley PetscBool monitorIntegral = PETSC_FALSE; 110649ef022SMatthew Knepley 111649ef022SMatthew Knepley PetscFunctionBegin; 1129566063dSJacob Faibussowitsch PetscCall(SNESConvergedDefault(snes, it, xnorm, gnorm, f, reason, ctx)); 113649ef022SMatthew Knepley if (monitorIntegral) { 114649ef022SMatthew Knepley Mat J; 115649ef022SMatthew Knepley Vec u; 116649ef022SMatthew Knepley MatNullSpace nullspace; 117649ef022SMatthew Knepley const Vec *nullvecs; 118649ef022SMatthew Knepley PetscScalar pintd; 119649ef022SMatthew Knepley 1209566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 1219566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 1229566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(J, &nullspace)); 1239566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(nullspace, NULL, NULL, &nullvecs)); 1249566063dSJacob Faibussowitsch PetscCall(VecDot(nullvecs[0], u, &pintd)); 1259566063dSJacob Faibussowitsch PetscCall(PetscInfo(snes, "SNES: Discrete integral of pressure: %g\n", (double)PetscRealPart(pintd))); 126649ef022SMatthew Knepley } 127649ef022SMatthew Knepley if (*reason > 0) { 128649ef022SMatthew Knepley Mat J; 129649ef022SMatthew Knepley Vec u; 130649ef022SMatthew Knepley MatNullSpace nullspace; 131649ef022SMatthew Knepley PetscInt pfield = 1; 132649ef022SMatthew Knepley 1339566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 1349566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 1359566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(J, &nullspace)); 1369566063dSJacob Faibussowitsch PetscCall(SNESCorrectDiscretePressure_Private(snes, pfield, nullspace, u, ctx)); 137649ef022SMatthew Knepley } 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 139649ef022SMatthew Knepley } 140649ef022SMatthew Knepley 141d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) 142d71ae5a4SJacob Faibussowitsch { 1436da023fcSToby Isaac PetscBool isPlex; 1446da023fcSToby Isaac 1456da023fcSToby Isaac PetscFunctionBegin; 1469566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex)); 1476da023fcSToby Isaac if (isPlex) { 1486da023fcSToby Isaac *plex = dm; 1499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 150f7148743SMatthew G. Knepley } else { 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex)); 152f7148743SMatthew G. Knepley if (!*plex) { 1539566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, plex)); 1549566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex)); 155cbf8eb3cSStefano Zampini } else { 156cbf8eb3cSStefano Zampini PetscCall(PetscObjectReference((PetscObject)*plex)); 157cbf8eb3cSStefano Zampini } 1586da023fcSToby Isaac if (copy) { 1599566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(dm, *plex)); 1609566063dSJacob Faibussowitsch PetscCall(DMCopyAuxiliaryVec(dm, *plex)); 1616da023fcSToby Isaac } 1626da023fcSToby Isaac } 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1646da023fcSToby Isaac } 1656da023fcSToby Isaac 1664267b1a3SMatthew G. Knepley /*@C 167cc0c4584SMatthew G. Knepley SNESMonitorFields - Monitors the residual for each field separately 168cc0c4584SMatthew G. Knepley 169c3339decSBarry Smith Collective 170cc0c4584SMatthew G. Knepley 171cc0c4584SMatthew G. Knepley Input Parameters: 172420bcc1bSBarry Smith + snes - the `SNES` context, must have an attached `DM` 173cc0c4584SMatthew G. Knepley . its - iteration number 174cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual 175f6dfbefdSBarry Smith - vf - `PetscViewerAndFormat` of `PetscViewerType` `PETSCVIEWERASCII` 176cc0c4584SMatthew G. Knepley 1772fe279fdSBarry Smith Level: intermediate 1782fe279fdSBarry Smith 179f6dfbefdSBarry Smith Note: 180cc0c4584SMatthew G. Knepley This routine prints the residual norm at each iteration. 181cc0c4584SMatthew G. Knepley 182420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()` 183cc0c4584SMatthew G. Knepley @*/ 184d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) 185d71ae5a4SJacob Faibussowitsch { 186d43b4f6eSBarry Smith PetscViewer viewer = vf->viewer; 187cc0c4584SMatthew G. Knepley Vec res; 188cc0c4584SMatthew G. Knepley DM dm; 189cc0c4584SMatthew G. Knepley PetscSection s; 190cc0c4584SMatthew G. Knepley const PetscScalar *r; 191cc0c4584SMatthew G. Knepley PetscReal *lnorms, *norms; 192cc0c4584SMatthew G. Knepley PetscInt numFields, f, pStart, pEnd, p; 193cc0c4584SMatthew G. Knepley 194cc0c4584SMatthew G. Knepley PetscFunctionBegin; 1954d4332d5SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 1969566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &res, NULL, NULL)); 1979566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1989566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 1999566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &numFields)); 2009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 2019566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms)); 2029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(res, &r)); 203cc0c4584SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 204cc0c4584SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 205cc0c4584SMatthew G. Knepley PetscInt fdof, foff, d; 206cc0c4584SMatthew G. Knepley 2079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 2089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff)); 209cc0c4584SMatthew G. Knepley for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d])); 210cc0c4584SMatthew G. Knepley } 211cc0c4584SMatthew G. Knepley } 2129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(res, &r)); 213e91c04dfSPierre Jolivet PetscCallMPI(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm))); 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 21663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm)); 217cc0c4584SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 2189566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 2199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f]))); 220cc0c4584SMatthew G. Knepley } 2219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "]\n")); 2229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 2239566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 2249566063dSJacob Faibussowitsch PetscCall(PetscFree2(lnorms, norms)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226cc0c4584SMatthew G. Knepley } 22724cdb843SMatthew G. Knepley 2286493148fSStefano Zampini /********************* SNES callbacks **************************/ 2296493148fSStefano Zampini 2306493148fSStefano Zampini /*@ 2316493148fSStefano Zampini DMPlexSNESComputeObjectiveFEM - Sums the local objectives from the local input X using pointwise functions specified by the user 2326493148fSStefano Zampini 2336493148fSStefano Zampini Input Parameters: 2346493148fSStefano Zampini + dm - The mesh 2356493148fSStefano Zampini . X - Local solution 2366493148fSStefano Zampini - user - The user context 2376493148fSStefano Zampini 2386493148fSStefano Zampini Output Parameter: 2396493148fSStefano Zampini . obj - Local objective value 2406493148fSStefano Zampini 2416493148fSStefano Zampini Level: developer 2426493148fSStefano Zampini 2436493148fSStefano Zampini .seealso: `DM`, `DMPlexSNESComputeResidualFEM()` 2446493148fSStefano Zampini @*/ 2456493148fSStefano Zampini PetscErrorCode DMPlexSNESComputeObjectiveFEM(DM dm, Vec X, PetscReal *obj, void *user) 2466493148fSStefano Zampini { 2476493148fSStefano Zampini PetscInt Nf, cellHeight, cStart, cEnd; 2486493148fSStefano Zampini PetscScalar *cintegral; 2496493148fSStefano Zampini 2506493148fSStefano Zampini PetscFunctionBegin; 2516493148fSStefano Zampini PetscCall(DMGetNumFields(dm, &Nf)); 2526493148fSStefano Zampini PetscCall(DMPlexGetVTKCellHeight(dm, &cellHeight)); 2536493148fSStefano Zampini PetscCall(DMPlexGetSimplexOrBoxCells(dm, cellHeight, &cStart, &cEnd)); 2546493148fSStefano Zampini PetscCall(PetscCalloc1((cEnd - cStart) * Nf, &cintegral)); 2556493148fSStefano Zampini PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0)); 2566493148fSStefano Zampini PetscCall(DMPlexComputeIntegral_Internal(dm, X, cStart, cEnd, cintegral, user)); 2576493148fSStefano Zampini /* Sum up values */ 2586493148fSStefano Zampini *obj = 0; 2596493148fSStefano Zampini for (PetscInt c = cStart; c < cEnd; ++c) 2606493148fSStefano Zampini for (PetscInt f = 0; f < Nf; ++f) *obj += PetscRealPart(cintegral[(c - cStart) * Nf + f]); 2616493148fSStefano Zampini PetscCall(PetscLogEventBegin(DMPLEX_IntegralFEM, dm, 0, 0, 0)); 2626493148fSStefano Zampini PetscCall(PetscFree(cintegral)); 2636493148fSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 2646493148fSStefano Zampini } 26524cdb843SMatthew G. Knepley 26624cdb843SMatthew G. Knepley /*@ 267420bcc1bSBarry Smith DMPlexSNESComputeResidualFEM - Sums the local residual into vector `F` from the local input `X` using pointwise functions specified by the user 26824cdb843SMatthew G. Knepley 26924cdb843SMatthew G. Knepley Input Parameters: 27024cdb843SMatthew G. Knepley + dm - The mesh 27124cdb843SMatthew G. Knepley . X - Local solution 27224cdb843SMatthew G. Knepley - user - The user context 27324cdb843SMatthew G. Knepley 27424cdb843SMatthew G. Knepley Output Parameter: 27524cdb843SMatthew G. Knepley . F - Local output vector 27624cdb843SMatthew G. Knepley 2772fe279fdSBarry Smith Level: developer 2782fe279fdSBarry Smith 279f6dfbefdSBarry Smith Note: 280420bcc1bSBarry Smith The residual is summed into `F`; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in `F` is intentional. 2818db23a53SJed Brown 2826493148fSStefano Zampini .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMSNESComputeJacobianAction()` 28324cdb843SMatthew G. Knepley @*/ 284d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 285d71ae5a4SJacob Faibussowitsch { 2866da023fcSToby Isaac DM plex; 287083401c6SMatthew G. Knepley IS allcellIS; 2886528b96dSMatthew G. Knepley PetscInt Nds, s; 28924cdb843SMatthew G. Knepley 29024cdb843SMatthew G. Knepley PetscFunctionBegin; 2919566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 2929566063dSJacob Faibussowitsch PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 2939566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 2946528b96dSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 2956528b96dSMatthew G. Knepley PetscDS ds; 2966528b96dSMatthew G. Knepley IS cellIS; 29706ad1575SMatthew G. Knepley PetscFormKey key; 2986528b96dSMatthew G. Knepley 29907218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL)); 3006528b96dSMatthew G. Knepley key.value = 0; 3016528b96dSMatthew G. Knepley key.field = 0; 30206ad1575SMatthew G. Knepley key.part = 0; 3036528b96dSMatthew G. Knepley if (!key.label) { 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)allcellIS)); 3056528b96dSMatthew G. Knepley cellIS = allcellIS; 3066528b96dSMatthew G. Knepley } else { 3076528b96dSMatthew G. Knepley IS pointIS; 3086528b96dSMatthew G. Knepley 3096528b96dSMatthew G. Knepley key.value = 1; 3109566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS)); 3119566063dSJacob Faibussowitsch PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 3129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 3136528b96dSMatthew G. Knepley } 314*754e4fbaSMatthew G. Knepley PetscCall(DMPlexComputeResidualByKey(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user)); 3159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 3166528b96dSMatthew G. Knepley } 3179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&allcellIS)); 3189566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 3193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3206528b96dSMatthew G. Knepley } 3216528b96dSMatthew G. Knepley 322bdd6f66aSToby Isaac /*@ 323420bcc1bSBarry Smith DMPlexSNESComputeResidualDS - Sums the local residual into vector `F` from the local input `X` using all pointwise functions with unique keys in the `PetscDS` 324d2b2dc1eSMatthew G. Knepley 325d2b2dc1eSMatthew G. Knepley Input Parameters: 326d2b2dc1eSMatthew G. Knepley + dm - The mesh 327d2b2dc1eSMatthew G. Knepley . X - Local solution 328d2b2dc1eSMatthew G. Knepley - user - The user context 329d2b2dc1eSMatthew G. Knepley 330d2b2dc1eSMatthew G. Knepley Output Parameter: 331d2b2dc1eSMatthew G. Knepley . F - Local output vector 332d2b2dc1eSMatthew G. Knepley 333d2b2dc1eSMatthew G. Knepley Level: developer 334d2b2dc1eSMatthew G. Knepley 335d2b2dc1eSMatthew G. Knepley Note: 336420bcc1bSBarry Smith The residual is summed into `F`; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in `F` is intentional. 337d2b2dc1eSMatthew G. Knepley 338420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()` 339d2b2dc1eSMatthew G. Knepley @*/ 340d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualDS(DM dm, Vec X, Vec F, void *user) 341d2b2dc1eSMatthew G. Knepley { 342d2b2dc1eSMatthew G. Knepley DM plex; 343d2b2dc1eSMatthew G. Knepley IS allcellIS; 344d2b2dc1eSMatthew G. Knepley PetscInt Nds, s; 345d2b2dc1eSMatthew G. Knepley 346d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 347d2b2dc1eSMatthew G. Knepley PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 348d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 349d2b2dc1eSMatthew G. Knepley PetscCall(DMGetNumDS(dm, &Nds)); 350d2b2dc1eSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 351d2b2dc1eSMatthew G. Knepley PetscDS ds; 352d2b2dc1eSMatthew G. Knepley DMLabel label; 353d2b2dc1eSMatthew G. Knepley IS cellIS; 354d2b2dc1eSMatthew G. Knepley 355d2b2dc1eSMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL)); 356d2b2dc1eSMatthew G. Knepley { 357d2b2dc1eSMatthew G. Knepley PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1}; 358d2b2dc1eSMatthew G. Knepley PetscWeakForm wf; 359d2b2dc1eSMatthew G. Knepley PetscInt Nm = 2, m, Nk = 0, k, kp, off = 0; 360d2b2dc1eSMatthew G. Knepley PetscFormKey *reskeys; 361d2b2dc1eSMatthew G. Knepley 362d2b2dc1eSMatthew G. Knepley /* Get unique residual keys */ 363d2b2dc1eSMatthew G. Knepley for (m = 0; m < Nm; ++m) { 364d2b2dc1eSMatthew G. Knepley PetscInt Nkm; 365d2b2dc1eSMatthew G. Knepley PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm)); 366d2b2dc1eSMatthew G. Knepley Nk += Nkm; 367d2b2dc1eSMatthew G. Knepley } 368d2b2dc1eSMatthew G. Knepley PetscCall(PetscMalloc1(Nk, &reskeys)); 369d2b2dc1eSMatthew G. Knepley for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys)); 370d2b2dc1eSMatthew G. Knepley PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk); 371d2b2dc1eSMatthew G. Knepley PetscCall(PetscFormKeySort(Nk, reskeys)); 372d2b2dc1eSMatthew G. Knepley for (k = 0, kp = 1; kp < Nk; ++kp) { 373d2b2dc1eSMatthew G. Knepley if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) { 374d2b2dc1eSMatthew G. Knepley ++k; 375d2b2dc1eSMatthew G. Knepley if (kp != k) reskeys[k] = reskeys[kp]; 376d2b2dc1eSMatthew G. Knepley } 377d2b2dc1eSMatthew G. Knepley } 378d2b2dc1eSMatthew G. Knepley Nk = k; 379d2b2dc1eSMatthew G. Knepley 380d2b2dc1eSMatthew G. Knepley PetscCall(PetscDSGetWeakForm(ds, &wf)); 381d2b2dc1eSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 382d2b2dc1eSMatthew G. Knepley DMLabel label = reskeys[k].label; 383d2b2dc1eSMatthew G. Knepley PetscInt val = reskeys[k].value; 384d2b2dc1eSMatthew G. Knepley 385d2b2dc1eSMatthew G. Knepley if (!label) { 386d2b2dc1eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)allcellIS)); 387d2b2dc1eSMatthew G. Knepley cellIS = allcellIS; 388d2b2dc1eSMatthew G. Knepley } else { 389d2b2dc1eSMatthew G. Knepley IS pointIS; 390d2b2dc1eSMatthew G. Knepley 391d2b2dc1eSMatthew G. Knepley PetscCall(DMLabelGetStratumIS(label, val, &pointIS)); 392d2b2dc1eSMatthew G. Knepley PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 393d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&pointIS)); 394d2b2dc1eSMatthew G. Knepley } 395*754e4fbaSMatthew G. Knepley PetscCall(DMPlexComputeResidualByKey(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user)); 396d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&cellIS)); 397d2b2dc1eSMatthew G. Knepley } 398d2b2dc1eSMatthew G. Knepley PetscCall(PetscFree(reskeys)); 399d2b2dc1eSMatthew G. Knepley } 400d2b2dc1eSMatthew G. Knepley } 401d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&allcellIS)); 402d2b2dc1eSMatthew G. Knepley PetscCall(DMDestroy(&plex)); 403d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 404d2b2dc1eSMatthew G. Knepley } 405d2b2dc1eSMatthew G. Knepley 406d2b2dc1eSMatthew G. Knepley /*@ 407420bcc1bSBarry Smith DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input `X` 408bdd6f66aSToby Isaac 409bdd6f66aSToby Isaac Input Parameters: 410bdd6f66aSToby Isaac + dm - The mesh 411bdd6f66aSToby Isaac - user - The user context 412bdd6f66aSToby Isaac 413bdd6f66aSToby Isaac Output Parameter: 414bdd6f66aSToby Isaac . X - Local solution 415bdd6f66aSToby Isaac 416bdd6f66aSToby Isaac Level: developer 417bdd6f66aSToby Isaac 418420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()` 419bdd6f66aSToby Isaac @*/ 420d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) 421d71ae5a4SJacob Faibussowitsch { 422bdd6f66aSToby Isaac DM plex; 423bdd6f66aSToby Isaac 424bdd6f66aSToby Isaac PetscFunctionBegin; 4259566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 4269566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL)); 4279566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 429bdd6f66aSToby Isaac } 430bdd6f66aSToby Isaac 4317a73cf09SMatthew G. Knepley /*@ 432c118d280SBarry Smith DMSNESComputeJacobianAction - Compute the action of the Jacobian J(`X`) on `Y` 4337a73cf09SMatthew G. Knepley 4347a73cf09SMatthew G. Knepley Input Parameters: 435f6dfbefdSBarry Smith + dm - The `DM` 4367a73cf09SMatthew G. Knepley . X - Local solution vector 4377a73cf09SMatthew G. Knepley . Y - Local input vector 4387a73cf09SMatthew G. Knepley - user - The user context 4397a73cf09SMatthew G. Knepley 4407a73cf09SMatthew G. Knepley Output Parameter: 44169d47153SPierre Jolivet . F - local output vector 4427a73cf09SMatthew G. Knepley 4437a73cf09SMatthew G. Knepley Level: developer 4447a73cf09SMatthew G. Knepley 445c118d280SBarry Smith Note: 446f6dfbefdSBarry Smith Users will typically use `DMSNESCreateJacobianMF()` followed by `MatMult()` instead of calling this routine directly. 4478e3a2eefSMatthew G. Knepley 448c118d280SBarry Smith This only works with `DMPLEX` 449c118d280SBarry Smith 450c118d280SBarry Smith Developer Note: 451c118d280SBarry Smith This should be called `DMPlexSNESComputeJacobianAction()` 452c118d280SBarry Smith 453a94f484eSPierre Jolivet .seealso: [](ch_snes), `DM`, `DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()` 4547a73cf09SMatthew G. Knepley @*/ 455d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user) 456d71ae5a4SJacob Faibussowitsch { 4578e3a2eefSMatthew G. Knepley DM plex; 4588e3a2eefSMatthew G. Knepley IS allcellIS; 4598e3a2eefSMatthew G. Knepley PetscInt Nds, s; 460a925c78cSMatthew G. Knepley 461a925c78cSMatthew G. Knepley PetscFunctionBegin; 4629566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 4639566063dSJacob Faibussowitsch PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 4649566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 4658e3a2eefSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 4668e3a2eefSMatthew G. Knepley PetscDS ds; 4678e3a2eefSMatthew G. Knepley DMLabel label; 4688e3a2eefSMatthew G. Knepley IS cellIS; 4697a73cf09SMatthew G. Knepley 47007218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL)); 4718e3a2eefSMatthew G. Knepley { 47206ad1575SMatthew G. Knepley PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3}; 4738e3a2eefSMatthew G. Knepley PetscWeakForm wf; 4748e3a2eefSMatthew G. Knepley PetscInt Nm = 4, m, Nk = 0, k, kp, off = 0; 47506ad1575SMatthew G. Knepley PetscFormKey *jackeys; 4768e3a2eefSMatthew G. Knepley 4778e3a2eefSMatthew G. Knepley /* Get unique Jacobian keys */ 4788e3a2eefSMatthew G. Knepley for (m = 0; m < Nm; ++m) { 4798e3a2eefSMatthew G. Knepley PetscInt Nkm; 4809566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm)); 4818e3a2eefSMatthew G. Knepley Nk += Nkm; 4828e3a2eefSMatthew G. Knepley } 4839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk, &jackeys)); 48448a46eb9SPierre Jolivet for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys)); 48563a3b9bcSJacob Faibussowitsch PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk); 4869566063dSJacob Faibussowitsch PetscCall(PetscFormKeySort(Nk, jackeys)); 4878e3a2eefSMatthew G. Knepley for (k = 0, kp = 1; kp < Nk; ++kp) { 4888e3a2eefSMatthew G. Knepley if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) { 4898e3a2eefSMatthew G. Knepley ++k; 4908e3a2eefSMatthew G. Knepley if (kp != k) jackeys[k] = jackeys[kp]; 4918e3a2eefSMatthew G. Knepley } 4928e3a2eefSMatthew G. Knepley } 4938e3a2eefSMatthew G. Knepley Nk = k; 4948e3a2eefSMatthew G. Knepley 4959566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 4968e3a2eefSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 4978e3a2eefSMatthew G. Knepley DMLabel label = jackeys[k].label; 4988e3a2eefSMatthew G. Knepley PetscInt val = jackeys[k].value; 4998e3a2eefSMatthew G. Knepley 5008e3a2eefSMatthew G. Knepley if (!label) { 5019566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)allcellIS)); 5028e3a2eefSMatthew G. Knepley cellIS = allcellIS; 5037a73cf09SMatthew G. Knepley } else { 5048e3a2eefSMatthew G. Knepley IS pointIS; 505a925c78cSMatthew G. Knepley 5069566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, val, &pointIS)); 5079566063dSJacob Faibussowitsch PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 5089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 509a925c78cSMatthew G. Knepley } 510*754e4fbaSMatthew G. Knepley PetscCall(DMPlexComputeJacobianActionByKey(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user)); 5119566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 5128e3a2eefSMatthew G. Knepley } 5139566063dSJacob Faibussowitsch PetscCall(PetscFree(jackeys)); 5148e3a2eefSMatthew G. Knepley } 5158e3a2eefSMatthew G. Knepley } 5169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&allcellIS)); 5179566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 5183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 519a925c78cSMatthew G. Knepley } 520a925c78cSMatthew G. Knepley 52124cdb843SMatthew G. Knepley /*@ 5222fe279fdSBarry Smith DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix `Jac` at the local solution `X` using pointwise functions specified by the user. 52324cdb843SMatthew G. Knepley 52424cdb843SMatthew G. Knepley Input Parameters: 5252fe279fdSBarry Smith + dm - The `DM` 52624cdb843SMatthew G. Knepley . X - Local input vector 52724cdb843SMatthew G. Knepley - user - The user context 52824cdb843SMatthew G. Knepley 5292fe279fdSBarry Smith Output Parameters: 5302fe279fdSBarry Smith + Jac - Jacobian matrix 5312fe279fdSBarry Smith - JacP - approximate Jacobian from which the preconditioner will be built, often `Jac` 5322fe279fdSBarry Smith 5332fe279fdSBarry Smith Level: developer 53424cdb843SMatthew G. Knepley 53524cdb843SMatthew G. Knepley Note: 53624cdb843SMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 53724cdb843SMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 53824cdb843SMatthew G. Knepley 539420bcc1bSBarry Smith .seealso: [](ch_snes), `DMPLEX`, `Mat` 54024cdb843SMatthew G. Knepley @*/ 541d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user) 542d71ae5a4SJacob Faibussowitsch { 5436da023fcSToby Isaac DM plex; 544083401c6SMatthew G. Knepley IS allcellIS; 545f04eb4edSMatthew G. Knepley PetscBool hasJac, hasPrec; 5466528b96dSMatthew G. Knepley PetscInt Nds, s; 54724cdb843SMatthew G. Knepley 54824cdb843SMatthew G. Knepley PetscFunctionBegin; 5499566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 5509566063dSJacob Faibussowitsch PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 5519566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 552083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 553083401c6SMatthew G. Knepley PetscDS ds; 554083401c6SMatthew G. Knepley IS cellIS; 55506ad1575SMatthew G. Knepley PetscFormKey key; 556083401c6SMatthew G. Knepley 55707218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL)); 5586528b96dSMatthew G. Knepley key.value = 0; 5596528b96dSMatthew G. Knepley key.field = 0; 56006ad1575SMatthew G. Knepley key.part = 0; 5616528b96dSMatthew G. Knepley if (!key.label) { 5629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)allcellIS)); 563083401c6SMatthew G. Knepley cellIS = allcellIS; 564083401c6SMatthew G. Knepley } else { 565083401c6SMatthew G. Knepley IS pointIS; 566083401c6SMatthew G. Knepley 5676528b96dSMatthew G. Knepley key.value = 1; 5689566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS)); 5699566063dSJacob Faibussowitsch PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 5709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 571083401c6SMatthew G. Knepley } 572083401c6SMatthew G. Knepley if (!s) { 5739566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobian(ds, &hasJac)); 5749566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec)); 5759566063dSJacob Faibussowitsch if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac)); 5769566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(JacP)); 577083401c6SMatthew G. Knepley } 578*754e4fbaSMatthew G. Knepley PetscCall(DMPlexComputeJacobianByKey(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user)); 5799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 580083401c6SMatthew G. Knepley } 5819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&allcellIS)); 5829566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58424cdb843SMatthew G. Knepley } 5851878804aSMatthew G. Knepley 5869371c9d4SSatish Balay struct _DMSNESJacobianMFCtx { 5878e3a2eefSMatthew G. Knepley DM dm; 5888e3a2eefSMatthew G. Knepley Vec X; 5898e3a2eefSMatthew G. Knepley void *ctx; 5908e3a2eefSMatthew G. Knepley }; 5918e3a2eefSMatthew G. Knepley 592d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A) 593d71ae5a4SJacob Faibussowitsch { 5948e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx *ctx; 5958e3a2eefSMatthew G. Knepley 5968e3a2eefSMatthew G. Knepley PetscFunctionBegin; 5979566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 5989566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, NULL)); 5999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 6009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->X)); 6019566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 6023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6038e3a2eefSMatthew G. Knepley } 6048e3a2eefSMatthew G. Knepley 605d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z) 606d71ae5a4SJacob Faibussowitsch { 6078e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx *ctx; 6088e3a2eefSMatthew G. Knepley 6098e3a2eefSMatthew G. Knepley PetscFunctionBegin; 6109566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 6119566063dSJacob Faibussowitsch PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx)); 6123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6138e3a2eefSMatthew G. Knepley } 6148e3a2eefSMatthew G. Knepley 6158e3a2eefSMatthew G. Knepley /*@ 616f6dfbefdSBarry Smith DMSNESCreateJacobianMF - Create a `Mat` which computes the action of the Jacobian matrix-free 6178e3a2eefSMatthew G. Knepley 618c3339decSBarry Smith Collective 6198e3a2eefSMatthew G. Knepley 6208e3a2eefSMatthew G. Knepley Input Parameters: 621f6dfbefdSBarry Smith + dm - The `DM` 6228e3a2eefSMatthew G. Knepley . X - The evaluation point for the Jacobian 6232fe279fdSBarry Smith - user - A user context, or `NULL` 6248e3a2eefSMatthew G. Knepley 6258e3a2eefSMatthew G. Knepley Output Parameter: 626f6dfbefdSBarry Smith . J - The `Mat` 6278e3a2eefSMatthew G. Knepley 6288e3a2eefSMatthew G. Knepley Level: advanced 6298e3a2eefSMatthew G. Knepley 630c118d280SBarry Smith Notes: 6312fe279fdSBarry Smith Vec `X` is kept in `J`, so updating `X` then updates the evaluation point. 6328e3a2eefSMatthew G. Knepley 633c118d280SBarry Smith This only works for `DMPLEX` 634c118d280SBarry Smith 635420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `SNES`, `DMSNESComputeJacobianAction()` 6368e3a2eefSMatthew G. Knepley @*/ 637d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J) 638d71ae5a4SJacob Faibussowitsch { 6398e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx *ctx; 6408e3a2eefSMatthew G. Knepley PetscInt n, N; 6418e3a2eefSMatthew G. Knepley 6428e3a2eefSMatthew G. Knepley PetscFunctionBegin; 6439566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 6449566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATSHELL)); 6459566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 6469566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N)); 6479566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, n, n, N, N)); 6489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 6499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)X)); 6509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx)); 6518e3a2eefSMatthew G. Knepley ctx->dm = dm; 6528e3a2eefSMatthew G. Knepley ctx->X = X; 6538e3a2eefSMatthew G. Knepley ctx->ctx = user; 6549566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(*J, ctx)); 6559566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private)); 6569566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private)); 6573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6588e3a2eefSMatthew G. Knepley } 6598e3a2eefSMatthew G. Knepley 660d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx) 661d71ae5a4SJacob Faibussowitsch { 66238cfc46eSPierre Jolivet SNES snes; 66338cfc46eSPierre Jolivet Mat pJ; 66438cfc46eSPierre Jolivet DM ovldm, origdm; 66538cfc46eSPierre Jolivet DMSNES sdm; 66638cfc46eSPierre Jolivet PetscErrorCode (*bfun)(DM, Vec, void *); 66738cfc46eSPierre Jolivet PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *); 66838cfc46eSPierre Jolivet void *bctx, *jctx; 66938cfc46eSPierre Jolivet 67038cfc46eSPierre Jolivet PetscFunctionBegin; 6719566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ)); 67228b400f6SJacob Faibussowitsch PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat"); 6739566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm)); 67428b400f6SJacob Faibussowitsch PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM"); 6759566063dSJacob Faibussowitsch PetscCall(MatGetDM(pJ, &ovldm)); 6769566063dSJacob Faibussowitsch PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx)); 6779566063dSJacob Faibussowitsch PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx)); 6789566063dSJacob Faibussowitsch PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx)); 6799566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes)); 68138cfc46eSPierre Jolivet if (!snes) { 6829566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes)); 6839566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, ovldm)); 6849566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes)); 6859566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)snes)); 68638cfc46eSPierre Jolivet } 6879566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(ovldm, &sdm)); 6889566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X)); 689800f99ffSJeremy L Thompson { 690800f99ffSJeremy L Thompson void *ctx; 691800f99ffSJeremy L Thompson PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *); 692800f99ffSJeremy L Thompson PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx)); 693800f99ffSJeremy L Thompson PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx)); 694800f99ffSJeremy L Thompson } 6959566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X)); 69638cfc46eSPierre Jolivet /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */ 69738cfc46eSPierre Jolivet { 69838cfc46eSPierre Jolivet Mat locpJ; 69938cfc46eSPierre Jolivet 7009566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(pJ, &locpJ)); 7019566063dSJacob Faibussowitsch PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN)); 70238cfc46eSPierre Jolivet } 7033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70438cfc46eSPierre Jolivet } 70538cfc46eSPierre Jolivet 706a925c78cSMatthew G. Knepley /*@ 7076493148fSStefano Zampini DMPlexSetSNESLocalFEM - Use `DMPLEX`'s internal FEM routines to compute `SNES` boundary values, objective, residual, and Jacobian. 7089f520fc2SToby Isaac 7099f520fc2SToby Isaac Input Parameters: 710f6dfbefdSBarry Smith + dm - The `DM` object 7116493148fSStefano Zampini . use_obj - Use the objective function callback 7126493148fSStefano Zampini - ctx - The user context that will be passed to pointwise evaluation routines 7131a244344SSatish Balay 7141a244344SSatish Balay Level: developer 715f6dfbefdSBarry Smith 7166493148fSStefano Zampini .seealso: [](ch_snes),`DMPLEX`, `SNES`, `PetscDSAddBoundary()`, `PetscDSSetObjective()`, `PetscDSSetResidual()`, `PetscDSSetJacobian()` 7179f520fc2SToby Isaac @*/ 7186493148fSStefano Zampini PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, PetscBool use_obj, void *ctx) 719d71ae5a4SJacob Faibussowitsch { 720d2b2dc1eSMatthew G. Knepley PetscBool useCeed; 721d2b2dc1eSMatthew G. Knepley 7229f520fc2SToby Isaac PetscFunctionBegin; 723d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 7246493148fSStefano Zampini PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, ctx)); 7256493148fSStefano Zampini if (use_obj) PetscCall(DMSNESSetObjectiveLocal(dm, DMPlexSNESComputeObjectiveFEM, ctx)); 726d2b2dc1eSMatthew G. Knepley if (useCeed) { 727d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 7286493148fSStefano Zampini PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualCEED, ctx)); 729d2b2dc1eSMatthew G. Knepley #else 730d2b2dc1eSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot use CEED traversals without LibCEED. Rerun configure with --download-ceed"); 731d2b2dc1eSMatthew G. Knepley #endif 7326493148fSStefano Zampini } else PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, ctx)); 7336493148fSStefano Zampini PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, ctx)); 7349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex)); 7353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7369f520fc2SToby Isaac } 7379f520fc2SToby Isaac 738cc4c1da9SBarry Smith /*@ 739f017bcb6SMatthew G. Knepley DMSNESCheckDiscretization - Check the discretization error of the exact solution 740f017bcb6SMatthew G. Knepley 741f017bcb6SMatthew G. Knepley Input Parameters: 742f6dfbefdSBarry Smith + snes - the `SNES` object 743f6dfbefdSBarry Smith . dm - the `DM` 744f2cacb80SMatthew G. Knepley . t - the time 745f6dfbefdSBarry Smith . u - a `DM` vector 746f0fc11ceSJed Brown - tol - A tolerance for the check, or -1 to print the results instead 747f017bcb6SMatthew G. Knepley 7482fe279fdSBarry Smith Output Parameter: 7492fe279fdSBarry Smith . error - An array which holds the discretization error in each field, or `NULL` 7502fe279fdSBarry Smith 7512fe279fdSBarry Smith Level: developer 752f3c94560SMatthew G. Knepley 753f6dfbefdSBarry Smith Note: 754f6dfbefdSBarry Smith The user must call `PetscDSSetExactSolution()` beforehand 7557f96f943SMatthew G. Knepley 756420bcc1bSBarry Smith Developer Note: 757420bcc1bSBarry Smith How is this related to `PetscConvEst`? 758420bcc1bSBarry Smith 759420bcc1bSBarry Smith .seealso: [](ch_snes), `PetscDSSetExactSolution()`, `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()` 760f017bcb6SMatthew G. Knepley @*/ 761d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[]) 762d71ae5a4SJacob Faibussowitsch { 763f017bcb6SMatthew G. Knepley PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 764f017bcb6SMatthew G. Knepley void **ectxs; 765f3c94560SMatthew G. Knepley PetscReal *err; 7667f96f943SMatthew G. Knepley MPI_Comm comm; 7677f96f943SMatthew G. Knepley PetscInt Nf, f; 7681878804aSMatthew G. Knepley 7691878804aSMatthew G. Knepley PetscFunctionBegin; 770f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 771f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 772064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 7734f572ea9SToby Isaac if (error) PetscAssertPointer(error, 6); 7747f96f943SMatthew G. Knepley 7759566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 7769566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-vec_view")); 7777f96f943SMatthew G. Knepley 7789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 7799566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 7809566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err)); 7817f96f943SMatthew G. Knepley { 7827f96f943SMatthew G. Knepley PetscInt Nds, s; 7837f96f943SMatthew G. Knepley 7849566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 785083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 7867f96f943SMatthew G. Knepley PetscDS ds; 787083401c6SMatthew G. Knepley DMLabel label; 788083401c6SMatthew G. Knepley IS fieldIS; 7897f96f943SMatthew G. Knepley const PetscInt *fields; 7907f96f943SMatthew G. Knepley PetscInt dsNf, f; 791083401c6SMatthew G. Knepley 79207218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL)); 7939566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &dsNf)); 7949566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fieldIS, &fields)); 795083401c6SMatthew G. Knepley for (f = 0; f < dsNf; ++f) { 796083401c6SMatthew G. Knepley const PetscInt field = fields[f]; 7979566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field])); 798083401c6SMatthew G. Knepley } 7999566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fieldIS, &fields)); 800083401c6SMatthew G. Knepley } 801083401c6SMatthew G. Knepley } 802f017bcb6SMatthew G. Knepley if (Nf > 1) { 8039566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err)); 804f017bcb6SMatthew G. Knepley if (tol >= 0.0) { 805ad540459SPierre Jolivet for (f = 0; f < Nf; ++f) PetscCheck(err[f] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g for field %" PetscInt_FMT " exceeds tolerance %g", (double)err[f], f, (double)tol); 806b878532fSMatthew G. Knepley } else if (error) { 807b878532fSMatthew G. Knepley for (f = 0; f < Nf; ++f) error[f] = err[f]; 8081878804aSMatthew G. Knepley } else { 8099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Error: [")); 810f017bcb6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 8119566063dSJacob Faibussowitsch if (f) PetscCall(PetscPrintf(comm, ", ")); 8129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "%g", (double)err[f])); 8131878804aSMatthew G. Knepley } 8149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "]\n")); 815f017bcb6SMatthew G. Knepley } 816f017bcb6SMatthew G. Knepley } else { 8179566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0])); 818f017bcb6SMatthew G. Knepley if (tol >= 0.0) { 81908401ef6SPierre Jolivet PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol); 820b878532fSMatthew G. Knepley } else if (error) { 821b878532fSMatthew G. Knepley error[0] = err[0]; 822f017bcb6SMatthew G. Knepley } else { 8239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0])); 824f017bcb6SMatthew G. Knepley } 825f017bcb6SMatthew G. Knepley } 8269566063dSJacob Faibussowitsch PetscCall(PetscFree3(exacts, ectxs, err)); 8273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 828f017bcb6SMatthew G. Knepley } 829f017bcb6SMatthew G. Knepley 830cc4c1da9SBarry Smith /*@ 831f017bcb6SMatthew G. Knepley DMSNESCheckResidual - Check the residual of the exact solution 832f017bcb6SMatthew G. Knepley 833f017bcb6SMatthew G. Knepley Input Parameters: 834f6dfbefdSBarry Smith + snes - the `SNES` object 835f6dfbefdSBarry Smith . dm - the `DM` 836f6dfbefdSBarry Smith . u - a `DM` vector 837f0fc11ceSJed Brown - tol - A tolerance for the check, or -1 to print the results instead 838f017bcb6SMatthew G. Knepley 839f6dfbefdSBarry Smith Output Parameter: 8402fe279fdSBarry Smith . residual - The residual norm of the exact solution, or `NULL` 841f3c94560SMatthew G. Knepley 842f017bcb6SMatthew G. Knepley Level: developer 843f017bcb6SMatthew G. Knepley 844420bcc1bSBarry Smith .seealso: [](ch_snes), `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()` 845f017bcb6SMatthew G. Knepley @*/ 846d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual) 847d71ae5a4SJacob Faibussowitsch { 848f017bcb6SMatthew G. Knepley MPI_Comm comm; 849f017bcb6SMatthew G. Knepley Vec r; 850f017bcb6SMatthew G. Knepley PetscReal res; 851f017bcb6SMatthew G. Knepley 852f017bcb6SMatthew G. Knepley PetscFunctionBegin; 853f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 854f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 855f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 8564f572ea9SToby Isaac if (residual) PetscAssertPointer(residual, 5); 8579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 8589566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL)); 8599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 8609566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r)); 8619566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res)); 862f017bcb6SMatthew G. Knepley if (tol >= 0.0) { 86308401ef6SPierre Jolivet PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol); 864b878532fSMatthew G. Knepley } else if (residual) { 865b878532fSMatthew G. Knepley *residual = res; 866f017bcb6SMatthew G. Knepley } else { 8679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res)); 86893ec0da9SPierre Jolivet PetscCall(VecFilter(r, 1.0e-10)); 8699566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual")); 8709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_")); 87164bbb635SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", (PetscObject)snes)); 8729566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(r, NULL, "-vec_view")); 87364bbb635SMatthew G. Knepley PetscCall(PetscObjectCompose((PetscObject)r, "__Vec_bc_zero__", NULL)); 874f017bcb6SMatthew G. Knepley } 8759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 8763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 877f017bcb6SMatthew G. Knepley } 878f017bcb6SMatthew G. Knepley 879cc4c1da9SBarry Smith /*@ 880f3c94560SMatthew G. Knepley DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 881f017bcb6SMatthew G. Knepley 882f017bcb6SMatthew G. Knepley Input Parameters: 883f6dfbefdSBarry Smith + snes - the `SNES` object 884f6dfbefdSBarry Smith . dm - the `DM` 885f6dfbefdSBarry Smith . u - a `DM` vector 886f0fc11ceSJed Brown - tol - A tolerance for the check, or -1 to print the results instead 887f017bcb6SMatthew G. Knepley 888f3c94560SMatthew G. Knepley Output Parameters: 8892fe279fdSBarry Smith + isLinear - Flag indicaing that the function looks linear, or `NULL` 8902fe279fdSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL` 891f3c94560SMatthew G. Knepley 892f017bcb6SMatthew G. Knepley Level: developer 893f017bcb6SMatthew G. Knepley 894420bcc1bSBarry Smith .seealso: [](ch_snes), `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()` 895f017bcb6SMatthew G. Knepley @*/ 896d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) 897d71ae5a4SJacob Faibussowitsch { 898f017bcb6SMatthew G. Knepley MPI_Comm comm; 899f017bcb6SMatthew G. Knepley PetscDS ds; 900f017bcb6SMatthew G. Knepley Mat J, M; 901f017bcb6SMatthew G. Knepley MatNullSpace nullspace; 902f3c94560SMatthew G. Knepley PetscReal slope, intercept; 9037f96f943SMatthew G. Knepley PetscBool hasJac, hasPrec, isLin = PETSC_FALSE; 904f017bcb6SMatthew G. Knepley 905f017bcb6SMatthew G. Knepley PetscFunctionBegin; 906f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 907d5d53b22SStefano Zampini if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 908d5d53b22SStefano Zampini if (u) PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 9094f572ea9SToby Isaac if (isLinear) PetscAssertPointer(isLinear, 5); 9104f572ea9SToby Isaac if (convRate) PetscAssertPointer(convRate, 6); 9119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 912d5d53b22SStefano Zampini if (!dm) PetscCall(SNESGetDM(snes, &dm)); 913d5d53b22SStefano Zampini if (u) PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL)); 914d5d53b22SStefano Zampini else PetscCall(SNESGetSolution(snes, &u)); 915f017bcb6SMatthew G. Knepley /* Create and view matrices */ 9169566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J)); 9179566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 9189566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobian(ds, &hasJac)); 9199566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec)); 920282e7bb4SMatthew G. Knepley if (hasJac && hasPrec) { 9219566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 9229566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, J, M)); 9239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix")); 9249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_")); 9259566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-mat_view")); 9269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 927282e7bb4SMatthew G. Knepley } else { 9289566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, J, J)); 929282e7bb4SMatthew G. Knepley } 9309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 9319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_")); 9329566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(J, NULL, "-mat_view")); 933f017bcb6SMatthew G. Knepley /* Check nullspace */ 9349566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(J, &nullspace)); 935f017bcb6SMatthew G. Knepley if (nullspace) { 9361878804aSMatthew G. Knepley PetscBool isNull; 9379566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nullspace, J, &isNull)); 93828b400f6SJacob Faibussowitsch PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 9391878804aSMatthew G. Knepley } 940f017bcb6SMatthew G. Knepley /* Taylor test */ 941f017bcb6SMatthew G. Knepley { 942f017bcb6SMatthew G. Knepley PetscRandom rand; 943f017bcb6SMatthew G. Knepley Vec du, uhat, r, rhat, df; 944f3c94560SMatthew G. Knepley PetscReal h; 945f017bcb6SMatthew G. Knepley PetscReal *es, *hs, *errors; 946f017bcb6SMatthew G. Knepley PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 947f017bcb6SMatthew G. Knepley PetscInt Nv, v; 948f017bcb6SMatthew G. Knepley 949f017bcb6SMatthew G. Knepley /* Choose a perturbation direction */ 9509566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand)); 9519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &du)); 9529566063dSJacob Faibussowitsch PetscCall(VecSetRandom(du, rand)); 9539566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 9549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &df)); 9559566063dSJacob Faibussowitsch PetscCall(MatMult(J, du, df)); 956f017bcb6SMatthew G. Knepley /* Evaluate residual at u, F(u), save in vector r */ 9579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 9589566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r)); 959f017bcb6SMatthew G. Knepley /* Look at the convergence of our Taylor approximation as we approach u */ 960fbccb6d4SPierre Jolivet for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv); 9619566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors)); 9629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &uhat)); 9639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &rhat)); 964f017bcb6SMatthew G. Knepley for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 9659566063dSJacob Faibussowitsch PetscCall(VecWAXPY(uhat, h, du, u)); 966f017bcb6SMatthew G. Knepley /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */ 9679566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, uhat, rhat)); 9689566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df)); 9699566063dSJacob Faibussowitsch PetscCall(VecNorm(rhat, NORM_2, &errors[Nv])); 970f017bcb6SMatthew G. Knepley 97138d69901SMatthew G. Knepley es[Nv] = errors[Nv] == 0 ? -16.0 : PetscLog10Real(errors[Nv]); 972f017bcb6SMatthew G. Knepley hs[Nv] = PetscLog10Real(h); 973f017bcb6SMatthew G. Knepley } 9749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uhat)); 9759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rhat)); 9769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&df)); 9779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 9789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&du)); 979f017bcb6SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 980f017bcb6SMatthew G. Knepley if ((tol >= 0) && (errors[v] > tol)) break; 981f017bcb6SMatthew G. Knepley else if (errors[v] > PETSC_SMALL) break; 982f017bcb6SMatthew G. Knepley } 983f3c94560SMatthew G. Knepley if (v == Nv) isLin = PETSC_TRUE; 9849566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept)); 9859566063dSJacob Faibussowitsch PetscCall(PetscFree3(es, hs, errors)); 986f017bcb6SMatthew G. Knepley /* Slope should be about 2 */ 987f017bcb6SMatthew G. Knepley if (tol >= 0) { 9880b121fc5SBarry Smith PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope); 989b878532fSMatthew G. Knepley } else if (isLinear || convRate) { 990b878532fSMatthew G. Knepley if (isLinear) *isLinear = isLin; 991b878532fSMatthew G. Knepley if (convRate) *convRate = slope; 992f017bcb6SMatthew G. Knepley } else { 9939566063dSJacob Faibussowitsch if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope)); 9949566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function appears to be linear\n")); 995f017bcb6SMatthew G. Knepley } 996f017bcb6SMatthew G. Knepley } 9979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9991878804aSMatthew G. Knepley } 10001878804aSMatthew G. Knepley 1001d6acfc2dSPierre Jolivet static PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u) 1002d71ae5a4SJacob Faibussowitsch { 1003f017bcb6SMatthew G. Knepley PetscFunctionBegin; 10049566063dSJacob Faibussowitsch PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL)); 10059566063dSJacob Faibussowitsch PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL)); 10069566063dSJacob Faibussowitsch PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL)); 10073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1008f017bcb6SMatthew G. Knepley } 1009f017bcb6SMatthew G. Knepley 1010cc4c1da9SBarry Smith /*@ 1011bee9a294SMatthew G. Knepley DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 1012bee9a294SMatthew G. Knepley 1013bee9a294SMatthew G. Knepley Input Parameters: 1014f6dfbefdSBarry Smith + snes - the `SNES` object 1015f6dfbefdSBarry Smith - u - representative `SNES` vector 10167f96f943SMatthew G. Knepley 10172fe279fdSBarry Smith Level: developer 10182fe279fdSBarry Smith 1019f6dfbefdSBarry Smith Note: 1020420bcc1bSBarry Smith The user must call `PetscDSSetExactSolution()` before this call 1021bee9a294SMatthew G. Knepley 1022420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DM` 1023bee9a294SMatthew G. Knepley @*/ 1024d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u) 1025d71ae5a4SJacob Faibussowitsch { 10261878804aSMatthew G. Knepley DM dm; 10271878804aSMatthew G. Knepley Vec sol; 10281878804aSMatthew G. Knepley PetscBool check; 10291878804aSMatthew G. Knepley 10301878804aSMatthew G. Knepley PetscFunctionBegin; 10319566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check)); 10323ba16761SJacob Faibussowitsch if (!check) PetscFunctionReturn(PETSC_SUCCESS); 10339566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 10349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &sol)); 10359566063dSJacob Faibussowitsch PetscCall(SNESSetSolution(snes, sol)); 10369566063dSJacob Faibussowitsch PetscCall(DMSNESCheck_Internal(snes, dm, sol)); 10379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sol)); 10383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1039552f7358SJed Brown } 10402eace19aSMatthew G. Knepley 10412eace19aSMatthew G. Knepley static PetscErrorCode lower_inf(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 10422eace19aSMatthew G. Knepley { 10432eace19aSMatthew G. Knepley PetscFunctionBegin; 10442eace19aSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) u[c] = PETSC_NINFINITY; 10452eace19aSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 10462eace19aSMatthew G. Knepley } 10472eace19aSMatthew G. Knepley 10482eace19aSMatthew G. Knepley static PetscErrorCode upper_inf(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 10492eace19aSMatthew G. Knepley { 10502eace19aSMatthew G. Knepley PetscFunctionBegin; 10512eace19aSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) u[c] = PETSC_INFINITY; 10522eace19aSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 10532eace19aSMatthew G. Knepley } 10542eace19aSMatthew G. Knepley 10552eace19aSMatthew G. Knepley PetscErrorCode DMPlexSetSNESVariableBounds(DM dm, SNES snes) 10562eace19aSMatthew G. Knepley { 10572eace19aSMatthew G. Knepley PetscDS ds; 10582eace19aSMatthew G. Knepley Vec lb, ub; 10592eace19aSMatthew G. Knepley PetscSimplePointFunc *lfuncs, *ufuncs; 10602eace19aSMatthew G. Knepley void **lctxs, **uctxs; 10612eace19aSMatthew G. Knepley PetscBool hasBound = PETSC_FALSE; 10622eace19aSMatthew G. Knepley PetscInt Nf; 10632eace19aSMatthew G. Knepley 10642eace19aSMatthew G. Knepley PetscFunctionBegin; 10652eace19aSMatthew G. Knepley // TODO Generalize for multiple DSes 10662eace19aSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 10672eace19aSMatthew G. Knepley PetscCall(PetscDSGetNumFields(ds, &Nf)); 10682eace19aSMatthew G. Knepley PetscCall(PetscMalloc4(Nf, &lfuncs, Nf, &lctxs, Nf, &ufuncs, Nf, &uctxs)); 10692eace19aSMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) { 10702eace19aSMatthew G. Knepley PetscCall(PetscDSGetLowerBound(ds, f, &lfuncs[f], &lctxs[f])); 10712eace19aSMatthew G. Knepley PetscCall(PetscDSGetUpperBound(ds, f, &ufuncs[f], &uctxs[f])); 10722eace19aSMatthew G. Knepley if (lfuncs[f] || ufuncs[f]) hasBound = PETSC_TRUE; 10732eace19aSMatthew G. Knepley if (!lfuncs[f]) lfuncs[f] = lower_inf; 10742eace19aSMatthew G. Knepley if (!ufuncs[f]) ufuncs[f] = upper_inf; 10752eace19aSMatthew G. Knepley } 10762eace19aSMatthew G. Knepley if (!hasBound) goto end; 10772eace19aSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &lb)); 10782eace19aSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &ub)); 10792eace19aSMatthew G. Knepley PetscCall(DMProjectFunction(dm, 0., lfuncs, lctxs, INSERT_VALUES, lb)); 10802eace19aSMatthew G. Knepley PetscCall(DMProjectFunction(dm, 0., ufuncs, uctxs, INSERT_VALUES, ub)); 10812eace19aSMatthew G. Knepley PetscCall(VecViewFromOptions(lb, NULL, "-dm_plex_snes_lb_view")); 10822eace19aSMatthew G. Knepley PetscCall(VecViewFromOptions(ub, NULL, "-dm_plex_snes_ub_view")); 10832eace19aSMatthew G. Knepley PetscCall(SNESVISetVariableBounds(snes, lb, ub)); 10842eace19aSMatthew G. Knepley PetscCall(VecDestroy(&lb)); 10852eace19aSMatthew G. Knepley PetscCall(VecDestroy(&ub)); 10862eace19aSMatthew G. Knepley end: 10872eace19aSMatthew G. Knepley PetscCall(PetscFree4(lfuncs, lctxs, ufuncs, uctxs)); 10882eace19aSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 10892eace19aSMatthew G. Knepley } 1090