/[lmdze]/trunk/IOIPSL/Histcom/histvert.f
ViewVC logotype

Contents of /trunk/IOIPSL/Histcom/histvert.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 82 - (show annotations)
Wed Mar 5 14:57:53 2014 UTC (10 years, 3 months ago) by guez
File size: 4668 byte(s)
Changed all ".f90" suffixes to ".f".
1 module histvert_m
2
3 implicit none
4
5 contains
6
7 SUBROUTINE histvert(pfileid, pzaxname, pzaxtitle, pzaxunit, pzvalues, &
8 pzaxid, pdirect)
9
10 ! This subroutine defines a vertical axis and returns its id. It
11 ! gives the user the possibility to define different vertical
12 ! axes. For each variable defined with histdef a vertical axis can
13 ! be specified by its ID.
14
15 USE errioipsl, ONLY: histerr
16 USE find_str_m, ONLY: find_str
17 USE histcom_var, ONLY: nb_zax, nb_zax_max, ncdf_ids, zax_ids, &
18 zax_name, zax_name_length, zax_size
19 USE netcdf, ONLY: nf90_def_dim, nf90_def_var, nf90_enddef, &
20 nf90_float, nf90_put_att, nf90_put_var, nf90_redef
21 USE strlowercase_m, ONLY: strlowercase
22
23 INTEGER, INTENT(IN):: pfileid
24 ! ID of the file the variable should be archived in
25
26 CHARACTER(len=*), INTENT(IN):: pzaxname ! name of the vertical axis
27 CHARACTER(len=*), INTENT(IN):: pzaxtitle ! title of the vertical axis
28 CHARACTER(len=*), INTENT(IN):: pzaxunit ! units of the vertical axis
29
30 REAL, INTENT(IN):: pzvalues(:) ! (pzsize)
31 ! coordinate values of the vertical axis
32
33 INTEGER, INTENT(OUT):: pzaxid ! ID of the axis (not the netCDF ID)
34
35 CHARACTER(len=*), INTENT(IN), OPTIONAL:: pdirect
36 ! positive direction of the axis: up or down
37
38 ! Local:
39 INTEGER pzsize ! size of the vertical axis
40 INTEGER:: pos, iv, nb, zdimid, zaxid_tmp
41 CHARACTER(len=20):: str20, tab_str20(nb_zax_max)
42 INTEGER:: tab_str20_length(nb_zax_max)
43 CHARACTER(len=70):: str70, str71, str72
44 CHARACTER(len=80):: str80
45 CHARACTER(len=20):: direction
46 INTEGER:: iret, leng, ncid
47
48 !---------------------------------------------------------------------
49
50 pzsize = size(pzvalues)
51
52 ! 1.0 Verifications:
53 ! Do we have enough space for an extra axis ?
54 ! Is the name already in use ?
55
56 ! - Direction of axis. Can we get if from the user.
57 ! If not we put unknown.
58
59 IF (present(pdirect)) THEN
60 direction = trim(pdirect)
61 CALL strlowercase(direction)
62 ELSE
63 direction = 'unknown'
64 END IF
65
66 ! Check the consistency of the attribute
67
68 IF ((direction/='unknown') .AND. (direction/='up') .AND. &
69 (direction/='down')) THEN
70 direction = 'unknown'
71 str80 = 'The specified axis was: ' // trim(direction)
72 CALL histerr(2, 'histvert', &
73 'The specified direction for the vertical axis is not possible.', &
74 'it is replaced by: unknown', str80)
75 END IF
76
77 IF (nb_zax(pfileid)+1>nb_zax_max) THEN
78 CALL histerr(3, 'histvert', &
79 'Table of vertical axes too small. You should increase ', &
80 'nb_zax_max in M_HISTCOM.f90 in order to accomodate all ', &
81 'these variables ')
82 END IF
83
84 iv = nb_zax(pfileid)
85 IF (iv>1) THEN
86 str20 = pzaxname
87 nb = iv - 1
88 tab_str20(1:nb) = zax_name(pfileid, 1:nb)
89 tab_str20_length(1:nb) = zax_name_length(pfileid, 1:nb)
90 CALL find_str(nb, tab_str20, tab_str20_length, str20, pos)
91 ELSE
92 pos = 0
93 END IF
94
95 IF (pos>0) THEN
96 str70 = 'Vertical axis already exists'
97 WRITE(str71, '("Check variable ", a, " in file", I3)') str20, &
98 pfileid
99 str72 = 'Can also be a wrong file ID in another declaration'
100 CALL histerr(3, 'histvert', str70, str71, str72)
101 END IF
102
103 iv = nb_zax(pfileid) + 1
104
105 ! 2.0 Add the information to the file
106
107 ncid = ncdf_ids(pfileid)
108
109 leng = min(len_trim(pzaxname), 20)
110 iret = nf90_def_dim(ncid, pzaxname(1:leng), pzsize, zaxid_tmp)
111 iret = nf90_def_var(ncid, pzaxname(1:leng), nf90_float, zaxid_tmp, zdimid)
112
113 leng = min(len_trim(pzaxunit), 20)
114 iret = nf90_put_att(ncid, zdimid, 'units', pzaxunit(1:leng))
115 iret = nf90_put_att(ncid, zdimid, 'positive', trim(direction))
116
117 iret = nf90_put_att(ncid, zdimid, 'valid_min', real(minval( &
118 pzvalues(1:pzsize))))
119 iret = nf90_put_att(ncid, zdimid, 'valid_max', real(maxval( &
120 pzvalues(1:pzsize))))
121
122 leng = min(len_trim(pzaxname), 20)
123 iret = nf90_put_att(ncid, zdimid, 'title', pzaxname(1:leng))
124 leng = min(len_trim(pzaxtitle), 80)
125 iret = nf90_put_att(ncid, zdimid, 'long_name', pzaxtitle(1:leng))
126
127 iret = nf90_enddef(ncid)
128
129 iret = nf90_put_var(ncid, zdimid, pzvalues(1:pzsize))
130
131 iret = nf90_redef(ncid)
132
133 ! 3.0 add the information to the common
134
135 nb_zax(pfileid) = iv
136 zax_size(pfileid, iv) = pzsize
137 zax_name(pfileid, iv) = pzaxname
138 zax_name_length(pfileid, iv) = len_trim(pzaxname)
139 zax_ids(pfileid, iv) = zaxid_tmp
140 pzaxid = iv
141
142 END SUBROUTINE histvert
143
144 end module histvert_m

  ViewVC Help
Powered by ViewVC 1.1.21