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