1c4762a1bSJed Brown static char help[] = "Test the PetscDTAltV interface for k-forms (alternating k-linear maps).\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscviewer.h> 4c4762a1bSJed Brown #include <petscdt.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *x, PetscBool verbose, PetscViewer viewer) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown PetscInt Nk, Mk, i, j, l; 9c4762a1bSJed Brown PetscReal *Lstarw, *Lx, *Lstar, *Lstarwcheck, wLx, Lstarwx; 10c4762a1bSJed Brown PetscReal diff, diffMat, normMat; 11c4762a1bSJed Brown PetscReal *walloc = NULL; 12c4762a1bSJed Brown const PetscReal *ww = NULL; 13c4762a1bSJed Brown PetscBool negative = (PetscBool)(k < 0); 14c4762a1bSJed Brown 15c4762a1bSJed Brown PetscFunctionBegin; 16c4762a1bSJed Brown k = PetscAbsInt(k); 179566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 189566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, k, &Mk)); 19c4762a1bSJed Brown if (negative) { 209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Mk, &walloc)); 219566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc)); 22c4762a1bSJed Brown ww = walloc; 23c4762a1bSJed Brown } else { 24c4762a1bSJed Brown ww = w; 25c4762a1bSJed Brown } 26*57508eceSPierre Jolivet PetscCall(PetscMalloc2(Nk, &Lstarw, M * k, &Lx)); 279566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nk * Mk, &Lstar, Nk, &Lstarwcheck)); 289566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullback(N, M, L, negative ? -k : k, w, Lstarw)); 299566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(N, M, L, negative ? -k : k, Lstar)); 30c4762a1bSJed Brown if (negative) { 31c4762a1bSJed Brown PetscReal *sLsw; 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk, &sLsw)); 349566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(N, N - k, 1, Lstarw, sLsw)); 359566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k, sLsw, x, &Lstarwx)); 369566063dSJacob Faibussowitsch PetscCall(PetscFree(sLsw)); 37c4762a1bSJed Brown } else { 389566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k, Lstarw, x, &Lstarwx)); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown for (l = 0; l < k; l++) { 41c4762a1bSJed Brown for (i = 0; i < M; i++) { 42c4762a1bSJed Brown PetscReal sum = 0.; 43c4762a1bSJed Brown 44c4762a1bSJed Brown for (j = 0; j < N; j++) sum += L[i * N + j] * x[l * N + j]; 45c4762a1bSJed Brown Lx[l * M + i] = sum; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown } 48c4762a1bSJed Brown diffMat = 0.; 49c4762a1bSJed Brown normMat = 0.; 50c4762a1bSJed Brown for (i = 0; i < Nk; i++) { 51c4762a1bSJed Brown PetscReal sum = 0.; 52ad540459SPierre Jolivet for (j = 0; j < Mk; j++) sum += Lstar[i * Mk + j] * w[j]; 53c4762a1bSJed Brown Lstarwcheck[i] = sum; 54c4762a1bSJed Brown diffMat += PetscSqr(PetscAbsReal(Lstarwcheck[i] - Lstarw[i])); 55c4762a1bSJed Brown normMat += PetscSqr(Lstarwcheck[i]) + PetscSqr(Lstarw[i]); 56c4762a1bSJed Brown } 57c4762a1bSJed Brown diffMat = PetscSqrtReal(diffMat); 58c4762a1bSJed Brown normMat = PetscSqrtReal(normMat); 59c4762a1bSJed Brown if (verbose) { 609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "L:\n")); 619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 629566063dSJacob Faibussowitsch if (M * N > 0) PetscCall(PetscRealView(M * N, L, viewer)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "L*:\n")); 669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 679566063dSJacob Faibussowitsch if (Nk * Mk > 0) PetscCall(PetscRealView(Nk * Mk, Lstar, viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "L*w:\n")); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 729566063dSJacob Faibussowitsch if (Nk > 0) PetscCall(PetscRealView(Nk, Lstarw, viewer)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 74c4762a1bSJed Brown } 759566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(M, k, ww, Lx, &wLx)); 76c4762a1bSJed Brown diff = PetscAbsReal(wLx - Lstarwx); 771dca8a05SBarry Smith PetscCheck(diff <= 10. * PETSC_SMALL * (PetscAbsReal(wLx) + PetscAbsReal(Lstarwx)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback does not commute with application: w(Lx)(%g) != (L* w)(x)(%g)", (double)wLx, (double)Lstarwx); 7801c1178eSBarry Smith PetscCheck(diffMat <= PETSC_SMALL * normMat, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback matrix does match matrix-free result"); 799566063dSJacob Faibussowitsch PetscCall(PetscFree2(Lstar, Lstarwcheck)); 809566063dSJacob Faibussowitsch PetscCall(PetscFree2(Lstarw, Lx)); 819566063dSJacob Faibussowitsch PetscCall(PetscFree(walloc)); 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 85d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 86d71ae5a4SJacob Faibussowitsch { 87c4762a1bSJed Brown PetscInt i, numTests = 5, n[5] = {0, 1, 2, 3, 4}; 88c4762a1bSJed Brown PetscBool verbose = PETSC_FALSE; 89c4762a1bSJed Brown PetscRandom rand; 90c4762a1bSJed Brown PetscViewer viewer; 91c4762a1bSJed Brown 92327415f7SBarry Smith PetscFunctionBeginUser; 939566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 94d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for exterior algebra tests", "none"); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-N", "Up to 5 vector space dimensions to test", "ex7.c", n, &numTests, NULL)); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-verbose", "Verbose test output", "ex7.c", verbose, &verbose, NULL)); 97d0609cedSBarry Smith PetscOptionsEnd(); 989566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 999566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1., 1.)); 1009566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 101c4762a1bSJed Brown if (!numTests) numTests = 5; 102c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD); 103c4762a1bSJed Brown for (i = 0; i < numTests; i++) { 104c4762a1bSJed Brown PetscInt k, N = n[i]; 105c4762a1bSJed Brown 10663a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "N = %" PetscInt_FMT ":\n", N)); 1079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown if (verbose) { 110c4762a1bSJed Brown PetscInt *perm; 111c4762a1bSJed Brown PetscInt fac = 1; 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &perm)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown for (k = 1; k <= N; k++) fac *= k; 11663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Permutations of %" PetscInt_FMT ":\n", N)); 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 118c4762a1bSJed Brown for (k = 0; k < fac; k++) { 119c4762a1bSJed Brown PetscBool isOdd, isOddCheck; 120c4762a1bSJed Brown PetscInt j, kCheck; 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(PetscDTEnumPerm(N, k, perm, &isOdd)); 12363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ":", k)); 12448a46eb9SPierre Jolivet for (j = 0; j < N; j++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, perm[j])); 1259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even")); 1269566063dSJacob Faibussowitsch PetscCall(PetscDTPermIndex(N, perm, &kCheck, &isOddCheck)); 1271dca8a05SBarry Smith PetscCheck(kCheck == k && isOddCheck == isOdd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumPerm / PetscDTPermIndex mismatch for (%" PetscInt_FMT ", %" PetscInt_FMT ")", N, k); 128c4762a1bSJed Brown } 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(perm)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown for (k = 0; k <= N; k++) { 133c4762a1bSJed Brown PetscInt j, Nk, M; 134c4762a1bSJed Brown PetscReal *w, *v, wv; 135c4762a1bSJed Brown PetscInt *subset; 136c4762a1bSJed Brown 1379566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k, &Nk)); 13863a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "k = %" PetscInt_FMT ":\n", k)); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 14063a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT " choose %" PetscInt_FMT "): %" PetscInt_FMT "\n", N, k, Nk)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Test subset and complement enumeration */ 1439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N, &subset)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 145c4762a1bSJed Brown for (j = 0; j < Nk; j++) { 146c4762a1bSJed Brown PetscBool isOdd, isOddCheck; 147c4762a1bSJed Brown PetscInt jCheck, kCheck; 148c4762a1bSJed Brown 1499566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(N, k, j, subset, &isOdd)); 1509566063dSJacob Faibussowitsch PetscCall(PetscDTPermIndex(N, subset, &kCheck, &isOddCheck)); 15108401ef6SPierre Jolivet PetscCheck(isOddCheck == isOdd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumSplit sign does not mmatch PetscDTPermIndex sign"); 152c4762a1bSJed Brown if (verbose) { 153c4762a1bSJed Brown PetscInt l; 154c4762a1bSJed Brown 15563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "subset %" PetscInt_FMT ":", j)); 15648a46eb9SPierre Jolivet for (l = 0; l < k; l++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, subset[l])); 1579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " |")); 15848a46eb9SPierre Jolivet for (l = k; l < N; l++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, subset[l])); 1599566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even")); 160c4762a1bSJed Brown } 1619566063dSJacob Faibussowitsch PetscCall(PetscDTSubsetIndex(N, k, subset, &jCheck)); 16263a3b9bcSJacob Faibussowitsch PetscCheck(jCheck == j, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "jCheck (%" PetscInt_FMT ") != j (%" PetscInt_FMT ")", jCheck, j); 163c4762a1bSJed Brown } 1649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1659566063dSJacob Faibussowitsch PetscCall(PetscFree(subset)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* Make a random k form */ 1689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk, &w)); 1699566063dSJacob Faibussowitsch for (j = 0; j < Nk; j++) PetscCall(PetscRandomGetValueReal(rand, &w[j])); 170c4762a1bSJed Brown /* Make a set of random vectors */ 1719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * k, &v)); 1729566063dSJacob Faibussowitsch for (j = 0; j < N * k; j++) PetscCall(PetscRandomGetValueReal(rand, &v[j])); 173c4762a1bSJed Brown 1749566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k, w, v, &wv)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown if (verbose) { 1779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "w:\n")); 1789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1799566063dSJacob Faibussowitsch if (Nk) PetscCall(PetscRealView(Nk, w, viewer)); 1809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "v:\n")); 1829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1839566063dSJacob Faibussowitsch if (N * k > 0) PetscCall(PetscRealView(N * k, v, viewer)); 1849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "w(v): %g\n", (double)wv)); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* sanity checks */ 189c4762a1bSJed Brown if (k == 1) { /* 1-forms are functionals (dot products) */ 190c4762a1bSJed Brown PetscInt l; 191c4762a1bSJed Brown PetscReal wvcheck = 0.; 192c4762a1bSJed Brown PetscReal diff; 193c4762a1bSJed Brown 194c4762a1bSJed Brown for (l = 0; l < N; l++) wvcheck += w[l] * v[l]; 195c4762a1bSJed Brown diff = PetscSqrtReal(PetscSqr(wvcheck - wv)); 1961dca8a05SBarry Smith PetscCheck(diff < PETSC_SMALL * (PetscAbsReal(wv) + PetscAbsReal(wvcheck)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "1-form / dot product equivalence: wvcheck (%g) != wv (%g)", (double)wvcheck, (double)wv); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown if (k == N && N < 5) { /* n-forms are scaled determinants */ 199c4762a1bSJed Brown PetscReal det, wvcheck, diff; 200c4762a1bSJed Brown 201c4762a1bSJed Brown switch (k) { 202d71ae5a4SJacob Faibussowitsch case 0: 203d71ae5a4SJacob Faibussowitsch det = 1.; 204d71ae5a4SJacob Faibussowitsch break; 205d71ae5a4SJacob Faibussowitsch case 1: 206d71ae5a4SJacob Faibussowitsch det = v[0]; 207d71ae5a4SJacob Faibussowitsch break; 208d71ae5a4SJacob Faibussowitsch case 2: 209d71ae5a4SJacob Faibussowitsch det = v[0] * v[3] - v[1] * v[2]; 210d71ae5a4SJacob Faibussowitsch break; 211d71ae5a4SJacob Faibussowitsch case 3: 212d71ae5a4SJacob Faibussowitsch det = 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]); 213d71ae5a4SJacob Faibussowitsch break; 214c4762a1bSJed Brown case 4: 2159371c9d4SSatish Balay det = v[0] * (v[5] * (v[10] * v[15] - v[11] * v[14]) + v[6] * (v[11] * v[13] - v[9] * v[15]) + v[7] * (v[9] * v[14] - v[10] * v[13])) - v[1] * (v[4] * (v[10] * v[15] - v[11] * v[14]) + v[6] * (v[11] * v[12] - v[8] * v[15]) + v[7] * (v[8] * v[14] - v[10] * v[12])) + v[2] * (v[4] * (v[9] * v[15] - v[11] * v[13]) + v[5] * (v[11] * v[12] - v[8] * v[15]) + v[7] * (v[8] * v[13] - v[9] * v[12])) - v[3] * (v[4] * (v[9] * v[14] - v[10] * v[13]) + v[5] * (v[10] * v[12] - v[8] * v[14]) + v[6] * (v[8] * v[13] - v[9] * v[12])); 216c4762a1bSJed Brown break; 217d71ae5a4SJacob Faibussowitsch default: 218d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "invalid k"); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown wvcheck = det * w[0]; 221c4762a1bSJed Brown diff = PetscSqrtReal(PetscSqr(wvcheck - wv)); 2221dca8a05SBarry Smith PetscCheck(diff < PETSC_SMALL * (PetscAbsReal(wv) + PetscAbsReal(wvcheck)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "n-form / determinant equivalence: wvcheck (%g) != wv (%g) %g", (double)wvcheck, (double)wv, (double)diff); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown if (k > 0) { /* k-forms are linear in each component */ 225c4762a1bSJed Brown PetscReal alpha; 226c4762a1bSJed Brown PetscReal *x, *axv, wx, waxv, waxvcheck; 227c4762a1bSJed Brown PetscReal diff; 228c4762a1bSJed Brown PetscReal rj; 229c4762a1bSJed Brown PetscInt l; 230c4762a1bSJed Brown 2319566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(N * k, &x, N * k, &axv)); 2329566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &alpha)); 2339566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, 0, k)); 2349566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rj)); 235c4762a1bSJed Brown j = (PetscInt)rj; 2369566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1., 1.)); 237c4762a1bSJed Brown for (l = 0; l < N * k; l++) x[l] = v[l]; 238c4762a1bSJed Brown for (l = 0; l < N * k; l++) axv[l] = v[l]; 239c4762a1bSJed Brown for (l = 0; l < N; l++) { 240c4762a1bSJed Brown PetscReal val; 241c4762a1bSJed Brown 2429566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &val)); 243c4762a1bSJed Brown x[j * N + l] = val; 244c4762a1bSJed Brown axv[j * N + l] += alpha * val; 245c4762a1bSJed Brown } 2469566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k, w, x, &wx)); 2479566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k, w, axv, &waxv)); 248c4762a1bSJed Brown waxvcheck = alpha * wx + wv; 249c4762a1bSJed Brown diff = waxv - waxvcheck; 2501dca8a05SBarry Smith PetscCheck(PetscAbsReal(diff) <= 10. * PETSC_SMALL * (PetscAbsReal(waxv) + PetscAbsReal(waxvcheck)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "linearity check: component %" PetscInt_FMT ", waxvcheck (%g) != waxv (%g)", j, (double)waxvcheck, (double)waxv); 2519566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, axv)); 252c4762a1bSJed Brown } 253c4762a1bSJed Brown if (k > 1) { /* k-forms are antisymmetric */ 254c4762a1bSJed Brown PetscReal rj, rl, *swapv, wswapv, diff; 255c4762a1bSJed Brown PetscInt l, m; 256c4762a1bSJed Brown 2579566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, 0, k)); 2589566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rj)); 259c4762a1bSJed Brown j = (PetscInt)rj; 260c4762a1bSJed Brown l = j; 261c4762a1bSJed Brown while (l == j) { 2629566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &rl)); 263c4762a1bSJed Brown l = (PetscInt)rl; 264c4762a1bSJed Brown } 2659566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rand, -1., 1.)); 2669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(N * k, &swapv)); 267c4762a1bSJed Brown for (m = 0; m < N * k; m++) swapv[m] = v[m]; 268c4762a1bSJed Brown for (m = 0; m < N; m++) { 269c4762a1bSJed Brown swapv[j * N + m] = v[l * N + m]; 270c4762a1bSJed Brown swapv[l * N + m] = v[j * N + m]; 271c4762a1bSJed Brown } 2729566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k, w, swapv, &wswapv)); 273c4762a1bSJed Brown diff = PetscAbsReal(wswapv + wv); 2741dca8a05SBarry Smith PetscCheck(diff <= PETSC_SMALL * (PetscAbsReal(wswapv) + PetscAbsReal(wv)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "antisymmetry check: components %" PetscInt_FMT " & %" PetscInt_FMT ", wswapv (%g) != -wv (%g)", j, l, (double)wswapv, (double)wv); 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(swapv)); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown for (j = 0; j <= k && j + k <= N; j++) { /* wedge product */ 278c4762a1bSJed Brown PetscInt Nj, Njk, l, JKj; 279c4762a1bSJed Brown PetscReal *u, *uWw, *uWwcheck, *uWwmat, *x, *xsplit, uWwx, uWwxcheck, diff, norm; 280c4762a1bSJed Brown PetscInt *split; 281c4762a1bSJed Brown 28263a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "wedge j = %" PetscInt_FMT ":\n", j)); 2839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2849566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j, &Nj)); 2859566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, j + k, &Njk)); 2869566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(Nj, &u, Njk, &uWw, N * (j + k), &x, N * (j + k), &xsplit)); 2879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(j + k, &split)); 2889566063dSJacob Faibussowitsch for (l = 0; l < Nj; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l])); 2899566063dSJacob Faibussowitsch for (l = 0; l < N * (j + k); l++) PetscCall(PetscRandomGetValueReal(rand, &x[l])); 2909566063dSJacob Faibussowitsch PetscCall(PetscDTAltVWedge(N, j, k, u, w, uWw)); 2919566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, j + k, uWw, x, &uWwx)); 292c4762a1bSJed Brown if (verbose) { 2939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "u:\n")); 2949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2959566063dSJacob Faibussowitsch if (Nj) PetscCall(PetscRealView(Nj, u, viewer)); 2969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "u wedge w:\n")); 2989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2999566063dSJacob Faibussowitsch if (Njk) PetscCall(PetscRealView(Njk, uWw, viewer)); 3009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "x:\n")); 3029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3039566063dSJacob Faibussowitsch if (N * (j + k) > 0) PetscCall(PetscRealView(N * (j + k), x, viewer)); 3049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "u wedge w(x): %g\n", (double)uWwx)); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown /* verify wedge formula */ 308c4762a1bSJed Brown uWwxcheck = 0.; 3099566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(j + k, j, &JKj)); 310c4762a1bSJed Brown for (l = 0; l < JKj; l++) { 311c4762a1bSJed Brown PetscBool isOdd; 312c4762a1bSJed Brown PetscReal ux, wx; 313c4762a1bSJed Brown PetscInt m, p; 314c4762a1bSJed Brown 3159566063dSJacob Faibussowitsch PetscCall(PetscDTEnumSplit(j + k, j, l, split, &isOdd)); 3169371c9d4SSatish Balay for (m = 0; m < j + k; m++) { 317ad540459SPierre Jolivet for (p = 0; p < N; p++) xsplit[m * N + p] = x[split[m] * N + p]; 3189371c9d4SSatish Balay } 3199566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, j, u, xsplit, &ux)); 3208e3a54c0SPierre Jolivet PetscCall(PetscDTAltVApply(N, k, w, PetscSafePointerPlusOffset(xsplit, j * N), &wx)); 321c4762a1bSJed Brown uWwxcheck += isOdd ? -(ux * wx) : (ux * wx); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown diff = PetscAbsReal(uWwx - uWwxcheck); 3241dca8a05SBarry Smith PetscCheck(diff <= 10. * PETSC_SMALL * (PetscAbsReal(uWwx) + PetscAbsReal(uWwxcheck)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge check: forms %" PetscInt_FMT " & %" PetscInt_FMT ", uWwxcheck (%g) != uWwx (%g)", j, k, (double)uWwxcheck, (double)uWwx); 3259566063dSJacob Faibussowitsch PetscCall(PetscFree(split)); 3269566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nk * Njk, &uWwmat, Njk, &uWwcheck)); 3279566063dSJacob Faibussowitsch PetscCall(PetscDTAltVWedgeMatrix(N, j, k, u, uWwmat)); 328c4762a1bSJed Brown if (verbose) { 3299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "(u wedge):\n")); 3309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3319566063dSJacob Faibussowitsch if ((Nk * Njk) > 0) PetscCall(PetscRealView(Nk * Njk, uWwmat, viewer)); 3329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown diff = 0.; 335c4762a1bSJed Brown norm = 0.; 336c4762a1bSJed Brown for (l = 0; l < Njk; l++) { 337c4762a1bSJed Brown PetscInt m; 338c4762a1bSJed Brown PetscReal sum = 0.; 339c4762a1bSJed Brown 340c4762a1bSJed Brown for (m = 0; m < Nk; m++) sum += uWwmat[l * Nk + m] * w[m]; 341c4762a1bSJed Brown uWwcheck[l] = sum; 342c4762a1bSJed Brown diff += PetscSqr(uWwcheck[l] - uWw[l]); 343c4762a1bSJed Brown norm += PetscSqr(uWwcheck[l]) + PetscSqr(uWw[l]); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown diff = PetscSqrtReal(diff); 346c4762a1bSJed Brown norm = PetscSqrtReal(norm); 3471dca8a05SBarry Smith PetscCheck(diff <= PETSC_SMALL * norm, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge matrix check: wedge matrix application does not match wedge direct application"); 3489566063dSJacob Faibussowitsch PetscCall(PetscFree2(uWwmat, uWwcheck)); 3499566063dSJacob Faibussowitsch PetscCall(PetscFree4(u, uWw, x, xsplit)); 3509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown for (M = PetscMax(1, k); M <= N; M++) { /* pullback */ 353c4762a1bSJed Brown PetscReal *L, *u, *x; 354c4762a1bSJed Brown PetscInt Mk, l; 355c4762a1bSJed Brown 3569566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(M, k, &Mk)); 3579566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(M * N, &L, Mk, &u, M * k, &x)); 3589566063dSJacob Faibussowitsch for (l = 0; l < M * N; l++) PetscCall(PetscRandomGetValueReal(rand, &L[l])); 3599566063dSJacob Faibussowitsch for (l = 0; l < Mk; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l])); 3609566063dSJacob Faibussowitsch for (l = 0; l < M * k; l++) PetscCall(PetscRandomGetValueReal(rand, &x[l])); 36163a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "pullback M = %" PetscInt_FMT ":\n", M)); 3629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3639566063dSJacob Faibussowitsch PetscCall(CheckPullback(M, N, L, k, w, x, verbose, viewer)); 3649566063dSJacob Faibussowitsch if (M != N) PetscCall(CheckPullback(N, M, L, k, u, v, PETSC_FALSE, viewer)); 3659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 366c4762a1bSJed Brown if ((k % N) && (N > 1)) { 36763a3b9bcSJacob Faibussowitsch if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "negative pullback M = %" PetscInt_FMT ":\n", M)); 3689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3699566063dSJacob Faibussowitsch PetscCall(CheckPullback(M, N, L, -k, w, x, verbose, viewer)); 3709566063dSJacob Faibussowitsch if (M != N) PetscCall(CheckPullback(N, M, L, -k, u, v, PETSC_FALSE, viewer)); 3719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 372c4762a1bSJed Brown } 3739566063dSJacob Faibussowitsch PetscCall(PetscFree3(L, u, x)); 374c4762a1bSJed Brown } 375c4762a1bSJed Brown if (k > 0) { /* Interior */ 376c4762a1bSJed Brown PetscInt Nkm, l, m; 377c4762a1bSJed Brown PetscReal *wIntv0, *wIntv0check, wvcheck, diff, diffMat, normMat; 378c4762a1bSJed Brown PetscReal *intv0mat, *matcheck; 379c4762a1bSJed Brown PetscInt(*indices)[3]; 380c4762a1bSJed Brown 3819566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm)); 3829566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(Nkm, &wIntv0, Nkm, &wIntv0check, Nk * Nkm, &intv0mat, Nk * Nkm, &matcheck, Nk * k, &indices)); 3839566063dSJacob Faibussowitsch PetscCall(PetscDTAltVInterior(N, k, w, v, wIntv0)); 3849566063dSJacob Faibussowitsch PetscCall(PetscDTAltVInteriorMatrix(N, k, v, intv0mat)); 3859566063dSJacob Faibussowitsch PetscCall(PetscDTAltVInteriorPattern(N, k, indices)); 386c4762a1bSJed Brown if (verbose) { 3879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "interior product matrix pattern:\n")); 3889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 389c4762a1bSJed Brown for (l = 0; l < Nk * k; l++) { 390c4762a1bSJed Brown PetscInt row = indices[l][0]; 391c4762a1bSJed Brown PetscInt col = indices[l][1]; 392c4762a1bSJed Brown PetscInt x = indices[l][2]; 393c4762a1bSJed Brown 39463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "intV[%" PetscInt_FMT ",%" PetscInt_FMT "] = %sV[%" PetscInt_FMT "]\n", row, col, x < 0 ? "-" : " ", x < 0 ? -(x + 1) : x)); 395c4762a1bSJed Brown } 3969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 397c4762a1bSJed Brown } 398c4762a1bSJed Brown for (l = 0; l < Nkm * Nk; l++) matcheck[l] = 0.; 399c4762a1bSJed Brown for (l = 0; l < Nk * k; l++) { 400c4762a1bSJed Brown PetscInt row = indices[l][0]; 401c4762a1bSJed Brown PetscInt col = indices[l][1]; 402c4762a1bSJed Brown PetscInt x = indices[l][2]; 403c4762a1bSJed Brown 404c4762a1bSJed Brown if (x < 0) { 405c4762a1bSJed Brown matcheck[row * Nk + col] = -v[-(x + 1)]; 406c4762a1bSJed Brown } else { 407c4762a1bSJed Brown matcheck[row * Nk + col] = v[x]; 408c4762a1bSJed Brown } 409c4762a1bSJed Brown } 410c4762a1bSJed Brown diffMat = 0.; 411c4762a1bSJed Brown normMat = 0.; 412c4762a1bSJed Brown for (l = 0; l < Nkm * Nk; l++) { 413c4762a1bSJed Brown diffMat += PetscSqr(PetscAbsReal(matcheck[l] - intv0mat[l])); 414c4762a1bSJed Brown normMat += PetscSqr(matcheck[l]) + PetscSqr(intv0mat[l]); 415c4762a1bSJed Brown } 416c4762a1bSJed Brown diffMat = PetscSqrtReal(diffMat); 417c4762a1bSJed Brown normMat = PetscSqrtReal(normMat); 4181dca8a05SBarry Smith PetscCheck(diffMat <= PETSC_SMALL * normMat, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: matrix pattern does not match matrix"); 419c4762a1bSJed Brown diffMat = 0.; 420c4762a1bSJed Brown normMat = 0.; 421c4762a1bSJed Brown for (l = 0; l < Nkm; l++) { 422c4762a1bSJed Brown PetscReal sum = 0.; 423c4762a1bSJed Brown 424c4762a1bSJed Brown for (m = 0; m < Nk; m++) sum += intv0mat[l * Nk + m] * w[m]; 425c4762a1bSJed Brown wIntv0check[l] = sum; 426c4762a1bSJed Brown 427c4762a1bSJed Brown diffMat += PetscSqr(PetscAbsReal(wIntv0check[l] - wIntv0[l])); 428c4762a1bSJed Brown normMat += PetscSqr(wIntv0check[l]) + PetscSqr(wIntv0[l]); 429c4762a1bSJed Brown } 430c4762a1bSJed Brown diffMat = PetscSqrtReal(diffMat); 431c4762a1bSJed Brown normMat = PetscSqrtReal(normMat); 4321dca8a05SBarry Smith PetscCheck(diffMat <= PETSC_SMALL * normMat, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: application does not match matrix"); 433c4762a1bSJed Brown if (verbose) { 4349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "(w int v_0):\n")); 4359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 4369566063dSJacob Faibussowitsch if (Nkm) PetscCall(PetscRealView(Nkm, wIntv0, viewer)); 4379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 438c4762a1bSJed Brown 4399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "(int v_0):\n")); 4409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 4419566063dSJacob Faibussowitsch if (Nk * Nkm > 0) PetscCall(PetscRealView(Nk * Nkm, intv0mat, viewer)); 4429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 443c4762a1bSJed Brown } 4449566063dSJacob Faibussowitsch PetscCall(PetscDTAltVApply(N, k - 1, wIntv0, &v[N], &wvcheck)); 445c4762a1bSJed Brown diff = PetscSqrtReal(PetscSqr(wvcheck - wv)); 4461dca8a05SBarry Smith PetscCheck(diff < PETSC_SMALL * (PetscAbsReal(wv) + PetscAbsReal(wvcheck)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: (w Int v0)(v_rem) (%g) != w(v) (%g)", (double)wvcheck, (double)wv); 4479566063dSJacob Faibussowitsch PetscCall(PetscFree5(wIntv0, wIntv0check, intv0mat, matcheck, indices)); 448c4762a1bSJed Brown } 449c4762a1bSJed Brown if (k >= N - k) { /* Hodge star */ 450c4762a1bSJed Brown PetscReal *u, *starw, *starstarw, wu, starwdotu; 451c4762a1bSJed Brown PetscReal diff, norm; 452c4762a1bSJed Brown PetscBool isOdd; 453c4762a1bSJed Brown PetscInt l; 454c4762a1bSJed Brown 455c4762a1bSJed Brown isOdd = (PetscBool)((k * (N - k)) & 1); 4569566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Nk, &u, Nk, &starw, Nk, &starstarw)); 4579566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(N, k, 1, w, starw)); 4589566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(N, N - k, 1, starw, starstarw)); 459c4762a1bSJed Brown if (verbose) { 4609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "star w:\n")); 4619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 4629566063dSJacob Faibussowitsch if (Nk) PetscCall(PetscRealView(Nk, starw, viewer)); 4639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 464c4762a1bSJed Brown 4659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "star star w:\n")); 4669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 4679566063dSJacob Faibussowitsch if (Nk) PetscCall(PetscRealView(Nk, starstarw, viewer)); 4689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 469c4762a1bSJed Brown } 4709566063dSJacob Faibussowitsch for (l = 0; l < Nk; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l])); 4719566063dSJacob Faibussowitsch PetscCall(PetscDTAltVWedge(N, k, N - k, w, u, &wu)); 472c4762a1bSJed Brown starwdotu = 0.; 473c4762a1bSJed Brown for (l = 0; l < Nk; l++) starwdotu += starw[l] * u[l]; 474c4762a1bSJed Brown diff = PetscAbsReal(wu - starwdotu); 4751dca8a05SBarry Smith PetscCheck(diff <= PETSC_SMALL * (PetscAbsReal(wu) + PetscAbsReal(starwdotu)), PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: (star w, u) (%g) != (w wedge u) (%g)", (double)starwdotu, (double)wu); 476c4762a1bSJed Brown 477c4762a1bSJed Brown diff = 0.; 478c4762a1bSJed Brown norm = 0.; 479c4762a1bSJed Brown for (l = 0; l < Nk; l++) { 480c4762a1bSJed Brown diff += PetscSqr(w[l] - (isOdd ? -starstarw[l] : starstarw[l])); 481c4762a1bSJed Brown norm += PetscSqr(w[l]) + PetscSqr(starstarw[l]); 482c4762a1bSJed Brown } 483c4762a1bSJed Brown diff = PetscSqrtReal(diff); 484c4762a1bSJed Brown norm = PetscSqrtReal(norm); 4851dca8a05SBarry Smith PetscCheck(diff <= PETSC_SMALL * norm, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: star(star(w)) != (-1)^(N*(N-k)) w"); 4869566063dSJacob Faibussowitsch PetscCall(PetscFree3(u, starw, starstarw)); 487c4762a1bSJed Brown } 4889566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 4899566063dSJacob Faibussowitsch PetscCall(PetscFree(w)); 4909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 491c4762a1bSJed Brown } 4929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 493c4762a1bSJed Brown } 4949566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 4959566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 496b122ec5aSJacob Faibussowitsch return 0; 497c4762a1bSJed Brown } 498c4762a1bSJed Brown 499c4762a1bSJed Brown /*TEST 500c4762a1bSJed Brown test: 501c4762a1bSJed Brown suffix: 1234 502c4762a1bSJed Brown args: -verbose 503c4762a1bSJed Brown test: 504c4762a1bSJed Brown suffix: 56 505c4762a1bSJed Brown args: -N 5,6 506c4762a1bSJed Brown TEST*/ 507