Changeset 764 for LMDZ4/trunk/libf/dyn3dpar/parallel.F90
- Timestamp:
- Jun 4, 2007, 4:13:10 PM (18 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ4/trunk/libf/dyn3dpar/parallel.F90
r630 r764 3 3 integer, save :: mpi_size 4 4 integer, save :: mpi_rank 5 integer, save :: COMM_LMDZ 5 6 integer, save :: jj_begin 6 7 integer, save :: jj_end … … 14 15 integer, allocatable, save, dimension(:) :: jj_end_para 15 16 integer, allocatable, save, dimension(:) :: jj_nb_para 17 integer, save :: OMP_CHUNK 16 18 17 19 contains … … 19 21 subroutine init_parallel 20 22 USE vampir 23 #ifdef CPP_COUPLE 24 #ifdef CPP_PSMILE 25 USE mod_prism_proto 26 #endif 27 #endif 21 28 implicit none 22 29 … … 25 32 integer :: type_size 26 33 integer, dimension(3) :: blocklen,type 34 integer :: comp_id 27 35 28 36 … … 30 38 #include "dimensions90.h" 31 39 #include "paramet90.h" 32 40 41 #ifdef CPP_COUPLE 42 #ifdef CPP_PSMILE 43 call prism_init_comp_proto (comp_id, 'lmdz.x', ierr) 44 call prism_get_localcomm_proto(COMM_LMDZ,ierr) 45 #endif 46 #else 33 47 call MPI_INIT(ierr) 48 COMM_LMDZ=MPI_COMM_WORLD 49 #endif 34 50 call InitVampir 35 call MPI_COMM_SIZE( MPI_COMM_WORLD,mpi_size,ierr)36 call MPI_COMM_RANK( MPI_COMM_WORLD,mpi_rank,ierr)51 call MPI_COMM_SIZE(COMM_LMDZ,mpi_size,ierr) 52 call MPI_COMM_RANK(COMM_LMDZ,mpi_rank,ierr) 37 53 38 54 … … 50 66 print *," ---> diminuez le nombre de CPU ou augmentez la taille en lattitude" 51 67 52 call MPI_ABORT( MPI_COMM_WORLD,-1, ierr)68 call MPI_ABORT(COMM_LMDZ,-1, ierr) 53 69 54 70 endif … … 133 149 134 150 subroutine Finalize_parallel 151 #ifdef CPP_COUPLE 152 #ifdef CPP_PSMILE 153 use mod_prism_proto 154 #endif 155 #endif 135 156 implicit none 136 157 … … 144 165 deallocate(jj_end_para) 145 166 deallocate(jj_nb_para) 146 167 168 #ifdef CPP_COUPLE 169 #ifdef CPP_PSMILE 170 call prism_terminate_proto(ierr) 171 IF (ierr .ne. PRISM_Ok) THEN 172 call abort_gcm('Finalize_parallel',' Probleme dans prism_terminate_proto ',1) 173 endif 174 #endif 175 #else 147 176 call MPI_FINALIZE(ierr) 177 #endif 148 178 149 179 end subroutine Finalize_parallel 150 180 151 181 subroutine Pack_Data(Field,ij,ll,row,Buffer) 152 182 implicit none … … 217 247 INTEGER :: Buffer_size 218 248 219 call MPI_Barrier( MPI_COMM_WORLD,ierr)249 call MPI_Barrier(COMM_LMDZ,ierr) 220 250 call VTb(VThallo) 221 251 … … 253 283 call PACK_Data(Field(ij_begin,1),ij,ll,down,Buffer_Send_up) 254 284 call MPI_ISSEND(Buffer_send_up,Buffer_Size,MPI_REAL8,MPI_Rank-1,1, & 255 MPI_COMM_WORLD,Request(NbRequest),ierr)285 COMM_LMDZ,Request(NbRequest),ierr) 256 286 ENDIF 257 287 … … 264 294 265 295 call MPI_ISSEND(Buffer_send_down,Buffer_Size,MPI_REAL8,MPI_Rank+1,1, & 266 MPI_COMM_WORLD,Request(NbRequest),ierr)296 COMM_LMDZ,Request(NbRequest),ierr) 267 297 ENDIF 268 298 … … 274 304 275 305 call MPI_IRECV(Buffer_recv_up,Buffer_size,MPI_REAL8,MPI_Rank-1,1, & 276 MPI_COMM_WORLD,Request(NbRequest),ierr)306 COMM_LMDZ,Request(NbRequest),ierr) 277 307 278 308 … … 285 315 286 316 call MPI_IRECV(Buffer_recv_down,Buffer_size,MPI_REAL8,MPI_Rank+1,1, & 287 MPI_COMM_WORLD,Request(NbRequest),ierr)317 COMM_LMDZ,Request(NbRequest),ierr) 288 318 289 319 … … 295 325 296 326 call VTe(VThallo) 297 call MPI_Barrier( MPI_COMM_WORLD,ierr)327 call MPI_Barrier(COMM_LMDZ,ierr) 298 328 RETURN 299 329 … … 349 379 350 380 call MPI_GATHERV(Buffer_send,(min(ij_end,ij)-ij_begin+1)*ll,MPI_REAL8, & 351 Buffer_Recv,Recv_count,displ,MPI_REAL8,rank, MPI_COMM_WORLD,ierr)381 Buffer_Recv,Recv_count,displ,MPI_REAL8,rank,COMM_LMDZ,ierr) 352 382 353 383 if (MPI_Rank==rank) then … … 380 410 381 411 call Gather_Field(Field,ij,ll,0) 382 call MPI_BCAST(Field,ij*ll,MPI_REAL8,0, MPI_COMM_WORLD)412 call MPI_BCAST(Field,ij*ll,MPI_REAL8,0,COMM_LMDZ,ierr) 383 413 384 414 end subroutine AllGather_Field … … 395 425 INTEGER :: ierr 396 426 397 call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank, MPI_COMM_WORLD)427 call MPI_BCAST(Field,ij*ll,MPI_REAL8,rank,COMM_LMDZ,ierr) 398 428 399 429 end subroutine Broadcast_Field
Note: See TracChangeset
for help on using the changeset viewer.