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

Last change on this file since 2616 was 2311, checked in by emillour, 5 years ago

Mars GCM:
Code tidying: use getin_p() instead of getin() and use "call abort_physic"
instead of "stop" or "call abort".
EM

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