xref: /petsc/src/snes/tutorials/ex7.c (revision 7b197a28a09f6917309a5286b58f91736188ce0c)
1*7b197a28SMatthew G. Knepley static char help[] = "Fermions on a hypercubic lattice.\n\n";
2*7b197a28SMatthew G. Knepley 
3*7b197a28SMatthew G. Knepley #include <petscdmplex.h>
4*7b197a28SMatthew G. Knepley #include <petscsnes.h>
5*7b197a28SMatthew G. Knepley 
6*7b197a28SMatthew G. Knepley /* Common operations:
7*7b197a28SMatthew G. Knepley 
8*7b197a28SMatthew G. Knepley  - View the input \psi as ASCII in lexicographic order: -psi_view
9*7b197a28SMatthew G. Knepley */
10*7b197a28SMatthew G. Knepley 
11*7b197a28SMatthew G. Knepley const PetscReal M = 1.0;
12*7b197a28SMatthew G. Knepley 
13*7b197a28SMatthew G. Knepley typedef struct {
14*7b197a28SMatthew G. Knepley   PetscBool usePV; /* Use Pauli-Villars preconditioning */
15*7b197a28SMatthew G. Knepley } AppCtx;
16*7b197a28SMatthew G. Knepley 
17*7b197a28SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
18*7b197a28SMatthew G. Knepley {
19*7b197a28SMatthew G. Knepley   PetscFunctionBegin;
20*7b197a28SMatthew G. Knepley   options->usePV = PETSC_TRUE;
21*7b197a28SMatthew G. Knepley 
22*7b197a28SMatthew G. Knepley   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
23*7b197a28SMatthew G. Knepley   PetscCall(PetscOptionsBool("-use_pv", "Use Pauli-Villars preconditioning", "ex1.c", options->usePV, &options->usePV, NULL));
24*7b197a28SMatthew G. Knepley   PetscOptionsEnd();
25*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
26*7b197a28SMatthew G. Knepley }
27*7b197a28SMatthew G. Knepley 
28*7b197a28SMatthew G. Knepley PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
29*7b197a28SMatthew G. Knepley {
30*7b197a28SMatthew G. Knepley   PetscFunctionBegin;
31*7b197a28SMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
32*7b197a28SMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
33*7b197a28SMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
34*7b197a28SMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
35*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
36*7b197a28SMatthew G. Knepley }
37*7b197a28SMatthew G. Knepley 
38*7b197a28SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
39*7b197a28SMatthew G. Knepley {
40*7b197a28SMatthew G. Knepley   PetscSection s;
41*7b197a28SMatthew G. Knepley   PetscInt     vStart, vEnd, v;
42*7b197a28SMatthew G. Knepley 
43*7b197a28SMatthew G. Knepley   PetscFunctionBegin;
44*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
45*7b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
46*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetChart(s, vStart, vEnd));
47*7b197a28SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
48*7b197a28SMatthew G. Knepley     PetscCall(PetscSectionSetDof(s, v, 12));
49*7b197a28SMatthew G. Knepley     /* TODO Divide the values into fields/components */
50*7b197a28SMatthew G. Knepley   }
51*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetUp(s));
52*7b197a28SMatthew G. Knepley   PetscCall(DMSetLocalSection(dm, s));
53*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&s));
54*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
55*7b197a28SMatthew G. Knepley }
56*7b197a28SMatthew G. Knepley 
57*7b197a28SMatthew G. Knepley static PetscErrorCode SetupAuxDiscretization(DM dm, AppCtx *user)
58*7b197a28SMatthew G. Knepley {
59*7b197a28SMatthew G. Knepley   DM           dmAux, coordDM;
60*7b197a28SMatthew G. Knepley   PetscSection s;
61*7b197a28SMatthew G. Knepley   Vec          gauge;
62*7b197a28SMatthew G. Knepley   PetscInt     eStart, eEnd, e;
63*7b197a28SMatthew G. Knepley 
64*7b197a28SMatthew G. Knepley   PetscFunctionBegin;
65*7b197a28SMatthew G. Knepley   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
66*7b197a28SMatthew G. Knepley   PetscCall(DMGetCoordinateDM(dm, &coordDM));
67*7b197a28SMatthew G. Knepley   PetscCall(DMClone(dm, &dmAux));
68*7b197a28SMatthew G. Knepley   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
69*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s));
70*7b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
71*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetChart(s, eStart, eEnd));
72*7b197a28SMatthew G. Knepley   for (e = eStart; e < eEnd; ++e) {
73*7b197a28SMatthew G. Knepley     /* TODO Should we store the whole SU(3) matrix, or the symmetric part? */
74*7b197a28SMatthew G. Knepley     PetscCall(PetscSectionSetDof(s, e, 9));
75*7b197a28SMatthew G. Knepley   }
76*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionSetUp(s));
77*7b197a28SMatthew G. Knepley   PetscCall(DMSetLocalSection(dmAux, s));
78*7b197a28SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&s));
79*7b197a28SMatthew G. Knepley   PetscCall(DMCreateLocalVector(dmAux, &gauge));
80*7b197a28SMatthew G. Knepley   PetscCall(DMDestroy(&dmAux));
81*7b197a28SMatthew G. Knepley   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, gauge));
82*7b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&gauge));
83*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
84*7b197a28SMatthew G. Knepley }
85*7b197a28SMatthew G. Knepley 
86*7b197a28SMatthew G. Knepley static PetscErrorCode PrintVertex(DM dm, PetscInt v)
87*7b197a28SMatthew G. Knepley {
88*7b197a28SMatthew G. Knepley   MPI_Comm       comm;
89*7b197a28SMatthew G. Knepley   PetscContainer c;
90*7b197a28SMatthew G. Knepley   PetscInt      *extent;
91*7b197a28SMatthew G. Knepley   PetscInt       dim, cStart, cEnd, sum;
92*7b197a28SMatthew G. Knepley 
93*7b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
94*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
95*7b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
96*7b197a28SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
97*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectQuery((PetscObject)dm, "_extent", (PetscObject *)&c));
98*7b197a28SMatthew G. Knepley   PetscCall(PetscContainerGetPointer(c, (void **)&extent));
99*7b197a28SMatthew G. Knepley   sum = 1;
100*7b197a28SMatthew G. Knepley   PetscCall(PetscPrintf(comm, "Vertex %" PetscInt_FMT ":", v));
101*7b197a28SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
102*7b197a28SMatthew G. Knepley     PetscCall(PetscPrintf(comm, " %" PetscInt_FMT, (v / sum) % extent[d]));
103*7b197a28SMatthew G. Knepley     if (d < dim) sum *= extent[d];
104*7b197a28SMatthew G. Knepley   }
105*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
106*7b197a28SMatthew G. Knepley }
107*7b197a28SMatthew G. Knepley 
108*7b197a28SMatthew G. Knepley // Apply \gamma_\mu
109*7b197a28SMatthew G. Knepley static PetscErrorCode ComputeGamma(PetscInt d, PetscInt ldx, PetscScalar f[])
110*7b197a28SMatthew G. Knepley {
111*7b197a28SMatthew G. Knepley   const PetscScalar fin[4] = {f[0 * ldx], f[1 * ldx], f[2 * ldx], f[3 * ldx]};
112*7b197a28SMatthew G. Knepley 
113*7b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
114*7b197a28SMatthew G. Knepley   switch (d) {
115*7b197a28SMatthew G. Knepley   case 0:
116*7b197a28SMatthew G. Knepley     f[0 * ldx] = PETSC_i * fin[3];
117*7b197a28SMatthew G. Knepley     f[1 * ldx] = PETSC_i * fin[2];
118*7b197a28SMatthew G. Knepley     f[2 * ldx] = -PETSC_i * fin[1];
119*7b197a28SMatthew G. Knepley     f[3 * ldx] = -PETSC_i * fin[0];
120*7b197a28SMatthew G. Knepley     break;
121*7b197a28SMatthew G. Knepley   case 1:
122*7b197a28SMatthew G. Knepley     f[0 * ldx] = -fin[3];
123*7b197a28SMatthew G. Knepley     f[1 * ldx] = fin[2];
124*7b197a28SMatthew G. Knepley     f[2 * ldx] = fin[1];
125*7b197a28SMatthew G. Knepley     f[3 * ldx] = -fin[0];
126*7b197a28SMatthew G. Knepley     break;
127*7b197a28SMatthew G. Knepley   case 2:
128*7b197a28SMatthew G. Knepley     f[0 * ldx] = PETSC_i * fin[2];
129*7b197a28SMatthew G. Knepley     f[1 * ldx] = -PETSC_i * fin[3];
130*7b197a28SMatthew G. Knepley     f[2 * ldx] = -PETSC_i * fin[0];
131*7b197a28SMatthew G. Knepley     f[3 * ldx] = PETSC_i * fin[1];
132*7b197a28SMatthew G. Knepley     break;
133*7b197a28SMatthew G. Knepley   case 3:
134*7b197a28SMatthew G. Knepley     f[0 * ldx] = fin[2];
135*7b197a28SMatthew G. Knepley     f[1 * ldx] = fin[3];
136*7b197a28SMatthew G. Knepley     f[2 * ldx] = fin[0];
137*7b197a28SMatthew G. Knepley     f[3 * ldx] = fin[1];
138*7b197a28SMatthew G. Knepley     break;
139*7b197a28SMatthew G. Knepley   default:
140*7b197a28SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Direction for gamma %" PetscInt_FMT " not in [0, 4)", d);
141*7b197a28SMatthew G. Knepley   }
142*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
143*7b197a28SMatthew G. Knepley }
144*7b197a28SMatthew G. Knepley 
145*7b197a28SMatthew G. Knepley // Apply (1 \pm \gamma_\mu)/2
146*7b197a28SMatthew G. Knepley static PetscErrorCode ComputeGammaFactor(PetscInt d, PetscBool forward, PetscInt ldx, PetscScalar f[])
147*7b197a28SMatthew G. Knepley {
148*7b197a28SMatthew G. Knepley   const PetscReal   sign   = forward ? -1. : 1.;
149*7b197a28SMatthew G. Knepley   const PetscScalar fin[4] = {f[0 * ldx], f[1 * ldx], f[2 * ldx], f[3 * ldx]};
150*7b197a28SMatthew G. Knepley 
151*7b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
152*7b197a28SMatthew G. Knepley   switch (d) {
153*7b197a28SMatthew G. Knepley   case 0:
154*7b197a28SMatthew G. Knepley     f[0 * ldx] += sign * PETSC_i * fin[3];
155*7b197a28SMatthew G. Knepley     f[1 * ldx] += sign * PETSC_i * fin[2];
156*7b197a28SMatthew G. Knepley     f[2 * ldx] += sign * -PETSC_i * fin[1];
157*7b197a28SMatthew G. Knepley     f[3 * ldx] += sign * -PETSC_i * fin[0];
158*7b197a28SMatthew G. Knepley     break;
159*7b197a28SMatthew G. Knepley   case 1:
160*7b197a28SMatthew G. Knepley     f[0 * ldx] += -sign * fin[3];
161*7b197a28SMatthew G. Knepley     f[1 * ldx] += sign * fin[2];
162*7b197a28SMatthew G. Knepley     f[2 * ldx] += sign * fin[1];
163*7b197a28SMatthew G. Knepley     f[3 * ldx] += -sign * fin[0];
164*7b197a28SMatthew G. Knepley     break;
165*7b197a28SMatthew G. Knepley   case 2:
166*7b197a28SMatthew G. Knepley     f[0 * ldx] += sign * PETSC_i * fin[2];
167*7b197a28SMatthew G. Knepley     f[1 * ldx] += sign * -PETSC_i * fin[3];
168*7b197a28SMatthew G. Knepley     f[2 * ldx] += sign * -PETSC_i * fin[0];
169*7b197a28SMatthew G. Knepley     f[3 * ldx] += sign * PETSC_i * fin[1];
170*7b197a28SMatthew G. Knepley     break;
171*7b197a28SMatthew G. Knepley   case 3:
172*7b197a28SMatthew G. Knepley     f[0 * ldx] += sign * fin[2];
173*7b197a28SMatthew G. Knepley     f[1 * ldx] += sign * fin[3];
174*7b197a28SMatthew G. Knepley     f[2 * ldx] += sign * fin[0];
175*7b197a28SMatthew G. Knepley     f[3 * ldx] += sign * fin[1];
176*7b197a28SMatthew G. Knepley     break;
177*7b197a28SMatthew G. Knepley   default:
178*7b197a28SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Direction for gamma %" PetscInt_FMT " not in [0, 4)", d);
179*7b197a28SMatthew G. Knepley   }
180*7b197a28SMatthew G. Knepley   f[0 * ldx] *= 0.5;
181*7b197a28SMatthew G. Knepley   f[1 * ldx] *= 0.5;
182*7b197a28SMatthew G. Knepley   f[2 * ldx] *= 0.5;
183*7b197a28SMatthew G. Knepley   f[3 * ldx] *= 0.5;
184*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
185*7b197a28SMatthew G. Knepley }
186*7b197a28SMatthew G. Knepley 
187*7b197a28SMatthew G. Knepley #include <petsc/private/dmpleximpl.h>
188*7b197a28SMatthew G. Knepley 
189*7b197a28SMatthew G. Knepley // ComputeAction() sums the action of 1/2 (1 \pm \gamma_\mu) U \psi into f
190*7b197a28SMatthew G. Knepley static PetscErrorCode ComputeAction(PetscInt d, PetscBool forward, const PetscScalar U[], const PetscScalar psi[], PetscScalar f[])
191*7b197a28SMatthew G. Knepley {
192*7b197a28SMatthew G. Knepley   PetscScalar tmp[12];
193*7b197a28SMatthew G. Knepley 
194*7b197a28SMatthew G. Knepley   PetscFunctionBeginHot;
195*7b197a28SMatthew G. Knepley   // Apply U
196*7b197a28SMatthew G. Knepley   for (PetscInt beta = 0; beta < 4; ++beta) {
197*7b197a28SMatthew G. Knepley     if (forward) DMPlex_Mult3D_Internal(U, 1, &psi[beta * 3], &tmp[beta * 3]);
198*7b197a28SMatthew G. Knepley     else DMPlex_MultTranspose3D_Internal(U, 1, &psi[beta * 3], &tmp[beta * 3]);
199*7b197a28SMatthew G. Knepley   }
200*7b197a28SMatthew G. Knepley   // Apply (1 \pm \gamma_\mu)/2 to each color
201*7b197a28SMatthew G. Knepley   for (PetscInt c = 0; c < 3; ++c) PetscCall(ComputeGammaFactor(d, forward, 3, &tmp[c]));
202*7b197a28SMatthew G. Knepley   // Note that we are subtracting this contribution
203*7b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < 12; ++i) f[i] -= tmp[i];
204*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
205*7b197a28SMatthew G. Knepley }
206*7b197a28SMatthew G. Knepley 
207*7b197a28SMatthew G. Knepley /*
208*7b197a28SMatthew 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.
209*7b197a28SMatthew G. Knepley */
210*7b197a28SMatthew G. Knepley static PetscErrorCode ComputeResidual(DM dm, Vec u, Vec f)
211*7b197a28SMatthew G. Knepley {
212*7b197a28SMatthew G. Knepley   PetscSection       s;
213*7b197a28SMatthew G. Knepley   const PetscScalar *ua;
214*7b197a28SMatthew G. Knepley   PetscScalar       *fa;
215*7b197a28SMatthew G. Knepley   PetscScalar        id[9] = {1., 0., 0., 0., 1., 0., 0., 0., 1.};
216*7b197a28SMatthew G. Knepley   PetscInt           dim, vStart, vEnd;
217*7b197a28SMatthew G. Knepley 
218*7b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
219*7b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
220*7b197a28SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
221*7b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
222*7b197a28SMatthew G. Knepley   PetscCall(VecGetArrayRead(u, &ua));
223*7b197a28SMatthew G. Knepley   PetscCall(VecGetArray(f, &fa));
224*7b197a28SMatthew G. Knepley   // Loop over y
225*7b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
226*7b197a28SMatthew G. Knepley     const PetscInt *supp;
227*7b197a28SMatthew G. Knepley     PetscInt        xdof, xoff;
228*7b197a28SMatthew G. Knepley 
229*7b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, v, &supp));
230*7b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetDof(s, v, &xdof));
231*7b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetOffset(s, v, &xoff));
232*7b197a28SMatthew G. Knepley     // Diagonal
233*7b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < xdof; ++i) fa[xoff + i] += (M + 4) * ua[xoff + i];
234*7b197a28SMatthew G. Knepley     // Loop over mu
235*7b197a28SMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) {
236*7b197a28SMatthew G. Knepley       const PetscInt *cone;
237*7b197a28SMatthew G. Knepley       PetscInt        yoff;
238*7b197a28SMatthew G. Knepley 
239*7b197a28SMatthew G. Knepley       // Left action -(1 + \gamma_\mu)/2 \otimes U^\dagger_\mu(y) \delta_{x - \mu,y} \psi(y)
240*7b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[2 * d + 0], &cone));
241*7b197a28SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(s, cone[0], &yoff));
242*7b197a28SMatthew G. Knepley       PetscCall(ComputeAction(d, PETSC_FALSE, id, &ua[yoff], &fa[xoff]));
243*7b197a28SMatthew G. Knepley       // Right edge -(1 - \gamma_\mu)/2 \otimes U_\mu(x) \delta_{x + \mu,y} \psi(y)
244*7b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[2 * d + 1], &cone));
245*7b197a28SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(s, cone[1], &yoff));
246*7b197a28SMatthew G. Knepley       PetscCall(ComputeAction(d, PETSC_TRUE, id, &ua[yoff], &fa[xoff]));
247*7b197a28SMatthew G. Knepley     }
248*7b197a28SMatthew G. Knepley   }
249*7b197a28SMatthew G. Knepley   PetscCall(VecRestoreArray(f, &fa));
250*7b197a28SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(u, &ua));
251*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
252*7b197a28SMatthew G. Knepley }
253*7b197a28SMatthew G. Knepley 
254*7b197a28SMatthew G. Knepley static PetscErrorCode CreateJacobian(DM dm, Mat *J)
255*7b197a28SMatthew G. Knepley {
256*7b197a28SMatthew G. Knepley   PetscFunctionBegin;
257*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
258*7b197a28SMatthew G. Knepley }
259*7b197a28SMatthew G. Knepley 
260*7b197a28SMatthew G. Knepley static PetscErrorCode ComputeJacobian(DM dm, Vec u, Mat J)
261*7b197a28SMatthew G. Knepley {
262*7b197a28SMatthew G. Knepley   PetscFunctionBegin;
263*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
264*7b197a28SMatthew G. Knepley }
265*7b197a28SMatthew G. Knepley 
266*7b197a28SMatthew G. Knepley static PetscErrorCode PrintTraversal(DM dm)
267*7b197a28SMatthew G. Knepley {
268*7b197a28SMatthew G. Knepley   MPI_Comm comm;
269*7b197a28SMatthew G. Knepley   PetscInt vStart, vEnd;
270*7b197a28SMatthew G. Knepley 
271*7b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
272*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
273*7b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
274*7b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
275*7b197a28SMatthew G. Knepley     const PetscInt *supp;
276*7b197a28SMatthew G. Knepley     PetscInt        Ns;
277*7b197a28SMatthew G. Knepley 
278*7b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupportSize(dm, v, &Ns));
279*7b197a28SMatthew G. Knepley     PetscCall(DMPlexGetSupport(dm, v, &supp));
280*7b197a28SMatthew G. Knepley     PetscCall(PrintVertex(dm, v));
281*7b197a28SMatthew G. Knepley     PetscCall(PetscPrintf(comm, "\n"));
282*7b197a28SMatthew G. Knepley     for (PetscInt s = 0; s < Ns; ++s) {
283*7b197a28SMatthew G. Knepley       const PetscInt *cone;
284*7b197a28SMatthew G. Knepley 
285*7b197a28SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, supp[s], &cone));
286*7b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "  Edge %" PetscInt_FMT ": ", supp[s]));
287*7b197a28SMatthew G. Knepley       PetscCall(PrintVertex(dm, cone[0]));
288*7b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, " -- "));
289*7b197a28SMatthew G. Knepley       PetscCall(PrintVertex(dm, cone[1]));
290*7b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(comm, "\n"));
291*7b197a28SMatthew G. Knepley     }
292*7b197a28SMatthew G. Knepley   }
293*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
294*7b197a28SMatthew G. Knepley }
295*7b197a28SMatthew G. Knepley 
296*7b197a28SMatthew G. Knepley static PetscErrorCode ComputeFFT(Mat FT, PetscInt Nc, Vec x, Vec p)
297*7b197a28SMatthew G. Knepley {
298*7b197a28SMatthew G. Knepley   Vec     *xComp, *pComp;
299*7b197a28SMatthew G. Knepley   PetscInt n, N;
300*7b197a28SMatthew G. Knepley 
301*7b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
302*7b197a28SMatthew G. Knepley   PetscCall(PetscMalloc2(Nc, &xComp, Nc, &pComp));
303*7b197a28SMatthew G. Knepley   PetscCall(VecGetLocalSize(x, &n));
304*7b197a28SMatthew G. Knepley   PetscCall(VecGetSize(x, &N));
305*7b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
306*7b197a28SMatthew G. Knepley     const char *vtype;
307*7b197a28SMatthew G. Knepley 
308*7b197a28SMatthew G. Knepley     // HACK: Make these from another DM up front
309*7b197a28SMatthew G. Knepley     PetscCall(VecCreate(PetscObjectComm((PetscObject)x), &xComp[i]));
310*7b197a28SMatthew G. Knepley     PetscCall(VecGetType(x, &vtype));
311*7b197a28SMatthew G. Knepley     PetscCall(VecSetType(xComp[i], vtype));
312*7b197a28SMatthew G. Knepley     PetscCall(VecSetSizes(xComp[i], n / Nc, N / Nc));
313*7b197a28SMatthew G. Knepley     PetscCall(VecDuplicate(xComp[i], &pComp[i]));
314*7b197a28SMatthew G. Knepley   }
315*7b197a28SMatthew G. Knepley   PetscCall(VecStrideGatherAll(x, xComp, INSERT_VALUES));
316*7b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) PetscCall(MatMult(FT, xComp[i], pComp[i]));
317*7b197a28SMatthew G. Knepley   PetscCall(VecStrideScatterAll(pComp, p, INSERT_VALUES));
318*7b197a28SMatthew G. Knepley   for (PetscInt i = 0; i < Nc; ++i) {
319*7b197a28SMatthew G. Knepley     PetscCall(VecDestroy(&xComp[i]));
320*7b197a28SMatthew G. Knepley     PetscCall(VecDestroy(&pComp[i]));
321*7b197a28SMatthew G. Knepley   }
322*7b197a28SMatthew G. Knepley   PetscCall(PetscFree2(xComp, pComp));
323*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
324*7b197a28SMatthew G. Knepley }
325*7b197a28SMatthew G. Knepley 
326*7b197a28SMatthew G. Knepley /*
327*7b197a28SMatthew G. Knepley   Test the action of the Wilson operator in the free field case U = I,
328*7b197a28SMatthew G. Knepley 
329*7b197a28SMatthew G. Knepley     \eta(x) = D_W(x - y) \psi(y)
330*7b197a28SMatthew G. Knepley 
331*7b197a28SMatthew G. Knepley   The Wilson operator is a convolution for the free field, so we can check that by the convolution theorem
332*7b197a28SMatthew G. Knepley 
333*7b197a28SMatthew G. Knepley     \hat\eta(x) = \mathcal{F}(D_W(x - y) \psi(y))
334*7b197a28SMatthew G. Knepley                 = \hat D_W(p) \mathcal{F}\psi(p)
335*7b197a28SMatthew G. Knepley 
336*7b197a28SMatthew G. Knepley   The Fourier series for the Wilson operator is
337*7b197a28SMatthew G. Knepley 
338*7b197a28SMatthew G. Knepley     M + \sum_\mu 2 \sin^2(p_\mu / 2) + i \gamma_\mu \sin(p_\mu)
339*7b197a28SMatthew G. Knepley */
340*7b197a28SMatthew G. Knepley static PetscErrorCode TestFreeField(DM dm)
341*7b197a28SMatthew G. Knepley {
342*7b197a28SMatthew G. Knepley   PetscSection       s;
343*7b197a28SMatthew G. Knepley   Mat                FT;
344*7b197a28SMatthew G. Knepley   Vec                psi, psiHat;
345*7b197a28SMatthew G. Knepley   Vec                eta, etaHat;
346*7b197a28SMatthew G. Knepley   Vec                DHat; // The product \hat D_w \hat psi
347*7b197a28SMatthew G. Knepley   PetscRandom        r;
348*7b197a28SMatthew G. Knepley   const PetscScalar *psih;
349*7b197a28SMatthew G. Knepley   PetscScalar       *dh;
350*7b197a28SMatthew G. Knepley   PetscReal         *coef, nrm;
351*7b197a28SMatthew G. Knepley   const PetscInt    *extent, Nc = 12;
352*7b197a28SMatthew G. Knepley   PetscInt           dim, V     = 1, vStart, vEnd;
353*7b197a28SMatthew G. Knepley   PetscContainer     c;
354*7b197a28SMatthew G. Knepley   PetscBool          constRhs = PETSC_FALSE;
355*7b197a28SMatthew G. Knepley 
356*7b197a28SMatthew G. Knepley   PetscFunctionBeginUser;
357*7b197a28SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-const_rhs", &constRhs, NULL));
358*7b197a28SMatthew G. Knepley 
359*7b197a28SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &s));
360*7b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &psi));
361*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)psi, "psi"));
362*7b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &psiHat));
363*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)psiHat, "psihat"));
364*7b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &eta));
365*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)eta, "eta"));
366*7b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &etaHat));
367*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)etaHat, "etahat"));
368*7b197a28SMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &DHat));
369*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)DHat, "Dhat"));
370*7b197a28SMatthew G. Knepley   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
371*7b197a28SMatthew G. Knepley   PetscCall(PetscRandomSetType(r, PETSCRAND48));
372*7b197a28SMatthew G. Knepley   if (constRhs) PetscCall(VecSet(psi, 1.));
373*7b197a28SMatthew G. Knepley   else PetscCall(VecSetRandom(psi, r));
374*7b197a28SMatthew G. Knepley   PetscCall(PetscRandomDestroy(&r));
375*7b197a28SMatthew G. Knepley 
376*7b197a28SMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
377*7b197a28SMatthew G. Knepley   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
378*7b197a28SMatthew G. Knepley   PetscCall(PetscObjectQuery((PetscObject)dm, "_extent", (PetscObject *)&c));
379*7b197a28SMatthew G. Knepley   PetscCall(PetscContainerGetPointer(c, (void **)&extent));
380*7b197a28SMatthew G. Knepley   PetscCall(MatCreateFFT(PetscObjectComm((PetscObject)dm), dim, extent, MATFFTW, &FT));
381*7b197a28SMatthew G. Knepley 
382*7b197a28SMatthew G. Knepley   PetscCall(PetscMalloc1(dim, &coef));
383*7b197a28SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
384*7b197a28SMatthew G. Knepley     coef[d] = 2. * PETSC_PI / extent[d];
385*7b197a28SMatthew G. Knepley     V *= extent[d];
386*7b197a28SMatthew G. Knepley   }
387*7b197a28SMatthew G. Knepley   PetscCall(ComputeResidual(dm, psi, eta));
388*7b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(eta, NULL, "-psi_view"));
389*7b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(eta, NULL, "-eta_view"));
390*7b197a28SMatthew G. Knepley   PetscCall(ComputeFFT(FT, Nc, psi, psiHat));
391*7b197a28SMatthew G. Knepley   PetscCall(VecScale(psiHat, 1. / V));
392*7b197a28SMatthew G. Knepley   PetscCall(ComputeFFT(FT, Nc, eta, etaHat));
393*7b197a28SMatthew G. Knepley   PetscCall(VecScale(etaHat, 1. / V));
394*7b197a28SMatthew G. Knepley   PetscCall(VecGetArrayRead(psiHat, &psih));
395*7b197a28SMatthew G. Knepley   PetscCall(VecGetArray(DHat, &dh));
396*7b197a28SMatthew G. Knepley   for (PetscInt v = vStart; v < vEnd; ++v) {
397*7b197a28SMatthew G. Knepley     PetscScalar tmp[12], tmp1 = 0.;
398*7b197a28SMatthew G. Knepley     PetscInt    dof, off;
399*7b197a28SMatthew G. Knepley 
400*7b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetDof(s, v, &dof));
401*7b197a28SMatthew G. Knepley     PetscCall(PetscSectionGetOffset(s, v, &off));
402*7b197a28SMatthew G. Knepley     for (PetscInt d = 0, prod = 1; d < dim; prod *= extent[d], ++d) {
403*7b197a28SMatthew G. Knepley       const PetscInt idx = (v / prod) % extent[d];
404*7b197a28SMatthew G. Knepley 
405*7b197a28SMatthew G. Knepley       tmp1 += 2. * PetscSqr(PetscSinReal(0.5 * coef[d] * idx));
406*7b197a28SMatthew G. Knepley       for (PetscInt i = 0; i < dof; ++i) tmp[i] = PETSC_i * PetscSinReal(coef[d] * idx) * psih[off + i];
407*7b197a28SMatthew G. Knepley       for (PetscInt c = 0; c < 3; ++c) ComputeGamma(d, 3, &tmp[c]);
408*7b197a28SMatthew G. Knepley       for (PetscInt i = 0; i < dof; ++i) dh[off + i] += tmp[i];
409*7b197a28SMatthew G. Knepley     }
410*7b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < dof; ++i) dh[off + i] += (M + tmp1) * psih[off + i];
411*7b197a28SMatthew G. Knepley   }
412*7b197a28SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(psiHat, &psih));
413*7b197a28SMatthew G. Knepley   PetscCall(VecRestoreArray(DHat, &dh));
414*7b197a28SMatthew G. Knepley 
415*7b197a28SMatthew G. Knepley   {
416*7b197a28SMatthew G. Knepley     Vec     *etaComp, *DComp;
417*7b197a28SMatthew G. Knepley     PetscInt n, N;
418*7b197a28SMatthew G. Knepley 
419*7b197a28SMatthew G. Knepley     PetscCall(PetscMalloc2(Nc, &etaComp, Nc, &DComp));
420*7b197a28SMatthew G. Knepley     PetscCall(VecGetLocalSize(etaHat, &n));
421*7b197a28SMatthew G. Knepley     PetscCall(VecGetSize(etaHat, &N));
422*7b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
423*7b197a28SMatthew G. Knepley       const char *vtype;
424*7b197a28SMatthew G. Knepley 
425*7b197a28SMatthew G. Knepley       // HACK: Make these from another DM up front
426*7b197a28SMatthew G. Knepley       PetscCall(VecCreate(PetscObjectComm((PetscObject)etaHat), &etaComp[i]));
427*7b197a28SMatthew G. Knepley       PetscCall(VecGetType(etaHat, &vtype));
428*7b197a28SMatthew G. Knepley       PetscCall(VecSetType(etaComp[i], vtype));
429*7b197a28SMatthew G. Knepley       PetscCall(VecSetSizes(etaComp[i], n / Nc, N / Nc));
430*7b197a28SMatthew G. Knepley       PetscCall(VecDuplicate(etaComp[i], &DComp[i]));
431*7b197a28SMatthew G. Knepley     }
432*7b197a28SMatthew G. Knepley     PetscCall(VecStrideGatherAll(etaHat, etaComp, INSERT_VALUES));
433*7b197a28SMatthew G. Knepley     PetscCall(VecStrideGatherAll(DHat, DComp, INSERT_VALUES));
434*7b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
435*7b197a28SMatthew G. Knepley       if (!i) {
436*7b197a28SMatthew G. Knepley         PetscCall(VecViewFromOptions(etaComp[i], NULL, "-etahat_view"));
437*7b197a28SMatthew G. Knepley         PetscCall(VecViewFromOptions(DComp[i], NULL, "-dhat_view"));
438*7b197a28SMatthew G. Knepley       }
439*7b197a28SMatthew G. Knepley       PetscCall(VecAXPY(etaComp[i], -1., DComp[i]));
440*7b197a28SMatthew G. Knepley       PetscCall(VecNorm(etaComp[i], NORM_INFINITY, &nrm));
441*7b197a28SMatthew G. Knepley       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Slice %" PetscInt_FMT ": %g\n", i, (double)nrm));
442*7b197a28SMatthew G. Knepley     }
443*7b197a28SMatthew G. Knepley     PetscCall(VecStrideScatterAll(etaComp, etaHat, INSERT_VALUES));
444*7b197a28SMatthew G. Knepley     for (PetscInt i = 0; i < Nc; ++i) {
445*7b197a28SMatthew G. Knepley       PetscCall(VecDestroy(&etaComp[i]));
446*7b197a28SMatthew G. Knepley       PetscCall(VecDestroy(&DComp[i]));
447*7b197a28SMatthew G. Knepley     }
448*7b197a28SMatthew G. Knepley     PetscCall(PetscFree2(etaComp, DComp));
449*7b197a28SMatthew G. Knepley     PetscCall(VecNorm(etaHat, NORM_INFINITY, &nrm));
450*7b197a28SMatthew G. Knepley     PetscCheck(nrm < PETSC_SMALL, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Free field test failed: %g", (double)nrm);
451*7b197a28SMatthew G. Knepley   }
452*7b197a28SMatthew G. Knepley 
453*7b197a28SMatthew G. Knepley   PetscCall(PetscFree(coef));
454*7b197a28SMatthew G. Knepley   PetscCall(MatDestroy(&FT));
455*7b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &psi));
456*7b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &psiHat));
457*7b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &eta));
458*7b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &etaHat));
459*7b197a28SMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &DHat));
460*7b197a28SMatthew G. Knepley   PetscFunctionReturn(0);
461*7b197a28SMatthew G. Knepley }
462*7b197a28SMatthew G. Knepley 
463*7b197a28SMatthew G. Knepley int main(int argc, char **argv)
464*7b197a28SMatthew G. Knepley {
465*7b197a28SMatthew G. Knepley   DM     dm;
466*7b197a28SMatthew G. Knepley   Vec    u, f;
467*7b197a28SMatthew G. Knepley   Mat    J;
468*7b197a28SMatthew G. Knepley   AppCtx user;
469*7b197a28SMatthew G. Knepley 
470*7b197a28SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
471*7b197a28SMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
472*7b197a28SMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
473*7b197a28SMatthew G. Knepley   PetscCall(SetupDiscretization(dm, &user));
474*7b197a28SMatthew G. Knepley   PetscCall(SetupAuxDiscretization(dm, &user));
475*7b197a28SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &u));
476*7b197a28SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &f));
477*7b197a28SMatthew G. Knepley   PetscCall(PrintTraversal(dm));
478*7b197a28SMatthew G. Knepley   PetscCall(ComputeResidual(dm, u, f));
479*7b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(f, NULL, "-res_view"));
480*7b197a28SMatthew G. Knepley   PetscCall(CreateJacobian(dm, &J));
481*7b197a28SMatthew G. Knepley   PetscCall(ComputeJacobian(dm, u, J));
482*7b197a28SMatthew G. Knepley   PetscCall(VecViewFromOptions(u, NULL, "-sol_view"));
483*7b197a28SMatthew G. Knepley   PetscCall(TestFreeField(dm));
484*7b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&f));
485*7b197a28SMatthew G. Knepley   PetscCall(VecDestroy(&u));
486*7b197a28SMatthew G. Knepley   PetscCall(DMDestroy(&dm));
487*7b197a28SMatthew G. Knepley   PetscCall(PetscFinalize());
488*7b197a28SMatthew G. Knepley   return 0;
489*7b197a28SMatthew G. Knepley }
490*7b197a28SMatthew G. Knepley 
491*7b197a28SMatthew G. Knepley /*TEST
492*7b197a28SMatthew G. Knepley 
493*7b197a28SMatthew G. Knepley   build:
494*7b197a28SMatthew G. Knepley     requires: complex
495*7b197a28SMatthew G. Knepley 
496*7b197a28SMatthew G. Knepley   test:
497*7b197a28SMatthew G. Knepley     requires: fftw
498*7b197a28SMatthew G. Knepley     suffix: dirac_free_field
499*7b197a28SMatthew G. Knepley     args: -dm_plex_dim 4 -dm_plex_shape hypercubic -dm_plex_box_faces 3,3,3,3 -dm_view -sol_view \
500*7b197a28SMatthew G. Knepley           -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces -dm_plex_check_pointsf
501*7b197a28SMatthew G. Knepley 
502*7b197a28SMatthew G. Knepley TEST*/
503