xref: /petsc/src/dm/dt/tests/ex7.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
69371c9d4SSatish Balay static PetscErrorCode CheckPullback(PetscInt N, PetscInt M, const PetscReal *L, PetscInt k, const PetscReal *w, PetscReal *x, PetscBool verbose, PetscViewer viewer) {
7c4762a1bSJed Brown   PetscInt         Nk, Mk, i, j, l;
8c4762a1bSJed Brown   PetscReal       *Lstarw, *Lx, *Lstar, *Lstarwcheck, wLx, Lstarwx;
9c4762a1bSJed Brown   PetscReal        diff, diffMat, normMat;
10c4762a1bSJed Brown   PetscReal       *walloc   = NULL;
11c4762a1bSJed Brown   const PetscReal *ww       = NULL;
12c4762a1bSJed Brown   PetscBool        negative = (PetscBool)(k < 0);
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   PetscFunctionBegin;
15c4762a1bSJed Brown   k = PetscAbsInt(k);
169566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
179566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(M, k, &Mk));
18c4762a1bSJed Brown   if (negative) {
199566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Mk, &walloc));
209566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc));
21c4762a1bSJed Brown     ww = walloc;
22c4762a1bSJed Brown   } else {
23c4762a1bSJed Brown     ww = w;
24c4762a1bSJed Brown   }
259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nk, &Lstarw, (M * k), &Lx));
269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nk * Mk, &Lstar, Nk, &Lstarwcheck));
279566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullback(N, M, L, negative ? -k : k, w, Lstarw));
289566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(N, M, L, negative ? -k : k, Lstar));
29c4762a1bSJed Brown   if (negative) {
30c4762a1bSJed Brown     PetscReal *sLsw;
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nk, &sLsw));
339566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVStar(N, N - k, 1, Lstarw, sLsw));
349566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVApply(N, k, sLsw, x, &Lstarwx));
359566063dSJacob Faibussowitsch     PetscCall(PetscFree(sLsw));
36c4762a1bSJed Brown   } else {
379566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVApply(N, k, Lstarw, x, &Lstarwx));
38c4762a1bSJed Brown   }
39c4762a1bSJed Brown   for (l = 0; l < k; l++) {
40c4762a1bSJed Brown     for (i = 0; i < M; i++) {
41c4762a1bSJed Brown       PetscReal sum = 0.;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown       for (j = 0; j < N; j++) sum += L[i * N + j] * x[l * N + j];
44c4762a1bSJed Brown       Lx[l * M + i] = sum;
45c4762a1bSJed Brown     }
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown   diffMat = 0.;
48c4762a1bSJed Brown   normMat = 0.;
49c4762a1bSJed Brown   for (i = 0; i < Nk; i++) {
50c4762a1bSJed Brown     PetscReal sum = 0.;
519371c9d4SSatish Balay     for (j = 0; j < Mk; j++) { sum += Lstar[i * Mk + j] * w[j]; }
52c4762a1bSJed Brown     Lstarwcheck[i] = sum;
53c4762a1bSJed Brown     diffMat += PetscSqr(PetscAbsReal(Lstarwcheck[i] - Lstarw[i]));
54c4762a1bSJed Brown     normMat += PetscSqr(Lstarwcheck[i]) + PetscSqr(Lstarw[i]);
55c4762a1bSJed Brown   }
56c4762a1bSJed Brown   diffMat = PetscSqrtReal(diffMat);
57c4762a1bSJed Brown   normMat = PetscSqrtReal(normMat);
58c4762a1bSJed Brown   if (verbose) {
599566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L:\n"));
609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
619566063dSJacob Faibussowitsch     if (M * N > 0) PetscCall(PetscRealView(M * N, L, viewer));
629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L*:\n"));
659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
669566063dSJacob Faibussowitsch     if (Nk * Mk > 0) PetscCall(PetscRealView(Nk * Mk, Lstar, viewer));
679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L*w:\n"));
709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
719566063dSJacob Faibussowitsch     if (Nk > 0) PetscCall(PetscRealView(Nk, Lstarw, viewer));
729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
73c4762a1bSJed Brown   }
749566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVApply(M, k, ww, Lx, &wLx));
75c4762a1bSJed Brown   diff = PetscAbsReal(wLx - Lstarwx);
761dca8a05SBarry 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);
771dca8a05SBarry Smith   PetscCheck(diffMat <= PETSC_SMALL * normMat, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback matrix does match matrix free result");
789566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Lstar, Lstarwcheck));
799566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Lstarw, Lx));
809566063dSJacob Faibussowitsch   PetscCall(PetscFree(walloc));
81c4762a1bSJed Brown   PetscFunctionReturn(0);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
849371c9d4SSatish Balay int main(int argc, char **argv) {
85c4762a1bSJed Brown   PetscInt    i, numTests = 5, n[5] = {0, 1, 2, 3, 4};
86c4762a1bSJed Brown   PetscBool   verbose = PETSC_FALSE;
87c4762a1bSJed Brown   PetscRandom rand;
88c4762a1bSJed Brown   PetscViewer viewer;
89c4762a1bSJed Brown 
90327415f7SBarry Smith   PetscFunctionBeginUser;
919566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
92d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for exterior algebra tests", "none");
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-N", "Up to 5 vector space dimensions to test", "ex7.c", n, &numTests, NULL));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-verbose", "Verbose test output", "ex7.c", verbose, &verbose, NULL));
95d0609cedSBarry Smith   PetscOptionsEnd();
969566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
979566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand, -1., 1.));
989566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
99c4762a1bSJed Brown   if (!numTests) numTests = 5;
100c4762a1bSJed Brown   viewer = PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD);
101c4762a1bSJed Brown   for (i = 0; i < numTests; i++) {
102c4762a1bSJed Brown     PetscInt k, N = n[i];
103c4762a1bSJed Brown 
10463a3b9bcSJacob Faibussowitsch     if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "N = %" PetscInt_FMT ":\n", N));
1059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown     if (verbose) {
108c4762a1bSJed Brown       PetscInt *perm;
109c4762a1bSJed Brown       PetscInt  fac = 1;
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(N, &perm));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown       for (k = 1; k <= N; k++) fac *= k;
11463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Permutations of %" PetscInt_FMT ":\n", N));
1159566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
116c4762a1bSJed Brown       for (k = 0; k < fac; k++) {
117c4762a1bSJed Brown         PetscBool isOdd, isOddCheck;
118c4762a1bSJed Brown         PetscInt  j, kCheck;
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumPerm(N, k, perm, &isOdd));
12163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ":", k));
122*48a46eb9SPierre Jolivet         for (j = 0; j < N; j++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, perm[j]));
1239566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even"));
1249566063dSJacob Faibussowitsch         PetscCall(PetscDTPermIndex(N, perm, &kCheck, &isOddCheck));
1251dca8a05SBarry Smith         PetscCheck(kCheck == k && isOddCheck == isOdd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumPerm / PetscDTPermIndex mismatch for (%" PetscInt_FMT ", %" PetscInt_FMT ")", N, k);
126c4762a1bSJed Brown       }
1279566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
1289566063dSJacob Faibussowitsch       PetscCall(PetscFree(perm));
129c4762a1bSJed Brown     }
130c4762a1bSJed Brown     for (k = 0; k <= N; k++) {
131c4762a1bSJed Brown       PetscInt   j, Nk, M;
132c4762a1bSJed Brown       PetscReal *w, *v, wv;
133c4762a1bSJed Brown       PetscInt  *subset;
134c4762a1bSJed Brown 
1359566063dSJacob Faibussowitsch       PetscCall(PetscDTBinomialInt(N, k, &Nk));
13663a3b9bcSJacob Faibussowitsch       if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "k = %" PetscInt_FMT ":\n", k));
1379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
13863a3b9bcSJacob Faibussowitsch       if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT " choose %" PetscInt_FMT "): %" PetscInt_FMT "\n", N, k, Nk));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown       /* Test subset and complement enumeration */
1419566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(N, &subset));
1429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
143c4762a1bSJed Brown       for (j = 0; j < Nk; j++) {
144c4762a1bSJed Brown         PetscBool isOdd, isOddCheck;
145c4762a1bSJed Brown         PetscInt  jCheck, kCheck;
146c4762a1bSJed Brown 
1479566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, k, j, subset, &isOdd));
1489566063dSJacob Faibussowitsch         PetscCall(PetscDTPermIndex(N, subset, &kCheck, &isOddCheck));
14908401ef6SPierre Jolivet         PetscCheck(isOddCheck == isOdd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumSplit sign does not mmatch PetscDTPermIndex sign");
150c4762a1bSJed Brown         if (verbose) {
151c4762a1bSJed Brown           PetscInt l;
152c4762a1bSJed Brown 
15363a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "subset %" PetscInt_FMT ":", j));
154*48a46eb9SPierre Jolivet           for (l = 0; l < k; l++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, subset[l]));
1559566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, " |"));
156*48a46eb9SPierre Jolivet           for (l = k; l < N; l++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, subset[l]));
1579566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even"));
158c4762a1bSJed Brown         }
1599566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subset, &jCheck));
16063a3b9bcSJacob Faibussowitsch         PetscCheck(jCheck == j, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "jCheck (%" PetscInt_FMT ") != j (%" PetscInt_FMT ")", jCheck, j);
161c4762a1bSJed Brown       }
1629566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
1639566063dSJacob Faibussowitsch       PetscCall(PetscFree(subset));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown       /* Make a random k form */
1669566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &w));
1679566063dSJacob Faibussowitsch       for (j = 0; j < Nk; j++) PetscCall(PetscRandomGetValueReal(rand, &w[j]));
168c4762a1bSJed Brown       /* Make a set of random vectors */
1699566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(N * k, &v));
1709566063dSJacob Faibussowitsch       for (j = 0; j < N * k; j++) PetscCall(PetscRandomGetValueReal(rand, &v[j]));
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVApply(N, k, w, v, &wv));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown       if (verbose) {
1759566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "w:\n"));
1769566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
1779566063dSJacob Faibussowitsch         if (Nk) PetscCall(PetscRealView(Nk, w, viewer));
1789566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
1799566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "v:\n"));
1809566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
1819566063dSJacob Faibussowitsch         if (N * k > 0) PetscCall(PetscRealView(N * k, v, viewer));
1829566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
1839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "w(v): %g\n", (double)wv));
184c4762a1bSJed Brown       }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown       /* sanity checks */
187c4762a1bSJed Brown       if (k == 1) { /* 1-forms are functionals (dot products) */
188c4762a1bSJed Brown         PetscInt  l;
189c4762a1bSJed Brown         PetscReal wvcheck = 0.;
190c4762a1bSJed Brown         PetscReal diff;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown         for (l = 0; l < N; l++) wvcheck += w[l] * v[l];
193c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
1941dca8a05SBarry 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);
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown       if (k == N && N < 5) { /* n-forms are scaled determinants */
197c4762a1bSJed Brown         PetscReal det, wvcheck, diff;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown         switch (k) {
2009371c9d4SSatish Balay         case 0: det = 1.; break;
2019371c9d4SSatish Balay         case 1: det = v[0]; break;
2029371c9d4SSatish Balay         case 2: det = v[0] * v[3] - v[1] * v[2]; break;
2039371c9d4SSatish Balay         case 3: 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]); break;
204c4762a1bSJed Brown         case 4:
2059371c9d4SSatish 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]));
206c4762a1bSJed Brown           break;
2079371c9d4SSatish Balay         default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "invalid k");
208c4762a1bSJed Brown         }
209c4762a1bSJed Brown         wvcheck = det * w[0];
210c4762a1bSJed Brown         diff    = PetscSqrtReal(PetscSqr(wvcheck - wv));
2111dca8a05SBarry 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);
212c4762a1bSJed Brown       }
213c4762a1bSJed Brown       if (k > 0) { /* k-forms are linear in each component */
214c4762a1bSJed Brown         PetscReal  alpha;
215c4762a1bSJed Brown         PetscReal *x, *axv, wx, waxv, waxvcheck;
216c4762a1bSJed Brown         PetscReal  diff;
217c4762a1bSJed Brown         PetscReal  rj;
218c4762a1bSJed Brown         PetscInt   l;
219c4762a1bSJed Brown 
2209566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(N * k, &x, N * k, &axv));
2219566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rand, &alpha));
2229566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, 0, k));
2239566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rand, &rj));
224c4762a1bSJed Brown         j = (PetscInt)rj;
2259566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, -1., 1.));
226c4762a1bSJed Brown         for (l = 0; l < N * k; l++) x[l] = v[l];
227c4762a1bSJed Brown         for (l = 0; l < N * k; l++) axv[l] = v[l];
228c4762a1bSJed Brown         for (l = 0; l < N; l++) {
229c4762a1bSJed Brown           PetscReal val;
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rand, &val));
232c4762a1bSJed Brown           x[j * N + l] = val;
233c4762a1bSJed Brown           axv[j * N + l] += alpha * val;
234c4762a1bSJed Brown         }
2359566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k, w, x, &wx));
2369566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k, w, axv, &waxv));
237c4762a1bSJed Brown         waxvcheck = alpha * wx + wv;
238c4762a1bSJed Brown         diff      = waxv - waxvcheck;
2391dca8a05SBarry 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);
2409566063dSJacob Faibussowitsch         PetscCall(PetscFree2(x, axv));
241c4762a1bSJed Brown       }
242c4762a1bSJed Brown       if (k > 1) { /* k-forms are antisymmetric */
243c4762a1bSJed Brown         PetscReal rj, rl, *swapv, wswapv, diff;
244c4762a1bSJed Brown         PetscInt  l, m;
245c4762a1bSJed Brown 
2469566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, 0, k));
2479566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rand, &rj));
248c4762a1bSJed Brown         j = (PetscInt)rj;
249c4762a1bSJed Brown         l = j;
250c4762a1bSJed Brown         while (l == j) {
2519566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rand, &rl));
252c4762a1bSJed Brown           l = (PetscInt)rl;
253c4762a1bSJed Brown         }
2549566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, -1., 1.));
2559566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(N * k, &swapv));
256c4762a1bSJed Brown         for (m = 0; m < N * k; m++) swapv[m] = v[m];
257c4762a1bSJed Brown         for (m = 0; m < N; m++) {
258c4762a1bSJed Brown           swapv[j * N + m] = v[l * N + m];
259c4762a1bSJed Brown           swapv[l * N + m] = v[j * N + m];
260c4762a1bSJed Brown         }
2619566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k, w, swapv, &wswapv));
262c4762a1bSJed Brown         diff = PetscAbsReal(wswapv + wv);
2631dca8a05SBarry 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);
2649566063dSJacob Faibussowitsch         PetscCall(PetscFree(swapv));
265c4762a1bSJed Brown       }
266c4762a1bSJed Brown       for (j = 0; j <= k && j + k <= N; j++) { /* wedge product */
267c4762a1bSJed Brown         PetscInt   Nj, Njk, l, JKj;
268c4762a1bSJed Brown         PetscReal *u, *uWw, *uWwcheck, *uWwmat, *x, *xsplit, uWwx, uWwxcheck, diff, norm;
269c4762a1bSJed Brown         PetscInt  *split;
270c4762a1bSJed Brown 
27163a3b9bcSJacob Faibussowitsch         if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "wedge j = %" PetscInt_FMT ":\n", j));
2729566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
2739566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(N, j, &Nj));
2749566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(N, j + k, &Njk));
2759566063dSJacob Faibussowitsch         PetscCall(PetscMalloc4(Nj, &u, Njk, &uWw, N * (j + k), &x, N * (j + k), &xsplit));
2769566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(j + k, &split));
2779566063dSJacob Faibussowitsch         for (l = 0; l < Nj; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l]));
2789566063dSJacob Faibussowitsch         for (l = 0; l < N * (j + k); l++) PetscCall(PetscRandomGetValueReal(rand, &x[l]));
2799566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedge(N, j, k, u, w, uWw));
2809566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, j + k, uWw, x, &uWwx));
281c4762a1bSJed Brown         if (verbose) {
2829566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "u:\n"));
2839566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2849566063dSJacob Faibussowitsch           if (Nj) PetscCall(PetscRealView(Nj, u, viewer));
2859566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2869566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "u wedge w:\n"));
2879566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2889566063dSJacob Faibussowitsch           if (Njk) PetscCall(PetscRealView(Njk, uWw, viewer));
2899566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2909566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "x:\n"));
2919566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
2929566063dSJacob Faibussowitsch           if (N * (j + k) > 0) PetscCall(PetscRealView(N * (j + k), x, viewer));
2939566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
2949566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "u wedge w(x): %g\n", (double)uWwx));
295c4762a1bSJed Brown         }
296c4762a1bSJed Brown         /* verify wedge formula */
297c4762a1bSJed Brown         uWwxcheck = 0.;
2989566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(j + k, j, &JKj));
299c4762a1bSJed Brown         for (l = 0; l < JKj; l++) {
300c4762a1bSJed Brown           PetscBool isOdd;
301c4762a1bSJed Brown           PetscReal ux, wx;
302c4762a1bSJed Brown           PetscInt  m, p;
303c4762a1bSJed Brown 
3049566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumSplit(j + k, j, l, split, &isOdd));
3059371c9d4SSatish Balay           for (m = 0; m < j + k; m++) {
3069371c9d4SSatish Balay             for (p = 0; p < N; p++) { xsplit[m * N + p] = x[split[m] * N + p]; }
3079371c9d4SSatish Balay           }
3089566063dSJacob Faibussowitsch           PetscCall(PetscDTAltVApply(N, j, u, xsplit, &ux));
3099566063dSJacob Faibussowitsch           PetscCall(PetscDTAltVApply(N, k, w, &xsplit[j * N], &wx));
310c4762a1bSJed Brown           uWwxcheck += isOdd ? -(ux * wx) : (ux * wx);
311c4762a1bSJed Brown         }
312c4762a1bSJed Brown         diff = PetscAbsReal(uWwx - uWwxcheck);
3131dca8a05SBarry 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);
3149566063dSJacob Faibussowitsch         PetscCall(PetscFree(split));
3159566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(Nk * Njk, &uWwmat, Njk, &uWwcheck));
3169566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedgeMatrix(N, j, k, u, uWwmat));
317c4762a1bSJed Brown         if (verbose) {
3189566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(u wedge):\n"));
3199566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3209566063dSJacob Faibussowitsch           if ((Nk * Njk) > 0) PetscCall(PetscRealView(Nk * Njk, uWwmat, viewer));
3219566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
322c4762a1bSJed Brown         }
323c4762a1bSJed Brown         diff = 0.;
324c4762a1bSJed Brown         norm = 0.;
325c4762a1bSJed Brown         for (l = 0; l < Njk; l++) {
326c4762a1bSJed Brown           PetscInt  m;
327c4762a1bSJed Brown           PetscReal sum = 0.;
328c4762a1bSJed Brown 
329c4762a1bSJed Brown           for (m = 0; m < Nk; m++) sum += uWwmat[l * Nk + m] * w[m];
330c4762a1bSJed Brown           uWwcheck[l] = sum;
331c4762a1bSJed Brown           diff += PetscSqr(uWwcheck[l] - uWw[l]);
332c4762a1bSJed Brown           norm += PetscSqr(uWwcheck[l]) + PetscSqr(uWw[l]);
333c4762a1bSJed Brown         }
334c4762a1bSJed Brown         diff = PetscSqrtReal(diff);
335c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
3361dca8a05SBarry Smith         PetscCheck(diff <= PETSC_SMALL * norm, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge matrix check: wedge matrix application does not match wedge direct application");
3379566063dSJacob Faibussowitsch         PetscCall(PetscFree2(uWwmat, uWwcheck));
3389566063dSJacob Faibussowitsch         PetscCall(PetscFree4(u, uWw, x, xsplit));
3399566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
340c4762a1bSJed Brown       }
341c4762a1bSJed Brown       for (M = PetscMax(1, k); M <= N; M++) { /* pullback */
342c4762a1bSJed Brown         PetscReal *L, *u, *x;
343c4762a1bSJed Brown         PetscInt   Mk, l;
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(M, k, &Mk));
3469566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(M * N, &L, Mk, &u, M * k, &x));
3479566063dSJacob Faibussowitsch         for (l = 0; l < M * N; l++) PetscCall(PetscRandomGetValueReal(rand, &L[l]));
3489566063dSJacob Faibussowitsch         for (l = 0; l < Mk; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l]));
3499566063dSJacob Faibussowitsch         for (l = 0; l < M * k; l++) PetscCall(PetscRandomGetValueReal(rand, &x[l]));
35063a3b9bcSJacob Faibussowitsch         if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "pullback M = %" PetscInt_FMT ":\n", M));
3519566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3529566063dSJacob Faibussowitsch         PetscCall(CheckPullback(M, N, L, k, w, x, verbose, viewer));
3539566063dSJacob Faibussowitsch         if (M != N) PetscCall(CheckPullback(N, M, L, k, u, v, PETSC_FALSE, viewer));
3549566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
355c4762a1bSJed Brown         if ((k % N) && (N > 1)) {
35663a3b9bcSJacob Faibussowitsch           if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "negative pullback M = %" PetscInt_FMT ":\n", M));
3579566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3589566063dSJacob Faibussowitsch           PetscCall(CheckPullback(M, N, L, -k, w, x, verbose, viewer));
3599566063dSJacob Faibussowitsch           if (M != N) PetscCall(CheckPullback(N, M, L, -k, u, v, PETSC_FALSE, viewer));
3609566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
361c4762a1bSJed Brown         }
3629566063dSJacob Faibussowitsch         PetscCall(PetscFree3(L, u, x));
363c4762a1bSJed Brown       }
364c4762a1bSJed Brown       if (k > 0) { /* Interior */
365c4762a1bSJed Brown         PetscInt   Nkm, l, m;
366c4762a1bSJed Brown         PetscReal *wIntv0, *wIntv0check, wvcheck, diff, diffMat, normMat;
367c4762a1bSJed Brown         PetscReal *intv0mat, *matcheck;
368c4762a1bSJed Brown         PetscInt(*indices)[3];
369c4762a1bSJed Brown 
3709566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(N, k - 1, &Nkm));
3719566063dSJacob Faibussowitsch         PetscCall(PetscMalloc5(Nkm, &wIntv0, Nkm, &wIntv0check, Nk * Nkm, &intv0mat, Nk * Nkm, &matcheck, Nk * k, &indices));
3729566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVInterior(N, k, w, v, wIntv0));
3739566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVInteriorMatrix(N, k, v, intv0mat));
3749566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVInteriorPattern(N, k, indices));
375c4762a1bSJed Brown         if (verbose) {
3769566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "interior product matrix pattern:\n"));
3779566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
378c4762a1bSJed Brown           for (l = 0; l < Nk * k; l++) {
379c4762a1bSJed Brown             PetscInt row = indices[l][0];
380c4762a1bSJed Brown             PetscInt col = indices[l][1];
381c4762a1bSJed Brown             PetscInt x   = indices[l][2];
382c4762a1bSJed Brown 
38363a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "intV[%" PetscInt_FMT ",%" PetscInt_FMT "] = %sV[%" PetscInt_FMT "]\n", row, col, x < 0 ? "-" : " ", x < 0 ? -(x + 1) : x));
384c4762a1bSJed Brown           }
3859566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
386c4762a1bSJed Brown         }
387c4762a1bSJed Brown         for (l = 0; l < Nkm * Nk; l++) matcheck[l] = 0.;
388c4762a1bSJed Brown         for (l = 0; l < Nk * k; l++) {
389c4762a1bSJed Brown           PetscInt row = indices[l][0];
390c4762a1bSJed Brown           PetscInt col = indices[l][1];
391c4762a1bSJed Brown           PetscInt x   = indices[l][2];
392c4762a1bSJed Brown 
393c4762a1bSJed Brown           if (x < 0) {
394c4762a1bSJed Brown             matcheck[row * Nk + col] = -v[-(x + 1)];
395c4762a1bSJed Brown           } else {
396c4762a1bSJed Brown             matcheck[row * Nk + col] = v[x];
397c4762a1bSJed Brown           }
398c4762a1bSJed Brown         }
399c4762a1bSJed Brown         diffMat = 0.;
400c4762a1bSJed Brown         normMat = 0.;
401c4762a1bSJed Brown         for (l = 0; l < Nkm * Nk; l++) {
402c4762a1bSJed Brown           diffMat += PetscSqr(PetscAbsReal(matcheck[l] - intv0mat[l]));
403c4762a1bSJed Brown           normMat += PetscSqr(matcheck[l]) + PetscSqr(intv0mat[l]);
404c4762a1bSJed Brown         }
405c4762a1bSJed Brown         diffMat = PetscSqrtReal(diffMat);
406c4762a1bSJed Brown         normMat = PetscSqrtReal(normMat);
4071dca8a05SBarry Smith         PetscCheck(diffMat <= PETSC_SMALL * normMat, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: matrix pattern does not match matrix");
408c4762a1bSJed Brown         diffMat = 0.;
409c4762a1bSJed Brown         normMat = 0.;
410c4762a1bSJed Brown         for (l = 0; l < Nkm; l++) {
411c4762a1bSJed Brown           PetscReal sum = 0.;
412c4762a1bSJed Brown 
413c4762a1bSJed Brown           for (m = 0; m < Nk; m++) sum += intv0mat[l * Nk + m] * w[m];
414c4762a1bSJed Brown           wIntv0check[l] = sum;
415c4762a1bSJed Brown 
416c4762a1bSJed Brown           diffMat += PetscSqr(PetscAbsReal(wIntv0check[l] - wIntv0[l]));
417c4762a1bSJed Brown           normMat += PetscSqr(wIntv0check[l]) + PetscSqr(wIntv0[l]);
418c4762a1bSJed Brown         }
419c4762a1bSJed Brown         diffMat = PetscSqrtReal(diffMat);
420c4762a1bSJed Brown         normMat = PetscSqrtReal(normMat);
4211dca8a05SBarry Smith         PetscCheck(diffMat <= PETSC_SMALL * normMat, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: application does not match matrix");
422c4762a1bSJed Brown         if (verbose) {
4239566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(w int v_0):\n"));
4249566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4259566063dSJacob Faibussowitsch           if (Nkm) PetscCall(PetscRealView(Nkm, wIntv0, viewer));
4269566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
427c4762a1bSJed Brown 
4289566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(int v_0):\n"));
4299566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4309566063dSJacob Faibussowitsch           if (Nk * Nkm > 0) PetscCall(PetscRealView(Nk * Nkm, intv0mat, viewer));
4319566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
432c4762a1bSJed Brown         }
4339566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k - 1, wIntv0, &v[N], &wvcheck));
434c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
4351dca8a05SBarry 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);
4369566063dSJacob Faibussowitsch         PetscCall(PetscFree5(wIntv0, wIntv0check, intv0mat, matcheck, indices));
437c4762a1bSJed Brown       }
438c4762a1bSJed Brown       if (k >= N - k) { /* Hodge star */
439c4762a1bSJed Brown         PetscReal *u, *starw, *starstarw, wu, starwdotu;
440c4762a1bSJed Brown         PetscReal  diff, norm;
441c4762a1bSJed Brown         PetscBool  isOdd;
442c4762a1bSJed Brown         PetscInt   l;
443c4762a1bSJed Brown 
444c4762a1bSJed Brown         isOdd = (PetscBool)((k * (N - k)) & 1);
4459566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(Nk, &u, Nk, &starw, Nk, &starstarw));
4469566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVStar(N, k, 1, w, starw));
4479566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVStar(N, N - k, 1, starw, starstarw));
448c4762a1bSJed Brown         if (verbose) {
4499566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "star w:\n"));
4509566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4519566063dSJacob Faibussowitsch           if (Nk) PetscCall(PetscRealView(Nk, starw, viewer));
4529566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
453c4762a1bSJed Brown 
4549566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "star star w:\n"));
4559566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4569566063dSJacob Faibussowitsch           if (Nk) PetscCall(PetscRealView(Nk, starstarw, viewer));
4579566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
458c4762a1bSJed Brown         }
4599566063dSJacob Faibussowitsch         for (l = 0; l < Nk; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l]));
4609566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedge(N, k, N - k, w, u, &wu));
461c4762a1bSJed Brown         starwdotu = 0.;
462c4762a1bSJed Brown         for (l = 0; l < Nk; l++) starwdotu += starw[l] * u[l];
463c4762a1bSJed Brown         diff = PetscAbsReal(wu - starwdotu);
4641dca8a05SBarry 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);
465c4762a1bSJed Brown 
466c4762a1bSJed Brown         diff = 0.;
467c4762a1bSJed Brown         norm = 0.;
468c4762a1bSJed Brown         for (l = 0; l < Nk; l++) {
469c4762a1bSJed Brown           diff += PetscSqr(w[l] - (isOdd ? -starstarw[l] : starstarw[l]));
470c4762a1bSJed Brown           norm += PetscSqr(w[l]) + PetscSqr(starstarw[l]);
471c4762a1bSJed Brown         }
472c4762a1bSJed Brown         diff = PetscSqrtReal(diff);
473c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
4741dca8a05SBarry Smith         PetscCheck(diff <= PETSC_SMALL * norm, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: star(star(w)) != (-1)^(N*(N-k)) w");
4759566063dSJacob Faibussowitsch         PetscCall(PetscFree3(u, starw, starstarw));
476c4762a1bSJed Brown       }
4779566063dSJacob Faibussowitsch       PetscCall(PetscFree(v));
4789566063dSJacob Faibussowitsch       PetscCall(PetscFree(w));
4799566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
480c4762a1bSJed Brown     }
4819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
482c4762a1bSJed Brown   }
4839566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
4849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
485b122ec5aSJacob Faibussowitsch   return 0;
486c4762a1bSJed Brown }
487c4762a1bSJed Brown 
488c4762a1bSJed Brown /*TEST
489c4762a1bSJed Brown   test:
490c4762a1bSJed Brown     suffix: 1234
491c4762a1bSJed Brown     args: -verbose
492c4762a1bSJed Brown   test:
493c4762a1bSJed Brown     suffix: 56
494c4762a1bSJed Brown     args: -N 5,6
495c4762a1bSJed Brown TEST*/
496