/[lmdze]/trunk/libf/IOIPSL/Flincom/flinfindcood.f90
ViewVC logotype

Annotation of /trunk/libf/IOIPSL/Flincom/flinfindcood.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 68 - (hide annotations)
Wed Nov 14 16:59:30 2012 UTC (11 years, 7 months ago) by guez
File size: 4179 byte(s)
Split "flincom.f90" into "flinclo.f90", "flinfindcood.f90",
"flininfo.f90" and "flinopen_nozoom.f90", in directory
"IOIPSL/Flincom".

Renamed "etat0_lim" to "ce0l", as in LMDZ.

Split "readsulfate.f" into "readsulfate.f90", "readsulfate_preind.f90"
and "getso4fromfile.f90".

In etat0, renamed variable q3d to q, as in "dynredem1". Replaced calls
to Flicom procedures by calls to NetCDF95.

In leapfrog, added call to writehist.

Extracted ASCII art from "grid_noro" into a file
"grid_noro.txt". Transformed explicit-shape local arrays into
automatic arrays, so that test on values of iim and jjm is no longer
needed. Test on weight:
          IF (weight(ii, jj) /= 0.) THEN
is useless. There is already a test before:
    if (any(weight == 0.)) stop "zero weight in grid_noro"

In "aeropt", replaced duplicated lines with different values of inu by
a loop on inu.

Removed arguments of "conf_phys". Corresponding variables are now
defined in "physiq", in a namelist. In "conf_phys", read a namelist
instead of using getin.

1 guez 68 module flinfindcood_m
2    
3     implicit none
4    
5     contains
6    
7     SUBROUTINE flinfindcood (fid_in, axtype, vid, ndim)
8    
9     ! This subroutine explores the file in order to find
10     ! the coordinate according to a number of rules
11    
12     USE errioipsl, ONLY: histerr
13     use flininfo_m, only: NCIDS, NCNBVA
14     USE netcdf, ONLY: nf90_get_att, nf90_inquire_dimension, &
15     nf90_inquire_variable, nf90_noerr
16     USE strlowercase_m, ONLY: strlowercase
17    
18     ! ARGUMENTS
19    
20     INTEGER, intent(in):: fid_in
21     integer vid, ndim
22     CHARACTER(LEN=3):: axtype
23    
24     ! LOCAL
25    
26     INTEGER:: iv, iret, dimnb
27     CHARACTER(LEN=40):: dimname, dimuni1, dimuni2, dimuni3
28     CHARACTER(LEN=30):: str1
29     LOGICAL:: found_rule = .FALSE.
30     !---------------------------------------------------------------------
31     vid = -1
32    
33     ! Make sure all strings are invalid
34    
35     dimname = '?-?'
36     dimuni1 = '?-?'
37     dimuni2 = '?-?'
38     dimuni3 = '?-?'
39    
40     ! First rule: we look for the correct units
41     ! lon: east
42     ! lat: north
43     ! We make an exact check as it would be too easy to mistake
44     ! some units by just comparing the substrings.
45    
46     SELECTCASE(axtype)
47     CASE ('lon')
48     dimuni1 = 'degree_e'
49     dimuni2 = 'degrees_e'
50     found_rule = .TRUE.
51     CASE('lat')
52     dimuni1 = 'degree_n'
53     dimuni2 = 'degrees_n'
54     found_rule = .TRUE.
55     CASE('lev')
56     dimuni1 = 'm'
57     dimuni2 = 'km'
58     dimuni3 = 'hpa'
59     found_rule = .TRUE.
60     CASE DEFAULT
61     found_rule = .FALSE.
62     END SELECT
63    
64     IF (found_rule) THEN
65     iv = 0
66     DO WHILE ((vid < 0).AND.(iv < ncnbva(fid_in)))
67     iv = iv+1
68     str1 = ''
69     iret = NF90_GET_ATT (ncids(fid_in), iv, 'units', str1)
70     IF (iret == NF90_NOERR) THEN
71     CALL strlowercase (str1)
72     IF ((INDEX(str1, TRIM(dimuni1)) == 1) &
73     .OR.(INDEX(str1, TRIM(dimuni2)) == 1) &
74     .OR.(INDEX(str1, TRIM(dimuni3)) == 1)) THEN
75     vid = iv
76     iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, ndims=ndim)
77     ENDIF
78     ENDIF
79     ENDDO
80     ENDIF
81    
82     ! Second rule: we find specific names:
83     ! lon: nav_lon
84     ! lat: nav_lat
85     ! Here we can check if we find the substring as the
86     ! names are more specific.
87    
88     SELECTCASE(axtype)
89     CASE ('lon')
90     dimname = 'nav_lon lon longitude'
91     found_rule = .TRUE.
92     CASE('lat')
93     dimname = 'nav_lat lat latitude'
94     found_rule = .TRUE.
95     CASE('lev')
96     dimname = 'plev level depth deptht'
97     found_rule = .TRUE.
98     CASE DEFAULT
99     found_rule = .FALSE.
100     END SELECT
101    
102     IF (found_rule) THEN
103     iv = 0
104     DO WHILE ((vid < 0).AND.(iv < ncnbva(fid_in)))
105     iv = iv+1
106     str1=''
107     iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
108     name=str1, ndims=ndim)
109     IF (INDEX(dimname, TRIM(str1)) >= 1) THEN
110     vid = iv
111     ENDIF
112     ENDDO
113     ENDIF
114    
115     ! Third rule: we find a variable with the same name as the dimension
116     ! lon = 1
117     ! lat = 2
118     ! lev = 3
119    
120     IF (vid < 0) THEN
121     SELECTCASE(axtype)
122     CASE ('lon')
123     dimnb = 1
124     found_rule = .TRUE.
125     CASE('lat')
126     dimnb = 2
127     found_rule = .TRUE.
128     CASE('lev')
129     dimnb = 3
130     found_rule = .TRUE.
131     CASE DEFAULT
132     found_rule = .FALSE.
133     END SELECT
134    
135     IF (found_rule) THEN
136     iret = NF90_INQUIRE_DIMENSION (ncids(fid_in), dimnb, name=dimname)
137     iv = 0
138     DO WHILE ((vid < 0).AND.(iv < ncnbva(fid_in)))
139     iv = iv+1
140     str1=''
141     iret = NF90_INQUIRE_VARIABLE (ncids(fid_in), iv, &
142     name=str1, ndims=ndim)
143     IF (INDEX(dimname, TRIM(str1)) == 1) THEN
144     vid = iv
145     ENDIF
146     ENDDO
147     ENDIF
148     ENDIF
149    
150     ! Stop the program if no coordinate was found
151    
152     IF (vid < 0) THEN
153     CALL histerr (3, 'flinfindcood', &
154     'No coordinate axis was found in the file', &
155     'The data in this file can not be used', axtype)
156     ENDIF
157    
158     END SUBROUTINE flinfindcood
159    
160     end module flinfindcood_m

  ViewVC Help
Powered by ViewVC 1.1.21