1 | subroutine soil_settings(nid,ngrid,nsoil,nqsoil,tsurf,tsoil, |
---|
2 | & qsoil,indextime) |
---|
3 | |
---|
4 | ! use netcdf |
---|
5 | use comsoil_h, only: layer, mlayer, inertiedat, inertiesoil, |
---|
6 | & volcapa,flux_geo,adsorption_soil, |
---|
7 | & igcm_h2o_vap_soil,igcm_h2o_ice_soil, |
---|
8 | & igcm_h2o_vap_ads |
---|
9 | use iostart, only: inquire_field_ndims, get_var, get_field, |
---|
10 | & inquire_field, inquire_dimension_length |
---|
11 | USE comslope_mod, ONLY: nslope |
---|
12 | implicit none |
---|
13 | |
---|
14 | !====================================================================== |
---|
15 | ! Author: Ehouarn Millour (07/2006) |
---|
16 | ! |
---|
17 | ! Purpose: Read and/or initialise soil depths and properties |
---|
18 | ! |
---|
19 | ! Modifications: Aug.2010 EM : use NetCDF90 to load variables (enables using |
---|
20 | ! r4 or r8 restarts independently of having compiled |
---|
21 | ! the GCM in r4 or r8) |
---|
22 | ! June 2013 TN : Possibility to read files with a time axis |
---|
23 | ! |
---|
24 | ! |
---|
25 | ! This subroutine reads from a NetCDF file (opened by the caller) |
---|
26 | ! of "startfi.nc" format. |
---|
27 | ! The various actions and variable read/initialized are: |
---|
28 | ! 1. Check out the number of soil layers (in datafile); if it isn't equal |
---|
29 | ! to nsoil, then some interpolation will be required |
---|
30 | ! Also check if data in file "startfi.nc" is in older format (ie: |
---|
31 | ! thermal inertia was depth-independent; and there was no "depth" |
---|
32 | ! coordinate. |
---|
33 | ! Read/build layer (and midlayer) depths |
---|
34 | ! 2. Read volumetric specific heat (or initialise it to default value) |
---|
35 | ! 3. Read Thermal inertia |
---|
36 | ! 4. Read soil temperatures |
---|
37 | ! 5. Interpolate thermal inertia and temperature on the new grid, if |
---|
38 | ! necessary |
---|
39 | !====================================================================== |
---|
40 | |
---|
41 | !====================================================================== |
---|
42 | ! arguments |
---|
43 | ! --------- |
---|
44 | ! inputs: |
---|
45 | integer,intent(in) :: nid ! Input Netcdf file ID |
---|
46 | integer,intent(in) :: ngrid ! # of horizontal grid points |
---|
47 | integer,intent(in) :: nsoil ! # of soil layers |
---|
48 | integer,intent(in) :: nqsoil ! # of tracers in the soil |
---|
49 | real,intent(in) :: tsurf(ngrid,nslope) ! surface temperature [K] |
---|
50 | integer,intent(in) :: indextime ! position on time axis |
---|
51 | ! output: |
---|
52 | real,intent(out) :: tsoil(ngrid,nsoil,nslope) ! soil temperature [K] |
---|
53 | real,intent(out) :: qsoil(ngrid,nsoil,nqsoil,nslope) ! Tracers in the subsurface [kg/kg] |
---|
54 | |
---|
55 | !====================================================================== |
---|
56 | ! local variables: |
---|
57 | integer ierr ! status (returned by NetCDF functions) |
---|
58 | integer nvarid ! ID of NetCDF variable |
---|
59 | integer dimid ! ID of NetCDF dimension |
---|
60 | integer dimlen ! length along the "depth" dimension |
---|
61 | integer ndims ! # of dimensions of read <inertiedat> data |
---|
62 | integer ig,iloop ! loop counters |
---|
63 | integer islope |
---|
64 | |
---|
65 | integer edges(3),corner(3) ! to read a specific time |
---|
66 | |
---|
67 | logical :: olddepthdef=.false. ! flag |
---|
68 | logical :: interpol=.false. ! flag: true if interpolation will be requiered |
---|
69 | |
---|
70 | ! to store "old" values |
---|
71 | real,dimension(:),allocatable :: surfinertia !surface thermal inertia |
---|
72 | real,dimension(:),allocatable :: oldmlayer |
---|
73 | real,dimension(:,:),allocatable :: oldinertiedat |
---|
74 | real,dimension(:,:,:),allocatable :: oldinertiesoil |
---|
75 | real,dimension(:,:,:),allocatable :: oldtsoil |
---|
76 | |
---|
77 | ! for interpolation |
---|
78 | real,dimension(:),allocatable :: oldgrid |
---|
79 | real,dimension(:),allocatable :: oldval |
---|
80 | real,dimension(:),allocatable :: newval |
---|
81 | |
---|
82 | real alpha,lay1 ! coefficients for building layers |
---|
83 | real xmin,xmax ! to display min and max of a field |
---|
84 | |
---|
85 | real,parameter :: default_volcapa=1.e6 |
---|
86 | |
---|
87 | logical :: found,ok |
---|
88 | |
---|
89 | character (len=30):: txt |
---|
90 | |
---|
91 | !====================================================================== |
---|
92 | |
---|
93 | ! 1. Depth coordinate |
---|
94 | ! ------------------- |
---|
95 | |
---|
96 | |
---|
97 | ! 1.1 Start by reading how many layers of soil there are |
---|
98 | |
---|
99 | dimlen=inquire_dimension_length("subsurface_layers") |
---|
100 | |
---|
101 | if (dimlen.ne.nsoil) then |
---|
102 | write(*,*)'soil_settings: Interpolation of soil temperature ', |
---|
103 | & 'and thermal inertia will be required!' |
---|
104 | ! if dimlen doesn't match nsoil, then interpolation of |
---|
105 | ! soil temperatures and thermal inertia will be requiered |
---|
106 | interpol=.true. |
---|
107 | ! allocate oldmlayer |
---|
108 | if (.not.allocated(oldmlayer)) then |
---|
109 | allocate(oldmlayer(dimlen),stat=ierr) |
---|
110 | if (ierr.ne.0) then |
---|
111 | write(*,*) 'soil_settings: failed allocation of oldmlayer!' |
---|
112 | call abort_physic("soil_settings", |
---|
113 | & "failed oldmlayer allocation",1) |
---|
114 | endif |
---|
115 | endif |
---|
116 | endif |
---|
117 | |
---|
118 | ! 1.2 Find out the # of dimensions <inertiedat> was defined as using |
---|
119 | ! (in ye old days, thermal inertia was only given at the "surface") |
---|
120 | |
---|
121 | ndims=inquire_field_ndims("inertiedat") |
---|
122 | |
---|
123 | ! 1.3 Read depths values or set olddepthdef flag and values |
---|
124 | if (ndims.eq.1) then ! we know that there is none |
---|
125 | write(*,*)'soil_settings: no <soildepth> field expected' |
---|
126 | write(*,*)'building one mimicking old definitions' |
---|
127 | olddepthdef=.true. |
---|
128 | interpol=.true. |
---|
129 | |
---|
130 | do iloop=1,dimlen |
---|
131 | oldmlayer(iloop)=sqrt(887.75/3.14)*((2.**(iloop-0.5))-1.) |
---|
132 | enddo |
---|
133 | else ! Look for depth |
---|
134 | ! read <depth> coordinate |
---|
135 | if (interpol) then !put values in oldmlayer |
---|
136 | call get_var("soildepth",oldmlayer,found) |
---|
137 | if (.not.found) then |
---|
138 | write(*,*)'soil_settings: Problem while reading <soildepth>' |
---|
139 | endif |
---|
140 | else ! put values in mlayer |
---|
141 | call get_var("soildepth",mlayer,found) |
---|
142 | if (.not.found) then |
---|
143 | write(*,*)'soil_settings: Problem while reading <soildepth>' |
---|
144 | endif |
---|
145 | endif !of if (interpol) |
---|
146 | endif !of if (ndims.eq.1) |
---|
147 | |
---|
148 | ! 1.4 Build mlayer(), if necessary |
---|
149 | if (interpol) then |
---|
150 | ! default mlayer distribution, following a power law: |
---|
151 | ! mlayer(k)=lay1*alpha**(k-1/2) |
---|
152 | lay1=2.e-4 |
---|
153 | alpha=2 |
---|
154 | do iloop=0,nsoil-1 |
---|
155 | mlayer(iloop)=lay1*(alpha**(iloop-0.5)) |
---|
156 | enddo |
---|
157 | endif |
---|
158 | |
---|
159 | ! 1.5 Build layer(); following the same law as mlayer() |
---|
160 | ! Assuming layer distribution follows mid-layer law: |
---|
161 | ! layer(k)=lay1*alpha**(k-1) |
---|
162 | lay1=sqrt(mlayer(0)*mlayer(1)) |
---|
163 | alpha=mlayer(1)/mlayer(0) |
---|
164 | do iloop=1,nsoil |
---|
165 | layer(iloop)=lay1*(alpha**(iloop-1)) |
---|
166 | enddo |
---|
167 | |
---|
168 | |
---|
169 | ! 2. Volumetric heat capacity (note: it is declared in comsoil_h) |
---|
170 | ! --------------------------- |
---|
171 | ! "volcapa" is (so far) 0D and written in "controle" table of startfi file |
---|
172 | ! volcapa is read or set when "controle" is read (see tabfi.F) |
---|
173 | ! Just in case, we check here that it is not zero. If it is, we |
---|
174 | ! set it to "default_volcapa" |
---|
175 | |
---|
176 | if (volcapa.le.0.0) then |
---|
177 | write(*,*)'soil_settings: Warning, volcapa = ',volcapa |
---|
178 | write(*,*)' That doesn t seem right' |
---|
179 | write(*,*)' Initializing Volumetric heat capacity to ', |
---|
180 | & default_volcapa |
---|
181 | volcapa=default_volcapa |
---|
182 | endif |
---|
183 | |
---|
184 | ! 3. Thermal inertia for present day climate -inertiedat- & the one given by the pem -ineritesoil-(note: it is declared in comsoil_h) |
---|
185 | ! ------------------ |
---|
186 | ! Knowing the # of dimensions <inertidat> was defined as using, |
---|
187 | ! read/build thermal inertia |
---|
188 | |
---|
189 | ! 3.1 Present day - inertiedat |
---|
190 | |
---|
191 | if (ndims.eq.1) then ! "old 2D-surface" format |
---|
192 | write(*,*)'soil_settings: Thermal inertia is only given as surfac |
---|
193 | &e data!' |
---|
194 | ! Read Surface thermal inertia |
---|
195 | allocate(surfinertia(ngrid)) |
---|
196 | ! ierr=nf90_get_var(nid,nvarid,surfinertia) |
---|
197 | ! if (ierr.NE.nf90_noerr) then |
---|
198 | ! write(*,*)'soil_settings: Problem while reading <inertiedat>' |
---|
199 | ! write(*,*)trim(nf90_strerror(ierr)) |
---|
200 | ! call abort |
---|
201 | ! endif |
---|
202 | call get_field("inertiedat",surfinertia,found) |
---|
203 | if (.not.found) then |
---|
204 | write(*,*) "soil_settings: Failed loading <inertiedat>" |
---|
205 | call abort_physic("soil_settings", |
---|
206 | & "failed loading <inertiedat>",1) |
---|
207 | endif |
---|
208 | |
---|
209 | write(*,*)' => Building soil thermal inertia (using reference sur |
---|
210 | &face thermal inertia)' |
---|
211 | do iloop=1,nsoil |
---|
212 | inertiedat(:,iloop)=surfinertia(:) |
---|
213 | enddo |
---|
214 | deallocate(surfinertia) |
---|
215 | |
---|
216 | else ! "3D surface+depth" format |
---|
217 | if (interpol) then ! put values in oldinertiedat |
---|
218 | if (.not.allocated(oldinertiedat)) then |
---|
219 | allocate(oldinertiedat(ngrid,dimlen),stat=ierr) |
---|
220 | if (ierr.ne.0) then |
---|
221 | write(*,*) 'soil_settings: failed allocation of ', |
---|
222 | & 'oldinertiedat!' |
---|
223 | call abort_physic("soil_settings", |
---|
224 | & "failed allocation of oldinertiedat",1) |
---|
225 | endif |
---|
226 | endif ! of if (.not.allocated(oldinertiedat)) |
---|
227 | |
---|
228 | call get_field("inertiedat",oldinertiedat,found) |
---|
229 | if (.not.found) then |
---|
230 | write(*,*) "soil_settings: Failed loading <inertiedat>" |
---|
231 | call abort_physic("soil_settings", |
---|
232 | & "failed loading <inertiedat>",1) |
---|
233 | endif |
---|
234 | else ! put values in therm_i |
---|
235 | call get_field("inertiedat",inertiedat,found) |
---|
236 | if (.not.found) then |
---|
237 | write(*,*) "soil_settings: Failed loading <inertiedat>" |
---|
238 | call abort_physic("soil_settings", |
---|
239 | & "failed loading <inertiedat>",1) |
---|
240 | endif |
---|
241 | endif ! of if (interpol) |
---|
242 | endif ! of if (ndims.eq.1) |
---|
243 | |
---|
244 | ! 3.2 Inertie from the PEM |
---|
245 | |
---|
246 | ok=inquire_field("inertiesoil") |
---|
247 | |
---|
248 | if (.not.ok) then |
---|
249 | write(*,*)'soil_settings: Field <inertiesoil> not found!' |
---|
250 | write(*,*)' => Building <inertiesoil> from surface values'// |
---|
251 | & ' <inertiedat>' |
---|
252 | ! Case 1: No interpolation needed, we just copy past inertiedat |
---|
253 | if(.not.(interpol)) then |
---|
254 | do islope=1,nslope |
---|
255 | inertiesoil(:,:,islope)=inertiedat(:,:) |
---|
256 | enddo |
---|
257 | else |
---|
258 | ! Case 2: Interpolation needed: we copy past old value from inertiedat |
---|
259 | if (.not.allocated(oldinertiesoil)) then |
---|
260 | allocate(oldinertiesoil(ngrid,dimlen,nslope),stat=ierr) |
---|
261 | endif |
---|
262 | do islope=1,nslope |
---|
263 | oldinertiesoil(:,:,islope)=oldinertiedat(:,:) |
---|
264 | enddo |
---|
265 | endif |
---|
266 | else ! <inertiesoil> found |
---|
267 | if (interpol) then ! put values in oldinertiesoil |
---|
268 | if (.not.allocated(oldinertiesoil)) then |
---|
269 | allocate(oldinertiesoil(ngrid,dimlen,nslope),stat=ierr) |
---|
270 | if (ierr.ne.0) then |
---|
271 | write(*,*) 'soil_settings: failed allocation of ', |
---|
272 | & 'oldinertiesoil!' |
---|
273 | call abort_physic("soil_settings", |
---|
274 | & "failed allocation of oldinertiesoil",1) |
---|
275 | endif |
---|
276 | endif |
---|
277 | |
---|
278 | call get_field("inertiesoil",oldinertiesoil,found) |
---|
279 | if (.not.found) then |
---|
280 | write(*,*) "soil_settings: Failed loading <inertiesoil>" |
---|
281 | call abort_physic("soil_settings", |
---|
282 | & "failed loading <inertiesoil>",1) |
---|
283 | endif |
---|
284 | else ! put values in inertiesoil |
---|
285 | call get_field("inertiesoil",inertiesoil,found, |
---|
286 | & timeindex=indextime) |
---|
287 | if (.not.found) then |
---|
288 | write(*,*) "soil_settings: Failed loading <inertiesoil>" |
---|
289 | call abort_physic("soil_settings", |
---|
290 | & "failed loading <inertiesoil>",1) |
---|
291 | endif |
---|
292 | endif ! of if (interpol) |
---|
293 | endif! of if (.not.ok) |
---|
294 | |
---|
295 | |
---|
296 | |
---|
297 | ! 4. Read soil temperatures |
---|
298 | ! ------------------------- |
---|
299 | |
---|
300 | ok=inquire_field("tsoil") |
---|
301 | |
---|
302 | if (.not.ok) then |
---|
303 | write(*,*)'soil_settings: Field <tsoil> not found!' |
---|
304 | write(*,*)' => Building <tsoil> from surface values <tsurf>' |
---|
305 | ! Case 1: No interpolation needed, we just copy past inertiedat |
---|
306 | if(.not.(interpol)) then |
---|
307 | do iloop=1,nsoil |
---|
308 | do islope=1,nslope |
---|
309 | tsoil(:,iloop,islope)=tsurf(:,islope) |
---|
310 | enddo |
---|
311 | enddo |
---|
312 | else |
---|
313 | ! Case 2: Interpolation needed: we copy past old value from inertiedat |
---|
314 | if (.not.allocated(oldtsoil)) then |
---|
315 | allocate(oldtsoil(ngrid,dimlen,nslope),stat=ierr) |
---|
316 | endif |
---|
317 | do iloop=1,dimlen |
---|
318 | do islope=1,nslope |
---|
319 | oldtsoil(:,iloop,islope)=tsurf(:,islope) |
---|
320 | enddo |
---|
321 | enddo |
---|
322 | endif |
---|
323 | else ! <tsoil> found |
---|
324 | if (interpol) then ! put values in oldtsoil |
---|
325 | if (.not.allocated(oldtsoil)) then |
---|
326 | allocate(oldtsoil(ngrid,dimlen,nslope),stat=ierr) |
---|
327 | if (ierr.ne.0) then |
---|
328 | write(*,*) 'soil_settings: failed allocation of ', |
---|
329 | & 'oldtsoil!' |
---|
330 | call abort_physic("soil_settings", |
---|
331 | & "failed allocation of oldtsoil",1) |
---|
332 | endif |
---|
333 | endif |
---|
334 | call get_field("tsoil",oldtsoil,found) |
---|
335 | if (.not.found) then |
---|
336 | write(*,*) "soil_settings: Failed loading <tsoil>" |
---|
337 | call abort_physic("soil_settings", |
---|
338 | & "failed loading <tsoil>",1) |
---|
339 | endif |
---|
340 | else ! put values in tsoil |
---|
341 | call get_field("tsoil",tsoil,found,timeindex=indextime) |
---|
342 | if (.not.found) then |
---|
343 | write(*,*) "soil_settings: Failed loading <tsoil>" |
---|
344 | call abort_physic("soil_settings", |
---|
345 | & "failed loading <tsoil>",1) |
---|
346 | endif |
---|
347 | endif ! of if (interpol) |
---|
348 | endif! of if (.not.ok) |
---|
349 | |
---|
350 | |
---|
351 | ! 5. If necessary, interpolate soil temperatures and thermal inertias |
---|
352 | ! ------------------------------------------------------------------- |
---|
353 | |
---|
354 | if (olddepthdef) then |
---|
355 | ! tsoil needs to be interpolated, but not therm_i |
---|
356 | allocate(oldgrid(dimlen+1)) |
---|
357 | allocate(oldval(dimlen+1)) |
---|
358 | allocate(newval(nsoil)) |
---|
359 | do ig=1,ngrid |
---|
360 | ! copy values |
---|
361 | do islope=1,nslope |
---|
362 | oldval(1)=tsurf(ig,islope) |
---|
363 | oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen,islope) |
---|
364 | ! build vertical coordinate |
---|
365 | oldgrid(1)=0. ! ground |
---|
366 | oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)* |
---|
367 | & (inertiesoil(ig,1,islope)/volcapa) |
---|
368 | ! interpolate |
---|
369 | call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil) |
---|
370 | ! copy result in tsoil |
---|
371 | tsoil(ig,:,islope)=newval(:) |
---|
372 | enddo |
---|
373 | enddo |
---|
374 | ! cleanup |
---|
375 | deallocate(oldgrid) |
---|
376 | deallocate(oldval) |
---|
377 | deallocate(newval) |
---|
378 | interpol=.false. ! no need for interpolation any more |
---|
379 | endif !of if (olddepthdef) |
---|
380 | if (interpol) then |
---|
381 | write(*,*)'soil_settings: Vertical interpolation along new grid' |
---|
382 | ! interpolate soil temperatures and thermal inertias |
---|
383 | if (.not.allocated(oldgrid)) then |
---|
384 | allocate(oldgrid(dimlen+1)) |
---|
385 | allocate(oldval(dimlen+1)) |
---|
386 | allocate(newval(nsoil)) |
---|
387 | endif |
---|
388 | |
---|
389 | ! thermal inertia - present day (inertiedat) |
---|
390 | do ig=1,ngrid |
---|
391 | ! copy data |
---|
392 | oldval(1:dimlen)=oldinertiedat(ig,1:dimlen) |
---|
393 | ! interpolate |
---|
394 | call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil) |
---|
395 | !copy result in inertiedat |
---|
396 | inertiedat(ig,:)=newval(:) |
---|
397 | ! thermal inertia - general case with PEM (inertie soil) |
---|
398 | do islope=1,nslope |
---|
399 | ! copy data |
---|
400 | oldval(1:dimlen)=oldinertiesoil(ig,1:dimlen,islope) |
---|
401 | ! interpolate |
---|
402 | call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil) |
---|
403 | !copy result in inertiedat |
---|
404 | inertiesoil(ig,:,islope)=newval(:) |
---|
405 | enddo |
---|
406 | enddo ! ig |
---|
407 | |
---|
408 | ! soil temperature |
---|
409 | ! vertical coordinate |
---|
410 | oldgrid(1)=0.0 |
---|
411 | oldgrid(2:dimlen+1)=oldmlayer(1:dimlen) |
---|
412 | do ig=1,ngrid |
---|
413 | do islope=1,nslope |
---|
414 | ! data |
---|
415 | oldval(1)=tsurf(ig,islope) |
---|
416 | oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen,islope) |
---|
417 | ! interpolate |
---|
418 | call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil) |
---|
419 | ! copy result in inertiedat |
---|
420 | tsoil(ig,:,islope)=newval(:) |
---|
421 | enddo |
---|
422 | enddo |
---|
423 | |
---|
424 | !cleanup |
---|
425 | deallocate(oldgrid) |
---|
426 | deallocate(oldval) |
---|
427 | deallocate(newval) |
---|
428 | deallocate(oldinertiedat) |
---|
429 | deallocate(oldinertiesoil) |
---|
430 | deallocate(oldtsoil) |
---|
431 | endif ! of if (interpol) |
---|
432 | |
---|
433 | |
---|
434 | |
---|
435 | ! 6. Load Geothermal Flux |
---|
436 | ! ---------------------------------------------------------------------- |
---|
437 | |
---|
438 | |
---|
439 | call get_field("flux_geo",flux_geo,found,timeindex=indextime) |
---|
440 | if (.not.found) then |
---|
441 | write(*,*) "soil_settings: Failed loading <flux_geo>" |
---|
442 | write(*,*) "soil_settings: Will put <flux_geo> to 0." |
---|
443 | flux_geo(:,:) = 0. |
---|
444 | endif |
---|
445 | |
---|
446 | |
---|
447 | ! 7. Adsorption |
---|
448 | ! ---------------------------------------------------------------------- |
---|
449 | |
---|
450 | |
---|
451 | if (adsorption_soil) then |
---|
452 | ! Subsurface water vapor |
---|
453 | txt="h2o_vap_soil" |
---|
454 | write(*,*) 'phyetat0: loading subsurface tracer', |
---|
455 | & ' h2o_vap_soil' |
---|
456 | call get_field(txt,qsoil(:,:,igcm_h2o_vap_soil,:),found, |
---|
457 | & indextime) |
---|
458 | write(*,*) 'found',found |
---|
459 | write(*,*) 'igcm_',igcm_h2o_vap_soil |
---|
460 | write(*,*) 'q=',qsoil(:,:,:,:) |
---|
461 | if (.not.found) then |
---|
462 | write(*,*) "phyetat0: Failed loading <",trim(txt),">" |
---|
463 | write(*,*) " ",trim(txt)," is set to zero" |
---|
464 | qsoil(:,:,igcm_h2o_vap_soil,:)= 0. |
---|
465 | else |
---|
466 | write(*,*) "phyetat0: suburface tracer <",trim(txt), |
---|
467 | & "> range:", minval(qsoil(:,:,igcm_h2o_vap_soil,:)), |
---|
468 | & maxval(qsoil(:,:,igcm_h2o_vap_soil,:)) |
---|
469 | endif |
---|
470 | ! Subsurface ice |
---|
471 | txt="h2o_ice_soil" |
---|
472 | write(*,*) 'phyetat0: loading subsurface tracer', |
---|
473 | & ' h2o_ice_soil' |
---|
474 | call get_field(txt,qsoil(:,:,igcm_h2o_ice_soil,:),found, |
---|
475 | & indextime) |
---|
476 | if (.not.found) then |
---|
477 | write(*,*) "phyetat0: Failed loading <",trim(txt),">" |
---|
478 | write(*,*) " ",trim(txt)," is set to zero" |
---|
479 | qsoil(:,:,igcm_h2o_ice_soil,:)= 0. |
---|
480 | else |
---|
481 | write(*,*) "phyetat0: suburface tracer <",trim(txt), |
---|
482 | & "> range:", |
---|
483 | & minval(qsoil(:,:,igcm_h2o_ice_soil,:)), |
---|
484 | & maxval(qsoil(:,:,igcm_h2o_ice_soil,:)) |
---|
485 | endif |
---|
486 | ! Adsorbed water |
---|
487 | txt="h2o_vap_ads" |
---|
488 | write(*,*) 'phyetat0: loading subsurface tracer', |
---|
489 | & ' h2o_vap_ads' |
---|
490 | call get_field(txt,qsoil(:,:,igcm_h2o_vap_ads,:),found, |
---|
491 | & indextime) |
---|
492 | if (.not.found) then |
---|
493 | write(*,*) "phyetat0: Failed loading <",trim(txt),">" |
---|
494 | write(*,*) " ",trim(txt)," is set to zero" |
---|
495 | qsoil(:,:,igcm_h2o_vap_ads,:)= 0. |
---|
496 | else |
---|
497 | write(*,*) "phyetat0: suburface tracer <",trim(txt),"> |
---|
498 | & range:", |
---|
499 | & minval(qsoil(:,:,igcm_h2o_vap_ads,:)), |
---|
500 | & maxval(qsoil(:,:,igcm_h2o_vap_ads,:)) |
---|
501 | endif |
---|
502 | |
---|
503 | endif ! of adsorption_soil |
---|
504 | |
---|
505 | |
---|
506 | |
---|
507 | ! 8. Report min and max values of soil temperatures and thermal inertias |
---|
508 | ! ---------------------------------------------------------------------- |
---|
509 | |
---|
510 | write(*,*) |
---|
511 | write(*,*)'Soil volumetric heat capacity:',volcapa |
---|
512 | |
---|
513 | xmin = MINVAL(inertiedat) |
---|
514 | xmax = MAXVAL(inertiedat) |
---|
515 | write(*,*)'Soil thermal inertia <inertiedat>:',xmin,xmax |
---|
516 | |
---|
517 | xmin = MINVAL(inertiesoil) |
---|
518 | xmax = MAXVAL(inertiesoil) |
---|
519 | write(*,*)'Soil thermal inertia <inertiesoil>:',xmin,xmax |
---|
520 | |
---|
521 | xmin = 1.0E+20 |
---|
522 | xmax = -1.0E+20 |
---|
523 | xmin = MINVAL(tsoil) |
---|
524 | xmax = MAXVAL(tsoil) |
---|
525 | write(*,*)'Soil temperature <tsoil>:',xmin,xmax |
---|
526 | |
---|
527 | end |
---|