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