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