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

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

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


Revision 1.1 - (show annotations) (download)
Mon Dec 15 02:35:29 2003 UTC (20 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint57t_post, checkpoint57o_post, checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint58e_post, checkpoint57v_post, checkpoint52n_post, checkpoint52j_post, checkpoint53d_post, checkpoint58u_post, checkpoint58w_post, checkpoint54a_pre, checkpoint57m_post, checkpoint55c_post, checkpoint54e_post, checkpoint52e_post, checkpoint57s_post, checkpoint54a_post, checkpoint53c_post, checkpoint57k_post, checkpoint55d_pre, checkpoint57d_post, checkpoint57g_post, checkpoint57b_post, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint56b_post, checkpoint57i_post, checkpoint57y_post, checkpoint57e_post, checkpoint52l_post, checkpoint55h_post, checkpoint58n_post, checkpoint58x_post, checkpoint52k_post, checkpoint57g_pre, checkpoint54b_post, checkpoint53b_pre, checkpoint55b_post, checkpoint58t_post, checkpoint58h_post, checkpoint54d_post, checkpoint56c_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55, checkpoint53a_post, checkpoint57f_pre, checkpoint57a_post, checkpoint54, checkpoint58q_post, checkpoint54f_post, checkpoint53b_post, checkpoint55g_post, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint55f_post, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, checkpoint53, checkpoint52d_post, eckpoint57e_pre, checkpoint57h_done, checkpoint58f_post, checkpoint53g_post, checkpoint52f_post, checkpoint57x_post, checkpoint57n_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint57f_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, hrcube5, checkpoint58o_post, checkpoint57z_post, checkpoint57c_post, checkpoint58y_post, checkpoint55e_post, checkpoint58k_post, checkpoint52i_post, checkpoint52j_pre, checkpoint58v_post, checkpoint53f_post, checkpoint55a_post, checkpoint53d_pre, checkpoint54c_post, checkpoint58s_post, checkpoint58p_post, checkpoint57j_post, checkpoint58b_post, checkpoint57h_pre, checkpoint58m_post, checkpoint57l_post, checkpoint52i_pre, checkpoint52h_pre, checkpoint52f_pre, checkpoint57h_post, hrcube_2, hrcube_3, checkpoint56a_post, checkpoint55d_post
coupler - multi-components interface: low-level communication S/R package.

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

  ViewVC Help
Powered by ViewVC 1.1.22