1*c4762a1bSJed Brown static char help[] = "An example of writing a global Vec from a DMPlex with HDF5 format.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown #include <petscviewerhdf5.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc, char **argv) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown MPI_Comm comm; 9*c4762a1bSJed Brown DM dm; 10*c4762a1bSJed Brown Vec v, nv, rv, coord; 11*c4762a1bSJed Brown PetscBool test_read = PETSC_FALSE, verbose = PETSC_FALSE, flg; 12*c4762a1bSJed Brown PetscViewer hdf5Viewer; 13*c4762a1bSJed Brown PetscInt dim = 2; 14*c4762a1bSJed Brown PetscInt numFields = 1; 15*c4762a1bSJed Brown PetscInt numBC = 0; 16*c4762a1bSJed Brown PetscInt numComp[1] = {2}; 17*c4762a1bSJed Brown PetscInt numDof[3] = {2, 0, 0}; 18*c4762a1bSJed Brown PetscInt bcFields[1] = {0}; 19*c4762a1bSJed Brown IS bcPoints[1] = {NULL}; 20*c4762a1bSJed Brown PetscSection section; 21*c4762a1bSJed Brown PetscReal norm; 22*c4762a1bSJed Brown PetscErrorCode ierr; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, (char *) 0, help);if (ierr) return ierr; 25*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Test Options","none");CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = PetscOptionsBool("-test_read","Test reading from the HDF5 file","",PETSC_FALSE,&test_read,NULL);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = PetscOptionsBool("-verbose","print the Vecs","",PETSC_FALSE,&verbose,NULL);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim","the dimension of the problem","",2,&dim,NULL,1,3);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 35*c4762a1bSJed Brown numDof[0] = dim; 36*c4762a1bSJed Brown ierr = DMSetNumFields(dm, numFields);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcFields, bcPoints, NULL, NULL, §ion);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = DMSetUseNatural(dm, PETSC_TRUE);CHKERRQ(ierr); 41*c4762a1bSJed Brown { 42*c4762a1bSJed Brown PetscPartitioner part; 43*c4762a1bSJed Brown DM dmDist; 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr); 48*c4762a1bSJed Brown if (dmDist) { 49*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 50*c4762a1bSJed Brown dm = dmDist; 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &v);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) v, "V");CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = DMGetCoordinates(dm, &coord);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = VecCopy(coord, v);CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown if (verbose) { 60*c4762a1bSJed Brown PetscInt size, bs; 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = VecGetSize(v, &size);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecGetBlockSize(v, &bs);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = VecView(v, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &nv);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) nv, "NV");CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = DMPlexGlobalToNaturalBegin(dm, v, nv);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = DMPlexGlobalToNaturalEnd(dm, v, nv);CHKERRQ(ierr); 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown if (verbose) { 74*c4762a1bSJed Brown PetscInt size, bs; 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown ierr = VecGetSize(nv, &size);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = VecGetBlockSize(nv, &bs);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "==== V in natural ordering. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = VecView(nv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown ierr = VecViewFromOptions(v, NULL, "-global_vec_view");CHKERRQ(ierr); 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown if (test_read) { 85*c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &rv);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) rv, "V");CHKERRQ(ierr); 87*c4762a1bSJed Brown /* Test native read */ 88*c4762a1bSJed Brown ierr = PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = VecLoad(rv, hdf5Viewer);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = PetscViewerPopFormat(hdf5Viewer);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = PetscViewerDestroy(&hdf5Viewer);CHKERRQ(ierr); 93*c4762a1bSJed Brown if (verbose) { 94*c4762a1bSJed Brown PetscInt size, bs; 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown ierr = VecGetSize(rv, &size);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = VecGetBlockSize(rv, &bs);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown ierr = VecEqual(rv, v, &flg);CHKERRQ(ierr); 102*c4762a1bSJed Brown if (flg) { 103*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n");CHKERRQ(ierr); 104*c4762a1bSJed Brown } else { 105*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n");CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = VecAXPY(rv, -1.0, v);CHKERRQ(ierr); 107*c4762a1bSJed Brown ierr = VecNorm(rv, NORM_INFINITY, &norm);CHKERRQ(ierr); 108*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 110*c4762a1bSJed Brown } 111*c4762a1bSJed Brown /* Test raw read */ 112*c4762a1bSJed Brown ierr = PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = VecLoad(rv, hdf5Viewer);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscViewerDestroy(&hdf5Viewer);CHKERRQ(ierr); 115*c4762a1bSJed Brown if (verbose) { 116*c4762a1bSJed Brown PetscInt size, bs; 117*c4762a1bSJed Brown 118*c4762a1bSJed Brown ierr = VecGetSize(rv, &size);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = VecGetBlockSize(rv, &bs);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 122*c4762a1bSJed Brown } 123*c4762a1bSJed Brown ierr = VecEqual(rv, nv, &flg);CHKERRQ(ierr); 124*c4762a1bSJed Brown if (flg) { 125*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n");CHKERRQ(ierr); 126*c4762a1bSJed Brown } else { 127*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n");CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = VecAXPY(rv, -1.0, v);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = VecNorm(rv, NORM_INFINITY, &norm);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = VecView(rv, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 132*c4762a1bSJed Brown } 133*c4762a1bSJed Brown ierr = VecDestroy(&rv);CHKERRQ(ierr); 134*c4762a1bSJed Brown } 135*c4762a1bSJed Brown ierr = VecDestroy(&nv);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = VecDestroy(&v);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = PetscFinalize(); 139*c4762a1bSJed Brown return ierr; 140*c4762a1bSJed Brown } 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown /*TEST 143*c4762a1bSJed Brown build: 144*c4762a1bSJed Brown requires: triangle hdf5 145*c4762a1bSJed Brown test: 146*c4762a1bSJed Brown suffix: 0 147*c4762a1bSJed Brown requires: triangle hdf5 148*c4762a1bSJed Brown nsize: 2 149*c4762a1bSJed Brown args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view 150*c4762a1bSJed Brown test: 151*c4762a1bSJed Brown suffix: 1 152*c4762a1bSJed Brown requires: triangle hdf5 153*c4762a1bSJed Brown nsize: 2 154*c4762a1bSJed Brown args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown TEST*/ 157