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