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

Last change on this file since 601 was 135, checked in by aslmd, 14 years ago

CHANGEMENT ARBORESCENCE ETAPE 2 -- NON COMPLET

File size: 13.0 KB
RevLine 
[135]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#ifdef NC_DOUBLE
136           ierr = NF_GET_VAR_DOUBLE(nid,nvarid,oldmlayer)
137#else
138           ierr = NF_GET_VAR_REAL(nid,nvarid,oldmlayer)
139#endif
140          if (ierr.ne.NF_NOERR) then
141           write(*,*)'soil_settings: Problem while reading <soildepth>'
142           call abort
143          endif
144        else ! put values in mlayer
145#ifdef NC_DOUBLE
146           ierr = NF_GET_VAR_DOUBLE(nid,nvarid,mlayer)
147#else
148           ierr = NF_GET_VAR_REAL(nid,nvarid,mlayer)
149#endif
150          if (ierr.ne.NF_NOERR) then
151           write(*,*)'soil_settings: Problem while reading <soildepth>'
152           call abort
153          endif
154        endif !of if (interpol)
155      endif !of if (ndims.eq.1)
156
157! 1.4 Build mlayer(), if necessary
158      if (interpol) then
159      ! default mlayer distribution, following a power law:
160      !  mlayer(k)=lay1*alpha**(k-1/2)
161        lay1=2.e-4
162        alpha=2
163        do iloop=0,nsoil-1
164          mlayer(iloop)=lay1*(alpha**(iloop-0.5))
165        enddo
166      endif
167
168! 1.5 Build layer(); following the same law as mlayer()
169      ! Assuming layer distribution follows mid-layer law:
170      ! layer(k)=lay1*alpha**(k-1)
171      lay1=sqrt(mlayer(0)*mlayer(1))
172      alpha=mlayer(1)/mlayer(0)
173      do iloop=1,nsoil
174        layer(iloop)=lay1*(alpha**(iloop-1))
175      enddo
176
177! 2. Volumetric heat capacity (note: it is declared in comsoil.h)
178! ---------------------------
179! "volcapa" is (so far) 0D and written in "controle" table of startfi file
180! volcapa is read or set when "controle" is read (see tabfi.F)
181! Just in case, we check here that it is not zero. If it is, we
182! set it to "default_volcapa"
183
184      if (volcapa.le.0.0) then
185        write(*,*)'soil_settings: Warning, volcapa = ',volcapa
186        write(*,*)'               That doesn t seem right'
187        write(*,*)'        Initializing Volumetric heat capacity to ',
188     &             default_volcapa
189        volcapa=default_volcapa
190      endif
191! Look for it
192!      ierr=NF_INQ_VARID(nid,"volcapa",nvarid)
193!      if (ierr.NE.NF_NOERR) then
194!         write(*,*)'soil_settings: Field <volcapa> not found!'
195!         write(*,*)'Initializing Volumetric heat capacity to ',
196!     &             default_volcapa
197!         volcapa=default_volcapa
198!      else
199!#ifdef NC_DOUBLE
200!       ierr = NF_GET_VAR_DOUBLE(nid,nvarid,volcapa)
201!#else
202!       ierr = NF_GET_VAR_REAL(nid,nvarid,volcapa)
203!#endif
204!        if (ierr.ne.NF_NOERR) then
205!         write(*,*)'soil_settings: Problem while reading <volcapa>'
206!         call abort
207!       endif
208!      endif
209
210! 3. Thermal inertia (note: it is declared in comsoil.h)
211! ------------------
212
213! 3.1 Look (again) for thermal inertia data (to reset nvarid)
214      ierr=NF_INQ_VARID(nid, "inertiedat", nvarid)
215      if (ierr.NE.NF_NOERR) then
216         write(*,*)'soil_settings: Field <inertiedat> not found!'
217         call abort
218      endif
219
220! 3.2 Knowing the # of dimensions <inertidat> was defined as using,
221!     read/build thermal inertia
222
223      if (ndims.eq.1) then ! "old 2D-surface" format
224       write(*,*)'soil_settings: Thermal inertia is only given as surfac
225     &e data!'
226       ! Read Surface thermal inertia
227       allocate(surfinertia(ngrid))
228#ifdef NC_DOUBLE
229       ierr = NF_GET_VAR_DOUBLE(nid, nvarid, surfinertia)
230#else
231       ierr = NF_GET_VAR_REAL(nid, nvarid, surfinertia)
232#endif
233        if (ierr.NE.NF_NOERR) then
234         write(*,*)'soil_settings: Problem while reading <inertiedat>'
235         call abort
236        endif
237       
238       write(*,*)' => Building soil thermal inertia (using reference sur
239     &face thermal inertia)'
240       do iloop=1,nsoil
241         inertiedat(:,iloop)=surfinertia(:)
242       enddo
243       deallocate(surfinertia)
244
245      else ! "3D surface+depth" format
246       if (interpol) then ! put values in oldinertiedat
247         if (.not.allocated(oldinertiedat)) then
248           allocate(oldinertiedat(ngrid,dimlen),stat=ierr)
249           if (ierr.ne.0) then
250            write(*,*) 'soil_settings: failed allocation of ',
251     &                 'oldinertiedat!'
252            stop
253           endif
254         endif
255#ifdef NC_DOUBLE
256        ierr = NF_GET_VAR_DOUBLE(nid,nvarid,oldinertiedat)
257#else
258        ierr = NF_GET_VAR_REAL(nid,nvarid,oldinertiedat)
259#endif
260        if (ierr.NE.NF_NOERR) then
261         write(*,*)'soil_settings: Problem while reading <inertiedat>'
262         call abort
263        endif
264       else ! put values in therm_i
265#ifdef NC_DOUBLE
266        ierr = NF_GET_VAR_DOUBLE(nid,nvarid,inertiedat)
267#else
268        ierr = NF_GET_VAR_REAL(nid,nvarid,inertiedat)
269#endif
270        if (ierr.NE.NF_NOERR) then
271         write(*,*)'soil_settings: Problem while reading <inertiedat>'
272         call abort
273        endif
274       endif ! of if (interpol)
275      endif ! of if (ndims.eq.1)
276     
277! 4. Read soil temperatures
278! -------------------------
279
280      ierr=NF_INQ_VARID(nid,"tsoil",nvarid)
281      if (ierr.ne.NF_NOERR) then
282        write(*,*)'soil_settings: Field <tsoil> not found!'
283        write(*,*)' => Building <tsoil> from surface values <tsurf>'
284        do iloop=1,nsoil
285          tsoil(:,iloop)=tsurf(:)
286        enddo
287      else ! <tsoil> found
288       if (interpol) then ! put values in oldtsoil
289         if (.not.allocated(oldtsoil)) then
290           allocate(oldtsoil(ngrid,dimlen),stat=ierr)
291           if (ierr.ne.0) then
292             write(*,*) 'soil_settings: failed allocation of ',
293     &                  'oldtsoil!'
294             stop
295           endif
296         endif
297#ifdef NC_DOUBLE
298        ierr=NF_GET_VAR_DOUBLE(nid,nvarid,oldtsoil)
299#else
300        ierr=NF_GET_VAR_REAL(nid,nvarid,oldtsoil)
301#endif
302        if (ierr.ne.NF_NOERR) then
303         write(*,*)'soil_settings: Problem while reading <tsoil>'
304         call abort
305        endif
306       else ! put values in tsoil
307#ifdef NC_DOUBLE
308        ierr=NF_GET_VAR_DOUBLE(nid,nvarid,tsoil)
309#else
310        ierr=NF_GET_VAR_REAL(nid,nvarid,tsoil)
311#endif
312        if (ierr.ne.NF_NOERR) then
313         write(*,*)'soil_settings: Problem while reading <tsoil>'
314         call abort
315        endif
316       endif ! of if (interpol)
317      endif
318
319! 5. If necessary, interpolate soil temperatures and thermal inertias
320! -------------------------------------------------------------------
321
322      if (olddepthdef) then
323      ! tsoil needs to be interpolated, but not therm_i
324        allocate(oldgrid(dimlen+1))
325        allocate(oldval(dimlen+1))
326        allocate(newval(nsoil))
327
328        do ig=1,ngrid
329          ! copy values
330          oldval(1)=tsurf(ig)
331          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen)
332          ! build vertical coordinate
333          oldgrid(1)=0. ! ground
334          oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)*
335     &                (inertiedat(ig,1)/volcapa)
336          ! interpolate
337          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
338          ! copy result in tsoil
339          tsoil(ig,:)=newval(:)
340        enddo
341
342        ! cleanup
343        deallocate(oldgrid)
344        deallocate(oldval)
345        deallocate(newval)
346        interpol=.false. ! no need for interpolation any more     
347      endif !of if (olddepthdef)
348
349      if (interpol) then
350      write(*,*)'soil_settings: Vertical interpolation along new grid'
351      ! interpolate soil temperatures and thermal inertias
352        if (.not.allocated(oldgrid)) then
353          allocate(oldgrid(dimlen+1))
354          allocate(oldval(dimlen+1))
355          allocate(newval(nsoil))
356        endif
357
358      ! thermal inertia
359        do ig=1,ngrid
360          ! copy data
361          oldval(1:dimlen)=oldinertiedat(ig,dimlen)
362          ! interpolate
363          call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil)
364          !copy result in inertiedat
365          inertiedat(ig,:)=newval(:)
366        enddo
367       
368      ! soil temperature
369        ! vertical coordinate
370        oldgrid(1)=0.0
371        oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)
372        do ig=1,ngrid
373          ! data
374          oldval(1)=tsurf(ig)
375          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen)
376          ! interpolate
377          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
378          ! copy result in inertiedat
379          tsoil(ig,:)=newval(:)
380        enddo
381       
382        !cleanup
383        deallocate(oldgrid)
384        deallocate(oldval)
385        deallocate(newval)
386        deallocate(oldinertiedat)
387        deallocate(oldtsoil)
388      endif ! of if (interpol)
389     
390! 6. Report min and max values of soil temperatures and thermal inertias
391! ----------------------------------------------------------------------
392
393      write(*,*)
394      write(*,*)'Soil volumetric heat capacity:',volcapa
395
396      xmin = MINVAL(inertiedat)
397      xmax = MAXVAL(inertiedat)
398      write(*,*)'Soil thermal inertia <inertiedat>:',xmin,xmax
399
400      xmin = 1.0E+20
401      xmax = -1.0E+20
402      xmin = MINVAL(tsoil)
403      xmax = MAXVAL(tsoil)
404      write(*,*)'Soil temperature <tsoil>:',xmin,xmax
405
406      end
Note: See TracBrowser for help on using the repository browser.