00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_mg_coars_subgrid_3d_real ( &
00012 chfmin, chfmax, midfp, &
00013 levdim1_fine, levdim2_fine, levdim3_fine, &
00014 chcmin, chcmax, midcp, &
00015 levdim1, levdim2, levdim3, icoarse, ierror)
00016
00017
00018
00019 use PRISM_constants
00020
00021 use PSMILe, dummy_interface => PSMILe_MG_coars_subgrid_3d_real
00022
00023 implicit none
00024
00025
00026
00027 Integer, Intent (In) :: levdim1_fine, levdim2_fine,
00028 levdim3_fine
00029
00030 Real, Intent (In) :: chfmin (0:levdim1_fine,
00031 0:levdim2_fine,
00032 0:levdim3_fine)
00033 Real, Intent (In) :: chfmax (0:levdim1_fine,
00034 0:levdim2_fine,
00035 0:levdim3_fine)
00036 Real, Intent (In) :: midfp (0:levdim1_fine,
00037 0:levdim2_fine,
00038 0:levdim3_fine)
00039
00040
00041 integer, Intent (In) :: levdim1, levdim2, levdim3
00042
00043
00044
00045 Integer, Intent (In) :: icoarse (ndim_3d)
00046
00047
00048
00049
00050
00051 Real, Intent (Out) :: chcmin (0:levdim1, 0:levdim2,
00052 0:levdim3)
00053 Real, Intent (Out) :: chcmax (0:levdim1, 0:levdim2,
00054 0:levdim3)
00055 Real, Intent (Out) :: midcp (0:levdim1, 0:levdim2,
00056 0:levdim3)
00057
00058 integer, Intent (Out) :: ierror
00059
00060
00061
00062
00063
00064
00065
00066 Integer :: i, if, if1, j, jf, jf1
00067 Integer :: k, kf, kf1, m, n1
00068
00069 Integer :: inc (ndim_3d)
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089 Character(len=len_cvs_string), save :: mycvs =
00090 '$Id: psmile_mg_coars_subgrid_3d_real.F90 2325 2010-04-21 15:00:07Z valcke $'
00091
00092
00093
00094 #ifdef VERBOSE
00095 print *, trim(ch_id), ': PSMILe_MG_coars_subgrid_3d_real'
00096
00097
00098
00099
00100 call psmile_flushstd
00101 #endif /* VERBOSE */
00102
00103
00104
00105 ierror = 0
00106
00107 #ifdef PRISM_ASSERTION
00108 if (MINVAL(icoarse(:)) < 1 .or. MAXVAL (icoarse(:)) > 2) then
00109 call psmile_assert ( __FILE__, __LINE__, &
00110 'Bad coarsening: Must be in range [1:2]')
00111 endif
00112 #endif /* PRISM_ASSERTION */
00113
00114
00115
00116 n1 = 0
00117 do m = 1, ndim_3d
00118 if (icoarse(m) == 1) n1 = n1 + 1
00119 end do
00120
00121
00122
00123
00124
00125 if (n1 == 0) then
00126
00127 do k = 0, levdim3
00128 kf = 2*k
00129 kf1 = min (kf+1, levdim3_fine)
00130
00131 do j = 0, levdim2
00132 jf = 2*j
00133 jf1 = min (jf+1, levdim2_fine)
00134
00135
00136 do i = 0, levdim1
00137 if = 2*i
00138 if1 = min (if+1, levdim1_fine)
00139
00140
00141
00142 chcmax (i,j,k) = max (chfmax (if, jf, kf), &
00143 chfmax (if1, jf, kf), &
00144 chfmax (if, jf1, kf), &
00145 chfmax (if1, jf1, kf), &
00146 chfmax (if, jf, kf1), &
00147 chfmax (if1, jf, kf1), &
00148 chfmax (if, jf1, kf1), &
00149 chfmax (if1, jf1, kf1) )
00150
00151
00152
00153 chcmin (i,j,k) = min (chfmin (if, jf, kf), &
00154 chfmin (if1, jf, kf), &
00155 chfmin (if, jf1, kf), &
00156 chfmin (if1, jf1, kf), &
00157 chfmin (if, jf, kf1), &
00158 chfmin (if1, jf, kf1), &
00159 chfmin (if, jf1, kf1), &
00160 chfmin (if1, jf1, kf1) )
00161
00162
00163
00164 midcp (i,j,k) = (midfp (if, jf, kf) + &
00165 midfp (if1, jf, kf) + &
00166 midfp (if, jf1, kf) + &
00167 midfp (if1, jf1, kf) + &
00168 midfp (if, jf, kf1) + &
00169 midfp (if1, jf, kf1) + &
00170 midfp (if, jf1, kf1) + &
00171 midfp (if1, jf1, kf1)) * 0.125
00172
00173 end do
00174 end do
00175 end do
00176
00177 else if (n1 == 2) then
00178
00179
00180
00181
00182
00183 inc (:) = icoarse(:) - 1
00184
00185 do k = 0, levdim3
00186 kf = icoarse(3)*k
00187 kf1 = min (kf+inc(3), levdim3_fine)
00188
00189 do j = 0, levdim2
00190 jf = icoarse(2)*j
00191 jf1 = min (jf+inc(2), levdim2_fine)
00192
00193
00194 do i = 0, levdim1
00195 if = icoarse(1)*i
00196 if1 = min (if+inc(1), levdim1_fine)
00197
00198
00199
00200 chcmax (i,j,k) = max (chfmax (if, jf, kf ), &
00201 chfmax (if1, jf, kf ), &
00202 chfmax (if, jf1, kf ), &
00203 chfmax (if1, jf1, kf ), &
00204 chfmax (if, jf, kf1), &
00205 chfmax (if1, jf, kf1), &
00206 chfmax (if, jf1, kf1), &
00207 chfmax (if1, jf1, kf1) )
00208
00209
00210
00211 chcmin (i,j,k) = min (chfmin (if, jf, kf ), &
00212 chfmin (if1, jf, kf ), &
00213 chfmin (if, jf1, kf ), &
00214 chfmin (if1, jf1, kf ), &
00215 chfmin (if, jf, kf1), &
00216 chfmin (if1, jf, kf1), &
00217 chfmin (if, jf1, kf1), &
00218 chfmin (if1, jf1, kf1) )
00219
00220
00221
00222 midcp (i,j,k) = (midfp (if, jf, kf ) + &
00223 midfp (if1, jf, kf ) + &
00224 midfp (if, jf1, kf ) + &
00225 midfp (if1, jf1, kf ) + &
00226 midfp (if, jf, kf1) + &
00227 midfp (if1, jf, kf1) + &
00228 midfp (if, jf1, kf1) + &
00229 midfp (if1, jf1, kf1)) * 0.125
00230
00231 end do
00232 end do
00233 end do
00234
00235
00236
00237
00238
00239 else if (n1 == ndim_3d) then
00240
00241 chcmin (0:levdim1,0:levdim2,0:levdim3) = &
00242 chfmin (0:levdim1,0:levdim2,0:levdim3)
00243
00244 chcmax (0:levdim1,0:levdim2,0:levdim3) = &
00245 chfmax (0:levdim1,0:levdim2,0:levdim3)
00246
00247 midcp (0:levdim1,0:levdim2,0:levdim3) = &
00248 midfp (0:levdim1,0:levdim2,0:levdim3)
00249
00250 #ifdef UNSINN
00251 ierror = PRISM_Error_Internal
00252 call psmile_error ( ierror, 'all coarsening == 1', &
00253 icoarse, ndim_3d, __FILE__, __LINE__ )
00254 #endif /* UNSINN */
00255
00256
00257
00258
00259
00260
00261 else if (icoarse(1) == 1) then
00262
00263 do k = 0, levdim3
00264 kf = icoarse(3)*k
00265 kf1 = min (kf+1, levdim3_fine)
00266
00267 do j = 0, levdim2
00268 jf = icoarse(2)*j
00269 jf1 = min (jf+1, levdim2_fine)
00270
00271
00272 do i = 0, levdim1
00273 if = i
00274
00275
00276
00277 chcmax (i,j,k) = max (chfmax (if, jf, kf), &
00278 chfmax (if, jf1, kf), &
00279 chfmax (if, jf, kf1), &
00280 chfmax (if, jf1, kf1) )
00281
00282
00283
00284 chcmin (i,j,k) = min (chfmin (if, jf, kf), &
00285 chfmin (if, jf1, kf), &
00286 chfmin (if, jf, kf1), &
00287 chfmin (if, jf1, kf1) )
00288
00289
00290
00291 midcp (i,j,k) = (midfp (if, jf, kf) + &
00292 midfp (if, jf1, kf) + &
00293 midfp (if, jf, kf1) + &
00294 midfp (if, jf1, kf1)) * 0.25
00295
00296 end do
00297 end do
00298 end do
00299
00300
00301
00302
00303
00304 else if (icoarse(2) == 1) then
00305
00306 do k = 0, levdim3
00307 kf = 2*k
00308 kf1 = min (kf+1, levdim3_fine)
00309
00310 do j = 0, levdim2
00311 jf = j
00312
00313
00314 do i = 0, levdim1
00315 if = 2*i
00316 if1 = min (if+1, levdim1_fine)
00317
00318
00319
00320 chcmax (i,j,k) = max (chfmax (if, jf, kf ), &
00321 chfmax (if1, jf, kf ), &
00322 chfmax (if, jf, kf1), &
00323 chfmax (if1, jf, kf1) )
00324
00325
00326
00327 chcmin (i,j,k) = min (chfmin (if, jf, kf ), &
00328 chfmin (if1, jf, kf ), &
00329 chfmin (if, jf, kf1), &
00330 chfmin (if1, jf, kf1) )
00331
00332
00333
00334 midcp (i,j,k) = (midfp (if, jf, kf ) + &
00335 midfp (if1, jf, kf ) + &
00336 midfp (if, jf, kf1) + &
00337 midfp (if1, jf, kf1) ) * 0.25
00338
00339 end do
00340 end do
00341 end do
00342
00343
00344
00345
00346
00347 else
00348
00349 #ifdef PRISM_ASSERTION
00350 if (.not. (icoarse(3) == 1)) then
00351 call psmile_assert ( __FILE__, __LINE__, &
00352 'else icoarse(3)')
00353 endif
00354 #endif /* PRISM_ASSERTION */
00355
00356 do k = 0, levdim3
00357 kf = k
00358
00359 do j = 0, levdim2
00360 jf = 2*j
00361 jf1 = min (jf+1, levdim2_fine)
00362
00363
00364 do i = 0, levdim1
00365 if = 2*i
00366 if1 = min (if+1, levdim1_fine)
00367
00368
00369
00370 chcmax (i,j,k) = max (chfmax (if, jf, kf), &
00371 chfmax (if1, jf, kf), &
00372 chfmax (if, jf1, kf), &
00373 chfmax (if1, jf1, kf) )
00374
00375
00376
00377 chcmin (i,j,k) = min (chfmin (if, jf, kf), &
00378 chfmin (if1, jf, kf), &
00379 chfmin (if, jf1, kf), &
00380 chfmin (if1, jf1, kf) )
00381
00382
00383
00384 midcp (i,j,k) = (midfp (if, jf, kf) + &
00385 midfp (if1, jf, kf) + &
00386 midfp (if, jf1, kf) + &
00387 midfp (if1, jf1, kf)) * 0.25
00388
00389 end do
00390 end do
00391 end do
00392
00393 endif
00394
00395
00396
00397 #ifdef DEBUGX
00398 do k = 0, levdim3
00399 do j = 0, levdim2
00400 do i = 0, levdim1
00401 print *, 'i,j,k, chmin, chmax', i,j,k, &
00402 chcmin(i,j,k), chcmax (i,j,k)
00403 end do
00404 end do
00405 end do
00406 #endif
00407 #ifdef VERBOSE
00408
00409
00410
00411
00412
00413
00414
00415
00416 print *, trim(ch_id), ': PSMILe_MG_coars_subgrid_3d_real eof', &
00417 ': ierror =', ierror
00418
00419 call psmile_flushstd
00420 #endif /* VERBOSE */
00421
00422 end subroutine PSMILe_MG_coars_subgrid_3d_real