xref: /petsc/src/snes/tutorials/ex7.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 {
147b197a28SMatthew G. Knepley   PetscBool usePV; /* Use Pauli-Villars preconditioning */
157b197a28SMatthew G. Knepley } AppCtx;
167b197a28SMatthew G. Knepley 
177b197a28SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
187b197a28SMatthew G. Knepley {
197b197a28SMatthew G. Knepley   PetscFunctionBegin;
207b197a28SMatthew G. Knepley   options->usePV = PETSC_TRUE;
217b197a28SMatthew G. Knepley 
227b197a28SMatthew G. Knepley   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
237b197a28SMatthew G. Knepley   PetscCall(PetscOptionsBool("-use_pv", "Use Pauli-Villars preconditioning", "ex1.c", options->usePV, &options->usePV, NULL));
247b197a28SMatthew G. Knepley   PetscOptionsEnd();
25*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
267b197a28SMatthew G. Knepley }
277b197a28SMatthew G. Knepley 
287b197a28SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
297b197a28SMatthew G. Knepley {
307b197a28SMatthew G. Knepley   PetscFunctionBegin;
317b197a28SMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
327b197a28SMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
337b197a28SMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
347b197a28SMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
35*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
367b197a28SMatthew G. Knepley }
377b197a28SMatthew G. Knepley 
387b197a28SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
397b197a28SMatthew G. Knepley {
407b197a28SMatthew G. Knepley   PetscSection s;
417b197a28SMatthew G. Knepley   PetscInt     vStart, vEnd, v;
427b197a28SMatthew G. Knepley 
437b197a28SMatthew G. Knepley   PetscFunctionBegin;
447b197a28SMatthew G. Knepley   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
457b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
467b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
477b197a28SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
487b197a28SMatthew G. Knepley     PetscCall(PetscSectionSetDof(s, v, 12));
497b197a28SMatthew G. Knepley     /* TODO Divide the values into fields/components */
507b197a28SMatthew G. Knepley   }
517b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetUp(s));
527b197a28SMatthew G. Knepley   PetscCall(DMSetLocalSection(dm, s));
537b197a28SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&s));
54*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557b197a28SMatthew G. Knepley }
567b197a28SMatthew G. Knepley 
577b197a28SMatthew G. Knepley static PetscErrorCode SetupAuxDiscretization(DM dm, AppCtx *user)
587b197a28SMatthew G. Knepley {
597b197a28SMatthew G. Knepley   DM           dmAux, coordDM;
607b197a28SMatthew G. Knepley   PetscSection s;
617b197a28SMatthew G. Knepley   Vec          gauge;
627b197a28SMatthew G. Knepley   PetscInt     eStart, eEnd, e;
637b197a28SMatthew G. Knepley 
647b197a28SMatthew G. Knepley   PetscFunctionBegin;
657b197a28SMatthew G. Knepley   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
667b197a28SMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &coordDM));
677b197a28SMatthew G. Knepley   PetscCall(DMClone(dm, &dmAux));
687b197a28SMatthew G. Knepley   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
697b197a28SMatthew G. Knepley   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
707b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
717b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetChart(s, eStart, eEnd));
727b197a28SMatthew G. Knepley   for (e = eStart; e < eEnd; ++e) {
737b197a28SMatthew G. Knepley     /* TODO Should we store the whole SU(3) matrix, or the symmetric part? */
747b197a28SMatthew G. Knepley     PetscCall(PetscSectionSetDof(s, e, 9));
757b197a28SMatthew G. Knepley   }
767b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetUp(s));
777b197a28SMatthew G. Knepley   PetscCall(DMSetLocalSection(dmAux, s));
787b197a28SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&s));
797b197a28SMatthew G. Knepley   PetscCall(DMCreateLocalVector(dmAux, &gauge));
807b197a28SMatthew G. Knepley   PetscCall(DMDestroy(&dmAux));
817b197a28SMatthew G. Knepley   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, gauge));
827b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&gauge));
83*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
847b197a28SMatthew G. Knepley }
857b197a28SMatthew G. Knepley 
867b197a28SMatthew G. Knepley static PetscErrorCode PrintVertex(DM dm, PetscInt v)
877b197a28SMatthew G. Knepley {
887b197a28SMatthew G. Knepley   MPI_Comm       comm;
897b197a28SMatthew G. Knepley   PetscContainer c;
907b197a28SMatthew G. Knepley   PetscInt      *extent;
917b197a28SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, sum;
927b197a28SMatthew G. Knepley 
937b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
947b197a28SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
957b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
967b197a28SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
977b197a28SMatthew G. Knepley   PetscCall(PetscObjectQuery((PetscObject)dm, "_extent", (PetscObject *)&c));
987b197a28SMatthew G. Knepley   PetscCall(PetscContainerGetPointer(c, (void **)&extent));
997b197a28SMatthew G. Knepley   sum = 1;
1007b197a28SMatthew G. Knepley   PetscCall(PetscPrintf(comm, "Vertex %" PetscInt_FMT ":", v));
1017b197a28SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
1027b197a28SMatthew G. Knepley     PetscCall(PetscPrintf(comm, " %" PetscInt_FMT, (v / sum) % extent[d]));
1037b197a28SMatthew G. Knepley     if (d < dim) sum *= extent[d];
1047b197a28SMatthew G. Knepley   }
105*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1067b197a28SMatthew G. Knepley }
1077b197a28SMatthew G. Knepley 
1087b197a28SMatthew G. Knepley // Apply \gamma_\mu
1097b197a28SMatthew G. Knepley static PetscErrorCode ComputeGamma(PetscInt d, PetscInt ldx, PetscScalar f[])
1107b197a28SMatthew G. Knepley {
1117b197a28SMatthew G. Knepley   const PetscScalar fin[4] = {f[0 * ldx], f[1 * ldx], f[2 * ldx], f[3 * ldx]};
1127b197a28SMatthew G. Knepley 
1137b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
1147b197a28SMatthew G. Knepley   switch (d) {
1157b197a28SMatthew G. Knepley   case 0:
1167b197a28SMatthew G. Knepley     f[0 * ldx] = PETSC_i * fin[3];
1177b197a28SMatthew G. Knepley     f[1 * ldx] = PETSC_i * fin[2];
1187b197a28SMatthew G. Knepley     f[2 * ldx] = -PETSC_i * fin[1];
1197b197a28SMatthew G. Knepley     f[3 * ldx] = -PETSC_i * fin[0];
1207b197a28SMatthew G. Knepley     break;
1217b197a28SMatthew G. Knepley   case 1:
1227b197a28SMatthew G. Knepley     f[0 * ldx] = -fin[3];
1237b197a28SMatthew G. Knepley     f[1 * ldx] = fin[2];
1247b197a28SMatthew G. Knepley     f[2 * ldx] = fin[1];
1257b197a28SMatthew G. Knepley     f[3 * ldx] = -fin[0];
1267b197a28SMatthew G. Knepley     break;
1277b197a28SMatthew G. Knepley   case 2:
1287b197a28SMatthew G. Knepley     f[0 * ldx] = PETSC_i * fin[2];
1297b197a28SMatthew G. Knepley     f[1 * ldx] = -PETSC_i * fin[3];
1307b197a28SMatthew G. Knepley     f[2 * ldx] = -PETSC_i * fin[0];
1317b197a28SMatthew G. Knepley     f[3 * ldx] = PETSC_i * fin[1];
1327b197a28SMatthew G. Knepley     break;
1337b197a28SMatthew G. Knepley   case 3:
1347b197a28SMatthew G. Knepley     f[0 * ldx] = fin[2];
1357b197a28SMatthew G. Knepley     f[1 * ldx] = fin[3];
1367b197a28SMatthew G. Knepley     f[2 * ldx] = fin[0];
1377b197a28SMatthew G. Knepley     f[3 * ldx] = fin[1];
1387b197a28SMatthew G. Knepley     break;
1397b197a28SMatthew G. Knepley   default:
1407b197a28SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Direction for gamma %" PetscInt_FMT " not in [0, 4)", d);
1417b197a28SMatthew G. Knepley   }
142*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1437b197a28SMatthew G. Knepley }
1447b197a28SMatthew G. Knepley 
1457b197a28SMatthew G. Knepley // Apply (1 \pm \gamma_\mu)/2
1467b197a28SMatthew G. Knepley static PetscErrorCode ComputeGammaFactor(PetscInt d, PetscBool forward, PetscInt ldx, PetscScalar f[])
1477b197a28SMatthew G. Knepley {
1487b197a28SMatthew G. Knepley   const PetscReal   sign   = forward ? -1. : 1.;
1497b197a28SMatthew G. Knepley   const PetscScalar fin[4] = {f[0 * ldx], f[1 * ldx], f[2 * ldx], f[3 * ldx]};
1507b197a28SMatthew G. Knepley 
1517b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
1527b197a28SMatthew G. Knepley   switch (d) {
1537b197a28SMatthew G. Knepley   case 0:
1547b197a28SMatthew G. Knepley     f[0 * ldx] += sign * PETSC_i * fin[3];
1557b197a28SMatthew G. Knepley     f[1 * ldx] += sign * PETSC_i * fin[2];
1567b197a28SMatthew G. Knepley     f[2 * ldx] += sign * -PETSC_i * fin[1];
1577b197a28SMatthew G. Knepley     f[3 * ldx] += sign * -PETSC_i * fin[0];
1587b197a28SMatthew G. Knepley     break;
1597b197a28SMatthew G. Knepley   case 1:
1607b197a28SMatthew G. Knepley     f[0 * ldx] += -sign * fin[3];
1617b197a28SMatthew G. Knepley     f[1 * ldx] += sign * fin[2];
1627b197a28SMatthew G. Knepley     f[2 * ldx] += sign * fin[1];
1637b197a28SMatthew G. Knepley     f[3 * ldx] += -sign * fin[0];
1647b197a28SMatthew G. Knepley     break;
1657b197a28SMatthew G. Knepley   case 2:
1667b197a28SMatthew G. Knepley     f[0 * ldx] += sign * PETSC_i * fin[2];
1677b197a28SMatthew G. Knepley     f[1 * ldx] += sign * -PETSC_i * fin[3];
1687b197a28SMatthew G. Knepley     f[2 * ldx] += sign * -PETSC_i * fin[0];
1697b197a28SMatthew G. Knepley     f[3 * ldx] += sign * PETSC_i * fin[1];
1707b197a28SMatthew G. Knepley     break;
1717b197a28SMatthew G. Knepley   case 3:
1727b197a28SMatthew G. Knepley     f[0 * ldx] += sign * fin[2];
1737b197a28SMatthew G. Knepley     f[1 * ldx] += sign * fin[3];
1747b197a28SMatthew G. Knepley     f[2 * ldx] += sign * fin[0];
1757b197a28SMatthew G. Knepley     f[3 * ldx] += sign * fin[1];
1767b197a28SMatthew G. Knepley     break;
1777b197a28SMatthew G. Knepley   default:
1787b197a28SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Direction for gamma %" PetscInt_FMT " not in [0, 4)", d);
1797b197a28SMatthew G. Knepley   }
1807b197a28SMatthew G. Knepley   f[0 * ldx] *= 0.5;
1817b197a28SMatthew G. Knepley   f[1 * ldx] *= 0.5;
1827b197a28SMatthew G. Knepley   f[2 * ldx] *= 0.5;
1837b197a28SMatthew G. Knepley   f[3 * ldx] *= 0.5;
184*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1857b197a28SMatthew G. Knepley }
1867b197a28SMatthew G. Knepley 
1877b197a28SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>
1887b197a28SMatthew G. Knepley 
1897b197a28SMatthew G. Knepley // ComputeAction() sums the action of 1/2 (1 \pm \gamma_\mu) U \psi into f
1907b197a28SMatthew G. Knepley static PetscErrorCode ComputeAction(PetscInt d, PetscBool forward, const PetscScalar U[], const PetscScalar psi[], PetscScalar f[])
1917b197a28SMatthew G. Knepley {
1927b197a28SMatthew G. Knepley   PetscScalar tmp[12];
1937b197a28SMatthew G. Knepley 
1947b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
1957b197a28SMatthew G. Knepley   // Apply U
1967b197a28SMatthew G. Knepley   for (PetscInt beta = 0; beta < 4; ++beta) {
1977b197a28SMatthew G. Knepley     if (forward) DMPlex_Mult3D_Internal(U, 1, &psi[beta * 3], &tmp[beta * 3]);
1987b197a28SMatthew G. Knepley     else DMPlex_MultTranspose3D_Internal(U, 1, &psi[beta * 3], &tmp[beta * 3]);
1997b197a28SMatthew G. Knepley   }
2007b197a28SMatthew G. Knepley   // Apply (1 \pm \gamma_\mu)/2 to each color
2017b197a28SMatthew G. Knepley   for (PetscInt c = 0; c < 3; ++c) PetscCall(ComputeGammaFactor(d, forward, 3, &tmp[c]));
2027b197a28SMatthew G. Knepley   // Note that we are subtracting this contribution
2037b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < 12; ++i) f[i] -= tmp[i];
204*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2057b197a28SMatthew G. Knepley }
2067b197a28SMatthew G. Knepley 
2077b197a28SMatthew G. Knepley /*
2087b197a28SMatthew 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.
2097b197a28SMatthew G. Knepley */
2107b197a28SMatthew G. Knepley static PetscErrorCode ComputeResidual(DM dm, Vec u, Vec f)
2117b197a28SMatthew G. Knepley {
212de19e7feSJoseph Pusztay   DM                 dmAux;
213de19e7feSJoseph Pusztay   Vec                gauge;
214de19e7feSJoseph Pusztay   PetscSection       s, sGauge;
2157b197a28SMatthew G. Knepley   const PetscScalar *ua;
216de19e7feSJoseph Pusztay   PetscScalar       *fa, *link;
2177b197a28SMatthew G. Knepley   PetscInt           dim, vStart, vEnd;
2187b197a28SMatthew G. Knepley 
2197b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
2207b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
2217b197a28SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
2227b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2237b197a28SMatthew G. Knepley   PetscCall(VecGetArrayRead(u, &ua));
2247b197a28SMatthew G. Knepley   PetscCall(VecGetArray(f, &fa));
225de19e7feSJoseph Pusztay 
226de19e7feSJoseph Pusztay   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &gauge));
227de19e7feSJoseph Pusztay   PetscCall(VecGetDM(gauge, &dmAux));
228de19e7feSJoseph Pusztay   PetscCall(DMGetLocalSection(dmAux, &sGauge));
229de19e7feSJoseph Pusztay   PetscCall(VecGetArray(gauge, &link));
2307b197a28SMatthew G. Knepley   // Loop over y
2317b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
2327b197a28SMatthew G. Knepley     const PetscInt *supp;
2337b197a28SMatthew G. Knepley     PetscInt        xdof, xoff;
2347b197a28SMatthew G. Knepley 
2357b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, v, &supp));
2367b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetDof(s, v, &xdof));
2377b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetOffset(s, v, &xoff));
2387b197a28SMatthew G. Knepley     // Diagonal
2397b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < xdof; ++i) fa[xoff + i] += (M + 4) * ua[xoff + i];
2407b197a28SMatthew G. Knepley     // Loop over mu
2417b197a28SMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) {
2427b197a28SMatthew G. Knepley       const PetscInt *cone;
243de19e7feSJoseph Pusztay       PetscInt        yoff, goff;
2447b197a28SMatthew G. Knepley 
2457b197a28SMatthew G. Knepley       // Left action -(1 + \gamma_\mu)/2 \otimes U^\dagger_\mu(y) \delta_{x - \mu,y} \psi(y)
2467b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[2 * d + 0], &cone));
2477b197a28SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(s, cone[0], &yoff));
248de19e7feSJoseph Pusztay       PetscCall(PetscSectionGetOffset(sGauge, supp[2 * d + 0], &goff));
249de19e7feSJoseph Pusztay       PetscCall(ComputeAction(d, PETSC_FALSE, &link[goff], &ua[yoff], &fa[xoff]));
2507b197a28SMatthew G. Knepley       // Right edge -(1 - \gamma_\mu)/2 \otimes U_\mu(x) \delta_{x + \mu,y} \psi(y)
2517b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[2 * d + 1], &cone));
2527b197a28SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(s, cone[1], &yoff));
253de19e7feSJoseph Pusztay       PetscCall(PetscSectionGetOffset(sGauge, supp[2 * d + 1], &goff));
254de19e7feSJoseph Pusztay       PetscCall(ComputeAction(d, PETSC_TRUE, &link[goff], &ua[yoff], &fa[xoff]));
2557b197a28SMatthew G. Knepley     }
2567b197a28SMatthew G. Knepley   }
2577b197a28SMatthew G. Knepley   PetscCall(VecRestoreArray(f, &fa));
258de19e7feSJoseph Pusztay   PetscCall(VecRestoreArray(gauge, &link));
2597b197a28SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(u, &ua));
260*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2617b197a28SMatthew G. Knepley }
2627b197a28SMatthew G. Knepley 
2637b197a28SMatthew G. Knepley static PetscErrorCode CreateJacobian(DM dm, Mat *J)
2647b197a28SMatthew G. Knepley {
2657b197a28SMatthew G. Knepley   PetscFunctionBegin;
266*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2677b197a28SMatthew G. Knepley }
2687b197a28SMatthew G. Knepley 
2697b197a28SMatthew G. Knepley static PetscErrorCode ComputeJacobian(DM dm, Vec u, Mat J)
2707b197a28SMatthew G. Knepley {
2717b197a28SMatthew G. Knepley   PetscFunctionBegin;
272*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2737b197a28SMatthew G. Knepley }
2747b197a28SMatthew G. Knepley 
2757b197a28SMatthew G. Knepley static PetscErrorCode PrintTraversal(DM dm)
2767b197a28SMatthew G. Knepley {
2777b197a28SMatthew G. Knepley   MPI_Comm comm;
2787b197a28SMatthew G. Knepley   PetscInt vStart, vEnd;
2797b197a28SMatthew G. Knepley 
2807b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
2817b197a28SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2827b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
2837b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
2847b197a28SMatthew G. Knepley     const PetscInt *supp;
2857b197a28SMatthew G. Knepley     PetscInt        Ns;
2867b197a28SMatthew G. Knepley 
2877b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupportSize(dm, v, &Ns));
2887b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, v, &supp));
2897b197a28SMatthew G. Knepley     PetscCall(PrintVertex(dm, v));
2907b197a28SMatthew G. Knepley     PetscCall(PetscPrintf(comm, "\n"));
2917b197a28SMatthew G. Knepley     for (PetscInt s = 0; s < Ns; ++s) {
2927b197a28SMatthew G. Knepley       const PetscInt *cone;
2937b197a28SMatthew G. Knepley 
2947b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[s], &cone));
2957b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "  Edge %" PetscInt_FMT ": ", supp[s]));
2967b197a28SMatthew G. Knepley       PetscCall(PrintVertex(dm, cone[0]));
2977b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, " -- "));
2987b197a28SMatthew G. Knepley       PetscCall(PrintVertex(dm, cone[1]));
2997b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "\n"));
3007b197a28SMatthew G. Knepley     }
3017b197a28SMatthew G. Knepley   }
302*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3037b197a28SMatthew G. Knepley }
3047b197a28SMatthew G. Knepley 
3057b197a28SMatthew G. Knepley static PetscErrorCode ComputeFFT(Mat FT, PetscInt Nc, Vec x, Vec p)
3067b197a28SMatthew G. Knepley {
3077b197a28SMatthew G. Knepley   Vec     *xComp, *pComp;
3087b197a28SMatthew G. Knepley   PetscInt n, N;
3097b197a28SMatthew G. Knepley 
3107b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
3117b197a28SMatthew G. Knepley   PetscCall(PetscMalloc2(Nc, &xComp, Nc, &pComp));
3127b197a28SMatthew G. Knepley   PetscCall(VecGetLocalSize(x, &n));
3137b197a28SMatthew G. Knepley   PetscCall(VecGetSize(x, &N));
3147b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
3157b197a28SMatthew G. Knepley     const char *vtype;
3167b197a28SMatthew G. Knepley 
3177b197a28SMatthew G. Knepley     // HACK: Make these from another DM up front
3187b197a28SMatthew G. Knepley     PetscCall(VecCreate(PetscObjectComm((PetscObject)x), &xComp[i]));
3197b197a28SMatthew G. Knepley     PetscCall(VecGetType(x, &vtype));
3207b197a28SMatthew G. Knepley     PetscCall(VecSetType(xComp[i], vtype));
3217b197a28SMatthew G. Knepley     PetscCall(VecSetSizes(xComp[i], n / Nc, N / Nc));
3227b197a28SMatthew G. Knepley     PetscCall(VecDuplicate(xComp[i], &pComp[i]));
3237b197a28SMatthew G. Knepley   }
3247b197a28SMatthew G. Knepley   PetscCall(VecStrideGatherAll(x, xComp, INSERT_VALUES));
3257b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) PetscCall(MatMult(FT, xComp[i], pComp[i]));
3267b197a28SMatthew G. Knepley   PetscCall(VecStrideScatterAll(pComp, p, INSERT_VALUES));
3277b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
3287b197a28SMatthew G. Knepley     PetscCall(VecDestroy(&xComp[i]));
3297b197a28SMatthew G. Knepley     PetscCall(VecDestroy(&pComp[i]));
3307b197a28SMatthew G. Knepley   }
3317b197a28SMatthew G. Knepley   PetscCall(PetscFree2(xComp, pComp));
332*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3337b197a28SMatthew G. Knepley }
3347b197a28SMatthew G. Knepley 
335de19e7feSJoseph Pusztay // Sets each link to be the identity for the free field test
336de19e7feSJoseph Pusztay static PetscErrorCode SetGauge_Identity(DM dm)
337de19e7feSJoseph Pusztay {
338de19e7feSJoseph Pusztay   DM           auxDM;
339de19e7feSJoseph Pusztay   Vec          auxVec;
340de19e7feSJoseph Pusztay   PetscSection s;
341de19e7feSJoseph Pusztay   PetscScalar  id[9] = {1., 0., 0., 0., 1., 0., 0., 0., 1.};
342de19e7feSJoseph Pusztay   PetscInt     eStart, eEnd;
343de19e7feSJoseph Pusztay 
344de19e7feSJoseph Pusztay   PetscFunctionBegin;
345de19e7feSJoseph Pusztay   PetscCall(DMGetAuxiliaryVec(dm, NULL, 0, 0, &auxVec));
346de19e7feSJoseph Pusztay   PetscCall(VecGetDM(auxVec, &auxDM));
347de19e7feSJoseph Pusztay   PetscCall(DMGetLocalSection(auxDM, &s));
348de19e7feSJoseph Pusztay   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
349de19e7feSJoseph Pusztay   for (PetscInt i = eStart; i < eEnd; ++i) { PetscCall(VecSetValuesSection(auxVec, s, i, id, INSERT_VALUES)); }
350de19e7feSJoseph Pusztay   PetscCall(VecViewFromOptions(auxVec, NULL, "-gauge_view"));
351*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
352de19e7feSJoseph Pusztay }
353de19e7feSJoseph Pusztay 
3547b197a28SMatthew G. Knepley /*
3557b197a28SMatthew G. Knepley   Test the action of the Wilson operator in the free field case U = I,
3567b197a28SMatthew G. Knepley 
3577b197a28SMatthew G. Knepley     \eta(x) = D_W(x - y) \psi(y)
3587b197a28SMatthew G. Knepley 
3597b197a28SMatthew G. Knepley   The Wilson operator is a convolution for the free field, so we can check that by the convolution theorem
3607b197a28SMatthew G. Knepley 
3617b197a28SMatthew G. Knepley     \hat\eta(x) = \mathcal{F}(D_W(x - y) \psi(y))
3627b197a28SMatthew G. Knepley                 = \hat D_W(p) \mathcal{F}\psi(p)
3637b197a28SMatthew G. Knepley 
3647b197a28SMatthew G. Knepley   The Fourier series for the Wilson operator is
3657b197a28SMatthew G. Knepley 
3667b197a28SMatthew G. Knepley     M + \sum_\mu 2 \sin^2(p_\mu / 2) + i \gamma_\mu \sin(p_\mu)
3677b197a28SMatthew G. Knepley */
3687b197a28SMatthew G. Knepley static PetscErrorCode TestFreeField(DM dm)
3697b197a28SMatthew G. Knepley {
3707b197a28SMatthew G. Knepley   PetscSection       s;
3717b197a28SMatthew G. Knepley   Mat                FT;
3727b197a28SMatthew G. Knepley   Vec                psi, psiHat;
3737b197a28SMatthew G. Knepley   Vec                eta, etaHat;
3747b197a28SMatthew G. Knepley   Vec                DHat; // The product \hat D_w \hat psi
3757b197a28SMatthew G. Knepley   PetscRandom        r;
3767b197a28SMatthew G. Knepley   const PetscScalar *psih;
3777b197a28SMatthew G. Knepley   PetscScalar       *dh;
3787b197a28SMatthew G. Knepley   PetscReal         *coef, nrm;
3797b197a28SMatthew G. Knepley   const PetscInt    *extent, Nc = 12;
3807b197a28SMatthew G. Knepley   PetscInt           dim, V     = 1, vStart, vEnd;
3817b197a28SMatthew G. Knepley   PetscContainer     c;
3827b197a28SMatthew G. Knepley   PetscBool          constRhs = PETSC_FALSE;
3837b197a28SMatthew G. Knepley 
3847b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
3857b197a28SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-const_rhs", &constRhs, NULL));
3867b197a28SMatthew G. Knepley 
387de19e7feSJoseph Pusztay   PetscCall(SetGauge_Identity(dm));
3887b197a28SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
3897b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &psi));
3907b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)psi, "psi"));
3917b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &psiHat));
3927b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)psiHat, "psihat"));
3937b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &eta));
3947b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)eta, "eta"));
3957b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &etaHat));
3967b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)etaHat, "etahat"));
3977b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &DHat));
3987b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)DHat, "Dhat"));
3997b197a28SMatthew G. Knepley   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
4007b197a28SMatthew G. Knepley   PetscCall(PetscRandomSetType(r, PETSCRAND48));
4017b197a28SMatthew G. Knepley   if (constRhs) PetscCall(VecSet(psi, 1.));
4027b197a28SMatthew G. Knepley   else PetscCall(VecSetRandom(psi, r));
4037b197a28SMatthew G. Knepley   PetscCall(PetscRandomDestroy(&r));
4047b197a28SMatthew G. Knepley 
4057b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
4067b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
4077b197a28SMatthew G. Knepley   PetscCall(PetscObjectQuery((PetscObject)dm, "_extent", (PetscObject *)&c));
4087b197a28SMatthew G. Knepley   PetscCall(PetscContainerGetPointer(c, (void **)&extent));
4097b197a28SMatthew G. Knepley   PetscCall(MatCreateFFT(PetscObjectComm((PetscObject)dm), dim, extent, MATFFTW, &FT));
4107b197a28SMatthew G. Knepley 
4117b197a28SMatthew G. Knepley   PetscCall(PetscMalloc1(dim, &coef));
4127b197a28SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
4137b197a28SMatthew G. Knepley     coef[d] = 2. * PETSC_PI / extent[d];
4147b197a28SMatthew G. Knepley     V *= extent[d];
4157b197a28SMatthew G. Knepley   }
4167b197a28SMatthew G. Knepley   PetscCall(ComputeResidual(dm, psi, eta));
4177b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(eta, NULL, "-psi_view"));
4187b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(eta, NULL, "-eta_view"));
4197b197a28SMatthew G. Knepley   PetscCall(ComputeFFT(FT, Nc, psi, psiHat));
4207b197a28SMatthew G. Knepley   PetscCall(VecScale(psiHat, 1. / V));
4217b197a28SMatthew G. Knepley   PetscCall(ComputeFFT(FT, Nc, eta, etaHat));
4227b197a28SMatthew G. Knepley   PetscCall(VecScale(etaHat, 1. / V));
4237b197a28SMatthew G. Knepley   PetscCall(VecGetArrayRead(psiHat, &psih));
4247b197a28SMatthew G. Knepley   PetscCall(VecGetArray(DHat, &dh));
4257b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
4267b197a28SMatthew G. Knepley     PetscScalar tmp[12], tmp1 = 0.;
4277b197a28SMatthew G. Knepley     PetscInt    dof, off;
4287b197a28SMatthew G. Knepley 
4297b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetDof(s, v, &dof));
4307b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetOffset(s, v, &off));
4317b197a28SMatthew G. Knepley     for (PetscInt d = 0, prod = 1; d < dim; prod *= extent[d], ++d) {
4327b197a28SMatthew G. Knepley       const PetscInt idx = (v / prod) % extent[d];
4337b197a28SMatthew G. Knepley 
4347b197a28SMatthew G. Knepley       tmp1 += 2. * PetscSqr(PetscSinReal(0.5 * coef[d] * idx));
4357b197a28SMatthew G. Knepley       for (PetscInt i = 0; i < dof; ++i) tmp[i] = PETSC_i * PetscSinReal(coef[d] * idx) * psih[off + i];
436*3ba16761SJacob Faibussowitsch       for (PetscInt c = 0; c < 3; ++c) PetscCall(ComputeGamma(d, 3, &tmp[c]));
4377b197a28SMatthew G. Knepley       for (PetscInt i = 0; i < dof; ++i) dh[off + i] += tmp[i];
4387b197a28SMatthew G. Knepley     }
4397b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < dof; ++i) dh[off + i] += (M + tmp1) * psih[off + i];
4407b197a28SMatthew G. Knepley   }
4417b197a28SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(psiHat, &psih));
4427b197a28SMatthew G. Knepley   PetscCall(VecRestoreArray(DHat, &dh));
4437b197a28SMatthew G. Knepley 
4447b197a28SMatthew G. Knepley   {
4457b197a28SMatthew G. Knepley     Vec     *etaComp, *DComp;
4467b197a28SMatthew G. Knepley     PetscInt n, N;
4477b197a28SMatthew G. Knepley 
4487b197a28SMatthew G. Knepley     PetscCall(PetscMalloc2(Nc, &etaComp, Nc, &DComp));
4497b197a28SMatthew G. Knepley     PetscCall(VecGetLocalSize(etaHat, &n));
4507b197a28SMatthew G. Knepley     PetscCall(VecGetSize(etaHat, &N));
4517b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
4527b197a28SMatthew G. Knepley       const char *vtype;
4537b197a28SMatthew G. Knepley 
4547b197a28SMatthew G. Knepley       // HACK: Make these from another DM up front
4557b197a28SMatthew G. Knepley       PetscCall(VecCreate(PetscObjectComm((PetscObject)etaHat), &etaComp[i]));
4567b197a28SMatthew G. Knepley       PetscCall(VecGetType(etaHat, &vtype));
4577b197a28SMatthew G. Knepley       PetscCall(VecSetType(etaComp[i], vtype));
4587b197a28SMatthew G. Knepley       PetscCall(VecSetSizes(etaComp[i], n / Nc, N / Nc));
4597b197a28SMatthew G. Knepley       PetscCall(VecDuplicate(etaComp[i], &DComp[i]));
4607b197a28SMatthew G. Knepley     }
4617b197a28SMatthew G. Knepley     PetscCall(VecStrideGatherAll(etaHat, etaComp, INSERT_VALUES));
4627b197a28SMatthew G. Knepley     PetscCall(VecStrideGatherAll(DHat, DComp, INSERT_VALUES));
4637b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
4647b197a28SMatthew G. Knepley       if (!i) {
4657b197a28SMatthew G. Knepley         PetscCall(VecViewFromOptions(etaComp[i], NULL, "-etahat_view"));
4667b197a28SMatthew G. Knepley         PetscCall(VecViewFromOptions(DComp[i], NULL, "-dhat_view"));
4677b197a28SMatthew G. Knepley       }
4687b197a28SMatthew G. Knepley       PetscCall(VecAXPY(etaComp[i], -1., DComp[i]));
4697b197a28SMatthew G. Knepley       PetscCall(VecNorm(etaComp[i], NORM_INFINITY, &nrm));
4707b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Slice %" PetscInt_FMT ": %g\n", i, (double)nrm));
4717b197a28SMatthew G. Knepley     }
4727b197a28SMatthew G. Knepley     PetscCall(VecStrideScatterAll(etaComp, etaHat, INSERT_VALUES));
4737b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
4747b197a28SMatthew G. Knepley       PetscCall(VecDestroy(&etaComp[i]));
4757b197a28SMatthew G. Knepley       PetscCall(VecDestroy(&DComp[i]));
4767b197a28SMatthew G. Knepley     }
4777b197a28SMatthew G. Knepley     PetscCall(PetscFree2(etaComp, DComp));
4787b197a28SMatthew G. Knepley     PetscCall(VecNorm(etaHat, NORM_INFINITY, &nrm));
4797b197a28SMatthew G. Knepley     PetscCheck(nrm < PETSC_SMALL, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Free field test failed: %g", (double)nrm);
4807b197a28SMatthew G. Knepley   }
4817b197a28SMatthew G. Knepley 
4827b197a28SMatthew G. Knepley   PetscCall(PetscFree(coef));
4837b197a28SMatthew G. Knepley   PetscCall(MatDestroy(&FT));
4847b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &psi));
4857b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &psiHat));
4867b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &eta));
4877b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &etaHat));
4887b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &DHat));
489*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4907b197a28SMatthew G. Knepley }
4917b197a28SMatthew G. Knepley 
4927b197a28SMatthew G. Knepley int main(int argc, char **argv)
4937b197a28SMatthew G. Knepley {
4947b197a28SMatthew G. Knepley   DM     dm;
4957b197a28SMatthew G. Knepley   Vec    u, f;
4967b197a28SMatthew G. Knepley   Mat    J;
4977b197a28SMatthew G. Knepley   AppCtx user;
4987b197a28SMatthew G. Knepley 
4997b197a28SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
5007b197a28SMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
5017b197a28SMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
5027b197a28SMatthew G. Knepley   PetscCall(SetupDiscretization(dm, &user));
5037b197a28SMatthew G. Knepley   PetscCall(SetupAuxDiscretization(dm, &user));
5047b197a28SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &u));
5057b197a28SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &f));
5067b197a28SMatthew G. Knepley   PetscCall(PrintTraversal(dm));
5077b197a28SMatthew G. Knepley   PetscCall(ComputeResidual(dm, u, f));
5087b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(f, NULL, "-res_view"));
5097b197a28SMatthew G. Knepley   PetscCall(CreateJacobian(dm, &J));
5107b197a28SMatthew G. Knepley   PetscCall(ComputeJacobian(dm, u, J));
5117b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
5127b197a28SMatthew G. Knepley   PetscCall(TestFreeField(dm));
5137b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&f));
5147b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&u));
5157b197a28SMatthew G. Knepley   PetscCall(DMDestroy(&dm));
5167b197a28SMatthew G. Knepley   PetscCall(PetscFinalize());
5177b197a28SMatthew G. Knepley   return 0;
5187b197a28SMatthew G. Knepley }
5197b197a28SMatthew G. Knepley 
5207b197a28SMatthew G. Knepley /*TEST
5217b197a28SMatthew G. Knepley 
5227b197a28SMatthew G. Knepley   build:
5237b197a28SMatthew G. Knepley     requires: complex
5247b197a28SMatthew G. Knepley 
5257b197a28SMatthew G. Knepley   test:
5267b197a28SMatthew G. Knepley     requires: fftw
5277b197a28SMatthew G. Knepley     suffix: dirac_free_field
5287b197a28SMatthew G. Knepley     args: -dm_plex_dim 4 -dm_plex_shape hypercubic -dm_plex_box_faces 3,3,3,3 -dm_view -sol_view \
5297b197a28SMatthew G. Knepley           -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces -dm_plex_check_pointsf
5307b197a28SMatthew G. Knepley 
5317b197a28SMatthew G. Knepley TEST*/
532