xref: /petsc/src/dm/dt/interface/dtaltv.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
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
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 
71db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
72fad4db65SToby Isaac @*/
73*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv)
74*d71ae5a4SJacob Faibussowitsch {
751a989b97SToby Isaac   PetscFunctionBegin;
7608401ef6SPierre Jolivet   PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
771dca8a05SBarry Smith   PetscCheck(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 {
829371c9d4SSatish Balay       if (N == 1) {
839371c9d4SSatish Balay         *wv = w[0] * v[0];
849371c9d4SSatish Balay       } else if (N == 2) {
859371c9d4SSatish Balay         if (k == 1) {
869371c9d4SSatish Balay           *wv = w[0] * v[0] + w[1] * v[1];
871a989b97SToby Isaac         } else {
889371c9d4SSatish Balay           *wv = w[0] * (v[0] * v[3] - v[1] * v[2]);
899371c9d4SSatish Balay         }
901a989b97SToby Isaac       } else {
919371c9d4SSatish Balay         if (k == 1) {
929371c9d4SSatish Balay           *wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];
939371c9d4SSatish Balay         } else if (k == 2) {
949371c9d4SSatish 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]);
959371c9d4SSatish Balay         } else {
969371c9d4SSatish 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]));
971a989b97SToby Isaac         }
981a989b97SToby Isaac       }
991a989b97SToby Isaac     }
1001a989b97SToby Isaac   } else {
1011a989b97SToby Isaac     PetscInt  Nk, Nf;
102fad4db65SToby Isaac     PetscInt *subset, *perm;
1031a989b97SToby Isaac     PetscInt  i, j, l;
1041a989b97SToby Isaac     PetscReal sum = 0.;
1051a989b97SToby Isaac 
1069566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(k, &Nf));
1079566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, k, &Nk));
1089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &perm));
1091a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
1101a989b97SToby Isaac       PetscReal subsum = 0.;
1111a989b97SToby Isaac 
1129566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
1131a989b97SToby Isaac       for (j = 0; j < Nf; j++) {
1141a989b97SToby Isaac         PetscBool permOdd;
1151a989b97SToby Isaac         PetscReal prod;
1161a989b97SToby Isaac 
1179566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumPerm(k, j, perm, &permOdd));
1181a989b97SToby Isaac         prod = permOdd ? -1. : 1.;
119ad540459SPierre Jolivet         for (l = 0; l < k; l++) prod *= v[perm[l] * N + subset[l]];
1201a989b97SToby Isaac         subsum += prod;
1211a989b97SToby Isaac       }
1221a989b97SToby Isaac       sum += w[i] * subsum;
1231a989b97SToby Isaac     }
1249566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, perm));
1251a989b97SToby Isaac     *wv = sum;
1261a989b97SToby Isaac   }
1271a989b97SToby Isaac   PetscFunctionReturn(0);
1281a989b97SToby Isaac }
1291a989b97SToby Isaac 
130fad4db65SToby Isaac /*@
13128222859SToby Isaac    PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form
132fad4db65SToby Isaac 
1334165533cSJose E. Roman    Input Parameters:
13428222859SToby Isaac +  N - the dimension of the vector space, N >= 0
13528222859SToby Isaac .  j - the degree j of the j-form a, 0 <= j <= N
13628222859SToby Isaac .  k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N
13728222859SToby Isaac .  a - a j-form, size [N choose j]
13828222859SToby Isaac -  b - a k-form, size [N choose k]
139fad4db65SToby Isaac 
1404165533cSJose E. Roman    Output Parameter:
14128222859SToby 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}}),
14228222859SToby 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}.
143fad4db65SToby Isaac 
144fad4db65SToby Isaac    Level: intermediate
145fad4db65SToby Isaac 
146db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVWedgeMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
147fad4db65SToby Isaac @*/
148*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb)
149*d71ae5a4SJacob Faibussowitsch {
1501a989b97SToby Isaac   PetscInt i;
1511a989b97SToby Isaac 
1521a989b97SToby Isaac   PetscFunctionBegin;
15308401ef6SPierre Jolivet   PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
1541dca8a05SBarry Smith   PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
1551dca8a05SBarry Smith   PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
1561a989b97SToby Isaac   if (N <= 3) {
1571a989b97SToby Isaac     PetscInt Njk;
1581a989b97SToby Isaac 
1599566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
1609371c9d4SSatish Balay     if (!j) {
161ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeb[i] = a[0] * b[i];
1629371c9d4SSatish Balay     } else if (!k) {
163ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeb[i] = a[i] * b[0];
1649371c9d4SSatish Balay     } else {
1659371c9d4SSatish Balay       if (N == 2) {
1669371c9d4SSatish Balay         awedgeb[0] = a[0] * b[1] - a[1] * b[0];
1679371c9d4SSatish Balay       } else {
1681a989b97SToby Isaac         if (j + k == 2) {
1691a989b97SToby Isaac           awedgeb[0] = a[0] * b[1] - a[1] * b[0];
1701a989b97SToby Isaac           awedgeb[1] = a[0] * b[2] - a[2] * b[0];
1711a989b97SToby Isaac           awedgeb[2] = a[1] * b[2] - a[2] * b[1];
1721a989b97SToby Isaac         } else {
1731a989b97SToby Isaac           awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0];
1741a989b97SToby Isaac         }
1751a989b97SToby Isaac       }
1761a989b97SToby Isaac     }
1771a989b97SToby Isaac   } else {
1781a989b97SToby Isaac     PetscInt  Njk;
1791a989b97SToby Isaac     PetscInt  JKj;
1801a989b97SToby Isaac     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
1811a989b97SToby Isaac     PetscInt  i;
1821a989b97SToby Isaac 
1839566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
1849566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(j + k, j, &JKj));
1859566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk));
1861a989b97SToby Isaac     for (i = 0; i < Njk; i++) {
1871a989b97SToby Isaac       PetscReal sum = 0.;
1881a989b97SToby Isaac       PetscInt  l;
1891a989b97SToby Isaac 
1909566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, j + k, i, subset));
1911a989b97SToby Isaac       for (l = 0; l < JKj; l++) {
1921a989b97SToby Isaac         PetscBool jkOdd;
1931a989b97SToby Isaac         PetscInt  m, jInd, kInd;
1941a989b97SToby Isaac 
1959566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd));
196ad540459SPierre Jolivet         for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]];
197ad540459SPierre Jolivet         for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]];
1989566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd));
1999566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd));
2001a989b97SToby Isaac         sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]);
2011a989b97SToby Isaac       }
2021a989b97SToby Isaac       awedgeb[i] = sum;
2031a989b97SToby Isaac     }
2049566063dSJacob Faibussowitsch     PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk));
2051a989b97SToby Isaac   }
2061a989b97SToby Isaac   PetscFunctionReturn(0);
2071a989b97SToby Isaac }
2081a989b97SToby Isaac 
209fad4db65SToby Isaac /*@
21028222859SToby Isaac    PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms
211fad4db65SToby Isaac 
2124165533cSJose E. Roman    Input Parameters:
21328222859SToby Isaac +  N - the dimension of the vector space, N >= 0
21428222859SToby Isaac .  j - the degree j of the j-form a, 0 <= j <= N
21528222859SToby Isaac .  k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N
21628222859SToby Isaac -  a - a j-form, size [N choose j]
217fad4db65SToby Isaac 
2184165533cSJose E. Roman    Output Parameter:
21928222859SToby 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
220fad4db65SToby Isaac 
221fad4db65SToby Isaac    Level: intermediate
222fad4db65SToby Isaac 
223db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
224fad4db65SToby Isaac @*/
225*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat)
226*d71ae5a4SJacob Faibussowitsch {
2271a989b97SToby Isaac   PetscInt i;
2281a989b97SToby Isaac 
2291a989b97SToby Isaac   PetscFunctionBegin;
23008401ef6SPierre Jolivet   PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension");
2311dca8a05SBarry Smith   PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree");
2321dca8a05SBarry Smith   PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension");
2331a989b97SToby Isaac   if (N <= 3) {
2341a989b97SToby Isaac     PetscInt Njk;
2351a989b97SToby Isaac 
2369566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
2371a989b97SToby Isaac     if (!j) {
238ad540459SPierre Jolivet       for (i = 0; i < Njk * Njk; i++) awedgeMat[i] = 0.;
239ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeMat[i * (Njk + 1)] = a[0];
2401a989b97SToby Isaac     } else if (!k) {
241ad540459SPierre Jolivet       for (i = 0; i < Njk; i++) awedgeMat[i] = a[i];
2421a989b97SToby Isaac     } else {
2431a989b97SToby Isaac       if (N == 2) {
2449371c9d4SSatish Balay         awedgeMat[0] = -a[1];
2459371c9d4SSatish Balay         awedgeMat[1] = a[0];
2461a989b97SToby Isaac       } else {
2471a989b97SToby Isaac         if (j + k == 2) {
2489371c9d4SSatish Balay           awedgeMat[0] = -a[1];
2499371c9d4SSatish Balay           awedgeMat[1] = a[0];
2509371c9d4SSatish Balay           awedgeMat[2] = 0.;
2519371c9d4SSatish Balay           awedgeMat[3] = -a[2];
2529371c9d4SSatish Balay           awedgeMat[4] = 0.;
2539371c9d4SSatish Balay           awedgeMat[5] = a[0];
2549371c9d4SSatish Balay           awedgeMat[6] = 0.;
2559371c9d4SSatish Balay           awedgeMat[7] = -a[2];
2569371c9d4SSatish Balay           awedgeMat[8] = a[1];
2571a989b97SToby Isaac         } else {
2589371c9d4SSatish Balay           awedgeMat[0] = a[2];
2599371c9d4SSatish Balay           awedgeMat[1] = -a[1];
2609371c9d4SSatish Balay           awedgeMat[2] = a[0];
2611a989b97SToby Isaac         }
2621a989b97SToby Isaac       }
2631a989b97SToby Isaac     }
2641a989b97SToby Isaac   } else {
2651a989b97SToby Isaac     PetscInt  Njk;
2661a989b97SToby Isaac     PetscInt  Nk;
2671a989b97SToby Isaac     PetscInt  JKj, i;
2681a989b97SToby Isaac     PetscInt *subset, *subsetjk, *subsetj, *subsetk;
2691a989b97SToby Isaac 
2709566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, k, &Nk));
2719566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
2729566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(j + k, j, &JKj));
2739566063dSJacob Faibussowitsch     PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk));
2741a989b97SToby Isaac     for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.;
2751a989b97SToby Isaac     for (i = 0; i < Njk; i++) {
2761a989b97SToby Isaac       PetscInt l;
2771a989b97SToby Isaac 
2789566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, j + k, i, subset));
2791a989b97SToby Isaac       for (l = 0; l < JKj; l++) {
2801a989b97SToby Isaac         PetscBool jkOdd;
2811a989b97SToby Isaac         PetscInt  m, jInd, kInd;
2821a989b97SToby Isaac 
2839566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd));
284ad540459SPierre Jolivet         for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]];
285ad540459SPierre Jolivet         for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]];
2869566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd));
2879566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd));
2881a989b97SToby Isaac         awedgeMat[i * Nk + kInd] += jkOdd ? -a[jInd] : a[jInd];
2891a989b97SToby Isaac       }
2901a989b97SToby Isaac     }
2919566063dSJacob Faibussowitsch     PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk));
2921a989b97SToby Isaac   }
2931a989b97SToby Isaac   PetscFunctionReturn(0);
2941a989b97SToby Isaac }
2951a989b97SToby Isaac 
296fad4db65SToby Isaac /*@
29728222859SToby Isaac    PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space
298fad4db65SToby Isaac 
2994165533cSJose E. Roman    Input Parameters:
30028222859SToby Isaac +  N - the dimension of the origin vector space of the linear transformation, M >= 0
30128222859SToby Isaac .  M - the dimension of the image vector space of the linear transformation, N >= 0
30228222859SToby Isaac .  L - a linear transformation, an [M x N] matrix in row-major format
30328222859SToby 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).
30428222859SToby Isaac -  w - a |k|-form in the image space, size [M choose |k|]
305fad4db65SToby Isaac 
3064165533cSJose E. Roman    Output Parameter:
30728222859SToby 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).
308fad4db65SToby Isaac 
309fad4db65SToby Isaac    Level: intermediate
310fad4db65SToby Isaac 
311a5b23f4aSJose 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,
31228222859SToby 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,
313fad4db65SToby Isaac    then the inverse Hodge star transformation.
314fad4db65SToby Isaac 
315db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullbackMatrix()`, `PetscDTAltVStar()`
316fad4db65SToby Isaac @*/
317*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw)
318*d71ae5a4SJacob Faibussowitsch {
3191a989b97SToby Isaac   PetscInt i, j, Nk, Mk;
3201a989b97SToby Isaac 
3211a989b97SToby Isaac   PetscFunctionBegin;
3221dca8a05SBarry Smith   PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
3231dca8a05SBarry Smith   PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
3241a989b97SToby Isaac   if (N <= 3 && M <= 3) {
3259566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
3269566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
3271a989b97SToby Isaac     if (!k) {
3281a989b97SToby Isaac       Lstarw[0] = w[0];
3291a989b97SToby Isaac     } else if (k == 1) {
3301a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3311a989b97SToby Isaac         PetscReal sum = 0.;
3321a989b97SToby Isaac 
333ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) sum += L[j * Nk + i] * w[j];
3341a989b97SToby Isaac         Lstarw[i] = sum;
3351a989b97SToby Isaac       }
3361a989b97SToby Isaac     } else if (k == -1) {
3371a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
3381a989b97SToby Isaac 
3391a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3401a989b97SToby Isaac         PetscReal sum = 0.;
3411a989b97SToby Isaac 
342ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j];
3431a989b97SToby Isaac         Lstarw[i] = mult[i] * sum;
3441a989b97SToby Isaac       }
3451a989b97SToby Isaac     } else if (k == 2) {
3469371c9d4SSatish Balay       PetscInt pairs[3][2] = {
3479371c9d4SSatish Balay         {0, 1},
3489371c9d4SSatish Balay         {0, 2},
3499371c9d4SSatish Balay         {1, 2}
3509371c9d4SSatish Balay       };
3511a989b97SToby Isaac 
3521a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3531a989b97SToby Isaac         PetscReal sum = 0.;
354ad540459SPierre 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];
3551a989b97SToby Isaac         Lstarw[i] = sum;
3561a989b97SToby Isaac       }
3571a989b97SToby Isaac     } else if (k == -2) {
3589371c9d4SSatish Balay       PetscInt pairs[3][2] = {
3599371c9d4SSatish Balay         {1, 2},
3609371c9d4SSatish Balay         {2, 0},
3619371c9d4SSatish Balay         {0, 1}
3629371c9d4SSatish Balay       };
3631a989b97SToby Isaac       PetscInt offi = (N == 2) ? 2 : 0;
3641a989b97SToby Isaac       PetscInt offj = (M == 2) ? 2 : 0;
3651a989b97SToby Isaac 
3661a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
3671a989b97SToby Isaac         PetscReal sum = 0.;
3681a989b97SToby Isaac 
369ad540459SPierre 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];
3701a989b97SToby Isaac         Lstarw[i] = sum;
3711a989b97SToby Isaac       }
3721a989b97SToby Isaac     } else {
3739371c9d4SSatish 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]);
3741a989b97SToby Isaac 
375ad540459SPierre Jolivet       for (i = 0; i < Nk; i++) Lstarw[i] = detL * w[i];
3761a989b97SToby Isaac     }
3771a989b97SToby Isaac   } else {
3781a989b97SToby Isaac     PetscInt         Nf, l, p;
3791a989b97SToby Isaac     PetscReal       *Lw, *Lwv;
3801a989b97SToby Isaac     PetscInt        *subsetw, *subsetv;
381fad4db65SToby Isaac     PetscInt        *perm;
3821a989b97SToby Isaac     PetscReal       *walloc   = NULL;
3831a989b97SToby Isaac     const PetscReal *ww       = NULL;
3841a989b97SToby Isaac     PetscBool        negative = PETSC_FALSE;
3851a989b97SToby Isaac 
3869566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
3879566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
3889566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf));
3891a989b97SToby Isaac     if (k < 0) {
3901a989b97SToby Isaac       negative = PETSC_TRUE;
3911a989b97SToby Isaac       k        = -k;
3929566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Mk, &walloc));
3939566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc));
3941a989b97SToby Isaac       ww = walloc;
3951a989b97SToby Isaac     } else {
3961a989b97SToby Isaac       ww = w;
3971a989b97SToby Isaac     }
3989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv));
3991a989b97SToby Isaac     for (i = 0; i < Nk; i++) Lstarw[i] = 0.;
4001a989b97SToby Isaac     for (i = 0; i < Mk; i++) {
4019566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(M, k, i, subsetw));
4021a989b97SToby Isaac       for (j = 0; j < Nk; j++) {
4039566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSubset(N, k, j, subsetv));
4041a989b97SToby Isaac         for (p = 0; p < Nf; p++) {
4051a989b97SToby Isaac           PetscReal prod;
4061a989b97SToby Isaac           PetscBool isOdd;
4071a989b97SToby Isaac 
4089566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd));
4091a989b97SToby Isaac           prod = isOdd ? -ww[i] : ww[i];
410ad540459SPierre Jolivet           for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]];
4111a989b97SToby Isaac           Lstarw[j] += prod;
4121a989b97SToby Isaac         }
4131a989b97SToby Isaac       }
4141a989b97SToby Isaac     }
4151a989b97SToby Isaac     if (negative) {
4161a989b97SToby Isaac       PetscReal *sLsw;
4171a989b97SToby Isaac 
4189566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &sLsw));
4199566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw));
4201a989b97SToby Isaac       for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i];
4219566063dSJacob Faibussowitsch       PetscCall(PetscFree(sLsw));
4221a989b97SToby Isaac     }
4239566063dSJacob Faibussowitsch     PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv));
4249566063dSJacob Faibussowitsch     PetscCall(PetscFree(walloc));
4251a989b97SToby Isaac   }
4261a989b97SToby Isaac   PetscFunctionReturn(0);
4271a989b97SToby Isaac }
4281a989b97SToby Isaac 
429fad4db65SToby Isaac /*@
430fad4db65SToby Isaac    PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation
431fad4db65SToby Isaac 
4324165533cSJose E. Roman    Input Parameters:
43328222859SToby Isaac +  N - the dimension of the origin vector space of the linear transformation, N >= 0
43428222859SToby Isaac .  M - the dimension of the image vector space of the linear transformation, M >= 0
43528222859SToby Isaac .  L - a linear transformation, an [M x N] matrix in row-major format
43628222859SToby 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())
437fad4db65SToby Isaac 
4384165533cSJose E. Roman    Output Parameter:
43928222859SToby Isaac .  Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w
440fad4db65SToby Isaac 
441fad4db65SToby Isaac    Level: intermediate
442fad4db65SToby Isaac 
443db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVStar()`
444fad4db65SToby Isaac @*/
445*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar)
446*d71ae5a4SJacob Faibussowitsch {
4471a989b97SToby Isaac   PetscInt   Nk, Mk, Nf, i, j, l, p;
4481a989b97SToby Isaac   PetscReal *Lw, *Lwv;
4491a989b97SToby Isaac   PetscInt  *subsetw, *subsetv;
450fad4db65SToby Isaac   PetscInt  *perm;
4511a989b97SToby Isaac   PetscBool  negative = PETSC_FALSE;
4521a989b97SToby Isaac 
4531a989b97SToby Isaac   PetscFunctionBegin;
4541dca8a05SBarry Smith   PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions");
4551dca8a05SBarry Smith   PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
4561a989b97SToby Isaac   if (N <= 3 && M <= 3) {
4571a989b97SToby Isaac     PetscReal mult[3] = {1., -1., 1.};
4581a989b97SToby Isaac 
4599566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
4609566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
4611a989b97SToby Isaac     if (!k) {
4621a989b97SToby Isaac       Lstar[0] = 1.;
4631a989b97SToby Isaac     } else if (k == 1) {
4649371c9d4SSatish Balay       for (i = 0; i < Nk; i++) {
465ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[j * Nk + i];
4669371c9d4SSatish Balay       }
4671a989b97SToby Isaac     } else if (k == -1) {
4681a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
469ad540459SPierre Jolivet         for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j];
4701a989b97SToby Isaac       }
4711a989b97SToby Isaac     } else if (k == 2) {
4729371c9d4SSatish Balay       PetscInt pairs[3][2] = {
4739371c9d4SSatish Balay         {0, 1},
4749371c9d4SSatish Balay         {0, 2},
4759371c9d4SSatish Balay         {1, 2}
4769371c9d4SSatish Balay       };
4771a989b97SToby Isaac 
4781a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
479ad540459SPierre 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]];
4801a989b97SToby Isaac       }
4811a989b97SToby Isaac     } else if (k == -2) {
4829371c9d4SSatish Balay       PetscInt pairs[3][2] = {
4839371c9d4SSatish Balay         {1, 2},
4849371c9d4SSatish Balay         {2, 0},
4859371c9d4SSatish Balay         {0, 1}
4869371c9d4SSatish Balay       };
4871a989b97SToby Isaac       PetscInt offi = (N == 2) ? 2 : 0;
4881a989b97SToby Isaac       PetscInt offj = (M == 2) ? 2 : 0;
4891a989b97SToby Isaac 
4901a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
4911a989b97SToby Isaac         for (j = 0; j < Mk; j++) {
4929371c9d4SSatish 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]];
4931a989b97SToby Isaac         }
4941a989b97SToby Isaac       }
4951a989b97SToby Isaac     } else {
4969371c9d4SSatish 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]);
4971a989b97SToby Isaac 
498ad540459SPierre Jolivet       for (i = 0; i < Nk; i++) Lstar[i] = detL;
4991a989b97SToby Isaac     }
5001a989b97SToby Isaac   } else {
5011a989b97SToby Isaac     if (k < 0) {
5021a989b97SToby Isaac       negative = PETSC_TRUE;
5031a989b97SToby Isaac       k        = -k;
5041a989b97SToby Isaac     }
5059566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk));
5069566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk));
5079566063dSJacob Faibussowitsch     PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf));
5089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv));
5091a989b97SToby Isaac     for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.;
5101a989b97SToby Isaac     for (i = 0; i < Mk; i++) {
5111a989b97SToby Isaac       PetscBool iOdd;
5121a989b97SToby Isaac       PetscInt  iidx, jidx;
5131a989b97SToby Isaac 
5149566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSplit(M, k, i, subsetw, &iOdd));
5151a989b97SToby Isaac       iidx = negative ? Mk - 1 - i : i;
51628222859SToby Isaac       iOdd = negative ? (PetscBool)(iOdd ^ ((k * (M - k)) & 1)) : PETSC_FALSE;
5171a989b97SToby Isaac       for (j = 0; j < Nk; j++) {
5181a989b97SToby Isaac         PetscBool jOdd;
5191a989b97SToby Isaac 
5209566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, k, j, subsetv, &jOdd));
5211a989b97SToby Isaac         jidx = negative ? Nk - 1 - j : j;
52228222859SToby Isaac         jOdd = negative ? (PetscBool)(iOdd ^ jOdd ^ ((k * (N - k)) & 1)) : PETSC_FALSE;
5231a989b97SToby Isaac         for (p = 0; p < Nf; p++) {
5241a989b97SToby Isaac           PetscReal prod;
5251a989b97SToby Isaac           PetscBool isOdd;
5261a989b97SToby Isaac 
5279566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd));
52828222859SToby Isaac           isOdd = (PetscBool)(isOdd ^ jOdd);
5291a989b97SToby Isaac           prod  = isOdd ? -1. : 1.;
530ad540459SPierre Jolivet           for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]];
5311a989b97SToby Isaac           Lstar[jidx * Mk + iidx] += prod;
5321a989b97SToby Isaac         }
5331a989b97SToby Isaac       }
5341a989b97SToby Isaac     }
5359566063dSJacob Faibussowitsch     PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv));
5361a989b97SToby Isaac   }
5371a989b97SToby Isaac   PetscFunctionReturn(0);
5381a989b97SToby Isaac }
5391a989b97SToby Isaac 
540fad4db65SToby Isaac /*@
54128222859SToby Isaac    PetscDTAltVInterior - Compute the interior product of a k-form with a vector
542fad4db65SToby Isaac 
5434165533cSJose E. Roman    Input Parameters:
54428222859SToby Isaac +  N - the dimension of the vector space, N >= 0
54528222859SToby Isaac .  k - the degree k of the k-form w, 0 <= k <= N
54628222859SToby Isaac .  w - a k-form, size [N choose k]
54728222859SToby Isaac -  v - an N dimensional vector
548fad4db65SToby Isaac 
5494165533cSJose E. Roman    Output Parameter:
55028222859SToby 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}).
551fad4db65SToby Isaac 
552fad4db65SToby Isaac    Level: intermediate
553fad4db65SToby Isaac 
554db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
555fad4db65SToby Isaac @*/
556*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv)
557*d71ae5a4SJacob Faibussowitsch {
5581a989b97SToby Isaac   PetscInt i, Nk, Nkm;
5591a989b97SToby Isaac 
5601a989b97SToby Isaac   PetscFunctionBegin;
5611dca8a05SBarry Smith   PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
5629566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
5639566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
5641a989b97SToby Isaac   if (N <= 3) {
5651a989b97SToby Isaac     if (k == 1) {
5661a989b97SToby Isaac       PetscReal sum = 0.;
5671a989b97SToby Isaac 
568ad540459SPierre Jolivet       for (i = 0; i < N; i++) sum += w[i] * v[i];
5691a989b97SToby Isaac       wIntv[0] = sum;
5701a989b97SToby Isaac     } else if (k == N) {
5711a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
5721a989b97SToby Isaac 
573ad540459SPierre Jolivet       for (i = 0; i < N; i++) wIntv[N - 1 - i] = w[0] * v[i] * mult[i];
5741a989b97SToby Isaac     } else {
5751a989b97SToby Isaac       wIntv[0] = -w[0] * v[1] - w[1] * v[2];
5761a989b97SToby Isaac       wIntv[1] = w[0] * v[0] - w[2] * v[2];
5771a989b97SToby Isaac       wIntv[2] = w[1] * v[0] + w[2] * v[1];
5781a989b97SToby Isaac     }
5791a989b97SToby Isaac   } else {
5801a989b97SToby Isaac     PetscInt *subset, *work;
5811a989b97SToby Isaac 
5829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
5831a989b97SToby Isaac     for (i = 0; i < Nkm; i++) wIntv[i] = 0.;
5841a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
5851a989b97SToby Isaac       PetscInt j, l, m;
5861a989b97SToby Isaac 
5879566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
5881a989b97SToby Isaac       for (j = 0; j < k; j++) {
5891a989b97SToby Isaac         PetscInt  idx;
59028222859SToby Isaac         PetscBool flip = (PetscBool)(j & 1);
5911a989b97SToby Isaac 
5921a989b97SToby Isaac         for (l = 0, m = 0; l < k; l++) {
5931a989b97SToby Isaac           if (l != j) work[m++] = subset[l];
5941a989b97SToby Isaac         }
5959566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
5961a989b97SToby Isaac         wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]);
5971a989b97SToby Isaac       }
5981a989b97SToby Isaac     }
5999566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
6001a989b97SToby Isaac   }
6011a989b97SToby Isaac   PetscFunctionReturn(0);
6021a989b97SToby Isaac }
6031a989b97SToby Isaac 
604fad4db65SToby Isaac /*@
60528222859SToby Isaac    PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector
606fad4db65SToby Isaac 
6074165533cSJose E. Roman    Input Parameters:
60828222859SToby Isaac +  N - the dimension of the vector space, N >= 0
60928222859SToby Isaac .  k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N
61028222859SToby Isaac -  v - an N dimensional vector
611fad4db65SToby Isaac 
6124165533cSJose E. Roman    Output Parameter:
613fad4db65SToby Isaac .  intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v)
614fad4db65SToby Isaac 
615fad4db65SToby Isaac    Level: intermediate
616fad4db65SToby Isaac 
617db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
618fad4db65SToby Isaac @*/
619*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat)
620*d71ae5a4SJacob Faibussowitsch {
6211a989b97SToby Isaac   PetscInt i, Nk, Nkm;
6221a989b97SToby Isaac 
6231a989b97SToby Isaac   PetscFunctionBegin;
6241dca8a05SBarry Smith   PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
6259566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
6269566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
6271a989b97SToby Isaac   if (N <= 3) {
6281a989b97SToby Isaac     if (k == 1) {
6291a989b97SToby Isaac       for (i = 0; i < N; i++) intvMat[i] = v[i];
6301a989b97SToby Isaac     } else if (k == N) {
6311a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
6321a989b97SToby Isaac 
6331a989b97SToby Isaac       for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i];
6341a989b97SToby Isaac     } else {
6359371c9d4SSatish Balay       intvMat[0] = -v[1];
6369371c9d4SSatish Balay       intvMat[1] = -v[2];
6379371c9d4SSatish Balay       intvMat[2] = 0.;
6389371c9d4SSatish Balay       intvMat[3] = v[0];
6399371c9d4SSatish Balay       intvMat[4] = 0.;
6409371c9d4SSatish Balay       intvMat[5] = -v[2];
6419371c9d4SSatish Balay       intvMat[6] = 0.;
6429371c9d4SSatish Balay       intvMat[7] = v[0];
6439371c9d4SSatish Balay       intvMat[8] = v[1];
6441a989b97SToby Isaac     }
6451a989b97SToby Isaac   } else {
6461a989b97SToby Isaac     PetscInt *subset, *work;
6471a989b97SToby Isaac 
6489566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
6491a989b97SToby Isaac     for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.;
6501a989b97SToby Isaac     for (i = 0; i < Nk; i++) {
6511a989b97SToby Isaac       PetscInt j, l, m;
6521a989b97SToby Isaac 
6539566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
6541a989b97SToby Isaac       for (j = 0; j < k; j++) {
6551a989b97SToby Isaac         PetscInt  idx;
65628222859SToby Isaac         PetscBool flip = (PetscBool)(j & 1);
6571a989b97SToby Isaac 
6581a989b97SToby Isaac         for (l = 0, m = 0; l < k; l++) {
6591a989b97SToby Isaac           if (l != j) work[m++] = subset[l];
6601a989b97SToby Isaac         }
6619566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
6621a989b97SToby Isaac         intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]];
6631a989b97SToby Isaac       }
6641a989b97SToby Isaac     }
6659566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
6661a989b97SToby Isaac   }
6671a989b97SToby Isaac   PetscFunctionReturn(0);
6681a989b97SToby Isaac }
6691a989b97SToby Isaac 
670fad4db65SToby Isaac /*@
671fad4db65SToby Isaac    PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix()
672fad4db65SToby Isaac 
6734165533cSJose E. Roman    Input Parameters:
67428222859SToby Isaac +  N - the dimension of the vector space, N >= 0
67528222859SToby Isaac -  k - the degree of the k-forms on which intvMat from PetscDTAltVInteriorMatrix() acts, 0 <= k <= N.
676fad4db65SToby Isaac 
6774165533cSJose E. Roman    Output Parameter:
67828222859SToby Isaac .  indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k
67928222859SToby 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
68028222859SToby Isaac              coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1)
681fad4db65SToby Isaac 
682fad4db65SToby Isaac    Level: intermediate
683fad4db65SToby Isaac 
68428222859SToby Isaac    Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential
685fad4db65SToby Isaac 
686db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
687fad4db65SToby Isaac @*/
688*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3])
689*d71ae5a4SJacob Faibussowitsch {
690dda711d0SToby Isaac   PetscInt i, Nk, Nkm;
691dda711d0SToby Isaac 
692dda711d0SToby Isaac   PetscFunctionBegin;
6931dca8a05SBarry Smith   PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
6949566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
6959566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
696dda711d0SToby Isaac   if (N <= 3) {
697dda711d0SToby Isaac     if (k == 1) {
698dda711d0SToby Isaac       for (i = 0; i < N; i++) {
699dda711d0SToby Isaac         indices[i][0] = 0;
700dda711d0SToby Isaac         indices[i][1] = i;
701dda711d0SToby Isaac         indices[i][2] = i;
702dda711d0SToby Isaac       }
703dda711d0SToby Isaac     } else if (k == N) {
704dda711d0SToby Isaac       PetscInt val[3] = {0, -2, 2};
705dda711d0SToby Isaac 
706dda711d0SToby Isaac       for (i = 0; i < N; i++) {
707dda711d0SToby Isaac         indices[i][0] = N - 1 - i;
708dda711d0SToby Isaac         indices[i][1] = 0;
709dda711d0SToby Isaac         indices[i][2] = val[i];
710dda711d0SToby Isaac       }
711dda711d0SToby Isaac     } else {
7129371c9d4SSatish Balay       indices[0][0] = 0;
7139371c9d4SSatish Balay       indices[0][1] = 0;
7149371c9d4SSatish Balay       indices[0][2] = -(1 + 1);
7159371c9d4SSatish Balay       indices[1][0] = 0;
7169371c9d4SSatish Balay       indices[1][1] = 1;
7179371c9d4SSatish Balay       indices[1][2] = -(2 + 1);
7189371c9d4SSatish Balay       indices[2][0] = 1;
7199371c9d4SSatish Balay       indices[2][1] = 0;
7209371c9d4SSatish Balay       indices[2][2] = 0;
7219371c9d4SSatish Balay       indices[3][0] = 1;
7229371c9d4SSatish Balay       indices[3][1] = 2;
7239371c9d4SSatish Balay       indices[3][2] = -(2 + 1);
7249371c9d4SSatish Balay       indices[4][0] = 2;
7259371c9d4SSatish Balay       indices[4][1] = 1;
7269371c9d4SSatish Balay       indices[4][2] = 0;
7279371c9d4SSatish Balay       indices[5][0] = 2;
7289371c9d4SSatish Balay       indices[5][1] = 2;
7299371c9d4SSatish Balay       indices[5][2] = 1;
730dda711d0SToby Isaac     }
731dda711d0SToby Isaac   } else {
732dda711d0SToby Isaac     PetscInt *subset, *work;
733dda711d0SToby Isaac 
7349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(k, &subset, k, &work));
735dda711d0SToby Isaac     for (i = 0; i < Nk; i++) {
736dda711d0SToby Isaac       PetscInt j, l, m;
737dda711d0SToby Isaac 
7389566063dSJacob Faibussowitsch       PetscCall(PetscDTEnumSubset(N, k, i, subset));
739dda711d0SToby Isaac       for (j = 0; j < k; j++) {
740dda711d0SToby Isaac         PetscInt  idx;
74128222859SToby Isaac         PetscBool flip = (PetscBool)(j & 1);
742dda711d0SToby Isaac 
743dda711d0SToby Isaac         for (l = 0, m = 0; l < k; l++) {
744dda711d0SToby Isaac           if (l != j) work[m++] = subset[l];
745dda711d0SToby Isaac         }
7469566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx));
747dda711d0SToby Isaac         indices[i * k + j][0] = idx;
748dda711d0SToby Isaac         indices[i * k + j][1] = i;
749dda711d0SToby Isaac         indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j];
750dda711d0SToby Isaac       }
751dda711d0SToby Isaac     }
7529566063dSJacob Faibussowitsch     PetscCall(PetscFree2(subset, work));
753dda711d0SToby Isaac   }
754dda711d0SToby Isaac   PetscFunctionReturn(0);
755dda711d0SToby Isaac }
756dda711d0SToby Isaac 
757fad4db65SToby Isaac /*@
75828222859SToby Isaac    PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form
759fad4db65SToby Isaac 
7604165533cSJose E. Roman    Input Parameters:
76128222859SToby Isaac +  N - the dimension of the vector space, N >= 0
76228222859SToby Isaac .  k - the degree k of the k-form w, 0 <= k <= N
76328222859SToby 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.
76428222859SToby Isaac -  w - a k-form, size [N choose k]
765fad4db65SToby Isaac 
7664165533cSJose E. Roman    Output Parameter:
76728222859SToby 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.
768fad4db65SToby Isaac 
769fad4db65SToby Isaac    Level: intermediate
770fad4db65SToby Isaac 
771db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()`
772fad4db65SToby Isaac @*/
773*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw)
774*d71ae5a4SJacob Faibussowitsch {
7751a989b97SToby Isaac   PetscInt Nk, i;
7761a989b97SToby Isaac 
7771a989b97SToby Isaac   PetscFunctionBegin;
7781dca8a05SBarry Smith   PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree");
7799566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
7801a989b97SToby Isaac   pow = pow % 4;
7811a989b97SToby Isaac   pow = (pow + 4) % 4; /* make non-negative */
7821a989b97SToby Isaac   /* pow is now 0, 1, 2, 3 */
7831a989b97SToby Isaac   if (N <= 3) {
7841a989b97SToby Isaac     if (pow & 1) {
7851a989b97SToby Isaac       PetscReal mult[3] = {1., -1., 1.};
7861a989b97SToby Isaac 
7871a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i];
7881a989b97SToby Isaac     } else {
7891a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = w[i];
7901a989b97SToby Isaac     }
7911a989b97SToby Isaac     if (pow > 1 && ((k * (N - k)) & 1)) {
7921a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
7931a989b97SToby Isaac     }
7941a989b97SToby Isaac   } else {
7951a989b97SToby Isaac     PetscInt *subset;
7961a989b97SToby Isaac 
7979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(N, &subset));
7981a989b97SToby Isaac     if (pow % 2) {
7991a989b97SToby Isaac       PetscInt l = (pow == 1) ? k : N - k;
8001a989b97SToby Isaac       for (i = 0; i < Nk; i++) {
8011a989b97SToby Isaac         PetscBool sOdd;
8021a989b97SToby Isaac         PetscInt  j, idx;
8031a989b97SToby Isaac 
8049566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, l, i, subset, &sOdd));
8059566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, l, subset, &idx));
8069566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, N - l, &subset[l], &j));
8071a989b97SToby Isaac         starw[j] = sOdd ? -w[idx] : w[idx];
8081a989b97SToby Isaac       }
8091a989b97SToby Isaac     } else {
8101a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = w[i];
8111a989b97SToby Isaac     }
8121a989b97SToby Isaac     /* star^2 = -1^(k * (N - k)) */
8131a989b97SToby Isaac     if (pow > 1 && (k * (N - k)) % 2) {
8141a989b97SToby Isaac       for (i = 0; i < Nk; i++) starw[i] = -starw[i];
8151a989b97SToby Isaac     }
8169566063dSJacob Faibussowitsch     PetscCall(PetscFree(subset));
8171a989b97SToby Isaac   }
8181a989b97SToby Isaac   PetscFunctionReturn(0);
8191a989b97SToby Isaac }
820