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

Diff of /trunk/libf/IOIPSL/Histcom/histend.f90

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 61 by guez, Fri Apr 20 14:58:43 2012 UTC revision 62 by guez, Thu Jul 26 14:37:37 2012 UTC
# Line 4  module histend_m Line 4  module histend_m
4    
5  contains  contains
6    
7    SUBROUTINE histend(pfileid)    SUBROUTINE histend(fileid)
8    
9      ! This subroutine ends the declaration of variables and sets the      ! This subroutine ends the declaration of variables, sets the time
10      ! time axes in the netcdf file and puts it into write mode.      ! axes in the NetCDF file and puts it into write mode.
   
     ! INPUT  
   
     ! pfileid: ID of the file to be worked on  
11    
12      USE ioipslmpp, ONLY: ioipslmpp_addatt      USE ioipslmpp, ONLY: ioipslmpp_addatt
13      USE errioipsl, ONLY: histerr      USE errioipsl, ONLY: histerr
# Line 19  contains Line 15  contains
15           missing_val, name, nb_tax, nb_var, ncdf_ids, ncvar_ids, regular, &           missing_val, name, nb_tax, nb_var, ncdf_ids, ncvar_ids, regular, &
16           tax_name, tdimid, tid, title, topp, unit_name, var_axid, var_zaxid, &           tax_name, tdimid, tid, title, topp, unit_name, var_axid, var_zaxid, &
17           xid, yid, zax_ids, zax_name           xid, yid, zax_ids, zax_name
18      USE calendar, ONLY: ioget_calendar, ju2ymds      USE ioget_calendar_m, ONLY: ioget_calendar
19      USE netcdf, ONLY: nf90_def_dim, nf90_def_var, nf90_enddef, &      USE calendar, ONLY: ju2ymds
20           nf90_float, nf90_put_att, nf90_unlimited      USE netcdf, ONLY: nf90_float, nf90_unlimited
21        use netcdf95, only: nf95_def_dim, nf95_def_var, nf95_put_att, nf95_enddef
22      INTEGER, INTENT (IN):: pfileid  
23        INTEGER, INTENT(IN):: fileid ! ID of the file to be worked on
24      INTEGER:: ncid, ncvarid  
25      INTEGER:: iret, ndim, iv, itx, ziv      ! Local:
26      INTEGER:: itax      INTEGER ncid, varid
27      INTEGER:: dims(4), dim_cnt      INTEGER iret, ndim, iv, itx, ziv
28      INTEGER:: year, month, day, hours, minutes      INTEGER itax
29      REAL:: sec      INTEGER dims(4), dim_cnt
30      REAL:: rtime0      INTEGER year, month, day, hours, minutes
31      CHARACTER (len=20):: tname, tunit      REAL sec
32      CHARACTER (len=30):: str30      REAL rtime0
33      CHARACTER (len=80):: ttitle      CHARACTER(len=20) tname, tunit
34      CHARACTER (len=120):: assoc      CHARACTER(len=30) str30
35      CHARACTER (len=70):: str70      CHARACTER(len=80) ttitle
36      CHARACTER (len=3), DIMENSION (12):: cal = (/ 'JAN', 'FEB', 'MAR', &      CHARACTER(len=120) assoc
37           'APR', 'MAY', 'JUN', 'JUL', 'AUG', 'SEP', 'OCT', 'NOV', 'DEC'/)      CHARACTER(len=70) str70
38      CHARACTER (len=7):: tmp_opp      CHARACTER(len=3):: cal(12) = (/ 'JAN', 'FEB', 'MAR', 'APR', 'MAY', 'JUN', &
39             'JUL', 'AUG', 'SEP', 'OCT', 'NOV', 'DEC'/)
40        CHARACTER(len=7) tmp_opp
41    
42      !---------------------------------------------------------------------      !---------------------------------------------------------------------
43      ncid = ncdf_ids(pfileid)  
44        ncid = ncdf_ids(fileid)
45    
46      ! 1.0 Create the time axes      ! 1.0 Create the time axes
47    
48      iret = nf90_def_dim(ncid, 'time_counter', nf90_unlimited, tid(pfileid))      call nf95_def_dim(ncid, 'time_counter', nf90_unlimited, tid(fileid))
49    
50      ! 1.1 Define all the time axes needed for this file      ! 1.1 Define all the time axes needed for this file
51    
52      DO itx = 1, nb_tax(pfileid)      DO itx = 1, nb_tax(fileid)
53         dims(1) = tid(pfileid)         IF (nb_tax(fileid)>1) THEN
54         IF (nb_tax(pfileid)>1) THEN            str30 = 't_' // tax_name(fileid, itx)
           str30 = 't_' // tax_name(pfileid, itx)  
