/[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.3 - (hide annotations) (download)
Mon Sep 14 16:22:01 2009 UTC (14 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint62, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61v, checkpoint61w, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.2: +2 -2 lines
change conversion S/R calls to use the right type (real*8).

1 jmc 1.3 C $Header: /u/gcmpack/MITgcm/pkg/compon_communic/comprecv_r8tiles.F,v 1.2 2007/10/08 23:58:20 jmc Exp $
2 jmc 1.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 jmc 1.3 call mitcplr_dbl2char( r8buf(9), recvdname )
64 jmc 1.1
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