11a989b97SToby Isaac #include <petsc/private/petscimpl.h> 21a989b97SToby Isaac #include <petsc/private/dtimpl.h> 31a989b97SToby Isaac 4*fad4db65SToby Isaac /*@ 5*fad4db65SToby Isaac PetscDTAltVApply - Apply a k-form to a set of k N-dimensional vectors 6*fad4db65SToby Isaac 7*fad4db65SToby Isaac Input Arguments: 8*fad4db65SToby Isaac + N - the dimension of the space 9*fad4db65SToby Isaac . k - the index of the alternating form 10*fad4db65SToby Isaac . w - the alternating form 11*fad4db65SToby Isaac - v - the set of k vectors of size N, size [k x N], row-major 12*fad4db65SToby Isaac 13*fad4db65SToby Isaac Output Arguments: 14*fad4db65SToby Isaac . wv - w[v[0],...,v[k-1]] 15*fad4db65SToby Isaac 16*fad4db65SToby Isaac Level: intermediate 17*fad4db65SToby Isaac 18*fad4db65SToby Isaac .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 19*fad4db65SToby Isaac @*/ 201a989b97SToby Isaac PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv) 211a989b97SToby Isaac { 221a989b97SToby Isaac PetscErrorCode ierr; 231a989b97SToby Isaac 241a989b97SToby Isaac PetscFunctionBegin; 251a989b97SToby Isaac if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 261a989b97SToby Isaac if (k < 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 271a989b97SToby Isaac if (N <= 3) { 281a989b97SToby Isaac if (!k) { 291a989b97SToby Isaac *wv = w[0]; 301a989b97SToby Isaac } else { 311a989b97SToby Isaac if (N == 1) {*wv = w[0] * v[0];} 321a989b97SToby Isaac else if (N == 2) { 331a989b97SToby Isaac if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1];} 341a989b97SToby Isaac else {*wv = w[0] * (v[0] * v[3] - v[1] * v[2]);} 351a989b97SToby Isaac } else { 361a989b97SToby Isaac if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];} 371a989b97SToby Isaac else if (k == 2) { 381a989b97SToby Isaac *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) + 391a989b97SToby Isaac w[1] * (v[0] * v[5] - v[2] * v[3]) + 401a989b97SToby Isaac w[2] * (v[1] * v[5] - v[2] * v[4]); 411a989b97SToby Isaac } else { 421a989b97SToby Isaac *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) + 431a989b97SToby Isaac v[1] * (v[5] * v[6] - v[3] * v[8]) + 441a989b97SToby Isaac v[2] * (v[3] * v[7] - v[4] * v[6])); 451a989b97SToby Isaac } 461a989b97SToby Isaac } 471a989b97SToby Isaac } 481a989b97SToby Isaac } else { 491a989b97SToby Isaac PetscInt Nk, Nf; 50*fad4db65SToby Isaac PetscInt *subset, *perm; 511a989b97SToby Isaac PetscInt i, j, l; 521a989b97SToby Isaac PetscReal sum = 0.; 531a989b97SToby Isaac 54*fad4db65SToby Isaac ierr = PetscDTFactorialInt(k, &Nf);CHKERRQ(ierr); 55*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 56*fad4db65SToby Isaac ierr = PetscMalloc2(k, &subset, k, &perm);CHKERRQ(ierr); 571a989b97SToby Isaac for (i = 0; i < Nk; i++) { 581a989b97SToby Isaac PetscReal subsum = 0.; 591a989b97SToby Isaac 601a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 611a989b97SToby Isaac for (j = 0; j < Nf; j++) { 621a989b97SToby Isaac PetscBool permOdd; 631a989b97SToby Isaac PetscReal prod; 641a989b97SToby Isaac 65*fad4db65SToby Isaac ierr = PetscDTEnumPerm(k, j, perm, &permOdd);CHKERRQ(ierr); 661a989b97SToby Isaac prod = permOdd ? -1. : 1.; 671a989b97SToby Isaac for (l = 0; l < k; l++) { 681a989b97SToby Isaac prod *= v[perm[l] * N + subset[l]]; 691a989b97SToby Isaac } 701a989b97SToby Isaac subsum += prod; 711a989b97SToby Isaac } 721a989b97SToby Isaac sum += w[i] * subsum; 731a989b97SToby Isaac } 74*fad4db65SToby Isaac ierr = PetscFree2(subset, perm);CHKERRQ(ierr); 751a989b97SToby Isaac *wv = sum; 761a989b97SToby Isaac } 771a989b97SToby Isaac PetscFunctionReturn(0); 781a989b97SToby Isaac } 791a989b97SToby Isaac 80*fad4db65SToby Isaac /*@ 81*fad4db65SToby Isaac PetscDTAltVWedge - Compute the wedge product of two alternating forms 82*fad4db65SToby Isaac 83*fad4db65SToby Isaac Input Arguments: 84*fad4db65SToby Isaac + N - the dimension of the space 85*fad4db65SToby Isaac . j - the index of the form a 86*fad4db65SToby Isaac . k - the index of the form b 87*fad4db65SToby Isaac . a - the j-form 88*fad4db65SToby Isaac - b - the k-form 89*fad4db65SToby Isaac 90*fad4db65SToby Isaac Output Arguments: 91*fad4db65SToby Isaac . awedgeb - a wedge b, a (j+k)-form 92*fad4db65SToby Isaac 93*fad4db65SToby Isaac Level: intermediate 94*fad4db65SToby Isaac 95*fad4db65SToby Isaac .seealso: PetscDTAltVWedgeMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 96*fad4db65SToby Isaac @*/ 971a989b97SToby Isaac PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) 981a989b97SToby Isaac { 991a989b97SToby Isaac PetscInt i; 1001a989b97SToby Isaac PetscErrorCode ierr; 1011a989b97SToby Isaac 1021a989b97SToby Isaac PetscFunctionBegin; 1031a989b97SToby Isaac if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 1041a989b97SToby Isaac if (j < 0 || k < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 1051a989b97SToby Isaac if (j + k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 1061a989b97SToby Isaac if (N <= 3) { 1071a989b97SToby Isaac PetscInt Njk; 1081a989b97SToby Isaac 109*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 1101a989b97SToby Isaac if (!j) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[0] * b[i];}} 1111a989b97SToby Isaac else if (!k) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[i] * b[0];}} 1121a989b97SToby Isaac else { 1131a989b97SToby Isaac if (N == 2) {awedgeb[0] = a[0] * b[1] - a[1] * b[0];} 1141a989b97SToby Isaac else { 1151a989b97SToby Isaac if (j+k == 2) { 1161a989b97SToby Isaac awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 1171a989b97SToby Isaac awedgeb[1] = a[0] * b[2] - a[2] * b[0]; 1181a989b97SToby Isaac awedgeb[2] = a[1] * b[2] - a[2] * b[1]; 1191a989b97SToby Isaac } else { 1201a989b97SToby Isaac awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0]; 1211a989b97SToby Isaac } 1221a989b97SToby Isaac } 1231a989b97SToby Isaac } 1241a989b97SToby Isaac } else { 1251a989b97SToby Isaac PetscInt Njk; 1261a989b97SToby Isaac PetscInt JKj; 1271a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 1281a989b97SToby Isaac PetscInt i; 1291a989b97SToby Isaac 130*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 131*fad4db65SToby Isaac ierr = PetscDTBinomialInt(j+k, j, &JKj);CHKERRQ(ierr); 1321a989b97SToby Isaac ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr); 1331a989b97SToby Isaac for (i = 0; i < Njk; i++) { 1341a989b97SToby Isaac PetscReal sum = 0.; 1351a989b97SToby Isaac PetscInt l; 1361a989b97SToby Isaac 1371a989b97SToby Isaac ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr); 1381a989b97SToby Isaac for (l = 0; l < JKj; l++) { 1391a989b97SToby Isaac PetscBool jkOdd; 1401a989b97SToby Isaac PetscInt m, jInd, kInd; 1411a989b97SToby Isaac 1421a989b97SToby Isaac ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr); 1431a989b97SToby Isaac for (m = 0; m < j; m++) { 1441a989b97SToby Isaac subsetj[m] = subset[subsetjk[m]]; 1451a989b97SToby Isaac } 1461a989b97SToby Isaac for (m = 0; m < k; m++) { 1471a989b97SToby Isaac subsetk[m] = subset[subsetjk[j+m]]; 1481a989b97SToby Isaac } 1491a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr); 1501a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr); 1511a989b97SToby Isaac sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]); 1521a989b97SToby Isaac } 1531a989b97SToby Isaac awedgeb[i] = sum; 1541a989b97SToby Isaac } 1551a989b97SToby Isaac ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr); 1561a989b97SToby Isaac } 1571a989b97SToby Isaac PetscFunctionReturn(0); 1581a989b97SToby Isaac } 1591a989b97SToby Isaac 160*fad4db65SToby Isaac /*@ 161*fad4db65SToby Isaac PetscDTAltVWedgeMatrix - Compute the matrix defined by the wedge product with an alternating form 162*fad4db65SToby Isaac 163*fad4db65SToby Isaac Input Arguments: 164*fad4db65SToby Isaac + N - the dimension of the space 165*fad4db65SToby Isaac . j - the index of the form a 166*fad4db65SToby Isaac . k - the index of the form that (a wedge) will be applied to 167*fad4db65SToby Isaac - a - the j-form 168*fad4db65SToby Isaac 169*fad4db65SToby Isaac Output Arguments: 170*fad4db65SToby 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 171*fad4db65SToby Isaac 172*fad4db65SToby Isaac Level: intermediate 173*fad4db65SToby Isaac 174*fad4db65SToby Isaac .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 175*fad4db65SToby Isaac @*/ 1761a989b97SToby Isaac PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat) 1771a989b97SToby Isaac { 1781a989b97SToby Isaac PetscInt i; 1791a989b97SToby Isaac PetscErrorCode ierr; 1801a989b97SToby Isaac 1811a989b97SToby Isaac PetscFunctionBegin; 1821a989b97SToby Isaac if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 1831a989b97SToby Isaac if (j < 0 || k < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 1841a989b97SToby Isaac if (j + k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 1851a989b97SToby Isaac if (N <= 3) { 1861a989b97SToby Isaac PetscInt Njk; 1871a989b97SToby Isaac 188*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 1891a989b97SToby Isaac if (!j) { 1901a989b97SToby Isaac for (i = 0; i < Njk * Njk; i++) {awedgeMat[i] = 0.;} 1911a989b97SToby Isaac for (i = 0; i < Njk; i++) {awedgeMat[i * (Njk + 1)] = a[0];} 1921a989b97SToby Isaac } else if (!k) { 1931a989b97SToby Isaac for (i = 0; i < Njk; i++) {awedgeMat[i] = a[i];} 1941a989b97SToby Isaac } else { 1951a989b97SToby Isaac if (N == 2) { 1961a989b97SToby Isaac awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; 1971a989b97SToby Isaac } else { 1981a989b97SToby Isaac if (j+k == 2) { 1991a989b97SToby Isaac awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; awedgeMat[2] = 0.; 2001a989b97SToby Isaac awedgeMat[3] = -a[2]; awedgeMat[4] = 0.; awedgeMat[5] = a[0]; 2011a989b97SToby Isaac awedgeMat[6] = 0.; awedgeMat[7] = -a[2]; awedgeMat[8] = a[1]; 2021a989b97SToby Isaac } else { 2031a989b97SToby Isaac awedgeMat[0] = a[2]; awedgeMat[1] = -a[1]; awedgeMat[2] = a[0]; 2041a989b97SToby Isaac } 2051a989b97SToby Isaac } 2061a989b97SToby Isaac } 2071a989b97SToby Isaac } else { 2081a989b97SToby Isaac PetscInt Njk; 2091a989b97SToby Isaac PetscInt Nk; 2101a989b97SToby Isaac PetscInt JKj, i; 2111a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 2121a989b97SToby Isaac 213*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 214*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr); 215*fad4db65SToby Isaac ierr = PetscDTBinomialInt(j+k, j, &JKj);CHKERRQ(ierr); 2161a989b97SToby Isaac ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr); 2171a989b97SToby Isaac for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.; 2181a989b97SToby Isaac for (i = 0; i < Njk; i++) { 2191a989b97SToby Isaac PetscInt l; 2201a989b97SToby Isaac 2211a989b97SToby Isaac ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr); 2221a989b97SToby Isaac for (l = 0; l < JKj; l++) { 2231a989b97SToby Isaac PetscBool jkOdd; 2241a989b97SToby Isaac PetscInt m, jInd, kInd; 2251a989b97SToby Isaac 2261a989b97SToby Isaac ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr); 2271a989b97SToby Isaac for (m = 0; m < j; m++) { 2281a989b97SToby Isaac subsetj[m] = subset[subsetjk[m]]; 2291a989b97SToby Isaac } 2301a989b97SToby Isaac for (m = 0; m < k; m++) { 2311a989b97SToby Isaac subsetk[m] = subset[subsetjk[j+m]]; 2321a989b97SToby Isaac } 2331a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr); 2341a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr); 2351a989b97SToby Isaac awedgeMat[i * Nk + kInd] += jkOdd ? - a[jInd] : a[jInd]; 2361a989b97SToby Isaac } 2371a989b97SToby Isaac } 2381a989b97SToby Isaac ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr); 2391a989b97SToby Isaac } 2401a989b97SToby Isaac PetscFunctionReturn(0); 2411a989b97SToby Isaac } 2421a989b97SToby Isaac 243*fad4db65SToby Isaac /*@ 244*fad4db65SToby Isaac PetscDTAltVPullback - Compute the pullback of an alternating form under a linear transformation 245*fad4db65SToby Isaac 246*fad4db65SToby Isaac Input Arguments: 247*fad4db65SToby Isaac + N - the dimension of the origin space 248*fad4db65SToby Isaac . M - the dimension of the image space 249*fad4db65SToby Isaac . L - the linear transformation, an [M x N] matrix in row-major format 250*fad4db65SToby Isaac . k - the index of the form. A negative form degree indicates that Pullback should be conjugated by the Hodge star operator (see note) 251*fad4db65SToby Isaac - w - the k-form in the image space 252*fad4db65SToby Isaac 253*fad4db65SToby Isaac Output Arguments: 254*fad4db65SToby Isaac . Lstarw - the pullback of w to a k-form in the origin space 255*fad4db65SToby Isaac 256*fad4db65SToby Isaac Level: intermediate 257*fad4db65SToby Isaac 258*fad4db65SToby Isaac Note: negative form degrees accomodate, 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, 259*fad4db65SToby 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 1-form to 2-form, pullback as a 2-form, 260*fad4db65SToby Isaac then the inverse Hodge star transformation. 261*fad4db65SToby Isaac 262*fad4db65SToby Isaac .seealso: PetscDTAltVPullbackMatrix(), PetscDTAltVStar() 263*fad4db65SToby Isaac @*/ 2641a989b97SToby Isaac PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw) 2651a989b97SToby Isaac { 2661a989b97SToby Isaac PetscInt i, j, Nk, Mk; 2671a989b97SToby Isaac PetscErrorCode ierr; 2681a989b97SToby Isaac 2691a989b97SToby Isaac PetscFunctionBegin; 2701a989b97SToby Isaac if (N < 0 || M < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 2711a989b97SToby Isaac if (PetscAbsInt(k) > N || PetscAbsInt(k) > M) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 2721a989b97SToby Isaac if (N <= 3 && M <= 3) { 2731a989b97SToby Isaac 274*fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 275*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 2761a989b97SToby Isaac if (!k) { 2771a989b97SToby Isaac Lstarw[0] = w[0]; 2781a989b97SToby Isaac } else if (k == 1) { 2791a989b97SToby Isaac for (i = 0; i < Nk; i++) { 2801a989b97SToby Isaac PetscReal sum = 0.; 2811a989b97SToby Isaac 2821a989b97SToby Isaac for (j = 0; j < Mk; j++) {sum += L[j * Nk + i] * w[j];} 2831a989b97SToby Isaac Lstarw[i] = sum; 2841a989b97SToby Isaac } 2851a989b97SToby Isaac } else if (k == -1) { 2861a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 2871a989b97SToby Isaac 2881a989b97SToby Isaac for (i = 0; i < Nk; i++) { 2891a989b97SToby Isaac PetscReal sum = 0.; 2901a989b97SToby Isaac 2911a989b97SToby Isaac for (j = 0; j < Mk; j++) { 2921a989b97SToby Isaac sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j]; 2931a989b97SToby Isaac } 2941a989b97SToby Isaac Lstarw[i] = mult[i] * sum; 2951a989b97SToby Isaac } 2961a989b97SToby Isaac } else if (k == 2) { 2971a989b97SToby Isaac PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}}; 2981a989b97SToby Isaac 2991a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3001a989b97SToby Isaac PetscReal sum = 0.; 3011a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3021a989b97SToby Isaac sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] - 3031a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j]; 3041a989b97SToby Isaac } 3051a989b97SToby Isaac Lstarw[i] = sum; 3061a989b97SToby Isaac } 3071a989b97SToby Isaac } else if (k == -2) { 3081a989b97SToby Isaac PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}}; 3091a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 3101a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 3111a989b97SToby Isaac 3121a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3131a989b97SToby Isaac PetscReal sum = 0.; 3141a989b97SToby Isaac 3151a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3161a989b97SToby Isaac sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] * 3171a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][1]] - 3181a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][0]] * 3191a989b97SToby Isaac L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j]; 3201a989b97SToby Isaac 3211a989b97SToby Isaac } 3221a989b97SToby Isaac Lstarw[i] = sum; 3231a989b97SToby Isaac } 3241a989b97SToby Isaac } else { 3251a989b97SToby Isaac PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + 3261a989b97SToby Isaac L[1] * (L[5] * L[6] - L[3] * L[8]) + 3271a989b97SToby Isaac L[2] * (L[3] * L[7] - L[4] * L[6]); 3281a989b97SToby Isaac 3291a989b97SToby Isaac for (i = 0; i < Nk; i++) {Lstarw[i] = detL * w[i];} 3301a989b97SToby Isaac } 3311a989b97SToby Isaac } else { 3321a989b97SToby Isaac PetscInt Nf, l, p; 3331a989b97SToby Isaac PetscReal *Lw, *Lwv; 3341a989b97SToby Isaac PetscInt *subsetw, *subsetv; 335*fad4db65SToby Isaac PetscInt *perm; 3361a989b97SToby Isaac PetscReal *walloc = NULL; 3371a989b97SToby Isaac const PetscReal *ww = NULL; 3381a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 3391a989b97SToby Isaac 340*fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 341*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 342*fad4db65SToby Isaac ierr = PetscDTFactorialInt(PetscAbsInt(k), &Nf);CHKERRQ(ierr); 3431a989b97SToby Isaac if (k < 0) { 3441a989b97SToby Isaac negative = PETSC_TRUE; 3451a989b97SToby Isaac k = -k; 3461a989b97SToby Isaac ierr = PetscMalloc1(Mk, &walloc);CHKERRQ(ierr); 3471a989b97SToby Isaac ierr = PetscDTAltVStar(M, M - k, 1, w, walloc);CHKERRQ(ierr); 3481a989b97SToby Isaac ww = walloc; 3491a989b97SToby Isaac } else { 3501a989b97SToby Isaac ww = w; 3511a989b97SToby Isaac } 352*fad4db65SToby Isaac ierr = PetscMalloc5(k, &subsetw, k, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr); 3531a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = 0.; 3541a989b97SToby Isaac for (i = 0; i < Mk; i++) { 3551a989b97SToby Isaac ierr = PetscDTEnumSubset(M, k, i, subsetw);CHKERRQ(ierr); 3561a989b97SToby Isaac for (j = 0; j < Nk; j++) { 3571a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, j, subsetv);CHKERRQ(ierr); 3581a989b97SToby Isaac for (p = 0; p < Nf; p++) { 3591a989b97SToby Isaac PetscReal prod; 3601a989b97SToby Isaac PetscBool isOdd; 3611a989b97SToby Isaac 362*fad4db65SToby Isaac ierr = PetscDTEnumPerm(k, p, perm, &isOdd);CHKERRQ(ierr); 3631a989b97SToby Isaac prod = isOdd ? -ww[i] : ww[i]; 3641a989b97SToby Isaac for (l = 0; l < k; l++) { 3651a989b97SToby Isaac prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 3661a989b97SToby Isaac } 3671a989b97SToby Isaac Lstarw[j] += prod; 3681a989b97SToby Isaac } 3691a989b97SToby Isaac } 3701a989b97SToby Isaac } 3711a989b97SToby Isaac if (negative) { 3721a989b97SToby Isaac PetscReal *sLsw; 3731a989b97SToby Isaac 3741a989b97SToby Isaac ierr = PetscMalloc1(Nk, &sLsw);CHKERRQ(ierr); 3751a989b97SToby Isaac ierr = PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw);CHKERRQ(ierr); 3761a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i]; 3771a989b97SToby Isaac ierr = PetscFree(sLsw);CHKERRQ(ierr); 3781a989b97SToby Isaac } 379*fad4db65SToby Isaac ierr = PetscFree5(subsetw, subsetv, perm, Lw, Lwv);CHKERRQ(ierr); 3801a989b97SToby Isaac ierr = PetscFree(walloc);CHKERRQ(ierr); 3811a989b97SToby Isaac } 3821a989b97SToby Isaac PetscFunctionReturn(0); 3831a989b97SToby Isaac } 3841a989b97SToby Isaac 385*fad4db65SToby Isaac /*@ 386*fad4db65SToby Isaac PetscDTAltVPullbackMatrix - Compute the pullback matrix for k-forms under a linear transformation 387*fad4db65SToby Isaac 388*fad4db65SToby Isaac Input Arguments: 389*fad4db65SToby Isaac + N - the dimension of the origin space 390*fad4db65SToby Isaac . M - the dimension of the image space 391*fad4db65SToby Isaac . L - the linear transformation, an [M x N] matrix in row-major format 392*fad4db65SToby Isaac - k - the index of the alternating forms. A negative form degree indicates that the pullback should be conjugated by the Hodge star operator (see note in PetscDTAltvPullback()) 393*fad4db65SToby Isaac 394*fad4db65SToby Isaac Output Arguments: 395*fad4db65SToby Isaac . Lstar - the pullback matrix, an [(N choose k) x (M choose k)] matrix in row-major format such that Lstar * w = L^* w 396*fad4db65SToby Isaac 397*fad4db65SToby Isaac Level: intermediate 398*fad4db65SToby Isaac 399*fad4db65SToby Isaac .seealso: PetscDTAltVPullback(), PetscDTAltVStar() 400*fad4db65SToby Isaac @*/ 4011a989b97SToby Isaac PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) 4021a989b97SToby Isaac { 4031a989b97SToby Isaac PetscInt Nk, Mk, Nf, i, j, l, p; 4041a989b97SToby Isaac PetscReal *Lw, *Lwv; 4051a989b97SToby Isaac PetscInt *subsetw, *subsetv; 406*fad4db65SToby Isaac PetscInt *perm; 4071a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 4081a989b97SToby Isaac PetscErrorCode ierr; 4091a989b97SToby Isaac 4101a989b97SToby Isaac PetscFunctionBegin; 4111a989b97SToby Isaac if (N < 0 || M < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 4121a989b97SToby Isaac if (PetscAbsInt(k) > N || PetscAbsInt(k) > M) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 4131a989b97SToby Isaac if (N <= 3 && M <= 3) { 4141a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 4151a989b97SToby Isaac 416*fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 417*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 4181a989b97SToby Isaac if (!k) { 4191a989b97SToby Isaac Lstar[0] = 1.; 4201a989b97SToby Isaac } else if (k == 1) { 4211a989b97SToby Isaac for (i = 0; i < Nk; i++) {for (j = 0; j < Mk; j++) {Lstar[i * Mk + j] = L[j * Nk + i];}} 4221a989b97SToby Isaac } else if (k == -1) { 4231a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4241a989b97SToby Isaac for (j = 0; j < Mk; j++) { 4251a989b97SToby Isaac Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j]; 4261a989b97SToby Isaac } 4271a989b97SToby Isaac } 4281a989b97SToby Isaac } else if (k == 2) { 4291a989b97SToby Isaac PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}}; 4301a989b97SToby Isaac 4311a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4321a989b97SToby Isaac for (j = 0; j < Mk; j++) { 4331a989b97SToby Isaac Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] * 4341a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][1]] - 4351a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][0]] * 4361a989b97SToby Isaac L[pairs[j][0] * N + pairs[i][1]]; 4371a989b97SToby Isaac } 4381a989b97SToby Isaac } 4391a989b97SToby Isaac } else if (k == -2) { 4401a989b97SToby Isaac PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}}; 4411a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 4421a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 4431a989b97SToby Isaac 4441a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4451a989b97SToby Isaac for (j = 0; j < Mk; j++) { 4461a989b97SToby Isaac Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] * 4471a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][1]] - 4481a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][0]] * 4491a989b97SToby Isaac L[pairs[offj + j][0] * N + pairs[offi + i][1]]; 4501a989b97SToby Isaac } 4511a989b97SToby Isaac } 4521a989b97SToby Isaac } else { 4531a989b97SToby Isaac PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + 4541a989b97SToby Isaac L[1] * (L[5] * L[6] - L[3] * L[8]) + 4551a989b97SToby Isaac L[2] * (L[3] * L[7] - L[4] * L[6]); 4561a989b97SToby Isaac 4571a989b97SToby Isaac for (i = 0; i < Nk; i++) {Lstar[i] = detL;} 4581a989b97SToby Isaac } 4591a989b97SToby Isaac } else { 4601a989b97SToby Isaac if (k < 0) { 4611a989b97SToby Isaac negative = PETSC_TRUE; 4621a989b97SToby Isaac k = -k; 4631a989b97SToby Isaac } 464*fad4db65SToby Isaac ierr = PetscDTBinomialInt(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 465*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 466*fad4db65SToby Isaac ierr = PetscDTFactorialInt(PetscAbsInt(k), &Nf);CHKERRQ(ierr); 467*fad4db65SToby Isaac ierr = PetscMalloc5(M, &subsetw, N, &subsetv, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr); 4681a989b97SToby Isaac for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.; 4691a989b97SToby Isaac for (i = 0; i < Mk; i++) { 4701a989b97SToby Isaac PetscBool iOdd; 4711a989b97SToby Isaac PetscInt iidx, jidx; 4721a989b97SToby Isaac 4731a989b97SToby Isaac ierr = PetscDTEnumSplit(M, k, i, subsetw, &iOdd);CHKERRQ(ierr); 4741a989b97SToby Isaac iidx = negative ? Mk - 1 - i : i; 4751a989b97SToby Isaac iOdd = negative ? iOdd ^ ((k * (M-k)) & 1) : PETSC_FALSE; 4761a989b97SToby Isaac for (j = 0; j < Nk; j++) { 4771a989b97SToby Isaac PetscBool jOdd; 4781a989b97SToby Isaac 4791a989b97SToby Isaac ierr = PetscDTEnumSplit(N, k, j, subsetv, &jOdd);CHKERRQ(ierr); 4801a989b97SToby Isaac jidx = negative ? Nk - 1 - j : j; 4811a989b97SToby Isaac jOdd = negative ? iOdd ^ jOdd ^ ((k * (N-k)) & 1) : PETSC_FALSE; 4821a989b97SToby Isaac for (p = 0; p < Nf; p++) { 4831a989b97SToby Isaac PetscReal prod; 4841a989b97SToby Isaac PetscBool isOdd; 4851a989b97SToby Isaac 486*fad4db65SToby Isaac ierr = PetscDTEnumPerm(k, p, perm, &isOdd);CHKERRQ(ierr); 4871a989b97SToby Isaac isOdd ^= jOdd; 4881a989b97SToby Isaac prod = isOdd ? -1. : 1.; 4891a989b97SToby Isaac for (l = 0; l < k; l++) { 4901a989b97SToby Isaac prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 4911a989b97SToby Isaac } 4921a989b97SToby Isaac Lstar[jidx * Mk + iidx] += prod; 4931a989b97SToby Isaac } 4941a989b97SToby Isaac } 4951a989b97SToby Isaac } 496*fad4db65SToby Isaac ierr = PetscFree5(subsetw, subsetv, perm, Lw, Lwv);CHKERRQ(ierr); 4971a989b97SToby Isaac } 4981a989b97SToby Isaac PetscFunctionReturn(0); 4991a989b97SToby Isaac } 5001a989b97SToby Isaac 501*fad4db65SToby Isaac /*@ 502*fad4db65SToby Isaac PetscDTAltVInterior - Compute the interior product of an alternating form with a vector 503*fad4db65SToby Isaac 504*fad4db65SToby Isaac Input Arguments: 505*fad4db65SToby Isaac + N - the dimension of the origin space 506*fad4db65SToby Isaac . k - the index of the alternating forms 507*fad4db65SToby Isaac . w - the k-form 508*fad4db65SToby Isaac - v - the N dimensional vector 509*fad4db65SToby Isaac 510*fad4db65SToby Isaac Output Arguments: 511*fad4db65SToby Isaac . wIntv - the (k-1) form (w int v). 512*fad4db65SToby Isaac 513*fad4db65SToby Isaac Level: intermediate 514*fad4db65SToby Isaac 515*fad4db65SToby Isaac .seealso: PetscDTAltVInteriorMatrix(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 516*fad4db65SToby Isaac @*/ 5171a989b97SToby Isaac PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv) 5181a989b97SToby Isaac { 5191a989b97SToby Isaac PetscInt i, Nk, Nkm; 5201a989b97SToby Isaac PetscErrorCode ierr; 5211a989b97SToby Isaac 5221a989b97SToby Isaac PetscFunctionBegin; 5231a989b97SToby Isaac if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 524*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 525*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr); 5261a989b97SToby Isaac if (N <= 3) { 5271a989b97SToby Isaac if (k == 1) { 5281a989b97SToby Isaac PetscReal sum = 0.; 5291a989b97SToby Isaac 5301a989b97SToby Isaac for (i = 0; i < N; i++) { 5311a989b97SToby Isaac sum += w[i] * v[i]; 5321a989b97SToby Isaac } 5331a989b97SToby Isaac wIntv[0] = sum; 5341a989b97SToby Isaac } else if (k == N) { 5351a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 5361a989b97SToby Isaac 5371a989b97SToby Isaac for (i = 0; i < N; i++) { 5381a989b97SToby Isaac wIntv[N - 1 - i] = w[0] * v[i] * mult[i]; 5391a989b97SToby Isaac } 5401a989b97SToby Isaac } else { 5411a989b97SToby Isaac wIntv[0] = - w[0]*v[1] - w[1]*v[2]; 5421a989b97SToby Isaac wIntv[1] = w[0]*v[0] - w[2]*v[2]; 5431a989b97SToby Isaac wIntv[2] = w[1]*v[0] + w[2]*v[1]; 5441a989b97SToby Isaac } 5451a989b97SToby Isaac } else { 5461a989b97SToby Isaac PetscInt *subset, *work; 5471a989b97SToby Isaac 5481a989b97SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 5491a989b97SToby Isaac for (i = 0; i < Nkm; i++) wIntv[i] = 0.; 5501a989b97SToby Isaac for (i = 0; i < Nk; i++) { 5511a989b97SToby Isaac PetscInt j, l, m; 5521a989b97SToby Isaac 5531a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 5541a989b97SToby Isaac for (j = 0; j < k; j++) { 5551a989b97SToby Isaac PetscInt idx; 5561a989b97SToby Isaac PetscBool flip = (j & 1); 5571a989b97SToby Isaac 5581a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 5591a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 5601a989b97SToby Isaac } 5611a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 5621a989b97SToby Isaac wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]); 5631a989b97SToby Isaac } 5641a989b97SToby Isaac } 5651a989b97SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 5661a989b97SToby Isaac } 5671a989b97SToby Isaac PetscFunctionReturn(0); 5681a989b97SToby Isaac } 5691a989b97SToby Isaac 570*fad4db65SToby Isaac /*@ 571*fad4db65SToby Isaac PetscDTAltVInteriorMatrix - Compute the matrix of the linear transformation induced on an alternating form by the interior product with a vector 572*fad4db65SToby Isaac 573*fad4db65SToby Isaac Input Arguments: 574*fad4db65SToby Isaac + N - the dimension of the origin space 575*fad4db65SToby Isaac . k - the index of the alternating forms 576*fad4db65SToby Isaac - v - the N dimensional vector 577*fad4db65SToby Isaac 578*fad4db65SToby Isaac Output Arguments: 579*fad4db65SToby Isaac . intvMat - an [(N choose (k-1)) x (N choose k)] matrix, row-major: (intvMat) * w = (w int v) 580*fad4db65SToby Isaac 581*fad4db65SToby Isaac Level: intermediate 582*fad4db65SToby Isaac 583*fad4db65SToby Isaac .seealso: PetscDTAltVInterior(), PetscDTAltVInteriorPattern(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 584*fad4db65SToby Isaac @*/ 5851a989b97SToby Isaac PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat) 5861a989b97SToby Isaac { 5871a989b97SToby Isaac PetscInt i, Nk, Nkm; 5881a989b97SToby Isaac PetscErrorCode ierr; 5891a989b97SToby Isaac 5901a989b97SToby Isaac PetscFunctionBegin; 5911a989b97SToby Isaac if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 592*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 593*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr); 5941a989b97SToby Isaac if (N <= 3) { 5951a989b97SToby Isaac if (k == 1) { 5961a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[i] = v[i]; 5971a989b97SToby Isaac } else if (k == N) { 5981a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 5991a989b97SToby Isaac 6001a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i]; 6011a989b97SToby Isaac } else { 6021a989b97SToby Isaac intvMat[0] = -v[1]; intvMat[1] = -v[2]; intvMat[2] = 0.; 6031a989b97SToby Isaac intvMat[3] = v[0]; intvMat[4] = 0.; intvMat[5] = -v[2]; 6041a989b97SToby Isaac intvMat[6] = 0.; intvMat[7] = v[0]; intvMat[8] = v[1]; 6051a989b97SToby Isaac } 6061a989b97SToby Isaac } else { 6071a989b97SToby Isaac PetscInt *subset, *work; 6081a989b97SToby Isaac 6091a989b97SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 6101a989b97SToby Isaac for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.; 6111a989b97SToby Isaac for (i = 0; i < Nk; i++) { 6121a989b97SToby Isaac PetscInt j, l, m; 6131a989b97SToby Isaac 6141a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 6151a989b97SToby Isaac for (j = 0; j < k; j++) { 6161a989b97SToby Isaac PetscInt idx; 6171a989b97SToby Isaac PetscBool flip = (j & 1); 6181a989b97SToby Isaac 6191a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 6201a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 6211a989b97SToby Isaac } 6221a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 6231a989b97SToby Isaac intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]]; 6241a989b97SToby Isaac } 6251a989b97SToby Isaac } 6261a989b97SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 6271a989b97SToby Isaac } 6281a989b97SToby Isaac PetscFunctionReturn(0); 6291a989b97SToby Isaac } 6301a989b97SToby Isaac 631*fad4db65SToby Isaac /*@ 632*fad4db65SToby Isaac PetscDTAltVInteriorPattern - compute the sparsity and sign pattern of the interior product matrix computed in PetscDTAltVInteriorMatrix() 633*fad4db65SToby Isaac 634*fad4db65SToby Isaac Input Arguments: 635*fad4db65SToby Isaac + N - the dimension of the origin space 636*fad4db65SToby Isaac - k - the index of the alternating forms 637*fad4db65SToby Isaac 638*fad4db65SToby Isaac Output Arguments: 639*fad4db65SToby Isaac . indices - The interior product matrix has (N choose k) * k non-zeros. indices[i][0] and indices[i][1] are the row and column of a non-zero, 640*fad4db65SToby Isaac and its value is equal to the vector coordinate v[j] if indices[i][2] = j, or -v[j] if indices[i][2] = -(j+1) 641*fad4db65SToby Isaac 642*fad4db65SToby Isaac Level: intermediate 643*fad4db65SToby Isaac 644*fad4db65SToby Isaac Note: this form is useful when the interior product needs to be computed at multiple locations, as when computing the Koszul differential 645*fad4db65SToby Isaac 646*fad4db65SToby Isaac .seealso: PetscDTAltVInterior(), PetscDTAltVInteriorMatrix(), PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 647*fad4db65SToby Isaac @*/ 648dda711d0SToby Isaac PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3]) 649dda711d0SToby Isaac { 650dda711d0SToby Isaac PetscInt i, Nk, Nkm; 651dda711d0SToby Isaac PetscErrorCode ierr; 652dda711d0SToby Isaac 653dda711d0SToby Isaac PetscFunctionBegin; 654dda711d0SToby Isaac if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 655*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 656*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr); 657dda711d0SToby Isaac if (N <= 3) { 658dda711d0SToby Isaac if (k == 1) { 659dda711d0SToby Isaac for (i = 0; i < N; i++) { 660dda711d0SToby Isaac indices[i][0] = 0; 661dda711d0SToby Isaac indices[i][1] = i; 662dda711d0SToby Isaac indices[i][2] = i; 663dda711d0SToby Isaac } 664dda711d0SToby Isaac } else if (k == N) { 665dda711d0SToby Isaac PetscInt val[3] = {0, -2, 2}; 666dda711d0SToby Isaac 667dda711d0SToby Isaac for (i = 0; i < N; i++) { 668dda711d0SToby Isaac indices[i][0] = N - 1 - i; 669dda711d0SToby Isaac indices[i][1] = 0; 670dda711d0SToby Isaac indices[i][2] = val[i]; 671dda711d0SToby Isaac } 672dda711d0SToby Isaac } else { 673dda711d0SToby Isaac indices[0][0] = 0; indices[0][1] = 0; indices[0][2] = -(1 + 1); 674dda711d0SToby Isaac indices[1][0] = 0; indices[1][1] = 1; indices[1][2] = -(2 + 1); 675dda711d0SToby Isaac indices[2][0] = 1; indices[2][1] = 0; indices[2][2] = 0; 676dda711d0SToby Isaac indices[3][0] = 1; indices[3][1] = 2; indices[3][2] = -(2 + 1); 677dda711d0SToby Isaac indices[4][0] = 2; indices[4][1] = 1; indices[4][2] = 0; 678dda711d0SToby Isaac indices[5][0] = 2; indices[5][1] = 2; indices[5][2] = 1; 679dda711d0SToby Isaac } 680dda711d0SToby Isaac } else { 681dda711d0SToby Isaac PetscInt *subset, *work; 682dda711d0SToby Isaac 683dda711d0SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 684dda711d0SToby Isaac for (i = 0; i < Nk; i++) { 685dda711d0SToby Isaac PetscInt j, l, m; 686dda711d0SToby Isaac 687dda711d0SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 688dda711d0SToby Isaac for (j = 0; j < k; j++) { 689dda711d0SToby Isaac PetscInt idx; 690dda711d0SToby Isaac PetscBool flip = (j & 1); 691dda711d0SToby Isaac 692dda711d0SToby Isaac for (l = 0, m = 0; l < k; l++) { 693dda711d0SToby Isaac if (l != j) work[m++] = subset[l]; 694dda711d0SToby Isaac } 695dda711d0SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 696dda711d0SToby Isaac indices[i * k + j][0] = idx; 697dda711d0SToby Isaac indices[i * k + j][1] = i; 698dda711d0SToby Isaac indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j]; 699dda711d0SToby Isaac } 700dda711d0SToby Isaac } 701dda711d0SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 702dda711d0SToby Isaac } 703dda711d0SToby Isaac PetscFunctionReturn(0); 704dda711d0SToby Isaac } 705dda711d0SToby Isaac 706*fad4db65SToby Isaac /*@ 707*fad4db65SToby Isaac PetscDTAltVStar - Apply a power of the Hodge star transformation to an alternating form 708*fad4db65SToby Isaac 709*fad4db65SToby Isaac Input Arguments: 710*fad4db65SToby Isaac + N - the dimension of the space 711*fad4db65SToby Isaac . k - the index of the alternating form 712*fad4db65SToby Isaac . pow - the number of times to apply the Hodge star operator 713*fad4db65SToby Isaac - w - the alternating form 714*fad4db65SToby Isaac 715*fad4db65SToby Isaac Output Arguments: 716*fad4db65SToby Isaac . starw = (star)^pow w 717*fad4db65SToby Isaac 718*fad4db65SToby Isaac Level: intermediate 719*fad4db65SToby Isaac 720*fad4db65SToby Isaac .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix() 721*fad4db65SToby Isaac @*/ 7221a989b97SToby Isaac PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw) 7231a989b97SToby Isaac { 7241a989b97SToby Isaac PetscInt Nk, i; 7251a989b97SToby Isaac PetscErrorCode ierr; 7261a989b97SToby Isaac 7271a989b97SToby Isaac PetscFunctionBegin; 7281a989b97SToby Isaac if (k < 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 729*fad4db65SToby Isaac ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr); 7301a989b97SToby Isaac pow = pow % 4; 7311a989b97SToby Isaac pow = (pow + 4) % 4; /* make non-negative */ 7321a989b97SToby Isaac /* pow is now 0, 1, 2, 3 */ 7331a989b97SToby Isaac if (N <= 3) { 7341a989b97SToby Isaac if (pow & 1) { 7351a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 7361a989b97SToby Isaac 7371a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i]; 7381a989b97SToby Isaac } else { 7391a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 7401a989b97SToby Isaac } 7411a989b97SToby Isaac if (pow > 1 && ((k * (N - k)) & 1)) { 7421a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 7431a989b97SToby Isaac } 7441a989b97SToby Isaac } else { 7451a989b97SToby Isaac PetscInt *subset; 7461a989b97SToby Isaac 7471a989b97SToby Isaac ierr = PetscMalloc1(N, &subset);CHKERRQ(ierr); 7481a989b97SToby Isaac if (pow % 2) { 7491a989b97SToby Isaac PetscInt l = (pow == 1) ? k : N - k; 7501a989b97SToby Isaac for (i = 0; i < Nk; i++) { 7511a989b97SToby Isaac PetscBool sOdd; 7521a989b97SToby Isaac PetscInt j, idx; 7531a989b97SToby Isaac 7541a989b97SToby Isaac ierr = PetscDTEnumSplit(N, l, i, subset, &sOdd);CHKERRQ(ierr); 7551a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, l, subset, &idx);CHKERRQ(ierr); 7561a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, N-l, &subset[l], &j);CHKERRQ(ierr); 7571a989b97SToby Isaac starw[j] = sOdd ? -w[idx] : w[idx]; 7581a989b97SToby Isaac } 7591a989b97SToby Isaac } else { 7601a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 7611a989b97SToby Isaac } 7621a989b97SToby Isaac /* star^2 = -1^(k * (N - k)) */ 7631a989b97SToby Isaac if (pow > 1 && (k * (N - k)) % 2) { 7641a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 7651a989b97SToby Isaac } 7661a989b97SToby Isaac ierr = PetscFree(subset);CHKERRQ(ierr); 7671a989b97SToby Isaac } 7681a989b97SToby Isaac PetscFunctionReturn(0); 7691a989b97SToby Isaac } 770