1c4762a1bSJed Brown static char help[] = "Tests for DMLabel\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 40fdc7489SMatthew Knepley #include <petsc/private/dmimpl.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown static PetscErrorCode TestInsertion() 7c4762a1bSJed Brown { 8c4762a1bSJed Brown DMLabel label, label2; 9c4762a1bSJed Brown const PetscInt values[5] = {0, 3, 4, -1, 176}, N = 10000; 10c4762a1bSJed Brown PetscInt i, v; 11c4762a1bSJed Brown 12c4762a1bSJed Brown PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Test Label", &label)); 149566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(label, -100)); 15c4762a1bSJed Brown for (i = 0; i < N; ++i) { 169566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(label, i, values[i%5])); 17c4762a1bSJed Brown } 18c4762a1bSJed Brown /* Test get in hash mode */ 19c4762a1bSJed Brown for (i = 0; i < N; ++i) { 20c4762a1bSJed Brown PetscInt val; 21c4762a1bSJed Brown 229566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, i, &val)); 2363a3b9bcSJacob Faibussowitsch PetscCheck(val == values[i%5],PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Value %" PetscInt_FMT " for point %" PetscInt_FMT " should be %" PetscInt_FMT, val, i, values[i%5]); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown /* Test stratum */ 26c4762a1bSJed Brown for (v = 0; v < 5; ++v) { 27c4762a1bSJed Brown IS stratum; 28c4762a1bSJed Brown const PetscInt *points; 29c4762a1bSJed Brown PetscInt n; 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(label, values[v], &stratum)); 3263a3b9bcSJacob Faibussowitsch PetscCheck(stratum,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Stratum %" PetscInt_FMT " is empty!", v); 339566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratum, &points)); 349566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(stratum, &n)); 35c4762a1bSJed Brown for (i = 0; i < n; ++i) { 3663a3b9bcSJacob Faibussowitsch PetscCheck(points[i] == i*5+v,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %" PetscInt_FMT " should be %" PetscInt_FMT, points[i], i*5+v); 37c4762a1bSJed Brown } 389566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratum, &points)); 399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratum)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown /* Test get in array mode */ 42c4762a1bSJed Brown for (i = 0; i < N; ++i) { 43c4762a1bSJed Brown PetscInt val; 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, i, &val)); 4663a3b9bcSJacob Faibussowitsch PetscCheck(val == values[i%5],PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Value %" PetscInt_FMT " should be %" PetscInt_FMT, val, values[i%5]); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown /* Test Duplicate */ 499566063dSJacob Faibussowitsch PetscCall(DMLabelDuplicate(label, &label2)); 50c4762a1bSJed Brown for (i = 0; i < N; ++i) { 51c4762a1bSJed Brown PetscInt val; 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label2, i, &val)); 5463a3b9bcSJacob Faibussowitsch PetscCheck(val == values[i%5],PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Value %" PetscInt_FMT " should be %" PetscInt_FMT, val, values[i%5]); 55c4762a1bSJed Brown } 569566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label2)); 579566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&label)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown static PetscErrorCode TestEmptyStrata(MPI_Comm comm) 62c4762a1bSJed Brown { 63c4762a1bSJed Brown DM dm, dmDist; 64c4762a1bSJed Brown PetscPartitioner part; 65c4762a1bSJed Brown PetscInt c0[6] = {2,3,6,7,9,11}; 66c4762a1bSJed Brown PetscInt c1[6] = {4,5,7,8,10,12}; 67c4762a1bSJed Brown PetscInt c2[4] = {13,15,19,21}; 68c4762a1bSJed Brown PetscInt c3[4] = {14,16,20,22}; 69c4762a1bSJed Brown PetscInt c4[4] = {15,17,21,23}; 70c4762a1bSJed Brown PetscInt c5[4] = {16,18,22,24}; 71c4762a1bSJed Brown PetscInt c6[4] = {13,14,19,20}; 72c4762a1bSJed Brown PetscInt c7[4] = {15,16,21,22}; 73c4762a1bSJed Brown PetscInt c8[4] = {17,18,23,24}; 74c4762a1bSJed Brown PetscInt c9[4] = {13,14,15,16}; 75c4762a1bSJed Brown PetscInt c10[4] = {15,16,17,18}; 76c4762a1bSJed Brown PetscInt c11[4] = {19,20,21,22}; 77c4762a1bSJed Brown PetscInt c12[4] = {21,22,23,24}; 78c4762a1bSJed Brown PetscInt dim = 3; 79c4762a1bSJed Brown PetscMPIInt rank; 80c4762a1bSJed Brown 81c4762a1bSJed Brown PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 83c4762a1bSJed Brown /* A 3D box with two adjacent cells, sharing one face and four vertices */ 849566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm)); 859566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX)); 869566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm, dim)); 87dd400576SPatrick Sanan if (rank == 0) { 889566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(dm, 0, 25)); 899566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 0, 6)); 909566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 1, 6)); 919566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 2, 4)); 929566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 3, 4)); 939566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 4, 4)); 949566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 5, 4)); 959566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 6, 4)); 969566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 7, 4)); 979566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 8, 4)); 989566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 9, 4)); 999566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 10, 4)); 1009566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 11, 4)); 1019566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(dm, 12, 4)); 102c4762a1bSJed Brown } 1039566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 104dd400576SPatrick Sanan if (rank == 0) { 1059566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 0, c0)); 1069566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 1, c1)); 1079566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 2, c2)); 1089566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 3, c3)); 1099566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 4, c4)); 1109566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 5, c5)); 1119566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 6, c6)); 1129566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 7, c7)); 1139566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 8, c8)); 1149566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 9, c9)); 1159566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 10, c10)); 1169566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 11, c11)); 1179566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(dm, 12, c12)); 118c4762a1bSJed Brown } 1199566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(dm)); 120c4762a1bSJed Brown /* Create a user managed depth label, so that we can leave out edges */ 121c4762a1bSJed Brown { 122c4762a1bSJed Brown DMLabel label; 123c4762a1bSJed Brown PetscInt numValues, maxValues = 0, v; 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "depth")); 1269566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &label)); 127dd400576SPatrick Sanan if (rank == 0) { 128c4762a1bSJed Brown PetscInt i; 129c4762a1bSJed Brown 130c4762a1bSJed Brown for (i = 0; i < 25; ++i) { 1319566063dSJacob Faibussowitsch if (i < 2) PetscCall(DMLabelSetValue(label, i, 3)); 1329566063dSJacob Faibussowitsch else if (i < 13) PetscCall(DMLabelSetValue(label, i, 2)); 133c4762a1bSJed Brown else { 1349566063dSJacob Faibussowitsch if (i==13) PetscCall(DMLabelAddStratum(label, 1)); 1359566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(label, i, 0)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown } 138c4762a1bSJed Brown } 1399566063dSJacob Faibussowitsch PetscCall(DMLabelGetNumValues(label, &numValues)); 1409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&numValues, &maxValues, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) dm))); 1419566063dSJacob Faibussowitsch for (v = numValues; v < maxValues; ++v) PetscCall(DMLabelAddStratum(label,v)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown { 144c4762a1bSJed Brown DMLabel label; 1459566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &label)); 1469566063dSJacob Faibussowitsch PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_(comm))); 147c4762a1bSJed Brown } 1489566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm,&part)); 1499566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 1509566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 1, NULL, &dmDist)); 151c4762a1bSJed Brown if (dmDist) { 1529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 153c4762a1bSJed Brown dm = dmDist; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown { 156c4762a1bSJed Brown DMLabel label; 1579566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthLabel(dm, &label)); 1589566063dSJacob Faibussowitsch PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_(comm))); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown /* Create a cell vector */ 161c4762a1bSJed Brown { 162c4762a1bSJed Brown Vec v; 163c4762a1bSJed Brown PetscSection s; 164c4762a1bSJed Brown PetscInt numComp[] = {1}; 165c4762a1bSJed Brown PetscInt dof[] = {0,0,0,1}; 166c4762a1bSJed Brown PetscInt N; 167c4762a1bSJed Brown 1689566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm, 1)); 1699566063dSJacob Faibussowitsch PetscCall(DMPlexCreateSection(dm, NULL, numComp, dof, 0, NULL, NULL, NULL, NULL, &s)); 1709566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, s)); 1719566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s)); 1729566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &v)); 1739566063dSJacob Faibussowitsch PetscCall(VecGetSize(v, &N)); 174c4762a1bSJed Brown if (N != 2) { 1759566063dSJacob Faibussowitsch PetscCall(DMView(dm, PETSC_VIEWER_STDOUT_(comm))); 17663a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "FAIL: Vector size %" PetscInt_FMT " != 2", N); 177c4762a1bSJed Brown } 1789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 179c4762a1bSJed Brown } 1809566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 181c4762a1bSJed Brown PetscFunctionReturn(0); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown static PetscErrorCode TestDistribution(MPI_Comm comm) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown DM dm, dmDist; 187c4762a1bSJed Brown PetscPartitioner part; 188c4762a1bSJed Brown DMLabel label; 1890fdc7489SMatthew Knepley char filename[PETSC_MAX_PATH_LEN]; 190c4762a1bSJed Brown const char *name = "test label"; 191c4762a1bSJed Brown PetscInt overlap = 0, cStart, cEnd, c; 192c4762a1bSJed Brown PetscMPIInt rank; 193c4762a1bSJed Brown PetscBool flg; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBegin; 1969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-filename", filename, sizeof(filename), &flg)); 198c4762a1bSJed Brown if (!flg) PetscFunctionReturn(0); 1999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-overlap", &overlap, NULL)); 2009566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex11_plex", PETSC_TRUE, &dm)); 2019566063dSJacob Faibussowitsch PetscCall(DMSetBasicAdjacency(dm, PETSC_TRUE, PETSC_FALSE)); 2029566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, name)); 2039566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 2049566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 205c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 2069566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(label, c, c)); 207c4762a1bSJed Brown } 2089566063dSJacob Faibussowitsch PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD)); 2099566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm,&part)); 2109566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 2119566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, overlap, NULL, &dmDist)); 212c4762a1bSJed Brown if (dmDist) { 2139566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 214c4762a1bSJed Brown dm = dmDist; 215c4762a1bSJed Brown } 2169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm, "Mesh")); 2179566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 2189566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, name, &label)); 2199566063dSJacob Faibussowitsch PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD)); 2209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 221c4762a1bSJed Brown PetscFunctionReturn(0); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 2240fdc7489SMatthew Knepley static PetscErrorCode TestUniversalLabel(MPI_Comm comm) 2250fdc7489SMatthew Knepley { 2260fdc7489SMatthew Knepley DM dm1, dm2; 2270fdc7489SMatthew Knepley DMLabel bd1, bd2, ulabel; 2280fdc7489SMatthew Knepley DMUniversalLabel universal; 2290fdc7489SMatthew Knepley PetscInt pStart, pEnd, p; 23030602db0SMatthew G. Knepley PetscBool run = PETSC_FALSE, notFile; 2310fdc7489SMatthew Knepley 2320fdc7489SMatthew Knepley PetscFunctionBeginUser; 2339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-universal", &run, NULL)); 2340fdc7489SMatthew Knepley if (!run) PetscFunctionReturn(0); 23530602db0SMatthew G. Knepley 23630602db0SMatthew G. Knepley char filename[PETSC_MAX_PATH_LEN]; 23730602db0SMatthew G. Knepley PetscBool flg; 23830602db0SMatthew G. Knepley 2399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-filename", filename, sizeof(filename), &flg)); 2400fdc7489SMatthew Knepley if (flg) { 2419566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(comm, filename, "ex11_plex", PETSC_TRUE, &dm1)); 2420fdc7489SMatthew Knepley } else { 2439566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm1)); 2449566063dSJacob Faibussowitsch PetscCall(DMSetType(dm1, DMPLEX)); 2459566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm1)); 24630602db0SMatthew G. Knepley } 2479566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dm1, "marker", ¬File)); 24830602db0SMatthew G. Knepley if (notFile) { 2499566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm1, "Boundary Faces")); 2509566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm1, "Boundary Faces", &bd1)); 2519566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm1, 13, bd1)); 2529566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm1, "Boundary")); 2539566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm1, "Boundary", &bd2)); 2549566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm1, 121, bd2)); 2559566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(dm1, bd2)); 2560fdc7489SMatthew Knepley } 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm1, "First Mesh")); 2589566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm1, NULL, "-dm_view")); 2590fdc7489SMatthew Knepley 2609566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreate(dm1, &universal)); 2619566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelGetLabel(universal, &ulabel)); 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject) ulabel, NULL, "-universal_view")); 2630fdc7489SMatthew Knepley 26430602db0SMatthew G. Knepley if (!notFile) { 2650fdc7489SMatthew Knepley PetscInt Nl, l; 2660fdc7489SMatthew Knepley 2679566063dSJacob Faibussowitsch PetscCall(DMClone(dm1, &dm2)); 2689566063dSJacob Faibussowitsch PetscCall(DMGetNumLabels(dm2, &Nl)); 2690fdc7489SMatthew Knepley for (l = Nl-1; l >= 0; --l) { 2700fdc7489SMatthew Knepley PetscBool isdepth, iscelltype; 2710fdc7489SMatthew Knepley const char *name; 2720fdc7489SMatthew Knepley 2739566063dSJacob Faibussowitsch PetscCall(DMGetLabelName(dm2, l, &name)); 2749566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "depth", 6, &isdepth)); 2759566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "celltype", 9, &iscelltype)); 2769566063dSJacob Faibussowitsch if (!isdepth && !iscelltype) PetscCall(DMRemoveLabel(dm2, name, NULL)); 2770fdc7489SMatthew Knepley } 2780fdc7489SMatthew Knepley } else { 2799566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, &dm2)); 2809566063dSJacob Faibussowitsch PetscCall(DMSetType(dm2, DMPLEX)); 2819566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm2)); 2820fdc7489SMatthew Knepley } 2839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dm2, "Second Mesh")); 2849566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelCreateLabels(universal, PETSC_TRUE, dm2)); 2859566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm2, &pStart, &pEnd)); 2860fdc7489SMatthew Knepley for (p = pStart; p < pEnd; ++p) { 2870fdc7489SMatthew Knepley PetscInt val; 2880fdc7489SMatthew Knepley 2899566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(ulabel, p, &val)); 2900fdc7489SMatthew Knepley if (val < 0) continue; 2919566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelSetLabelValue(universal, dm2, PETSC_TRUE, p, val)); 2920fdc7489SMatthew Knepley } 2939566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm2, NULL, "-dm_view")); 2940fdc7489SMatthew Knepley 2959566063dSJacob Faibussowitsch PetscCall(DMUniversalLabelDestroy(&universal)); 2969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm1)); 2979566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm2)); 2980fdc7489SMatthew Knepley PetscFunctionReturn(0); 2990fdc7489SMatthew Knepley } 3000fdc7489SMatthew Knepley 301c4762a1bSJed Brown int main(int argc, char **argv) 302c4762a1bSJed Brown { 303c4762a1bSJed Brown 304*327415f7SBarry Smith PetscFunctionBeginUser; 3059566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3069566063dSJacob Faibussowitsch /*PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));*/ 3079566063dSJacob Faibussowitsch PetscCall(TestInsertion()); 3089566063dSJacob Faibussowitsch PetscCall(TestEmptyStrata(PETSC_COMM_WORLD)); 3099566063dSJacob Faibussowitsch PetscCall(TestDistribution(PETSC_COMM_WORLD)); 3109566063dSJacob Faibussowitsch PetscCall(TestUniversalLabel(PETSC_COMM_WORLD)); 3119566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 312b122ec5aSJacob Faibussowitsch return 0; 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 315c4762a1bSJed Brown /*TEST 316c4762a1bSJed Brown 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 0 3190fdc7489SMatthew Knepley requires: triangle 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 1 3220fdc7489SMatthew Knepley requires: triangle 323c4762a1bSJed Brown nsize: 2 324c4762a1bSJed Brown args: -petscpartitioner_type simple 325c4762a1bSJed Brown 326c4762a1bSJed Brown testset: 327c4762a1bSJed Brown suffix: gmsh 328c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/square.msh -petscpartitioner_type simple 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: 1 331c4762a1bSJed Brown nsize: 1 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 2 334c4762a1bSJed Brown nsize: 2 335c4762a1bSJed Brown 336c4762a1bSJed Brown testset: 337c4762a1bSJed Brown suffix: exodusii 338c4762a1bSJed Brown requires: exodusii 339c4762a1bSJed Brown args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/2Dgrd.exo -petscpartitioner_type simple 340c4762a1bSJed Brown test: 341c4762a1bSJed Brown suffix: 1 342c4762a1bSJed Brown nsize: 1 343c4762a1bSJed Brown test: 344c4762a1bSJed Brown suffix: 2 345c4762a1bSJed Brown nsize: 2 346c4762a1bSJed Brown 3470fdc7489SMatthew Knepley test: 3480fdc7489SMatthew Knepley suffix: univ 3490fdc7489SMatthew Knepley requires: triangle 3500fdc7489SMatthew Knepley args: -universal -dm_view -universal_view 3510fdc7489SMatthew Knepley 3520fdc7489SMatthew Knepley test: 3530fdc7489SMatthew Knepley # Note that the labels differ because we have multiply-marked some points during EGADS creation 3540fdc7489SMatthew Knepley suffix: univ_egads_sphere 3550fdc7489SMatthew Knepley requires: egads 35630602db0SMatthew G. Knepley args: -universal -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_view -universal_view 3570fdc7489SMatthew Knepley 3580fdc7489SMatthew Knepley test: 3590fdc7489SMatthew Knepley # Note that the labels differ because we have multiply-marked some points during EGADS creation 3600fdc7489SMatthew Knepley suffix: univ_egads_ball 3610fdc7489SMatthew Knepley requires: egads ctetgen 36230602db0SMatthew G. Knepley args: -universal -dm_plex_boundary_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_view -universal_view 3630fdc7489SMatthew Knepley 364c4762a1bSJed Brown TEST*/ 365