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

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

Mars GCM:
Series of changes to enable running in parallel (using LMDZ.COMMON dynamics);
Current LMDZ.MARS can still notheless be compiled and run in serial mode
"as previously".
Summary of main changes:

  • Main programs (newstart, start2archive, xvik) that used to be in dyn3d have been moved to phymars.
  • dyn3d/control.h is now module control_mod.F90
  • 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 parallel case)
  • created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the min and max of a field over the whole planet.

EM

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