[1262] | 1 | ! (c) 2008, Lawrence Livermore National Security Limited Liability Corporation. |
---|
| 2 | ! All rights reserved. |
---|
[2428] | 3 | ! $Revision: 23 $, $Date: 2011-03-31 15:41:37 +0200 (jeu. 31 mars 2011) $ |
---|
| 4 | ! $URL: http://cfmip-obs-sim.googlecode.com/svn/stable/v1.4.0/llnl/prec_scops.f $ |
---|
[1262] | 5 | ! |
---|
| 6 | ! Redistribution and use in source and binary forms, with or without modification, are permitted |
---|
| 7 | ! provided that the following conditions are met: |
---|
| 8 | ! |
---|
| 9 | ! * Redistributions of source code must retain the above copyright notice, this list |
---|
| 10 | ! of conditions and the following disclaimer. |
---|
| 11 | ! * 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 materials |
---|
| 13 | ! provided with the distribution. |
---|
| 14 | ! * Neither the name of the Lawrence Livermore National Security Limited Liability Corporation |
---|
| 15 | ! nor the names of its contributors may be used to endorse or promote products derived from |
---|
| 16 | ! this software without specific prior written permission. |
---|
| 17 | ! |
---|
| 18 | ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR |
---|
| 19 | ! IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND |
---|
| 20 | ! FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR |
---|
| 21 | ! CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL |
---|
| 22 | ! DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, |
---|
| 23 | ! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER |
---|
| 24 | ! IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT |
---|
| 25 | ! OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. |
---|
| 26 | |
---|
| 27 | subroutine prec_scops(npoints,nlev,ncol,ls_p_rate,cv_p_rate, |
---|
| 28 | & frac_out,prec_frac) |
---|
| 29 | |
---|
| 30 | |
---|
| 31 | implicit none |
---|
| 32 | |
---|
| 33 | INTEGER npoints ! number of model points in the horizontal |
---|
| 34 | INTEGER nlev ! number of model levels in column |
---|
| 35 | INTEGER ncol ! number of subcolumns |
---|
| 36 | |
---|
| 37 | INTEGER i,j,ilev,ibox,cv_col |
---|
| 38 | |
---|
| 39 | REAL ls_p_rate(npoints,nlev),cv_p_rate(npoints,nlev) |
---|
| 40 | |
---|
| 41 | REAL frac_out(npoints,ncol,nlev) ! boxes gridbox divided up into |
---|
| 42 | ! Equivalent of BOX in original version, but |
---|
| 43 | ! indexed by column then row, rather than |
---|
| 44 | ! by row then column |
---|
| 45 | !TOA to SURFACE!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 46 | REAL prec_frac(npoints,ncol,nlev) ! 0 -> clear sky |
---|
| 47 | ! 1 -> LS precipitation |
---|
| 48 | ! 2 -> CONV precipitation |
---|
[2428] | 49 | ! 3 -> both |
---|
[1262] | 50 | !TOA to SURFACE!!!!!!!!!!!!!!!!!! |
---|
[2428] | 51 | |
---|
[1262] | 52 | INTEGER flag_ls, flag_cv |
---|
| 53 | INTEGER frac_out_ls(npoints,ncol),frac_out_cv(npoints,ncol) !flag variables for |
---|
| 54 | ! stratiform cloud and convective cloud in the vertical column |
---|
| 55 | |
---|
| 56 | cv_col = 0.05*ncol |
---|
| 57 | if (cv_col .eq. 0) cv_col=1 |
---|
| 58 | |
---|
| 59 | do ilev=1,nlev |
---|
[2428] | 60 | do ibox=1,ncol |
---|
| 61 | do j=1,npoints |
---|
| 62 | prec_frac(j,ibox,ilev) = 0 |
---|
| 63 | enddo |
---|
[1262] | 64 | enddo |
---|
| 65 | enddo |
---|
| 66 | |
---|
| 67 | do j=1,npoints |
---|
| 68 | do ibox=1,ncol |
---|
[2428] | 69 | frac_out_ls(j,ibox)=0 |
---|
| 70 | frac_out_cv(j,ibox)=0 |
---|
| 71 | flag_ls=0 |
---|
| 72 | flag_cv=0 |
---|
[1262] | 73 | do ilev=1,nlev |
---|
[2428] | 74 | if (frac_out(j,ibox,ilev) .eq. 1) then |
---|
| 75 | flag_ls=1 |
---|
| 76 | endif |
---|
| 77 | if (frac_out(j,ibox,ilev) .eq. 2) then |
---|
| 78 | flag_cv=1 |
---|
| 79 | endif |
---|
| 80 | enddo !loop over nlev |
---|
| 81 | if (flag_ls .eq. 1) then |
---|
| 82 | frac_out_ls(j,ibox)=1 |
---|
| 83 | endif |
---|
| 84 | if (flag_cv .eq. 1) then |
---|
| 85 | frac_out_cv(j,ibox)=1 |
---|
| 86 | endif |
---|
[1262] | 87 | enddo ! loop over ncol |
---|
| 88 | enddo ! loop over npoints |
---|
| 89 | |
---|
| 90 | ! initialize the top layer |
---|
| 91 | do j=1,npoints |
---|
| 92 | flag_ls=0 |
---|
[2428] | 93 | flag_cv=0 |
---|
| 94 | |
---|
[1262] | 95 | if (ls_p_rate(j,1) .gt. 0.) then |
---|
[2428] | 96 | do ibox=1,ncol ! possibility ONE |
---|
| 97 | if (frac_out(j,ibox,1) .eq. 1) then |
---|
| 98 | prec_frac(j,ibox,1) = 1 |
---|
| 99 | flag_ls=1 |
---|
| 100 | endif |
---|
| 101 | enddo ! loop over ncol |
---|
| 102 | if (flag_ls .eq. 0) then ! possibility THREE |
---|
| 103 | do ibox=1,ncol |
---|
| 104 | if (frac_out(j,ibox,2) .eq. 1) then |
---|
| 105 | prec_frac(j,ibox,1) = 1 |
---|
| 106 | flag_ls=1 |
---|
| 107 | endif |
---|
| 108 | enddo ! loop over ncol |
---|
| 109 | endif |
---|
| 110 | if (flag_ls .eq. 0) then ! possibility Four |
---|
| 111 | do ibox=1,ncol |
---|
| 112 | if (frac_out_ls(j,ibox) .eq. 1) then |
---|
| 113 | prec_frac(j,ibox,1) = 1 |
---|
| 114 | flag_ls=1 |
---|
| 115 | endif |
---|
| 116 | enddo ! loop over ncol |
---|
| 117 | endif |
---|
| 118 | if (flag_ls .eq. 0) then ! possibility Five |
---|
| 119 | do ibox=1,ncol |
---|
| 120 | ! prec_frac(j,1:ncol,1) = 1 |
---|
| 121 | prec_frac(j,ibox,1) = 1 |
---|
| 122 | enddo ! loop over ncol |
---|
| 123 | endif |
---|
| 124 | endif |
---|
[1262] | 125 | ! There is large scale precipitation |
---|
[2428] | 126 | |
---|
[1262] | 127 | if (cv_p_rate(j,1) .gt. 0.) then |
---|
| 128 | do ibox=1,ncol ! possibility ONE |
---|
| 129 | if (frac_out(j,ibox,1) .eq. 2) then |
---|
| 130 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
[2428] | 131 | prec_frac(j,ibox,1) = 2 |
---|
| 132 | else |
---|
| 133 | prec_frac(j,ibox,1) = 3 |
---|
| 134 | endif |
---|
| 135 | flag_cv=1 |
---|
| 136 | endif |
---|
| 137 | enddo ! loop over ncol |
---|
| 138 | if (flag_cv .eq. 0) then ! possibility THREE |
---|
| 139 | do ibox=1,ncol |
---|
| 140 | if (frac_out(j,ibox,2) .eq. 2) then |
---|
| 141 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
| 142 | prec_frac(j,ibox,1) = 2 |
---|
| 143 | else |
---|
| 144 | prec_frac(j,ibox,1) = 3 |
---|
| 145 | endif |
---|
| 146 | flag_cv=1 |
---|
| 147 | endif |
---|
| 148 | enddo ! loop over ncol |
---|
| 149 | endif |
---|
| 150 | if (flag_cv .eq. 0) then ! possibility Four |
---|
| 151 | do ibox=1,ncol |
---|
| 152 | if (frac_out_cv(j,ibox) .eq. 1) then |
---|
| 153 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
| 154 | prec_frac(j,ibox,1) = 2 |
---|
| 155 | else |
---|
| 156 | prec_frac(j,ibox,1) = 3 |
---|
| 157 | endif |
---|
| 158 | flag_cv=1 |
---|
| 159 | endif |
---|
| 160 | enddo ! loop over ncol |
---|
| 161 | endif |
---|
| 162 | if (flag_cv .eq. 0) then ! possibility Five |
---|
| 163 | do ibox=1,cv_col |
---|
| 164 | if (prec_frac(j,ibox,1) .eq. 0) then |
---|
| 165 | prec_frac(j,ibox,1) = 2 |
---|
| 166 | else |
---|
| 167 | prec_frac(j,ibox,1) = 3 |
---|
| 168 | endif |
---|
| 169 | enddo !loop over cv_col |
---|
| 170 | endif |
---|
| 171 | endif |
---|
| 172 | ! There is convective precipitation |
---|
| 173 | |
---|
| 174 | enddo ! loop over npoints |
---|
[1262] | 175 | ! end of initializing the top layer |
---|
| 176 | |
---|
| 177 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 178 | |
---|
| 179 | ! working on the levels from top to surface |
---|
| 180 | do ilev=2,nlev |
---|
| 181 | do j=1,npoints |
---|
| 182 | flag_ls=0 |
---|
[2428] | 183 | flag_cv=0 |
---|
| 184 | |
---|
[1262] | 185 | if (ls_p_rate(j,ilev) .gt. 0.) then |
---|
| 186 | do ibox=1,ncol ! possibility ONE&TWO |
---|
| 187 | if ((frac_out(j,ibox,ilev) .eq. 1) .or. |
---|
| 188 | & ((prec_frac(j,ibox,ilev-1) .eq. 1) |
---|
| 189 | & .or. (prec_frac(j,ibox,ilev-1) .eq. 3))) then |
---|
| 190 | prec_frac(j,ibox,ilev) = 1 |
---|
[2428] | 191 | flag_ls=1 |
---|
[1262] | 192 | endif |
---|
[2428] | 193 | enddo ! loop over ncol |
---|
| 194 | if ((flag_ls .eq. 0) .and. (ilev .lt. nlev)) then ! possibility THREE |
---|
| 195 | do ibox=1,ncol |
---|
| 196 | if (frac_out(j,ibox,ilev+1) .eq. 1) then |
---|
| 197 | prec_frac(j,ibox,ilev) = 1 |
---|
| 198 | flag_ls=1 |
---|
| 199 | endif |
---|
| 200 | enddo ! loop over ncol |
---|
| 201 | endif |
---|
| 202 | if (flag_ls .eq. 0) then ! possibility Four |
---|
| 203 | do ibox=1,ncol |
---|
| 204 | if (frac_out_ls(j,ibox) .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 Five |
---|
| 211 | do ibox=1,ncol |
---|
| 212 | ! prec_frac(j,1:ncol,ilev) = 1 |
---|
| 213 | prec_frac(j,ibox,ilev) = 1 |
---|
| 214 | enddo ! loop over ncol |
---|
| 215 | endif |
---|
| 216 | endif ! There is large scale precipitation |
---|
| 217 | |
---|
[1262] | 218 | if (cv_p_rate(j,ilev) .gt. 0.) then |
---|
| 219 | do ibox=1,ncol ! possibility ONE&TWO |
---|
| 220 | if ((frac_out(j,ibox,ilev) .eq. 2) .or. |
---|
| 221 | & ((prec_frac(j,ibox,ilev-1) .eq. 2) |
---|
| 222 | & .or. (prec_frac(j,ibox,ilev-1) .eq. 3))) then |
---|
| 223 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
[2428] | 224 | prec_frac(j,ibox,ilev) = 2 |
---|
| 225 | else |
---|
| 226 | prec_frac(j,ibox,ilev) = 3 |
---|
| 227 | endif |
---|
| 228 | flag_cv=1 |
---|
| 229 | endif |
---|
| 230 | enddo ! loop over ncol |
---|
| 231 | if ((flag_cv .eq. 0) .and. (ilev .lt. nlev)) then ! possibility THREE |
---|
| 232 | do ibox=1,ncol |
---|
| 233 | if (frac_out(j,ibox,ilev+1) .eq. 2) then |
---|
| 234 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
| 235 | prec_frac(j,ibox,ilev) = 2 |
---|
| 236 | else |
---|
| 237 | prec_frac(j,ibox,ilev) = 3 |
---|
| 238 | endif |
---|
| 239 | flag_cv=1 |
---|
| 240 | endif |
---|
| 241 | enddo ! loop over ncol |
---|
| 242 | endif |
---|
| 243 | if (flag_cv .eq. 0) then ! possibility Four |
---|
| 244 | do ibox=1,ncol |
---|
| 245 | if (frac_out_cv(j,ibox) .eq. 1) then |
---|
| 246 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
| 247 | prec_frac(j,ibox,ilev) = 2 |
---|
| 248 | else |
---|
| 249 | prec_frac(j,ibox,ilev) = 3 |
---|
| 250 | endif |
---|
| 251 | flag_cv=1 |
---|
| 252 | endif |
---|
| 253 | enddo ! loop over ncol |
---|
| 254 | endif |
---|
| 255 | if (flag_cv .eq. 0) then ! possibility Five |
---|
| 256 | do ibox=1,cv_col |
---|
| 257 | if (prec_frac(j,ibox,ilev) .eq. 0) then |
---|
| 258 | prec_frac(j,ibox,ilev) = 2 |
---|
| 259 | else |
---|
| 260 | prec_frac(j,ibox,ilev) = 3 |
---|
| 261 | endif |
---|
| 262 | enddo !loop over cv_col |
---|
| 263 | endif |
---|
| 264 | endif ! There is convective precipitation |
---|
| 265 | |
---|
| 266 | enddo ! loop over npoints |
---|
| 267 | enddo ! loop over nlev |
---|
[1262] | 268 | |
---|
| 269 | end |
---|
| 270 | |
---|