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

Last change on this file since 2829 was 2818, checked in by llange, 3 years ago

Log:
Mars GCM
Minor fix in soil_settings: oldmlayer was not allocated if dim(nlayer)
in the start != nsoil creating a seg.fault. The number of layers in the
soil can now be changed, as well as the power law.
LL

File size: 15.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      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!       ierr=nf90_inq_dimid(nid,"subsurface_layers",dimid)
88!       if (ierr.ne.nf90_noerr) then
89!        write(*,*)'soil_settings: Problem reading <subsurface_layers>'
90!         write(*,*)trim(nf90_strerror(ierr))
91!        call abort
92!       endif
93
94!       ierr=nf90_inquire_dimension(nid,dimid,len=dimlen)
95!       if (ierr.ne.nf90_noerr) then
96!        write(*,*)'soil_settings: Problem getting <subsurface_layers>',
97!     &             'length'
98!        write(*,*)trim(nf90_strerror(ierr))
99!         call abort
100!       endif
101        dimlen=inquire_dimension_length("subsurface_layers")
102
103        if (dimlen.ne.nsoil) then
104          write(*,*)'soil_settings: Interpolation of soil temperature ',
105     &              'and thermal inertia will be required!'
106        ! if dimlen doesn't match nsoil, then interpolation of
107        ! soil temperatures and thermal inertia will be requiered
108          interpol=.true.
109        ! allocate oldmlayer
110        if (.not.allocated(oldmlayer)) then
111          allocate(oldmlayer(dimlen),stat=ierr)
112          if (ierr.ne.0) then
113            write(*,*) 'soil_settings: failed allocation of oldmlayer!'
114            call abort_physic("soil_settings",
115     &           "failed oldmlayer allocation",1)
116          endif
117        endif
118        endif
119
120! 1.2 Find out the # of dimensions <inertiedat> was defined as using
121!     (in ye old days, thermal inertia was only given at the "surface")
122      ! Look for thermal inertia data
123!      ierr=nf90_inq_varid(nid,"inertiedat",nvarid)
124!      if (ierr.NE.nf90_noerr) then
125!         write(*,*)'soil_settings: Field <inertiedat> not found!'
126!         write(*,*)trim(nf90_strerror(ierr))
127!         call abort
128!      endif
129!
130!      ! Read the # of dimensions <inertidat> was defined as using
131!      ierr=nf90_inquire_variable(nid,nvarid,ndims=ndims)
132!      ! if (ndims.eq.1) then we have the "old 2D-surface" format
133      ndims=inquire_field_ndims("inertiedat")
134
135! 1.3 Read depths values or set olddepthdef flag and values
136      if (ndims.eq.1) then ! we know that there is none
137        write(*,*)'soil_settings: no <soildepth> field expected'
138        write(*,*)'building one mimicking old definitions'
139        olddepthdef=.true.
140        interpol=.true.
141
142        do iloop=1,dimlen
143          oldmlayer(iloop)=sqrt(887.75/3.14)*((2.**(iloop-0.5))-1.)
144        enddo
145      else ! Look for depth
146!        ierr=nf90_inq_varid(nid,"soildepth",nvarid)
147!        if (ierr.ne.nf90_noerr) then
148!          write(*,*)'soil_settings: Field <soildepth> not found!'
149!          write(*,*)trim(nf90_strerror(ierr))
150!         call abort
151!        endif
152        ! read <depth> coordinate
153        if (interpol) then !put values in oldmlayer
154!          ierr=nf90_get_var(nid,nvarid,oldmlayer)
155!          if (ierr.ne.nf90_noerr) then
156!           write(*,*)'soil_settings: Problem while reading <soildepth>'
157!           write(*,*)trim(nf90_strerror(ierr))
158!           call abort
159!          endif
160          call get_var("soildepth",oldmlayer,found)
161          if (.not.found) then
162            write(*,*)'soil_settings: Problem while reading <soildepth>'
163          endif
164        else ! put values in mlayer
165!          ierr=nf90_get_var(nid,nvarid,mlayer)
166!          if (ierr.ne.nf90_noerr) then
167!           write(*,*)'soil_settings: Problem while reading <soildepth>'
168!           write(*,*)trim(nf90_strerror(ierr))
169!           call abort
170!          endif
171          call get_var("soildepth",mlayer,found)
172          if (.not.found) then
173            write(*,*)'soil_settings: Problem while reading <soildepth>'
174          endif
175        endif !of if (interpol)
176      endif !of if (ndims.eq.1)
177
178! 1.4 Build mlayer(), if necessary
179      if (interpol) then
180      ! default mlayer distribution, following a power law:
181      !  mlayer(k)=lay1*alpha**(k-1/2)
182        lay1=2.e-4
183        alpha=2
184        do iloop=0,nsoil-1
185          mlayer(iloop)=lay1*(alpha**(iloop-0.5))
186        enddo
187      endif
188
189! 1.5 Build layer(); following the same law as mlayer()
190      ! Assuming layer distribution follows mid-layer law:
191      ! layer(k)=lay1*alpha**(k-1)
192      lay1=sqrt(mlayer(0)*mlayer(1))
193      alpha=mlayer(1)/mlayer(0)
194      do iloop=1,nsoil
195        layer(iloop)=lay1*(alpha**(iloop-1))
196      enddo
197
198! 2. Volumetric heat capacity (note: it is declared in comsoil_h)
199! ---------------------------
200! "volcapa" is (so far) 0D and written in "controle" table of startfi file
201! volcapa is read or set when "controle" is read (see tabfi.F)
202! Just in case, we check here that it is not zero. If it is, we
203! set it to "default_volcapa"
204
205      if (volcapa.le.0.0) then
206        write(*,*)'soil_settings: Warning, volcapa = ',volcapa
207        write(*,*)'               That doesn t seem right'
208        write(*,*)'        Initializing Volumetric heat capacity to ',
209     &             default_volcapa
210        volcapa=default_volcapa
211      endif
212! Look for it
213!      ierr=NF_INQ_VARID(nid,"volcapa",nvarid)
214!      if (ierr.NE.nf90_noerr) then
215!         write(*,*)'soil_settings: Field <volcapa> not found!'
216!         write(*,*)'Initializing Volumetric heat capacity to ',
217!     &             default_volcapa
218!         volcapa=default_volcapa
219!      else
220!#ifdef NC_DOUBLE
221!       ierr = NF_GET_VAR_DOUBLE(nid,nvarid,volcapa)
222!#else
223!       ierr = NF_GET_VAR_REAL(nid,nvarid,volcapa)
224!#endif
225!        if (ierr.ne.nf90_noerr) then
226!         write(*,*)'soil_settings: Problem while reading <volcapa>'
227!         call abort
228!       endif
229!      endif
230
231! 3. Thermal inertia (note: it is declared in comsoil_h)
232! ------------------
233
234! 3.1 Look (again) for thermal inertia data (to reset nvarid)
235!      ierr=nf90_inq_varid(nid,"inertiedat",nvarid)
236!      if (ierr.NE.nf90_noerr) then
237!         write(*,*)'soil_settings: Field <inertiedat> not found!'
238!         write(*,*)trim(nf90_strerror(ierr))
239!         call abort
240!      endif
241
242! 3.2 Knowing the # of dimensions <inertidat> was defined as using,
243!     read/build thermal inertia
244
245      if (ndims.eq.1) then ! "old 2D-surface" format
246       write(*,*)'soil_settings: Thermal inertia is only given as surfac
247     &e data!'
248       ! Read Surface thermal inertia
249       allocate(surfinertia(ngrid))
250!       ierr=nf90_get_var(nid,nvarid,surfinertia)
251!        if (ierr.NE.nf90_noerr) then
252!         write(*,*)'soil_settings: Problem while reading <inertiedat>'
253!         write(*,*)trim(nf90_strerror(ierr))
254!         call abort
255!        endif
256       call get_field("inertiedat",surfinertia,found)
257       if (.not.found) then
258         write(*,*) "soil_settings: Failed loading <inertiedat>"
259         call abort_physic("soil_settings",
260     &        "failed loading <inertiedat>",1)
261       endif
262       
263       write(*,*)' => Building soil thermal inertia (using reference sur
264     &face thermal inertia)'
265       do iloop=1,nsoil
266         inertiedat(:,iloop)=surfinertia(:)
267       enddo
268       deallocate(surfinertia)
269
270      else ! "3D surface+depth" format
271       if (interpol) then ! put values in oldinertiedat
272         if (.not.allocated(oldinertiedat)) then
273           allocate(oldinertiedat(ngrid,dimlen),stat=ierr)
274           if (ierr.ne.0) then
275            write(*,*) 'soil_settings: failed allocation of ',
276     &                 'oldinertiedat!'
277            call abort_physic("soil_settings",
278     &        "failed allocation of oldinertiedat",1)
279           endif
280         endif ! of if (.not.allocated(oldinertiedat))
281!         ierr=nf90_get_var(nid,nvarid,oldinertiedat)
282!        if (ierr.NE.nf90_noerr) then
283!         write(*,*)'soil_settings: Problem while reading <inertiedat>'
284!         write(*,*)trim(nf90_strerror(ierr))
285!         call abort
286!        endif
287        call get_field("inertiedat",oldinertiedat,found)
288        if (.not.found) then
289          write(*,*) "soil_settings: Failed loading <inertiedat>"
290         call abort_physic("soil_settings",
291     &        "failed loading <inertiedat>",1)
292        endif
293       else ! put values in therm_i
294!        ierr=nf90_get_var(nid,nvarid,inertiedat)
295!        if (ierr.NE.nf90_noerr) then
296!         write(*,*)'soil_settings: Problem while reading <inertiedat>'
297!         write(*,*)trim(nf90_strerror(ierr))
298!         call abort
299         call get_field("inertiedat",inertiedat,found)
300         if (.not.found) then
301           write(*,*) "soil_settings: Failed loading <inertiedat>"
302           call abort_physic("soil_settings",
303     &        "failed loading <inertiedat>",1)
304         endif
305!        endif
306       endif ! of if (interpol)
307      endif ! of if (ndims.eq.1)
308     
309! 4. Read soil temperatures
310! -------------------------
311
312!      ierr=nf90_inq_varid(nid,"tsoil",nvarid)
313      ok=inquire_field("tsoil")
314!      if (ierr.ne.nf90_noerr) then
315      if (.not.ok) then
316        write(*,*)'soil_settings: Field <tsoil> not found!'
317        write(*,*)' => Building <tsoil> from surface values <tsurf>'
318        do iloop=1,nsoil
319          tsoil(:,iloop)=tsurf(:)
320        enddo
321      else ! <tsoil> found
322       if (interpol) then ! put values in oldtsoil
323         if (.not.allocated(oldtsoil)) then
324           allocate(oldtsoil(ngrid,dimlen),stat=ierr)
325           if (ierr.ne.0) then
326             write(*,*) 'soil_settings: failed allocation of ',
327     &                  'oldtsoil!'
328             call abort_physic("soil_settings",
329     &        "failed allocation of oldtsoil",1)
330           endif
331         endif
332!        ierr=nf90_get_var(nid,nvarid,oldtsoil)
333!        if (ierr.ne.nf90_noerr) then
334!        write(*,*)'soil_settings: Problem while reading <tsoil>'
335!         write(*,*)trim(nf90_strerror(ierr))
336!        call abort
337!       endif
338         call get_field("tsoil",oldtsoil,found)
339         if (.not.found) then
340           write(*,*) "soil_settings: Failed loading <tsoil>"
341           call abort_physic("soil_settings",
342     &          "failed loading <tsoil>",1)
343         endif
344       else ! put values in tsoil
345!        corner(1)=1
346!        corner(2)=1
347!        corner(3)=indextime
348!        edges(1)=ngrid
349!        edges(2)=nsoil
350!        edges(3)=1
351!        !ierr=nf90_get_var(nid,nvarid,tsoil,corner,edges)
352!        ierr=nf90_get_var(nid,nvarid,tsoil)
353!        if (ierr.ne.nf90_noerr) then
354!        write(*,*)'soil_settings: Problem while reading <tsoil>'
355!         write(*,*)trim(nf90_strerror(ierr))
356!        call abort
357!       endif
358         call get_field("tsoil",tsoil,found,timeindex=indextime)
359         if (.not.found) then
360           write(*,*) "soil_settings: Failed loading <tsoil>"
361           call abort_physic("soil_settings",
362     &          "failed loading <tsoil>",1)
363         endif
364       endif ! of if (interpol)
365      endif! of if (.not.ok)
366
367! 5. If necessary, interpolate soil temperatures and thermal inertias
368! -------------------------------------------------------------------
369
370      if (olddepthdef) then
371      ! tsoil needs to be interpolated, but not therm_i
372        allocate(oldgrid(dimlen+1))
373        allocate(oldval(dimlen+1))
374        allocate(newval(nsoil))
375
376        do ig=1,ngrid
377          ! copy values
378          oldval(1)=tsurf(ig)
379          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen)
380          ! build vertical coordinate
381          oldgrid(1)=0. ! ground
382          oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)*
383     &                (inertiedat(ig,1)/volcapa)
384          ! interpolate
385          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
386          ! copy result in tsoil
387          tsoil(ig,:)=newval(:)
388        enddo
389
390        ! cleanup
391        deallocate(oldgrid)
392        deallocate(oldval)
393        deallocate(newval)
394        interpol=.false. ! no need for interpolation any more     
395      endif !of if (olddepthdef)
396
397      if (interpol) then
398      write(*,*)'soil_settings: Vertical interpolation along new grid'
399      ! interpolate soil temperatures and thermal inertias
400        if (.not.allocated(oldgrid)) then
401          allocate(oldgrid(dimlen+1))
402          allocate(oldval(dimlen+1))
403          allocate(newval(nsoil))
404        endif
405
406      ! thermal inertia
407        do ig=1,ngrid
408          ! copy data
409          oldval(1:dimlen)=oldinertiedat(ig,dimlen)
410          ! interpolate
411          call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil)
412          !copy result in inertiedat
413          inertiedat(ig,:)=newval(:)
414        enddo
415       
416      ! soil temperature
417        ! vertical coordinate
418        oldgrid(1)=0.0
419        oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)
420        do ig=1,ngrid
421          ! data
422          oldval(1)=tsurf(ig)
423          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen)
424          ! interpolate
425          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
426          ! copy result in inertiedat
427          tsoil(ig,:)=newval(:)
428        enddo
429       
430        !cleanup
431        deallocate(oldgrid)
432        deallocate(oldval)
433        deallocate(newval)
434        deallocate(oldinertiedat)
435        deallocate(oldtsoil)
436      endif ! of if (interpol)
437     
438! 6. Report min and max values of soil temperatures and thermal inertias
439! ----------------------------------------------------------------------
440
441      write(*,*)
442      write(*,*)'Soil volumetric heat capacity:',volcapa
443
444      xmin = MINVAL(inertiedat)
445      xmax = MAXVAL(inertiedat)
446      write(*,*)'Soil thermal inertia <inertiedat>:',xmin,xmax
447
448      xmin = 1.0E+20
449      xmax = -1.0E+20
450      xmin = MINVAL(tsoil)
451      xmax = MAXVAL(tsoil)
452      write(*,*)'Soil temperature <tsoil>:',xmin,xmax
453
454      end
Note: See TracBrowser for help on using the repository browser.