1
2!=============================================================================
3!
4! Note that this code was adapted from fanc/pntf_test.F, which was written by:
5! John Tannahill, LLNL
6!
7! This a multi-variable write test on Fortran.
8!
9! This code writes the array, tt(k)(j)(i), into the file 'testfile.nc'. It
10! then reads the array from the file, and compares it with the original
11! values.
12!
13! i=longitude, j=latitude, k=level
14!
15! $Id: mcoll_testf.f90 2512 2016-09-29 01:29:37Z wkliao $
16!
17!=============================================================================
18
19      program Mcoll_Testf
20
21      use mpi
22      use pnetcdf
23      implicit none
24
25!     -----------------------
26!     Parameter declarations.
27!     -----------------------
28
29      integer NREADS, NWRITES
30      parameter (NREADS = 5, NWRITES = 5 )
31      ! number of read samples
32      ! number of write samples
33      INTEGER(KIND=MPI_OFFSET_KIND) TOTSIZ_3D(3) ! global sizes of 3D field
34
35!     ----------------------
36!     Variable declarations.
37!     ----------------------
38
39      logical reorder
40      logical isperiodic(3)
41      integer comm_cart                   ! Cartesian communicator
42      integer err, ierr, get_args
43      INTEGER(KIND=MPI_OFFSET_KIND) istart, jstart, kstart      ! offsets of 3D field
44      INTEGER(KIND=MPI_OFFSET_KIND) locsiz
45      integer mype                        ! rank in comm_cart
46      integer totpes                      ! total number of PEs
47
48      INTEGER(KIND=MPI_OFFSET_KIND) locsiz_3d(3)  ! local sizes of 3D fields
49      integer pe_coords(3)                ! Cartesian PE coords
50      integer numpes(3)                   ! number of PEs along axes;
51                                          !   determined by MPI where a
52                                          !   zero is specified
53      integer rank, Write_File
54      character(len=256) :: filename, cmd, msg
55
56      real*4  filsiz
57      real*4  rdt_l(2)
58      real*4  wrt_g(2)
59      real*4  wrt_l(2)
60
61      real*4  wrates_g(2)
62      real*4  wrates_l(2)
63
64      data reorder / .false. /
65      data isperiodic / .false., .false., .false. /
66      data numpes / 1, 1, 0 /
67!      data TOTSIZ_3D / 256, 256, 256 /
68      data TOTSIZ_3D / 8, 8, 8 /
69
70!     ----------------
71!     Begin execution.
72!     ----------------
73
74      call MPI_Init (ierr)
75      call MPI_Comm_Size(MPI_COMM_WORLD, totpes, ierr)
76      call MPI_Comm_rank(MPI_COMM_WORLD, rank, ierr)
77
78      if (rank .EQ. 0) then
79          filename = 'testfile.nc'
80          err = get_args(cmd, filename)
81      endif
82      call MPI_Bcast(err, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
83      if (err .EQ. 0) goto 999
84
85      call MPI_Bcast(filename, 256, MPI_CHARACTER, 0, MPI_COMM_WORLD, ierr)
86
87      call MPI_Dims_Create (totpes, 3, numpes, ierr)
88
89      call MPI_Cart_Create(MPI_COMM_WORLD, 3, numpes, isperiodic, &
90                           reorder, comm_cart, ierr)
91
92      call MPI_Comm_Rank (comm_cart, mype, ierr)
93
94      call MPI_Cart_Coords (comm_cart, mype, 3, pe_coords, ierr)
95
96      rdt_l(1) = 1.0e38
97      rdt_l(2) = 1.0e38
98      wrt_l(1) = 1.0e38
99      wrt_l(2) = 1.0e38
100!      rdt_l(:) = Huge (rdt_l)   ! initialize for timing
101!      wrt_l(:) = Huge (wrt_l)
102
103!     ----------------------------------------
104!     Determine local size for tt (locsiz_3d).
105!     ----------------------------------------
106
107!     ===============
108      call Find_Locnx(TOTSIZ_3D(1), pe_coords(1), numpes(1), &
109                      locsiz_3d(1), istart)
110      call Find_Locnx(TOTSIZ_3D(2), pe_coords(2), numpes(2), &
111                      locsiz_3d(2), jstart)
112      call Find_Locnx(TOTSIZ_3D(3), pe_coords(3), numpes(3), &
113                      locsiz_3d(3), kstart)
114!     ===============
115
116!     -------------------------------
117!     Compute file size in 1d6 bytes.
118!     -------------------------------
119
120      filsiz = real((TOTSIZ_3D(1) * TOTSIZ_3D(2) * TOTSIZ_3D(3)) * &
121               1.0d-6 * 4.0d0)
122
123!     -------------------------------------
124!     Print data decomposition information.
125!     -------------------------------------
126
127!      if (mype == 0) Write (6,900)
128
129      call MPI_Barrier (comm_cart, ierr)
130
131!      Write (6, 902)
132!     &  mype, pe_coords(1), pe_coords(2), pe_coords(3),
133!     &  TOTSIZ_3D(1), TOTSIZ_3D(2), TOTSIZ_3D(3),
134!     &  locsiz_3d(1), locsiz_3d(2), locsiz_3d(3),
135!     &  kstart, jstart, istart
136
137! 900  format ("mype  pe_coords    totsiz_3d         locsiz_3d       ",&
138!              "kstart,jstart,istart")
139! 902  format (i3,3x,i2,1x,i2,1x,i2,2x,i4,1x,i4,1x,i4,4x,i4,1x,i4,1x, &
140!              i4,3x,i6,1x,i6,1x,i6)
141
142!     -------------------------
143!     Write and then read back.
144!     -------------------------
145
146      locsiz = locsiz_3d(1) * locsiz_3d(2) * locsiz_3d(3)
147
148!     ===============
149      ierr = Write_File(filename, NWRITES, comm_cart, &
150                      istart, jstart, kstart, locsiz, locsiz_3d,  &
151                      TOTSIZ_3D, wrt_l)
152      if (ierr .NE. NF90_NOERR) then
153          write(6,*) trim(nf90mpi_strerror(ierr))
154          goto 999
155      endif
156!!!   Write (6,*) wrt_l(1), wrt_l(2)
157
158!     ----------------------------
159!     Compute and print I/O rates.
160!     ----------------------------
161
162      wrates_l(1) = filsiz / wrt_l(2)               ! write rate
163      wrates_l(2) = filsiz / (wrt_l(1) + wrt_l(2))  ! effective write rate
164
165      call MPI_Allreduce(wrates_l, wrates_g, 2, MPI_REAL, MPI_MIN, &
166                         comm_cart, ierr)
167      call MPI_Allreduce(wrt_l,    wrt_g,    2, MPI_REAL, MPI_MAX, &
168                         comm_cart, ierr)
169
170!      if (mype == 0) then
171!        Write (6,905) filsiz
172!        Write (6,910) wrates_g(1), wrates_g(2)
173!      end if
174
175! 905  format ("File size: ", e10.3, " MB")
176! 910  format ("    Write: ", f9.3, " MB/s  (eff., ", f9.3, " MB/s)")
177! 915  format ("    Read : ", f9.3, " MB/s  (eff., ", f9.3, " MB/s)")
178! 920  format ("Total number PEs: ", i4)
179! 922  format (e11.3, e11.3, f9.3, e11.3, e11.3, f9.3)
180
181      call MPI_Comm_Free (comm_cart, ierr)
182
183      if (rank .EQ. 0) then
184          msg='*** TESTING F90 '//trim(cmd)//' for nf90mpi_iput_var API'
185          call pass_fail(0, msg)
186      endif
187
188 999  call MPI_Finalize(ierr)
189
190      end program Mcoll_Testf
191
192!     ------------
193
194
195      integer function Write_File(filename, nwrites, comm_cart, &
196                            istart, jstart, kstart, locsiz,  &
197                            locsiz_3d, totsiz_3d, wrt_l)
198
199      use mpi
200      use pnetcdf
201      implicit none
202
203!     ----------------------
204!     Argument declarations.
205!     ----------------------
206
207      character (len=*) filename
208      integer nwrites
209      integer comm_cart
210      INTEGER(KIND=MPI_OFFSET_KIND) istart, jstart, kstart
211      INTEGER(KIND=MPI_OFFSET_KIND) locsiz
212      INTEGER(KIND=MPI_OFFSET_KIND) locsiz_3d(3)
213      INTEGER(KIND=MPI_OFFSET_KIND) totsiz_3d(3)
214      real*4  wrt_l(2)
215
216!     ----------------------------
217!     Local variable declarations.
218!     ----------------------------
219
220      integer ierr, info
221      integer lon_id, lat_id, lev_id
222      integer ncid
223      integer nw
224      integer tt1_id
225      integer req(nwrites)
226      integer stat(nwrites)
227      INTEGER(KIND=MPI_OFFSET_KIND) start_3d(3)
228      INTEGER(KIND=MPI_OFFSET_KIND) count_3d(3)
229      integer dim_id(3)
230      double precision  t1, t2, t3
231      integer max_loc_size
232      parameter( max_loc_size = 20000000 )
233      ! real*4  tt1(max_loc_size)   ! Need tt(locsiz)
234      real*4, dimension(locsiz_3d(1), locsiz_3d(2), locsiz_3d(3)) :: tt1
235
236      if (locsiz .gt. MAX_LOC_SIZE) then
237         print *, 'locsiz = ', locsiz, ' larger than MAX_LOC_SIZE'
238         stop
239      endif
240!     ----------------
241!     Begin execution.
242!     ----------------
243
244!      start_3d(1:3) = (/ kstart, jstart, istart /)
245!      count_3d(:)   = locsiz_3d(:)
246      start_3d(1) = istart
247      start_3d(2) = jstart
248      start_3d(3) = kstart
249      count_3d(1) = locsiz_3d(1)
250      count_3d(2) = locsiz_3d(2)
251      count_3d(3) = locsiz_3d(3)
252
253!       ==============
254        call Get_Field(istart, jstart, kstart, locsiz_3d, &
255                       totsiz_3d, tt1)
256
257        call MPI_Barrier (comm_cart, ierr)
258        t1 = MPI_Wtime ( )
259
260!       =================
261        call MPI_Info_create(info, ierr)
262        ! call MPI_Info_set(info, "romio_pvfs2_posix_write","enable",ierr)
263
264        Write_File = nf90mpi_create(comm_cart, filename, NF90_CLOBBER, &
265                                    info, ncid)
266        if (Write_File .NE. NF90_NOERR) return
267
268        call MPI_Info_free(info, ierr)
269
270!       ==================
271        Write_File = nf90mpi_def_dim(ncid, "level",     totsiz_3d(1)*nwrites, &
272                             lon_id)
273        if (Write_File .NE. NF90_NOERR) return
274        Write_File = nf90mpi_def_dim(ncid, "latitude",  totsiz_3d(2), lat_id)
275        if (Write_File .NE. NF90_NOERR) return
276        Write_File = nf90mpi_def_dim(ncid, "longitude", totsiz_3d(3), lev_id)
277        if (Write_File .NE. NF90_NOERR) return
278!       ==================
279
280        dim_id(1) = lon_id
281        dim_id(2) = lat_id
282        dim_id(3) = lev_id
283
284!       ==================
285        Write_File = nf90mpi_def_var(ncid, "tt1", NF90_REAL, dim_id, tt1_id)
286        if (Write_File .NE. NF90_NOERR) return
287
288!       =================
289        Write_File = nf90mpi_enddef (ncid)
290        if (Write_File .NE. NF90_NOERR) return
291!       =================
292
293        t2 = MPI_Wtime ( )
294
295      do nw = 1, nwrites
296         Write_File = nf90mpi_iput_var(ncid, tt1_id, tt1, req(nw), &
297                                 start_3d, count_3d)
298        if (Write_File .NE. NF90_NOERR) return
299
300         start_3d(1) = start_3d(1) + count_3d(1)
301      end do
302
303      Write_File = nf90mpi_wait_all(ncid, nwrites, req, stat)
304      if (Write_File .NE. NF90_NOERR) return
305
306!       ================
307      Write_File = nf90mpi_close (ncid)
308      if (Write_File .NE. NF90_NOERR) return
309!       ================
310
311! 900  format ("mynod:", i1, " reqid : ", i1)
312
313      call MPI_Barrier (comm_cart, ierr)
314      t3 = MPI_Wtime ( )
315
316      if (t2 - t1 < wrt_l(1)) wrt_l(1) = real(t2 - t1)
317      if (t3 - t2 < wrt_l(2)) wrt_l(2) = real(t3 - t2)
318
319      end
320
321!     ------------
322
323      subroutine Find_Locnx(nx, mype, totpes, locnx, ibegin)
324
325      use mpi
326      implicit none
327
328!     ----------------------
329!     Argument declarations.
330!     ----------------------
331
332      INTEGER(KIND=MPI_OFFSET_KIND) nx
333      integer mype
334      integer totpes
335      INTEGER(KIND=MPI_OFFSET_KIND) locnx
336      INTEGER(KIND=MPI_OFFSET_KIND) ibegin
337
338!     ----------------------------
339!     Local variable declarations.
340!     ----------------------------
341
342      INTEGER(KIND=MPI_OFFSET_KIND) iremain
343
344!     ----------------
345!     Begin execution.
346!     ----------------
347
348      locnx = nx / totpes
349
350      iremain = nx - (totpes * locnx)
351
352      if (mype < iremain) locnx = locnx + 1
353
354      ibegin = mype * (nx / totpes) + iremain + 1
355
356      if (mype < iremain) ibegin = ibegin + (mype - iremain)
357
358      end
359
360!     ------------
361
362      subroutine Get_Field(istart, jstart, kstart, locsiz_3d, &
363                           totsiz_3d, tt)
364
365      use mpi
366      implicit none
367
368!     ----------------------
369!     Argument declarations.
370!     ----------------------
371
372      INTEGER(KIND=MPI_OFFSET_KIND) istart, jstart, kstart
373      INTEGER(KIND=MPI_OFFSET_KIND) locsiz_3d(3)
374      INTEGER(KIND=MPI_OFFSET_KIND) totsiz_3d(3)
375      real*4, dimension(locsiz_3d(1),locsiz_3d(2),locsiz_3d(3)) ::tt
376
377!     ----------------------------
378!     Local variable declarations.
379!     ----------------------------
380
381      integer ii, jj, kk
382      integer ind
383
384!     ----------------
385!     Begin execution.
386!     ----------------
387
388      ind = 1
389
390      do kk = 1, int(locsiz_3d(3))
391        do jj = 1, int(locsiz_3d(2))
392          do ii = 1, int(locsiz_3d(1))
393
394             tt(ii,jj,kk) = real( &
395               (istart-1 +(ii - 1) + 1 + totsiz_3d(3)*(jstart-1 +  &
396                       (jj - 1) + totsiz_3d(2)*(kstart-1 +  &
397                       (kk-1)))) * 1.0d-3)
398             ind = ind + 1
399
400          end do
401        end do
402      end do
403
404      end
405
406!     ------------
407
408      subroutine Compare_Vec(comm_cart, locsiz, tt, buf)
409
410      use mpi
411      implicit none
412
413!     ----------------------
414!     Argument declarations.
415!     ----------------------
416
417      integer comm_cart
418      INTEGER(KIND=MPI_OFFSET_KIND) locsiz
419      real*4  tt (locsiz)
420      real*4  buf(locsiz)
421
422!     ----------------------------
423!     Local variable declarations.
424!     ----------------------------
425
426      integer ierr
427      integer ii
428      real*4  delmax(1), delmin(1), delta
429      real*4  diff
430      real*4  wr(5)
431      real*4  ws(5)
432
433!     ----------------
434!     Begin execution.
435!     ----------------
436
437      ws(1) = 0.0d0      ! diff
438      ws(2) = 0.0d0      ! sumsq
439      ws(3) = real(locsiz)     ! locsiz
440      ws(4) = 0.0d0      ! delmax
441      ws(5) = 1.0d38     ! Huge (ws)  ! delmin
442
443      do ii = 1, int(locsiz)
444        delta = (tt(ii) - buf(ii)) * (tt(ii) - buf(ii))
445        ws(1) = ws(1) + delta
446        ws(2) = ws(2) + tt(ii) * tt(ii)
447        if (delta > ws(4)) ws(4) = delta
448        if (delta < ws(5)) ws(5) = delta
449      end do
450
451      call MPI_Allreduce(ws,    wr,     3, MPI_REAL, MPI_SUM, &
452                         comm_cart, ierr)
453      call MPI_Allreduce(ws(4), delmax, 1, MPI_REAL, MPI_MAX, &
454                         comm_cart, ierr)
455      call MPI_Allreduce(ws(5), delmin(1), 1, MPI_REAL, MPI_MIN, &
456                         comm_cart, ierr)
457
458      diff      = Sqrt (wr(1) / wr(2))         ! normalized error
459      delmax(1) = Sqrt (wr(3) * delmax(1)/wr(2))  ! normalized max difference
460      delmin(1) = Sqrt (wr(3) * delmin(1)/wr(2))  ! normalized min difference
461
462      end
463
464