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

Last change on this file since 1217 was 1216, checked in by emillour, 11 years ago

Generic model:
Major cleanup, in order to ease the use of LMDZ.GENERIC with (parallel) dynamics
in LMDZ.COMMON: (NB: this will break LMDZ.UNIVERSAL, which should be thrashed
in the near future)

  • Updated makegcm_* scripts (and makdim) and added the "-full" (to enforce full recomputation of the model) option
  • In dyn3d: converted control.h to module control_mod.F90 and converted iniadvtrac.F to module infotrac.F90
  • Added module mod_const_mpi.F90 in dyn3d (not used in serial mode)
  • Rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallelism)
  • added created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the max and min of a field over the whole planet. This should be further imroved with computation of means (possibly area weighed), etc.

EM

File size: 11.5 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     
68      ! for interpolation
69      real,dimension(:),allocatable :: oldgrid
70      real,dimension(:),allocatable :: oldval
71      real,dimension(:),allocatable :: newval
72
73      real alpha,lay1 ! 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!======================================================================
81
82! 1. Depth coordinate
83! -------------------
84
85! 1.1 Start by reading how many layers of soil there are
86
87        dimlen=inquire_dimension_length("subsurface_layers")
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
99      ndims=inquire_field_ndims("inertiedat")
100
101! 1.3 Read depths values or set olddepthdef flag and values
102      if (ndims.eq.1) then ! we know that there is none
103        write(*,*)'soil_settings: no <soildepth> field expected'
104        write(*,*)'building one mimicking old definitions'
105        olddepthdef=.true.
106        interpol=.true.
107        ! allocate oldmlayer
108        if (.not.allocated(oldmlayer)) then
109          allocate(oldmlayer(dimlen),stat=ierr)
110          if (ierr.ne.0) then
111            write(*,*) 'soil_settings: failed allocation of oldmlayer!'
112            stop
113          endif
114        endif
115        do iloop=1,dimlen
116          oldmlayer(iloop)=sqrt(887.75/3.14)*((2.**(iloop-0.5))-1.)
117        enddo
118      else ! Look for depth
119        ! read <depth> coordinate
120        if (interpol) then !put values in oldmlayer
121          call get_var("soildepth",oldmlayer,found)
122          if (.not.found) then
123            write(*,*)'soil_settings: Problem while reading <soildepth>'
124          endif
125        else ! put values in mlayer
126          call get_var("soildepth",mlayer,found)
127          if (.not.found) then
128            write(*,*)'soil_settings: Problem while reading <soildepth>'
129          endif
130        endif !of if (interpol)
131      endif !of if (ndims.eq.1)
132
133! 1.4 Build mlayer(), if necessary
134      if (interpol) then
135      ! default mlayer distribution, following a power law:
136      !  mlayer(k)=lay1*alpha**(k-1/2)
137        lay1=2.e-4
138        alpha=2
139        do iloop=0,nsoil-1
140          mlayer(iloop)=lay1*(alpha**(iloop-0.5))
141        enddo
142      endif
143
144! 1.5 Build layer(); following the same law as mlayer()
145      ! Assuming layer distribution follows mid-layer law:
146      ! layer(k)=lay1*alpha**(k-1)
147      lay1=sqrt(mlayer(0)*mlayer(1))
148      alpha=mlayer(1)/mlayer(0)
149      do iloop=1,nsoil
150        layer(iloop)=lay1*(alpha**(iloop-1))
151      enddo
152
153! 2. Volumetric heat capacity (note: it is declared in comsoil_h)
154! ---------------------------
155! "volcapa" is (so far) 0D and written in "controle" table of startfi file
156! volcapa is read or set when "controle" is read (see tabfi.F)
157! Just in case, we check here that it is not zero. If it is, we
158! set it to "default_volcapa"
159
160      if (volcapa.le.0.0) then
161        write(*,*)'soil_settings: Warning, volcapa = ',volcapa
162        write(*,*)'               That doesn t seem right'
163        write(*,*)'        Initializing Volumetric heat capacity to ',
164     &             default_volcapa
165        volcapa=default_volcapa
166      endif
167
168! 3. Thermal inertia (note: it is declared in comsoil_h)
169! ------------------
170
171! 3.1 Look (again) for thermal inertia data (to reset nvarid)
172
173! 3.2 Knowing the # of dimensions <inertidat> was defined as using,
174!     read/build thermal inertia
175
176      if (ndims.eq.1) then ! "old 2D-surface" format
177       write(*,*)'soil_settings: Thermal inertia is only given as surfac
178     &e data!'
179       ! Read Surface thermal inertia
180       allocate(surfinertia(ngrid))
181       call get_field("inertiedat",surfinertia,found)
182       if (.not.found) then
183         write(*,*) "soil_settings: Failed loading <inertiedat>"
184         call abort
185       endif
186       
187       write(*,*)' => Building soil thermal inertia (using reference sur
188     &face thermal inertia)'
189       do iloop=1,nsoil
190         inertiedat(:,iloop)=surfinertia(:)
191       enddo
192       deallocate(surfinertia)
193
194      else ! "3D surface+depth" format
195       if (interpol) then ! put values in oldinertiedat
196         if (.not.allocated(oldinertiedat)) then
197           allocate(oldinertiedat(ngrid,dimlen),stat=ierr)
198           if (ierr.ne.0) then
199            write(*,*) 'soil_settings: failed allocation of ',
200     &                 'oldinertiedat!'
201            stop
202           endif
203         endif ! of if (.not.allocated(oldinertiedat))
204        call get_field("inertiedat",oldinertiedat,found)
205        if (.not.found) then
206          write(*,*) "soil_settings: Failed loading <inertiedat>"
207          call abort
208        endif
209       else ! put values in therm_i
210         call get_field("inertiedat",inertiedat,found)
211         if (.not.found) then
212           write(*,*) "soil_settings: Failed loading <inertiedat>"
213           call abort
214         endif
215!        endif
216       endif ! of if (interpol)
217      endif ! of if (ndims.eq.1)
218     
219! 4. Read soil temperatures
220! -------------------------
221
222!      ierr=nf90_inq_varid(nid,"tsoil",nvarid)
223      ok=inquire_field("tsoil")
224!      if (ierr.ne.nf90_noerr) then
225      if (.not.ok) then
226        write(*,*)'soil_settings: Field <tsoil> not found!'
227        write(*,*)' => Building <tsoil> from surface values <tsurf>'
228        do iloop=1,nsoil
229          tsoil(:,iloop)=tsurf(:)
230        enddo
231      else ! <tsoil> found
232       if (interpol) then ! put values in oldtsoil
233         if (.not.allocated(oldtsoil)) then
234           allocate(oldtsoil(ngrid,dimlen),stat=ierr)
235           if (ierr.ne.0) then
236             write(*,*) 'soil_settings: failed allocation of ',
237     &                  'oldtsoil!'
238             stop
239           endif
240         endif
241         call get_field("tsoil",oldtsoil,found)
242         if (.not.found) then
243           write(*,*) "soil_settings: Failed loading <tsoil>"
244           call abort
245         endif
246       else ! put values in tsoil
247         call get_field("tsoil",tsoil,found,timeindex=indextime)
248         if (.not.found) then
249           write(*,*) "soil_settings: Failed loading <tsoil>"
250           call abort
251         endif
252       endif ! of if (interpol)
253      endif! of if (.not.ok)
254
255! 5. If necessary, interpolate soil temperatures and thermal inertias
256! -------------------------------------------------------------------
257
258      if (olddepthdef) then
259      ! tsoil needs to be interpolated, but not therm_i
260        allocate(oldgrid(dimlen+1))
261        allocate(oldval(dimlen+1))
262        allocate(newval(nsoil))
263
264        do ig=1,ngrid
265          ! copy values
266          oldval(1)=tsurf(ig)
267          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen)
268          ! build vertical coordinate
269          oldgrid(1)=0. ! ground
270          oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)*
271     &                (inertiedat(ig,1)/volcapa)
272          ! interpolate
273          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
274          ! copy result in tsoil
275          tsoil(ig,:)=newval(:)
276        enddo
277
278        ! cleanup
279        deallocate(oldgrid)
280        deallocate(oldval)
281        deallocate(newval)
282        interpol=.false. ! no need for interpolation any more     
283      endif !of if (olddepthdef)
284
285      if (interpol) then
286      write(*,*)'soil_settings: Vertical interpolation along new grid'
287      ! interpolate soil temperatures and thermal inertias
288        if (.not.allocated(oldgrid)) then
289          allocate(oldgrid(dimlen+1))
290          allocate(oldval(dimlen+1))
291          allocate(newval(nsoil))
292        endif
293
294      ! thermal inertia
295        do ig=1,ngrid
296          ! copy data
297          oldval(1:dimlen)=oldinertiedat(ig,dimlen)
298          ! interpolate
299          call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil)
300          !copy result in inertiedat
301          inertiedat(ig,:)=newval(:)
302        enddo
303       
304      ! soil temperature
305        ! vertical coordinate
306        oldgrid(1)=0.0
307        oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)
308        do ig=1,ngrid
309          ! data
310          oldval(1)=tsurf(ig)
311          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen)
312          ! interpolate
313          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
314          ! copy result in inertiedat
315          tsoil(ig,:)=newval(:)
316        enddo
317       
318        !cleanup
319        deallocate(oldgrid)
320        deallocate(oldval)
321        deallocate(newval)
322        deallocate(oldinertiedat)
323        deallocate(oldtsoil)
324      endif ! of if (interpol)
325     
326! 6. Report min and max values of soil temperatures and thermal inertias
327! ----------------------------------------------------------------------
328
329      write(*,*)
330      write(*,*)'Soil volumetric heat capacity:',volcapa
331
332      xmin = MINVAL(inertiedat)
333      xmax = MAXVAL(inertiedat)
334      write(*,*)'Soil thermal inertia <inertiedat>:',xmin,xmax
335
336      xmin = 1.0E+20
337      xmax = -1.0E+20
338      xmin = MINVAL(tsoil)
339      xmax = MAXVAL(tsoil)
340      write(*,*)'Soil temperature <tsoil>:',xmin,xmax
341
342      end
Note: See TracBrowser for help on using the repository browser.