source: trunk/LMDZ.GENERIC/libf/phystd/soil_settings.F @ 786

Last change on this file since 786 was 786, checked in by jleconte, 13 years ago

19/09/2012 == JL

  • Correction in largescale to improve robustness when large water vapor amount
  • Correction in soil_setting to allow change of the number of subsurface layers


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