00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine prismtrs_set_neighbors3d(id_process_global_rank, &
00012 id_epio_id, &
00013 id_epio_tgt_size, &
00014 id_nb_neighcorn, &
00015 id_nbtot_links, &
00016 id_source_grid_type, &
00017 id_err)
00018
00019
00020
00021 USE PRISMDrv, dummy_interface => PRISMTrs_Set_neighbors3d
00022
00023 IMPLICIT NONE
00024
00025
00026
00027
00028
00029 INTEGER, INTENT (IN) :: id_process_global_rank
00030
00031
00032 INTEGER, INTENT (IN) :: id_epio_id
00033
00034
00035 INTEGER, INTENT (IN) :: id_epio_tgt_size
00036
00037
00038 INTEGER, INTENT (IN) :: id_nb_neighcorn
00039
00040
00041 INTEGER, INTENT (IN) :: id_nbtot_links
00042
00043
00044 INTEGER, INTENT (IN) :: id_source_grid_type
00045
00046
00047
00048 INTEGER, INTENT (Out) :: id_err
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068 CHARACTER(LEN=len_cvs_string), SAVE :: mycvs =
00069 '$Id: prismtrs_set_neighbors3d.F90 2685 2010-10-28 14:05:10Z coquart $'
00070
00071 INTEGER :: il_status(MPI_STATUS_SIZE)
00072
00073 INTEGER, PARAMETER :: nerrp=3
00074 INTEGER :: ierrp(nerrp)
00075 INTEGER :: id_loop, id_loopint, il_count
00076
00077 #ifdef DEBUGX
00078 DOUBLE PRECISION :: dbl_rad2deg
00079 DOUBLE PRECISION :: dl_latmin, dl_latmax, dl_lonmin, dl_lonmax
00080 REAL :: rea_rad2deg
00081 LOGICAL :: ll_gaussred
00082 INTEGER :: il_stride, il_unit, ii
00083 #endif
00084
00085
00086
00087 #ifdef VERBOSE
00088 PRINT *, '| | | Enter PRISMTrs_Set_neighbors3d'
00089 call psmile_flushstd
00090 #endif
00091
00092
00093
00094
00095
00096 Drv_Epios(id_epio_id)%src_grid_type = id_source_grid_type
00097
00098 IF ( id_nbtot_links == PSMILe_Undef ) THEN
00099
00100
00101
00102
00103
00104 #ifdef VERBOSE
00105 PRINT *, '| | | | | | Reception of the neighbors from ', &
00106 id_process_global_rank, ': size(',id_epio_tgt_size, &
00107 ',',id_nb_neighcorn, ')'
00108 call psmile_flushstd
00109 #endif
00110
00111 ALLOCATE( &
00112 Drv_Epios(id_epio_id)%index_array(id_epio_tgt_size,id_nb_neighcorn), &
00113 stat=id_err)
00114
00115 IF ( id_err > 0 ) THEN
00116 ierrp (1) = id_err
00117 ierrp (2) = id_epio_tgt_size
00118 ierrp (3) = id_nb_neighcorn
00119 call psmile_error_common ( PRISM_Error_Alloc, 'Epio neighbors', &
00120 ierrp, 3, __FILE__, __LINE__ )
00121 ENDIF
00122
00123 CALL MPI_Recv (Drv_Epios(id_epio_id)%index_array(1,1), &
00124 id_epio_tgt_size*id_nb_neighcorn, MPI_Integer, &
00125 id_process_global_rank, 4, comm_drv_trans, &
00126 il_status, id_err)
00127
00128 ALLOCATE( &
00129 Drv_Epios(id_epio_id)%same_lat(id_epio_tgt_size), stat=id_err)
00130
00131 IF ( id_err > 0 ) THEN
00132 ierrp (1) = id_err
00133 ierrp (2) = id_epio_tgt_size
00134 call psmile_error_common ( PRISM_Error_Alloc, 'Same lat', &
00135 ierrp, 3, __FILE__, __LINE__ )
00136 ENDIF
00137
00138 Drv_Epios(id_epio_id)%same_lat = PSMILe_Undef
00139
00140 IF (id_source_grid_type == PRISM_Gaussreduced_regvrt) THEN
00141
00142 CALL MPI_Recv (Drv_Epios(id_epio_id)%same_lat(1), &
00143 id_epio_tgt_size, MPI_Integer, &
00144 id_process_global_rank, 5, comm_drv_trans, &
00145 il_status, id_err)
00146 ENDIF
00147
00148 #ifdef DEBUGX
00149
00150 dbl_rad2deg = 360.0/6.2831853
00151 il_unit = 70+id_epio_id
00152 dl_latmin = 85.8
00153 dl_latmax = 86.8
00154 dl_lonmin = 79.8
00155 dl_lonmax = 80.8
00156 DO ii = 1, id_epio_tgt_size
00157 IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Double_Precision) THEN
00158 IF (dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(ii) <= dl_latmax &
00159 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(ii) >= dl_latmin &
00160 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(ii) <= dl_lonmax &
00161 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(ii) >= dl_lonmin) THEN
00162 WRITE(il_unit, 117) dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(ii), &
00163 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(ii)
00164 WRITE(il_unit, 118) dbl_rad2deg*Drv_Epios(id_epio_id)%src_lon_pointer_dble(Drv_Epios(id_epio_id)%index_array(ii,:))
00165 WRITE(il_unit, 119) dbl_rad2deg*Drv_Epios(id_epio_id)%src_lat_pointer_dble(Drv_Epios(id_epio_id)%index_array(ii,:))
00166 ENDIF
00167 ENDIF
00168 ENDDO
00169
00170 117 FORMAT ('tgt lon and lat', 2X, F12.4, 3X, F12.4)
00171 118 FORMAT ('src lon =', 2X, F12.4, 3X, F12.4,3X, F12.4, 3X, F12.4)
00172 119 FORMAT ('src lat =', 2X, F12.4, 3X, F12.4,3X, F12.4, 3X, F12.4)
00173 #endif
00174
00175 #ifdef VERBOSE
00176 PRINT *, '| | | | | | Reception of the neighbors done'
00177 call psmile_flushstd
00178 #endif
00179
00180 ELSE
00181
00182
00183
00184
00185
00186 #ifdef VERBOSE
00187 PRINT *, '| | | | | | Reception of the neighbour cells from ', &
00188 id_process_global_rank, ': size(',id_nbtot_links, &
00189 ',',id_nb_neighcorn, ')'
00190 call psmile_flushstd
00191 #endif
00192
00193 Drv_Epios(id_epio_id)%src_nbr_corner = id_nb_neighcorn
00194
00195 ALLOCATE( &
00196 Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_epio_tgt_size), &
00197 stat=id_err)
00198
00199 IF ( id_err > 0 ) THEN
00200 ierrp (1) = id_err
00201 ierrp (2) = id_epio_tgt_size
00202 call psmile_error_common ( PRISM_Error_Alloc, 'Epio neighbors', &
00203 ierrp, 2, __FILE__, __LINE__ )
00204 ENDIF
00205
00206 CALL MPI_Recv (Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(1), &
00207 id_epio_tgt_size, &
00208 MPI_Integer, id_process_global_rank, 5, comm_drv_trans, &
00209 il_status, id_err)
00210
00211 ALLOCATE( &
00212 Drv_Epios(id_epio_id)%index_array(id_nbtot_links,id_nb_neighcorn), &
00213 stat=id_err)
00214
00215 IF ( id_err > 0 ) THEN
00216 ierrp (1) = id_err
00217 ierrp (2) = id_nbtot_links
00218 ierrp (3) = id_nb_neighcorn
00219 call psmile_error_common ( PRISM_Error_Alloc, 'Epio neighbors', &
00220 ierrp, 3, __FILE__, __LINE__ )
00221 ENDIF
00222
00223 CALL MPI_Recv (Drv_Epios(id_epio_id)%index_array(1,1), &
00224 id_nbtot_links*id_nb_neighcorn, &
00225 MPI_Integer, id_process_global_rank, 6, comm_drv_trans, &
00226 il_status, id_err)
00227
00228 ALLOCATE( &
00229 Drv_Epios(id_epio_id)%srcepio_add(id_nbtot_links), &
00230 stat=id_err)
00231
00232 IF ( id_err > 0 ) THEN
00233 ierrp (1) = id_err
00234 ierrp (2) = id_nbtot_links
00235 call psmile_error_common ( PRISM_Error_Alloc, 'Epio compact rank', &
00236 ierrp, 2, __FILE__, __LINE__ )
00237 ENDIF
00238
00239 CALL MPI_Recv (Drv_Epios(id_epio_id)%srcepio_add(1), &
00240 id_nbtot_links, &
00241 MPI_Integer, id_process_global_rank, 7, comm_drv_trans, &
00242 il_status, id_err)
00243
00244
00245 #ifdef VERBOSE
00246 PRINT *, '| | | | | | Reception of the neighbour cell done'
00247 call psmile_flushstd
00248 #endif
00249
00250 #ifdef DEBUGX
00251 ll_gaussred = .FALSE.
00252 il_stride = 0
00253 IF (Drv_Epios(id_epio_id)%src_grid_type == PRISM_Gaussreduced_regvrt .OR. &
00254 Drv_Epios(id_epio_id)%src_grid_type == PRISM_Gaussreduced_sigmavrt) &
00255 ll_gaussred = .true.
00256 IF (ll_gaussred) il_stride = Drv_Epios(id_epio_id)%gaussred_stride
00257
00258 WRITE(89,*)'id_epio_id', id_epio_id
00259 WRITE(89,*)'id_epio_tgt_size',id_epio_tgt_size
00260 WRITE(89,*)'Drv_Epios(id_epio_id)%index_array(:,:)', Drv_Epios(id_epio_id)%index_array(:,:)
00261 WRITE(89,*)'nbsrccells_pertgtpt', Drv_Epios(id_epio_id)%nbsrccells_pertgtpt
00262 CALL flush(89)
00263 dbl_rad2deg = 360.0/6.2831853
00264 rea_rad2deg = 360.0/6.2831853
00265 il_unit=105+id_epio_id
00266 WRITE(il_unit,*)'id_epio_id', id_epio_id
00267 WRITE(il_unit,*)'id_nb_neighcorn',id_nb_neighcorn
00268 WRITE(il_unit,*)'id_nbtot_links', id_nbtot_links
00269 WRITE(il_unit,*)'il_stride = ', il_stride
00270 il_count = 0
00271 dl_latmin = 87
00272 dl_latmax = 91.
00273 dl_lonmin = 0.
00274 dl_lonmax = 7.
00275 DO id_loop = 1, id_epio_tgt_size
00276 IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Real) THEN
00277 IF (dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_real(id_loop) <= dl_latmax &
00278 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_real(id_loop) >= dl_latmin &
00279 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_real(id_loop) <= dl_lonmax &
00280 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_real(id_loop) >= dl_lonmin) THEN
00281 WRITE(il_unit,*) ' '
00282 WRITE(il_unit,*) 'id_loop (ind epiotgt) =', id_loop
00283 WRITE(il_unit,111) &
00284 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_real(id_loop), &
00285 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_real(id_loop+id_epio_tgt_size), &
00286 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_real(id_loop+2*id_epio_tgt_size), &
00287 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_real(id_loop+3*id_epio_tgt_size)
00288 WRITE(il_unit,116) &
00289 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_real(id_loop), &
00290 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_real(id_loop+id_epio_tgt_size), &
00291 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_real(id_loop+2*id_epio_tgt_size), &
00292 rea_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_real(id_loop+3*id_epio_tgt_size)
00293 WRITE(il_unit,*) 'Z= ', Drv_Epios(id_epio_id)%tgt_z_pointer_real(id_loop)
00294 WRITE(il_unit,*) 'Nbsrccells=', Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_loop)
00295 CALL flush(il_unit)
00296
00297 DO id_loopint = 1, Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_loop)
00298 WRITE(il_unit,112) id_loopint
00299 il_count = il_count + 1
00300 WRITE(il_unit,114) rea_rad2deg*Drv_Epios(id_epio_id)%src_lon_pointer_real(Drv_Epios(id_epio_id) &
00301 %index_array(il_count,:))
00302 WRITE(il_unit,114) rea_rad2deg*Drv_Epios(id_epio_id)%src_lat_pointer_real(Drv_Epios(id_epio_id) &
00303 %index_array(il_count,:))
00304 WRITE(il_unit,114) Drv_Epios(id_epio_id)%src_z_pointer_real(Drv_Epios(id_epio_id)%index_array(il_count,:))
00305 CALL flush(il_unit)
00306 ENDDO
00307 ELSE
00308 il_count=il_count+Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_loop)
00309 ENDIF
00310 ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Double_Precision) THEN
00311 IF (dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(id_loop) <= dl_latmax &
00312 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(id_loop) >= dl_latmin &
00313 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(id_loop) <= dl_lonmax &
00314 .AND. dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(id_loop) >= dl_lonmin) THEN
00315 WRITE(il_unit,*) ' '
00316 WRITE(il_unit,*) 'id_loop =', id_loop
00317 WRITE(il_unit,121) &
00318 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(id_loop), &
00319 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(id_loop+id_epio_tgt_size), &
00320 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(id_loop+2*id_epio_tgt_size), &
00321 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lon_pointer_dble(id_loop+3*id_epio_tgt_size)
00322 WRITE(il_unit,126) &
00323 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(id_loop), &
00324 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(id_loop+id_epio_tgt_size), &
00325 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(id_loop+2*id_epio_tgt_size), &
00326 dbl_rad2deg*Drv_Epios(id_epio_id)%tgt_lat_pointer_dble(id_loop+3*id_epio_tgt_size)
00327 WRITE(il_unit,*) 'Z= ', Drv_Epios(id_epio_id)%tgt_z_pointer_dble(id_loop)
00328 WRITE(il_unit,*) 'Nbsrccells=', Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_loop)
00329 CALL flush(il_unit)
00330
00331 DO id_loopint = 1, Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_loop)
00332 WRITE(il_unit,112) id_loopint
00333 112 FORMAT ('associated source cell number', 1X, I2)
00334 CALL flush(il_unit)
00335 il_count = il_count + 1
00336 WRITE(il_unit,*) 'il_count= ', il_count
00337 WRITE(il_unit,*) 'EPIO source number :',Drv_Epios(id_epio_id)%srcepio_add(il_count)
00338 CALL flush(il_unit)
00339 IF (ll_gaussred) THEN
00340 WRITE(il_unit,125) dbl_rad2deg*Drv_Epios(id_epio_id)%src_lon_pointer_dble(Drv_Epios(id_epio_id) &
00341 %index_array(il_count,1)), dbl_rad2deg*Drv_Epios(id_epio_id) &
00342 %src_lon_pointer_dble(Drv_Epios(id_epio_id)%index_array(il_count,1)+il_stride)
00343 CALL flush(il_unit)
00344 WRITE(il_unit,125) dbl_rad2deg*Drv_Epios(id_epio_id)%src_lat_pointer_dble(Drv_Epios(id_epio_id) &
00345 %index_array(il_count,1)), dbl_rad2deg*Drv_Epios(id_epio_id) &
00346 %src_lat_pointer_dble(Drv_Epios(id_epio_id)%index_array(il_count,1)+il_stride)
00347 CALL flush(il_unit)
00348 WRITE(il_unit,125) dbl_rad2deg*Drv_Epios(id_epio_id)%src_z_pointer_dble(Drv_Epios(id_epio_id) &
00349 %index_array(il_count,1)), dbl_rad2deg*Drv_Epios(id_epio_id) &
00350 %src_z_pointer_dble(Drv_Epios(id_epio_id)%index_array(il_count,1)+il_stride)
00351 CALL flush(il_unit)
00352 ELSE
00353 WRITE(il_unit,114) dbl_rad2deg*Drv_Epios(id_epio_id)%src_lon_pointer_dble(Drv_Epios(id_epio_id) &
00354 %index_array(il_count,:))
00355 CALL flush(il_unit)
00356 WRITE(il_unit,114) dbl_rad2deg*Drv_Epios(id_epio_id)%src_lat_pointer_dble(Drv_Epios(id_epio_id) &
00357 %index_array(il_count,:))
00358 CALL flush(il_unit)
00359 WRITE(il_unit,114) dbl_rad2deg*Drv_Epios(id_epio_id)%src_z_pointer_dble(Drv_Epios(id_epio_id) &
00360 %index_array(il_count,:))
00361 CALL flush(il_unit)
00362 ENDIF
00363 ENDDO
00364 ELSE
00365 il_count=il_count+Drv_Epios(id_epio_id)%nbsrccells_pertgtpt(id_loop)
00366 ENDIF
00367 ENDIF
00368 ENDDO
00369 #endif
00370
00371
00372 #ifdef VERBOSE
00373 PRINT *, '| | | | | | Reception of the source epio compact rank'
00374 call psmile_flushstd
00375 #endif
00376 #ifdef DEBUGX
00377 WRITE(89,*)'srcepio_add', Drv_Epios(id_epio_id)%srcepio_add
00378 #endif
00379 ENDIF
00380 #ifdef DEBUGX
00381 111 FORMAT ('tgt lon =', 2X, F12.4, 3X, F12.4,3X, F12.4, 3X, F12.4)
00382 121 FORMAT ('tgt lon =', 2X, D12.4, 3X, D12.4,3X, D12.4, 3X, D12.4)
00383 116 FORMAT ('tgt lat =', 2X, F12.4, 3X, F12.4,3X, F12.4, 3X, F12.4)
00384 126 FORMAT ('tgt lat =', 2X, D12.4, 3X, D12.4,3X, D12.4, 3X, D12.4)
00385 114 FORMAT (2X, F12.4, 3X, F12.4, 3X, F12.4, 3X, F12.4)
00386 115 FORMAT (2X, F12.4, 3X, F12.4)
00387 125 FORMAT (2X, D12.4, 3X, D12.4)
00388 #endif
00389
00390 #ifdef VERBOSE
00391 PRINT *, '| | | Quit PRISMTrs_Set_neighbors3d'
00392 call psmile_flushstd
00393 #endif
00394
00395 END SUBROUTINE PRISMTrs_Set_neighbors3d
00396
00397
00398
00399