xref: /petsc/src/dm/dt/tests/ex7.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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   PetscErrorCode  ierr;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   PetscFunctionBegin;
18c4762a1bSJed Brown   k = PetscAbsInt(k);
19c4762a1bSJed Brown   ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr);
20c4762a1bSJed Brown   ierr = PetscDTBinomialInt(M, k, &Mk);CHKERRQ(ierr);
21c4762a1bSJed Brown   if (negative) {
22c4762a1bSJed Brown     ierr = PetscMalloc1(Mk, &walloc);CHKERRQ(ierr);
23c4762a1bSJed Brown     ierr = PetscDTAltVStar(M, M - k, 1, w, walloc);CHKERRQ(ierr);
24c4762a1bSJed Brown     ww = walloc;
25c4762a1bSJed Brown   } else {
26c4762a1bSJed Brown     ww = w;
27c4762a1bSJed Brown   }
28c4762a1bSJed Brown   ierr = PetscMalloc2(Nk, &Lstarw, (M*k), &Lx);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = PetscMalloc2(Nk * Mk, &Lstar, Nk, &Lstarwcheck);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = PetscDTAltVPullback(N, M, L, negative ? -k : k, w, Lstarw);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = PetscDTAltVPullbackMatrix(N, M, L, negative ? -k : k, Lstar);CHKERRQ(ierr);
32c4762a1bSJed Brown   if (negative) {
33c4762a1bSJed Brown     PetscReal *sLsw;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown     ierr = PetscMalloc1(Nk, &sLsw);CHKERRQ(ierr);
36c4762a1bSJed Brown     ierr = PetscDTAltVStar(N, N - k, 1, Lstarw, sLsw);CHKERRQ(ierr);
37c4762a1bSJed Brown     ierr = PetscDTAltVApply(N, k, sLsw, x, &Lstarwx);CHKERRQ(ierr);
38c4762a1bSJed Brown     ierr = PetscFree(sLsw);CHKERRQ(ierr);
39c4762a1bSJed Brown   } else {
40c4762a1bSJed Brown     ierr = PetscDTAltVApply(N, k, Lstarw, x, &Lstarwx);CHKERRQ(ierr);
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown   for (l = 0; l < k; l++) {
43c4762a1bSJed Brown     for (i = 0; i < M; i++) {
44c4762a1bSJed Brown       PetscReal sum = 0.;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown       for (j = 0; j < N; j++) sum += L[i * N + j] * x[l * N + j];
47c4762a1bSJed Brown       Lx[l * M + i] = sum;
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown   diffMat = 0.;
51c4762a1bSJed Brown   normMat = 0.;
52c4762a1bSJed Brown   for (i = 0; i < Nk; i++) {
53c4762a1bSJed Brown     PetscReal sum = 0.;
54c4762a1bSJed Brown     for (j = 0; j < Mk; j++) {
55c4762a1bSJed Brown       sum += Lstar[i * Mk + j] * w[j];
56c4762a1bSJed Brown     }
57c4762a1bSJed Brown     Lstarwcheck[i] = sum;
58c4762a1bSJed Brown     diffMat += PetscSqr(PetscAbsReal(Lstarwcheck[i] - Lstarw[i]));
59c4762a1bSJed Brown     normMat += PetscSqr(Lstarwcheck[i]) +  PetscSqr(Lstarw[i]);
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown   diffMat = PetscSqrtReal(diffMat);
62c4762a1bSJed Brown   normMat = PetscSqrtReal(normMat);
63c4762a1bSJed Brown   if (verbose) {
64c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(viewer, "L:\n");CHKERRQ(ierr);
65c4762a1bSJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
66c4762a1bSJed Brown     if (M*N > 0) {ierr = PetscRealView(M*N, L, viewer);CHKERRQ(ierr);}
67c4762a1bSJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(viewer, "L*:\n");CHKERRQ(ierr);
70c4762a1bSJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
71c4762a1bSJed Brown     if (Nk * Mk > 0) {ierr = PetscRealView(Nk * Mk, Lstar, viewer);CHKERRQ(ierr);}
72c4762a1bSJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
73c4762a1bSJed Brown 
74c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(viewer, "L*w:\n");CHKERRQ(ierr);
75c4762a1bSJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
76c4762a1bSJed Brown     if (Nk > 0) {ierr = PetscRealView(Nk, Lstarw, viewer);CHKERRQ(ierr);}
77c4762a1bSJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown   ierr = PetscDTAltVApply(M, k, ww, Lx, &wLx);CHKERRQ(ierr);
80c4762a1bSJed Brown   diff = PetscAbsReal(wLx - Lstarwx);
81c4762a1bSJed Brown   if (diff > 10. * PETSC_SMALL * (PetscAbsReal(wLx) + PetscAbsReal(Lstarwx))) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback does not commute with application: w(Lx)(%g) != (L* w)(x)(%g)", wLx, Lstarwx);
82c4762a1bSJed Brown   if (diffMat > PETSC_SMALL * normMat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "pullback check: pullback matrix does match matrix free result");
83c4762a1bSJed Brown   ierr = PetscFree2(Lstar, Lstarwcheck);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = PetscFree2(Lstarw, Lx);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = PetscFree(walloc);CHKERRQ(ierr);
86c4762a1bSJed Brown   PetscFunctionReturn(0);
87c4762a1bSJed Brown }
88c4762a1bSJed Brown 
89c4762a1bSJed Brown int main(int argc, char **argv)
90c4762a1bSJed Brown {
91c4762a1bSJed Brown   PetscInt       i, numTests = 5, n[5] = {0, 1, 2, 3, 4};
92c4762a1bSJed Brown   PetscBool      verbose = PETSC_FALSE;
93c4762a1bSJed Brown   PetscRandom    rand;
94c4762a1bSJed Brown   PetscViewer    viewer;
95c4762a1bSJed Brown   PetscErrorCode ierr;
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
98c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for exterior algebra tests","none");CHKERRQ(ierr);
99fbdc3dfeSToby Isaac   ierr = PetscOptionsIntArray("-N", "Up to 5 vector space dimensions to test","ex7.c",n,&numTests,NULL);CHKERRQ(ierr);
100fbdc3dfeSToby Isaac   ierr = PetscOptionsBool("-verbose", "Verbose test output","ex7.c",verbose,&verbose,NULL);CHKERRQ(ierr);
101*1e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_SELF, &rand);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = PetscRandomSetInterval(rand, -1., 1.);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
105c4762a1bSJed Brown   if (!numTests) numTests = 5;
106c4762a1bSJed Brown   viewer = PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD);
107c4762a1bSJed Brown   for (i = 0; i < numTests; i++) {
108c4762a1bSJed Brown     PetscInt       k, N = n[i];
109c4762a1bSJed Brown 
110c4762a1bSJed Brown     if (verbose) {ierr = PetscViewerASCIIPrintf(viewer, "N = %D:\n", N);CHKERRQ(ierr);}
111c4762a1bSJed Brown     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
112c4762a1bSJed Brown 
113c4762a1bSJed Brown     if (verbose) {
114c4762a1bSJed Brown       PetscInt *perm;
115c4762a1bSJed Brown       PetscInt fac = 1;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown       ierr = PetscMalloc1(N, &perm);CHKERRQ(ierr);
118c4762a1bSJed Brown 
119c4762a1bSJed Brown       for (k = 1; k <= N; k++) fac *= k;
120c4762a1bSJed Brown       ierr = PetscViewerASCIIPrintf(viewer, "Permutations of %D:\n", N);CHKERRQ(ierr);
121c4762a1bSJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
122c4762a1bSJed Brown       for (k = 0; k < fac; k++) {
123c4762a1bSJed Brown         PetscBool isOdd, isOddCheck;
124c4762a1bSJed Brown         PetscInt  j, kCheck;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown         ierr = PetscDTEnumPerm(N, k, perm, &isOdd);CHKERRQ(ierr);
127c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer, "%D:", k);CHKERRQ(ierr);
128c4762a1bSJed Brown         for (j = 0; j < N; j++) {
129c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_WORLD, " %D", perm[j]);CHKERRQ(ierr);
130c4762a1bSJed Brown         }
131c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even");CHKERRQ(ierr);
132c4762a1bSJed Brown         ierr = PetscDTPermIndex(N, perm, &kCheck, &isOddCheck);CHKERRQ(ierr);
133c4762a1bSJed Brown         if (kCheck != k || isOddCheck != isOdd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumPerm / PetscDTPermIndex mismatch for (%D, %D)\n", N, k);
134c4762a1bSJed Brown       }
135c4762a1bSJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
136c4762a1bSJed Brown       ierr = PetscFree(perm);CHKERRQ(ierr);
137c4762a1bSJed Brown     }
138c4762a1bSJed Brown     for (k = 0; k <= N; k++) {
139c4762a1bSJed Brown       PetscInt   j, Nk, M;
140c4762a1bSJed Brown       PetscReal *w, *v, wv;
141c4762a1bSJed Brown       PetscInt  *subset;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown       ierr = PetscDTBinomialInt(N, k, &Nk);CHKERRQ(ierr);
144c4762a1bSJed Brown       if (verbose) {ierr = PetscViewerASCIIPrintf(viewer, "k = %D:\n", k);CHKERRQ(ierr);}
145c4762a1bSJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
146c4762a1bSJed Brown       if (verbose) {ierr = PetscViewerASCIIPrintf(viewer, "(%D choose %D): %D\n", N, k, Nk);CHKERRQ(ierr);}
147c4762a1bSJed Brown 
148c4762a1bSJed Brown       /* Test subset and complement enumeration */
149c4762a1bSJed Brown       ierr = PetscMalloc1(N, &subset);CHKERRQ(ierr);
150c4762a1bSJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
151c4762a1bSJed Brown       for (j = 0; j < Nk; j++) {
152c4762a1bSJed Brown         PetscBool isOdd, isOddCheck;
153c4762a1bSJed Brown         PetscInt  jCheck, kCheck;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown         ierr = PetscDTEnumSplit(N, k, j, subset, &isOdd);CHKERRQ(ierr);
156c4762a1bSJed Brown         ierr = PetscDTPermIndex(N, subset, &kCheck, &isOddCheck);CHKERRQ(ierr);
157c4762a1bSJed Brown         if (isOddCheck != isOdd) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTEnumSplit sign does not mmatch PetscDTPermIndex sign");
158c4762a1bSJed Brown         if (verbose) {
159c4762a1bSJed Brown           PetscInt l;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "subset %D:", j);CHKERRQ(ierr);
162c4762a1bSJed Brown           for (l = 0; l < k; l++) {
163c4762a1bSJed Brown             ierr = PetscPrintf(PETSC_COMM_WORLD, " %D", subset[l]);CHKERRQ(ierr);
164c4762a1bSJed Brown           }
165c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_WORLD, " |");CHKERRQ(ierr);
166c4762a1bSJed Brown           for (l = k; l < N; l++) {
167c4762a1bSJed Brown             ierr = PetscPrintf(PETSC_COMM_WORLD, " %D", subset[l]);CHKERRQ(ierr);
168c4762a1bSJed Brown           }
169c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_WORLD, ", %s\n", isOdd ? "odd" : "even");CHKERRQ(ierr);
170c4762a1bSJed Brown         }
171c4762a1bSJed Brown         ierr = PetscDTSubsetIndex(N, k, subset, &jCheck);CHKERRQ(ierr);
172c4762a1bSJed Brown         if (jCheck != j) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "jCheck (%D) != j (%D)", jCheck, j);
173c4762a1bSJed Brown       }
174c4762a1bSJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
175c4762a1bSJed Brown       ierr = PetscFree(subset);CHKERRQ(ierr);
176c4762a1bSJed Brown 
177c4762a1bSJed Brown       /* Make a random k form */
178c4762a1bSJed Brown       ierr = PetscMalloc1(Nk, &w);CHKERRQ(ierr);
179c4762a1bSJed Brown       for (j = 0; j < Nk; j++) {ierr = PetscRandomGetValueReal(rand, &w[j]);CHKERRQ(ierr);}
180c4762a1bSJed Brown       /* Make a set of random vectors */
181c4762a1bSJed Brown       ierr = PetscMalloc1(N*k, &v);CHKERRQ(ierr);
182c4762a1bSJed Brown       for (j = 0; j < N*k; j++) {ierr = PetscRandomGetValueReal(rand, &v[j]);CHKERRQ(ierr);}
183c4762a1bSJed Brown 
184c4762a1bSJed Brown       ierr = PetscDTAltVApply(N, k, w, v, &wv);CHKERRQ(ierr);
185c4762a1bSJed Brown 
186c4762a1bSJed Brown       if (verbose) {
187c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer, "w:\n");CHKERRQ(ierr);
188c4762a1bSJed Brown         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
189c4762a1bSJed Brown         if (Nk) {ierr = PetscRealView(Nk, w, viewer);CHKERRQ(ierr);}
190c4762a1bSJed Brown         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
191c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer, "v:\n");CHKERRQ(ierr);
192c4762a1bSJed Brown         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
193c4762a1bSJed Brown         if (N*k > 0) {ierr = PetscRealView(N*k, v, viewer);CHKERRQ(ierr);}
194c4762a1bSJed Brown         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
195c4762a1bSJed Brown         ierr = PetscViewerASCIIPrintf(viewer, "w(v): %g\n", (double) wv);CHKERRQ(ierr);
196c4762a1bSJed Brown       }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown       /* sanity checks */
199c4762a1bSJed Brown       if (k == 1) { /* 1-forms are functionals (dot products) */
200c4762a1bSJed Brown         PetscInt  l;
201c4762a1bSJed Brown         PetscReal wvcheck = 0.;
202c4762a1bSJed Brown         PetscReal diff;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown         for (l = 0; l < N; l++) wvcheck += w[l] * v[l];
205c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
206c4762a1bSJed Brown         if (diff >= PETSC_SMALL * (PetscAbsReal(wv) + PetscAbsReal(wvcheck))) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "1-form / dot product equivalence: wvcheck (%g) != wv (%g)", (double) wvcheck, (double) wv);
207c4762a1bSJed Brown       }
208c4762a1bSJed Brown       if (k == N && N < 5) { /* n-forms are scaled determinants */
209c4762a1bSJed Brown         PetscReal det, wvcheck, diff;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown         switch (k) {
212c4762a1bSJed Brown         case 0:
213c4762a1bSJed Brown           det = 1.;
214c4762a1bSJed Brown           break;
215c4762a1bSJed Brown         case 1:
216c4762a1bSJed Brown           det = v[0];
217c4762a1bSJed Brown           break;
218c4762a1bSJed Brown         case 2:
219c4762a1bSJed Brown           det = v[0] * v[3] - v[1] * v[2];
220c4762a1bSJed Brown           break;
221c4762a1bSJed Brown         case 3:
222c4762a1bSJed Brown           det = v[0] * (v[4] * v[8] - v[5] * v[7]) +
223c4762a1bSJed Brown                 v[1] * (v[5] * v[6] - v[3] * v[8]) +
224c4762a1bSJed Brown                 v[2] * (v[3] * v[7] - v[4] * v[6]);
225c4762a1bSJed Brown           break;
226c4762a1bSJed Brown         case 4:
227c4762a1bSJed Brown           det = v[0] * (v[5] * (v[10] * v[15] - v[11] * v[14]) +
228c4762a1bSJed Brown                         v[6] * (v[11] * v[13] - v[ 9] * v[15]) +
229c4762a1bSJed Brown                         v[7] * (v[ 9] * v[14] - v[10] * v[13])) -
230c4762a1bSJed Brown                 v[1] * (v[4] * (v[10] * v[15] - v[11] * v[14]) +
231c4762a1bSJed Brown                         v[6] * (v[11] * v[12] - v[ 8] * v[15]) +
232c4762a1bSJed Brown                         v[7] * (v[ 8] * v[14] - v[10] * v[12])) +
233c4762a1bSJed Brown                 v[2] * (v[4] * (v[ 9] * v[15] - v[11] * v[13]) +
234c4762a1bSJed Brown                         v[5] * (v[11] * v[12] - v[ 8] * v[15]) +
235c4762a1bSJed Brown                         v[7] * (v[ 8] * v[13] - v[ 9] * v[12])) -
236c4762a1bSJed Brown                 v[3] * (v[4] * (v[ 9] * v[14] - v[10] * v[13]) +
237c4762a1bSJed Brown                         v[5] * (v[10] * v[12] - v[ 8] * v[14]) +
238c4762a1bSJed Brown                         v[6] * (v[ 8] * v[13] - v[ 9] * v[12]));
239c4762a1bSJed Brown           break;
240c4762a1bSJed Brown         default:
241c4762a1bSJed Brown           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "invalid k");
242c4762a1bSJed Brown         }
243c4762a1bSJed Brown         wvcheck = det * w[0];
244c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
245c4762a1bSJed Brown         if (diff >= PETSC_SMALL * (PetscAbsReal(wv) + PetscAbsReal(wvcheck))) SETERRQ3(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "n-form / determinant equivalence: wvcheck (%g) != wv (%g) %g", (double) wvcheck, (double) wv, (double) diff);
246c4762a1bSJed Brown       }
247c4762a1bSJed Brown       if (k > 0) { /* k-forms are linear in each component */
248c4762a1bSJed Brown         PetscReal alpha;
249c4762a1bSJed Brown         PetscReal *x, *axv, wx, waxv, waxvcheck;
250c4762a1bSJed Brown         PetscReal diff;
251c4762a1bSJed Brown         PetscReal rj;
252c4762a1bSJed Brown         PetscInt  l;
253c4762a1bSJed Brown 
254c4762a1bSJed Brown         ierr = PetscMalloc2(N * k, &x, N * k, &axv);CHKERRQ(ierr);
255c4762a1bSJed Brown         ierr = PetscRandomGetValueReal(rand, &alpha);CHKERRQ(ierr);
256c4762a1bSJed Brown         ierr = PetscRandomSetInterval(rand, 0, k);CHKERRQ(ierr);
257c4762a1bSJed Brown         ierr = PetscRandomGetValueReal(rand, &rj);CHKERRQ(ierr);
258c4762a1bSJed Brown         j = (PetscInt) rj;
259c4762a1bSJed Brown         ierr = PetscRandomSetInterval(rand, -1., 1.);CHKERRQ(ierr);
260c4762a1bSJed Brown         for (l = 0; l < N*k; l++) x[l] = v[l];
261c4762a1bSJed Brown         for (l = 0; l < N*k; l++) axv[l] = v[l];
262c4762a1bSJed Brown         for (l = 0; l < N; l++) {
263c4762a1bSJed Brown           PetscReal val;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown           ierr = PetscRandomGetValueReal(rand, &val);CHKERRQ(ierr);
266c4762a1bSJed Brown           x[j * N + l] = val;
267c4762a1bSJed Brown           axv[j * N + l] += alpha * val;
268c4762a1bSJed Brown         }
269c4762a1bSJed Brown         ierr = PetscDTAltVApply(N, k, w, x, &wx);CHKERRQ(ierr);
270c4762a1bSJed Brown         ierr = PetscDTAltVApply(N, k, w, axv, &waxv);CHKERRQ(ierr);
271c4762a1bSJed Brown         waxvcheck = alpha * wx + wv;
272c4762a1bSJed Brown         diff = waxv - waxvcheck;
273c4762a1bSJed Brown         if (PetscAbsReal(diff) > 10. * PETSC_SMALL * (PetscAbsReal(waxv) + PetscAbsReal(waxvcheck))) SETERRQ3(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "linearity check: component %D, waxvcheck (%g) != waxv (%g)", j, (double) waxvcheck, (double) waxv);
274c4762a1bSJed Brown         ierr = PetscFree2(x,axv);CHKERRQ(ierr);
275c4762a1bSJed Brown       }
276c4762a1bSJed Brown       if (k > 1) { /* k-forms are antisymmetric */
277c4762a1bSJed Brown         PetscReal rj, rl, *swapv, wswapv, diff;
278c4762a1bSJed Brown         PetscInt  l, m;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown         ierr = PetscRandomSetInterval(rand, 0, k);CHKERRQ(ierr);
281c4762a1bSJed Brown         ierr = PetscRandomGetValueReal(rand, &rj);CHKERRQ(ierr);
282c4762a1bSJed Brown         j = (PetscInt) rj;
283c4762a1bSJed Brown         l = j;
284c4762a1bSJed Brown         while (l == j) {
285c4762a1bSJed Brown           ierr = PetscRandomGetValueReal(rand, &rl);CHKERRQ(ierr);
286c4762a1bSJed Brown           l = (PetscInt) rl;
287c4762a1bSJed Brown         }
288c4762a1bSJed Brown         ierr = PetscRandomSetInterval(rand, -1., 1.);CHKERRQ(ierr);
289c4762a1bSJed Brown         ierr = PetscMalloc1(N * k, &swapv);CHKERRQ(ierr);
290c4762a1bSJed Brown         for (m = 0; m < N * k; m++) swapv[m] = v[m];
291c4762a1bSJed Brown         for (m = 0; m < N; m++) {
292c4762a1bSJed Brown           swapv[j * N + m] = v[l * N + m];
293c4762a1bSJed Brown           swapv[l * N + m] = v[j * N + m];
294c4762a1bSJed Brown         }
295c4762a1bSJed Brown         ierr = PetscDTAltVApply(N, k, w, swapv, &wswapv);CHKERRQ(ierr);
296c4762a1bSJed Brown         diff = PetscAbsReal(wswapv + wv);
297c4762a1bSJed Brown         if (diff > PETSC_SMALL * (PetscAbsReal(wswapv) + PetscAbsReal(wv))) SETERRQ4(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "antisymmetry check: components %D & %D, wswapv (%g) != -wv (%g)", j, l, (double) wswapv, (double) wv);
298c4762a1bSJed Brown         ierr = PetscFree(swapv);CHKERRQ(ierr);
299c4762a1bSJed Brown       }
300c4762a1bSJed Brown       for (j = 0; j <= k && j + k <= N; j++) { /* wedge product */
301c4762a1bSJed Brown         PetscInt   Nj, Njk, l, JKj;
302c4762a1bSJed Brown         PetscReal *u, *uWw, *uWwcheck, *uWwmat, *x, *xsplit, uWwx, uWwxcheck, diff, norm;
303c4762a1bSJed Brown         PetscInt  *split;
304c4762a1bSJed Brown 
305c4762a1bSJed Brown         if (verbose) {ierr = PetscViewerASCIIPrintf(viewer, "wedge j = %D:\n", j);CHKERRQ(ierr);}
306c4762a1bSJed Brown         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
307c4762a1bSJed Brown         ierr = PetscDTBinomialInt(N, j,   &Nj);CHKERRQ(ierr);
308c4762a1bSJed Brown         ierr = PetscDTBinomialInt(N, j+k, &Njk);CHKERRQ(ierr);
309c4762a1bSJed Brown         ierr = PetscMalloc4(Nj, &u, Njk, &uWw, N*(j+k), &x, N*(j+k), &xsplit);CHKERRQ(ierr);
310c4762a1bSJed Brown         ierr = PetscMalloc1(j+k,&split);CHKERRQ(ierr);
311c4762a1bSJed Brown         for (l = 0; l < Nj; l++) {ierr = PetscRandomGetValueReal(rand, &u[l]);CHKERRQ(ierr);}
312c4762a1bSJed Brown         for (l = 0; l < N*(j+k); l++) {ierr = PetscRandomGetValueReal(rand, &x[l]);CHKERRQ(ierr);}
313c4762a1bSJed Brown         ierr = PetscDTAltVWedge(N, j, k, u, w, uWw);CHKERRQ(ierr);
314c4762a1bSJed Brown         ierr = PetscDTAltVApply(N, j+k, uWw, x, &uWwx);CHKERRQ(ierr);
315c4762a1bSJed Brown         if (verbose) {
316c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "u:\n");CHKERRQ(ierr);
317c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
318c4762a1bSJed Brown           if (Nj) {ierr = PetscRealView(Nj, u, viewer);CHKERRQ(ierr);}
319c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
320c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "u wedge w:\n");CHKERRQ(ierr);
321c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
322c4762a1bSJed Brown           if (Njk) {ierr = PetscRealView(Njk, uWw, viewer);CHKERRQ(ierr);}
323c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
324c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "x:\n");CHKERRQ(ierr);
325c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
326c4762a1bSJed Brown           if (N*(j+k) > 0) {ierr = PetscRealView(N*(j+k), x, viewer);CHKERRQ(ierr);}
327c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
328c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "u wedge w(x): %g\n", (double) uWwx);CHKERRQ(ierr);
329c4762a1bSJed Brown         }
330c4762a1bSJed Brown         /* verify wedge formula */
331c4762a1bSJed Brown         uWwxcheck = 0.;
332c4762a1bSJed Brown         ierr = PetscDTBinomialInt(j+k, j, &JKj);CHKERRQ(ierr);
333c4762a1bSJed Brown         for (l = 0; l < JKj; l++) {
334c4762a1bSJed Brown           PetscBool isOdd;
335c4762a1bSJed Brown           PetscReal ux, wx;
336c4762a1bSJed Brown           PetscInt  m, p;
337c4762a1bSJed Brown 
338c4762a1bSJed Brown           ierr = PetscDTEnumSplit(j+k, j, l, split, &isOdd);CHKERRQ(ierr);
339c4762a1bSJed Brown           for (m = 0; m < j+k; m++) {for (p = 0; p < N; p++) {xsplit[m * N + p] = x[split[m] * N + p];}}
340c4762a1bSJed Brown           ierr = PetscDTAltVApply(N, j, u, xsplit, &ux);CHKERRQ(ierr);
341c4762a1bSJed Brown           ierr = PetscDTAltVApply(N, k, w, &xsplit[j*N], &wx);CHKERRQ(ierr);
342c4762a1bSJed Brown           uWwxcheck += isOdd ? -(ux * wx) : (ux * wx);
343c4762a1bSJed Brown         }
344c4762a1bSJed Brown         diff = PetscAbsReal(uWwx - uWwxcheck);
345c4762a1bSJed Brown         if (diff > 10. * PETSC_SMALL * (PetscAbsReal(uWwx) + PetscAbsReal(uWwxcheck))) SETERRQ4(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge check: forms %D & %D, uWwxcheck (%g) != uWwx (%g)", j, k, (double) uWwxcheck, (double) uWwx);
346c4762a1bSJed Brown         ierr = PetscFree(split);CHKERRQ(ierr);
347c4762a1bSJed Brown         ierr = PetscMalloc2(Nk * Njk, &uWwmat, Njk, &uWwcheck);CHKERRQ(ierr);
348c4762a1bSJed Brown         ierr = PetscDTAltVWedgeMatrix(N, j, k, u, uWwmat);CHKERRQ(ierr);
349c4762a1bSJed Brown         if (verbose) {
350c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "(u wedge):\n");CHKERRQ(ierr);
351c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
352c4762a1bSJed Brown           if ((Nk * Njk) > 0) {ierr = PetscRealView(Nk * Njk, uWwmat, viewer);CHKERRQ(ierr);}
353c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
354c4762a1bSJed Brown         }
355c4762a1bSJed Brown         diff = 0.;
356c4762a1bSJed Brown         norm = 0.;
357c4762a1bSJed Brown         for (l = 0; l < Njk; l++) {
358c4762a1bSJed Brown           PetscInt  m;
359c4762a1bSJed Brown           PetscReal sum = 0.;
360c4762a1bSJed Brown 
361c4762a1bSJed Brown           for (m = 0; m < Nk; m++) sum += uWwmat[l * Nk + m] * w[m];
362c4762a1bSJed Brown           uWwcheck[l] = sum;
363c4762a1bSJed Brown           diff += PetscSqr(uWwcheck[l] - uWw[l]);
364c4762a1bSJed Brown           norm += PetscSqr(uWwcheck[l]) + PetscSqr(uWw[l]);
365c4762a1bSJed Brown         }
366c4762a1bSJed Brown         diff = PetscSqrtReal(diff);
367c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
368c4762a1bSJed Brown         if (diff > PETSC_SMALL * norm) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "wedge matrix check: wedge matrix application does not match wedge direct application");
369c4762a1bSJed Brown         ierr = PetscFree2(uWwmat, uWwcheck);CHKERRQ(ierr);
370c4762a1bSJed Brown         ierr = PetscFree4(u, uWw, x, xsplit);CHKERRQ(ierr);
371c4762a1bSJed Brown         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
372c4762a1bSJed Brown       }
373c4762a1bSJed Brown       for (M = PetscMax(1,k); M <= N; M++) { /* pullback */
374c4762a1bSJed Brown         PetscReal   *L, *u, *x;
375c4762a1bSJed Brown         PetscInt     Mk, l;
376c4762a1bSJed Brown 
377c4762a1bSJed Brown         ierr = PetscDTBinomialInt(M, k, &Mk);CHKERRQ(ierr);
378c4762a1bSJed Brown         ierr = PetscMalloc3(M*N, &L, Mk, &u, M*k, &x);CHKERRQ(ierr);
379c4762a1bSJed Brown         for (l = 0; l < M*N; l++) {ierr = PetscRandomGetValueReal(rand, &L[l]);CHKERRQ(ierr);}
380c4762a1bSJed Brown         for (l = 0; l < Mk; l++) {ierr = PetscRandomGetValueReal(rand, &u[l]);CHKERRQ(ierr);}
381c4762a1bSJed Brown         for (l = 0; l < M*k; l++) {ierr = PetscRandomGetValueReal(rand, &x[l]);CHKERRQ(ierr);}
382c4762a1bSJed Brown         if (verbose) {ierr = PetscViewerASCIIPrintf(viewer, "pullback M = %D:\n", M);CHKERRQ(ierr);}
383c4762a1bSJed Brown         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
384c4762a1bSJed Brown         ierr = CheckPullback(M, N, L, k, w, x, verbose, viewer);CHKERRQ(ierr);
385c4762a1bSJed Brown         if (M != N) {ierr = CheckPullback(N, M, L, k, u, v, PETSC_FALSE, viewer);CHKERRQ(ierr);}
386c4762a1bSJed Brown         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
387c4762a1bSJed Brown         if ((k % N) && (N > 1)) {
388c4762a1bSJed Brown           if (verbose) {ierr = PetscViewerASCIIPrintf(viewer, "negative pullback M = %D:\n", M);CHKERRQ(ierr);}
389c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
390c4762a1bSJed Brown           ierr = CheckPullback(M, N, L, -k, w, x, verbose, viewer);CHKERRQ(ierr);
391c4762a1bSJed Brown           if (M != N) {ierr = CheckPullback(N, M, L, -k, u, v, PETSC_FALSE, viewer);CHKERRQ(ierr);}
392c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
393c4762a1bSJed Brown         }
394c4762a1bSJed Brown         ierr = PetscFree3(L, u, x);CHKERRQ(ierr);
395c4762a1bSJed Brown       }
396c4762a1bSJed Brown       if (k > 0) { /* Interior */
397c4762a1bSJed Brown         PetscInt    Nkm, l, m;
398c4762a1bSJed Brown         PetscReal  *wIntv0, *wIntv0check, wvcheck, diff, diffMat, normMat;
399c4762a1bSJed Brown         PetscReal  *intv0mat, *matcheck;
400c4762a1bSJed Brown         PetscInt  (*indices)[3];
401c4762a1bSJed Brown 
402c4762a1bSJed Brown         ierr = PetscDTBinomialInt(N, k-1, &Nkm);CHKERRQ(ierr);
403c4762a1bSJed Brown         ierr = PetscMalloc5(Nkm, &wIntv0, Nkm, &wIntv0check, Nk * Nkm, &intv0mat, Nk * Nkm, &matcheck, Nk * k, &indices);CHKERRQ(ierr);
404c4762a1bSJed Brown         ierr = PetscDTAltVInterior(N, k, w, v, wIntv0);CHKERRQ(ierr);
405c4762a1bSJed Brown         ierr = PetscDTAltVInteriorMatrix(N, k, v, intv0mat);CHKERRQ(ierr);
406c4762a1bSJed Brown         ierr = PetscDTAltVInteriorPattern(N, k, indices);CHKERRQ(ierr);
407c4762a1bSJed Brown         if (verbose) {
408c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "interior product matrix pattern:\n");CHKERRQ(ierr);
409c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
410c4762a1bSJed Brown           for (l = 0; l < Nk * k; l++) {
411c4762a1bSJed Brown             PetscInt row = indices[l][0];
412c4762a1bSJed Brown             PetscInt col = indices[l][1];
413c4762a1bSJed Brown             PetscInt x   = indices[l][2];
414c4762a1bSJed Brown 
415c4762a1bSJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"intV[%D,%D] = %sV[%D]\n", row, col, x < 0 ? "-" : " ", x < 0 ? -(x + 1) : x);CHKERRQ(ierr);
416c4762a1bSJed Brown           }
417c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
418c4762a1bSJed Brown         }
419c4762a1bSJed Brown         for (l = 0; l < Nkm * Nk; l++) matcheck[l] = 0.;
420c4762a1bSJed Brown         for (l = 0; l < Nk * k; l++) {
421c4762a1bSJed Brown           PetscInt row = indices[l][0];
422c4762a1bSJed Brown           PetscInt col = indices[l][1];
423c4762a1bSJed Brown           PetscInt x   = indices[l][2];
424c4762a1bSJed Brown 
425c4762a1bSJed Brown           if (x < 0) {
426c4762a1bSJed Brown             matcheck[row * Nk + col] = -v[-(x+1)];
427c4762a1bSJed Brown           } else {
428c4762a1bSJed Brown             matcheck[row * Nk + col] = v[x];
429c4762a1bSJed Brown           }
430c4762a1bSJed Brown         }
431c4762a1bSJed Brown         diffMat = 0.;
432c4762a1bSJed Brown         normMat = 0.;
433c4762a1bSJed Brown         for (l = 0; l < Nkm * Nk; l++) {
434c4762a1bSJed Brown           diffMat += PetscSqr(PetscAbsReal(matcheck[l] - intv0mat[l]));
435c4762a1bSJed Brown           normMat += PetscSqr(matcheck[l]) + PetscSqr(intv0mat[l]);
436c4762a1bSJed Brown         }
437c4762a1bSJed Brown         diffMat = PetscSqrtReal(diffMat);
438c4762a1bSJed Brown         normMat = PetscSqrtReal(normMat);
439c4762a1bSJed Brown         if (diffMat > PETSC_SMALL * normMat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: matrix pattern does not match matrix");
440c4762a1bSJed Brown         diffMat = 0.;
441c4762a1bSJed Brown         normMat = 0.;
442c4762a1bSJed Brown         for (l = 0; l < Nkm; l++) {
443c4762a1bSJed Brown           PetscReal sum = 0.;
444c4762a1bSJed Brown 
445c4762a1bSJed Brown           for (m = 0; m < Nk; m++) sum += intv0mat[l * Nk + m] * w[m];
446c4762a1bSJed Brown           wIntv0check[l] = sum;
447c4762a1bSJed Brown 
448c4762a1bSJed Brown           diffMat += PetscSqr(PetscAbsReal(wIntv0check[l] - wIntv0[l]));
449c4762a1bSJed Brown           normMat += PetscSqr(wIntv0check[l]) + PetscSqr(wIntv0[l]);
450c4762a1bSJed Brown         }
451c4762a1bSJed Brown         diffMat = PetscSqrtReal(diffMat);
452c4762a1bSJed Brown         normMat = PetscSqrtReal(normMat);
453c4762a1bSJed Brown         if (diffMat > PETSC_SMALL * normMat) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: application does not match matrix");
454c4762a1bSJed Brown         if (verbose) {
455c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "(w int v_0):\n");CHKERRQ(ierr);
456c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
457c4762a1bSJed Brown           if (Nkm) {ierr = PetscRealView(Nkm, wIntv0, viewer);CHKERRQ(ierr);}
458c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
459c4762a1bSJed Brown 
460c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "(int v_0):\n");CHKERRQ(ierr);
461c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
462c4762a1bSJed Brown           if (Nk * Nkm > 0) {ierr = PetscRealView(Nk * Nkm, intv0mat, viewer);CHKERRQ(ierr);}
463c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
464c4762a1bSJed Brown         }
465c4762a1bSJed Brown         ierr = PetscDTAltVApply(N, k - 1, wIntv0, &v[N], &wvcheck);CHKERRQ(ierr);
466c4762a1bSJed Brown         diff = PetscSqrtReal(PetscSqr(wvcheck - wv));
467c4762a1bSJed Brown         if (diff >= PETSC_SMALL * (PetscAbsReal(wv) + PetscAbsReal(wvcheck))) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Interior product check: (w Int v0)(v_rem) (%g) != w(v) (%g)", (double) wvcheck, (double) wv);
468c4762a1bSJed Brown         ierr = PetscFree5(wIntv0,wIntv0check,intv0mat,matcheck,indices);CHKERRQ(ierr);
469c4762a1bSJed Brown       }
470c4762a1bSJed Brown       if (k >= N - k) { /* Hodge star */
471c4762a1bSJed Brown         PetscReal *u, *starw, *starstarw, wu, starwdotu;
472c4762a1bSJed Brown         PetscReal diff, norm;
473c4762a1bSJed Brown         PetscBool isOdd;
474c4762a1bSJed Brown         PetscInt l;
475c4762a1bSJed Brown 
476c4762a1bSJed Brown         isOdd = (PetscBool) ((k * (N - k)) & 1);
477c4762a1bSJed Brown         ierr = PetscMalloc3(Nk, &u, Nk, &starw, Nk, &starstarw);CHKERRQ(ierr);
478c4762a1bSJed Brown         ierr = PetscDTAltVStar(N, k, 1, w, starw);CHKERRQ(ierr);
479c4762a1bSJed Brown         ierr = PetscDTAltVStar(N, N-k, 1, starw, starstarw);CHKERRQ(ierr);
480c4762a1bSJed Brown         if (verbose) {
481c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "star w:\n");CHKERRQ(ierr);
482c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
483c4762a1bSJed Brown           if (Nk) {ierr = PetscRealView(Nk, starw, viewer);CHKERRQ(ierr);}
484c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
485c4762a1bSJed Brown 
486c4762a1bSJed Brown           ierr = PetscViewerASCIIPrintf(viewer, "star star w:\n");CHKERRQ(ierr);
487c4762a1bSJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
488c4762a1bSJed Brown           if (Nk) {ierr = PetscRealView(Nk, starstarw, viewer);CHKERRQ(ierr);}
489c4762a1bSJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
490c4762a1bSJed Brown         }
491c4762a1bSJed Brown         for (l = 0; l < Nk; l++) {ierr = PetscRandomGetValueReal(rand,&u[l]);CHKERRQ(ierr);}
492c4762a1bSJed Brown         ierr = PetscDTAltVWedge(N, k, N - k, w, u, &wu);CHKERRQ(ierr);
493c4762a1bSJed Brown         starwdotu = 0.;
494c4762a1bSJed Brown         for (l = 0; l < Nk; l++) starwdotu += starw[l] * u[l];
495c4762a1bSJed Brown         diff = PetscAbsReal(wu - starwdotu);
496c4762a1bSJed Brown         if (diff > PETSC_SMALL * (PetscAbsReal(wu) + PetscAbsReal(starwdotu))) SETERRQ2(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: (star w, u) (%g) != (w wedge u) (%g)", (double) starwdotu, (double) wu);
497c4762a1bSJed Brown 
498c4762a1bSJed Brown         diff = 0.;
499c4762a1bSJed Brown         norm = 0.;
500c4762a1bSJed Brown         for (l = 0; l < Nk; l++) {
501c4762a1bSJed Brown           diff += PetscSqr(w[l] - (isOdd ? -starstarw[l] : starstarw[l]));
502c4762a1bSJed Brown           norm += PetscSqr(w[l]) + PetscSqr(starstarw[l]);
503c4762a1bSJed Brown         }
504c4762a1bSJed Brown         diff = PetscSqrtReal(diff);
505c4762a1bSJed Brown         norm = PetscSqrtReal(norm);
506c4762a1bSJed Brown         if (diff > PETSC_SMALL * norm) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Hodge star check: star(star(w)) != (-1)^(N*(N-k)) w");
507c4762a1bSJed Brown         ierr = PetscFree3(u, starw, starstarw);CHKERRQ(ierr);
508c4762a1bSJed Brown       }
509c4762a1bSJed Brown       ierr = PetscFree(v);CHKERRQ(ierr);
510c4762a1bSJed Brown       ierr = PetscFree(w);CHKERRQ(ierr);
511c4762a1bSJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
512c4762a1bSJed Brown     }
513c4762a1bSJed Brown     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
514c4762a1bSJed Brown   }
515c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
516c4762a1bSJed Brown   ierr = PetscFinalize();
517c4762a1bSJed Brown   return ierr;
518c4762a1bSJed Brown }
519c4762a1bSJed Brown 
520c4762a1bSJed Brown /*TEST
521c4762a1bSJed Brown   test:
522c4762a1bSJed Brown     suffix: 1234
523c4762a1bSJed Brown     args: -verbose
524c4762a1bSJed Brown   test:
525c4762a1bSJed Brown     suffix: 56
526c4762a1bSJed Brown     args: -N 5,6
527c4762a1bSJed Brown TEST*/
528