xref: /petsc/src/dm/impls/plex/tests/ex6.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests for DMLabel lookup\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmplex.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown typedef struct {
6*c4762a1bSJed Brown   PetscInt  debug;        /* The debugging level */
7*c4762a1bSJed Brown   PetscInt  pStart, pEnd; /* The label chart */
8*c4762a1bSJed Brown   PetscInt  numStrata;    /* The number of label strata */
9*c4762a1bSJed Brown   PetscReal fill;         /* Percentage of label to fill */
10*c4762a1bSJed Brown   PetscInt  size;         /* The number of set values */
11*c4762a1bSJed Brown } AppCtx;
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
14*c4762a1bSJed Brown {
15*c4762a1bSJed Brown   PetscErrorCode ierr;
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown   PetscFunctionBegin;
18*c4762a1bSJed Brown   options->debug     = 0;
19*c4762a1bSJed Brown   options->pStart    = 0;
20*c4762a1bSJed Brown   options->pEnd      = 1000;
21*c4762a1bSJed Brown   options->numStrata = 5;
22*c4762a1bSJed Brown   options->fill      = 0.10;
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex6.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-num_strata", "The number of label values", "ex6.c", options->numStrata, &options->numStrata, NULL,0);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-pend", "The label point limit", "ex6.c", options->pEnd, &options->pEnd, NULL,0);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = PetscOptionsReal("-fill", "The percentage of label chart to set", "ex6.c", options->fill, &options->fill, NULL);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
30*c4762a1bSJed Brown   PetscFunctionReturn(0);
31*c4762a1bSJed Brown }
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown PetscErrorCode TestSetup(DMLabel label, AppCtx *user)
34*c4762a1bSJed Brown {
35*c4762a1bSJed Brown   PetscRandom    r;
36*c4762a1bSJed Brown   PetscInt       n = (PetscInt) (user->fill*(user->pEnd - user->pStart)), i;
37*c4762a1bSJed Brown   PetscErrorCode ierr;
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown   PetscFunctionBegin;
40*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);/* -random_type <> */
42*c4762a1bSJed Brown   ierr = PetscRandomSetInterval(r, user->pStart, user->pEnd);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = PetscRandomSetSeed(r, 123456789L);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = PetscRandomSeed(r);CHKERRQ(ierr);
45*c4762a1bSJed Brown   user->size = 0;
46*c4762a1bSJed Brown   for(i = 0; i < n; ++i) {
47*c4762a1bSJed Brown     PetscReal p;
48*c4762a1bSJed Brown     PetscInt  val;
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown     ierr = PetscRandomGetValueReal(r, &p);CHKERRQ(ierr);
51*c4762a1bSJed Brown     ierr = DMLabelGetValue(label, (PetscInt) p, &val);CHKERRQ(ierr);
52*c4762a1bSJed Brown     if (val < 0) {
53*c4762a1bSJed Brown       ++user->size;
54*c4762a1bSJed Brown       ierr = DMLabelSetValue(label, (PetscInt) p, i % user->numStrata);CHKERRQ(ierr);
55*c4762a1bSJed Brown     }
56*c4762a1bSJed Brown   }
57*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = DMLabelCreateIndex(label, user->pStart, user->pEnd);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Created label with chart [%D, %D) and set %D values\n", user->pStart, user->pEnd, user->size);CHKERRQ(ierr);
60*c4762a1bSJed Brown   PetscFunctionReturn(0);
61*c4762a1bSJed Brown }
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown PetscErrorCode TestLookup(DMLabel label, AppCtx *user)
64*c4762a1bSJed Brown {
65*c4762a1bSJed Brown   const PetscInt pStart = user->pStart;
66*c4762a1bSJed Brown   const PetscInt pEnd   = user->pEnd;
67*c4762a1bSJed Brown   PetscInt       p, n = 0;
68*c4762a1bSJed Brown   PetscErrorCode ierr;
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   PetscFunctionBegin;
71*c4762a1bSJed Brown   for(p = pStart; p < pEnd; ++p) {
72*c4762a1bSJed Brown     PetscInt  val;
73*c4762a1bSJed Brown     PetscBool has;
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown     ierr = DMLabelGetValue(label, p, &val);CHKERRQ(ierr);
76*c4762a1bSJed Brown     ierr = DMLabelHasPoint(label, p, &has);CHKERRQ(ierr);
77*c4762a1bSJed Brown     if (((val >= 0) && !has) || ((val < 0) && has)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Label value %D does not match contains check %D for point %D", val, (PetscInt) has, p);
78*c4762a1bSJed Brown     if (has) ++n;
79*c4762a1bSJed Brown   }
80*c4762a1bSJed Brown   if (n != user->size) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of label points detected %D does not match number set %D", n, user->size);
81*c4762a1bSJed Brown   /* Also put in timing code */
82*c4762a1bSJed Brown   PetscFunctionReturn(0);
83*c4762a1bSJed Brown }
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown PetscErrorCode TestClear(DMLabel label, AppCtx *user)
86*c4762a1bSJed Brown {
87*c4762a1bSJed Brown   PetscInt       pStart = user->pStart, pEnd = user->pEnd, p;
88*c4762a1bSJed Brown   PetscInt       defaultValue;
89*c4762a1bSJed Brown   PetscErrorCode ierr;
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown   PetscFunctionBegin;
92*c4762a1bSJed Brown   ierr = DMLabelGetDefaultValue(label,&defaultValue);CHKERRQ(ierr);
93*c4762a1bSJed Brown   for (p = pStart; p < pEnd; p++) {
94*c4762a1bSJed Brown     PetscInt  val;
95*c4762a1bSJed Brown     PetscBool hasPoint;
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown     ierr = DMLabelGetValue(label,p,&val);CHKERRQ(ierr);
98*c4762a1bSJed Brown     if (val != defaultValue) {
99*c4762a1bSJed Brown       ierr = DMLabelClearValue(label,p,val);CHKERRQ(ierr);
100*c4762a1bSJed Brown     }
101*c4762a1bSJed Brown     ierr = DMLabelGetValue(label,p,&val);CHKERRQ(ierr);
102*c4762a1bSJed Brown     ierr = DMLabelHasPoint(label,p,&hasPoint);CHKERRQ(ierr);
103*c4762a1bSJed Brown     if (val != defaultValue) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected default value %D after clearing point %D, got %D",defaultValue,p,val);
104*c4762a1bSJed Brown     if (hasPoint) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Label contains %D after clearing",p);
105*c4762a1bSJed Brown   }
106*c4762a1bSJed Brown   PetscFunctionReturn(0);
107*c4762a1bSJed Brown }
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown int main(int argc, char **argv)
110*c4762a1bSJed Brown {
111*c4762a1bSJed Brown   DMLabel        label;
112*c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
113*c4762a1bSJed Brown   PetscErrorCode ierr;
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
116*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF, "Test Label", &label);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = TestSetup(label, &user);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = TestLookup(label, &user);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = TestClear(label,&user);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscFinalize();
123*c4762a1bSJed Brown   return ierr;
124*c4762a1bSJed Brown }
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown /*TEST
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown   test:
129*c4762a1bSJed Brown     suffix: 0
130*c4762a1bSJed Brown     args: -malloc_dump
131*c4762a1bSJed Brown   test:
132*c4762a1bSJed Brown     suffix: 1
133*c4762a1bSJed Brown     args: -malloc_dump -pend 10000
134*c4762a1bSJed Brown   test:
135*c4762a1bSJed Brown     suffix: 2
136*c4762a1bSJed Brown     args: -malloc_dump -pend 10000 -fill 0.05
137*c4762a1bSJed Brown   test:
138*c4762a1bSJed Brown     suffix: 3
139*c4762a1bSJed Brown     args: -malloc_dump -pend 10000 -fill 0.25
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown TEST*/
142