00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine prismtrs_loop(id_err)
00012
00013
00014
00015
00016 USE PRISMDrv, dummy_interface => PRISMTrs_Loop
00017
00018 IMPLICIT NONE
00019
00020
00021
00022
00023
00024
00025
00026 INTEGER, INTENT (Out) :: id_err
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 CHARACTER(LEN=len_cvs_string), SAVE :: mycvs =
00047 '$Id: prismtrs_loop.F90 2937 2011-02-03 10:18:38Z redler $'
00048
00049
00050 INTEGER :: il_stop, il_stop_trans, ib
00051 LOGICAL :: ll_loop
00052
00053 INTEGER, PARAMETER :: nerrp=3
00054 INTEGER :: ierrp(nerrp)
00055
00056 INTEGER, DIMENSION(PSMILe_trans_Header_length) :: ila_loop
00057
00058 INTEGER :: il_status(MPI_STATUS_SIZE)
00059 INTEGER :: il_index, il_rank
00060
00061
00062
00063 #ifdef VERBOSE
00064 PRINT *, '| Enter PRISMTrs_Loop'
00065 call psmile_flushstd
00066 #endif
00067
00068
00069
00070
00071 il_stop = 0
00072 il_stop_trans = 0
00073 ila_loop = 280177
00074 ll_loop = .true.
00075 DO WHILE (ll_loop)
00076
00077 #ifdef VERBOSE
00078 PRINT *, '| | '
00079 PRINT *, '| | Trs ready to receive '
00080 call psmile_flushstd
00081 #endif
00082
00083
00084 CALL MPI_Recv (ila_loop, PSMILe_trans_Header_length, MPI_Integer, &
00085 MPI_ANY_SOURCE, MPI_ANY_TAG, comm_drv_trans, &
00086 il_status, id_err)
00087
00088 IF ( id_err /= MPI_SUCCESS ) THEN
00089 ierrp (1) = id_err
00090 ierrp (2) = PRISM_root
00091 ierrp (3) = PSMILe_Init_tag
00092 id_err = PRISM_Error_Recv
00093 call psmile_error_common ( id_err, 'MPI_Recv', &
00094 ierrp, 3, __FILE__, __LINE__ )
00095 RETURN
00096 ENDIF
00097
00098
00099 #ifdef VERBOSE
00100 PRINT *, '| | Trs receives : ',ila_loop(1), &
00101 ' from global rank ', ila_loop(2)
00102 call psmile_flushstd
00103 #endif
00104
00105 SELECT CASE (ila_loop(1))
00106
00107
00108
00109 CASE (PSMILe_trans_Abort)
00110 PRINT *, '| | Trs received abort message from global rank', ila_loop(2)
00111 PRINT *, '| | Application sequence number is ', ila_loop(4)
00112 PRINT *, '| | Rank within application is ', ila_loop(5)
00113
00114 IF ( ila_loop(3) == 1 ) THEN
00115 PRINT *, '| | invoked by user via call to PRISM_Abort'
00116 ELSE
00117 PRINT *, '| | invoked by PSMILe via call to PSMILe_Abort.'
00118 ENDIF
00119
00120 CALL MPI_ABORT (MPI_COMM_WORLD, MPI_ERR_UNKNOWN, id_err)
00121
00122
00123
00124 CASE (PSMILe_trans_Set_rank_trans)
00125
00126 call prismtrs_get_trans_rank(il_rank, id_err)
00127
00128 CALL MPI_Send(il_rank, 1, MPI_Integer, ila_loop(2), 5, &
00129 comm_drv_trans, id_err)
00130
00131 CASE (PSMILe_trans_Set_epio_trans)
00132
00133 call prismtrs_get_epio_handle(il_index, id_err)
00134
00135 CALL MPI_Send(il_index, 1, MPI_Integer, ila_loop(2), 6, &
00136 comm_drv_trans, id_err)
00137
00138
00139 CASE (PSMILe_trans_Set_src_epio_info)
00140
00141 IF (ila_loop(9) .eq. PSMILe_3D) THEN
00142 #ifdef VERBOSE
00143 PRINT *, &
00144 '| | Trs updates 3d src epio ', ila_loop(6)
00145 call psmile_flushstd
00146 #endif
00147 ELSE IF (ila_loop(9) .eq. PSMILe_2D1D) THEN
00148 #ifdef VERBOSE
00149 PRINT *, &
00150 '| | Trs updates 2d1d src epio ', &
00151 '. Size latlon ', ila_loop(7), ' and size z ', ila_loop(8)
00152 call psmile_flushstd
00153 #endif
00154 ELSE
00155 #ifdef VERBOSE
00156 PRINT *, &
00157 '| | Trs updates 1d1d1d src epio ', &
00158 '. Size ', ila_loop(7), ' x ', ila_loop(8), ' x ', ila_loop(9)
00159 call psmile_flushstd
00160 #endif
00161 END IF
00162
00163 IF (ila_loop(3) .eq. PRISM_Double_Precision) THEN
00164 call prismtrs_set_src_epio_dble(ila_loop, id_err)
00165
00166 ELSE IF (ila_loop(3) .eq. PRISM_Real) THEN
00167 call prismtrs_set_src_epio_real(ila_loop, id_err)
00168 END IF
00169
00170
00171 CASE (PSMILe_trans_Set_tgt_epio_info)
00172
00173 #ifdef VERBOSE
00174 IF (ila_loop(9) .eq. PSMILe_3D) THEN
00175 PRINT *, &
00176 '| | Trs updates 3d tgt epio ', ila_loop(6)
00177 call psmile_flushstd
00178 ELSE IF (ila_loop(9) .eq. PSMILe_2D1D) THEN
00179 PRINT *, &
00180 '| | Trs updates 2d1d tgt epio ', &
00181 '. Size latlon ', ila_loop(7), ' and size z ', ila_loop(8)
00182 call psmile_flushstd
00183 ELSE
00184 PRINT *, &
00185 '| | Trs updates 1d1d1d tgt epio ', &
00186 '. Size ', ila_loop(7), ' x ', ila_loop(8), ' x ', ila_loop(9)
00187 call psmile_flushstd
00188 END IF
00189 #endif
00190
00191 IF (ila_loop(3) .eq. PRISM_Double_Precision) THEN
00192 call prismtrs_set_tgt_epio_dble(ila_loop, id_err)
00193
00194 ELSE IF (ila_loop(3) .eq. PRISM_Real) THEN
00195 call prismtrs_set_tgt_epio_real(ila_loop, id_err)
00196 END IF
00197
00198
00199 CASE (PSMILe_trans_Set_triple_links)
00200
00201 call prismtrs_set_triple_links(ila_loop(3), &
00202 ila_loop(4), &
00203 ila_loop(5), &
00204 id_err)
00205
00206
00207 CASE (PSMILe_trans_Set_neighbors_info)
00208
00209 IF (ila_loop(6) .eq. PSMILe_3D) THEN
00210
00211 #ifdef VERBOSE
00212 PRINT *, &
00213 '| | Trs updates 3d neighbors info for epio ', ila_loop(3)
00214 call psmile_flushstd
00215 #endif
00216 call prismtrs_set_neighbors3d(ila_loop(2), &
00217 ila_loop(3), &
00218 ila_loop(4), &
00219 ila_loop(5), &
00220 ila_loop(7), &
00221 ila_loop(11), &
00222 id_err)
00223 ELSE
00224
00225 #ifdef VERBOSE
00226 PRINT *, &
00227 '| | Trs does not know how to update such neighbors info for epio ', ila_loop(3)
00228 call psmile_flushstd
00229 #endif
00230
00231 END IF
00232
00233
00234 CASE (PSMILe_trans_Put)
00235
00236 #ifdef VERBOSE
00237 PRINT *, '| | Trs receives the transient_out ', &
00238 ila_loop(3),' and epio ', ila_loop(4)
00239 call psmile_flushstd
00240 #endif
00241 IF (ila_loop(5) .eq. PRISM_Integer) THEN
00242 call prismtrs_mind_int(ila_loop(2), &
00243 ila_loop(3), &
00244 ila_loop(4), &
00245 ila_loop(6), &
00246 ila_loop(7), &
00247 id_err)
00248 ELSE IF (ila_loop(5) .eq. PRISM_Real) THEN
00249 call prismtrs_mind_real(ila_loop(2), &
00250 ila_loop(3), &
00251 ila_loop(4), &
00252 ila_loop(6), &
00253 ila_loop(7), &
00254 id_err)
00255 ELSE IF (ila_loop(5) .eq. PRISM_Double_Precision) THEN
00256 call prismtrs_mind_dble(ila_loop(2), &
00257 ila_loop(3), &
00258 ila_loop(4), &
00259 ila_loop(6), &
00260 ila_loop(7), &
00261 id_err)
00262 END IF
00263
00264
00265 CASE (PSMILe_trans_Get)
00266
00267 #ifdef VERBOSE
00268 PRINT *, '| | Trs is aked to send the transient_in ', &
00269 ila_loop(3),' and epio ', ila_loop(4)
00270 call psmile_flushstd
00271 #endif
00272
00273 IF (ila_loop(5) .eq. PRISM_Integer) THEN
00274 call prismtrs_target_int(ila_loop(2), &
00275 ila_loop(3), &
00276 ila_loop(4), &
00277 ila_loop(6), &
00278 id_err)
00279 ELSE IF (ila_loop(5) .eq. PRISM_Real) THEN
00280 call prismtrs_target_real(ila_loop(2), &
00281 ila_loop(3), &
00282 ila_loop(4), &
00283 ila_loop(6), &
00284 id_err)
00285 ELSE IF (ila_loop(5) .eq. PRISM_Double_Precision) THEN
00286 call prismtrs_target_dble(ila_loop(2), &
00287 ila_loop(3), &
00288 ila_loop(4), &
00289 ila_loop(6), &
00290 id_err)
00291 END IF
00292
00293
00294 CASE (PSMILe_trans_Finalize)
00295
00296 IF (appl%rank .eq. PRISM_Root) THEN
00297 #ifdef VERBOSE
00298 PRINT *, '| | Trs is asked to finalize by process ', &
00299 ila_loop(2)
00300 #endif
00301 il_stop = il_stop + 1
00302 #ifdef VERBOSE
00303 PRINT *, '| | ', il_stop, ' pes on ',ig_nb_tot_pes, &
00304 ' have already asked to finalize'
00305 call psmile_flushstd
00306 #endif
00307
00308 IF (il_stop .eq. ig_nb_tot_pes) THEN
00309
00310
00311
00312 IF (ig_driver_nb_pes .gt. 1) THEN
00313 IF (il_stop_trans .eq. 0) THEN
00314 ila_loop(1) = PSMILe_trans_End_trans
00315 DO ib = 1, ig_driver_nb_pes-1
00316 CALL MPI_Send(ila_loop(1), &
00317 PSMILe_trans_Header_length, &
00318 MPI_Integer, ib, 0, comm_drv_trans, id_err)
00319 END DO
00320 END IF
00321 END IF
00322
00323 ll_loop = .false.
00324 END IF
00325 ELSE
00326 ila_loop(1) = PSMILe_trans_Finalize_trans
00327 CALL MPI_Send(ila_loop(1), PSMILe_trans_Header_length, &
00328 MPI_Integer, PRISM_root, 0, comm_drv_trans, id_err)
00329 END IF
00330
00331
00332 CASE (PSMILe_trans_End_trans)
00333
00334 ll_loop = .false.
00335
00336 CASE DEFAULT
00337
00338 PRINT *, &
00339 '| | Trs is not inteligent enough to understand'
00340 call psmile_flushstd
00341
00342 END SELECT
00343
00344 END DO
00345
00346
00347 #ifdef VERBOSE
00348 PRINT *, '| Quit PRISMTrs_Loop'
00349 PRINT *, '|'
00350 call psmile_flushstd
00351 #endif
00352
00353 END SUBROUTINE PRISMTrs_Loop
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370