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 14124cdb843SMatthew G. Knepley /************************** Interpolation *******************************/ 14224cdb843SMatthew G. Knepley 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESConvertPlex(DM dm, DM *plex, PetscBool copy) 144d71ae5a4SJacob Faibussowitsch { 1456da023fcSToby Isaac PetscBool isPlex; 1466da023fcSToby Isaac 1476da023fcSToby Isaac PetscFunctionBegin; 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &isPlex)); 1496da023fcSToby Isaac if (isPlex) { 1506da023fcSToby Isaac *plex = dm; 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 152f7148743SMatthew G. Knepley } else { 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "dm_plex", (PetscObject *)plex)); 154f7148743SMatthew G. Knepley if (!*plex) { 1559566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, DMPLEX, plex)); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "dm_plex", (PetscObject)*plex)); 1576da023fcSToby Isaac if (copy) { 1589566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(dm, *plex)); 1599566063dSJacob Faibussowitsch PetscCall(DMCopyAuxiliaryVec(dm, *plex)); 1606da023fcSToby Isaac } 161f7148743SMatthew G. Knepley } else { 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)*plex)); 163f7148743SMatthew G. Knepley } 1646da023fcSToby Isaac } 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1666da023fcSToby Isaac } 1676da023fcSToby Isaac 1684267b1a3SMatthew G. Knepley /*@C 169f6dfbefdSBarry Smith DMInterpolationCreate - Creates a `DMInterpolationInfo` context 1704267b1a3SMatthew G. Knepley 171d083f849SBarry Smith Collective 1724267b1a3SMatthew G. Knepley 1734267b1a3SMatthew G. Knepley Input Parameter: 1744267b1a3SMatthew G. Knepley . comm - the communicator 1754267b1a3SMatthew G. Knepley 1764267b1a3SMatthew G. Knepley Output Parameter: 1774267b1a3SMatthew G. Knepley . ctx - the context 1784267b1a3SMatthew G. Knepley 1794267b1a3SMatthew G. Knepley Level: beginner 1804267b1a3SMatthew G. Knepley 181420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()` 1824267b1a3SMatthew G. Knepley @*/ 183d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 184d71ae5a4SJacob Faibussowitsch { 185552f7358SJed Brown PetscFunctionBegin; 1864f572ea9SToby Isaac PetscAssertPointer(ctx, 2); 1879566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 1881aa26658SKarl Rupp 189552f7358SJed Brown (*ctx)->comm = comm; 190552f7358SJed Brown (*ctx)->dim = -1; 191552f7358SJed Brown (*ctx)->nInput = 0; 1920298fd71SBarry Smith (*ctx)->points = NULL; 1930298fd71SBarry Smith (*ctx)->cells = NULL; 194552f7358SJed Brown (*ctx)->n = -1; 1950298fd71SBarry Smith (*ctx)->coords = NULL; 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197552f7358SJed Brown } 198552f7358SJed Brown 1994267b1a3SMatthew G. Knepley /*@C 2004267b1a3SMatthew G. Knepley DMInterpolationSetDim - Sets the spatial dimension for the interpolation context 2014267b1a3SMatthew G. Knepley 20220f4b53cSBarry Smith Not Collective 2034267b1a3SMatthew G. Knepley 2044267b1a3SMatthew G. Knepley Input Parameters: 2054267b1a3SMatthew G. Knepley + ctx - the context 2064267b1a3SMatthew G. Knepley - dim - the spatial dimension 2074267b1a3SMatthew G. Knepley 2084267b1a3SMatthew G. Knepley Level: intermediate 2094267b1a3SMatthew G. Knepley 210420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 2114267b1a3SMatthew G. Knepley @*/ 212d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 213d71ae5a4SJacob Faibussowitsch { 214552f7358SJed Brown PetscFunctionBegin; 21563a3b9bcSJacob Faibussowitsch PetscCheck(!(dim < 1) && !(dim > 3), ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim); 216552f7358SJed Brown ctx->dim = dim; 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218552f7358SJed Brown } 219552f7358SJed Brown 2204267b1a3SMatthew G. Knepley /*@C 2214267b1a3SMatthew G. Knepley DMInterpolationGetDim - Gets the spatial dimension for the interpolation context 2224267b1a3SMatthew G. Knepley 22320f4b53cSBarry Smith Not Collective 2244267b1a3SMatthew G. Knepley 2254267b1a3SMatthew G. Knepley Input Parameter: 2264267b1a3SMatthew G. Knepley . ctx - the context 2274267b1a3SMatthew G. Knepley 2284267b1a3SMatthew G. Knepley Output Parameter: 2294267b1a3SMatthew G. Knepley . dim - the spatial dimension 2304267b1a3SMatthew G. Knepley 2314267b1a3SMatthew G. Knepley Level: intermediate 2324267b1a3SMatthew G. Knepley 233420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 2344267b1a3SMatthew G. Knepley @*/ 235d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 236d71ae5a4SJacob Faibussowitsch { 237552f7358SJed Brown PetscFunctionBegin; 2384f572ea9SToby Isaac PetscAssertPointer(dim, 2); 239552f7358SJed Brown *dim = ctx->dim; 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 241552f7358SJed Brown } 242552f7358SJed Brown 2434267b1a3SMatthew G. Knepley /*@C 2444267b1a3SMatthew G. Knepley DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context 2454267b1a3SMatthew G. Knepley 24620f4b53cSBarry Smith Not Collective 2474267b1a3SMatthew G. Knepley 2484267b1a3SMatthew G. Knepley Input Parameters: 2494267b1a3SMatthew G. Knepley + ctx - the context 2504267b1a3SMatthew G. Knepley - dof - the number of fields 2514267b1a3SMatthew G. Knepley 2524267b1a3SMatthew G. Knepley Level: intermediate 2534267b1a3SMatthew G. Knepley 254420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 2554267b1a3SMatthew G. Knepley @*/ 256d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 257d71ae5a4SJacob Faibussowitsch { 258552f7358SJed Brown PetscFunctionBegin; 25963a3b9bcSJacob Faibussowitsch PetscCheck(dof >= 1, ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof); 260552f7358SJed Brown ctx->dof = dof; 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 262552f7358SJed Brown } 263552f7358SJed Brown 2644267b1a3SMatthew G. Knepley /*@C 2654267b1a3SMatthew G. Knepley DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context 2664267b1a3SMatthew G. Knepley 26720f4b53cSBarry Smith Not Collective 2684267b1a3SMatthew G. Knepley 2694267b1a3SMatthew G. Knepley Input Parameter: 2704267b1a3SMatthew G. Knepley . ctx - the context 2714267b1a3SMatthew G. Knepley 2724267b1a3SMatthew G. Knepley Output Parameter: 2734267b1a3SMatthew G. Knepley . dof - the number of fields 2744267b1a3SMatthew G. Knepley 2754267b1a3SMatthew G. Knepley Level: intermediate 2764267b1a3SMatthew G. Knepley 277420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 2784267b1a3SMatthew G. Knepley @*/ 279d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 280d71ae5a4SJacob Faibussowitsch { 281552f7358SJed Brown PetscFunctionBegin; 2824f572ea9SToby Isaac PetscAssertPointer(dof, 2); 283552f7358SJed Brown *dof = ctx->dof; 2843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 285552f7358SJed Brown } 286552f7358SJed Brown 2874267b1a3SMatthew G. Knepley /*@C 2884267b1a3SMatthew G. Knepley DMInterpolationAddPoints - Add points at which we will interpolate the fields 2894267b1a3SMatthew G. Knepley 29020f4b53cSBarry Smith Not Collective 2914267b1a3SMatthew G. Knepley 2924267b1a3SMatthew G. Knepley Input Parameters: 2934267b1a3SMatthew G. Knepley + ctx - the context 2944267b1a3SMatthew G. Knepley . n - the number of points 2954267b1a3SMatthew G. Knepley - points - the coordinates for each point, an array of size n * dim 2964267b1a3SMatthew G. Knepley 2972fe279fdSBarry Smith Level: intermediate 2982fe279fdSBarry Smith 299f6dfbefdSBarry Smith Note: 300f6dfbefdSBarry Smith The coordinate information is copied. 3014267b1a3SMatthew G. Knepley 302420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()` 3034267b1a3SMatthew G. Knepley @*/ 304d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 305d71ae5a4SJacob Faibussowitsch { 306552f7358SJed Brown PetscFunctionBegin; 30708401ef6SPierre Jolivet PetscCheck(ctx->dim >= 0, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 30828b400f6SJacob Faibussowitsch PetscCheck(!ctx->points, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 309552f7358SJed Brown ctx->nInput = n; 3101aa26658SKarl Rupp 3119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n * ctx->dim, &ctx->points)); 3129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ctx->points, points, n * ctx->dim)); 3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 314552f7358SJed Brown } 315552f7358SJed Brown 3164267b1a3SMatthew G. Knepley /*@C 31752aa1562SMatthew G. Knepley DMInterpolationSetUp - Compute spatial indices for point location during interpolation 3184267b1a3SMatthew G. Knepley 319c3339decSBarry Smith Collective 3204267b1a3SMatthew G. Knepley 3214267b1a3SMatthew G. Knepley Input Parameters: 3224267b1a3SMatthew G. Knepley + ctx - the context 323f6dfbefdSBarry Smith . dm - the `DM` for the function space used for interpolation 324f6dfbefdSBarry Smith . redundantPoints - If `PETSC_TRUE`, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes. 325f6dfbefdSBarry Smith - ignoreOutsideDomain - If `PETSC_TRUE`, ignore points outside the domain, otherwise return an error 3264267b1a3SMatthew G. Knepley 3274267b1a3SMatthew G. Knepley Level: intermediate 3284267b1a3SMatthew G. Knepley 329420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 3304267b1a3SMatthew G. Knepley @*/ 331d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain) 332d71ae5a4SJacob Faibussowitsch { 333552f7358SJed Brown MPI_Comm comm = ctx->comm; 334552f7358SJed Brown PetscScalar *a; 335552f7358SJed Brown PetscInt p, q, i; 336552f7358SJed Brown PetscMPIInt rank, size; 337552f7358SJed Brown Vec pointVec; 3383a93e3b7SToby Isaac PetscSF cellSF; 339552f7358SJed Brown PetscLayout layout; 340552f7358SJed Brown PetscReal *globalPoints; 341cb313848SJed Brown PetscScalar *globalPointsScalar; 342552f7358SJed Brown const PetscInt *ranges; 343552f7358SJed Brown PetscMPIInt *counts, *displs; 3443a93e3b7SToby Isaac const PetscSFNode *foundCells; 3453a93e3b7SToby Isaac const PetscInt *foundPoints; 346552f7358SJed Brown PetscMPIInt *foundProcs, *globalProcs; 3473a93e3b7SToby Isaac PetscInt n, N, numFound; 348552f7358SJed Brown 34919436ca2SJed Brown PetscFunctionBegin; 350064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 3519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 3529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 35308401ef6SPierre Jolivet PetscCheck(ctx->dim >= 0, comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 35419436ca2SJed Brown /* Locate points */ 35519436ca2SJed Brown n = ctx->nInput; 356552f7358SJed Brown if (!redundantPoints) { 3579566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &layout)); 3589566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(layout, 1)); 3599566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(layout, n)); 3609566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(layout)); 3619566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetSize(layout, &N)); 362552f7358SJed Brown /* Communicate all points to all processes */ 3639566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(N * ctx->dim, &globalPoints, size, &counts, size, &displs)); 3649566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRanges(layout, &ranges)); 365552f7358SJed Brown for (p = 0; p < size; ++p) { 366552f7358SJed Brown counts[p] = (ranges[p + 1] - ranges[p]) * ctx->dim; 367552f7358SJed Brown displs[p] = ranges[p] * ctx->dim; 368552f7358SJed Brown } 3699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv(ctx->points, n * ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm)); 370552f7358SJed Brown } else { 371552f7358SJed Brown N = n; 372552f7358SJed Brown globalPoints = ctx->points; 37338ea73c8SJed Brown counts = displs = NULL; 37438ea73c8SJed Brown layout = NULL; 375552f7358SJed Brown } 376552f7358SJed Brown #if 0 3779566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs)); 37819436ca2SJed Brown /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 379552f7358SJed Brown #else 380cb313848SJed Brown #if defined(PETSC_USE_COMPLEX) 3819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * ctx->dim, &globalPointsScalar)); 38246b3086cSToby Isaac for (i = 0; i < N * ctx->dim; i++) globalPointsScalar[i] = globalPoints[i]; 383cb313848SJed Brown #else 384cb313848SJed Brown globalPointsScalar = globalPoints; 385cb313848SJed Brown #endif 3869566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N * ctx->dim, globalPointsScalar, &pointVec)); 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N, &foundProcs, N, &globalProcs)); 388ad540459SPierre Jolivet for (p = 0; p < N; ++p) foundProcs[p] = size; 3893a93e3b7SToby Isaac cellSF = NULL; 3909566063dSJacob Faibussowitsch PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF)); 3919566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(cellSF, NULL, &numFound, &foundPoints, &foundCells)); 392552f7358SJed Brown #endif 3933a93e3b7SToby Isaac for (p = 0; p < numFound; ++p) { 3943a93e3b7SToby Isaac if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank; 395552f7358SJed Brown } 396552f7358SJed Brown /* Let the lowest rank process own each point */ 3971c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm)); 398552f7358SJed Brown ctx->n = 0; 399552f7358SJed Brown for (p = 0; p < N; ++p) { 40052aa1562SMatthew G. Knepley if (globalProcs[p] == size) { 4019371c9d4SSatish Balay PetscCheck(ignoreOutsideDomain, comm, PETSC_ERR_PLIB, "Point %" PetscInt_FMT ": %g %g %g not located in mesh", p, (double)globalPoints[p * ctx->dim + 0], (double)(ctx->dim > 1 ? globalPoints[p * ctx->dim + 1] : 0.0), 4029371c9d4SSatish Balay (double)(ctx->dim > 2 ? globalPoints[p * ctx->dim + 2] : 0.0)); 403f7d195e4SLawrence Mitchell if (rank == 0) ++ctx->n; 40452aa1562SMatthew G. Knepley } else if (globalProcs[p] == rank) ++ctx->n; 405552f7358SJed Brown } 406552f7358SJed Brown /* Create coordinates vector and array of owned cells */ 4079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ctx->n, &ctx->cells)); 4089566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &ctx->coords)); 4099566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ctx->coords, ctx->n * ctx->dim, PETSC_DECIDE)); 4109566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(ctx->coords, ctx->dim)); 4119566063dSJacob Faibussowitsch PetscCall(VecSetType(ctx->coords, VECSTANDARD)); 4129566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->coords, &a)); 413552f7358SJed Brown for (p = 0, q = 0, i = 0; p < N; ++p) { 414552f7358SJed Brown if (globalProcs[p] == rank) { 415552f7358SJed Brown PetscInt d; 416552f7358SJed Brown 4171aa26658SKarl Rupp for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p * ctx->dim + d]; 418f317b747SMatthew G. Knepley ctx->cells[q] = foundCells[q].index; 419f317b747SMatthew G. Knepley ++q; 420552f7358SJed Brown } 421dd400576SPatrick Sanan if (globalProcs[p] == size && rank == 0) { 42252aa1562SMatthew G. Knepley PetscInt d; 42352aa1562SMatthew G. Knepley 42452aa1562SMatthew G. Knepley for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.; 42552aa1562SMatthew G. Knepley ctx->cells[q] = -1; 42652aa1562SMatthew G. Knepley ++q; 42752aa1562SMatthew G. Knepley } 428552f7358SJed Brown } 4299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->coords, &a)); 430552f7358SJed Brown #if 0 4319566063dSJacob Faibussowitsch PetscCall(PetscFree3(foundCells,foundProcs,globalProcs)); 432552f7358SJed Brown #else 4339566063dSJacob Faibussowitsch PetscCall(PetscFree2(foundProcs, globalProcs)); 4349566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&cellSF)); 4359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pointVec)); 436552f7358SJed Brown #endif 4379566063dSJacob Faibussowitsch if ((void *)globalPointsScalar != (void *)globalPoints) PetscCall(PetscFree(globalPointsScalar)); 4389566063dSJacob Faibussowitsch if (!redundantPoints) PetscCall(PetscFree3(globalPoints, counts, displs)); 4399566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&layout)); 4403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 441552f7358SJed Brown } 442552f7358SJed Brown 4434267b1a3SMatthew G. Knepley /*@C 444f6dfbefdSBarry Smith DMInterpolationGetCoordinates - Gets a `Vec` with the coordinates of each interpolation point 4454267b1a3SMatthew G. Knepley 446c3339decSBarry Smith Collective 4474267b1a3SMatthew G. Knepley 4484267b1a3SMatthew G. Knepley Input Parameter: 4494267b1a3SMatthew G. Knepley . ctx - the context 4504267b1a3SMatthew G. Knepley 4514267b1a3SMatthew G. Knepley Output Parameter: 4524267b1a3SMatthew G. Knepley . coordinates - the coordinates of interpolation points 4534267b1a3SMatthew G. Knepley 45420f4b53cSBarry Smith Level: intermediate 45520f4b53cSBarry Smith 456f6dfbefdSBarry Smith Note: 457f6dfbefdSBarry Smith The local vector entries correspond to interpolation points lying on this process, according to the associated `DM`. 458f6dfbefdSBarry Smith This is a borrowed vector that the user should not destroy. 4594267b1a3SMatthew G. Knepley 460420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 4614267b1a3SMatthew G. Knepley @*/ 462d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 463d71ae5a4SJacob Faibussowitsch { 464552f7358SJed Brown PetscFunctionBegin; 4654f572ea9SToby Isaac PetscAssertPointer(coordinates, 2); 46628b400f6SJacob Faibussowitsch PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 467552f7358SJed Brown *coordinates = ctx->coords; 4683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 469552f7358SJed Brown } 470552f7358SJed Brown 4714267b1a3SMatthew G. Knepley /*@C 472f6dfbefdSBarry Smith DMInterpolationGetVector - Gets a `Vec` which can hold all the interpolated field values 4734267b1a3SMatthew G. Knepley 474c3339decSBarry Smith Collective 4754267b1a3SMatthew G. Knepley 4764267b1a3SMatthew G. Knepley Input Parameter: 4774267b1a3SMatthew G. Knepley . ctx - the context 4784267b1a3SMatthew G. Knepley 4794267b1a3SMatthew G. Knepley Output Parameter: 4804267b1a3SMatthew G. Knepley . v - a vector capable of holding the interpolated field values 4814267b1a3SMatthew G. Knepley 4822fe279fdSBarry Smith Level: intermediate 4832fe279fdSBarry Smith 484f6dfbefdSBarry Smith Note: 485f6dfbefdSBarry Smith This vector should be returned using `DMInterpolationRestoreVector()`. 4864267b1a3SMatthew G. Knepley 487420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 4884267b1a3SMatthew G. Knepley @*/ 489d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 490d71ae5a4SJacob Faibussowitsch { 491552f7358SJed Brown PetscFunctionBegin; 4924f572ea9SToby Isaac PetscAssertPointer(v, 2); 49328b400f6SJacob Faibussowitsch PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 4949566063dSJacob Faibussowitsch PetscCall(VecCreate(ctx->comm, v)); 4959566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*v, ctx->n * ctx->dof, PETSC_DECIDE)); 4969566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*v, ctx->dof)); 4979566063dSJacob Faibussowitsch PetscCall(VecSetType(*v, VECSTANDARD)); 4983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 499552f7358SJed Brown } 500552f7358SJed Brown 5014267b1a3SMatthew G. Knepley /*@C 502f6dfbefdSBarry Smith DMInterpolationRestoreVector - Returns a `Vec` which can hold all the interpolated field values 5034267b1a3SMatthew G. Knepley 504c3339decSBarry Smith Collective 5054267b1a3SMatthew G. Knepley 5064267b1a3SMatthew G. Knepley Input Parameters: 5074267b1a3SMatthew G. Knepley + ctx - the context 5084267b1a3SMatthew G. Knepley - v - a vector capable of holding the interpolated field values 5094267b1a3SMatthew G. Knepley 5104267b1a3SMatthew G. Knepley Level: intermediate 5114267b1a3SMatthew G. Knepley 512420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 5134267b1a3SMatthew G. Knepley @*/ 514d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 515d71ae5a4SJacob Faibussowitsch { 516552f7358SJed Brown PetscFunctionBegin; 5174f572ea9SToby Isaac PetscAssertPointer(v, 2); 51828b400f6SJacob Faibussowitsch PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 5199566063dSJacob Faibussowitsch PetscCall(VecDestroy(v)); 5203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 521552f7358SJed Brown } 522552f7358SJed Brown 523d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 524d71ae5a4SJacob Faibussowitsch { 525cfe5744fSMatthew G. Knepley PetscReal v0, J, invJ, detJ; 526cfe5744fSMatthew G. Knepley const PetscInt dof = ctx->dof; 527cfe5744fSMatthew G. Knepley const PetscScalar *coords; 528cfe5744fSMatthew G. Knepley PetscScalar *a; 529cfe5744fSMatthew G. Knepley PetscInt p; 530cfe5744fSMatthew G. Knepley 531cfe5744fSMatthew G. Knepley PetscFunctionBegin; 5329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->coords, &coords)); 5339566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &a)); 534cfe5744fSMatthew G. Knepley for (p = 0; p < ctx->n; ++p) { 535cfe5744fSMatthew G. Knepley PetscInt c = ctx->cells[p]; 536cfe5744fSMatthew G. Knepley PetscScalar *x = NULL; 537cfe5744fSMatthew G. Knepley PetscReal xir[1]; 538cfe5744fSMatthew G. Knepley PetscInt xSize, comp; 539cfe5744fSMatthew G. Knepley 5409566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ)); 541cfe5744fSMatthew G. Knepley PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 542cfe5744fSMatthew G. Knepley xir[0] = invJ * PetscRealPart(coords[p] - v0); 5439566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 544cfe5744fSMatthew G. Knepley if (2 * dof == xSize) { 545cfe5744fSMatthew G. Knepley for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) + x[1 * dof + comp] * xir[0]; 546cfe5744fSMatthew G. Knepley } else if (dof == xSize) { 547cfe5744fSMatthew G. Knepley for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp]; 548cfe5744fSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input closure size %" PetscInt_FMT " must be either %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 2 * dof, dof); 5499566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 550cfe5744fSMatthew G. Knepley } 5519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &a)); 5529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 5533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 554cfe5744fSMatthew G. Knepley } 555cfe5744fSMatthew G. Knepley 556d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 557d71ae5a4SJacob Faibussowitsch { 558552f7358SJed Brown PetscReal *v0, *J, *invJ, detJ; 55956044e6dSMatthew G. Knepley const PetscScalar *coords; 56056044e6dSMatthew G. Knepley PetscScalar *a; 561552f7358SJed Brown PetscInt p; 562552f7358SJed Brown 563552f7358SJed Brown PetscFunctionBegin; 5649566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ)); 5659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->coords, &coords)); 5669566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &a)); 567552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 568552f7358SJed Brown PetscInt c = ctx->cells[p]; 569a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 570552f7358SJed Brown PetscReal xi[4]; 571552f7358SJed Brown PetscInt d, f, comp; 572552f7358SJed Brown 5739566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 57463a3b9bcSJacob Faibussowitsch PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 5759566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x)); 5761aa26658SKarl Rupp for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 5771aa26658SKarl Rupp 578552f7358SJed Brown for (d = 0; d < ctx->dim; ++d) { 579552f7358SJed Brown xi[d] = 0.0; 5801aa26658SKarl Rupp for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]); 5811aa26658SKarl Rupp for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] += PetscRealPart(x[(d + 1) * ctx->dof + comp] - x[0 * ctx->dof + comp]) * xi[d]; 582552f7358SJed Brown } 5839566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x)); 584552f7358SJed Brown } 5859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &a)); 5869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 5879566063dSJacob Faibussowitsch PetscCall(PetscFree3(v0, J, invJ)); 5883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 589552f7358SJed Brown } 590552f7358SJed Brown 591d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 592d71ae5a4SJacob Faibussowitsch { 5937a1931ceSMatthew G. Knepley PetscReal *v0, *J, *invJ, detJ; 59456044e6dSMatthew G. Knepley const PetscScalar *coords; 59556044e6dSMatthew G. Knepley PetscScalar *a; 5967a1931ceSMatthew G. Knepley PetscInt p; 5977a1931ceSMatthew G. Knepley 5987a1931ceSMatthew G. Knepley PetscFunctionBegin; 5999566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ)); 6009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->coords, &coords)); 6019566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &a)); 6027a1931ceSMatthew G. Knepley for (p = 0; p < ctx->n; ++p) { 6037a1931ceSMatthew G. Knepley PetscInt c = ctx->cells[p]; 6047a1931ceSMatthew G. Knepley const PetscInt order[3] = {2, 1, 3}; 6052584bbe8SMatthew G. Knepley PetscScalar *x = NULL; 6067a1931ceSMatthew G. Knepley PetscReal xi[4]; 6077a1931ceSMatthew G. Knepley PetscInt d, f, comp; 6087a1931ceSMatthew G. Knepley 6099566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 61063a3b9bcSJacob Faibussowitsch PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 6119566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x)); 6127a1931ceSMatthew G. Knepley for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 6137a1931ceSMatthew G. Knepley 6147a1931ceSMatthew G. Knepley for (d = 0; d < ctx->dim; ++d) { 6157a1931ceSMatthew G. Knepley xi[d] = 0.0; 6167a1931ceSMatthew G. Knepley for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]); 6177a1931ceSMatthew G. Knepley for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] += PetscRealPart(x[order[d] * ctx->dof + comp] - x[0 * ctx->dof + comp]) * xi[d]; 6187a1931ceSMatthew G. Knepley } 6199566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x)); 6207a1931ceSMatthew G. Knepley } 6219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &a)); 6229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 6239566063dSJacob Faibussowitsch PetscCall(PetscFree3(v0, J, invJ)); 6243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6257a1931ceSMatthew G. Knepley } 6267a1931ceSMatthew G. Knepley 627d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 628d71ae5a4SJacob Faibussowitsch { 629552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar *)ctx; 630552f7358SJed Brown const PetscScalar x0 = vertices[0]; 631552f7358SJed Brown const PetscScalar y0 = vertices[1]; 632552f7358SJed Brown const PetscScalar x1 = vertices[2]; 633552f7358SJed Brown const PetscScalar y1 = vertices[3]; 634552f7358SJed Brown const PetscScalar x2 = vertices[4]; 635552f7358SJed Brown const PetscScalar y2 = vertices[5]; 636552f7358SJed Brown const PetscScalar x3 = vertices[6]; 637552f7358SJed Brown const PetscScalar y3 = vertices[7]; 638552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 639552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 640552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 641552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 642552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 643552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 64456044e6dSMatthew G. Knepley const PetscScalar *ref; 64556044e6dSMatthew G. Knepley PetscScalar *real; 646552f7358SJed Brown 647552f7358SJed Brown PetscFunctionBegin; 6489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xref, &ref)); 6499566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xreal, &real)); 650552f7358SJed Brown { 651552f7358SJed Brown const PetscScalar p0 = ref[0]; 652552f7358SJed Brown const PetscScalar p1 = ref[1]; 653552f7358SJed Brown 654552f7358SJed Brown real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 655552f7358SJed Brown real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 656552f7358SJed Brown } 6579566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(28)); 6589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xref, &ref)); 6599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xreal, &real)); 6603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 661552f7358SJed Brown } 662552f7358SJed Brown 663af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 664d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 665d71ae5a4SJacob Faibussowitsch { 666552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar *)ctx; 667552f7358SJed Brown const PetscScalar x0 = vertices[0]; 668552f7358SJed Brown const PetscScalar y0 = vertices[1]; 669552f7358SJed Brown const PetscScalar x1 = vertices[2]; 670552f7358SJed Brown const PetscScalar y1 = vertices[3]; 671552f7358SJed Brown const PetscScalar x2 = vertices[4]; 672552f7358SJed Brown const PetscScalar y2 = vertices[5]; 673552f7358SJed Brown const PetscScalar x3 = vertices[6]; 674552f7358SJed Brown const PetscScalar y3 = vertices[7]; 675552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 676552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 67756044e6dSMatthew G. Knepley const PetscScalar *ref; 678552f7358SJed Brown 679552f7358SJed Brown PetscFunctionBegin; 6809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xref, &ref)); 681552f7358SJed Brown { 682552f7358SJed Brown const PetscScalar x = ref[0]; 683552f7358SJed Brown const PetscScalar y = ref[1]; 684552f7358SJed Brown const PetscInt rows[2] = {0, 1}; 685da80777bSKarl Rupp PetscScalar values[4]; 686da80777bSKarl Rupp 6879371c9d4SSatish Balay values[0] = (x1 - x0 + f_01 * y) * 0.5; 6889371c9d4SSatish Balay values[1] = (x3 - x0 + f_01 * x) * 0.5; 6899371c9d4SSatish Balay values[2] = (y1 - y0 + g_01 * y) * 0.5; 6909371c9d4SSatish Balay values[3] = (y3 - y0 + g_01 * x) * 0.5; 6919566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES)); 692552f7358SJed Brown } 6939566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(30)); 6949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xref, &ref)); 6959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 6969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 6973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 698552f7358SJed Brown } 699552f7358SJed Brown 700d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 701d71ae5a4SJacob Faibussowitsch { 702fafc0619SMatthew G Knepley DM dmCoord; 703de73a395SMatthew G. Knepley PetscFE fem = NULL; 704552f7358SJed Brown SNES snes; 705552f7358SJed Brown KSP ksp; 706552f7358SJed Brown PC pc; 707552f7358SJed Brown Vec coordsLocal, r, ref, real; 708552f7358SJed Brown Mat J; 709716009bfSMatthew G. Knepley PetscTabulation T = NULL; 71056044e6dSMatthew G. Knepley const PetscScalar *coords; 71156044e6dSMatthew G. Knepley PetscScalar *a; 712716009bfSMatthew G. Knepley PetscReal xir[2] = {0., 0.}; 713de73a395SMatthew G. Knepley PetscInt Nf, p; 7145509d985SMatthew G. Knepley const PetscInt dof = ctx->dof; 715552f7358SJed Brown 716552f7358SJed Brown PetscFunctionBegin; 7179566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 718716009bfSMatthew G. Knepley if (Nf) { 719cfe5744fSMatthew G. Knepley PetscObject obj; 720cfe5744fSMatthew G. Knepley PetscClassId id; 721cfe5744fSMatthew G. Knepley 7229566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &obj)); 7239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 7249371c9d4SSatish Balay if (id == PETSCFE_CLASSID) { 7259371c9d4SSatish Balay fem = (PetscFE)obj; 7269371c9d4SSatish Balay PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T)); 7279371c9d4SSatish Balay } 728716009bfSMatthew G. Knepley } 7299566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 7309566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 7319566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 7329566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_")); 7339566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &r)); 7349566063dSJacob Faibussowitsch PetscCall(VecSetSizes(r, 2, 2)); 7359566063dSJacob Faibussowitsch PetscCall(VecSetType(r, dm->vectype)); 7369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r, &ref)); 7379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r, &real)); 7389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 7399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, 2, 2, 2, 2)); 7409566063dSJacob Faibussowitsch PetscCall(MatSetType(J, MATSEQDENSE)); 7419566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 7429566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL)); 7439566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL)); 7449566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 7459566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 7469566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCLU)); 7479566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 748552f7358SJed Brown 7499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->coords, &coords)); 7509566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &a)); 751552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 752a1e44745SMatthew G. Knepley PetscScalar *x = NULL, *vertices = NULL; 753552f7358SJed Brown PetscScalar *xi; 754552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 755552f7358SJed Brown 756552f7358SJed Brown /* Can make this do all points at once */ 7579566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 75863a3b9bcSJacob Faibussowitsch PetscCheck(4 * 2 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4 * 2); 7599566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 7609566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, vertices)); 7619566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices)); 7629566063dSJacob Faibussowitsch PetscCall(VecGetArray(real, &xi)); 763552f7358SJed Brown xi[0] = coords[p * ctx->dim + 0]; 764552f7358SJed Brown xi[1] = coords[p * ctx->dim + 1]; 7659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(real, &xi)); 7669566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, real, ref)); 7679566063dSJacob Faibussowitsch PetscCall(VecGetArray(ref, &xi)); 768cb313848SJed Brown xir[0] = PetscRealPart(xi[0]); 769cb313848SJed Brown xir[1] = PetscRealPart(xi[1]); 770cfe5744fSMatthew G. Knepley if (4 * dof == xSize) { 7719371c9d4SSatish Balay for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) * (1 - xir[1]) + x[1 * dof + comp] * xir[0] * (1 - xir[1]) + x[2 * dof + comp] * xir[0] * xir[1] + x[3 * dof + comp] * (1 - xir[0]) * xir[1]; 772cfe5744fSMatthew G. Knepley } else if (dof == xSize) { 773cfe5744fSMatthew G. Knepley for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp]; 774cfe5744fSMatthew G. Knepley } else { 7755509d985SMatthew G. Knepley PetscInt d; 7761aa26658SKarl Rupp 7772c71b3e2SJacob Faibussowitsch PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE"); 7789371c9d4SSatish Balay xir[0] = 2.0 * xir[0] - 1.0; 7799371c9d4SSatish Balay xir[1] = 2.0 * xir[1] - 1.0; 7809566063dSJacob Faibussowitsch PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T)); 7815509d985SMatthew G. Knepley for (comp = 0; comp < dof; ++comp) { 7825509d985SMatthew G. Knepley a[p * dof + comp] = 0.0; 783ad540459SPierre Jolivet for (d = 0; d < xSize / dof; ++d) a[p * dof + comp] += x[d * dof + comp] * T->T[0][d * dof + comp]; 7845509d985SMatthew G. Knepley } 7855509d985SMatthew G. Knepley } 7869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ref, &xi)); 7879566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 7889566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 789552f7358SJed Brown } 7909566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&T)); 7919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &a)); 7929566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 793552f7358SJed Brown 7949566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 7959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 7969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ref)); 7979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&real)); 7989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 7993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 800552f7358SJed Brown } 801552f7358SJed Brown 802d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 803d71ae5a4SJacob Faibussowitsch { 804552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar *)ctx; 805552f7358SJed Brown const PetscScalar x0 = vertices[0]; 806552f7358SJed Brown const PetscScalar y0 = vertices[1]; 807552f7358SJed Brown const PetscScalar z0 = vertices[2]; 8087a1931ceSMatthew G. Knepley const PetscScalar x1 = vertices[9]; 8097a1931ceSMatthew G. Knepley const PetscScalar y1 = vertices[10]; 8107a1931ceSMatthew G. Knepley const PetscScalar z1 = vertices[11]; 811552f7358SJed Brown const PetscScalar x2 = vertices[6]; 812552f7358SJed Brown const PetscScalar y2 = vertices[7]; 813552f7358SJed Brown const PetscScalar z2 = vertices[8]; 8147a1931ceSMatthew G. Knepley const PetscScalar x3 = vertices[3]; 8157a1931ceSMatthew G. Knepley const PetscScalar y3 = vertices[4]; 8167a1931ceSMatthew G. Knepley const PetscScalar z3 = vertices[5]; 817552f7358SJed Brown const PetscScalar x4 = vertices[12]; 818552f7358SJed Brown const PetscScalar y4 = vertices[13]; 819552f7358SJed Brown const PetscScalar z4 = vertices[14]; 820552f7358SJed Brown const PetscScalar x5 = vertices[15]; 821552f7358SJed Brown const PetscScalar y5 = vertices[16]; 822552f7358SJed Brown const PetscScalar z5 = vertices[17]; 823552f7358SJed Brown const PetscScalar x6 = vertices[18]; 824552f7358SJed Brown const PetscScalar y6 = vertices[19]; 825552f7358SJed Brown const PetscScalar z6 = vertices[20]; 826552f7358SJed Brown const PetscScalar x7 = vertices[21]; 827552f7358SJed Brown const PetscScalar y7 = vertices[22]; 828552f7358SJed Brown const PetscScalar z7 = vertices[23]; 829552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 830552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 831552f7358SJed Brown const PetscScalar h_1 = z1 - z0; 832552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 833552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 834552f7358SJed Brown const PetscScalar h_3 = z3 - z0; 835552f7358SJed Brown const PetscScalar f_4 = x4 - x0; 836552f7358SJed Brown const PetscScalar g_4 = y4 - y0; 837552f7358SJed Brown const PetscScalar h_4 = z4 - z0; 838552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 839552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 840552f7358SJed Brown const PetscScalar h_01 = z2 - z1 - z3 + z0; 841552f7358SJed Brown const PetscScalar f_12 = x7 - x3 - x4 + x0; 842552f7358SJed Brown const PetscScalar g_12 = y7 - y3 - y4 + y0; 843552f7358SJed Brown const PetscScalar h_12 = z7 - z3 - z4 + z0; 844552f7358SJed Brown const PetscScalar f_02 = x5 - x1 - x4 + x0; 845552f7358SJed Brown const PetscScalar g_02 = y5 - y1 - y4 + y0; 846552f7358SJed Brown const PetscScalar h_02 = z5 - z1 - z4 + z0; 847552f7358SJed Brown const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 848552f7358SJed Brown const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 849552f7358SJed Brown const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 85056044e6dSMatthew G. Knepley const PetscScalar *ref; 85156044e6dSMatthew G. Knepley PetscScalar *real; 852552f7358SJed Brown 853552f7358SJed Brown PetscFunctionBegin; 8549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xref, &ref)); 8559566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xreal, &real)); 856552f7358SJed Brown { 857552f7358SJed Brown const PetscScalar p0 = ref[0]; 858552f7358SJed Brown const PetscScalar p1 = ref[1]; 859552f7358SJed Brown const PetscScalar p2 = ref[2]; 860552f7358SJed Brown 861552f7358SJed Brown real[0] = x0 + f_1 * p0 + f_3 * p1 + f_4 * p2 + f_01 * p0 * p1 + f_12 * p1 * p2 + f_02 * p0 * p2 + f_012 * p0 * p1 * p2; 862552f7358SJed Brown real[1] = y0 + g_1 * p0 + g_3 * p1 + g_4 * p2 + g_01 * p0 * p1 + g_01 * p0 * p1 + g_12 * p1 * p2 + g_02 * p0 * p2 + g_012 * p0 * p1 * p2; 863552f7358SJed Brown real[2] = z0 + h_1 * p0 + h_3 * p1 + h_4 * p2 + h_01 * p0 * p1 + h_01 * p0 * p1 + h_12 * p1 * p2 + h_02 * p0 * p2 + h_012 * p0 * p1 * p2; 864552f7358SJed Brown } 8659566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(114)); 8669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xref, &ref)); 8679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xreal, &real)); 8683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 869552f7358SJed Brown } 870552f7358SJed Brown 871d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 872d71ae5a4SJacob Faibussowitsch { 873552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar *)ctx; 874552f7358SJed Brown const PetscScalar x0 = vertices[0]; 875552f7358SJed Brown const PetscScalar y0 = vertices[1]; 876552f7358SJed Brown const PetscScalar z0 = vertices[2]; 8777a1931ceSMatthew G. Knepley const PetscScalar x1 = vertices[9]; 8787a1931ceSMatthew G. Knepley const PetscScalar y1 = vertices[10]; 8797a1931ceSMatthew G. Knepley const PetscScalar z1 = vertices[11]; 880552f7358SJed Brown const PetscScalar x2 = vertices[6]; 881552f7358SJed Brown const PetscScalar y2 = vertices[7]; 882552f7358SJed Brown const PetscScalar z2 = vertices[8]; 8837a1931ceSMatthew G. Knepley const PetscScalar x3 = vertices[3]; 8847a1931ceSMatthew G. Knepley const PetscScalar y3 = vertices[4]; 8857a1931ceSMatthew G. Knepley const PetscScalar z3 = vertices[5]; 886552f7358SJed Brown const PetscScalar x4 = vertices[12]; 887552f7358SJed Brown const PetscScalar y4 = vertices[13]; 888552f7358SJed Brown const PetscScalar z4 = vertices[14]; 889552f7358SJed Brown const PetscScalar x5 = vertices[15]; 890552f7358SJed Brown const PetscScalar y5 = vertices[16]; 891552f7358SJed Brown const PetscScalar z5 = vertices[17]; 892552f7358SJed Brown const PetscScalar x6 = vertices[18]; 893552f7358SJed Brown const PetscScalar y6 = vertices[19]; 894552f7358SJed Brown const PetscScalar z6 = vertices[20]; 895552f7358SJed Brown const PetscScalar x7 = vertices[21]; 896552f7358SJed Brown const PetscScalar y7 = vertices[22]; 897552f7358SJed Brown const PetscScalar z7 = vertices[23]; 898552f7358SJed Brown const PetscScalar f_xy = x2 - x1 - x3 + x0; 899552f7358SJed Brown const PetscScalar g_xy = y2 - y1 - y3 + y0; 900552f7358SJed Brown const PetscScalar h_xy = z2 - z1 - z3 + z0; 901552f7358SJed Brown const PetscScalar f_yz = x7 - x3 - x4 + x0; 902552f7358SJed Brown const PetscScalar g_yz = y7 - y3 - y4 + y0; 903552f7358SJed Brown const PetscScalar h_yz = z7 - z3 - z4 + z0; 904552f7358SJed Brown const PetscScalar f_xz = x5 - x1 - x4 + x0; 905552f7358SJed Brown const PetscScalar g_xz = y5 - y1 - y4 + y0; 906552f7358SJed Brown const PetscScalar h_xz = z5 - z1 - z4 + z0; 907552f7358SJed Brown const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 908552f7358SJed Brown const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 909552f7358SJed Brown const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 91056044e6dSMatthew G. Knepley const PetscScalar *ref; 911552f7358SJed Brown 912552f7358SJed Brown PetscFunctionBegin; 9139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xref, &ref)); 914552f7358SJed Brown { 915552f7358SJed Brown const PetscScalar x = ref[0]; 916552f7358SJed Brown const PetscScalar y = ref[1]; 917552f7358SJed Brown const PetscScalar z = ref[2]; 918552f7358SJed Brown const PetscInt rows[3] = {0, 1, 2}; 919da80777bSKarl Rupp PetscScalar values[9]; 920da80777bSKarl Rupp 921da80777bSKarl Rupp values[0] = (x1 - x0 + f_xy * y + f_xz * z + f_xyz * y * z) / 2.0; 922da80777bSKarl Rupp values[1] = (x3 - x0 + f_xy * x + f_yz * z + f_xyz * x * z) / 2.0; 923da80777bSKarl Rupp values[2] = (x4 - x0 + f_yz * y + f_xz * x + f_xyz * x * y) / 2.0; 924da80777bSKarl Rupp values[3] = (y1 - y0 + g_xy * y + g_xz * z + g_xyz * y * z) / 2.0; 925da80777bSKarl Rupp values[4] = (y3 - y0 + g_xy * x + g_yz * z + g_xyz * x * z) / 2.0; 926da80777bSKarl Rupp values[5] = (y4 - y0 + g_yz * y + g_xz * x + g_xyz * x * y) / 2.0; 927da80777bSKarl Rupp values[6] = (z1 - z0 + h_xy * y + h_xz * z + h_xyz * y * z) / 2.0; 928da80777bSKarl Rupp values[7] = (z3 - z0 + h_xy * x + h_yz * z + h_xyz * x * z) / 2.0; 929da80777bSKarl Rupp values[8] = (z4 - z0 + h_yz * y + h_xz * x + h_xyz * x * y) / 2.0; 9301aa26658SKarl Rupp 9319566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES)); 932552f7358SJed Brown } 9339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(152)); 9349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xref, &ref)); 9359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 9369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 9373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 938552f7358SJed Brown } 939552f7358SJed Brown 940d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 941d71ae5a4SJacob Faibussowitsch { 942fafc0619SMatthew G Knepley DM dmCoord; 943552f7358SJed Brown SNES snes; 944552f7358SJed Brown KSP ksp; 945552f7358SJed Brown PC pc; 946552f7358SJed Brown Vec coordsLocal, r, ref, real; 947552f7358SJed Brown Mat J; 94856044e6dSMatthew G. Knepley const PetscScalar *coords; 94956044e6dSMatthew G. Knepley PetscScalar *a; 950552f7358SJed Brown PetscInt p; 951552f7358SJed Brown 952552f7358SJed Brown PetscFunctionBegin; 9539566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 9549566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 9559566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 9569566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_")); 9579566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &r)); 9589566063dSJacob Faibussowitsch PetscCall(VecSetSizes(r, 3, 3)); 9599566063dSJacob Faibussowitsch PetscCall(VecSetType(r, dm->vectype)); 9609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r, &ref)); 9619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(r, &real)); 9629566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 9639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, 3, 3, 3, 3)); 9649566063dSJacob Faibussowitsch PetscCall(MatSetType(J, MATSEQDENSE)); 9659566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 9669566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL)); 9679566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL)); 9689566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 9699566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 9709566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCLU)); 9719566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 972552f7358SJed Brown 9739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->coords, &coords)); 9749566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &a)); 975552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 976a1e44745SMatthew G. Knepley PetscScalar *x = NULL, *vertices = NULL; 977552f7358SJed Brown PetscScalar *xi; 978cb313848SJed Brown PetscReal xir[3]; 979552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 980552f7358SJed Brown 981552f7358SJed Brown /* Can make this do all points at once */ 9829566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 983cfe5744fSMatthew G. Knepley PetscCheck(8 * 3 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8 * 3); 9849566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 985cfe5744fSMatthew G. Knepley PetscCheck((8 * ctx->dof == xSize) || (ctx->dof == xSize), ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input closure size %" PetscInt_FMT " should be %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 8 * ctx->dof, ctx->dof); 9869566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, NULL, NULL, vertices)); 9879566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices)); 9889566063dSJacob Faibussowitsch PetscCall(VecGetArray(real, &xi)); 989552f7358SJed Brown xi[0] = coords[p * ctx->dim + 0]; 990552f7358SJed Brown xi[1] = coords[p * ctx->dim + 1]; 991552f7358SJed Brown xi[2] = coords[p * ctx->dim + 2]; 9929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(real, &xi)); 9939566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, real, ref)); 9949566063dSJacob Faibussowitsch PetscCall(VecGetArray(ref, &xi)); 995cb313848SJed Brown xir[0] = PetscRealPart(xi[0]); 996cb313848SJed Brown xir[1] = PetscRealPart(xi[1]); 997cb313848SJed Brown xir[2] = PetscRealPart(xi[2]); 998cfe5744fSMatthew G. Knepley if (8 * ctx->dof == xSize) { 999552f7358SJed Brown for (comp = 0; comp < ctx->dof; ++comp) { 10009371c9d4SSatish Balay a[p * ctx->dof + comp] = x[0 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * (1 - xir[2]) + x[3 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * (1 - xir[2]) + x[2 * ctx->dof + comp] * xir[0] * xir[1] * (1 - xir[2]) + x[1 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * (1 - xir[2]) + 10019371c9d4SSatish Balay x[4 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * xir[2] + x[5 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * xir[2] + x[6 * ctx->dof + comp] * xir[0] * xir[1] * xir[2] + x[7 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * xir[2]; 1002552f7358SJed Brown } 1003cfe5744fSMatthew G. Knepley } else { 1004cfe5744fSMatthew G. Knepley for (comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 1005cfe5744fSMatthew G. Knepley } 10069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ref, &xi)); 10079566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 10089566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 1009552f7358SJed Brown } 10109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &a)); 10119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 1012552f7358SJed Brown 10139566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 10149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 10159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ref)); 10169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&real)); 10179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 10183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1019552f7358SJed Brown } 1020552f7358SJed Brown 10214267b1a3SMatthew G. Knepley /*@C 10224267b1a3SMatthew G. Knepley DMInterpolationEvaluate - Using the input from dm and x, calculates interpolated field values at the interpolation points. 10234267b1a3SMatthew G. Knepley 1024552f7358SJed Brown Input Parameters: 1025f6dfbefdSBarry Smith + ctx - The `DMInterpolationInfo` context 1026f6dfbefdSBarry Smith . dm - The `DM` 1027552f7358SJed Brown - x - The local vector containing the field to be interpolated 1028552f7358SJed Brown 10292fe279fdSBarry Smith Output Parameter: 1030552f7358SJed Brown . v - The vector containing the interpolated values 10314267b1a3SMatthew G. Knepley 10324267b1a3SMatthew G. Knepley Level: beginner 10334267b1a3SMatthew G. Knepley 10342fe279fdSBarry Smith Note: 10352fe279fdSBarry Smith A suitable `v` can be obtained using `DMInterpolationGetVector()`. 10362fe279fdSBarry Smith 1037420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 10384267b1a3SMatthew G. Knepley @*/ 1039d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 1040d71ae5a4SJacob Faibussowitsch { 104166a0a883SMatthew G. Knepley PetscDS ds; 104266a0a883SMatthew G. Knepley PetscInt n, p, Nf, field; 104366a0a883SMatthew G. Knepley PetscBool useDS = PETSC_FALSE; 1044552f7358SJed Brown 1045552f7358SJed Brown PetscFunctionBegin; 1046552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1047552f7358SJed Brown PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 1048552f7358SJed Brown PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 10499566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 105063a3b9bcSJacob Faibussowitsch PetscCheck(n == ctx->n * ctx->dof, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %" PetscInt_FMT " should be %" PetscInt_FMT, n, ctx->n * ctx->dof); 10513ba16761SJacob Faibussowitsch if (!n) PetscFunctionReturn(PETSC_SUCCESS); 10529566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 1053680d18d5SMatthew G. Knepley if (ds) { 105466a0a883SMatthew G. Knepley useDS = PETSC_TRUE; 10559566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &Nf)); 1056680d18d5SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 105766a0a883SMatthew G. Knepley PetscObject obj; 1058680d18d5SMatthew G. Knepley PetscClassId id; 1059680d18d5SMatthew G. Knepley 10609566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, &obj)); 10619566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(obj, &id)); 10629371c9d4SSatish Balay if (id != PETSCFE_CLASSID) { 10639371c9d4SSatish Balay useDS = PETSC_FALSE; 10649371c9d4SSatish Balay break; 10659371c9d4SSatish Balay } 106666a0a883SMatthew G. Knepley } 106766a0a883SMatthew G. Knepley } 106866a0a883SMatthew G. Knepley if (useDS) { 106966a0a883SMatthew G. Knepley const PetscScalar *coords; 107066a0a883SMatthew G. Knepley PetscScalar *interpolant; 107166a0a883SMatthew G. Knepley PetscInt cdim, d; 107266a0a883SMatthew G. Knepley 10739566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 10749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->coords, &coords)); 10759566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(v, &interpolant)); 1076680d18d5SMatthew G. Knepley for (p = 0; p < ctx->n; ++p) { 107766a0a883SMatthew G. Knepley PetscReal pcoords[3], xi[3]; 1078680d18d5SMatthew G. Knepley PetscScalar *xa = NULL; 107966a0a883SMatthew G. Knepley PetscInt coff = 0, foff = 0, clSize; 1080680d18d5SMatthew G. Knepley 108152aa1562SMatthew G. Knepley if (ctx->cells[p] < 0) continue; 1082680d18d5SMatthew G. Knepley for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p * cdim + d]); 10839566063dSJacob Faibussowitsch PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi)); 10849566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa)); 108566a0a883SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 108666a0a883SMatthew G. Knepley PetscTabulation T; 108766a0a883SMatthew G. Knepley PetscFE fe; 108866a0a883SMatthew G. Knepley 10899566063dSJacob Faibussowitsch PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 10909566063dSJacob Faibussowitsch PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T)); 1091680d18d5SMatthew G. Knepley { 1092680d18d5SMatthew G. Knepley const PetscReal *basis = T->T[0]; 1093680d18d5SMatthew G. Knepley const PetscInt Nb = T->Nb; 1094680d18d5SMatthew G. Knepley const PetscInt Nc = T->Nc; 109566a0a883SMatthew G. Knepley PetscInt f, fc; 109666a0a883SMatthew G. Knepley 1097680d18d5SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 109866a0a883SMatthew G. Knepley interpolant[p * ctx->dof + coff + fc] = 0.0; 1099ad540459SPierre Jolivet for (f = 0; f < Nb; ++f) interpolant[p * ctx->dof + coff + fc] += xa[foff + f] * basis[(0 * Nb + f) * Nc + fc]; 1100680d18d5SMatthew G. Knepley } 110166a0a883SMatthew G. Knepley coff += Nc; 110266a0a883SMatthew G. Knepley foff += Nb; 1103680d18d5SMatthew G. Knepley } 11049566063dSJacob Faibussowitsch PetscCall(PetscTabulationDestroy(&T)); 1105680d18d5SMatthew G. Knepley } 11069566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa)); 110763a3b9bcSJacob Faibussowitsch PetscCheck(coff == ctx->dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof); 110863a3b9bcSJacob Faibussowitsch PetscCheck(foff == clSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize); 110966a0a883SMatthew G. Knepley } 11109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 11119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(v, &interpolant)); 111266a0a883SMatthew G. Knepley } else { 111366a0a883SMatthew G. Knepley DMPolytopeType ct; 111466a0a883SMatthew G. Knepley 1115680d18d5SMatthew G. Knepley /* TODO Check each cell individually */ 11169566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, ctx->cells[0], &ct)); 1117680d18d5SMatthew G. Knepley switch (ct) { 1118d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_SEGMENT: 1119d71ae5a4SJacob Faibussowitsch PetscCall(DMInterpolate_Segment_Private(ctx, dm, x, v)); 1120d71ae5a4SJacob Faibussowitsch break; 1121d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TRIANGLE: 1122d71ae5a4SJacob Faibussowitsch PetscCall(DMInterpolate_Triangle_Private(ctx, dm, x, v)); 1123d71ae5a4SJacob Faibussowitsch break; 1124d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_QUADRILATERAL: 1125d71ae5a4SJacob Faibussowitsch PetscCall(DMInterpolate_Quad_Private(ctx, dm, x, v)); 1126d71ae5a4SJacob Faibussowitsch break; 1127d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_TETRAHEDRON: 1128d71ae5a4SJacob Faibussowitsch PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, x, v)); 1129d71ae5a4SJacob Faibussowitsch break; 1130d71ae5a4SJacob Faibussowitsch case DM_POLYTOPE_HEXAHEDRON: 1131d71ae5a4SJacob Faibussowitsch PetscCall(DMInterpolate_Hex_Private(ctx, dm, x, v)); 1132d71ae5a4SJacob Faibussowitsch break; 1133d71ae5a4SJacob Faibussowitsch default: 1134d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]); 1135680d18d5SMatthew G. Knepley } 1136552f7358SJed Brown } 11373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1138552f7358SJed Brown } 1139552f7358SJed Brown 11404267b1a3SMatthew G. Knepley /*@C 1141f6dfbefdSBarry Smith DMInterpolationDestroy - Destroys a `DMInterpolationInfo` context 11424267b1a3SMatthew G. Knepley 1143c3339decSBarry Smith Collective 11444267b1a3SMatthew G. Knepley 11454267b1a3SMatthew G. Knepley Input Parameter: 11464267b1a3SMatthew G. Knepley . ctx - the context 11474267b1a3SMatthew G. Knepley 11484267b1a3SMatthew G. Knepley Level: beginner 11494267b1a3SMatthew G. Knepley 1150420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 11514267b1a3SMatthew G. Knepley @*/ 1152d71ae5a4SJacob Faibussowitsch PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 1153d71ae5a4SJacob Faibussowitsch { 1154552f7358SJed Brown PetscFunctionBegin; 11554f572ea9SToby Isaac PetscAssertPointer(ctx, 1); 11569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*ctx)->coords)); 11579566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->points)); 11589566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->cells)); 11599566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 11600298fd71SBarry Smith *ctx = NULL; 11613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1162552f7358SJed Brown } 1163cc0c4584SMatthew G. Knepley 1164cc0c4584SMatthew G. Knepley /*@C 1165cc0c4584SMatthew G. Knepley SNESMonitorFields - Monitors the residual for each field separately 1166cc0c4584SMatthew G. Knepley 1167c3339decSBarry Smith Collective 1168cc0c4584SMatthew G. Knepley 1169cc0c4584SMatthew G. Knepley Input Parameters: 1170420bcc1bSBarry Smith + snes - the `SNES` context, must have an attached `DM` 1171cc0c4584SMatthew G. Knepley . its - iteration number 1172cc0c4584SMatthew G. Knepley . fgnorm - 2-norm of residual 1173f6dfbefdSBarry Smith - vf - `PetscViewerAndFormat` of `PetscViewerType` `PETSCVIEWERASCII` 1174cc0c4584SMatthew G. Knepley 11752fe279fdSBarry Smith Level: intermediate 11762fe279fdSBarry Smith 1177f6dfbefdSBarry Smith Note: 1178cc0c4584SMatthew G. Knepley This routine prints the residual norm at each iteration. 1179cc0c4584SMatthew G. Knepley 1180420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESMonitorSet()`, `SNESMonitorDefault()` 1181cc0c4584SMatthew G. Knepley @*/ 1182d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, PetscViewerAndFormat *vf) 1183d71ae5a4SJacob Faibussowitsch { 1184d43b4f6eSBarry Smith PetscViewer viewer = vf->viewer; 1185cc0c4584SMatthew G. Knepley Vec res; 1186cc0c4584SMatthew G. Knepley DM dm; 1187cc0c4584SMatthew G. Knepley PetscSection s; 1188cc0c4584SMatthew G. Knepley const PetscScalar *r; 1189cc0c4584SMatthew G. Knepley PetscReal *lnorms, *norms; 1190cc0c4584SMatthew G. Knepley PetscInt numFields, f, pStart, pEnd, p; 1191cc0c4584SMatthew G. Knepley 1192cc0c4584SMatthew G. Knepley PetscFunctionBegin; 11934d4332d5SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 4); 11949566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &res, NULL, NULL)); 11959566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 11969566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, &s)); 11979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(s, &numFields)); 11989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 11999566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numFields, &lnorms, numFields, &norms)); 12009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(res, &r)); 1201cc0c4584SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 1202cc0c4584SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 1203cc0c4584SMatthew G. Knepley PetscInt fdof, foff, d; 1204cc0c4584SMatthew G. Knepley 12059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 12069566063dSJacob Faibussowitsch PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff)); 1207cc0c4584SMatthew G. Knepley for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff + d])); 1208cc0c4584SMatthew G. Knepley } 1209cc0c4584SMatthew G. Knepley } 12109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(res, &r)); 12111c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm))); 12129566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 12139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)snes)->tablevel)); 121463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " SNES Function norm %14.12e [", its, (double)fgnorm)); 1215cc0c4584SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 12169566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 12179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%14.12e", (double)PetscSqrtReal(norms[f]))); 1218cc0c4584SMatthew G. Knepley } 12199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "]\n")); 12209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)snes)->tablevel)); 12219566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 12229566063dSJacob Faibussowitsch PetscCall(PetscFree2(lnorms, norms)); 12233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1224cc0c4584SMatthew G. Knepley } 122524cdb843SMatthew G. Knepley 122624cdb843SMatthew G. Knepley /********************* Residual Computation **************************/ 122724cdb843SMatthew G. Knepley 122824cdb843SMatthew G. Knepley /*@ 1229420bcc1bSBarry Smith DMPlexSNESComputeResidualFEM - Sums the local residual into vector `F` from the local input `X` using pointwise functions specified by the user 123024cdb843SMatthew G. Knepley 123124cdb843SMatthew G. Knepley Input Parameters: 123224cdb843SMatthew G. Knepley + dm - The mesh 123324cdb843SMatthew G. Knepley . X - Local solution 123424cdb843SMatthew G. Knepley - user - The user context 123524cdb843SMatthew G. Knepley 123624cdb843SMatthew G. Knepley Output Parameter: 123724cdb843SMatthew G. Knepley . F - Local output vector 123824cdb843SMatthew G. Knepley 12392fe279fdSBarry Smith Level: developer 12402fe279fdSBarry Smith 1241f6dfbefdSBarry Smith Note: 1242420bcc1bSBarry Smith The residual is summed into `F`; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in `F` is intentional. 12438db23a53SJed Brown 1244420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()` 124524cdb843SMatthew G. Knepley @*/ 1246d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 1247d71ae5a4SJacob Faibussowitsch { 12486da023fcSToby Isaac DM plex; 1249083401c6SMatthew G. Knepley IS allcellIS; 12506528b96dSMatthew G. Knepley PetscInt Nds, s; 125124cdb843SMatthew G. Knepley 125224cdb843SMatthew G. Knepley PetscFunctionBegin; 12539566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 12549566063dSJacob Faibussowitsch PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 12559566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 12566528b96dSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 12576528b96dSMatthew G. Knepley PetscDS ds; 12586528b96dSMatthew G. Knepley IS cellIS; 125906ad1575SMatthew G. Knepley PetscFormKey key; 12606528b96dSMatthew G. Knepley 126107218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL)); 12626528b96dSMatthew G. Knepley key.value = 0; 12636528b96dSMatthew G. Knepley key.field = 0; 126406ad1575SMatthew G. Knepley key.part = 0; 12656528b96dSMatthew G. Knepley if (!key.label) { 12669566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)allcellIS)); 12676528b96dSMatthew G. Knepley cellIS = allcellIS; 12686528b96dSMatthew G. Knepley } else { 12696528b96dSMatthew G. Knepley IS pointIS; 12706528b96dSMatthew G. Knepley 12716528b96dSMatthew G. Knepley key.value = 1; 12729566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS)); 12739566063dSJacob Faibussowitsch PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 12749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 12756528b96dSMatthew G. Knepley } 12769566063dSJacob Faibussowitsch PetscCall(DMPlexComputeResidual_Internal(plex, key, cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user)); 12779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 12786528b96dSMatthew G. Knepley } 12799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&allcellIS)); 12809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 12813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12826528b96dSMatthew G. Knepley } 12836528b96dSMatthew G. Knepley 1284bdd6f66aSToby Isaac /*@ 1285420bcc1bSBarry Smith DMPlexSNESComputeResidualDS - Sums the local residual into vector `F` from the local input `X` using all pointwise functions with unique keys in the `PetscDS` 1286d2b2dc1eSMatthew G. Knepley 1287d2b2dc1eSMatthew G. Knepley Input Parameters: 1288d2b2dc1eSMatthew G. Knepley + dm - The mesh 1289d2b2dc1eSMatthew G. Knepley . X - Local solution 1290d2b2dc1eSMatthew G. Knepley - user - The user context 1291d2b2dc1eSMatthew G. Knepley 1292d2b2dc1eSMatthew G. Knepley Output Parameter: 1293d2b2dc1eSMatthew G. Knepley . F - Local output vector 1294d2b2dc1eSMatthew G. Knepley 1295d2b2dc1eSMatthew G. Knepley Level: developer 1296d2b2dc1eSMatthew G. Knepley 1297d2b2dc1eSMatthew G. Knepley Note: 1298420bcc1bSBarry Smith The residual is summed into `F`; the caller is responsible for using `VecZeroEntries()` or otherwise ensuring that any data in `F` is intentional. 1299d2b2dc1eSMatthew G. Knepley 1300420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()` 1301d2b2dc1eSMatthew G. Knepley @*/ 1302d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualDS(DM dm, Vec X, Vec F, void *user) 1303d2b2dc1eSMatthew G. Knepley { 1304d2b2dc1eSMatthew G. Knepley DM plex; 1305d2b2dc1eSMatthew G. Knepley IS allcellIS; 1306d2b2dc1eSMatthew G. Knepley PetscInt Nds, s; 1307d2b2dc1eSMatthew G. Knepley 1308d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 1309d2b2dc1eSMatthew G. Knepley PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 1310d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 1311d2b2dc1eSMatthew G. Knepley PetscCall(DMGetNumDS(dm, &Nds)); 1312d2b2dc1eSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 1313d2b2dc1eSMatthew G. Knepley PetscDS ds; 1314d2b2dc1eSMatthew G. Knepley DMLabel label; 1315d2b2dc1eSMatthew G. Knepley IS cellIS; 1316d2b2dc1eSMatthew G. Knepley 1317d2b2dc1eSMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL)); 1318d2b2dc1eSMatthew G. Knepley { 1319d2b2dc1eSMatthew G. Knepley PetscWeakFormKind resmap[2] = {PETSC_WF_F0, PETSC_WF_F1}; 1320d2b2dc1eSMatthew G. Knepley PetscWeakForm wf; 1321d2b2dc1eSMatthew G. Knepley PetscInt Nm = 2, m, Nk = 0, k, kp, off = 0; 1322d2b2dc1eSMatthew G. Knepley PetscFormKey *reskeys; 1323d2b2dc1eSMatthew G. Knepley 1324d2b2dc1eSMatthew G. Knepley /* Get unique residual keys */ 1325d2b2dc1eSMatthew G. Knepley for (m = 0; m < Nm; ++m) { 1326d2b2dc1eSMatthew G. Knepley PetscInt Nkm; 1327d2b2dc1eSMatthew G. Knepley PetscCall(PetscHMapFormGetSize(ds->wf->form[resmap[m]], &Nkm)); 1328d2b2dc1eSMatthew G. Knepley Nk += Nkm; 1329d2b2dc1eSMatthew G. Knepley } 1330d2b2dc1eSMatthew G. Knepley PetscCall(PetscMalloc1(Nk, &reskeys)); 1331d2b2dc1eSMatthew G. Knepley for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[resmap[m]], &off, reskeys)); 1332d2b2dc1eSMatthew G. Knepley PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk); 1333d2b2dc1eSMatthew G. Knepley PetscCall(PetscFormKeySort(Nk, reskeys)); 1334d2b2dc1eSMatthew G. Knepley for (k = 0, kp = 1; kp < Nk; ++kp) { 1335d2b2dc1eSMatthew G. Knepley if ((reskeys[k].label != reskeys[kp].label) || (reskeys[k].value != reskeys[kp].value)) { 1336d2b2dc1eSMatthew G. Knepley ++k; 1337d2b2dc1eSMatthew G. Knepley if (kp != k) reskeys[k] = reskeys[kp]; 1338d2b2dc1eSMatthew G. Knepley } 1339d2b2dc1eSMatthew G. Knepley } 1340d2b2dc1eSMatthew G. Knepley Nk = k; 1341d2b2dc1eSMatthew G. Knepley 1342d2b2dc1eSMatthew G. Knepley PetscCall(PetscDSGetWeakForm(ds, &wf)); 1343d2b2dc1eSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 1344d2b2dc1eSMatthew G. Knepley DMLabel label = reskeys[k].label; 1345d2b2dc1eSMatthew G. Knepley PetscInt val = reskeys[k].value; 1346d2b2dc1eSMatthew G. Knepley 1347d2b2dc1eSMatthew G. Knepley if (!label) { 1348d2b2dc1eSMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)allcellIS)); 1349d2b2dc1eSMatthew G. Knepley cellIS = allcellIS; 1350d2b2dc1eSMatthew G. Knepley } else { 1351d2b2dc1eSMatthew G. Knepley IS pointIS; 1352d2b2dc1eSMatthew G. Knepley 1353d2b2dc1eSMatthew G. Knepley PetscCall(DMLabelGetStratumIS(label, val, &pointIS)); 1354d2b2dc1eSMatthew G. Knepley PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 1355d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&pointIS)); 1356d2b2dc1eSMatthew G. Knepley } 1357d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexComputeResidual_Internal(plex, reskeys[k], cellIS, PETSC_MIN_REAL, X, NULL, 0.0, F, user)); 1358d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&cellIS)); 1359d2b2dc1eSMatthew G. Knepley } 1360d2b2dc1eSMatthew G. Knepley PetscCall(PetscFree(reskeys)); 1361d2b2dc1eSMatthew G. Knepley } 1362d2b2dc1eSMatthew G. Knepley } 1363d2b2dc1eSMatthew G. Knepley PetscCall(ISDestroy(&allcellIS)); 1364d2b2dc1eSMatthew G. Knepley PetscCall(DMDestroy(&plex)); 1365d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 1366d2b2dc1eSMatthew G. Knepley } 1367d2b2dc1eSMatthew G. Knepley 1368d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 1369d2b2dc1eSMatthew G. Knepley PetscErrorCode DMPlexSNESComputeResidualCEED(DM dm, Vec locX, Vec locF, void *user) 1370d2b2dc1eSMatthew G. Knepley { 1371d2b2dc1eSMatthew G. Knepley Ceed ceed; 1372d2b2dc1eSMatthew G. Knepley DMCeed sd = dm->dmceed; 1373d2b2dc1eSMatthew G. Knepley CeedVector clocX, clocF; 1374d2b2dc1eSMatthew G. Knepley 1375d2b2dc1eSMatthew G. Knepley PetscFunctionBegin; 1376d2b2dc1eSMatthew G. Knepley PetscCall(DMGetCeed(dm, &ceed)); 1377d2b2dc1eSMatthew G. Knepley PetscCheck(sd, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "This DM has no CEED data. Call DMCeedCreate() before computing the residual."); 1378d2b2dc1eSMatthew G. Knepley PetscCall(DMCeedComputeGeometry(dm, sd)); 1379d2b2dc1eSMatthew G. Knepley 1380d2b2dc1eSMatthew G. Knepley PetscCall(VecGetCeedVectorRead(locX, ceed, &clocX)); 1381d2b2dc1eSMatthew G. Knepley PetscCall(VecGetCeedVector(locF, ceed, &clocF)); 1382d2b2dc1eSMatthew G. Knepley PetscCallCEED(CeedOperatorApplyAdd(sd->op, clocX, clocF, CEED_REQUEST_IMMEDIATE)); 1383d2b2dc1eSMatthew G. Knepley PetscCall(VecRestoreCeedVectorRead(locX, &clocX)); 1384d2b2dc1eSMatthew G. Knepley PetscCall(VecRestoreCeedVector(locF, &clocF)); 1385d2b2dc1eSMatthew G. Knepley 1386d2b2dc1eSMatthew G. Knepley { 1387d2b2dc1eSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *)dm->data; 1388d2b2dc1eSMatthew G. Knepley 1389d2b2dc1eSMatthew G. Knepley if (mesh->printFEM) { 1390d2b2dc1eSMatthew G. Knepley PetscSection section; 1391d2b2dc1eSMatthew G. Knepley Vec locFbc; 1392d2b2dc1eSMatthew G. Knepley PetscInt pStart, pEnd, p, maxDof; 1393d2b2dc1eSMatthew G. Knepley PetscScalar *zeroes; 1394d2b2dc1eSMatthew G. Knepley 1395d2b2dc1eSMatthew G. Knepley PetscCall(DMGetLocalSection(dm, §ion)); 1396d2b2dc1eSMatthew G. Knepley PetscCall(VecDuplicate(locF, &locFbc)); 1397d2b2dc1eSMatthew G. Knepley PetscCall(VecCopy(locF, locFbc)); 1398d2b2dc1eSMatthew G. Knepley PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 1399d2b2dc1eSMatthew G. Knepley PetscCall(PetscSectionGetMaxDof(section, &maxDof)); 1400d2b2dc1eSMatthew G. Knepley PetscCall(PetscCalloc1(maxDof, &zeroes)); 1401d2b2dc1eSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) PetscCall(VecSetValuesSection(locFbc, section, p, zeroes, INSERT_BC_VALUES)); 1402d2b2dc1eSMatthew G. Knepley PetscCall(PetscFree(zeroes)); 1403d2b2dc1eSMatthew G. Knepley PetscCall(DMPrintLocalVec(dm, "Residual", mesh->printTol, locFbc)); 1404d2b2dc1eSMatthew G. Knepley PetscCall(VecDestroy(&locFbc)); 1405d2b2dc1eSMatthew G. Knepley } 1406d2b2dc1eSMatthew G. Knepley } 1407d2b2dc1eSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 1408d2b2dc1eSMatthew G. Knepley } 1409d2b2dc1eSMatthew G. Knepley #endif 1410d2b2dc1eSMatthew G. Knepley 1411d2b2dc1eSMatthew G. Knepley /*@ 1412420bcc1bSBarry Smith DMPlexSNESComputeBoundaryFEM - Form the boundary values for the local input `X` 1413bdd6f66aSToby Isaac 1414bdd6f66aSToby Isaac Input Parameters: 1415bdd6f66aSToby Isaac + dm - The mesh 1416bdd6f66aSToby Isaac - user - The user context 1417bdd6f66aSToby Isaac 1418bdd6f66aSToby Isaac Output Parameter: 1419bdd6f66aSToby Isaac . X - Local solution 1420bdd6f66aSToby Isaac 1421bdd6f66aSToby Isaac Level: developer 1422bdd6f66aSToby Isaac 1423420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `DMPLEX`, `DMPlexComputeJacobianAction()` 1424bdd6f66aSToby Isaac @*/ 1425d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeBoundaryFEM(DM dm, Vec X, void *user) 1426d71ae5a4SJacob Faibussowitsch { 1427bdd6f66aSToby Isaac DM plex; 1428bdd6f66aSToby Isaac 1429bdd6f66aSToby Isaac PetscFunctionBegin; 14309566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 14319566063dSJacob Faibussowitsch PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, X, PETSC_MIN_REAL, NULL, NULL, NULL)); 14329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 14333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1434bdd6f66aSToby Isaac } 1435bdd6f66aSToby Isaac 14367a73cf09SMatthew G. Knepley /*@ 1437*c118d280SBarry Smith DMSNESComputeJacobianAction - Compute the action of the Jacobian J(`X`) on `Y` 14387a73cf09SMatthew G. Knepley 14397a73cf09SMatthew G. Knepley Input Parameters: 1440f6dfbefdSBarry Smith + dm - The `DM` 14417a73cf09SMatthew G. Knepley . X - Local solution vector 14427a73cf09SMatthew G. Knepley . Y - Local input vector 14437a73cf09SMatthew G. Knepley - user - The user context 14447a73cf09SMatthew G. Knepley 14457a73cf09SMatthew G. Knepley Output Parameter: 144669d47153SPierre Jolivet . F - local output vector 14477a73cf09SMatthew G. Knepley 14487a73cf09SMatthew G. Knepley Level: developer 14497a73cf09SMatthew G. Knepley 1450*c118d280SBarry Smith Note: 1451f6dfbefdSBarry Smith Users will typically use `DMSNESCreateJacobianMF()` followed by `MatMult()` instead of calling this routine directly. 14528e3a2eefSMatthew G. Knepley 1453*c118d280SBarry Smith This only works with `DMPLEX` 1454*c118d280SBarry Smith 1455*c118d280SBarry Smith Developer Note: 1456*c118d280SBarry Smith This should be called `DMPlexSNESComputeJacobianAction()` 1457*c118d280SBarry Smith 1458420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, ``DMSNESCreateJacobianMF()`, `DMPlexSNESComputeResidualFEM()` 14597a73cf09SMatthew G. Knepley @*/ 1460d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESComputeJacobianAction(DM dm, Vec X, Vec Y, Vec F, void *user) 1461d71ae5a4SJacob Faibussowitsch { 14628e3a2eefSMatthew G. Knepley DM plex; 14638e3a2eefSMatthew G. Knepley IS allcellIS; 14648e3a2eefSMatthew G. Knepley PetscInt Nds, s; 1465a925c78cSMatthew G. Knepley 1466a925c78cSMatthew G. Knepley PetscFunctionBegin; 14679566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 14689566063dSJacob Faibussowitsch PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 14699566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 14708e3a2eefSMatthew G. Knepley for (s = 0; s < Nds; ++s) { 14718e3a2eefSMatthew G. Knepley PetscDS ds; 14728e3a2eefSMatthew G. Knepley DMLabel label; 14738e3a2eefSMatthew G. Knepley IS cellIS; 14747a73cf09SMatthew G. Knepley 147507218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &label, NULL, &ds, NULL)); 14768e3a2eefSMatthew G. Knepley { 147706ad1575SMatthew G. Knepley PetscWeakFormKind jacmap[4] = {PETSC_WF_G0, PETSC_WF_G1, PETSC_WF_G2, PETSC_WF_G3}; 14788e3a2eefSMatthew G. Knepley PetscWeakForm wf; 14798e3a2eefSMatthew G. Knepley PetscInt Nm = 4, m, Nk = 0, k, kp, off = 0; 148006ad1575SMatthew G. Knepley PetscFormKey *jackeys; 14818e3a2eefSMatthew G. Knepley 14828e3a2eefSMatthew G. Knepley /* Get unique Jacobian keys */ 14838e3a2eefSMatthew G. Knepley for (m = 0; m < Nm; ++m) { 14848e3a2eefSMatthew G. Knepley PetscInt Nkm; 14859566063dSJacob Faibussowitsch PetscCall(PetscHMapFormGetSize(ds->wf->form[jacmap[m]], &Nkm)); 14868e3a2eefSMatthew G. Knepley Nk += Nkm; 14878e3a2eefSMatthew G. Knepley } 14889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk, &jackeys)); 148948a46eb9SPierre Jolivet for (m = 0; m < Nm; ++m) PetscCall(PetscHMapFormGetKeys(ds->wf->form[jacmap[m]], &off, jackeys)); 149063a3b9bcSJacob Faibussowitsch PetscCheck(off == Nk, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of keys %" PetscInt_FMT " should be %" PetscInt_FMT, off, Nk); 14919566063dSJacob Faibussowitsch PetscCall(PetscFormKeySort(Nk, jackeys)); 14928e3a2eefSMatthew G. Knepley for (k = 0, kp = 1; kp < Nk; ++kp) { 14938e3a2eefSMatthew G. Knepley if ((jackeys[k].label != jackeys[kp].label) || (jackeys[k].value != jackeys[kp].value)) { 14948e3a2eefSMatthew G. Knepley ++k; 14958e3a2eefSMatthew G. Knepley if (kp != k) jackeys[k] = jackeys[kp]; 14968e3a2eefSMatthew G. Knepley } 14978e3a2eefSMatthew G. Knepley } 14988e3a2eefSMatthew G. Knepley Nk = k; 14998e3a2eefSMatthew G. Knepley 15009566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 15018e3a2eefSMatthew G. Knepley for (k = 0; k < Nk; ++k) { 15028e3a2eefSMatthew G. Knepley DMLabel label = jackeys[k].label; 15038e3a2eefSMatthew G. Knepley PetscInt val = jackeys[k].value; 15048e3a2eefSMatthew G. Knepley 15058e3a2eefSMatthew G. Knepley if (!label) { 15069566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)allcellIS)); 15078e3a2eefSMatthew G. Knepley cellIS = allcellIS; 15087a73cf09SMatthew G. Knepley } else { 15098e3a2eefSMatthew G. Knepley IS pointIS; 1510a925c78cSMatthew G. Knepley 15119566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, val, &pointIS)); 15129566063dSJacob Faibussowitsch PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 15139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 1514a925c78cSMatthew G. Knepley } 15159566063dSJacob Faibussowitsch PetscCall(DMPlexComputeJacobian_Action_Internal(plex, jackeys[k], cellIS, 0.0, 0.0, X, NULL, Y, F, user)); 15169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 15178e3a2eefSMatthew G. Knepley } 15189566063dSJacob Faibussowitsch PetscCall(PetscFree(jackeys)); 15198e3a2eefSMatthew G. Knepley } 15208e3a2eefSMatthew G. Knepley } 15219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&allcellIS)); 15229566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 15233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1524a925c78cSMatthew G. Knepley } 1525a925c78cSMatthew G. Knepley 152624cdb843SMatthew G. Knepley /*@ 15272fe279fdSBarry Smith DMPlexSNESComputeJacobianFEM - Form the local portion of the Jacobian matrix `Jac` at the local solution `X` using pointwise functions specified by the user. 152824cdb843SMatthew G. Knepley 152924cdb843SMatthew G. Knepley Input Parameters: 15302fe279fdSBarry Smith + dm - The `DM` 153124cdb843SMatthew G. Knepley . X - Local input vector 153224cdb843SMatthew G. Knepley - user - The user context 153324cdb843SMatthew G. Knepley 15342fe279fdSBarry Smith Output Parameters: 15352fe279fdSBarry Smith + Jac - Jacobian matrix 15362fe279fdSBarry Smith - JacP - approximate Jacobian from which the preconditioner will be built, often `Jac` 15372fe279fdSBarry Smith 15382fe279fdSBarry Smith Level: developer 153924cdb843SMatthew G. Knepley 154024cdb843SMatthew G. Knepley Note: 154124cdb843SMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 154224cdb843SMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 154324cdb843SMatthew G. Knepley 1544420bcc1bSBarry Smith .seealso: [](ch_snes), `DMPLEX`, `Mat` 154524cdb843SMatthew G. Knepley @*/ 1546d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSNESComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, void *user) 1547d71ae5a4SJacob Faibussowitsch { 15486da023fcSToby Isaac DM plex; 1549083401c6SMatthew G. Knepley IS allcellIS; 1550f04eb4edSMatthew G. Knepley PetscBool hasJac, hasPrec; 15516528b96dSMatthew G. Knepley PetscInt Nds, s; 155224cdb843SMatthew G. Knepley 155324cdb843SMatthew G. Knepley PetscFunctionBegin; 15549566063dSJacob Faibussowitsch PetscCall(DMSNESConvertPlex(dm, &plex, PETSC_TRUE)); 15559566063dSJacob Faibussowitsch PetscCall(DMPlexGetAllCells_Internal(plex, &allcellIS)); 15569566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 1557083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 1558083401c6SMatthew G. Knepley PetscDS ds; 1559083401c6SMatthew G. Knepley IS cellIS; 156006ad1575SMatthew G. Knepley PetscFormKey key; 1561083401c6SMatthew G. Knepley 156207218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &key.label, NULL, &ds, NULL)); 15636528b96dSMatthew G. Knepley key.value = 0; 15646528b96dSMatthew G. Knepley key.field = 0; 156506ad1575SMatthew G. Knepley key.part = 0; 15666528b96dSMatthew G. Knepley if (!key.label) { 15679566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)allcellIS)); 1568083401c6SMatthew G. Knepley cellIS = allcellIS; 1569083401c6SMatthew G. Knepley } else { 1570083401c6SMatthew G. Knepley IS pointIS; 1571083401c6SMatthew G. Knepley 15726528b96dSMatthew G. Knepley key.value = 1; 15739566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(key.label, key.value, &pointIS)); 15749566063dSJacob Faibussowitsch PetscCall(ISIntersect_Caching_Internal(allcellIS, pointIS, &cellIS)); 15759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&pointIS)); 1576083401c6SMatthew G. Knepley } 1577083401c6SMatthew G. Knepley if (!s) { 15789566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobian(ds, &hasJac)); 15799566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec)); 15809566063dSJacob Faibussowitsch if (hasJac && hasPrec) PetscCall(MatZeroEntries(Jac)); 15819566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(JacP)); 1582083401c6SMatthew G. Knepley } 15839566063dSJacob Faibussowitsch PetscCall(DMPlexComputeJacobian_Internal(plex, key, cellIS, 0.0, 0.0, X, NULL, Jac, JacP, user)); 15849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 1585083401c6SMatthew G. Knepley } 15869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&allcellIS)); 15879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&plex)); 15883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158924cdb843SMatthew G. Knepley } 15901878804aSMatthew G. Knepley 15919371c9d4SSatish Balay struct _DMSNESJacobianMFCtx { 15928e3a2eefSMatthew G. Knepley DM dm; 15938e3a2eefSMatthew G. Knepley Vec X; 15948e3a2eefSMatthew G. Knepley void *ctx; 15958e3a2eefSMatthew G. Knepley }; 15968e3a2eefSMatthew G. Knepley 1597d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Destroy_Private(Mat A) 1598d71ae5a4SJacob Faibussowitsch { 15998e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx *ctx; 16008e3a2eefSMatthew G. Knepley 16018e3a2eefSMatthew G. Knepley PetscFunctionBegin; 16029566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 16039566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, NULL)); 16049566063dSJacob Faibussowitsch PetscCall(DMDestroy(&ctx->dm)); 16059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->X)); 16069566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx)); 16073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16088e3a2eefSMatthew G. Knepley } 16098e3a2eefSMatthew G. Knepley 1610d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSNESJacobianMF_Mult_Private(Mat A, Vec Y, Vec Z) 1611d71ae5a4SJacob Faibussowitsch { 16128e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx *ctx; 16138e3a2eefSMatthew G. Knepley 16148e3a2eefSMatthew G. Knepley PetscFunctionBegin; 16159566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &ctx)); 16169566063dSJacob Faibussowitsch PetscCall(DMSNESComputeJacobianAction(ctx->dm, ctx->X, Y, Z, ctx->ctx)); 16173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16188e3a2eefSMatthew G. Knepley } 16198e3a2eefSMatthew G. Knepley 16208e3a2eefSMatthew G. Knepley /*@ 1621f6dfbefdSBarry Smith DMSNESCreateJacobianMF - Create a `Mat` which computes the action of the Jacobian matrix-free 16228e3a2eefSMatthew G. Knepley 1623c3339decSBarry Smith Collective 16248e3a2eefSMatthew G. Knepley 16258e3a2eefSMatthew G. Knepley Input Parameters: 1626f6dfbefdSBarry Smith + dm - The `DM` 16278e3a2eefSMatthew G. Knepley . X - The evaluation point for the Jacobian 16282fe279fdSBarry Smith - user - A user context, or `NULL` 16298e3a2eefSMatthew G. Knepley 16308e3a2eefSMatthew G. Knepley Output Parameter: 1631f6dfbefdSBarry Smith . J - The `Mat` 16328e3a2eefSMatthew G. Knepley 16338e3a2eefSMatthew G. Knepley Level: advanced 16348e3a2eefSMatthew G. Knepley 1635*c118d280SBarry Smith Notes: 16362fe279fdSBarry Smith Vec `X` is kept in `J`, so updating `X` then updates the evaluation point. 16378e3a2eefSMatthew G. Knepley 1638*c118d280SBarry Smith This only works for `DMPLEX` 1639*c118d280SBarry Smith 1640420bcc1bSBarry Smith .seealso: [](ch_snes), `DM`, `SNES`, `DMSNESComputeJacobianAction()` 16418e3a2eefSMatthew G. Knepley @*/ 1642d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCreateJacobianMF(DM dm, Vec X, void *user, Mat *J) 1643d71ae5a4SJacob Faibussowitsch { 16448e3a2eefSMatthew G. Knepley struct _DMSNESJacobianMFCtx *ctx; 16458e3a2eefSMatthew G. Knepley PetscInt n, N; 16468e3a2eefSMatthew G. Knepley 16478e3a2eefSMatthew G. Knepley PetscFunctionBegin; 16489566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 16499566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATSHELL)); 16509566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(X, &n)); 16519566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N)); 16529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, n, n, N, N)); 16539566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 16549566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)X)); 16559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, &ctx)); 16568e3a2eefSMatthew G. Knepley ctx->dm = dm; 16578e3a2eefSMatthew G. Knepley ctx->X = X; 16588e3a2eefSMatthew G. Knepley ctx->ctx = user; 16599566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(*J, ctx)); 16609566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_DESTROY, (void (*)(void))DMSNESJacobianMF_Destroy_Private)); 16619566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(*J, MATOP_MULT, (void (*)(void))DMSNESJacobianMF_Mult_Private)); 16623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16638e3a2eefSMatthew G. Knepley } 16648e3a2eefSMatthew G. Knepley 166538cfc46eSPierre Jolivet /* 166638cfc46eSPierre Jolivet MatComputeNeumannOverlap - Computes an unassembled (Neumann) local overlapping Mat in nonlinear context. 166738cfc46eSPierre Jolivet 166838cfc46eSPierre Jolivet Input Parameters: 16692fe279fdSBarry Smith + X - `SNES` linearization point 167038cfc46eSPierre Jolivet . ovl - index set of overlapping subdomains 167138cfc46eSPierre Jolivet 167238cfc46eSPierre Jolivet Output Parameter: 167338cfc46eSPierre Jolivet . J - unassembled (Neumann) local matrix 167438cfc46eSPierre Jolivet 167538cfc46eSPierre Jolivet Level: intermediate 167638cfc46eSPierre Jolivet 1677420bcc1bSBarry Smith .seealso: [](ch_snes), `DMCreateNeumannOverlap()`, `MATIS`, `PCHPDDMSetAuxiliaryMat()` 167838cfc46eSPierre Jolivet */ 1679d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatComputeNeumannOverlap_Plex(Mat J, PetscReal t, Vec X, Vec X_t, PetscReal s, IS ovl, void *ctx) 1680d71ae5a4SJacob Faibussowitsch { 168138cfc46eSPierre Jolivet SNES snes; 168238cfc46eSPierre Jolivet Mat pJ; 168338cfc46eSPierre Jolivet DM ovldm, origdm; 168438cfc46eSPierre Jolivet DMSNES sdm; 168538cfc46eSPierre Jolivet PetscErrorCode (*bfun)(DM, Vec, void *); 168638cfc46eSPierre Jolivet PetscErrorCode (*jfun)(DM, Vec, Mat, Mat, void *); 168738cfc46eSPierre Jolivet void *bctx, *jctx; 168838cfc46eSPierre Jolivet 168938cfc46eSPierre Jolivet PetscFunctionBegin; 16909566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_MATIS", (PetscObject *)&pJ)); 169128b400f6SJacob Faibussowitsch PetscCheck(pJ, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing overlapping Mat"); 16929566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Original_HPDDM", (PetscObject *)&origdm)); 169328b400f6SJacob Faibussowitsch PetscCheck(origdm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing original DM"); 16949566063dSJacob Faibussowitsch PetscCall(MatGetDM(pJ, &ovldm)); 16959566063dSJacob Faibussowitsch PetscCall(DMSNESGetBoundaryLocal(origdm, &bfun, &bctx)); 16969566063dSJacob Faibussowitsch PetscCall(DMSNESSetBoundaryLocal(ovldm, bfun, bctx)); 16979566063dSJacob Faibussowitsch PetscCall(DMSNESGetJacobianLocal(origdm, &jfun, &jctx)); 16989566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobianLocal(ovldm, jfun, jctx)); 16999566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject *)&snes)); 170038cfc46eSPierre Jolivet if (!snes) { 17019566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)ovl), &snes)); 17029566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, ovldm)); 17039566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)ovl, "_DM_Overlap_HPDDM_SNES", (PetscObject)snes)); 17049566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)snes)); 170538cfc46eSPierre Jolivet } 17069566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(ovldm, &sdm)); 17079566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(X)); 1708800f99ffSJeremy L Thompson { 1709800f99ffSJeremy L Thompson void *ctx; 1710800f99ffSJeremy L Thompson PetscErrorCode (*J)(SNES, Vec, Mat, Mat, void *); 1711800f99ffSJeremy L Thompson PetscCall(DMSNESGetJacobian(ovldm, &J, &ctx)); 1712800f99ffSJeremy L Thompson PetscCallBack("SNES callback Jacobian", (*J)(snes, X, pJ, pJ, ctx)); 1713800f99ffSJeremy L Thompson } 17149566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(X)); 171538cfc46eSPierre Jolivet /* this is a no-hop, just in case we decide to change the placeholder for the local Neumann matrix */ 171638cfc46eSPierre Jolivet { 171738cfc46eSPierre Jolivet Mat locpJ; 171838cfc46eSPierre Jolivet 17199566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(pJ, &locpJ)); 17209566063dSJacob Faibussowitsch PetscCall(MatCopy(locpJ, J, SAME_NONZERO_PATTERN)); 172138cfc46eSPierre Jolivet } 17223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172338cfc46eSPierre Jolivet } 172438cfc46eSPierre Jolivet 1725a925c78cSMatthew G. Knepley /*@ 1726f6dfbefdSBarry Smith DMPlexSetSNESLocalFEM - Use `DMPLEX`'s internal FEM routines to compute `SNES` boundary values, residual, and Jacobian. 17279f520fc2SToby Isaac 17289f520fc2SToby Isaac Input Parameters: 1729f6dfbefdSBarry Smith + dm - The `DM` object 1730f6dfbefdSBarry Smith . boundaryctx - the user context that will be passed to pointwise evaluation of boundary values (see `PetscDSAddBoundary()`) 1731f6dfbefdSBarry Smith . residualctx - the user context that will be passed to pointwise evaluation of finite element residual computations (see `PetscDSSetResidual()`) 1732f6dfbefdSBarry Smith - jacobianctx - the user context that will be passed to pointwise evaluation of finite element Jacobian construction (see `PetscDSSetJacobian()`) 17331a244344SSatish Balay 17341a244344SSatish Balay Level: developer 1735f6dfbefdSBarry Smith 1736420bcc1bSBarry Smith .seealso: [](ch_snes), `DMPLEX`, `SNES` 17379f520fc2SToby Isaac @*/ 1738d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetSNESLocalFEM(DM dm, void *boundaryctx, void *residualctx, void *jacobianctx) 1739d71ae5a4SJacob Faibussowitsch { 1740d2b2dc1eSMatthew G. Knepley PetscBool useCeed; 1741d2b2dc1eSMatthew G. Knepley 17429f520fc2SToby Isaac PetscFunctionBegin; 1743d2b2dc1eSMatthew G. Knepley PetscCall(DMPlexGetUseCeed(dm, &useCeed)); 17449566063dSJacob Faibussowitsch PetscCall(DMSNESSetBoundaryLocal(dm, DMPlexSNESComputeBoundaryFEM, boundaryctx)); 1745d2b2dc1eSMatthew G. Knepley if (useCeed) { 1746d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 1747d2b2dc1eSMatthew G. Knepley PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualCEED, residualctx)); 1748d2b2dc1eSMatthew G. Knepley #else 1749d2b2dc1eSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Cannot use CEED traversals without LibCEED. Rerun configure with --download-ceed"); 1750d2b2dc1eSMatthew G. Knepley #endif 1751d2b2dc1eSMatthew G. Knepley } else PetscCall(DMSNESSetFunctionLocal(dm, DMPlexSNESComputeResidualFEM, residualctx)); 17529566063dSJacob Faibussowitsch PetscCall(DMSNESSetJacobianLocal(dm, DMPlexSNESComputeJacobianFEM, jacobianctx)); 17539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", MatComputeNeumannOverlap_Plex)); 17543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17559f520fc2SToby Isaac } 17569f520fc2SToby Isaac 1757f017bcb6SMatthew G. Knepley /*@C 1758f017bcb6SMatthew G. Knepley DMSNESCheckDiscretization - Check the discretization error of the exact solution 1759f017bcb6SMatthew G. Knepley 1760f017bcb6SMatthew G. Knepley Input Parameters: 1761f6dfbefdSBarry Smith + snes - the `SNES` object 1762f6dfbefdSBarry Smith . dm - the `DM` 1763f2cacb80SMatthew G. Knepley . t - the time 1764f6dfbefdSBarry Smith . u - a `DM` vector 1765f0fc11ceSJed Brown - tol - A tolerance for the check, or -1 to print the results instead 1766f017bcb6SMatthew G. Knepley 17672fe279fdSBarry Smith Output Parameter: 17682fe279fdSBarry Smith . error - An array which holds the discretization error in each field, or `NULL` 17692fe279fdSBarry Smith 17702fe279fdSBarry Smith Level: developer 1771f3c94560SMatthew G. Knepley 1772f6dfbefdSBarry Smith Note: 1773f6dfbefdSBarry Smith The user must call `PetscDSSetExactSolution()` beforehand 17747f96f943SMatthew G. Knepley 1775420bcc1bSBarry Smith Developer Note: 1776420bcc1bSBarry Smith How is this related to `PetscConvEst`? 1777420bcc1bSBarry Smith 1778420bcc1bSBarry Smith .seealso: [](ch_snes), `PetscDSSetExactSolution()`, `DNSNESCheckFromOptions()`, `DMSNESCheckResidual()`, `DMSNESCheckJacobian()` 1779f017bcb6SMatthew G. Knepley @*/ 1780d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckDiscretization(SNES snes, DM dm, PetscReal t, Vec u, PetscReal tol, PetscReal error[]) 1781d71ae5a4SJacob Faibussowitsch { 1782f017bcb6SMatthew G. Knepley PetscErrorCode (**exacts)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 1783f017bcb6SMatthew G. Knepley void **ectxs; 1784f3c94560SMatthew G. Knepley PetscReal *err; 17857f96f943SMatthew G. Knepley MPI_Comm comm; 17867f96f943SMatthew G. Knepley PetscInt Nf, f; 17871878804aSMatthew G. Knepley 17881878804aSMatthew G. Knepley PetscFunctionBegin; 1789f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1790f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1791064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 17924f572ea9SToby Isaac if (error) PetscAssertPointer(error, 6); 17937f96f943SMatthew G. Knepley 17949566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, t, u, NULL)); 17959566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-vec_view")); 17967f96f943SMatthew G. Knepley 17979566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 17989566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(dm, &Nf)); 17999566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(Nf, &exacts, Nf, &ectxs, PetscMax(1, Nf), &err)); 18007f96f943SMatthew G. Knepley { 18017f96f943SMatthew G. Knepley PetscInt Nds, s; 18027f96f943SMatthew G. Knepley 18039566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(dm, &Nds)); 1804083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 18057f96f943SMatthew G. Knepley PetscDS ds; 1806083401c6SMatthew G. Knepley DMLabel label; 1807083401c6SMatthew G. Knepley IS fieldIS; 18087f96f943SMatthew G. Knepley const PetscInt *fields; 18097f96f943SMatthew G. Knepley PetscInt dsNf, f; 1810083401c6SMatthew G. Knepley 181107218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(dm, s, &label, &fieldIS, &ds, NULL)); 18129566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &dsNf)); 18139566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fieldIS, &fields)); 1814083401c6SMatthew G. Knepley for (f = 0; f < dsNf; ++f) { 1815083401c6SMatthew G. Knepley const PetscInt field = fields[f]; 18169566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, field, &exacts[field], &ectxs[field])); 1817083401c6SMatthew G. Knepley } 18189566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fieldIS, &fields)); 1819083401c6SMatthew G. Knepley } 1820083401c6SMatthew G. Knepley } 1821f017bcb6SMatthew G. Knepley if (Nf > 1) { 18229566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, t, exacts, ectxs, u, err)); 1823f017bcb6SMatthew G. Knepley if (tol >= 0.0) { 1824ad540459SPierre 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); 1825b878532fSMatthew G. Knepley } else if (error) { 1826b878532fSMatthew G. Knepley for (f = 0; f < Nf; ++f) error[f] = err[f]; 18271878804aSMatthew G. Knepley } else { 18289566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Error: [")); 1829f017bcb6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 18309566063dSJacob Faibussowitsch if (f) PetscCall(PetscPrintf(comm, ", ")); 18319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "%g", (double)err[f])); 18321878804aSMatthew G. Knepley } 18339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "]\n")); 1834f017bcb6SMatthew G. Knepley } 1835f017bcb6SMatthew G. Knepley } else { 18369566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, t, exacts, ectxs, u, &err[0])); 1837f017bcb6SMatthew G. Knepley if (tol >= 0.0) { 183808401ef6SPierre Jolivet PetscCheck(err[0] <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Error %g exceeds tolerance %g", (double)err[0], (double)tol); 1839b878532fSMatthew G. Knepley } else if (error) { 1840b878532fSMatthew G. Knepley error[0] = err[0]; 1841f017bcb6SMatthew G. Knepley } else { 18429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Error: %g\n", (double)err[0])); 1843f017bcb6SMatthew G. Knepley } 1844f017bcb6SMatthew G. Knepley } 18459566063dSJacob Faibussowitsch PetscCall(PetscFree3(exacts, ectxs, err)); 18463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1847f017bcb6SMatthew G. Knepley } 1848f017bcb6SMatthew G. Knepley 1849f017bcb6SMatthew G. Knepley /*@C 1850f017bcb6SMatthew G. Knepley DMSNESCheckResidual - Check the residual of the exact solution 1851f017bcb6SMatthew G. Knepley 1852f017bcb6SMatthew G. Knepley Input Parameters: 1853f6dfbefdSBarry Smith + snes - the `SNES` object 1854f6dfbefdSBarry Smith . dm - the `DM` 1855f6dfbefdSBarry Smith . u - a `DM` vector 1856f0fc11ceSJed Brown - tol - A tolerance for the check, or -1 to print the results instead 1857f017bcb6SMatthew G. Knepley 1858f6dfbefdSBarry Smith Output Parameter: 18592fe279fdSBarry Smith . residual - The residual norm of the exact solution, or `NULL` 1860f3c94560SMatthew G. Knepley 1861f017bcb6SMatthew G. Knepley Level: developer 1862f017bcb6SMatthew G. Knepley 1863420bcc1bSBarry Smith .seealso: [](ch_snes), `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckJacobian()` 1864f017bcb6SMatthew G. Knepley @*/ 1865d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckResidual(SNES snes, DM dm, Vec u, PetscReal tol, PetscReal *residual) 1866d71ae5a4SJacob Faibussowitsch { 1867f017bcb6SMatthew G. Knepley MPI_Comm comm; 1868f017bcb6SMatthew G. Knepley Vec r; 1869f017bcb6SMatthew G. Knepley PetscReal res; 1870f017bcb6SMatthew G. Knepley 1871f017bcb6SMatthew G. Knepley PetscFunctionBegin; 1872f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1873f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1874f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 18754f572ea9SToby Isaac if (residual) PetscAssertPointer(residual, 5); 18769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 18779566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL)); 18789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 18799566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r)); 18809566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res)); 1881f017bcb6SMatthew G. Knepley if (tol >= 0.0) { 188208401ef6SPierre Jolivet PetscCheck(res <= tol, comm, PETSC_ERR_ARG_WRONG, "L_2 Residual %g exceeds tolerance %g", (double)res, (double)tol); 1883b878532fSMatthew G. Knepley } else if (residual) { 1884b878532fSMatthew G. Knepley *residual = res; 1885f017bcb6SMatthew G. Knepley } else { 18869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Residual: %g\n", (double)res)); 188793ec0da9SPierre Jolivet PetscCall(VecFilter(r, 1.0e-10)); 18889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r, "Initial Residual")); 18899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)r, "res_")); 18909566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(r, NULL, "-vec_view")); 1891f017bcb6SMatthew G. Knepley } 18929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 18933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1894f017bcb6SMatthew G. Knepley } 1895f017bcb6SMatthew G. Knepley 1896f017bcb6SMatthew G. Knepley /*@C 1897f3c94560SMatthew G. Knepley DMSNESCheckJacobian - Check the Jacobian of the exact solution against the residual using the Taylor Test 1898f017bcb6SMatthew G. Knepley 1899f017bcb6SMatthew G. Knepley Input Parameters: 1900f6dfbefdSBarry Smith + snes - the `SNES` object 1901f6dfbefdSBarry Smith . dm - the `DM` 1902f6dfbefdSBarry Smith . u - a `DM` vector 1903f0fc11ceSJed Brown - tol - A tolerance for the check, or -1 to print the results instead 1904f017bcb6SMatthew G. Knepley 1905f3c94560SMatthew G. Knepley Output Parameters: 19062fe279fdSBarry Smith + isLinear - Flag indicaing that the function looks linear, or `NULL` 19072fe279fdSBarry Smith - convRate - The rate of convergence of the linear model, or `NULL` 1908f3c94560SMatthew G. Knepley 1909f017bcb6SMatthew G. Knepley Level: developer 1910f017bcb6SMatthew G. Knepley 1911420bcc1bSBarry Smith .seealso: [](ch_snes), `DNSNESCheckFromOptions()`, `DMSNESCheckDiscretization()`, `DMSNESCheckResidual()` 1912f017bcb6SMatthew G. Knepley @*/ 1913d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckJacobian(SNES snes, DM dm, Vec u, PetscReal tol, PetscBool *isLinear, PetscReal *convRate) 1914d71ae5a4SJacob Faibussowitsch { 1915f017bcb6SMatthew G. Knepley MPI_Comm comm; 1916f017bcb6SMatthew G. Knepley PetscDS ds; 1917f017bcb6SMatthew G. Knepley Mat J, M; 1918f017bcb6SMatthew G. Knepley MatNullSpace nullspace; 1919f3c94560SMatthew G. Knepley PetscReal slope, intercept; 19207f96f943SMatthew G. Knepley PetscBool hasJac, hasPrec, isLin = PETSC_FALSE; 1921f017bcb6SMatthew G. Knepley 1922f017bcb6SMatthew G. Knepley PetscFunctionBegin; 1923f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1924f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 1925f3c94560SMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 3); 19264f572ea9SToby Isaac if (isLinear) PetscAssertPointer(isLinear, 5); 19274f572ea9SToby Isaac if (convRate) PetscAssertPointer(convRate, 6); 19289566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 19299566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, 0.0, u, NULL)); 1930f017bcb6SMatthew G. Knepley /* Create and view matrices */ 19319566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J)); 19329566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 19339566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobian(ds, &hasJac)); 19349566063dSJacob Faibussowitsch PetscCall(PetscDSHasJacobianPreconditioner(ds, &hasPrec)); 1935282e7bb4SMatthew G. Knepley if (hasJac && hasPrec) { 19369566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &M)); 19379566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, J, M)); 19389566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)M, "Preconditioning Matrix")); 19399566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "jacpre_")); 19409566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(M, NULL, "-mat_view")); 19419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&M)); 1942282e7bb4SMatthew G. Knepley } else { 19439566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, J, J)); 1944282e7bb4SMatthew G. Knepley } 19459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)J, "Jacobian")); 19469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)J, "jac_")); 19479566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(J, NULL, "-mat_view")); 1948f017bcb6SMatthew G. Knepley /* Check nullspace */ 19499566063dSJacob Faibussowitsch PetscCall(MatGetNullSpace(J, &nullspace)); 1950f017bcb6SMatthew G. Knepley if (nullspace) { 19511878804aSMatthew G. Knepley PetscBool isNull; 19529566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(nullspace, J, &isNull)); 195328b400f6SJacob Faibussowitsch PetscCheck(isNull, comm, PETSC_ERR_PLIB, "The null space calculated for the system operator is invalid."); 19541878804aSMatthew G. Knepley } 1955f017bcb6SMatthew G. Knepley /* Taylor test */ 1956f017bcb6SMatthew G. Knepley { 1957f017bcb6SMatthew G. Knepley PetscRandom rand; 1958f017bcb6SMatthew G. Knepley Vec du, uhat, r, rhat, df; 1959f3c94560SMatthew G. Knepley PetscReal h; 1960f017bcb6SMatthew G. Knepley PetscReal *es, *hs, *errors; 1961f017bcb6SMatthew G. Knepley PetscReal hMax = 1.0, hMin = 1e-6, hMult = 0.1; 1962f017bcb6SMatthew G. Knepley PetscInt Nv, v; 1963f017bcb6SMatthew G. Knepley 1964f017bcb6SMatthew G. Knepley /* Choose a perturbation direction */ 19659566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand)); 19669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &du)); 19679566063dSJacob Faibussowitsch PetscCall(VecSetRandom(du, rand)); 19689566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 19699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &df)); 19709566063dSJacob Faibussowitsch PetscCall(MatMult(J, du, df)); 1971f017bcb6SMatthew G. Knepley /* Evaluate residual at u, F(u), save in vector r */ 19729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 19739566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r)); 1974f017bcb6SMatthew G. Knepley /* Look at the convergence of our Taylor approximation as we approach u */ 19759371c9d4SSatish Balay for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) 19769371c9d4SSatish Balay ; 19779566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(Nv, &es, Nv, &hs, Nv, &errors)); 19789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &uhat)); 19799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &rhat)); 1980f017bcb6SMatthew G. Knepley for (h = hMax, Nv = 0; h >= hMin; h *= hMult, ++Nv) { 19819566063dSJacob Faibussowitsch PetscCall(VecWAXPY(uhat, h, du, u)); 1982f017bcb6SMatthew G. Knepley /* F(\hat u) \approx F(u) + J(u) (uhat - u) = F(u) + h * J(u) du */ 19839566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, uhat, rhat)); 19849566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(rhat, -1.0, -h, 1.0, r, df)); 19859566063dSJacob Faibussowitsch PetscCall(VecNorm(rhat, NORM_2, &errors[Nv])); 1986f017bcb6SMatthew G. Knepley 198738d69901SMatthew G. Knepley es[Nv] = errors[Nv] == 0 ? -16.0 : PetscLog10Real(errors[Nv]); 1988f017bcb6SMatthew G. Knepley hs[Nv] = PetscLog10Real(h); 1989f017bcb6SMatthew G. Knepley } 19909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uhat)); 19919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rhat)); 19929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&df)); 19939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 19949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&du)); 1995f017bcb6SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 1996f017bcb6SMatthew G. Knepley if ((tol >= 0) && (errors[v] > tol)) break; 1997f017bcb6SMatthew G. Knepley else if (errors[v] > PETSC_SMALL) break; 1998f017bcb6SMatthew G. Knepley } 1999f3c94560SMatthew G. Knepley if (v == Nv) isLin = PETSC_TRUE; 20009566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nv, hs, es, &slope, &intercept)); 20019566063dSJacob Faibussowitsch PetscCall(PetscFree3(es, hs, errors)); 2002f017bcb6SMatthew G. Knepley /* Slope should be about 2 */ 2003f017bcb6SMatthew G. Knepley if (tol >= 0) { 20040b121fc5SBarry Smith PetscCheck(isLin || PetscAbsReal(2 - slope) <= tol, comm, PETSC_ERR_ARG_WRONG, "Taylor approximation convergence rate should be 2, not %0.2f", (double)slope); 2005b878532fSMatthew G. Knepley } else if (isLinear || convRate) { 2006b878532fSMatthew G. Knepley if (isLinear) *isLinear = isLin; 2007b878532fSMatthew G. Knepley if (convRate) *convRate = slope; 2008f017bcb6SMatthew G. Knepley } else { 20099566063dSJacob Faibussowitsch if (!isLin) PetscCall(PetscPrintf(comm, "Taylor approximation converging at order %3.2f\n", (double)slope)); 20109566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function appears to be linear\n")); 2011f017bcb6SMatthew G. Knepley } 2012f017bcb6SMatthew G. Knepley } 20139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 20143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20151878804aSMatthew G. Knepley } 20161878804aSMatthew G. Knepley 2017d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheck_Internal(SNES snes, DM dm, Vec u) 2018d71ae5a4SJacob Faibussowitsch { 2019f017bcb6SMatthew G. Knepley PetscFunctionBegin; 20209566063dSJacob Faibussowitsch PetscCall(DMSNESCheckDiscretization(snes, dm, 0.0, u, -1.0, NULL)); 20219566063dSJacob Faibussowitsch PetscCall(DMSNESCheckResidual(snes, dm, u, -1.0, NULL)); 20229566063dSJacob Faibussowitsch PetscCall(DMSNESCheckJacobian(snes, dm, u, -1.0, NULL, NULL)); 20233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2024f017bcb6SMatthew G. Knepley } 2025f017bcb6SMatthew G. Knepley 2026bee9a294SMatthew G. Knepley /*@C 2027bee9a294SMatthew G. Knepley DMSNESCheckFromOptions - Check the residual and Jacobian functions using the exact solution by outputting some diagnostic information 2028bee9a294SMatthew G. Knepley 2029bee9a294SMatthew G. Knepley Input Parameters: 2030f6dfbefdSBarry Smith + snes - the `SNES` object 2031f6dfbefdSBarry Smith - u - representative `SNES` vector 20327f96f943SMatthew G. Knepley 20332fe279fdSBarry Smith Level: developer 20342fe279fdSBarry Smith 2035f6dfbefdSBarry Smith Note: 2036420bcc1bSBarry Smith The user must call `PetscDSSetExactSolution()` before this call 2037bee9a294SMatthew G. Knepley 2038420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `DM` 2039bee9a294SMatthew G. Knepley @*/ 2040d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSNESCheckFromOptions(SNES snes, Vec u) 2041d71ae5a4SJacob Faibussowitsch { 20421878804aSMatthew G. Knepley DM dm; 20431878804aSMatthew G. Knepley Vec sol; 20441878804aSMatthew G. Knepley PetscBool check; 20451878804aSMatthew G. Knepley 20461878804aSMatthew G. Knepley PetscFunctionBegin; 20479566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-dmsnes_check", &check)); 20483ba16761SJacob Faibussowitsch if (!check) PetscFunctionReturn(PETSC_SUCCESS); 20499566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 20509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &sol)); 20519566063dSJacob Faibussowitsch PetscCall(SNESSetSolution(snes, sol)); 20529566063dSJacob Faibussowitsch PetscCall(DMSNESCheck_Internal(snes, dm, sol)); 20539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sol)); 20543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2055552f7358SJed Brown } 2056