Ignore:
Timestamp:
03/12/13 16:34:45 (11 years ago)
Author:
ymipsl
Message:

Set constant sign for wind way :
ne(ij,right)==ne_right=1
ne(ij,rup)==ne_rup=-1
ne(ij,lup)==ne_lup=1
ne(ij,left)==ne_left=-1
ne(ij,ldown)==ne_ldown=1
ne(ij,rdown)==ne_rdown=-1

Modified transfert function to be compliant for this convention.

YM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • codes/icosagcm/trunk/src/transfert_mpi.f90

    r74 r146  
    44  TYPE array 
    55    INTEGER,POINTER :: value(:) 
     6    INTEGER,POINTER :: sign(:) 
    67    INTEGER         :: domain 
    78    INTEGER         :: rank 
     
    1718    INTEGER :: max_size 
    1819    INTEGER :: size 
     20    LOGICAL :: vector 
    1921    INTEGER,POINTER :: src_domain(:) 
    2022    INTEGER,POINTER :: src_i(:) 
     
    2426    INTEGER,POINTER :: target_i(:) 
    2527    INTEGER,POINTER :: target_j(:) 
     28    INTEGER,POINTER :: target_sign(:) 
    2629    INTEGER :: nrecv 
    2730    TYPE(ARRAY),POINTER :: recv(:) 
     
    3134   
    3235  TYPE(t_request),POINTER :: req_i1(:) 
    33   TYPE(t_request),POINTER :: req_e1(:) 
     36  TYPE(t_request),POINTER :: req_e1_scal(:) 
     37  TYPE(t_request),POINTER :: req_e1_vect(:) 
    3438   
    3539   
     
    8488    CALL finalize_request(req_i1) 
    8589   
    86     CALL create_request(field_u,req_e1) 
     90    CALL create_request(field_u,req_e1_scal) 
    8791    DO ind=1,ndomain 
    8892      CALL swap_dimensions(ind) 
    8993      DO i=ii_begin,ii_end 
    90         CALL request_add_point(ind,i,jj_begin-1,req_e1,rup) 
    91         CALL request_add_point(ind,i+1,jj_begin-1,req_e1,lup) 
     94        CALL request_add_point(ind,i,jj_begin-1,req_e1_scal,rup) 
     95        CALL request_add_point(ind,i+1,jj_begin-1,req_e1_scal,lup) 
    9296      ENDDO 
    9397 
    9498      DO j=jj_begin,jj_end 
    95         CALL request_add_point(ind,ii_end+1,j,req_e1,left) 
    96         CALL request_add_point(ind,ii_end+1,j-1,req_e1,lup) 
     99        CALL request_add_point(ind,ii_end+1,j,req_e1_scal,left) 
     100        CALL request_add_point(ind,ii_end+1,j-1,req_e1_scal,lup) 
    97101      ENDDO     
    98102     
    99103      DO i=ii_begin,ii_end 
    100         CALL request_add_point(ind,i,jj_end+1,req_e1,ldown) 
    101         CALL request_add_point(ind,i-1,jj_end+1,req_e1,rdown) 
     104        CALL request_add_point(ind,i,jj_end+1,req_e1_scal,ldown) 
     105        CALL request_add_point(ind,i-1,jj_end+1,req_e1_scal,rdown) 
    102106      ENDDO     
    103107 
    104108      DO j=jj_begin,jj_end 
    105         CALL request_add_point(ind,ii_begin-1,j,req_e1,right) 
    106         CALL request_add_point(ind,ii_begin-1,j+1,req_e1,rdown) 
     109        CALL request_add_point(ind,ii_begin-1,j,req_e1_scal,right) 
     110        CALL request_add_point(ind,ii_begin-1,j+1,req_e1_scal,rdown) 
    107111      ENDDO    
    108112 
    109113      DO i=ii_begin+1,ii_end-1 
    110         CALL request_add_point(ind,i,jj_begin,req_e1,right) 
    111         CALL request_add_point(ind,i,jj_end,req_e1,right) 
     114        CALL request_add_point(ind,i,jj_begin,req_e1_scal,right) 
     115        CALL request_add_point(ind,i,jj_end,req_e1_scal,right) 
    112116      ENDDO 
    113117     
    114118      DO j=jj_begin+1,jj_end-1 
    115         CALL request_add_point(ind,ii_begin,j,req_e1,rup) 
    116         CALL request_add_point(ind,ii_end,j,req_e1,rup) 
     119        CALL request_add_point(ind,ii_begin,j,req_e1_scal,rup) 
     120        CALL request_add_point(ind,ii_end,j,req_e1_scal,rup) 
    117121      ENDDO    
    118122 
    119       CALL request_add_point(ind,ii_begin+1,jj_begin,req_e1,left) 
    120       CALL request_add_point(ind,ii_begin,jj_begin+1,req_e1,ldown) 
    121       CALL request_add_point(ind,ii_begin+1,jj_end,req_e1,left) 
    122       CALL request_add_point(ind,ii_end,jj_begin+1,req_e1,ldown) 
    123  
    124       CALL finalize_request(req_e1) 
    125       
    126     ENDDO 
    127    
     123      CALL request_add_point(ind,ii_begin+1,jj_begin,req_e1_scal,left) 
     124      CALL request_add_point(ind,ii_begin,jj_begin+1,req_e1_scal,ldown) 
     125      CALL request_add_point(ind,ii_begin+1,jj_end,req_e1_scal,left) 
     126      CALL request_add_point(ind,ii_end,jj_begin+1,req_e1_scal,ldown) 
     127 
     128    ENDDO 
     129 
     130    CALL finalize_request(req_e1_scal) 
     131     
     132    CALL create_request(field_u,req_e1_vect,.TRUE.) 
     133    DO ind=1,ndomain 
     134      CALL swap_dimensions(ind) 
     135      DO i=ii_begin,ii_end 
     136        CALL request_add_point(ind,i,jj_begin-1,req_e1_vect,rup) 
     137        CALL request_add_point(ind,i+1,jj_begin-1,req_e1_vect,lup) 
     138      ENDDO 
     139 
     140      DO j=jj_begin,jj_end 
     141        CALL request_add_point(ind,ii_end+1,j,req_e1_vect,left) 
     142        CALL request_add_point(ind,ii_end+1,j-1,req_e1_vect,lup) 
     143      ENDDO     
     144     
     145      DO i=ii_begin,ii_end 
     146        CALL request_add_point(ind,i,jj_end+1,req_e1_vect,ldown) 
     147        CALL request_add_point(ind,i-1,jj_end+1,req_e1_vect,rdown) 
     148      ENDDO     
     149 
     150      DO j=jj_begin,jj_end 
     151        CALL request_add_point(ind,ii_begin-1,j,req_e1_vect,right) 
     152        CALL request_add_point(ind,ii_begin-1,j+1,req_e1_vect,rdown) 
     153      ENDDO    
     154 
     155      DO i=ii_begin+1,ii_end-1 
     156        CALL request_add_point(ind,i,jj_begin,req_e1_vect,right) 
     157        CALL request_add_point(ind,i,jj_end,req_e1_vect,right) 
     158      ENDDO 
     159     
     160      DO j=jj_begin+1,jj_end-1 
     161        CALL request_add_point(ind,ii_begin,j,req_e1_vect,rup) 
     162        CALL request_add_point(ind,ii_end,j,req_e1_vect,rup) 
     163      ENDDO    
     164 
     165      CALL request_add_point(ind,ii_begin+1,jj_begin,req_e1_vect,left) 
     166      CALL request_add_point(ind,ii_begin,jj_begin+1,req_e1_vect,ldown) 
     167      CALL request_add_point(ind,ii_begin+1,jj_end,req_e1_vect,left) 
     168      CALL request_add_point(ind,ii_end,jj_begin+1,req_e1_vect,ldown) 
     169 
     170   
     171    ENDDO   
     172 
     173    CALL finalize_request(req_e1_vect) 
     174 
    128175  END SUBROUTINE init_transfert 
    129176   
    130   SUBROUTINE create_request(type_field,request) 
     177  SUBROUTINE create_request(type_field,request,vector) 
    131178  USE domain_mod 
    132179  USE field_mod 
     
    134181    INTEGER :: type_field 
    135182    TYPE(t_request),POINTER :: request(:) 
     183    LOGICAL,OPTIONAL :: vector 
     184     
    136185    TYPE(t_request),POINTER :: req 
    137186    TYPE(t_domain),POINTER :: d 
    138187    INTEGER :: ind 
    139188    INTEGER :: max_size 
    140          
     189        
    141190    ALLOCATE(request(ndomain)) 
    142191 
     
    155204      req%max_size=max_size*2 
    156205      req%size=0 
     206      req%vector=.FALSE. 
     207      IF (PRESENT(vector)) req%vector=vector 
    157208      ALLOCATE(req%src_domain(req%max_size)) 
    158209      ALLOCATE(req%src_ind(req%max_size)) 
     
    162213      ALLOCATE(req%target_i(req%max_size)) 
    163214      ALLOCATE(req%target_j(req%max_size)) 
     215      ALLOCATE(req%target_sign(req%max_size)) 
    164216    ENDDO 
    165217   
     
    177229    INTEGER,POINTER :: target_i(:) 
    178230    INTEGER,POINTER :: target_j(:) 
     231    INTEGER,POINTER :: target_sign(:) 
    179232 
    180233    PRINT *,"REALLOCATE_REQUEST" 
     
    186239    target_i=>req%target_i 
    187240    target_j=>req%target_j 
     241    target_sign=>req%target_sign 
    188242!    req%max_size=req%max_size*2 
    189243    ALLOCATE(req%src_domain(req%max_size*2)) 
     
    194248    ALLOCATE(req%target_i(req%max_size*2)) 
    195249    ALLOCATE(req%target_j(req%max_size*2)) 
     250    ALLOCATE(req%target_sign(req%max_size*2)) 
    196251     
    197252    req%src_domain(1:req%max_size)=src_domain(:) 
     
    202257    req%target_i(1:req%max_size)=target_i(:) 
    203258    req%target_j(1:req%max_size)=target_j(:) 
     259    req%target_sign(1:req%max_size)=target_sign(:) 
    204260     
    205261    req%max_size=req%max_size*2 
     
    212268    DEALLOCATE(target_i) 
    213269    DEALLOCATE(target_j) 
     270    DEALLOCATE(target_sign) 
    214271 
    215272  END SUBROUTINE reallocate_request 
     
    243300 
    244301        req%target_ind(req%size)=(j-1)*d%iim+i 
     302        req%target_sign(req%size)=1 
    245303        req%src_domain(req%size)=src_domain 
    246304        req%src_ind(req%size)=(src_j-1)*src_iim+src_i 
     
    254312        src_n=(src_j-1)*src_iim+src_i 
    255313        src_delta=domain(ind)%delta(i,j) 
    256          
    257 !        src_pos=MOD(pos-1+src_delta+6,6)+1 
    258314        src_pos=domain(ind)%edge_assign_pos(pos-1,i,j)+1 
    259315                 
    260316        req%target_ind(req%size)=(j-1)*d%iim+i+d%u_pos(pos) 
     317 
     318        req%target_sign(req%size)= 1 
     319        IF (req%vector) req%target_sign(req%size)= domain(ind)%edge_assign_sign(pos-1,i,j) 
     320 
    261321        req%src_domain(req%size)=src_domain 
    262322        req%src_ind(req%size)=src_n+domain_glo(src_domain)%u_pos(src_pos) 
    263  
    264 !        req%target_i(req%size)=i 
    265 !        req%target_j(req%size)=j 
    266 !        req%src_i(req%size)=domain(ind)%assign_i(i,j) 
    267 !        req%src_j(req%size)=domain(ind)%assign_j(i,j) 
    268          
    269 !        PRINT *,"1--->",ind,i,j,pos 
    270 !        PRINT *,"2--->",src_domain,src_i,src_j,src_pos 
    271323 
    272324      ELSE IF (req%type_field==field_z) THEN 
     
    283335         
    284336        req%target_ind(req%size)=(j-1)*d%iim+i+d%z_pos(pos) 
     337        req%target_sign(req%size)=1 
    285338        req%src_domain(req%size)=src_domain 
    286339        req%src_ind(req%size)=src_n+domain_glo(src_domain)%z_pos(src_pos) 
     
    338391          req%recv(irecv)%domain=domglo_loc_ind(ind_glo) 
    339392          ALLOCATE(req%recv(irecv)%value(req%recv(irecv)%size)) 
     393          ALLOCATE(req%recv(irecv)%sign(req%recv(irecv)%size)) 
    340394          ALLOCATE(req%recv(irecv)%buffer(req%recv(irecv)%size)) 
    341395        ENDIF 
     
    350404        req%recv(irecv)%value(size)=req%src_ind(i) 
    351405        req%recv(irecv)%buffer(size)=req%target_ind(i) 
     406        req%recv(irecv)%sign(size)=req%target_sign(i) 
    352407      ENDDO 
    353408    ENDDO 
     
    362417      ENDDO 
    363418    ENDDO 
    364      
    365419     
    366420    CALL MPI_Alltoall(nb_domain_recv,1,MPI_INTEGER,nb_domain_send,1,MPI_INTEGER,comm_icosa,ierr)      
     
    465519     DO irecv=1,req%nrecv 
    466520       req%recv(irecv)%value(:)=req%recv(irecv)%buffer(:) 
     521       req%recv(irecv)%sign(:) =req%recv(irecv)%sign(:) 
    467522       DEALLOCATE(req%recv(irecv)%buffer) 
    468523     ENDDO 
     
    470525    
    471526     
    472 !   nb_domain_recv(:)=0 
    473 !   nb_data_domain_recv(:)=0 
    474 !    
    475 !   DO ind_loc=1,ndomain 
    476 !    
    477 !     DO i=1,req%size 
    478 !       ind_glo=req%src_domain(i) 
    479 !       nb_data_domain_recv(ind_glo)=nb_data_domain_recv(ind_glo)+1 
    480 !     ENDDO 
    481 !    
    482 !     DO ind_glo=1,ndomain_glo 
    483 !       IF ( nb_data_domain_recv(ind_glo) > 0 )  nb_domain_recv(domglo_rank(ind_glo))=nb_domain_recv(domglo_rank(ind_glo))+1 
    484 !     ENDDO 
    485 !        
    486 !     CALL MPI_Alltoall(nb_domain_recv,1,MPI_INTEGER,nb_domain_send,1,MPI_INTEGER,comm_icosa,ierr) 
    487 !   ENDDO 
    488 !    
    489 !   DO  
    490 !   recv=sum(nb_domain_recv(:)) 
    491 !   send=sum(nb_domain_send(:)) 
    492  
    493 !   ALLOCATE(req%recv(recv)) 
    494 !   ALLOCATE(req%send(send)) 
    495  
    496 !   ALLOCATE(mpi_req(2*(send+recv))) 
    497 !   ALLOCATE(status(MPI_STATUS_SIZE,2*(send+recv))) 
    498 !    
    499 !   recv=0 
    500 !   ireq=0 
    501 !   DO ind_glo=1,ndomain_glo 
    502 !     IF (nb_data_domain_recv(ind_glo)>0) THEN 
    503 !       recv=recv+1 
    504 !       list_domain_recv(ind_glo)=recv 
    505 !       req%recv(recv)%rank=domglo_rank(ind_glo) 
    506 !       req%recv(recv)%size=nb_data_domain_recv(ind_glo) 
    507 !       req%recv(recv)%domain=domglo_loc_ind(ind_glo) 
    508 !       ALLOCATE(req%recv(recv)%value(req%recv(recv)%size)) 
    509 !       ireq=ireq+1 
    510 !       CALL MPI_ISEND(req%recv(recv)%domain,1,MPI_INTEGER,req%recv(recv)%rank,0,comm_icosa, mpi_req(ireq),ierr) 
    511 !       ireq=ireq+1 
    512 !       CALL MPI_ISEND(req%recv(recv)%size,1,MPI_INTEGER,req%recv(recv)%rank,0,comm_icosa, mpi_req(ireq),ierr) 
    513 !     ENDIF 
    514 !   ENDDO 
    515 !    
    516 !     
    517 !   send=0 
    518 !   DO rank=0,mpi_size-1 
    519 !     DO j=1,nb_domain_send(rank) 
    520 !       send=send+1 
    521 !       req%send(send)%rank=rank 
    522 !       ireq=ireq+1 
    523 !       CALL MPI_IRECV(req%send(send)%domain,1,MPI_INTEGER,rank,0,comm_icosa, mpi_req(ireq),ierr) 
    524 !       ireq=ireq+1 
    525 !       CALL MPI_IRECV(req%send(send)%size,1,MPI_INTEGER,rank,0,comm_icosa, mpi_req(ireq),ierr) 
    526 !     ENDDO 
    527 !   ENDDO 
    528 !    
    529 !   CALL MPI_WAITALL(2*(send+recv),mpi_req,status,ierr)  
    530  
    531 !   req%recv(:)%size=0 
    532 !    
    533 !   DO i=1,req%size 
    534 !     j=list_domain_recv(req%src_domain(i)) 
    535 !     req%recv(j)%size=req%recv(j)%size+1 
    536 !     size=req%recv(j)%size 
    537 !     req%recv(j)%value(size)=req%src_ind(i) 
    538 !   ENDDO 
    539 !        
    540 !   ireq=0 
    541 !   DO i=1,recv 
    542 !     ireq=ireq+1 
    543 !     CALL MPI_ISEND(req%recv(i)%value,req%recv(i)%size,MPI_INTEGER,req%recv(i)%rank,req%recv(i)%domain,comm_icosa, mpi_req(ireq),ierr)    
    544 !   ENDDO 
    545  
    546 !   DO i=1,send 
    547 !     ireq=ireq+1 
    548 !     ALLOCATE(req%send(i)%value(req%send(i)%size)) 
    549 !     CALL MPI_IRECV(req%send(i)%value,req%send(i)%size,MPI_INTEGER,req%send(i)%rank,req%send(i)%domain,comm_icosa, mpi_req(ireq),ierr)    
    550 !   ENDDO 
    551 !    
    552 !   CALL MPI_WAITALL(send+recv,mpi_req,status,ierr)  
    553      
    554  
    555    END SUBROUTINE Finalize_request  
     527  END SUBROUTINE Finalize_request  
    556528 
    557529 
     
    571543    REAL(rstd),POINTER :: buffer_r4(:,:,:)  
    572544    INTEGER,POINTER :: value(:)  
     545    INTEGER,POINTER :: sgn(:)  
    573546    TYPE(ARRAY),POINTER :: recv,send  
    574547    TYPE(t_request),POINTER :: req 
     
    634607            buffer_r2=>recv%buffer_r2 
    635608            value=>recv%value 
     609            sgn=>recv%sign 
    636610            DO n=1,recv%size 
    637               rval2d(value(n))=buffer_r2(n)   
     611              rval2d(value(n))=buffer_r2(n)*sgn(n)   
    638612            ENDDO         
    639613            DEALLOCATE(recv%buffer_r2) 
     
    697671            buffer_r3=>recv%buffer_r3 
    698672            value=>recv%value 
     673            sgn=>recv%sign 
    699674            DO n=1,recv%size 
    700               rval3d(value(n),:)=buffer_r3(n,:)   
     675              rval3d(value(n),:)=buffer_r3(n,:)*sgn(n)   
    701676            ENDDO         
    702677            DEALLOCATE(recv%buffer_r3) 
     
    760735            buffer_r4=>recv%buffer_r4 
    761736            value=>recv%value 
     737            sgn=>recv%sign 
    762738            DO n=1,recv%size 
    763               rval4d(value(n),:,:)=buffer_r4(n,:,:)   
     739              rval4d(value(n),:,:)=buffer_r4(n,:,:)*sgn(n)  
    764740            ENDDO         
    765741            DEALLOCATE(recv%buffer_r4) 
     
    797773        IF (field(ind)%ndim==2) THEN 
    798774          DO n=1,req%size 
    799             rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n)) 
     775            rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n))*req%target_sign(n) 
    800776          ENDDO 
    801777        ELSE IF (field(ind)%ndim==3) THEN 
    802778          DO n=1,req%size 
    803             rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:) 
     779            rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:)*req%target_sign(n) 
    804780          ENDDO 
    805781        ELSE IF (field(ind)%ndim==4) THEN 
    806782          DO n=1,req%size 
    807             rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:) 
     783            rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:)*req%target_sign(n) 
    808784          ENDDO 
    809785        ENDIF 
Note: See TracChangeset for help on using the changeset viewer.