source: trunk/LMDZ.MARS/libf/phymars/nirco2abs.F @ 2398

Last change on this file since 2398 was 2398, checked in by emillour, 4 years ago

Mars GCM:
Some code cleanup: use "call abort_physic()" instead of "stop" or "call abort"
EM

File size: 9.0 KB
Line 
1      SUBROUTINE nirco2abs(ngrid,nlayer,pplay,dist_sol,nq,pq,
2     $     mu0,fract,declin,pdtnirco2)
3                                                   
4       use tracer_mod, only: igcm_co2, igcm_o
5       use comgeomfi_h, only: sinlon, coslon, sinlat, coslat
6       USE comcstfi_h, ONLY: pi
7       USE time_phylmdz_mod, ONLY: daysec
8       IMPLICIT NONE
9c=======================================================================
10c   subject:
11c   --------
12c   Computing heating rate due to
13c   absorption by CO2 in the near-infrared
14c   This version includes NLTE effects
15c
16c   (Scheme to be described in Forget et al., JGR, 2003)
17c   (old Scheme described in Forget et al., JGR, 1999)
18c
19c   This version updated with a new functional fit,
20c   see NLTE correction-factor of Lopez-Valverde et al (1998)
21c   Stephen Lewis 2000
22c
23c   jul 2011 malv+fgg    New corrections for NLTE implemented
24c   08/2002 : correction for bug when running with diurnal=F
25c
26c   author:  Frederic Hourdin 1996
27c   ------
28c            Francois Forget 1999
29c
30c   input:
31c   -----
32c   ngrid                 number of gridpoint of horizontal grid
33c   nlayer                Number of layer
34c   dist_sol              sun-Mars distance (AU)
35c   mu0(ngrid)         
36c   fract(ngrid)          day fraction of the time interval
37c   declin                latitude of subslar point
38c
39c   output:
40c   -------
41c
42c   pdtnirco2(ngrid,nlayer)      Heating rate (K/s)
43c
44c
45c=======================================================================
46c
47c    0.  Declarations :
48c    ------------------
49c
50      include "callkeys.h"
51      include "nirdata.h"
52
53c-----------------------------------------------------------------------
54c    Input/Output
55c    ------------
56      integer,intent(in) :: ngrid ! number of (horizontal) grid points
57      integer,intent(in) :: nlayer ! number of atmospheric layers
58      real,intent(in) :: pplay(ngrid,nlayer) ! Pressure
59      real,intent(in) :: dist_sol ! Sun-Mars distance (in AU)
60      integer,intent(in) :: nq ! number of tracers
61      real,intent(in) :: pq(ngrid,nlayer,nq) ! tracers
62      real,intent(in) :: mu0(ngrid) ! solar angle
63      real,intent(in) :: fract(ngrid) ! day fraction of the time interval
64      real,intent(in) :: declin ! latitude of sub-solar point
65     
66      real,intent(out) :: pdtnirco2(ngrid,nlayer) ! heating rate (K/s)
67c
68c    Local variables :
69c    -----------------
70      INTEGER l,ig, n, nstep,i
71      REAL co2heat0, zmu(ngrid)
72
73c     special diurnal=F
74      real mu0_int(ngrid),fract_int(ngrid),zday_int
75      real ztim1,ztim2,ztim3,step
76
77c
78c   local saved variables
79c   ---------------------
80      logical,save :: firstcall=.true.
81      integer,save :: ico2=0 ! index of "co2" tracer
82      integer,save :: io=0 ! index of "o" tracer
83c     p0noonlte is a pressure below which non LTE effects are significant.
84c     REAL p0nonlte
85c     DATA p0nonlte/7.5e-3/
86c     SAVE p0nonlte
87
88c     parameters for CO2 heating fit
89      real,parameter :: n_a=1.1956475
90      real,parameter :: n_p0=0.0015888279
91      real,parameter :: n_b=1.9628251
92
93c     Variables added to implement NLTE correction factor (feb 2011)
94      real    pyy(nlayer)
95      real    cor1(nlayer),oldoco2(nlayer),alfa2(nlayer)
96      real    p2011,cociente1,merge
97      real    cor0,oco2gcm
98
99c----------------------------------------------------------------------
100
101c     Initialisation
102c     --------------
103      ! AS: OK firstcall absolute
104      if (firstcall) then
105        if (nircorr.eq.1) then
106          ! we will need co2 and o tracers
107          ico2=igcm_co2
108          if (ico2==0) then
109            write(*,*) "nirco2abs error: I need a CO2 tracer"
110            write(*,*) "     when running with nircorr==1"
111            call abort_physic("nirco2abs","need a CO2 tracer",1)
112          endif
113          io=igcm_o
114          if (io==0) then
115            write(*,*) "nirco2abs error: I need an O tracer"
116            write(*,*) "     when running with nircorr==1"
117            call abort_physic("nirco2abs","need an O tracer",1)
118          endif
119        endif
120        firstcall=.false.
121      endif
122
123
124c     co2heat is the heating by CO2 at 700Pa for a zero zenithal angle.
125      co2heat0=n_a*(1.52/dist_sol)**2/daysec
126
127c     Simple calcul for a given sun incident angle (if diurnal=T)
128c     --------------------------------------------
129
130      IF (diurnal) THEN 
131         do ig=1,ngrid
132            zmu(ig)=sqrt(1224.*mu0(ig)*mu0(ig)+1.)/35.
133
134            if(nircorr.eq.1) then
135               do l=1,nlayer
136                  pyy(l)=pplay(ig,l)
137               enddo
138
139               call interpnir(cor1,pyy,nlayer,corgcm,pres1d,npres)
140               call interpnir(oldoco2,pyy,nlayer,oco21d,pres1d,npres)
141               call interpnir(alfa2,pyy,nlayer,alfa,pres1d,npres)
142            endif
143
144            do l=1,nlayer
145!           Calculations for the O/CO2 correction
146               if(nircorr.eq.1) then
147                  cor0=1./(1.+n_p0/pplay(ig,l))**n_b
148                  if(pq(ig,l,ico2).gt.1.e-6) then                   
149                     oco2gcm=pq(ig,l,io)/pq(ig,l,ico2)
150                     ! handle the rare cases when pq(ig,l,io)<0
151                     if (pq(ig,l,io).lt.0) then
152                       write(*,*) "nirco2abs: warning ig=",ig," l=",l,
153     &                            " pq(ig,l,io)=",pq(ig,l,io)
154                       oco2gcm=1.e6
155                     endif
156                  else
157                     oco2gcm=1.e6
158                  endif
159                  cociente1=oco2gcm/oldoco2(l)
160                  merge=alog10(cociente1)*alfa2(l)+alog10(cor0)*
161     $                 (1.-alfa2(l))
162                  merge=10**merge
163                  p2011=sqrt(merge)*cor0
164               else if (nircorr.eq.0) then
165                  p2011=1.
166                  cor1(l)=1.
167               endif
168
169               if(fract(ig).gt.0.) pdtnirco2(ig,l)=
170     &             co2heat0*sqrt((700.*zmu(ig))/pplay(ig,l))
171     &             /(1.+n_p0/pplay(ig,l))**n_b
172!           Corrections from tabulation
173     $              * cor1(l) * p2011
174c           OLD SCHEME (forget et al. 1999)
175c    s           co2heat0*sqrt((700.*zmu(ig))/pplay(ig,l))
176c    s          / (1+p0nonlte/pplay(ig,l))
177           enddo
178         enddo
179         
180
181c     Averaging over diurnal cycle (if diurnal=F)
182c     -------------------------------------------
183c     NIR CO2 abs is slightly non linear. To remove the diurnal
184c     cycle, it is better to average the heating rate over 1 day rather
185c     than using the mean mu0 computed by mucorr in physiq.F (FF, 1998)
186
187      ELSE      ! if (.not.diurnal) then
188
189         nstep = 20   ! number of integration step /sol
190         mu0_int(1:ngrid) = 0.
191         ztim1 = 0.
192         do n=1,nstep
193            zday_int = (n-1)/float(nstep)
194            ztim2=COS(declin)*COS(2.*pi*(zday_int-.5))
195            ztim3=-COS(declin)*SIN(2.*pi*(zday_int-.5))
196            CALL solang(ngrid,sinlon,coslon,sinlat,coslat,
197     s             ztim1,ztim2,ztim3,
198     s             mu0_int,fract_int)
199            do ig=1,ngrid               
200               zmu(ig)=sqrt(1224.*mu0_int(ig)*mu0_int(ig)+1.)/35.
201
202               if(nircorr.eq.1) then
203                  do l=1,nlayer
204                     pyy(l)=pplay(ig,l)
205                  enddo
206               call interpnir(cor1,pyy,nlayer,corgcm,pres1d,npres)
207               call interpnir(oldoco2,pyy,nlayer,oco21d,pres1d,npres)
208               call interpnir(alfa2,pyy,nlayer,alfa,pres1d,npres)
209               endif
210
211               do l=1,nlayer
212                  if(nircorr.eq.1) then
213                     cor0=1./(1.+n_p0/pplay(ig,l))**n_b
214                     oco2gcm=pq(ig,l,io)/pq(ig,l,ico2)
215                     cociente1=oco2gcm/oldoco2(l)
216                     merge=alog10(cociente1)*alfa2(l)+alog10(cor0)*
217     $                    (1.-alfa2(l))
218                     merge=10**merge
219                     p2011=sqrt(merge)*cor0
220                  else if (nircorr.eq.0) then
221                     p2011=1.
222                     cor1(l)=1.
223                  endif
224
225                  if(fract_int(ig).gt.0.) pdtnirco2(ig,l)=
226     &                 pdtnirco2(ig,l) + (1/float(nstep))*
227     &                 co2heat0*sqrt((700.*zmu(ig))/pplay(ig,l))
228     &                 /(1.+n_p0/pplay(ig,l))**n_b
229!                      Corrections from tabulation
230     $                 * cor1(l) * p2011
231               enddo
232            enddo
233         end do
234      END IF
235
236      end
237
238
239     
240      subroutine interpnir(escout,p,nlayer,escin,pin,nl)
241C
242C subroutine to perform linear interpolation in pressure from 1D profile
243C escin(nl) sampled on pressure grid pin(nl) to profile
244C escout(nlayer) on pressure grid p(nlayer).
245C
246      real escout(nlayer),p(nlayer)
247      real escin(nl),pin(nl),wm,wp
248      integer nl,nlayer,n1,n,nm,np
249      do n1=1,nlayer
250         if(p(n1) .gt. 1500. .or. p(n1) .lt. 1.0e-13) then
251            escout(n1) = 0.0
252         else
253            do n = 1,nl-1
254               if (p(n1).le.pin(n).and.p(n1).ge.pin(n+1)) then
255                  nm=n
256                  np=n+1
257                  wm=abs(pin(np)-p(n1))/(pin(nm)-pin(np))
258                  wp=1.0 - wm
259               endif
260            enddo
261            escout(n1) = escin(nm)*wm + escin(np)*wp
262         endif
263      enddo
264      end
Note: See TracBrowser for help on using the repository browser.