1 | subroutine phyetat0 (ngrid,fichnom,tab0,Lmodif,nsoil,nq, & |
---|
2 | day_ini,time,tsurf,tsoil, & |
---|
3 | emis,q2,qsurf,cloudfrac,totcloudfrac,hice) |
---|
4 | |
---|
5 | USE infotrac, ONLY: tname |
---|
6 | USE surfdat_h, only: phisfi, albedodat, zmea, zstd, zsig, zgam, zthe |
---|
7 | use iostart, only: nid_start, open_startphy, close_startphy, & |
---|
8 | get_field, get_var, inquire_field, & |
---|
9 | inquire_dimension, inquire_dimension_length |
---|
10 | |
---|
11 | implicit none |
---|
12 | |
---|
13 | !====================================================================== |
---|
14 | ! Auteur(s) Z.X. Li (LMD/CNRS) date: 19930818 |
---|
15 | ! Adaptation à Mars : Yann Wanherdrick |
---|
16 | ! Objet: Lecture de l etat initial pour la physique |
---|
17 | !====================================================================== |
---|
18 | #include "netcdf.inc" |
---|
19 | #include "dimensions.h" |
---|
20 | #include "dimphys.h" |
---|
21 | #include "planete.h" |
---|
22 | #include "comcstfi.h" |
---|
23 | |
---|
24 | !====================================================================== |
---|
25 | ! INTEGER nbsrf !Mars nbsrf a 1 au lieu de 4 |
---|
26 | ! PARAMETER (nbsrf=1) ! nombre de sous-fractions pour une maille |
---|
27 | !====================================================================== |
---|
28 | ! Arguments: |
---|
29 | ! --------- |
---|
30 | ! inputs: |
---|
31 | integer,intent(in) :: ngrid |
---|
32 | character*(*),intent(in) :: fichnom ! "startfi.nc" file |
---|
33 | integer,intent(in) :: tab0 |
---|
34 | integer,intent(in) :: Lmodif |
---|
35 | integer,intent(in) :: nsoil ! # of soil layers |
---|
36 | integer,intent(in) :: nq |
---|
37 | integer,intent(in) :: day_ini |
---|
38 | real,intent(in) :: time |
---|
39 | |
---|
40 | ! outputs: |
---|
41 | real,intent(out) :: tsurf(ngrid) ! surface temperature |
---|
42 | real,intent(out) :: tsoil(ngrid,nsoil) ! soil temperature |
---|
43 | real,intent(out) :: emis(ngrid) ! surface emissivity |
---|
44 | real,intent(out) :: q2(ngrid, llm+1) ! |
---|
45 | real,intent(out) :: qsurf(ngrid,nq) ! tracers on surface |
---|
46 | ! real co2ice(ngrid) ! co2 ice cover |
---|
47 | real,intent(out) :: cloudfrac(ngrid,nlayermx) |
---|
48 | real,intent(out) :: hice(ngrid), totcloudfrac(ngrid) |
---|
49 | |
---|
50 | !====================================================================== |
---|
51 | ! Local variables: |
---|
52 | |
---|
53 | ! INTEGER radpas |
---|
54 | ! REAL co2_ppm |
---|
55 | ! REAL solaire |
---|
56 | |
---|
57 | real xmin,xmax ! to display min and max of a field |
---|
58 | ! |
---|
59 | INTEGER ig,iq,lmax |
---|
60 | INTEGER nid, nvarid |
---|
61 | INTEGER ierr, i, nsrf |
---|
62 | ! integer isoil |
---|
63 | ! INTEGER length |
---|
64 | ! PARAMETER (length=100) |
---|
65 | CHARACTER*7 str7 |
---|
66 | CHARACTER*2 str2 |
---|
67 | CHARACTER*1 yes |
---|
68 | ! |
---|
69 | REAL p_rad,p_omeg,p_g,p_cpp,p_mugaz,p_daysec |
---|
70 | INTEGER nqold |
---|
71 | |
---|
72 | ! flag which identifies if 'startfi.nc' file is using old names (qsurf01,...) |
---|
73 | ! logical :: oldtracernames=.false. |
---|
74 | integer :: count |
---|
75 | character(len=30) :: txt ! to store some text |
---|
76 | |
---|
77 | INTEGER :: indextime=1 ! index of selected time, default value=1 |
---|
78 | logical :: found |
---|
79 | |
---|
80 | !! added variables by AS to allow to read only slices of startfi |
---|
81 | ! real :: toto(ngrid) |
---|
82 | ! integer :: sta !! subscript in starti where we start to read |
---|
83 | ! integer, dimension(2) :: sta2d |
---|
84 | ! integer, dimension(2) :: siz2d |
---|
85 | |
---|
86 | ! AS: get the cursor that is stored in dimphys.h |
---|
87 | ! ... this allows to read only a part of startfi horiz grid |
---|
88 | !sta = cursor |
---|
89 | |
---|
90 | ! |
---|
91 | ! ALLOCATE ARRAYS IN surfdat_h |
---|
92 | ! |
---|
93 | IF (.not. ALLOCATED(albedodat)) ALLOCATE(albedodat(ngrid)) |
---|
94 | IF (.not. ALLOCATED(phisfi)) ALLOCATE(phisfi(ngrid)) |
---|
95 | IF (.not. ALLOCATED(zmea)) ALLOCATE(zmea(ngrid)) |
---|
96 | IF (.not. ALLOCATED(zstd)) ALLOCATE(zstd(ngrid)) |
---|
97 | IF (.not. ALLOCATED(zsig)) ALLOCATE(zsig(ngrid)) |
---|
98 | IF (.not. ALLOCATED(zgam)) ALLOCATE(zgam(ngrid)) |
---|
99 | IF (.not. ALLOCATED(zthe)) ALLOCATE(zthe(ngrid)) |
---|
100 | |
---|
101 | ! open physics initial state file: |
---|
102 | call open_startphy(fichnom) |
---|
103 | |
---|
104 | |
---|
105 | ! possibility to modify tab_cntrl in tabfi |
---|
106 | write(*,*) |
---|
107 | write(*,*) 'TABFI in phyeta0: Lmodif=',Lmodif," tab0=",tab0 |
---|
108 | call tabfi (ngrid,nid_start,Lmodif,tab0,day_ini,lmax,p_rad, & |
---|
109 | p_omeg,p_g,p_cpp,p_mugaz,p_daysec,time) |
---|
110 | |
---|
111 | !c |
---|
112 | !c Lecture des latitudes (coordonnees): |
---|
113 | !c |
---|
114 | ! ierr = NF_INQ_VARID (nid, "latitude", nvarid) |
---|
115 | ! IF (ierr.NE.NF_NOERR) THEN |
---|
116 | ! PRINT*, 'phyetat0: Le champ <latitude> est absent' |
---|
117 | ! CALL abort |
---|
118 | ! ENDIF |
---|
119 | !#ifdef NC_DOUBLE |
---|
120 | ! ierr = NF_GET_VARA_DOUBLE(nid,nvarid,sta,ngrid,lati) |
---|
121 | !#else |
---|
122 | ! ierr = NF_GET_VARA_REAL(nid,nvarid,sta,ngrid,lati) |
---|
123 | !#endif |
---|
124 | ! IF (ierr.NE.NF_NOERR) THEN |
---|
125 | ! PRINT*, 'phyetat0: Lecture echouee pour <latitude>' |
---|
126 | ! CALL abort |
---|
127 | ! ENDIF |
---|
128 | !c |
---|
129 | !c Lecture des longitudes (coordonnees): |
---|
130 | !c |
---|
131 | ! ierr = NF_INQ_VARID (nid, "longitude", nvarid) |
---|
132 | ! IF (ierr.NE.NF_NOERR) THEN |
---|
133 | ! PRINT*, 'phyetat0: Le champ <longitude> est absent' |
---|
134 | ! CALL abort |
---|
135 | ! ENDIF |
---|
136 | !#ifdef NC_DOUBLE |
---|
137 | ! ierr = NF_GET_VARA_DOUBLE(nid,nvarid,sta,ngrid,long) |
---|
138 | !#else |
---|
139 | ! ierr = NF_GET_VARA_REAL(nid,nvarid,sta,ngrid,long) |
---|
140 | !#endif |
---|
141 | ! IF (ierr.NE.NF_NOERR) THEN |
---|
142 | ! PRINT*, 'phyetat0: Lecture echouee pour <longitude>' |
---|
143 | ! CALL abort |
---|
144 | ! ENDIF |
---|
145 | !c |
---|
146 | !c Lecture des aires des mailles: |
---|
147 | !c |
---|
148 | ! ierr = NF_INQ_VARID (nid, "area", nvarid) |
---|
149 | ! IF (ierr.NE.NF_NOERR) THEN |
---|
150 | ! PRINT*, 'phyetat0: Le champ <area> est absent' |
---|
151 | ! CALL abort |
---|
152 | ! ENDIF |
---|
153 | !#ifdef NC_DOUBLE |
---|
154 | ! ierr = NF_GET_VARA_DOUBLE(nid,nvarid,sta,ngrid,area) |
---|
155 | !#else |
---|
156 | ! ierr = NF_GET_VARA_REAL(nid,nvarid,sta,ngrid,area) |
---|
157 | !#endif |
---|
158 | ! IF (ierr.NE.NF_NOERR) THEN |
---|
159 | ! PRINT*, 'phyetat0: Lecture echouee pour <area>' |
---|
160 | ! CALL abort |
---|
161 | ! ENDIF |
---|
162 | ! xmin = 1.0E+20 |
---|
163 | ! xmax = -1.0E+20 |
---|
164 | ! xmin = MINVAL(area) |
---|
165 | ! xmax = MAXVAL(area) |
---|
166 | ! PRINT*,'Aires des mailles <area>:', xmin, xmax |
---|
167 | |
---|
168 | ! Load surface geopotential: |
---|
169 | call get_field("phisfi",phisfi,found) |
---|
170 | if (.not.found) then |
---|
171 | write(*,*) "phyetat0: Failed loading <phisfi>" |
---|
172 | call abort |
---|
173 | else |
---|
174 | write(*,*) "phyetat0: surface geopotential <phisfi> range:", & |
---|
175 | minval(phisfi), maxval(phisfi) |
---|
176 | endif |
---|
177 | |
---|
178 | ! Load bare ground albedo: |
---|
179 | call get_field("albedodat",albedodat,found) |
---|
180 | if (.not.found) then |
---|
181 | write(*,*) "phyetat0: Failed loading <albedodat>" |
---|
182 | call abort |
---|
183 | else |
---|
184 | write(*,*) "phyetat0: Bare ground albedo <albedodat> range:", & |
---|
185 | minval(albedodat), maxval(albedodat) |
---|
186 | endif |
---|
187 | |
---|
188 | ! ZMEA |
---|
189 | call get_field("ZMEA",zmea,found) |
---|
190 | if (.not.found) then |
---|
191 | write(*,*) "phyetat0: Failed loading <ZMEA>" |
---|
192 | call abort |
---|
193 | else |
---|
194 | write(*,*) "phyetat0: <ZMEA> range:", & |
---|
195 | minval(zmea), maxval(zmea) |
---|
196 | endif |
---|
197 | |
---|
198 | ! ZSTD |
---|
199 | call get_field("ZSTD",zstd,found) |
---|
200 | if (.not.found) then |
---|
201 | write(*,*) "phyetat0: Failed loading <ZSTD>" |
---|
202 | call abort |
---|
203 | else |
---|
204 | write(*,*) "phyetat0: <ZSTD> range:", & |
---|
205 | minval(zstd), maxval(zstd) |
---|
206 | endif |
---|
207 | |
---|
208 | ! ZSIG |
---|
209 | call get_field("ZSIG",zsig,found) |
---|
210 | if (.not.found) then |
---|
211 | write(*,*) "phyetat0: Failed loading <ZSIG>" |
---|
212 | call abort |
---|
213 | else |
---|
214 | write(*,*) "phyetat0: <ZSIG> range:", & |
---|
215 | minval(zsig), maxval(zsig) |
---|
216 | endif |
---|
217 | |
---|
218 | ! ZGAM |
---|
219 | call get_field("ZGAM",zgam,found) |
---|
220 | if (.not.found) then |
---|
221 | write(*,*) "phyetat0: Failed loading <ZGAM>" |
---|
222 | call abort |
---|
223 | else |
---|
224 | write(*,*) "phyetat0: <ZGAM> range:", & |
---|
225 | minval(zgam), maxval(zgam) |
---|
226 | endif |
---|
227 | |
---|
228 | ! ZTHE |
---|
229 | call get_field("ZTHE",zthe,found) |
---|
230 | if (.not.found) then |
---|
231 | write(*,*) "phyetat0: Failed loading <ZTHE>" |
---|
232 | call abort |
---|
233 | else |
---|
234 | write(*,*) "phyetat0: <ZTHE> range:", & |
---|
235 | minval(zthe), maxval(zthe) |
---|
236 | endif |
---|
237 | |
---|
238 | ! Surface temperature : |
---|
239 | call get_field("tsurf",tsurf,found,indextime) |
---|
240 | if (.not.found) then |
---|
241 | write(*,*) "phyetat0: Failed loading <tsurf>" |
---|
242 | call abort |
---|
243 | else |
---|
244 | write(*,*) "phyetat0: Surface temperature <tsurf> range:", & |
---|
245 | minval(tsurf), maxval(tsurf) |
---|
246 | endif |
---|
247 | |
---|
248 | ! Surface emissivity |
---|
249 | call get_field("emis",emis,found,indextime) |
---|
250 | if (.not.found) then |
---|
251 | write(*,*) "phyetat0: Failed loading <emis>" |
---|
252 | call abort |
---|
253 | else |
---|
254 | write(*,*) "phyetat0: Surface emissivity <emis> range:", & |
---|
255 | minval(emis), maxval(emis) |
---|
256 | endif |
---|
257 | |
---|
258 | ! Cloud fraction (added by BC 2010) |
---|
259 | call get_field("cloudfrac",cloudfrac,found,indextime) |
---|
260 | if (.not.found) then |
---|
261 | write(*,*) "phyetat0: Failed loading <cloudfrac>" |
---|
262 | call abort |
---|
263 | else |
---|
264 | write(*,*) "phyetat0: Cloud fraction <cloudfrac> range:", & |
---|
265 | minval(cloudfrac), maxval(cloudfrac) |
---|
266 | endif |
---|
267 | |
---|
268 | ! Total cloud fraction (added by BC 2010) |
---|
269 | call get_field("totcloudfrac",totcloudfrac,found,indextime) |
---|
270 | if (.not.found) then |
---|
271 | write(*,*) "phyetat0: Failed loading <totcloudfrac>" |
---|
272 | call abort |
---|
273 | else |
---|
274 | write(*,*) "phyetat0: Total cloud fraction <totcloudfrac> range:", & |
---|
275 | minval(totcloudfrac), maxval(totcloudfrac) |
---|
276 | endif |
---|
277 | |
---|
278 | ! Height of oceanic ice (added by BC 2010) |
---|
279 | call get_field("hice",hice,found,indextime) |
---|
280 | if (.not.found) then |
---|
281 | write(*,*) "phyetat0: Failed loading <hice>" |
---|
282 | call abort |
---|
283 | else |
---|
284 | write(*,*) "phyetat0: Height of oceanic ice <hice> range:", & |
---|
285 | minval(hice), maxval(hice) |
---|
286 | endif |
---|
287 | |
---|
288 | ! pbl wind variance |
---|
289 | call get_field("q2",q2,found,indextime) |
---|
290 | if (.not.found) then |
---|
291 | write(*,*) "phyetat0: Failed loading <q2>" |
---|
292 | call abort |
---|
293 | else |
---|
294 | write(*,*) "phyetat0: PBL wind variance <q2> range:", & |
---|
295 | minval(q2), maxval(q2) |
---|
296 | endif |
---|
297 | |
---|
298 | ! tracer on surface |
---|
299 | if (nq.ge.1) then |
---|
300 | do iq=1,nq |
---|
301 | txt=tname(iq) |
---|
302 | if (txt.eq."h2o_vap") then |
---|
303 | ! There is no surface tracer for h2o_vap; |
---|
304 | ! "h2o_ice" should be loaded instead |
---|
305 | txt="h2o_ice" |
---|
306 | write(*,*) 'phyetat0: loading surface tracer', & |
---|
307 | ' h2o_ice instead of h2o_vap' |
---|
308 | endif |
---|
309 | call get_field(txt,qsurf(:,iq),found,indextime) |
---|
310 | if (.not.found) then |
---|
311 | write(*,*) "phyetat0: Failed loading <",trim(txt),">" |
---|
312 | write(*,*) " ",trim(txt)," is set to zero" |
---|
313 | else |
---|
314 | write(*,*) "phyetat0: Surface tracer <",trim(txt),"> range:", & |
---|
315 | minval(qsurf(:,iq)), maxval(qsurf(:,iq)) |
---|
316 | endif |
---|
317 | enddo |
---|
318 | endif ! of if (nq.ge.1) |
---|
319 | |
---|
320 | ! Call to soil_settings, in order to read soil temperatures, |
---|
321 | ! as well as thermal inertia and volumetric heat capacity |
---|
322 | |
---|
323 | call soil_settings(nid_start,ngrid,nsoil,tsurf,tsoil,indextime) |
---|
324 | |
---|
325 | ! |
---|
326 | ! close file: |
---|
327 | ! |
---|
328 | call close_startphy |
---|
329 | |
---|
330 | END SUBROUTINE phyetat0 |
---|