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