1 | subroutine soil_settings(nid,ngrid,nsoil,tsurf,tsoil) |
---|
2 | |
---|
3 | use comsoil_h |
---|
4 | |
---|
5 | implicit none |
---|
6 | |
---|
7 | !====================================================================== |
---|
8 | ! Author: Ehouarn Millour (07/2006) |
---|
9 | ! |
---|
10 | ! Purpose: Read and/or initialise soil depths and properties |
---|
11 | ! |
---|
12 | ! This subroutine reads from a NetCDF file (opened by the caller) |
---|
13 | ! of "startfi.nc" format. |
---|
14 | ! The various actions and variable read/initialized are: |
---|
15 | ! 1. Check out the number of soil layers (in datafile); if it isn't equal |
---|
16 | ! to nsoil, then some interpolation will be required |
---|
17 | ! Also check if data in file "startfi.nc" is in older format (ie: |
---|
18 | ! thermal inertia was depth-independent; and there was no "depth" |
---|
19 | ! coordinate. |
---|
20 | ! Read/build layer (and midlayer) depths |
---|
21 | ! 2. Read volumetric specific heat (or initialise it to default value) |
---|
22 | ! 3. Read Thermal inertia |
---|
23 | ! 4. Read soil temperatures |
---|
24 | ! 5. Interpolate thermal inertia and temperature on the new grid, if |
---|
25 | ! necessary |
---|
26 | !====================================================================== |
---|
27 | |
---|
28 | #include "dimensions.h" |
---|
29 | #include "dimphys.h" |
---|
30 | |
---|
31 | #include "netcdf.inc" |
---|
32 | !====================================================================== |
---|
33 | ! arguments |
---|
34 | ! --------- |
---|
35 | ! inputs: |
---|
36 | integer nid ! Input Netcdf file ID |
---|
37 | integer ngrid ! # of horizontal grid points |
---|
38 | integer nsoil ! # of soil layers |
---|
39 | integer sta ! # at which reading starts |
---|
40 | real tsurf(ngrid) ! surface temperature |
---|
41 | ! output: |
---|
42 | real tsoil(ngrid,nsoilmx) ! soil temperature |
---|
43 | |
---|
44 | !====================================================================== |
---|
45 | ! local variables: |
---|
46 | integer ierr ! status (returned by NetCDF functions) |
---|
47 | integer nvarid ! ID of NetCDF variable |
---|
48 | integer dimid ! ID of NetCDF dimension |
---|
49 | integer dimlen ! length along the "depth" dimension |
---|
50 | integer ndims ! # of dimensions of read <inertiedat> data |
---|
51 | integer ig,iloop ! loop counters |
---|
52 | |
---|
53 | logical :: olddepthdef=.false. ! flag |
---|
54 | logical :: interpol=.false. ! flag: true if interpolation will be requiered |
---|
55 | |
---|
56 | ! to store "old" values |
---|
57 | real,dimension(:),allocatable :: surfinertia !surface thermal inertia |
---|
58 | real,dimension(:),allocatable :: oldmlayer |
---|
59 | real,dimension(:,:),allocatable :: oldinertiedat |
---|
60 | real,dimension(:,:),allocatable :: oldtsoil |
---|
61 | |
---|
62 | ! for interpolation |
---|
63 | real,dimension(:),allocatable :: oldgrid |
---|
64 | real,dimension(:),allocatable :: oldval |
---|
65 | real,dimension(:),allocatable :: newval |
---|
66 | |
---|
67 | real alpha,lay1 ! coefficients for building layers |
---|
68 | real xmin,xmax ! to display min and max of a field |
---|
69 | |
---|
70 | real,parameter :: default_volcapa=1.e6 |
---|
71 | |
---|
72 | integer, dimension(2) :: sta2d |
---|
73 | integer, dimension(2) :: siz2d |
---|
74 | |
---|
75 | !====================================================================== |
---|
76 | |
---|
77 | !! this is defined in comsoil_h |
---|
78 | IF (.not.ALLOCATED(layer)) |
---|
79 | . ALLOCATE(layer(nsoilmx)) |
---|
80 | IF (.not.ALLOCATED(mlayer)) |
---|
81 | . ALLOCATE(mlayer(0:nsoilmx-1)) |
---|
82 | IF (.not.ALLOCATED(inertiedat)) |
---|
83 | . ALLOCATE(inertiedat(ngrid,nsoilmx)) |
---|
84 | |
---|
85 | !! this is defined in dimphys.h |
---|
86 | sta = cursor |
---|
87 | |
---|
88 | ! 1. Depth coordinate |
---|
89 | ! ------------------- |
---|
90 | |
---|
91 | ! 1.1 Start by reading how many layers of soil there are |
---|
92 | |
---|
93 | ierr=NF_INQ_DIMID(nid,"subsurface_layers",dimid) |
---|
94 | if (ierr.ne.NF_NOERR) then |
---|
95 | write(*,*)'soil_settings: Problem reading <subsurface_layers>' |
---|
96 | call abort |
---|
97 | endif |
---|
98 | |
---|
99 | ierr=NF_INQ_DIMLEN(nid,dimid,dimlen) |
---|
100 | if (ierr.ne.NF_NOERR) then |
---|
101 | write(*,*)'soil_settings: Problem getting <subsurface_layers>', |
---|
102 | & 'length' |
---|
103 | call abort |
---|
104 | endif |
---|
105 | |
---|
106 | if (dimlen.ne.nsoil) then |
---|
107 | write(*,*)'soil_settings: Interpolation of soil temperature ', |
---|
108 | & 'and thermal inertia will be required!' |
---|
109 | ! if dimlen doesn't match nsoil, then interpolation of |
---|
110 | ! soil temperatures and thermal inertia will be requiered |
---|
111 | interpol=.true. |
---|
112 | endif |
---|
113 | |
---|
114 | sta2d = (/sta,1/) |
---|
115 | siz2d = (/ngrid,dimlen/) |
---|
116 | |
---|
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 | ! Look for thermal inertia data |
---|
121 | ierr=NF_INQ_VARID(nid, "inertiedat", nvarid) |
---|
122 | if (ierr.NE.NF_NOERR) then |
---|
123 | write(*,*)'soil_settings: Field <inertiedat> not found!' |
---|
124 | call abort |
---|
125 | endif |
---|
126 | |
---|
127 | ! Read the # of dimensions <inertidat> was defined as using |
---|
128 | ierr=NF_INQ_VARNDIMS(nid,nvarid,ndims) |
---|
129 | ! if (ndims.eq.1) then we have the "old 2D-surface" format |
---|
130 | |
---|
131 | ! 1.3 Read depths values or set olddepthdef flag and values |
---|
132 | if (ndims.eq.1) then ! we know that there is none |
---|
133 | write(*,*)'soil_settings: no <soildepth> field expected' |
---|
134 | write(*,*)'building one mimicking old definitions' |
---|
135 | olddepthdef=.true. |
---|
136 | interpol=.true. |
---|
137 | ! allocate oldmlayer |
---|
138 | if (.not.allocated(oldmlayer)) then |
---|
139 | allocate(oldmlayer(dimlen),stat=ierr) |
---|
140 | if (ierr.ne.0) then |
---|
141 | write(*,*) 'soil_settings: failed allocation of oldmlayer!' |
---|
142 | stop |
---|
143 | endif |
---|
144 | endif |
---|
145 | do iloop=1,dimlen |
---|
146 | oldmlayer(iloop)=sqrt(887.75/3.14)*((2.**(iloop-0.5))-1.) |
---|
147 | enddo |
---|
148 | else ! Look for depth |
---|
149 | ierr=NF_INQ_VARID(nid,"soildepth",nvarid) |
---|
150 | if (ierr.ne.NF_NOERR) then |
---|
151 | write(*,*)'soil_settings: Field <soildepth> not found!' |
---|
152 | call abort |
---|
153 | endif |
---|
154 | ! read <depth> coordinate |
---|
155 | if (interpol) then !put values in oldmlayer |
---|
156 | if (.not.allocated(oldmlayer)) then |
---|
157 | allocate(oldmlayer(dimlen),stat=ierr) |
---|
158 | endif |
---|
159 | #ifdef NC_DOUBLE |
---|
160 | ierr = NF_GET_VAR_DOUBLE(nid,nvarid,oldmlayer) |
---|
161 | #else |
---|
162 | ierr = NF_GET_VAR_REAL(nid,nvarid,oldmlayer) |
---|
163 | #endif |
---|
164 | if (ierr.ne.NF_NOERR) then |
---|
165 | write(*,*)'soil_settings: Problem while reading <soildepth>' |
---|
166 | call abort |
---|
167 | endif |
---|
168 | else ! put values in mlayer |
---|
169 | #ifdef NC_DOUBLE |
---|
170 | ierr = NF_GET_VAR_DOUBLE(nid,nvarid,mlayer) |
---|
171 | #else |
---|
172 | ierr = NF_GET_VAR_REAL(nid,nvarid,mlayer) |
---|
173 | #endif |
---|
174 | if (ierr.ne.NF_NOERR) then |
---|
175 | write(*,*)'soil_settings: Problem while reading <soildepth>' |
---|
176 | call abort |
---|
177 | endif |
---|
178 | endif !of if (interpol) |
---|
179 | endif !of if (ndims.eq.1) |
---|
180 | |
---|
181 | ! 1.4 Build mlayer(), if necessary |
---|
182 | if (interpol) then |
---|
183 | ! default mlayer distribution, following a power law: |
---|
184 | ! mlayer(k)=lay1*alpha**(k-1/2) |
---|
185 | lay1=2.e-4 !mars case (nsoilmx=18) |
---|
186 | ! lay1=3.e-2 !earth case (nsoilmx=13) |
---|
187 | alpha=2 |
---|
188 | do iloop=0,nsoil-1 |
---|
189 | mlayer(iloop)=lay1*(alpha**(iloop-0.5)) |
---|
190 | enddo |
---|
191 | endif |
---|
192 | |
---|
193 | ! 1.5 Build layer(); following the same law as mlayer() |
---|
194 | ! Assuming layer distribution follows mid-layer law: |
---|
195 | ! layer(k)=lay1*alpha**(k-1) |
---|
196 | lay1=sqrt(mlayer(0)*mlayer(1)) |
---|
197 | alpha=mlayer(1)/mlayer(0) |
---|
198 | do iloop=1,nsoil |
---|
199 | layer(iloop)=lay1*(alpha**(iloop-1)) |
---|
200 | enddo |
---|
201 | |
---|
202 | ! 2. Volumetric heat capacity (note: it is declared in comsoil.h) |
---|
203 | ! --------------------------- |
---|
204 | ! "volcapa" is (so far) 0D and written in "controle" table of startfi file |
---|
205 | ! volcapa is read or set when "controle" is read (see tabfi.F) |
---|
206 | ! Just in case, we check here that it is not zero. If it is, we |
---|
207 | ! set it to "default_volcapa" |
---|
208 | |
---|
209 | if (volcapa.le.0.0) then |
---|
210 | write(*,*)'soil_settings: Warning, volcapa = ',volcapa |
---|
211 | write(*,*)' That doesn t seem right' |
---|
212 | write(*,*)' Initializing Volumetric heat capacity to ', |
---|
213 | & default_volcapa |
---|
214 | volcapa=default_volcapa |
---|
215 | endif |
---|
216 | ! Look for it |
---|
217 | ! ierr=NF_INQ_VARID(nid,"volcapa",nvarid) |
---|
218 | ! if (ierr.NE.NF_NOERR) then |
---|
219 | ! write(*,*)'soil_settings: Field <volcapa> not found!' |
---|
220 | ! write(*,*)'Initializing Volumetric heat capacity to ', |
---|
221 | ! & default_volcapa |
---|
222 | ! volcapa=default_volcapa |
---|
223 | ! else |
---|
224 | !#ifdef NC_DOUBLE |
---|
225 | ! ierr = NF_GET_VAR_DOUBLE(nid,nvarid,volcapa) |
---|
226 | !#else |
---|
227 | ! ierr = NF_GET_VAR_REAL(nid,nvarid,volcapa) |
---|
228 | !#endif |
---|
229 | ! if (ierr.ne.NF_NOERR) then |
---|
230 | ! write(*,*)'soil_settings: Problem while reading <volcapa>' |
---|
231 | ! call abort |
---|
232 | ! endif |
---|
233 | ! endif |
---|
234 | |
---|
235 | ! 3. Thermal inertia (note: it is declared in comsoil.h) |
---|
236 | ! ------------------ |
---|
237 | |
---|
238 | ! 3.1 Look (again) for thermal inertia data (to reset nvarid) |
---|
239 | ierr=NF_INQ_VARID(nid, "inertiedat", nvarid) |
---|
240 | if (ierr.NE.NF_NOERR) then |
---|
241 | write(*,*)'soil_settings: Field <inertiedat> not found!' |
---|
242 | call abort |
---|
243 | endif |
---|
244 | |
---|
245 | ! 3.2 Knowing the # of dimensions <inertidat> was defined as using, |
---|
246 | ! read/build thermal inertia |
---|
247 | |
---|
248 | if (ndims.eq.1) then ! "old 2D-surface" format |
---|
249 | write(*,*)'soil_settings: Thermal inertia is only given as surfac |
---|
250 | &e data!' |
---|
251 | ! Read Surface thermal inertia |
---|
252 | allocate(surfinertia(ngrid)) |
---|
253 | #ifdef NC_DOUBLE |
---|
254 | ierr = NF_GET_VARA_DOUBLE(nid, nvarid, sta, ngrid, surfinertia) |
---|
255 | #else |
---|
256 | ierr = NF_GET_VARA_REAL(nid, nvarid, sta, ngrid, surfinertia) |
---|
257 | #endif |
---|
258 | if (ierr.NE.NF_NOERR) then |
---|
259 | write(*,*)'soil_settings: Problem while reading <inertiedat>' |
---|
260 | call abort |
---|
261 | endif |
---|
262 | |
---|
263 | write(*,*)' => Building soil thermal inertia (using reference sur |
---|
264 | &face thermal inertia)' |
---|
265 | do iloop=1,nsoil |
---|
266 | inertiedat(:,iloop)=surfinertia(:) |
---|
267 | enddo |
---|
268 | deallocate(surfinertia) |
---|
269 | |
---|
270 | else ! "3D surface+depth" format |
---|
271 | if (interpol) then ! put values in oldinertiedat |
---|
272 | if (.not.allocated(oldinertiedat)) then |
---|
273 | allocate(oldinertiedat(ngrid,dimlen),stat=ierr) |
---|
274 | if (ierr.ne.0) then |
---|
275 | write(*,*) 'soil_settings: failed allocation of ', |
---|
276 | & 'oldinertiedat!' |
---|
277 | stop |
---|
278 | endif |
---|
279 | endif |
---|
280 | #ifdef NC_DOUBLE |
---|
281 | ierr = NF_GET_VARA_DOUBLE(nid,nvarid,sta2d,siz2d,oldinertiedat) |
---|
282 | #else |
---|
283 | ierr = NF_GET_VARA_REAL(nid,nvarid,sta2d,siz2d,oldinertiedat) |
---|
284 | #endif |
---|
285 | if (ierr.NE.NF_NOERR) then |
---|
286 | write(*,*)'soil_settings: Problem while reading <inertiedat>' |
---|
287 | call abort |
---|
288 | endif |
---|
289 | else ! put values in therm_i |
---|
290 | #ifdef NC_DOUBLE |
---|
291 | ierr = NF_GET_VARA_DOUBLE(nid,nvarid,sta2d,siz2d,inertiedat) |
---|
292 | #else |
---|
293 | ierr = NF_GET_VARA_REAL(nid,nvarid,sta2d,siz2d,inertiedat) |
---|
294 | #endif |
---|
295 | if (ierr.NE.NF_NOERR) then |
---|
296 | write(*,*)'soil_settings: Problem while reading <inertiedat>' |
---|
297 | call abort |
---|
298 | endif |
---|
299 | endif ! of if (interpol) |
---|
300 | endif ! of if (ndims.eq.1) |
---|
301 | |
---|
302 | ! 4. Read soil temperatures |
---|
303 | ! ------------------------- |
---|
304 | |
---|
305 | ierr=NF_INQ_VARID(nid,"tsoil",nvarid) |
---|
306 | if (ierr.ne.NF_NOERR) then |
---|
307 | write(*,*)'soil_settings: Field <tsoil> not found!' |
---|
308 | write(*,*)' => Building <tsoil> from surface values <tsurf>' |
---|
309 | do iloop=1,nsoil |
---|
310 | tsoil(:,iloop)=tsurf(:) |
---|
311 | enddo |
---|
312 | else ! <tsoil> found |
---|
313 | if (interpol) then ! put values in oldtsoil |
---|
314 | if (.not.allocated(oldtsoil)) then |
---|
315 | allocate(oldtsoil(ngrid,dimlen),stat=ierr) |
---|
316 | if (ierr.ne.0) then |
---|
317 | write(*,*) 'soil_settings: failed allocation of ', |
---|
318 | & 'oldtsoil!' |
---|
319 | stop |
---|
320 | endif |
---|
321 | endif |
---|
322 | #ifdef NC_DOUBLE |
---|
323 | ierr=NF_GET_VARA_DOUBLE(nid,nvarid,sta2d,siz2d,oldtsoil) |
---|
324 | #else |
---|
325 | ierr=NF_GET_VARA_REAL(nid,nvarid,sta2d,siz2d,oldtsoil) |
---|
326 | #endif |
---|
327 | if (ierr.ne.NF_NOERR) then |
---|
328 | write(*,*)'soil_settings: Problem while reading <tsoil>' |
---|
329 | call abort |
---|
330 | endif |
---|
331 | else ! put values in tsoil |
---|
332 | #ifdef NC_DOUBLE |
---|
333 | ierr=NF_GET_VARA_DOUBLE(nid,nvarid,sta2d,siz2d,tsoil) |
---|
334 | #else |
---|
335 | ierr=NF_GET_VARA_REAL(nid,nvarid,sta2d,siz2d,tsoil) |
---|
336 | #endif |
---|
337 | if (ierr.ne.NF_NOERR) then |
---|
338 | write(*,*)'soil_settings: Problem while reading <tsoil>' |
---|
339 | call abort |
---|
340 | endif |
---|
341 | endif ! of if (interpol) |
---|
342 | endif |
---|
343 | |
---|
344 | ! 5. If necessary, interpolate soil temperatures and thermal inertias |
---|
345 | ! ------------------------------------------------------------------- |
---|
346 | |
---|
347 | if (olddepthdef) then |
---|
348 | ! tsoil needs to be interpolated, but not therm_i |
---|
349 | allocate(oldgrid(dimlen+1)) |
---|
350 | allocate(oldval(dimlen+1)) |
---|
351 | allocate(newval(nsoil)) |
---|
352 | |
---|
353 | do ig=1,ngrid |
---|
354 | ! copy values |
---|
355 | oldval(1)=tsurf(ig) |
---|
356 | oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen) |
---|
357 | ! build vertical coordinate |
---|
358 | oldgrid(1)=0. ! ground |
---|
359 | oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)* |
---|
360 | & (inertiedat(ig,1)/volcapa) |
---|
361 | ! interpolate |
---|
362 | call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil) |
---|
363 | ! copy result in tsoil |
---|
364 | tsoil(ig,:)=newval(:) |
---|
365 | enddo |
---|
366 | |
---|
367 | ! cleanup |
---|
368 | deallocate(oldgrid) |
---|
369 | deallocate(oldval) |
---|
370 | deallocate(newval) |
---|
371 | interpol=.false. ! no need for interpolation any more |
---|
372 | endif !of if (olddepthdef) |
---|
373 | |
---|
374 | if (interpol) then |
---|
375 | write(*,*)'soil_settings: Vertical interpolation along new grid' |
---|
376 | ! interpolate soil temperatures and thermal inertias |
---|
377 | if (.not.allocated(oldgrid)) then |
---|
378 | allocate(oldgrid(dimlen+1)) |
---|
379 | allocate(oldval(dimlen+1)) |
---|
380 | allocate(newval(nsoil)) |
---|
381 | endif |
---|
382 | |
---|
383 | ! thermal inertia |
---|
384 | do ig=1,ngrid |
---|
385 | ! copy data |
---|
386 | oldval(1:dimlen)=oldinertiedat(ig,dimlen) |
---|
387 | ! interpolate |
---|
388 | call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil) |
---|
389 | !copy result in inertiedat |
---|
390 | inertiedat(ig,:)=newval(:) |
---|
391 | enddo |
---|
392 | |
---|
393 | ! soil temperature |
---|
394 | ! vertical coordinate |
---|
395 | oldgrid(1)=0.0 |
---|
396 | oldgrid(2:dimlen+1)=oldmlayer(1:dimlen) |
---|
397 | do ig=1,ngrid |
---|
398 | ! data |
---|
399 | oldval(1)=tsurf(ig) |
---|
400 | oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen) |
---|
401 | ! interpolate |
---|
402 | call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil) |
---|
403 | ! copy result in inertiedat |
---|
404 | tsoil(ig,:)=newval(:) |
---|
405 | enddo |
---|
406 | |
---|
407 | !cleanup |
---|
408 | deallocate(oldgrid) |
---|
409 | deallocate(oldval) |
---|
410 | deallocate(newval) |
---|
411 | deallocate(oldinertiedat) |
---|
412 | deallocate(oldtsoil) |
---|
413 | endif ! of if (interpol) |
---|
414 | |
---|
415 | ! 6. Report min and max values of soil temperatures and thermal inertias |
---|
416 | ! ---------------------------------------------------------------------- |
---|
417 | |
---|
418 | write(*,*) |
---|
419 | write(*,*)'Soil volumetric heat capacity:',volcapa |
---|
420 | |
---|
421 | xmin = MINVAL(inertiedat) |
---|
422 | xmax = MAXVAL(inertiedat) |
---|
423 | write(*,*)'Soil thermal inertia <inertiedat>:',xmin,xmax |
---|
424 | |
---|
425 | xmin = 1.0E+20 |
---|
426 | xmax = -1.0E+20 |
---|
427 | xmin = MINVAL(tsoil) |
---|
428 | xmax = MAXVAL(tsoil) |
---|
429 | write(*,*)'Soil temperature <tsoil>:',xmin,xmax |
---|
430 | |
---|
431 | end |
---|