source: trunk/LMDZ.MARS/libf/phymars/soil_settings.F @ 1047

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

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

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