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) 15*dce8aebaSBarry Smith This action is implemented as `PetscDTAltVApply()`. 1629a920c6SToby Isaac 1729a920c6SToby Isaac The k-forms on a vector space form a vector space themselves, Alt^k V. The dimension of Alt^k V, if V is N dimensional, is N choose k. (This 1829a920c6SToby Isaac shows that for an N dimensional space, only 0 <= k <= N are valid form degrees.) 1929a920c6SToby Isaac The standard basis for Alt^k V, used in PetscDTAltV, has one basis k-form for each ordered subset of k coordinates of the N dimensional space: 2029a920c6SToby Isaac For example, if the coordinate directions of a four dimensional space are (t, x, y, z), then there are 4 choose 2 = 6 ordered subsets of two coordinates. 2129a920c6SToby Isaac They are, in lexicographic order, (t, x), (t, y), (t, z), (x, y), (x, z) and (y, z). PetscDTAltV also orders the basis of Alt^k V lexicographically 2229a920c6SToby Isaac by the associated subsets. 2329a920c6SToby Isaac 2429a920c6SToby Isaac The unit basis k-form associated with coordinates (c_1, ..., c_k) acts on a set of k vectors (v_1, ..., v_k) by creating a square matrix V where 2529a920c6SToby Isaac V[i,j] = v_i[c_j] and taking the determinant of V. 2629a920c6SToby Isaac 2729a920c6SToby Isaac If j + k <= N, then a j-form f and a k-form g can be multiplied to create a (j+k)-form using the wedge or exterior product, (f wedge g). 2829a920c6SToby Isaac This is an anticommutative product, (f wedge g) = -(g wedge f). It is sufficient to describe the wedge product of two basis forms. 2929a920c6SToby Isaac Let f be the basis j-form associated with coordinates (f_1,...,f_j) and g be the basis k-form associated with coordinates (g_1,...,g_k): 3029a920c6SToby Isaac - If there is any coordinate in both sets, then (f wedge g) = 0. 3129a920c6SToby Isaac - Otherwise, (f wedge g) is a multiple of the basis (j+k)-form h associated with (f_1,...,f_j,g_1,...,g_k). 3229a920c6SToby Isaac - In fact it is equal to either h or -h depending on how (f_1,...,f_j,g_1,...,g_k) compares to the same list of coordinates given in ascending order: if it is an even permutation of that list, then (f wedge g) = h, otherwise (f wedge g) = -h. 33*dce8aebaSBarry Smith The wedge product is implemented for either two inputs (f and g) in `PetscDTAltVWedge()`, or for one (just f, giving a 34*dce8aebaSBarry Smith matrix to multiply against multiple choices of g) in `PetscDTAltVWedgeMatrix()`. 3529a920c6SToby Isaac 3629a920c6SToby Isaac If k > 0, a k-form w and a vector v can combine to make a (k-1)-formm through the interior product, (w int v), 3729a920c6SToby Isaac defined by (w int v)(v_1,...,v_{k-1}) = w(v,v_1,...,v_{k-1}). 3829a920c6SToby Isaac 3929a920c6SToby Isaac The interior product is implemented for either two inputs (w and v) in PetscDTAltVInterior, for one (just v, giving a 40*dce8aebaSBarry Smith matrix to multiply against multiple choices of w) in `PetscDTAltVInteriorMatrix()`, 41*dce8aebaSBarry Smith or for no inputs (giving the sparsity pattern of `PetscDTAltVInteriorMatrix()`) in `PetscDTAltVInteriorPattern()`. 4229a920c6SToby Isaac 4329a920c6SToby Isaac When there is a linear map L: V -> W from an N dimensional vector space to an M dimensional vector space, 4429a920c6SToby Isaac it induces the linear pullback map L^* : Alt^k W -> Alt^k V, defined by L^* w(v_1,...,v_k) = w(L v_1, ..., L v_k). 45*dce8aebaSBarry Smith The pullback is implemented as `PetscDTAltVPullback()` (acting on a known w) and `PetscDTAltVPullbackMatrix()` (creating a matrix that computes the actin of L^*). 4629a920c6SToby Isaac 4729a920c6SToby Isaac Alt^k V and Alt^(N-k) V have the same dimension, and the Hodge star operator maps between them. We note that Alt^N V is a one dimensional space, and its 4829a920c6SToby Isaac basis vector is sometime called vol. The Hodge star operator has the property that (f wedge (star g)) = (f,g) vol, where (f,g) is the simple inner product 4929a920c6SToby Isaac of the basis coefficients of f and g. 5029a920c6SToby Isaac Powers of the Hodge star operator can be applied with PetscDTAltVStar 5129a920c6SToby Isaac 526c877ef6SSatish Balay Level: intermediate 5329a920c6SToby Isaac 54db781477SPatrick Sanan .seealso: `PetscDTAltVApply()`, `PetscDTAltVWedge()`, `PetscDTAltVInterior()`, `PetscDTAltVPullback()`, `PetscDTAltVStar()` 5529a920c6SToby Isaac M*/ 5629a920c6SToby Isaac 57fad4db65SToby Isaac /*@ 5828222859SToby Isaac PetscDTAltVApply - Apply an a k-form (an alternating k-linear map) to a set of k N-dimensional vectors 59fad4db65SToby Isaac 604165533cSJose E. Roman Input Parameters: 6128222859SToby Isaac + N - the dimension of the vector space, N >= 0 6228222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 63*dce8aebaSBarry Smith . w - a k-form, size [N choose k] (each degree of freedom of a k-form is associated with a subset of k coordinates of the N-dimensional vectors. 64*dce8aebaSBarry Smith The degrees of freedom are ordered lexicographically by their associated subsets) 6528222859SToby Isaac - v - a set of k vectors of size N, size [k x N], each vector stored contiguously 66fad4db65SToby Isaac 674165533cSJose E. Roman Output Parameter: 68*dce8aebaSBarry Smith . wv - w(v_1,...,v_k) = \sum_i w_i * det(V_i): the degree of freedom w_i is associated with coordinates [s_{i,1},...,s_{i,k}], and the square matrix V_i has 69*dce8aebaSBarry Smith entry (j,k) given by the s_{i,k}'th coordinate of v_j 70fad4db65SToby Isaac 71fad4db65SToby Isaac Level: intermediate 72fad4db65SToby Isaac 73db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 74fad4db65SToby Isaac @*/ 75d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv) 76d71ae5a4SJacob Faibussowitsch { 771a989b97SToby Isaac PetscFunctionBegin; 7808401ef6SPierre Jolivet PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 791dca8a05SBarry Smith PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 801a989b97SToby Isaac if (N <= 3) { 811a989b97SToby Isaac if (!k) { 821a989b97SToby Isaac *wv = w[0]; 831a989b97SToby Isaac } else { 849371c9d4SSatish Balay if (N == 1) { 859371c9d4SSatish Balay *wv = w[0] * v[0]; 869371c9d4SSatish Balay } else if (N == 2) { 879371c9d4SSatish Balay if (k == 1) { 889371c9d4SSatish Balay *wv = w[0] * v[0] + w[1] * v[1]; 891a989b97SToby Isaac } else { 909371c9d4SSatish Balay *wv = w[0] * (v[0] * v[3] - v[1] * v[2]); 919371c9d4SSatish Balay } 921a989b97SToby Isaac } else { 939371c9d4SSatish Balay if (k == 1) { 949371c9d4SSatish Balay *wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2]; 959371c9d4SSatish Balay } else if (k == 2) { 969371c9d4SSatish Balay *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) + w[1] * (v[0] * v[5] - v[2] * v[3]) + w[2] * (v[1] * v[5] - v[2] * v[4]); 979371c9d4SSatish Balay } else { 989371c9d4SSatish Balay *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) + v[1] * (v[5] * v[6] - v[3] * v[8]) + v[2] * (v[3] * v[7] - v[4] * v[6])); 991a989b97SToby Isaac } 1001a989b97SToby Isaac } 1011a989b97SToby Isaac } 1021a989b97SToby Isaac } else { 1031a989b97SToby Isaac PetscInt Nk, Nf; 104fad4db65SToby Isaac PetscInt *subset, *perm; 1051a989b97SToby Isaac PetscInt i, j, l; 1061a989b97SToby Isaac PetscReal sum = 0.; 1071a989b97SToby Isaac 1089566063dSJacob Faibussowitsch PetscCall(PetscDTFactorialInt(k, &Nf)); 1099566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &perm)); 1111a989b97SToby Isaac for (i = 0; i < Nk; i++) { 1121a989b97SToby Isaac PetscReal subsum = 0.; 1131a989b97SToby Isaac 1149566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 1151a989b97SToby Isaac for (j = 0; j < Nf; j++) { 1161a989b97SToby Isaac PetscBool permOdd; 1171a989b97SToby Isaac PetscReal prod; 1181a989b97SToby Isaac 1199566063dSJacob Faibussowitsch PetscCall(PetscDTEnumPerm(k, j, perm, &permOdd)); 1201a989b97SToby Isaac prod = permOdd ? -1. : 1.; 121ad540459SPierre Jolivet for (l = 0; l < k; l++) prod *= v[perm[l] * N + subset[l]]; 1221a989b97SToby Isaac subsum += prod; 1231a989b97SToby Isaac } 1241a989b97SToby Isaac sum += w[i] * subsum; 1251a989b97SToby Isaac } 1269566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, perm)); 1271a989b97SToby Isaac *wv = sum; 1281a989b97SToby Isaac } 1291a989b97SToby Isaac PetscFunctionReturn(0); 1301a989b97SToby Isaac } 1311a989b97SToby Isaac 132fad4db65SToby Isaac /*@ 13328222859SToby Isaac PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form 134fad4db65SToby Isaac 1354165533cSJose E. Roman Input Parameters: 13628222859SToby Isaac + N - the dimension of the vector space, N >= 0 13728222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N 13828222859SToby Isaac . k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N 13928222859SToby Isaac . a - a j-form, size [N choose j] 14028222859SToby Isaac - b - a k-form, size [N choose k] 141fad4db65SToby Isaac 1424165533cSJose E. Roman Output Parameter: 14328222859SToby Isaac . awedgeb - the (j+k)-form a wedge b, size [N choose (j+k)]: (a wedge b)(v_1,...,v_{j+k}) = \sum_{s} sign(s) a(v_{s_1},...,v_{s_j}) b(v_{s_{j+1}},...,v_{s_{j+k}}), 14428222859SToby Isaac where the sum is over permutations s such that s_1 < s_2 < ... < s_j and s_{j+1} < s_{j+2} < ... < s_{j+k}. 145fad4db65SToby Isaac 146fad4db65SToby Isaac Level: intermediate 147fad4db65SToby Isaac 148db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVWedgeMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 149fad4db65SToby Isaac @*/ 150d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) 151d71ae5a4SJacob Faibussowitsch { 1521a989b97SToby Isaac PetscInt i; 1531a989b97SToby Isaac 1541a989b97SToby Isaac PetscFunctionBegin; 15508401ef6SPierre Jolivet PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 1561dca8a05SBarry Smith PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 1571dca8a05SBarry Smith PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 1581a989b97SToby Isaac if (N <= 3) { 1591a989b97SToby Isaac PetscInt Njk; 1601a989b97SToby Isaac 1619566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 1629371c9d4SSatish Balay if (!j) { 163ad540459SPierre Jolivet for (i = 0; i < Njk; i++) awedgeb[i] = a[0] * b[i]; 1649371c9d4SSatish Balay } else if (!k) { 165ad540459SPierre Jolivet for (i = 0; i < Njk; i++) awedgeb[i] = a[i] * b[0]; 1669371c9d4SSatish Balay } else { 1679371c9d4SSatish Balay if (N == 2) { 1689371c9d4SSatish Balay awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 1699371c9d4SSatish Balay } else { 1701a989b97SToby Isaac if (j + k == 2) { 1711a989b97SToby Isaac awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 1721a989b97SToby Isaac awedgeb[1] = a[0] * b[2] - a[2] * b[0]; 1731a989b97SToby Isaac awedgeb[2] = a[1] * b[2] - a[2] * b[1]; 1741a989b97SToby Isaac } else { 1751a989b97SToby Isaac awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0]; 1761a989b97SToby Isaac } 1771a989b97SToby Isaac } 1781a989b97SToby Isaac } 1791a989b97SToby Isaac } else { 1801a989b97SToby Isaac PetscInt Njk; 1811a989b97SToby Isaac PetscInt JKj; 1821a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 1831a989b97SToby Isaac PetscInt i; 1841a989b97SToby Isaac 1859566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 1869566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 1879566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk)); 1881a989b97SToby Isaac for (i = 0; i < Njk; i++) { 1891a989b97SToby Isaac PetscReal sum = 0.; 1901a989b97SToby Isaac PetscInt l; 1911a989b97SToby Isaac 1929566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, j + k, i, subset)); 1931a989b97SToby Isaac for (l = 0; l < JKj; l++) { 1941a989b97SToby Isaac PetscBool jkOdd; 1951a989b97SToby Isaac PetscInt m, jInd, kInd; 1961a989b97SToby Isaac 1979566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd)); 198ad540459SPierre Jolivet for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]]; 199ad540459SPierre Jolivet for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]]; 2009566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd)); 2019566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd)); 2021a989b97SToby Isaac sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]); 2031a989b97SToby Isaac } 2041a989b97SToby Isaac awedgeb[i] = sum; 2051a989b97SToby Isaac } 2069566063dSJacob Faibussowitsch PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk)); 2071a989b97SToby Isaac } 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 225db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 226fad4db65SToby Isaac @*/ 227d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat) 228d71ae5a4SJacob Faibussowitsch { 2291a989b97SToby Isaac PetscInt i; 2301a989b97SToby Isaac 2311a989b97SToby Isaac PetscFunctionBegin; 23208401ef6SPierre Jolivet PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 2331dca8a05SBarry Smith PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 2341dca8a05SBarry Smith PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 2351a989b97SToby Isaac if (N <= 3) { 2361a989b97SToby Isaac PetscInt Njk; 2371a989b97SToby Isaac 2389566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 2391a989b97SToby Isaac if (!j) { 240ad540459SPierre Jolivet for (i = 0; i < Njk * Njk; i++) awedgeMat[i] = 0.; 241ad540459SPierre Jolivet for (i = 0; i < Njk; i++) awedgeMat[i * (Njk + 1)] = a[0]; 2421a989b97SToby Isaac } else if (!k) { 243ad540459SPierre Jolivet for (i = 0; i < Njk; i++) awedgeMat[i] = a[i]; 2441a989b97SToby Isaac } else { 2451a989b97SToby Isaac if (N == 2) { 2469371c9d4SSatish Balay awedgeMat[0] = -a[1]; 2479371c9d4SSatish Balay awedgeMat[1] = a[0]; 2481a989b97SToby Isaac } else { 2491a989b97SToby Isaac if (j + k == 2) { 2509371c9d4SSatish Balay awedgeMat[0] = -a[1]; 2519371c9d4SSatish Balay awedgeMat[1] = a[0]; 2529371c9d4SSatish Balay awedgeMat[2] = 0.; 2539371c9d4SSatish Balay awedgeMat[3] = -a[2]; 2549371c9d4SSatish Balay awedgeMat[4] = 0.; 2559371c9d4SSatish Balay awedgeMat[5] = a[0]; 2569371c9d4SSatish Balay awedgeMat[6] = 0.; 2579371c9d4SSatish Balay awedgeMat[7] = -a[2]; 2589371c9d4SSatish Balay awedgeMat[8] = a[1]; 2591a989b97SToby Isaac } else { 2609371c9d4SSatish Balay awedgeMat[0] = a[2]; 2619371c9d4SSatish Balay awedgeMat[1] = -a[1]; 2629371c9d4SSatish Balay awedgeMat[2] = a[0]; 2631a989b97SToby Isaac } 2641a989b97SToby Isaac } 2651a989b97SToby Isaac } 2661a989b97SToby Isaac } else { 2671a989b97SToby Isaac PetscInt Njk; 2681a989b97SToby Isaac PetscInt Nk; 2691a989b97SToby Isaac PetscInt JKj, i; 2701a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 2711a989b97SToby Isaac 2729566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 2739566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 2749566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 2759566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk)); 2761a989b97SToby Isaac for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.; 2771a989b97SToby Isaac for (i = 0; i < Njk; i++) { 2781a989b97SToby Isaac PetscInt l; 2791a989b97SToby Isaac 2809566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, j + k, i, subset)); 2811a989b97SToby Isaac for (l = 0; l < JKj; l++) { 2821a989b97SToby Isaac PetscBool jkOdd; 2831a989b97SToby Isaac PetscInt m, jInd, kInd; 2841a989b97SToby Isaac 2859566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd)); 286ad540459SPierre Jolivet for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]]; 287ad540459SPierre Jolivet for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]]; 2889566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd)); 2899566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd)); 2901a989b97SToby Isaac awedgeMat[i * Nk + kInd] += jkOdd ? -a[jInd] : a[jInd]; 2911a989b97SToby Isaac } 2921a989b97SToby Isaac } 2939566063dSJacob Faibussowitsch PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk)); 2941a989b97SToby Isaac } 2951a989b97SToby Isaac PetscFunctionReturn(0); 2961a989b97SToby Isaac } 2971a989b97SToby Isaac 298fad4db65SToby Isaac /*@ 29928222859SToby Isaac PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space 300fad4db65SToby Isaac 3014165533cSJose E. Roman Input Parameters: 30228222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, M >= 0 30328222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, N >= 0 30428222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format 305*dce8aebaSBarry Smith . k - the *signed* degree k of the |k|-form w, -(min(M,N)) <= k <= min(M,N). A negative form degree indicates that the pullback should be conjugated by 306*dce8aebaSBarry Smith the Hodge star operator (see note). 30728222859SToby Isaac - w - a |k|-form in the image space, size [M choose |k|] 308fad4db65SToby Isaac 3094165533cSJose E. Roman Output Parameter: 31028222859SToby Isaac . Lstarw - the pullback of w to a |k|-form in the origin space, size [N choose |k|]: (Lstarw)(v_1,...v_k) = w(L*v_1,...,L*v_k). 311fad4db65SToby Isaac 312fad4db65SToby Isaac Level: intermediate 313fad4db65SToby Isaac 314a5b23f4aSJose 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, 31528222859SToby 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, 316fad4db65SToby Isaac then the inverse Hodge star transformation. 317fad4db65SToby Isaac 318db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullbackMatrix()`, `PetscDTAltVStar()` 319fad4db65SToby Isaac @*/ 320d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw) 321d71ae5a4SJacob Faibussowitsch { 3221a989b97SToby Isaac PetscInt i, j, Nk, Mk; 3231a989b97SToby Isaac 3241a989b97SToby Isaac PetscFunctionBegin; 3251dca8a05SBarry Smith PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 3261dca8a05SBarry Smith PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 3271a989b97SToby Isaac if (N <= 3 && M <= 3) { 3289566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 3299566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 3301a989b97SToby Isaac if (!k) { 3311a989b97SToby Isaac Lstarw[0] = w[0]; 3321a989b97SToby Isaac } else if (k == 1) { 3331a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3341a989b97SToby Isaac PetscReal sum = 0.; 3351a989b97SToby Isaac 336ad540459SPierre Jolivet for (j = 0; j < Mk; j++) sum += L[j * Nk + i] * w[j]; 3371a989b97SToby Isaac Lstarw[i] = sum; 3381a989b97SToby Isaac } 3391a989b97SToby Isaac } else if (k == -1) { 3401a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 3411a989b97SToby Isaac 3421a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3431a989b97SToby Isaac PetscReal sum = 0.; 3441a989b97SToby Isaac 345ad540459SPierre Jolivet for (j = 0; j < Mk; j++) sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j]; 3461a989b97SToby Isaac Lstarw[i] = mult[i] * sum; 3471a989b97SToby Isaac } 3481a989b97SToby Isaac } else if (k == 2) { 3499371c9d4SSatish Balay PetscInt pairs[3][2] = { 3509371c9d4SSatish Balay {0, 1}, 3519371c9d4SSatish Balay {0, 2}, 3529371c9d4SSatish Balay {1, 2} 3539371c9d4SSatish Balay }; 3541a989b97SToby Isaac 3551a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3561a989b97SToby Isaac PetscReal sum = 0.; 357ad540459SPierre 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]; 3581a989b97SToby Isaac Lstarw[i] = sum; 3591a989b97SToby Isaac } 3601a989b97SToby Isaac } else if (k == -2) { 3619371c9d4SSatish Balay PetscInt pairs[3][2] = { 3629371c9d4SSatish Balay {1, 2}, 3639371c9d4SSatish Balay {2, 0}, 3649371c9d4SSatish Balay {0, 1} 3659371c9d4SSatish Balay }; 3661a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 3671a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 3681a989b97SToby Isaac 3691a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3701a989b97SToby Isaac PetscReal sum = 0.; 3711a989b97SToby Isaac 372ad540459SPierre 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]; 3731a989b97SToby Isaac Lstarw[i] = sum; 3741a989b97SToby Isaac } 3751a989b97SToby Isaac } else { 3769371c9d4SSatish 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]); 3771a989b97SToby Isaac 378ad540459SPierre Jolivet 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]; 413ad540459SPierre Jolivet for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 4141a989b97SToby Isaac Lstarw[j] += prod; 4151a989b97SToby Isaac } 4161a989b97SToby Isaac } 4171a989b97SToby Isaac } 4181a989b97SToby Isaac if (negative) { 4191a989b97SToby Isaac PetscReal *sLsw; 4201a989b97SToby Isaac 4219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk, &sLsw)); 4229566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw)); 4231a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i]; 4249566063dSJacob Faibussowitsch PetscCall(PetscFree(sLsw)); 4251a989b97SToby Isaac } 4269566063dSJacob Faibussowitsch PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv)); 4279566063dSJacob Faibussowitsch PetscCall(PetscFree(walloc)); 4281a989b97SToby Isaac } 4291a989b97SToby Isaac PetscFunctionReturn(0); 4301a989b97SToby Isaac } 4311a989b97SToby Isaac 432fad4db65SToby Isaac /*@ 433fad4db65SToby Isaac PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation 434fad4db65SToby Isaac 4354165533cSJose E. Roman Input Parameters: 43628222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, N >= 0 43728222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, M >= 0 43828222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format 439*dce8aebaSBarry Smith - k - the *signed* degree k of the |k|-forms on which Lstar acts, -(min(M,N)) <= k <= min(M,N). 440*dce8aebaSBarry Smith A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note in `PetscDTAltvPullback()`) 441fad4db65SToby Isaac 4424165533cSJose E. Roman Output Parameter: 44328222859SToby Isaac . Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w 444fad4db65SToby Isaac 445fad4db65SToby Isaac Level: intermediate 446fad4db65SToby Isaac 447db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVStar()` 448fad4db65SToby Isaac @*/ 449d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) 450d71ae5a4SJacob Faibussowitsch { 4511a989b97SToby Isaac PetscInt Nk, Mk, Nf, i, j, l, p; 4521a989b97SToby Isaac PetscReal *Lw, *Lwv; 4531a989b97SToby Isaac PetscInt *subsetw, *subsetv; 454fad4db65SToby Isaac PetscInt *perm; 4551a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 4561a989b97SToby Isaac 4571a989b97SToby Isaac PetscFunctionBegin; 4581dca8a05SBarry Smith PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 4591dca8a05SBarry Smith PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 4601a989b97SToby Isaac if (N <= 3 && M <= 3) { 4611a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 4621a989b97SToby Isaac 4639566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 4649566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 4651a989b97SToby Isaac if (!k) { 4661a989b97SToby Isaac Lstar[0] = 1.; 4671a989b97SToby Isaac } else if (k == 1) { 4689371c9d4SSatish Balay for (i = 0; i < Nk; i++) { 469ad540459SPierre Jolivet for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[j * Nk + i]; 4709371c9d4SSatish Balay } 4711a989b97SToby Isaac } else if (k == -1) { 4721a989b97SToby Isaac for (i = 0; i < Nk; i++) { 473ad540459SPierre Jolivet for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j]; 4741a989b97SToby Isaac } 4751a989b97SToby Isaac } else if (k == 2) { 4769371c9d4SSatish Balay PetscInt pairs[3][2] = { 4779371c9d4SSatish Balay {0, 1}, 4789371c9d4SSatish Balay {0, 2}, 4799371c9d4SSatish Balay {1, 2} 4809371c9d4SSatish Balay }; 4811a989b97SToby Isaac 4821a989b97SToby Isaac for (i = 0; i < Nk; i++) { 483ad540459SPierre 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]]; 4841a989b97SToby Isaac } 4851a989b97SToby Isaac } else if (k == -2) { 4869371c9d4SSatish Balay PetscInt pairs[3][2] = { 4879371c9d4SSatish Balay {1, 2}, 4889371c9d4SSatish Balay {2, 0}, 4899371c9d4SSatish Balay {0, 1} 4909371c9d4SSatish Balay }; 4911a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 4921a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 4931a989b97SToby Isaac 4941a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4951a989b97SToby Isaac for (j = 0; j < Mk; j++) { 4969371c9d4SSatish 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]]; 4971a989b97SToby Isaac } 4981a989b97SToby Isaac } 4991a989b97SToby Isaac } else { 5009371c9d4SSatish 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]); 5011a989b97SToby Isaac 502ad540459SPierre Jolivet for (i = 0; i < Nk; i++) Lstar[i] = detL; 5031a989b97SToby Isaac } 5041a989b97SToby Isaac } else { 5051a989b97SToby Isaac if (k < 0) { 5061a989b97SToby Isaac negative = PETSC_TRUE; 5071a989b97SToby Isaac k = -k; 5081a989b97SToby Isaac } 5099566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 5109566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 5119566063dSJacob Faibussowitsch PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf)); 5129566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv)); 5131a989b97SToby Isaac for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.; 5141a989b97SToby Isaac for (i = 0; i < Mk; i++) { 5151a989b97SToby Isaac PetscBool iOdd; 5161a989b97SToby Isaac PetscInt iidx, jidx; 5171a989b97SToby Isaac 5189566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(M, k, i, subsetw, &iOdd)); 5191a989b97SToby Isaac iidx = negative ? Mk - 1 - i : i; 52028222859SToby Isaac iOdd = negative ? (PetscBool)(iOdd ^ ((k * (M - k)) & 1)) : PETSC_FALSE; 5211a989b97SToby Isaac for (j = 0; j < Nk; j++) { 5221a989b97SToby Isaac PetscBool jOdd; 5231a989b97SToby Isaac 5249566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(N, k, j, subsetv, &jOdd)); 5251a989b97SToby Isaac jidx = negative ? Nk - 1 - j : j; 52628222859SToby Isaac jOdd = negative ? (PetscBool)(iOdd ^ jOdd ^ ((k * (N - k)) & 1)) : PETSC_FALSE; 5271a989b97SToby Isaac for (p = 0; p < Nf; p++) { 5281a989b97SToby Isaac PetscReal prod; 5291a989b97SToby Isaac PetscBool isOdd; 5301a989b97SToby Isaac 5319566063dSJacob Faibussowitsch PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd)); 53228222859SToby Isaac isOdd = (PetscBool)(isOdd ^ jOdd); 5331a989b97SToby Isaac prod = isOdd ? -1. : 1.; 534ad540459SPierre Jolivet for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 5351a989b97SToby Isaac Lstar[jidx * Mk + iidx] += prod; 5361a989b97SToby Isaac } 5371a989b97SToby Isaac } 5381a989b97SToby Isaac } 5399566063dSJacob Faibussowitsch PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv)); 5401a989b97SToby Isaac } 5411a989b97SToby Isaac PetscFunctionReturn(0); 5421a989b97SToby Isaac } 5431a989b97SToby Isaac 544fad4db65SToby Isaac /*@ 54528222859SToby Isaac PetscDTAltVInterior - Compute the interior product of a k-form with a vector 546fad4db65SToby Isaac 5474165533cSJose E. Roman Input Parameters: 54828222859SToby Isaac + N - the dimension of the vector space, N >= 0 54928222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 55028222859SToby Isaac . w - a k-form, size [N choose k] 55128222859SToby Isaac - v - an N dimensional vector 552fad4db65SToby Isaac 5534165533cSJose E. Roman Output Parameter: 55428222859SToby 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}). 555fad4db65SToby Isaac 556fad4db65SToby Isaac Level: intermediate 557fad4db65SToby Isaac 558db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 559fad4db65SToby Isaac @*/ 560d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv) 561d71ae5a4SJacob Faibussowitsch { 5621a989b97SToby Isaac PetscInt i, Nk, Nkm; 5631a989b97SToby Isaac 5641a989b97SToby Isaac PetscFunctionBegin; 5651dca8a05SBarry Smith PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 5669566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 5679566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 5681a989b97SToby Isaac if (N <= 3) { 5691a989b97SToby Isaac if (k == 1) { 5701a989b97SToby Isaac PetscReal sum = 0.; 5711a989b97SToby Isaac 572ad540459SPierre Jolivet for (i = 0; i < N; i++) sum += w[i] * v[i]; 5731a989b97SToby Isaac wIntv[0] = sum; 5741a989b97SToby Isaac } else if (k == N) { 5751a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 5761a989b97SToby Isaac 577ad540459SPierre Jolivet for (i = 0; i < N; i++) wIntv[N - 1 - i] = w[0] * v[i] * mult[i]; 5781a989b97SToby Isaac } else { 5791a989b97SToby Isaac wIntv[0] = -w[0] * v[1] - w[1] * v[2]; 5801a989b97SToby Isaac wIntv[1] = w[0] * v[0] - w[2] * v[2]; 5811a989b97SToby Isaac wIntv[2] = w[1] * v[0] + w[2] * v[1]; 5821a989b97SToby Isaac } 5831a989b97SToby Isaac } else { 5841a989b97SToby Isaac PetscInt *subset, *work; 5851a989b97SToby Isaac 5869566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &work)); 5871a989b97SToby Isaac for (i = 0; i < Nkm; i++) wIntv[i] = 0.; 5881a989b97SToby Isaac for (i = 0; i < Nk; i++) { 5891a989b97SToby Isaac PetscInt j, l, m; 5901a989b97SToby Isaac 5919566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 5921a989b97SToby Isaac for (j = 0; j < k; j++) { 5931a989b97SToby Isaac PetscInt idx; 59428222859SToby Isaac PetscBool flip = (PetscBool)(j & 1); 5951a989b97SToby Isaac 5961a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 5971a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 5981a989b97SToby Isaac } 5999566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 6001a989b97SToby Isaac wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]); 6011a989b97SToby Isaac } 6021a989b97SToby Isaac } 6039566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, work)); 6041a989b97SToby Isaac } 6051a989b97SToby Isaac PetscFunctionReturn(0); 6061a989b97SToby Isaac } 6071a989b97SToby Isaac 608fad4db65SToby Isaac /*@ 60928222859SToby Isaac PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector 610fad4db65SToby Isaac 6114165533cSJose E. Roman Input Parameters: 61228222859SToby Isaac + N - the dimension of the vector space, N >= 0 61328222859SToby Isaac . k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N 61428222859SToby Isaac - v - an N dimensional vector 615fad4db65SToby Isaac 6164165533cSJose E. Roman Output Parameter: 617fad4db65SToby Isaac . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v) 618fad4db65SToby Isaac 619fad4db65SToby Isaac Level: intermediate 620fad4db65SToby Isaac 621db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 622fad4db65SToby Isaac @*/ 623d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat) 624d71ae5a4SJacob Faibussowitsch { 6251a989b97SToby Isaac PetscInt i, Nk, Nkm; 6261a989b97SToby Isaac 6271a989b97SToby Isaac PetscFunctionBegin; 6281dca8a05SBarry Smith PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 6299566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 6309566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 6311a989b97SToby Isaac if (N <= 3) { 6321a989b97SToby Isaac if (k == 1) { 6331a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[i] = v[i]; 6341a989b97SToby Isaac } else if (k == N) { 6351a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 6361a989b97SToby Isaac 6371a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i]; 6381a989b97SToby Isaac } else { 6399371c9d4SSatish Balay intvMat[0] = -v[1]; 6409371c9d4SSatish Balay intvMat[1] = -v[2]; 6419371c9d4SSatish Balay intvMat[2] = 0.; 6429371c9d4SSatish Balay intvMat[3] = v[0]; 6439371c9d4SSatish Balay intvMat[4] = 0.; 6449371c9d4SSatish Balay intvMat[5] = -v[2]; 6459371c9d4SSatish Balay intvMat[6] = 0.; 6469371c9d4SSatish Balay intvMat[7] = v[0]; 6479371c9d4SSatish Balay intvMat[8] = v[1]; 6481a989b97SToby Isaac } 6491a989b97SToby Isaac } else { 6501a989b97SToby Isaac PetscInt *subset, *work; 6511a989b97SToby Isaac 6529566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &work)); 6531a989b97SToby Isaac for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.; 6541a989b97SToby Isaac for (i = 0; i < Nk; i++) { 6551a989b97SToby Isaac PetscInt j, l, m; 6561a989b97SToby Isaac 6579566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 6581a989b97SToby Isaac for (j = 0; j < k; j++) { 6591a989b97SToby Isaac PetscInt idx; 66028222859SToby Isaac PetscBool flip = (PetscBool)(j & 1); 6611a989b97SToby Isaac 6621a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 6631a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 6641a989b97SToby Isaac } 6659566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 6661a989b97SToby Isaac intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]]; 6671a989b97SToby Isaac } 6681a989b97SToby Isaac } 6699566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, work)); 6701a989b97SToby Isaac } 6711a989b97SToby Isaac PetscFunctionReturn(0); 6721a989b97SToby Isaac } 6731a989b97SToby Isaac 674fad4db65SToby Isaac /*@ 675*dce8aebaSBarry Smith PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in `PetscDTAltVInteriorMatrix()` 676fad4db65SToby Isaac 6774165533cSJose E. Roman Input Parameters: 67828222859SToby Isaac + N - the dimension of the vector space, N >= 0 679*dce8aebaSBarry Smith - k - the degree of the k-forms on which intvMat from `PetscDTAltVInteriorMatrix()` acts, 0 <= k <= N. 680fad4db65SToby Isaac 6814165533cSJose E. Roman Output Parameter: 68228222859SToby Isaac . indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k 68328222859SToby 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 68428222859SToby Isaac coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1) 685fad4db65SToby Isaac 686fad4db65SToby Isaac Level: intermediate 687fad4db65SToby Isaac 688*dce8aebaSBarry Smith Note: 689*dce8aebaSBarry Smith This function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential 690fad4db65SToby Isaac 691db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 692fad4db65SToby Isaac @*/ 693d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3]) 694d71ae5a4SJacob Faibussowitsch { 695dda711d0SToby Isaac PetscInt i, Nk, Nkm; 696dda711d0SToby Isaac 697dda711d0SToby Isaac PetscFunctionBegin; 6981dca8a05SBarry Smith PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 6999566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 7009566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 701dda711d0SToby Isaac if (N <= 3) { 702dda711d0SToby Isaac if (k == 1) { 703dda711d0SToby Isaac for (i = 0; i < N; i++) { 704dda711d0SToby Isaac indices[i][0] = 0; 705dda711d0SToby Isaac indices[i][1] = i; 706dda711d0SToby Isaac indices[i][2] = i; 707dda711d0SToby Isaac } 708dda711d0SToby Isaac } else if (k == N) { 709dda711d0SToby Isaac PetscInt val[3] = {0, -2, 2}; 710dda711d0SToby Isaac 711dda711d0SToby Isaac for (i = 0; i < N; i++) { 712dda711d0SToby Isaac indices[i][0] = N - 1 - i; 713dda711d0SToby Isaac indices[i][1] = 0; 714dda711d0SToby Isaac indices[i][2] = val[i]; 715dda711d0SToby Isaac } 716dda711d0SToby Isaac } else { 7179371c9d4SSatish Balay indices[0][0] = 0; 7189371c9d4SSatish Balay indices[0][1] = 0; 7199371c9d4SSatish Balay indices[0][2] = -(1 + 1); 7209371c9d4SSatish Balay indices[1][0] = 0; 7219371c9d4SSatish Balay indices[1][1] = 1; 7229371c9d4SSatish Balay indices[1][2] = -(2 + 1); 7239371c9d4SSatish Balay indices[2][0] = 1; 7249371c9d4SSatish Balay indices[2][1] = 0; 7259371c9d4SSatish Balay indices[2][2] = 0; 7269371c9d4SSatish Balay indices[3][0] = 1; 7279371c9d4SSatish Balay indices[3][1] = 2; 7289371c9d4SSatish Balay indices[3][2] = -(2 + 1); 7299371c9d4SSatish Balay indices[4][0] = 2; 7309371c9d4SSatish Balay indices[4][1] = 1; 7319371c9d4SSatish Balay indices[4][2] = 0; 7329371c9d4SSatish Balay indices[5][0] = 2; 7339371c9d4SSatish Balay indices[5][1] = 2; 7349371c9d4SSatish Balay indices[5][2] = 1; 735dda711d0SToby Isaac } 736dda711d0SToby Isaac } else { 737dda711d0SToby Isaac PetscInt *subset, *work; 738dda711d0SToby Isaac 7399566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &work)); 740dda711d0SToby Isaac for (i = 0; i < Nk; i++) { 741dda711d0SToby Isaac PetscInt j, l, m; 742dda711d0SToby Isaac 7439566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 744dda711d0SToby Isaac for (j = 0; j < k; j++) { 745dda711d0SToby Isaac PetscInt idx; 74628222859SToby Isaac PetscBool flip = (PetscBool)(j & 1); 747dda711d0SToby Isaac 748dda711d0SToby Isaac for (l = 0, m = 0; l < k; l++) { 749dda711d0SToby Isaac if (l != j) work[m++] = subset[l]; 750dda711d0SToby Isaac } 7519566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 752dda711d0SToby Isaac indices[i * k + j][0] = idx; 753dda711d0SToby Isaac indices[i * k + j][1] = i; 754dda711d0SToby Isaac indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j]; 755dda711d0SToby Isaac } 756dda711d0SToby Isaac } 7579566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, work)); 758dda711d0SToby Isaac } 759dda711d0SToby Isaac PetscFunctionReturn(0); 760dda711d0SToby Isaac } 761dda711d0SToby Isaac 762fad4db65SToby Isaac /*@ 76328222859SToby Isaac PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form 764fad4db65SToby Isaac 7654165533cSJose E. Roman Input Parameters: 76628222859SToby Isaac + N - the dimension of the vector space, N >= 0 76728222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 76828222859SToby 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. 76928222859SToby Isaac - w - a k-form, size [N choose k] 770fad4db65SToby Isaac 7714165533cSJose E. Roman Output Parameter: 77228222859SToby 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. 773fad4db65SToby Isaac 774fad4db65SToby Isaac Level: intermediate 775fad4db65SToby Isaac 776db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 777fad4db65SToby Isaac @*/ 778d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw) 779d71ae5a4SJacob Faibussowitsch { 7801a989b97SToby Isaac PetscInt Nk, i; 7811a989b97SToby Isaac 7821a989b97SToby Isaac PetscFunctionBegin; 7831dca8a05SBarry Smith PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 7849566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 7851a989b97SToby Isaac pow = pow % 4; 7861a989b97SToby Isaac pow = (pow + 4) % 4; /* make non-negative */ 7871a989b97SToby Isaac /* pow is now 0, 1, 2, 3 */ 7881a989b97SToby Isaac if (N <= 3) { 7891a989b97SToby Isaac if (pow & 1) { 7901a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 7911a989b97SToby Isaac 7921a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i]; 7931a989b97SToby Isaac } else { 7941a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 7951a989b97SToby Isaac } 7961a989b97SToby Isaac if (pow > 1 && ((k * (N - k)) & 1)) { 7971a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 7981a989b97SToby Isaac } 7991a989b97SToby Isaac } else { 8001a989b97SToby Isaac PetscInt *subset; 8011a989b97SToby Isaac 8029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &subset)); 8031a989b97SToby Isaac if (pow % 2) { 8041a989b97SToby Isaac PetscInt l = (pow == 1) ? k : N - k; 8051a989b97SToby Isaac for (i = 0; i < Nk; i++) { 8061a989b97SToby Isaac PetscBool sOdd; 8071a989b97SToby Isaac PetscInt j, idx; 8081a989b97SToby Isaac 8099566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(N, l, i, subset, &sOdd)); 8109566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, l, subset, &idx)); 8119566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, N - l, &subset[l], &j)); 8121a989b97SToby Isaac starw[j] = sOdd ? -w[idx] : w[idx]; 8131a989b97SToby Isaac } 8141a989b97SToby Isaac } else { 8151a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 8161a989b97SToby Isaac } 8171a989b97SToby Isaac /* star^2 = -1^(k * (N - k)) */ 8181a989b97SToby Isaac if (pow > 1 && (k * (N - k)) % 2) { 8191a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 8201a989b97SToby Isaac } 8219566063dSJacob Faibussowitsch PetscCall(PetscFree(subset)); 8221a989b97SToby Isaac } 8231a989b97SToby Isaac PetscFunctionReturn(0); 8241a989b97SToby Isaac } 825