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