xref: /petsc/src/dm/dt/tests/ex7.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
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);
189566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(N, k, &Nk));
199566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(M, k, &Mk));
20c4762a1bSJed Brown   if (negative) {
219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Mk, &walloc));
229566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVStar(M, M - k, 1, w, walloc));
23c4762a1bSJed Brown     ww = walloc;
24c4762a1bSJed Brown   } else {
25c4762a1bSJed Brown     ww = w;
26c4762a1bSJed Brown   }
279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nk, &Lstarw, (M*k), &Lx));
289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nk * Mk, &Lstar, Nk, &Lstarwcheck));
299566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullback(N, M, L, negative ? -k : k, w, Lstarw));
309566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(N, M, L, negative ? -k : k, Lstar));
31c4762a1bSJed Brown   if (negative) {
32c4762a1bSJed Brown     PetscReal *sLsw;
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nk, &sLsw));
359566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVStar(N, N - k, 1, Lstarw, sLsw));
369566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVApply(N, k, sLsw, x, &Lstarwx));
379566063dSJacob Faibussowitsch     PetscCall(PetscFree(sLsw));
38c4762a1bSJed Brown   } else {
399566063dSJacob Faibussowitsch     PetscCall(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) {
639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L:\n"));
649566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
659566063dSJacob Faibussowitsch     if (M*N > 0) PetscCall(PetscRealView(M*N, L, viewer));
669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L*:\n"));
699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
709566063dSJacob Faibussowitsch     if (Nk * Mk > 0) PetscCall(PetscRealView(Nk * Mk, Lstar, viewer));
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "L*w:\n"));
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
759566063dSJacob Faibussowitsch     if (Nk > 0) PetscCall(PetscRealView(Nk, Lstarw, viewer));
769566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
77c4762a1bSJed Brown   }
789566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVApply(M, k, ww, Lx, &wLx));
79c4762a1bSJed Brown   diff = PetscAbsReal(wLx - Lstarwx);
80*1dca8a05SBarry 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);
81*1dca8a05SBarry Smith   PetscCheck(diffMat <= PETSC_SMALL * normMat,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback matrix does match matrix free result");
829566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Lstar, Lstarwcheck));
839566063dSJacob Faibussowitsch   PetscCall(PetscFree2(Lstarw, Lx));
849566063dSJacob Faibussowitsch   PetscCall(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 
959566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
96d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for exterior algebra tests","none");
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-N", "Up to 5 vector space dimensions to test","ex7.c",n,&numTests,NULL));
989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-verbose", "Verbose test output","ex7.c",verbose,&verbose,NULL));
99d0609cedSBarry Smith   PetscOptionsEnd();
1009566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
1019566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(rand, -1., 1.));
1029566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
103c4762a1bSJed Brown   if (!numTests) numTests = 5;
104c4762a1bSJed Brown   viewer = PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD);
105c4762a1bSJed Brown   for (i = 0; i < numTests; i++) {
106c4762a1bSJed Brown     PetscInt       k, N = n[i];
107c4762a1bSJed Brown 
10863a3b9bcSJacob Faibussowitsch     if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "N = %" PetscInt_FMT ":\n", N));
1099566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown     if (verbose) {
112c4762a1bSJed Brown       PetscInt *perm;
113c4762a1bSJed Brown       PetscInt fac = 1;
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(N, &perm));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown       for (k = 1; k <= N; k++) fac *= k;
11863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Permutations of %" PetscInt_FMT ":\n", N));
1199566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
120c4762a1bSJed Brown       for (k = 0; k < fac; k++) {
121c4762a1bSJed Brown         PetscBool isOdd, isOddCheck;
122c4762a1bSJed Brown         PetscInt  j, kCheck;
123c4762a1bSJed Brown 
1249566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumPerm(N, k, perm, &isOdd));
12563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ":", k));
126c4762a1bSJed Brown         for (j = 0; j < N; j++) {
12763a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, perm[j]));
128c4762a1bSJed Brown         }
1299566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even"));
1309566063dSJacob Faibussowitsch         PetscCall(PetscDTPermIndex(N, perm, &kCheck, &isOddCheck));
131*1dca8a05SBarry Smith         PetscCheck(kCheck == k && isOddCheck == isOdd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumPerm / PetscDTPermIndex mismatch for (%" PetscInt_FMT ", %" PetscInt_FMT ")", N, k);
132c4762a1bSJed Brown       }
1339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
1349566063dSJacob Faibussowitsch       PetscCall(PetscFree(perm));
135c4762a1bSJed Brown     }
136c4762a1bSJed Brown     for (k = 0; k <= N; k++) {
137c4762a1bSJed Brown       PetscInt   j, Nk, M;
138c4762a1bSJed Brown       PetscReal *w, *v, wv;
139c4762a1bSJed Brown       PetscInt  *subset;
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch       PetscCall(PetscDTBinomialInt(N, k, &Nk));
14263a3b9bcSJacob Faibussowitsch       if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "k = %" PetscInt_FMT ":\n", k));
1439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
14463a3b9bcSJacob Faibussowitsch       if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT " choose %" PetscInt_FMT "): %" PetscInt_FMT "\n", N, k, Nk));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown       /* Test subset and complement enumeration */
1479566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(N, &subset));
1489566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
149c4762a1bSJed Brown       for (j = 0; j < Nk; j++) {
150c4762a1bSJed Brown         PetscBool isOdd, isOddCheck;
151c4762a1bSJed Brown         PetscInt  jCheck, kCheck;
152c4762a1bSJed Brown 
1539566063dSJacob Faibussowitsch         PetscCall(PetscDTEnumSplit(N, k, j, subset, &isOdd));
1549566063dSJacob Faibussowitsch         PetscCall(PetscDTPermIndex(N, subset, &kCheck, &isOddCheck));
15508401ef6SPierre Jolivet         PetscCheck(isOddCheck == isOdd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumSplit sign does not mmatch PetscDTPermIndex sign");
156c4762a1bSJed Brown         if (verbose) {
157c4762a1bSJed Brown           PetscInt l;
158c4762a1bSJed Brown 
15963a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "subset %" PetscInt_FMT ":", j));
160c4762a1bSJed Brown           for (l = 0; l < k; l++) {
16163a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, subset[l]));
162c4762a1bSJed Brown           }
1639566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, " |"));
164c4762a1bSJed Brown           for (l = k; l < N; l++) {
16563a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %" PetscInt_FMT, subset[l]));
166c4762a1bSJed Brown           }
1679566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even"));
168c4762a1bSJed Brown         }
1699566063dSJacob Faibussowitsch         PetscCall(PetscDTSubsetIndex(N, k, subset, &jCheck));
17063a3b9bcSJacob Faibussowitsch         PetscCheck(jCheck == j,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "jCheck (%" PetscInt_FMT ") != j (%" PetscInt_FMT ")", jCheck, j);
171c4762a1bSJed Brown       }
1729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
1739566063dSJacob Faibussowitsch       PetscCall(PetscFree(subset));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown       /* Make a random k form */
1769566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(Nk, &w));
1779566063dSJacob Faibussowitsch       for (j = 0; j < Nk; j++) PetscCall(PetscRandomGetValueReal(rand, &w[j]));
178c4762a1bSJed Brown       /* Make a set of random vectors */
1799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(N*k, &v));
1809566063dSJacob Faibussowitsch       for (j = 0; j < N*k; j++) PetscCall(PetscRandomGetValueReal(rand, &v[j]));
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVApply(N, k, w, v, &wv));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown       if (verbose) {
1859566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "w:\n"));
1869566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
1879566063dSJacob Faibussowitsch         if (Nk) PetscCall(PetscRealView(Nk, w, viewer));
1889566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
1899566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "v:\n"));
1909566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
1919566063dSJacob Faibussowitsch         if (N*k > 0) PetscCall(PetscRealView(N*k, v, viewer));
1929566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
1939566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "w(v): %g\n", (double) wv));
194c4762a1bSJed Brown       }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown       /* sanity checks */
197c4762a1bSJed Brown       if (k == 1) { /* 1-forms are functionals (dot products) */
198c4762a1bSJed Brown         PetscInt  l;
199c4762a1bSJed Brown         PetscReal wvcheck = 0.;
200c4762a1bSJed Brown         PetscReal diff;
201c4762a1bSJed Brown 
202c4762a1bSJed Brown         for (l = 0; l < N; l++) wvcheck += w[l] * v[l];
203c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
204*1dca8a05SBarry 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);
205c4762a1bSJed Brown       }
206c4762a1bSJed Brown       if (k == N && N < 5) { /* n-forms are scaled determinants */
207c4762a1bSJed Brown         PetscReal det, wvcheck, diff;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown         switch (k) {
210c4762a1bSJed Brown         case 0:
211c4762a1bSJed Brown           det = 1.;
212c4762a1bSJed Brown           break;
213c4762a1bSJed Brown         case 1:
214c4762a1bSJed Brown           det = v[0];
215c4762a1bSJed Brown           break;
216c4762a1bSJed Brown         case 2:
217c4762a1bSJed Brown           det = v[0] * v[3] - v[1] * v[2];
218c4762a1bSJed Brown           break;
219c4762a1bSJed Brown         case 3:
220c4762a1bSJed Brown           det = v[0] * (v[4] * v[8] - v[5] * v[7]) +
221c4762a1bSJed Brown                 v[1] * (v[5] * v[6] - v[3] * v[8]) +
222c4762a1bSJed Brown                 v[2] * (v[3] * v[7] - v[4] * v[6]);
223c4762a1bSJed Brown           break;
224c4762a1bSJed Brown         case 4:
225c4762a1bSJed Brown           det = v[0] * (v[5] * (v[10] * v[15] - v[11] * v[14]) +
226c4762a1bSJed Brown                         v[6] * (v[11] * v[13] - v[ 9] * v[15]) +
227c4762a1bSJed Brown                         v[7] * (v[ 9] * v[14] - v[10] * v[13])) -
228c4762a1bSJed Brown                 v[1] * (v[4] * (v[10] * v[15] - v[11] * v[14]) +
229c4762a1bSJed Brown                         v[6] * (v[11] * v[12] - v[ 8] * v[15]) +
230c4762a1bSJed Brown                         v[7] * (v[ 8] * v[14] - v[10] * v[12])) +
231c4762a1bSJed Brown                 v[2] * (v[4] * (v[ 9] * v[15] - v[11] * v[13]) +
232c4762a1bSJed Brown                         v[5] * (v[11] * v[12] - v[ 8] * v[15]) +
233c4762a1bSJed Brown                         v[7] * (v[ 8] * v[13] - v[ 9] * v[12])) -
234c4762a1bSJed Brown                 v[3] * (v[4] * (v[ 9] * v[14] - v[10] * v[13]) +
235c4762a1bSJed Brown                         v[5] * (v[10] * v[12] - v[ 8] * v[14]) +
236c4762a1bSJed Brown                         v[6] * (v[ 8] * v[13] - v[ 9] * v[12]));
237c4762a1bSJed Brown           break;
238c4762a1bSJed Brown         default:
239c4762a1bSJed Brown           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "invalid k");
240c4762a1bSJed Brown         }
241c4762a1bSJed Brown         wvcheck = det * w[0];
242c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
243*1dca8a05SBarry 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);
244c4762a1bSJed Brown       }
245c4762a1bSJed Brown       if (k > 0) { /* k-forms are linear in each component */
246c4762a1bSJed Brown         PetscReal alpha;
247c4762a1bSJed Brown         PetscReal *x, *axv, wx, waxv, waxvcheck;
248c4762a1bSJed Brown         PetscReal diff;
249c4762a1bSJed Brown         PetscReal rj;
250c4762a1bSJed Brown         PetscInt  l;
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(N * k, &x, N * k, &axv));
2539566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rand, &alpha));
2549566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, 0, k));
2559566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rand, &rj));
256c4762a1bSJed Brown         j = (PetscInt) rj;
2579566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, -1., 1.));
258c4762a1bSJed Brown         for (l = 0; l < N*k; l++) x[l] = v[l];
259c4762a1bSJed Brown         for (l = 0; l < N*k; l++) axv[l] = v[l];
260c4762a1bSJed Brown         for (l = 0; l < N; l++) {
261c4762a1bSJed Brown           PetscReal val;
262c4762a1bSJed Brown 
2639566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rand, &val));
264c4762a1bSJed Brown           x[j * N + l] = val;
265c4762a1bSJed Brown           axv[j * N + l] += alpha * val;
266c4762a1bSJed Brown         }
2679566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k, w, x, &wx));
2689566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k, w, axv, &waxv));
269c4762a1bSJed Brown         waxvcheck = alpha * wx + wv;
270c4762a1bSJed Brown         diff = waxv - waxvcheck;
271*1dca8a05SBarry 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);
2729566063dSJacob Faibussowitsch         PetscCall(PetscFree2(x,axv));
273c4762a1bSJed Brown       }
274c4762a1bSJed Brown       if (k > 1) { /* k-forms are antisymmetric */
275c4762a1bSJed Brown         PetscReal rj, rl, *swapv, wswapv, diff;
276c4762a1bSJed Brown         PetscInt  l, m;
277c4762a1bSJed Brown 
2789566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, 0, k));
2799566063dSJacob Faibussowitsch         PetscCall(PetscRandomGetValueReal(rand, &rj));
280c4762a1bSJed Brown         j = (PetscInt) rj;
281c4762a1bSJed Brown         l = j;
282c4762a1bSJed Brown         while (l == j) {
2839566063dSJacob Faibussowitsch           PetscCall(PetscRandomGetValueReal(rand, &rl));
284c4762a1bSJed Brown           l = (PetscInt) rl;
285c4762a1bSJed Brown         }
2869566063dSJacob Faibussowitsch         PetscCall(PetscRandomSetInterval(rand, -1., 1.));
2879566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(N * k, &swapv));
288c4762a1bSJed Brown         for (m = 0; m < N * k; m++) swapv[m] = v[m];
289c4762a1bSJed Brown         for (m = 0; m < N; m++) {
290c4762a1bSJed Brown           swapv[j * N + m] = v[l * N + m];
291c4762a1bSJed Brown           swapv[l * N + m] = v[j * N + m];
292c4762a1bSJed Brown         }
2939566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k, w, swapv, &wswapv));
294c4762a1bSJed Brown         diff = PetscAbsReal(wswapv + wv);
295*1dca8a05SBarry 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);
2969566063dSJacob Faibussowitsch         PetscCall(PetscFree(swapv));
297c4762a1bSJed Brown       }
298c4762a1bSJed Brown       for (j = 0; j <= k && j + k <= N; j++) { /* wedge product */
299c4762a1bSJed Brown         PetscInt   Nj, Njk, l, JKj;
300c4762a1bSJed Brown         PetscReal *u, *uWw, *uWwcheck, *uWwmat, *x, *xsplit, uWwx, uWwxcheck, diff, norm;
301c4762a1bSJed Brown         PetscInt  *split;
302c4762a1bSJed Brown 
30363a3b9bcSJacob Faibussowitsch         if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "wedge j = %" PetscInt_FMT ":\n", j));
3049566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3059566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(N, j,   &Nj));
3069566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(N, j+k, &Njk));
3079566063dSJacob Faibussowitsch         PetscCall(PetscMalloc4(Nj, &u, Njk, &uWw, N*(j+k), &x, N*(j+k), &xsplit));
3089566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(j+k,&split));
3099566063dSJacob Faibussowitsch         for (l = 0; l < Nj; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l]));
3109566063dSJacob Faibussowitsch         for (l = 0; l < N*(j+k); l++) PetscCall(PetscRandomGetValueReal(rand, &x[l]));
3119566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedge(N, j, k, u, w, uWw));
3129566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, j+k, uWw, x, &uWwx));
313c4762a1bSJed Brown         if (verbose) {
3149566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "u:\n"));
3159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3169566063dSJacob Faibussowitsch           if (Nj) PetscCall(PetscRealView(Nj, u, viewer));
3179566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
3189566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "u wedge w:\n"));
3199566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3209566063dSJacob Faibussowitsch           if (Njk) PetscCall(PetscRealView(Njk, uWw, viewer));
3219566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
3229566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "x:\n"));
3239566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3249566063dSJacob Faibussowitsch           if (N*(j+k) > 0) PetscCall(PetscRealView(N*(j+k), x, viewer));
3259566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
3269566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "u wedge w(x): %g\n", (double) uWwx));
327c4762a1bSJed Brown         }
328c4762a1bSJed Brown         /* verify wedge formula */
329c4762a1bSJed Brown         uWwxcheck = 0.;
3309566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(j+k, j, &JKj));
331c4762a1bSJed Brown         for (l = 0; l < JKj; l++) {
332c4762a1bSJed Brown           PetscBool isOdd;
333c4762a1bSJed Brown           PetscReal ux, wx;
334c4762a1bSJed Brown           PetscInt  m, p;
335c4762a1bSJed Brown 
3369566063dSJacob Faibussowitsch           PetscCall(PetscDTEnumSplit(j+k, j, l, split, &isOdd));
337c4762a1bSJed Brown           for (m = 0; m < j+k; m++) {for (p = 0; p < N; p++) {xsplit[m * N + p] = x[split[m] * N + p];}}
3389566063dSJacob Faibussowitsch           PetscCall(PetscDTAltVApply(N, j, u, xsplit, &ux));
3399566063dSJacob Faibussowitsch           PetscCall(PetscDTAltVApply(N, k, w, &xsplit[j*N], &wx));
340c4762a1bSJed Brown           uWwxcheck += isOdd ? -(ux * wx) : (ux * wx);
341c4762a1bSJed Brown         }
342c4762a1bSJed Brown         diff = PetscAbsReal(uWwx - uWwxcheck);
343*1dca8a05SBarry 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);
3449566063dSJacob Faibussowitsch         PetscCall(PetscFree(split));
3459566063dSJacob Faibussowitsch         PetscCall(PetscMalloc2(Nk * Njk, &uWwmat, Njk, &uWwcheck));
3469566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedgeMatrix(N, j, k, u, uWwmat));
347c4762a1bSJed Brown         if (verbose) {
3489566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(u wedge):\n"));
3499566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3509566063dSJacob Faibussowitsch           if ((Nk * Njk) > 0) PetscCall(PetscRealView(Nk * Njk, uWwmat, viewer));
3519566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
352c4762a1bSJed Brown         }
353c4762a1bSJed Brown         diff = 0.;
354c4762a1bSJed Brown         norm = 0.;
355c4762a1bSJed Brown         for (l = 0; l < Njk; l++) {
356c4762a1bSJed Brown           PetscInt  m;
357c4762a1bSJed Brown           PetscReal sum = 0.;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown           for (m = 0; m < Nk; m++) sum += uWwmat[l * Nk + m] * w[m];
360c4762a1bSJed Brown           uWwcheck[l] = sum;
361c4762a1bSJed Brown           diff += PetscSqr(uWwcheck[l] - uWw[l]);
362c4762a1bSJed Brown           norm += PetscSqr(uWwcheck[l]) + PetscSqr(uWw[l]);
363c4762a1bSJed Brown         }
364c4762a1bSJed Brown         diff = PetscSqrtReal(diff);
365c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
366*1dca8a05SBarry Smith         PetscCheck(diff <= PETSC_SMALL * norm,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge matrix check: wedge matrix application does not match wedge direct application");
3679566063dSJacob Faibussowitsch         PetscCall(PetscFree2(uWwmat, uWwcheck));
3689566063dSJacob Faibussowitsch         PetscCall(PetscFree4(u, uWw, x, xsplit));
3699566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
370c4762a1bSJed Brown       }
371c4762a1bSJed Brown       for (M = PetscMax(1,k); M <= N; M++) { /* pullback */
372c4762a1bSJed Brown         PetscReal   *L, *u, *x;
373c4762a1bSJed Brown         PetscInt     Mk, l;
374c4762a1bSJed Brown 
3759566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(M, k, &Mk));
3769566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(M*N, &L, Mk, &u, M*k, &x));
3779566063dSJacob Faibussowitsch         for (l = 0; l < M*N; l++) PetscCall(PetscRandomGetValueReal(rand, &L[l]));
3789566063dSJacob Faibussowitsch         for (l = 0; l < Mk; l++) PetscCall(PetscRandomGetValueReal(rand, &u[l]));
3799566063dSJacob Faibussowitsch         for (l = 0; l < M*k; l++) PetscCall(PetscRandomGetValueReal(rand, &x[l]));
38063a3b9bcSJacob Faibussowitsch         if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "pullback M = %" PetscInt_FMT ":\n", M));
3819566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3829566063dSJacob Faibussowitsch         PetscCall(CheckPullback(M, N, L, k, w, x, verbose, viewer));
3839566063dSJacob Faibussowitsch         if (M != N) PetscCall(CheckPullback(N, M, L, k, u, v, PETSC_FALSE, viewer));
3849566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
385c4762a1bSJed Brown         if ((k % N) && (N > 1)) {
38663a3b9bcSJacob Faibussowitsch           if (verbose) PetscCall(PetscViewerASCIIPrintf(viewer, "negative pullback M = %" PetscInt_FMT ":\n", M));
3879566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
3889566063dSJacob Faibussowitsch           PetscCall(CheckPullback(M, N, L, -k, w, x, verbose, viewer));
3899566063dSJacob Faibussowitsch           if (M != N) PetscCall(CheckPullback(N, M, L, -k, u, v, PETSC_FALSE, viewer));
3909566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
391c4762a1bSJed Brown         }
3929566063dSJacob Faibussowitsch         PetscCall(PetscFree3(L, u, x));
393c4762a1bSJed Brown       }
394c4762a1bSJed Brown       if (k > 0) { /* Interior */
395c4762a1bSJed Brown         PetscInt    Nkm, l, m;
396c4762a1bSJed Brown         PetscReal  *wIntv0, *wIntv0check, wvcheck, diff, diffMat, normMat;
397c4762a1bSJed Brown         PetscReal  *intv0mat, *matcheck;
398c4762a1bSJed Brown         PetscInt  (*indices)[3];
399c4762a1bSJed Brown 
4009566063dSJacob Faibussowitsch         PetscCall(PetscDTBinomialInt(N, k-1, &Nkm));
4019566063dSJacob Faibussowitsch         PetscCall(PetscMalloc5(Nkm, &wIntv0, Nkm, &wIntv0check, Nk * Nkm, &intv0mat, Nk * Nkm, &matcheck, Nk * k, &indices));
4029566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVInterior(N, k, w, v, wIntv0));
4039566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVInteriorMatrix(N, k, v, intv0mat));
4049566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVInteriorPattern(N, k, indices));
405c4762a1bSJed Brown         if (verbose) {
4069566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "interior product matrix pattern:\n"));
4079566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
408c4762a1bSJed Brown           for (l = 0; l < Nk * k; l++) {
409c4762a1bSJed Brown             PetscInt row = indices[l][0];
410c4762a1bSJed Brown             PetscInt col = indices[l][1];
411c4762a1bSJed Brown             PetscInt x   = indices[l][2];
412c4762a1bSJed Brown 
41363a3b9bcSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer,"intV[%" PetscInt_FMT ",%" PetscInt_FMT "] = %sV[%" PetscInt_FMT "]\n", row, col, x < 0 ? "-" : " ", x < 0 ? -(x + 1) : x));
414c4762a1bSJed Brown           }
4159566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
416c4762a1bSJed Brown         }
417c4762a1bSJed Brown         for (l = 0; l < Nkm * Nk; l++) matcheck[l] = 0.;
418c4762a1bSJed Brown         for (l = 0; l < Nk * k; l++) {
419c4762a1bSJed Brown           PetscInt row = indices[l][0];
420c4762a1bSJed Brown           PetscInt col = indices[l][1];
421c4762a1bSJed Brown           PetscInt x   = indices[l][2];
422c4762a1bSJed Brown 
423c4762a1bSJed Brown           if (x < 0) {
424c4762a1bSJed Brown             matcheck[row * Nk + col] = -v[-(x+1)];
425c4762a1bSJed Brown           } else {
426c4762a1bSJed Brown             matcheck[row * Nk + col] = v[x];
427c4762a1bSJed Brown           }
428c4762a1bSJed Brown         }
429c4762a1bSJed Brown         diffMat = 0.;
430c4762a1bSJed Brown         normMat = 0.;
431c4762a1bSJed Brown         for (l = 0; l < Nkm * Nk; l++) {
432c4762a1bSJed Brown           diffMat += PetscSqr(PetscAbsReal(matcheck[l] - intv0mat[l]));
433c4762a1bSJed Brown           normMat += PetscSqr(matcheck[l]) + PetscSqr(intv0mat[l]);
434c4762a1bSJed Brown         }
435c4762a1bSJed Brown         diffMat = PetscSqrtReal(diffMat);
436c4762a1bSJed Brown         normMat = PetscSqrtReal(normMat);
437*1dca8a05SBarry Smith         PetscCheck(diffMat <= PETSC_SMALL * normMat,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: matrix pattern does not match matrix");
438c4762a1bSJed Brown         diffMat = 0.;
439c4762a1bSJed Brown         normMat = 0.;
440c4762a1bSJed Brown         for (l = 0; l < Nkm; l++) {
441c4762a1bSJed Brown           PetscReal sum = 0.;
442c4762a1bSJed Brown 
443c4762a1bSJed Brown           for (m = 0; m < Nk; m++) sum += intv0mat[l * Nk + m] * w[m];
444c4762a1bSJed Brown           wIntv0check[l] = sum;
445c4762a1bSJed Brown 
446c4762a1bSJed Brown           diffMat += PetscSqr(PetscAbsReal(wIntv0check[l] - wIntv0[l]));
447c4762a1bSJed Brown           normMat += PetscSqr(wIntv0check[l]) + PetscSqr(wIntv0[l]);
448c4762a1bSJed Brown         }
449c4762a1bSJed Brown         diffMat = PetscSqrtReal(diffMat);
450c4762a1bSJed Brown         normMat = PetscSqrtReal(normMat);
451*1dca8a05SBarry Smith         PetscCheck(diffMat <= PETSC_SMALL * normMat,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: application does not match matrix");
452c4762a1bSJed Brown         if (verbose) {
4539566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(w int v_0):\n"));
4549566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4559566063dSJacob Faibussowitsch           if (Nkm) PetscCall(PetscRealView(Nkm, wIntv0, viewer));
4569566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
457c4762a1bSJed Brown 
4589566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(int v_0):\n"));
4599566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4609566063dSJacob Faibussowitsch           if (Nk * Nkm > 0) PetscCall(PetscRealView(Nk * Nkm, intv0mat, viewer));
4619566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
462c4762a1bSJed Brown         }
4639566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVApply(N, k - 1, wIntv0, &v[N], &wvcheck));
464c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
465*1dca8a05SBarry 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);
4669566063dSJacob Faibussowitsch         PetscCall(PetscFree5(wIntv0,wIntv0check,intv0mat,matcheck,indices));
467c4762a1bSJed Brown       }
468c4762a1bSJed Brown       if (k >= N - k) { /* Hodge star */
469c4762a1bSJed Brown         PetscReal *u, *starw, *starstarw, wu, starwdotu;
470c4762a1bSJed Brown         PetscReal diff, norm;
471c4762a1bSJed Brown         PetscBool isOdd;
472c4762a1bSJed Brown         PetscInt l;
473c4762a1bSJed Brown 
474c4762a1bSJed Brown         isOdd = (PetscBool) ((k * (N - k)) & 1);
4759566063dSJacob Faibussowitsch         PetscCall(PetscMalloc3(Nk, &u, Nk, &starw, Nk, &starstarw));
4769566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVStar(N, k, 1, w, starw));
4779566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVStar(N, N-k, 1, starw, starstarw));
478c4762a1bSJed Brown         if (verbose) {
4799566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "star w:\n"));
4809566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4819566063dSJacob Faibussowitsch           if (Nk) PetscCall(PetscRealView(Nk, starw, viewer));
4829566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
483c4762a1bSJed Brown 
4849566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "star star w:\n"));
4859566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
4869566063dSJacob Faibussowitsch           if (Nk) PetscCall(PetscRealView(Nk, starstarw, viewer));
4879566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
488c4762a1bSJed Brown         }
4899566063dSJacob Faibussowitsch         for (l = 0; l < Nk; l++) PetscCall(PetscRandomGetValueReal(rand,&u[l]));
4909566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedge(N, k, N - k, w, u, &wu));
491c4762a1bSJed Brown         starwdotu = 0.;
492c4762a1bSJed Brown         for (l = 0; l < Nk; l++) starwdotu += starw[l] * u[l];
493c4762a1bSJed Brown         diff = PetscAbsReal(wu - starwdotu);
494*1dca8a05SBarry 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);
495c4762a1bSJed Brown 
496c4762a1bSJed Brown         diff = 0.;
497c4762a1bSJed Brown         norm = 0.;
498c4762a1bSJed Brown         for (l = 0; l < Nk; l++) {
499c4762a1bSJed Brown           diff += PetscSqr(w[l] - (isOdd ? -starstarw[l] : starstarw[l]));
500c4762a1bSJed Brown           norm += PetscSqr(w[l]) + PetscSqr(starstarw[l]);
501c4762a1bSJed Brown         }
502c4762a1bSJed Brown         diff = PetscSqrtReal(diff);
503c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
504*1dca8a05SBarry Smith         PetscCheck(diff <= PETSC_SMALL * norm,PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: star(star(w)) != (-1)^(N*(N-k)) w");
5059566063dSJacob Faibussowitsch         PetscCall(PetscFree3(u, starw, starstarw));
506c4762a1bSJed Brown       }
5079566063dSJacob Faibussowitsch       PetscCall(PetscFree(v));
5089566063dSJacob Faibussowitsch       PetscCall(PetscFree(w));
5099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
510c4762a1bSJed Brown     }
5119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
512c4762a1bSJed Brown   }
5139566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
5149566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
515b122ec5aSJacob Faibussowitsch   return 0;
516c4762a1bSJed Brown }
517c4762a1bSJed Brown 
518c4762a1bSJed Brown /*TEST
519c4762a1bSJed Brown   test:
520c4762a1bSJed Brown     suffix: 1234
521c4762a1bSJed Brown     args: -verbose
522c4762a1bSJed Brown   test:
523c4762a1bSJed Brown     suffix: 56
524c4762a1bSJed Brown     args: -N 5,6
525c4762a1bSJed Brown TEST*/
526