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

Last change on this file since 3026 was 2951, checked in by emillour, 20 months ago

Mars PCM:
Fix OpenMP bug on inertiesoil, unnecessary loop and a line
too long for fixed fortran format.
Minor fix for picky gfortran compiler in testphys1d (no .eq. for logicals!)
EM

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