1 #include <petsc/private/petscimpl.h> 2 #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/ 3 4 /*MC 5 PetscDTAltV - An interface for common operations on k-forms, also known as alternating algebraic forms or alternating k-linear maps. 6 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*. 7 8 A recommended reference for this material is Section 2 "Exterior algebra and exterior calculus" in "Finite element 9 exterior calculus, homological techniques, and applications", by Arnold, Falk, & Winther (2006, doi:10.1017/S0962492906210018). 10 11 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 12 vectors from a vector space V and producing a real number: 13 - 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) 14 - k-linear: w acts linear in each vector separately, e.g. w(a*v + b*y, v_2, ..., v_k) = a*w(v,v_2,...,v_k) + b*w(y,v_2,...,v_k) 15 This action is implemented as PetscDTAltVApply. 16 17 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 18 shows that for an N dimensional space, only 0 <= k <= N are valid form degrees.) 19 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: 20 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. 21 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 22 by the associated subsets. 23 24 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 25 V[i,j] = v_i[c_j] and taking the determinant of V. 26 27 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). 28 This is an anticommutative product, (f wedge g) = -(g wedge f). It is sufficient to describe the wedge product of two basis forms. 29 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): 30 - If there is any coordinate in both sets, then (f wedge g) = 0. 31 - Otherwise, (f wedge g) is a multiple of the basis (j+k)-form h associated with (f_1,...,f_j,g_1,...,g_k). 32 - In fact it is equal to either h or -h depending on how (f_1,...,f_j,g_1,...,g_k) compares to the same list of coordinates given in ascending order: if it is an even permutation of that list, then (f wedge g) = h, otherwise (f wedge g) = -h. 33 The wedge product is implemented for either two inputs (f and g) in PetscDTAltVWedge, or for one (just f, giving a 34 matrix to multiply against multiple choices of g) in PetscDTAltVWedgeMatrix. 35 36 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), 37 defined by (w int v)(v_1,...,v_{k-1}) = w(v,v_1,...,v_{k-1}). 38 39 The interior product is implemented for either two inputs (w and v) in PetscDTAltVInterior, for one (just v, giving a 40 matrix to multiply against multiple choices of w) in PetscDTAltVInteriorMatrix, 41 or for no inputs (giving the sparsity pattern of PetscDTAltVInteriorMatrix) in PetscDTAltVInteriorPattern. 42 43 When there is a linear map L: V -> W from an N dimensional vector space to an M dimensional vector space, 44 it induces the linear pullback map L^* : Alt^k W -> Alt^k V, defined by L^* w(v_1,...,v_k) = w(L v_1, ..., L v_k). 45 The pullback is implemented as PetscDTAltVPullback (acting on a known w) and PetscDTAltVPullbackMatrix (creating a matrix that computes the actin of L^*). 46 47 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 48 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 49 of the basis coefficients of f and g. 50 Powers of the Hodge star operator can be applied with PetscDTAltVStar 51 52 Level: intermediate 53 54 .seealso: `PetscDTAltVApply()`, `PetscDTAltVWedge()`, `PetscDTAltVInterior()`, `PetscDTAltVPullback()`, `PetscDTAltVStar()` 55 M*/ 56 57 /*@ 58 PetscDTAltVApply - Apply an a k-form (an alternating k-linear map) to a set of k N-dimensional vectors 59 60 Input Parameters: 61 + N - the dimension of the vector space, N >= 0 62 . k - the degree k of the k-form w, 0 <= k <= N 63 . 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) 64 - v - a set of k vectors of size N, size [k x N], each vector stored contiguously 65 66 Output Parameter: 67 . 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 68 69 Level: intermediate 70 71 .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 72 @*/ 73 PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv) 74 { 75 PetscFunctionBegin; 76 PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 77 PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 78 if (N <= 3) { 79 if (!k) { 80 *wv = w[0]; 81 } else { 82 if (N == 1) { 83 *wv = w[0] * v[0]; 84 } else if (N == 2) { 85 if (k == 1) { 86 *wv = w[0] * v[0] + w[1] * v[1]; 87 } else { 88 *wv = w[0] * (v[0] * v[3] - v[1] * v[2]); 89 } 90 } else { 91 if (k == 1) { 92 *wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2]; 93 } else if (k == 2) { 94 *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]); 95 } else { 96 *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])); 97 } 98 } 99 } 100 } else { 101 PetscInt Nk, Nf; 102 PetscInt *subset, *perm; 103 PetscInt i, j, l; 104 PetscReal sum = 0.; 105 106 PetscCall(PetscDTFactorialInt(k, &Nf)); 107 PetscCall(PetscDTBinomialInt(N, k, &Nk)); 108 PetscCall(PetscMalloc2(k, &subset, k, &perm)); 109 for (i = 0; i < Nk; i++) { 110 PetscReal subsum = 0.; 111 112 PetscCall(PetscDTEnumSubset(N, k, i, subset)); 113 for (j = 0; j < Nf; j++) { 114 PetscBool permOdd; 115 PetscReal prod; 116 117 PetscCall(PetscDTEnumPerm(k, j, perm, &permOdd)); 118 prod = permOdd ? -1. : 1.; 119 for (l = 0; l < k; l++) prod *= v[perm[l] * N + subset[l]]; 120 subsum += prod; 121 } 122 sum += w[i] * subsum; 123 } 124 PetscCall(PetscFree2(subset, perm)); 125 *wv = sum; 126 } 127 PetscFunctionReturn(0); 128 } 129 130 /*@ 131 PetscDTAltVWedge - Compute the wedge product of a j-form and a k-form, giving a (j+k) form 132 133 Input Parameters: 134 + N - the dimension of the vector space, N >= 0 135 . j - the degree j of the j-form a, 0 <= j <= N 136 . k - the degree k of the k-form b, 0 <= k <= N and 0 <= j+k <= N 137 . a - a j-form, size [N choose j] 138 - b - a k-form, size [N choose k] 139 140 Output Parameter: 141 . 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}}), 142 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}. 143 144 Level: intermediate 145 146 .seealso: `PetscDTAltV`, `PetscDTAltVWedgeMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 147 @*/ 148 PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) 149 { 150 PetscInt i; 151 152 PetscFunctionBegin; 153 PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 154 PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 155 PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 156 if (N <= 3) { 157 PetscInt Njk; 158 159 PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 160 if (!j) { 161 for (i = 0; i < Njk; i++) awedgeb[i] = a[0] * b[i]; 162 } else if (!k) { 163 for (i = 0; i < Njk; i++) awedgeb[i] = a[i] * b[0]; 164 } else { 165 if (N == 2) { 166 awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 167 } else { 168 if (j + k == 2) { 169 awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 170 awedgeb[1] = a[0] * b[2] - a[2] * b[0]; 171 awedgeb[2] = a[1] * b[2] - a[2] * b[1]; 172 } else { 173 awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0]; 174 } 175 } 176 } 177 } else { 178 PetscInt Njk; 179 PetscInt JKj; 180 PetscInt *subset, *subsetjk, *subsetj, *subsetk; 181 PetscInt i; 182 183 PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 184 PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 185 PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk)); 186 for (i = 0; i < Njk; i++) { 187 PetscReal sum = 0.; 188 PetscInt l; 189 190 PetscCall(PetscDTEnumSubset(N, j + k, i, subset)); 191 for (l = 0; l < JKj; l++) { 192 PetscBool jkOdd; 193 PetscInt m, jInd, kInd; 194 195 PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd)); 196 for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]]; 197 for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]]; 198 PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd)); 199 PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd)); 200 sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]); 201 } 202 awedgeb[i] = sum; 203 } 204 PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk)); 205 } 206 PetscFunctionReturn(0); 207 } 208 209 /*@ 210 PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with a given j-form that maps k-forms to (j+k)-forms 211 212 Input Parameters: 213 + N - the dimension of the vector space, N >= 0 214 . j - the degree j of the j-form a, 0 <= j <= N 215 . k - the degree k of the k-forms that (a wedge) will be applied to, 0 <= k <= N and 0 <= j+k <= N 216 - a - a j-form, size [N choose j] 217 218 Output Parameter: 219 . 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 220 221 Level: intermediate 222 223 .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 224 @*/ 225 PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat) 226 { 227 PetscInt i; 228 229 PetscFunctionBegin; 230 PetscCheck(N >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 231 PetscCheck(j >= 0 && k >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 232 PetscCheck(j + k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 233 if (N <= 3) { 234 PetscInt Njk; 235 236 PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 237 if (!j) { 238 for (i = 0; i < Njk * Njk; i++) awedgeMat[i] = 0.; 239 for (i = 0; i < Njk; i++) awedgeMat[i * (Njk + 1)] = a[0]; 240 } else if (!k) { 241 for (i = 0; i < Njk; i++) awedgeMat[i] = a[i]; 242 } else { 243 if (N == 2) { 244 awedgeMat[0] = -a[1]; 245 awedgeMat[1] = a[0]; 246 } else { 247 if (j + k == 2) { 248 awedgeMat[0] = -a[1]; 249 awedgeMat[1] = a[0]; 250 awedgeMat[2] = 0.; 251 awedgeMat[3] = -a[2]; 252 awedgeMat[4] = 0.; 253 awedgeMat[5] = a[0]; 254 awedgeMat[6] = 0.; 255 awedgeMat[7] = -a[2]; 256 awedgeMat[8] = a[1]; 257 } else { 258 awedgeMat[0] = a[2]; 259 awedgeMat[1] = -a[1]; 260 awedgeMat[2] = a[0]; 261 } 262 } 263 } 264 } else { 265 PetscInt Njk; 266 PetscInt Nk; 267 PetscInt JKj, i; 268 PetscInt *subset, *subsetjk, *subsetj, *subsetk; 269 270 PetscCall(PetscDTBinomialInt(N, k, &Nk)); 271 PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 272 PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 273 PetscCall(PetscMalloc4(j + k, &subset, j + k, &subsetjk, j, &subsetj, k, &subsetk)); 274 for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.; 275 for (i = 0; i < Njk; i++) { 276 PetscInt l; 277 278 PetscCall(PetscDTEnumSubset(N, j + k, i, subset)); 279 for (l = 0; l < JKj; l++) { 280 PetscBool jkOdd; 281 PetscInt m, jInd, kInd; 282 283 PetscCall(PetscDTEnumSplit(j + k, j, l, subsetjk, &jkOdd)); 284 for (m = 0; m < j; m++) subsetj[m] = subset[subsetjk[m]]; 285 for (m = 0; m < k; m++) subsetk[m] = subset[subsetjk[j + m]]; 286 PetscCall(PetscDTSubsetIndex(N, j, subsetj, &jInd)); 287 PetscCall(PetscDTSubsetIndex(N, k, subsetk, &kInd)); 288 awedgeMat[i * Nk + kInd] += jkOdd ? -a[jInd] : a[jInd]; 289 } 290 } 291 PetscCall(PetscFree4(subset, subsetjk, subsetj, subsetk)); 292 } 293 PetscFunctionReturn(0); 294 } 295 296 /*@ 297 PetscDTAltVPullback - Compute the pullback of a k-form under a linear transformation of the coordinate space 298 299 Input Parameters: 300 + N - the dimension of the origin vector space of the linear transformation, M >= 0 301 . M - the dimension of the image vector space of the linear transformation, N >= 0 302 . L - a linear transformation, an [M x N] matrix in row-major format 303 . 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). 304 - w - a |k|-form in the image space, size [M choose |k|] 305 306 Output Parameter: 307 . 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). 308 309 Level: intermediate 310 311 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, 312 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, 313 then the inverse Hodge star transformation. 314 315 .seealso: `PetscDTAltV`, `PetscDTAltVPullbackMatrix()`, `PetscDTAltVStar()` 316 @*/ 317 PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw) 318 { 319 PetscInt i, j, Nk, Mk; 320 321 PetscFunctionBegin; 322 PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 323 PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 324 if (N <= 3 && M <= 3) { 325 PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 326 PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 327 if (!k) { 328 Lstarw[0] = w[0]; 329 } else if (k == 1) { 330 for (i = 0; i < Nk; i++) { 331 PetscReal sum = 0.; 332 333 for (j = 0; j < Mk; j++) sum += L[j * Nk + i] * w[j]; 334 Lstarw[i] = sum; 335 } 336 } else if (k == -1) { 337 PetscReal mult[3] = {1., -1., 1.}; 338 339 for (i = 0; i < Nk; i++) { 340 PetscReal sum = 0.; 341 342 for (j = 0; j < Mk; j++) sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j]; 343 Lstarw[i] = mult[i] * sum; 344 } 345 } else if (k == 2) { 346 PetscInt pairs[3][2] = { 347 {0, 1}, 348 {0, 2}, 349 {1, 2} 350 }; 351 352 for (i = 0; i < Nk; i++) { 353 PetscReal sum = 0.; 354 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]; 355 Lstarw[i] = sum; 356 } 357 } else if (k == -2) { 358 PetscInt pairs[3][2] = { 359 {1, 2}, 360 {2, 0}, 361 {0, 1} 362 }; 363 PetscInt offi = (N == 2) ? 2 : 0; 364 PetscInt offj = (M == 2) ? 2 : 0; 365 366 for (i = 0; i < Nk; i++) { 367 PetscReal sum = 0.; 368 369 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]; 370 Lstarw[i] = sum; 371 } 372 } else { 373 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]); 374 375 for (i = 0; i < Nk; i++) Lstarw[i] = detL * w[i]; 376 } 377 } else { 378 PetscInt Nf, l, p; 379 PetscReal *Lw, *Lwv; 380 PetscInt *subsetw, *subsetv; 381 PetscInt *perm; 382 PetscReal *walloc = NULL; 383 const PetscReal *ww = NULL; 384 PetscBool negative = PETSC_FALSE; 385 386 PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 387 PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 388 PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf)); 389 if (k < 0) { 390 negative = PETSC_TRUE; 391 k = -k; 392 PetscCall(PetscMalloc1(Mk, &walloc)); 393 PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc)); 394 ww = walloc; 395 } else { 396 ww = w; 397 } 398 PetscCall(PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv)); 399 for (i = 0; i < Nk; i++) Lstarw[i] = 0.; 400 for (i = 0; i < Mk; i++) { 401 PetscCall(PetscDTEnumSubset(M, k, i, subsetw)); 402 for (j = 0; j < Nk; j++) { 403 PetscCall(PetscDTEnumSubset(N, k, j, subsetv)); 404 for (p = 0; p < Nf; p++) { 405 PetscReal prod; 406 PetscBool isOdd; 407 408 PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd)); 409 prod = isOdd ? -ww[i] : ww[i]; 410 for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 411 Lstarw[j] += prod; 412 } 413 } 414 } 415 if (negative) { 416 PetscReal *sLsw; 417 418 PetscCall(PetscMalloc1(Nk, &sLsw)); 419 PetscCall(PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw)); 420 for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i]; 421 PetscCall(PetscFree(sLsw)); 422 } 423 PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv)); 424 PetscCall(PetscFree(walloc)); 425 } 426 PetscFunctionReturn(0); 427 } 428 429 /*@ 430 PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation 431 432 Input Parameters: 433 + N - the dimension of the origin vector space of the linear transformation, N >= 0 434 . M - the dimension of the image vector space of the linear transformation, M >= 0 435 . L - a linear transformation, an [M x N] matrix in row-major format 436 - 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()) 437 438 Output Parameter: 439 . Lstar - the pullback matrix, an [(N choose |k|) x (M choose |k|)] matrix in row-major format such that Lstar * w = L^* w 440 441 Level: intermediate 442 443 .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVStar()` 444 @*/ 445 PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) 446 { 447 PetscInt Nk, Mk, Nf, i, j, l, p; 448 PetscReal *Lw, *Lwv; 449 PetscInt *subsetw, *subsetv; 450 PetscInt *perm; 451 PetscBool negative = PETSC_FALSE; 452 453 PetscFunctionBegin; 454 PetscCheck(N >= 0 && M >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 455 PetscCheck(PetscAbsInt(k) <= N && PetscAbsInt(k) <= M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 456 if (N <= 3 && M <= 3) { 457 PetscReal mult[3] = {1., -1., 1.}; 458 459 PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 460 PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 461 if (!k) { 462 Lstar[0] = 1.; 463 } else if (k == 1) { 464 for (i = 0; i < Nk; i++) { 465 for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[j * Nk + i]; 466 } 467 } else if (k == -1) { 468 for (i = 0; i < Nk; i++) { 469 for (j = 0; j < Mk; j++) Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j]; 470 } 471 } else if (k == 2) { 472 PetscInt pairs[3][2] = { 473 {0, 1}, 474 {0, 2}, 475 {1, 2} 476 }; 477 478 for (i = 0; i < Nk; i++) { 479 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]]; 480 } 481 } else if (k == -2) { 482 PetscInt pairs[3][2] = { 483 {1, 2}, 484 {2, 0}, 485 {0, 1} 486 }; 487 PetscInt offi = (N == 2) ? 2 : 0; 488 PetscInt offj = (M == 2) ? 2 : 0; 489 490 for (i = 0; i < Nk; i++) { 491 for (j = 0; j < Mk; j++) { 492 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]]; 493 } 494 } 495 } else { 496 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]); 497 498 for (i = 0; i < Nk; i++) Lstar[i] = detL; 499 } 500 } else { 501 if (k < 0) { 502 negative = PETSC_TRUE; 503 k = -k; 504 } 505 PetscCall(PetscDTBinomialInt(M, PetscAbsInt(k), &Mk)); 506 PetscCall(PetscDTBinomialInt(N, PetscAbsInt(k), &Nk)); 507 PetscCall(PetscDTFactorialInt(PetscAbsInt(k), &Nf)); 508 PetscCall(PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv)); 509 for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.; 510 for (i = 0; i < Mk; i++) { 511 PetscBool iOdd; 512 PetscInt iidx, jidx; 513 514 PetscCall(PetscDTEnumSplit(M, k, i, subsetw, &iOdd)); 515 iidx = negative ? Mk - 1 - i : i; 516 iOdd = negative ? (PetscBool)(iOdd ^ ((k * (M - k)) & 1)) : PETSC_FALSE; 517 for (j = 0; j < Nk; j++) { 518 PetscBool jOdd; 519 520 PetscCall(PetscDTEnumSplit(N, k, j, subsetv, &jOdd)); 521 jidx = negative ? Nk - 1 - j : j; 522 jOdd = negative ? (PetscBool)(iOdd ^ jOdd ^ ((k * (N - k)) & 1)) : PETSC_FALSE; 523 for (p = 0; p < Nf; p++) { 524 PetscReal prod; 525 PetscBool isOdd; 526 527 PetscCall(PetscDTEnumPerm(k, p, perm, &isOdd)); 528 isOdd = (PetscBool)(isOdd ^ jOdd); 529 prod = isOdd ? -1. : 1.; 530 for (l = 0; l < k; l++) prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 531 Lstar[jidx * Mk + iidx] += prod; 532 } 533 } 534 } 535 PetscCall(PetscFree5(subsetw, subsetv, perm, Lw, Lwv)); 536 } 537 PetscFunctionReturn(0); 538 } 539 540 /*@ 541 PetscDTAltVInterior - Compute the interior product of a k-form with a vector 542 543 Input Parameters: 544 + N - the dimension of the vector space, N >= 0 545 . k - the degree k of the k-form w, 0 <= k <= N 546 . w - a k-form, size [N choose k] 547 - v - an N dimensional vector 548 549 Output Parameter: 550 . 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}). 551 552 Level: intermediate 553 554 .seealso: `PetscDTAltV`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 555 @*/ 556 PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv) 557 { 558 PetscInt i, Nk, Nkm; 559 560 PetscFunctionBegin; 561 PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 562 PetscCall(PetscDTBinomialInt(N, k, &Nk)); 563 PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 564 if (N <= 3) { 565 if (k == 1) { 566 PetscReal sum = 0.; 567 568 for (i = 0; i < N; i++) sum += w[i] * v[i]; 569 wIntv[0] = sum; 570 } else if (k == N) { 571 PetscReal mult[3] = {1., -1., 1.}; 572 573 for (i = 0; i < N; i++) wIntv[N - 1 - i] = w[0] * v[i] * mult[i]; 574 } else { 575 wIntv[0] = -w[0] * v[1] - w[1] * v[2]; 576 wIntv[1] = w[0] * v[0] - w[2] * v[2]; 577 wIntv[2] = w[1] * v[0] + w[2] * v[1]; 578 } 579 } else { 580 PetscInt *subset, *work; 581 582 PetscCall(PetscMalloc2(k, &subset, k, &work)); 583 for (i = 0; i < Nkm; i++) wIntv[i] = 0.; 584 for (i = 0; i < Nk; i++) { 585 PetscInt j, l, m; 586 587 PetscCall(PetscDTEnumSubset(N, k, i, subset)); 588 for (j = 0; j < k; j++) { 589 PetscInt idx; 590 PetscBool flip = (PetscBool)(j & 1); 591 592 for (l = 0, m = 0; l < k; l++) { 593 if (l != j) work[m++] = subset[l]; 594 } 595 PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 596 wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]); 597 } 598 } 599 PetscCall(PetscFree2(subset, work)); 600 } 601 PetscFunctionReturn(0); 602 } 603 604 /*@ 605 PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on a k-form by the interior product with a vector 606 607 Input Parameters: 608 + N - the dimension of the vector space, N >= 0 609 . k - the degree k of the k-forms on which intvMat acts, 0 <= k <= N 610 - v - an N dimensional vector 611 612 Output Parameter: 613 . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v) 614 615 Level: intermediate 616 617 .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorPattern()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 618 @*/ 619 PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat) 620 { 621 PetscInt i, Nk, Nkm; 622 623 PetscFunctionBegin; 624 PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 625 PetscCall(PetscDTBinomialInt(N, k, &Nk)); 626 PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 627 if (N <= 3) { 628 if (k == 1) { 629 for (i = 0; i < N; i++) intvMat[i] = v[i]; 630 } else if (k == N) { 631 PetscReal mult[3] = {1., -1., 1.}; 632 633 for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i]; 634 } else { 635 intvMat[0] = -v[1]; 636 intvMat[1] = -v[2]; 637 intvMat[2] = 0.; 638 intvMat[3] = v[0]; 639 intvMat[4] = 0.; 640 intvMat[5] = -v[2]; 641 intvMat[6] = 0.; 642 intvMat[7] = v[0]; 643 intvMat[8] = v[1]; 644 } 645 } else { 646 PetscInt *subset, *work; 647 648 PetscCall(PetscMalloc2(k, &subset, k, &work)); 649 for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.; 650 for (i = 0; i < Nk; i++) { 651 PetscInt j, l, m; 652 653 PetscCall(PetscDTEnumSubset(N, k, i, subset)); 654 for (j = 0; j < k; j++) { 655 PetscInt idx; 656 PetscBool flip = (PetscBool)(j & 1); 657 658 for (l = 0, m = 0; l < k; l++) { 659 if (l != j) work[m++] = subset[l]; 660 } 661 PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 662 intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]]; 663 } 664 } 665 PetscCall(PetscFree2(subset, work)); 666 } 667 PetscFunctionReturn(0); 668 } 669 670 /*@ 671 PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix() 672 673 Input Parameters: 674 + N - the dimension of the vector space, N >= 0 675 - k - the degree of the k-forms on which intvMat from PetscDTAltVInteriorMatrix() acts, 0 <= k <= N. 676 677 Output Parameter: 678 . indices - The interior product matrix intvMat has size [(N choose (k-1)) x (N choose k)] and has (N choose k) * k 679 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 680 coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1) 681 682 Level: intermediate 683 684 Note: this function is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential 685 686 .seealso: `PetscDTAltV`, `PetscDTAltVInterior()`, `PetscDTAltVInteriorMatrix()`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 687 @*/ 688 PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3]) 689 { 690 PetscInt i, Nk, Nkm; 691 692 PetscFunctionBegin; 693 PetscCheck(k > 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 694 PetscCall(PetscDTBinomialInt(N, k, &Nk)); 695 PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 696 if (N <= 3) { 697 if (k == 1) { 698 for (i = 0; i < N; i++) { 699 indices[i][0] = 0; 700 indices[i][1] = i; 701 indices[i][2] = i; 702 } 703 } else if (k == N) { 704 PetscInt val[3] = {0, -2, 2}; 705 706 for (i = 0; i < N; i++) { 707 indices[i][0] = N - 1 - i; 708 indices[i][1] = 0; 709 indices[i][2] = val[i]; 710 } 711 } else { 712 indices[0][0] = 0; 713 indices[0][1] = 0; 714 indices[0][2] = -(1 + 1); 715 indices[1][0] = 0; 716 indices[1][1] = 1; 717 indices[1][2] = -(2 + 1); 718 indices[2][0] = 1; 719 indices[2][1] = 0; 720 indices[2][2] = 0; 721 indices[3][0] = 1; 722 indices[3][1] = 2; 723 indices[3][2] = -(2 + 1); 724 indices[4][0] = 2; 725 indices[4][1] = 1; 726 indices[4][2] = 0; 727 indices[5][0] = 2; 728 indices[5][1] = 2; 729 indices[5][2] = 1; 730 } 731 } else { 732 PetscInt *subset, *work; 733 734 PetscCall(PetscMalloc2(k, &subset, k, &work)); 735 for (i = 0; i < Nk; i++) { 736 PetscInt j, l, m; 737 738 PetscCall(PetscDTEnumSubset(N, k, i, subset)); 739 for (j = 0; j < k; j++) { 740 PetscInt idx; 741 PetscBool flip = (PetscBool)(j & 1); 742 743 for (l = 0, m = 0; l < k; l++) { 744 if (l != j) work[m++] = subset[l]; 745 } 746 PetscCall(PetscDTSubsetIndex(N, k - 1, work, &idx)); 747 indices[i * k + j][0] = idx; 748 indices[i * k + j][1] = i; 749 indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j]; 750 } 751 } 752 PetscCall(PetscFree2(subset, work)); 753 } 754 PetscFunctionReturn(0); 755 } 756 757 /*@ 758 PetscDTAltVStar - Apply a power of the Hodge star operator, which maps k-forms to (N-k) forms, to a k-form 759 760 Input Parameters: 761 + N - the dimension of the vector space, N >= 0 762 . k - the degree k of the k-form w, 0 <= k <= N 763 . 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. 764 - w - a k-form, size [N choose k] 765 766 Output Parameter: 767 . 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. 768 769 Level: intermediate 770 771 .seealso: `PetscDTAltV`, `PetscDTAltVPullback()`, `PetscDTAltVPullbackMatrix()` 772 @*/ 773 PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw) 774 { 775 PetscInt Nk, i; 776 777 PetscFunctionBegin; 778 PetscCheck(k >= 0 && k <= N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 779 PetscCall(PetscDTBinomialInt(N, k, &Nk)); 780 pow = pow % 4; 781 pow = (pow + 4) % 4; /* make non-negative */ 782 /* pow is now 0, 1, 2, 3 */ 783 if (N <= 3) { 784 if (pow & 1) { 785 PetscReal mult[3] = {1., -1., 1.}; 786 787 for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i]; 788 } else { 789 for (i = 0; i < Nk; i++) starw[i] = w[i]; 790 } 791 if (pow > 1 && ((k * (N - k)) & 1)) { 792 for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 793 } 794 } else { 795 PetscInt *subset; 796 797 PetscCall(PetscMalloc1(N, &subset)); 798 if (pow % 2) { 799 PetscInt l = (pow == 1) ? k : N - k; 800 for (i = 0; i < Nk; i++) { 801 PetscBool sOdd; 802 PetscInt j, idx; 803 804 PetscCall(PetscDTEnumSplit(N, l, i, subset, &sOdd)); 805 PetscCall(PetscDTSubsetIndex(N, l, subset, &idx)); 806 PetscCall(PetscDTSubsetIndex(N, N - l, &subset[l], &j)); 807 starw[j] = sOdd ? -w[idx] : w[idx]; 808 } 809 } else { 810 for (i = 0; i < Nk; i++) starw[i] = w[i]; 811 } 812 /* star^2 = -1^(k * (N - k)) */ 813 if (pow > 1 && (k * (N - k)) % 2) { 814 for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 815 } 816 PetscCall(PetscFree(subset)); 817 } 818 PetscFunctionReturn(0); 819 } 820