Actual source code: ex8.c
1: const char help[] = "Tests PetscDTBaryToIndex(), PetscDTIndexToBary(), PetscDTIndexToGradedOrder() and PetscDTGradedOrderToIndex()";
3: #include <petsc/private/petscimpl.h>
4: #include <petsc/private/dtimpl.h>
5: #include <petsc/private/petscfeimpl.h>
7: int main(int argc, char **argv)
8: {
9: PetscInt d, n, maxdim = 4;
10: PetscInt *btupprev, *btup;
11: PetscInt *gtup;
13: PetscFunctionBeginUser;
14: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
15: PetscCall(PetscMalloc3(maxdim + 1, &btup, maxdim + 1, &btupprev, maxdim, >up));
16: for (d = 0; d <= maxdim; d++) {
17: for (n = 0; n <= d + 2; n++) {
18: PetscInt j, k, Nk, kchk;
20: PetscCall(PetscDTBinomialInt(d + n, d, &Nk));
21: for (k = 0; k < Nk; k++) {
22: PetscInt sum;
24: PetscCall(PetscDTIndexToBary(d + 1, n, k, btup));
25: for (j = 0, sum = 0; j < d + 1; j++) {
26: PetscCheck(btup[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " negative entry", d, n, k);
27: sum += btup[j];
28: }
29: PetscCheck(sum == n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " incorrect sum", d, n, k);
30: PetscCall(PetscDTBaryToIndex(d + 1, n, btup, &kchk));
31: PetscCheck(kchk == k, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTBaryToIndex, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " mismatch", d, n, k);
32: if (k) {
33: j = d;
34: while (j >= 0 && btup[j] == btupprev[j]) j--;
35: PetscCheck(j >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " equal to previous", d, n, k);
36: PetscCheck(btup[j] >= btupprev[j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " less to previous", d, n, k);
37: } else PetscCall(PetscArraycpy(btupprev, btup, d + 1));
38: PetscCall(PetscDTIndexToGradedOrder(d, Nk - 1 - k, gtup));
39: PetscCall(PetscDTGradedOrderToIndex(d, gtup, &kchk));
40: PetscCheck(kchk == Nk - 1 - k, PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTGradedOrderToIndex, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " mismatch", d, n, Nk - 1 - k);
41: for (j = 0; j < d; j++) PetscCheck(gtup[j] == btup[d - 1 - j], PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToGradedOrder, d = %" PetscInt_FMT ", n = %" PetscInt_FMT ", k = %" PetscInt_FMT " incorrect", d, n, Nk - 1 - k);
42: }
43: }
44: }
45: PetscCall(PetscFree3(btup, btupprev, gtup));
46: PetscCall(PetscFinalize());
47: return 0;
48: }
50: /*TEST
52: test:
54: TEST*/