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

Last change on this file was 3126, checked in by llange, 12 months ago

Mars PCM

  • Update in soilwater: adding the possibility to run without adsorption, but with the possibility to run with seasonal frost forming in the subsurface
  • THe choice of the isotherm for adsorption can be now done by setting the integer choice_ads in the callphys.def choice_ads = 1 adsorption rate is computed with the H2O thermal speed; choice_ads = 2 adsorption rate is computed based on exeperimental resutls, choice_ads =3 no adsorption

LL

File size: 18.7 KB
Line 
1      subroutine soil_settings(nid,ngrid,nsoil,nqsoil,tsurf,tsoil,
2     &                         qsoil,indextime)
3
4!      use netcdf
5      use comsoil_h, only: layer, mlayer, inertiedat, inertiesoil,
6     &                     volcapa,flux_geo,adsorption_soil,
7     &                     igcm_h2o_vap_soil,igcm_h2o_ice_soil,
8     &                     igcm_h2o_vap_ads
9      use iostart, only: inquire_field_ndims, get_var, get_field,
10     &                   inquire_field, inquire_dimension_length
11      USE comslope_mod, ONLY: nslope
12      implicit none
13
14!======================================================================
15!  Author: Ehouarn Millour (07/2006)
16!
17!  Purpose: Read and/or initialise soil depths and properties
18!
19! Modifications: Aug.2010 EM : use NetCDF90 to load variables (enables using
20!                      r4 or r8 restarts independently of having compiled
21!                      the GCM in r4 or r8)
22!                June 2013 TN : Possibility to read files with a time axis
23!
24!
25!  This subroutine reads from a NetCDF file (opened by the caller)
26!  of "startfi.nc" format.
27!  The various actions and variable read/initialized are:
28!  1. Check out the number of soil layers (in datafile); if it isn't equal
29!     to nsoil, then some interpolation will be required
30!     Also check if data in file "startfi.nc" is in older format (ie:
31!     thermal inertia was depth-independent; and there was no "depth"
32!     coordinate.
33!     Read/build layer (and midlayer) depths
34!  2. Read volumetric specific heat (or initialise it to default value)
35!  3. Read Thermal inertia
36!  4. Read soil temperatures
37!  5. Interpolate thermal inertia and temperature on the new grid, if
38!     necessary
39!======================================================================
40
41!======================================================================
42!  arguments
43!  ---------
44!  inputs:
45      integer,intent(in) :: nid ! Input Netcdf file ID
46      integer,intent(in) :: ngrid       ! # of horizontal grid points
47      integer,intent(in) :: nsoil       ! # of soil layers
48      integer,intent(in) :: nqsoil      ! # of tracers in the soil
49      real,intent(in) :: tsurf(ngrid,nslope)   ! surface temperature [K]
50      integer,intent(in) :: indextime   ! position on time axis
51!  output:
52      real,intent(out) :: tsoil(ngrid,nsoil,nslope)     ! soil temperature [K]
53      real,intent(out) :: qsoil(ngrid,nsoil,nqsoil,nslope) ! Tracers in the subsurface [kg/kg]
54
55!======================================================================
56! local variables:
57      integer ierr      ! status (returned by NetCDF functions)
58      integer nvarid    ! ID of NetCDF variable
59      integer dimid     ! ID of NetCDF dimension
60      integer dimlen    ! length along the "depth" dimension
61      integer ndims     ! # of dimensions of read <inertiedat> data
62      integer ig,iloop  ! loop counters
63      integer islope
64     
65      integer edges(3),corner(3) ! to read a specific time
66
67      logical :: olddepthdef=.false. ! flag
68      logical :: interpol=.false. ! flag: true if interpolation will be requiered
69
70      ! to store "old" values
71      real,dimension(:),allocatable :: surfinertia !surface thermal inertia
72      real,dimension(:),allocatable :: oldmlayer
73      real,dimension(:,:),allocatable :: oldinertiedat
74      real,dimension(:,:,:),allocatable :: oldinertiesoil
75      real,dimension(:,:,:),allocatable :: oldtsoil
76     
77      ! for interpolation
78      real,dimension(:),allocatable :: oldgrid
79      real,dimension(:),allocatable :: oldval
80      real,dimension(:),allocatable :: newval
81
82      real alpha,lay1 ! coefficients for building layers
83      real xmin,xmax ! to display min and max of a field
84
85      real,parameter :: default_volcapa=1.e6
86
87      logical :: found,ok
88
89      character (len=30):: txt
90     
91!======================================================================
92
93! 1. Depth coordinate
94! -------------------
95
96
97! 1.1 Start by reading how many layers of soil there are
98
99        dimlen=inquire_dimension_length("subsurface_layers")
100
101        if (dimlen.ne.nsoil) then
102          write(*,*)'soil_settings: Interpolation of soil temperature ',
103     &              'and thermal inertia will be required!'
104        ! if dimlen doesn't match nsoil, then interpolation of
105        ! soil temperatures and thermal inertia will be requiered
106          interpol=.true.
107        ! allocate oldmlayer
108        if (.not.allocated(oldmlayer)) then
109          allocate(oldmlayer(dimlen),stat=ierr)
110          if (ierr.ne.0) then
111            write(*,*) 'soil_settings: failed allocation of oldmlayer!'
112            call abort_physic("soil_settings",
113     &           "failed oldmlayer allocation",1)
114          endif
115        endif
116        endif
117
118! 1.2 Find out the # of dimensions <inertiedat> was defined as using
119!     (in ye old days, thermal inertia was only given at the "surface")
120
121      ndims=inquire_field_ndims("inertiedat")
122
123! 1.3 Read depths values or set olddepthdef flag and values
124      if (ndims.eq.1) then ! we know that there is none
125        write(*,*)'soil_settings: no <soildepth> field expected'
126        write(*,*)'building one mimicking old definitions'
127        olddepthdef=.true.
128        interpol=.true.
129
130        do iloop=1,dimlen
131          oldmlayer(iloop)=sqrt(887.75/3.14)*((2.**(iloop-0.5))-1.)
132        enddo
133      else ! Look for depth
134        ! read <depth> coordinate
135        if (interpol) then !put values in oldmlayer
136          call get_var("soildepth",oldmlayer,found)
137          if (.not.found) then
138            write(*,*)'soil_settings: Problem while reading <soildepth>'
139          endif
140        else ! put values in mlayer
141          call get_var("soildepth",mlayer,found)
142          if (.not.found) then
143            write(*,*)'soil_settings: Problem while reading <soildepth>'
144          endif
145        endif !of if (interpol)
146      endif !of if (ndims.eq.1)
147
148! 1.4 Build mlayer(), if necessary
149      if (interpol) then
150      !  mlayer(k)=lay1*(1+k^(2.9)*(1-exp(-k/20)), PYM grid
151        lay1=2.e-4
152        do iloop=0,nsoil-1
153          mlayer(iloop) = lay1*(1+iloop**2.9*(1-exp(-real(iloop)/20.)))
154        enddo
155      endif
156
157! 1.5 Build layer();
158      do iloop=1,nsoil-1
159        layer(iloop)=(mlayer(iloop)+mlayer(iloop-1))/2
160      enddo
161      layer(nsoil)=2*mlayer(nsoil-1) - mlayer(nsoil-2)
162
163! 2. Volumetric heat capacity (note: it is declared in comsoil_h)
164! ---------------------------
165! "volcapa" is (so far) 0D and written in "controle" table of startfi file
166! volcapa is read or set when "controle" is read (see tabfi.F)
167! Just in case, we check here that it is not zero. If it is, we
168! set it to "default_volcapa"
169
170      if (volcapa.le.0.0) then
171        write(*,*)'soil_settings: Warning, volcapa = ',volcapa
172        write(*,*)'               That doesn t seem right'
173        write(*,*)'        Initializing Volumetric heat capacity to ',
174     &             default_volcapa
175        volcapa=default_volcapa
176      endif
177
178! 3. Thermal inertia for present day climate -inertiedat- & the one given by the pem -ineritesoil-(note: it is declared in comsoil_h)
179! ------------------
180!     Knowing the # of dimensions <inertidat> was defined as using,
181!     read/build thermal inertia
182
183! 3.1 Present day - inertiedat
184
185      if (ndims.eq.1) then ! "old 2D-surface" format
186       write(*,*)'soil_settings: Thermal inertia is only given as surfac
187     &e data!'
188       ! Read Surface thermal inertia
189       allocate(surfinertia(ngrid))
190!       ierr=nf90_get_var(nid,nvarid,surfinertia)
191!        if (ierr.NE.nf90_noerr) then
192!         write(*,*)'soil_settings: Problem while reading <inertiedat>'
193!         write(*,*)trim(nf90_strerror(ierr))
194!         call abort
195!        endif
196       call get_field("inertiedat",surfinertia,found)
197       if (.not.found) then
198         write(*,*) "soil_settings: Failed loading <inertiedat>"
199         call abort_physic("soil_settings",
200     &        "failed loading <inertiedat>",1)
201       endif
202       
203       write(*,*)' => Building soil thermal inertia (using reference sur
204     &face thermal inertia)'
205       do iloop=1,nsoil
206         inertiedat(:,iloop)=surfinertia(:)
207       enddo
208       deallocate(surfinertia)
209
210      else ! "3D surface+depth" format
211       if (interpol) then ! put values in oldinertiedat
212         if (.not.allocated(oldinertiedat)) then
213           allocate(oldinertiedat(ngrid,dimlen),stat=ierr)
214           if (ierr.ne.0) then
215            write(*,*) 'soil_settings: failed allocation of ',
216     &                 'oldinertiedat!'
217            call abort_physic("soil_settings",
218     &        "failed allocation of oldinertiedat",1)
219           endif
220         endif ! of if (.not.allocated(oldinertiedat))
221
222        call get_field("inertiedat",oldinertiedat,found)
223        if (.not.found) then
224          write(*,*) "soil_settings: Failed loading <inertiedat>"
225         call abort_physic("soil_settings",
226     &        "failed loading <inertiedat>",1)
227        endif
228       else ! put values in therm_i
229         call get_field("inertiedat",inertiedat,found)
230         if (.not.found) then
231           write(*,*) "soil_settings: Failed loading <inertiedat>"
232           call abort_physic("soil_settings",
233     &        "failed loading <inertiedat>",1)
234         endif
235       endif ! of if (interpol)
236      endif ! of if (ndims.eq.1)
237
238! 3.2 Inertie from the PEM
239
240      ok=inquire_field("inertiesoil")
241
242      if (.not.ok) then
243        write(*,*)'soil_settings: Field <inertiesoil> not found!'
244        write(*,*)' => Building <inertiesoil> from surface values'//
245     & ' <inertiedat>'
246! Case 1: No interpolation needed, we just copy past inertiedat
247        if(.not.(interpol)) then
248           do islope=1,nslope
249              inertiesoil(:,:,islope)=inertiedat(:,:)
250           enddo
251        else
252! Case 2: Interpolation needed: we copy past old value from inertiedat
253           if (.not.allocated(oldinertiesoil)) then
254              allocate(oldinertiesoil(ngrid,dimlen,nslope),stat=ierr)
255           endif
256           do islope=1,nslope
257              oldinertiesoil(:,:,islope)=oldinertiedat(:,:)
258           enddo
259        endif
260      else ! <inertiesoil> found
261       if (interpol) then ! put values in oldinertiesoil
262         if (.not.allocated(oldinertiesoil)) then
263           allocate(oldinertiesoil(ngrid,dimlen,nslope),stat=ierr)
264           if (ierr.ne.0) then
265             write(*,*) 'soil_settings: failed allocation of ',
266     &                  'oldinertiesoil!'
267             call abort_physic("soil_settings",
268     &        "failed allocation of oldinertiesoil",1)
269           endif
270         endif
271
272         call get_field("inertiesoil",oldinertiesoil,found)
273         if (.not.found) then
274           write(*,*) "soil_settings: Failed loading <inertiesoil>"
275           call abort_physic("soil_settings",
276     &          "failed loading <inertiesoil>",1)
277         endif
278       else ! put values in inertiesoil
279         call get_field("inertiesoil",inertiesoil,found,
280     &                 timeindex=indextime)
281         if (.not.found) then
282           write(*,*) "soil_settings: Failed loading <inertiesoil>"
283           call abort_physic("soil_settings",
284     &          "failed loading <inertiesoil>",1)
285         endif
286       endif ! of if (interpol)
287      endif! of if (.not.ok)
288
289
290     
291! 4. Read soil temperatures
292! -------------------------
293
294      ok=inquire_field("tsoil")
295
296      if (.not.ok) then
297        write(*,*)'soil_settings: Field <tsoil> not found!'
298        write(*,*)' => Building <tsoil> from surface values <tsurf>'
299! Case 1: No interpolation needed, we just copy past inertiedat
300        if(.not.(interpol)) then
301          do iloop=1,nsoil
302            do islope=1,nslope
303              tsoil(:,iloop,islope)=tsurf(:,islope)
304            enddo
305          enddo
306        else
307! Case 2: Interpolation needed: we copy past old value from inertiedat
308          if (.not.allocated(oldtsoil)) then
309             allocate(oldtsoil(ngrid,dimlen,nslope),stat=ierr)
310          endif
311          do iloop=1,dimlen
312            do islope=1,nslope
313              oldtsoil(:,iloop,islope)=tsurf(:,islope)
314            enddo
315          enddo
316        endif
317      else ! <tsoil> found
318       if (interpol) then ! put values in oldtsoil
319         if (.not.allocated(oldtsoil)) then
320           allocate(oldtsoil(ngrid,dimlen,nslope),stat=ierr)
321           if (ierr.ne.0) then
322             write(*,*) 'soil_settings: failed allocation of ',
323     &                  'oldtsoil!'
324             call abort_physic("soil_settings",
325     &        "failed allocation of oldtsoil",1)
326           endif
327         endif
328         call get_field("tsoil",oldtsoil,found)
329         if (.not.found) then
330           write(*,*) "soil_settings: Failed loading <tsoil>"
331           call abort_physic("soil_settings",
332     &          "failed loading <tsoil>",1)
333         endif
334       else ! put values in tsoil
335         call get_field("tsoil",tsoil,found,timeindex=indextime)
336         if (.not.found) then
337           write(*,*) "soil_settings: Failed loading <tsoil>"
338           call abort_physic("soil_settings",
339     &          "failed loading <tsoil>",1)
340         endif
341       endif ! of if (interpol)
342      endif! of if (.not.ok)
343
344
345! 5. If necessary, interpolate soil temperatures and thermal inertias
346! -------------------------------------------------------------------
347
348      if (olddepthdef) then
349      ! tsoil needs to be interpolated, but not therm_i
350        allocate(oldgrid(dimlen+1))
351        allocate(oldval(dimlen+1))
352            allocate(newval(nsoil))
353        do ig=1,ngrid
354          ! copy values
355          do islope=1,nslope
356          oldval(1)=tsurf(ig,islope)
357          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen,islope)
358          ! build vertical coordinate
359          oldgrid(1)=0. ! ground
360          oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)*
361     &                (inertiesoil(ig,1,islope)/volcapa)
362          ! interpolate
363          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
364          ! copy result in tsoil
365          tsoil(ig,:,islope)=newval(:)
366          enddo
367          enddo
368        ! cleanup
369          deallocate(oldgrid)
370          deallocate(oldval)
371          deallocate(newval)
372          interpol=.false. ! no need for interpolation any more     
373      endif !of if (olddepthdef)
374      if (interpol) then
375      write(*,*)'soil_settings: Vertical interpolation along new grid'
376      ! interpolate soil temperatures and thermal inertias
377      if (.not.allocated(oldgrid)) then
378          allocate(oldgrid(dimlen+1))
379          allocate(oldval(dimlen+1))
380          allocate(newval(nsoil))
381      endif
382
383      ! thermal inertia - present day (inertiedat)
384      do ig=1,ngrid
385          ! copy data
386          oldval(1:dimlen)=oldinertiedat(ig,1:dimlen)
387          ! interpolate
388          call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil)
389          !copy result in inertiedat
390          inertiedat(ig,:)=newval(:)
391      ! thermal inertia - general case with PEM (inertie soil)
392          do islope=1,nslope
393          ! copy data
394             oldval(1:dimlen)=oldinertiesoil(ig,1:dimlen,islope)
395          ! interpolate
396             call interp_line(oldmlayer,oldval,dimlen,mlayer,newval,nsoil)
397          !copy result in inertiedat
398             inertiesoil(ig,:,islope)=newval(:)
399           enddo
400        enddo ! ig
401       
402      ! soil temperature
403        ! vertical coordinate
404          oldgrid(1)=0.0
405          oldgrid(2:dimlen+1)=oldmlayer(1:dimlen)
406        do ig=1,ngrid
407          do islope=1,nslope
408          ! data
409          oldval(1)=tsurf(ig,islope)
410          oldval(2:dimlen+1)=oldtsoil(ig,1:dimlen,islope)
411          ! interpolate
412          call interp_line(oldgrid,oldval,dimlen+1,mlayer,newval,nsoil)
413          ! copy result in inertiedat
414          tsoil(ig,:,islope)=newval(:)
415          enddo
416         enddo
417       
418        !cleanup
419          deallocate(oldgrid)
420          deallocate(oldval)
421          deallocate(newval)
422          deallocate(oldinertiedat)
423          deallocate(oldinertiesoil)
424          deallocate(oldtsoil)
425      endif ! of if (interpol)
426
427
428     
429! 6. Load Geothermal Flux
430! ----------------------------------------------------------------------
431
432
433      call get_field("flux_geo",flux_geo,found,timeindex=indextime)
434      if (.not.found) then
435           write(*,*) "soil_settings: Failed loading <flux_geo>"
436           write(*,*) "soil_settings: Will put  <flux_geo> to 0."
437           flux_geo(:,:) = 0.
438      endif
439
440
441! 7. Adsorption
442! ----------------------------------------------------------------------
443
444
445      if (adsorption_soil) then
446! Subsurface water vapor
447             txt="h2o_vap_soil"
448             write(*,*) 'phyetat0: loading subsurface tracer',
449     &                               ' h2o_vap_soil'
450             call get_field(txt,qsoil(:,:,igcm_h2o_vap_soil,:),found,
451     &                                                 indextime)
452             if (.not.found) then
453               write(*,*) "phyetat0: Failed loading <",trim(txt),">"
454               write(*,*) "         ",trim(txt)," is set to zero"
455               qsoil(:,:,igcm_h2o_vap_soil,:)= 0.
456             else
457               write(*,*) "phyetat0: suburface tracer <",trim(txt),
458     &              "> range:", minval(qsoil(:,:,igcm_h2o_vap_soil,:)),
459     &                          maxval(qsoil(:,:,igcm_h2o_vap_soil,:))
460             endif
461! Subsurface ice
462               txt="h2o_ice_soil"
463               write(*,*) 'phyetat0: loading subsurface tracer',
464     &                                ' h2o_ice_soil'
465             call get_field(txt,qsoil(:,:,igcm_h2o_ice_soil,:),found,
466     &                                                  indextime)
467             if (.not.found) then
468               write(*,*) "phyetat0: Failed loading <",trim(txt),">"
469               write(*,*) "         ",trim(txt)," is set to zero"
470               qsoil(:,:,igcm_h2o_ice_soil,:)= 0.
471             else
472               write(*,*) "phyetat0: suburface tracer <",trim(txt),
473     &                     "> range:",
474     &                     minval(qsoil(:,:,igcm_h2o_ice_soil,:)),
475     &                     maxval(qsoil(:,:,igcm_h2o_ice_soil,:))
476             endif
477! Adsorbed water
478             txt="h2o_vap_ads"
479               write(*,*) 'phyetat0: loading subsurface tracer',
480     &                                ' h2o_vap_ads'
481             call get_field(txt,qsoil(:,:,igcm_h2o_vap_ads,:),found,
482     &                                                indextime)
483             if (.not.found) then
484               write(*,*) "phyetat0: Failed loading <",trim(txt),">"
485               write(*,*) "         ",trim(txt)," is set to zero"
486               qsoil(:,:,igcm_h2o_vap_ads,:)= 0.
487             else
488               write(*,*) "phyetat0: suburface tracer <",trim(txt),">
489     &                     range:",
490     &                     minval(qsoil(:,:,igcm_h2o_vap_ads,:)),
491     &                     maxval(qsoil(:,:,igcm_h2o_vap_ads,:))
492             endif
493
494      endif ! of adsorption_soil
495
496
497     
498! 8. Report min and max values of soil temperatures and thermal inertias
499! ----------------------------------------------------------------------
500
501      write(*,*)
502      write(*,*)'Soil volumetric heat capacity:',volcapa
503
504      xmin = MINVAL(inertiedat)
505      xmax = MAXVAL(inertiedat)
506      write(*,*)'Soil thermal inertia <inertiedat>:',xmin,xmax
507
508      xmin = MINVAL(inertiesoil)
509      xmax = MAXVAL(inertiesoil)
510      write(*,*)'Soil thermal inertia <inertiesoil>:',xmin,xmax
511
512      xmin = 1.0E+20
513      xmax = -1.0E+20
514      xmin = MINVAL(tsoil)
515      xmax = MAXVAL(tsoil)
516      write(*,*)'Soil temperature <tsoil>:',xmin,xmax
517
518      end
Note: See TracBrowser for help on using the repository browser.