1!
2!  Copyright (C) 2013, Northwestern University and Argonne National Laboratory
3!  See COPYRIGHT notice in top-level directory.
4!
5! $Id: test_write.F90 2011 2015-02-14 19:48:35Z wkliao $
6!
7
8! Test nf90mpi_create
9!    For mode in NF90_NOCLOBBER, NF90_CLOBBER do:
10!       create netcdf file 'scratch.nc' with no data, close it
11!       test that it can be opened, do nf90mpi_inq to check nvars = 0, etc.
12!    Try again in NF90_NOCLOBBER mode, check error return
13! On exit, delete this file
14        subroutine test_nf90mpi_create()
15        use pnetcdf
16        implicit        none
17#include "tests.inc"
18
19        integer clobber         !/* 0 for NF90_NOCLOBBER, 1 for NF90_CLOBBER */
20        integer err
21        integer ncid
22        integer ndims           !/* number of dimensions */
23        integer nvars           !/* number of variables */
24        integer ngatts          !/* number of global attributes */
25        integer recdim          !/* id of unlimited dimension */
26        integer flags
27        integer nok
28
29        flags = IOR(NF90_NOCLOBBER, extra_flags)
30        nok = 0
31        do 1, clobber = 0, 1
32            err = nf90mpi_create(comm, scratch, flags,  info, &
33                               ncid)
34            if (err .ne. NF90_NOERR) then
35                call errore('nf90mpi_create: ', err)
36            end if
37            nok = nok + 1
38            err = nf90mpi_close(ncid)
39            if (err .ne. NF90_NOERR) then
40                call errore('nf90mpi_close: ', err)
41            end if
42            err = nf90mpi_open(comm, scratch, NF90_NOWRITE, info,  &
43                             ncid)
44            if (err .ne. NF90_NOERR) then
45                call errore('nf90mpi_open: ', err)
46            end if
47            err = nf90mpi_inquire(ncid, ndims, nvars, ngatts, recdim)
48            if (err .ne. NF90_NOERR) then
49                call errore('nf90mpi_inquire: ', err)
50            else if (ndims .ne. 0) then
51                call errori( &
52                    'nf90mpi_inquire: wrong number of dimensions returned, ', &
53                    ndims)
54            else if (nvars .ne. 0) then
55                call errori( &
56                    'nf90mpi_inquire: wrong number of variables returned, ', &
57                    nvars)
58            else if (ngatts .ne. 0) then
59                call errori( &
60                    'nf90mpi_inquire: wrong number of global atts returned, ', &
61                    ngatts)
62            else if (recdim .ge. 1) then
63                call errori( &
64                    'nf90mpi_inquire: wrong record dimension ID returned, ', &
65                    recdim)
66            end if
67            err = nf90mpi_close(ncid)
68            if (err .ne. NF90_NOERR) then
69                call errore('nf90mpi_close: ', err)
70            end if
71
72            flags = IOR(NF90_CLOBBER, extra_flags)
731       continue
74
75        flags = IOR(NF90_NOCLOBBER, extra_flags)
76        err = nf90mpi_create(comm, scratch, flags,  info, &
77                           ncid)
78        if (err .ne. NF90_EEXIST) then
79            call errore('attempt to overwrite file: ', err)
80        end if
81        nok = nok + 1
82        err = nf90mpi_delete(scratch, info)
83        if (err .ne. NF90_NOERR) then
84            call errori('delete of scratch file failed: ', err)
85        end if
86        call print_nok(nok)
87        end
88
89
90! Test nf90mpi_redef
91! (In fact also tests nf90mpi_enddef - called from test_nf90mpi_enddef)
92!    BAD_ID
93!    attempt redef (error) & enddef on read-only file
94!    create file, define dims & vars.
95!    attempt put var (error)
96!    attempt redef (error) & enddef.
97!    put vars
98!    attempt def new dims (error)
99!    redef
100!    def new dims, vars.
101!    put atts
102!    enddef
103!    put vars
104!    close
105!    check file: vars & atts
106        subroutine test_nf90mpi_redef()
107        use pnetcdf
108        implicit        none
109#include "tests.inc"
110        integer         title_len
111        parameter       (title_len = 9)
112
113        integer                 ncid            !/* netcdf id */
114        integer                 dimid           !/* dimension id */
115        integer                 vid             !/* variable id */
116        integer                 err, flags
117        character*(title_len)   title
118        doubleprecision         var
119        character*(NF90_MAX_NAME) name
120        integer(kind=MPI_OFFSET_KIND)                 start(1)
121        integer(kind=MPI_OFFSET_KIND)                 length
122        integer                 intlen
123        integer                 dimids(1)
124        integer                 nok
125
126        nok = 0
127
128        title = 'Not funny'
129
130        ! BAD_ID tests
131        err = nf90mpi_redef(BAD_ID)
132        if (err .ne. NF90_EBADID) then
133            call errore('bad ncid: ', err)
134        endif
135        nok = nok + 1
136        err = nf90mpi_enddef(BAD_ID)
137        if (err .ne. NF90_EBADID) then
138            call errore('bad ncid: ', err)
139        endif
140        nok = nok + 1
141
142        ! read-only tests
143        err = nf90mpi_open(comm, testfile, NF90_NOWRITE, info, &
144                         ncid)
145        if (err .ne. NF90_NOERR) &
146            call errore('nf90mpi_open: ', err)
147        err = nf90mpi_redef(ncid)
148        if (err .ne. NF90_EPERM) then
149            call errore('nf90mpi_redef in NF90_NOWRITE mode: ', err)
150        endif
151        nok = nok + 1
152        err = nf90mpi_enddef(ncid)
153        if (err .ne. NF90_ENOTINDEFINE) then
154            call errore('nf90mpi_redef in NF90_NOWRITE mode: ', err)
155        endif
156        nok = nok + 1
157        err = nf90mpi_close(ncid)
158        if (err .ne. NF90_NOERR)  &
159            call errore('nf90mpi_close: ', err)
160
161!           /* tests using scratch file */
162        flags = IOR(NF90_NOCLOBBER, extra_flags)
163        err = nf90mpi_create(comm, scratch, flags, info, &
164                           ncid)
165        if (err .ne. NF90_NOERR) then
166            call errore('nf90mpi_create: ', err)
167            return
168        end if
169        call def_dims(ncid)
170        call def_vars(ncid)
171        call put_atts(ncid)
172        err = nf90mpi_inq_varid(ncid, 'd', vid)
173        if (err .ne. NF90_NOERR)  &
174            call errore('nf90mpi_inq_varid: ', err)
175        var = 1.0
176!       should not enter indep mode in define mode
177        err = nf90mpi_begin_indep_data(ncid)
178        if (err .ne. NF90_EINDEFINE) &
179          call errore('nf90mpi_begin_indep_data... in define mode: ', err)
180        start = 0
181        err = nf90mpi_put_var(ncid, vid, var, start)
182        if (err .ne. NF90_EINDEFINE) &
183            call errore('nf90mpi_put_var... in define mode: ', err)
184        err = nf90mpi_end_indep_data(ncid)
185        if (err .ne. NF90_ENOTINDEP) &
186          call errore('nf90mpi_end_indep_data... not in indep mode: ', err)
187        err = nf90mpi_redef(ncid)
188        if (err .ne. NF90_EINDEFINE) then
189            call errore('nf90mpi_redef in define mode: ', err)
190        endif
191        nok = nok + 1
192        err = nf90mpi_enddef(ncid)
193        if (err .ne. NF90_NOERR) &
194            call errore('nf90mpi_enddef: ', err)
195        call put_vars(ncid)
196        length = 8
197        err = nf90mpi_def_dim(ncid, 'abc', length, dimid)
198        if (err .ne. NF90_ENOTINDEFINE) &
199            call errore('nf90mpi_def_dim in define mode: ', err)
200        err = nf90mpi_redef(ncid)
201        if (err .ne. NF90_NOERR) then
202            call errore('nf90mpi_redef: ', err)
203        endif
204        nok = nok + 1
205        length = 8
206        err = nf90mpi_def_dim(ncid, 'abc', length, dimid)
207        if (err .ne. NF90_NOERR) &
208            call errore('nf90mpi_def_dim: ', err)
209        dimids(1) = 0
210        err = nf90mpi_def_var(ncid, 'abc', NF90_INT, varid=vid)
211        if (err .ne. NF90_NOERR) &
212            call errore('nf90mpi_def_var: ', err)
213        length = len(title)
214        err = nf90mpi_put_att(ncid, NF90_GLOBAL, 'title', &
215                              title)
216        if (err .ne. NF90_NOERR) &
217            call errore('nf90mpi_put_att: ', err)
218        err = nf90mpi_enddef(ncid)
219        if (err .ne. NF90_NOERR) &
220            call errore('nf90mpi_enddef: ', err)
221        var = 1.0
222        err = nf90mpi_end_indep_data(ncid)
223        if (err .ne. NF90_ENOTINDEP) &
224          call errore('nf90mpi_end_indep_data: in collective mode: ',err)
225        err = nf90mpi_begin_indep_data(ncid)
226        if (err .ne. NF90_NOERR) &
227            call errore('nf90mpi_begin_indep_data: ', err)
228        err = nf90mpi_put_var(ncid, vid, var, start)
229        if (err .ne. NF90_NOERR) &
230            call errore('nf90mpi_put_var: ', err)
231        err = nf90mpi_end_indep_data(ncid)
232        if (err .ne. NF90_NOERR) &
233            call errore('nf90mpi_end_indep_data: ', err)
234        err = nf90mpi_close(ncid)
235        if (err .ne. NF90_NOERR)  &
236            call errore('nf90mpi_close: ', err)
237
238!           /* check scratch file written as expected */
239        call check_file(scratch)
240        err = nf90mpi_open(comm, scratch, NF90_NOWRITE, &
241              info, ncid)
242        if (err .ne. NF90_NOERR) &
243            call errore('nf90mpi_open: ', err)
244        err = nf90mpi_inquire_dimension(ncid, dimid, name, length)
245        if (err .ne. NF90_NOERR)  &
246            call errore('nf90mpi_inquire_dimension: ', err)
247        if (name .ne. "abc") &
248            call errori('Unexpected dim name in netCDF ', ncid)
249        if (length .ne. 8) then
250            intlen = INT(length)
251            call errori('Unexpected dim length: ', intlen)
252        end if
253        err = nf90mpi_begin_indep_data(ncid)
254        err = nf90mpi_get_var(ncid, vid, var, start)
255        if (err .ne. NF90_NOERR) &
256            call errore('nf90mpi_get_var: ', err)
257        if (var .ne. 1.0) &
258            call errori( &
259                'nf90mpi_get_var: unexpected value in netCDF ' &
260                , ncid)
261        err = nf90mpi_end_indep_data(ncid)
262        err = nf90mpi_close(ncid)
263        if (err .ne. NF90_NOERR) &
264            call errore('nf90mpi_close: ', err)
265
266        err = nf90mpi_delete(scratch, info)
267        if (err .ne. NF90_NOERR) &
268            call errori('delete failed for netCDF: ', err)
269        call print_nok(nok)
270        end
271
272! Test nf90mpi_enddef
273! Simply calls test_nf90mpi_redef which tests both nf90mpi_redef & nf90mpi_enddef
274
275        subroutine test_nf90mpi_enddef()
276        use pnetcdf
277        implicit        none
278#include "tests.inc"
279
280        call test_nf90mpi_redef
281        end
282
283
284! Test nf90mpi_sync
285!    try with bad handle, check error
286!    try in define mode, check error
287!    try writing with one handle, reading with another on same netCDF
288        subroutine test_nf90mpi_sync()
289        use pnetcdf
290        implicit        none
291#include "tests.inc"
292
293        integer ncidw         !/* netcdf id for writing */
294        integer ncidr         !/* netcdf id for reading */
295        integer err, flags
296        integer nok
297
298        nok = 0
299!           /* BAD_ID test */
300        err = nf90mpi_sync(BAD_ID)
301        if (err .ne. NF90_EBADID) then
302            call errore('bad ncid: ', err)
303        else
304            nok = nok + 1
305        endif
306
307!           /* create scratch file & try nf90mpi_sync in define mode */
308        flags = IOR(NF90_NOCLOBBER, extra_flags)
309        err = nf90mpi_create(comm, scratch, flags, info, &
310                           ncidw)
311        if (err .ne. NF90_NOERR) then
312            call errore('nf90mpi_create: ', err)
313            return
314        end if
315        err = nf90mpi_sync(ncidw)
316        if (err .ne. NF90_EINDEFINE) then
317            call errore('nf90mpi_sync called in define mode: ', err)
318        else
319            nok = nok + 1
320        endif
321
322!           /* write using same handle */
323        call def_dims(ncidw)
324        call def_vars(ncidw)
325        call put_atts(ncidw)
326        err = nf90mpi_enddef(ncidw)
327        if (err .ne. NF90_NOERR) &
328            call errore('nf90mpi_enddef: ', err)
329        call put_vars(ncidw)
330        err = nf90mpi_sync(ncidw)
331        if (err .ne. NF90_NOERR) then
332            call errore('nf90mpi_sync of ncidw failed: ', err)
333        else
334            nok = nok + 1
335        endif
336
337!           /* open another handle, nf90mpi_sync, read (check) */
338        err = nf90mpi_open(comm, scratch, NF90_NOWRITE, info, &
339                         ncidr)
340        if (err .ne. NF90_NOERR) &
341            call errore('nf90mpi_open: ', err)
342        err = nf90mpi_sync(ncidr)
343        if (err .ne. NF90_NOERR) then
344            call errore('nf90mpi_sync of ncidr failed: ', err)
345        else
346            nok = nok + 1
347        endif
348        call check_dims(ncidr)
349        call check_atts(ncidr)
350        call check_vars(ncidr)
351
352!           /* close both handles */
353        err = nf90mpi_close(ncidr)
354        if (err .ne. NF90_NOERR) &
355            call errore('nf90mpi_close: ', err)
356        err = nf90mpi_close(ncidw)
357        if (err .ne. NF90_NOERR) &
358            call errore('nf90mpi_close: ', err)
359
360        err = nf90mpi_delete(scratch, info)
361        if (err .ne. NF90_NOERR) &
362            call errori('delete of scratch file failed: ', err)
363        call print_nok(nok)
364        end
365
366
367! Test nf90mpi_abort
368!    try with bad handle, check error
369!    try in define mode before anything written, check that file was deleted
370!    try after nf90mpi_enddef, nf90mpi_redef, define new dims, vars, atts
371!    try after writing variable
372        subroutine test_nf90mpi_abort()
373        use pnetcdf
374        implicit        none
375#include "tests.inc"
376
377        integer ncid          !/* netcdf id */
378        integer err, flags
379        integer ndims
380        integer nvars
381        integer ngatts
382        integer recdim
383        integer nok
384
385        nok = 0
386
387!           /* BAD_ID test */
388        err = nf90mpi_abort(BAD_ID)
389        if (err .ne. NF90_EBADID) then
390            call errore('bad ncid: status = ', err)
391        else
392            nok = nok + 1
393        endif
394
395!           /* create scratch file & try nf90mpi_abort in define mode */
396        flags = IOR(NF90_NOCLOBBER, extra_flags)
397        err = nf90mpi_create(comm, scratch, flags, info, &
398                           ncid)
399        if (err .ne. NF90_NOERR) then
400            call errore('nf90mpi_create: ', err)
401            return
402        end if
403        call def_dims(ncid)
404        call def_vars(ncid)
405        call put_atts(ncid)
406        err = nf90mpi_abort(ncid)
407        if (err .ne. NF90_NOERR) then
408            call errore('nf90mpi_abort of ncid failed: ', err)
409        else
410            nok = nok + 1
411        endif
412        err = nf90mpi_close(ncid)    !/* should already be closed */
413        if (err .ne. NF90_EBADID) &
414            call errore('bad ncid: ', err)
415        err = nf90mpi_delete(scratch, info)    !/* should already be deleted */
416        if (err .eq. NF90_NOERR) &
417            call errori('scratch file should not exist: ', err)
418
419!            create scratch file
420!            do nf90mpi_enddef & nf90mpi_redef
421!            define new dims, vars, atts
422!            try nf90mpi_abort: should restore previous state (no dims, vars, atts)
423        flags = IOR(NF90_NOCLOBBER, extra_flags)
424        err = nf90mpi_create(comm, scratch, flags, info, &
425                           ncid)
426        if (err .ne. NF90_NOERR) then
427            call errore('nf90mpi_create: ', err)
428            return
429        end if
430        err = nf90mpi_enddef(ncid)
431        if (err .ne. NF90_NOERR) &
432            call errore('nf90mpi_enddef: ', err)
433        err = nf90mpi_redef(ncid)
434        if (err .ne. NF90_NOERR) &
435            call errore('nf90mpi_redef: ', err)
436        call def_dims(ncid)
437        call def_vars(ncid)
438        call put_atts(ncid)
439        err = nf90mpi_abort(ncid)
440        if (err .ne. NF90_NOERR) then
441            call errore('nf90mpi_abort of ncid failed: ', err)
442        else
443            nok = nok + 1
444        endif
445        err = nf90mpi_close(ncid)    !/* should already be closed */
446        if (err .ne. NF90_EBADID) &
447            call errore('bad ncid: ', err)
448        err = nf90mpi_open(comm, scratch, NF90_NOWRITE, info, &
449                         ncid)
450        if (err .ne. NF90_NOERR) &
451            call errore('nf90mpi_open: ', err)
452        err = nf90mpi_inquire(ncid, ndims, nvars, ngatts, recdim)
453        if (err .ne. NF90_NOERR) &
454            call errore('nf90mpi_inquire: ', err)
455        if (ndims .ne. 0) &
456            call errori('ndims should be ', 0)
457        if (nvars .ne. 0) &
458            call errori('nvars should be ', 0)
459        if (ngatts .ne. 0) &
460            call errori('ngatts should be ', 0)
461        err = nf90mpi_close (ncid)
462        if (err .ne. NF90_NOERR) &
463            call errore('nf90mpi_close: ', err)
464
465!           /* try nf90mpi_abort in data mode - should just close */
466        flags = IOR(NF90_CLOBBER, extra_flags)
467        err = nf90mpi_create(comm, scratch, flags, info, &
468                           ncid)
469        if (err .ne. NF90_NOERR) then
470            call errore('nf90mpi_create: ', err)
471            return
472        end if
473        call def_dims(ncid)
474        call def_vars(ncid)
475        call put_atts(ncid)
476        err = nf90mpi_enddef(ncid)
477        if (err .ne. NF90_NOERR) &
478            call errore('nf90mpi_enddef: ', err)
479        call put_vars(ncid)
480        err = nf90mpi_abort(ncid)
481        if (err .ne. NF90_NOERR) then
482            call errore('nf90mpi_abort of ncid failed: ', err)
483        else
484            nok = nok + 1
485        endif
486        err = nf90mpi_close(ncid)       !/* should already be closed */
487        if (err .ne. NF90_EBADID) &
488            call errore('bad ncid: ', err)
489        call check_file(scratch)
490        err = nf90mpi_delete(scratch, info)
491        if (err .ne. NF90_NOERR) &
492            call errori('delete of scratch file failed: ', err)
493        call print_nok(nok)
494        end
495
496
497! Test nf90mpi_def_dim
498!    try with bad netCDF handle, check error
499!    try in data mode, check error
500!    check that returned id is one more than previous id
501!    try adding same dimension twice, check error
502!    try with illegal sizes, check error
503!    make sure unlimited size works, shows up in nf90mpi_inq_unlimdim
504!    try to define a second unlimited dimension, check error
505        subroutine test_nf90mpi_def_dim()
506        use pnetcdf
507        implicit        none
508#include "tests.inc"
509
510        integer ncid
511        integer err             !/* status */
512        integer i
513        integer dimid         !/* dimension id */
514        integer(kind=MPI_OFFSET_KIND) length
515        integer nok, flags
516
517        nok = 0
518
519!           /* BAD_ID test */
520        length = 8
521        err = nf90mpi_def_dim(BAD_ID, 'abc', length, dimid)
522        if (err .ne. NF90_EBADID) then
523            call errore('bad ncid: ', err)
524        else
525            nok = nok + 1
526        endif
527
528!           /* data mode test */
529        flags = IOR(NF90_CLOBBER, extra_flags)
530        err = nf90mpi_create(comm, scratch, flags, info, &
531                           ncid)
532        if (err .ne. NF90_NOERR) then
533            call errore('nf90mpi_create: ', err)
534            return
535        end if
536        err = nf90mpi_enddef(ncid)
537        if (err .ne. NF90_NOERR) &
538            call errore('nf90mpi_enddef: ', err)
539        length = 8
540        err = nf90mpi_def_dim(ncid, 'abc', length, dimid)
541        if (err .ne. NF90_ENOTINDEFINE) then
542            call errore('bad ncid: ', err)
543        else
544            nok = nok + 1
545        endif
546
547!           /* define-mode tests: unlimited dim */
548        err = nf90mpi_redef(ncid)
549        if (err .ne. NF90_NOERR) &
550            call errore('nf90mpi_redef: ', err)
551        err = nf90mpi_def_dim(ncid, dim_name(1), NF90MPI_UNLIMITED, dimid)
552        if (err .ne. NF90_NOERR)  then
553            call errore('nf90mpi_def_dim: ', err)
554        else
555            nok = nok + 1
556        endif
557        if (dimid .ne. 1)  &
558            call errori('Unexpected dimid: ', dimid)
559        err = nf90mpi_inquire(ncid, unlimitedDimId=dimid)
560        if (err .ne. NF90_NOERR)  &
561            call errore('nf90mpi_inquire: ', err)
562        if (dimid .ne. RECDIM)  &
563            call error('Unexpected recdim: ')
564        err = nf90mpi_inquire_dimension(ncid, dimid, len=length)
565        if (length .ne. 0)  &
566            call errori('Unexpected length: ', 0)
567        err = nf90mpi_def_dim(ncid, 'abc', NF90MPI_UNLIMITED, dimid)
568        if (err .ne. NF90_EUNLIMIT) then
569            call errore('2nd unlimited dimension: ', err)
570        else
571            nok = nok + 1
572        endif
573
574!           /* define-mode tests: remaining dims */
575        do 1, i = 2, NDIMS
576            err = nf90mpi_def_dim(ncid, dim_name(i-1), dim_len(i),  &
577                             dimid)
578            if (err .ne. NF90_ENAMEINUSE) then
579                call errore('duplicate name: ', err)
580            else
581                nok = nok + 1
582            endif
583            err = nf90mpi_def_dim(ncid, BAD_NAME, dim_len(i), dimid)
584            if (err .ne. NF90_EBADNAME) then
585                call errore('bad name: ', err)
586            else
587                nok = nok + 1
588            endif
589            length = NF90MPI_UNLIMITED - 1
590            err = nf90mpi_def_dim(ncid, dim_name(i), length, &
591                             dimid)
592            if (err .ne. NF90_EDIMSIZE) then
593                call errore('bad size: ', err)
594            else
595                nok = nok + 1
596            endif
597            err = nf90mpi_def_dim(ncid, dim_name(i), dim_len(i), dimid)
598            if (err .ne. NF90_NOERR)  then
599                call errore('nf90mpi_def_dim: ', err)
600            else
601                nok = nok + 1
602            endif
603            if (dimid .ne. i)  &
604                call errori('Unexpected dimid: ', 0)
6051       continue
606
607!           /* Following just to expand unlimited dim */
608        call def_vars(ncid)
609        err = nf90mpi_enddef(ncid)
610        if (err .ne. NF90_NOERR) &
611            call errore('nf90mpi_enddef: ', err)
612        call put_vars(ncid)
613
614!           /* Check all dims */
615        call check_dims(ncid)
616
617        err = nf90mpi_close(ncid)
618        if (err .ne. NF90_NOERR) &
619            call errore('nf90mpi_close: ', err)
620        err = nf90mpi_delete(scratch, info)
621        if (err .ne. NF90_NOERR) &
622            call errori('delete of scratch file failed: ', err)
623        call print_nok(nok)
624        end
625
626
627! Test nf90mpi_rename_dim
628!    try with bad netCDF handle, check error
629!    check that proper rename worked with nf90mpi_inquire_dimension
630!    try renaming to existing dimension name, check error
631!    try with bad dimension handle, check error
632        subroutine test_nf90mpi_rename_dim()
633        use pnetcdf
634        implicit        none
635#include "tests.inc"
636
637        integer ncid
638        integer err             !/* status */
639        character*(NF90_MAX_NAME) name
640        integer nok, flags
641
642        nok = 0
643
644!           /* BAD_ID test */
645        err = nf90mpi_rename_dim(BAD_ID, 1, 'abc')
646        if (err .ne. NF90_EBADID) then
647            call errore('bad ncid: ', err)
648        else
649            nok = nok + 1
650        endif
651
652!           /* main tests */
653        flags = IOR(NF90_NOCLOBBER, extra_flags)
654        err = nf90mpi_create(comm, scratch, flags, info, &
655                           ncid)
656        if (err .ne. NF90_NOERR) then
657            call errore('nf90mpi_create: ', err)
658            return
659        end if
660        call def_dims(ncid)
661        err = nf90mpi_rename_dim(ncid, BAD_DIMID, 'abc')
662        if (err .ne. NF90_EBADDIM) then
663            call errore('bad dimid: ', err)
664        else
665            nok = nok + 1
666        endif
667        err = nf90mpi_rename_dim(ncid, 3, 'abc')
668        if (err .ne. NF90_NOERR) then
669            call errore('nf90mpi_rename_dim: ', err)
670        else
671            nok = nok + 1
672        endif
673        err = nf90mpi_inquire_dimension(ncid, 3, name)
674        if (name .ne. 'abc') &
675            call errorc('Unexpected name: ', name)
676        err = nf90mpi_rename_dim(ncid, 1, 'abc')
677        if (err .ne. NF90_ENAMEINUSE) then
678            call errore('duplicate name: ', err)
679        else
680            nok = nok + 1
681        endif
682
683        err = nf90mpi_close(ncid)
684        if (err .ne. NF90_NOERR) &
685            call errore('nf90mpi_close: ', err)
686        err = nf90mpi_delete(scratch, info)
687        if (err .ne. NF90_NOERR) &
688            call errori('delete of scratch file failed: ', err)
689        call print_nok(nok)
690        end
691
692
693! Test nf90mpi_def_var
694!    try with bad netCDF handle, check error
695!    try with bad name, check error
696!    scalar tests:
697!      check that proper define worked with nf90mpi_inq_var
698!      try redefining an existing variable, check error
699!      try with bad datatype, check error
700!      try with bad number of dimensions, check error
701!      try in data mode, check error
702!    check that returned id is one more than previous id
703!    try with bad dimension ids, check error
704        subroutine test_nf90mpi_def_var()
705        use pnetcdf
706        implicit        none
707#include "tests.inc"
708
709        integer ncid
710        integer vid
711        integer err             !/* status */
712        integer i
713        integer ndims
714        integer na
715        character*(NF90_MAX_NAME) name
716        integer dimids(MAX_RANK)
717        integer datatype
718        integer nok, flags
719
720        nok = 0
721
722!           /* BAD_ID test */
723        err = nf90mpi_def_var(BAD_ID, 'abc', NF90_SHORT, varid=vid)
724        if (err .ne. NF90_EBADID) then
725            call errore('bad ncid: status = ', err)
726        else
727            nok = nok + 1
728        endif
729
730!       scalar tests
731        flags = IOR(NF90_NOCLOBBER, extra_flags)
732        err = nf90mpi_create(comm, scratch, flags, info, &
733                           ncid)
734        if (err .ne. NF90_NOERR) then
735            call errore('nf90mpi_create: ', err)
736            return
737        end if
738        err = nf90mpi_def_var(ncid, 'abc', NF90_SHORT, varid=vid)
739        if (err .ne. NF90_NOERR) then
740            call errore('nf90mpi_def_var: ', err)
741        else
742            nok = nok + 1
743        endif
744        err = nf90mpi_inquire_variable(ncid, vid, name, datatype, ndims, dimids,  &
745                         na)
746        if (err .ne. NF90_NOERR) &
747            call errore('nf90mpi_inquire_variable: ', err)
748        if (name .ne. 'abc') &
749            call errorc('Unexpected name: ', name)
750        if (datatype .ne. NF90_SHORT) &
751            call error('Unexpected datatype')
752        if (ndims .ne. 0) &
753            call error('Unexpected rank')
754        err = nf90mpi_def_var(ncid, BAD_NAME, NF90_SHORT, varid=vid)
755        if (err .ne. NF90_EBADNAME) then
756            call errore('bad name: ', err)
757        else
758            nok = nok + 1
759        endif
760        err = nf90mpi_def_var(ncid, 'abc', NF90_SHORT, varid=vid)
761        if (err .ne. NF90_ENAMEINUSE) then
762            call errore('duplicate name: ', err)
763        else
764            nok = nok + 1
765        endif
766        err = nf90mpi_def_var(ncid, 'ABC', BAD_TYPE, varid=vid)
767        if (err .ne. NF90_EBADTYPE) then
768            call errore('bad type: ', err)
769        else
770            nok = nok + 1
771        endif
772        err = nf90mpi_enddef(ncid)
773        if (err .ne. NF90_NOERR) &
774            call errore('nf90mpi_enddef: ', err)
775        err = nf90mpi_def_var(ncid, 'ABC', NF90_SHORT, varid=vid)
776        if (err .ne. NF90_ENOTINDEFINE) then
777            call errore('nf90mpi_def_var called in data mode: ', err)
778        else
779            nok = nok + 1
780        endif
781        err = nf90mpi_close(ncid)
782        if (err .ne. NF90_NOERR) &
783            call errore('nf90mpi_close: ', err)
784        err = nf90mpi_delete(scratch, info)
785        if (err .ne. NF90_NOERR) &
786            call errorc('delete of scratch file failed: ', scratch)
787
788!           /* general tests using global vars */
789        flags = IOR(NF90_CLOBBER, extra_flags)
790        err = nf90mpi_create(comm, scratch, flags, info, &
791                           ncid)
792        if (err .ne. NF90_NOERR) then
793            call errore('nf90mpi_create: ', err)
794            return
795        end if
796        call def_dims(ncid)
797        do 1, i = 1, numVars
798            err = nf90mpi_def_var(ncid, var_name(i), var_type(i),  &
799                             var_dimid(1:var_rank(i),i), vid)
800            if (err .ne. NF90_NOERR)  then
801                call errore('nf90mpi_def_var: ', err)
802            else
803                nok = nok + 1
804            endif
805            if (vid .ne. i) &
806                call error('Unexpected varid')
8071       continue
808
809!           /* try bad dim ids */
810        dimids(1) = BAD_DIMID
811        err = nf90mpi_def_var(ncid, 'abc', NF90_SHORT, dimids(1:1), vid)
812        if (err .ne. NF90_EBADDIM) then
813            call errore('bad dim ids: ', err)
814        else
815            nok = nok + 1
816        endif
817        err = nf90mpi_close(ncid)
818        if (err .ne. NF90_NOERR) &
819            call errore('nf90mpi_close: ', err)
820
821        err = nf90mpi_delete(scratch, info)
822        if (err .ne. NF90_NOERR) &
823            call errorc('delete of scratch file failed: ', scratch)
824        call print_nok(nok)
825        end
826
827
828! Test nf90mpi_rename_var
829!    try with bad netCDF handle, check error
830!    try with bad variable handle, check error
831!    try renaming to existing variable name, check error
832!    check that proper rename worked with nf90mpi_inq_varid
833!    try in data mode, check error
834        subroutine test_nf90mpi_rename_var()
835        use pnetcdf
836        implicit        none
837#include "tests.inc"
838
839        integer ncid
840        integer vid
841        integer err
842        integer i
843        character*(NF90_MAX_NAME) name
844        integer nok, flags
845
846        nok = 0
847
848        flags = IOR(NF90_NOCLOBBER, extra_flags)
849        err = nf90mpi_create(comm, scratch, flags, info, &
850                           ncid)
851        if (err .ne. NF90_NOERR) then
852            call errore('nf90mpi_create: ', err)
853            return
854        end if
855        err = nf90mpi_rename_var(ncid, BAD_VARID, 'newName')
856        if (err .ne. NF90_ENOTVAR) then
857            call errore('bad var id: ', err)
858        else
859            nok = nok + 1
860        endif
861        call def_dims(ncid)
862        call def_vars(ncid)
863
864!           /* Prefix "new_" to each name */
865        do 1, i = 1, numVars
866            err = nf90mpi_rename_var(BAD_ID, i, 'newName')
867            if (err .ne. NF90_EBADID) then
868                call errore('bad ncid: ', err)
869            else
870                nok = nok + 1
871            endif
872            err = nf90mpi_rename_var(ncid, i, var_name(numVars))
873            if (err .ne. NF90_ENAMEINUSE) then
874                call errore('duplicate name: ', err)
875            else
876                nok = nok + 1
877            endif
878            name = 'new_' // var_name(i)
879            err = nf90mpi_rename_var(ncid, i, name)
880            if (err .ne. NF90_NOERR) then
881                call errore('nf90mpi_rename_var: ', err)
882            else
883                nok = nok + 1
884            endif
885            err = nf90mpi_inq_varid(ncid, name, vid)
886            if (err .ne. NF90_NOERR) &
887                call errore('nf90mpi_inq_varid: ', err)
888            if (vid .ne. i) &
889                call error('Unexpected varid')
8901       continue
891
892!           /* Change to data mode */
893!           /* Try making names even longer. Then restore original names */
894        err = nf90mpi_enddef(ncid)
895        if (err .ne. NF90_NOERR) &
896            call errore('nf90mpi_enddef: ', err)
897        do 2, i = 1, numVars
898            name = 'even_longer_' // var_name(i)
899            err = nf90mpi_rename_var(ncid, i, name)
900            if (err .ne. NF90_ENOTINDEFINE) then
901                call errore('longer name in data mode: ', err)
902            else
903                nok = nok + 1
904            endif
905            err = nf90mpi_rename_var(ncid, i, var_name(i))
906            if (err .ne. NF90_NOERR) then
907                call errore('nf90mpi_rename_var: ', err)
908            else
909                nok = nok + 1
910            endif
911            err = nf90mpi_inq_varid(ncid, var_name(i), vid)
912            if (err .ne. NF90_NOERR) &
913                call errore('nf90mpi_inq_varid: ', err)
914            if (vid .ne. i) &
915                call error('Unexpected varid')
9162       continue
917
918        call put_vars(ncid)
919        call check_vars(ncid)
920
921        err = nf90mpi_close(ncid)
922        if (err .ne. NF90_NOERR) &
923            call errore('nf90mpi_close: ', err)
924
925        err = nf90mpi_delete(scratch, info)
926        if (err .ne. NF90_NOERR) &
927            call errorc('delete of scratch file failed: ', scratch)
928        call print_nok(nok)
929        end
930
931
932! Test nf90mpi_copy_att
933!    try with bad source or target netCDF handles, check error
934!    try with bad source or target variable handle, check error
935!    try with nonexisting attribute, check error
936!    check that NF90_GLOBAL variable for source or target works
937!    check that new attribute put works with target in define mode
938!    check that old attribute put works with target in data mode
939!    check that changing type and length of an attribute work OK
940!    try with same ncid for source and target, different variables
941!    try with same ncid for source and target, same variable
942        subroutine test_nf90mpi_copy_att()
943        use pnetcdf
944        implicit        none
945#include "tests.inc"
946        character*2 ATT_NAME
947        integer VARID, NATTS, ATT_LEN
948
949        integer ncid_in
950        integer ncid_out
951        integer vid
952        integer err
953        integer i
954        integer j
955        character*(NF90_MAX_NAME) name    !/* of att */
956        integer datatype                !/* of att */
957        integer(kind=MPI_OFFSET_KIND) length                  !/* of att */
958        character*1     value
959        integer nok, flags
960
961        nok = 0
962        err = nf90mpi_open(comm, testfile, NF90_NOWRITE, info, &
963                         ncid_in)
964        if (err .ne. NF90_NOERR) &
965            call errore('nf90mpi_open: ', err)
966        flags = IOR(NF90_NOCLOBBER, extra_flags)
967        err = nf90mpi_create(comm, scratch, flags, info, &
968                           ncid_out)
969        if (err .ne. NF90_NOERR) then
970            call errore('nf90mpi_create: ', err)
971            return
972        end if
973        call def_dims(ncid_out)
974        call def_vars(ncid_out)
975
976        do 1, i = 0, numVars
977            vid = VARID(i)
978            do 2, j = 1, NATTS(i)
979                name = ATT_NAME(j,i)
980                err = nf90mpi_copy_att(ncid_in, BAD_VARID, name, &
981                                      ncid_out, vid)
982                if (err .ne. NF90_ENOTVAR) &
983                    call errore('bad var id: ', err)
984                nok = nok + 1
985                err = nf90mpi_copy_att(ncid_in, vid, name, ncid_out,  &
986                                  BAD_VARID)
987                if (err .ne. NF90_ENOTVAR) &
988                    call errore('bad var id: ', err)
989                nok = nok + 1
990                err = nf90mpi_copy_att(BAD_ID, vid, name,  &
991                      ncid_out, vid)
992                if (err .ne. NF90_EBADID) &
993                    call errore('bad ncid: ', err)
994                nok = nok + 1
995                err = nf90mpi_copy_att(ncid_in, vid, name,  &
996                      BAD_ID, vid)
997                if (err .ne. NF90_EBADID) &
998                    call errore('bad ncid: ', err)
999                nok = nok + 1
1000                err = nf90mpi_copy_att(ncid_in, vid, 'noSuch', &
1001                                      ncid_out, vid)
1002                if (err .ne. NF90_ENOTATT) &
1003                    call errore('bad attname: ', err)
1004                nok = nok + 1
1005                err = nf90mpi_copy_att(ncid_in, vid, name,  &
1006                       ncid_out, vid)
1007                if (err .ne. NF90_NOERR) &
1008                    call errore('nf90mpi_copy_att: ', err)
1009                nok = nok + 1
1010                err = nf90mpi_copy_att(ncid_out, vid, name, &
1011                                      ncid_out, vid)
1012                if (err .ne. NF90_NOERR) &
1013                    call errore('source = target: ', err)
1014                nok = nok + 1
10152           continue
10161       continue
1017
1018        err = nf90mpi_close(ncid_in)
1019        if (err .ne. NF90_NOERR) &
1020            call errore('nf90mpi_close: ', err)
1021
1022!           /* Close scratch. Reopen & check attributes */
1023        err = nf90mpi_close(ncid_out)
1024        if (err .ne. NF90_NOERR) &
1025            call errore('nf90mpi_close: ', err)
1026        err = nf90mpi_open(comm, scratch, NF90_WRITE, info, &
1027                         ncid_out)
1028        if (err .ne. NF90_NOERR) &
1029            call errore('nf90mpi_open: ', err)
1030        call check_atts(ncid_out)
1031
1032!           change to define mode
1033!           define single char. global att. ':a' with value 'A'
1034!           This will be used as source for following copies
1035        err = nf90mpi_redef(ncid_out)
1036        if (err .ne. NF90_NOERR) &
1037            call errore('nf90mpi_redef: ', err)
1038        length = 1
1039        err = nf90mpi_put_att(ncid_out, NF90_GLOBAL, 'a', 'A')
1040        if (err .ne. NF90_NOERR) &
1041            call errore('nf90mpi_put_att: ', err)
1042
1043!           change to data mode
1044!           Use scratch as both source & dest.
1045!           try copy to existing att. change type & decrease length
1046!           rename 1st existing att of each var (if any) 'a'
1047!           if this att. exists them copy ':a' to it
1048        err = nf90mpi_enddef(ncid_out)
1049        if (err .ne. NF90_NOERR) &
1050            call errore('nf90mpi_enddef: ', err)
1051        do 3, i = 1, numVars
1052            if (NATTS(i) .gt. 0 .and. ATT_LEN(1,i) .gt. 0) then
1053                err = nf90mpi_rename_att(ncid_out, i,  &
1054                      att_name(1,i), 'a')
1055                if (err .ne. NF90_NOERR) &
1056                    call errore('nf90mpi_rename_att: ', err)
1057                err = nf90mpi_copy_att(ncid_out, NF90_GLOBAL, 'a', &
1058                                      ncid_out, i)
1059                if (err .ne. NF90_NOERR) &
1060                    call errore('nf90mpi_copy_att: ', err)
1061                nok = nok + 1
1062            end if
10633       continue
1064        err = nf90mpi_close(ncid_out)
1065        if (err .ne. NF90_NOERR) &
1066            call errore('nf90mpi_close: ', err)
1067
1068!           /* Reopen & check */
1069        err = nf90mpi_open(comm, scratch, NF90_WRITE, info, &
1070                         ncid_out)
1071        if (err .ne. NF90_NOERR) &
1072            call errore('nf90mpi_open: ', err)
1073        do 4, i = 1, numVars
1074            if (NATTS(i) .gt. 0 .and. ATT_LEN(1,i) .gt. 0) then
1075                err = nf90mpi_inquire_attribute(ncid_out, i, 'a',  &
1076                      datatype, length)
1077                if (err .ne. NF90_NOERR) &
1078                    call errore('nf90mpi_inquire_attribute: ', err)
1079                if (datatype .ne. NF90_CHAR) &
1080                    call error('Unexpected type')
1081                if (length .ne. 1) &
1082                    call error('Unexpected length')
1083                err = nf90mpi_get_att(ncid_out, i, 'a', value)
1084                if (err .ne. NF90_NOERR) &
1085                    call errore('nf90mpi_get_att: ', err)
1086                if (value .ne. 'A') &
1087                    call error('Unexpected value')
1088            end if
10894       continue
1090
1091        err = nf90mpi_close(ncid_out)
1092        if (err .ne. NF90_NOERR) &
1093            call errore('nf90mpi_close: ', err)
1094        err = nf90mpi_delete(scratch, info)
1095        if (err .ne. NF90_NOERR) &
1096            call errorc('delete of scratch file failed', scratch)
1097        call print_nok(nok)
1098        end
1099
1100
1101! Test nf90mpi_rename_att
1102!    try with bad netCDF handle, check error
1103!    try with bad variable handle, check error
1104!    try with nonexisting att name, check error
1105!    try renaming to existing att name, check error
1106!    check that proper rename worked with nf90mpi_inq_attid
1107!    try in data mode, check error
1108        subroutine test_nf90mpi_rename_att()
1109        use pnetcdf
1110        implicit        none
1111#include "tests.inc"
1112        character*2 ATT_NAME
1113        integer VARID, ATT_TYPE, NATTS, ATT_LEN
1114        double precision hash
1115        logical equal, inrange
1116
1117        integer ncid
1118        integer vid
1119        integer err, flags
1120        integer i
1121        integer j
1122        integer  k
1123        integer attnum
1124        character*(NF90_MAX_NAME) atnam
1125        character*(NF90_MAX_NAME) name
1126        character*(NF90_MAX_NAME) oldname
1127        character*(NF90_MAX_NAME) newname
1128        integer nok             !/* count of valid comparisons */
1129        integer datatype
1130        integer attyp
1131        integer(kind=MPI_OFFSET_KIND) length
1132        integer(kind=MPI_OFFSET_KIND) attlength
1133        integer(kind=MPI_OFFSET_KIND) ndx(1)
1134        character*(MAX_NELS)    text
1135        doubleprecision value(MAX_NELS)
1136        doubleprecision expect
1137
1138        nok = 0
1139
1140        flags = IOR(NF90_NOCLOBBER, extra_flags)
1141        err = nf90mpi_create(comm, scratch, flags, info, &
1142                           ncid)
1143        if (err .ne. NF90_NOERR) then
1144            call errore('nf90mpi_create: ', err)
1145            return
1146        end if
1147        err = nf90mpi_rename_att(ncid, BAD_VARID, 'abc', 'newName')
1148        if (err .ne. NF90_ENOTVAR) &
1149            call errore('bad var id: ', err)
1150        call def_dims(ncid)
1151        call def_vars(ncid)
1152        call put_atts(ncid)
1153
1154        do 1, i = 0, numVars
1155            vid = VARID(i)
1156            do 2, j = 1, NATTS(i)
1157                atnam = ATT_NAME(j,i)
1158                err = nf90mpi_rename_att(BAD_ID, vid, atnam,  &
1159                       'newName')
1160                if (err .ne. NF90_EBADID) &
1161                    call errore('bad ncid: ', err)
1162                err = nf90mpi_rename_att(ncid, vid, 'noSuch',  &
1163                        'newName')
1164                if (err .ne. NF90_ENOTATT) &
1165                    call errore('bad attname: ', err)
1166                newname = 'new_' // trim(atnam)
1167                err = nf90mpi_rename_att(ncid, vid, atnam, newname)
1168                if (err .ne. NF90_NOERR) &
1169                    call errore('nf90mpi_rename_att: ', err)
1170                err = nf90mpi_inquire_attribute(ncid, vid, newname, attnum=attnum)
1171                if (err .ne. NF90_NOERR) &
1172                    call errore('nf90mpi_inquire_attribute: ', err)
1173                if (attnum .ne. j) &
1174                    call error('Unexpected attnum')
11752           continue
11761       continue
1177
1178!           /* Close. Reopen & check */
1179        err = nf90mpi_close(ncid)
1180        if (err .ne. NF90_NOERR) &
1181            call errore('nf90mpi_close: ', err)
1182        err = nf90mpi_open(comm, scratch, NF90_WRITE, info,  &
1183            ncid)
1184        if (err .ne. NF90_NOERR) &
1185            call errore('nf90mpi_open: ', err)
1186
1187        do 3, i = 0, numVars
1188            vid = VARID(i)
1189            do 4, j = 1, NATTS(i)
1190                atnam = ATT_NAME(j,i)
1191                attyp = ATT_TYPE(j,i)
1192                attlength = ATT_LEN(j,i)
1193                newname = 'new_' // trim(atnam)
1194                err = nf90mpi_inq_attname(ncid, vid, j, name)
1195                if (err .ne. NF90_NOERR) &
1196                    call errore('nf90mpi_inq_attname: ', err)
1197                if (name .ne. newname) &
1198                    call error('nf90mpi_inq_attname: unexpected name')
1199                err = nf90mpi_inquire_attribute(ncid, vid, name,  &
1200                    datatype, length)
1201                if (err .ne. NF90_NOERR) &
1202                    call errore('nf90mpi_inquire_attribute: ', err)
1203                if (datatype .ne. attyp) &
1204                    call error('nf90mpi_inquire_attribute: unexpected type')
1205                if (length .ne. attlength) &
1206                    call error('nf90mpi_inquire_attribute: unexpected length')
1207                if (datatype .eq. NF90_CHAR) then
1208                    err = nf90mpi_get_att(ncid, vid, name, text)
1209                    if (err .ne. NF90_NOERR) &
1210                        call errore('nf90mpi_get_att: ', err)
1211                    do 5, k = 1, INT(attlength)
1212                        ndx(1) = k
1213                        expect = hash(datatype, -1, ndx)
1214                        if (ichar(text(k:k)) .ne. expect) then
1215                            call error( &
1216                                'nf90mpi_get_att: unexpected value')
1217                        else
1218                            nok = nok + 1
1219                        end if
12205                   continue
1221                else
1222                    err = nf90mpi_get_att(ncid, vid, name,  &
1223                           value)
1224                    if (err .ne. NF90_NOERR) &
1225                        call errore('nf90mpi_get_att: ', err)
1226                    do 6, k = 1, INT(attlength)
1227                        ndx(1) = k
1228                        expect = hash(datatype, -1, ndx)
1229                        if (inRange(expect, datatype)) then
1230                            if (.not. equal(value(k),expect,datatype, &
1231                                            NF90_DOUBLE)) then
1232                                call error( &
1233                              'nf90mpi_get_att: unexpected value')
1234                            else
1235                                nok = nok + 1
1236                            end if
1237                        end if
12386                   continue
1239                end if
12404           continue
12413       continue
1242        call print_nok(nok)
1243
1244!           /* Now in data mode */
1245!           /* Try making names even longer. Then restore original names */
1246
1247        do 7, i = 0, numVars
1248            vid = VARID(i)
1249            do 8, j = 1, NATTS(i)
1250                atnam = ATT_NAME(j,i)
1251                oldname = 'new_' // trim(atnam)
1252                newname = 'even_longer_' // trim(atnam)
1253                err = nf90mpi_rename_att(ncid, vid, oldname, newname)
1254                if (err .ne. NF90_ENOTINDEFINE) &
1255                    call errore('longer name in data mode: ', err)
1256                err = nf90mpi_rename_att(ncid, vid, oldname, atnam)
1257                if (err .ne. NF90_NOERR) &
1258                    call errore('nf90mpi_rename_att: ', err)
1259                err = nf90mpi_inquire_attribute(ncid, vid, atnam, attnum=attnum)
1260                if (err .ne. NF90_NOERR) &
1261                    call errore('nf90mpi_inquire_attribute: ', err)
1262                if (attnum .ne. j) &
1263                    call error('Unexpected attnum')
12648           continue
12657       continue
1266
1267        err = nf90mpi_close(ncid)
1268        if (err .ne. NF90_NOERR) &
1269            call errore('nf90mpi_close: ', err)
1270
1271        err = nf90mpi_delete(scratch, info)
1272        if (err .ne. NF90_NOERR) &
1273            call errori('delete of scratch file failed: ', err)
1274        end
1275
1276
1277! Test nf90mpi_del_att
1278!    try with bad netCDF handle, check error
1279!    try with bad variable handle, check error
1280!    try with nonexisting att name, check error
1281!    check that proper delete worked using:
1282!      nf90mpi_inq_attid, nf90mpi_inq_natts, nf90mpi_inq_varnatts
1283        subroutine test_nf90mpi_del_att()
1284        use pnetcdf
1285        implicit        none
1286#include "tests.inc"
1287        character*2 ATT_NAME
1288        integer VARID, NATTS
1289
1290        integer ncid
1291        integer err, flags
1292        integer i
1293        integer j
1294        integer attnum
1295        integer na
1296        integer numatts
1297        integer vid
1298        character*(NF90_MAX_NAME)  name           !/* of att */
1299        integer nok             !/* count of valid comparisons */
1300
1301        nok = 0
1302
1303        flags = IOR(NF90_NOCLOBBER, extra_flags)
1304        err = nf90mpi_create(comm, scratch, flags, info, &
1305                           ncid)
1306        if (err .ne. NF90_NOERR) then
1307            call errore('nf90mpi_create: ', err)
1308            return
1309        end if
1310        err = nf90mpi_del_att(ncid, BAD_VARID, 'abc')
1311        if (err .ne. NF90_ENOTVAR) then
1312            call errore('bad var id: ', err)
1313        else
1314            nok = nok + 1
1315        endif
1316        call def_dims(ncid)
1317        call def_vars(ncid)
1318        call put_atts(ncid)
1319
1320        do 1, i = 0, numVars
1321            vid = VARID(i)
1322            numatts = NATTS(i)
1323            do 2, j = 1, numatts
1324                name = ATT_NAME(j,i)
1325                err = nf90mpi_del_att(BAD_ID, vid, name)
1326                if (err .ne. NF90_EBADID) then
1327                    call errore('bad ncid: ', err)
1328                else
1329                    nok = nok + 1
1330                endif
1331                err = nf90mpi_del_att(ncid, vid, 'noSuch')
1332                if (err .ne. NF90_ENOTATT) then
1333                    call errore('bad attname: ', err)
1334                else
1335                    nok = nok + 1
1336                endif
1337                err = nf90mpi_del_att(ncid, vid, name)
1338                if (err .ne. NF90_NOERR) then
1339                    call errore('nf90mpi_del_att: ', err)
1340                else
1341                    nok = nok + 1
1342                endif
1343                err = nf90mpi_inquire_attribute(ncid, vid, name, attnum=attnum)
1344                if (err .ne. NF90_ENOTATT) &
1345                    call errore('bad attname: ', err)
1346                if (i .lt. 1) then
1347                    err = nf90mpi_inquire(ncid, nAttributes=na)
1348                    if (err .ne. NF90_NOERR) &
1349                        call errore('nf90mpi_inquire: ', err)
1350                    if (na .ne. numatts-j) then
1351                        call errori('natts: expected: ', numatts-j)
1352                        call errori('natts: got:      ', na)
1353                    endif
1354                else
1355                    err = nf90mpi_inquire_variable(ncid, vid, nAtts=na)
1356                    if (err .ne. NF90_NOERR) &
1357                        call errore('nf90mpi_inquire_variable: ', err)
1358                    if (na .ne. numatts-j) then
1359                        call errori('natts: expected: ', numatts-j)
1360                        call errori('natts: got:      ', na)
1361                    endif
1362                endif
13632           continue
13641       continue
1365
1366!           /* Close. Reopen & check no attributes left */
1367        err = nf90mpi_close(ncid)
1368        if (err .ne. NF90_NOERR) &
1369            call errore('nf90mpi_close: ', err)
1370        err = nf90mpi_open(comm, scratch, NF90_WRITE,  &
1371           info, ncid)
1372        if (err .ne. NF90_NOERR) &
1373            call errore('nf90mpi_open: ', err)
1374        err = nf90mpi_inquire(ncid, nAttributes=na)
1375        if (err .ne. NF90_NOERR) &
1376            call errore('nf90mpi_inquire: ', err)
1377        if (na .ne. 0) &
1378            call errori('natts: expected 0, got ', na)
1379        do 3, i = 0, numVars
1380            vid = VARID(i)
1381            err = nf90mpi_inquire(ncid, nAttributes=na)
1382            if (err .ne. NF90_NOERR) &
1383                call errore('nf90mpi_inquire: ', err)
1384            if (na .ne. 0) &
1385                call errori('natts: expected 0, got ', na)
13863       continue
1387
1388!           /* restore attributes. change to data mode. try to delete */
1389        err = nf90mpi_redef(ncid)
1390        if (err .ne. NF90_NOERR) &
1391            call errore('nf90mpi_redef: ', err)
1392        call put_atts(ncid)
1393        err = nf90mpi_enddef(ncid)
1394        if (err .ne. NF90_NOERR) &
1395            call errore('nf90mpi_enddef: ', err)
1396
1397        do 4, i = 0, numVars
1398            vid = VARID(i)
1399            numatts = NATTS(i)
1400            do 5, j = 1, numatts
1401                name = ATT_NAME(j,i)
1402                err = nf90mpi_del_att(ncid, vid, name)
1403                if (err .ne. NF90_ENOTINDEFINE) then
1404                    call errore('in data mode: ', err)
1405                else
1406                    nok = nok + 1
1407                endif
14085           continue
14094       continue
1410
1411        err = nf90mpi_close(ncid)
1412        if (err .ne. NF90_NOERR) &
1413            call errore('nf90mpi_close: ', err)
1414        err = nf90mpi_delete(scratch, info)
1415        if (err .ne. NF90_NOERR) &
1416            call errori('delete of scratch file failed: ', err)
1417        call print_nok(nok)
1418        end
1419
1420! parallel-netcdf doesn't implement set_fill, so i have not
1421!    parallel-netcdfified this subroutine
1422! Test nf90mpi_set_fill
1423!    try with bad netCDF handle, check error
1424!    try in read-only mode, check error
1425!    try with bad new_fillmode, check error
1426!    try in data mode, check error
1427!    check that proper set to NF90_FILL works for record & non-record variables
1428!    (note that it is not possible to test NF90_NOFILL mode!)
1429!    close file & create again for test using attribute _FillValue
1430        subroutine test_nf90mpi_set_fill()
1431        use pnetcdf
1432        implicit none
1433#include "tests.inc"
1434        ! character*2 ATT_NAME
1435        ! integer VARID, ATT_TYPE, NATTS
1436
1437        integer ncid
1438        integer vid
1439        integer err, flags
1440        integer i
1441        integer j
1442        integer old_fillmode
1443        character*1 text
1444        doubleprecision value
1445        doubleprecision fill, fills(1)
1446        integer(kind=MPI_OFFSET_KIND) index(MAX_RANK)
1447        integer(kind=MPI_OFFSET_KIND) length
1448        integer index2indexes
1449        integer nok             !/* count of valid comparisons */
1450
1451        value = 0
1452        nok = 0
1453
1454!           /* bad ncid */
1455        err = nf90mpi_set_fill(BAD_ID, NF90_NOFILL, old_fillmode)
1456        if (err .ne. NF90_EBADID) &
1457            call errore('bad ncid: ', err)
1458
1459!           /* try in read-only mode */
1460        err = nf90mpi_open(comm, testfile, NF90_NOWRITE, info, &
1461                         ncid)
1462        if (err .ne. NF90_NOERR) &
1463            call errore('nf90mpi_open: ', err)
1464        err = nf90mpi_set_fill(ncid, NF90_NOFILL, old_fillmode)
1465        if (err .ne. NF90_EPERM) &
1466            call errore('read-only: ', err)
1467        err = nf90mpi_close(ncid)
1468        if (err .ne. NF90_NOERR) &
1469            call errore('nf90mpi_close: ', err)
1470
1471!           /* create scratch */
1472        flags = IOR(NF90_NOCLOBBER, extra_flags)
1473        err = nf90mpi_create(comm, scratch, flags, info, &
1474                           ncid)
1475        if (err .ne. NF90_NOERR) then
1476            call errore('nf90mpi_create: ', err)
1477            return
1478        end if
1479
1480!           /* BAD_FILLMODE */
1481        err = nf90mpi_set_fill(ncid, BAD_FILLMODE, old_fillmode)
1482        if (err .ne. NF90_EINVAL) &
1483            call errore('bad fillmode: ', err)
1484
1485!           /* proper calls */
1486        err = nf90mpi_set_fill(ncid, NF90_FILL, old_fillmode)
1487        if (err .ne. NF90_NOERR) &
1488            call errore('nf90mpi_set_fill: ', err)
1489        if (old_fillmode .ne. NF90_NOFILL) &
1490            call errori('Unexpected old fill mode: ', old_fillmode)
1491        err = nf90mpi_set_fill(ncid, NF90_NOFILL, old_fillmode)
1492        if (err .ne. NF90_NOERR) &
1493            call errore('nf90mpi_set_fill: ', err)
1494        if (old_fillmode .ne. NF90_FILL) &
1495            call errori('Unexpected old fill mode: ', old_fillmode)
1496        err = nf90mpi_set_fill(ncid, NF90_FILL, old_fillmode)
1497        if (err .ne. NF90_NOERR) &
1498            call errore('nf90mpi_set_fill: ', err)
1499
1500!           /* define dims & vars */
1501        call def_dims(ncid)
1502        call def_vars(ncid)
1503
1504!           /* Change to data mode. Set fillmode again */
1505        err = nf90mpi_enddef(ncid)
1506        if (err .ne. NF90_NOERR) &
1507            call errore('nf90mpi_enddef: ', err)
1508        err = nf90mpi_set_fill(ncid, NF90_FILL, old_fillmode)
1509        if (err .ne. NF90_ENOTINDEFINE) &
1510            call errore('nf90mpi_set_fill: ', err)
1511
1512!       /* Write record number NRECS to force writing of preceding records */
1513!       /* Assumes variable cr is char vector with UNLIMITED dimension */
1514        err = nf90mpi_inq_varid(ncid, 'cr', vid)
1515        if (err .ne. NF90_NOERR) &
1516            call errore('nf90mpi_inq_varid: ', err)
1517        index(1) = NRECS
1518        text = char(NF90_FILL_CHAR)
1519        err = nf90mpi_put_var_all(ncid, vid, text, index)
1520        if (err .ne. NF90_NOERR) &
1521            call errore('nf90mpi_put_var_all: ', err)
1522
1523!           /* get all variables & check all values equal default fill */
1524        do 1, i = 1, numVars
1525            if (var_dimid(var_rank(i),i) .eq. RECDIM) cycle ! skip record variables
1526
1527            if (var_type(i) .eq. NF90_CHAR) then
1528                fill = NF90_FILL_CHAR
1529            else if (var_type(i) .eq. NF90_BYTE) then
1530                fill = NF90_FILL_BYTE
1531            else if (var_type(i) .eq. NF90_SHORT) then
1532                fill = NF90_FILL_SHORT
1533            else if (var_type(i) .eq. NF90_INT) then
1534                fill = NF90_FILL_INT
1535            else if (var_type(i) .eq. NF90_FLOAT) then
1536                fill = NF90_FILL_FLOAT
1537            else if (var_type(i) .eq. NF90_DOUBLE) then
1538                fill = NF90_FILL_DOUBLE
1539            else if (var_type(i) .eq. NF90_UBYTE) then
1540                fill = NF90_FILL_UBYTE
1541            else if (var_type(i) .eq. NF90_USHORT) then
1542                fill = NF90_FILL_USHORT
1543            else if (var_type(i) .eq. NF90_UINT) then
1544                fill = NF90_FILL_UINT
1545            else if (var_type(i) .eq. NF90_INT64) then
1546                fill = NF90_FILL_INT64
1547            else if (var_type(i) .eq. NF90_UINT64) then
1548                fill = NF90_FILL_UINT64
1549            else
1550                stop 'test_nf90mpi_set_fill(): impossible var_type(i)'
1551            end if
1552
1553            do 2, j = 1, var_nels(i)
1554                err = index2indexes(j, var_rank(i), var_shape(1,i),  &
1555                                    index)
1556                if (err .ne. NF90_NOERR) &
1557                    call error('error in index2indexes()')
1558                if (var_type(i) .eq. NF90_CHAR) then
1559                    err = nf90mpi_get_var_all(ncid, i, text, index)
1560                    if (err .ne. NF90_NOERR) &
1561                        call errore('nf90mpi_get_var_all failed: ',err)
1562                    value = ichar(text)
1563                else
1564                    err = nf90mpi_get_var_all(ncid, i, value, index)
1565                    if (err .ne. NF90_NOERR) &
1566                        call errore &
1567                                 ('nf90mpi_get_var_all failed: ',err)
1568                end if
1569                if (value .ne. fill .and.  &
1570                    abs((fill - value)/fill) .gt. 1.0e-9) then
1571                    call errord('Unexpected fill value: ', value)
1572                else
1573                    nok = nok + 1
1574                end if
15752           continue
15761       continue
1577
1578!       /* close scratch & create again for test using attribute _FillValue */
1579        err = nf90mpi_close(ncid)
1580        if (err .ne. NF90_NOERR) &
1581            call errore('nf90mpi_close: ', err)
1582        flags = IOR(NF90_CLOBBER, extra_flags)
1583        err = nf90mpi_create(comm, scratch, flags, info, &
1584                           ncid)
1585        if (err .ne. NF90_NOERR) then
1586            call errore('nf90mpi_create: ', err)
1587            return
1588        end if
1589        call def_dims(ncid)
1590        call def_vars(ncid)
1591
1592!           /* set _FillValue = 42 for all vars */
1593        fill = 42
1594        text = char(int(fill))
1595        length = 1
1596        do 3, i = 1, numVars
1597            if (var_dimid(var_rank(i),i) .eq. RECDIM) cycle ! skip record variables
1598
1599            if (var_type(i) .eq. NF90_CHAR) then
1600                err = nf90mpi_put_att(ncid, i, '_FillValue', &
1601                   text)
1602                if (err .ne. NF90_NOERR) &
1603                    call errore('nf90mpi_put_att: ', err)
1604            else
1605                ! cannot use overloading, as fill's type must match with variable's
1606                ! err = nf90mpi_put_att(ncid, i, '_FillValue', fill)
1607                fills(1) = fill
1608                err = nfmpi_put_att_double(ncid, i, '_FillValue', &
1609                                           var_type(i),length,fills(1))
1610                if (err .ne. NF90_NOERR) &
1611                    call errore('nf90mpi_put_att: ', err)
1612            end if
16133       continue
1614
1615!           /* data mode. write records */
1616        err = nf90mpi_enddef(ncid)
1617        if (err .ne. NF90_NOERR) &
1618            call errore('nf90mpi_enddef: ', err)
1619        index(1) = NRECS
1620        err = nf90mpi_put_var_all(ncid, vid, text, index)
1621        if (err .ne. NF90_NOERR) &
1622            call errore('nf90mpi_put_var_all: ', err)
1623
1624!           /* get all variables & check all values equal 42 */
1625        do 4, i = 1, numVars
1626            if (var_dimid(var_rank(i),i) .eq. RECDIM) cycle ! skip record variables
1627
1628            do 5, j = 1, var_nels(i)
1629                err = index2indexes(j, var_rank(i), var_shape(1,i),  &
1630                                    index)
1631                if (err .ne. NF90_NOERR) &
1632                    call error('error in index2indexes')
1633                if (var_type(i) .eq. NF90_CHAR) then
1634                    err = nf90mpi_get_var_all(ncid, i, text, index)
1635                    if (err .ne. NF90_NOERR) &
1636                        call errore('nf90mpi_get_var_all failed: ',err)
1637                    value = ichar(text)
1638                else
1639                    err = nf90mpi_get_var_all(ncid, i, value, index)
1640                    if (err .ne. NF90_NOERR) &
1641                        call errore &
1642                                ('nf90mpi_get_var_all failed: ', err)
1643                end if
1644                if (value .ne. fill) then
1645                    call errord(' Value expected: ', fill)
1646                    call errord(' Value read:     ', value)
1647                else
1648                    nok = nok + 1
1649                end if
16505           continue
16514       continue
1652
1653        err = nf90mpi_close(ncid)
1654        if (err .ne. NF90_NOERR) &
1655            call errore('nf90mpi_close: ', err)
1656        err = nf90mpi_delete(scratch, info)
1657        if (err .ne. NF90_NOERR) &
1658            call errori('delete of scratch file failed: ', err)
1659
1660        call print_nok(nok)
1661        end
1662
1663#if 0
1664! * Test nc_set_default_format
1665! *    try with bad default format
1666! *    try with NULL old_formatp
1667! *    try in data mode, check error
1668! *    check that proper set to NC_FILL works for record & non-record variables
1669! *    (note that it is not possible to test NC_NOFILL mode!)
1670! *    close file & create again for test using attribute _FillValue
1671      subroutine test_nf90mpi_set_default_format()
1672        use pnetcdf
1673      implicit none
1674#include "tests.inc"
1675
1676      integer ncid
1677      integer err, flags
1678      integer i
1679      integer version
1680      integer old_format
1681      integer nf90mpi_get_file_version
1682
1683!     /* bad format */
1684      err = nf90mpi_set_default_format(3, old_format)
1685      IF (err .ne. NF90_EINVAL) &
1686           call errore("bad default format: status = %d", err)
1687
1688!    /* Cycle through available formats. */
1689      do 1 i=1, 2
1690         err = nf90mpi_set_default_format(i, old_format)
1691         if (err .ne. NF90_NOERR)  &
1692               call errore("setting classic format: status = %d", err)
1693         flags = IOR(NF90_CLOBBER, extra_flags)
1694         err = nf90mpi_create(comm, scratch, flags,  &
1695                      info, ncid)
1696         if (err .ne. NF90_NOERR)  &
1697              call errore("bad nf90mpi_create: status = %d", err)
1698         err = nf90mpi_put_att_text(ncid, NF90_GLOBAL, "testatt",  &
1699              "blah")
1700         if (err .ne. NF90_NOERR) call errore("bad put_att: status = %d", err)
1701         err = nf90mpi_close(ncid)
1702         if (err .ne. NF90_NOERR) call errore("bad close: status = %d", err)
1703         err = nf90mpi_get_file_version(scratch, version)
1704         if (err .ne. NF90_NOERR) call errore("bad file version = %d", err)
1705         if (version .ne. i) &
1706              call errore("bad file version = %d", err)
1707 1    continue
1708
1709!    /* Remove the left-over file. */
1710      err = nf90mpi_delete(scratch)
1711      if (err .ne. NF90_NOERR) call errore("remove failed", err)
1712      end
1713
1714#endif
1715
1716!     This function looks in a file for the netCDF magic number.
1717      integer function nf90mpi_get_file_version(path, version)
1718        use pnetcdf
1719      implicit none
1720#include "tests.inc"
1721
1722      character*(*) path
1723      integer version
1724      character magic*4
1725      integer ver
1726      integer f
1727      parameter (f = 10)
1728
1729      open(f, file=path, status='OLD', form='UNFORMATTED', &
1730           access='DIRECT', recl=4)
1731
1732!     Assume this is not a netcdf file.
1733      nf90mpi_get_file_version = NF90_ENOTNC
1734      version = 0
1735
1736!     Read the magic number, the first 4 bytes of the file.
1737      read(f, rec=1, err = 1) magic
1738
1739!     If the first three characters are not "CDF" we're done.
1740      if (index(magic, 'CDF') .eq. 1) then
1741         ver = ichar(magic(4:4))
1742         if (ver .eq. 1) then
1743            version = 1
1744            nf90mpi_get_file_version = NF90_NOERR
1745         elseif (ver .eq. 2) then
1746            version = 2
1747            nf90mpi_get_file_version = NF90_NOERR
1748         endif
1749      endif
1750
1751 1    close(f)
1752      return
1753      end
1754
1755
1756