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

Last change on this file since 1477 was 1350, checked in by emillour, 10 years ago

Generic model physics:

  • Added missing allocation in soil_settings, needed when changing number of soil layers.

EM

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