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

Annotation of /trunk/libf/IOIPSL/Histcom/histvert.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (hide annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years, 2 months ago) by guez
File size: 4832 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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

  ViewVC Help
Powered by ViewVC 1.1.21