read_setnparticles.c 2.02 KB
Newer Older
gsell's avatar
gsell committed
1
/*
gsell's avatar
gsell committed
2
  Copyright (c) 2006-2015, The Regents of the University of California,
gsell's avatar
gsell committed
3 4 5 6 7 8 9 10 11
  through Lawrence Berkeley National Laboratory (subject to receipt of any
  required approvals from the U.S. Dept. of Energy) and the Paul Scherrer
  Institut (Switzerland).  All rights reserved.

  License: see file COPYING in top level of source distribution.
*/

#include "H5hut.h"

gsell's avatar
gsell committed
12 13
#include <stdlib.h>

gsell's avatar
gsell committed
14 15 16 17 18
// name of input file
const char* fname = "example_setnparticles.h5";

// H5hut verbosity level
const h5_int64_t h5_verbosity = H5_VERBOSE_DEFAULT;
gsell's avatar
gsell committed
19 20 21 22 23

int
main (
        int argc, char* argv[]
        ){
gsell's avatar
gsell committed
24

gsell's avatar
gsell committed
25 26 27
        // initialize MPI & H5hut
        MPI_Init (&argc, &argv);
        MPI_Comm comm = MPI_COMM_WORLD;
gsell's avatar
gsell committed
28 29
        int comm_size = 1;
        MPI_Comm_size (comm, &comm_size);
gsell's avatar
gsell committed
30 31 32
        int comm_rank = 0;
        MPI_Comm_rank (comm, &comm_rank);
        H5AbortOnError ();
gsell's avatar
gsell committed
33
        H5SetVerbosityLevel (h5_verbosity);
gsell's avatar
gsell committed
34

gsell's avatar
gsell committed
35
        // open file and go to first step
gsell's avatar
gsell committed
36
        h5_file_t file = H5OpenFile (fname, H5_O_RDONLY, H5_PROP_DEFAULT);
gsell's avatar
gsell committed
37 38
        H5SetStep (file, 0);

39
        // compute number of particles this process has to read
gsell's avatar
gsell committed
40 41
        h5_ssize_t num_particles_total = H5PartGetNumParticles (file);
        h5_ssize_t num_particles = num_particles_total / comm_size;
gsell's avatar
gsell committed
42
        if (comm_rank+1 == comm_size)
gsell's avatar
gsell committed
43 44
                num_particles += num_particles_total % comm_size;

gsell's avatar
gsell committed
45
	printf ("[proc %d]: particles in view: %lld\n", comm_rank, (long long)num_particles);
gsell's avatar
gsell committed
46 47
	printf ("[proc %d]: total number of particles: %lld\n",
		comm_rank, (long long unsigned)num_particles_total);
gsell's avatar
gsell committed
48

49
	// set number of particles
gsell's avatar
gsell committed
50
        H5PartSetNumParticles (file, num_particles);
gsell's avatar
gsell committed
51

gsell's avatar
gsell committed
52 53
        // read and print data
        h5_int32_t* data = calloc (num_particles, sizeof (*data));
gsell's avatar
gsell committed
54
        H5PartReadDataInt32 (file, "data", data);
gsell's avatar
gsell committed
55 56 57 58
        for (int i = 0; i < num_particles; i++) {
                printf ("[proc %d]: local index = %d, value = %d\n",
                        comm_rank, i, data[i]);
        }
gsell's avatar
gsell committed
59 60

        // cleanup
gsell's avatar
gsell committed
61
	free (data);
gsell's avatar
gsell committed
62
        H5CloseFile (file);
gsell's avatar
gsell committed
63 64
	MPI_Finalize ();
        return 0;
gsell's avatar
gsell committed
65
}