ROMS
Loading...
Searching...
No Matches
def_var.F
Go to the documentation of this file.
1#include "cppdefs.h"
3!
4!git $Id$
5!================================================== Hernan G. Arango ===
6! Copyright (c) 2002-2025 The ROMS Group !
7! Licensed under a MIT/X style license !
8! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine defines the requested NetCDF variable. !
12! !
13! On Input: !
14! !
15! ng Nested grid number (integer) !
16! model Calling model identifier (integer) !
17! ncid NetCDF file ID (integer) !
18#if defined PIO_LIB && defined DISTRIBUTE
19!or pioFile PIO file descriptor structure, TYPE(file_desc_t) !
20! pioFile%fh file handler !
21! pioFile%iosystem IO system descriptor (struct) !
22#endif
23! Vtype NetCDF variable type (integer) !
24! nVdim Number of variable dimensions (integer; 0=scalar) !
25! Vdim Dimensions IDs for this variable (integer vector) !
26! Aval Attribute values (real vector): !
27! Aval(1) => Add offset value !
28! Aval(2) => Valid minimum value !
29! Aval(3) => Valid maximum value !
30! Aval(4) => Missing value !
31! Aval(5) => C-grid variable type !
32! Aval(6) => Fill value !
33! Vinfo Variable information (character array): !
34! Vinfo( 1) => Variable name !
35! Vinfo( 2) => Variable "longname" attribute !
36! Vinfo( 3) => Variable "units" attribute !
37! Vinfo( 4) => Variable "calendar" attribute !
38! Vinfo( 5) => Variable "valid_min" attribute !
39! Vinfo( 6) => Variable "valid_max" attribute !
40! Vinfo( 7) => Variable "option_T" attribute !
41! Vinfo( 8) => Variable "option_F" attribute !
42! Vinfo( 9) => Variable "option_0" attribute !
43! Vinfo(10) => Variable "option_1" attribute !
44! Vinfo(11) => Variable "negative_value" attribute !
45! Vinfo(12) => Variable "positive_value" attribute !
46! Vinfo(13) => Variable "cycle" attribute !
47! Vinfo(14) => Variable "field" attribute !
48! Vinfo(15) => Variable "positions" attribute !
49! Vinfo(16) => Variable "time" attribute !
50! Vinfo(17) => Variable "missing_value" attribute !
51! Vinfo(18) => Variable "add_offset" attribute !
52! Vinfo(19) => Variable "size_class" attribute !
53! Vinfo(20) => Variable "water_points" attribute !
54! Vinfo(21) => Variable "standard_name" attribute !
55! Vinfo(22) => Variable "coordinates" attribute !
56! Vinfo(23) => Variable "formula_terms" attribute !
57! Vinfo(24) => Variable "_FillValue" attribute !
58! Vinfo(25) => Variable "positive" attribute !
59! ncname NetCDF file name !
60! SetFillVal Logical switch to set fill value in land areas !
61! (optional) !
62! SetParAccess Logical switch to set parallel I/O access flag to !
63! either collective or independent (optional) !
64! !
65! On Output: !
66! !
67! def_var Error flag (integer) !
68! Vid NetCDF variable ID (integer) !
69#if defined PIO_LIB && defined DISTRIBUTE
70!or pioVar PIO variable descriptor structure, TYPE(Var_desc_t) !
71! pioVar%varID Variable ID !
72! pioVar%ncid File ID !
73#endif
74! !
75! Notice that arrays "Aval" and "Vinfo" is destroyed on output to !
76! facilitate the definition of the next variable. !
77! !
78!=======================================================================
79!
80 USE mod_param
81 USE mod_parallel
82 USE mod_iounits
83 USE mod_ncparam
84 USE mod_scalars
85!
86 implicit none
87!
88 INTERFACE def_var
89 MODULE PROCEDURE def_var_nf90
90#if defined PIO_LIB && defined DISTRIBUTE
91 MODULE PROCEDURE def_var_pio
92#endif
93 END INTERFACE def_var
94!
95 CONTAINS
96!
97!***********************************************************************
98 FUNCTION def_var_nf90 (ng, model, ncid, Vid, &
99 & Vtype, nVdim, Vdim, &
100 & Aval, Vinfo, ncname, &
101 & SetFillVal, SetParAccess) RESULT (status)
102!***********************************************************************
103!
104 USE mod_netcdf
105!
106#if !defined PARALLEL_IO && defined DISTRIBUTE
107 USE distribute_mod, ONLY : mp_bcasti
108#endif
109 USE strings_mod, ONLY : founderror
110!
111! Imported variable declarations.
112!
113 logical, intent(in), optional :: setfillval
114 logical, intent(in), optional :: setparaccess
115!
116 integer, intent (in) :: ng, model, ncid, vtype, nvdim
117
118 integer, dimension(:), intent (in) :: vdim
119
120 integer, intent (out) :: vid
121!
122 real(r8), dimension(:), intent(inout) :: aval
123!
124 character (len=*), intent(in) :: ncname
125 character (len=*), intent(inout) :: vinfo(25)
126!
127! Local variable declarations.
128!
129#ifdef MASKING
130 logical :: landfill
131#endif
132#if defined PARALLEL_IO && defined DISTRIBUTE
133 logical :: ltiled
134!
135#endif
136 integer :: i, j, latt
137 integer :: status
138!
139 character (len= 5) location
140 character (len=160) text
141
142 character (len=*), parameter :: myfile = &
143 & __FILE__//", def_var_nf90"
144!
145!-----------------------------------------------------------------------
146! Define requested variable and its attributes.
147!-----------------------------------------------------------------------
148!
149 status=nf90_noerr
150!
151 IF (outthread) THEN
152!
153! Define variable.
154!
155 IF (exit_flag.eq.noerror) THEN
156 IF (len_trim(vinfo(1)).gt.0) THEN
157 IF ((nvdim.eq.1).and.(vdim(1).eq.0)) THEN
158 status=nf90_def_var(ncid, trim(vinfo(1)), vtype, &
159 & varid = vid)
160 ELSE
161 status=nf90_def_var(ncid, trim(vinfo(1)), vtype, &
162 & vdim(1:nvdim), vid)
163 END IF
164 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
165 IF (master) WRITE (stdout,10) ng, trim(vinfo(1)), &
166 & trim(ncname)
167 exit_flag=3
168 ioerror=status
169 END IF
170 END IF
171 END IF
172
173#if !defined PARALLEL_IO && (defined HDF5 && defined DEFLATE)
174!
175! Define deflate (file compresion) parameters. Notice that deflation
176! cannot be used in parallel I/O for writing data. This is because
177! the compression makes it impossible for the HDF5 library to exactly
178! map the data to the disk location. However, deflated data can be
179! read with parallel I/O
180!
181 IF (exit_flag.eq.noerror) THEN
182 IF (len_trim(vinfo(1)).gt.0) THEN
183 IF ((nvdim.gt.1).and.(vdim(1).ne.0)) THEN
184 status=nf90_def_var_deflate(ncid, vid, shuffle, deflate, &
186 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
187 IF (master) WRITE (stdout,20) trim(vinfo(1)), &
188 & trim(ncname)
189 exit_flag=3
190 ioerror=status
191 END IF
192 END IF
193 END IF
194 END IF
195#endif
196!
197! Define special attributes for SGRID conventions variable "grid".
198!
199 IF (exit_flag.eq.noerror) THEN
200 IF (trim(vinfo(1)).eq.'grid') THEN
201 status=nf90_put_att(ncid, vid, 'cf_role', &
202 & 'grid_topology')
203 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
204 IF (master) WRITE (stdout,30) 'cf_role', &
205 & trim(vinfo(1)), &
206 & trim(ncname)
207 exit_flag=3
208 ioerror=status
209 END IF
210!
211 status=nf90_put_att(ncid, vid, 'topology_dimension', &
212 & (/2/))
213 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
214 IF (master) WRITE (stdout,30) 'topology_dimension', &
215 & trim(vinfo(1)), &
216 & trim(ncname)
217 exit_flag=3
218 ioerror=status
219 END IF
220!
221 status=nf90_put_att(ncid, vid, 'node_dimensions', &
222 & 'xi_psi eta_psi')
223 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
224 IF (master) WRITE (stdout,30) 'node_dimensions', &
225 & trim(vinfo(1)), &
226 & trim(ncname)
227 exit_flag=3
228 ioerror=status
229 END IF
230!
231 text='xi_rho: xi_psi (padding: both) '// &
232 & 'eta_rho: eta_psi (padding: both)'
233 status=nf90_put_att(ncid, vid, 'face_dimensions', &
234 & trim(text))
235 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
236 IF (master) WRITE (stdout,30) 'face_dimensions', &
237 & trim(vinfo(1)), &
238 & trim(ncname)
239 exit_flag=3
240 ioerror=status
241 END IF
242!
243 text='xi_u: xi_psi eta_u: eta_psi (padding: both)'
244 status=nf90_put_att(ncid, vid, 'edge1_dimensions', &
245 & trim(text))
246 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
247 IF (master) WRITE (stdout,30) 'edge1_dimensions', &
248 & trim(vinfo(1)), &
249 & trim(ncname)
250 exit_flag=3
251 ioerror=status
252 END IF
253!
254 text='xi_v: xi_psi (padding: both) eta_v: eta_psi'
255 status=nf90_put_att(ncid, vid, 'edge2_dimensions', &
256 & trim(text))
257 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
258 IF (master) WRITE (stdout,30) 'edge2_dimensions', &
259 & trim(vinfo(1)), &
260 & trim(ncname)
261 exit_flag=3
262 ioerror=status
263 END IF
264!
265 IF (spherical) THEN
266 status=nf90_put_att(ncid, vid, 'node_coordinates', &
267 & 'lon_psi lat_psi')
268 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
269 IF (master) WRITE (stdout,30) 'node_coordinates', &
270 & trim(vinfo(1)), &
271 & trim(ncname)
272 exit_flag=3
273 ioerror=status
274 END IF
275!
276 status=nf90_put_att(ncid, vid, 'face_coordinates', &
277 & 'lon_rho lat_rho')
278 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
279 IF (master) WRITE (stdout,30) 'face_coordinates', &
280 & trim(vinfo(1)), &
281 & trim(ncname)
282 exit_flag=3
283 ioerror=status
284 END IF
285!
286 status=nf90_put_att(ncid, vid, 'edge1_coordinates', &
287 & 'lon_u lat_u')
288 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
289 IF (master) WRITE (stdout,30) 'edge1_coordinates', &
290 & trim(vinfo(1)), &
291 & trim(ncname)
292 exit_flag=3
293 ioerror=status
294 END IF
295!
296 status=nf90_put_att(ncid, vid, 'edge2_coordinates', &
297 & 'lon_v lat_v')
298 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
299 IF (master) WRITE (stdout,30) 'edge2_coordinates', &
300 & trim(vinfo(1)), &
301 & trim(ncname)
302 exit_flag=3
303 ioerror=status
304 END IF
305 ELSE
306 status=nf90_put_att(ncid, vid, 'node_coordinates', &
307 & 'x_psi y_psi')
308 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
309 IF (master) WRITE (stdout,30) 'node_coordinates', &
310 & trim(vinfo(1)), &
311 & trim(ncname)
312 exit_flag=3
313 ioerror=status
314 END IF
315!
316 status=nf90_put_att(ncid, vid, 'face_coordinates', &
317 & 'x_rho y_rho')
318 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
319 IF (master) WRITE (stdout,30) 'face_coordinates', &
320 & trim(vinfo(1)), &
321 & trim(ncname)
322 exit_flag=3
323 ioerror=status
324 END IF
325!
326 status=nf90_put_att(ncid, vid, 'edge1_coordinates', &
327 & 'x_u y_u')
328 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
329 IF (master) WRITE (stdout,30) 'edge1_coordinates', &
330 & trim(vinfo(1)), &
331 & trim(ncname)
332 exit_flag=3
333 ioerror=status
334 END IF
335!
336 status=nf90_put_att(ncid, vid, 'edge2_coordinates', &
337 & 'x_v y_v')
338 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
339 IF (master) WRITE (stdout,30) 'edge2_coordinates', &
340 & trim(vinfo(1)), &
341 & trim(ncname)
342 exit_flag=3
343 ioerror=status
344 END IF
345 END IF
346#ifdef SOLVE3D
347!
348 status=nf90_put_att(ncid, vid, 'vertical_dimensions', &
349 & 's_rho: s_w (padding: none)')
350 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
351 IF (master) WRITE (stdout,30) 'vertical_dimensions', &
352 & trim(vinfo(1)), &
353 & trim(ncname)
354 exit_flag=3
355 ioerror=status
356 END IF
357#endif
358 END IF
359 END IF
360!
361! If applicable, define "standard_name" attribute.
362!
363 IF (exit_flag.eq.noerror) THEN
364 latt=len_trim(vinfo(21))
365 IF (latt.gt.0.and.(vinfo(21)(1:6).ne.'nulval')) THEN
366 status=nf90_put_att(ncid, vid, 'standard_name', &
367 & vinfo(21)(1:latt))
368 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
369 IF (master) WRITE (stdout,30) 'standard_name', &
370 & trim(vinfo(1)), &
371 & trim(ncname)
372 exit_flag=3
373 ioerror=status
374 END IF
375 END IF
376 END IF
377!
378! Define "long_name" attribute.
379!
380 IF (exit_flag.eq.noerror) THEN
381 IF (len_trim(vinfo(2)).gt.0) THEN
382 vinfo(2)=trim(adjustl(vinfo(2)))
383 latt=len_trim(vinfo(2))
384 status=nf90_put_att(ncid, vid, 'long_name', &
385 & vinfo(2)(1:latt))
386 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
387 IF (master) WRITE (stdout,30) 'long_name', &
388 & trim(vinfo(1)), &
389 & trim(ncname)
390 exit_flag=3
391 ioerror=status
392 END IF
393 END IF
394 END IF
395!
396! Define "size_class" attribute.
397!
398 IF (exit_flag.eq.noerror) THEN
399 latt=len_trim(vinfo(19))
400 IF (latt.gt.0) THEN
401 status=nf90_put_att(ncid, vid, 'size_class', &
402 & vinfo(19)(1:latt))
403 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
404 IF (master) WRITE (stdout,30) 'size_class', &
405 & trim(vinfo(1)), &
406 & trim(ncname)
407 exit_flag=3
408 ioerror=status
409 END IF
410 END IF
411 END IF
412!
413! If applicable, define "units" attribute.
414!
415 IF (exit_flag.eq.noerror) THEN
416 latt=len_trim(vinfo(3))
417 IF (latt.gt.0) THEN
418 IF (trim(vinfo(3)).ne.'nondimensional') THEN
419 status=nf90_put_att(ncid, vid, 'units', &
420 & vinfo(3)(1:latt))
421 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
422 IF (master) WRITE (stdout,30) 'units', &
423 & trim(vinfo(1)), &
424 & trim(ncname)
425 exit_flag=3
426 ioerror=status
427 END IF
428 END IF
429 END IF
430 END IF
431!
432! If applicable, define "calendar" attribute.
433!
434 IF (exit_flag.eq.noerror) THEN
435 latt=len_trim(vinfo(4))
436 IF (latt.gt.0) THEN
437 status=nf90_put_att(ncid, vid, 'calendar', &
438 & vinfo(4)(1:latt))
439 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
440 IF (master) WRITE (stdout,30) 'calendar', &
441 & trim(vinfo(1)), &
442 & trim(ncname)
443 exit_flag=3
444 ioerror=status
445 END IF
446 END IF
447 END IF
448!
449! If applicable, define "valid_min" attribute.
450!
451 IF (exit_flag.eq.noerror) THEN
452 latt=len_trim(vinfo(5))
453 IF (latt.gt.0) THEN
454 IF (vtype.eq.nf90_int) THEN
455 status=nf90_put_att(ncid, vid, trim(vinfo(5)), &
456 & int(aval(2)))
457#ifndef NO_4BYTE_REALS
458 ELSE IF (vtype.eq.nf90_float) THEN
459 status=nf90_put_att(ncid, vid, trim(vinfo(5)), &
460 & real(aval(2),r4))
461#endif
462 ELSE
463 status=nf90_put_att(ncid, vid, trim(vinfo(5)), &
464 & aval(2))
465 END IF
466 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
467 IF (master) WRITE (stdout,30) trim(vinfo(5)), &
468 & trim(vinfo(1)), &
469 & trim(ncname)
470 exit_flag=3
471 ioerror=status
472 END IF
473 aval(2)=0.0_r8
474 END IF
475 END IF
476!
477! If applicable, define "valid_max" attribute.
478!
479 IF (exit_flag.eq.noerror) THEN
480 latt=len_trim(vinfo(6))
481 IF (latt.gt.0) THEN
482 IF (vtype.eq.nf90_int) THEN
483 status=nf90_put_att(ncid, vid, trim(vinfo(6)), &
484 & int(aval(3)))
485#ifndef NO_4BYTE_REALS
486 ELSE IF (vtype.eq.nf90_float) THEN
487 status=nf90_put_att(ncid, vid, trim(vinfo(6)), &
488 & real(aval(3),r4))
489#endif
490 ELSE
491 status=nf90_put_att(ncid, vid, trim(vinfo(6)), &
492 & aval(3))
493 END IF
494 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
495 IF (master) WRITE (stdout,30) trim(vinfo(6)), &
496 & trim(vinfo(1)), &
497 & trim(ncname)
498 exit_flag=3
499 ioerror=status
500 END IF
501 aval(3)=0.0_r8
502 END IF
503 END IF
504!
505! If applicable, define "flag_values" and "flag_meanings" attributes
506! for logical variables.
507!
508 IF (exit_flag.eq.noerror) THEN
509 IF ((len_trim(vinfo(7)).gt.0).and. &
510 & (len_trim(vinfo(8)).gt.0)) THEN
511 text='T, F'
512 latt=len_trim(text)
513 status=nf90_put_att(ncid, vid, 'flag_values', &
514 & text(1:latt))
515 IF (status.eq.nf90_noerr) THEN
516 text=trim(vinfo(7))//' '//trim(vinfo(8))
517 latt=len_trim(text)
518 status=nf90_put_att(ncid, vid, 'flag_meanings', &
519 & text(1:latt))
520 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
521 IF (master) WRITE (stdout,30) 'flag_meanings (T/F)', &
522 & trim(vinfo(1)), &
523 & trim(ncname)
524 exit_flag=3
525 ioerror=status
526 END IF
527 ELSE
528 IF (master) WRITE (stdout,30) 'flag_values (T/F)', &
529 & trim(vinfo(1)), &
530 & trim(ncname)
531 exit_flag=3
532 ioerror=status
533 END IF
534 END IF
535 END IF
536!
537! If applicable, define "flag_values" and "flag_meanings" attributes
538! for integer and floating point variables.
539!
540 IF (exit_flag.eq.noerror) THEN
541 IF ((len_trim(vinfo( 9)).gt.0).and. &
542 & (len_trim(vinfo(10)).gt.0)) THEN
543 IF (vtype.eq.nf90_int) THEN
544 status=nf90_put_att(ncid, vid, 'flag_values', &
545 & (/0,1/))
546#ifndef NO_4BYTE_REALS
547 ELSE IF (vtype.eq.nf90_float) THEN
548 status=nf90_put_att(ncid, vid, 'flag_values', &
549 & (/0.0_r4, 1.0_r4/))
550#endif
551 ELSE
552 status=nf90_put_att(ncid, vid, 'flag_values', &
553 & (/0.0_r8, 1.0_r8/))
554 END IF
555 IF (status.eq.nf90_noerr) THEN
556 text=trim(vinfo(9))//' '//trim(vinfo(10))
557 latt=len_trim(text)
558 status=nf90_put_att(ncid, vid, 'flag_meanings', &
559 & text(1:latt))
560 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
561 IF (master) WRITE (stdout,30) 'flag_meanings', &
562 & trim(vinfo(1)), &
563 & trim(ncname)
564 exit_flag=3
565 ioerror=status
566 END IF
567 ELSE
568 IF (master) WRITE (stdout,30) 'flag_values', &
569 & trim(vinfo(1)), &
570 & trim(ncname)
571 exit_flag=3
572 ioerror=status
573 END IF
574 END IF
575 END IF
576!
577! If applicable, define "negative_value" attribute.
578!
579 IF (exit_flag.eq.noerror) THEN
580 latt=len_trim(vinfo(11))
581 IF (latt.gt.0) THEN
582 status=nf90_put_att(ncid, vid, 'negative_value', &
583 & vinfo(11)(1:latt))
584 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
585 IF (master) WRITE (stdout,30) 'negative_value', &
586 & trim(vinfo(1)), &
587 & trim(ncname)
588 exit_flag=3
589 ioerror=status
590 END IF
591 END IF
592 END IF
593!
594! If applicable, define "positive_value" attribute.
595!
596 IF (exit_flag.eq.noerror) THEN
597 latt=len_trim(vinfo(12))
598 IF (latt.gt.0) THEN
599 status=nf90_put_att(ncid, vid, 'positive_value', &
600 & vinfo(12)(1:latt))
601 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
602 IF (master) WRITE (stdout,30) 'positive_value', &
603 & trim(vinfo(1)), &
604 & trim(ncname)
605 exit_flag=3
606 ioerror=status
607 END IF
608 END IF
609 END IF
610!
611! If applicable, define CF "positive" attribute.
612!
613 IF (exit_flag.eq.noerror) THEN
614 latt=len_trim(vinfo(25))
615 IF (latt.gt.0) THEN
616 status=nf90_put_att(ncid, vid, 'positive', &
617 & vinfo(25)(1:latt))
618 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
619 IF (master) WRITE (stdout,30) 'positive', &
620 & trim(vinfo(1)), &
621 & trim(ncname)
622 exit_flag=3
623 ioerror=status
624 END IF
625 END IF
626 END IF
627!
628! If applicable, define "cycle" attribute.
629!
630 IF (exit_flag.eq.noerror) THEN
631 latt=len_trim(vinfo(13))
632 IF (latt.gt.0) THEN
633 status=nf90_put_att(ncid, vid, 'cycle', &
634 & vinfo(13)(1:latt))
635 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
636 IF (master) WRITE (stdout,30) 'cycle', &
637 & trim(vinfo(1)), &
638 & trim(ncname)
639 exit_flag=3
640 ioerror=status
641 END IF
642 END IF
643 END IF
644!
645! If applicable, define "positions" attribute.
646!
647 IF (exit_flag.eq.noerror) THEN
648 latt=len_trim(vinfo(15))
649 IF (latt.gt.0) THEN
650 status=nf90_put_att(ncid, vid, 'positions', &
651 & vinfo(15)(1:latt))
652 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
653 IF (master) WRITE (stdout,30) 'positions', &
654 & trim(vinfo(1)), &
655 & trim(ncname)
656 exit_flag=3
657 ioerror=status
658 END IF
659 END IF
660 END IF
661!
662! If applicable, define "time" attribute.
663!
664 IF (exit_flag.eq.noerror) THEN
665 latt=len_trim(vinfo(16))
666 IF (latt.gt.0) THEN
667 status=nf90_put_att(ncid, vid, 'time', &
668 & vinfo(16)(1:latt))
669 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
670 IF (master) WRITE (stdout,30) 'time', &
671 & trim(vinfo(1)), &
672 & trim(ncname)
673 exit_flag=3
674 ioerror=status
675 END IF
676 END IF
677 END IF
678!
679! If applicable, define "cell_methods" attribute.
680!
681 IF (exit_flag.eq.noerror) THEN
682 IF (((abs(int(aval(5))).le.4).and.(nvdim.gt.2)).or. &
683 & ((abs(int(aval(5))).gt.4).and.(nvdim.gt.3))) THEN
684 IF ((ncid.eq.avg(ng)%ncid).or. &
685 & (ncid.eq.dia(ng)%ncid)) THEN
686 text='ocean_time: mean'
687 ELSE
688 text='ocean_time: point'
689 END IF
690 latt=len_trim(text)
691 status=nf90_put_att(ncid, vid, 'cell_methods', &
692 & text(1:latt))
693 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
694 IF (master) WRITE (stdout,30) 'cell_methods', &
695 & text(1:latt), &
696 & trim(ncname)
697 exit_flag=3
698 ioerror=status
699 END IF
700 END IF
701 END IF
702!
703! If applicable, define "missing_value" attribute.
704!
705 IF (exit_flag.eq.noerror) THEN
706 latt=len_trim(vinfo(17))
707 IF (latt.gt.0) THEN
708 IF (vtype.eq.nf90_int) THEN
709 status=nf90_put_att(ncid, vid, trim(vinfo(17)), &
710 & int(aval(4)))
711#ifndef NO_4BYTE_REALS
712 ELSE IF (vtype.eq.nf90_float) THEN
713 status=nf90_put_att(ncid, vid, trim(vinfo(17)), &
714 & real(aval(4),r4))
715#endif
716 ELSE
717 status=nf90_put_att(ncid, vid, trim(vinfo(17)), &
718 & aval(4))
719 END IF
720 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
721 IF (master) WRITE (stdout,30) trim(vinfo(17)), &
722 & trim(vinfo(1)), &
723 & trim(ncname)
724 exit_flag=3
725 ioerror=status
726 END IF
727 aval(4)=0.0_r8
728 END IF
729 END IF
730!
731! If applicable, define "add_offset" attribute.
732!
733 IF (exit_flag.eq.noerror) THEN
734 latt=len_trim(vinfo(18))
735 IF (latt.gt.0) THEN
736 IF (vtype.eq.nf90_int) THEN
737 status=nf90_put_att(ncid, vid, trim(vinfo(18)), &
738 & int(aval(1)))
739#ifndef NO_4BYTE_REALS
740 ELSE IF (vtype.eq.nf90_float) THEN
741 status=nf90_put_att(ncid, vid, trim(vinfo(18)), &
742 & real(aval(1),r4))
743#endif
744 ELSE
745 status=nf90_put_att(ncid, vid, trim(vinfo(18)), &
746 & aval(1))
747 END IF
748 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
749 IF (master) WRITE (stdout,30) trim(vinfo(18)), &
750 & trim(vinfo(1)), &
751 & trim(ncname)
752 exit_flag=3
753 ioerror=status
754 END IF
755 aval(1)=0.0_r8
756 END IF
757 END IF
758!
759! If applicable, define "_FillValue" attribute.
760!
761 IF (exit_flag.eq.noerror) THEN
762 latt=len_trim(vinfo(24))
763 IF (latt.gt.0) THEN
764 IF (vtype.eq.nf90_int) THEN
765 status=nf90_put_att(ncid, vid, trim(vinfo(24)), &
766 & int(aval(6)))
767#ifndef NO_4BYTE_REALS
768 ELSE IF (vtype.eq.nf90_float) THEN
769 status=nf90_put_att(ncid, vid, trim(vinfo(24)), &
770 & real(aval(6),r4))
771#endif
772 ELSE
773 status=nf90_put_att(ncid, vid, trim(vinfo(24)), &
774 & aval(6))
775 END IF
776 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
777 IF (master) WRITE (stdout,30) trim(vinfo(24)), &
778 & trim(vinfo(1)), &
779 & trim(ncname)
780 exit_flag=3
781 ioerror=status
782 END IF
783 aval(6)=0.0_r8
784 END IF
785 END IF
786!
787! If applicable, define "water_points" attribute.
788!
789 IF (exit_flag.eq.noerror) THEN
790 latt=len_trim(vinfo(20))
791 IF (latt.gt.0) THEN
792 status=nf90_put_att(ncid, vid, 'water_points', &
793 & vinfo(20)(1:latt))
794 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
795 IF (master) WRITE (stdout,30) 'water_points', &
796 & trim(vinfo(1)), &
797 & trim(ncname)
798 exit_flag=3
799 ioerror=status
800 END IF
801 END IF
802 END IF
803!
804! If applicable, define "coordinates" attribute.
805!
806 IF (exit_flag.eq.noerror) THEN
807 latt=len_trim(vinfo(22))
808 IF (latt.gt.0) THEN
809 IF (spherical) THEN
810 IF (int(aval(5)).eq.r2dvar) THEN
811 text='lon_rho lat_rho'
812 location='face'
813 ELSE IF (int(aval(5)).eq.r3dvar) THEN
814 text='lon_rho lat_rho s_rho'
815 location='face'
816 ELSE IF (int(aval(5)).eq.w3dvar) THEN
817 text='lon_rho lat_rho s_w'
818 location='face'
819 ELSE IF (int(aval(5)).eq.b3dvar) THEN
820 text='lon_rho lat_rho'
821 location='face'
822 ELSE IF (int(aval(5)).eq.u2dvar) THEN
823 text='lon_u lat_u'
824 location='edge1'
825 ELSE IF (int(aval(5)).eq.u3dvar) THEN
826 text='lon_u lat_u s_rho'
827 location='edge1'
828 ELSE IF (int(aval(5)).eq.-u3dvar) THEN
829 text='lon_u lat_u s_w'
830 location='edge1'
831 ELSE IF (int(aval(5)).eq.v2dvar) THEN
832 text='lon_v lat_v'
833 location='edge2'
834 ELSE IF (int(aval(5)).eq.v3dvar) THEN
835 text='lon_v lat_v s_rho'
836 location='edge2'
837 ELSE IF (int(aval(5)).eq.-v3dvar) THEN
838 text='lon_v lat_v s_w'
839 location='edge2'
840 ELSE IF (int(aval(5)).eq.p2dvar) THEN
841 text='lon_psi lat_psi'
842 location='node'
843 ELSE IF (int(aval(5)).eq.p3dvar) THEN
844 text='lon_psi lat_psi s_rho'
845 location='node'
846 ELSE IF (int(aval(5)).eq.l3dvar) THEN
847 text='lon_rho lat_rho light'
848 location='face'
849 ELSE IF (int(aval(5)).eq.l4dvar) THEN
850 text='lon_rho lat_rho s_rho light'
851 location='face'
852 END IF
853 ELSE
854 IF (int(aval(5)).eq.r2dvar) THEN
855 text='x_rho y_rho'
856 location='face'
857 ELSE IF (int(aval(5)).eq.r3dvar) THEN
858 text='x_rho y_rho s_rho'
859 location='face'
860 ELSE IF (int(aval(5)).eq.w3dvar) THEN
861 text='x_rho y_rho s_w'
862 location='face'
863 ELSE IF (int(aval(5)).eq.b3dvar) THEN
864 text='x_rho y_rho'
865 location='face'
866 ELSE IF (int(aval(5)).eq.u2dvar) THEN
867 text='x_u y_u'
868 location='edge1'
869 ELSE IF (int(aval(5)).eq.u3dvar) THEN
870 text='x_u y_u s_rho'
871 location='edge1'
872 ELSE IF (int(aval(5)).eq.-u3dvar) THEN
873 text='x_u y_u s_w'
874 location='edge1'
875 ELSE IF (int(aval(5)).eq.v2dvar) THEN
876 text='x_v y_v'
877 location='edge2'
878 ELSE IF (int(aval(5)).eq.v3dvar) THEN
879 text='x_v y_v s_rho'
880 location='edge2'
881 ELSE IF (int(aval(5)).eq.-v3dvar) THEN
882 text='x_v y_v s_w'
883 location='edge2'
884 ELSE IF (int(aval(5)).eq.p2dvar) THEN
885 text='x_psi y_psi'
886 location='node'
887 ELSE IF (int(aval(5)).eq.p3dvar) THEN
888 text='x_psi y_psi s_rho'
889 location='node'
890 ELSE IF (int(aval(5)).eq.l3dvar) THEN
891 text='x_rho y_rho Nbands'
892 location='face'
893 ELSE IF (int(aval(5)).eq.l4dvar) THEN
894 text='x_rho y_rho s_rho Nbands'
895 location='face'
896 END IF
897 END IF
898
899 status=nf90_put_att(ncid, vid, 'grid', &
900 & 'grid')
901 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
902 IF (master) WRITE (stdout,30) 'grid', &
903 & trim(vinfo(1)), &
904 & trim(ncname)
905 exit_flag=3
906 ioerror=status
907 END IF
908!
909 status=nf90_put_att(ncid, vid, 'location', &
910 & trim(location))
911 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
912 IF (master) WRITE (stdout,30) 'location', &
913 & trim(vinfo(1)), &
914 & trim(ncname)
915 exit_flag=3
916 ioerror=status
917 END IF
918!
919 latt=len_trim(text)
920 IF (((abs(int(aval(5))).le.4).and.(nvdim.gt.2)).or. &
921 & ((abs(int(aval(5))).gt.4).and.(nvdim.gt.3))) THEN
922 text=text(1:latt)//' ocean_time'
923 latt=len_trim(text)
924 END IF
925 status=nf90_put_att(ncid, vid, trim(vinfo(22)), &
926 & text(1:latt))
927 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
928 IF (master) WRITE (stdout,30) 'coordinates', &
929 & trim(vinfo(1)), &
930 & trim(ncname)
931 exit_flag=3
932 ioerror=status
933 END IF
934 aval(5)=0.0_r8
935 END IF
936 END IF
937!
938! If applicable, define "formula_terms" attribute.
939!
940 IF (exit_flag.eq.noerror) THEN
941 latt=len_trim(vinfo(23))
942 IF (latt.gt.0) THEN
943 status=nf90_put_att(ncid, vid, 'formula_terms', &
944 & vinfo(23)(1:latt))
945 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
946 IF (master) WRITE (stdout,30) 'formula_terms', &
947 & trim(vinfo(1)), &
948 & trim(ncname)
949 exit_flag=3
950 ioerror=status
951 END IF
952 END IF
953 END IF
954!
955! If applicable, define "field" attribute (always last).
956!
957 IF (exit_flag.eq.noerror) THEN
958 latt=len_trim(vinfo(14))
959 IF (latt.gt.0) THEN
960 status=nf90_put_att(ncid, vid, 'field', &
961 & vinfo(14)(1:latt))
962 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
963 IF (master) WRITE (stdout,30) 'field', &
964 & trim(vinfo(1)), &
965 & trim(ncname)
966 exit_flag=3
967 ioerror=status
968 END IF
969 END IF
970 END IF
971
972#if defined MASKING && !defined WRITE_WATER
973!
974! If land/sea masking, define "_FillValue" attribute since masking
975! areas are overwritten with special value (spval) during output.
976! Notice that the coordinate attribute is used to check which
977! variables need the "_FillValue" attribute.
978!
979 IF (exit_flag.eq.noerror) THEN
980 latt=len_trim(vinfo(22))
981 IF (PRESENT(setfillval)) THEN
982 landfill=setfillval
983 ELSE
984 landfill=(latt.gt.0).and.(nvdim.gt.2)
985 END IF
986 IF (landfill) THEN
987 IF (vtype.eq.nf90_double) THEN
988 status=nf90_put_att(ncid, vid, '_FillValue', &
989 & spval)
990# ifndef NO_4BYTE_REALS
991 ELSE IF (vtype.eq.nf90_float) THEN
992 status=nf90_put_att(ncid, vid, '_FillValue', &
993 & real(spval,r4))
994# endif
995 END IF
996 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
997 IF (master) WRITE (stdout,30) '_FillValue', &
998 & trim(vinfo(1)), &
999 & trim(ncname)
1000 exit_flag=3
1001 ioerror=status
1002 END IF
1003 END IF
1004 END IF
1005#endif
1006#if defined PARALLEL_IO && defined DISTRIBUTE
1007!
1008! Set parallel access for tiled and non-tiled variables. The optional
1009! argument "SetParAccess" is only available for non-tiled variables.
1010! Otherwise, it is assumed that defined variables are tiled. That is,
1011! they are variables (2D/3D) with parallel domain decomposition.
1012!
1013 IF (exit_flag.eq.noerror) THEN
1014 IF (PRESENT(setparaccess)) THEN
1015 ltiled=setparaccess
1016 ELSE
1017 ltiled=.true.
1018 END IF
1019 IF (ltiled) THEN
1020 status=nf90_var_par_access(ncid, vid, io_tiled_access)
1021 ELSE
1022 status=nf90_var_par_access(ncid, vid, io_nontiled_access)
1023 END IF
1024 IF (founderror(status, nf90_noerr, __line__, myfile)) THEN
1025 IF (master) WRITE (stdout,40) trim(vinfo(1)), &
1026 & trim(ncname)
1027 exit_flag=3
1028 ioerror=status
1029 END IF
1030 END IF
1031#endif
1032 END IF
1033!
1034! Clean information variables.
1035!
1036 DO i=1,SIZE(vinfo)
1037 DO j=1,len(vinfo(1))
1038 vinfo(i)(j:j)=' '
1039 END DO
1040 END DO
1041
1042#if !defined PARALLEL_IO && defined DISTRIBUTE
1043!
1044! Broadcast error flag.
1045!
1046 CALL mp_bcasti (ng, model, exit_flag)
1047#endif
1048!
1049 10 FORMAT (/,' DEF_VAR_NF90 - Grid ',i2.2, &
1050 & ', unable to define variable: ',a,/, &
1051 & 16x,'in NetCDF file: ',a)
1052 20 FORMAT (/,' DEF_VAR_NF90 - error while setting deflate', &
1053 & ' parameters for variable: ',a,/,16x,'in NetCDF file: ',a)
1054 30 FORMAT (/,'DEF_VAR_NF90 - error while defining attribute: ',a, &
1055 & ' for variable: ',a,/,16x,'in NetCDF file: ',a)
1056#if defined PARALLEL_IO && defined DISTRIBUTE
1057 40 FORMAT (/,'DEF_VAR - error while setting parallel access flag', &
1058 & ' for variable: ',a,/,16x,'in NetCDF file: ',a)
1059#endif
1060!
1061 RETURN
1062 END FUNCTION def_var_nf90
1063
1064#if defined PIO_LIB && defined DISTRIBUTE
1065!
1066!***********************************************************************
1067 FUNCTION def_var_pio (ng, model, pioFile, pioVar, &
1068 & Vtype, nVdim, Vdim, &
1069 & Aval, Vinfo, ncname, &
1070 & SetFillVal, SetParAccess) RESULT (status)
1071!***********************************************************************
1072!
1073 USE pio
1074!
1075 USE strings_mod, ONLY : founderror
1076!
1077! Imported variable declarations.
1078!
1079 logical, intent(in), optional :: setfillval
1080 logical, intent(in), optional :: setparaccess
1081!
1082 integer, intent(in) :: ng, model, vtype, nvdim
1083
1084 integer, dimension(:), intent (in) :: vdim
1085!
1086 real(r8), dimension(:), intent(inout) :: aval
1087!
1088 character (len=*), intent(in) :: ncname
1089 character (len=*), intent(inout) :: vinfo(25)
1090!
1091 TYPE (file_desc_t) :: piofile
1092 TYPE (var_desc_t), intent(out) :: piovar
1093!
1094! Local variable declarations.
1095!
1096#ifdef MASKING
1097 logical :: landfill
1098#endif
1099 integer :: fileh
1100 integer :: i, j, latt
1101 integer :: status
1102!
1103 character (len= 5) location
1104 character (len=160) text
1105
1106 character (len=*), parameter :: myfile = &
1107 & __FILE__//", def_var_pio"
1108!
1109!-----------------------------------------------------------------------
1110! Define requested variable and its attributes.
1111!-----------------------------------------------------------------------
1112!
1113 status=pio_noerr
1114!
1115! Get NetCDF file handle from descriptor.
1116!
1117 fileh=abs(piofile%fh)
1118!
1119! Define variable.
1120!
1121 IF (exit_flag.eq.noerror) THEN
1122 IF (len_trim(vinfo(1)).gt.0) THEN
1123 IF ((nvdim.eq.1).and.(vdim(1).eq.0)) THEN
1124 status=pio_def_var(piofile, &
1125 & name = trim(vinfo(1)), &
1126 & type = vtype, &
1127 & vardesc = piovar)
1128 ELSE
1129 status=pio_def_var(piofile, &
1130 & name = trim(vinfo(1)), &
1131 & type = vtype, &
1132 & dimids = vdim(1:nvdim), &
1133 & vardesc = piovar)
1134 END IF
1135 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1136 IF (master) WRITE (stdout,10) ng, trim(vinfo(1)), &
1137 & trim(ncname)
1138 exit_flag=3
1139 ioerror=status
1140 END IF
1141 END IF
1142 END IF
1143!
1144! Define special attributes for SGRID conventions variable "grid".
1145!
1146 IF (exit_flag.eq.noerror) THEN
1147 IF (trim(vinfo(1)).eq.'grid') THEN
1148 status=pio_put_att(piofile, piovar, 'cf_role', &
1149 & 'grid_topology')
1150 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1151 IF (master) WRITE (stdout,30) 'cf_role', &
1152 & trim(vinfo(1)), &
1153 & trim(ncname)
1154 exit_flag=3
1155 ioerror=status
1156 END IF
1157!
1158 status=pio_put_att(piofile, piovar, 'topology_dimension', &
1159 & (/2/))
1160 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1161 IF (master) WRITE (stdout,30) 'topology_dimension', &
1162 & trim(vinfo(1)), &
1163 & trim(ncname)
1164 exit_flag=3
1165 ioerror=status
1166 END IF
1167!
1168 status=pio_put_att(piofile, piovar, 'node_dimensions', &
1169 & 'xi_psi eta_psi')
1170 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1171 IF (master) WRITE (stdout,30) 'node_dimensions', &
1172 & trim(vinfo(1)), &
1173 & trim(ncname)
1174 exit_flag=3
1175 ioerror=status
1176 END IF
1177!
1178 text='xi_rho: xi_psi (padding: both) '// &
1179 & 'eta_rho: eta_psi (padding: both)'
1180 status=pio_put_att(piofile, piovar, 'face_dimensions', &
1181 & trim(text))
1182 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1183 IF (master) WRITE (stdout,30) 'face_dimensions', &
1184 & trim(vinfo(1)), &
1185 & trim(ncname)
1186 exit_flag=3
1187 ioerror=status
1188 END IF
1189!
1190 text='xi_u: xi_psi eta_u: eta_psi (padding: both)'
1191 status=pio_put_att(piofile, piovar, 'edge1_dimensions', &
1192 & trim(text))
1193 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1194 IF (master) WRITE (stdout,30) 'edge1_dimensions', &
1195 & trim(vinfo(1)), &
1196 & trim(ncname)
1197 exit_flag=3
1198 ioerror=status
1199 END IF
1200!
1201 text='xi_v: xi_psi (padding: both) eta_v: eta_psi'
1202 status=pio_put_att(piofile, piovar, 'edge2_dimensions', &
1203 & trim(text))
1204 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1205 IF (master) WRITE (stdout,30) 'edge2_dimensions', &
1206 & trim(vinfo(1)), &
1207 & trim(ncname)
1208 exit_flag=3
1209 ioerror=status
1210 END IF
1211!
1212 IF (spherical) THEN
1213 status=pio_put_att(piofile, piovar, 'node_coordinates', &
1214 & 'lon_psi lat_psi')
1215 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1216 IF (master) WRITE (stdout,30) 'node_coordinates', &
1217 & trim(vinfo(1)), &
1218 & trim(ncname)
1219 exit_flag=3
1220 ioerror=status
1221 END IF
1222!
1223 status=pio_put_att(piofile, piovar, 'face_coordinates', &
1224 & 'lon_rho lat_rho')
1225 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1226 IF (master) WRITE (stdout,30) 'face_coordinates', &
1227 & trim(vinfo(1)), &
1228 & trim(ncname)
1229 exit_flag=3
1230 ioerror=status
1231 END IF
1232!
1233 status=pio_put_att(piofile, piovar, 'edge1_coordinates', &
1234 & 'lon_u lat_u')
1235 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1236 IF (master) WRITE (stdout,30) 'edge1_coordinates', &
1237 & trim(vinfo(1)), &
1238 & trim(ncname)
1239 exit_flag=3
1240 ioerror=status
1241 END IF
1242!
1243 status=pio_put_att(piofile, piovar, 'edge2_coordinates', &
1244 & 'lon_v lat_v')
1245 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1246 IF (master) WRITE (stdout,30) 'edge2_coordinates', &
1247 & trim(vinfo(1)), &
1248 & trim(ncname)
1249 exit_flag=3
1250 ioerror=status
1251 END IF
1252 ELSE
1253 status=pio_put_att(piofile, piovar, 'node_coordinates', &
1254 & 'x_psi y_psi')
1255 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1256 IF (master) WRITE (stdout,30) 'node_coordinates', &
1257 & trim(vinfo(1)), &
1258 & trim(ncname)
1259 exit_flag=3
1260 ioerror=status
1261 END IF
1262!
1263 status=pio_put_att(piofile, piovar, 'face_coordinates', &
1264 & 'x_rho y_rho')
1265 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1266 IF (master) WRITE (stdout,30) 'face_coordinates', &
1267 & trim(vinfo(1)), &
1268 & trim(ncname)
1269 exit_flag=3
1270 ioerror=status
1271 END IF
1272!
1273 status=pio_put_att(piofile, piovar, 'edge1_coordinates', &
1274 & 'x_u y_u')
1275 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1276 IF (master) WRITE (stdout,30) 'edge1_coordinates', &
1277 & trim(vinfo(1)), &
1278 & trim(ncname)
1279 exit_flag=3
1280 ioerror=status
1281 END IF
1282!
1283 status=pio_put_att(piofile, piovar, 'edge2_coordinates', &
1284 & 'x_v y_v')
1285 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1286 IF (master) WRITE (stdout,30) 'edge2_coordinates', &
1287 & trim(vinfo(1)), &
1288 & trim(ncname)
1289 exit_flag=3
1290 ioerror=status
1291 END IF
1292 END IF
1293#ifdef SOLVE3D
1294!
1295 status=pio_put_att(piofile, piovar, 'vertical_dimensions', &
1296 & 's_rho: s_w (padding: none)')
1297 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1298 IF (master) WRITE (stdout,30) 'vertical_dimensions', &
1299 & trim(vinfo(1)), &
1300 & trim(ncname)
1301 exit_flag=3
1302 ioerror=status
1303 END IF
1304#endif
1305 END IF
1306 END IF
1307!
1308! Define "long_name" attribute.
1309!
1310 IF (exit_flag.eq.noerror) THEN
1311 IF (len_trim(vinfo(2)).gt.0) THEN
1312 vinfo(2)=trim(adjustl(vinfo(2)))
1313 latt=len_trim(vinfo(2))
1314 status=pio_put_att(piofile, piovar, 'long_name', &
1315 & vinfo(2)(1:latt))
1316 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1317 IF (master) WRITE (stdout,30) 'long_name', &
1318 & trim(vinfo(1)), &
1319 & trim(ncname)
1320 exit_flag=3
1321 ioerror=status
1322 END IF
1323 END IF
1324 END IF
1325!
1326! Define "size_class" attribute.
1327!
1328 IF (exit_flag.eq.noerror) THEN
1329 latt=len_trim(vinfo(19))
1330 IF (latt.gt.0) THEN
1331 status=pio_put_att(piofile, piovar, 'size_class', &
1332 & vinfo(19)(1:latt))
1333 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1334 IF (master) WRITE (stdout,30) 'size_class', &
1335 & trim(vinfo(1)), &
1336 & trim(ncname)
1337 exit_flag=3
1338 ioerror=status
1339 END IF
1340 END IF
1341 END IF
1342!
1343! If applicable, define "units" attribute.
1344!
1345 IF (exit_flag.eq.noerror) THEN
1346 latt=len_trim(vinfo(3))
1347 IF (latt.gt.0) THEN
1348 IF (trim(vinfo(3)).ne.'nondimensional') THEN
1349 status=pio_put_att(piofile, piovar, 'units', &
1350 & vinfo(3)(1:latt))
1351 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1352 IF (master) WRITE (stdout,30) 'units', &
1353 & trim(vinfo(1)), &
1354 & trim(ncname)
1355 exit_flag=3
1356 ioerror=status
1357 END IF
1358 END IF
1359 END IF
1360 END IF
1361!
1362! If applicable, define "calendar" attribute.
1363!
1364 IF (exit_flag.eq.noerror) THEN
1365 latt=len_trim(vinfo(4))
1366 IF (latt.gt.0) THEN
1367 status=pio_put_att(piofile, piovar, 'calendar', &
1368 & vinfo(4)(1:latt))
1369 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1370 IF (master) WRITE (stdout,30) 'calendar', &
1371 & trim(vinfo(1)), &
1372 & trim(ncname)
1373 exit_flag=3
1374 ioerror=status
1375 END IF
1376 END IF
1377 END IF
1378!
1379! If applicable, define "valid_min" attribute.
1380!
1381 IF (exit_flag.eq.noerror) THEN
1382 latt=len_trim(vinfo(5))
1383 IF (latt.gt.0) THEN
1384 IF (vtype.eq.pio_int) THEN
1385 status=pio_put_att(piofile, piovar, trim(vinfo(5)), &
1386 & int(aval(2)))
1387#ifndef NO_4BYTE_REALS
1388 ELSE IF (vtype.eq.pio_real) THEN
1389 status=pio_put_att(piofile, piovar, trim(vinfo(5)), &
1390 & real(aval(2),r4))
1391#endif
1392 ELSE
1393 status=pio_put_att(piofile, piovar, trim(vinfo(5)), &
1394 & aval(2))
1395 END IF
1396 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1397 IF (master) WRITE (stdout,30) trim(vinfo(5)), &
1398 & trim(vinfo(1)), &
1399 & trim(ncname)
1400 exit_flag=3
1401 ioerror=status
1402 END IF
1403 aval(2)=0.0_r8
1404 END IF
1405 END IF
1406!
1407! If applicable, define "valid_max" attribute.
1408!
1409 IF (exit_flag.eq.noerror) THEN
1410 latt=len_trim(vinfo(6))
1411 IF (latt.gt.0) THEN
1412 IF (vtype.eq.pio_int) THEN
1413 status=pio_put_att(piofile, piovar, trim(vinfo(6)), &
1414 & int(aval(3)))
1415#ifndef NO_4BYTE_REALS
1416 ELSE IF (vtype.eq.pio_real) THEN
1417 status=pio_put_att(piofile, piovar, trim(vinfo(6)), &
1418 & real(aval(3),r4))
1419#endif
1420 ELSE
1421 status=pio_put_att(piofile, piovar, trim(vinfo(6)), &
1422 & aval(3))
1423 END IF
1424 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1425 IF (master) WRITE (stdout,30) trim(vinfo(6)), &
1426 & trim(vinfo(1)), &
1427 & trim(ncname)
1428 exit_flag=3
1429 ioerror=status
1430 END IF
1431 aval(3)=0.0_r8
1432 END IF
1433 END IF
1434!
1435! If applicable, define "flag_values" and "flag_meanings" attributes
1436! for logical variables.
1437!
1438 IF (exit_flag.eq.noerror) THEN
1439 IF ((len_trim(vinfo(7)).gt.0).and. &
1440 & (len_trim(vinfo(8)).gt.0)) THEN
1441 text='T, F'
1442 latt=len_trim(text)
1443 status=pio_put_att(piofile, piovar, 'flag_values', &
1444 & text(1:latt))
1445 IF (status.eq.pio_noerr) THEN
1446 text=trim(vinfo(7))//' '//trim(vinfo(8))
1447 latt=len_trim(text)
1448 status=pio_put_att(piofile, piovar, 'flag_meanings', &
1449 & text(1:latt))
1450 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1451 IF (master) WRITE (stdout,30) 'flag_meanings (T/F)', &
1452 & trim(vinfo(1)), &
1453 & trim(ncname)
1454 exit_flag=3
1455 ioerror=status
1456 END IF
1457 ELSE
1458 IF (master) WRITE (stdout,30) 'flag_values (T/F)', &
1459 & trim(vinfo(1)), &
1460 & trim(ncname)
1461 exit_flag=3
1462 ioerror=status
1463 END IF
1464 END IF
1465 END IF
1466!
1467! If applicable, define "flag_values" and "flag_meanings" attributes
1468! for integer and floating point variables.
1469!
1470 IF (exit_flag.eq.noerror) THEN
1471 IF ((len_trim(vinfo( 9)).gt.0).and. &
1472 & (len_trim(vinfo(10)).gt.0)) THEN
1473 IF (vtype.eq.pio_int) THEN
1474 status=pio_put_att(piofile, piovar, 'flag_values', &
1475 & (/0,1/))
1476#ifndef NO_4BYTE_REALS
1477 ELSE IF (vtype.eq.pio_real) THEN
1478 status=pio_put_att(piofile, piovar, 'flag_values', &
1479 & (/0.0_r4, 1.0_r4/))
1480#endif
1481 ELSE
1482 status=pio_put_att(piofile, piovar, 'flag_values', &
1483 & (/0.0_r8, 1.0_r8/))
1484 END IF
1485 IF (status.eq.pio_noerr) THEN
1486 text=trim(vinfo(9))//' '//trim(vinfo(10))
1487 latt=len_trim(text)
1488 status=pio_put_att(piofile, piovar, 'flag_meanings', &
1489 & text(1:latt))
1490 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1491 IF (master) WRITE (stdout,30) 'flag_meanings', &
1492 & trim(vinfo(1)), &
1493 & trim(ncname)
1494 exit_flag=3
1495 ioerror=status
1496 END IF
1497 ELSE
1498 IF (master) WRITE (stdout,30) 'flag_values', &
1499 & trim(vinfo(1)), &
1500 & trim(ncname)
1501 exit_flag=3
1502 ioerror=status
1503 END IF
1504 END IF
1505 END IF
1506!
1507! If applicable, define "negative_value" attribute.
1508!
1509 IF (exit_flag.eq.noerror) THEN
1510 latt=len_trim(vinfo(11))
1511 IF (latt.gt.0) THEN
1512 status=pio_put_att(piofile, piovar, 'negative_value', &
1513 & vinfo(11)(1:latt))
1514 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1515 IF (master) WRITE (stdout,30) 'negative_value', &
1516 & trim(vinfo(1)), &
1517 & trim(ncname)
1518 exit_flag=3
1519 ioerror=status
1520 END IF
1521 END IF
1522 END IF
1523!
1524! If applicable, define "positive_value" attribute.
1525!
1526 IF (exit_flag.eq.noerror) THEN
1527 latt=len_trim(vinfo(12))
1528 IF (latt.gt.0) THEN
1529 status=pio_put_att(piofile, piovar, 'positive_value', &
1530 & vinfo(12)(1:latt))
1531 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1532 IF (master) WRITE (stdout,30) 'positive_value', &
1533 & trim(vinfo(1)), &
1534 & trim(ncname)
1535 exit_flag=3
1536 ioerror=status
1537 END IF
1538 END IF
1539 END IF
1540!
1541! If applicable, define CF "positive" attribute.
1542!
1543 IF (exit_flag.eq.noerror) THEN
1544 latt=len_trim(vinfo(25))
1545 IF (latt.gt.0) THEN
1546 status=pio_put_att(piofile, piovar, 'positive', &
1547 & vinfo(25)(1:latt))
1548 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1549 IF (master) WRITE (stdout,30) 'positive', &
1550 & trim(vinfo(1)), &
1551 & trim(ncname)
1552 exit_flag=3
1553 ioerror=status
1554 END IF
1555 END IF
1556 END IF
1557!
1558! If applicable, define "cycle" attribute.
1559!
1560 IF (exit_flag.eq.noerror) THEN
1561 latt=len_trim(vinfo(13))
1562 IF (latt.gt.0) THEN
1563 status=pio_put_att(piofile, piovar, 'cycle', &
1564 & vinfo(13)(1:latt))
1565 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1566 IF (master) WRITE (stdout,30) 'cycle', &
1567 & trim(vinfo(1)), &
1568 & trim(ncname)
1569 exit_flag=3
1570 ioerror=status
1571 END IF
1572 END IF
1573 END IF
1574!
1575! If applicable, define "positions" attribute.
1576!
1577 IF (exit_flag.eq.noerror) THEN
1578 latt=len_trim(vinfo(15))
1579 IF (latt.gt.0) THEN
1580 status=pio_put_att(piofile, piovar, 'positions', &
1581 & vinfo(15)(1:latt))
1582 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1583 IF (master) WRITE (stdout,30) 'positions', &
1584 & trim(vinfo(1)), &
1585 & trim(ncname)
1586 exit_flag=3
1587 ioerror=status
1588 END IF
1589 END IF
1590 END IF
1591!
1592! If applicable, define "time" attribute.
1593!
1594 IF (exit_flag.eq.noerror) THEN
1595 latt=len_trim(vinfo(16))
1596 IF (latt.gt.0) THEN
1597 status=pio_put_att(piofile, piovar, 'time', &
1598 & vinfo(16)(1:latt))
1599 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1600 IF (master) WRITE (stdout,30) 'time', &
1601 & trim(vinfo(1)), &
1602 & trim(ncname)
1603 exit_flag=3
1604 ioerror=status
1605 END IF
1606 END IF
1607 END IF
1608!
1609! If applicable, define "cell_methods" attribute.
1610!
1611 IF (exit_flag.eq.noerror) THEN
1612 IF (((abs(int(aval(5))).le.4).and.(nvdim.gt.2)).or. &
1613 & ((abs(int(aval(5))).gt.4).and.(nvdim.gt.3))) THEN
1614 IF ((fileh.eq.abs(avg(ng)%pioFile%fh)).or. &
1615 & (fileh.eq.abs(dia(ng)%pioFile%fh))) THEN
1616 text='ocean_time: mean'
1617 ELSE
1618 text='ocean_time: point'
1619 END IF
1620 latt=len_trim(text)
1621 status=pio_put_att(piofile, piovar, 'cell_methods', &
1622 & text(1:latt))
1623 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1624 IF (master) WRITE (stdout,30) 'cell_methods', &
1625 & text(1:latt), &
1626 & trim(ncname)
1627 exit_flag=3
1628 ioerror=status
1629 END IF
1630 END IF
1631 END IF
1632!
1633! If applicable, define "missing_value" attribute.
1634!
1635 IF (exit_flag.eq.noerror) THEN
1636 latt=len_trim(vinfo(17))
1637 IF (latt.gt.0) THEN
1638 IF (vtype.eq.pio_int) THEN
1639 status=pio_put_att(piofile, piovar, trim(vinfo(17)), &
1640 & int(aval(4)))
1641#ifndef NO_4BYTE_REALS
1642 ELSE IF (vtype.eq.pio_real) THEN
1643 status=pio_put_att(piofile, piovar, trim(vinfo(17)), &
1644 & real(aval(4),r4))
1645#endif
1646 ELSE
1647 status=pio_put_att(piofile, piovar, trim(vinfo(17)), &
1648 & aval(4))
1649 END IF
1650 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1651 IF (master) WRITE (stdout,30) trim(vinfo(17)), &
1652 & trim(vinfo(1)), &
1653 & trim(ncname)
1654 exit_flag=3
1655 ioerror=status
1656 END IF
1657 aval(4)=0.0_r8
1658 END IF
1659 END IF
1660!
1661! If applicable, define "add_offset" attribute.
1662!
1663 IF (exit_flag.eq.noerror) THEN
1664 latt=len_trim(vinfo(18))
1665 IF (latt.gt.0) THEN
1666 IF (vtype.eq.pio_int) THEN
1667 status=pio_put_att(piofile, piovar, trim(vinfo(18)), &
1668 & int(aval(1)))
1669#ifndef NO_4BYTE_REALS
1670 ELSE IF (vtype.eq.pio_real) THEN
1671 status=pio_put_att(piofile, piovar, trim(vinfo(18)), &
1672 & real(aval(1),r4))
1673#endif
1674 ELSE
1675 status=pio_put_att(piofile, piovar, trim(vinfo(18)), &
1676 & aval(1))
1677 END IF
1678 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1679 IF (master) WRITE (stdout,30) trim(vinfo(18)), &
1680 & trim(vinfo(1)), &
1681 & trim(ncname)
1682 exit_flag=3
1683 ioerror=status
1684 END IF
1685 aval(1)=0.0_r8
1686 END IF
1687 END IF
1688!
1689! If applicable, define "_FillValue" attribute.
1690!
1691 IF (exit_flag.eq.noerror) THEN
1692 latt=len_trim(vinfo(24))
1693 IF (latt.gt.0) THEN
1694 IF (vtype.eq.pio_int) THEN
1695 status=pio_put_att(piofile, piovar, trim(vinfo(24)), &
1696 & int(aval(6)))
1697#ifndef NO_4BYTE_REALS
1698 ELSE IF (vtype.eq.pio_real) THEN
1699 status=pio_put_att(piofile, piovar, trim(vinfo(24)), &
1700 & real(aval(6),r4))
1701#endif
1702 ELSE
1703 status=pio_put_att(piofile, piovar, trim(vinfo(24)), &
1704 & aval(6))
1705 END IF
1706 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1707 IF (master) WRITE (stdout,30) trim(vinfo(24)), &
1708 & trim(vinfo(1)), &
1709 & trim(ncname)
1710 exit_flag=3
1711 ioerror=status
1712 END IF
1713 aval(6)=0.0_r8
1714 END IF
1715 END IF
1716!
1717! If applicable, define "water_points" attribute.
1718!
1719 IF (exit_flag.eq.noerror) THEN
1720 latt=len_trim(vinfo(20))
1721 IF (latt.gt.0) THEN
1722 status=pio_put_att(piofile, piovar, 'water_points', &
1723 & vinfo(20)(1:latt))
1724 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1725 IF (master) WRITE (stdout,30) 'water_points', &
1726 & trim(vinfo(1)), &
1727 & trim(ncname)
1728 exit_flag=3
1729 ioerror=status
1730 END IF
1731 END IF
1732 END IF
1733!
1734! If applicable, define "standard_name" attribute.
1735!
1736 IF (exit_flag.eq.noerror) THEN
1737 latt=len_trim(vinfo(21))
1738 IF (latt.gt.0) THEN
1739 status=pio_put_att(piofile, piovar, 'standard_name', &
1740 & vinfo(21)(1:latt))
1741 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1742 IF (master) WRITE (stdout,30) 'standard_name', &
1743 & trim(vinfo(1)), &
1744 & trim(ncname)
1745 exit_flag=3
1746 ioerror=status
1747 END IF
1748 END IF
1749 END IF
1750!
1751! If applicable, define "coordinates" attribute.
1752!
1753 IF (exit_flag.eq.noerror) THEN
1754 latt=len_trim(vinfo(22))
1755 IF (latt.gt.0) THEN
1756 IF (spherical) THEN
1757 IF (int(aval(5)).eq.r2dvar) THEN
1758 text='lon_rho lat_rho'
1759 location='face'
1760 ELSE IF (int(aval(5)).eq.r3dvar) THEN
1761 text='lon_rho lat_rho s_rho'
1762 location='face'
1763 ELSE IF (int(aval(5)).eq.w3dvar) THEN
1764 text='lon_rho lat_rho s_w'
1765 location='face'
1766 ELSE IF (int(aval(5)).eq.b3dvar) THEN
1767 text='lon_rho lat_rho'
1768 location='face'
1769 ELSE IF (int(aval(5)).eq.u2dvar) THEN
1770 text='lon_u lat_u'
1771 location='edge1'
1772 ELSE IF (int(aval(5)).eq.u3dvar) THEN
1773 text='lon_u lat_u s_rho'
1774 location='edge1'
1775 ELSE IF (int(aval(5)).eq.-u3dvar) THEN
1776 text='lon_u lat_u s_w'
1777 location='edge1'
1778 ELSE IF (int(aval(5)).eq.v2dvar) THEN
1779 text='lon_v lat_v'
1780 location='edge2'
1781 ELSE IF (int(aval(5)).eq.v3dvar) THEN
1782 text='lon_v lat_v s_rho'
1783 location='edge2'
1784 ELSE IF (int(aval(5)).eq.-v3dvar) THEN
1785 text='lon_v lat_v s_w'
1786 location='edge2'
1787 ELSE IF (int(aval(5)).eq.p2dvar) THEN
1788 text='lon_psi lat_psi'
1789 location='node'
1790 ELSE IF (int(aval(5)).eq.p3dvar) THEN
1791 text='lon_psi lat_psi s_rho'
1792 location='node'
1793 ELSE IF (int(aval(5)).eq.l3dvar) THEN
1794 text='lon_rho lat_rho light'
1795 location='face'
1796 ELSE IF (int(aval(5)).eq.l4dvar) THEN
1797 text='lon_rho lat_rho s_rho light'
1798 location='face'
1799 END IF
1800 ELSE
1801 IF (int(aval(5)).eq.r2dvar) THEN
1802 text='x_rho y_rho'
1803 location='face'
1804 ELSE IF (int(aval(5)).eq.r3dvar) THEN
1805 text='x_rho y_rho s_rho'
1806 location='face'
1807 ELSE IF (int(aval(5)).eq.w3dvar) THEN
1808 text='x_rho y_rho s_w'
1809 location='face'
1810 ELSE IF (int(aval(5)).eq.b3dvar) THEN
1811 text='x_rho y_rho'
1812 location='face'
1813 ELSE IF (int(aval(5)).eq.u2dvar) THEN
1814 text='x_u y_u'
1815 location='edge1'
1816 ELSE IF (int(aval(5)).eq.u3dvar) THEN
1817 text='x_u y_u s_rho'
1818 location='edge1'
1819 ELSE IF (int(aval(5)).eq.-u3dvar) THEN
1820 text='x_u y_u s_w'
1821 location='edge1'
1822 ELSE IF (int(aval(5)).eq.v2dvar) THEN
1823 text='x_v y_v'
1824 location='edge2'
1825 ELSE IF (int(aval(5)).eq.v3dvar) THEN
1826 text='x_v y_v s_rho'
1827 location='edge2'
1828 ELSE IF (int(aval(5)).eq.-v3dvar) THEN
1829 text='x_v y_v s_w'
1830 location='edge2'
1831 ELSE IF (int(aval(5)).eq.p2dvar) THEN
1832 text='x_psi y_psi'
1833 location='node'
1834 ELSE IF (int(aval(5)).eq.p3dvar) THEN
1835 text='x_psi y_psi s_rho'
1836 location='node'
1837 ELSE IF (int(aval(5)).eq.l3dvar) THEN
1838 text='x_rho y_rho Nbands'
1839 location='face'
1840 ELSE IF (int(aval(5)).eq.l4dvar) THEN
1841 text='x_rho y_rho s_rho Nbands'
1842 location='face'
1843 END IF
1844 END IF
1845
1846 status=pio_put_att(piofile, piovar, 'grid', &
1847 & 'grid')
1848 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1849 IF (master) WRITE (stdout,30) 'grid', &
1850 & trim(vinfo(1)), &
1851 & trim(ncname)
1852 exit_flag=3
1853 ioerror=status
1854 END IF
1855!
1856 status=pio_put_att(piofile, piovar, 'location', &
1857 & trim(location))
1858 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1859 IF (master) WRITE (stdout,30) 'location', &
1860 & trim(vinfo(1)), &
1861 & trim(ncname)
1862 exit_flag=3
1863 ioerror=status
1864 END IF
1865!
1866 latt=len_trim(text)
1867 IF (((abs(int(aval(5))).le.4).and.(nvdim.gt.2)).or. &
1868 & ((abs(int(aval(5))).gt.4).and.(nvdim.gt.3))) THEN
1869 text=text(1:latt)//' ocean_time'
1870 latt=len_trim(text)
1871 END IF
1872 status=pio_put_att(piofile, piovar, trim(vinfo(22)), &
1873 & text(1:latt))
1874 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1875 IF (master) WRITE (stdout,30) 'coordinates', &
1876 & trim(vinfo(1)), &
1877 & trim(ncname)
1878 exit_flag=3
1879 ioerror=status
1880 END IF
1881 aval(5)=0.0_r8
1882 END IF
1883 END IF
1884!
1885! If applicable, define "formula_terms" attribute.
1886!
1887 IF (exit_flag.eq.noerror) THEN
1888 latt=len_trim(vinfo(23))
1889 IF (latt.gt.0) THEN
1890 status=pio_put_att(piofile, piovar, 'formula_terms', &
1891 & vinfo(23)(1:latt))
1892 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1893 IF (master) WRITE (stdout,30) 'formula_terms', &
1894 & trim(vinfo(1)), &
1895 & trim(ncname)
1896 exit_flag=3
1897 ioerror=status
1898 END IF
1899 END IF
1900 END IF
1901!
1902! If applicable, define "field" attribute (always last).
1903!
1904 IF (exit_flag.eq.noerror) THEN
1905 latt=len_trim(vinfo(14))
1906 IF (latt.gt.0) THEN
1907 status=pio_put_att(piofile, piovar, 'field', &
1908 & vinfo(14)(1:latt))
1909 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1910 IF (master) WRITE (stdout,30) 'field', &
1911 & trim(vinfo(1)), &
1912 & trim(ncname)
1913 exit_flag=3
1914 ioerror=status
1915 END IF
1916 END IF
1917 END IF
1918
1919#if defined MASKING && !defined WRITE_WATER
1920!
1921! If land/sea masking, define "_FillValue" attribute since masking
1922! areas are overwritten with special value (spval) during output.
1923! Notice that the coordinate attribute is used to check which
1924! variables need the "_FillValue" attribute.
1925!
1926 IF (exit_flag.eq.noerror) THEN
1927 latt=len_trim(vinfo(22))
1928 IF (PRESENT(setfillval)) THEN
1929 landfill=setfillval
1930 ELSE
1931 landfill=(latt.gt.0).and.(nvdim.gt.2)
1932 END IF
1933 IF (landfill) THEN
1934 IF (vtype.eq.pio_double) THEN
1935 status=pio_put_att(piofile, piovar, '_FillValue', &
1936 & spval)
1937# ifndef NO_4BYTE_REALS
1938 ELSE IF (vtype.eq.pio_real) THEN
1939 status=pio_put_att(piofile, piovar, '_FillValue', &
1940 & real(spval,r4))
1941# endif
1942 END IF
1943 IF (founderror(status, pio_noerr, __line__, myfile)) THEN
1944 IF (master) WRITE (stdout,30) '_FillValue', &
1945 & trim(vinfo(1)), &
1946 & trim(ncname)
1947 exit_flag=3
1948 ioerror=status
1949 END IF
1950 END IF
1951 END IF
1952#endif
1953!
1954! Clean information variables.
1955!
1956 DO i=1,SIZE(vinfo)
1957 DO j=1,len(vinfo(1))
1958 vinfo(i)(j:j)=' '
1959 END DO
1960 END DO
1961!
1962 10 FORMAT (/,' DEF_VAR_PIO - Grid ',i2.2, &
1963 & ', unable to define variable: ',a,/, &
1964 & 15x,'in NetCDF file: ',a)
1965 20 FORMAT (/,' DEF_VAR_PIO - error while setting deflate', &
1966 & ' parameters for variable: ',a,/,15x,'in NetCDF file: ',a)
1967 30 FORMAT (/,'DEF_VAR_NF90 - error while defining attribute: ',a, &
1968 & ' for variable: ',a,/,15x,'in NetCDF file: ',a)
1969!
1970 RETURN
1971 END FUNCTION def_var_pio
1972#endif
1973
1974 END MODULE def_var_mod
integer function def_var_pio(ng, model, piofile, piovar, vtype, nvdim, vdim, aval, vinfo, ncname, setfillval, setparaccess)
Definition def_var.F:1071
integer function def_var_nf90(ng, model, ncid, vid, vtype, nvdim, vdim, aval, vinfo, ncname, setfillval, setparaccess)
Definition def_var.F:102
integer ioerror
type(t_io), dimension(:), allocatable avg
integer stdout
type(t_io), dimension(:), allocatable dia
integer shuffle
Definition mod_netcdf.F:222
integer deflate
Definition mod_netcdf.F:223
integer, parameter io_tiled_access
Definition mod_netcdf.F:253
integer, parameter io_nontiled_access
Definition mod_netcdf.F:252
integer deflate_level
Definition mod_netcdf.F:224
logical master
logical outthread
integer, parameter b3dvar
Definition mod_param.F:725
integer, parameter r3dvar
Definition mod_param.F:721
integer, parameter l4dvar
Definition mod_param.F:727
integer, parameter u3dvar
Definition mod_param.F:722
integer, parameter u2dvar
Definition mod_param.F:718
integer, parameter w3dvar
Definition mod_param.F:724
integer, parameter p2dvar
Definition mod_param.F:716
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter l3dvar
Definition mod_param.F:726
integer, parameter v2dvar
Definition mod_param.F:719
integer, parameter p3dvar
Definition mod_param.F:720
integer, parameter v3dvar
Definition mod_param.F:723
logical spherical
real(dp), parameter spval
integer exit_flag
integer noerror
logical function, public founderror(flag, noerr, line, routine)
Definition strings.F:52