source: trunk/LMDZ.TITAN/libf/phytitan/aeropacity.F90 @ 1647

Last change on this file since 1647 was 1647, checked in by jvatant, 8 years ago

+ Major clean of the new LMDZ.TITAN from too-generic options and routines (water, co2, ocean, surface type ...)
+ From this revision LMDZ.TITAN begins to be really separated from LMDZ.GENERIC
+ Partial desactivation of aerosols, only the dummy case is still enabled to keep the code running ( new aerosol routines to come in followings commits )

JVO

File size: 6.3 KB
Line 
1      Subroutine aeropacity(ngrid,nlayer,nq,pplay,pplev,pq, &
2         aerosol,reffrad,QREFvis3d,QREFir3d,tau_col)
3
4       use radinc_h, only : L_TAUMAX,naerkind
5       use aerosol_mod
6       USE tracer_h, only: noms
7       use comcstfi_mod, only: g
8       use callkeys_mod, only: pres_bottom_tropo,pres_top_tropo,obs_tau_col_tropo,  &
9                pres_bottom_strato,pres_top_strato,obs_tau_col_strato
10                 
11       implicit none
12
13!==================================================================
14!     
15!     Purpose
16!     -------
17!     Compute aerosol optical depth in each gridbox.
18!     
19!     Authors
20!     -------
21!     F. Forget
22!     F. Montmessin (water ice scheme)
23!     update J.-B. Madeleine (2008)
24!     dust removal, simplification by Robin Wordsworth (2009)
25!
26!     Input
27!     -----
28!     ngrid             Number of horizontal gridpoints
29!     nlayer            Number of layers
30!     nq                Number of tracers
31!     pplev             Pressure (Pa) at each layer boundary
32!     pq                Aerosol mixing ratio
33!     reffrad(ngrid,nlayer,naerkind)         Aerosol effective radius
34!     QREFvis3d(ngrid,nlayer,naerkind) \ 3d extinction coefficients
35!     QREFir3d(ngrid,nlayer,naerkind)  / at reference wavelengths
36!
37!     Output
38!     ------
39!     aerosol            Aerosol optical depth in layer l, grid point ig
40!     tau_col            Total column optical depth at grid point ig
41!
42!=======================================================================
43
44      INTEGER,INTENT(IN) :: ngrid  ! number of atmospheric columns
45      INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers
46      INTEGER,INTENT(IN) :: nq     ! number of tracers
47      REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa)
48      REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa)
49      REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers (.../kg_of_air)
50      REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! aerosol optical depth
51      REAL,INTENT(IN) :: reffrad(ngrid,nlayer,naerkind) ! aerosol effective radius
52      REAL,INTENT(IN) :: QREFvis3d(ngrid,nlayer,naerkind) ! extinction coefficient in the visible
53      REAL,INTENT(IN) :: QREFir3d(ngrid,nlayer,naerkind)
54      REAL,INTENT(OUT):: tau_col(ngrid) !column integrated visible optical depth
55
56      real aerosol0
57
58      INTEGER l,ig,iq,iaer
59
60      LOGICAL,SAVE :: firstcall=.true.
61!$OMP THREADPRIVATE(firstcall)
62      REAL CBRT
63      EXTERNAL CBRT
64
65      ! for fixed dust profiles
66      real expfactor, zp
67
68      ! identify tracers
69      IF (firstcall) THEN
70
71        write(*,*) "Tracers found in aeropacity:"
72
73        if (noaero) then
74          print*, "No active aerosols found in aeropacity"
75        else
76          print*, "If you would like to use aerosols, make sure any old"
77          print*, "start files are updated in newstart using the option"
78          print*, "q=0"
79          write(*,*) "Active aerosols found in aeropacity:"
80        endif
81
82        if (iaero_back2lay.ne.0) then
83          print*,'iaero_back2lay= ',iaero_back2lay
84        endif
85
86        firstcall=.false.
87      ENDIF ! of IF (firstcall)
88
89
90!     ---------------------------------------------------------   
91
92!==================================================================
93!    Two-layer aerosols (unknown composition)
94!    S. Guerlet (2013)
95!==================================================================
96
97      if (iaero_back2lay .ne.0) then
98           iaer=iaero_back2lay
99!       1. Initialization
100            aerosol(1:ngrid,1:nlayer,iaer)=0.0
101!       2. Opacity calculation
102          DO ig=1,ngrid
103           DO l=1,nlayer-1
104             aerosol(ig,l,iaer) = ( pplev(ig,l) - pplev(ig,l+1) )
105             !! 1. below tropospheric layer: no aerosols
106             IF (pplev(ig,l) .gt. pres_bottom_tropo) THEN
107               aerosol(ig,l,iaer) = 0.*aerosol(ig,l,iaer)
108             !! 2. tropo layer
109             ELSEIF (pplev(ig,l) .le. pres_bottom_tropo .and. pplev(ig,l) .ge. pres_top_tropo) THEN
110               aerosol(ig,l,iaer) = obs_tau_col_tropo*aerosol(ig,l,iaer)
111             !! 3. linear transition
112             ELSEIF (pplev(ig,l) .lt. pres_top_tropo .and. pplev(ig,l) .gt. pres_bottom_strato) THEN
113               expfactor=log(obs_tau_col_strato/obs_tau_col_tropo)/log(pres_bottom_strato/pres_top_tropo)
114               aerosol(ig,l,iaer)= obs_tau_col_tropo*((pplev(ig,l)/pres_top_tropo)**expfactor)*aerosol(ig,l,iaer)/1.5
115             !! 4. strato layer
116             ELSEIF (pplev(ig,l) .le. pres_bottom_strato .and. pplev(ig,l) .gt. pres_top_strato) THEN
117               aerosol(ig,l,iaer)= obs_tau_col_strato*aerosol(ig,l,iaer)
118             !! 5. above strato layer: no aerosols
119             ELSEIF (pplev(ig,l) .lt. pres_top_strato) THEN
120               aerosol(ig,l,iaer) = 0.*aerosol(ig,l,iaer)
121             ENDIF
122           ENDDO
123          ENDDO
124
125 !       3. Re-normalize to observed total column
126         tau_col(:)=0.0
127         DO l=1,nlayer
128          DO ig=1,ngrid
129               tau_col(ig) = tau_col(ig) &
130                     + aerosol(ig,l,iaer)/(obs_tau_col_tropo+obs_tau_col_strato)
131            ENDDO
132         ENDDO
133
134         DO ig=1,ngrid
135           DO l=1,nlayer-1
136                aerosol(ig,l,iaer)=aerosol(ig,l,iaer)/tau_col(ig)
137           ENDDO
138         ENDDO
139
140
141      end if ! if Two-layer aerosols 
142
143
144! --------------------------------------------------------------------------
145! Column integrated visible optical depth in each point (used for diagnostic)
146
147      tau_col(:)=0.0
148      do iaer = 1, naerkind
149         do l=1,nlayer
150            do ig=1,ngrid
151               tau_col(ig) = tau_col(ig) + aerosol(ig,l,iaer)
152            end do
153         end do
154      end do
155
156      do ig=1,ngrid
157         do l=1,nlayer
158            do iaer = 1, naerkind
159               if(aerosol(ig,l,iaer).gt.1.e3)then
160                  print*,'WARNING: aerosol=',aerosol(ig,l,iaer)
161                  print*,'at ig=',ig,',  l=',l,', iaer=',iaer
162                  print*,'QREFvis3d=',QREFvis3d(ig,l,iaer)
163                  print*,'reffrad=',reffrad(ig,l,iaer)
164               endif
165            end do
166         end do
167      end do
168
169      do ig=1,ngrid
170         if(tau_col(ig).gt.1.e3)then
171            print*,'WARNING: tau_col=',tau_col(ig)
172            print*,'at ig=',ig
173            print*,'aerosol=',aerosol(ig,:,:)
174            print*,'QREFvis3d=',QREFvis3d(ig,:,:)
175            print*,'reffrad=',reffrad(ig,:,:)
176         endif
177      end do
178      return
179    end subroutine aeropacity
180     
Note: See TracBrowser for help on using the repository browser.