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*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv) { 741a989b97SToby Isaac PetscFunctionBegin; 7508401ef6SPierre Jolivet PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 761dca8a05SBarry Smith PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 771a989b97SToby Isaac if (N <= 3) { 781a989b97SToby Isaac if (!k) { 791a989b97SToby Isaac *wv = w[0]; 801a989b97SToby Isaac } else { 81*9371c9d4SSatish Balay if (N == 1) { 82*9371c9d4SSatish Balay *wv = w[0] * v[0]; 83*9371c9d4SSatish Balay } else if (N == 2) { 84*9371c9d4SSatish Balay if (k == 1) { 85*9371c9d4SSatish Balay *wv = w[0] * v[0] + w[1] * v[1]; 861a989b97SToby Isaac } else { 87*9371c9d4SSatish Balay *wv = w[0] * (v[0] * v[3] - v[1] * v[2]); 88*9371c9d4SSatish Balay } 891a989b97SToby Isaac } else { 90*9371c9d4SSatish Balay if (k == 1) { 91*9371c9d4SSatish Balay *wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2]; 92*9371c9d4SSatish Balay } else if (k == 2) { 93*9371c9d4SSatish 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]); 94*9371c9d4SSatish Balay } else { 95*9371c9d4SSatish 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])); 961a989b97SToby Isaac } 971a989b97SToby Isaac } 981a989b97SToby Isaac } 991a989b97SToby Isaac } else { 1001a989b97SToby Isaac PetscInt Nk, Nf; 101fad4db65SToby Isaac PetscInt *subset, *perm; 1021a989b97SToby Isaac PetscInt i, j, l; 1031a989b97SToby Isaac PetscReal sum = 0.; 1041a989b97SToby Isaac 1059566063dSJacob Faibussowitsch PetscCall(PetscDTFactorialInt(k, &Nf)); 1069566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 1079566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &perm)); 1081a989b97SToby Isaac for (i = 0; i < Nk; i++) { 1091a989b97SToby Isaac PetscReal subsum = 0.; 1101a989b97SToby Isaac 1119566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 1121a989b97SToby Isaac for (j = 0; j < Nf; j++) { 1131a989b97SToby Isaac PetscBool permOdd; 1141a989b97SToby Isaac PetscReal prod; 1151a989b97SToby Isaac 1169566063dSJacob Faibussowitsch PetscCall(PetscDTEnumPerm(k, j, perm, &permOdd)); 1171a989b97SToby Isaac prod = permOdd ? -1. : 1.; 118*9371c9d4SSatish Balay for (l = 0; l < k; l++) { prod *= v[perm[l] * N + subset[l]]; } 1191a989b97SToby Isaac subsum += prod; 1201a989b97SToby Isaac } 1211a989b97SToby Isaac sum += w[i] * subsum; 1221a989b97SToby Isaac } 1239566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, perm)); 1241a989b97SToby Isaac *wv = sum; 1251a989b97SToby Isaac } 1261a989b97SToby Isaac PetscFunctionReturn(0); 1271a989b97SToby Isaac } 1281a989b97SToby Isaac 129fad4db65SToby Isaac /*@ 13028222859SToby Isaac PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form 131fad4db65SToby Isaac 1324165533cSJose E. Roman Input Parameters: 13328222859SToby Isaac + N - the dimension of the vector space, N >= 0 13428222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N 13528222859SToby Isaac . k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N 13628222859SToby Isaac . a - a j-form, size [N choose j] 13728222859SToby Isaac - b - a k-form, size [N choose k] 138fad4db65SToby Isaac 1394165533cSJose E. Roman Output Parameter: 14028222859SToby 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}}), 14128222859SToby 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}. 142fad4db65SToby Isaac 143fad4db65SToby Isaac Level: intermediate 144fad4db65SToby Isaac 145db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVWedgeMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 146fad4db65SToby Isaac @*/ 147*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) { 1481a989b97SToby Isaac PetscInt i; 1491a989b97SToby Isaac 1501a989b97SToby Isaac PetscFunctionBegin; 15108401ef6SPierre Jolivet PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 1521dca8a05SBarry Smith PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 1531dca8a05SBarry Smith PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 1541a989b97SToby Isaac if (N <= 3) { 1551a989b97SToby Isaac PetscInt Njk; 1561a989b97SToby Isaac 1579566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 158*9371c9d4SSatish Balay if (!j) { 159*9371c9d4SSatish Balay for (i = 0; i < Njk; i++) { awedgeb[i] = a[0] * b[i]; } 160*9371c9d4SSatish Balay } else if (!k) { 161*9371c9d4SSatish Balay for (i = 0; i < Njk; i++) { awedgeb[i] = a[i] * b[0]; } 162*9371c9d4SSatish Balay } else { 163*9371c9d4SSatish Balay if (N == 2) { 164*9371c9d4SSatish Balay awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 165*9371c9d4SSatish Balay } else { 1661a989b97SToby Isaac if (j + k == 2) { 1671a989b97SToby Isaac awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 1681a989b97SToby Isaac awedgeb[1] = a[0] * b[2] - a[2] * b[0]; 1691a989b97SToby Isaac awedgeb[2] = a[1] * b[2] - a[2] * b[1]; 1701a989b97SToby Isaac } else { 1711a989b97SToby Isaac awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0]; 1721a989b97SToby Isaac } 1731a989b97SToby Isaac } 1741a989b97SToby Isaac } 1751a989b97SToby Isaac } else { 1761a989b97SToby Isaac PetscInt Njk; 1771a989b97SToby Isaac PetscInt JKj; 1781a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 1791a989b97SToby Isaac PetscInt i; 1801a989b97SToby Isaac 1819566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 1829566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 1839566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk)); 1841a989b97SToby Isaac for (i = 0; i < Njk; i++) { 1851a989b97SToby Isaac PetscReal sum = 0.; 1861a989b97SToby Isaac PetscInt l; 1871a989b97SToby Isaac 1889566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, j + k, i, subset)); 1891a989b97SToby Isaac for (l = 0; l < JKj; l++) { 1901a989b97SToby Isaac PetscBool jkOdd; 1911a989b97SToby Isaac PetscInt m, jInd, kInd; 1921a989b97SToby Isaac 1939566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd)); 194*9371c9d4SSatish Balay for (m = 0; m < j; m++) { subsetj[m] = subset[subsetjk[m]]; } 195*9371c9d4SSatish Balay for (m = 0; m < k; m++) { subsetk[m] = subset[subsetjk[j + m]]; } 1969566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd)); 1979566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd)); 1981a989b97SToby Isaac sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]); 1991a989b97SToby Isaac } 2001a989b97SToby Isaac awedgeb[i] = sum; 2011a989b97SToby Isaac } 2029566063dSJacob Faibussowitsch PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk)); 2031a989b97SToby Isaac } 2041a989b97SToby Isaac PetscFunctionReturn(0); 2051a989b97SToby Isaac } 2061a989b97SToby Isaac 207fad4db65SToby Isaac /*@ 20828222859SToby Isaac PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms 209fad4db65SToby Isaac 2104165533cSJose E. Roman Input Parameters: 21128222859SToby Isaac + N - the dimension of the vector space, N >= 0 21228222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N 21328222859SToby Isaac . k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N 21428222859SToby Isaac - a - a j-form, size [N choose j] 215fad4db65SToby Isaac 2164165533cSJose E. Roman Output Parameter: 21728222859SToby 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 218fad4db65SToby Isaac 219fad4db65SToby Isaac Level: intermediate 220fad4db65SToby Isaac 221db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 222fad4db65SToby Isaac @*/ 223*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat) { 2241a989b97SToby Isaac PetscInt i; 2251a989b97SToby Isaac 2261a989b97SToby Isaac PetscFunctionBegin; 22708401ef6SPierre Jolivet PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 2281dca8a05SBarry Smith PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 2291dca8a05SBarry Smith PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 2301a989b97SToby Isaac if (N <= 3) { 2311a989b97SToby Isaac PetscInt Njk; 2321a989b97SToby Isaac 2339566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 2341a989b97SToby Isaac if (!j) { 2351a989b97SToby Isaac for (i = 0; i < Njk * Njk; i++) { awedgeMat[i] = 0.; } 2361a989b97SToby Isaac for (i = 0; i < Njk; i++) { awedgeMat[i * (Njk + 1)] = a[0]; } 2371a989b97SToby Isaac } else if (!k) { 2381a989b97SToby Isaac for (i = 0; i < Njk; i++) { awedgeMat[i] = a[i]; } 2391a989b97SToby Isaac } else { 2401a989b97SToby Isaac if (N == 2) { 241*9371c9d4SSatish Balay awedgeMat[0] = -a[1]; 242*9371c9d4SSatish Balay awedgeMat[1] = a[0]; 2431a989b97SToby Isaac } else { 2441a989b97SToby Isaac if (j + k == 2) { 245*9371c9d4SSatish Balay awedgeMat[0] = -a[1]; 246*9371c9d4SSatish Balay awedgeMat[1] = a[0]; 247*9371c9d4SSatish Balay awedgeMat[2] = 0.; 248*9371c9d4SSatish Balay awedgeMat[3] = -a[2]; 249*9371c9d4SSatish Balay awedgeMat[4] = 0.; 250*9371c9d4SSatish Balay awedgeMat[5] = a[0]; 251*9371c9d4SSatish Balay awedgeMat[6] = 0.; 252*9371c9d4SSatish Balay awedgeMat[7] = -a[2]; 253*9371c9d4SSatish Balay awedgeMat[8] = a[1]; 2541a989b97SToby Isaac } else { 255*9371c9d4SSatish Balay awedgeMat[0] = a[2]; 256*9371c9d4SSatish Balay awedgeMat[1] = -a[1]; 257*9371c9d4SSatish Balay awedgeMat[2] = a[0]; 2581a989b97SToby Isaac } 2591a989b97SToby Isaac } 2601a989b97SToby Isaac } 2611a989b97SToby Isaac } else { 2621a989b97SToby Isaac PetscInt Njk; 2631a989b97SToby Isaac PetscInt Nk; 2641a989b97SToby Isaac PetscInt JKj, i; 2651a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 2661a989b97SToby Isaac 2679566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 2689566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 2699566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 2709566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk)); 2711a989b97SToby Isaac for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.; 2721a989b97SToby Isaac for (i = 0; i < Njk; i++) { 2731a989b97SToby Isaac PetscInt l; 2741a989b97SToby Isaac 2759566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, j + k, i, subset)); 2761a989b97SToby Isaac for (l = 0; l < JKj; l++) { 2771a989b97SToby Isaac PetscBool jkOdd; 2781a989b97SToby Isaac PetscInt m, jInd, kInd; 2791a989b97SToby Isaac 2809566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd)); 281*9371c9d4SSatish Balay for (m = 0; m < j; m++) { subsetj[m] = subset[subsetjk[m]]; } 282*9371c9d4SSatish Balay for (m = 0; m < k; m++) { subsetk[m] = subset[subsetjk[j + m]]; } 2839566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd)); 2851a989b97SToby Isaac awedgeMat[i * Nk + kInd] += jkOdd ? -a[jInd] : a[jInd]; 2861a989b97SToby Isaac } 2871a989b97SToby Isaac } 2889566063dSJacob Faibussowitsch PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk)); 2891a989b97SToby Isaac } 2901a989b97SToby Isaac PetscFunctionReturn(0); 2911a989b97SToby Isaac } 2921a989b97SToby Isaac 293fad4db65SToby Isaac /*@ 29428222859SToby Isaac PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space 295fad4db65SToby Isaac 2964165533cSJose E. Roman Input Parameters: 29728222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, M >= 0 29828222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, N >= 0 29928222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format 30028222859SToby Isaac . k - the *signed* degree k of the |k|-form w, -(min(M,N)) <= k <= min(M,N). A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note). 30128222859SToby Isaac - w - a |k|-form in the image space, size [M choose |k|] 302fad4db65SToby Isaac 3034165533cSJose E. Roman Output Parameter: 30428222859SToby Isaac . Lstarw - the pullback of w to a |k|-form in the origin space, size [N choose |k|]: (Lstarw)(v_1,...v_k) = w(L*v_1,...,L*v_k). 305fad4db65SToby Isaac 306fad4db65SToby Isaac Level: intermediate 307fad4db65SToby Isaac 308a5b23f4aSJose E. Roman Note: negative form degrees accommodate, e.g., H-div conforming vector fields. An H-div conforming vector field stores its degrees of freedom as (dx, dy, dz), like a 1-form, 30928222859SToby Isaac but its normal trace is integrated on faces, like a 2-form. The correct pullback then is to apply the Hodge star transformation from (M-2)-form to 2-form, pullback as a 2-form, 310fad4db65SToby Isaac then the inverse Hodge star transformation. 311fad4db65SToby Isaac 312db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullbackMatrix()`, `PetscDTAltVStar()` 313fad4db65SToby Isaac @*/ 314*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw) { 3151a989b97SToby Isaac PetscInt i, j, Nk, Mk; 3161a989b97SToby Isaac 3171a989b97SToby Isaac PetscFunctionBegin; 3181dca8a05SBarry Smith PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 3191dca8a05SBarry Smith PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 3201a989b97SToby Isaac if (N <= 3 && M <= 3) { 3219566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 3229566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 3231a989b97SToby Isaac if (!k) { 3241a989b97SToby Isaac Lstarw[0] = w[0]; 3251a989b97SToby Isaac } else if (k == 1) { 3261a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3271a989b97SToby Isaac PetscReal sum = 0.; 3281a989b97SToby Isaac 3291a989b97SToby Isaac for (j = 0; j < Mk; j++) { sum += L[j * Nk + i] * w[j]; } 3301a989b97SToby Isaac Lstarw[i] = sum; 3311a989b97SToby Isaac } 3321a989b97SToby Isaac } else if (k == -1) { 3331a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 3341a989b97SToby Isaac 3351a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3361a989b97SToby Isaac PetscReal sum = 0.; 3371a989b97SToby Isaac 338*9371c9d4SSatish Balay for (j = 0; j < Mk; j++) { sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j]; } 3391a989b97SToby Isaac Lstarw[i] = mult[i] * sum; 3401a989b97SToby Isaac } 3411a989b97SToby Isaac } else if (k == 2) { 342*9371c9d4SSatish Balay PetscInt pairs[3][2] = { 343*9371c9d4SSatish Balay {0, 1}, 344*9371c9d4SSatish Balay {0, 2}, 345*9371c9d4SSatish Balay {1, 2} 346*9371c9d4SSatish Balay }; 3471a989b97SToby Isaac 3481a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3491a989b97SToby Isaac PetscReal sum = 0.; 350*9371c9d4SSatish Balay 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]; } 3511a989b97SToby Isaac Lstarw[i] = sum; 3521a989b97SToby Isaac } 3531a989b97SToby Isaac } else if (k == -2) { 354*9371c9d4SSatish Balay PetscInt pairs[3][2] = { 355*9371c9d4SSatish Balay {1, 2}, 356*9371c9d4SSatish Balay {2, 0}, 357*9371c9d4SSatish Balay {0, 1} 358*9371c9d4SSatish Balay }; 3591a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 3601a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 3611a989b97SToby Isaac 3621a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3631a989b97SToby Isaac PetscReal sum = 0.; 3641a989b97SToby Isaac 365*9371c9d4SSatish Balay 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]; } 3661a989b97SToby Isaac Lstarw[i] = sum; 3671a989b97SToby Isaac } 3681a989b97SToby Isaac } else { 369*9371c9d4SSatish 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]); 3701a989b97SToby Isaac 3711a989b97SToby Isaac for (i = 0; i < Nk; i++) { Lstarw[i] = detL * w[i]; } 3721a989b97SToby Isaac } 3731a989b97SToby Isaac } else { 3741a989b97SToby Isaac PetscInt Nf, l, p; 3751a989b97SToby Isaac PetscReal *Lw, *Lwv; 3761a989b97SToby Isaac PetscInt *subsetw, *subsetv; 377fad4db65SToby Isaac PetscInt *perm; 3781a989b97SToby Isaac PetscReal *walloc = NULL; 3791a989b97SToby Isaac const PetscReal *ww = NULL; 3801a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 3811a989b97SToby Isaac 3829566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 3839566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 3849566063dSJacob Faibussowitsch PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf)); 3851a989b97SToby Isaac if (k < 0) { 3861a989b97SToby Isaac negative = PETSC_TRUE; 3871a989b97SToby Isaac k = -k; 3889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Mk, &walloc)); 3899566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc)); 3901a989b97SToby Isaac ww = walloc; 3911a989b97SToby Isaac } else { 3921a989b97SToby Isaac ww = w; 3931a989b97SToby Isaac } 3949566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv)); 3951a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = 0.; 3961a989b97SToby Isaac for (i = 0; i < Mk; i++) { 3979566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(M, k, i, subsetw)); 3981a989b97SToby Isaac for (j = 0; j < Nk; j++) { 3999566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, j, subsetv)); 4001a989b97SToby Isaac for (p = 0; p < Nf; p++) { 4011a989b97SToby Isaac PetscReal prod; 4021a989b97SToby Isaac PetscBool isOdd; 4031a989b97SToby Isaac 4049566063dSJacob Faibussowitsch PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd)); 4051a989b97SToby Isaac prod = isOdd ? -ww[i] : ww[i]; 406*9371c9d4SSatish Balay for (l = 0; l < k; l++) { prod *= L[subsetw[perm[l]] * N + subsetv[l]]; } 4071a989b97SToby Isaac Lstarw[j] += prod; 4081a989b97SToby Isaac } 4091a989b97SToby Isaac } 4101a989b97SToby Isaac } 4111a989b97SToby Isaac if (negative) { 4121a989b97SToby Isaac PetscReal *sLsw; 4131a989b97SToby Isaac 4149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk, &sLsw)); 4159566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw)); 4161a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i]; 4179566063dSJacob Faibussowitsch PetscCall(PetscFree(sLsw)); 4181a989b97SToby Isaac } 4199566063dSJacob Faibussowitsch PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv)); 4209566063dSJacob Faibussowitsch PetscCall(PetscFree(walloc)); 4211a989b97SToby Isaac } 4221a989b97SToby Isaac PetscFunctionReturn(0); 4231a989b97SToby Isaac } 4241a989b97SToby Isaac 425fad4db65SToby Isaac /*@ 426fad4db65SToby Isaac PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation 427fad4db65SToby Isaac 4284165533cSJose E. Roman Input Parameters: 42928222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, N >= 0 43028222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, M >= 0 43128222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format 43228222859SToby 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()) 433fad4db65SToby Isaac 4344165533cSJose E. Roman Output Parameter: 43528222859SToby Isaac . Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w 436fad4db65SToby Isaac 437fad4db65SToby Isaac Level: intermediate 438fad4db65SToby Isaac 439db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVStar()` 440fad4db65SToby Isaac @*/ 441*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) { 4421a989b97SToby Isaac PetscInt Nk, Mk, Nf, i, j, l, p; 4431a989b97SToby Isaac PetscReal *Lw, *Lwv; 4441a989b97SToby Isaac PetscInt *subsetw, *subsetv; 445fad4db65SToby Isaac PetscInt *perm; 4461a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 4471a989b97SToby Isaac 4481a989b97SToby Isaac PetscFunctionBegin; 4491dca8a05SBarry Smith PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 4501dca8a05SBarry Smith PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 4511a989b97SToby Isaac if (N <= 3 && M <= 3) { 4521a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 4531a989b97SToby Isaac 4549566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 4559566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 4561a989b97SToby Isaac if (!k) { 4571a989b97SToby Isaac Lstar[0] = 1.; 4581a989b97SToby Isaac } else if (k == 1) { 459*9371c9d4SSatish Balay for (i = 0; i < Nk; i++) { 460*9371c9d4SSatish Balay for (j = 0; j < Mk; j++) { Lstar[i * Mk + j] = L[j * Nk + i]; } 461*9371c9d4SSatish Balay } 4621a989b97SToby Isaac } else if (k == -1) { 4631a989b97SToby Isaac for (i = 0; i < Nk; i++) { 464*9371c9d4SSatish Balay for (j = 0; j < Mk; j++) { Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j]; } 4651a989b97SToby Isaac } 4661a989b97SToby Isaac } else if (k == 2) { 467*9371c9d4SSatish Balay PetscInt pairs[3][2] = { 468*9371c9d4SSatish Balay {0, 1}, 469*9371c9d4SSatish Balay {0, 2}, 470*9371c9d4SSatish Balay {1, 2} 471*9371c9d4SSatish Balay }; 4721a989b97SToby Isaac 4731a989b97SToby Isaac for (i = 0; i < Nk; i++) { 474*9371c9d4SSatish Balay 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]]; } 4751a989b97SToby Isaac } 4761a989b97SToby Isaac } else if (k == -2) { 477*9371c9d4SSatish Balay PetscInt pairs[3][2] = { 478*9371c9d4SSatish Balay {1, 2}, 479*9371c9d4SSatish Balay {2, 0}, 480*9371c9d4SSatish Balay {0, 1} 481*9371c9d4SSatish Balay }; 4821a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 4831a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 4841a989b97SToby Isaac 4851a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4861a989b97SToby Isaac for (j = 0; j < Mk; j++) { 487*9371c9d4SSatish 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]]; 4881a989b97SToby Isaac } 4891a989b97SToby Isaac } 4901a989b97SToby Isaac } else { 491*9371c9d4SSatish 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]); 4921a989b97SToby Isaac 4931a989b97SToby Isaac for (i = 0; i < Nk; i++) { Lstar[i] = detL; } 4941a989b97SToby Isaac } 4951a989b97SToby Isaac } else { 4961a989b97SToby Isaac if (k < 0) { 4971a989b97SToby Isaac negative = PETSC_TRUE; 4981a989b97SToby Isaac k = -k; 4991a989b97SToby Isaac } 5009566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 5019566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 5029566063dSJacob Faibussowitsch PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf)); 5039566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv)); 5041a989b97SToby Isaac for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.; 5051a989b97SToby Isaac for (i = 0; i < Mk; i++) { 5061a989b97SToby Isaac PetscBool iOdd; 5071a989b97SToby Isaac PetscInt iidx, jidx; 5081a989b97SToby Isaac 5099566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(M, k, i, subsetw, &iOdd)); 5101a989b97SToby Isaac iidx = negative ? Mk - 1 - i : i; 51128222859SToby Isaac iOdd = negative ? (PetscBool)(iOdd ^ ((k * (M - k)) & 1)) : PETSC_FALSE; 5121a989b97SToby Isaac for (j = 0; j < Nk; j++) { 5131a989b97SToby Isaac PetscBool jOdd; 5141a989b97SToby Isaac 5159566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(N, k, j, subsetv, &jOdd)); 5161a989b97SToby Isaac jidx = negative ? Nk - 1 - j : j; 51728222859SToby Isaac jOdd = negative ? (PetscBool)(iOdd ^ jOdd ^ ((k * (N - k)) & 1)) : PETSC_FALSE; 5181a989b97SToby Isaac for (p = 0; p < Nf; p++) { 5191a989b97SToby Isaac PetscReal prod; 5201a989b97SToby Isaac PetscBool isOdd; 5211a989b97SToby Isaac 5229566063dSJacob Faibussowitsch PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd)); 52328222859SToby Isaac isOdd = (PetscBool)(isOdd ^ jOdd); 5241a989b97SToby Isaac prod = isOdd ? -1. : 1.; 525*9371c9d4SSatish Balay for (l = 0; l < k; l++) { prod *= L[subsetw[perm[l]] * N + subsetv[l]]; } 5261a989b97SToby Isaac Lstar[jidx * Mk + iidx] += prod; 5271a989b97SToby Isaac } 5281a989b97SToby Isaac } 5291a989b97SToby Isaac } 5309566063dSJacob Faibussowitsch PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv)); 5311a989b97SToby Isaac } 5321a989b97SToby Isaac PetscFunctionReturn(0); 5331a989b97SToby Isaac } 5341a989b97SToby Isaac 535fad4db65SToby Isaac /*@ 53628222859SToby Isaac PetscDTAltVInterior - Compute the interior product of a k-form with a vector 537fad4db65SToby Isaac 5384165533cSJose E. Roman Input Parameters: 53928222859SToby Isaac + N - the dimension of the vector space, N >= 0 54028222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 54128222859SToby Isaac . w - a k-form, size [N choose k] 54228222859SToby Isaac - v - an N dimensional vector 543fad4db65SToby Isaac 5444165533cSJose E. Roman Output Parameter: 54528222859SToby 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}). 546fad4db65SToby Isaac 547fad4db65SToby Isaac Level: intermediate 548fad4db65SToby Isaac 549db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 550fad4db65SToby Isaac @*/ 551*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv) { 5521a989b97SToby Isaac PetscInt i, Nk, Nkm; 5531a989b97SToby Isaac 5541a989b97SToby Isaac PetscFunctionBegin; 5551dca8a05SBarry Smith PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 5569566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 5579566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 5581a989b97SToby Isaac if (N <= 3) { 5591a989b97SToby Isaac if (k == 1) { 5601a989b97SToby Isaac PetscReal sum = 0.; 5611a989b97SToby Isaac 562*9371c9d4SSatish Balay for (i = 0; i < N; i++) { sum += w[i] * v[i]; } 5631a989b97SToby Isaac wIntv[0] = sum; 5641a989b97SToby Isaac } else if (k == N) { 5651a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 5661a989b97SToby Isaac 567*9371c9d4SSatish Balay for (i = 0; i < N; i++) { wIntv[N - 1 - i] = w[0] * v[i] * mult[i]; } 5681a989b97SToby Isaac } else { 5691a989b97SToby Isaac wIntv[0] = -w[0] * v[1] - w[1] * v[2]; 5701a989b97SToby Isaac wIntv[1] = w[0] * v[0] - w[2] * v[2]; 5711a989b97SToby Isaac wIntv[2] = w[1] * v[0] + w[2] * v[1]; 5721a989b97SToby Isaac } 5731a989b97SToby Isaac } else { 5741a989b97SToby Isaac PetscInt *subset, *work; 5751a989b97SToby Isaac 5769566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &work)); 5771a989b97SToby Isaac for (i = 0; i < Nkm; i++) wIntv[i] = 0.; 5781a989b97SToby Isaac for (i = 0; i < Nk; i++) { 5791a989b97SToby Isaac PetscInt j, l, m; 5801a989b97SToby Isaac 5819566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 5821a989b97SToby Isaac for (j = 0; j < k; j++) { 5831a989b97SToby Isaac PetscInt idx; 58428222859SToby Isaac PetscBool flip = (PetscBool)(j & 1); 5851a989b97SToby Isaac 5861a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 5871a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 5881a989b97SToby Isaac } 5899566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 5901a989b97SToby Isaac wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]); 5911a989b97SToby Isaac } 5921a989b97SToby Isaac } 5939566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, work)); 5941a989b97SToby Isaac } 5951a989b97SToby Isaac PetscFunctionReturn(0); 5961a989b97SToby Isaac } 5971a989b97SToby Isaac 598fad4db65SToby Isaac /*@ 59928222859SToby Isaac PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector 600fad4db65SToby Isaac 6014165533cSJose E. Roman Input Parameters: 60228222859SToby Isaac + N - the dimension of the vector space, N >= 0 60328222859SToby Isaac . k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N 60428222859SToby Isaac - v - an N dimensional vector 605fad4db65SToby Isaac 6064165533cSJose E. Roman Output Parameter: 607fad4db65SToby Isaac . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v) 608fad4db65SToby Isaac 609fad4db65SToby Isaac Level: intermediate 610fad4db65SToby Isaac 611db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 612fad4db65SToby Isaac @*/ 613*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat) { 6141a989b97SToby Isaac PetscInt i, Nk, Nkm; 6151a989b97SToby Isaac 6161a989b97SToby Isaac PetscFunctionBegin; 6171dca8a05SBarry Smith PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 6189566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 6199566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 6201a989b97SToby Isaac if (N <= 3) { 6211a989b97SToby Isaac if (k == 1) { 6221a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[i] = v[i]; 6231a989b97SToby Isaac } else if (k == N) { 6241a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 6251a989b97SToby Isaac 6261a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i]; 6271a989b97SToby Isaac } else { 628*9371c9d4SSatish Balay intvMat[0] = -v[1]; 629*9371c9d4SSatish Balay intvMat[1] = -v[2]; 630*9371c9d4SSatish Balay intvMat[2] = 0.; 631*9371c9d4SSatish Balay intvMat[3] = v[0]; 632*9371c9d4SSatish Balay intvMat[4] = 0.; 633*9371c9d4SSatish Balay intvMat[5] = -v[2]; 634*9371c9d4SSatish Balay intvMat[6] = 0.; 635*9371c9d4SSatish Balay intvMat[7] = v[0]; 636*9371c9d4SSatish Balay intvMat[8] = v[1]; 6371a989b97SToby Isaac } 6381a989b97SToby Isaac } else { 6391a989b97SToby Isaac PetscInt *subset, *work; 6401a989b97SToby Isaac 6419566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &work)); 6421a989b97SToby Isaac for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.; 6431a989b97SToby Isaac for (i = 0; i < Nk; i++) { 6441a989b97SToby Isaac PetscInt j, l, m; 6451a989b97SToby Isaac 6469566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 6471a989b97SToby Isaac for (j = 0; j < k; j++) { 6481a989b97SToby Isaac PetscInt idx; 64928222859SToby Isaac PetscBool flip = (PetscBool)(j & 1); 6501a989b97SToby Isaac 6511a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 6521a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 6531a989b97SToby Isaac } 6549566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 6551a989b97SToby Isaac intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]]; 6561a989b97SToby Isaac } 6571a989b97SToby Isaac } 6589566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, work)); 6591a989b97SToby Isaac } 6601a989b97SToby Isaac PetscFunctionReturn(0); 6611a989b97SToby Isaac } 6621a989b97SToby Isaac 663fad4db65SToby Isaac /*@ 664fad4db65SToby Isaac PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix() 665fad4db65SToby Isaac 6664165533cSJose E. Roman Input Parameters: 66728222859SToby Isaac + N - the dimension of the vector space, N >= 0 66828222859SToby Isaac - k - the degree of the k-forms on which intvMat from PetscDTAltVInteriorMatrix() acts, 0 <= k <= N. 669fad4db65SToby Isaac 6704165533cSJose E. Roman Output Parameter: 67128222859SToby Isaac . indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k 67228222859SToby 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 67328222859SToby Isaac coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1) 674fad4db65SToby Isaac 675fad4db65SToby Isaac Level: intermediate 676fad4db65SToby Isaac 67728222859SToby Isaac Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential 678fad4db65SToby Isaac 679db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 680fad4db65SToby Isaac @*/ 681*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3]) { 682dda711d0SToby Isaac PetscInt i, Nk, Nkm; 683dda711d0SToby Isaac 684dda711d0SToby Isaac PetscFunctionBegin; 6851dca8a05SBarry Smith PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 6869566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 6879566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 688dda711d0SToby Isaac if (N <= 3) { 689dda711d0SToby Isaac if (k == 1) { 690dda711d0SToby Isaac for (i = 0; i < N; i++) { 691dda711d0SToby Isaac indices[i][0] = 0; 692dda711d0SToby Isaac indices[i][1] = i; 693dda711d0SToby Isaac indices[i][2] = i; 694dda711d0SToby Isaac } 695dda711d0SToby Isaac } else if (k == N) { 696dda711d0SToby Isaac PetscInt val[3] = {0, -2, 2}; 697dda711d0SToby Isaac 698dda711d0SToby Isaac for (i = 0; i < N; i++) { 699dda711d0SToby Isaac indices[i][0] = N - 1 - i; 700dda711d0SToby Isaac indices[i][1] = 0; 701dda711d0SToby Isaac indices[i][2] = val[i]; 702dda711d0SToby Isaac } 703dda711d0SToby Isaac } else { 704*9371c9d4SSatish Balay indices[0][0] = 0; 705*9371c9d4SSatish Balay indices[0][1] = 0; 706*9371c9d4SSatish Balay indices[0][2] = -(1 + 1); 707*9371c9d4SSatish Balay indices[1][0] = 0; 708*9371c9d4SSatish Balay indices[1][1] = 1; 709*9371c9d4SSatish Balay indices[1][2] = -(2 + 1); 710*9371c9d4SSatish Balay indices[2][0] = 1; 711*9371c9d4SSatish Balay indices[2][1] = 0; 712*9371c9d4SSatish Balay indices[2][2] = 0; 713*9371c9d4SSatish Balay indices[3][0] = 1; 714*9371c9d4SSatish Balay indices[3][1] = 2; 715*9371c9d4SSatish Balay indices[3][2] = -(2 + 1); 716*9371c9d4SSatish Balay indices[4][0] = 2; 717*9371c9d4SSatish Balay indices[4][1] = 1; 718*9371c9d4SSatish Balay indices[4][2] = 0; 719*9371c9d4SSatish Balay indices[5][0] = 2; 720*9371c9d4SSatish Balay indices[5][1] = 2; 721*9371c9d4SSatish Balay indices[5][2] = 1; 722dda711d0SToby Isaac } 723dda711d0SToby Isaac } else { 724dda711d0SToby Isaac PetscInt *subset, *work; 725dda711d0SToby Isaac 7269566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(k, &subset, k, &work)); 727dda711d0SToby Isaac for (i = 0; i < Nk; i++) { 728dda711d0SToby Isaac PetscInt j, l, m; 729dda711d0SToby Isaac 7309566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSubset(N, k, i, subset)); 731dda711d0SToby Isaac for (j = 0; j < k; j++) { 732dda711d0SToby Isaac PetscInt idx; 73328222859SToby Isaac PetscBool flip = (PetscBool)(j & 1); 734dda711d0SToby Isaac 735dda711d0SToby Isaac for (l = 0, m = 0; l < k; l++) { 736dda711d0SToby Isaac if (l != j) work[m++] = subset[l]; 737dda711d0SToby Isaac } 7389566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 739dda711d0SToby Isaac indices[i * k + j][0] = idx; 740dda711d0SToby Isaac indices[i * k + j][1] = i; 741dda711d0SToby Isaac indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j]; 742dda711d0SToby Isaac } 743dda711d0SToby Isaac } 7449566063dSJacob Faibussowitsch PetscCall(PetscFree2(subset, work)); 745dda711d0SToby Isaac } 746dda711d0SToby Isaac PetscFunctionReturn(0); 747dda711d0SToby Isaac } 748dda711d0SToby Isaac 749fad4db65SToby Isaac /*@ 75028222859SToby Isaac PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form 751fad4db65SToby Isaac 7524165533cSJose E. Roman Input Parameters: 75328222859SToby Isaac + N - the dimension of the vector space, N >= 0 75428222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 75528222859SToby 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. 75628222859SToby Isaac - w - a k-form, size [N choose k] 757fad4db65SToby Isaac 7584165533cSJose E. Roman Output Parameter: 75928222859SToby 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. 760fad4db65SToby Isaac 761fad4db65SToby Isaac Level: intermediate 762fad4db65SToby Isaac 763db781477SPatrick Sanan .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 764fad4db65SToby Isaac @*/ 765*9371c9d4SSatish Balay PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw) { 7661a989b97SToby Isaac PetscInt Nk, i; 7671a989b97SToby Isaac 7681a989b97SToby Isaac PetscFunctionBegin; 7691dca8a05SBarry Smith PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 7709566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 7711a989b97SToby Isaac pow = pow % 4; 7721a989b97SToby Isaac pow = (pow + 4) % 4; /* make non-negative */ 7731a989b97SToby Isaac /* pow is now 0, 1, 2, 3 */ 7741a989b97SToby Isaac if (N <= 3) { 7751a989b97SToby Isaac if (pow & 1) { 7761a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 7771a989b97SToby Isaac 7781a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i]; 7791a989b97SToby Isaac } else { 7801a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 7811a989b97SToby Isaac } 7821a989b97SToby Isaac if (pow > 1 && ((k * (N - k)) & 1)) { 7831a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 7841a989b97SToby Isaac } 7851a989b97SToby Isaac } else { 7861a989b97SToby Isaac PetscInt *subset; 7871a989b97SToby Isaac 7889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &subset)); 7891a989b97SToby Isaac if (pow % 2) { 7901a989b97SToby Isaac PetscInt l = (pow == 1) ? k : N - k; 7911a989b97SToby Isaac for (i = 0; i < Nk; i++) { 7921a989b97SToby Isaac PetscBool sOdd; 7931a989b97SToby Isaac PetscInt j, idx; 7941a989b97SToby Isaac 7959566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(N, l, i, subset, &sOdd)); 7969566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, l, subset, &idx)); 7979566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, N - l, &subset[l], &j)); 7981a989b97SToby Isaac starw[j] = sOdd ? -w[idx] : w[idx]; 7991a989b97SToby Isaac } 8001a989b97SToby Isaac } else { 8011a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 8021a989b97SToby Isaac } 8031a989b97SToby Isaac /* star^2 = -1^(k * (N - k)) */ 8041a989b97SToby Isaac if (pow > 1 && (k * (N - k)) % 2) { 8051a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 8061a989b97SToby Isaac } 8079566063dSJacob Faibussowitsch PetscCall(PetscFree(subset)); 8081a989b97SToby Isaac } 8091a989b97SToby Isaac PetscFunctionReturn(0); 8101a989b97SToby Isaac } 811