source: trunk/LMDZ.PLUTO/libf/phypluto/evolch4.F90 @ 3945

Last change on this file since 3945 was 3932, checked in by tbertrand, 8 weeks ago

PLUTO PCM: Allow 1D model to be run with empirical seasonal cycles of N2 and CH4 (change in pressure and CH4 abundance with time)
TB

File size: 2.4 KB
Line 
1      subroutine evolch4(ngrid,nlayer,pls,pplev,pdpsrf,zqch4evol)
2      use datafile_mod
3      use comcstfi_mod, only: pi
4      use mod_phys_lmdz_para, only : is_master
5      use tracer_h, only: igcm_ch4_gas, igcm_n2, mmol
6
7      implicit none
8
9!==================================================================
10!     Purpose
11!     -------
12!     Get tracer fields according to the solar longitude (for 1D purpose for now)
13!
14!     Inputs
15!     ------
16!     pls                 solar longitude
17!     pplev               pressure
18!     pdpsrf              tendency on pressure
19!
20!     Outputs
21!     -------
22!     zdqch4evol          tendency on zqch4
23
24!     Both
25!     ----
26!
27!     Authors
28!     -------
29!     Tanguy Bertrand
30!==================================================================
31
32!-----------------------------------------------------------------------
33!     Arguments
34
35      INTEGER ngrid, nlayer
36      REAL,INTENT(IN) :: pls
37      REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! pressure levels/
38      REAL,INTENT(IN) :: pdpsrf(ngrid)
39      REAL,INTENT(OUT) :: zqch4evol(nlayer) ! Final tendancy
40!-----------------------------------------------------------------------
41!     Local variables
42      LOGICAL,SAVE :: firstcall=.true.
43!$OMP THREADPRIVATE(firstcall)
44
45      integer Nfine
46      parameter(Nfine=143)
47      integer ifine
48      character(len=100) :: file_path
49      real,save,allocatable :: lssdat(:)
50      real,save,allocatable :: vmrdat(:)
51      real :: vmrsrf
52
53
54!---------------- INPUT ------------------------------------------------
55
56      IF (firstcall) then
57         firstcall=.false.
58
59!$OMP MASTER
60         file_path=trim(datadir)//'/cycle_vmrch4.txt'
61         if (is_master) print*,file_path
62
63         open(222,file=file_path,form='formatted')
64
65            if(.not.allocated(lssdat)) then
66               allocate(lssdat(Nfine))
67            endif
68            if(.not.allocated(vmrdat)) then
69               allocate(vmrdat(Nfine))
70            endif
71
72            do ifine=1,Nfine
73               read(222,*) lssdat(ifine), vmrdat(ifine)
74            enddo
75            close(222)
76
77!$OMP END MASTER
78!$OMP BARRIER
79      ENDIF
80
81      CALL interp_line(lssdat,vmrdat,Nfine,pls*180./pi,vmrsrf,1)
82
83      ! Reconstruct CH4 profile in percent
84      zqch4evol(:)=min(vmrsrf*(pplev(1,1:nlayer)/pplev(1,1))**((16.-28.)/28.),50.)
85      ! Convert in MMR
86      zqch4evol(:)=zqch4evol(:)/100.*mmol(igcm_n2)/mmol(igcm_ch4_gas)
87
88      end subroutine evolch4
89
Note: See TracBrowser for help on using the repository browser.