xref: /petsc/src/dm/dt/interface/dtaltv.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
11a989b97SToby Isaac #include <petsc/private/petscimpl.h>
228222859SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
31a989b97SToby Isaac 
429a920c6SToby Isaac /*MC
529a920c6SToby Isaac    PetscDTAltV - An interface for common operations on k-forms, also known as alternating algebraic forms or alternating k-linear maps.
629a920c6SToby Isaac    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*.
729a920c6SToby Isaac 
829a920c6SToby Isaac    A recommended reference for this material is Section 2 "Exterior algebra and exterior calculus" in "Finite element
929a920c6SToby Isaac    exterior calculus, homological techniques, and applications", by Arnold, Falk, & Winther (2006, doi:10.1017/S0962492906210018).
1029a920c6SToby Isaac 
1129a920c6SToby Isaac    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
1229a920c6SToby Isaac    vectors from a vector space V and producing a real number:
1329a920c6SToby Isaac    - 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)
1429a920c6SToby Isaac    - 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)
15dce8aebaSBarry Smith    This action is implemented as `PetscDTAltVApply()`.
1629a920c6SToby Isaac 
1729a920c6SToby Isaac    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
1829a920c6SToby Isaac    shows that for an N dimensional space, only 0 <= k <= N are valid form degrees.)
1929a920c6SToby Isaac    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:
2029a920c6SToby Isaac    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.
2129a920c6SToby Isaac    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
2229a920c6SToby Isaac    by the associated subsets.
2329a920c6SToby Isaac 
2429a920c6SToby Isaac    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
2529a920c6SToby Isaac    V[i,j] = v_i[c_j] and taking the determinant of V.
2629a920c6SToby Isaac 
2729a920c6SToby Isaac    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).
2829a920c6SToby Isaac    This is an anticommutative product, (f wedge g) = -(g wedge f).  It is sufficient to describe the wedge product of two basis forms.
2929a920c6SToby Isaac    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):
3029a920c6SToby Isaac    - If there is any coordinate in both sets, then (f wedge g) = 0.
3129a920c6SToby Isaac    - Otherwise, (f wedge g) is a multiple of the basis (j+k)-form h associated with (f_1,...,f_j,g_1,...,g_k).
3229a920c6SToby Isaac    - 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.
33dce8aebaSBarry Smith    The wedge product is implemented for either two inputs (f and g) in `PetscDTAltVWedge()`, or for one (just f, giving a
34dce8aebaSBarry Smith    matrix to multiply against multiple choices of g) in `PetscDTAltVWedgeMatrix()`.
3529a920c6SToby Isaac 
3629a920c6SToby Isaac    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),
3729a920c6SToby Isaac    defined by (w int v)(v_1,...,v_{k-1}) = w(v,v_1,...,v_{k-1}).
3829a920c6SToby Isaac 
3929a920c6SToby Isaac    The interior product is implemented for either two inputs (w and v) in PetscDTAltVInterior, for one (just v, giving a
40dce8aebaSBarry Smith    matrix to multiply against multiple choices of w) in `PetscDTAltVInteriorMatrix()`,
41dce8aebaSBarry Smith    or for no inputs (giving the sparsity pattern of `PetscDTAltVInteriorMatrix()`) in `PetscDTAltVInteriorPattern()`.
4229a920c6SToby Isaac 
4329a920c6SToby Isaac    When there is a linear map L: V -> W from an N dimensional vector space to an M dimensional vector space,
4429a920c6SToby Isaac    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).
45dce8aebaSBarry Smith    The pullback is implemented as `PetscDTAltVPullback()` (acting on a known w) and `PetscDTAltVPullbackMatrix()` (creating a matrix that computes the actin of L^*).
4629a920c6SToby Isaac 
4729a920c6SToby Isaac    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
4829a920c6SToby Isaac    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
4929a920c6SToby Isaac    of the basis coefficients of f and g.
5029a920c6SToby Isaac    Powers of the Hodge star operator can be applied with PetscDTAltVStar
5129a920c6SToby Isaac 
526c877ef6SSatish Balay    Level: intermediate
5329a920c6SToby Isaac 
54db781477SPatrick Sanan .seealso: `PetscDTAltVApply()`, `PetscDTAltVWedge()`, `PetscDTAltVInterior()`, `PetscDTAltVPullback()`, `PetscDTAltVStar()`
5529a920c6SToby Isaac M*/
5629a920c6SToby Isaac 
57fad4db65SToby Isaac /*@
5828222859SToby Isaac   PetscDTAltVApply - Apply an a k-form (an alternating k-linear map) to a set of k N-dimensional vectors
59fad4db65SToby Isaac 
604165533cSJose E. Roman   Input Parameters:
6128222859SToby Isaac + N - the dimension of the vector space, N >= 0
6228222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N
63dce8aebaSBarry Smith . 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.
64dce8aebaSBarry Smith        The degrees of freedom are ordered lexicographically by their associated subsets)
6528222859SToby Isaac - v - a set of k vectors of size N, size [k x N], each vector stored contiguously
66fad4db65SToby Isaac 
674165533cSJose E. Roman   Output Parameter:
68dce8aebaSBarry Smith . 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
69dce8aebaSBarry Smith    entry (j,k) given by the s_{i,k}'th coordinate of v_j
70fad4db65SToby Isaac 
71fad4db65SToby Isaac   Level: intermediate
72fad4db65SToby Isaac 
73db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
74fad4db65SToby Isaac @*/
75d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv)
76d71ae5a4SJacob Faibussowitsch {
771a989b97SToby Isaac   PetscFunctionBegin;
7808401ef6SPierre Jolivet   PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
791dca8a05SBarry Smith   PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
801a989b97SToby Isaac   if (N <= 3) {
811a989b97SToby Isaac     if (!k) {
821a989b97SToby Isaac       *wv = w[0];
831a989b97SToby Isaac     } else {
849371c9d4SSatish Balay       if (N == 1) {
859371c9d4SSatish Balay         *wv = w[0] * v[0];
869371c9d4SSatish Balay       } else if (N == 2) {
879371c9d4SSatish Balay         if (k == 1) {
889371c9d4SSatish Balay           *wv = w[0] * v[0] + w[1] * v[1];
891a989b97SToby Isaac         } else {
909371c9d4SSatish Balay           *wv = w[0] * (v[0] * v[3] - v[1] * v[2]);
919371c9d4SSatish Balay         }
921a989b97SToby Isaac       } else {
939371c9d4SSatish Balay         if (k == 1) {
949371c9d4SSatish Balay           *wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];
959371c9d4SSatish Balay         } else if (k == 2) {
969371c9d4SSatish Balay           *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) + w[1] * (v[0] * v[5] - v[2] * v[3]) + w[2] * (v[1] * v[5] - v[2] * v[4]);
979371c9d4SSatish Balay         } else {
989371c9d4SSatish Balay           *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) + v[1] * (v[5] * v[6] - v[3] * v[8]) + v[2] * (v[3] * v[7] - v[4] * v[6]));
991a989b97SToby Isaac         }
1001a989b97SToby Isaac       }
1011a989b97SToby Isaac     }
1021a989b97SToby Isaac   } else {
1031a989b97SToby Isaac     PetscInt  Nk, Nf;
104fad4db65SToby Isaac     PetscInt *subset, *perm;
1051a989b97SToby Isaac     PetscInt  i, j, l;
1061a989b97SToby Isaac     PetscReal sum = 0.;
1071a989b97SToby Isaac 
1089566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(k, &Nf));
1099566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, k, &Nk));
1109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &perm));
1111a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
1121a989b97SToby Isaac       PetscReal subsum = 0.;
1131a989b97SToby Isaac 
1149566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
1151a989b97SToby Isaac       for (j = 0; j < Nf; j++) {
1161a989b97SToby Isaac         PetscBool permOdd;
1171a989b97SToby Isaac         PetscReal prod;
1181a989b97SToby Isaac 
1199566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumPerm(k, j, perm, &permOdd));
1201a989b97SToby Isaac         prod = permOdd ? -1. : 1.;
121ad540459SPierre Jolivet         for (l = 0; l < k; l++) prod *= v[perm[l] * N + subset[l]];
1221a989b97SToby Isaac         subsum += prod;
1231a989b97SToby Isaac       }
1241a989b97SToby Isaac       sum += w[i] * subsum;
1251a989b97SToby Isaac     }
1269566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, perm));
1271a989b97SToby Isaac     *wv = sum;
1281a989b97SToby Isaac   }
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1301a989b97SToby Isaac }
1311a989b97SToby Isaac 
132fad4db65SToby Isaac /*@
13328222859SToby Isaac   PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form
134fad4db65SToby Isaac 
1354165533cSJose E. Roman   Input Parameters:
13628222859SToby Isaac + N - the dimension of the vector space, N >= 0
13728222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N
13828222859SToby Isaac . k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N
13928222859SToby Isaac . a - a j-form, size [N choose j]
14028222859SToby Isaac - b - a k-form, size [N choose k]
141fad4db65SToby Isaac 
1424165533cSJose E. Roman   Output Parameter:
14328222859SToby Isaac . 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}}),
14428222859SToby Isaac              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}.
145fad4db65SToby Isaac 
146fad4db65SToby Isaac   Level: intermediate
147fad4db65SToby Isaac 
148db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVWedgeMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
149fad4db65SToby Isaac @*/
150d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb)
151d71ae5a4SJacob Faibussowitsch {
1521a989b97SToby Isaac   PetscInt i;
1531a989b97SToby Isaac 
1541a989b97SToby Isaac   PetscFunctionBegin;
15508401ef6SPierre Jolivet   PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
1561dca8a05SBarry Smith   PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
1571dca8a05SBarry Smith   PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
1581a989b97SToby Isaac   if (N <= 3) {
1591a989b97SToby Isaac     PetscInt Njk;
1601a989b97SToby Isaac 
1619566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
1629371c9d4SSatish Balay     if (!j) {
163ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeb[i] = a[0] * b[i];
1649371c9d4SSatish Balay     } else if (!k) {
165ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeb[i] = a[i] * b[0];
1669371c9d4SSatish Balay     } else {
1679371c9d4SSatish Balay       if (N == 2) {
1689371c9d4SSatish Balay         awedgeb[0] = a[0] * b[1] - a[1] * b[0];
1699371c9d4SSatish Balay       } else {
1701a989b97SToby Isaac         if (j + k == 2) {
1711a989b97SToby Isaac           awedgeb[0] = a[0] * b[1] - a[1] * b[0];
1721a989b97SToby Isaac           awedgeb[1] = a[0] * b[2] - a[2] * b[0];
1731a989b97SToby Isaac           awedgeb[2] = a[1] * b[2] - a[2] * b[1];
1741a989b97SToby Isaac         } else {
1751a989b97SToby Isaac           awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0];
1761a989b97SToby Isaac         }
1771a989b97SToby Isaac       }
1781a989b97SToby Isaac     }
1791a989b97SToby Isaac   } else {
1801a989b97SToby Isaac     PetscInt  Njk;
1811a989b97SToby Isaac     PetscInt  JKj;
1821a989b97SToby Isaac     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
1831a989b97SToby Isaac     PetscInt  i;
1841a989b97SToby Isaac 
1859566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
1869566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(j + k, j, &JKj));
1879566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk));
1881a989b97SToby Isaac     for (i = 0; i < Njk; i++) {
1891a989b97SToby Isaac       PetscReal sum = 0.;
1901a989b97SToby Isaac       PetscInt  l;
1911a989b97SToby Isaac 
1929566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, j + k, i, subset));
1931a989b97SToby Isaac       for (l = 0; l < JKj; l++) {
1941a989b97SToby Isaac         PetscBool jkOdd;
1951a989b97SToby Isaac         PetscInt  m, jInd, kInd;
1961a989b97SToby Isaac 
1979566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd));
198ad540459SPierre Jolivet         for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]];
199ad540459SPierre Jolivet         for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]];
2009566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd));
2019566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd));
2021a989b97SToby Isaac         sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]);
2031a989b97SToby Isaac       }
2041a989b97SToby Isaac       awedgeb[i] = sum;
2051a989b97SToby Isaac     }
2069566063dSJacob Faibussowitsch     PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk));
2071a989b97SToby Isaac   }
2083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2091a989b97SToby Isaac }
2101a989b97SToby Isaac 
211fad4db65SToby Isaac /*@
21228222859SToby Isaac   PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms
213fad4db65SToby Isaac 
2144165533cSJose E. Roman   Input Parameters:
21528222859SToby Isaac + N - the dimension of the vector space, N >= 0
21628222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N
21728222859SToby Isaac . k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N
21828222859SToby Isaac - a - a j-form, size [N choose j]
219fad4db65SToby Isaac 
2204165533cSJose E. Roman   Output Parameter:
22160225df5SJacob Faibussowitsch . awedgeMat - (a wedge), an [(N choose j+k) x (N choose k)] matrix in row-major order, such that (a wedge) * b = a wedge b
222fad4db65SToby Isaac 
223fad4db65SToby Isaac   Level: intermediate
224fad4db65SToby Isaac 
225db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
226fad4db65SToby Isaac @*/
227d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat)
228d71ae5a4SJacob Faibussowitsch {
2291a989b97SToby Isaac   PetscInt i;
2301a989b97SToby Isaac 
2311a989b97SToby Isaac   PetscFunctionBegin;
23208401ef6SPierre Jolivet   PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
2331dca8a05SBarry Smith   PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
2341dca8a05SBarry Smith   PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
2351a989b97SToby Isaac   if (N <= 3) {
2361a989b97SToby Isaac     PetscInt Njk;
2371a989b97SToby Isaac 
2389566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
2391a989b97SToby Isaac     if (!j) {
240ad540459SPierre Jolivet       for (i = 0; i < Njk * Njk; i++) awedgeMat[i] = 0.;
241ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeMat[i * (Njk + 1)] = a[0];
2421a989b97SToby Isaac     } else if (!k) {
243ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeMat[i] = a[i];
2441a989b97SToby Isaac     } else {
2451a989b97SToby Isaac       if (N == 2) {
2469371c9d4SSatish Balay         awedgeMat[0] = -a[1];
2479371c9d4SSatish Balay         awedgeMat[1] = a[0];
2481a989b97SToby Isaac       } else {
2491a989b97SToby Isaac         if (j + k == 2) {
2509371c9d4SSatish Balay           awedgeMat[0] = -a[1];
2519371c9d4SSatish Balay           awedgeMat[1] = a[0];
2529371c9d4SSatish Balay           awedgeMat[2] = 0.;
2539371c9d4SSatish Balay           awedgeMat[3] = -a[2];
2549371c9d4SSatish Balay           awedgeMat[4] = 0.;
2559371c9d4SSatish Balay           awedgeMat[5] = a[0];
2569371c9d4SSatish Balay           awedgeMat[6] = 0.;
2579371c9d4SSatish Balay           awedgeMat[7] = -a[2];
2589371c9d4SSatish Balay           awedgeMat[8] = a[1];
2591a989b97SToby Isaac         } else {
2609371c9d4SSatish Balay           awedgeMat[0] = a[2];
2619371c9d4SSatish Balay           awedgeMat[1] = -a[1];
2629371c9d4SSatish Balay           awedgeMat[2] = a[0];
2631a989b97SToby Isaac         }
2641a989b97SToby Isaac       }
2651a989b97SToby Isaac     }
2661a989b97SToby Isaac   } else {
2671a989b97SToby Isaac     PetscInt  Njk;
2681a989b97SToby Isaac     PetscInt  Nk;
2691a989b97SToby Isaac     PetscInt  JKj, i;
2701a989b97SToby Isaac     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
2711a989b97SToby Isaac 
2729566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, k, &Nk));
2739566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
2749566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(j + k, j, &JKj));
2759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk));
2761a989b97SToby Isaac     for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.;
2771a989b97SToby Isaac     for (i = 0; i < Njk; i++) {
2781a989b97SToby Isaac       PetscInt l;
2791a989b97SToby Isaac 
2809566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, j + k, i, subset));
2811a989b97SToby Isaac       for (l = 0; l < JKj; l++) {
2821a989b97SToby Isaac         PetscBool jkOdd;
2831a989b97SToby Isaac         PetscInt  m, jInd, kInd;
2841a989b97SToby Isaac 
2859566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd));
286ad540459SPierre Jolivet         for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]];
287ad540459SPierre Jolivet         for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]];
2889566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd));
2899566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd));
2901a989b97SToby Isaac         awedgeMat[i * Nk + kInd] += jkOdd ? -a[jInd] : a[jInd];
2911a989b97SToby Isaac       }
2921a989b97SToby Isaac     }
2939566063dSJacob Faibussowitsch     PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk));
2941a989b97SToby Isaac   }
2953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2961a989b97SToby Isaac }
2971a989b97SToby Isaac 
298fad4db65SToby Isaac /*@
29928222859SToby Isaac   PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space
300fad4db65SToby Isaac 
3014165533cSJose E. Roman   Input Parameters:
30228222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, M >= 0
30328222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, N >= 0
30428222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format
305dce8aebaSBarry Smith . 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
306dce8aebaSBarry Smith        the Hodge star operator (see note).
30728222859SToby Isaac - w - a |k|-form in the image space, size [M choose |k|]
308fad4db65SToby Isaac 
3094165533cSJose E. Roman   Output Parameter:
31028222859SToby Isaac . 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).
311fad4db65SToby Isaac 
312fad4db65SToby Isaac   Level: intermediate
313fad4db65SToby Isaac 
314*a4e35b19SJacob Faibussowitsch   Note:
315*a4e35b19SJacob Faibussowitsch   Negative form degrees accommodate, e.g., H-div conforming vector fields.  An H-div conforming
316*a4e35b19SJacob Faibussowitsch   vector field stores its degrees of freedom as (dx, dy, dz), like a 1-form, but its normal
317*a4e35b19SJacob Faibussowitsch   trace is integrated on faces, like a 2-form.  The correct pullback then is to apply the Hodge
318*a4e35b19SJacob Faibussowitsch   star transformation from (M-2)-form to 2-form, pullback as a 2-form, then invert the Hodge
319*a4e35b19SJacob Faibussowitsch   star transformation.
320fad4db65SToby Isaac 
321db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullbackMatrix()`, `PetscDTAltVStar()`
322fad4db65SToby Isaac @*/
323d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw)
324d71ae5a4SJacob Faibussowitsch {
3251a989b97SToby Isaac   PetscInt i, j, Nk, Mk;
3261a989b97SToby Isaac 
3271a989b97SToby Isaac   PetscFunctionBegin;
3281dca8a05SBarry Smith   PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
3291dca8a05SBarry Smith   PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
3301a989b97SToby Isaac   if (N <= 3 && M <= 3) {
3319566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
3329566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
3331a989b97SToby Isaac     if (!k) {
3341a989b97SToby Isaac       Lstarw[0] = w[0];
3351a989b97SToby Isaac     } else if (k == 1) {
3361a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3371a989b97SToby Isaac         PetscReal sum = 0.;
3381a989b97SToby Isaac 
339ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) sum += L[j * Nk + i] * w[j];
3401a989b97SToby Isaac         Lstarw[i] = sum;
3411a989b97SToby Isaac       }
3421a989b97SToby Isaac     } else if (k == -1) {
3431a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
3441a989b97SToby Isaac 
3451a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3461a989b97SToby Isaac         PetscReal sum = 0.;
3471a989b97SToby Isaac 
348ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j];
3491a989b97SToby Isaac         Lstarw[i] = mult[i] * sum;
3501a989b97SToby Isaac       }
3511a989b97SToby Isaac     } else if (k == 2) {
3529371c9d4SSatish Balay       PetscInt pairs[3][2] = {
3539371c9d4SSatish Balay         {0, 1},
3549371c9d4SSatish Balay         {0, 2},
3559371c9d4SSatish Balay         {1, 2}
3569371c9d4SSatish Balay       };
3571a989b97SToby Isaac 
3581a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3591a989b97SToby Isaac         PetscReal sum = 0.;
360ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] - L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j];
3611a989b97SToby Isaac         Lstarw[i] = sum;
3621a989b97SToby Isaac       }
3631a989b97SToby Isaac     } else if (k == -2) {
3649371c9d4SSatish Balay       PetscInt pairs[3][2] = {
3659371c9d4SSatish Balay         {1, 2},
3669371c9d4SSatish Balay         {2, 0},
3679371c9d4SSatish Balay         {0, 1}
3689371c9d4SSatish Balay       };
3691a989b97SToby Isaac       PetscInt offi = (N == 2) ? 2 : 0;
3701a989b97SToby Isaac       PetscInt offj = (M == 2) ? 2 : 0;
3711a989b97SToby Isaac 
3721a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3731a989b97SToby Isaac         PetscReal sum = 0.;
3741a989b97SToby Isaac 
375ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] * L[pairs[offj + j][1] * N + pairs[offi + i][1]] - L[pairs[offj + j][1] * N + pairs[offi + i][0]] * L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j];
3761a989b97SToby Isaac         Lstarw[i] = sum;
3771a989b97SToby Isaac       }
3781a989b97SToby Isaac     } else {
3799371c9d4SSatish Balay       PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + L[1] * (L[5] * L[6] - L[3] * L[8]) + L[2] * (L[3] * L[7] - L[4] * L[6]);
3801a989b97SToby Isaac 
381ad540459SPierre Jolivet       for (i = 0; i < Nk; i++) Lstarw[i] = detL * w[i];
3821a989b97SToby Isaac     }
3831a989b97SToby Isaac   } else {
3841a989b97SToby Isaac     PetscInt         Nf, l, p;
3851a989b97SToby Isaac     PetscReal       *Lw, *Lwv;
3861a989b97SToby Isaac     PetscInt        *subsetw, *subsetv;
387fad4db65SToby Isaac     PetscInt        *perm;
3881a989b97SToby Isaac     PetscReal       *walloc   = NULL;
3891a989b97SToby Isaac     const PetscReal *ww       = NULL;
3901a989b97SToby Isaac     PetscBool        negative = PETSC_FALSE;
3911a989b97SToby Isaac 
3929566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
3939566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
3949566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf));
3951a989b97SToby Isaac     if (k < 0) {
3961a989b97SToby Isaac       negative = PETSC_TRUE;
3971a989b97SToby Isaac       k        = -k;
3989566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Mk, &walloc));
3999566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc));
4001a989b97SToby Isaac       ww = walloc;
4011a989b97SToby Isaac     } else {
4021a989b97SToby Isaac       ww = w;
4031a989b97SToby Isaac     }
4049566063dSJacob Faibussowitsch     PetscCall(PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv));
4051a989b97SToby Isaac     for (i = 0; i < Nk; i++) Lstarw[i] = 0.;
4061a989b97SToby Isaac     for (i = 0; i < Mk; i++) {
4079566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(M, k, i, subsetw));
4081a989b97SToby Isaac       for (j = 0; j < Nk; j++) {
4099566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSubset(N, k, j, subsetv));
4101a989b97SToby Isaac         for (p = 0; p < Nf; p++) {
4111a989b97SToby Isaac           PetscReal prod;
4121a989b97SToby Isaac           PetscBool isOdd;
4131a989b97SToby Isaac 
4149566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd));
4151a989b97SToby Isaac           prod = isOdd ? -ww[i] : ww[i];
416ad540459SPierre Jolivet           for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]];
4171a989b97SToby Isaac           Lstarw[j] += prod;
4181a989b97SToby Isaac         }
4191a989b97SToby Isaac       }
4201a989b97SToby Isaac     }
4211a989b97SToby Isaac     if (negative) {
4221a989b97SToby Isaac       PetscReal *sLsw;
4231a989b97SToby Isaac 
4249566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &sLsw));
4259566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw));
4261a989b97SToby Isaac       for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i];
4279566063dSJacob Faibussowitsch       PetscCall(PetscFree(sLsw));
4281a989b97SToby Isaac     }
4299566063dSJacob Faibussowitsch     PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv));
4309566063dSJacob Faibussowitsch     PetscCall(PetscFree(walloc));
4311a989b97SToby Isaac   }
4323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4331a989b97SToby Isaac }
4341a989b97SToby Isaac 
435fad4db65SToby Isaac /*@
436fad4db65SToby Isaac   PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation
437fad4db65SToby Isaac 
4384165533cSJose E. Roman   Input Parameters:
43928222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, N >= 0
44028222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, M >= 0
44128222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format
442dce8aebaSBarry Smith - k - the *signed* degree k of the |k|-forms on which Lstar acts, -(min(M,N)) <= k <= min(M,N).
443dce8aebaSBarry Smith        A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note in `PetscDTAltvPullback()`)
444fad4db65SToby Isaac 
4454165533cSJose E. Roman   Output Parameter:
44628222859SToby Isaac . Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w
447fad4db65SToby Isaac 
448fad4db65SToby Isaac   Level: intermediate
449fad4db65SToby Isaac 
450db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVStar()`
451fad4db65SToby Isaac @*/
452d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar)
453d71ae5a4SJacob Faibussowitsch {
4541a989b97SToby Isaac   PetscInt   Nk, Mk, Nf, i, j, l, p;
4551a989b97SToby Isaac   PetscReal *Lw, *Lwv;
4561a989b97SToby Isaac   PetscInt  *subsetw, *subsetv;
457fad4db65SToby Isaac   PetscInt  *perm;
4581a989b97SToby Isaac   PetscBool  negative = PETSC_FALSE;
4591a989b97SToby Isaac 
4601a989b97SToby Isaac   PetscFunctionBegin;
4611dca8a05SBarry Smith   PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
4621dca8a05SBarry Smith   PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
4631a989b97SToby Isaac   if (N <= 3 && M <= 3) {
4641a989b97SToby Isaac     PetscReal mult[3] = {1., -1., 1.};
4651a989b97SToby Isaac 
4669566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
4679566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
4681a989b97SToby Isaac     if (!k) {
4691a989b97SToby Isaac       Lstar[0] = 1.;
4701a989b97SToby Isaac     } else if (k == 1) {
4719371c9d4SSatish Balay       for (i = 0; i < Nk; i++) {
472ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[j * Nk + i];
4739371c9d4SSatish Balay       }
4741a989b97SToby Isaac     } else if (k == -1) {
4751a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
476ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j];
4771a989b97SToby Isaac       }
4781a989b97SToby Isaac     } else if (k == 2) {
4799371c9d4SSatish Balay       PetscInt pairs[3][2] = {
4809371c9d4SSatish Balay         {0, 1},
4819371c9d4SSatish Balay         {0, 2},
4829371c9d4SSatish Balay         {1, 2}
4839371c9d4SSatish Balay       };
4841a989b97SToby Isaac 
4851a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
486ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] - L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]];
4871a989b97SToby Isaac       }
4881a989b97SToby Isaac     } else if (k == -2) {
4899371c9d4SSatish Balay       PetscInt pairs[3][2] = {
4909371c9d4SSatish Balay         {1, 2},
4919371c9d4SSatish Balay         {2, 0},
4929371c9d4SSatish Balay         {0, 1}
4939371c9d4SSatish Balay       };
4941a989b97SToby Isaac       PetscInt offi = (N == 2) ? 2 : 0;
4951a989b97SToby Isaac       PetscInt offj = (M == 2) ? 2 : 0;
4961a989b97SToby Isaac 
4971a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
4981a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
4999371c9d4SSatish Balay           Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] * L[pairs[offj + j][1] * N + pairs[offi + i][1]] - L[pairs[offj + j][1] * N + pairs[offi + i][0]] * L[pairs[offj + j][0] * N + pairs[offi + i][1]];
5001a989b97SToby Isaac         }
5011a989b97SToby Isaac       }
5021a989b97SToby Isaac     } else {
5039371c9d4SSatish Balay       PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + L[1] * (L[5] * L[6] - L[3] * L[8]) + L[2] * (L[3] * L[7] - L[4] * L[6]);
5041a989b97SToby Isaac 
505ad540459SPierre Jolivet       for (i = 0; i < Nk; i++) Lstar[i] = detL;
5061a989b97SToby Isaac     }
5071a989b97SToby Isaac   } else {
5081a989b97SToby Isaac     if (k < 0) {
5091a989b97SToby Isaac       negative = PETSC_TRUE;
5101a989b97SToby Isaac       k        = -k;
5111a989b97SToby Isaac     }
5129566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
5139566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
5149566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf));
5159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv));
5161a989b97SToby Isaac     for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.;
5171a989b97SToby Isaac     for (i = 0; i < Mk; i++) {
5181a989b97SToby Isaac       PetscBool iOdd;
5191a989b97SToby Isaac       PetscInt  iidx, jidx;
5201a989b97SToby Isaac 
5219566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSplit(M, k, i, subsetw, &iOdd));
5221a989b97SToby Isaac       iidx = negative ? Mk - 1 - i : i;
52328222859SToby Isaac       iOdd = negative ? (PetscBool)(iOdd ^ ((k * (M - k)) & 1)) : PETSC_FALSE;
5241a989b97SToby Isaac       for (j = 0; j < Nk; j++) {
5251a989b97SToby Isaac         PetscBool jOdd;
5261a989b97SToby Isaac 
5279566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, k, j, subsetv, &jOdd));
5281a989b97SToby Isaac         jidx = negative ? Nk - 1 - j : j;
52928222859SToby Isaac         jOdd = negative ? (PetscBool)(iOdd ^ jOdd ^ ((k * (N - k)) & 1)) : PETSC_FALSE;
5301a989b97SToby Isaac         for (p = 0; p < Nf; p++) {
5311a989b97SToby Isaac           PetscReal prod;
5321a989b97SToby Isaac           PetscBool isOdd;
5331a989b97SToby Isaac 
5349566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd));
53528222859SToby Isaac           isOdd = (PetscBool)(isOdd ^ jOdd);
5361a989b97SToby Isaac           prod  = isOdd ? -1. : 1.;
537ad540459SPierre Jolivet           for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]];
5381a989b97SToby Isaac           Lstar[jidx * Mk + iidx] += prod;
5391a989b97SToby Isaac         }
5401a989b97SToby Isaac       }
5411a989b97SToby Isaac     }
5429566063dSJacob Faibussowitsch     PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv));
5431a989b97SToby Isaac   }
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5451a989b97SToby Isaac }
5461a989b97SToby Isaac 
547fad4db65SToby Isaac /*@
54828222859SToby Isaac   PetscDTAltVInterior - Compute the interior product of a k-form with a vector
549fad4db65SToby Isaac 
5504165533cSJose E. Roman   Input Parameters:
55128222859SToby Isaac + N - the dimension of the vector space, N >= 0
55228222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N
55328222859SToby Isaac . w - a k-form, size [N choose k]
55428222859SToby Isaac - v - an N dimensional vector
555fad4db65SToby Isaac 
5564165533cSJose E. Roman   Output Parameter:
55728222859SToby Isaac . 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}).
558fad4db65SToby Isaac 
559fad4db65SToby Isaac   Level: intermediate
560fad4db65SToby Isaac 
561db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
562fad4db65SToby Isaac @*/
563d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv)
564d71ae5a4SJacob Faibussowitsch {
5651a989b97SToby Isaac   PetscInt i, Nk, Nkm;
5661a989b97SToby Isaac 
5671a989b97SToby Isaac   PetscFunctionBegin;
5681dca8a05SBarry Smith   PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
5699566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
5709566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
5711a989b97SToby Isaac   if (N <= 3) {
5721a989b97SToby Isaac     if (k == 1) {
5731a989b97SToby Isaac       PetscReal sum = 0.;
5741a989b97SToby Isaac 
575ad540459SPierre Jolivet       for (i = 0; i < N; i++) sum += w[i] * v[i];
5761a989b97SToby Isaac       wIntv[0] = sum;
5771a989b97SToby Isaac     } else if (k == N) {
5781a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
5791a989b97SToby Isaac 
580ad540459SPierre Jolivet       for (i = 0; i < N; i++) wIntv[N - 1 - i] = w[0] * v[i] * mult[i];
5811a989b97SToby Isaac     } else {
5821a989b97SToby Isaac       wIntv[0] = -w[0] * v[1] - w[1] * v[2];
5831a989b97SToby Isaac       wIntv[1] = w[0] * v[0] - w[2] * v[2];
5841a989b97SToby Isaac       wIntv[2] = w[1] * v[0] + w[2] * v[1];
5851a989b97SToby Isaac     }
5861a989b97SToby Isaac   } else {
5871a989b97SToby Isaac     PetscInt *subset, *work;
5881a989b97SToby Isaac 
5899566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
5901a989b97SToby Isaac     for (i = 0; i < Nkm; i++) wIntv[i] = 0.;
5911a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
5921a989b97SToby Isaac       PetscInt j, l, m;
5931a989b97SToby Isaac 
5949566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
5951a989b97SToby Isaac       for (j = 0; j < k; j++) {
5961a989b97SToby Isaac         PetscInt  idx;
59728222859SToby Isaac         PetscBool flip = (PetscBool)(j & 1);
5981a989b97SToby Isaac 
5991a989b97SToby Isaac         for (l = 0, m = 0; l < k; l++) {
6001a989b97SToby Isaac           if (l != j) work[m++] = subset[l];
6011a989b97SToby Isaac         }
6029566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
6031a989b97SToby Isaac         wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]);
6041a989b97SToby Isaac       }
6051a989b97SToby Isaac     }
6069566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
6071a989b97SToby Isaac   }
6083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6091a989b97SToby Isaac }
6101a989b97SToby Isaac 
611fad4db65SToby Isaac /*@
61228222859SToby Isaac   PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector
613fad4db65SToby Isaac 
6144165533cSJose E. Roman   Input Parameters:
61528222859SToby Isaac + N - the dimension of the vector space, N >= 0
61628222859SToby Isaac . k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N
61728222859SToby Isaac - v - an N dimensional vector
618fad4db65SToby Isaac 
6194165533cSJose E. Roman   Output Parameter:
620fad4db65SToby Isaac . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v)
621fad4db65SToby Isaac 
622fad4db65SToby Isaac   Level: intermediate
623fad4db65SToby Isaac 
624db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
625fad4db65SToby Isaac @*/
626d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat)
627d71ae5a4SJacob Faibussowitsch {
6281a989b97SToby Isaac   PetscInt i, Nk, Nkm;
6291a989b97SToby Isaac 
6301a989b97SToby Isaac   PetscFunctionBegin;
6311dca8a05SBarry Smith   PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
6329566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
6339566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
6341a989b97SToby Isaac   if (N <= 3) {
6351a989b97SToby Isaac     if (k == 1) {
6361a989b97SToby Isaac       for (i = 0; i < N; i++) intvMat[i] = v[i];
6371a989b97SToby Isaac     } else if (k == N) {
6381a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
6391a989b97SToby Isaac 
6401a989b97SToby Isaac       for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i];
6411a989b97SToby Isaac     } else {
6429371c9d4SSatish Balay       intvMat[0] = -v[1];
6439371c9d4SSatish Balay       intvMat[1] = -v[2];
6449371c9d4SSatish Balay       intvMat[2] = 0.;
6459371c9d4SSatish Balay       intvMat[3] = v[0];
6469371c9d4SSatish Balay       intvMat[4] = 0.;
6479371c9d4SSatish Balay       intvMat[5] = -v[2];
6489371c9d4SSatish Balay       intvMat[6] = 0.;
6499371c9d4SSatish Balay       intvMat[7] = v[0];
6509371c9d4SSatish Balay       intvMat[8] = v[1];
6511a989b97SToby Isaac     }
6521a989b97SToby Isaac   } else {
6531a989b97SToby Isaac     PetscInt *subset, *work;
6541a989b97SToby Isaac 
6559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
6561a989b97SToby Isaac     for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.;
6571a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
6581a989b97SToby Isaac       PetscInt j, l, m;
6591a989b97SToby Isaac 
6609566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
6611a989b97SToby Isaac       for (j = 0; j < k; j++) {
6621a989b97SToby Isaac         PetscInt  idx;
66328222859SToby Isaac         PetscBool flip = (PetscBool)(j & 1);
6641a989b97SToby Isaac 
6651a989b97SToby Isaac         for (l = 0, m = 0; l < k; l++) {
6661a989b97SToby Isaac           if (l != j) work[m++] = subset[l];
6671a989b97SToby Isaac         }
6689566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
6691a989b97SToby Isaac         intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]];
6701a989b97SToby Isaac       }
6711a989b97SToby Isaac     }
6729566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
6731a989b97SToby Isaac   }
6743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6751a989b97SToby Isaac }
6761a989b97SToby Isaac 
677fad4db65SToby Isaac /*@
678dce8aebaSBarry Smith   PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in `PetscDTAltVInteriorMatrix()`
679fad4db65SToby Isaac 
6804165533cSJose E. Roman   Input Parameters:
68128222859SToby Isaac + N - the dimension of the vector space, N >= 0
682dce8aebaSBarry Smith - k - the degree of the k-forms on which intvMat from `PetscDTAltVInteriorMatrix()` acts, 0 <= k <= N.
683fad4db65SToby Isaac 
6844165533cSJose E. Roman   Output Parameter:
68528222859SToby Isaac . indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k
68628222859SToby Isaac              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
68728222859SToby Isaac              coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1)
688fad4db65SToby Isaac 
689fad4db65SToby Isaac   Level: intermediate
690fad4db65SToby Isaac 
691dce8aebaSBarry Smith   Note:
692dce8aebaSBarry Smith   This function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential
693fad4db65SToby Isaac 
694db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
695fad4db65SToby Isaac @*/
696d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3])
697d71ae5a4SJacob Faibussowitsch {
698dda711d0SToby Isaac   PetscInt i, Nk, Nkm;
699dda711d0SToby Isaac 
700dda711d0SToby Isaac   PetscFunctionBegin;
7011dca8a05SBarry Smith   PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
7029566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
7039566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
704dda711d0SToby Isaac   if (N <= 3) {
705dda711d0SToby Isaac     if (k == 1) {
706dda711d0SToby Isaac       for (i = 0; i < N; i++) {
707dda711d0SToby Isaac         indices[i][0] = 0;
708dda711d0SToby Isaac         indices[i][1] = i;
709dda711d0SToby Isaac         indices[i][2] = i;
710dda711d0SToby Isaac       }
711dda711d0SToby Isaac     } else if (k == N) {
712dda711d0SToby Isaac       PetscInt val[3] = {0, -2, 2};
713dda711d0SToby Isaac 
714dda711d0SToby Isaac       for (i = 0; i < N; i++) {
715dda711d0SToby Isaac         indices[i][0] = N - 1 - i;
716dda711d0SToby Isaac         indices[i][1] = 0;
717dda711d0SToby Isaac         indices[i][2] = val[i];
718dda711d0SToby Isaac       }
719dda711d0SToby Isaac     } else {
7209371c9d4SSatish Balay       indices[0][0] = 0;
7219371c9d4SSatish Balay       indices[0][1] = 0;
7229371c9d4SSatish Balay       indices[0][2] = -(1 + 1);
7239371c9d4SSatish Balay       indices[1][0] = 0;
7249371c9d4SSatish Balay       indices[1][1] = 1;
7259371c9d4SSatish Balay       indices[1][2] = -(2 + 1);
7269371c9d4SSatish Balay       indices[2][0] = 1;
7279371c9d4SSatish Balay       indices[2][1] = 0;
7289371c9d4SSatish Balay       indices[2][2] = 0;
7299371c9d4SSatish Balay       indices[3][0] = 1;
7309371c9d4SSatish Balay       indices[3][1] = 2;
7319371c9d4SSatish Balay       indices[3][2] = -(2 + 1);
7329371c9d4SSatish Balay       indices[4][0] = 2;
7339371c9d4SSatish Balay       indices[4][1] = 1;
7349371c9d4SSatish Balay       indices[4][2] = 0;
7359371c9d4SSatish Balay       indices[5][0] = 2;
7369371c9d4SSatish Balay       indices[5][1] = 2;
7379371c9d4SSatish Balay       indices[5][2] = 1;
738dda711d0SToby Isaac     }
739dda711d0SToby Isaac   } else {
740dda711d0SToby Isaac     PetscInt *subset, *work;
741dda711d0SToby Isaac 
7429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
743dda711d0SToby Isaac     for (i = 0; i < Nk; i++) {
744dda711d0SToby Isaac       PetscInt j, l, m;
745dda711d0SToby Isaac 
7469566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
747dda711d0SToby Isaac       for (j = 0; j < k; j++) {
748dda711d0SToby Isaac         PetscInt  idx;
74928222859SToby Isaac         PetscBool flip = (PetscBool)(j & 1);
750dda711d0SToby Isaac 
751dda711d0SToby Isaac         for (l = 0, m = 0; l < k; l++) {
752dda711d0SToby Isaac           if (l != j) work[m++] = subset[l];
753dda711d0SToby Isaac         }
7549566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
755dda711d0SToby Isaac         indices[i * k + j][0] = idx;
756dda711d0SToby Isaac         indices[i * k + j][1] = i;
757dda711d0SToby Isaac         indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j];
758dda711d0SToby Isaac       }
759dda711d0SToby Isaac     }
7609566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
761dda711d0SToby Isaac   }
7623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
763dda711d0SToby Isaac }
764dda711d0SToby Isaac 
765fad4db65SToby Isaac /*@
76628222859SToby Isaac   PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form
767fad4db65SToby Isaac 
7684165533cSJose E. Roman   Input Parameters:
76928222859SToby Isaac + N   - the dimension of the vector space, N >= 0
77028222859SToby Isaac . k   - the degree k of the k-form w, 0 <= k <= N
77128222859SToby Isaac . 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.
77228222859SToby Isaac - w   - a k-form, size [N choose k]
773fad4db65SToby Isaac 
7744165533cSJose E. Roman   Output Parameter:
775*a4e35b19SJacob Faibussowitsch . starw - (star)^pow w
776fad4db65SToby Isaac 
777fad4db65SToby Isaac   Level: intermediate
778fad4db65SToby Isaac 
779*a4e35b19SJacob Faibussowitsch   Notes:
780*a4e35b19SJacob Faibussowitsch   Each degree of freedom of a k-form is associated with a subset S of k coordinates of the N
781*a4e35b19SJacob Faibussowitsch   dimensional vector space: the Hodge start operator (star) maps that degree of freedom to the
782*a4e35b19SJacob Faibussowitsch   degree of freedom associated with S', the complement of S, with a sign change if the
783*a4e35b19SJacob Faibussowitsch   permutation of coordinates {S[0], ... S[k-1], S'[0], starw- 1]} is an odd permutation.  This
784*a4e35b19SJacob Faibussowitsch   implies (star)^2 w = (-1)^{k(N-k)} w, and (star)^4 w = w.
785*a4e35b19SJacob Faibussowitsch 
786db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
787fad4db65SToby Isaac @*/
788d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw)
789d71ae5a4SJacob Faibussowitsch {
7901a989b97SToby Isaac   PetscInt Nk, i;
7911a989b97SToby Isaac 
7921a989b97SToby Isaac   PetscFunctionBegin;
7931dca8a05SBarry Smith   PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
7949566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
7951a989b97SToby Isaac   pow = pow % 4;
7961a989b97SToby Isaac   pow = (pow + 4) % 4; /* make non-negative */
7971a989b97SToby Isaac   /* pow is now 0, 1, 2, 3 */
7981a989b97SToby Isaac   if (N <= 3) {
7991a989b97SToby Isaac     if (pow & 1) {
8001a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
8011a989b97SToby Isaac 
8021a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i];
8031a989b97SToby Isaac     } else {
8041a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = w[i];
8051a989b97SToby Isaac     }
8061a989b97SToby Isaac     if (pow > 1 && ((k * (N - k)) & 1)) {
8071a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
8081a989b97SToby Isaac     }
8091a989b97SToby Isaac   } else {
8101a989b97SToby Isaac     PetscInt *subset;
8111a989b97SToby Isaac 
8129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N, &subset));
8131a989b97SToby Isaac     if (pow % 2) {
8141a989b97SToby Isaac       PetscInt l = (pow == 1) ? k : N - k;
8151a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
8161a989b97SToby Isaac         PetscBool sOdd;
8171a989b97SToby Isaac         PetscInt  j, idx;
8181a989b97SToby Isaac 
8199566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, l, i, subset, &sOdd));
8209566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, l, subset, &idx));
8219566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, N - l, &subset[l], &j));
8221a989b97SToby Isaac         starw[j] = sOdd ? -w[idx] : w[idx];
8231a989b97SToby Isaac       }
8241a989b97SToby Isaac     } else {
8251a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = w[i];
8261a989b97SToby Isaac     }
8271a989b97SToby Isaac     /* star^2 = -1^(k * (N - k)) */
8281a989b97SToby Isaac     if (pow > 1 && (k * (N - k)) % 2) {
8291a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
8301a989b97SToby Isaac     }
8319566063dSJacob Faibussowitsch     PetscCall(PetscFree(subset));
8321a989b97SToby Isaac   }
8333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8341a989b97SToby Isaac }
835