55         ELSE         ELSE
56            str30 = 'time_counter'            str30 = 'time_counter'
57         END IF         END IF
58         iret = nf90_def_var(ncid, str30, nf90_float, dims(1), &         call nf95_def_var(ncid, str30, nf90_float, tid(fileid), &
59              tdimid(pfileid, itx))              tdimid(fileid, itx))
   
        !   To transform the current itau into a real date and take it  
        !   as the origin of the file requires the time counter to change.  
        !   Thus it is an operation the user has to ask for.  
        !   This function should thus only be re-instated  
        !   if there is a ioconf routine to control it.  
   
        ! rtime0 = itau2date(itau0(pfileid), date0(pfileid), deltat(pfileid))  
        rtime0 = date0(pfileid)  
60    
61           rtime0 = date0(fileid)
62         CALL ju2ymds(rtime0, year, month, day, sec)         CALL ju2ymds(rtime0, year, month, day, sec)
63    
64         !   Catch any error induced by a change in calendar !         ! Catch any error induced by a change in calendar
65    
66         IF (year<0) THEN         IF (year <  0) THEN
67            year = 2000 + year            year = 2000 + year
68         END IF         END IF
69    
# Line 81  contains Line 71  contains
71         minutes = int((sec-hours*60.*60.)/60.)         minutes = int((sec-hours*60.*60.)/60.)
72         sec = sec - (hours*60.*60.+minutes*60.)         sec = sec - (hours*60.*60.+minutes*60.)
73    
74         WRITE (str70, 7000) year, month, day, hours, minutes, int(sec)         WRITE(str70, 7000) year, month, day, hours, minutes, int(sec)
75         iret = nf90_put_att(ncid, tdimid(pfileid, itx), 'units', trim(str70))         call nf95_put_att(ncid, tdimid(fileid, itx), 'units', trim(str70))
76    
77         CALL ioget_calendar(str30)         CALL ioget_calendar(str30)
78         iret = nf90_put_att(ncid, tdimid(pfileid, itx), 'calendar', trim(str30))         call nf95_put_att(ncid, tdimid(fileid, itx), 'calendar', trim(str30))
79    
80         iret = nf90_put_att(ncid, tdimid(pfileid, itx), 'title', 'Time')         call nf95_put_att(ncid, tdimid(fileid, itx), 'title', 'Time')
81    
82         iret = nf90_put_att(ncid, tdimid(pfileid, itx), 'long_name', &         call nf95_put_att(ncid, tdimid(fileid, itx), 'long_name', &
83              'Time axis')              'Time axis')
84    
85         WRITE (str70, 7001) year, cal(month), day, hours, minutes, int(sec)         WRITE(str70, 7001) year, cal(month), day, hours, minutes, int(sec)
86         iret = nf90_put_att(ncid, tdimid(pfileid, itx), 'time_origin', &         call nf95_put_att(ncid, tdimid(fileid, itx), 'time_origin', &
87              trim(str70))              trim(str70))
88      END DO      END DO
89    
     ! The formats we need  
   
 7000 FORMAT ('seconds since ', I4.4, '-', I2.2, '-', I2.2, ' ', I2.2, ':', I2.2, ':', &  
          I2.2)  
 7001 FORMAT (' ', I4.4, '-', A3, '-', I2.2, ' ', I2.2, ':', I2.2, ':', I2.2)  
   
