1*c0888a1bSJDBetteridge 2*c0888a1bSJDBetteridge static char help[] = "Tests `PetscGarbageKeySortedIntersect()`\n\n"; 3*c0888a1bSJDBetteridge 4*c0888a1bSJDBetteridge #include <petscsys.h> 5*c0888a1bSJDBetteridge #include <petsc/private/garbagecollector.h> 6*c0888a1bSJDBetteridge 7*c0888a1bSJDBetteridge /* This program tests `PetscGarbageKeySortedIntersect(), which is the 8*c0888a1bSJDBetteridge public (MPI) interface to 9*c0888a1bSJDBetteridge `PetscErrorCode GarbageKeySortedIntersect_Private()`. 10*c0888a1bSJDBetteridge Sets are sent packed in arrays, with the first entry as the number of 11*c0888a1bSJDBetteridge set elements and the sets the remaining elements. This is because the 12*c0888a1bSJDBetteridge MPI reduction operation must have the call signature: 13*c0888a1bSJDBetteridge void PetscGarbageKeySortedIntersect(void *inset, void *inoutset, PetscMPIInt *length, MPI_Datatype *dtype) 14*c0888a1bSJDBetteridge This is a thin wrapper for the private routine: 15*c0888a1bSJDBetteridge PetscErrorCode GarbageKeySortedIntersect_Private(PetscInt64 seta[], PetscInt *lena, PetscInt64 setb[], PetscInt lenb) 16*c0888a1bSJDBetteridge Where 17*c0888a1bSJDBetteridge seta = (PetscInt64 *)inoutset; 18*c0888a1bSJDBetteridge setb = (PetscInt64 *)inset; 19*c0888a1bSJDBetteridge And the arguments are passed as: 20*c0888a1bSJDBetteridge &seta[1], (PetscInt *)&seta[0], &setb[1], (PetscInt)setb[0] 21*c0888a1bSJDBetteridge */ 22*c0888a1bSJDBetteridge 23*c0888a1bSJDBetteridge /* Populate a set with upto the first 49 unique Fibonnaci numbers */ 24*c0888a1bSJDBetteridge PetscErrorCode Fibonnaci(PetscInt64 **set, PetscInt n) 25*c0888a1bSJDBetteridge { 26*c0888a1bSJDBetteridge PetscInt ii; 27*c0888a1bSJDBetteridge PetscInt64 fib[] = {1, 2, 3, 5, 8, 13, 21, 34, 55, 89, 144, 233, 377, 610, 987, 1597, 2584, 28*c0888a1bSJDBetteridge 4181, 6765, 10946, 17711, 28657, 46368, 75025, 121393, 196418, 317811, 514229, 832040, 1346269, 2178309, 3524578, 5702887, 9227465, 29*c0888a1bSJDBetteridge 14930352, 24157817, 39088169, 63245986, 102334155, 165580141, 267914296, 433494437, 701408733, 1134903170, 1836311903, 2971215073, 4807526976, 7778742049, 12586269025}; 30*c0888a1bSJDBetteridge 31*c0888a1bSJDBetteridge PetscFunctionBeginUser; 32*c0888a1bSJDBetteridge PetscAssert((n < 50), PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONGSTATE, "n must be less than 50\n"); 33*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(n + 1, set)); 34*c0888a1bSJDBetteridge (*set)[0] = (PetscInt64)n; 35*c0888a1bSJDBetteridge for (ii = 0; ii < n; ii++) { (*set)[ii + 1] = fib[ii]; } 36*c0888a1bSJDBetteridge PetscFunctionReturn(0); 37*c0888a1bSJDBetteridge } 38*c0888a1bSJDBetteridge 39*c0888a1bSJDBetteridge /* Populate a set with Square numbers */ 40*c0888a1bSJDBetteridge PetscErrorCode Square(PetscInt64 **set, PetscInt n) 41*c0888a1bSJDBetteridge { 42*c0888a1bSJDBetteridge PetscInt64 ii; 43*c0888a1bSJDBetteridge 44*c0888a1bSJDBetteridge PetscFunctionBeginUser; 45*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(n + 1, set)); 46*c0888a1bSJDBetteridge (*set)[0] = (PetscInt64)n; 47*c0888a1bSJDBetteridge for (ii = 1; ii < n + 1; ii++) { (*set)[ii] = ii * ii; } 48*c0888a1bSJDBetteridge PetscFunctionReturn(0); 49*c0888a1bSJDBetteridge } 50*c0888a1bSJDBetteridge 51*c0888a1bSJDBetteridge /* Populate a set with Cube numbers */ 52*c0888a1bSJDBetteridge PetscErrorCode Cube(PetscInt64 **set, PetscInt n) 53*c0888a1bSJDBetteridge { 54*c0888a1bSJDBetteridge PetscInt64 ii; 55*c0888a1bSJDBetteridge 56*c0888a1bSJDBetteridge PetscFunctionBeginUser; 57*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(n + 1, set)); 58*c0888a1bSJDBetteridge (*set)[0] = (PetscInt64)n; 59*c0888a1bSJDBetteridge for (ii = 1; ii < n + 1; ii++) { (*set)[ii] = ii * ii * ii; } 60*c0888a1bSJDBetteridge PetscFunctionReturn(0); 61*c0888a1bSJDBetteridge } 62*c0888a1bSJDBetteridge 63*c0888a1bSJDBetteridge /* Populate a set with numbers to sixth power */ 64*c0888a1bSJDBetteridge PetscErrorCode Sixth(PetscInt64 **set, PetscInt n) 65*c0888a1bSJDBetteridge { 66*c0888a1bSJDBetteridge PetscInt64 ii; 67*c0888a1bSJDBetteridge 68*c0888a1bSJDBetteridge PetscFunctionBeginUser; 69*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(n + 1, set)); 70*c0888a1bSJDBetteridge (*set)[0] = (PetscInt64)n; 71*c0888a1bSJDBetteridge for (ii = 1; ii < n + 1; ii++) { (*set)[ii] = ii * ii * ii * ii * ii * ii; } 72*c0888a1bSJDBetteridge PetscFunctionReturn(0); 73*c0888a1bSJDBetteridge } 74*c0888a1bSJDBetteridge 75*c0888a1bSJDBetteridge /* Print out the contents of a set */ 76*c0888a1bSJDBetteridge PetscErrorCode PrintSet(PetscInt64 *set) 77*c0888a1bSJDBetteridge { 78*c0888a1bSJDBetteridge char text[64]; 79*c0888a1bSJDBetteridge PetscInt ii; 80*c0888a1bSJDBetteridge 81*c0888a1bSJDBetteridge PetscFunctionBeginUser; 82*c0888a1bSJDBetteridge PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[")); 83*c0888a1bSJDBetteridge for (ii = 1; ii <= (PetscInt)set[0]; ii++) { 84*c0888a1bSJDBetteridge PetscCall(PetscFormatConvert(" %" PetscInt64_FMT ",", text)); 85*c0888a1bSJDBetteridge PetscCall(PetscPrintf(PETSC_COMM_WORLD, text, set[ii])); 86*c0888a1bSJDBetteridge } 87*c0888a1bSJDBetteridge PetscCall(PetscPrintf(PETSC_COMM_WORLD, "]\n")); 88*c0888a1bSJDBetteridge PetscFunctionReturn(0); 89*c0888a1bSJDBetteridge } 90*c0888a1bSJDBetteridge 91*c0888a1bSJDBetteridge /* Check set equality */ 92*c0888a1bSJDBetteridge PetscErrorCode AssertSetsEqual(PetscInt64 *set, PetscInt64 *true_set) 93*c0888a1bSJDBetteridge { 94*c0888a1bSJDBetteridge PetscInt ii; 95*c0888a1bSJDBetteridge 96*c0888a1bSJDBetteridge PetscFunctionBeginUser; 97*c0888a1bSJDBetteridge PetscAssert((set[0] == true_set[0]), PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets of different sizes"); 98*c0888a1bSJDBetteridge for (ii = 1; ii < set[0] + 1; ii++) { PetscAssert((set[ii] == true_set[ii]), PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Sets are different"); } 99*c0888a1bSJDBetteridge PetscFunctionReturn(0); 100*c0888a1bSJDBetteridge } 101*c0888a1bSJDBetteridge 102*c0888a1bSJDBetteridge /* Tests functionality when two enpty sets are passed */ 103*c0888a1bSJDBetteridge PetscErrorCode test_empty_empty() 104*c0888a1bSJDBetteridge { 105*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 106*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 107*c0888a1bSJDBetteridge PetscMPIInt length = 1; 108*c0888a1bSJDBetteridge 109*c0888a1bSJDBetteridge PetscFunctionBeginUser; 110*c0888a1bSJDBetteridge 111*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(1, &set_a)); 112*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(1, &set_b)); 113*c0888a1bSJDBetteridge 114*c0888a1bSJDBetteridge set_a[0] = 0; 115*c0888a1bSJDBetteridge 116*c0888a1bSJDBetteridge set_b[0] = 0; 117*c0888a1bSJDBetteridge 118*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 119*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 120*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 121*c0888a1bSJDBetteridge 122*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 123*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 124*c0888a1bSJDBetteridge 125*c0888a1bSJDBetteridge PetscFunctionReturn(0); 126*c0888a1bSJDBetteridge } 127*c0888a1bSJDBetteridge 128*c0888a1bSJDBetteridge /* Tests functionality when seta is empty */ 129*c0888a1bSJDBetteridge PetscErrorCode test_a_empty() 130*c0888a1bSJDBetteridge { 131*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 132*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 133*c0888a1bSJDBetteridge PetscMPIInt length = 1; 134*c0888a1bSJDBetteridge 135*c0888a1bSJDBetteridge PetscFunctionBeginUser; 136*c0888a1bSJDBetteridge 137*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(1, &set_a)); 138*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(2, &set_b)); 139*c0888a1bSJDBetteridge 140*c0888a1bSJDBetteridge set_a[0] = 0; 141*c0888a1bSJDBetteridge 142*c0888a1bSJDBetteridge set_b[0] = 1; 143*c0888a1bSJDBetteridge set_b[1] = 1; 144*c0888a1bSJDBetteridge 145*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 146*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 147*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 148*c0888a1bSJDBetteridge 149*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 150*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 151*c0888a1bSJDBetteridge 152*c0888a1bSJDBetteridge PetscFunctionReturn(0); 153*c0888a1bSJDBetteridge } 154*c0888a1bSJDBetteridge 155*c0888a1bSJDBetteridge /* Tests functionality when setb is empty */ 156*c0888a1bSJDBetteridge PetscErrorCode test_b_empty() 157*c0888a1bSJDBetteridge { 158*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 159*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 160*c0888a1bSJDBetteridge PetscMPIInt length = 1; 161*c0888a1bSJDBetteridge 162*c0888a1bSJDBetteridge PetscFunctionBeginUser; 163*c0888a1bSJDBetteridge 164*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(2, &set_a)); 165*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(1, &set_b)); 166*c0888a1bSJDBetteridge 167*c0888a1bSJDBetteridge set_a[0] = 1; 168*c0888a1bSJDBetteridge set_a[1] = 1; 169*c0888a1bSJDBetteridge 170*c0888a1bSJDBetteridge set_b[0] = 0; 171*c0888a1bSJDBetteridge 172*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 173*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 174*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 175*c0888a1bSJDBetteridge 176*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 177*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 178*c0888a1bSJDBetteridge 179*c0888a1bSJDBetteridge PetscFunctionReturn(0); 180*c0888a1bSJDBetteridge } 181*c0888a1bSJDBetteridge 182*c0888a1bSJDBetteridge /* Tests functionality when both sets are identical */ 183*c0888a1bSJDBetteridge PetscErrorCode test_identical() 184*c0888a1bSJDBetteridge { 185*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 186*c0888a1bSJDBetteridge PetscInt64 truth[] = {3, 1, 4, 9}; 187*c0888a1bSJDBetteridge PetscMPIInt length = 4; 188*c0888a1bSJDBetteridge 189*c0888a1bSJDBetteridge PetscFunctionBeginUser; 190*c0888a1bSJDBetteridge 191*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(4, &set_a)); 192*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(4, &set_b)); 193*c0888a1bSJDBetteridge 194*c0888a1bSJDBetteridge set_a[0] = 3; 195*c0888a1bSJDBetteridge set_a[1] = 1; 196*c0888a1bSJDBetteridge set_a[2] = 4; 197*c0888a1bSJDBetteridge set_a[3] = 9; 198*c0888a1bSJDBetteridge 199*c0888a1bSJDBetteridge set_b[0] = 3; 200*c0888a1bSJDBetteridge set_b[1] = 1; 201*c0888a1bSJDBetteridge set_b[2] = 4; 202*c0888a1bSJDBetteridge set_b[3] = 9; 203*c0888a1bSJDBetteridge 204*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 205*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 206*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 207*c0888a1bSJDBetteridge 208*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 209*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 210*c0888a1bSJDBetteridge 211*c0888a1bSJDBetteridge PetscFunctionReturn(0); 212*c0888a1bSJDBetteridge } 213*c0888a1bSJDBetteridge 214*c0888a1bSJDBetteridge /* Tests functionality when sets have no elements in common */ 215*c0888a1bSJDBetteridge PetscErrorCode test_disjoint() 216*c0888a1bSJDBetteridge { 217*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 218*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 219*c0888a1bSJDBetteridge PetscMPIInt length = 1; 220*c0888a1bSJDBetteridge 221*c0888a1bSJDBetteridge PetscFunctionBeginUser; 222*c0888a1bSJDBetteridge 223*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(4, &set_a)); 224*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(4, &set_b)); 225*c0888a1bSJDBetteridge 226*c0888a1bSJDBetteridge set_a[0] = 3; 227*c0888a1bSJDBetteridge set_a[1] = 1; 228*c0888a1bSJDBetteridge set_a[2] = 4; 229*c0888a1bSJDBetteridge set_a[3] = 9; 230*c0888a1bSJDBetteridge 231*c0888a1bSJDBetteridge set_b[0] = 3; 232*c0888a1bSJDBetteridge set_b[1] = 2; 233*c0888a1bSJDBetteridge set_b[2] = 6; 234*c0888a1bSJDBetteridge set_b[3] = 8; 235*c0888a1bSJDBetteridge 236*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 237*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 238*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 239*c0888a1bSJDBetteridge 240*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 241*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 242*c0888a1bSJDBetteridge 243*c0888a1bSJDBetteridge PetscFunctionReturn(0); 244*c0888a1bSJDBetteridge } 245*c0888a1bSJDBetteridge 246*c0888a1bSJDBetteridge /* Tests functionality when sets only have one element in common */ 247*c0888a1bSJDBetteridge PetscErrorCode test_single_common() 248*c0888a1bSJDBetteridge { 249*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 250*c0888a1bSJDBetteridge PetscInt64 truth[] = {1, 4}; 251*c0888a1bSJDBetteridge PetscMPIInt length = 1; 252*c0888a1bSJDBetteridge 253*c0888a1bSJDBetteridge PetscFunctionBeginUser; 254*c0888a1bSJDBetteridge 255*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(4, &set_a)); 256*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(5, &set_b)); 257*c0888a1bSJDBetteridge 258*c0888a1bSJDBetteridge set_a[0] = 3; 259*c0888a1bSJDBetteridge set_a[1] = 1; 260*c0888a1bSJDBetteridge set_a[2] = 4; 261*c0888a1bSJDBetteridge set_a[3] = 9; 262*c0888a1bSJDBetteridge 263*c0888a1bSJDBetteridge set_b[0] = 3; 264*c0888a1bSJDBetteridge set_b[1] = 2; 265*c0888a1bSJDBetteridge set_b[2] = 4; 266*c0888a1bSJDBetteridge set_b[3] = 6; 267*c0888a1bSJDBetteridge set_b[4] = 8; 268*c0888a1bSJDBetteridge 269*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 270*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 271*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 272*c0888a1bSJDBetteridge 273*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 274*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 275*c0888a1bSJDBetteridge 276*c0888a1bSJDBetteridge PetscFunctionReturn(0); 277*c0888a1bSJDBetteridge } 278*c0888a1bSJDBetteridge 279*c0888a1bSJDBetteridge /* Specific test case flagged by PETSc issue #1247 */ 280*c0888a1bSJDBetteridge PetscErrorCode test_issue_1247() 281*c0888a1bSJDBetteridge { 282*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 283*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 284*c0888a1bSJDBetteridge PetscMPIInt length = 1; 285*c0888a1bSJDBetteridge 286*c0888a1bSJDBetteridge PetscFunctionBeginUser; 287*c0888a1bSJDBetteridge 288*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(3, &set_a)); 289*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(2, &set_b)); 290*c0888a1bSJDBetteridge 291*c0888a1bSJDBetteridge set_a[0] = 2; 292*c0888a1bSJDBetteridge set_a[1] = 2; 293*c0888a1bSJDBetteridge set_a[2] = 3; 294*c0888a1bSJDBetteridge 295*c0888a1bSJDBetteridge set_b[0] = 1; 296*c0888a1bSJDBetteridge set_b[1] = 1; 297*c0888a1bSJDBetteridge 298*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 299*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 300*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 301*c0888a1bSJDBetteridge 302*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 303*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 304*c0888a1bSJDBetteridge 305*c0888a1bSJDBetteridge PetscFunctionReturn(0); 306*c0888a1bSJDBetteridge } 307*c0888a1bSJDBetteridge 308*c0888a1bSJDBetteridge /* Tests functionality when seta is empty and setb is large */ 309*c0888a1bSJDBetteridge PetscErrorCode test_empty_big() 310*c0888a1bSJDBetteridge { 311*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 312*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 313*c0888a1bSJDBetteridge PetscMPIInt length = 1; 314*c0888a1bSJDBetteridge 315*c0888a1bSJDBetteridge PetscFunctionBeginUser; 316*c0888a1bSJDBetteridge 317*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(1, &set_a)); 318*c0888a1bSJDBetteridge PetscCall(Square(&set_b, 999)); 319*c0888a1bSJDBetteridge 320*c0888a1bSJDBetteridge set_a[0] = 0; 321*c0888a1bSJDBetteridge 322*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 323*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 324*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 325*c0888a1bSJDBetteridge 326*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 327*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 328*c0888a1bSJDBetteridge 329*c0888a1bSJDBetteridge PetscFunctionReturn(0); 330*c0888a1bSJDBetteridge } 331*c0888a1bSJDBetteridge 332*c0888a1bSJDBetteridge /* Tests functionality when seta is small and setb is large */ 333*c0888a1bSJDBetteridge PetscErrorCode test_small_big() 334*c0888a1bSJDBetteridge { 335*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 336*c0888a1bSJDBetteridge PetscInt64 truth[] = {3, 1, 4, 9}; 337*c0888a1bSJDBetteridge PetscMPIInt length = 1; 338*c0888a1bSJDBetteridge 339*c0888a1bSJDBetteridge PetscFunctionBeginUser; 340*c0888a1bSJDBetteridge 341*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(5, &set_a)); 342*c0888a1bSJDBetteridge PetscCall(Square(&set_b, 999)); 343*c0888a1bSJDBetteridge 344*c0888a1bSJDBetteridge set_a[0] = 4; 345*c0888a1bSJDBetteridge set_a[1] = 1; 346*c0888a1bSJDBetteridge set_a[2] = 4; 347*c0888a1bSJDBetteridge set_a[3] = 8; 348*c0888a1bSJDBetteridge set_a[4] = 9; 349*c0888a1bSJDBetteridge 350*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 351*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 352*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 353*c0888a1bSJDBetteridge 354*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 355*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 356*c0888a1bSJDBetteridge 357*c0888a1bSJDBetteridge PetscFunctionReturn(0); 358*c0888a1bSJDBetteridge } 359*c0888a1bSJDBetteridge 360*c0888a1bSJDBetteridge /* Tests functionality when seta is medium sized and setb is large */ 361*c0888a1bSJDBetteridge PetscErrorCode test_moderate_big() 362*c0888a1bSJDBetteridge { 363*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 364*c0888a1bSJDBetteridge PetscInt64 truth[] = {2, 1, 144}; 365*c0888a1bSJDBetteridge PetscMPIInt length = 1; 366*c0888a1bSJDBetteridge 367*c0888a1bSJDBetteridge PetscFunctionBeginUser; 368*c0888a1bSJDBetteridge 369*c0888a1bSJDBetteridge PetscCall(Fibonnaci(&set_a, 49)); 370*c0888a1bSJDBetteridge PetscCall(Square(&set_b, 999)); 371*c0888a1bSJDBetteridge 372*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 373*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 374*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 375*c0888a1bSJDBetteridge 376*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 377*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 378*c0888a1bSJDBetteridge 379*c0888a1bSJDBetteridge PetscFunctionReturn(0); 380*c0888a1bSJDBetteridge } 381*c0888a1bSJDBetteridge 382*c0888a1bSJDBetteridge /* Tests functionality when seta and setb are large */ 383*c0888a1bSJDBetteridge PetscErrorCode test_big_big() 384*c0888a1bSJDBetteridge { 385*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 386*c0888a1bSJDBetteridge PetscInt64 *truth; 387*c0888a1bSJDBetteridge PetscMPIInt length = 1; 388*c0888a1bSJDBetteridge 389*c0888a1bSJDBetteridge PetscFunctionBeginUser; 390*c0888a1bSJDBetteridge 391*c0888a1bSJDBetteridge PetscCall(Cube(&set_a, 999)); 392*c0888a1bSJDBetteridge PetscCall(Square(&set_b, 999)); 393*c0888a1bSJDBetteridge 394*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 395*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 396*c0888a1bSJDBetteridge 397*c0888a1bSJDBetteridge PetscCall(Sixth(&truth, 9)); 398*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 399*c0888a1bSJDBetteridge 400*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 401*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 402*c0888a1bSJDBetteridge PetscCall(PetscFree(truth)); 403*c0888a1bSJDBetteridge 404*c0888a1bSJDBetteridge PetscFunctionReturn(0); 405*c0888a1bSJDBetteridge } 406*c0888a1bSJDBetteridge 407*c0888a1bSJDBetteridge /* Tests functionality when setb is empty and setb is large */ 408*c0888a1bSJDBetteridge PetscErrorCode test_big_empty() 409*c0888a1bSJDBetteridge { 410*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 411*c0888a1bSJDBetteridge PetscInt64 truth[] = {0}; 412*c0888a1bSJDBetteridge PetscMPIInt length = 1; 413*c0888a1bSJDBetteridge 414*c0888a1bSJDBetteridge PetscFunctionBeginUser; 415*c0888a1bSJDBetteridge 416*c0888a1bSJDBetteridge PetscCall(Cube(&set_a, 999)); 417*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(1, &set_b)); 418*c0888a1bSJDBetteridge 419*c0888a1bSJDBetteridge set_b[0] = 0; 420*c0888a1bSJDBetteridge 421*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 422*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 423*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 424*c0888a1bSJDBetteridge 425*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 426*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 427*c0888a1bSJDBetteridge 428*c0888a1bSJDBetteridge PetscFunctionReturn(0); 429*c0888a1bSJDBetteridge } 430*c0888a1bSJDBetteridge 431*c0888a1bSJDBetteridge /* Tests functionality when setb is small and setb is large */ 432*c0888a1bSJDBetteridge PetscErrorCode test_big_small() 433*c0888a1bSJDBetteridge { 434*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 435*c0888a1bSJDBetteridge PetscInt64 truth[] = {2, 1, 8}; 436*c0888a1bSJDBetteridge PetscMPIInt length = 1; 437*c0888a1bSJDBetteridge 438*c0888a1bSJDBetteridge PetscFunctionBeginUser; 439*c0888a1bSJDBetteridge 440*c0888a1bSJDBetteridge PetscCall(Cube(&set_a, 999)); 441*c0888a1bSJDBetteridge PetscCall(PetscMalloc1(5, &set_b)); 442*c0888a1bSJDBetteridge 443*c0888a1bSJDBetteridge set_b[0] = 4; 444*c0888a1bSJDBetteridge set_b[1] = 1; 445*c0888a1bSJDBetteridge set_b[2] = 4; 446*c0888a1bSJDBetteridge set_b[3] = 8; 447*c0888a1bSJDBetteridge set_b[4] = 9; 448*c0888a1bSJDBetteridge 449*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 450*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 451*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 452*c0888a1bSJDBetteridge 453*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 454*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 455*c0888a1bSJDBetteridge 456*c0888a1bSJDBetteridge PetscFunctionReturn(0); 457*c0888a1bSJDBetteridge } 458*c0888a1bSJDBetteridge 459*c0888a1bSJDBetteridge /* Tests functionality when setb is medium sized and setb is large */ 460*c0888a1bSJDBetteridge PetscErrorCode test_big_moderate() 461*c0888a1bSJDBetteridge { 462*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 463*c0888a1bSJDBetteridge PetscInt64 truth[] = {2, 1, 8}; 464*c0888a1bSJDBetteridge PetscMPIInt length = 1; 465*c0888a1bSJDBetteridge 466*c0888a1bSJDBetteridge PetscFunctionBeginUser; 467*c0888a1bSJDBetteridge 468*c0888a1bSJDBetteridge PetscCall(Cube(&set_a, 999)); 469*c0888a1bSJDBetteridge PetscCall(Fibonnaci(&set_b, 49)); 470*c0888a1bSJDBetteridge 471*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 472*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 473*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 474*c0888a1bSJDBetteridge 475*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 476*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 477*c0888a1bSJDBetteridge 478*c0888a1bSJDBetteridge PetscFunctionReturn(0); 479*c0888a1bSJDBetteridge } 480*c0888a1bSJDBetteridge 481*c0888a1bSJDBetteridge /* Tests functionality when seta and setb are large, in the opposite 482*c0888a1bSJDBetteridge order to test_big_big() */ 483*c0888a1bSJDBetteridge PetscErrorCode test_big_big_reversed() 484*c0888a1bSJDBetteridge { 485*c0888a1bSJDBetteridge PetscInt64 *set_a, *set_b; 486*c0888a1bSJDBetteridge PetscInt64 *truth; 487*c0888a1bSJDBetteridge PetscMPIInt length = 1; 488*c0888a1bSJDBetteridge 489*c0888a1bSJDBetteridge PetscFunctionBeginUser; 490*c0888a1bSJDBetteridge 491*c0888a1bSJDBetteridge PetscCall(Cube(&set_a, 999)); 492*c0888a1bSJDBetteridge PetscCall(Square(&set_b, 999)); 493*c0888a1bSJDBetteridge 494*c0888a1bSJDBetteridge PetscGarbageKeySortedIntersect((void *)set_b, (void *)set_a, &length, NULL); 495*c0888a1bSJDBetteridge PetscCall(PrintSet(set_a)); 496*c0888a1bSJDBetteridge 497*c0888a1bSJDBetteridge PetscCall(Sixth(&truth, 9)); 498*c0888a1bSJDBetteridge PetscCall(AssertSetsEqual(set_a, truth)); 499*c0888a1bSJDBetteridge 500*c0888a1bSJDBetteridge PetscCall(PetscFree(set_a)); 501*c0888a1bSJDBetteridge PetscCall(PetscFree(set_b)); 502*c0888a1bSJDBetteridge PetscCall(PetscFree(truth)); 503*c0888a1bSJDBetteridge 504*c0888a1bSJDBetteridge PetscFunctionReturn(0); 505*c0888a1bSJDBetteridge } 506*c0888a1bSJDBetteridge 507*c0888a1bSJDBetteridge /* Main executes the individual tests in a predefined order */ 508*c0888a1bSJDBetteridge int main(int argc, char **argv) 509*c0888a1bSJDBetteridge { 510*c0888a1bSJDBetteridge PetscFunctionBeginUser; 511*c0888a1bSJDBetteridge PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 512*c0888a1bSJDBetteridge 513*c0888a1bSJDBetteridge /* Small tests */ 514*c0888a1bSJDBetteridge /* Test different edge cases with small sets */ 515*c0888a1bSJDBetteridge PetscCall(test_empty_empty()); 516*c0888a1bSJDBetteridge PetscCall(test_a_empty()); 517*c0888a1bSJDBetteridge PetscCall(test_b_empty()); 518*c0888a1bSJDBetteridge PetscCall(test_identical()); 519*c0888a1bSJDBetteridge PetscCall(test_disjoint()); 520*c0888a1bSJDBetteridge PetscCall(test_single_common()); 521*c0888a1bSJDBetteridge PetscCall(test_issue_1247()); 522*c0888a1bSJDBetteridge 523*c0888a1bSJDBetteridge /* Big tests */ 524*c0888a1bSJDBetteridge /* Test different edge cases with big sets */ 525*c0888a1bSJDBetteridge PetscCall(test_empty_big()); 526*c0888a1bSJDBetteridge PetscCall(test_small_big()); 527*c0888a1bSJDBetteridge PetscCall(test_moderate_big()); 528*c0888a1bSJDBetteridge PetscCall(test_big_big()); 529*c0888a1bSJDBetteridge PetscCall(test_big_empty()); 530*c0888a1bSJDBetteridge PetscCall(test_big_small()); 531*c0888a1bSJDBetteridge PetscCall(test_big_moderate()); 532*c0888a1bSJDBetteridge PetscCall(test_big_big_reversed()); 533*c0888a1bSJDBetteridge 534*c0888a1bSJDBetteridge PetscCall(PetscPrintf(PETSC_COMM_WORLD, "ALL PASSED\n")); 535*c0888a1bSJDBetteridge PetscCall(PetscFinalize()); 536*c0888a1bSJDBetteridge return 0; 537*c0888a1bSJDBetteridge } 538*c0888a1bSJDBetteridge 539*c0888a1bSJDBetteridge /*TEST 540*c0888a1bSJDBetteridge 541*c0888a1bSJDBetteridge test: 542*c0888a1bSJDBetteridge output_file: output/ex62.out 543*c0888a1bSJDBetteridge 544*c0888a1bSJDBetteridge TEST*/ 545