source: LMDZ6/branches/Amaury_dev/libf/phylmd/ecrad.v1.5.1/srtm_taumol17.F90 @ 5157

Last change on this file since 5157 was 3908, checked in by idelkadi, 3 years ago

Online implementation of the radiative transfer code ECRAD in the LMDZ model.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 8.5 KB
Line 
1SUBROUTINE SRTM_TAUMOL17 &
2 & ( KIDIA   , KFDIA    , KLEV,&
3 & P_FAC00   , P_FAC01  , P_FAC10   , P_FAC11,&
4 & K_JP      , K_JT     , K_JT1     , P_ONEMINUS,&
5 & P_COLH2O  , P_COLCO2 , P_COLMOL,&
6 & K_LAYTROP , P_SELFFAC, P_SELFFRAC, K_INDSELF  , P_FORFAC, P_FORFRAC, K_INDFOR,&
7 & P_SFLUXZEN, P_TAUG   , P_TAUR    , PRMU0   &
8 & ) 
9
10!     Written by Eli J. Mlawer, Atmospheric & Environmental Research.
11
12!     BAND 17:  3250-4000 cm-1 (low - H2O,CO2; high - H2O,CO2)
13
14! Modifications
15!        M.Hamrud      01-Oct-2003 CY28 Cleaning
16
17!     JJMorcrette 2003-02-24 adapted to ECMWF environment
18!        D.Salmond  31-Oct-2007 Vector version in the style of RRTM from Meteo France & NEC
19!     JJMorcrette 20010610 Flexible configuration for number of g-points
20
21USE PARKIND1 , ONLY : JPIM, JPRB
22USE YOMHOOK  , ONLY : LHOOK, DR_HOOK
23USE PARSRTM  , ONLY : JPG
24USE YOESRTM  , ONLY : NG17
25USE YOESRTA17, ONLY : ABSA, ABSB, FORREFC, SELFREFC, SFLUXREFC, RAYL, LAYREFFR, STRRAT 
26USE YOESRTWN , ONLY : NSPA, NSPB
27
28IMPLICIT NONE
29
30!-- Output
31INTEGER(KIND=JPIM),INTENT(IN)    :: KIDIA, KFDIA
32INTEGER(KIND=JPIM),INTENT(IN)    :: KLEV
33REAL(KIND=JPRB)   ,INTENT(IN)    :: P_FAC00(KIDIA:KFDIA,KLEV)
34REAL(KIND=JPRB)   ,INTENT(IN)    :: P_FAC01(KIDIA:KFDIA,KLEV)
35REAL(KIND=JPRB)   ,INTENT(IN)    :: P_FAC10(KIDIA:KFDIA,KLEV)
36REAL(KIND=JPRB)   ,INTENT(IN)    :: P_FAC11(KIDIA:KFDIA,KLEV)
37INTEGER(KIND=JPIM),INTENT(IN)    :: K_JP(KIDIA:KFDIA,KLEV)
38INTEGER(KIND=JPIM),INTENT(IN)    :: K_JT(KIDIA:KFDIA,KLEV)
39INTEGER(KIND=JPIM),INTENT(IN)    :: K_JT1(KIDIA:KFDIA,KLEV)
40REAL(KIND=JPRB)   ,INTENT(IN)    :: P_ONEMINUS(KIDIA:KFDIA)
41REAL(KIND=JPRB)   ,INTENT(IN)    :: P_COLH2O(KIDIA:KFDIA,KLEV)
42REAL(KIND=JPRB)   ,INTENT(IN)    :: P_COLCO2(KIDIA:KFDIA,KLEV)
43REAL(KIND=JPRB)   ,INTENT(IN)    :: P_COLMOL(KIDIA:KFDIA,KLEV)
44INTEGER(KIND=JPIM),INTENT(IN)    :: K_LAYTROP(KIDIA:KFDIA)
45REAL(KIND=JPRB)   ,INTENT(IN)    :: P_SELFFAC(KIDIA:KFDIA,KLEV)
46REAL(KIND=JPRB)   ,INTENT(IN)    :: P_SELFFRAC(KIDIA:KFDIA,KLEV)
47INTEGER(KIND=JPIM),INTENT(IN)    :: K_INDSELF(KIDIA:KFDIA,KLEV)
48REAL(KIND=JPRB)   ,INTENT(IN)    :: P_FORFAC(KIDIA:KFDIA,KLEV)
49REAL(KIND=JPRB)   ,INTENT(IN)    :: P_FORFRAC(KIDIA:KFDIA,KLEV)
50INTEGER(KIND=JPIM),INTENT(IN)    :: K_INDFOR(KIDIA:KFDIA,KLEV)
51
52REAL(KIND=JPRB)   ,INTENT(OUT)   :: P_SFLUXZEN(KIDIA:KFDIA,JPG)
53REAL(KIND=JPRB)   ,INTENT(OUT)   :: P_TAUG(KIDIA:KFDIA,KLEV,JPG)
54REAL(KIND=JPRB)   ,INTENT(OUT)   :: P_TAUR(KIDIA:KFDIA,KLEV,JPG)
55REAL(KIND=JPRB)   ,INTENT(IN)    :: PRMU0(KIDIA:KFDIA)
56!- from INTFAC     
57!- from INTIND
58!- from PRECISE             
59!- from PROFDATA             
60!- from SELF             
61INTEGER(KIND=JPIM) :: IG, IND0, IND1, INDS, INDF, JS, I_LAY, I_LAYSOLFR(KIDIA:KFDIA), I_NLAYERS, IPLON
62
63! REAL(KIND=JPRB) :: Z_FAC000, Z_FAC001, Z_FAC010, Z_FAC011, Z_FAC100, Z_FAC101,&
64!  & Z_FAC110, Z_FAC111
65REAL(KIND=JPRB) :: Z_FS, Z_SPECCOMB, Z_SPECMULT, Z_SPECPARM, Z_TAURAY 
66REAL(KIND=JPRB) :: ZHOOK_HANDLE
67
68ASSOCIATE(NFLEVG=>KLEV)
69IF (LHOOK) CALL DR_HOOK('SRTM_TAUMOL17',0,ZHOOK_HANDLE)
70
71I_NLAYERS = KLEV
72
73!     Compute the optical depth by interpolating in ln(pressure),
74!     temperature, and appropriate species.  Below LAYTROP, the water
75!     vapor self-continuum is interpolated (in temperature) separately. 
76
77DO I_LAY = 1, I_NLAYERS
78  DO IPLON = KIDIA, KFDIA
79    IF (PRMU0(IPLON) > 0.0_JPRB) THEN
80      IF (I_LAY <= K_LAYTROP(IPLON)) THEN
81        Z_SPECCOMB = P_COLH2O(IPLON,I_LAY) + STRRAT*P_COLCO2(IPLON,I_LAY)
82        Z_SPECPARM = P_COLH2O(IPLON,I_LAY)/Z_SPECCOMB
83        IF (Z_SPECPARM >= P_ONEMINUS(IPLON)) Z_SPECPARM = P_ONEMINUS(IPLON)
84        Z_SPECMULT = 8.*(Z_SPECPARM)
85        JS = 1 + INT(Z_SPECMULT)
86        Z_FS = MOD(Z_SPECMULT, 1.0_JPRB )
87        ! Z_FAC000 = (1. - Z_FS) * P_FAC00(I_LAY)
88        ! Z_FAC010 = (1. - Z_FS) * P_FAC10(I_LAY)
89        ! Z_FAC100 = Z_FS * P_FAC00(I_LAY)
90        ! Z_FAC110 = Z_FS * P_FAC10(I_LAY)
91        ! Z_FAC001 = (1. - Z_FS) * P_FAC01(I_LAY)
92        ! Z_FAC011 = (1. - Z_FS) * P_FAC11(I_LAY)
93        ! Z_FAC101 = Z_FS * P_FAC01(I_LAY)
94        ! Z_FAC111 = Z_FS * P_FAC11(I_LAY)
95        IND0 = ((K_JP(IPLON,I_LAY)-1)*5+(K_JT(IPLON,I_LAY)-1))*NSPA(17) + JS
96        IND1 = (K_JP(IPLON,I_LAY)*5+(K_JT1(IPLON,I_LAY)-1))*NSPA(17) + JS
97        INDS = K_INDSELF(IPLON,I_LAY)
98        INDF = K_INDFOR(IPLON,I_LAY)
99        Z_TAURAY = P_COLMOL(IPLON,I_LAY) * RAYL
100
101        !  DO IG = 1, NG(17)
102!CDIR UNROLL=NG17
103        DO IG = 1, NG17
104          P_TAUG(IPLON,I_LAY,IG) = Z_SPECCOMB * &
105           & (&
106           & (1. - Z_FS) * ( ABSA(IND0,IG) * P_FAC00(IPLON,I_LAY) + &
107           &                 ABSA(IND0+9,IG) * P_FAC10(IPLON,I_LAY) + &
108           &                 ABSA(IND1,IG) * P_FAC01(IPLON,I_LAY) + &
109           &                 ABSA(IND1+9,IG) * P_FAC11(IPLON,I_LAY) ) + &
110           & Z_FS        * ( ABSA(IND0+1,IG) * P_FAC00(IPLON,I_LAY) + &
111           &                 ABSA(IND0+10,IG) * P_FAC10(IPLON,I_LAY) + &
112           &                 ABSA(IND1+1,IG) * P_FAC01(IPLON,I_LAY) + &
113           &                 ABSA(IND1+10,IG) * P_FAC11(IPLON,I_LAY) ) &
114           & )  + &
115           & P_COLH2O(IPLON,I_LAY) * &
116           & (P_SELFFAC(IPLON,I_LAY) * (SELFREFC(INDS,IG) + &
117           & P_SELFFRAC(IPLON,I_LAY) * &
118           & (SELFREFC(INDS+1,IG) - SELFREFC(INDS,IG))) + &
119           & P_FORFAC(IPLON,I_LAY) * (FORREFC(INDF,IG) + &
120           & P_FORFRAC(IPLON,I_LAY) * &
121           & (FORREFC(INDF+1,IG) - FORREFC(INDF,IG))))   
122
123          !     &           + TAURAY
124          !    SSA(LAY,IG) = TAURAY/TAUG(LAY,IG)
125          P_TAUR(IPLON,I_LAY,IG) = Z_TAURAY
126        ENDDO
127      ENDIF
128    ENDIF
129  ENDDO
130ENDDO
131
132I_LAYSOLFR(:) = I_NLAYERS
133
134DO I_LAY = 1, I_NLAYERS
135  DO IPLON = KIDIA, KFDIA
136    IF (PRMU0(IPLON) > 0.0_JPRB) THEN
137      IF (I_LAY >= K_LAYTROP(IPLON)+1) THEN
138        IF (K_JP(IPLON,I_LAY-1) < LAYREFFR .AND. K_JP(IPLON,I_LAY) >= LAYREFFR) &
139         & I_LAYSOLFR(IPLON) = I_LAY 
140        Z_SPECCOMB = P_COLH2O(IPLON,I_LAY) + STRRAT*P_COLCO2(IPLON,I_LAY)
141        Z_SPECPARM = P_COLH2O(IPLON,I_LAY)/Z_SPECCOMB
142        IF (Z_SPECPARM >= P_ONEMINUS(IPLON)) Z_SPECPARM = P_ONEMINUS(IPLON)
143        Z_SPECMULT = 4.*(Z_SPECPARM)
144        JS = 1 + INT(Z_SPECMULT)
145        Z_FS = MOD(Z_SPECMULT, 1.0_JPRB )
146        ! Z_FAC000 = (1. - Z_FS) * P_FAC00(I_LAY)
147        ! Z_FAC010 = (1. - Z_FS) * P_FAC10(I_LAY)
148        ! Z_FAC100 = Z_FS * P_FAC00(I_LAY)
149        ! Z_FAC110 = Z_FS * P_FAC10(I_LAY)
150        ! Z_FAC001 = (1. - Z_FS) * P_FAC01(I_LAY)
151        ! Z_FAC011 = (1. - Z_FS) * P_FAC11(I_LAY)
152        ! Z_FAC101 = Z_FS * P_FAC01(I_LAY)
153        ! Z_FAC111 = Z_FS * P_FAC11(I_LAY)
154        IND0 = ((K_JP(IPLON,I_LAY)-13)*5+(K_JT(IPLON,I_LAY)-1))*NSPB(17) + JS
155        IND1 = ((K_JP(IPLON,I_LAY)-12)*5+(K_JT1(IPLON,I_LAY)-1))*NSPB(17) + JS
156        INDF = K_INDFOR(IPLON,I_LAY)
157        Z_TAURAY = P_COLMOL(IPLON,I_LAY) * RAYL
158
159        !  DO IG = 1, NG(17)
160!CDIR UNROLL=NG17
161        DO IG = 1, NG17
162          P_TAUG(IPLON,I_LAY,IG) = Z_SPECCOMB * &
163           !    & (Z_FAC000 * ABSB(IND0,IG) + &
164           !    & Z_FAC100 * ABSB(IND0+1,IG) + &
165           !    & Z_FAC010 * ABSB(IND0+5,IG) + &
166           !    & Z_FAC110 * ABSB(IND0+6,IG) + &
167           !    & Z_FAC001 * ABSB(IND1,IG) + &
168           !    & Z_FAC101 * ABSB(IND1+1,IG) + &
169           !    & Z_FAC011 * ABSB(IND1+5,IG) + &
170           !    & Z_FAC111 * ABSB(IND1+6,IG)) + &
171           & (&
172           & (1. - Z_FS) * ( ABSB(IND0,IG) * P_FAC00(IPLON,I_LAY) + &
173           &                 ABSB(IND0+5,IG) * P_FAC10(IPLON,I_LAY) + &
174           &                 ABSB(IND1,IG) * P_FAC01(IPLON,I_LAY) + &
175           &                 ABSB(IND1+5,IG) * P_FAC11(IPLON,I_LAY) ) + &
176           & Z_FS        * ( ABSB(IND0+1,IG) * P_FAC00(IPLON,I_LAY) + &
177           &                 ABSB(IND0+6,IG) * P_FAC10(IPLON,I_LAY) + &
178           &                 ABSB(IND1+1,IG) * P_FAC01(IPLON,I_LAY) + &
179           &                 ABSB(IND1+6,IG) * P_FAC11(IPLON,I_LAY) ) &
180           & ) + &
181           & P_COLH2O(IPLON,I_LAY) * &
182           & P_FORFAC(IPLON,I_LAY) * (FORREFC(INDF,IG) + &
183           & P_FORFRAC(IPLON,I_LAY) * &
184           & (FORREFC(INDF+1,IG) - FORREFC(INDF,IG)))   
185          !     &           + TAURAY
186          !    SSA(LAY,IG) = TAURAY/TAUG(LAY,IG)
187          IF (I_LAY == I_LAYSOLFR(IPLON)) P_SFLUXZEN(IPLON,IG) = SFLUXREFC(IG,JS) &
188           & + Z_FS * (SFLUXREFC(IG,JS+1) - SFLUXREFC(IG,JS)) 
189          P_TAUR(IPLON,I_LAY,IG) = Z_TAURAY
190        ENDDO
191      ENDIF
192    ENDIF
193  ENDDO
194ENDDO
195
196!-----------------------------------------------------------------------
197IF (LHOOK) CALL DR_HOOK('SRTM_TAUMOL17',1,ZHOOK_HANDLE)
198END ASSOCIATE
199END SUBROUTINE SRTM_TAUMOL17
Note: See TracBrowser for help on using the repository browser.