/[MITgcm]/MITgcm/pkg/compon_communic/comprecv_r8tiles.F
ViewVC logotype

Annotation of /MITgcm/pkg/compon_communic/comprecv_r8tiles.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (hide annotations) (download)
Mon Oct 8 23:58:20 2007 UTC (16 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59k, checkpoint59j, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q
Changes since 1.1: +3 -0 lines
add missing cvs $Header:$ or $Name:$

1 jmc 1.2 C $Header: $
2     C $Name: $
3    
4 jmc 1.1 !=======================================================================
5     subroutine comprecv_r8tiles( dataname, Ni,Oi,Nj,Oj,Nk,Tx,Ty, arr )
6     implicit none
7     ! Arguments
8     character*(*) dataname
9     integer Ni,Oi,Nj,Oj,Io,Jo,Nk,Tx,Ty
10     real*8 arr(1-Oi:Ni+Oi,1-Oj:Nj+Oj,Nk,Tx,Ty)
11     ! Predefined constants/arrays
12     #include "CPLR_SIG.h"
13     ! MPI variables
14     #include "mpif.h"
15     integer count,dtype,rank,tag,comm,ierr
16     integer stat(MPI_STATUS_SIZE)
17     ! Functions
18     integer generate_tag
19     ! Local
20     integer i,j,ij,nx,ny,k,bibj,bi,bj
21     character*(MAXLEN_COMP_NAME) recvdname
22     ! ------------------------------------------------------------------
23    
24     if (HEADER_SIZE+Ni*Nj.gt.MAX_R8_BUFLEN)
25     & stop 'comprecv_r8tiles: Nx*Ny too big'
26    
27     ! Foreach tile which is non-blank
28     do bibj=1,my_num_tiles
29    
30     bi=my_tile_bi(bibj)
31     bj=my_tile_bj(bibj)
32    
33     ! Receive message
34     count=MAX_R8_BUFLEN
35     dtype=MPI_DOUBLE_PRECISION
36     tag=generate_tag(123,bibj,dataname)
37     rank=my_coupler_rank
38     comm=MPI_COMM_myglobal
39    
40     if (VERB) then
41     write(LogUnit,*) 'comprecv_r8tiles: calling MPI_Recv rank=',rank
42     write(LogUnit,*) 'comprecv_r8tiles: dataname=',dataname
43     call flush(LogUnit)
44     endif
45     call MPI_Recv(r8buf, count, dtype, rank, tag, comm, stat, ierr)
46     if (VERB) then
47     write(LogUnit,*) 'comprecv_r8tiles: returned ierr=',ierr
48     call flush(LogUnit)
49     endif
50    
51     if (ierr.ne.0) then
52     write(LogUnit,*) 'comprecv_r8tiles: rank(W,G)=',
53     & my_rank_in_world,my_rank_in_global,
54     & ' ierr=',ierr
55     stop 'comprecv_r8tiles: MPI_Recv failed'
56     endif
57    
58     ! Extract buffer
59     Io=int(0.5+r8buf(1))
60     Jo=int(0.5+r8buf(2))
61     nx=int(0.5+r8buf(3))
62     ny=int(0.5+r8buf(4))
63     call mitcplr_real2char( r8buf(9), recvdname )
64    
65     if (Io.ne.my_tile_i0(bibj)) stop 'comprecv_r8tiles: bad Io'
66     if (Jo.ne.my_tile_j0(bibj)) stop 'comprecv_r8tiles: bad Jo'
67     if (nx.ne.my_tile_nx(bibj)) stop 'comprecv_r8tiles: bad nx'
68     if (ny.ne.my_tile_ny(bibj)) stop 'comprecv_r8tiles: bad ny'
69     if (recvdname .ne. dataname) then
70     write(LogUnit,*) 'comprecv_r8tiles: recvdname = ',recvdname
71     write(LogUnit,*) 'comprecv_r8tiles: dataname = ',dataname
72     stop 'comprecv_r8tiles: recvdname != dataname'
73     endif
74    
75     ! Copy buffer to interior of tile
76     k=1
77     do j=1,Nj
78     do i=1,Ni
79     ij=HEADER_SIZE+i+Ni*(j-1)
80     arr(i,j,k,bi,bj)=r8buf(ij)
81     enddo
82     enddo
83    
84     enddo ! bibj
85    
86     ! ------------------------------------------------------------------
87     return
88     end
89     !=======================================================================

  ViewVC Help
Powered by ViewVC 1.1.22