1 | ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
2 | ! Copyright (c) 2008, Lawrence Livermore National Security Limited Liability Corporation |
---|
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 | ! May 2015- D. Swales - Modified for COSPv2.0 |
---|
31 | |
---|
32 | ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
33 | module mod_prec_scops |
---|
34 | implicit none |
---|
35 | contains |
---|
36 | |
---|
37 | subroutine prec_scops(npoints,nlev,ncol,ls_p_rate,cv_p_rate,frac_out,prec_frac) |
---|
38 | |
---|
39 | USE COSP_KINDS, ONLY: wp |
---|
40 | use mod_cosp_config, ONLY: scops_ccfrac |
---|
41 | |
---|
42 | INTEGER npoints ! number of model points in the horizontal |
---|
43 | INTEGER nlev ! number of model levels in column |
---|
44 | INTEGER ncol ! number of subcolumns |
---|
45 | |
---|
46 | INTEGER j,ilev,ibox,cv_col |
---|
47 | |
---|
48 | REAL(WP) ls_p_rate(npoints,nlev),cv_p_rate(npoints,nlev) |
---|
49 | |
---|
50 | REAL(WP) frac_out(npoints,ncol,nlev) ! boxes gridbox divided up into |
---|
51 | ! Equivalent of BOX in original version, but |
---|
52 | ! indexed by column then row, rather than |
---|
53 | ! by row then column |
---|
54 | !TOA to SURFACE!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
55 | REAL(WP) prec_frac(npoints,ncol,nlev) ! 0 -> clear sky |
---|
56 | ! 1 -> LS precipitation |
---|
57 | ! 2 -> CONV precipitation |
---|
58 | ! 3 -> both |
---|
59 | !TOA to SURFACE!!!!!!!!!!!!!!!!!! |
---|
60 | |
---|
61 | INTEGER flag_ls, flag_cv |
---|
62 | INTEGER frac_out_ls(npoints,ncol),frac_out_cv(npoints,ncol) !flag variables for |
---|
63 | ! stratiform cloud and convective cloud in the vertical column |
---|
64 | |
---|
65 | cv_col = scops_ccfrac*ncol |
---|
66 | if (cv_col .eq. 0) cv_col=1 |
---|
67 | |
---|
68 | do ilev=1,nlev |
---|
69 | do ibox=1,ncol |
---|
70 | do j=1,npoints |
---|
71 | prec_frac(j,ibox,ilev) = 0 |
---|
72 | enddo |
---|
73 | enddo |
---|
74 | enddo |
---|
75 | |
---|
76 | do j=1,npoints |
---|
77 | do ibox=1,ncol |
---|
78 | frac_out_ls(j,ibox)=0 |
---|
79 | frac_out_cv(j,ibox)=0 |
---|
80 | flag_ls=0 |
---|
81 | flag_cv=0 |
---|
82 | do ilev=1,nlev |
---|
83 | if (frac_out(j,ibox,ilev) .eq. 1) then |
---|
84 | flag_ls=1 |
---|
85 | endif |
---|
86 | if (frac_out(j,ibox,ilev) .eq. 2) then |
---|
87 | flag_cv=1 |
---|
88 | endif |
---|
89 | enddo !loop over nlev |
---|
90 | if (flag_ls .eq. 1) then |
---|
91 | frac_out_ls(j,ibox)=1 |
---|
92 | endif |
---|
93 | if (flag_cv .eq. 1) then |
---|
94 | frac_out_cv(j,ibox)=1 |
---|
95 | endif |
---|
96 | enddo ! loop over ncol |
---|
97 | enddo ! loop over npoints |
---|
98 | |
---|
99 | ! initialize the top layer |
---|
100 | do j=1,npoints |
---|
101 | flag_ls=0 |
---|
102 | flag_cv=0 |
---|
103 | |
---|
104 | if (ls_p_rate(j,1) .gt. 0.) then |
---|
105 | do ibox=1,ncol ! possibility ONE |
---|
106 | if (frac_out(j,ibox,1) .eq. 1) then |
---|
107 | prec_frac(j,ibox,1) = 1 |
---|
108 | flag_ls=1 |
---|
109 | endif |
---|
110 | enddo ! loop over ncol |
---|
111 | if (flag_ls .eq. 0) then ! possibility THREE |
---|
112 | do ibox=1,ncol |
---|
113 | if (frac_out(j,ibox,2) .eq. 1) then |
---|
114 | prec_frac(j,ibox,1) = 1 |
---|
115 | flag_ls=1 |
---|
116 | endif |
---|
117 | enddo ! loop over ncol |
---|
118 | endif |
---|
119 | if (flag_ls .eq. 0) then ! possibility Four |
---|
120 | do ibox=1,ncol |
---|
121 | if (frac_out_ls(j,ibox) .eq. 1) then |
---|
122 | prec_frac(j,ibox,1) = 1 |
---|
123 | flag_ls=1 |
---|
124 | endif |
---|
125 | enddo ! loop over ncol |
---|
126 | endif |
---|
127 | if (flag_ls .eq. 0) then ! possibility Five |
---|
128 | do ibox=1,ncol |
---|
129 | ! prec_frac(j,1:ncol,1) = 1 |
---|
130 | prec_frac(j,ibox,1) = 1 |
---|
131 | enddo ! loop over ncol |
---|
132 | endif |
---|
133 | endif |
---|
134 | ! There is large scale precipitation |
---|
135 | |
---|
136 | if (cv_p_rate(j,1) .gt. 0.) then |
---|
137 | do ibox=1,ncol ! possibility ONE |
---|
138 | if (frac_out(j,ibox,1) .eq. 2) then |
---|
139 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
140 | prec_frac(j,ibox,1) = 2 |
---|
141 | else |
---|
142 | prec_frac(j,ibox,1) = 3 |
---|
143 | endif |
---|
144 | flag_cv=1 |
---|
145 | endif |
---|
146 | enddo ! loop over ncol |
---|
147 | if (flag_cv .eq. 0) then ! possibility THREE |
---|
148 | do ibox=1,ncol |
---|
149 | if (frac_out(j,ibox,2) .eq. 2) then |
---|
150 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
151 | prec_frac(j,ibox,1) = 2 |
---|
152 | else |
---|
153 | prec_frac(j,ibox,1) = 3 |
---|
154 | endif |
---|
155 | flag_cv=1 |
---|
156 | endif |
---|
157 | enddo ! loop over ncol |
---|
158 | endif |
---|
159 | if (flag_cv .eq. 0) then ! possibility Four |
---|
160 | do ibox=1,ncol |
---|
161 | if (frac_out_cv(j,ibox) .eq. 1) then |
---|
162 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
163 | prec_frac(j,ibox,1) = 2 |
---|
164 | else |
---|
165 | prec_frac(j,ibox,1) = 3 |
---|
166 | endif |
---|
167 | flag_cv=1 |
---|
168 | endif |
---|
169 | enddo ! loop over ncol |
---|
170 | endif |
---|
171 | if (flag_cv .eq. 0) then ! possibility Five |
---|
172 | do ibox=1,cv_col |
---|
173 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
174 | prec_frac(j,ibox,1) = 2 |
---|
175 | else |
---|
176 | prec_frac(j,ibox,1) = 3 |
---|
177 | endif |
---|
178 | enddo !loop over cv_col |
---|
179 | endif |
---|
180 | endif |
---|
181 | ! There is convective precipitation |
---|
182 | |
---|
183 | enddo ! loop over npoints |
---|
184 | ! end of initializing the top layer |
---|
185 | |
---|
186 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
187 | |
---|
188 | ! working on the levels from top to surface |
---|
189 | do ilev=2,nlev |
---|
190 | do j=1,npoints |
---|
191 | flag_ls=0 |
---|
192 | flag_cv=0 |
---|
193 | |
---|
194 | if (ls_p_rate(j,ilev) .gt. 0.) then |
---|
195 | do ibox=1,ncol ! possibility ONE&TWO |
---|
196 | if ((frac_out(j,ibox,ilev) .eq. 1) .or. ((prec_frac(j,ibox,ilev-1) .eq. 1) & |
---|
197 | .or. (prec_frac(j,ibox,ilev-1) .eq. 3))) then |
---|
198 | prec_frac(j,ibox,ilev) = 1 |
---|
199 | flag_ls=1 |
---|
200 | endif |
---|
201 | enddo ! loop over ncol |
---|
202 | if ((flag_ls .eq. 0) .and. (ilev .lt. nlev)) then ! possibility THREE |
---|
203 | do ibox=1,ncol |
---|
204 | if (frac_out(j,ibox,ilev+1) .eq. 1) then |
---|
205 | prec_frac(j,ibox,ilev) = 1 |
---|
206 | flag_ls=1 |
---|
207 | endif |
---|
208 | enddo ! loop over ncol |
---|
209 | endif |
---|
210 | if (flag_ls .eq. 0) then ! possibility Four |
---|
211 | do ibox=1,ncol |
---|
212 | if (frac_out_ls(j,ibox) .eq. 1) then |
---|
213 | prec_frac(j,ibox,ilev) = 1 |
---|
214 | flag_ls=1 |
---|
215 | endif |
---|
216 | enddo ! loop over ncol |
---|
217 | endif |
---|
218 | if (flag_ls .eq. 0) then ! possibility Five |
---|
219 | do ibox=1,ncol |
---|
220 | ! prec_frac(j,1:ncol,ilev) = 1 |
---|
221 | prec_frac(j,ibox,ilev) = 1 |
---|
222 | enddo ! loop over ncol |
---|
223 | endif |
---|
224 | endif ! There is large scale precipitation |
---|
225 | |
---|
226 | if (cv_p_rate(j,ilev) .gt. 0.) then |
---|
227 | do ibox=1,ncol ! possibility ONE&TWO |
---|
228 | if ((frac_out(j,ibox,ilev) .eq. 2) .or. ((prec_frac(j,ibox,ilev-1) .eq. 2) & |
---|
229 | .or. (prec_frac(j,ibox,ilev-1) .eq. 3))) then |
---|
230 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
231 | prec_frac(j,ibox,ilev) = 2 |
---|
232 | else |
---|
233 | prec_frac(j,ibox,ilev) = 3 |
---|
234 | endif |
---|
235 | flag_cv=1 |
---|
236 | endif |
---|
237 | enddo ! loop over ncol |
---|
238 | if ((flag_cv .eq. 0) .and. (ilev .lt. nlev)) then ! possibility THREE |
---|
239 | do ibox=1,ncol |
---|
240 | if (frac_out(j,ibox,ilev+1) .eq. 2) then |
---|
241 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
242 | prec_frac(j,ibox,ilev) = 2 |
---|
243 | else |
---|
244 | prec_frac(j,ibox,ilev) = 3 |
---|
245 | endif |
---|
246 | flag_cv=1 |
---|
247 | endif |
---|
248 | enddo ! loop over ncol |
---|
249 | endif |
---|
250 | if (flag_cv .eq. 0) then ! possibility Four |
---|
251 | do ibox=1,ncol |
---|
252 | if (frac_out_cv(j,ibox) .eq. 1) then |
---|
253 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
254 | prec_frac(j,ibox,ilev) = 2 |
---|
255 | else |
---|
256 | prec_frac(j,ibox,ilev) = 3 |
---|
257 | endif |
---|
258 | flag_cv=1 |
---|
259 | endif |
---|
260 | enddo ! loop over ncol |
---|
261 | endif |
---|
262 | if (flag_cv .eq. 0) then ! possibility Five |
---|
263 | do ibox=1,cv_col |
---|
264 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
265 | prec_frac(j,ibox,ilev) = 2 |
---|
266 | else |
---|
267 | prec_frac(j,ibox,ilev) = 3 |
---|
268 | endif |
---|
269 | enddo !loop over cv_col |
---|
270 | endif |
---|
271 | endif ! There is convective precipitation |
---|
272 | |
---|
273 | enddo ! loop over npoints |
---|
274 | enddo ! loop over nlev |
---|
275 | |
---|
276 | end subroutine prec_scops |
---|
277 | end module mod_prec_scops |
---|