Actual source code: plexhdf5xdmf.c
petsc-3.10.3 2018-12-18
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/isimpl.h>
3: #include <petsc/private/vecimpl.h>
4: #include <petscviewerhdf5.h>
6: #if defined(PETSC_HAVE_HDF5)
7: PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM dm, PetscViewer viewer)
8: {
9: Vec coordinates;
10: IS cells;
11: PetscInt dim, spatialDim, N, numCells, numVertices, numCorners, bs;
12: //hid_t fileId;
13: PetscMPIInt rank;
14: MPI_Comm comm;
15: PetscErrorCode ierr;
18: PetscObjectGetComm((PetscObject)dm, &comm);
19: MPI_Comm_rank(comm, &rank);
21: /* Read topology */
22: PetscViewerHDF5ReadAttribute(viewer, "/viz/topology/cells", "cell_dim", PETSC_INT, (void *) &dim);
23: PetscViewerHDF5ReadAttribute(viewer, "/viz/topology/cells", "cell_corners", PETSC_INT, (void *) &numCorners);
24: PetscViewerHDF5PushGroup(viewer, "/viz/topology");
25: ISCreate(comm, &cells);
26: PetscObjectSetName((PetscObject) cells, "cells");
27: ISLoad(cells, viewer);
28: ISGetLocalSize(cells, &numCells);
29: ISGetBlockSize(cells, &bs);
30: if (bs != numCorners) SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "cell_corners and 2-dimension of cells not equal");
31: PetscViewerHDF5PopGroup(viewer);
32: numCells /= numCorners;
34: /* Read geometry */
35: PetscViewerHDF5PushGroup(viewer, "/geometry");
36: VecCreate(comm, &coordinates);
37: PetscObjectSetName((PetscObject) coordinates, "vertices");
38: PetscViewerHDF5ReadSizes(viewer, "vertices", &spatialDim, &N);
39: VecSetSizes(coordinates, PETSC_DECIDE, N);
40: VecSetBlockSize(coordinates, spatialDim);
41: VecLoad(coordinates, viewer);
42: VecGetLocalSize(coordinates, &numVertices);
43: VecGetBlockSize(coordinates, &spatialDim);
44: PetscViewerHDF5PopGroup(viewer);
45: numVertices /= spatialDim;
47: #if defined(PETSC_USE_DEBUG)
48: /* Check that maximum index referred in cells is in line with global number of vertices */
49: {
50: PetscInt max1, max2;
51: ISGetMinMax(cells, NULL, &max1);
52: VecGetSize(coordinates, &max2);
53: max2 /= spatialDim; max2--;
54: if (max1 > max2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "maximum index in cells = %d > %d = total number of vertices - 1", max1, max2);
55: }
56: #endif
58: {
59: const PetscScalar *coordinates_arr;
60: PetscReal *coordinates_arr_real;
61: const PetscInt *cells_arr;
62: int *cells_arr_int;
63: PetscSF sfVert=NULL;
64: #if defined(PETSC_USE_64BIT_INDICES) || defined(PETSC_USE_COMPLEX)
65: PetscInt i;
66: #endif
68: VecGetArrayRead(coordinates, &coordinates_arr);
69: ISGetIndices(cells, &cells_arr);
71: #if defined(PETSC_USE_64BIT_INDICES)
72: /* convert to 32-bit integers if PetscInt is 64-bit */
73: /*TODO More systematic would be to change all the function arguments to PetscInt */
74: PetscMalloc1(numCells*numCorners, &cells_arr_int);
75: for (i = 0; i < numCells*numCorners; ++i) {
76: PetscMPIIntCast(cells_arr[i], &cells_arr_int[i]);
77: }
78: #else
79: cells_arr_int = (int*)cells_arr;
80: #endif
82: #if defined(PETSC_USE_COMPLEX)
83: /* convert to real numbers if PetscScalar is complex */
84: /*TODO More systematic would be to change all the function arguments to PetscScalar */
85: PetscMalloc1(numVertices*spatialDim, &coordinates_arr_real);
86: for (i = 0; i < numVertices*spatialDim; ++i) {
87: coordinates_arr_real[i] = PetscRealPart(coordinates_arr[i]);
88: #if defined(PETSC_USE_DEBUG)
89: if (PetscImaginaryPart(coordinates_arr[i])) {
90: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vector of coordinates contains complex numbers but only real vectors are currently supported.");
91: }
92: #endif /* defined(PETSC_USE_DEBUG) */
93: }
94: #else
95: coordinates_arr_real = (PetscReal*)coordinates_arr;
96: #endif
98: DMSetDimension(dm, dim);
99: DMPlexBuildFromCellList_Parallel_Internal(dm, spatialDim, numCells, numVertices, numCorners, cells_arr_int, PETSC_TRUE, &sfVert);
100: DMPlexBuildCoordinates_Parallel_Internal( dm, spatialDim, numCells, numVertices, sfVert, coordinates_arr_real);
101: VecRestoreArrayRead(coordinates, &coordinates_arr);
102: ISRestoreIndices(cells, &cells_arr);
103: PetscSFDestroy(&sfVert);
104: #if defined(PETSC_USE_64BIT_INDICES)
105: PetscFree(cells_arr_int);
106: #endif
107: #if defined(PETSC_USE_COMPLEX)
108: PetscFree(coordinates_arr_real);
109: #endif
110: }
111: ISDestroy(&cells);
112: VecDestroy(&coordinates);
114: /* scale coordinates - unlike in DMPlexLoad_HDF5_Internal, this can only be done after DM is populated */
115: {
116: PetscReal lengthScale;
118: DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale);
119: DMGetCoordinates(dm, &coordinates);
120: VecScale(coordinates, 1.0/lengthScale);
121: }
123: /* Read Labels */
124: DMPlexLoadLabels_HDF5_Internal(dm, viewer);
125: return(0);
126: }
127: #endif