1c4762a1bSJed Brown static char help[] = "An example of writing a global Vec from a DMPlex with HDF5 format.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscviewerhdf5.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc, char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown MPI_Comm comm; 9c4762a1bSJed Brown DM dm; 10c4762a1bSJed Brown Vec v, nv, rv, coord; 11c4762a1bSJed Brown PetscBool test_read = PETSC_FALSE, verbose = PETSC_FALSE, flg; 12c4762a1bSJed Brown PetscViewer hdf5Viewer; 13c4762a1bSJed Brown PetscInt numFields = 1; 14c4762a1bSJed Brown PetscInt numBC = 0; 15c4762a1bSJed Brown PetscInt numComp[1] = {2}; 16c4762a1bSJed Brown PetscInt numDof[3] = {2, 0, 0}; 17c4762a1bSJed Brown PetscInt bcFields[1] = {0}; 18c4762a1bSJed Brown IS bcPoints[1] = {NULL}; 19c4762a1bSJed Brown PetscSection section; 20c4762a1bSJed Brown PetscReal norm; 2130602db0SMatthew G. Knepley PetscInt dim; 22c4762a1bSJed Brown PetscErrorCode ierr; 23c4762a1bSJed Brown 24*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, (char *) 0, help)); 25c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Test Options","none");CHKERRQ(ierr); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_read","Test reading from the HDF5 file","",PETSC_FALSE,&test_read,NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-verbose","print the Vecs","",PETSC_FALSE,&verbose,NULL)); 301e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 31c4762a1bSJed Brown 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, &dm)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(dm, DMPLEX)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 375f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 38c4762a1bSJed Brown numDof[0] = dim; 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetNumFields(dm, numFields)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateSection(dm, NULL, numComp, numDof, numBC, bcFields, bcPoints, NULL, NULL, §ion)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(dm, section)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(§ion)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUseNatural(dm, PETSC_TRUE)); 44c4762a1bSJed Brown { 45c4762a1bSJed Brown PetscPartitioner part; 46c4762a1bSJed Brown DM dmDist; 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm,&part)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm, 0, NULL, &dmDist)); 51c4762a1bSJed Brown if (dmDist) { 525f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 53c4762a1bSJed Brown dm = dmDist; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 575f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &v)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) v, "V")); 595f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(dm, &coord)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(coord, v)); 61c4762a1bSJed Brown 62c4762a1bSJed Brown if (verbose) { 63c4762a1bSJed Brown PetscInt size, bs; 64c4762a1bSJed Brown 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(v, &size)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(v, &bs)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "==== original V in global ordering. size==%d\tblock size=%d\n", size, bs)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v, PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 715f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &nv)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) nv, "NV")); 735f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGlobalToNaturalBegin(dm, v, nv)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGlobalToNaturalEnd(dm, v, nv)); 75c4762a1bSJed Brown 76c4762a1bSJed Brown if (verbose) { 77c4762a1bSJed Brown PetscInt size, bs; 78c4762a1bSJed Brown 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(nv, &size)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(nv, &bs)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "==== V in natural ordering. size==%d\tblock size=%d\n", size, bs)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(nv, PETSC_VIEWER_STDOUT_WORLD)); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(v, NULL, "-global_vec_view")); 86c4762a1bSJed Brown 87c4762a1bSJed Brown if (test_read) { 885f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &rv)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) rv, "V")); 90c4762a1bSJed Brown /* Test native read */ 915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(hdf5Viewer, PETSC_VIEWER_NATIVE)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(rv, hdf5Viewer)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(hdf5Viewer)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&hdf5Viewer)); 96c4762a1bSJed Brown if (verbose) { 97c4762a1bSJed Brown PetscInt size, bs; 98c4762a1bSJed Brown 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(rv, &size)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(rv, &bs)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(rv, PETSC_VIEWER_STDOUT_WORLD)); 103c4762a1bSJed Brown } 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(rv, v, &flg)); 105c4762a1bSJed Brown if (flg) { 1065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "V and RV are equal\n")); 107c4762a1bSJed Brown } else { 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "V and RV are not equal\n\n")); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(rv, -1.0, v)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(rv, NORM_INFINITY, &norm)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(rv, PETSC_VIEWER_STDOUT_WORLD)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown /* Test raw read */ 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerHDF5Open(comm, "V.h5", FILE_MODE_READ, &hdf5Viewer)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(rv, hdf5Viewer)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&hdf5Viewer)); 118c4762a1bSJed Brown if (verbose) { 119c4762a1bSJed Brown PetscInt size, bs; 120c4762a1bSJed Brown 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(rv, &size)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(rv, &bs)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "==== Vector from file. size==%d\tblock size=%d\n", size, bs)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(rv, PETSC_VIEWER_STDOUT_WORLD)); 125c4762a1bSJed Brown } 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(rv, nv, &flg)); 127c4762a1bSJed Brown if (flg) { 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are equal\n")); 129c4762a1bSJed Brown } else { 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "NV and RV are not equal\n\n")); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(rv, -1.0, v)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(rv, NORM_INFINITY, &norm)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "diff norm is = %g\n", (double) norm)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(rv, PETSC_VIEWER_STDOUT_WORLD)); 135c4762a1bSJed Brown } 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&rv)); 137c4762a1bSJed Brown } 1385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&nv)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 141*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 142*b122ec5aSJacob Faibussowitsch return 0; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /*TEST 146c4762a1bSJed Brown build: 147c4762a1bSJed Brown requires: triangle hdf5 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: 0 150c4762a1bSJed Brown requires: triangle hdf5 151c4762a1bSJed Brown nsize: 2 152c4762a1bSJed Brown args: -petscpartitioner_type simple -verbose -globaltonatural_sf_view 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown suffix: 1 155c4762a1bSJed Brown requires: triangle hdf5 156c4762a1bSJed Brown nsize: 2 157c4762a1bSJed Brown args: -petscpartitioner_type simple -verbose -global_vec_view hdf5:V.h5:native -test_read 158c4762a1bSJed Brown 159c4762a1bSJed Brown TEST*/ 160