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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 328 - (show annotations)
Thu Jun 13 14:40:06 2019 UTC (4 years, 11 months ago) by guez
File size: 4663 byte(s)
Change all `.f` suffixes to `.f90`. (The opposite was done in revision
82.)  Because of change of philosopy in GNUmakefile: we already had a
rewritten rule for `.f`, so it does not make the makefile longer to
replace it by a rule for `.f90`. And it spares us options of
makedepf90 and of the compiler. Also we prepare the way for a simpler
`CMakeLists.txt`.

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 netcdf95, ONLY: nf95_def_dim, nf95_def_var, nf95_enddef, &
20 nf95_put_att, nf95_put_var, nf95_redef
21 use netcdf, only: nf90_float
22 USE strlowercase_m, ONLY: strlowercase
23
24 INTEGER, INTENT(IN):: pfileid
25 ! ID of the file the variable should be archived in
26
27 CHARACTER(len=*), INTENT(IN):: pzaxname ! name of the vertical axis
28 CHARACTER(len=*), INTENT(IN):: pzaxtitle ! title of the vertical axis
29 CHARACTER(len=*), INTENT(IN):: pzaxunit ! units of the vertical axis
30
31 REAL, INTENT(IN):: pzvalues(:) ! (pzsize)
32 ! coordinate values of the vertical axis
33
34 INTEGER, INTENT(OUT):: pzaxid ! ID of the axis (not the netCDF ID)
35
36 CHARACTER(len=*), INTENT(IN), OPTIONAL:: pdirect
37 ! positive direction of the axis: up or down
38
39 ! Local:
40 INTEGER pzsize ! size of the vertical axis
41 INTEGER:: pos, iv, nb, zdimid, zaxid_tmp
42 CHARACTER(len=20):: str20, tab_str20(nb_zax_max)
43 INTEGER:: tab_str20_length(nb_zax_max)
44 CHARACTER(len=70):: str70, str71, str72
45 CHARACTER(len=80):: str80
46 CHARACTER(len=20):: direction
47 INTEGER:: leng, ncid
48
49 !---------------------------------------------------------------------
50
51 pzsize = size(pzvalues)
52
53 ! 1.0 Verifications:
54 ! Do we have enough space for an extra axis ?
55 ! Is the name already in use ?
56
57 ! - Direction of axis. Can we get if from the user.
58 ! If not we put unknown.
59
60 IF (present(pdirect)) THEN
61 direction = trim(pdirect)
62 CALL strlowercase(direction)
63 ELSE
64 direction = 'unknown'
65 END IF
66
67 ! Check the consistency of the attribute
68
69 IF ((direction/='unknown') .AND. (direction/='up') .AND. &
70 (direction/='down')) THEN
71 direction = 'unknown'
72 str80 = 'The specified axis was: ' // trim(direction)
73 CALL histerr(2, 'histvert', &
74 'The specified direction for the vertical axis is not possible.', &
75 'it is replaced by: unknown', str80)
76 END IF
77
78 IF (nb_zax(pfileid)+1>nb_zax_max) THEN
79 CALL histerr(3, 'histvert', &
80 'Table of vertical axes too small. You should increase ', &
81 'nb_zax_max in M_HISTCOM.f90 in order to accomodate all ', &
82 'these variables ')
83 END IF
84
85 iv = nb_zax(pfileid)
86 IF (iv>1) THEN
87 str20 = pzaxname
88 nb = iv - 1
89 tab_str20(1:nb) = zax_name(pfileid, 1:nb)
90 tab_str20_length(1:nb) = zax_name_length(pfileid, 1:nb)
91 CALL find_str(nb, tab_str20, tab_str20_length, str20, pos)
92 ELSE
93 pos = 0
94 END IF
95
96 IF (pos>0) THEN
97 str70 = 'Vertical axis already exists'
98 WRITE(str71, '("Check variable ", a, " in file", I3)') str20, &
99 pfileid
100 str72 = 'Can also be a wrong file ID in another declaration'
101 CALL histerr(3, 'histvert', str70, str71, str72)
102 END IF
103
104 iv = nb_zax(pfileid) + 1
105
106 ! 2.0 Add the information to the file
107
108 ncid = ncdf_ids(pfileid)
109
110 leng = min(len_trim(pzaxname), 20)
111 call nf95_def_dim(ncid, pzaxname(1:leng), pzsize, zaxid_tmp)
112 call nf95_def_var(ncid, pzaxname(1:leng), nf90_float, zaxid_tmp, zdimid)
113
114 leng = min(len_trim(pzaxunit), 20)
115 call NF95_PUT_ATT(ncid, zdimid, 'units', pzaxunit(1:leng))
116 call NF95_PUT_ATT(ncid, zdimid, 'positive', trim(direction))
117
118 call NF95_PUT_ATT(ncid, zdimid, 'valid_min', real(minval( &
119 pzvalues(1:pzsize))))
120 call NF95_PUT_ATT(ncid, zdimid, 'valid_max', real(maxval( &
121 pzvalues(1:pzsize))))
122
123 leng = min(len_trim(pzaxname), 20)
124 call NF95_PUT_ATT(ncid, zdimid, 'title', pzaxname(1:leng))
125 leng = min(len_trim(pzaxtitle), 80)
126 call NF95_PUT_ATT(ncid, zdimid, 'long_name', pzaxtitle(1:leng))
127
128 call nf95_enddef(ncid)
129
130 call nf95_put_var(ncid, zdimid, pzvalues(1:pzsize))
131
132 call nf95_redef(ncid)
133
134 ! 3.0 add the information to the common
135
136 nb_zax(pfileid) = iv
137 zax_size(pfileid, iv) = pzsize
138 zax_name(pfileid, iv) = pzaxname
139 zax_name_length(pfileid, iv) = len_trim(pzaxname)
140 zax_ids(pfileid, iv) = zaxid_tmp
141 pzaxid = iv
142
143 END SUBROUTINE histvert
144
145 end module histvert_m

  ViewVC Help
Powered by ViewVC 1.1.21