Changeset 854 for LMDZ4/trunk


Ignore:
Timestamp:
Oct 23, 2007, 2:03:16 PM (17 years ago)
Author:
lsce
Message:

AC + JG + YM + ACo : - Corrected bugs found on 16 and 24 CPUs on ccrt scalar machine platine

  • Added use of lmdz communicator in inca model
Location:
LMDZ4/trunk/libf
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • LMDZ4/trunk/libf/dyn3dpar/gcm.F

    r807 r854  
    230230#ifdef INCA
    231231      call init_const_lmdz(nbtrac,anneeref,dayref,iphysiq,day_step,nday)
    232       call init_inca_para(iim,jjm+1,klon_glo,mpi_size,klon_mpi_para_nb)
     232      call init_inca_para(iim,jjm+1,klon_glo,mpi_size,klon_mpi_para_nb,
     233     $ COMM_LMDZ)
    233234#endif
    234235
  • LMDZ4/trunk/libf/dyn3dpar/vlspltgen_p.F

    r774 r854  
    184184       
    185185          call vlx_p(zq(1,1,iq),pente_max,zm(1,1,iq),mu,
    186      &               ij_begin,ij_begin+2*iip1-1)
    187           call vlx_p(zq(1,1,iq),pente_max,zm(1,1,iq),mu,
    188      &               ij_end-2*iip1+1,ij_end)
     186     &               ij_begin,ij_end)
    189187
    190188c$OMP MASTER
     
    197195     
    198196          call vlxqs_p(zq(1,1,iq),pente_max,zm(1,1,iq),mu,qsat,
    199      &                 ij_begin,ij_begin+2*iip1-1)
    200           call vlxqs_p(zq(1,1,iq),pente_max,zm(1,1,iq),mu,qsat,
    201      &                 ij_end-2*iip1+1,ij_end)
     197     &                 ij_begin,ij_end)
    202198c$OMP MASTER
    203199          call VTb(VTHallo)
     
    219215      call VTb(VTHallo)
    220216      call SendRequest(MyRequest1)
    221       call VTe(VTHallo)
    222 c$OMP END MASTER       
    223 c$OMP BARRIER
    224       do iq=1,nqmx
    225 
    226         if(iadv(iq) == 0) then
    227        
    228           cycle
    229        
    230         else if (iadv(iq)==10) then
    231 
    232           call vlx_p(zq(1,1,iq),pente_max,zm(1,1,iq),mu,
    233      &               ij_begin+2*iip1,ij_end-2*iip1)
    234        
    235         else if (iadv(iq)==14) then
    236 
    237           call vlxqs_p(zq(1,1,iq),pente_max,zm(1,1,iq),mu,qsat,
    238      &                 ij_begin+2*iip1,ij_end-2*iip1)
    239    
    240         else
    241        
    242           stop 'vlspltgen_p : schema non parallelise'
    243      
    244         endif
    245      
    246       enddo
    247 c$OMP BARRIER     
    248 c$OMP MASTER
    249       call VTb(VTHallo)
    250217      call WaitRecvRequest(MyRequest1)
    251218      call WaitSendRequest(MyRequest1)
     
    286253c$OMP BARRIER       
    287254          call vlz_p(zq(1,1,iq),pente_max,zm(1,1,iq),mw,
    288      &               ij_begin,ij_begin+2*iip1-1)
    289           call vlz_p(zq(1,1,iq),pente_max,zm(1,1,iq),mw,
    290      &               ij_end-2*iip1+1,ij_end)
     255     &               ij_begin,ij_end)
    291256c$OMP BARRIER
    292257
     
    309274      call VTb(VTHallo)
    310275      call SendRequest(MyRequest2)
    311       call VTe(VTHallo)
    312 c$OMP END MASTER       
    313 c$OMP BARRIER
    314       do iq=1,nqmx
    315 
    316         if(iadv(iq) == 0) then
    317          
    318           cycle
    319        
    320         else if (iadv(iq)==10 .or. iadv(iq)==14 ) then
    321 c$OMP BARRIER       
    322           call vlz_p(zq(1,1,iq),pente_max,zm(1,1,iq),mw,
    323      &               ij_begin+2*iip1,ij_end-2*iip1)
    324 c$OMP BARRIER       
    325         else
    326        
    327           stop 'vlspltgen_p : schema non parallelise'
    328      
    329         endif
    330      
    331       enddo
    332 
    333 c$OMP BARRIER
    334 c$OMP MASTER
    335       call VTb(VTHallo)
    336276      call WaitRecvRequest(MyRequest2)
    337277      call WaitSendRequest(MyRequest2)
  • LMDZ4/trunk/libf/phylmd/surf_land_orchidee_mod.F90

    r818 r854  
    485485    INTEGER                               :: MyLastPoint
    486486    INTEGER                               :: LastPoint
    487     INTEGER                               :: mpi_rank
    488     INTEGER                               :: mpi_size
    489     INTEGER                               :: ierr  
     487    INTEGER                               :: mpi_rank_orch
     488    INTEGER                               :: mpi_size_orch
     489    INTEGER                               :: ierr
    490490!
    491491! End definition
     
    496496    IF (is_parallel) THEN
    497497#ifdef CPP_PARA   
    498        CALL MPI_COMM_SIZE(orch_comm,mpi_size,ierr)
    499        CALL MPI_COMM_RANK(orch_comm,mpi_rank,ierr)
     498       CALL MPI_COMM_SIZE(orch_comm,mpi_size_orch,ierr)
     499       CALL MPI_COMM_RANK(orch_comm,mpi_rank_orch,ierr)
    500500#endif
    501501    ELSE
    502        mpi_rank=0
    503        mpi_size=1
    504     ENDIF
    505    
     502       mpi_rank_orch=0
     503       mpi_size_orch=1
     504    ENDIF
     505
    506506    IF (is_parallel) THEN
    507        IF (.NOT. is_north_pole) THEN
     507       IF (mpi_rank_orch /= 0) THEN
    508508#ifdef CPP_PARA
    509           CALL MPI_RECV(LastPoint,1,MPI_INTEGER,mpi_rank-1,1234,orch_comm,status,ierr)
    510 #endif
    511        ENDIF
    512        
    513        IF (.NOT. is_south_pole) THEN
     509          CALL MPI_RECV(LastPoint,1,MPI_INTEGER,mpi_rank_orch-1,1234,orch_comm,status,ierr)
     510#endif
     511       ENDIF
     512       
     513       IF (mpi_rank_orch /= mpi_size_orch-1) THEN
    514514#ifdef CPP_PARA
    515           CALL MPI_SEND(MyLastPoint,1,MPI_INTEGER,mpi_rank+1,1234,orch_comm,ierr) 
    516 #endif
    517        ENDIF
    518     ENDIF
    519    
    520     IF (is_north_pole) THEN
     515          CALL MPI_SEND(MyLastPoint,1,MPI_INTEGER,mpi_rank_orch+1,1234,orch_comm,ierr) 
     516#endif
     517       ENDIF
     518    ENDIF
     519   
     520    IF (mpi_rank_orch == 0) THEN
    521521       offset=0
    522522    ELSE
  • LMDZ4/trunk/libf/phytherm/surf_land_orchidee_mod.F90

    r837 r854  
    485485    INTEGER                               :: MyLastPoint
    486486    INTEGER                               :: LastPoint
    487     INTEGER                               :: mpi_rank
    488     INTEGER                               :: mpi_size
    489     INTEGER                               :: ierr  
     487    INTEGER                               :: mpi_rank_orch
     488    INTEGER                               :: mpi_size_orch
     489    INTEGER                               :: ierr
    490490!
    491491! End definition
     
    496496    IF (is_parallel) THEN
    497497#ifdef CPP_PARA   
    498        CALL MPI_COMM_SIZE(orch_comm,mpi_size,ierr)
    499        CALL MPI_COMM_RANK(orch_comm,mpi_rank,ierr)
     498       CALL MPI_COMM_SIZE(orch_comm,mpi_size_orch,ierr)
     499       CALL MPI_COMM_RANK(orch_comm,mpi_rank_orch,ierr)
    500500#endif
    501501    ELSE
    502        mpi_rank=0
    503        mpi_size=1
    504     ENDIF
    505    
     502       mpi_rank_orch=0
     503       mpi_size_orch=1
     504    ENDIF
     505
    506506    IF (is_parallel) THEN
    507        IF (.NOT. is_north_pole) THEN
     507       IF (mpi_rank_orch /= 0) THEN
    508508#ifdef CPP_PARA
    509           CALL MPI_RECV(LastPoint,1,MPI_INTEGER,mpi_rank-1,1234,orch_comm,status,ierr)
    510 #endif
    511        ENDIF
    512        
    513        IF (.NOT. is_south_pole) THEN
     509          CALL MPI_RECV(LastPoint,1,MPI_INTEGER,mpi_rank_orch-1,1234,orch_comm,status,ierr)
     510#endif
     511       ENDIF
     512       
     513       IF (mpi_rank_orch /= mpi_size_orch-1) THEN
    514514#ifdef CPP_PARA
    515           CALL MPI_SEND(MyLastPoint,1,MPI_INTEGER,mpi_rank+1,1234,orch_comm,ierr) 
    516 #endif
    517        ENDIF
    518     ENDIF
    519    
    520     IF (is_north_pole) THEN
     515          CALL MPI_SEND(MyLastPoint,1,MPI_INTEGER,mpi_rank_orch+1,1234,orch_comm,ierr) 
     516#endif
     517       ENDIF
     518    ENDIF
     519   
     520    IF (mpi_rank_orch == 0) THEN
    521521       offset=0
    522522    ELSE
Note: See TracChangeset for help on using the changeset viewer.