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

Last change on this file since 1009 was 999, checked in by tnavarro, 11 years ago

Possibility to store multiple initial states in one start/startfi. This is RETROCOMPATIBLE. New option ecrithist in run.def to write data in start/startfi every ecrithist dynamical timestep. New option timestart in run.def to initialize the GCM with the time timestart stored in start

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