/[lmdze]/trunk/libf/IOIPSL/flincom.f90
ViewVC logotype

Diff of /trunk/libf/IOIPSL/flincom.f90

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

revision 36 by guez, Thu Dec 2 17:11:04 2010 UTC revision 42 by guez, Thu Mar 24 11:52:41 2011 UTC
# Line 15  MODULE flincom Line 15  MODULE flincom
15    
16  CONTAINS  CONTAINS
17    
18    SUBROUTINE flinopen_nozoom(iim, jjm, llm, lon, lat, lev, &    SUBROUTINE flinopen_nozoom(iim, jjm, llm, lon, lat, lev, ttm, itaus, date0, &
19         ttm, itaus, date0, dt, fid_out)         dt, fid_out)
20    
21      ! The routine will open an input file      ! This procedure opens an input file. There is no test of the
22      ! INPUT      ! content of the file against the input from the model.
     ! There is no test of the content of the file against the input  
     ! from the model  
   
     ! iim: size in the x direction in the file (longitude)  
     ! jjm: size in the y direction  
     ! llm: number of levels  
     ! (llm = 0 means no axis to be expected)  
   
     ! WARNING:  
     ! It is for the user to check  
     ! that the dimensions of lon lat and lev are correct when passed to  
     ! flinopen. This can be done after the call when iim and jjm have  
     ! been retrieved from the netCDF file. In F90 this problem will  
     ! be solved with an internal assign  
     ! IF iim, jjm, llm or ttm are parameters in the calling program  
     ! it will create a segmentation fault  
   
     ! ttm: size of time axis  
   
     ! OUTPUT  
   
     ! lon: array of (iim, jjm),  
     ! that contains the longitude of each point  
     ! lat: same for latitude  
     ! lev: An array of llm for the latitude  
     ! itaus: Time steps within this file  
     ! date0: Julian date at which itau = 0  
     ! dt: length of the time steps of the data  
23    
24      !---------------------------------------------------------------------      ! The user should check that the dimensions of lon lat and lev are
25        ! correct when passed to flinopen. This can be done after the call
26        ! when iim and jjm have been retrieved from the netCDF file. In
27        ! Fortran 90 this problem will be solved with an internal assign
28        ! IF iim, jjm, llm or ttm are parameters in the calling program it
29        ! will create a segmentation fault
30    
31      USE calendar, ONLY: ymds2ju, ioconf_calendar      USE calendar, ONLY: ymds2ju, ioconf_calendar
32      USE errioipsl, ONLY: histerr      USE errioipsl, ONLY: histerr
33      USE netcdf, ONLY: nf90_get_att, nf90_get_var, nf90_global, &      USE netcdf, ONLY: nf90_get_att, nf90_get_var, nf90_global, &
34           nf90_inquire_variable           nf90_inquire_variable
35    
36      ! ARGUMENTS      INTEGER, intent(in):: iim ! size in the x direction in the file (longitude)
37        INTEGER, intent(in):: jjm ! size in the y direction
38    
39      INTEGER, intent(in):: iim, jjm, llm, ttm      INTEGER, intent(in):: llm ! number of levels
40      real, intent(out):: lon(iim, jjm), lat(iim, jjm), lev(llm)      ! (llm = 0 means no axis to be expected)
41      INTEGER, intent(out):: itaus(ttm)  
42      REAL, intent(out):: date0, dt      INTEGER, intent(in):: ttm ! size of time axis
43        real, intent(out):: lon(iim, jjm) ! longitude
44        real, intent(out):: lat(iim, jjm) ! latitude
45        real, intent(out):: lev(llm)
46        INTEGER, intent(out):: itaus(ttm) ! time steps within this file
47        REAL, intent(out):: date0 ! Julian date at which itau = 0
48        REAL, intent(out):: dt ! length of the time steps of the data
49    
50      INTEGER, intent(in):: fid_out      INTEGER, intent(in):: fid_out
51      ! (file ID which is later used to read the data)      ! (file ID which is later used to read the data)
52    
53      ! LOCAL      ! Variables local to the procedure:
54        INTEGER iret, vid, fid, nbdim, i
55      INTEGER:: iret, vid, fid, nbdim, i      INTEGER gdtt_id, old_id, iv, gdtmaf_id
56      INTEGER:: gdtt_id, old_id, iv, gdtmaf_id      CHARACTER(LEN=250) name
57      CHARACTER(LEN=250):: name      CHARACTER(LEN=80) units, my_calendar
58      CHARACTER(LEN=80):: units, my_calendar      INTEGER year, month, day
59      INTEGER:: year, month, day      REAL r_year, r_month, r_day
60      REAL:: r_year, r_month, r_day      INTEGER year0, month0, day0, hours0, minutes0, seci
61      INTEGER:: year0, month0, day0, hours0, minutes0, seci      REAL sec, sec0
62      REAL:: sec, sec0      CHARACTER strc
     CHARACTER:: strc  
   
