xref: /petsc/src/snes/tutorials/ex7.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
17b197a28SMatthew G. Knepley static char help[] = "Fermions on a hypercubic lattice.\n\n";
27b197a28SMatthew G. Knepley 
37b197a28SMatthew G. Knepley #include <petscdmplex.h>
47b197a28SMatthew G. Knepley #include <petscsnes.h>
57b197a28SMatthew G. Knepley 
67b197a28SMatthew G. Knepley /* Common operations:
77b197a28SMatthew G. Knepley 
87b197a28SMatthew G. Knepley  - View the input \psi as ASCII in lexicographic order: -psi_view
97b197a28SMatthew G. Knepley */
107b197a28SMatthew G. Knepley 
117b197a28SMatthew G. Knepley const PetscReal M = 1.0;
127b197a28SMatthew G. Knepley 
137b197a28SMatthew G. Knepley typedef struct {
14fc229b03SMatthew G. Knepley   PetscBool printTr; // Print the traversal
15fc229b03SMatthew G. Knepley   PetscBool usePV;   // Use Pauli-Villars preconditioning
167b197a28SMatthew G. Knepley } AppCtx;
177b197a28SMatthew G. Knepley 
ProcessOptions(MPI_Comm comm,AppCtx * options)187b197a28SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
197b197a28SMatthew G. Knepley {
207b197a28SMatthew G. Knepley   PetscFunctionBegin;
21fc229b03SMatthew G. Knepley   options->printTr = PETSC_FALSE;
227b197a28SMatthew G. Knepley   options->usePV   = PETSC_TRUE;
237b197a28SMatthew G. Knepley 
247b197a28SMatthew G. Knepley   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
25fc229b03SMatthew G. Knepley   PetscCall(PetscOptionsBool("-print_traversal", "Print the mesh traversal", "ex1.c", options->printTr, &options->printTr, NULL));
267b197a28SMatthew G. Knepley   PetscCall(PetscOptionsBool("-use_pv", "Use Pauli-Villars preconditioning", "ex1.c", options->usePV, &options->usePV, NULL));
277b197a28SMatthew G. Knepley   PetscOptionsEnd();
283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297b197a28SMatthew G. Knepley }
307b197a28SMatthew G. Knepley 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)317b197a28SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
327b197a28SMatthew G. Knepley {
337b197a28SMatthew G. Knepley   PetscFunctionBegin;
347b197a28SMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
357b197a28SMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
367b197a28SMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
377b197a28SMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
397b197a28SMatthew G. Knepley }
407b197a28SMatthew G. Knepley 
SetupDiscretization(DM dm,AppCtx * ctx)417b197a28SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
427b197a28SMatthew G. Knepley {
437b197a28SMatthew G. Knepley   PetscSection s;
447b197a28SMatthew G. Knepley   PetscInt     vStart, vEnd, v;
457b197a28SMatthew G. Knepley 
467b197a28SMatthew G. Knepley   PetscFunctionBegin;
477b197a28SMatthew G. Knepley   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
487b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
497b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
507b197a28SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
517b197a28SMatthew G. Knepley     PetscCall(PetscSectionSetDof(s, v, 12));
527b197a28SMatthew G. Knepley     /* TODO Divide the values into fields/components */
537b197a28SMatthew G. Knepley   }
547b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetUp(s));
557b197a28SMatthew G. Knepley   PetscCall(DMSetLocalSection(dm, s));
567b197a28SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&s));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
587b197a28SMatthew G. Knepley }
597b197a28SMatthew G. Knepley 
SetupAuxDiscretization(DM dm,AppCtx * user)607b197a28SMatthew G. Knepley static PetscErrorCode SetupAuxDiscretization(DM dm, AppCtx *user)
617b197a28SMatthew G. Knepley {
627b197a28SMatthew G. Knepley   DM           dmAux, coordDM;
637b197a28SMatthew G. Knepley   PetscSection s;
647b197a28SMatthew G. Knepley   Vec          gauge;
657b197a28SMatthew G. Knepley   PetscInt     eStart, eEnd, e;
667b197a28SMatthew G. Knepley 
677b197a28SMatthew G. Knepley   PetscFunctionBegin;
687b197a28SMatthew G. Knepley   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
697b197a28SMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &coordDM));
707b197a28SMatthew G. Knepley   PetscCall(DMClone(dm, &dmAux));
717b197a28SMatthew G. Knepley   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
727b197a28SMatthew G. Knepley   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
737b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
747b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetChart(s, eStart, eEnd));
757b197a28SMatthew G. Knepley   for (e = eStart; e < eEnd; ++e) {
767b197a28SMatthew G. Knepley     /* TODO Should we store the whole SU(3) matrix, or the symmetric part? */
777b197a28SMatthew G. Knepley     PetscCall(PetscSectionSetDof(s, e, 9));
787b197a28SMatthew G. Knepley   }
797b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetUp(s));
807b197a28SMatthew G. Knepley   PetscCall(DMSetLocalSection(dmAux, s));
817b197a28SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&s));
827b197a28SMatthew G. Knepley   PetscCall(DMCreateLocalVector(dmAux, &gauge));
837b197a28SMatthew G. Knepley   PetscCall(DMDestroy(&dmAux));
847b197a28SMatthew G. Knepley   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, gauge));
857b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&gauge));
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
877b197a28SMatthew G. Knepley }
887b197a28SMatthew G. Knepley 
PrintVertex(DM dm,PetscInt v)897b197a28SMatthew G. Knepley static PetscErrorCode PrintVertex(DM dm, PetscInt v)
907b197a28SMatthew G. Knepley {
917b197a28SMatthew G. Knepley   MPI_Comm       comm;
927b197a28SMatthew G. Knepley   PetscContainer c;
937b197a28SMatthew G. Knepley   PetscInt      *extent;
947b197a28SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, sum;
957b197a28SMatthew G. Knepley 
967b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
977b197a28SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
987b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
997b197a28SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
1007b197a28SMatthew G. Knepley   PetscCall(PetscObjectQuery((PetscObject)dm, "_extent", (PetscObject *)&c));
101*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(c, &extent));
1027b197a28SMatthew G. Knepley   sum = 1;
1037b197a28SMatthew G. Knepley   PetscCall(PetscPrintf(comm, "Vertex %" PetscInt_FMT ":", v));
1047b197a28SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
1057b197a28SMatthew G. Knepley     PetscCall(PetscPrintf(comm, " %" PetscInt_FMT, (v / sum) % extent[d]));
1067b197a28SMatthew G. Knepley     if (d < dim) sum *= extent[d];
1077b197a28SMatthew G. Knepley   }
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1097b197a28SMatthew G. Knepley }
1107b197a28SMatthew G. Knepley 
1117b197a28SMatthew G. Knepley // Apply \gamma_\mu
ComputeGamma(PetscInt d,PetscInt ldx,PetscScalar f[])1127b197a28SMatthew G. Knepley static PetscErrorCode ComputeGamma(PetscInt d, PetscInt ldx, PetscScalar f[])
1137b197a28SMatthew G. Knepley {
1147b197a28SMatthew G. Knepley   const PetscScalar fin[4] = {f[0 * ldx], f[1 * ldx], f[2 * ldx], f[3 * ldx]};
1157b197a28SMatthew G. Knepley 
1167b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
1177b197a28SMatthew G. Knepley   switch (d) {
1187b197a28SMatthew G. Knepley   case 0:
1197b197a28SMatthew G. Knepley     f[0 * ldx] = PETSC_i * fin[3];
1207b197a28SMatthew G. Knepley     f[1 * ldx] = PETSC_i * fin[2];
1217b197a28SMatthew G. Knepley     f[2 * ldx] = -PETSC_i * fin[1];
1227b197a28SMatthew G. Knepley     f[3 * ldx] = -PETSC_i * fin[0];
1237b197a28SMatthew G. Knepley     break;
1247b197a28SMatthew G. Knepley   case 1:
1257b197a28SMatthew G. Knepley     f[0 * ldx] = -fin[3];
1267b197a28SMatthew G. Knepley     f[1 * ldx] = fin[2];
1277b197a28SMatthew G. Knepley     f[2 * ldx] = fin[1];
1287b197a28SMatthew G. Knepley     f[3 * ldx] = -fin[0];
1297b197a28SMatthew G. Knepley     break;
1307b197a28SMatthew G. Knepley   case 2:
1317b197a28SMatthew G. Knepley     f[0 * ldx] = PETSC_i * fin[2];
1327b197a28SMatthew G. Knepley     f[1 * ldx] = -PETSC_i * fin[3];
1337b197a28SMatthew G. Knepley     f[2 * ldx] = -PETSC_i * fin[0];
1347b197a28SMatthew G. Knepley     f[3 * ldx] = PETSC_i * fin[1];
1357b197a28SMatthew G. Knepley     break;
1367b197a28SMatthew G. Knepley   case 3:
1377b197a28SMatthew G. Knepley     f[0 * ldx] = fin[2];
1387b197a28SMatthew G. Knepley     f[1 * ldx] = fin[3];
1397b197a28SMatthew G. Knepley     f[2 * ldx] = fin[0];
1407b197a28SMatthew G. Knepley     f[3 * ldx] = fin[1];
1417b197a28SMatthew G. Knepley     break;
1427b197a28SMatthew G. Knepley   default:
1437b197a28SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Direction for gamma %" PetscInt_FMT " not in [0, 4)", d);
1447b197a28SMatthew G. Knepley   }
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1467b197a28SMatthew G. Knepley }
1477b197a28SMatthew G. Knepley 
1487b197a28SMatthew G. Knepley // Apply (1 \pm \gamma_\mu)/2
ComputeGammaFactor(PetscInt d,PetscBool forward,PetscInt ldx,PetscScalar f[])1497b197a28SMatthew G. Knepley static PetscErrorCode ComputeGammaFactor(PetscInt d, PetscBool forward, PetscInt ldx, PetscScalar f[])
1507b197a28SMatthew G. Knepley {
1517b197a28SMatthew G. Knepley   const PetscReal   sign   = forward ? -1. : 1.;
1527b197a28SMatthew G. Knepley   const PetscScalar fin[4] = {f[0 * ldx], f[1 * ldx], f[2 * ldx], f[3 * ldx]};
1537b197a28SMatthew G. Knepley 
1547b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
1557b197a28SMatthew G. Knepley   switch (d) {
1567b197a28SMatthew G. Knepley   case 0:
1577b197a28SMatthew G. Knepley     f[0 * ldx] += sign * PETSC_i * fin[3];
1587b197a28SMatthew G. Knepley     f[1 * ldx] += sign * PETSC_i * fin[2];
1597b197a28SMatthew G. Knepley     f[2 * ldx] += sign * -PETSC_i * fin[1];
1607b197a28SMatthew G. Knepley     f[3 * ldx] += sign * -PETSC_i * fin[0];
1617b197a28SMatthew G. Knepley     break;
1627b197a28SMatthew G. Knepley   case 1:
1637b197a28SMatthew G. Knepley     f[0 * ldx] += -sign * fin[3];
1647b197a28SMatthew G. Knepley     f[1 * ldx] += sign * fin[2];
1657b197a28SMatthew G. Knepley     f[2 * ldx] += sign * fin[1];
1667b197a28SMatthew G. Knepley     f[3 * ldx] += -sign * fin[0];
1677b197a28SMatthew G. Knepley     break;
1687b197a28SMatthew G. Knepley   case 2:
1697b197a28SMatthew G. Knepley     f[0 * ldx] += sign * PETSC_i * fin[2];
1707b197a28SMatthew G. Knepley     f[1 * ldx] += sign * -PETSC_i * fin[3];
1717b197a28SMatthew G. Knepley     f[2 * ldx] += sign * -PETSC_i * fin[0];
1727b197a28SMatthew G. Knepley     f[3 * ldx] += sign * PETSC_i * fin[1];
1737b197a28SMatthew G. Knepley     break;
1747b197a28SMatthew G. Knepley   case 3:
1757b197a28SMatthew G. Knepley     f[0 * ldx] += sign * fin[2];
1767b197a28SMatthew G. Knepley     f[1 * ldx] += sign * fin[3];
1777b197a28SMatthew G. Knepley     f[2 * ldx] += sign * fin[0];
1787b197a28SMatthew G. Knepley     f[3 * ldx] += sign * fin[1];
1797b197a28SMatthew G. Knepley     break;
1807b197a28SMatthew G. Knepley   default:
1817b197a28SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Direction for gamma %" PetscInt_FMT " not in [0, 4)", d);
1827b197a28SMatthew G. Knepley   }
1837b197a28SMatthew G. Knepley   f[0 * ldx] *= 0.5;
1847b197a28SMatthew G. Knepley   f[1 * ldx] *= 0.5;
1857b197a28SMatthew G. Knepley   f[2 * ldx] *= 0.5;
1867b197a28SMatthew G. Knepley   f[3 * ldx] *= 0.5;
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1887b197a28SMatthew G. Knepley }
1897b197a28SMatthew G. Knepley 
1907b197a28SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>
1917b197a28SMatthew G. Knepley 
1927b197a28SMatthew G. Knepley // ComputeAction() sums the action of 1/2 (1 \pm \gamma_\mu) U \psi into f
ComputeAction(PetscInt d,PetscBool forward,const PetscScalar U[],const PetscScalar psi[],PetscScalar f[])1937b197a28SMatthew G. Knepley static PetscErrorCode ComputeAction(PetscInt d, PetscBool forward, const PetscScalar U[], const PetscScalar psi[], PetscScalar f[])
1947b197a28SMatthew G. Knepley {
1957b197a28SMatthew G. Knepley   PetscScalar tmp[12];
1967b197a28SMatthew G. Knepley 
1977b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
1987b197a28SMatthew G. Knepley   // Apply U
1997b197a28SMatthew G. Knepley   for (PetscInt beta = 0; beta < 4; ++beta) {
2007b197a28SMatthew G. Knepley     if (forward) DMPlex_Mult3D_Internal(U, 1, &psi[beta * 3], &tmp[beta * 3]);
2017b197a28SMatthew G. Knepley     else DMPlex_MultTranspose3D_Internal(U, 1, &psi[beta * 3], &tmp[beta * 3]);
2027b197a28SMatthew G. Knepley   }
2037b197a28SMatthew G. Knepley   // Apply (1 \pm \gamma_\mu)/2 to each color
2047b197a28SMatthew G. Knepley   for (PetscInt c = 0; c < 3; ++c) PetscCall(ComputeGammaFactor(d, forward, 3, &tmp[c]));
2057b197a28SMatthew G. Knepley   // Note that we are subtracting this contribution
2067b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < 12; ++i) f[i] -= tmp[i];
2073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2087b197a28SMatthew G. Knepley }
2097b197a28SMatthew G. Knepley 
2107b197a28SMatthew G. Knepley /*
2117b197a28SMatthew G. Knepley   The assembly loop runs over vertices. Each vertex has 2d edges in its support. The edges are ordered first by the dimension along which they run, and second from smaller to larger index, expect for edges which loop periodically. The vertices on edges are also ordered from smaller to larger index except for periodic edges.
2127b197a28SMatthew G. Knepley */
ComputeResidualLocal(DM dm,Vec u,Vec f)213fc229b03SMatthew G. Knepley static PetscErrorCode ComputeResidualLocal(DM dm, Vec u, Vec f)
2147b197a28SMatthew G. Knepley {
215de19e7feSJoseph Pusztay   DM                 dmAux;
216de19e7feSJoseph Pusztay   Vec                gauge;
217de19e7feSJoseph Pusztay   PetscSection       s, sGauge;
2187b197a28SMatthew G. Knepley   const PetscScalar *ua;
219de19e7feSJoseph Pusztay   PetscScalar       *fa, *link;
2207b197a28SMatthew G. Knepley   PetscInt           dim, vStart, vEnd;
2217b197a28SMatthew G. Knepley 
2227b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
2237b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
2247b197a28SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
2257b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2267b197a28SMatthew G. Knepley   PetscCall(VecGetArrayRead(u, &ua));
2277b197a28SMatthew G. Knepley   PetscCall(VecGetArray(f, &fa));
228de19e7feSJoseph Pusztay 
229de19e7feSJoseph Pusztay   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &gauge));
230de19e7feSJoseph Pusztay   PetscCall(VecGetDM(gauge, &dmAux));
231de19e7feSJoseph Pusztay   PetscCall(DMGetLocalSection(dmAux, &sGauge));
232de19e7feSJoseph Pusztay   PetscCall(VecGetArray(gauge, &link));
2337b197a28SMatthew G. Knepley   // Loop over y
2347b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
2357b197a28SMatthew G. Knepley     const PetscInt *supp;
2367b197a28SMatthew G. Knepley     PetscInt        xdof, xoff;
2377b197a28SMatthew G. Knepley 
2387b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, v, &supp));
2397b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetDof(s, v, &xdof));
2407b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetOffset(s, v, &xoff));
2417b197a28SMatthew G. Knepley     // Diagonal
2427b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < xdof; ++i) fa[xoff + i] += (M + 4) * ua[xoff + i];
2437b197a28SMatthew G. Knepley     // Loop over mu
2447b197a28SMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) {
2457b197a28SMatthew G. Knepley       const PetscInt *cone;
246de19e7feSJoseph Pusztay       PetscInt        yoff, goff;
2477b197a28SMatthew G. Knepley 
2487b197a28SMatthew G. Knepley       // Left action -(1 + \gamma_\mu)/2 \otimes U^\dagger_\mu(y) \delta_{x - \mu,y} \psi(y)
2497b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[2 * d + 0], &cone));
2507b197a28SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(s, cone[0], &yoff));
251de19e7feSJoseph Pusztay       PetscCall(PetscSectionGetOffset(sGauge, supp[2 * d + 0], &goff));
252de19e7feSJoseph Pusztay       PetscCall(ComputeAction(d, PETSC_FALSE, &link[goff], &ua[yoff], &fa[xoff]));
2537b197a28SMatthew G. Knepley       // Right edge -(1 - \gamma_\mu)/2 \otimes U_\mu(x) \delta_{x + \mu,y} \psi(y)
2547b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[2 * d + 1], &cone));
2557b197a28SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(s, cone[1], &yoff));
256de19e7feSJoseph Pusztay       PetscCall(PetscSectionGetOffset(sGauge, supp[2 * d + 1], &goff));
257de19e7feSJoseph Pusztay       PetscCall(ComputeAction(d, PETSC_TRUE, &link[goff], &ua[yoff], &fa[xoff]));
2587b197a28SMatthew G. Knepley     }
2597b197a28SMatthew G. Knepley   }
2607b197a28SMatthew G. Knepley   PetscCall(VecRestoreArray(f, &fa));
261de19e7feSJoseph Pusztay   PetscCall(VecRestoreArray(gauge, &link));
2627b197a28SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(u, &ua));
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2647b197a28SMatthew G. Knepley }
2657b197a28SMatthew G. Knepley 
ComputeResidual(DM dm,Vec u,Vec f)266fc229b03SMatthew G. Knepley static PetscErrorCode ComputeResidual(DM dm, Vec u, Vec f)
267fc229b03SMatthew G. Knepley {
268fc229b03SMatthew G. Knepley   Vec lu, lf;
269fc229b03SMatthew G. Knepley 
270fc229b03SMatthew G. Knepley   PetscFunctionBegin;
271fc229b03SMatthew G. Knepley   PetscCall(DMGetLocalVector(dm, &lu));
272fc229b03SMatthew G. Knepley   PetscCall(DMGetLocalVector(dm, &lf));
273fc229b03SMatthew G. Knepley   PetscCall(DMGlobalToLocalBegin(dm, u, INSERT_VALUES, lu));
274fc229b03SMatthew G. Knepley   PetscCall(DMGlobalToLocalEnd(dm, u, INSERT_VALUES, lu));
275fc229b03SMatthew G. Knepley   PetscCall(VecSet(lf, 0.));
276fc229b03SMatthew G. Knepley   PetscCall(ComputeResidualLocal(dm, lu, lf));
277fc229b03SMatthew G. Knepley   PetscCall(DMLocalToGlobalBegin(dm, lf, INSERT_VALUES, f));
278fc229b03SMatthew G. Knepley   PetscCall(DMLocalToGlobalEnd(dm, lf, INSERT_VALUES, f));
279fc229b03SMatthew G. Knepley   PetscCall(DMRestoreLocalVector(dm, &lu));
280fc229b03SMatthew G. Knepley   PetscCall(DMRestoreLocalVector(dm, &lf));
281fc229b03SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
282fc229b03SMatthew G. Knepley }
283fc229b03SMatthew G. Knepley 
CreateJacobian(DM dm,Mat * J)2847b197a28SMatthew G. Knepley static PetscErrorCode CreateJacobian(DM dm, Mat *J)
2857b197a28SMatthew G. Knepley {
2867b197a28SMatthew G. Knepley   PetscFunctionBegin;
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2887b197a28SMatthew G. Knepley }
2897b197a28SMatthew G. Knepley 
ComputeJacobian(DM dm,Vec u,Mat J)2907b197a28SMatthew G. Knepley static PetscErrorCode ComputeJacobian(DM dm, Vec u, Mat J)
2917b197a28SMatthew G. Knepley {
2927b197a28SMatthew G. Knepley   PetscFunctionBegin;
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2947b197a28SMatthew G. Knepley }
2957b197a28SMatthew G. Knepley 
PrintTraversal(DM dm)2967b197a28SMatthew G. Knepley static PetscErrorCode PrintTraversal(DM dm)
2977b197a28SMatthew G. Knepley {
2987b197a28SMatthew G. Knepley   MPI_Comm comm;
2997b197a28SMatthew G. Knepley   PetscInt vStart, vEnd;
3007b197a28SMatthew G. Knepley 
3017b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
3027b197a28SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3037b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3047b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
3057b197a28SMatthew G. Knepley     const PetscInt *supp;
3067b197a28SMatthew G. Knepley     PetscInt        Ns;
3077b197a28SMatthew G. Knepley 
3087b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupportSize(dm, v, &Ns));
3097b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, v, &supp));
3107b197a28SMatthew G. Knepley     PetscCall(PrintVertex(dm, v));
3117b197a28SMatthew G. Knepley     PetscCall(PetscPrintf(comm, "\n"));
3127b197a28SMatthew G. Knepley     for (PetscInt s = 0; s < Ns; ++s) {
3137b197a28SMatthew G. Knepley       const PetscInt *cone;
3147b197a28SMatthew G. Knepley 
3157b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[s], &cone));
3167b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "  Edge %" PetscInt_FMT ": ", supp[s]));
3177b197a28SMatthew G. Knepley       PetscCall(PrintVertex(dm, cone[0]));
3187b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, " -- "));
3197b197a28SMatthew G. Knepley       PetscCall(PrintVertex(dm, cone[1]));
3207b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "\n"));
3217b197a28SMatthew G. Knepley     }
3227b197a28SMatthew G. Knepley   }
3233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3247b197a28SMatthew G. Knepley }
3257b197a28SMatthew G. Knepley 
ComputeFFT(Mat FT,PetscInt Nc,Vec x,Vec p)3267b197a28SMatthew G. Knepley static PetscErrorCode ComputeFFT(Mat FT, PetscInt Nc, Vec x, Vec p)
3277b197a28SMatthew G. Knepley {
328fc229b03SMatthew G. Knepley   Vec     *xComp, *pComp, xTmp, pTmp;
3297b197a28SMatthew G. Knepley   PetscInt n, N;
3307b197a28SMatthew G. Knepley 
3317b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
3327b197a28SMatthew G. Knepley   PetscCall(PetscMalloc2(Nc, &xComp, Nc, &pComp));
3337b197a28SMatthew G. Knepley   PetscCall(VecGetLocalSize(x, &n));
3347b197a28SMatthew G. Knepley   PetscCall(VecGetSize(x, &N));
3357b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
3367b197a28SMatthew G. Knepley     const char *vtype;
3377b197a28SMatthew G. Knepley 
3387b197a28SMatthew G. Knepley     // HACK: Make these from another DM up front
3397b197a28SMatthew G. Knepley     PetscCall(VecCreate(PetscObjectComm((PetscObject)x), &xComp[i]));
3407b197a28SMatthew G. Knepley     PetscCall(VecGetType(x, &vtype));
3417b197a28SMatthew G. Knepley     PetscCall(VecSetType(xComp[i], vtype));
3427b197a28SMatthew G. Knepley     PetscCall(VecSetSizes(xComp[i], n / Nc, N / Nc));
3437b197a28SMatthew G. Knepley     PetscCall(VecDuplicate(xComp[i], &pComp[i]));
3447b197a28SMatthew G. Knepley   }
345fc229b03SMatthew G. Knepley   PetscCall(MatCreateVecsFFTW(FT, &xTmp, &pTmp, NULL));
3467b197a28SMatthew G. Knepley   PetscCall(VecStrideGatherAll(x, xComp, INSERT_VALUES));
347fc229b03SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
348fc229b03SMatthew G. Knepley     PetscCall(VecScatterPetscToFFTW(FT, xComp[i], xTmp));
349fc229b03SMatthew G. Knepley     PetscCall(MatMult(FT, xTmp, pTmp));
350fc229b03SMatthew G. Knepley     PetscCall(VecScatterFFTWToPetsc(FT, pTmp, pComp[i]));
351fc229b03SMatthew G. Knepley   }
3527b197a28SMatthew G. Knepley   PetscCall(VecStrideScatterAll(pComp, p, INSERT_VALUES));
3537b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
3547b197a28SMatthew G. Knepley     PetscCall(VecDestroy(&xComp[i]));
3557b197a28SMatthew G. Knepley     PetscCall(VecDestroy(&pComp[i]));
3567b197a28SMatthew G. Knepley   }
357fc229b03SMatthew G. Knepley   PetscCall(VecDestroy(&xTmp));
358fc229b03SMatthew G. Knepley   PetscCall(VecDestroy(&pTmp));
3597b197a28SMatthew G. Knepley   PetscCall(PetscFree2(xComp, pComp));
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3617b197a28SMatthew G. Knepley }
3627b197a28SMatthew G. Knepley 
363de19e7feSJoseph Pusztay // Sets each link to be the identity for the free field test
SetGauge_Identity(DM dm)364de19e7feSJoseph Pusztay static PetscErrorCode SetGauge_Identity(DM dm)
365de19e7feSJoseph Pusztay {
366de19e7feSJoseph Pusztay   DM           auxDM;
367de19e7feSJoseph Pusztay   Vec          auxVec;
368de19e7feSJoseph Pusztay   PetscSection s;
369de19e7feSJoseph Pusztay   PetscScalar  id[9] = {1., 0., 0., 0., 1., 0., 0., 0., 1.};
370de19e7feSJoseph Pusztay   PetscInt     eStart, eEnd;
371de19e7feSJoseph Pusztay 
372de19e7feSJoseph Pusztay   PetscFunctionBegin;
373de19e7feSJoseph Pusztay   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &auxVec));
374de19e7feSJoseph Pusztay   PetscCall(VecGetDM(auxVec, &auxDM));
375de19e7feSJoseph Pusztay   PetscCall(DMGetLocalSection(auxDM, &s));
376de19e7feSJoseph Pusztay   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
377aa624791SPierre Jolivet   for (PetscInt i = eStart; i < eEnd; ++i) PetscCall(VecSetValuesSection(auxVec, s, i, id, INSERT_VALUES));
378de19e7feSJoseph Pusztay   PetscCall(VecViewFromOptions(auxVec, NULL, "-gauge_view"));
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
380de19e7feSJoseph Pusztay }
381de19e7feSJoseph Pusztay 
3827b197a28SMatthew G. Knepley /*
3837b197a28SMatthew G. Knepley   Test the action of the Wilson operator in the free field case U = I,
3847b197a28SMatthew G. Knepley 
3857b197a28SMatthew G. Knepley     \eta(x) = D_W(x - y) \psi(y)
3867b197a28SMatthew G. Knepley 
3877b197a28SMatthew G. Knepley   The Wilson operator is a convolution for the free field, so we can check that by the convolution theorem
3887b197a28SMatthew G. Knepley 
3897b197a28SMatthew G. Knepley     \hat\eta(x) = \mathcal{F}(D_W(x - y) \psi(y))
3907b197a28SMatthew G. Knepley                 = \hat D_W(p) \mathcal{F}\psi(p)
3917b197a28SMatthew G. Knepley 
3927b197a28SMatthew G. Knepley   The Fourier series for the Wilson operator is
3937b197a28SMatthew G. Knepley 
3947b197a28SMatthew G. Knepley     M + \sum_\mu 2 \sin^2(p_\mu / 2) + i \gamma_\mu \sin(p_\mu)
3957b197a28SMatthew G. Knepley */
TestFreeField(DM dm)3967b197a28SMatthew G. Knepley static PetscErrorCode TestFreeField(DM dm)
3977b197a28SMatthew G. Knepley {
3987b197a28SMatthew G. Knepley   PetscSection       s;
3997b197a28SMatthew G. Knepley   Mat                FT;
4007b197a28SMatthew G. Knepley   Vec                psi, psiHat;
4017b197a28SMatthew G. Knepley   Vec                eta, etaHat;
4027b197a28SMatthew G. Knepley   Vec                DHat; // The product \hat D_w \hat psi
4037b197a28SMatthew G. Knepley   PetscRandom        r;
4047b197a28SMatthew G. Knepley   const PetscScalar *psih;
4057b197a28SMatthew G. Knepley   PetscScalar       *dh;
4067b197a28SMatthew G. Knepley   PetscReal         *coef, nrm;
4077b197a28SMatthew G. Knepley   const PetscInt    *extent, Nc = 12;
4087b197a28SMatthew G. Knepley   PetscInt           dim, V     = 1, vStart, vEnd;
4097b197a28SMatthew G. Knepley   PetscContainer     c;
4107b197a28SMatthew G. Knepley   PetscBool          constRhs = PETSC_FALSE;
411fc229b03SMatthew G. Knepley   PetscMPIInt        size;
4127b197a28SMatthew G. Knepley 
4137b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
414fc229b03SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
415fc229b03SMatthew G. Knepley   if (size > 1) PetscFunctionReturn(PETSC_SUCCESS);
4167b197a28SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-const_rhs", &constRhs, NULL));
4177b197a28SMatthew G. Knepley 
418de19e7feSJoseph Pusztay   PetscCall(SetGauge_Identity(dm));
4197b197a28SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
4207b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &psi));
4217b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)psi, "psi"));
4227b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &psiHat));
4237b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)psiHat, "psihat"));
4247b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &eta));
4257b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)eta, "eta"));
4267b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &etaHat));
4277b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)etaHat, "etahat"));
4287b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &DHat));
4297b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)DHat, "Dhat"));
4307b197a28SMatthew G. Knepley   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
4317b197a28SMatthew G. Knepley   PetscCall(PetscRandomSetType(r, PETSCRAND48));
4327b197a28SMatthew G. Knepley   if (constRhs) PetscCall(VecSet(psi, 1.));
4337b197a28SMatthew G. Knepley   else PetscCall(VecSetRandom(psi, r));
4347b197a28SMatthew G. Knepley   PetscCall(PetscRandomDestroy(&r));
4357b197a28SMatthew G. Knepley 
4367b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
4377b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4387b197a28SMatthew G. Knepley   PetscCall(PetscObjectQuery((PetscObject)dm, "_extent", (PetscObject *)&c));
439*2a8381b2SBarry Smith   PetscCall(PetscContainerGetPointer(c, &extent));
4407b197a28SMatthew G. Knepley   PetscCall(MatCreateFFT(PetscObjectComm((PetscObject)dm), dim, extent, MATFFTW, &FT));
4417b197a28SMatthew G. Knepley 
4427b197a28SMatthew G. Knepley   PetscCall(PetscMalloc1(dim, &coef));
4437b197a28SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
4447b197a28SMatthew G. Knepley     coef[d] = 2. * PETSC_PI / extent[d];
4457b197a28SMatthew G. Knepley     V *= extent[d];
4467b197a28SMatthew G. Knepley   }
4477b197a28SMatthew G. Knepley   PetscCall(ComputeResidual(dm, psi, eta));
448fc229b03SMatthew G. Knepley   PetscCall(VecViewFromOptions(psi, NULL, "-psi_view"));
4497b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(eta, NULL, "-eta_view"));
4507b197a28SMatthew G. Knepley   PetscCall(ComputeFFT(FT, Nc, psi, psiHat));
4517b197a28SMatthew G. Knepley   PetscCall(VecScale(psiHat, 1. / V));
4527b197a28SMatthew G. Knepley   PetscCall(ComputeFFT(FT, Nc, eta, etaHat));
4537b197a28SMatthew G. Knepley   PetscCall(VecScale(etaHat, 1. / V));
4547b197a28SMatthew G. Knepley   PetscCall(VecGetArrayRead(psiHat, &psih));
4557b197a28SMatthew G. Knepley   PetscCall(VecGetArray(DHat, &dh));
4567b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
4577b197a28SMatthew G. Knepley     PetscScalar tmp[12], tmp1 = 0.;
4587b197a28SMatthew G. Knepley     PetscInt    dof, off;
4597b197a28SMatthew G. Knepley 
4607b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetDof(s, v, &dof));
4617b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetOffset(s, v, &off));
4627b197a28SMatthew G. Knepley     for (PetscInt d = 0, prod = 1; d < dim; prod *= extent[d], ++d) {
4637b197a28SMatthew G. Knepley       const PetscInt idx = (v / prod) % extent[d];
4647b197a28SMatthew G. Knepley 
4657b197a28SMatthew G. Knepley       tmp1 += 2. * PetscSqr(PetscSinReal(0.5 * coef[d] * idx));
4667b197a28SMatthew G. Knepley       for (PetscInt i = 0; i < dof; ++i) tmp[i] = PETSC_i * PetscSinReal(coef[d] * idx) * psih[off + i];
4673ba16761SJacob Faibussowitsch       for (PetscInt c = 0; c < 3; ++c) PetscCall(ComputeGamma(d, 3, &tmp[c]));
4687b197a28SMatthew G. Knepley       for (PetscInt i = 0; i < dof; ++i) dh[off + i] += tmp[i];
4697b197a28SMatthew G. Knepley     }
4707b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < dof; ++i) dh[off + i] += (M + tmp1) * psih[off + i];
4717b197a28SMatthew G. Knepley   }
4727b197a28SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(psiHat, &psih));
4737b197a28SMatthew G. Knepley   PetscCall(VecRestoreArray(DHat, &dh));
4747b197a28SMatthew G. Knepley 
4757b197a28SMatthew G. Knepley   {
4767b197a28SMatthew G. Knepley     Vec     *etaComp, *DComp;
4777b197a28SMatthew G. Knepley     PetscInt n, N;
4787b197a28SMatthew G. Knepley 
4797b197a28SMatthew G. Knepley     PetscCall(PetscMalloc2(Nc, &etaComp, Nc, &DComp));
4807b197a28SMatthew G. Knepley     PetscCall(VecGetLocalSize(etaHat, &n));
4817b197a28SMatthew G. Knepley     PetscCall(VecGetSize(etaHat, &N));
4827b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
4837b197a28SMatthew G. Knepley       const char *vtype;
4847b197a28SMatthew G. Knepley 
4857b197a28SMatthew G. Knepley       // HACK: Make these from another DM up front
4867b197a28SMatthew G. Knepley       PetscCall(VecCreate(PetscObjectComm((PetscObject)etaHat), &etaComp[i]));
4877b197a28SMatthew G. Knepley       PetscCall(VecGetType(etaHat, &vtype));
4887b197a28SMatthew G. Knepley       PetscCall(VecSetType(etaComp[i], vtype));
4897b197a28SMatthew G. Knepley       PetscCall(VecSetSizes(etaComp[i], n / Nc, N / Nc));
4907b197a28SMatthew G. Knepley       PetscCall(VecDuplicate(etaComp[i], &DComp[i]));
4917b197a28SMatthew G. Knepley     }
4927b197a28SMatthew G. Knepley     PetscCall(VecStrideGatherAll(etaHat, etaComp, INSERT_VALUES));
4937b197a28SMatthew G. Knepley     PetscCall(VecStrideGatherAll(DHat, DComp, INSERT_VALUES));
4947b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
4957b197a28SMatthew G. Knepley       if (!i) {
4967b197a28SMatthew G. Knepley         PetscCall(VecViewFromOptions(etaComp[i], NULL, "-etahat_view"));
4977b197a28SMatthew G. Knepley         PetscCall(VecViewFromOptions(DComp[i], NULL, "-dhat_view"));
4987b197a28SMatthew G. Knepley       }
4997b197a28SMatthew G. Knepley       PetscCall(VecAXPY(etaComp[i], -1., DComp[i]));
5007b197a28SMatthew G. Knepley       PetscCall(VecNorm(etaComp[i], NORM_INFINITY, &nrm));
5017b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Slice %" PetscInt_FMT ": %g\n", i, (double)nrm));
5027b197a28SMatthew G. Knepley     }
5037b197a28SMatthew G. Knepley     PetscCall(VecStrideScatterAll(etaComp, etaHat, INSERT_VALUES));
5047b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
5057b197a28SMatthew G. Knepley       PetscCall(VecDestroy(&etaComp[i]));
5067b197a28SMatthew G. Knepley       PetscCall(VecDestroy(&DComp[i]));
5077b197a28SMatthew G. Knepley     }
5087b197a28SMatthew G. Knepley     PetscCall(PetscFree2(etaComp, DComp));
5097b197a28SMatthew G. Knepley     PetscCall(VecNorm(etaHat, NORM_INFINITY, &nrm));
5107b197a28SMatthew G. Knepley     PetscCheck(nrm < PETSC_SMALL, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Free field test failed: %g", (double)nrm);
5117b197a28SMatthew G. Knepley   }
5127b197a28SMatthew G. Knepley 
5137b197a28SMatthew G. Knepley   PetscCall(PetscFree(coef));
5147b197a28SMatthew G. Knepley   PetscCall(MatDestroy(&FT));
5157b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &psi));
5167b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &psiHat));
5177b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &eta));
5187b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &etaHat));
5197b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &DHat));
5203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5217b197a28SMatthew G. Knepley }
5227b197a28SMatthew G. Knepley 
main(int argc,char ** argv)5237b197a28SMatthew G. Knepley int main(int argc, char **argv)
5247b197a28SMatthew G. Knepley {
5257b197a28SMatthew G. Knepley   DM     dm;
5267b197a28SMatthew G. Knepley   Vec    u, f;
5277b197a28SMatthew G. Knepley   Mat    J;
5287b197a28SMatthew G. Knepley   AppCtx user;
5297b197a28SMatthew G. Knepley 
5307b197a28SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
5317b197a28SMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
5327b197a28SMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
5337b197a28SMatthew G. Knepley   PetscCall(SetupDiscretization(dm, &user));
5347b197a28SMatthew G. Knepley   PetscCall(SetupAuxDiscretization(dm, &user));
5357b197a28SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &u));
5367b197a28SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &f));
537fc229b03SMatthew G. Knepley   if (user.printTr) PetscCall(PrintTraversal(dm));
5387b197a28SMatthew G. Knepley   PetscCall(ComputeResidual(dm, u, f));
5397b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(f, NULL, "-res_view"));
5407b197a28SMatthew G. Knepley   PetscCall(CreateJacobian(dm, &J));
5417b197a28SMatthew G. Knepley   PetscCall(ComputeJacobian(dm, u, J));
5427b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
5437b197a28SMatthew G. Knepley   PetscCall(TestFreeField(dm));
5447b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&f));
5457b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&u));
5467b197a28SMatthew G. Knepley   PetscCall(DMDestroy(&dm));
5477b197a28SMatthew G. Knepley   PetscCall(PetscFinalize());
5487b197a28SMatthew G. Knepley   return 0;
5497b197a28SMatthew G. Knepley }
5507b197a28SMatthew G. Knepley 
5517b197a28SMatthew G. Knepley /*TEST
5527b197a28SMatthew G. Knepley 
5537b197a28SMatthew G. Knepley   build:
5547b197a28SMatthew G. Knepley     requires: complex
5557b197a28SMatthew G. Knepley 
5567b197a28SMatthew G. Knepley   test:
5577b197a28SMatthew G. Knepley     requires: fftw
5587b197a28SMatthew G. Knepley     suffix: dirac_free_field
5597b197a28SMatthew G. Knepley     args: -dm_plex_dim 4 -dm_plex_shape hypercubic -dm_plex_box_faces 3,3,3,3 -dm_view -sol_view \
5607b197a28SMatthew G. Knepley           -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces -dm_plex_check_pointsf
5617b197a28SMatthew G. Knepley 
562fc229b03SMatthew G. Knepley   # Problem on rank 3 on MacOSX
563fc229b03SMatthew G. Knepley   test:
564fc229b03SMatthew G. Knepley     requires: fftw
565fc229b03SMatthew G. Knepley     suffix: dirac_free_field_par
566fc229b03SMatthew G. Knepley     nsize: 16
567fc229b03SMatthew G. Knepley     args: -dm_plex_dim 4 -dm_plex_shape hypercubic -dm_plex_box_faces 4,4,4,4 -dm_view -sol_view \
568fc229b03SMatthew G. Knepley           -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces -dm_plex_check_pointsf
569fc229b03SMatthew G. Knepley 
5707b197a28SMatthew G. Knepley TEST*/
571