source: trunk/LMDZ.PLUTO/libf/phypluto/mugas_mod.F90 @ 3959

Last change on this file since 3959 was 3949, checked in by debatzbr, 7 weeks ago

Pluto PCM: Add condensable gas tracers through muphi
BBT

File size: 5.8 KB
Line 
1subroutine mugas_prof(ngrid,nlayer,nq,zzlay,zzlev,pplay,pt,pdqchem)
2      use datafile_mod
3      use comcstfi_mod, only: r
4      use tracer_h, only: noms,igcm_C2H2_mugas,igcm_C2H6_mugas,&
5                          igcm_C4H2_mugas,igcm_C6H6_mugas,igcm_HCN_mugas
6      use mod_phys_lmdz_para, only : is_master
7
8      implicit none
9
10!==================================================================
11!     Purpose
12!     -------
13!     Get fixed fluxes of microphysical condensable gas
14!     from txt file.
15!
16!     Inputs
17!     ------
18!     ngrid   | Number of atmospheric columns.
19!     nlayer  | Number of atmospheric layers.
20!     nq      | Number of tracers.
21!     zzlay   | Altitude at the middle of the atmospheric layers.
22!     zzlev   | Altitude at the atmospheric layer boundaries.
23!     pplay   | Mid-layer pressure (Pa)
24!     pt      | Temperature (K).
25!
26!     Outputs
27!     -------
28!
29!     Both
30!     ----
31!     pdqchem ! Tracer tendencies for condensable gases through muphi (kg/kg_of_air/s).
32!
33!     Authors
34!     -------
35!     Bruno de Batz de Trenquelléon (2025)
36!==================================================================
37     
38      ! Arguments:
39      !~~~~~~~~~~~
40      integer,intent(in) :: ngrid
41      integer,intent(in) :: nlayer
42      integer,intent(in) :: nq
43      real,intent(in)    :: zzlay(ngrid,nlayer)
44      real,intent(in)    :: zzlev(ngrid,nlayer+1)
45      real,intent(in)    :: pplay(ngrid,nlayer)
46      real,intent(in)    :: pt(ngrid,nlayer)
47
48      real,intent(inout) :: pdqchem(ngrid,nlayer,nq)
49
50      ! Local variables:
51      !~~~~~~~~~~~~~~~~~
52      integer :: ig, iq, l, ifine
53
54      real :: dz(ngrid,nlayer)
55      real :: zdq_HCN(ngrid,nlayer)
56      real :: zdq_C4H2(ngrid,nlayer)
57      real :: zdq_C6H6(ngrid,nlayer)
58      real :: zdq_C2H2(ngrid,nlayer)
59      real :: zdq_C2H6(ngrid,nlayer)
60     
61      integer Nfine
62      parameter(Nfine=999)
63      character(len=100) :: file_path
64      character(len=100) :: file_name
65      real,save,allocatable :: Altdata(:)
66      real,save,allocatable :: Pdata(:)
67      real,save,allocatable :: Tdata(:)
68      real,save,allocatable :: C2H2_fluxdata(:)
69      real,save,allocatable :: C2H6_fluxdata(:)
70      real,save,allocatable :: C4H2_fluxdata(:)
71      real,save,allocatable :: C6H6_fluxdata(:)
72      real,save,allocatable :: HCN_fluxdata(:)
73
74      LOGICAL,SAVE :: firstcall=.true.
75!$OMP THREADPRIVATE(firstcall)
76
77      ! 0. Initialization:
78      !~~~~~~~~~~~~~~~~~~~
79      dz(:,:) = zzlev(:,2:nlayer+1)-zzlev(:,1:nlayer)
80
81      zdq_C2H2(:,:) = 0.d0
82      zdq_C2H6(:,:) = 0.d0
83      zdq_C4H2(:,:) = 0.d0
84      zdq_C6H6(:,:) = 0.d0
85      zdq_HCN(:,:)  = 0.d0
86
87      ! 1. Read data:
88      !~~~~~~~~~~~~~~
89      if (firstcall) then
90        firstcall=.false.
91
92        if (mugasflux_file/='None') then
93            file_name = mugasflux_file
94        else
95            STOP "No filename given for condensable gas fluxes profile. Need to set mugasflux_file..."
96        endif
97!$OMP MASTER
98
99        if(.not.allocated(Altdata))       allocate(Altdata(Nfine))
100        if(.not.allocated(Pdata))         allocate(Pdata(Nfine))
101        if(.not.allocated(Tdata))         allocate(Tdata(Nfine))
102        if(.not.allocated(C2H2_fluxdata)) allocate(C2H2_fluxdata(Nfine))
103        if(.not.allocated(C2H6_fluxdata)) allocate(C2H6_fluxdata(Nfine))
104        if(.not.allocated(C4H2_fluxdata)) allocate(C4H2_fluxdata(Nfine))
105        if(.not.allocated(C6H6_fluxdata)) allocate(C6H6_fluxdata(Nfine))
106        if(.not.allocated(HCN_fluxdata))  allocate(HCN_fluxdata(Nfine))
107       
108        file_path=file_name
109        open(224,file=file_path,form='formatted')
110        do ifine = 1, Nfine
111            read(224,*) Altdata(ifine), Pdata(ifine), Tdata(ifine), HCN_fluxdata(ifine), &
112                        C4H2_fluxdata(ifine), C6H6_fluxdata(ifine), C2H2_fluxdata(ifine), C2H6_fluxdata(ifine)
113        enddo
114        close(224)
115!$OMP END MASTER
116!$OMP BARRIER
117      endif
118
119      ! 2. Interpolate on the model vertical grid:
120      !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
121      do iq = 1, nq
122       
123        ! C2H2
124        if (noms(iq).eq."C2H2_mugas") then
125          do ig = 1, ngrid
126            call interp_line(Altdata*1000.,C2H2_fluxdata,Nfine,zzlay(ig,:),zdq_C2H2(ig,:),nlayer)
127            do l = 1, nlayer
128              pdqchem(ig,l,igcm_C2H2_mugas) = zdq_C2H2(ig,l) / dz(ig,l) / (pplay(ig,l)/(r*pt(ig,l)))
129            enddo
130          enddo
131        endif
132
133        ! C2H6
134        if (noms(iq).eq."C2H6_mugas") then
135          do ig = 1, ngrid
136            call interp_line(Altdata*1000.,C2H6_fluxdata,Nfine,zzlay(ig,:),zdq_C2H6(ig,:),nlayer)
137            do l = 1, nlayer
138              pdqchem(ig,l,igcm_C2H6_mugas) = zdq_C2H6(ig,l) / dz(ig,l) / (pplay(ig,l)/(r*pt(ig,l)))
139            enddo
140          enddo
141        endif
142
143        ! C4H2
144        if (noms(iq).eq."C4H2_mugas") then
145          do ig = 1, ngrid
146            call interp_line(Altdata*1000.,C4H2_fluxdata,Nfine,zzlay(ig,:),zdq_C4H2(ig,:),nlayer)
147            do l = 1, nlayer
148              pdqchem(ig,l,igcm_C4H2_mugas) = zdq_C4H2(ig,l) / dz(ig,l) / (pplay(ig,l)/(r*pt(ig,l)))
149            enddo
150          enddo
151        endif
152       
153        ! C6H6
154        if (noms(iq).eq."C6H6_mugas") then
155          do ig = 1, ngrid
156            call interp_line(Altdata*1000.,C6H6_fluxdata,Nfine,zzlay(ig,:),zdq_C6H6(ig,:),nlayer)
157            do l = 1, nlayer
158              pdqchem(ig,l,igcm_C6H6_mugas) = zdq_C6H6(ig,l) / dz(ig,l) / (pplay(ig,l)/(r*pt(ig,l)))
159            enddo
160          enddo
161        endif
162
163        ! HCN
164        if (noms(iq).eq."HCN_mugas") then
165          do ig = 1, ngrid
166            call interp_line(Altdata*1000.,HCN_fluxdata,Nfine,zzlay(ig,:),zdq_HCN(ig,:),nlayer)
167            do l = 1, nlayer
168              pdqchem(ig,l,igcm_HCN_mugas) = zdq_HCN(ig,l) / dz(ig,l) / (pplay(ig,l)/(r*pt(ig,l)))
169            enddo
170          enddo
171        endif
172     
173      enddo ! end of nq
174
175end subroutine mugas_prof
Note: See TracBrowser for help on using the repository browser.