63      REAL, DIMENSION(:), ALLOCATABLE:: vec_tmp      REAL, DIMENSION(:), ALLOCATABLE:: vec_tmp
64    
65      !---------------------------------------------------------------------      !---------------------------------------------------------------------
66    
67      IF ( (fid_out < 1).OR.(fid_out > nbfile_max) ) THEN      IF ((fid_out < 1) .OR. (fid_out > nbfile_max)) THEN
68         ! Either the fid_out has not been initialized (0 or very large)         ! Either the fid_out has not been initialized (0 or very large)
69         ! then we have to open anyway. Else we only need to open the file         ! then we have to open anyway. Else we only need to open the file
70         ! if it has not been opened before.         ! if it has not been opened before.
71         print *, "Call flinfo before flinopen"         print *, "Call flinfo before flinopen"
72         stop 1         stop 1
73      end IF      end IF
74      IF (.NOT.ncfileopen(fid_out)) THEN      IF (.NOT. ncfileopen(fid_out)) THEN
75         print *, "Call flinfo before flinopen"         print *, "Call flinfo before flinopen"
76         stop 1         stop 1
77      end IF      end IF
# Line 99  CONTAINS Line 81  CONTAINS
81    
82      fid = ncids(fid_out)      fid = ncids(fid_out)
83    
     ! 2.0 get the sizes and names of the different coordinates  
     ! and do a first set of verification.  
   
     ! 3.0 Check if we are realy talking about the same coodinate system  
     ! if not then we get the lon, lat and lev variables from the file  
   
84      ! 4.0 extracting the coordinates      ! 4.0 extracting the coordinates
85    
86      CALL flinfindcood (fid_out, 'lon', vid, nbdim)      CALL flinfindcood (fid_out, 'lon', vid, nbdim)
# Line 123  CONTAINS Line 99  CONTAINS
99    
100      CALL flinfindcood (fid_out, 'lat', vid, nbdim)      CALL flinfindcood (fid_out, 'lat', vid, nbdim)
101      IF (nbdim == 2) THEN      IF (nbdim == 2) THEN
102         iret = NF90_GET_VAR (fid, vid, lat, &         iret = NF90_GET_VAR (fid, vid, lat, start=(/ 1, 1 /), &
103              start=(/ 1, 1 /), count=(/ iim, jjm /))              count=(/ iim, jjm /))
104      ELSE      ELSE
105         ALLOCATE(vec_tmp(jjm))         ALLOCATE(vec_tmp(jjm))
106         iret = NF90_GET_VAR (fid, vid, vec_tmp, &         iret = NF90_GET_VAR (fid, vid, vec_tmp, start=(/ 1 /), count=(/ jjm /))
             start=(/ 1 /), count=(/ jjm /))  
107         DO i=1, iim         DO i=1, iim
108            lat(i, :) = vec_tmp(:)            lat(i, :) = vec_tmp(:)
109         ENDDO         ENDDO
# Line 138  CONTAINS Line 113  CONTAINS
113      IF (llm > 0) THEN      IF (llm > 0) THEN
114         CALL flinfindcood (fid_out, 'lev', vid, nbdim)         CALL flinfindcood (fid_out, 'lev', vid, nbdim)
115         IF (nbdim == 1) THEN         IF (nbdim == 1) THEN
116            iret = NF90_GET_VAR (fid, vid, lev, &            iret = NF90_GET_VAR (fid, vid, lev, start=(/ 1 /), count=(/ llm /))
                start=(/ 1 /), count=(/ llm /))  