90      ! 2.0 declare the variables      ! 2.0 declare the variables
91    
92      DO iv = 1, nb_var(pfileid)      DO iv = 1, nb_var(fileid)
93    
94         itax = var_axid(pfileid, iv)         itax = var_axid(fileid, iv)
95    
96         tname = name(pfileid, iv)         tname = name(fileid, iv)
97         tunit = unit_name(pfileid, iv)         tunit = unit_name(fileid, iv)
98         ttitle = title(pfileid, iv)         ttitle = title(fileid, iv)
99    
100         IF (regular(pfileid)) THEN         IF (regular(fileid)) THEN
101            dims(1:2) = (/ xid(pfileid), yid(pfileid) /)            dims(1:2) = (/ xid(fileid), yid(fileid) /)
102            dim_cnt = 2            dim_cnt = 2
103         ELSE         ELSE
104            dims(1) = xid(pfileid)            dims(1) = xid(fileid)
105            dim_cnt = 1            dim_cnt = 1
106         END IF         END IF
107    
108         tmp_opp = topp(pfileid, iv)         tmp_opp = topp(fileid, iv)
109         ziv = var_zaxid(pfileid, iv)         ziv = var_zaxid(fileid, iv)
110    
111         !   2.1 dimension of field         ! 2.1 dimension of field
112    
113         IF ((trim(tmp_opp)/='never')) THEN         IF ((trim(tmp_opp)/='never')) THEN
114            IF ((trim(tmp_opp)/='once') .AND. (trim( &            IF ((trim(tmp_opp)/='once') .AND. (trim( &
115                 tmp_opp)/='l_max') .AND. (trim(tmp_opp)/='l_min')) THEN                 tmp_opp)/='l_max') .AND. (trim(tmp_opp)/='l_min')) THEN
116               IF (ziv==-99) THEN               IF (ziv==-99) THEN
117                  ndim = dim_cnt + 1                  ndim = dim_cnt + 1
118                  dims(dim_cnt+1:dim_cnt+2) = (/ tid(pfileid), 0 /)                  dims(dim_cnt+1:dim_cnt+2) = (/ tid(fileid), 0 /)
119               ELSE               ELSE
120                  ndim = dim_cnt + 2                  ndim = dim_cnt + 2
121                  dims(dim_cnt+1:dim_cnt+2) = (/ zax_ids(pfileid, ziv), &                  dims(dim_cnt+1:dim_cnt+2) = (/ zax_ids(fileid, ziv), &
122                       tid(pfileid) /)                       tid(fileid) /)
123               END IF               END IF
124            ELSE            ELSE
125               IF (ziv==-99) THEN               IF (ziv==-99) THEN
# Line 143  contains Line 127  contains
127                  dims(dim_cnt+1:dim_cnt+2) = (/ 0, 0 /)                  dims(dim_cnt+1:dim_cnt+2) = (/ 0, 0 /)
128               ELSE               ELSE
129                  ndim = dim_cnt + 1                  ndim = dim_cnt + 1
130                  dims(dim_cnt+1:dim_cnt+2) = (/ zax_ids(pfileid, ziv), 0 /)                  dims(dim_cnt+1:dim_cnt+2) = (/ zax_ids(fileid, ziv), 0 /)
131               END IF               END IF
132            END IF            END IF
133    
134            iret = nf90_def_var(ncid, trim(tname), nf90_float, dims(1:abs(ndim)), &            call nf95_def_var(ncid, trim(tname), nf90_float, dims(1:abs(ndim)), &
135                 ncvarid)                 varid)
   
           ncvar_ids(pfileid, iv) = ncvarid  
   
           iret = nf90_put_att(ncid, ncvarid, 'units', trim(tunit))  
   
           iret = nf90_put_att(ncid, ncvarid, 'missing_value', &  
                real(missing_val))  
           iret = nf90_put_att(ncid, ncvarid, 'long_name', trim(ttitle))  
