11a989b97SToby Isaac #include <petsc/private/petscimpl.h> 21a989b97SToby Isaac #include <petsc/private/dtimpl.h> 31a989b97SToby Isaac 41a989b97SToby Isaac PetscErrorCode PetscDTAltVApply(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wv) 51a989b97SToby Isaac { 61a989b97SToby Isaac PetscErrorCode ierr; 71a989b97SToby Isaac 81a989b97SToby Isaac PetscFunctionBegin; 91a989b97SToby Isaac if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 101a989b97SToby Isaac if (k < 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 111a989b97SToby Isaac if (N <= 3) { 121a989b97SToby Isaac if (!k) { 131a989b97SToby Isaac *wv = w[0]; 141a989b97SToby Isaac } else { 151a989b97SToby Isaac if (N == 1) {*wv = w[0] * v[0];} 161a989b97SToby Isaac else if (N == 2) { 171a989b97SToby Isaac if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1];} 181a989b97SToby Isaac else {*wv = w[0] * (v[0] * v[3] - v[1] * v[2]);} 191a989b97SToby Isaac } else { 201a989b97SToby Isaac if (k == 1) {*wv = w[0] * v[0] + w[1] * v[1] + w[2] * v[2];} 211a989b97SToby Isaac else if (k == 2) { 221a989b97SToby Isaac *wv = w[0] * (v[0] * v[4] - v[1] * v[3]) + 231a989b97SToby Isaac w[1] * (v[0] * v[5] - v[2] * v[3]) + 241a989b97SToby Isaac w[2] * (v[1] * v[5] - v[2] * v[4]); 251a989b97SToby Isaac } else { 261a989b97SToby Isaac *wv = w[0] * (v[0] * (v[4] * v[8] - v[5] * v[7]) + 271a989b97SToby Isaac v[1] * (v[5] * v[6] - v[3] * v[8]) + 281a989b97SToby Isaac v[2] * (v[3] * v[7] - v[4] * v[6])); 291a989b97SToby Isaac } 301a989b97SToby Isaac } 311a989b97SToby Isaac } 321a989b97SToby Isaac } else { 331a989b97SToby Isaac PetscInt Nk, Nf; 341a989b97SToby Isaac PetscInt *subset, *work, *perm; 351a989b97SToby Isaac PetscInt i, j, l; 361a989b97SToby Isaac PetscReal sum = 0.; 371a989b97SToby Isaac 381a989b97SToby Isaac ierr = PetscDTFactorialInt_Internal(k, &Nf);CHKERRQ(ierr); 391a989b97SToby Isaac ierr = PetscDTBinomial(N, k, &Nk);CHKERRQ(ierr); 401a989b97SToby Isaac ierr = PetscMalloc3(k, &subset, k, &work, k, &perm);CHKERRQ(ierr); 411a989b97SToby Isaac for (i = 0; i < Nk; i++) { 421a989b97SToby Isaac PetscReal subsum = 0.; 431a989b97SToby Isaac 441a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 451a989b97SToby Isaac for (j = 0; j < Nf; j++) { 461a989b97SToby Isaac PetscBool permOdd; 471a989b97SToby Isaac PetscReal prod; 481a989b97SToby Isaac 491a989b97SToby Isaac ierr = PetscDTEnumPerm(k, j, work, perm, &permOdd);CHKERRQ(ierr); 501a989b97SToby Isaac prod = permOdd ? -1. : 1.; 511a989b97SToby Isaac for (l = 0; l < k; l++) { 521a989b97SToby Isaac prod *= v[perm[l] * N + subset[l]]; 531a989b97SToby Isaac } 541a989b97SToby Isaac subsum += prod; 551a989b97SToby Isaac } 561a989b97SToby Isaac sum += w[i] * subsum; 571a989b97SToby Isaac } 581a989b97SToby Isaac ierr = PetscFree3(subset, work, perm);CHKERRQ(ierr); 591a989b97SToby Isaac *wv = sum; 601a989b97SToby Isaac } 611a989b97SToby Isaac PetscFunctionReturn(0); 621a989b97SToby Isaac } 631a989b97SToby Isaac 641a989b97SToby Isaac PetscErrorCode PetscDTAltVWedge(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, const PetscReal *b, PetscReal *awedgeb) 651a989b97SToby Isaac { 661a989b97SToby Isaac PetscInt i; 671a989b97SToby Isaac PetscErrorCode ierr; 681a989b97SToby Isaac 691a989b97SToby Isaac PetscFunctionBegin; 701a989b97SToby Isaac if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 711a989b97SToby Isaac if (j < 0 || k < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 721a989b97SToby Isaac if (j + k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 731a989b97SToby Isaac if (N <= 3) { 741a989b97SToby Isaac PetscInt Njk; 751a989b97SToby Isaac 761a989b97SToby Isaac ierr = PetscDTBinomial(N, j+k, &Njk);CHKERRQ(ierr); 771a989b97SToby Isaac if (!j) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[0] * b[i];}} 781a989b97SToby Isaac else if (!k) {for (i = 0; i < Njk; i++) {awedgeb[i] = a[i] * b[0];}} 791a989b97SToby Isaac else { 801a989b97SToby Isaac if (N == 2) {awedgeb[0] = a[0] * b[1] - a[1] * b[0];} 811a989b97SToby Isaac else { 821a989b97SToby Isaac if (j+k == 2) { 831a989b97SToby Isaac awedgeb[0] = a[0] * b[1] - a[1] * b[0]; 841a989b97SToby Isaac awedgeb[1] = a[0] * b[2] - a[2] * b[0]; 851a989b97SToby Isaac awedgeb[2] = a[1] * b[2] - a[2] * b[1]; 861a989b97SToby Isaac } else { 871a989b97SToby Isaac awedgeb[0] = a[0] * b[2] - a[1] * b[1] + a[2] * b[0]; 881a989b97SToby Isaac } 891a989b97SToby Isaac } 901a989b97SToby Isaac } 911a989b97SToby Isaac } else { 921a989b97SToby Isaac PetscInt Njk; 931a989b97SToby Isaac PetscInt JKj; 941a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 951a989b97SToby Isaac PetscInt i; 961a989b97SToby Isaac 971a989b97SToby Isaac ierr = PetscDTBinomial(N, j+k, &Njk);CHKERRQ(ierr); 981a989b97SToby Isaac ierr = PetscDTBinomial(j+k, j, &JKj);CHKERRQ(ierr); 991a989b97SToby Isaac ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr); 1001a989b97SToby Isaac for (i = 0; i < Njk; i++) { 1011a989b97SToby Isaac PetscReal sum = 0.; 1021a989b97SToby Isaac PetscInt l; 1031a989b97SToby Isaac 1041a989b97SToby Isaac ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr); 1051a989b97SToby Isaac for (l = 0; l < JKj; l++) { 1061a989b97SToby Isaac PetscBool jkOdd; 1071a989b97SToby Isaac PetscInt m, jInd, kInd; 1081a989b97SToby Isaac 1091a989b97SToby Isaac ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr); 1101a989b97SToby Isaac for (m = 0; m < j; m++) { 1111a989b97SToby Isaac subsetj[m] = subset[subsetjk[m]]; 1121a989b97SToby Isaac } 1131a989b97SToby Isaac for (m = 0; m < k; m++) { 1141a989b97SToby Isaac subsetk[m] = subset[subsetjk[j+m]]; 1151a989b97SToby Isaac } 1161a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr); 1171a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr); 1181a989b97SToby Isaac sum += jkOdd ? -(a[jInd] * b[kInd]) : (a[jInd] * b[kInd]); 1191a989b97SToby Isaac } 1201a989b97SToby Isaac awedgeb[i] = sum; 1211a989b97SToby Isaac } 1221a989b97SToby Isaac ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr); 1231a989b97SToby Isaac } 1241a989b97SToby Isaac PetscFunctionReturn(0); 1251a989b97SToby Isaac } 1261a989b97SToby Isaac 1271a989b97SToby Isaac PetscErrorCode PetscDTAltVWedgeMatrix(PetscInt N, PetscInt j, PetscInt k, const PetscReal *a, PetscReal *awedgeMat) 1281a989b97SToby Isaac { 1291a989b97SToby Isaac PetscInt i; 1301a989b97SToby Isaac PetscErrorCode ierr; 1311a989b97SToby Isaac 1321a989b97SToby Isaac PetscFunctionBegin; 1331a989b97SToby Isaac if (N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimension"); 1341a989b97SToby Isaac if (j < 0 || k < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "negative form degree"); 1351a989b97SToby Isaac if (j + k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wedge greater than dimension"); 1361a989b97SToby Isaac if (N <= 3) { 1371a989b97SToby Isaac PetscInt Njk; 1381a989b97SToby Isaac 1391a989b97SToby Isaac ierr = PetscDTBinomial(N, j+k, &Njk);CHKERRQ(ierr); 1401a989b97SToby Isaac if (!j) { 1411a989b97SToby Isaac for (i = 0; i < Njk * Njk; i++) {awedgeMat[i] = 0.;} 1421a989b97SToby Isaac for (i = 0; i < Njk; i++) {awedgeMat[i * (Njk + 1)] = a[0];} 1431a989b97SToby Isaac } else if (!k) { 1441a989b97SToby Isaac for (i = 0; i < Njk; i++) {awedgeMat[i] = a[i];} 1451a989b97SToby Isaac } else { 1461a989b97SToby Isaac if (N == 2) { 1471a989b97SToby Isaac awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; 1481a989b97SToby Isaac } else { 1491a989b97SToby Isaac if (j+k == 2) { 1501a989b97SToby Isaac awedgeMat[0] = -a[1]; awedgeMat[1] = a[0]; awedgeMat[2] = 0.; 1511a989b97SToby Isaac awedgeMat[3] = -a[2]; awedgeMat[4] = 0.; awedgeMat[5] = a[0]; 1521a989b97SToby Isaac awedgeMat[6] = 0.; awedgeMat[7] = -a[2]; awedgeMat[8] = a[1]; 1531a989b97SToby Isaac } else { 1541a989b97SToby Isaac awedgeMat[0] = a[2]; awedgeMat[1] = -a[1]; awedgeMat[2] = a[0]; 1551a989b97SToby Isaac } 1561a989b97SToby Isaac } 1571a989b97SToby Isaac } 1581a989b97SToby Isaac } else { 1591a989b97SToby Isaac PetscInt Njk; 1601a989b97SToby Isaac PetscInt Nk; 1611a989b97SToby Isaac PetscInt JKj, i; 1621a989b97SToby Isaac PetscInt *subset, *subsetjk, *subsetj, *subsetk; 1631a989b97SToby Isaac 1641a989b97SToby Isaac ierr = PetscDTBinomial(N, k, &Nk);CHKERRQ(ierr); 1651a989b97SToby Isaac ierr = PetscDTBinomial(N, j+k, &Njk);CHKERRQ(ierr); 1661a989b97SToby Isaac ierr = PetscDTBinomial(j+k, j, &JKj);CHKERRQ(ierr); 1671a989b97SToby Isaac ierr = PetscMalloc4(j+k, &subset, j+k, &subsetjk, j, &subsetj, k, &subsetk);CHKERRQ(ierr); 1681a989b97SToby Isaac for (i = 0; i < Njk * Nk; i++) awedgeMat[i] = 0.; 1691a989b97SToby Isaac for (i = 0; i < Njk; i++) { 1701a989b97SToby Isaac PetscInt l; 1711a989b97SToby Isaac 1721a989b97SToby Isaac ierr = PetscDTEnumSubset(N, j+k, i, subset);CHKERRQ(ierr); 1731a989b97SToby Isaac for (l = 0; l < JKj; l++) { 1741a989b97SToby Isaac PetscBool jkOdd; 1751a989b97SToby Isaac PetscInt m, jInd, kInd; 1761a989b97SToby Isaac 1771a989b97SToby Isaac ierr = PetscDTEnumSplit(j+k, j, l, subsetjk, &jkOdd);CHKERRQ(ierr); 1781a989b97SToby Isaac for (m = 0; m < j; m++) { 1791a989b97SToby Isaac subsetj[m] = subset[subsetjk[m]]; 1801a989b97SToby Isaac } 1811a989b97SToby Isaac for (m = 0; m < k; m++) { 1821a989b97SToby Isaac subsetk[m] = subset[subsetjk[j+m]]; 1831a989b97SToby Isaac } 1841a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, j, subsetj, &jInd);CHKERRQ(ierr); 1851a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k, subsetk, &kInd);CHKERRQ(ierr); 1861a989b97SToby Isaac awedgeMat[i * Nk + kInd] += jkOdd ? - a[jInd] : a[jInd]; 1871a989b97SToby Isaac } 1881a989b97SToby Isaac } 1891a989b97SToby Isaac ierr = PetscFree4(subset, subsetjk, subsetj, subsetk);CHKERRQ(ierr); 1901a989b97SToby Isaac } 1911a989b97SToby Isaac PetscFunctionReturn(0); 1921a989b97SToby Isaac } 1931a989b97SToby Isaac 1941a989b97SToby Isaac /* L: V -> W [|W| by |V| array], L*: altW -> altV */ 1951a989b97SToby Isaac PetscErrorCode PetscDTAltVPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *Lstarw) 1961a989b97SToby Isaac { 1971a989b97SToby Isaac PetscInt i, j, Nk, Mk; 1981a989b97SToby Isaac PetscErrorCode ierr; 1991a989b97SToby Isaac 2001a989b97SToby Isaac PetscFunctionBegin; 2011a989b97SToby Isaac if (N < 0 || M < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 2021a989b97SToby Isaac if (PetscAbsInt(k) > N || PetscAbsInt(k) > M) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 2031a989b97SToby Isaac if (N <= 3 && M <= 3) { 2041a989b97SToby Isaac 2051a989b97SToby Isaac ierr = PetscDTBinomial(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 2061a989b97SToby Isaac ierr = PetscDTBinomial(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 2071a989b97SToby Isaac if (!k) { 2081a989b97SToby Isaac Lstarw[0] = w[0]; 2091a989b97SToby Isaac } else if (k == 1) { 2101a989b97SToby Isaac for (i = 0; i < Nk; i++) { 2111a989b97SToby Isaac PetscReal sum = 0.; 2121a989b97SToby Isaac 2131a989b97SToby Isaac for (j = 0; j < Mk; j++) {sum += L[j * Nk + i] * w[j];} 2141a989b97SToby Isaac Lstarw[i] = sum; 2151a989b97SToby Isaac } 2161a989b97SToby Isaac } else if (k == -1) { 2171a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 2181a989b97SToby Isaac 2191a989b97SToby Isaac for (i = 0; i < Nk; i++) { 2201a989b97SToby Isaac PetscReal sum = 0.; 2211a989b97SToby Isaac 2221a989b97SToby Isaac for (j = 0; j < Mk; j++) { 2231a989b97SToby Isaac sum += L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * w[j] * mult[j]; 2241a989b97SToby Isaac } 2251a989b97SToby Isaac Lstarw[i] = mult[i] * sum; 2261a989b97SToby Isaac } 2271a989b97SToby Isaac } else if (k == 2) { 2281a989b97SToby Isaac PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}}; 2291a989b97SToby Isaac 2301a989b97SToby Isaac for (i = 0; i < Nk; i++) { 2311a989b97SToby Isaac PetscReal sum = 0.; 2321a989b97SToby Isaac for (j = 0; j < Mk; j++) { 2331a989b97SToby Isaac sum += (L[pairs[j][0] * N + pairs[i][0]] * L[pairs[j][1] * N + pairs[i][1]] - 2341a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][0]] * L[pairs[j][0] * N + pairs[i][1]]) * w[j]; 2351a989b97SToby Isaac } 2361a989b97SToby Isaac Lstarw[i] = sum; 2371a989b97SToby Isaac } 2381a989b97SToby Isaac } else if (k == -2) { 2391a989b97SToby Isaac PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}}; 2401a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 2411a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 2421a989b97SToby Isaac 2431a989b97SToby Isaac for (i = 0; i < Nk; i++) { 2441a989b97SToby Isaac PetscReal sum = 0.; 2451a989b97SToby Isaac 2461a989b97SToby Isaac for (j = 0; j < Mk; j++) { 2471a989b97SToby Isaac sum += (L[pairs[offj + j][0] * N + pairs[offi + i][0]] * 2481a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][1]] - 2491a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][0]] * 2501a989b97SToby Isaac L[pairs[offj + j][0] * N + pairs[offi + i][1]]) * w[j]; 2511a989b97SToby Isaac 2521a989b97SToby Isaac } 2531a989b97SToby Isaac Lstarw[i] = sum; 2541a989b97SToby Isaac } 2551a989b97SToby Isaac } else { 2561a989b97SToby Isaac PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + 2571a989b97SToby Isaac L[1] * (L[5] * L[6] - L[3] * L[8]) + 2581a989b97SToby Isaac L[2] * (L[3] * L[7] - L[4] * L[6]); 2591a989b97SToby Isaac 2601a989b97SToby Isaac for (i = 0; i < Nk; i++) {Lstarw[i] = detL * w[i];} 2611a989b97SToby Isaac } 2621a989b97SToby Isaac } else { 2631a989b97SToby Isaac PetscInt Nf, l, p; 2641a989b97SToby Isaac PetscReal *Lw, *Lwv; 2651a989b97SToby Isaac PetscInt *subsetw, *subsetv; 2661a989b97SToby Isaac PetscInt *work, *perm; 2671a989b97SToby Isaac PetscReal *walloc = NULL; 2681a989b97SToby Isaac const PetscReal *ww = NULL; 2691a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 2701a989b97SToby Isaac 2711a989b97SToby Isaac ierr = PetscDTBinomial(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 2721a989b97SToby Isaac ierr = PetscDTBinomial(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 2731a989b97SToby Isaac ierr = PetscDTFactorialInt_Internal(PetscAbsInt(k), &Nf);CHKERRQ(ierr); 2741a989b97SToby Isaac if (k < 0) { 2751a989b97SToby Isaac negative = PETSC_TRUE; 2761a989b97SToby Isaac k = -k; 2771a989b97SToby Isaac ierr = PetscMalloc1(Mk, &walloc);CHKERRQ(ierr); 2781a989b97SToby Isaac ierr = PetscDTAltVStar(M, M - k, 1, w, walloc);CHKERRQ(ierr); 2791a989b97SToby Isaac ww = walloc; 2801a989b97SToby Isaac } else { 2811a989b97SToby Isaac ww = w; 2821a989b97SToby Isaac } 2831a989b97SToby Isaac ierr = PetscMalloc6(k, &subsetw, k, &subsetv, k, &work, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr); 2841a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = 0.; 2851a989b97SToby Isaac for (i = 0; i < Mk; i++) { 2861a989b97SToby Isaac ierr = PetscDTEnumSubset(M, k, i, subsetw);CHKERRQ(ierr); 2871a989b97SToby Isaac for (j = 0; j < Nk; j++) { 2881a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, j, subsetv);CHKERRQ(ierr); 2891a989b97SToby Isaac for (p = 0; p < Nf; p++) { 2901a989b97SToby Isaac PetscReal prod; 2911a989b97SToby Isaac PetscBool isOdd; 2921a989b97SToby Isaac 2931a989b97SToby Isaac ierr = PetscDTEnumPerm(k, p, work, perm, &isOdd);CHKERRQ(ierr); 2941a989b97SToby Isaac prod = isOdd ? -ww[i] : ww[i]; 2951a989b97SToby Isaac for (l = 0; l < k; l++) { 2961a989b97SToby Isaac prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 2971a989b97SToby Isaac } 2981a989b97SToby Isaac Lstarw[j] += prod; 2991a989b97SToby Isaac } 3001a989b97SToby Isaac } 3011a989b97SToby Isaac } 3021a989b97SToby Isaac if (negative) { 3031a989b97SToby Isaac PetscReal *sLsw; 3041a989b97SToby Isaac 3051a989b97SToby Isaac ierr = PetscMalloc1(Nk, &sLsw);CHKERRQ(ierr); 3061a989b97SToby Isaac ierr = PetscDTAltVStar(N, N - k, -1, Lstarw, sLsw);CHKERRQ(ierr); 3071a989b97SToby Isaac for (i = 0; i < Nk; i++) Lstarw[i] = sLsw[i]; 3081a989b97SToby Isaac ierr = PetscFree(sLsw);CHKERRQ(ierr); 3091a989b97SToby Isaac } 3101a989b97SToby Isaac ierr = PetscFree6(subsetw, subsetv, work, perm, Lw, Lwv);CHKERRQ(ierr); 3111a989b97SToby Isaac ierr = PetscFree(walloc);CHKERRQ(ierr); 3121a989b97SToby Isaac } 3131a989b97SToby Isaac PetscFunctionReturn(0); 3141a989b97SToby Isaac } 3151a989b97SToby Isaac 3161a989b97SToby Isaac PetscErrorCode PetscDTAltVPullbackMatrix(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, PetscReal *Lstar) 3171a989b97SToby Isaac { 3181a989b97SToby Isaac PetscInt Nk, Mk, Nf, i, j, l, p; 3191a989b97SToby Isaac PetscReal *Lw, *Lwv; 3201a989b97SToby Isaac PetscInt *subsetw, *subsetv; 3211a989b97SToby Isaac PetscInt *work, *perm; 3221a989b97SToby Isaac PetscBool negative = PETSC_FALSE; 3231a989b97SToby Isaac PetscErrorCode ierr; 3241a989b97SToby Isaac 3251a989b97SToby Isaac PetscFunctionBegin; 3261a989b97SToby Isaac if (N < 0 || M < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid dimensions"); 3271a989b97SToby Isaac if (PetscAbsInt(k) > N || PetscAbsInt(k) > M) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 3281a989b97SToby Isaac if (N <= 3 && M <= 3) { 3291a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 3301a989b97SToby Isaac 3311a989b97SToby Isaac ierr = PetscDTBinomial(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 3321a989b97SToby Isaac ierr = PetscDTBinomial(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 3331a989b97SToby Isaac if (!k) { 3341a989b97SToby Isaac Lstar[0] = 1.; 3351a989b97SToby Isaac } else if (k == 1) { 3361a989b97SToby Isaac for (i = 0; i < Nk; i++) {for (j = 0; j < Mk; j++) {Lstar[i * Mk + j] = L[j * Nk + i];}} 3371a989b97SToby Isaac } else if (k == -1) { 3381a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3391a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3401a989b97SToby Isaac Lstar[i * Mk + j] = L[(Mk - 1 - j) * Nk + (Nk - 1 - i)] * mult[i] * mult[j]; 3411a989b97SToby Isaac } 3421a989b97SToby Isaac } 3431a989b97SToby Isaac } else if (k == 2) { 3441a989b97SToby Isaac PetscInt pairs[3][2] = {{0,1},{0,2},{1,2}}; 3451a989b97SToby Isaac 3461a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3471a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3481a989b97SToby Isaac Lstar[i * Mk + j] = L[pairs[j][0] * N + pairs[i][0]] * 3491a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][1]] - 3501a989b97SToby Isaac L[pairs[j][1] * N + pairs[i][0]] * 3511a989b97SToby Isaac L[pairs[j][0] * N + pairs[i][1]]; 3521a989b97SToby Isaac } 3531a989b97SToby Isaac } 3541a989b97SToby Isaac } else if (k == -2) { 3551a989b97SToby Isaac PetscInt pairs[3][2] = {{1,2},{2,0},{0,1}}; 3561a989b97SToby Isaac PetscInt offi = (N == 2) ? 2 : 0; 3571a989b97SToby Isaac PetscInt offj = (M == 2) ? 2 : 0; 3581a989b97SToby Isaac 3591a989b97SToby Isaac for (i = 0; i < Nk; i++) { 3601a989b97SToby Isaac for (j = 0; j < Mk; j++) { 3611a989b97SToby Isaac Lstar[i * Mk + j] = L[pairs[offj + j][0] * N + pairs[offi + i][0]] * 3621a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][1]] - 3631a989b97SToby Isaac L[pairs[offj + j][1] * N + pairs[offi + i][0]] * 3641a989b97SToby Isaac L[pairs[offj + j][0] * N + pairs[offi + i][1]]; 3651a989b97SToby Isaac } 3661a989b97SToby Isaac } 3671a989b97SToby Isaac } else { 3681a989b97SToby Isaac PetscReal detL = L[0] * (L[4] * L[8] - L[5] * L[7]) + 3691a989b97SToby Isaac L[1] * (L[5] * L[6] - L[3] * L[8]) + 3701a989b97SToby Isaac L[2] * (L[3] * L[7] - L[4] * L[6]); 3711a989b97SToby Isaac 3721a989b97SToby Isaac for (i = 0; i < Nk; i++) {Lstar[i] = detL;} 3731a989b97SToby Isaac } 3741a989b97SToby Isaac } else { 3751a989b97SToby Isaac if (k < 0) { 3761a989b97SToby Isaac negative = PETSC_TRUE; 3771a989b97SToby Isaac k = -k; 3781a989b97SToby Isaac } 3791a989b97SToby Isaac ierr = PetscDTBinomial(M, PetscAbsInt(k), &Mk);CHKERRQ(ierr); 3801a989b97SToby Isaac ierr = PetscDTBinomial(N, PetscAbsInt(k), &Nk);CHKERRQ(ierr); 3811a989b97SToby Isaac ierr = PetscDTFactorialInt_Internal(PetscAbsInt(k), &Nf);CHKERRQ(ierr); 3821a989b97SToby Isaac ierr = PetscMalloc6(M, &subsetw, N, &subsetv, k, &work, k, &perm, N * k, &Lw, k * k, &Lwv);CHKERRQ(ierr); 3831a989b97SToby Isaac for (i = 0; i < Nk * Mk; i++) Lstar[i] = 0.; 3841a989b97SToby Isaac for (i = 0; i < Mk; i++) { 3851a989b97SToby Isaac PetscBool iOdd; 3861a989b97SToby Isaac PetscInt iidx, jidx; 3871a989b97SToby Isaac 3881a989b97SToby Isaac ierr = PetscDTEnumSplit(M, k, i, subsetw, &iOdd);CHKERRQ(ierr); 3891a989b97SToby Isaac iidx = negative ? Mk - 1 - i : i; 3901a989b97SToby Isaac iOdd = negative ? iOdd ^ ((k * (M-k)) & 1) : PETSC_FALSE; 3911a989b97SToby Isaac for (j = 0; j < Nk; j++) { 3921a989b97SToby Isaac PetscBool jOdd; 3931a989b97SToby Isaac 3941a989b97SToby Isaac ierr = PetscDTEnumSplit(N, k, j, subsetv, &jOdd);CHKERRQ(ierr); 3951a989b97SToby Isaac jidx = negative ? Nk - 1 - j : j; 3961a989b97SToby Isaac jOdd = negative ? iOdd ^ jOdd ^ ((k * (N-k)) & 1) : PETSC_FALSE; 3971a989b97SToby Isaac for (p = 0; p < Nf; p++) { 3981a989b97SToby Isaac PetscReal prod; 3991a989b97SToby Isaac PetscBool isOdd; 4001a989b97SToby Isaac 4011a989b97SToby Isaac ierr = PetscDTEnumPerm(k, p, work, perm, &isOdd);CHKERRQ(ierr); 4021a989b97SToby Isaac isOdd ^= jOdd; 4031a989b97SToby Isaac prod = isOdd ? -1. : 1.; 4041a989b97SToby Isaac for (l = 0; l < k; l++) { 4051a989b97SToby Isaac prod *= L[subsetw[perm[l]] * N + subsetv[l]]; 4061a989b97SToby Isaac } 4071a989b97SToby Isaac Lstar[jidx * Mk + iidx] += prod; 4081a989b97SToby Isaac } 4091a989b97SToby Isaac } 4101a989b97SToby Isaac } 4111a989b97SToby Isaac ierr = PetscFree6(subsetw, subsetv, work, perm, Lw, Lwv);CHKERRQ(ierr); 4121a989b97SToby Isaac } 4131a989b97SToby Isaac PetscFunctionReturn(0); 4141a989b97SToby Isaac } 4151a989b97SToby Isaac 4161a989b97SToby Isaac PetscErrorCode PetscDTAltVInterior(PetscInt N, PetscInt k, const PetscReal *w, const PetscReal *v, PetscReal *wIntv) 4171a989b97SToby Isaac { 4181a989b97SToby Isaac PetscInt i, Nk, Nkm; 4191a989b97SToby Isaac PetscErrorCode ierr; 4201a989b97SToby Isaac 4211a989b97SToby Isaac PetscFunctionBegin; 4221a989b97SToby Isaac if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 4231a989b97SToby Isaac ierr = PetscDTBinomial(N, k, &Nk);CHKERRQ(ierr); 4241a989b97SToby Isaac ierr = PetscDTBinomial(N, k-1, &Nkm);CHKERRQ(ierr); 4251a989b97SToby Isaac if (N <= 3) { 4261a989b97SToby Isaac if (k == 1) { 4271a989b97SToby Isaac PetscReal sum = 0.; 4281a989b97SToby Isaac 4291a989b97SToby Isaac for (i = 0; i < N; i++) { 4301a989b97SToby Isaac sum += w[i] * v[i]; 4311a989b97SToby Isaac } 4321a989b97SToby Isaac wIntv[0] = sum; 4331a989b97SToby Isaac } else if (k == N) { 4341a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 4351a989b97SToby Isaac 4361a989b97SToby Isaac for (i = 0; i < N; i++) { 4371a989b97SToby Isaac wIntv[N - 1 - i] = w[0] * v[i] * mult[i]; 4381a989b97SToby Isaac } 4391a989b97SToby Isaac } else { 4401a989b97SToby Isaac wIntv[0] = - w[0]*v[1] - w[1]*v[2]; 4411a989b97SToby Isaac wIntv[1] = w[0]*v[0] - w[2]*v[2]; 4421a989b97SToby Isaac wIntv[2] = w[1]*v[0] + w[2]*v[1]; 4431a989b97SToby Isaac } 4441a989b97SToby Isaac } else { 4451a989b97SToby Isaac PetscInt *subset, *work; 4461a989b97SToby Isaac 4471a989b97SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 4481a989b97SToby Isaac for (i = 0; i < Nkm; i++) wIntv[i] = 0.; 4491a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4501a989b97SToby Isaac PetscInt j, l, m; 4511a989b97SToby Isaac 4521a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 4531a989b97SToby Isaac for (j = 0; j < k; j++) { 4541a989b97SToby Isaac PetscInt idx; 4551a989b97SToby Isaac PetscBool flip = (j & 1); 4561a989b97SToby Isaac 4571a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 4581a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 4591a989b97SToby Isaac } 4601a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 4611a989b97SToby Isaac wIntv[idx] += flip ? -(w[i] * v[subset[j]]) : (w[i] * v[subset[j]]); 4621a989b97SToby Isaac } 4631a989b97SToby Isaac } 4641a989b97SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 4651a989b97SToby Isaac } 4661a989b97SToby Isaac PetscFunctionReturn(0); 4671a989b97SToby Isaac } 4681a989b97SToby Isaac 4691a989b97SToby Isaac PetscErrorCode PetscDTAltVInteriorMatrix(PetscInt N, PetscInt k, const PetscReal *v, PetscReal *intvMat) 4701a989b97SToby Isaac { 4711a989b97SToby Isaac PetscInt i, Nk, Nkm; 4721a989b97SToby Isaac PetscErrorCode ierr; 4731a989b97SToby Isaac 4741a989b97SToby Isaac PetscFunctionBegin; 4751a989b97SToby Isaac if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 4761a989b97SToby Isaac ierr = PetscDTBinomial(N, k, &Nk);CHKERRQ(ierr); 4771a989b97SToby Isaac ierr = PetscDTBinomial(N, k-1, &Nkm);CHKERRQ(ierr); 4781a989b97SToby Isaac if (N <= 3) { 4791a989b97SToby Isaac if (k == 1) { 4801a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[i] = v[i]; 4811a989b97SToby Isaac } else if (k == N) { 4821a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 4831a989b97SToby Isaac 4841a989b97SToby Isaac for (i = 0; i < N; i++) intvMat[N - 1 - i] = v[i] * mult[i]; 4851a989b97SToby Isaac } else { 4861a989b97SToby Isaac intvMat[0] = -v[1]; intvMat[1] = -v[2]; intvMat[2] = 0.; 4871a989b97SToby Isaac intvMat[3] = v[0]; intvMat[4] = 0.; intvMat[5] = -v[2]; 4881a989b97SToby Isaac intvMat[6] = 0.; intvMat[7] = v[0]; intvMat[8] = v[1]; 4891a989b97SToby Isaac } 4901a989b97SToby Isaac } else { 4911a989b97SToby Isaac PetscInt *subset, *work; 4921a989b97SToby Isaac 4931a989b97SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 4941a989b97SToby Isaac for (i = 0; i < Nk * Nkm; i++) intvMat[i] = 0.; 4951a989b97SToby Isaac for (i = 0; i < Nk; i++) { 4961a989b97SToby Isaac PetscInt j, l, m; 4971a989b97SToby Isaac 4981a989b97SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 4991a989b97SToby Isaac for (j = 0; j < k; j++) { 5001a989b97SToby Isaac PetscInt idx; 5011a989b97SToby Isaac PetscBool flip = (j & 1); 5021a989b97SToby Isaac 5031a989b97SToby Isaac for (l = 0, m = 0; l < k; l++) { 5041a989b97SToby Isaac if (l != j) work[m++] = subset[l]; 5051a989b97SToby Isaac } 5061a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 5071a989b97SToby Isaac intvMat[idx * Nk + i] += flip ? -v[subset[j]] : v[subset[j]]; 5081a989b97SToby Isaac } 5091a989b97SToby Isaac } 5101a989b97SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 5111a989b97SToby Isaac } 5121a989b97SToby Isaac PetscFunctionReturn(0); 5131a989b97SToby Isaac } 5141a989b97SToby Isaac 515*dda711d0SToby Isaac PetscErrorCode PetscDTAltVInteriorPattern(PetscInt N, PetscInt k, PetscInt (*indices)[3]) 516*dda711d0SToby Isaac { 517*dda711d0SToby Isaac PetscInt i, Nk, Nkm; 518*dda711d0SToby Isaac PetscErrorCode ierr; 519*dda711d0SToby Isaac 520*dda711d0SToby Isaac PetscFunctionBegin; 521*dda711d0SToby Isaac if (k <= 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 522*dda711d0SToby Isaac ierr = PetscDTBinomial(N, k, &Nk);CHKERRQ(ierr); 523*dda711d0SToby Isaac ierr = PetscDTBinomial(N, k-1, &Nkm);CHKERRQ(ierr); 524*dda711d0SToby Isaac if (N <= 3) { 525*dda711d0SToby Isaac if (k == 1) { 526*dda711d0SToby Isaac for (i = 0; i < N; i++) { 527*dda711d0SToby Isaac indices[i][0] = 0; 528*dda711d0SToby Isaac indices[i][1] = i; 529*dda711d0SToby Isaac indices[i][2] = i; 530*dda711d0SToby Isaac } 531*dda711d0SToby Isaac } else if (k == N) { 532*dda711d0SToby Isaac PetscInt val[3] = {0, -2, 2}; 533*dda711d0SToby Isaac 534*dda711d0SToby Isaac for (i = 0; i < N; i++) { 535*dda711d0SToby Isaac indices[i][0] = N - 1 - i; 536*dda711d0SToby Isaac indices[i][1] = 0; 537*dda711d0SToby Isaac indices[i][2] = val[i]; 538*dda711d0SToby Isaac } 539*dda711d0SToby Isaac } else { 540*dda711d0SToby Isaac indices[0][0] = 0; indices[0][1] = 0; indices[0][2] = -(1 + 1); 541*dda711d0SToby Isaac indices[1][0] = 0; indices[1][1] = 1; indices[1][2] = -(2 + 1); 542*dda711d0SToby Isaac indices[2][0] = 1; indices[2][1] = 0; indices[2][2] = 0; 543*dda711d0SToby Isaac indices[3][0] = 1; indices[3][1] = 2; indices[3][2] = -(2 + 1); 544*dda711d0SToby Isaac indices[4][0] = 2; indices[4][1] = 1; indices[4][2] = 0; 545*dda711d0SToby Isaac indices[5][0] = 2; indices[5][1] = 2; indices[5][2] = 1; 546*dda711d0SToby Isaac } 547*dda711d0SToby Isaac } else { 548*dda711d0SToby Isaac PetscInt *subset, *work; 549*dda711d0SToby Isaac 550*dda711d0SToby Isaac ierr = PetscMalloc2(k, &subset, k, &work);CHKERRQ(ierr); 551*dda711d0SToby Isaac for (i = 0; i < Nk; i++) { 552*dda711d0SToby Isaac PetscInt j, l, m; 553*dda711d0SToby Isaac 554*dda711d0SToby Isaac ierr = PetscDTEnumSubset(N, k, i, subset);CHKERRQ(ierr); 555*dda711d0SToby Isaac for (j = 0; j < k; j++) { 556*dda711d0SToby Isaac PetscInt idx; 557*dda711d0SToby Isaac PetscBool flip = (j & 1); 558*dda711d0SToby Isaac 559*dda711d0SToby Isaac for (l = 0, m = 0; l < k; l++) { 560*dda711d0SToby Isaac if (l != j) work[m++] = subset[l]; 561*dda711d0SToby Isaac } 562*dda711d0SToby Isaac ierr = PetscDTSubsetIndex(N, k - 1, work, &idx);CHKERRQ(ierr); 563*dda711d0SToby Isaac indices[i * k + j][0] = idx; 564*dda711d0SToby Isaac indices[i * k + j][1] = i; 565*dda711d0SToby Isaac indices[i * k + j][2] = flip ? -(subset[j] + 1) : subset[j]; 566*dda711d0SToby Isaac } 567*dda711d0SToby Isaac } 568*dda711d0SToby Isaac ierr = PetscFree2(subset, work);CHKERRQ(ierr); 569*dda711d0SToby Isaac } 570*dda711d0SToby Isaac PetscFunctionReturn(0); 571*dda711d0SToby Isaac } 572*dda711d0SToby Isaac 5731a989b97SToby Isaac PetscErrorCode PetscDTAltVStar(PetscInt N, PetscInt k, PetscInt pow, const PetscReal *w, PetscReal *starw) 5741a989b97SToby Isaac { 5751a989b97SToby Isaac PetscInt Nk, i; 5761a989b97SToby Isaac PetscErrorCode ierr; 5771a989b97SToby Isaac 5781a989b97SToby Isaac PetscFunctionBegin; 5791a989b97SToby Isaac if (k < 0 || k > N) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid form degree"); 5801a989b97SToby Isaac ierr = PetscDTBinomial(N, k, &Nk);CHKERRQ(ierr); 5811a989b97SToby Isaac pow = pow % 4; 5821a989b97SToby Isaac pow = (pow + 4) % 4; /* make non-negative */ 5831a989b97SToby Isaac /* pow is now 0, 1, 2, 3 */ 5841a989b97SToby Isaac if (N <= 3) { 5851a989b97SToby Isaac if (pow & 1) { 5861a989b97SToby Isaac PetscReal mult[3] = {1., -1., 1.}; 5871a989b97SToby Isaac 5881a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[Nk - 1 - i] = w[i] * mult[i]; 5891a989b97SToby Isaac } else { 5901a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 5911a989b97SToby Isaac } 5921a989b97SToby Isaac if (pow > 1 && ((k * (N - k)) & 1)) { 5931a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 5941a989b97SToby Isaac } 5951a989b97SToby Isaac } else { 5961a989b97SToby Isaac PetscInt *subset; 5971a989b97SToby Isaac 5981a989b97SToby Isaac ierr = PetscMalloc1(N, &subset);CHKERRQ(ierr); 5991a989b97SToby Isaac if (pow % 2) { 6001a989b97SToby Isaac PetscInt l = (pow == 1) ? k : N - k; 6011a989b97SToby Isaac for (i = 0; i < Nk; i++) { 6021a989b97SToby Isaac PetscBool sOdd; 6031a989b97SToby Isaac PetscInt j, idx; 6041a989b97SToby Isaac 6051a989b97SToby Isaac ierr = PetscDTEnumSplit(N, l, i, subset, &sOdd);CHKERRQ(ierr); 6061a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, l, subset, &idx);CHKERRQ(ierr); 6071a989b97SToby Isaac ierr = PetscDTSubsetIndex(N, N-l, &subset[l], &j);CHKERRQ(ierr); 6081a989b97SToby Isaac starw[j] = sOdd ? -w[idx] : w[idx]; 6091a989b97SToby Isaac } 6101a989b97SToby Isaac } else { 6111a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = w[i]; 6121a989b97SToby Isaac } 6131a989b97SToby Isaac /* star^2 = -1^(k * (N - k)) */ 6141a989b97SToby Isaac if (pow > 1 && (k * (N - k)) % 2) { 6151a989b97SToby Isaac for (i = 0; i < Nk; i++) starw[i] = -starw[i]; 6161a989b97SToby Isaac } 6171a989b97SToby Isaac ierr = PetscFree(subset);CHKERRQ(ierr); 6181a989b97SToby Isaac } 6191a989b97SToby Isaac PetscFunctionReturn(0); 6201a989b97SToby Isaac } 621