11a989b97SToby Isaac #include <petsc/private/petscimpl.h> 228222859SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/ 31a989b97SToby Isaac 429a920c6SToby Isaac /*MC 529a920c6SToby Isaac PetscDTAltV - An interface for common operations on k-forms, also known as alternating algebraic forms or alternating k-linear maps. 629a920c6SToby Isaac The name of the interface comes from the notation "Alt V" for the algebra of all k-forms acting vectors in the space V, also known as the exterior algebra of V*. 729a920c6SToby Isaac 829a920c6SToby Isaac A recommended reference for this material is Section 2 "Exterior algebra and exterior calculus" in "Finite element 929a920c6SToby Isaac exterior calculus, homological techniques, and applications", by Arnold, Falk, & Winther (2006, doi:10.1017/S0962492906210018). 1029a920c6SToby Isaac 1129a920c6SToby Isaac A k-form w (k is called the "form degree" of w) is an alternating k-linear map acting on tuples (v_1, ..., v_k) of 1229a920c6SToby Isaac vectors from a vector space V and producing a real number: 1329a920c6SToby Isaac - alternating: swapping any two vectors in a tuple reverses the sign of the result, e.g. w(v_1, v_2, ..., v_k) = -w(v_2, v_1, ..., v_k) 1429a920c6SToby Isaac - k-linear: w acts linear in each vector separately, e.g. w(a*v + b*y, v_2, ..., v_k) = a*w(v,v_2,...,v_k) + b*w(y,v_2,...,v_k) 1529a920c6SToby Isaac This action is implemented as PetscDTAltVApply. 1629a920c6SToby Isaac 1729a920c6SToby Isaac The k-forms on a vector space form a vector space themselves, Alt^k V. The dimension of Alt^k V, if V is N dimensional, is N choose k. (This 1829a920c6SToby Isaac shows that for an N dimensional space, only 0 <= k <= N are valid form degrees.) 1929a920c6SToby Isaac The standard basis for Alt^k V, used in PetscDTAltV, has one basis k-form for each ordered subset of k coordinates of the N dimensional space: 2029a920c6SToby Isaac For example, if the coordinate directions of a four dimensional space are (t, x, y, z), then there are 4 choose 2 = 6 ordered subsets of two coordinates. 2129a920c6SToby Isaac They are, in lexicographic order, (t, x), (t, y), (t, z), (x, y), (x, z) and (y, z). PetscDTAltV also orders the basis of Alt^k V lexicographically 2229a920c6SToby Isaac by the associated subsets. 2329a920c6SToby Isaac 2429a920c6SToby Isaac The unit basis k-form associated with coordinates (c_1, ..., c_k) acts on a set of k vectors (v_1, ..., v_k) by creating a square matrix V where 2529a920c6SToby Isaac V[i,j] = v_i[c_j] and taking the determinant of V. 2629a920c6SToby Isaac 2729a920c6SToby Isaac If j + k <= N, then a j-form f and a k-form g can be multiplied to create a (j+k)-form using the wedge or exterior product, (f wedge g). 2829a920c6SToby Isaac This is an anticommutative product, (f wedge g) = -(g wedge f). It is sufficient to describe the wedge product of two basis forms. 2929a920c6SToby Isaac Let f be the basis j-form associated with coordinates (f_1,...,f_j) and g be the basis k-form associated with coordinates (g_1,...,g_k): 3029a920c6SToby Isaac - If there is any coordinate in both sets, then (f wedge g) = 0. 3129a920c6SToby Isaac - Otherwise, (f wedge g) is a multiple of the basis (j+k)-form h associated with (f_1,...,f_j,g_1,...,g_k). 3229a920c6SToby Isaac - In fact it is equal to either h or -h depending on how (f_1,...,f_j,g_1,...,g_k) compares to the same list of coordinates given in ascending order: if it is an even permutation of that list, then (f wedge g) = h, otherwise (f wedge g) = -h. 3329a920c6SToby Isaac The wedge product is implemented for either two inputs (f and g) in PetscDTAltVWedge, or for one (just f, giving a 3429a920c6SToby Isaac matrix to multiply against multiple choices of g) in PetscDTAltVWedgeMatrix. 3529a920c6SToby Isaac 3629a920c6SToby Isaac If k > 0, a k-form w and a vector v can combine to make a (k-1)-formm through the interior product, (w int v), 3729a920c6SToby Isaac defined by (w int v)(v_1,...,v_{k-1}) = w(v,v_1,...,v_{k-1}). 3829a920c6SToby Isaac 3929a920c6SToby Isaac The interior product is implemented for either two inputs (w and v) in PetscDTAltVInterior, for one (just v, giving a 4029a920c6SToby Isaac matrix to multiply against multiple choices of w) in PetscDTAltVInteriorMatrix, 4129a920c6SToby Isaac or for no inputs (giving the sparsity pattern of PetscDTAltVInteriorMatrix) in PetscDTAltVInteriorPattern. 4229a920c6SToby Isaac 4329a920c6SToby Isaac When there is a linear map L: V -> W from an N dimensional vector space to an M dimensional vector space, 4429a920c6SToby Isaac it induces the linear pullback map L^* : Alt^k W -> Alt^k V, defined by L^* w(v_1,...,v_k) = w(L v_1, ..., L v_k). 4529a920c6SToby Isaac The pullback is implemented as PetscDTAltVPullback (acting on a known w) and PetscDTAltVPullbackMatrix (creating a matrix that computes the actin of L^*). 4629a920c6SToby Isaac 4729a920c6SToby Isaac Alt^k V and Alt^(N-k) V have the same dimension, and the Hodge star operator maps between them. We note that Alt^N V is a one dimensional space, and its 4829a920c6SToby Isaac basis vector is sometime called vol. The Hodge star operator has the property that (f wedge (star g)) = (f,g) vol, where (f,g) is the simple inner product 4929a920c6SToby Isaac of the basis coefficients of f and g. 5029a920c6SToby Isaac Powers of the Hodge star operator can be applied with PetscDTAltVStar 5129a920c6SToby Isaac 526c877ef6SSatish Balay Level: intermediate 5329a920c6SToby Isaac 5429a920c6SToby Isaac .seealso: PetscDTAltVApply(), PetscDTAltVWedge(), PetscDTAltVInterior(), PetscDTAltVPullback(), PetscDTAltVStar() 5529a920c6SToby Isaac M*/ 5629a920c6SToby Isaac 57fad4db65SToby Isaac /*@ 5828222859SToby Isaac PetscDTAltVApply - Apply an a k-form (an alternating k-linear map) to a set of k N-dimensional vectors 59fad4db65SToby Isaac 604165533cSJose E. Roman Input Parameters: 6128222859SToby Isaac + N - the dimension of the vector space, N >= 0 6228222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 6328222859SToby Isaac . w - a k-form, size [N choose k] (each degree of freedom of a k-form is associated with a subset of k coordinates of the N-dimensional vectors: the degrees of freedom are ordered lexicographically by their associated subsets) 6428222859SToby Isaac - v - a set of k vectors of size N, size [k x N], each vector stored contiguously 65fad4db65SToby Isaac 664165533cSJose E. Roman Output Parameter: 6728222859SToby Isaac . wv - w(v_1,...,v_k) = \sum_i w_i * det(V_i): the degree of freedom w_i is associated with coordinates [s_{i,1},...,s_{i,k}], and the square matrix V_i has entry (j,k) given by the s_{i,k}'th coordinate of v_j 68fad4db65SToby Isaac 69fad4db65SToby Isaac Level: intermediate 70fad4db65SToby Isaac 7129a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 72fad4db65SToby Isaac @*/ 731a989b97SToby Isaac PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv) 741a989b97SToby Isaac { 751a989b97SToby Isaac PetscErrorCode ierr; 761a989b97SToby Isaac 771a989b97SToby Isaac PetscFunctionBegin; 78*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 79*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(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 { 841a989b97SToby Isaac if (N == 1) {*wv = w[0] * v[0];} 851a989b97SToby Isaac else if (N == 2) { 861a989b97SToby Isaac if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1];} 871a989b97SToby Isaac else {*wv = w[0] * (v[0] * v[3] - v[1] * v[2]);} 881a989b97SToby Isaac } else { 891a989b97SToby Isaac if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];} 901a989b97SToby Isaac else if (k == 2) { 911a989b97SToby Isaac *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) + 921a989b97SToby Isaac w[1] * (v[0] * v[5] - v[2] * v[3]) + 931a989b97SToby Isaac w[2] * (v[1] * v[5] - v[2] * v[4]); 941a989b97SToby Isaac } else { 951a989b97SToby Isaac *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) + 961a989b97SToby Isaac v[1] * (v[5] * v[6] - v[3] * v[8]) + 971a989b97SToby Isaac v[2] * (v[3] * v[7] - v[4] * v[6])); 981a989b97SToby Isaac } 991a989b97SToby Isaac } 1001a989b97SToby Isaac } 1011a989b97SToby Isaac } else { 1021a989b97SToby Isaac PetscInt Nk, Nf; 103fad4db65SToby Isaac PetscInt *subset, *perm; 1041a989b97SToby Isaac PetscInt i, j, l; 1051a989b97SToby Isaac PetscReal sum = 0.; 1061a989b97SToby Isaac 107fad4db65SToby Isaac ierr = PetscDTFactorialInt(k, &Nf);CHKERRQ(ierr); 108fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 109fad4db65SToby Isaac ierr = PetscMalloc2(k, &subset, k, &perm);CHKERRQ(ierr); 1101a989b97SToby Isaac for (i = 0; i < Nk; i++) { 1111a989b97SToby Isaac PetscReal subsum = 0.; 1121a989b97SToby Isaac 1131a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 1141a989b97SToby Isaac for (j = 0; j < Nf; j++) { 1151a989b97SToby Isaac PetscBool permOdd; 1161a989b97SToby Isaac PetscReal prod; 1171a989b97SToby Isaac 118fad4db65SToby Isaac ierr = PetscDTEnumPerm(k, j, perm, &permOdd);CHKERRQ(ierr); 1191a989b97SToby Isaac prod = permOdd ? -1. : 1.; 1201a989b97SToby Isaac for (l = 0; l < k; l++) { 1211a989b97SToby Isaac prod *= v[perm[l] * N + subset[l]]; 1221a989b97SToby Isaac } 1231a989b97SToby Isaac subsum += prod; 1241a989b97SToby Isaac } 1251a989b97SToby Isaac sum += w[i] * subsum; 1261a989b97SToby Isaac } 127fad4db65SToby Isaac ierr = PetscFree2(subset, perm);CHKERRQ(ierr); 1281a989b97SToby Isaac *wv = sum; 1291a989b97SToby Isaac } 1301a989b97SToby Isaac PetscFunctionReturn(0); 1311a989b97SToby Isaac } 1321a989b97SToby Isaac 133fad4db65SToby Isaac /*@ 13428222859SToby Isaac PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form 135fad4db65SToby Isaac 1364165533cSJose E. Roman Input Parameters: 13728222859SToby Isaac + N - the dimension of the vector space, N >= 0 13828222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N 13928222859SToby Isaac . k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N 14028222859SToby Isaac . a - a j-form, size [N choose j] 14128222859SToby Isaac - b - a k-form, size [N choose k] 142fad4db65SToby Isaac 1434165533cSJose E. Roman Output Parameter: 14428222859SToby 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}}), 14528222859SToby 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}. 146fad4db65SToby Isaac 147fad4db65SToby Isaac Level: intermediate 148fad4db65SToby Isaac 14929a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVWedgeMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 150fad4db65SToby Isaac @*/ 1511a989b97SToby Isaac PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) 1521a989b97SToby Isaac { 1531a989b97SToby Isaac PetscInt i; 1541a989b97SToby Isaac PetscErrorCode ierr; 1551a989b97SToby Isaac 1561a989b97SToby Isaac PetscFunctionBegin; 157*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 158*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(j < 0 || k < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 159*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(j + k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 1601a989b97SToby Isaac if (N <= 3) { 1611a989b97SToby Isaac PetscInt Njk; 1621a989b97SToby Isaac 163fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 1641a989b97SToby Isaac if (!j) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[0] * b[i];}} 1651a989b97SToby Isaac else if (!k) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[i] * b[0];}} 1661a989b97SToby Isaac else { 1671a989b97SToby Isaac if (N == 2) {awedgeb[0] = a[0] * b[1] - a[1] * b[0];} 1681a989b97SToby Isaac else { 1691a989b97SToby Isaac if (j+k == 2) { 1701a989b97SToby Isaac awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 1711a989b97SToby Isaac awedgeb[1] = a[0] * b[2] - a[2] * b[0]; 1721a989b97SToby Isaac awedgeb[2] = a[1] * b[2] - a[2] * b[1]; 1731a989b97SToby Isaac } else { 1741a989b97SToby Isaac awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0]; 1751a989b97SToby Isaac } 1761a989b97SToby Isaac } 1771a989b97SToby Isaac } 1781a989b97SToby Isaac } else { 1791a989b97SToby Isaac PetscInt Njk; 1801a989b97SToby Isaac PetscInt JKj; 1811a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 1821a989b97SToby Isaac PetscInt i; 1831a989b97SToby Isaac 184fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 185fad4db65SToby Isaac ierr = PetscDTBinomialInt(j+k, j, &JKj);CHKERRQ(ierr); 1861a989b97SToby Isaac ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr); 1871a989b97SToby Isaac for (i = 0; i < Njk; i++) { 1881a989b97SToby Isaac PetscReal sum = 0.; 1891a989b97SToby Isaac PetscInt l; 1901a989b97SToby Isaac 1911a989b97SToby Isaac ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr); 1921a989b97SToby Isaac for (l = 0; l < JKj; l++) { 1931a989b97SToby Isaac PetscBool jkOdd; 1941a989b97SToby Isaac PetscInt m, jInd, kInd; 1951a989b97SToby Isaac 1961a989b97SToby Isaac ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr); 1971a989b97SToby Isaac for (m = 0; m < j; m++) { 1981a989b97SToby Isaac subsetj[m] = subset[subsetjk[m]]; 1991a989b97SToby Isaac } 2001a989b97SToby Isaac for (m = 0; m < k; m++) { 2011a989b97SToby Isaac subsetk[m] = subset[subsetjk[j+m]]; 2021a989b97SToby Isaac } 2031a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr); 2041a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr); 2051a989b97SToby Isaac sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]); 2061a989b97SToby Isaac } 2071a989b97SToby Isaac awedgeb[i] = sum; 2081a989b97SToby Isaac } 2091a989b97SToby Isaac ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr); 2101a989b97SToby Isaac } 2111a989b97SToby Isaac PetscFunctionReturn(0); 2121a989b97SToby Isaac } 2131a989b97SToby Isaac 214fad4db65SToby Isaac /*@ 21528222859SToby Isaac PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms 216fad4db65SToby Isaac 2174165533cSJose E. Roman Input Parameters: 21828222859SToby Isaac + N - the dimension of the vector space, N >= 0 21928222859SToby Isaac . j - the degree j of the j-form a, 0 <= j <= N 22028222859SToby Isaac . k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N 22128222859SToby Isaac - a - a j-form, size [N choose j] 222fad4db65SToby Isaac 2234165533cSJose E. Roman Output Parameter: 22428222859SToby 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 225fad4db65SToby Isaac 226fad4db65SToby Isaac Level: intermediate 227fad4db65SToby Isaac 22829a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 229fad4db65SToby Isaac @*/ 2301a989b97SToby Isaac PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat) 2311a989b97SToby Isaac { 2321a989b97SToby Isaac PetscInt i; 2331a989b97SToby Isaac PetscErrorCode ierr; 2341a989b97SToby Isaac 2351a989b97SToby Isaac PetscFunctionBegin; 236*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 237*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(j < 0 || k < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 238*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(j + k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 2391a989b97SToby Isaac if (N <= 3) { 2401a989b97SToby Isaac PetscInt Njk; 2411a989b97SToby Isaac 242fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 2431a989b97SToby Isaac if (!j) { 2441a989b97SToby Isaac for (i = 0; i < Njk * Njk; i++) {awedgeMat[i] = 0.;} 2451a989b97SToby Isaac for (i = 0; i < Njk; i++) {awedgeMat[i * (Njk + 1)] = a[0];} 2461a989b97SToby Isaac } else if (!k) { 2471a989b97SToby Isaac for (i = 0; i < Njk; i++) {awedgeMat[i] = a[i];} 2481a989b97SToby Isaac } else { 2491a989b97SToby Isaac if (N == 2) { 2501a989b97SToby Isaac awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; 2511a989b97SToby Isaac } else { 2521a989b97SToby Isaac if (j+k == 2) { 2531a989b97SToby Isaac awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; awedgeMat[2] = 0.; 2541a989b97SToby Isaac awedgeMat[3] = -a[2]; awedgeMat[4] = 0.; awedgeMat[5] = a[0]; 2551a989b97SToby Isaac awedgeMat[6] = 0.; awedgeMat[7] = -a[2]; awedgeMat[8] = a[1]; 2561a989b97SToby Isaac } else { 2571a989b97SToby Isaac awedgeMat[0] = a[2]; awedgeMat[1] = -a[1]; 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 267fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 268fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 269fad4db65SToby Isaac ierr = PetscDTBinomialInt(j+k, j, &JKj);CHKERRQ(ierr); 2701a989b97SToby Isaac ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr); 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 2751a989b97SToby Isaac ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr); 2761a989b97SToby Isaac for (l = 0; l < JKj; l++) { 2771a989b97SToby Isaac PetscBool jkOdd; 2781a989b97SToby Isaac PetscInt m, jInd, kInd; 2791a989b97SToby Isaac 2801a989b97SToby Isaac ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr); 2811a989b97SToby Isaac for (m = 0; m < j; m++) { 2821a989b97SToby Isaac subsetj[m] = subset[subsetjk[m]]; 2831a989b97SToby Isaac } 2841a989b97SToby Isaac for (m = 0; m < k; m++) { 2851a989b97SToby Isaac subsetk[m] = subset[subsetjk[j+m]]; 2861a989b97SToby Isaac } 2871a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr); 2881a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr); 2891a989b97SToby Isaac awedgeMat[i * Nk + kInd] += jkOdd ? - a[jInd] : a[jInd]; 2901a989b97SToby Isaac } 2911a989b97SToby Isaac } 2921a989b97SToby Isaac ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr); 2931a989b97SToby Isaac } 2941a989b97SToby Isaac PetscFunctionReturn(0); 2951a989b97SToby Isaac } 2961a989b97SToby Isaac 297fad4db65SToby Isaac /*@ 29828222859SToby Isaac PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space 299fad4db65SToby Isaac 3004165533cSJose E. Roman Input Parameters: 30128222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, M >= 0 30228222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, N >= 0 30328222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format 30428222859SToby 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). 30528222859SToby Isaac - w - a |k|-form in the image space, size [M choose |k|] 306fad4db65SToby Isaac 3074165533cSJose E. Roman Output Parameter: 30828222859SToby 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). 309fad4db65SToby Isaac 310fad4db65SToby Isaac Level: intermediate 311fad4db65SToby Isaac 312a5b23f4aSJose 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, 31328222859SToby 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, 314fad4db65SToby Isaac then the inverse Hodge star transformation. 315fad4db65SToby Isaac 31629a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullbackMatrix(), PetscDTAltVStar() 317fad4db65SToby Isaac @*/ 3181a989b97SToby Isaac PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw) 3191a989b97SToby Isaac { 3201a989b97SToby Isaac PetscInt i, j, Nk, Mk; 3211a989b97SToby Isaac PetscErrorCode ierr; 3221a989b97SToby Isaac 3231a989b97SToby Isaac PetscFunctionBegin; 324*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 0 || M < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 325*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsInt(k) > N || PetscAbsInt(k) > M,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 3261a989b97SToby Isaac if (N <= 3 && M <= 3) { 3271a989b97SToby Isaac 328fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 329fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 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 3361a989b97SToby Isaac 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 3451a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3461a989b97SToby Isaac sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j]; 3471a989b97SToby Isaac } 3481a989b97SToby Isaac Lstarw[i] = mult[i] * sum; 3491a989b97SToby Isaac } 3501a989b97SToby Isaac } else if (k == 2) { 3511a989b97SToby Isaac PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}}; 3521a989b97SToby Isaac 3531a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3541a989b97SToby Isaac PetscReal sum = 0.; 3551a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3561a989b97SToby Isaac sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] - 3571a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j]; 3581a989b97SToby Isaac } 3591a989b97SToby Isaac Lstarw[i] = sum; 3601a989b97SToby Isaac } 3611a989b97SToby Isaac } else if (k == -2) { 3621a989b97SToby Isaac PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}}; 3631a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 3641a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 3651a989b97SToby Isaac 3661a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3671a989b97SToby Isaac PetscReal sum = 0.; 3681a989b97SToby Isaac 3691a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3701a989b97SToby Isaac sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] * 3711a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][1]] - 3721a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][0]] * 3731a989b97SToby Isaac L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j]; 3741a989b97SToby Isaac 3751a989b97SToby Isaac } 3761a989b97SToby Isaac Lstarw[i] = sum; 3771a989b97SToby Isaac } 3781a989b97SToby Isaac } else { 3791a989b97SToby Isaac PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + 3801a989b97SToby Isaac L[1] * (L[5] * L[6] - L[3] * L[8]) + 3811a989b97SToby Isaac L[2] * (L[3] * L[7] - L[4] * L[6]); 3821a989b97SToby Isaac 3831a989b97SToby Isaac for (i = 0; i < Nk; i++) {Lstarw[i] = detL * w[i];} 3841a989b97SToby Isaac } 3851a989b97SToby Isaac } else { 3861a989b97SToby Isaac PetscInt Nf, l, p; 3871a989b97SToby Isaac PetscReal *Lw, *Lwv; 3881a989b97SToby Isaac PetscInt *subsetw, *subsetv; 389fad4db65SToby Isaac PetscInt *perm; 3901a989b97SToby Isaac PetscReal *walloc = NULL; 3911a989b97SToby Isaac const PetscReal *ww = NULL; 3921a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 3931a989b97SToby Isaac 394fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 395fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 396fad4db65SToby Isaac ierr = PetscDTFactorialInt(PetscAbsInt(k), &Nf);CHKERRQ(ierr); 3971a989b97SToby Isaac if (k < 0) { 3981a989b97SToby Isaac negative = PETSC_TRUE; 3991a989b97SToby Isaac k = -k; 4001a989b97SToby Isaac ierr = PetscMalloc1(Mk, &walloc);CHKERRQ(ierr); 4011a989b97SToby Isaac ierr = PetscDTAltVStar(M, M - k, 1, w, walloc);CHKERRQ(ierr); 4021a989b97SToby Isaac ww = walloc; 4031a989b97SToby Isaac } else { 4041a989b97SToby Isaac ww = w; 4051a989b97SToby Isaac } 406fad4db65SToby Isaac ierr = PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr); 4071a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = 0.; 4081a989b97SToby Isaac for (i = 0; i < Mk; i++) { 4091a989b97SToby Isaac ierr = PetscDTEnumSubset(M, k, i, subsetw);CHKERRQ(ierr); 4101a989b97SToby Isaac for (j = 0; j < Nk; j++) { 4111a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, j, subsetv);CHKERRQ(ierr); 4121a989b97SToby Isaac for (p = 0; p < Nf; p++) { 4131a989b97SToby Isaac PetscReal prod; 4141a989b97SToby Isaac PetscBool isOdd; 4151a989b97SToby Isaac 416fad4db65SToby Isaac ierr = PetscDTEnumPerm(k, p, perm, &isOdd);CHKERRQ(ierr); 4171a989b97SToby Isaac prod = isOdd ? -ww[i] : ww[i]; 4181a989b97SToby Isaac for (l = 0; l < k; l++) { 4191a989b97SToby Isaac prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 4201a989b97SToby Isaac } 4211a989b97SToby Isaac Lstarw[j] += prod; 4221a989b97SToby Isaac } 4231a989b97SToby Isaac } 4241a989b97SToby Isaac } 4251a989b97SToby Isaac if (negative) { 4261a989b97SToby Isaac PetscReal *sLsw; 4271a989b97SToby Isaac 4281a989b97SToby Isaac ierr = PetscMalloc1(Nk, &sLsw);CHKERRQ(ierr); 4291a989b97SToby Isaac ierr = PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw);CHKERRQ(ierr); 4301a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i]; 4311a989b97SToby Isaac ierr = PetscFree(sLsw);CHKERRQ(ierr); 4321a989b97SToby Isaac } 433fad4db65SToby Isaac ierr = PetscFree5(subsetw, subsetv, perm, Lw, Lwv);CHKERRQ(ierr); 4341a989b97SToby Isaac ierr = PetscFree(walloc);CHKERRQ(ierr); 4351a989b97SToby Isaac } 4361a989b97SToby Isaac PetscFunctionReturn(0); 4371a989b97SToby Isaac } 4381a989b97SToby Isaac 439fad4db65SToby Isaac /*@ 440fad4db65SToby Isaac PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation 441fad4db65SToby Isaac 4424165533cSJose E. Roman Input Parameters: 44328222859SToby Isaac + N - the dimension of the origin vector space of the linear transformation, N >= 0 44428222859SToby Isaac . M - the dimension of the image vector space of the linear transformation, M >= 0 44528222859SToby Isaac . L - a linear transformation, an [M x N] matrix in row-major format 44628222859SToby 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()) 447fad4db65SToby Isaac 4484165533cSJose E. Roman Output Parameter: 44928222859SToby Isaac . Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w 450fad4db65SToby Isaac 451fad4db65SToby Isaac Level: intermediate 452fad4db65SToby Isaac 45329a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVStar() 454fad4db65SToby Isaac @*/ 4551a989b97SToby Isaac PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) 4561a989b97SToby Isaac { 4571a989b97SToby Isaac PetscInt Nk, Mk, Nf, i, j, l, p; 4581a989b97SToby Isaac PetscReal *Lw, *Lwv; 4591a989b97SToby Isaac PetscInt *subsetw, *subsetv; 460fad4db65SToby Isaac PetscInt *perm; 4611a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 4621a989b97SToby Isaac PetscErrorCode ierr; 4631a989b97SToby Isaac 4641a989b97SToby Isaac PetscFunctionBegin; 465*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 0 || M < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 466*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsInt(k) > N || PetscAbsInt(k) > M,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 4671a989b97SToby Isaac if (N <= 3 && M <= 3) { 4681a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 4691a989b97SToby Isaac 470fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 471fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 4721a989b97SToby Isaac if (!k) { 4731a989b97SToby Isaac Lstar[0] = 1.; 4741a989b97SToby Isaac } else if (k == 1) { 4751a989b97SToby Isaac for (i = 0; i < Nk; i++) {for (j = 0; j < Mk; j++) {Lstar[i * Mk + j] = L[j * Nk + i];}} 4761a989b97SToby Isaac } else if (k == -1) { 4771a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4781a989b97SToby Isaac for (j = 0; j < Mk; j++) { 4791a989b97SToby Isaac Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j]; 4801a989b97SToby Isaac } 4811a989b97SToby Isaac } 4821a989b97SToby Isaac } else if (k == 2) { 4831a989b97SToby Isaac PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}}; 4841a989b97SToby Isaac 4851a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4861a989b97SToby Isaac for (j = 0; j < Mk; j++) { 4871a989b97SToby Isaac Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] * 4881a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][1]] - 4891a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][0]] * 4901a989b97SToby Isaac L[pairs[j][0] * N + pairs[i][1]]; 4911a989b97SToby Isaac } 4921a989b97SToby Isaac } 4931a989b97SToby Isaac } else if (k == -2) { 4941a989b97SToby Isaac PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}}; 4951a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 4961a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 4971a989b97SToby Isaac 4981a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4991a989b97SToby Isaac for (j = 0; j < Mk; j++) { 5001a989b97SToby Isaac Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] * 5011a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][1]] - 5021a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][0]] * 5031a989b97SToby Isaac L[pairs[offj + j][0] * N + pairs[offi + i][1]]; 5041a989b97SToby Isaac } 5051a989b97SToby Isaac } 5061a989b97SToby Isaac } else { 5071a989b97SToby Isaac PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + 5081a989b97SToby Isaac L[1] * (L[5] * L[6] - L[3] * L[8]) + 5091a989b97SToby Isaac L[2] * (L[3] * L[7] - L[4] * L[6]); 5101a989b97SToby Isaac 5111a989b97SToby Isaac for (i = 0; i < Nk; i++) {Lstar[i] = detL;} 5121a989b97SToby Isaac } 5131a989b97SToby Isaac } else { 5141a989b97SToby Isaac if (k < 0) { 5151a989b97SToby Isaac negative = PETSC_TRUE; 5161a989b97SToby Isaac k = -k; 5171a989b97SToby Isaac } 518fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 519fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 520fad4db65SToby Isaac ierr = PetscDTFactorialInt(PetscAbsInt(k), &Nf);CHKERRQ(ierr); 521fad4db65SToby Isaac ierr = PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr); 5221a989b97SToby Isaac for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.; 5231a989b97SToby Isaac for (i = 0; i < Mk; i++) { 5241a989b97SToby Isaac PetscBool iOdd; 5251a989b97SToby Isaac PetscInt iidx, jidx; 5261a989b97SToby Isaac 5271a989b97SToby Isaac ierr = PetscDTEnumSplit(M, k, i, subsetw, &iOdd);CHKERRQ(ierr); 5281a989b97SToby Isaac iidx = negative ? Mk - 1 - i : i; 52928222859SToby Isaac iOdd = negative ? (PetscBool) (iOdd ^ ((k * (M-k)) & 1)) : PETSC_FALSE; 5301a989b97SToby Isaac for (j = 0; j < Nk; j++) { 5311a989b97SToby Isaac PetscBool jOdd; 5321a989b97SToby Isaac 5331a989b97SToby Isaac ierr = PetscDTEnumSplit(N, k, j, subsetv, &jOdd);CHKERRQ(ierr); 5341a989b97SToby Isaac jidx = negative ? Nk - 1 - j : j; 53528222859SToby Isaac jOdd = negative ? (PetscBool) (iOdd ^ jOdd ^ ((k * (N-k)) & 1)) : PETSC_FALSE; 5361a989b97SToby Isaac for (p = 0; p < Nf; p++) { 5371a989b97SToby Isaac PetscReal prod; 5381a989b97SToby Isaac PetscBool isOdd; 5391a989b97SToby Isaac 540fad4db65SToby Isaac ierr = PetscDTEnumPerm(k, p, perm, &isOdd);CHKERRQ(ierr); 54128222859SToby Isaac isOdd = (PetscBool) (isOdd ^ jOdd); 5421a989b97SToby Isaac prod = isOdd ? -1. : 1.; 5431a989b97SToby Isaac for (l = 0; l < k; l++) { 5441a989b97SToby Isaac prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 5451a989b97SToby Isaac } 5461a989b97SToby Isaac Lstar[jidx * Mk + iidx] += prod; 5471a989b97SToby Isaac } 5481a989b97SToby Isaac } 5491a989b97SToby Isaac } 550fad4db65SToby Isaac ierr = PetscFree5(subsetw, subsetv, perm, Lw, Lwv);CHKERRQ(ierr); 5511a989b97SToby Isaac } 5521a989b97SToby Isaac PetscFunctionReturn(0); 5531a989b97SToby Isaac } 5541a989b97SToby Isaac 555fad4db65SToby Isaac /*@ 55628222859SToby Isaac PetscDTAltVInterior - Compute the interior product of a k-form with a vector 557fad4db65SToby Isaac 5584165533cSJose E. Roman Input Parameters: 55928222859SToby Isaac + N - the dimension of the vector space, N >= 0 56028222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 56128222859SToby Isaac . w - a k-form, size [N choose k] 56228222859SToby Isaac - v - an N dimensional vector 563fad4db65SToby Isaac 5644165533cSJose E. Roman Output Parameter: 56528222859SToby 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}). 566fad4db65SToby Isaac 567fad4db65SToby Isaac Level: intermediate 568fad4db65SToby Isaac 56929a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVInteriorMatrix(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 570fad4db65SToby Isaac @*/ 5711a989b97SToby Isaac PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv) 5721a989b97SToby Isaac { 5731a989b97SToby Isaac PetscInt i, Nk, Nkm; 5741a989b97SToby Isaac PetscErrorCode ierr; 5751a989b97SToby Isaac 5761a989b97SToby Isaac PetscFunctionBegin; 577*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k <= 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 578fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 579fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr); 5801a989b97SToby Isaac if (N <= 3) { 5811a989b97SToby Isaac if (k == 1) { 5821a989b97SToby Isaac PetscReal sum = 0.; 5831a989b97SToby Isaac 5841a989b97SToby Isaac for (i = 0; i < N; i++) { 5851a989b97SToby Isaac sum += w[i] * v[i]; 5861a989b97SToby Isaac } 5871a989b97SToby Isaac wIntv[0] = sum; 5881a989b97SToby Isaac } else if (k == N) { 5891a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 5901a989b97SToby Isaac 5911a989b97SToby Isaac for (i = 0; i < N; i++) { 5921a989b97SToby Isaac wIntv[N - 1 - i] = w[0] * v[i] * mult[i]; 5931a989b97SToby Isaac } 5941a989b97SToby Isaac } else { 5951a989b97SToby Isaac wIntv[0] = - w[0]*v[1] - w[1]*v[2]; 5961a989b97SToby Isaac wIntv[1] = w[0]*v[0] - w[2]*v[2]; 5971a989b97SToby Isaac wIntv[2] = w[1]*v[0] + w[2]*v[1]; 5981a989b97SToby Isaac } 5991a989b97SToby Isaac } else { 6001a989b97SToby Isaac PetscInt *subset, *work; 6011a989b97SToby Isaac 6021a989b97SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 6031a989b97SToby Isaac for (i = 0; i < Nkm; i++) wIntv[i] = 0.; 6041a989b97SToby Isaac for (i = 0; i < Nk; i++) { 6051a989b97SToby Isaac PetscInt j, l, m; 6061a989b97SToby Isaac 6071a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 6081a989b97SToby Isaac for (j = 0; j < k; j++) { 6091a989b97SToby Isaac PetscInt idx; 61028222859SToby Isaac PetscBool flip = (PetscBool) (j & 1); 6111a989b97SToby Isaac 6121a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 6131a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 6141a989b97SToby Isaac } 6151a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 6161a989b97SToby Isaac wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]); 6171a989b97SToby Isaac } 6181a989b97SToby Isaac } 6191a989b97SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 6201a989b97SToby Isaac } 6211a989b97SToby Isaac PetscFunctionReturn(0); 6221a989b97SToby Isaac } 6231a989b97SToby Isaac 624fad4db65SToby Isaac /*@ 62528222859SToby Isaac PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector 626fad4db65SToby Isaac 6274165533cSJose E. Roman Input Parameters: 62828222859SToby Isaac + N - the dimension of the vector space, N >= 0 62928222859SToby Isaac . k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N 63028222859SToby Isaac - v - an N dimensional vector 631fad4db65SToby Isaac 6324165533cSJose E. Roman Output Parameter: 633fad4db65SToby Isaac . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v) 634fad4db65SToby Isaac 635fad4db65SToby Isaac Level: intermediate 636fad4db65SToby Isaac 63729a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVInterior(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 638fad4db65SToby Isaac @*/ 6391a989b97SToby Isaac PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat) 6401a989b97SToby Isaac { 6411a989b97SToby Isaac PetscInt i, Nk, Nkm; 6421a989b97SToby Isaac PetscErrorCode ierr; 6431a989b97SToby Isaac 6441a989b97SToby Isaac PetscFunctionBegin; 645*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k <= 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 646fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 647fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr); 6481a989b97SToby Isaac if (N <= 3) { 6491a989b97SToby Isaac if (k == 1) { 6501a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[i] = v[i]; 6511a989b97SToby Isaac } else if (k == N) { 6521a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 6531a989b97SToby Isaac 6541a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i]; 6551a989b97SToby Isaac } else { 6561a989b97SToby Isaac intvMat[0] = -v[1]; intvMat[1] = -v[2]; intvMat[2] = 0.; 6571a989b97SToby Isaac intvMat[3] = v[0]; intvMat[4] = 0.; intvMat[5] = -v[2]; 6581a989b97SToby Isaac intvMat[6] = 0.; intvMat[7] = v[0]; intvMat[8] = v[1]; 6591a989b97SToby Isaac } 6601a989b97SToby Isaac } else { 6611a989b97SToby Isaac PetscInt *subset, *work; 6621a989b97SToby Isaac 6631a989b97SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 6641a989b97SToby Isaac for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.; 6651a989b97SToby Isaac for (i = 0; i < Nk; i++) { 6661a989b97SToby Isaac PetscInt j, l, m; 6671a989b97SToby Isaac 6681a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 6691a989b97SToby Isaac for (j = 0; j < k; j++) { 6701a989b97SToby Isaac PetscInt idx; 67128222859SToby Isaac PetscBool flip = (PetscBool) (j & 1); 6721a989b97SToby Isaac 6731a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 6741a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 6751a989b97SToby Isaac } 6761a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 6771a989b97SToby Isaac intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]]; 6781a989b97SToby Isaac } 6791a989b97SToby Isaac } 6801a989b97SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 6811a989b97SToby Isaac } 6821a989b97SToby Isaac PetscFunctionReturn(0); 6831a989b97SToby Isaac } 6841a989b97SToby Isaac 685fad4db65SToby Isaac /*@ 686fad4db65SToby Isaac PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix() 687fad4db65SToby Isaac 6884165533cSJose E. Roman Input Parameters: 68928222859SToby Isaac + N - the dimension of the vector space, N >= 0 69028222859SToby Isaac - k - the degree of the k-forms on which intvMat from PetscDTAltVInteriorMatrix() acts, 0 <= k <= N. 691fad4db65SToby Isaac 6924165533cSJose E. Roman Output Parameter: 69328222859SToby Isaac . indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k 69428222859SToby 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 69528222859SToby Isaac coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1) 696fad4db65SToby Isaac 697fad4db65SToby Isaac Level: intermediate 698fad4db65SToby Isaac 69928222859SToby Isaac Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential 700fad4db65SToby Isaac 70129a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVInterior(), PetscDTAltVInteriorMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 702fad4db65SToby Isaac @*/ 703dda711d0SToby Isaac PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3]) 704dda711d0SToby Isaac { 705dda711d0SToby Isaac PetscInt i, Nk, Nkm; 706dda711d0SToby Isaac PetscErrorCode ierr; 707dda711d0SToby Isaac 708dda711d0SToby Isaac PetscFunctionBegin; 709*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k <= 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 710fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 711fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr); 712dda711d0SToby Isaac if (N <= 3) { 713dda711d0SToby Isaac if (k == 1) { 714dda711d0SToby Isaac for (i = 0; i < N; i++) { 715dda711d0SToby Isaac indices[i][0] = 0; 716dda711d0SToby Isaac indices[i][1] = i; 717dda711d0SToby Isaac indices[i][2] = i; 718dda711d0SToby Isaac } 719dda711d0SToby Isaac } else if (k == N) { 720dda711d0SToby Isaac PetscInt val[3] = {0, -2, 2}; 721dda711d0SToby Isaac 722dda711d0SToby Isaac for (i = 0; i < N; i++) { 723dda711d0SToby Isaac indices[i][0] = N - 1 - i; 724dda711d0SToby Isaac indices[i][1] = 0; 725dda711d0SToby Isaac indices[i][2] = val[i]; 726dda711d0SToby Isaac } 727dda711d0SToby Isaac } else { 728dda711d0SToby Isaac indices[0][0] = 0; indices[0][1] = 0; indices[0][2] = -(1 + 1); 729dda711d0SToby Isaac indices[1][0] = 0; indices[1][1] = 1; indices[1][2] = -(2 + 1); 730dda711d0SToby Isaac indices[2][0] = 1; indices[2][1] = 0; indices[2][2] = 0; 731dda711d0SToby Isaac indices[3][0] = 1; indices[3][1] = 2; indices[3][2] = -(2 + 1); 732dda711d0SToby Isaac indices[4][0] = 2; indices[4][1] = 1; indices[4][2] = 0; 733dda711d0SToby Isaac indices[5][0] = 2; indices[5][1] = 2; indices[5][2] = 1; 734dda711d0SToby Isaac } 735dda711d0SToby Isaac } else { 736dda711d0SToby Isaac PetscInt *subset, *work; 737dda711d0SToby Isaac 738dda711d0SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 739dda711d0SToby Isaac for (i = 0; i < Nk; i++) { 740dda711d0SToby Isaac PetscInt j, l, m; 741dda711d0SToby Isaac 742dda711d0SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 743dda711d0SToby Isaac for (j = 0; j < k; j++) { 744dda711d0SToby Isaac PetscInt idx; 74528222859SToby Isaac PetscBool flip = (PetscBool) (j & 1); 746dda711d0SToby Isaac 747dda711d0SToby Isaac for (l = 0, m = 0; l < k; l++) { 748dda711d0SToby Isaac if (l != j) work[m++] = subset[l]; 749dda711d0SToby Isaac } 750dda711d0SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 751dda711d0SToby Isaac indices[i * k + j][0] = idx; 752dda711d0SToby Isaac indices[i * k + j][1] = i; 753dda711d0SToby Isaac indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j]; 754dda711d0SToby Isaac } 755dda711d0SToby Isaac } 756dda711d0SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 757dda711d0SToby Isaac } 758dda711d0SToby Isaac PetscFunctionReturn(0); 759dda711d0SToby Isaac } 760dda711d0SToby Isaac 761fad4db65SToby Isaac /*@ 76228222859SToby Isaac PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form 763fad4db65SToby Isaac 7644165533cSJose E. Roman Input Parameters: 76528222859SToby Isaac + N - the dimension of the vector space, N >= 0 76628222859SToby Isaac . k - the degree k of the k-form w, 0 <= k <= N 76728222859SToby 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. 76828222859SToby Isaac - w - a k-form, size [N choose k] 769fad4db65SToby Isaac 7704165533cSJose E. Roman Output Parameter: 77128222859SToby 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. 772fad4db65SToby Isaac 773fad4db65SToby Isaac Level: intermediate 774fad4db65SToby Isaac 77529a920c6SToby Isaac .seealso: PetscDTAltV, PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 776fad4db65SToby Isaac @*/ 7771a989b97SToby Isaac PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw) 7781a989b97SToby Isaac { 7791a989b97SToby Isaac PetscInt Nk, i; 7801a989b97SToby Isaac PetscErrorCode ierr; 7811a989b97SToby Isaac 7821a989b97SToby Isaac PetscFunctionBegin; 783*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k < 0 || k > N,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 784fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 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 8021a989b97SToby Isaac ierr = PetscMalloc1(N, &subset);CHKERRQ(ierr); 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 8091a989b97SToby Isaac ierr = PetscDTEnumSplit(N, l, i, subset, &sOdd);CHKERRQ(ierr); 8101a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, l, subset, &idx);CHKERRQ(ierr); 8111a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, N-l, &subset[l], &j);CHKERRQ(ierr); 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 } 8211a989b97SToby Isaac ierr = PetscFree(subset);CHKERRQ(ierr); 8221a989b97SToby Isaac } 8231a989b97SToby Isaac PetscFunctionReturn(0); 8241a989b97SToby Isaac } 825