xref: /petsc/src/dm/dt/tests/ex8.c (revision 4e8afd12df985612b40e26aef947da2f566aee4e)
1 const char help[] = "Tests PetscDTBaryToIndex() and PetscDTIndexToBary()";
2 
3 #include <petsc/private/petscimpl.h>
4 #include <petsc/private/dtimpl.h>
5 #include <petsc/private/petscfeimpl.h>
6 
7 int main(int argc, char **argv)
8 {
9   PetscInt       d, n, maxdim = 4;
10   PetscInt       *btupprev, *btup;
11   PetscErrorCode ierr;
12 
13   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
14   ierr = PetscMalloc2(maxdim + 1, &btup, maxdim + 1, &btupprev);CHKERRQ(ierr);
15   for (d = 0; d <= maxdim; d++) {
16     for (n = 0; n <= d + 2; n++) {
17       PetscInt j, k, Nk, kchk;
18 
19       ierr = PetscDTBinomialInt(d + n, d, &Nk);CHKERRQ(ierr);
20       for (k = 0; k < Nk; k++) {
21         PetscInt sum;
22 
23         ierr = PetscDTIndexToBary(d + 1, n, k, btup);CHKERRQ(ierr);
24         for (j = 0, sum = 0; j < d + 1; j++) {
25           if (btup[j] < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %d, n = %d, k = %d negative entry\n", d, n, k);
26           sum += btup[j];
27         }
28         if (sum != n) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %d, n = %d, k = %d incorrect sum\n", d, n, k);
29         ierr = PetscDTBaryToIndex(d + 1, n, btup, &kchk);CHKERRQ(ierr);
30         if (kchk != k) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTBaryToIndex, d = %d, n = %d, k = %d mismatch\n", d, n, k);
31         if (k) {
32           j = d;
33           while (j >= 0 && btup[j] == btupprev[j]) j--;
34           if (j < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %d, n = %d, k = %d equal to previous\n", d, n, k);
35           if (btup[j] < btupprev[j]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PetscDTIndexToBary, d = %d, n = %d, k = %d less to previous\n", d, n, k);
36         } else {
37           ierr = PetscArraycpy(btupprev, btup, d + 1);CHKERRQ(ierr);
38         }
39       }
40     }
41   }
42   ierr = PetscFree2(btup, btupprev);CHKERRQ(ierr);
43   ierr = PetscFinalize();
44   return ierr;
45 }
46 
47 /*TEST
48 
49   test:
50 
51 TEST*/
52