136    
137            iret = nf90_put_att(ncid, ncvarid, 'short_name', trim(tname))            ncvar_ids(fileid, iv) = varid
138    
139            iret = nf90_put_att(ncid, ncvarid, 'online_operation', trim(fullop( &            call nf95_put_att(ncid, varid, 'units', trim(tunit))
140                 pfileid, iv)))            call nf95_put_att(ncid, varid, 'missing_value', missing_val)
141              call nf95_put_att(ncid, varid, 'long_name', trim(ttitle))
142              call nf95_put_att(ncid, varid, 'short_name', trim(tname))
143              call nf95_put_att(ncid, varid, 'online_operation', trim(fullop( &
144                   fileid, iv)))
145    
146            SELECT CASE (ndim)            SELECT CASE (ndim)
147            CASE (-3)            CASE (-3)
# Line 178  contains Line 158  contains
158                    'allowed at this stage', ' ')                    'allowed at this stage', ' ')
159            END SELECT            END SELECT
160    
161            iret = nf90_put_att(ncid, ncvarid, 'axis', trim(str30))            call nf95_put_att(ncid, varid, 'axis', trim(str30))
162    
163            assoc = 'nav_lat nav_lon'            assoc = 'nav_lat nav_lon'
164            ziv = var_zaxid(pfileid, iv)            ziv = var_zaxid(fileid, iv)
165            IF (ziv>0) THEN            IF (ziv>0) THEN
166               str30 = zax_name(pfileid, ziv)               str30 = zax_name(fileid, ziv)
167               assoc = trim(str30) // ' ' // trim(assoc)               assoc = trim(str30) // ' ' // trim(assoc)
168            END IF            END IF
169    
170            IF (itax>0) THEN            IF (itax>0) THEN
171               IF (nb_tax(pfileid)>1) THEN               IF (nb_tax(fileid)>1) THEN
172                  str30 = 't_' // tax_name(pfileid, itax)                  str30 = 't_' // tax_name(fileid, itax)
173               ELSE               ELSE
174                  str30 = 'time_counter'                  str30 = 'time_counter'
175               END IF               END IF
176               assoc = trim(str30) // ' ' // trim(assoc)               assoc = trim(str30) // ' ' // trim(assoc)
177    
178               iret = nf90_put_att(ncid, ncvarid, 'interval_operation', &               call nf95_put_att(ncid, varid, 'interval_operation', &
179                    real(freq_opp(pfileid, iv)))                    real(freq_opp(fileid, iv)))
180               iret = nf90_put_att(ncid, ncvarid, 'interval_write', real(freq_wrt( &               call nf95_put_att(ncid, varid, 'interval_write', &
181                    pfileid, iv)))                    real(freq_wrt(fileid, iv)))
182            END IF            END IF
183            iret = nf90_put_att(ncid, ncvarid, 'associate', trim(assoc))            call nf95_put_att(ncid, varid, 'associate', trim(assoc))
184         END IF         END IF
185      END DO      END DO
186    
187      !  Add MPP attributes      ! Add MPP attributes
   
188      CALL ioipslmpp_addatt(ncid)      CALL ioipslmpp_addatt(ncid)
189    
190      ! 3.0 Put the netcdf file into wrte mode      ! 3.0 Put the netcdf file into write mode
191        call nf95_enddef(ncid)
     iret = nf90_enddef(ncid)  
192    
193    7000 FORMAT ('seconds since ', I4.4, '-', I2.2, '-', I2.2, ' ', I2.2, ':', &
194              I2.2, ':', I2.2)
195    7001 FORMAT (' ', I4.4, '-', A3, '-', I2.2, ' ', I2.2, ':', I2.2, ':', I2.2)
196        
197    END SUBROUTINE histend    END SUBROUTINE histend
198    
199  end module histend_m  end module histend_m

Legend:
Removed from v.61  
changed lines
  Added in v.62

  ViewVC Help
Powered by ViewVC 1.1.21