1c4762a1bSJed Brown static char help[] = "Test for mesh reordering\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown typedef struct { 6c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 7c4762a1bSJed Brown PetscReal refinementLimit; /* Maximum volume of a refined cell */ 8c4762a1bSJed Brown PetscInt numFields; /* The number of section fields */ 9c4762a1bSJed Brown PetscInt *numComponents; /* The number of field components */ 10c4762a1bSJed Brown PetscInt *numDof; /* The dof signature for the section */ 11c4762a1bSJed Brown PetscInt numGroups; /* If greater than 1, use grouping in test */ 12c4762a1bSJed Brown } AppCtx; 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscErrorCode ProcessOptions(AppCtx *options) 15c4762a1bSJed Brown { 16c4762a1bSJed Brown PetscInt len; 17c4762a1bSJed Brown PetscBool flg; 18c4762a1bSJed Brown PetscErrorCode ierr; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 21c4762a1bSJed Brown options->dim = 2; 22c4762a1bSJed Brown options->refinementLimit = 0.0; 23c4762a1bSJed Brown options->numFields = 1; 24c4762a1bSJed Brown options->numComponents = NULL; 25c4762a1bSJed Brown options->numDof = NULL; 26c4762a1bSJed Brown options->numGroups = 0; 27c4762a1bSJed Brown 28c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex10.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex10.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex10.c", options->numFields, &options->numFields, NULL,1);CHKERRQ(ierr); 32c4762a1bSJed Brown if (options->numFields) { 33c4762a1bSJed Brown len = options->numFields; 34c4762a1bSJed Brown ierr = PetscCalloc1(len, &options->numComponents);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = PetscOptionsIntArray("-num_components", "The number of components per field", "ex10.c", options->numComponents, &len, &flg);CHKERRQ(ierr); 36c4762a1bSJed Brown if (flg && (len != options->numFields)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %D should be %D", len, options->numFields); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown len = (options->dim+1) * PetscMax(1, options->numFields); 39c4762a1bSJed Brown ierr = PetscMalloc1(len, &options->numDof);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex10.c", options->numDof, &len, &flg);CHKERRQ(ierr); 41c4762a1bSJed Brown if (flg && (len != (options->dim+1) * PetscMax(1, options->numFields))) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of dof array is %D should be %D", len, (options->dim+1) * PetscMax(1, options->numFields)); 42c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-num_groups", "Group permutation by this many label values", "ex10.c", options->numGroups, &options->numGroups, NULL,0);CHKERRQ(ierr); 43c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown PetscErrorCode CleanupContext(AppCtx *user) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown PetscErrorCode ierr; 50c4762a1bSJed Brown 51c4762a1bSJed Brown PetscFunctionBegin; 52c4762a1bSJed Brown ierr = PetscFree(user->numComponents);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscFree(user->numDof);CHKERRQ(ierr); 54c4762a1bSJed Brown PetscFunctionReturn(0); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* This mesh comes from~\cite{saad2003}, Fig. 2.10, p. 70. */ 58c4762a1bSJed Brown PetscErrorCode CreateTestMesh(MPI_Comm comm, DM *dm, AppCtx *options) 59c4762a1bSJed Brown { 60a4a685f2SJacob Faibussowitsch const PetscInt cells[16*3] = {6, 7, 8, 7, 9, 10, 10, 11, 12, 11, 13, 14, 0, 6, 8, 6, 2, 7, 1, 8, 7, 1, 7, 10, 61c4762a1bSJed Brown 2, 9, 7, 10, 9, 4, 1, 10, 12, 10, 4, 11, 12, 11, 3, 3, 11, 14, 11, 4, 13, 14, 13, 5}; 62a4a685f2SJacob Faibussowitsch const PetscReal coords[15*2] = {0, -3, 0, -1, 2, -1, 0, 1, 2, 1, 63c4762a1bSJed Brown 0, 3, 1, -2, 1, -1, 0, -2, 2, 0, 64c4762a1bSJed Brown 1, 0, 1, 1, 0, 0, 1, 2, 0, 2}; 65c4762a1bSJed Brown PetscErrorCode ierr; 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscFunctionBegin; 68*0fdc7489SMatthew Knepley ierr = DMPlexCreateFromCellListPetsc(comm, 2, 16, 15, 3, PETSC_FALSE, cells, 2, coords, dm);CHKERRQ(ierr); 69c4762a1bSJed Brown PetscFunctionReturn(0); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscErrorCode TestReordering(DM dm, AppCtx *user) 73c4762a1bSJed Brown { 74c4762a1bSJed Brown DM pdm; 75c4762a1bSJed Brown IS perm; 76c4762a1bSJed Brown Mat A, pA; 77c4762a1bSJed Brown PetscInt bw, pbw; 78c4762a1bSJed Brown MatOrderingType order = MATORDERINGRCM; 79c4762a1bSJed Brown PetscErrorCode ierr; 80c4762a1bSJed Brown 81c4762a1bSJed Brown PetscFunctionBegin; 82c4762a1bSJed Brown ierr = DMPlexGetOrdering(dm, order, NULL, &perm);CHKERRQ(ierr); 83c4762a1bSJed Brown ierr = DMPlexPermute(dm, perm, &pdm);CHKERRQ(ierr); 84*0fdc7489SMatthew Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) pdm, "perm_");CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = DMSetFromOptions(pdm);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 87*0fdc7489SMatthew Knepley ierr = DMViewFromOptions(dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 88*0fdc7489SMatthew Knepley ierr = DMViewFromOptions(pdm, NULL, "-dm_view");CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &A);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = DMCreateMatrix(pdm, &pA);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatComputeBandwidth(A, 0.0, &bw);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatComputeBandwidth(pA, 0.0, &pbw);CHKERRQ(ierr); 93*0fdc7489SMatthew Knepley ierr = MatViewFromOptions(A, NULL, "-orig_mat_view");CHKERRQ(ierr); 94*0fdc7489SMatthew Knepley ierr = MatViewFromOptions(pA, NULL, "-perm_mat_view");CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = MatDestroy(&pA);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = DMDestroy(&pdm);CHKERRQ(ierr); 98c4762a1bSJed Brown if (pbw > bw) { 99c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s increased bandwidth from %D to %D\n", order, bw, pbw);CHKERRQ(ierr); 100c4762a1bSJed Brown } else { 101c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s reduced bandwidth from %D to %D\n", order, bw, pbw);CHKERRQ(ierr); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown PetscFunctionReturn(0); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown PetscErrorCode CreateGroupLabel(DM dm, PetscInt numGroups, DMLabel *label, AppCtx *options) 107c4762a1bSJed Brown { 108c4762a1bSJed Brown const PetscInt groupA[10] = {15, 3, 13, 12, 2, 10, 7, 6, 0, 4}; 109c4762a1bSJed Brown const PetscInt groupB[6] = {14, 11, 9, 1, 8, 5}; 110c4762a1bSJed Brown PetscInt c; 111c4762a1bSJed Brown PetscErrorCode ierr; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 114c4762a1bSJed Brown if (numGroups < 2) {*label = NULL; PetscFunctionReturn(0);} 115c4762a1bSJed Brown if (numGroups != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Test only coded for 2 groups, not %D", numGroups); 116c4762a1bSJed Brown ierr = DMLabelCreate(PETSC_COMM_SELF, "groups", label);CHKERRQ(ierr); 117c4762a1bSJed Brown for (c = 0; c < 10; ++c) {ierr = DMLabelSetValue(*label, groupA[c], 101);CHKERRQ(ierr);} 118c4762a1bSJed Brown for (c = 0; c < 6; ++c) {ierr = DMLabelSetValue(*label, groupB[c], 1001);CHKERRQ(ierr);} 119c4762a1bSJed Brown PetscFunctionReturn(0); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscErrorCode TestReorderingByGroup(DM dm, AppCtx *user) 123c4762a1bSJed Brown { 124c4762a1bSJed Brown DM pdm; 125c4762a1bSJed Brown DMLabel label; 126c4762a1bSJed Brown Mat A, pA; 127c4762a1bSJed Brown MatOrderingType order = MATORDERINGRCM; 128c4762a1bSJed Brown IS perm; 129c4762a1bSJed Brown PetscErrorCode ierr; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBegin; 132c4762a1bSJed Brown ierr = CreateGroupLabel(dm, user->numGroups, &label, user);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = DMPlexGetOrdering(dm, order, label, &perm);CHKERRQ(ierr); 134c4762a1bSJed Brown ierr = DMLabelDestroy(&label);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = DMPlexPermute(dm, perm, &pdm);CHKERRQ(ierr); 136*0fdc7489SMatthew Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) pdm, "perm_");CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = DMSetFromOptions(pdm);CHKERRQ(ierr); 138*0fdc7489SMatthew Knepley ierr = DMViewFromOptions(dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 139*0fdc7489SMatthew Knepley ierr = DMViewFromOptions(pdm, NULL, "-perm_dm_view");CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &A);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = DMCreateMatrix(pdm, &pA);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatViewFromOptions(A, NULL, "-orig_mat_view");CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = MatViewFromOptions(pA, NULL, "-perm_mat_view");CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = MatDestroy(&pA);CHKERRQ(ierr); 147c4762a1bSJed Brown ierr = DMDestroy(&pdm);CHKERRQ(ierr); 148c4762a1bSJed Brown PetscFunctionReturn(0); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown int main(int argc, char **argv) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown DM dm; 154c4762a1bSJed Brown PetscSection s; 155c4762a1bSJed Brown AppCtx user; 156c4762a1bSJed Brown PetscErrorCode ierr; 157c4762a1bSJed Brown 158c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 159c4762a1bSJed Brown ierr = ProcessOptions(&user);CHKERRQ(ierr); 160c4762a1bSJed Brown if (user.numGroups < 1) { 161*0fdc7489SMatthew Knepley ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, user.dim, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr); 162*0fdc7489SMatthew Knepley if (user.refinementLimit > 0.0) { 163*0fdc7489SMatthew Knepley DM rdm; 164*0fdc7489SMatthew Knepley const char *name; 165*0fdc7489SMatthew Knepley 166*0fdc7489SMatthew Knepley ierr = DMPlexSetRefinementUniform(dm, PETSC_FALSE);CHKERRQ(ierr); 167*0fdc7489SMatthew Knepley ierr = DMPlexSetRefinementLimit(dm, user.refinementLimit);CHKERRQ(ierr); 168*0fdc7489SMatthew Knepley ierr = DMRefine(dm, PETSC_COMM_WORLD, &rdm);CHKERRQ(ierr); 169*0fdc7489SMatthew Knepley ierr = PetscObjectGetName((PetscObject) dm, &name);CHKERRQ(ierr); 170*0fdc7489SMatthew Knepley ierr = PetscObjectSetName((PetscObject) rdm, name);CHKERRQ(ierr); 171*0fdc7489SMatthew Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 172*0fdc7489SMatthew Knepley dm = rdm; 173*0fdc7489SMatthew Knepley } 174c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = DMSetNumFields(dm, user.numFields);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = PetscSectionDestroy(&s);CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = TestReordering(dm, &user);CHKERRQ(ierr); 181c4762a1bSJed Brown } else { 182c4762a1bSJed Brown ierr = CreateTestMesh(PETSC_COMM_WORLD, &dm, &user);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = DMSetNumFields(dm, user.numFields);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = PetscSectionDestroy(&s);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = TestReorderingByGroup(dm, &user);CHKERRQ(ierr); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = CleanupContext(&user);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = PetscFinalize(); 194c4762a1bSJed Brown return ierr; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown /*TEST 198c4762a1bSJed Brown 199c4762a1bSJed Brown # Two cell tests 0-3 200c4762a1bSJed Brown test: 201c4762a1bSJed Brown suffix: 0 202*0fdc7489SMatthew Knepley requires: triangle 203*0fdc7489SMatthew Knepley args: -dim 2 -dm_plex_box_simplex 1 -num_dof 1,0,0 -mat_view 204c4762a1bSJed Brown test: 205c4762a1bSJed Brown suffix: 1 206*0fdc7489SMatthew Knepley args: -dim 2 -dm_plex_box_simplex 0 -num_dof 1,0,0 -mat_view 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown suffix: 2 209*0fdc7489SMatthew Knepley requires: ctetgen 210*0fdc7489SMatthew Knepley args: -dim 3 -dm_plex_box_simplex 1 -num_dof 1,0,0,0 -mat_view 211c4762a1bSJed Brown test: 212c4762a1bSJed Brown suffix: 3 213*0fdc7489SMatthew Knepley args: -dim 3 -dm_plex_box_simplex 0 -num_dof 1,0,0,0 -mat_view 214c4762a1bSJed Brown # Refined tests 4-7 215c4762a1bSJed Brown test: 216c4762a1bSJed Brown suffix: 4 217c4762a1bSJed Brown requires: triangle 218*0fdc7489SMatthew Knepley args: -dim 2 -dm_plex_box_simplex 1 -refinement_limit 0.00625 -num_dof 1,0,0 219c4762a1bSJed Brown test: 220c4762a1bSJed Brown suffix: 5 221*0fdc7489SMatthew Knepley args: -dim 2 -dm_plex_box_simplex 0 -dm_refine 1 -num_dof 1,0,0 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 6 224c4762a1bSJed Brown requires: ctetgen 225*0fdc7489SMatthew Knepley args: -dim 3 -dm_plex_box_simplex 1 -refinement_limit 0.00625 -num_dof 1,0,0,0 226c4762a1bSJed Brown test: 227c4762a1bSJed Brown suffix: 7 228*0fdc7489SMatthew Knepley args: -dim 3 -dm_plex_box_simplex 0 -dm_refine 1 -num_dof 1,0,0,0 229c4762a1bSJed Brown # Parallel tests 230c4762a1bSJed Brown # Grouping tests 231c4762a1bSJed Brown test: 232c4762a1bSJed Brown suffix: group_1 233c4762a1bSJed Brown args: -num_groups 1 -num_dof 1,0,0 -is_view -orig_mat_view -perm_mat_view 234c4762a1bSJed Brown test: 235c4762a1bSJed Brown suffix: group_2 236c4762a1bSJed Brown args: -num_groups 2 -num_dof 1,0,0 -is_view -perm_mat_view 237c4762a1bSJed Brown 238c4762a1bSJed Brown TEST*/ 239