xref: /petsc/src/dm/dt/interface/dtaltv.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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)
1529a920c6SToby Isaac    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.
3329a920c6SToby Isaac    The wedge product is implemented for either two inputs (f and g) in PetscDTAltVWedge, or for one (just f, giving a
3429a920c6SToby Isaac    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
4029a920c6SToby Isaac    matrix to multiply against multiple choices of w) in PetscDTAltVInteriorMatrix,
4129a920c6SToby Isaac    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).
4529a920c6SToby Isaac    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 
5429a920c6SToby Isaac .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
6328222859SToby Isaac .  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)
6428222859SToby Isaac -  v - a set of k vectors of size N, size [k x N], each vector stored contiguously
65fad4db65SToby Isaac 
664165533cSJose E. Roman    Output Parameter:
6728222859SToby Isaac .  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
68fad4db65SToby Isaac 
69fad4db65SToby Isaac    Level: intermediate
70fad4db65SToby Isaac 
7129a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
72fad4db65SToby Isaac @*/
731a989b97SToby Isaac PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv)
741a989b97SToby Isaac {
751a989b97SToby Isaac   PetscFunctionBegin;
76*08401ef6SPierre Jolivet   PetscCheck(N >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
772c71b3e2SJacob Faibussowitsch   PetscCheckFalse(k < 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
781a989b97SToby Isaac   if (N <= 3) {
791a989b97SToby Isaac     if (!k) {
801a989b97SToby Isaac       *wv = w[0];
811a989b97SToby Isaac     } else {
821a989b97SToby Isaac       if (N == 1)        {*wv = w[0] * v[0];}
831a989b97SToby Isaac       else if (N == 2) {
841a989b97SToby Isaac         if (k == 1)      {*wv = w[0] * v[0] + w[1] * v[1];}
851a989b97SToby Isaac         else             {*wv = w[0] * (v[0] * v[3] - v[1] * v[2]);}
861a989b97SToby Isaac       } else {
871a989b97SToby Isaac         if (k == 1)      {*wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];}
881a989b97SToby Isaac         else if (k == 2) {
891a989b97SToby Isaac           *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) +
901a989b97SToby Isaac                 w[1] * (v[0] * v[5] - v[2] * v[3]) +
911a989b97SToby Isaac                 w[2] * (v[1] * v[5] - v[2] * v[4]);
921a989b97SToby Isaac         } else {
931a989b97SToby Isaac           *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) +
941a989b97SToby Isaac                         v[1] * (v[5] * v[6] - v[3] * v[8]) +
951a989b97SToby Isaac                         v[2] * (v[3] * v[7] - v[4] * v[6]));
961a989b97SToby Isaac         }
971a989b97SToby Isaac       }
981a989b97SToby Isaac     }
991a989b97SToby Isaac   } else {
1001a989b97SToby Isaac     PetscInt Nk, Nf;
101fad4db65SToby Isaac     PetscInt *subset, *perm;
1021a989b97SToby Isaac     PetscInt i, j, l;
1031a989b97SToby Isaac     PetscReal sum = 0.;
1041a989b97SToby Isaac 
1059566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(k, &Nf));
1069566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, k, &Nk));
1079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &perm));
1081a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
1091a989b97SToby Isaac       PetscReal subsum = 0.;
1101a989b97SToby Isaac 
1119566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
1121a989b97SToby Isaac       for (j = 0; j < Nf; j++) {
1131a989b97SToby Isaac         PetscBool permOdd;
1141a989b97SToby Isaac         PetscReal prod;
1151a989b97SToby Isaac 
1169566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumPerm(k, j, perm, &permOdd));
1171a989b97SToby Isaac         prod = permOdd ? -1. : 1.;
1181a989b97SToby Isaac         for (l = 0; l < k; l++) {
1191a989b97SToby Isaac           prod *= v[perm[l] * N + subset[l]];
1201a989b97SToby Isaac         }
1211a989b97SToby Isaac         subsum += prod;
1221a989b97SToby Isaac       }
1231a989b97SToby Isaac       sum += w[i] * subsum;
1241a989b97SToby Isaac     }
1259566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, perm));
1261a989b97SToby Isaac     *wv = sum;
1271a989b97SToby Isaac   }
1281a989b97SToby Isaac   PetscFunctionReturn(0);
1291a989b97SToby Isaac }
1301a989b97SToby Isaac 
131fad4db65SToby Isaac /*@
13228222859SToby Isaac    PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form
133fad4db65SToby Isaac 
1344165533cSJose E. Roman    Input Parameters:
13528222859SToby Isaac +  N - the dimension of the vector space, N >= 0
13628222859SToby Isaac .  j - the degree j of the j-form a, 0 <= j <= N
13728222859SToby Isaac .  k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N
13828222859SToby Isaac .  a - a j-form, size [N choose j]
13928222859SToby Isaac -  b - a k-form, size [N choose k]
140fad4db65SToby Isaac 
1414165533cSJose E. Roman    Output Parameter:
14228222859SToby 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}}),
14328222859SToby 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}.
144fad4db65SToby Isaac 
145fad4db65SToby Isaac    Level: intermediate
146fad4db65SToby Isaac 
14729a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVWedgeMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
148fad4db65SToby Isaac @*/
1491a989b97SToby Isaac PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb)
1501a989b97SToby Isaac {
1511a989b97SToby Isaac   PetscInt       i;
1521a989b97SToby Isaac 
1531a989b97SToby Isaac   PetscFunctionBegin;
154*08401ef6SPierre Jolivet   PetscCheck(N >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
1552c71b3e2SJacob Faibussowitsch   PetscCheckFalse(j < 0 || k < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
1562c71b3e2SJacob Faibussowitsch   PetscCheckFalse(j + k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
1571a989b97SToby Isaac   if (N <= 3) {
1581a989b97SToby Isaac     PetscInt Njk;
1591a989b97SToby Isaac 
1609566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j+k, &Njk));
1611a989b97SToby Isaac     if (!j)      {for (i = 0; i < Njk; i++) {awedgeb[i] = a[0] * b[i];}}
1621a989b97SToby Isaac     else if (!k) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[i] * b[0];}}
1631a989b97SToby Isaac     else {
1641a989b97SToby Isaac       if (N == 2) {awedgeb[0] = a[0] * b[1] - a[1] * b[0];}
1651a989b97SToby Isaac       else {
1661a989b97SToby Isaac         if (j+k == 2) {
1671a989b97SToby Isaac           awedgeb[0] = a[0] * b[1] - a[1] * b[0];
1681a989b97SToby Isaac           awedgeb[1] = a[0] * b[2] - a[2] * b[0];
1691a989b97SToby Isaac           awedgeb[2] = a[1] * b[2] - a[2] * b[1];
1701a989b97SToby Isaac         } else {
1711a989b97SToby Isaac           awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0];
1721a989b97SToby Isaac         }
1731a989b97SToby Isaac       }
1741a989b97SToby Isaac     }
1751a989b97SToby Isaac   } else {
1761a989b97SToby Isaac     PetscInt  Njk;
1771a989b97SToby Isaac     PetscInt  JKj;
1781a989b97SToby Isaac     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
1791a989b97SToby Isaac     PetscInt  i;
1801a989b97SToby Isaac 
1819566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j+k, &Njk));
1829566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(j+k, j, &JKj));
1839566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk));
1841a989b97SToby Isaac     for (i = 0; i < Njk; i++) {
1851a989b97SToby Isaac       PetscReal sum = 0.;
1861a989b97SToby Isaac       PetscInt  l;
1871a989b97SToby Isaac 
1889566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, j+k, i, subset));
1891a989b97SToby Isaac       for (l = 0; l < JKj; l++) {
1901a989b97SToby Isaac         PetscBool jkOdd;
1911a989b97SToby Isaac         PetscInt  m, jInd, kInd;
1921a989b97SToby Isaac 
1939566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd));
1941a989b97SToby Isaac         for (m = 0; m < j; m++) {
1951a989b97SToby Isaac           subsetj[m] = subset[subsetjk[m]];
1961a989b97SToby Isaac         }
1971a989b97SToby Isaac         for (m = 0; m < k; m++) {
1981a989b97SToby Isaac           subsetk[m] = subset[subsetjk[j+m]];
1991a989b97SToby Isaac         }
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   }
2081a989b97SToby Isaac   PetscFunctionReturn(0);
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:
22128222859SToby Isaac .  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
222fad4db65SToby Isaac 
223fad4db65SToby Isaac    Level: intermediate
224fad4db65SToby Isaac 
22529a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
226fad4db65SToby Isaac @*/
2271a989b97SToby Isaac PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat)
2281a989b97SToby Isaac {
2291a989b97SToby Isaac   PetscInt       i;
2301a989b97SToby Isaac 
2311a989b97SToby Isaac   PetscFunctionBegin;
232*08401ef6SPierre Jolivet   PetscCheck(N >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
2332c71b3e2SJacob Faibussowitsch   PetscCheckFalse(j < 0 || k < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
2342c71b3e2SJacob Faibussowitsch   PetscCheckFalse(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) {
2401a989b97SToby Isaac       for (i = 0; i < Njk * Njk; i++) {awedgeMat[i] = 0.;}
2411a989b97SToby Isaac       for (i = 0; i < Njk; i++) {awedgeMat[i * (Njk + 1)] = a[0];}
2421a989b97SToby Isaac     } else if (!k) {
2431a989b97SToby Isaac       for (i = 0; i < Njk; i++) {awedgeMat[i] = a[i];}
2441a989b97SToby Isaac     } else {
2451a989b97SToby Isaac       if (N == 2) {
2461a989b97SToby Isaac         awedgeMat[0] = -a[1]; awedgeMat[1] =  a[0];
2471a989b97SToby Isaac       } else {
2481a989b97SToby Isaac         if (j+k == 2) {
2491a989b97SToby Isaac           awedgeMat[0] = -a[1]; awedgeMat[1] =  a[0]; awedgeMat[2] =    0.;
2501a989b97SToby Isaac           awedgeMat[3] = -a[2]; awedgeMat[4] =    0.; awedgeMat[5] =  a[0];
2511a989b97SToby Isaac           awedgeMat[6] =    0.; awedgeMat[7] = -a[2]; awedgeMat[8] =  a[1];
2521a989b97SToby Isaac         } else {
2531a989b97SToby Isaac           awedgeMat[0] =  a[2]; awedgeMat[1] = -a[1]; awedgeMat[2] =  a[0];
2541a989b97SToby Isaac         }
2551a989b97SToby Isaac       }
2561a989b97SToby Isaac     }
2571a989b97SToby Isaac   } else {
2581a989b97SToby Isaac     PetscInt  Njk;
2591a989b97SToby Isaac     PetscInt  Nk;
2601a989b97SToby Isaac     PetscInt  JKj, i;
2611a989b97SToby Isaac     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
2621a989b97SToby Isaac 
2639566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N,   k,   &Nk));
2649566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N,   j+k, &Njk));
2659566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(j+k, j,   &JKj));
2669566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk));
2671a989b97SToby Isaac     for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.;
2681a989b97SToby Isaac     for (i = 0; i < Njk; i++) {
2691a989b97SToby Isaac       PetscInt  l;
2701a989b97SToby Isaac 
2719566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, j+k, i, subset));
2721a989b97SToby Isaac       for (l = 0; l < JKj; l++) {
2731a989b97SToby Isaac         PetscBool jkOdd;
2741a989b97SToby Isaac         PetscInt  m, jInd, kInd;
2751a989b97SToby Isaac 
2769566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd));
2771a989b97SToby Isaac         for (m = 0; m < j; m++) {
2781a989b97SToby Isaac           subsetj[m] = subset[subsetjk[m]];
2791a989b97SToby Isaac         }
2801a989b97SToby Isaac         for (m = 0; m < k; m++) {
2811a989b97SToby Isaac           subsetk[m] = subset[subsetjk[j+m]];
2821a989b97SToby Isaac         }
2839566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd));
2849566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd));
2851a989b97SToby Isaac         awedgeMat[i * Nk + kInd] += jkOdd ? - a[jInd] : a[jInd];
2861a989b97SToby Isaac       }
2871a989b97SToby Isaac     }
2889566063dSJacob Faibussowitsch     PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk));
2891a989b97SToby Isaac   }
2901a989b97SToby Isaac   PetscFunctionReturn(0);
2911a989b97SToby Isaac }
2921a989b97SToby Isaac 
293fad4db65SToby Isaac /*@
29428222859SToby Isaac    PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space
295fad4db65SToby Isaac 
2964165533cSJose E. Roman    Input Parameters:
29728222859SToby Isaac +  N - the dimension of the origin vector space of the linear transformation, M >= 0
29828222859SToby Isaac .  M - the dimension of the image vector space of the linear transformation, N >= 0
29928222859SToby Isaac .  L - a linear transformation, an [M x N] matrix in row-major format
30028222859SToby Isaac .  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).
30128222859SToby Isaac -  w - a |k|-form in the image space, size [M choose |k|]
302fad4db65SToby Isaac 
3034165533cSJose E. Roman    Output Parameter:
30428222859SToby 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).
305fad4db65SToby Isaac 
306fad4db65SToby Isaac    Level: intermediate
307fad4db65SToby Isaac 
308a5b23f4aSJose E. Roman    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,
30928222859SToby Isaac    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,
310fad4db65SToby Isaac    then the inverse Hodge star transformation.
311fad4db65SToby Isaac 
31229a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullbackMatrix(), PetscDTAltVStar()
313fad4db65SToby Isaac @*/
3141a989b97SToby Isaac PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw)
3151a989b97SToby Isaac {
3161a989b97SToby Isaac   PetscInt         i, j, Nk, Mk;
3171a989b97SToby Isaac 
3181a989b97SToby Isaac   PetscFunctionBegin;
3192c71b3e2SJacob Faibussowitsch   PetscCheckFalse(N < 0 || M < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
3202c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsInt(k) > N || PetscAbsInt(k) > M,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
3211a989b97SToby Isaac   if (N <= 3 && M <= 3) {
3221a989b97SToby Isaac 
3239566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
3249566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
3251a989b97SToby Isaac     if (!k) {
3261a989b97SToby Isaac       Lstarw[0] = w[0];
3271a989b97SToby Isaac     } else if (k == 1) {
3281a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3291a989b97SToby Isaac         PetscReal sum = 0.;
3301a989b97SToby Isaac 
3311a989b97SToby Isaac         for (j = 0; j < Mk; j++) {sum += L[j * Nk + i] * w[j];}
3321a989b97SToby Isaac         Lstarw[i] = sum;
3331a989b97SToby Isaac       }
3341a989b97SToby Isaac     } else if (k == -1) {
3351a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
3361a989b97SToby Isaac 
3371a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3381a989b97SToby Isaac         PetscReal sum = 0.;
3391a989b97SToby Isaac 
3401a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
3411a989b97SToby Isaac           sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j];
3421a989b97SToby Isaac         }
3431a989b97SToby Isaac         Lstarw[i] = mult[i] * sum;
3441a989b97SToby Isaac       }
3451a989b97SToby Isaac     } else if (k == 2) {
3461a989b97SToby Isaac       PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}};
3471a989b97SToby Isaac 
3481a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3491a989b97SToby Isaac         PetscReal sum = 0.;
3501a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
3511a989b97SToby Isaac           sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] -
3521a989b97SToby Isaac                   L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j];
3531a989b97SToby Isaac         }
3541a989b97SToby Isaac         Lstarw[i] = sum;
3551a989b97SToby Isaac       }
3561a989b97SToby Isaac     } else if (k == -2) {
3571a989b97SToby Isaac       PetscInt  pairs[3][2] = {{1,2},{2,0},{0,1}};
3581a989b97SToby Isaac       PetscInt  offi = (N == 2) ? 2 : 0;
3591a989b97SToby Isaac       PetscInt  offj = (M == 2) ? 2 : 0;
3601a989b97SToby Isaac 
3611a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3621a989b97SToby Isaac         PetscReal sum   = 0.;
3631a989b97SToby Isaac 
3641a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
3651a989b97SToby Isaac           sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] *
3661a989b97SToby Isaac                   L[pairs[offj + j][1] * N + pairs[offi + i][1]] -
3671a989b97SToby Isaac                   L[pairs[offj + j][1] * N + pairs[offi + i][0]] *
3681a989b97SToby Isaac                   L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j];
3691a989b97SToby Isaac 
3701a989b97SToby Isaac         }
3711a989b97SToby Isaac         Lstarw[i] = sum;
3721a989b97SToby Isaac       }
3731a989b97SToby Isaac     } else {
3741a989b97SToby Isaac       PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) +
3751a989b97SToby Isaac                        L[1] * (L[5] * L[6] - L[3] * L[8]) +
3761a989b97SToby Isaac                        L[2] * (L[3] * L[7] - L[4] * L[6]);
3771a989b97SToby Isaac 
3781a989b97SToby Isaac       for (i = 0; i < Nk; i++) {Lstarw[i] = detL * w[i];}
3791a989b97SToby Isaac     }
3801a989b97SToby Isaac   } else {
3811a989b97SToby Isaac     PetscInt         Nf, l, p;
3821a989b97SToby Isaac     PetscReal       *Lw, *Lwv;
3831a989b97SToby Isaac     PetscInt        *subsetw, *subsetv;
384fad4db65SToby Isaac     PetscInt        *perm;
3851a989b97SToby Isaac     PetscReal       *walloc = NULL;
3861a989b97SToby Isaac     const PetscReal *ww = NULL;
3871a989b97SToby Isaac     PetscBool        negative = PETSC_FALSE;
3881a989b97SToby Isaac 
3899566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
3909566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
3919566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf));
3921a989b97SToby Isaac     if (k < 0) {
3931a989b97SToby Isaac       negative = PETSC_TRUE;
3941a989b97SToby Isaac       k = -k;
3959566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Mk, &walloc));
3969566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc));
3971a989b97SToby Isaac       ww = walloc;
3981a989b97SToby Isaac     } else {
3991a989b97SToby Isaac       ww = w;
4001a989b97SToby Isaac     }
4019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv));
4021a989b97SToby Isaac     for (i = 0; i < Nk; i++) Lstarw[i] = 0.;
4031a989b97SToby Isaac     for (i = 0; i < Mk; i++) {
4049566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(M, k, i, subsetw));
4051a989b97SToby Isaac       for (j = 0; j < Nk; j++) {
4069566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSubset(N, k, j, subsetv));
4071a989b97SToby Isaac         for (p = 0; p < Nf; p++) {
4081a989b97SToby Isaac           PetscReal prod;
4091a989b97SToby Isaac           PetscBool isOdd;
4101a989b97SToby Isaac 
4119566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd));
4121a989b97SToby Isaac           prod = isOdd ? -ww[i] : ww[i];
4131a989b97SToby Isaac           for (l = 0; l < k; l++) {
4141a989b97SToby Isaac             prod *= L[subsetw[perm[l]] * N + subsetv[l]];
4151a989b97SToby Isaac           }
4161a989b97SToby Isaac           Lstarw[j] += prod;
4171a989b97SToby Isaac         }
4181a989b97SToby Isaac       }
4191a989b97SToby Isaac     }
4201a989b97SToby Isaac     if (negative) {
4211a989b97SToby Isaac       PetscReal *sLsw;
4221a989b97SToby Isaac 
4239566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &sLsw));
4249566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(N, N - k, -1,  Lstarw, sLsw));
4251a989b97SToby Isaac       for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i];
4269566063dSJacob Faibussowitsch       PetscCall(PetscFree(sLsw));
4271a989b97SToby Isaac     }
4289566063dSJacob Faibussowitsch     PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv));
4299566063dSJacob Faibussowitsch     PetscCall(PetscFree(walloc));
4301a989b97SToby Isaac   }
4311a989b97SToby Isaac   PetscFunctionReturn(0);
4321a989b97SToby Isaac }
4331a989b97SToby Isaac 
434fad4db65SToby Isaac /*@
435fad4db65SToby Isaac    PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation
436fad4db65SToby Isaac 
4374165533cSJose E. Roman    Input Parameters:
43828222859SToby Isaac +  N - the dimension of the origin vector space of the linear transformation, N >= 0
43928222859SToby Isaac .  M - the dimension of the image vector space of the linear transformation, M >= 0
44028222859SToby Isaac .  L - a linear transformation, an [M x N] matrix in row-major format
44128222859SToby Isaac -  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())
442fad4db65SToby Isaac 
4434165533cSJose E. Roman    Output Parameter:
44428222859SToby Isaac .  Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w
445fad4db65SToby Isaac 
446fad4db65SToby Isaac    Level: intermediate
447fad4db65SToby Isaac 
44829a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVStar()
449fad4db65SToby Isaac @*/
4501a989b97SToby Isaac PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar)
4511a989b97SToby Isaac {
4521a989b97SToby Isaac   PetscInt        Nk, Mk, Nf, i, j, l, p;
4531a989b97SToby Isaac   PetscReal      *Lw, *Lwv;
4541a989b97SToby Isaac   PetscInt       *subsetw, *subsetv;
455fad4db65SToby Isaac   PetscInt       *perm;
4561a989b97SToby Isaac   PetscBool       negative = PETSC_FALSE;
4571a989b97SToby Isaac 
4581a989b97SToby Isaac   PetscFunctionBegin;
4592c71b3e2SJacob Faibussowitsch   PetscCheckFalse(N < 0 || M < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
4602c71b3e2SJacob Faibussowitsch   PetscCheckFalse(PetscAbsInt(k) > N || PetscAbsInt(k) > M,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
4611a989b97SToby Isaac   if (N <= 3 && M <= 3) {
4621a989b97SToby Isaac     PetscReal mult[3] = {1., -1., 1.};
4631a989b97SToby Isaac 
4649566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
4659566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
4661a989b97SToby Isaac     if (!k) {
4671a989b97SToby Isaac       Lstar[0] = 1.;
4681a989b97SToby Isaac     } else if (k == 1) {
4691a989b97SToby Isaac       for (i = 0; i < Nk; i++) {for (j = 0; j < Mk; j++) {Lstar[i * Mk + j] = L[j * Nk + i];}}
4701a989b97SToby Isaac     } else if (k == -1) {
4711a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
4721a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
4731a989b97SToby Isaac           Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j];
4741a989b97SToby Isaac         }
4751a989b97SToby Isaac       }
4761a989b97SToby Isaac     } else if (k == 2) {
4771a989b97SToby Isaac       PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}};
4781a989b97SToby Isaac 
4791a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
4801a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
4811a989b97SToby Isaac           Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] *
4821a989b97SToby Isaac                               L[pairs[j][1] * N + pairs[i][1]] -
4831a989b97SToby Isaac                               L[pairs[j][1] * N + pairs[i][0]] *
4841a989b97SToby Isaac                               L[pairs[j][0] * N + pairs[i][1]];
4851a989b97SToby Isaac         }
4861a989b97SToby Isaac       }
4871a989b97SToby Isaac     } else if (k == -2) {
4881a989b97SToby Isaac       PetscInt  pairs[3][2] = {{1,2},{2,0},{0,1}};
4891a989b97SToby Isaac       PetscInt  offi = (N == 2) ? 2 : 0;
4901a989b97SToby Isaac       PetscInt  offj = (M == 2) ? 2 : 0;
4911a989b97SToby Isaac 
4921a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
4931a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
4941a989b97SToby Isaac           Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] *
4951a989b97SToby Isaac                               L[pairs[offj + j][1] * N + pairs[offi + i][1]] -
4961a989b97SToby Isaac                               L[pairs[offj + j][1] * N + pairs[offi + i][0]] *
4971a989b97SToby Isaac                               L[pairs[offj + j][0] * N + pairs[offi + i][1]];
4981a989b97SToby Isaac         }
4991a989b97SToby Isaac       }
5001a989b97SToby Isaac     } else {
5011a989b97SToby Isaac       PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) +
5021a989b97SToby Isaac                        L[1] * (L[5] * L[6] - L[3] * L[8]) +
5031a989b97SToby Isaac                        L[2] * (L[3] * L[7] - L[4] * L[6]);
5041a989b97SToby Isaac 
5051a989b97SToby Isaac       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.;
5371a989b97SToby Isaac           for (l = 0; l < k; l++) {
5381a989b97SToby Isaac             prod *= L[subsetw[perm[l]] * N + subsetv[l]];
5391a989b97SToby Isaac           }
5401a989b97SToby Isaac           Lstar[jidx * Mk + iidx] += prod;
5411a989b97SToby Isaac         }
5421a989b97SToby Isaac       }
5431a989b97SToby Isaac     }
5449566063dSJacob Faibussowitsch     PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv));
5451a989b97SToby Isaac   }
5461a989b97SToby Isaac   PetscFunctionReturn(0);
5471a989b97SToby Isaac }
5481a989b97SToby Isaac 
549fad4db65SToby Isaac /*@
55028222859SToby Isaac    PetscDTAltVInterior - Compute the interior product of a k-form with a vector
551fad4db65SToby Isaac 
5524165533cSJose E. Roman    Input Parameters:
55328222859SToby Isaac +  N - the dimension of the vector space, N >= 0
55428222859SToby Isaac .  k - the degree k of the k-form w, 0 <= k <= N
55528222859SToby Isaac .  w - a k-form, size [N choose k]
55628222859SToby Isaac -  v - an N dimensional vector
557fad4db65SToby Isaac 
5584165533cSJose E. Roman    Output Parameter:
55928222859SToby 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}).
560fad4db65SToby Isaac 
561fad4db65SToby Isaac    Level: intermediate
562fad4db65SToby Isaac 
56329a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVInteriorMatrix(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
564fad4db65SToby Isaac @*/
5651a989b97SToby Isaac PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv)
5661a989b97SToby Isaac {
5671a989b97SToby Isaac   PetscInt        i, Nk, Nkm;
5681a989b97SToby Isaac 
5691a989b97SToby Isaac   PetscFunctionBegin;
5702c71b3e2SJacob Faibussowitsch   PetscCheckFalse(k <= 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
5719566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k,   &Nk));
5729566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k-1, &Nkm));
5731a989b97SToby Isaac   if (N <= 3) {
5741a989b97SToby Isaac     if (k == 1) {
5751a989b97SToby Isaac       PetscReal sum = 0.;
5761a989b97SToby Isaac 
5771a989b97SToby Isaac       for (i = 0; i < N; i++) {
5781a989b97SToby Isaac         sum += w[i] * v[i];
5791a989b97SToby Isaac       }
5801a989b97SToby Isaac       wIntv[0] = sum;
5811a989b97SToby Isaac     } else if (k == N) {
5821a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
5831a989b97SToby Isaac 
5841a989b97SToby Isaac       for (i = 0; i < N; i++) {
5851a989b97SToby Isaac         wIntv[N - 1 - i] = w[0] * v[i] * mult[i];
5861a989b97SToby Isaac       }
5871a989b97SToby Isaac     } else {
5881a989b97SToby Isaac       wIntv[0] = - w[0]*v[1] - w[1]*v[2];
5891a989b97SToby Isaac       wIntv[1] =   w[0]*v[0] - w[2]*v[2];
5901a989b97SToby Isaac       wIntv[2] =   w[1]*v[0] + w[2]*v[1];
5911a989b97SToby Isaac     }
5921a989b97SToby Isaac   } else {
5931a989b97SToby Isaac     PetscInt       *subset, *work;
5941a989b97SToby Isaac 
5959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
5961a989b97SToby Isaac     for (i = 0; i < Nkm; i++) wIntv[i] = 0.;
5971a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
5981a989b97SToby Isaac       PetscInt  j, l, m;
5991a989b97SToby Isaac 
6009566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
6011a989b97SToby Isaac       for (j = 0; j < k; j++) {
6021a989b97SToby Isaac         PetscInt  idx;
60328222859SToby Isaac         PetscBool flip = (PetscBool) (j & 1);
6041a989b97SToby Isaac 
6051a989b97SToby Isaac         for (l = 0, m = 0; l < k; l++) {
6061a989b97SToby Isaac           if (l != j) work[m++] = subset[l];
6071a989b97SToby Isaac         }
6089566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
6091a989b97SToby Isaac         wIntv[idx] += flip ? -(w[i] * v[subset[j]]) :  (w[i] * v[subset[j]]);
6101a989b97SToby Isaac       }
6111a989b97SToby Isaac     }
6129566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
6131a989b97SToby Isaac   }
6141a989b97SToby Isaac   PetscFunctionReturn(0);
6151a989b97SToby Isaac }
6161a989b97SToby Isaac 
617fad4db65SToby Isaac /*@
61828222859SToby Isaac    PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector
619fad4db65SToby Isaac 
6204165533cSJose E. Roman    Input Parameters:
62128222859SToby Isaac +  N - the dimension of the vector space, N >= 0
62228222859SToby Isaac .  k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N
62328222859SToby Isaac -  v - an N dimensional vector
624fad4db65SToby Isaac 
6254165533cSJose E. Roman    Output Parameter:
626fad4db65SToby Isaac .  intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v)
627fad4db65SToby Isaac 
628fad4db65SToby Isaac    Level: intermediate
629fad4db65SToby Isaac 
63029a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVInterior(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
631fad4db65SToby Isaac @*/
6321a989b97SToby Isaac PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat)
6331a989b97SToby Isaac {
6341a989b97SToby Isaac   PetscInt        i, Nk, Nkm;
6351a989b97SToby Isaac 
6361a989b97SToby Isaac   PetscFunctionBegin;
6372c71b3e2SJacob Faibussowitsch   PetscCheckFalse(k <= 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
6389566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k,   &Nk));
6399566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k-1, &Nkm));
6401a989b97SToby Isaac   if (N <= 3) {
6411a989b97SToby Isaac     if (k == 1) {
6421a989b97SToby Isaac       for (i = 0; i < N; i++) intvMat[i] = v[i];
6431a989b97SToby Isaac     } else if (k == N) {
6441a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
6451a989b97SToby Isaac 
6461a989b97SToby Isaac       for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i];
6471a989b97SToby Isaac     } else {
6481a989b97SToby Isaac       intvMat[0] = -v[1]; intvMat[1] = -v[2]; intvMat[2] =    0.;
6491a989b97SToby Isaac       intvMat[3] =  v[0]; intvMat[4] =    0.; intvMat[5] = -v[2];
6501a989b97SToby Isaac       intvMat[6] =    0.; intvMat[7] =  v[0]; 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   }
6741a989b97SToby Isaac   PetscFunctionReturn(0);
6751a989b97SToby Isaac }
6761a989b97SToby Isaac 
677fad4db65SToby Isaac /*@
678fad4db65SToby Isaac    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
68228222859SToby Isaac -  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 
69128222859SToby Isaac    Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential
692fad4db65SToby Isaac 
69329a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVInterior(), PetscDTAltVInteriorMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
694fad4db65SToby Isaac @*/
695dda711d0SToby Isaac PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3])
696dda711d0SToby Isaac {
697dda711d0SToby Isaac   PetscInt        i, Nk, Nkm;
698dda711d0SToby Isaac 
699dda711d0SToby Isaac   PetscFunctionBegin;
7002c71b3e2SJacob Faibussowitsch   PetscCheckFalse(k <= 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
7019566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k,   &Nk));
7029566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k-1, &Nkm));
703dda711d0SToby Isaac   if (N <= 3) {
704dda711d0SToby Isaac     if (k == 1) {
705dda711d0SToby Isaac       for (i = 0; i < N; i++) {
706dda711d0SToby Isaac         indices[i][0] = 0;
707dda711d0SToby Isaac         indices[i][1] = i;
708dda711d0SToby Isaac         indices[i][2] = i;
709dda711d0SToby Isaac       }
710dda711d0SToby Isaac     } else if (k == N) {
711dda711d0SToby Isaac       PetscInt val[3] = {0, -2, 2};
712dda711d0SToby Isaac 
713dda711d0SToby Isaac       for (i = 0; i < N; i++) {
714dda711d0SToby Isaac         indices[i][0] = N - 1 - i;
715dda711d0SToby Isaac         indices[i][1] = 0;
716dda711d0SToby Isaac         indices[i][2] = val[i];
717dda711d0SToby Isaac       }
718dda711d0SToby Isaac     } else {
719dda711d0SToby Isaac       indices[0][0] = 0; indices[0][1] = 0; indices[0][2] = -(1 + 1);
720dda711d0SToby Isaac       indices[1][0] = 0; indices[1][1] = 1; indices[1][2] = -(2 + 1);
721dda711d0SToby Isaac       indices[2][0] = 1; indices[2][1] = 0; indices[2][2] = 0;
722dda711d0SToby Isaac       indices[3][0] = 1; indices[3][1] = 2; indices[3][2] = -(2 + 1);
723dda711d0SToby Isaac       indices[4][0] = 2; indices[4][1] = 1; indices[4][2] = 0;
724dda711d0SToby Isaac       indices[5][0] = 2; indices[5][1] = 2; indices[5][2] = 1;
725dda711d0SToby Isaac     }
726dda711d0SToby Isaac   } else {
727dda711d0SToby Isaac     PetscInt       *subset, *work;
728dda711d0SToby Isaac 
7299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
730dda711d0SToby Isaac     for (i = 0; i < Nk; i++) {
731dda711d0SToby Isaac       PetscInt  j, l, m;
732dda711d0SToby Isaac 
7339566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
734dda711d0SToby Isaac       for (j = 0; j < k; j++) {
735dda711d0SToby Isaac         PetscInt  idx;
73628222859SToby Isaac         PetscBool flip = (PetscBool) (j & 1);
737dda711d0SToby Isaac 
738dda711d0SToby Isaac         for (l = 0, m = 0; l < k; l++) {
739dda711d0SToby Isaac           if (l != j) work[m++] = subset[l];
740dda711d0SToby Isaac         }
7419566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
742dda711d0SToby Isaac         indices[i * k + j][0] = idx;
743dda711d0SToby Isaac         indices[i * k + j][1] = i;
744dda711d0SToby Isaac         indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j];
745dda711d0SToby Isaac       }
746dda711d0SToby Isaac     }
7479566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
748dda711d0SToby Isaac   }
749dda711d0SToby Isaac   PetscFunctionReturn(0);
750dda711d0SToby Isaac }
751dda711d0SToby Isaac 
752fad4db65SToby Isaac /*@
75328222859SToby Isaac    PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form
754fad4db65SToby Isaac 
7554165533cSJose E. Roman    Input Parameters:
75628222859SToby Isaac +  N - the dimension of the vector space, N >= 0
75728222859SToby Isaac .  k - the degree k of the k-form w, 0 <= k <= N
75828222859SToby 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.
75928222859SToby Isaac -  w - a k-form, size [N choose k]
760fad4db65SToby Isaac 
7614165533cSJose E. Roman    Output Parameter:
76228222859SToby Isaac .  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.
763fad4db65SToby Isaac 
764fad4db65SToby Isaac    Level: intermediate
765fad4db65SToby Isaac 
76629a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
767fad4db65SToby Isaac @*/
7681a989b97SToby Isaac PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw)
7691a989b97SToby Isaac {
7701a989b97SToby Isaac   PetscInt        Nk, i;
7711a989b97SToby Isaac 
7721a989b97SToby Isaac   PetscFunctionBegin;
7732c71b3e2SJacob Faibussowitsch   PetscCheckFalse(k < 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
7749566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
7751a989b97SToby Isaac   pow = pow % 4;
7761a989b97SToby Isaac   pow = (pow + 4) % 4; /* make non-negative */
7771a989b97SToby Isaac   /* pow is now 0, 1, 2, 3 */
7781a989b97SToby Isaac   if (N <= 3) {
7791a989b97SToby Isaac     if (pow & 1) {
7801a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
7811a989b97SToby Isaac 
7821a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i];
7831a989b97SToby Isaac     } else {
7841a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = w[i];
7851a989b97SToby Isaac     }
7861a989b97SToby Isaac     if (pow > 1 && ((k * (N - k)) & 1)) {
7871a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
7881a989b97SToby Isaac     }
7891a989b97SToby Isaac   } else {
7901a989b97SToby Isaac     PetscInt       *subset;
7911a989b97SToby Isaac 
7929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N, &subset));
7931a989b97SToby Isaac     if (pow % 2) {
7941a989b97SToby Isaac       PetscInt l = (pow == 1) ? k : N - k;
7951a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
7961a989b97SToby Isaac         PetscBool sOdd;
7971a989b97SToby Isaac         PetscInt  j, idx;
7981a989b97SToby Isaac 
7999566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, l, i, subset, &sOdd));
8009566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, l, subset, &idx));
8019566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, N-l, &subset[l], &j));
8021a989b97SToby Isaac         starw[j] = sOdd ? -w[idx] : w[idx];
8031a989b97SToby Isaac       }
8041a989b97SToby Isaac     } else {
8051a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = w[i];
8061a989b97SToby Isaac     }
8071a989b97SToby Isaac     /* star^2 = -1^(k * (N - k)) */
8081a989b97SToby Isaac     if (pow > 1 && (k * (N - k)) % 2) {
8091a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
8101a989b97SToby Isaac     }
8119566063dSJacob Faibussowitsch     PetscCall(PetscFree(subset));
8121a989b97SToby Isaac   }
8131a989b97SToby Isaac   PetscFunctionReturn(0);
8141a989b97SToby Isaac }
815