117         ELSE         ELSE
118            CALL histerr (3, 'flinopen', &            CALL histerr (3, 'flinopen', &
119                 'Can not handle vertical coordinates that have more', &                 'Can not handle vertical coordinates that have more', &
# Line 187  CONTAINS Line 161  CONTAINS
161         ! Find the calendar         ! Find the calendar
162         my_calendar='XXXX'         my_calendar='XXXX'
163         iret = NF90_GET_ATT (fid, gdtmaf_id, 'calendar', my_calendar)         iret = NF90_GET_ATT (fid, gdtmaf_id, 'calendar', my_calendar)
164         IF ( INDEX(my_calendar, 'XXXX') < 1 ) THEN         IF (INDEX(my_calendar, 'XXXX') < 1) THEN
165            CALL ioconf_calendar(my_calendar)            CALL ioconf_calendar(my_calendar)
166         ENDIF         ENDIF
167    
# Line 231  CONTAINS Line 205  CONTAINS
205    SUBROUTINE flininfo(filename, iim, jjm, llm, ttm, fid_out)    SUBROUTINE flininfo(filename, iim, jjm, llm, ttm, fid_out)
206    
207      ! This subroutine allows to get some information.      ! This subroutine allows to get some information.
208      ! It is usualy done within flinopen but the user may want to call      ! It is usually called by "flinopen" but the user may want to call
209      ! it before in order to allocate the space needed to extract the      ! it before in order to allocate the space needed to extract the
210      ! data from the file.      ! data from the file.
211    
# Line 244  CONTAINS Line 218  CONTAINS
218      CHARACTER(LEN=*), intent(in):: filename      CHARACTER(LEN=*), intent(in):: filename
219      INTEGER, intent(out):: iim, jjm, llm, ttm, fid_out      INTEGER, intent(out):: iim, jjm, llm, ttm, fid_out
220    
221      ! LOCAL      ! Variables local to the procedure:
   
222      INTEGER, SAVE:: nbfiles = 0      INTEGER, SAVE:: nbfiles = 0
223      INTEGER, SAVE:: ncdims(nbfile_max, 4)      INTEGER, SAVE:: ncdims(nbfile_max, 4)
224      INTEGER:: iret, fid, ndims, nvars, nb_atts, id_unlim      INTEGER iret, fid, ndims, nvars, nb_atts, id_unlim
225      INTEGER:: iv, lll      INTEGER iv, lll
226      CHARACTER(LEN=80):: name      CHARACTER(LEN=80) name
227      CHARACTER(LEN=30):: axname      CHARACTER(LEN=30) axname
228    
229      !---------------------------------------------------------------------      !---------------------------------------------------------------------
230    
# Line 373  CONTAINS Line 346  CONTAINS
346    
347      IF (found_rule) THEN      IF (found_rule) THEN
348         iv = 0         iv = 0
349         DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )         DO WHILE ((vid < 0).AND.(iv < ncnbva(fid_in)))
350            iv = iv+1            iv = iv+1
351            str1 = ''            str1 = ''
352            iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1)            iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1)
353            IF (iret == NF90_NOERR) THEN            IF (iret == NF90_NOERR) THEN
354               CALL strlowercase (str1)               CALL strlowercase (str1)
355               IF ( (INDEX(str1, TRIM(dimuni1)) == 1) &               IF ((INDEX(str1, TRIM(dimuni1)) == 1) &
356                    .OR.(INDEX(str1, TRIM(dimuni2)) == 1) &                    .OR.(INDEX(str1, TRIM(dimuni2)) == 1) &
357                    .OR.(INDEX(str1, TRIM(dimuni3)) == 1) ) THEN                    .OR.(INDEX(str1, TRIM(dimuni3)) == 1)) THEN
358                  vid = iv                  vid = iv
359                  iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim)                  iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim)
360               ENDIF               ENDIF
# Line 411  CONTAINS Line 384  CONTAINS
384    
385      IF (found_rule) THEN      IF (found_rule) THEN
386         iv = 0         iv = 0
387         DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )         DO WHILE ((vid < 0).AND.(iv < ncnbva(fid_in)))
388            iv = iv+1            iv = iv+1
389            str1=''            str1=''
390            iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &            iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
# Line 445  CONTAINS Line 418  CONTAINS
418         IF (found_rule) THEN         IF (found_rule) THEN
419            iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname)            iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname)
420            iv = 0            iv = 0
421            DO WHILE ( (vid < 0).AND.(iv < ncnbva(fid_in)) )            DO WHILE ((vid < 0).AND.(iv < ncnbva(fid_in)))
422               iv = iv+1               iv = iv+1
423               str1=''               str1=''
424               iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &               iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &

Legend:
Removed from v.36  
changed lines
  Added in v.42

  ViewVC Help
Powered by ViewVC 1.1.21