1 | ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
2 | ! Copyright (c) 2009, Centre National de la Recherche Scientifique |
---|
3 | ! All rights reserved. |
---|
4 | |
---|
5 | ! Redistribution and use in source and binary forms, with or without modification, are |
---|
6 | ! permitted provided that the following conditions are met: |
---|
7 | |
---|
8 | ! 1. Redistributions of source code must retain the above copyright notice, this list of |
---|
9 | ! conditions and the following disclaimer. |
---|
10 | |
---|
11 | ! 2. Redistributions in binary form must reproduce the above copyright notice, this list |
---|
12 | ! of conditions and the following disclaimer in the documentation and/or other |
---|
13 | ! materials provided with the distribution. |
---|
14 | |
---|
15 | ! 3. Neither the name of the copyright holder nor the names of its contributors may be |
---|
16 | ! used to endorse or promote products derived from this software without specific prior |
---|
17 | ! written permission. |
---|
18 | |
---|
19 | ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY |
---|
20 | ! EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF |
---|
21 | ! MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL |
---|
22 | ! THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, |
---|
23 | ! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT |
---|
24 | ! OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS |
---|
25 | ! INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT |
---|
26 | ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE |
---|
27 | ! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. |
---|
28 | |
---|
29 | ! History |
---|
30 | ! December 2008, S. Bony, H. Chepfer and J-L. Dufresne : |
---|
31 | ! - optimization for vectorization |
---|
32 | ! Version 2.0 (October 2008) |
---|
33 | ! Version 2.1 (December 2008) |
---|
34 | ! May 2015 - D. Swales - Modified for COSPv2.0 |
---|
35 | ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
36 | module mod_parasol |
---|
37 | USE COSP_KINDS, ONLY: wp |
---|
38 | USE COSP_MATH_CONSTANTS, ONLY: pi |
---|
39 | use mod_cosp_config, ONLY: R_UNDEF,PARASOL_NREFL,PARASOL_NTAU,PARASOL_TAU,PARASOL_SZA,rlumA,rlumB |
---|
40 | implicit none |
---|
41 | |
---|
42 | contains |
---|
43 | SUBROUTINE parasol_subcolumn(npoints,nrefl,tautot_S_liq,tautot_S_ice,refl) |
---|
44 | ! ########################################################################## |
---|
45 | ! Purpose: To compute Parasol reflectance signal from model-simulated profiles |
---|
46 | ! of cloud water and cloud fraction in each sub-column of each model |
---|
47 | ! gridbox. |
---|
48 | |
---|
49 | |
---|
50 | ! December 2008, S. Bony, H. Chepfer and J-L. Dufresne : |
---|
51 | ! - optimization for vectorization |
---|
52 | |
---|
53 | ! Version 2.0 (October 2008) |
---|
54 | ! Version 2.1 (December 2008) |
---|
55 | ! ########################################################################## |
---|
56 | |
---|
57 | ! INPUTS |
---|
58 | INTEGER,intent(in) :: & |
---|
59 | npoints, & ! Number of horizontal gridpoints |
---|
60 | nrefl ! Number of angles for which the reflectance is computed |
---|
61 | REAL(WP),intent(inout),dimension(npoints) :: & |
---|
62 | tautot_S_liq, & ! Liquid water optical thickness, from TOA to SFC |
---|
63 | tautot_S_ice ! Ice water optical thickness, from TOA to SFC |
---|
64 | ! OUTPUTS |
---|
65 | REAL(WP),intent(inout),dimension(npoints,nrefl) :: & |
---|
66 | refl ! Parasol reflectances |
---|
67 | |
---|
68 | ! LOCAL VARIABLES |
---|
69 | REAL(WP),dimension(npoints) :: & |
---|
70 | tautot_S, & ! Cloud optical thickness, from TOA to surface |
---|
71 | frac_taucol_liq, & ! |
---|
72 | frac_taucol_ice ! |
---|
73 | |
---|
74 | ! Look up table variables: |
---|
75 | INTEGER :: ny,it |
---|
76 | REAL(WP),dimension(PARASOL_NREFL) :: r_norm |
---|
77 | REAL(WP),dimension(PARASOL_NREFL,PARASOL_NTAU-1) :: aa,ab,ba,bb |
---|
78 | REAL(WP),dimension(npoints,5) :: rlumA_mod,rlumB_mod |
---|
79 | |
---|
80 | !-------------------------------------------------------------------------------- |
---|
81 | ! Lum_norm=f(PARASOL_SZA,tau_cloud) derived from adding-doubling calculations |
---|
82 | ! valid ONLY ABOVE OCEAN (albedo_sfce=5%) |
---|
83 | ! valid only in one viewing direction (theta_v=30�, phi_s-phi_v=320�) |
---|
84 | ! based on adding-doubling radiative transfer computation |
---|
85 | ! for PARASOL_TAU values (0 to 100) and for PARASOL_SZA values (0 to 80) |
---|
86 | ! for 2 scattering phase functions: liquid spherical, ice non spherical |
---|
87 | |
---|
88 | ! Initialize |
---|
89 | rlumA_mod(1:npoints,1:5) = 0._wp |
---|
90 | rlumB_mod(1:npoints,1:5) = 0._wp |
---|
91 | |
---|
92 | r_norm(1:PARASOL_NREFL)=1._wp/ cos(pi/180._wp*PARASOL_SZA(1:PARASOL_NREFL)) |
---|
93 | |
---|
94 | tautot_S_liq(1:npoints) = max(tautot_S_liq(1:npoints),PARASOL_TAU(1)) |
---|
95 | tautot_S_ice(1:npoints) = max(tautot_S_ice(1:npoints),PARASOL_TAU(1)) |
---|
96 | tautot_S(1:npoints) = tautot_S_ice(1:npoints) + tautot_S_liq(1:npoints) |
---|
97 | |
---|
98 | ! Relative fraction of the opt. thick due to liquid or ice clouds |
---|
99 | WHERE (tautot_S(1:npoints) .gt. 0.) |
---|
100 | frac_taucol_liq(1:npoints) = tautot_S_liq(1:npoints) / tautot_S(1:npoints) |
---|
101 | frac_taucol_ice(1:npoints) = tautot_S_ice(1:npoints) / tautot_S(1:npoints) |
---|
102 | ELSEWHERE |
---|
103 | frac_taucol_liq(1:npoints) = 1._wp |
---|
104 | frac_taucol_ice(1:npoints) = 0._wp |
---|
105 | END WHERE |
---|
106 | tautot_S(1:npoints)=MIN(tautot_S(1:npoints),PARASOL_TAU(PARASOL_NTAU)) |
---|
107 | |
---|
108 | ! Linear interpolation |
---|
109 | DO ny=1,PARASOL_NTAU-1 |
---|
110 | ! Microphysics A (liquid clouds) |
---|
111 | aA(1:PARASOL_NREFL,ny) = (rlumA(1:PARASOL_NREFL,ny+1)-rlumA(1:PARASOL_NREFL,ny))/(PARASOL_TAU(ny+1)-PARASOL_TAU(ny)) |
---|
112 | bA(1:PARASOL_NREFL,ny) = rlumA(1:PARASOL_NREFL,ny) - aA(1:PARASOL_NREFL,ny)*PARASOL_TAU(ny) |
---|
113 | ! Microphysics B (ice clouds) |
---|
114 | aB(1:PARASOL_NREFL,ny) = (rlumB(1:PARASOL_NREFL,ny+1)-rlumB(1:PARASOL_NREFL,ny))/(PARASOL_TAU(ny+1)-PARASOL_TAU(ny)) |
---|
115 | bB(1:PARASOL_NREFL,ny) = rlumB(1:PARASOL_NREFL,ny) - aB(1:PARASOL_NREFL,ny)*PARASOL_TAU(ny) |
---|
116 | ENDDO |
---|
117 | |
---|
118 | DO it=1,PARASOL_NREFL |
---|
119 | DO ny=1,PARASOL_NTAU-1 |
---|
120 | WHERE (tautot_S(1:npoints) .ge. PARASOL_TAU(ny).and. & |
---|
121 | tautot_S(1:npoints) .le. PARASOL_TAU(ny+1)) |
---|
122 | rlumA_mod(1:npoints,it) = aA(it,ny)*tautot_S(1:npoints) + bA(it,ny) |
---|
123 | rlumB_mod(1:npoints,it) = aB(it,ny)*tautot_S(1:npoints) + bB(it,ny) |
---|
124 | END WHERE |
---|
125 | END DO |
---|
126 | END DO |
---|
127 | |
---|
128 | DO it=1,PARASOL_NREFL |
---|
129 | refl(1:npoints,it) = frac_taucol_liq(1:npoints) * rlumA_mod(1:npoints,it) & |
---|
130 | + frac_taucol_ice(1:npoints) * rlumB_mod(1:npoints,it) |
---|
131 | ! Normalized radiance -> reflectance: |
---|
132 | refl(1:npoints,it) = refl(1:npoints,it) * r_norm(it) |
---|
133 | ENDDO |
---|
134 | |
---|
135 | RETURN |
---|
136 | END SUBROUTINE parasol_subcolumn |
---|
137 | ! ###################################################################################### |
---|
138 | ! SUBROUTINE parasol_gridbox |
---|
139 | ! ###################################################################################### |
---|
140 | subroutine parasol_column(npoints,nrefl,ncol,land,refl,parasolrefl) |
---|
141 | |
---|
142 | ! Inputs |
---|
143 | integer,intent(in) :: & |
---|
144 | npoints, & ! Number of horizontal grid points |
---|
145 | ncol, & ! Number of subcolumns |
---|
146 | nrefl ! Number of solar zenith angles for parasol reflectances |
---|
147 | real(wp),intent(in),dimension(npoints) :: & |
---|
148 | land ! Landmask [0 - Ocean, 1 - Land] |
---|
149 | real(wp),intent(in),dimension(npoints,ncol,nrefl) :: & |
---|
150 | refl ! Subgrid parasol reflectance ! parasol |
---|
151 | |
---|
152 | ! Outputs |
---|
153 | real(wp),intent(out),dimension(npoints,nrefl) :: & |
---|
154 | parasolrefl ! Grid-averaged parasol reflectance |
---|
155 | |
---|
156 | ! Local variables |
---|
157 | integer :: k,ic |
---|
158 | |
---|
159 | ! Compute grid-box averaged Parasol reflectances |
---|
160 | parasolrefl(:,:) = 0._wp |
---|
161 | do k = 1, nrefl |
---|
162 | do ic = 1, ncol |
---|
163 | parasolrefl(:,k) = parasolrefl(:,k) + refl(:,ic,k) |
---|
164 | enddo |
---|
165 | enddo |
---|
166 | |
---|
167 | do k = 1, nrefl |
---|
168 | parasolrefl(:,k) = parasolrefl(:,k) / float(ncol) |
---|
169 | ! if land=1 -> parasolrefl=R_UNDEF |
---|
170 | ! if land=0 -> parasolrefl=parasolrefl |
---|
171 | parasolrefl(:,k) = parasolrefl(:,k) * MAX(1._wp-land(:),0.0) & |
---|
172 | + (1._wp - MAX(1._wp-land(:),0.0))*R_UNDEF |
---|
173 | enddo |
---|
174 | end subroutine parasol_column |
---|
175 | |
---|
176 | end module mod_parasol |
---|