xref: /petsc/src/dm/dt/interface/dtaltv.c (revision 28222859993257eecc81314090ee1fbecec5d4dc)
1 #include <petsc/private/petscimpl.h>
2 #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
3 
4 /*@
5    PetscDTAltVApply - Apply an a k-form (an alternating k-linear map) to a set of k N-dimensional vectors
6 
7    Input Arguments:
8 +  N - the dimension of the vector space, N >= 0
9 .  k - the degree k of the k-form w, 0 <= k <= N
10 .  w - a k-form, size [N choose k] (each degree of freedom of a k-form is associated with a subset of k coordinates of the N-dimensional vectors: the degrees of freedom are ordered lexicographically by their associated subsets)
11 -  v - a set of k vectors of size N, size [k x N], each vector stored contiguously
12 
13    Output Arguments:
14 .  wv - w(v_1,...,v_k) = \sum_i w_i * det(V_i): the degree of freedom w_i is associated with coordinates [s_{i,1},...,s_{i,k}], and the square matrix V_i has entry (j,k) given by the s_{i,k}'th coordinate of v_j
15 
16    Level: intermediate
17 
18 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
19 @*/
20 PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv)
21 {
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
26   if (k < 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
27   if (N <= 3) {
28     if (!k) {
29       *wv = w[0];
30     } else {
31       if (N == 1)        {*wv = w[0] * v[0];}
32       else if (N == 2) {
33         if (k == 1)      {*wv = w[0] * v[0] + w[1] * v[1];}
34         else             {*wv = w[0] * (v[0] * v[3] - v[1] * v[2]);}
35       } else {
36         if (k == 1)      {*wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];}
37         else if (k == 2) {
38           *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) +
39                 w[1] * (v[0] * v[5] - v[2] * v[3]) +
40                 w[2] * (v[1] * v[5] - v[2] * v[4]);
41         } else {
42           *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) +
43                         v[1] * (v[5] * v[6] - v[3] * v[8]) +
44                         v[2] * (v[3] * v[7] - v[4] * v[6]));
45         }
46       }
47     }
48   } else {
49     PetscInt Nk, Nf;
50     PetscInt *subset, *perm;
51     PetscInt i, j, l;
52     PetscReal sum = 0.;
53 
54     ierr = PetscDTFactorialInt(k, &Nf);CHKERRQ(ierr);
55     ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr);
56     ierr = PetscMalloc2(k, &subset, k, &perm);CHKERRQ(ierr);
57     for (i = 0; i < Nk; i++) {
58       PetscReal subsum = 0.;
59 
60       ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr);
61       for (j = 0; j < Nf; j++) {
62         PetscBool permOdd;
63         PetscReal prod;
64 
65         ierr = PetscDTEnumPerm(k, j, perm, &permOdd);CHKERRQ(ierr);
66         prod = permOdd ? -1. : 1.;
67         for (l = 0; l < k; l++) {
68           prod *= v[perm[l] * N + subset[l]];
69         }
70         subsum += prod;
71       }
72       sum += w[i] * subsum;
73     }
74     ierr = PetscFree2(subset, perm);CHKERRQ(ierr);
75     *wv = sum;
76   }
77   PetscFunctionReturn(0);
78 }
79 
80 /*@
81    PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form
82 
83    Input Arguments:
84 +  N - the dimension of the vector space, N >= 0
85 .  j - the degree j of the j-form a, 0 <= j <= N
86 .  k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N
87 .  a - a j-form, size [N choose j]
88 -  b - a k-form, size [N choose k]
89 
90    Output Arguments:
91 .  awedgeb - the (j+k)-form a wedge b, size [N choose (j+k)]: (a wedge b)(v_1,...,v_{j+k}) = \sum_{s} sign(s) a(v_{s_1},...,v_{s_j}) b(v_{s_{j+1}},...,v_{s_{j+k}}),
92              where the sum is over permutations s such that s_1 < s_2 < ... < s_j and s_{j+1} < s_{j+2} < ... < s_{j+k}.
93 
94    Level: intermediate
95 
96 .seealso: PetscDTAltVWedgeMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
97 @*/
98 PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb)
99 {
100   PetscInt       i;
101   PetscErrorCode ierr;
102 
103   PetscFunctionBegin;
104   if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
105   if (j < 0 || k < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
106   if (j + k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
107   if (N <= 3) {
108     PetscInt Njk;
109 
110     ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr);
111     if (!j)      {for (i = 0; i < Njk; i++) {awedgeb[i] = a[0] * b[i];}}
112     else if (!k) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[i] * b[0];}}
113     else {
114       if (N == 2) {awedgeb[0] = a[0] * b[1] - a[1] * b[0];}
115       else {
116         if (j+k == 2) {
117           awedgeb[0] = a[0] * b[1] - a[1] * b[0];
118           awedgeb[1] = a[0] * b[2] - a[2] * b[0];
119           awedgeb[2] = a[1] * b[2] - a[2] * b[1];
120         } else {
121           awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0];
122         }
123       }
124     }
125   } else {
126     PetscInt  Njk;
127     PetscInt  JKj;
128     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
129     PetscInt  i;
130 
131     ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr);
132     ierr = PetscDTBinomialInt(j+k, j, &JKj);CHKERRQ(ierr);
133     ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr);
134     for (i = 0; i < Njk; i++) {
135       PetscReal sum = 0.;
136       PetscInt  l;
137 
138       ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr);
139       for (l = 0; l < JKj; l++) {
140         PetscBool jkOdd;
141         PetscInt  m, jInd, kInd;
142 
143         ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr);
144         for (m = 0; m < j; m++) {
145           subsetj[m] = subset[subsetjk[m]];
146         }
147         for (m = 0; m < k; m++) {
148           subsetk[m] = subset[subsetjk[j+m]];
149         }
150         ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr);
151         ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr);
152         sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]);
153       }
154       awedgeb[i] = sum;
155     }
156     ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr);
157   }
158   PetscFunctionReturn(0);
159 }
160 
161 /*@
162    PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms
163 
164    Input Arguments:
165 +  N - the dimension of the vector space, N >= 0
166 .  j - the degree j of the j-form a, 0 <= j <= N
167 .  k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N
168 -  a - a j-form, size [N choose j]
169 
170    Output Arguments:
171 .  awedge - (a wedge), an [(N choose j+k) x (N choose k)] matrix in row-major order, such that (a wedge) * b = a wedge b
172 
173    Level: intermediate
174 
175 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
176 @*/
177 PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat)
178 {
179   PetscInt       i;
180   PetscErrorCode ierr;
181 
182   PetscFunctionBegin;
183   if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
184   if (j < 0 || k < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
185   if (j + k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
186   if (N <= 3) {
187     PetscInt Njk;
188 
189     ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr);
190     if (!j) {
191       for (i = 0; i < Njk * Njk; i++) {awedgeMat[i] = 0.;}
192       for (i = 0; i < Njk; i++) {awedgeMat[i * (Njk + 1)] = a[0];}
193     } else if (!k) {
194       for (i = 0; i < Njk; i++) {awedgeMat[i] = a[i];}
195     } else {
196       if (N == 2) {
197         awedgeMat[0] = -a[1]; awedgeMat[1] =  a[0];
198       } else {
199         if (j+k == 2) {
200           awedgeMat[0] = -a[1]; awedgeMat[1] =  a[0]; awedgeMat[2] =    0.;
201           awedgeMat[3] = -a[2]; awedgeMat[4] =    0.; awedgeMat[5] =  a[0];
202           awedgeMat[6] =    0.; awedgeMat[7] = -a[2]; awedgeMat[8] =  a[1];
203         } else {
204           awedgeMat[0] =  a[2]; awedgeMat[1] = -a[1]; awedgeMat[2] =  a[0];
205         }
206       }
207     }
208   } else {
209     PetscInt  Njk;
210     PetscInt  Nk;
211     PetscInt  JKj, i;
212     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
213 
214     ierr = PetscDTBinomialInt(N,   k,   &Nk);CHKERRQ(ierr);
215     ierr = PetscDTBinomialInt(N,   j+k, &Njk);CHKERRQ(ierr);
216     ierr = PetscDTBinomialInt(j+k, j,   &JKj);CHKERRQ(ierr);
217     ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr);
218     for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.;
219     for (i = 0; i < Njk; i++) {
220       PetscInt  l;
221 
222       ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr);
223       for (l = 0; l < JKj; l++) {
224         PetscBool jkOdd;
225         PetscInt  m, jInd, kInd;
226 
227         ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr);
228         for (m = 0; m < j; m++) {
229           subsetj[m] = subset[subsetjk[m]];
230         }
231         for (m = 0; m < k; m++) {
232           subsetk[m] = subset[subsetjk[j+m]];
233         }
234         ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr);
235         ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr);
236         awedgeMat[i * Nk + kInd] += jkOdd ? - a[jInd] : a[jInd];
237       }
238     }
239     ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr);
240   }
241   PetscFunctionReturn(0);
242 }
243 
244 /*@
245    PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space
246 
247    Input Arguments:
248 +  N - the dimension of the origin vector space of the linear transformation, M >= 0
249 .  M - the dimension of the image vector space of the linear transformation, N >= 0
250 .  L - a linear transformation, an [M x N] matrix in row-major format
251 .  k - the *signed* degree k of the |k|-form w, -(min(M,N)) <= k <= min(M,N).  A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note).
252 -  w - a |k|-form in the image space, size [M choose |k|]
253 
254    Output Arguments:
255 .  Lstarw - the pullback of w to a |k|-form in the origin space, size [N choose |k|]: (Lstarw)(v_1,...v_k) = w(L*v_1,...,L*v_k).
256 
257    Level: intermediate
258 
259    Note: negative form degrees accomodate, e.g., H-div conforming vector fields.  An H-div conforming vector field stores its degrees of freedom as (dx, dy, dz), like a 1-form,
260    but its normal trace is integrated on faces, like a 2-form.  The correct pullback then is to apply the Hodge star transformation from (M-2)-form to 2-form, pullback as a 2-form,
261    then the inverse Hodge star transformation.
262 
263 .seealso: PetscDTAltVPullbackMatrix(), PetscDTAltVStar()
264 @*/
265 PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw)
266 {
267   PetscInt         i, j, Nk, Mk;
268   PetscErrorCode   ierr;
269 
270   PetscFunctionBegin;
271   if (N < 0 || M < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
272   if (PetscAbsInt(k) > N || PetscAbsInt(k) > M) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
273   if (N <= 3 && M <= 3) {
274 
275     ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr);
276     ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
277     if (!k) {
278       Lstarw[0] = w[0];
279     } else if (k == 1) {
280       for (i = 0; i < Nk; i++) {
281         PetscReal sum = 0.;
282 
283         for (j = 0; j < Mk; j++) {sum += L[j * Nk + i] * w[j];}
284         Lstarw[i] = sum;
285       }
286     } else if (k == -1) {
287       PetscReal mult[3] = {1., -1., 1.};
288 
289       for (i = 0; i < Nk; i++) {
290         PetscReal sum = 0.;
291 
292         for (j = 0; j < Mk; j++) {
293           sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j];
294         }
295         Lstarw[i] = mult[i] * sum;
296       }
297     } else if (k == 2) {
298       PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}};
299 
300       for (i = 0; i < Nk; i++) {
301         PetscReal sum = 0.;
302         for (j = 0; j < Mk; j++) {
303           sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] -
304                   L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j];
305         }
306         Lstarw[i] = sum;
307       }
308     } else if (k == -2) {
309       PetscInt  pairs[3][2] = {{1,2},{2,0},{0,1}};
310       PetscInt  offi = (N == 2) ? 2 : 0;
311       PetscInt  offj = (M == 2) ? 2 : 0;
312 
313       for (i = 0; i < Nk; i++) {
314         PetscReal sum   = 0.;
315 
316         for (j = 0; j < Mk; j++) {
317           sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] *
318                   L[pairs[offj + j][1] * N + pairs[offi + i][1]] -
319                   L[pairs[offj + j][1] * N + pairs[offi + i][0]] *
320                   L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j];
321 
322         }
323         Lstarw[i] = sum;
324       }
325     } else {
326       PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) +
327                        L[1] * (L[5] * L[6] - L[3] * L[8]) +
328                        L[2] * (L[3] * L[7] - L[4] * L[6]);
329 
330       for (i = 0; i < Nk; i++) {Lstarw[i] = detL * w[i];}
331     }
332   } else {
333     PetscInt         Nf, l, p;
334     PetscReal       *Lw, *Lwv;
335     PetscInt        *subsetw, *subsetv;
336     PetscInt        *perm;
337     PetscReal       *walloc = NULL;
338     const PetscReal *ww = NULL;
339     PetscBool        negative = PETSC_FALSE;
340 
341     ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr);
342     ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
343     ierr = PetscDTFactorialInt(PetscAbsInt(k), &Nf);CHKERRQ(ierr);
344     if (k < 0) {
345       negative = PETSC_TRUE;
346       k = -k;
347       ierr = PetscMalloc1(Mk, &walloc);CHKERRQ(ierr);
348       ierr = PetscDTAltVStar(M, M - k, 1, w, walloc);CHKERRQ(ierr);
349       ww = walloc;
350     } else {
351       ww = w;
352     }
353     ierr = PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr);
354     for (i = 0; i < Nk; i++) Lstarw[i] = 0.;
355     for (i = 0; i < Mk; i++) {
356       ierr = PetscDTEnumSubset(M, k, i, subsetw);CHKERRQ(ierr);
357       for (j = 0; j < Nk; j++) {
358         ierr = PetscDTEnumSubset(N, k, j, subsetv);CHKERRQ(ierr);
359         for (p = 0; p < Nf; p++) {
360           PetscReal prod;
361           PetscBool isOdd;
362 
363           ierr = PetscDTEnumPerm(k, p, perm, &isOdd);CHKERRQ(ierr);
364           prod = isOdd ? -ww[i] : ww[i];
365           for (l = 0; l < k; l++) {
366             prod *= L[subsetw[perm[l]] * N + subsetv[l]];
367           }
368           Lstarw[j] += prod;
369         }
370       }
371     }
372     if (negative) {
373       PetscReal *sLsw;
374 
375       ierr = PetscMalloc1(Nk, &sLsw);CHKERRQ(ierr);
376       ierr = PetscDTAltVStar(N, N - k, -1,  Lstarw, sLsw);CHKERRQ(ierr);
377       for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i];
378       ierr = PetscFree(sLsw);CHKERRQ(ierr);
379     }
380     ierr = PetscFree5(subsetw, subsetv, perm, Lw, Lwv);CHKERRQ(ierr);
381     ierr = PetscFree(walloc);CHKERRQ(ierr);
382   }
383   PetscFunctionReturn(0);
384 }
385 
386 /*@
387    PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation
388 
389    Input Arguments:
390 +  N - the dimension of the origin vector space of the linear transformation, N >= 0
391 .  M - the dimension of the image vector space of the linear transformation, M >= 0
392 .  L - a linear transformation, an [M x N] matrix in row-major format
393 -  k - the *signed* degree k of the |k|-forms on which Lstar acts, -(min(M,N)) <= k <= min(M,N).  A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note in PetscDTAltvPullback())
394 
395    Output Arguments:
396 .  Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w
397 
398    Level: intermediate
399 
400 .seealso: PetscDTAltVPullback(), PetscDTAltVStar()
401 @*/
402 PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar)
403 {
404   PetscInt        Nk, Mk, Nf, i, j, l, p;
405   PetscReal      *Lw, *Lwv;
406   PetscInt       *subsetw, *subsetv;
407   PetscInt       *perm;
408   PetscBool       negative = PETSC_FALSE;
409   PetscErrorCode  ierr;
410 
411   PetscFunctionBegin;
412   if (N < 0 || M < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
413   if (PetscAbsInt(k) > N || PetscAbsInt(k) > M) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
414   if (N <= 3 && M <= 3) {
415     PetscReal mult[3] = {1., -1., 1.};
416 
417     ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr);
418     ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
419     if (!k) {
420       Lstar[0] = 1.;
421     } else if (k == 1) {
422       for (i = 0; i < Nk; i++) {for (j = 0; j < Mk; j++) {Lstar[i * Mk + j] = L[j * Nk + i];}}
423     } else if (k == -1) {
424       for (i = 0; i < Nk; i++) {
425         for (j = 0; j < Mk; j++) {
426           Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j];
427         }
428       }
429     } else if (k == 2) {
430       PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}};
431 
432       for (i = 0; i < Nk; i++) {
433         for (j = 0; j < Mk; j++) {
434           Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] *
435                               L[pairs[j][1] * N + pairs[i][1]] -
436                               L[pairs[j][1] * N + pairs[i][0]] *
437                               L[pairs[j][0] * N + pairs[i][1]];
438         }
439       }
440     } else if (k == -2) {
441       PetscInt  pairs[3][2] = {{1,2},{2,0},{0,1}};
442       PetscInt  offi = (N == 2) ? 2 : 0;
443       PetscInt  offj = (M == 2) ? 2 : 0;
444 
445       for (i = 0; i < Nk; i++) {
446         for (j = 0; j < Mk; j++) {
447           Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] *
448                               L[pairs[offj + j][1] * N + pairs[offi + i][1]] -
449                               L[pairs[offj + j][1] * N + pairs[offi + i][0]] *
450                               L[pairs[offj + j][0] * N + pairs[offi + i][1]];
451         }
452       }
453     } else {
454       PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) +
455                        L[1] * (L[5] * L[6] - L[3] * L[8]) +
456                        L[2] * (L[3] * L[7] - L[4] * L[6]);
457 
458       for (i = 0; i < Nk; i++) {Lstar[i] = detL;}
459     }
460   } else {
461     if (k < 0) {
462       negative = PETSC_TRUE;
463       k = -k;
464     }
465     ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr);
466     ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
467     ierr = PetscDTFactorialInt(PetscAbsInt(k), &Nf);CHKERRQ(ierr);
468     ierr = PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr);
469     for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.;
470     for (i = 0; i < Mk; i++) {
471       PetscBool iOdd;
472       PetscInt  iidx, jidx;
473 
474       ierr = PetscDTEnumSplit(M, k, i, subsetw, &iOdd);CHKERRQ(ierr);
475       iidx = negative ? Mk - 1 - i : i;
476       iOdd = negative ? (PetscBool) (iOdd ^ ((k * (M-k)) & 1)) : PETSC_FALSE;
477       for (j = 0; j < Nk; j++) {
478         PetscBool jOdd;
479 
480         ierr = PetscDTEnumSplit(N, k, j, subsetv, &jOdd);CHKERRQ(ierr);
481         jidx = negative ? Nk - 1 - j : j;
482         jOdd = negative ? (PetscBool) (iOdd ^ jOdd ^ ((k * (N-k)) & 1)) : PETSC_FALSE;
483         for (p = 0; p < Nf; p++) {
484           PetscReal prod;
485           PetscBool isOdd;
486 
487           ierr = PetscDTEnumPerm(k, p, perm, &isOdd);CHKERRQ(ierr);
488           isOdd = (PetscBool) (isOdd ^ jOdd);
489           prod = isOdd ? -1. : 1.;
490           for (l = 0; l < k; l++) {
491             prod *= L[subsetw[perm[l]] * N + subsetv[l]];
492           }
493           Lstar[jidx * Mk + iidx] += prod;
494         }
495       }
496     }
497     ierr = PetscFree5(subsetw, subsetv, perm, Lw, Lwv);CHKERRQ(ierr);
498   }
499   PetscFunctionReturn(0);
500 }
501 
502 /*@
503    PetscDTAltVInterior - Compute the interior product of a k-form with a vector
504 
505    Input Arguments:
506 +  N - the dimension of the vector space, N >= 0
507 .  k - the degree k of the k-form w, 0 <= k <= N
508 .  w - a k-form, size [N choose k]
509 -  v - an N dimensional vector
510 
511    Output Arguments:
512 .  wIntv - the (k-1)-form (w int v), size [N choose (k-1)]: (w int v) is defined by its action on (k-1) vectors {v_1, ..., v_{k-1}} as (w inv v)(v_1, ..., v_{k-1}) = w(v, v_1, ..., v_{k-1}).
513 
514    Level: intermediate
515 
516 .seealso: PetscDTAltVInteriorMatrix(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
517 @*/
518 PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv)
519 {
520   PetscInt        i, Nk, Nkm;
521   PetscErrorCode  ierr;
522 
523   PetscFunctionBegin;
524   if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
525   ierr = PetscDTBinomialInt(N, k,   &Nk);CHKERRQ(ierr);
526   ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr);
527   if (N <= 3) {
528     if (k == 1) {
529       PetscReal sum = 0.;
530 
531       for (i = 0; i < N; i++) {
532         sum += w[i] * v[i];
533       }
534       wIntv[0] = sum;
535     } else if (k == N) {
536       PetscReal mult[3] = {1., -1., 1.};
537 
538       for (i = 0; i < N; i++) {
539         wIntv[N - 1 - i] = w[0] * v[i] * mult[i];
540       }
541     } else {
542       wIntv[0] = - w[0]*v[1] - w[1]*v[2];
543       wIntv[1] =   w[0]*v[0] - w[2]*v[2];
544       wIntv[2] =   w[1]*v[0] + w[2]*v[1];
545     }
546   } else {
547     PetscInt       *subset, *work;
548 
549     ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr);
550     for (i = 0; i < Nkm; i++) wIntv[i] = 0.;
551     for (i = 0; i < Nk; i++) {
552       PetscInt  j, l, m;
553 
554       ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr);
555       for (j = 0; j < k; j++) {
556         PetscInt  idx;
557         PetscBool flip = (PetscBool) (j & 1);
558 
559         for (l = 0, m = 0; l < k; l++) {
560           if (l != j) work[m++] = subset[l];
561         }
562         ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr);
563         wIntv[idx] += flip ? -(w[i] * v[subset[j]]) :  (w[i] * v[subset[j]]);
564       }
565     }
566     ierr = PetscFree2(subset, work);CHKERRQ(ierr);
567   }
568   PetscFunctionReturn(0);
569 }
570 
571 /*@
572    PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector
573 
574    Input Arguments:
575 +  N - the dimension of the vector space, N >= 0
576 .  k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N
577 -  v - an N dimensional vector
578 
579    Output Arguments:
580 .  intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v)
581 
582    Level: intermediate
583 
584 .seealso: PetscDTAltVInterior(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
585 @*/
586 PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat)
587 {
588   PetscInt        i, Nk, Nkm;
589   PetscErrorCode  ierr;
590 
591   PetscFunctionBegin;
592   if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
593   ierr = PetscDTBinomialInt(N, k,   &Nk);CHKERRQ(ierr);
594   ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr);
595   if (N <= 3) {
596     if (k == 1) {
597       for (i = 0; i < N; i++) intvMat[i] = v[i];
598     } else if (k == N) {
599       PetscReal mult[3] = {1., -1., 1.};
600 
601       for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i];
602     } else {
603       intvMat[0] = -v[1]; intvMat[1] = -v[2]; intvMat[2] =    0.;
604       intvMat[3] =  v[0]; intvMat[4] =    0.; intvMat[5] = -v[2];
605       intvMat[6] =    0.; intvMat[7] =  v[0]; intvMat[8] =  v[1];
606     }
607   } else {
608     PetscInt       *subset, *work;
609 
610     ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr);
611     for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.;
612     for (i = 0; i < Nk; i++) {
613       PetscInt  j, l, m;
614 
615       ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr);
616       for (j = 0; j < k; j++) {
617         PetscInt  idx;
618         PetscBool flip = (PetscBool) (j & 1);
619 
620         for (l = 0, m = 0; l < k; l++) {
621           if (l != j) work[m++] = subset[l];
622         }
623         ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr);
624         intvMat[idx * Nk + i] += flip ? -v[subset[j]] :  v[subset[j]];
625       }
626     }
627     ierr = PetscFree2(subset, work);CHKERRQ(ierr);
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 /*@
633    PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix()
634 
635    Input Arguments:
636 +  N - the dimension of the vector space, N >= 0
637 -  k - the degree of the k-forms on which intvMat from PetscDTAltVInteriorMatrix() acts, 0 <= k <= N.
638 
639    Output Arguments:
640 .  indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k
641              non-zeros.  indices[i][0] and indices[i][1] are the row and column of a non-zero, and its value is equal to the vector
642              coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1)
643 
644    Level: intermediate
645 
646    Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential
647 
648 .seealso: PetscDTAltVInterior(), PetscDTAltVInteriorMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
649 @*/
650 PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3])
651 {
652   PetscInt        i, Nk, Nkm;
653   PetscErrorCode  ierr;
654 
655   PetscFunctionBegin;
656   if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
657   ierr = PetscDTBinomialInt(N, k,   &Nk);CHKERRQ(ierr);
658   ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr);
659   if (N <= 3) {
660     if (k == 1) {
661       for (i = 0; i < N; i++) {
662         indices[i][0] = 0;
663         indices[i][1] = i;
664         indices[i][2] = i;
665       }
666     } else if (k == N) {
667       PetscInt val[3] = {0, -2, 2};
668 
669       for (i = 0; i < N; i++) {
670         indices[i][0] = N - 1 - i;
671         indices[i][1] = 0;
672         indices[i][2] = val[i];
673       }
674     } else {
675       indices[0][0] = 0; indices[0][1] = 0; indices[0][2] = -(1 + 1);
676       indices[1][0] = 0; indices[1][1] = 1; indices[1][2] = -(2 + 1);
677       indices[2][0] = 1; indices[2][1] = 0; indices[2][2] = 0;
678       indices[3][0] = 1; indices[3][1] = 2; indices[3][2] = -(2 + 1);
679       indices[4][0] = 2; indices[4][1] = 1; indices[4][2] = 0;
680       indices[5][0] = 2; indices[5][1] = 2; indices[5][2] = 1;
681     }
682   } else {
683     PetscInt       *subset, *work;
684 
685     ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr);
686     for (i = 0; i < Nk; i++) {
687       PetscInt  j, l, m;
688 
689       ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr);
690       for (j = 0; j < k; j++) {
691         PetscInt  idx;
692         PetscBool flip = (PetscBool) (j & 1);
693 
694         for (l = 0, m = 0; l < k; l++) {
695           if (l != j) work[m++] = subset[l];
696         }
697         ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr);
698         indices[i * k + j][0] = idx;
699         indices[i * k + j][1] = i;
700         indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j];
701       }
702     }
703     ierr = PetscFree2(subset, work);CHKERRQ(ierr);
704   }
705   PetscFunctionReturn(0);
706 }
707 
708 /*@
709    PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form
710 
711    Input Arguments:
712 +  N - the dimension of the vector space, N >= 0
713 .  k - the degree k of the k-form w, 0 <= k <= N
714 .  pow - the number of times to apply the Hodge star operator: pow < 0 indicates that the inverse of the Hodge star operator should be applied |pow| times.
715 -  w - a k-form, size [N choose k]
716 
717    Output Arguments:
718 .  starw = (star)^pow w.  Each degree of freedom of a k-form is associated with a subset S of k coordinates of the N dimensional vector space: the Hodge start operator (star) maps that degree of freedom to the degree of freedom associated with S', the complement of S, with a sign change if the permutation of coordinates {S[0], ... S[k-1], S'[0], ... S'[N-k- 1]} is an odd permutation.  This implies (star)^2 w = (-1)^{k(N-k)} w, and (star)^4 w = w.
719 
720    Level: intermediate
721 
722 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
723 @*/
724 PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw)
725 {
726   PetscInt        Nk, i;
727   PetscErrorCode  ierr;
728 
729   PetscFunctionBegin;
730   if (k < 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
731   ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr);
732   pow = pow % 4;
733   pow = (pow + 4) % 4; /* make non-negative */
734   /* pow is now 0, 1, 2, 3 */
735   if (N <= 3) {
736     if (pow & 1) {
737       PetscReal mult[3] = {1., -1., 1.};
738 
739       for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i];
740     } else {
741       for (i = 0; i < Nk; i++) starw[i] = w[i];
742     }
743     if (pow > 1 && ((k * (N - k)) & 1)) {
744       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
745     }
746   } else {
747     PetscInt       *subset;
748 
749     ierr = PetscMalloc1(N, &subset);CHKERRQ(ierr);
750     if (pow % 2) {
751       PetscInt l = (pow == 1) ? k : N - k;
752       for (i = 0; i < Nk; i++) {
753         PetscBool sOdd;
754         PetscInt  j, idx;
755 
756         ierr = PetscDTEnumSplit(N, l, i, subset, &sOdd);CHKERRQ(ierr);
757         ierr = PetscDTSubsetIndex(N, l, subset, &idx);CHKERRQ(ierr);
758         ierr = PetscDTSubsetIndex(N, N-l, &subset[l], &j);CHKERRQ(ierr);
759         starw[j] = sOdd ? -w[idx] : w[idx];
760       }
761     } else {
762       for (i = 0; i < Nk; i++) starw[i] = w[i];
763     }
764     /* star^2 = -1^(k * (N - k)) */
765     if (pow > 1 && (k * (N - k)) % 2) {
766       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
767     }
768     ierr = PetscFree(subset);CHKERRQ(ierr);
769   }
770   PetscFunctionReturn(0);
771 }
772