1!
2!  Copyright (C) 2012, Northwestern University and Argonne National Laboratory
3!  See COPYRIGHT notice in top-level directory.
4!
5! $Id: varn_real.f90 2638 2016-11-18 14:51:02Z wkliao $
6
7!
8! This example shows how to use a single call of nf90mpi_put_varn_all()
9! to write a sequence of one-element requests with arbitrary array indices.
10!
11! The compile and run commands are given below, together with an ncmpidump of
12! the output file.
13!
14!    % mpif90 -O2 -o varn_real varn_real.f90 -lpnetcdf
15!    % mpiexec -n 4 ./varn_real /pvfs2/wkliao/testfile.nc
16!    % ncmpidump /pvfs2/wkliao/testfile.nc
17!    netcdf testfile {
18!    // file format: CDF-5 (big variables)
19!    dimensions:
20!             Y = 4 ;
21!             X = 10 ;
22!    variables:
23!             float var(Y, X) ;
24!    data:
25!
26!     var =
27!       3, 3, 3, 1, 1, 0, 0, 2, 1, 1,
28!       0, 2, 2, 2, 3, 1, 1, 2, 2, 2,
29!       1, 1, 2, 3, 3, 3, 0, 0, 1, 1,
30!       0, 0, 0, 2, 1, 1, 1, 3, 3, 3 ;
31!    }
32!
33
34      subroutine check(err, message)
35          use mpi
36          use pnetcdf
37          implicit none
38          integer err
39          character(len=*) message
40          character(len=128) msg
41
42          ! It is a good idea to check returned value for possible error
43          if (err .NE. NF90_NOERR) then
44              write(6,*) message, trim(nf90mpi_strerror(err))
45              msg = '*** TESTING F90 varn_real.f90 for varn API '
46              call pass_fail(1, msg)
47              ! call MPI_Abort(MPI_COMM_WORLD, -1, err)
48          end if
49      end subroutine check
50
51      program main
52          use mpi
53          use pnetcdf
54          implicit none
55
56          integer NDIMS
57          PARAMETER(NDIMS=2)
58
59          character(LEN=256) filename, cmd, msg
60          integer rank, nprocs, err, ierr, num_reqs, get_args
61          integer ncid, cmode, varid, dimid(2), y, x, i, j, nerrs
62          integer, allocatable :: req_lens(:)
63          real, allocatable :: buffer(:), buffer2D(:,:)
64          real oneReal
65          integer(kind=MPI_OFFSET_KIND) NY, NX
66          integer(kind=MPI_OFFSET_KIND) w_len, w_req_len
67          integer(kind=MPI_OFFSET_KIND), allocatable :: starts(:,:)
68          integer(kind=MPI_OFFSET_KIND), allocatable :: counts(:,:)
69          integer(kind=MPI_OFFSET_KIND) malloc_size, sum_size
70
71          NY = 4
72          NX = 10
73
74          call MPI_Init(ierr)
75          call MPI_Comm_rank(MPI_COMM_WORLD, rank, ierr)
76          call MPI_Comm_size(MPI_COMM_WORLD, nprocs, ierr)
77
78          ! take filename from command-line argument if there is any
79          if (rank .EQ. 0) then
80              filename = "testfile.nc"
81              err = get_args(cmd, filename)
82          endif
83          call MPI_Bcast(err, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
84          if (err .EQ. 0) goto 999
85
86          call MPI_Bcast(filename, 256, MPI_CHARACTER, 0, MPI_COMM_WORLD, ierr)
87
88          nerrs = 0
89
90          if (.FALSE. .AND. nprocs .NE. 4 .AND. rank .EQ. 0) &
91              print*,'Warning: ',trim(cmd),' is intended to run on ', &
92                     '4 processes'
93
94          ! create file, truncate it if exists
95          cmode = IOR(NF90_CLOBBER, NF90_64BIT_DATA)
96          err = nf90mpi_create(MPI_COMM_WORLD, filename, cmode, &
97                             MPI_INFO_NULL, ncid)
98          call check(err, 'In nf90mpi_create: ')
99
100          ! create a global array of size NY * NX */
101          err = nf90mpi_def_dim(ncid, "Y", NY, dimid(2))
102          call check(err, 'In nf90mpi_def_dim Y: ')
103          err = nf90mpi_def_dim(ncid, "X", NX, dimid(1))
104          call check(err, 'In nf90mpi_def_dim X: ')
105          err = nf90mpi_def_var(ncid, "var", NF90_FLOAT, dimid, varid)
106          call check(err, 'In nf90mpi_def_var var: ')
107          err = nf90mpi_enddef(ncid)
108          call check(err, 'In nf90mpi_enddef: ')
109
110          ! pick arbitrary numbers of requests for 4 processes
111          num_reqs = 0
112          if (rank .EQ.  0) then
113              num_reqs = 8
114          elseif (rank .EQ. 1) then
115              num_reqs = 13
116          elseif (rank .EQ. 2) then
117              num_reqs = 9
118          elseif (rank .EQ. 3) then
119              num_reqs = 10
120          endif
121
122          ! Note that in Fortran, array indices start with 1
123          ALLOCATE(starts(NDIMS, num_reqs+1))
124          ALLOCATE(counts(NDIMS, num_reqs))
125
126          ! assign arbitrary starts
127          y=2
128          x=1
129          if (rank .EQ. 0) then
130              starts(y, 1) = 1; starts(x, 1) = 6
131              starts(y, 2) = 2; starts(x, 2) = 1
132              starts(y, 3) = 3; starts(x, 3) = 7
133              starts(y, 4) = 4; starts(x, 4) = 1
134              starts(y, 5) = 1; starts(x, 5) = 7
135              starts(y, 6) = 3; starts(x, 6) = 8
136              starts(y, 7) = 4; starts(x, 7) = 2
137              starts(y, 8) = 4; starts(x, 8) = 3
138              ! rank 0 is writing the following locations: ("-" means skip)
139              !   -  -  -  -  -  0  0  -  -  -
140              !   0  -  -  -  -  -  -  -  -  -
141              !   -  -  -  -  -  -  0  0  -  -
142              !   0  0  0  -  -  -  -  -  -  -
143          elseif (rank .EQ. 1) then
144              starts(y,  1) = 1; starts(x,  1) = 4
145              starts(y,  2) = 1; starts(x,  2) = 9
146              starts(y,  3) = 2; starts(x,  3) = 6
147              starts(y,  4) = 3; starts(x,  4) = 1
148              starts(y,  5) = 3; starts(x,  5) = 9
149              starts(y,  6) = 4; starts(x,  6) = 5
150              starts(y,  7) = 1; starts(x,  7) = 5
151              starts(y,  8) = 1; starts(x,  8) = 10
152              starts(y,  9) = 2; starts(x,  9) = 7
153              starts(y, 10) = 3; starts(x, 10) = 2
154              starts(y, 11) = 3; starts(x, 11) = 10
155              starts(y, 12) = 4; starts(x, 12) = 6
156              starts(y, 13) = 4; starts(x, 13) = 7
157              ! rank 1 is writing the following locations: ("-" means skip)
158              !   -  -  -  1  1  -  -  -  1  1
159              !   -  -  -  -  -  1  1  -  -  -
160              !   1  1  -  -  -  -  -  -  1  1
161              !   -  -  -  -  1  1  1  -  -  -
162          elseif (rank .EQ. 2) then
163              starts(y, 1) = 1; starts(x, 1) = 8
164              starts(y, 2) = 2; starts(x, 2) = 2
165              starts(y, 3) = 2; starts(x, 3) = 8
166              starts(y, 4) = 3; starts(x, 4) = 3
167              starts(y, 5) = 4; starts(x, 5) = 4
168              starts(y, 6) = 2; starts(x, 6) = 3
169              starts(y, 7) = 2; starts(x, 7) = 9
170              starts(y, 8) = 2; starts(x, 8) = 4
171              starts(y, 9) = 2; starts(x, 9) = 10
172              ! rank 2 is writing the following locations: ("-" means skip)
173              !   -  -  -  -  -  -  -  2  -  -
174              !   -  2  2  2  -  -  -  2  2  2
175              !   -  -  2  -  -  -  -  -  -  -
176              !   -  -  -  2  -  -  -  -  -  -
177          elseif (rank .EQ. 3) then
178              starts(y,  1) = 1; starts(x,  1) = 1
179              starts(y,  2) = 2; starts(x,  2) = 5
180              starts(y,  3) = 3; starts(x,  3) = 4
181              starts(y,  4) = 4; starts(x,  4) = 8
182              starts(y,  5) = 1; starts(x,  5) = 2
183              starts(y,  6) = 3; starts(x,  6) = 5
184              starts(y,  7) = 4; starts(x,  7) = 9
185              starts(y,  8) = 1; starts(x,  8) = 3
186              starts(y,  9) = 3; starts(x,  9) = 6
187              starts(y, 10) = 4; starts(x, 10) = 10
188              ! rank 3 is writing the following locations: ("-" means skip)
189              !   3  3  3  -  -  -  -  -  -  -
190              !   -  -  -  -  3  -  -  -  -  -
191              !   -  -  -  3  3  3  -  -  -  -
192              !   -  -  -  -  -  -  -  3  3  3
193          endif
194          counts = 1
195
196          ALLOCATE(req_lens(num_reqs))
197
198          ! allocate I/O buffer and initialize its contents
199          w_len = 0
200          do i=1, num_reqs
201             w_req_len = 1
202             do j=1, NDIMS
203                w_req_len = w_req_len * counts(j, i)
204             enddo
205             req_lens(i) = INT(w_req_len)
206             w_len = w_len + w_req_len
207          enddo
208          ALLOCATE(buffer(w_len))
209          ALLOCATE(buffer2D(3, w_len/3+1))
210          ! Note on 2D buffer, put_varn will fill in the first w_len
211          ! elements in buffer2D, no matter the shape of buffer2D is
212
213          buffer   = rank
214          buffer2D = rank
215
216          ! varn write a 2D buffer in memory to a 2D array in file
217          err = nf90mpi_put_varn_all(ncid, varid, buffer2D, num_reqs, &
218                                     starts, counts)
219          call check(err, 'In nf90mpi_put_varn_all: ')
220
221          ! varn write a 1D buffer in memory to a 2D array in file
222          err = nf90mpi_put_varn_all(ncid, varid, buffer, num_reqs, &
223                                     starts, counts)
224          call check(err, 'In nf90mpi_put_varn_all: ')
225
226          if (nprocs .GT. 4) call MPI_Barrier(MPI_COMM_WORLD, err)
227
228          ! read a scalar back and check the content
229          oneReal = -1.0  ! a scalar
230          if (rank .GE. 4) starts = 1_MPI_OFFSET_KIND
231          err = nf90mpi_get_varn_all(ncid, varid, oneReal, starts)
232          call check(err, ' 22 In nf90mpi_get_varn_all: ')
233
234 995      format(A,I2,A,F4.1)
235          if (rank .LT. 4 .AND. oneReal .NE. rank) then
236              print 995, "Error: expecting OneReal=",rank, &
237                         " but got", oneReal
238              nerrs = nerrs + 1
239          endif
240
241          ! varn read a 2D array in file to a 2D buffer in memory
242          buffer2D = -1.0;
243          err = nf90mpi_get_varn_all(ncid, varid, buffer2D, num_reqs, &
244                                     starts, counts)
245          call check(err, 'In nf90mpi_get_varn_all: ')
246
247 996      format(A,I2,A,I2,A,I2,A,F4.1)
248          do i=1, INT(w_len/3)
249             do j=1, 3
250                if (buffer2D(j,i) .NE. rank) then
251                    print 996, "Error: expecting buffer2D(",j,",",i,")=", &
252                               rank, " but got", buffer2D(j,i)
253                    nerrs = nerrs + 1
254                endif
255             enddo
256          enddo
257
258          ! varn read a 2D array in file to a 1D buffer in memory
259          buffer = -1.0;
260          err = nf90mpi_get_varn_all(ncid, varid, buffer, num_reqs, &
261                                     starts, counts)
262          call check(err, 'In nf90mpi_get_varn_all: ')
263
264 997      format(A,I2,A,I2,A,F4.1)
265          do i=1, INT(w_len)
266             if (buffer(i) .NE. rank) then
267                 print 997, "Error: expecting buffer(",i,")=",rank, &
268                            " but got", buffer(i)
269                 nerrs = nerrs + 1
270             endif
271          enddo
272
273          err = nf90mpi_close(ncid);
274          call check(err, 'In nf90mpi_close: ')
275
276          DEALLOCATE(req_lens);
277          DEALLOCATE(buffer);
278          DEALLOCATE(buffer2D);
279          DEALLOCATE(starts);
280          DEALLOCATE(counts);
281
282          ! check if there is any PnetCDF internal malloc residue
283 998      format(A,I13,A)
284          err = nf90mpi_inq_malloc_size(malloc_size)
285          if (err == NF90_NOERR) then
286              call MPI_Reduce(malloc_size, sum_size, 1, MPI_INTEGER8, &
287                              MPI_SUM, 0, MPI_COMM_WORLD, ierr)
288              if (rank .EQ. 0 .AND. sum_size .GT. 0_MPI_OFFSET_KIND) print 998, &
289                  'heap memory allocated by PnetCDF internally has ',  &
290                  sum_size, ' bytes yet to be freed'
291          endif
292
293          if (rank .eq. 0) then
294              msg = '*** TESTING F90 '//trim(cmd)//' for varn API '
295              call pass_fail(nerrs, msg)
296          endif
297
298 999      call MPI_Finalize(ierr)
299
300      end program
301
302