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