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 6c4762a1bSJed Brown static PetscErrorCode CheckPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, 7c4762a1bSJed Brown const PetscReal *w, PetscReal *x, PetscBool verbose, PetscViewer viewer) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscInt Nk, Mk, i, j, l; 10c4762a1bSJed Brown PetscReal *Lstarw, *Lx, *Lstar, *Lstarwcheck, wLx, Lstarwx; 11c4762a1bSJed Brown PetscReal diff, diffMat, normMat; 12c4762a1bSJed Brown PetscReal *walloc = NULL; 13c4762a1bSJed Brown const PetscReal *ww = NULL; 14c4762a1bSJed Brown PetscBool negative = (PetscBool) (k < 0); 15c4762a1bSJed Brown 16c4762a1bSJed Brown PetscFunctionBegin; 17c4762a1bSJed Brown k = PetscAbsInt(k); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(N, k, &Nk)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(M, k, &Mk)); 20c4762a1bSJed Brown if (negative) { 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Mk, &walloc)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVStar(M, M - k, 1, w, walloc)); 23c4762a1bSJed Brown ww = walloc; 24c4762a1bSJed Brown } else { 25c4762a1bSJed Brown ww = w; 26c4762a1bSJed Brown } 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nk, &Lstarw, (M*k), &Lx)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nk * Mk, &Lstar, Nk, &Lstarwcheck)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVPullback(N, M, L, negative ? -k : k, w, Lstarw)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVPullbackMatrix(N, M, L, negative ? -k : k, Lstar)); 31c4762a1bSJed Brown if (negative) { 32c4762a1bSJed Brown PetscReal *sLsw; 33c4762a1bSJed Brown 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nk, &sLsw)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVStar(N, N - k, 1, Lstarw, sLsw)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, sLsw, x, &Lstarwx)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(sLsw)); 38c4762a1bSJed Brown } else { 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, Lstarw, x, &Lstarwx)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown for (l = 0; l < k; l++) { 42c4762a1bSJed Brown for (i = 0; i < M; i++) { 43c4762a1bSJed Brown PetscReal sum = 0.; 44c4762a1bSJed Brown 45c4762a1bSJed Brown for (j = 0; j < N; j++) sum += L[i * N + j] * x[l * N + j]; 46c4762a1bSJed Brown Lx[l * M + i] = sum; 47c4762a1bSJed Brown } 48c4762a1bSJed Brown } 49c4762a1bSJed Brown diffMat = 0.; 50c4762a1bSJed Brown normMat = 0.; 51c4762a1bSJed Brown for (i = 0; i < Nk; i++) { 52c4762a1bSJed Brown PetscReal sum = 0.; 53c4762a1bSJed Brown for (j = 0; j < Mk; j++) { 54c4762a1bSJed Brown sum += Lstar[i * Mk + j] * w[j]; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown Lstarwcheck[i] = sum; 57c4762a1bSJed Brown diffMat += PetscSqr(PetscAbsReal(Lstarwcheck[i] - Lstarw[i])); 58c4762a1bSJed Brown normMat += PetscSqr(Lstarwcheck[i]) + PetscSqr(Lstarw[i]); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown diffMat = PetscSqrtReal(diffMat); 61c4762a1bSJed Brown normMat = PetscSqrtReal(normMat); 62c4762a1bSJed Brown if (verbose) { 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "L:\n")); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 655f80ce2aSJacob Faibussowitsch if (M*N > 0) CHKERRQ(PetscRealView(M*N, L, viewer)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "L*:\n")); 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 705f80ce2aSJacob Faibussowitsch if (Nk * Mk > 0) CHKERRQ(PetscRealView(Nk * Mk, Lstar, viewer)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "L*w:\n")); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 755f80ce2aSJacob Faibussowitsch if (Nk > 0) CHKERRQ(PetscRealView(Nk, Lstarw, viewer)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 77c4762a1bSJed Brown } 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(M, k, ww, Lx, &wLx)); 79c4762a1bSJed Brown diff = PetscAbsReal(wLx - Lstarwx); 802c71b3e2SJacob Faibussowitsch PetscCheckFalse(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)", wLx, Lstarwx); 812c71b3e2SJacob Faibussowitsch PetscCheckFalse(diffMat > PETSC_SMALL * normMat,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback matrix does match matrix free result"); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(Lstar, Lstarwcheck)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(Lstarw, Lx)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(walloc)); 85c4762a1bSJed Brown PetscFunctionReturn(0); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown int main(int argc, char **argv) 89c4762a1bSJed Brown { 90c4762a1bSJed Brown PetscInt i, numTests = 5, n[5] = {0, 1, 2, 3, 4}; 91c4762a1bSJed Brown PetscBool verbose = PETSC_FALSE; 92c4762a1bSJed Brown PetscRandom rand; 93c4762a1bSJed Brown PetscViewer viewer; 94c4762a1bSJed Brown PetscErrorCode ierr; 95c4762a1bSJed Brown 96*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 97c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for exterior algebra tests","none");CHKERRQ(ierr); 985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-N", "Up to 5 vector space dimensions to test","ex7.c",n,&numTests,NULL)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-verbose", "Verbose test output","ex7.c",verbose,&verbose,NULL)); 1001e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rand)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand, -1., 1.)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 104c4762a1bSJed Brown if (!numTests) numTests = 5; 105c4762a1bSJed Brown viewer = PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD); 106c4762a1bSJed Brown for (i = 0; i < numTests; i++) { 107c4762a1bSJed Brown PetscInt k, N = n[i]; 108c4762a1bSJed Brown 1095f80ce2aSJacob Faibussowitsch if (verbose) CHKERRQ(PetscViewerASCIIPrintf(viewer, "N = %D:\n", N)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown if (verbose) { 113c4762a1bSJed Brown PetscInt *perm; 114c4762a1bSJed Brown PetscInt fac = 1; 115c4762a1bSJed Brown 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(N, &perm)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown for (k = 1; k <= N; k++) fac *= k; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "Permutations of %D:\n", N)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 121c4762a1bSJed Brown for (k = 0; k < fac; k++) { 122c4762a1bSJed Brown PetscBool isOdd, isOddCheck; 123c4762a1bSJed Brown PetscInt j, kCheck; 124c4762a1bSJed Brown 1255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTEnumPerm(N, k, perm, &isOdd)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "%D:", k)); 127c4762a1bSJed Brown for (j = 0; j < N; j++) { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, " %D", perm[j])); 129c4762a1bSJed Brown } 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even")); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTPermIndex(N, perm, &kCheck, &isOddCheck)); 1322c71b3e2SJacob Faibussowitsch PetscCheckFalse(kCheck != k || isOddCheck != isOdd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumPerm / PetscDTPermIndex mismatch for (%D, %D)", N, k); 133c4762a1bSJed Brown } 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(perm)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown for (k = 0; k <= N; k++) { 138c4762a1bSJed Brown PetscInt j, Nk, M; 139c4762a1bSJed Brown PetscReal *w, *v, wv; 140c4762a1bSJed Brown PetscInt *subset; 141c4762a1bSJed Brown 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(N, k, &Nk)); 1435f80ce2aSJacob Faibussowitsch if (verbose) CHKERRQ(PetscViewerASCIIPrintf(viewer, "k = %D:\n", k)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 1455f80ce2aSJacob Faibussowitsch if (verbose) CHKERRQ(PetscViewerASCIIPrintf(viewer, "(%D choose %D): %D\n", N, k, Nk)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Test subset and complement enumeration */ 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(N, &subset)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 150c4762a1bSJed Brown for (j = 0; j < Nk; j++) { 151c4762a1bSJed Brown PetscBool isOdd, isOddCheck; 152c4762a1bSJed Brown PetscInt jCheck, kCheck; 153c4762a1bSJed Brown 1545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTEnumSplit(N, k, j, subset, &isOdd)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTPermIndex(N, subset, &kCheck, &isOddCheck)); 1562c71b3e2SJacob Faibussowitsch PetscCheckFalse(isOddCheck != isOdd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumSplit sign does not mmatch PetscDTPermIndex sign"); 157c4762a1bSJed Brown if (verbose) { 158c4762a1bSJed Brown PetscInt l; 159c4762a1bSJed Brown 1605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "subset %D:", j)); 161c4762a1bSJed Brown for (l = 0; l < k; l++) { 1625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, " %D", subset[l])); 163c4762a1bSJed Brown } 1645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, " |")); 165c4762a1bSJed Brown for (l = k; l < N; l++) { 1665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, " %D", subset[l])); 167c4762a1bSJed Brown } 1685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even")); 169c4762a1bSJed Brown } 1705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTSubsetIndex(N, k, subset, &jCheck)); 1712c71b3e2SJacob Faibussowitsch PetscCheckFalse(jCheck != j,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "jCheck (%D) != j (%D)", jCheck, j); 172c4762a1bSJed Brown } 1735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(subset)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* Make a random k form */ 1775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nk, &w)); 1785f80ce2aSJacob Faibussowitsch for (j = 0; j < Nk; j++) CHKERRQ(PetscRandomGetValueReal(rand, &w[j])); 179c4762a1bSJed Brown /* Make a set of random vectors */ 1805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(N*k, &v)); 1815f80ce2aSJacob Faibussowitsch for (j = 0; j < N*k; j++) CHKERRQ(PetscRandomGetValueReal(rand, &v[j])); 182c4762a1bSJed Brown 1835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, w, v, &wv)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown if (verbose) { 1865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "w:\n")); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 1885f80ce2aSJacob Faibussowitsch if (Nk) CHKERRQ(PetscRealView(Nk, w, viewer)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "v:\n")); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 1925f80ce2aSJacob Faibussowitsch if (N*k > 0) CHKERRQ(PetscRealView(N*k, v, viewer)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "w(v): %g\n", (double) wv)); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* sanity checks */ 198c4762a1bSJed Brown if (k == 1) { /* 1-forms are functionals (dot products) */ 199c4762a1bSJed Brown PetscInt l; 200c4762a1bSJed Brown PetscReal wvcheck = 0.; 201c4762a1bSJed Brown PetscReal diff; 202c4762a1bSJed Brown 203c4762a1bSJed Brown for (l = 0; l < N; l++) wvcheck += w[l] * v[l]; 204c4762a1bSJed Brown diff = PetscSqrtReal(PetscSqr(wvcheck - wv)); 2052c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown if (k == N && N < 5) { /* n-forms are scaled determinants */ 208c4762a1bSJed Brown PetscReal det, wvcheck, diff; 209c4762a1bSJed Brown 210c4762a1bSJed Brown switch (k) { 211c4762a1bSJed Brown case 0: 212c4762a1bSJed Brown det = 1.; 213c4762a1bSJed Brown break; 214c4762a1bSJed Brown case 1: 215c4762a1bSJed Brown det = v[0]; 216c4762a1bSJed Brown break; 217c4762a1bSJed Brown case 2: 218c4762a1bSJed Brown det = v[0] * v[3] - v[1] * v[2]; 219c4762a1bSJed Brown break; 220c4762a1bSJed Brown case 3: 221c4762a1bSJed Brown det = v[0] * (v[4] * v[8] - v[5] * v[7]) + 222c4762a1bSJed Brown v[1] * (v[5] * v[6] - v[3] * v[8]) + 223c4762a1bSJed Brown v[2] * (v[3] * v[7] - v[4] * v[6]); 224c4762a1bSJed Brown break; 225c4762a1bSJed Brown case 4: 226c4762a1bSJed Brown det = v[0] * (v[5] * (v[10] * v[15] - v[11] * v[14]) + 227c4762a1bSJed Brown v[6] * (v[11] * v[13] - v[ 9] * v[15]) + 228c4762a1bSJed Brown v[7] * (v[ 9] * v[14] - v[10] * v[13])) - 229c4762a1bSJed Brown v[1] * (v[4] * (v[10] * v[15] - v[11] * v[14]) + 230c4762a1bSJed Brown v[6] * (v[11] * v[12] - v[ 8] * v[15]) + 231c4762a1bSJed Brown v[7] * (v[ 8] * v[14] - v[10] * v[12])) + 232c4762a1bSJed Brown v[2] * (v[4] * (v[ 9] * v[15] - v[11] * v[13]) + 233c4762a1bSJed Brown v[5] * (v[11] * v[12] - v[ 8] * v[15]) + 234c4762a1bSJed Brown v[7] * (v[ 8] * v[13] - v[ 9] * v[12])) - 235c4762a1bSJed Brown v[3] * (v[4] * (v[ 9] * v[14] - v[10] * v[13]) + 236c4762a1bSJed Brown v[5] * (v[10] * v[12] - v[ 8] * v[14]) + 237c4762a1bSJed Brown v[6] * (v[ 8] * v[13] - v[ 9] * v[12])); 238c4762a1bSJed Brown break; 239c4762a1bSJed Brown default: 240c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "invalid k"); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown wvcheck = det * w[0]; 243c4762a1bSJed Brown diff = PetscSqrtReal(PetscSqr(wvcheck - wv)); 2442c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown if (k > 0) { /* k-forms are linear in each component */ 247c4762a1bSJed Brown PetscReal alpha; 248c4762a1bSJed Brown PetscReal *x, *axv, wx, waxv, waxvcheck; 249c4762a1bSJed Brown PetscReal diff; 250c4762a1bSJed Brown PetscReal rj; 251c4762a1bSJed Brown PetscInt l; 252c4762a1bSJed Brown 2535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(N * k, &x, N * k, &axv)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand, &alpha)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand, 0, k)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand, &rj)); 257c4762a1bSJed Brown j = (PetscInt) rj; 2585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand, -1., 1.)); 259c4762a1bSJed Brown for (l = 0; l < N*k; l++) x[l] = v[l]; 260c4762a1bSJed Brown for (l = 0; l < N*k; l++) axv[l] = v[l]; 261c4762a1bSJed Brown for (l = 0; l < N; l++) { 262c4762a1bSJed Brown PetscReal val; 263c4762a1bSJed Brown 2645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand, &val)); 265c4762a1bSJed Brown x[j * N + l] = val; 266c4762a1bSJed Brown axv[j * N + l] += alpha * val; 267c4762a1bSJed Brown } 2685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, w, x, &wx)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, w, axv, &waxv)); 270c4762a1bSJed Brown waxvcheck = alpha * wx + wv; 271c4762a1bSJed Brown diff = waxv - waxvcheck; 2722c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(diff) > 10. * PETSC_SMALL * (PetscAbsReal(waxv) + PetscAbsReal(waxvcheck)),PETSC_COMM_WORLD, PETSC_ERR_PLIB, "linearity check: component %D, waxvcheck (%g) != waxv (%g)", j, (double) waxvcheck, (double) waxv); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(x,axv)); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown if (k > 1) { /* k-forms are antisymmetric */ 276c4762a1bSJed Brown PetscReal rj, rl, *swapv, wswapv, diff; 277c4762a1bSJed Brown PetscInt l, m; 278c4762a1bSJed Brown 2795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand, 0, k)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand, &rj)); 281c4762a1bSJed Brown j = (PetscInt) rj; 282c4762a1bSJed Brown l = j; 283c4762a1bSJed Brown while (l == j) { 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValueReal(rand, &rl)); 285c4762a1bSJed Brown l = (PetscInt) rl; 286c4762a1bSJed Brown } 2875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetInterval(rand, -1., 1.)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(N * k, &swapv)); 289c4762a1bSJed Brown for (m = 0; m < N * k; m++) swapv[m] = v[m]; 290c4762a1bSJed Brown for (m = 0; m < N; m++) { 291c4762a1bSJed Brown swapv[j * N + m] = v[l * N + m]; 292c4762a1bSJed Brown swapv[l * N + m] = v[j * N + m]; 293c4762a1bSJed Brown } 2945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, w, swapv, &wswapv)); 295c4762a1bSJed Brown diff = PetscAbsReal(wswapv + wv); 2962c71b3e2SJacob Faibussowitsch PetscCheckFalse(diff > PETSC_SMALL * (PetscAbsReal(wswapv) + PetscAbsReal(wv)),PETSC_COMM_WORLD, PETSC_ERR_PLIB, "antisymmetry check: components %D & %D, wswapv (%g) != -wv (%g)", j, l, (double) wswapv, (double) wv); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(swapv)); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown for (j = 0; j <= k && j + k <= N; j++) { /* wedge product */ 300c4762a1bSJed Brown PetscInt Nj, Njk, l, JKj; 301c4762a1bSJed Brown PetscReal *u, *uWw, *uWwcheck, *uWwmat, *x, *xsplit, uWwx, uWwxcheck, diff, norm; 302c4762a1bSJed Brown PetscInt *split; 303c4762a1bSJed Brown 3045f80ce2aSJacob Faibussowitsch if (verbose) CHKERRQ(PetscViewerASCIIPrintf(viewer, "wedge j = %D:\n", j)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(N, j, &Nj)); 3075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(N, j+k, &Njk)); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc4(Nj, &u, Njk, &uWw, N*(j+k), &x, N*(j+k), &xsplit)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(j+k,&split)); 3105f80ce2aSJacob Faibussowitsch for (l = 0; l < Nj; l++) CHKERRQ(PetscRandomGetValueReal(rand, &u[l])); 3115f80ce2aSJacob Faibussowitsch for (l = 0; l < N*(j+k); l++) CHKERRQ(PetscRandomGetValueReal(rand, &x[l])); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVWedge(N, j, k, u, w, uWw)); 3135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, j+k, uWw, x, &uWwx)); 314c4762a1bSJed Brown if (verbose) { 3155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "u:\n")); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3175f80ce2aSJacob Faibussowitsch if (Nj) CHKERRQ(PetscRealView(Nj, u, viewer)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "u wedge w:\n")); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3215f80ce2aSJacob Faibussowitsch if (Njk) CHKERRQ(PetscRealView(Njk, uWw, viewer)); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "x:\n")); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3255f80ce2aSJacob Faibussowitsch if (N*(j+k) > 0) CHKERRQ(PetscRealView(N*(j+k), x, viewer)); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "u wedge w(x): %g\n", (double) uWwx)); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown /* verify wedge formula */ 330c4762a1bSJed Brown uWwxcheck = 0.; 3315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(j+k, j, &JKj)); 332c4762a1bSJed Brown for (l = 0; l < JKj; l++) { 333c4762a1bSJed Brown PetscBool isOdd; 334c4762a1bSJed Brown PetscReal ux, wx; 335c4762a1bSJed Brown PetscInt m, p; 336c4762a1bSJed Brown 3375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTEnumSplit(j+k, j, l, split, &isOdd)); 338c4762a1bSJed Brown for (m = 0; m < j+k; m++) {for (p = 0; p < N; p++) {xsplit[m * N + p] = x[split[m] * N + p];}} 3395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, j, u, xsplit, &ux)); 3405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k, w, &xsplit[j*N], &wx)); 341c4762a1bSJed Brown uWwxcheck += isOdd ? -(ux * wx) : (ux * wx); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown diff = PetscAbsReal(uWwx - uWwxcheck); 3442c71b3e2SJacob Faibussowitsch PetscCheckFalse(diff > 10. * PETSC_SMALL * (PetscAbsReal(uWwx) + PetscAbsReal(uWwxcheck)),PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge check: forms %D & %D, uWwxcheck (%g) != uWwx (%g)", j, k, (double) uWwxcheck, (double) uWwx); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(split)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(Nk * Njk, &uWwmat, Njk, &uWwcheck)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVWedgeMatrix(N, j, k, u, uWwmat)); 348c4762a1bSJed Brown if (verbose) { 3495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "(u wedge):\n")); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3515f80ce2aSJacob Faibussowitsch if ((Nk * Njk) > 0) CHKERRQ(PetscRealView(Nk * Njk, uWwmat, viewer)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown diff = 0.; 355c4762a1bSJed Brown norm = 0.; 356c4762a1bSJed Brown for (l = 0; l < Njk; l++) { 357c4762a1bSJed Brown PetscInt m; 358c4762a1bSJed Brown PetscReal sum = 0.; 359c4762a1bSJed Brown 360c4762a1bSJed Brown for (m = 0; m < Nk; m++) sum += uWwmat[l * Nk + m] * w[m]; 361c4762a1bSJed Brown uWwcheck[l] = sum; 362c4762a1bSJed Brown diff += PetscSqr(uWwcheck[l] - uWw[l]); 363c4762a1bSJed Brown norm += PetscSqr(uWwcheck[l]) + PetscSqr(uWw[l]); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown diff = PetscSqrtReal(diff); 366c4762a1bSJed Brown norm = PetscSqrtReal(norm); 3672c71b3e2SJacob Faibussowitsch PetscCheckFalse(diff > PETSC_SMALL * norm,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge matrix check: wedge matrix application does not match wedge direct application"); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(uWwmat, uWwcheck)); 3695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree4(u, uWw, x, xsplit)); 3705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 371c4762a1bSJed Brown } 372c4762a1bSJed Brown for (M = PetscMax(1,k); M <= N; M++) { /* pullback */ 373c4762a1bSJed Brown PetscReal *L, *u, *x; 374c4762a1bSJed Brown PetscInt Mk, l; 375c4762a1bSJed Brown 3765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(M, k, &Mk)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(M*N, &L, Mk, &u, M*k, &x)); 3785f80ce2aSJacob Faibussowitsch for (l = 0; l < M*N; l++) CHKERRQ(PetscRandomGetValueReal(rand, &L[l])); 3795f80ce2aSJacob Faibussowitsch for (l = 0; l < Mk; l++) CHKERRQ(PetscRandomGetValueReal(rand, &u[l])); 3805f80ce2aSJacob Faibussowitsch for (l = 0; l < M*k; l++) CHKERRQ(PetscRandomGetValueReal(rand, &x[l])); 3815f80ce2aSJacob Faibussowitsch if (verbose) CHKERRQ(PetscViewerASCIIPrintf(viewer, "pullback M = %D:\n", M)); 3825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3835f80ce2aSJacob Faibussowitsch CHKERRQ(CheckPullback(M, N, L, k, w, x, verbose, viewer)); 3845f80ce2aSJacob Faibussowitsch if (M != N) CHKERRQ(CheckPullback(N, M, L, k, u, v, PETSC_FALSE, viewer)); 3855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 386c4762a1bSJed Brown if ((k % N) && (N > 1)) { 3875f80ce2aSJacob Faibussowitsch if (verbose) CHKERRQ(PetscViewerASCIIPrintf(viewer, "negative pullback M = %D:\n", M)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(CheckPullback(M, N, L, -k, w, x, verbose, viewer)); 3905f80ce2aSJacob Faibussowitsch if (M != N) CHKERRQ(CheckPullback(N, M, L, -k, u, v, PETSC_FALSE, viewer)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 392c4762a1bSJed Brown } 3935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(L, u, x)); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown if (k > 0) { /* Interior */ 396c4762a1bSJed Brown PetscInt Nkm, l, m; 397c4762a1bSJed Brown PetscReal *wIntv0, *wIntv0check, wvcheck, diff, diffMat, normMat; 398c4762a1bSJed Brown PetscReal *intv0mat, *matcheck; 399c4762a1bSJed Brown PetscInt (*indices)[3]; 400c4762a1bSJed Brown 4015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(N, k-1, &Nkm)); 4025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc5(Nkm, &wIntv0, Nkm, &wIntv0check, Nk * Nkm, &intv0mat, Nk * Nkm, &matcheck, Nk * k, &indices)); 4035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVInterior(N, k, w, v, wIntv0)); 4045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVInteriorMatrix(N, k, v, intv0mat)); 4055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVInteriorPattern(N, k, indices)); 406c4762a1bSJed Brown if (verbose) { 4075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "interior product matrix pattern:\n")); 4085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 409c4762a1bSJed Brown for (l = 0; l < Nk * k; l++) { 410c4762a1bSJed Brown PetscInt row = indices[l][0]; 411c4762a1bSJed Brown PetscInt col = indices[l][1]; 412c4762a1bSJed Brown PetscInt x = indices[l][2]; 413c4762a1bSJed Brown 4145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"intV[%D,%D] = %sV[%D]\n", row, col, x < 0 ? "-" : " ", x < 0 ? -(x + 1) : x)); 415c4762a1bSJed Brown } 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 417c4762a1bSJed Brown } 418c4762a1bSJed Brown for (l = 0; l < Nkm * Nk; l++) matcheck[l] = 0.; 419c4762a1bSJed Brown for (l = 0; l < Nk * k; l++) { 420c4762a1bSJed Brown PetscInt row = indices[l][0]; 421c4762a1bSJed Brown PetscInt col = indices[l][1]; 422c4762a1bSJed Brown PetscInt x = indices[l][2]; 423c4762a1bSJed Brown 424c4762a1bSJed Brown if (x < 0) { 425c4762a1bSJed Brown matcheck[row * Nk + col] = -v[-(x+1)]; 426c4762a1bSJed Brown } else { 427c4762a1bSJed Brown matcheck[row * Nk + col] = v[x]; 428c4762a1bSJed Brown } 429c4762a1bSJed Brown } 430c4762a1bSJed Brown diffMat = 0.; 431c4762a1bSJed Brown normMat = 0.; 432c4762a1bSJed Brown for (l = 0; l < Nkm * Nk; l++) { 433c4762a1bSJed Brown diffMat += PetscSqr(PetscAbsReal(matcheck[l] - intv0mat[l])); 434c4762a1bSJed Brown normMat += PetscSqr(matcheck[l]) + PetscSqr(intv0mat[l]); 435c4762a1bSJed Brown } 436c4762a1bSJed Brown diffMat = PetscSqrtReal(diffMat); 437c4762a1bSJed Brown normMat = PetscSqrtReal(normMat); 4382c71b3e2SJacob Faibussowitsch PetscCheckFalse(diffMat > PETSC_SMALL * normMat,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: matrix pattern does not match matrix"); 439c4762a1bSJed Brown diffMat = 0.; 440c4762a1bSJed Brown normMat = 0.; 441c4762a1bSJed Brown for (l = 0; l < Nkm; l++) { 442c4762a1bSJed Brown PetscReal sum = 0.; 443c4762a1bSJed Brown 444c4762a1bSJed Brown for (m = 0; m < Nk; m++) sum += intv0mat[l * Nk + m] * w[m]; 445c4762a1bSJed Brown wIntv0check[l] = sum; 446c4762a1bSJed Brown 447c4762a1bSJed Brown diffMat += PetscSqr(PetscAbsReal(wIntv0check[l] - wIntv0[l])); 448c4762a1bSJed Brown normMat += PetscSqr(wIntv0check[l]) + PetscSqr(wIntv0[l]); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown diffMat = PetscSqrtReal(diffMat); 451c4762a1bSJed Brown normMat = PetscSqrtReal(normMat); 4522c71b3e2SJacob Faibussowitsch PetscCheckFalse(diffMat > PETSC_SMALL * normMat,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: application does not match matrix"); 453c4762a1bSJed Brown if (verbose) { 4545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "(w int v_0):\n")); 4555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 4565f80ce2aSJacob Faibussowitsch if (Nkm) CHKERRQ(PetscRealView(Nkm, wIntv0, viewer)); 4575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 458c4762a1bSJed Brown 4595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "(int v_0):\n")); 4605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 4615f80ce2aSJacob Faibussowitsch if (Nk * Nkm > 0) CHKERRQ(PetscRealView(Nk * Nkm, intv0mat, viewer)); 4625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 463c4762a1bSJed Brown } 4645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVApply(N, k - 1, wIntv0, &v[N], &wvcheck)); 465c4762a1bSJed Brown diff = PetscSqrtReal(PetscSqr(wvcheck - wv)); 4662c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree5(wIntv0,wIntv0check,intv0mat,matcheck,indices)); 468c4762a1bSJed Brown } 469c4762a1bSJed Brown if (k >= N - k) { /* Hodge star */ 470c4762a1bSJed Brown PetscReal *u, *starw, *starstarw, wu, starwdotu; 471c4762a1bSJed Brown PetscReal diff, norm; 472c4762a1bSJed Brown PetscBool isOdd; 473c4762a1bSJed Brown PetscInt l; 474c4762a1bSJed Brown 475c4762a1bSJed Brown isOdd = (PetscBool) ((k * (N - k)) & 1); 4765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(Nk, &u, Nk, &starw, Nk, &starstarw)); 4775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVStar(N, k, 1, w, starw)); 4785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVStar(N, N-k, 1, starw, starstarw)); 479c4762a1bSJed Brown if (verbose) { 4805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "star w:\n")); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 4825f80ce2aSJacob Faibussowitsch if (Nk) CHKERRQ(PetscRealView(Nk, starw, viewer)); 4835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 484c4762a1bSJed Brown 4855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "star star w:\n")); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 4875f80ce2aSJacob Faibussowitsch if (Nk) CHKERRQ(PetscRealView(Nk, starstarw, viewer)); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 489c4762a1bSJed Brown } 4905f80ce2aSJacob Faibussowitsch for (l = 0; l < Nk; l++) CHKERRQ(PetscRandomGetValueReal(rand,&u[l])); 4915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVWedge(N, k, N - k, w, u, &wu)); 492c4762a1bSJed Brown starwdotu = 0.; 493c4762a1bSJed Brown for (l = 0; l < Nk; l++) starwdotu += starw[l] * u[l]; 494c4762a1bSJed Brown diff = PetscAbsReal(wu - starwdotu); 4952c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 496c4762a1bSJed Brown 497c4762a1bSJed Brown diff = 0.; 498c4762a1bSJed Brown norm = 0.; 499c4762a1bSJed Brown for (l = 0; l < Nk; l++) { 500c4762a1bSJed Brown diff += PetscSqr(w[l] - (isOdd ? -starstarw[l] : starstarw[l])); 501c4762a1bSJed Brown norm += PetscSqr(w[l]) + PetscSqr(starstarw[l]); 502c4762a1bSJed Brown } 503c4762a1bSJed Brown diff = PetscSqrtReal(diff); 504c4762a1bSJed Brown norm = PetscSqrtReal(norm); 5052c71b3e2SJacob Faibussowitsch PetscCheckFalse(diff > PETSC_SMALL * norm,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: star(star(w)) != (-1)^(N*(N-k)) w"); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(u, starw, starstarw)); 507c4762a1bSJed Brown } 5085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 5095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(w)); 5105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 511c4762a1bSJed Brown } 5125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 513c4762a1bSJed Brown } 5145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 515*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 516*b122ec5aSJacob Faibussowitsch return 0; 517c4762a1bSJed Brown } 518c4762a1bSJed Brown 519c4762a1bSJed Brown /*TEST 520c4762a1bSJed Brown test: 521c4762a1bSJed Brown suffix: 1234 522c4762a1bSJed Brown args: -verbose 523c4762a1bSJed Brown test: 524c4762a1bSJed Brown suffix: 56 525c4762a1bSJed Brown args: -N 5,6 526c4762a1bSJed Brown TEST*/ 527