00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 subroutine prismtrs_bcast2trs(id_err)
00011
00012
00013
00014
00015 USE PRISMDrv, dummy_interface => PRISMTrs_bcast2trs
00016
00017 IMPLICIT NONE
00018
00019 Integer, Intent (Out) :: id_err
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038 CHARACTER(LEN=len_cvs_string), SAVE :: mycvs =
00039 '$Id: prismtrs_bcast2trs.F90 2399 2010-06-21 08:09:39Z coquart $'
00040
00041 DOUBLE PRECISION, ALLOCATABLE :: dla_buffer(:)
00042
00043 INTEGER, DIMENSION(:,:), ALLOCATABLE :: ila_info_exch
00044 INTEGER, PARAMETER :: ila_info_exch_size = 12
00045 INTEGER, DIMENSION(:,:), ALLOCATABLE :: ila_info_interp
00046 INTEGER :: il_size, il_rank, ib
00047 INTEGER, DIMENSION(2) :: ila_bcast
00048
00049 INTEGER, PARAMETER :: nerrp = 2
00050 INTEGER :: ierrp (nerrp)
00051
00052
00053
00054 #ifdef VERBOSE
00055 PRINT *, '| | Enter PRISMDrv_bcast2trs'
00056 call psmile_flushstd
00057 #endif
00058
00059 id_err = 0
00060
00061
00062
00063
00064 CALL MPI_Comm_size(comm_drv_local, il_size, id_err)
00065 CALL MPI_Comm_rank(comm_drv_local, il_rank, id_err)
00066
00067
00068
00069
00070 IF (il_size .gt. 1) THEN
00071
00072
00073 IF (il_rank == PRISM_root) THEN
00074 ila_bcast(1) = Number_of_Exchanges
00075 ila_bcast(2) = Number_of_Interps
00076 END IF
00077
00078 CALL MPI_Bcast( ila_bcast, 2, MPI_Integer, &
00079 PRISM_root, comm_drv_local, id_err )
00080
00081 IF ( id_err /= MPI_SUCCESS ) THEN
00082 ierrp (1) = id_err
00083 id_err = PRISM_Error_MPI
00084
00085 call psmile_error_common ( id_err, 'MPI_Bcast', &
00086 ierrp, 1, __FILE__, __LINE__ )
00087 RETURN
00088 ENDIF
00089
00090 IF (il_rank .ge. PRISM_root+1) THEN
00091 Number_of_Exchanges = ila_bcast(1)
00092 Number_of_Interps = ila_bcast(2)
00093 END IF
00094
00095
00096 IF (il_rank .ne. PRISM_root) THEN
00097 ALLOCATE(Drv_Exchanges(Number_of_Exchanges), stat = id_err)
00098 IF (id_err > 0) THEN
00099 ierrp (1) = id_err
00100 ierrp (2) = Number_of_Exchanges
00101 id_err = 13
00102
00103 call psmile_error_common ( id_err, 'Drv_Exchanges', &
00104 ierrp, 2, __FILE__, __LINE__ )
00105 RETURN
00106 ENDIF
00107 END IF
00108
00109
00110 IF (il_rank .ne. PRISM_root) THEN
00111 ALLOCATE(Drv_Interps(Number_of_Interps), stat = id_err)
00112 IF (id_err > 0) THEN
00113 ierrp (1) = id_err
00114 ierrp (2) = Number_of_Interps
00115 id_err = 13
00116
00117 call psmile_error_common ( id_err, 'Drv_Interps', &
00118 ierrp, 2, __FILE__, __LINE__ )
00119 RETURN
00120 ENDIF
00121 END IF
00122
00123
00124 ALLOCATE(ila_info_exch(Number_of_Exchanges,ila_info_exch_size), stat = id_err)
00125 IF (id_err > 0) THEN
00126 ierrp (1) = id_err
00127 ierrp (2) = Number_of_Exchanges
00128 id_err = 13
00129
00130 call psmile_error_common ( id_err, 'ila_info_exch', &
00131 ierrp, 2, __FILE__, __LINE__ )
00132 RETURN
00133 ENDIF
00134
00135 ALLOCATE(ila_info_interp(Number_of_Interps,27), stat = id_err)
00136 IF (id_err > 0) THEN
00137 ierrp (1) = id_err
00138 ierrp (2) = Number_of_Interps
00139 id_err = 13
00140
00141 call psmile_error_common ( id_err, 'ila_info_interp', &
00142 ierrp, 2, __FILE__, __LINE__ )
00143 RETURN
00144 ENDIF
00145
00146
00147 IF (il_rank .eq. PRISM_root) THEN
00148
00149 DO ib = 1, Number_of_Exchanges
00150
00151 ila_info_exch(ib,1) = Drv_Exchanges(ib)%trans_out_id
00152 ila_info_exch(ib,2) = Drv_Exchanges(ib)%trans_in_id
00153 ila_info_exch(ib,3) = Drv_Exchanges(ib)%epio_id
00154
00155 ila_info_exch(ib,4) = Drv_Exchanges(ib)%interp_id
00156 ila_info_exch(ib,5) = Drv_Exchanges(ib)%interp_status
00157
00158 ila_info_exch(ib,6) = Drv_Exchanges(ib)%transf_id
00159 ila_info_exch(ib,7) = Drv_Exchanges(ib)%transf_status
00160
00161 ila_info_exch(ib,8) = Drv_Exchanges(ib)%trans_in_field_size
00162 ila_info_exch(ib,9) = Drv_Exchanges(ib)%trans_in_request
00163 ila_info_exch(ib,10) = Drv_Exchanges(ib)%trans_in_status
00164
00165 ila_info_exch(ib,11) = Drv_Exchanges(ib)%trans_in_field_type
00166 ila_info_exch(ib,12) = Drv_Exchanges(ib)%conservation
00167
00168 END DO
00169
00170 DO ib = 1, Number_of_Interps
00171
00172 ila_info_interp(ib,1) = Drv_Interps(ib)%interp_id
00173 ila_info_interp(ib,2) = Drv_Interps(ib)%interp_type
00174 ila_info_interp(ib,3:5) = Drv_Interps(ib)%interp_method
00175 ila_info_interp(ib,6) = Drv_Interps(ib)%nb_neighbors
00176 ila_info_interp(ib,7:9) = Drv_Interps(ib)%arg1
00177 ila_info_interp(ib,10:12) = Drv_Interps(ib)%arg2
00178 ila_info_interp(ib,13:15) = Drv_Interps(ib)%arg3
00179 ila_info_interp(ib,16:18) = Drv_Interps(ib)%arg4
00180 ila_info_interp(ib,19:21) = Drv_Interps(ib)%arg5
00181 ila_info_interp(ib,22:24) = Drv_Interps(ib)%arg6
00182 ila_info_interp(ib,25:27) = Drv_Interps(ib)%arg7
00183
00184 END DO
00185
00186 END IF
00187
00188
00189
00190 CALL MPI_Bcast( ila_info_exch, ila_info_exch_size*Number_of_Exchanges, &
00191 MPI_Integer, PRISM_root, comm_drv_local, id_err )
00192
00193 IF ( id_err /= MPI_SUCCESS ) THEN
00194 ierrp (1) = id_err
00195 id_err = PRISM_Error_MPI
00196
00197 call psmile_error_common ( id_err, 'MPI_Bcast', &
00198 ierrp, 1, __FILE__, __LINE__ )
00199 RETURN
00200 ENDIF
00201
00202 CALL MPI_Bcast( ila_info_interp, 27*Number_of_Interps, MPI_Integer, &
00203 PRISM_root, comm_drv_local, id_err )
00204
00205 IF ( id_err /= MPI_SUCCESS ) THEN
00206 ierrp (1) = id_err
00207 id_err = PRISM_Error_MPI
00208
00209 call psmile_error_common ( id_err, 'MPI_Bcast', &
00210 ierrp, 1, __FILE__, __LINE__ )
00211 RETURN
00212 ENDIF
00213
00214
00215 IF (il_rank .ne. PRISM_root) THEN
00216
00217 DO ib =1, Number_of_Exchanges
00218 Drv_Exchanges(ib)%trans_out_id = ila_info_exch(ib,1)
00219 Drv_Exchanges(ib)%trans_in_id = ila_info_exch(ib,2)
00220 Drv_Exchanges(ib)%epio_id = ila_info_exch(ib,3)
00221
00222 Drv_Exchanges(ib)%interp_id = ila_info_exch(ib,4)
00223 Drv_Exchanges(ib)%interp_status = ila_info_exch(ib,5)
00224
00225 Drv_Exchanges(ib)%transf_id = ila_info_exch(ib,6)
00226 Drv_Exchanges(ib)%transf_status = ila_info_exch(ib,7)
00227
00228 Drv_Exchanges(ib)%trans_in_field_size = ila_info_exch(ib,8)
00229 Drv_Exchanges(ib)%trans_in_request = ila_info_exch(ib,9)
00230 Drv_Exchanges(ib)%trans_in_status = ila_info_exch(ib,10)
00231
00232 Drv_Exchanges(ib)%trans_in_field_type = ila_info_exch(ib,11)
00233
00234 Drv_Exchanges(ib)%conservation = ila_info_exch(ib,12)
00235
00236 Drv_Exchanges(ib)%trans_in_nbr_allocated_fields = 0
00237
00238 Nullify(Drv_Exchanges(ib)%trans_in_field_dble)
00239 Nullify(Drv_Exchanges(ib)%trans_in_field_real)
00240 Nullify(Drv_Exchanges(ib)%trans_in_field_int)
00241
00242 Nullify(Drv_Exchanges(ib)%global_sum_dble)
00243 Nullify(Drv_Exchanges(ib)%global_sum_int)
00244
00245
00246 END DO
00247
00248 DO ib = 1, Number_of_Interps
00249
00250 Drv_Interps(ib)%interp_id = ila_info_interp(ib,1)
00251 Drv_Interps(ib)%interp_type = ila_info_interp(ib,2)
00252 Drv_Interps(ib)%interp_method = ila_info_interp(ib,3:5)
00253 Drv_Interps(ib)%nb_neighbors = ila_info_interp(ib,6)
00254 Drv_Interps(ib)%arg1 = ila_info_interp(ib,7:9)
00255 Drv_Interps(ib)%arg2 = ila_info_interp(ib,10:12)
00256 Drv_Interps(ib)%arg3 = ila_info_interp(ib,13:15)
00257 Drv_Interps(ib)%arg4 = ila_info_interp(ib,16:18)
00258 Drv_Interps(ib)%arg5 = ila_info_interp(ib,19:21)
00259 Drv_Interps(ib)%arg6 = ila_info_interp(ib,22:24)
00260 Drv_Interps(ib)%arg7 = ila_info_interp(ib,25:27)
00261
00262 END DO
00263
00264 END IF
00265
00266
00267
00268
00269 ALLOCATE (dla_buffer(Number_of_Interps), stat = id_err)
00270 IF (id_err > 0) THEN
00271 ierrp (1) = id_err
00272 ierrp (2) = Number_of_Interps
00273 id_err = 13
00274
00275 call psmile_error_common ( id_err, 'Drv_Exchanges', &
00276 ierrp, 2, __FILE__, __LINE__ )
00277 RETURN
00278 ENDIF
00279
00280 IF (il_rank == PRISM_root) THEN
00281
00282 DO ib = 1, Number_of_Interps
00283 dla_buffer(ib) = Drv_Interps(ib)%arg8
00284 ENDDO
00285
00286 ENDIF
00287
00288 CALL MPI_Bcast( dla_buffer, Number_of_Interps, MPI_DOUBLE_PRECISION, &
00289 PRISM_root, comm_drv_local, id_err )
00290
00291 IF ( id_err /= MPI_SUCCESS ) THEN
00292 ierrp (1) = id_err
00293 id_err = PRISM_Error_MPI
00294
00295 call psmile_error_common ( id_err, 'MPI_Bcast', &
00296 ierrp, 1, __FILE__, __LINE__ )
00297 RETURN
00298 ENDIF
00299
00300 IF (il_rank /= PRISM_root) THEN
00301
00302 DO ib = 1, Number_of_Interps
00303 Drv_Interps(ib)%arg8 = dla_buffer(ib)
00304 ENDDO
00305
00306 ENDIF
00307
00308 END IF
00309
00310 #ifdef VERBOSE
00311 PRINT *, '| | Quit PRISMDrv_bcast2trs'
00312 call psmile_flushstd
00313 #endif
00314
00315 END SUBROUTINE PRISMTrs_bcast2trs
00316
00317
00318
00319
00320
00321
00322