|
H5Part 1.6.6
|
Functions/Subroutines | |
| integer *8 function | h5pt_getdatasetname (filehandle, index, name) |
| See H5PartGetDatasetName. | |
| integer *8 function | h5pt_getndatasets (filehandle) |
| See H5PartGetNumDatasets. | |
| integer *8 function | h5pt_getnpoints (filehandle) |
| See H5PartGetNumParticles. | |
| integer *8 function | h5pt_getnsteps (filehandle) |
| See H5PartGetNumSteps. | |
| integer *8 function | h5pt_getview (filehandle, start, end) |
| See H5PartGetView. | |
| integer *8 function | h5pt_hasview (filehandle) |
| See H5PartResetView. | |
| integer *8 function | h5pt_resetview (filehandle) |
| See H5PartResetView. | |
| integer *8 function | h5pt_setnpoints (filehandle, npoints) |
| See H5PartSetNumParticles. | |
| integer *8 function | h5pt_setnpoints_strided (filehandle, npoints, stride) |
| See H5PartSetNumParticlesStrided. | |
| integer *8 function | h5pt_setstep (filehandle, step) |
| See H5PartSetStep. | |
| integer *8 function | h5pt_setview (filehandle, start, end) |
| See H5PartSetView. | |
| integer *8 function | h5pt_setview_indices (filehandle, indices, nelem) |
| See H5PartSetViewIndices. | |
| integer*8 function h5pt_getdatasetname | ( | integer*8, intent(in) | filehandle, |
| integer*8, intent(in) | index, | ||
| character(len=*), intent(out) | name ) |
See H5PartGetDatasetName.
| [in] | filehandle | the handle returned during file open |
| [in] | index | index of dataset to query (starting from 0) |
| [out] | name | buffer to read the dataset name into |
| integer*8 function h5pt_getndatasets | ( | integer*8, intent(in) | filehandle | ) |
See H5PartGetNumDatasets.
| [in] | filehandle | the handle returned during file open |
| integer*8 function h5pt_getnpoints | ( | integer*8, intent(in) | filehandle | ) |
| [in] | filehandle | the handle returned during file open |
| integer*8 function h5pt_getnsteps | ( | integer*8, intent(in) | filehandle | ) |
See H5PartGetNumSteps.
| [in] | filehandle | the handle returned during file open |
| integer*8 function h5pt_getview | ( | integer*8, intent(in) | filehandle, |
| integer*8, intent(out) | start, | ||
| integer*8, intent(out) | end ) |
See H5PartGetView.
| [in] | filehandle | the handle returned during file open |
| [out] | start | buffer to store the offset of the first particle in the view |
| [out] | end | buffer to store the offset of the last particle in the view (inclusive) |
| integer*8 function h5pt_hasview | ( | integer*8, intent(in) | filehandle | ) |
See H5PartResetView.
| [in] | filehandle | the handle returned during file open |
| integer*8 function h5pt_resetview | ( | integer*8, intent(in) | filehandle | ) |
See H5PartResetView.
| [in] | filehandle | the handle returned during file open |
| integer*8 function h5pt_setnpoints | ( | integer*8, intent(in) | filehandle, |
| integer*8, intent(in) | npoints ) |
| [in] | filehandle | the handle returned during file open |
| [in] | npoints | the number of particles on this processor |
| integer*8 function h5pt_setnpoints_strided | ( | integer*8, intent(in) | filehandle, |
| integer*8, intent(in) | npoints, | ||
| integer*8, intent(in) | stride ) |
See H5PartSetNumParticlesStrided.
| [in] | filehandle | the handle returned during file open |
| [in] | npoints | the number of particles on this processor |
| [in] | stride | the stride value (e.g. the number of fields in the particle data array) |
| integer*8 function h5pt_setstep | ( | integer*8, intent(in) | filehandle, |
| integer*8, intent(in) | step ) |
See H5PartSetStep.
| [in] | filehandle | the handle returned during file open |
| [in] | step | a timestep value >= 1 |
| integer*8 function h5pt_setview | ( | integer*8, intent(in) | filehandle, |
| integer*8, intent(in) | start, | ||
| integer*8, intent(in) | end ) |
See H5PartSetView.
| [in] | filehandle | the handle returned during file open |
| [in] | start | offset of the first particle in the view |
| [in] | end | offset of the last particle in the view (inclusive) |
| integer*8 function h5pt_setview_indices | ( | integer*8, intent(in) | filehandle, |
| integer*8, dimension(*), intent(in) | indices, | ||
| integer*8, intent(in) | nelem ) |
See H5PartSetViewIndices.
| [in] | filehandle | the handle returned during file open |
| [in] | indices | list of indicies to select in this view |
| [in] | nelem | number of particles in the list |