source: trunk/LMDZ.TITAN/libf/phytitan/soil_settings.F @ 3000

Last change on this file since 3000 was 1538, checked in by emillour, 9 years ago

Generic GCM:

  • Made nsoilmx be no longer a "parameter" and thus added the possibility to define the number of subsurface layers nsoilmx, along with first layer thickness "lay1_soil" and companion coefficient "alpha_soil", in callphys.def at run time. As before (when these were hard-coded), these are such that the depth of soil mid-layers are: mlayer(k)=lay1_soil*alpha_soil(k-1/2), for k=0,nsoil-1

EM

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