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

Annotation of /MITgcm/pkg/compon_communic/compsend_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, 7 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 compsend_r8tiles( dataname,Ni,Oi,Nj,Oj,Nk,Ti,Tj, arr )
6     implicit none
7     ! Arguments
8     character*(*) dataname
9     integer Ni,Oi,Nj,Oj,Nk,Ti,Tj
10     real*8 arr(1-Oi:Ni+Oi,1-Oj:Nj+Oj,Nk,Ti,Tj)
11     ! Predefined constants/arrays
12     #include "CPLR_SIG.h"
13     ! MPI variables
14     #include "mpif.h"
15     integer count,datatype,dest,tag,comm,ierr
16     ! Functions
17     integer generate_tag
18     ! Local
19     integer i,j,ij,bi,bj,k,bibj
20     ! ------------------------------------------------------------------
21    
22     if (HEADER_SIZE+Ni*Nj.gt.MAX_R8_BUFLEN)
23     & stop 'compsend_r8tiles: Nx*Ny too big'
24    
25     ! Foreach tile which is non-blank
26     do bibj=1,my_num_tiles
27    
28     bi=my_tile_bi(bibj)
29     bj=my_tile_bj(bibj)
30    
31     if (Ni.ne.my_tile_nx(bibj))
32     & stop 'compsend_r8tiles: ni != my_tile_nx'
33     if (Nj.ne.my_tile_ny(bibj))
34     & stop 'compsend_r8tiles: nj != my_tile_ny'
35    
36     ! Set up buffer
37     r8buf(1)=float( my_tile_i0(bibj) )
38     r8buf(2)=float( my_tile_nx(bibj) )
39     r8buf(3)=float( my_tile_j0(bibj) )
40     r8buf(4)=float( my_tile_ny(bibj) )
41     call mitcplr_char2real( dataname, r8buf(9) )
42    
43     ! Copy interior of tile to buffer
44     k=1
45     do j=1,Nj
46     do i=1,Ni
47     ij=HEADER_SIZE+i+Ni*(j-1)
48     r8buf(ij)=arr(i,j,k,bi,bj)
49     enddo
50     enddo
51    
52     ! Send message
53     count=HEADER_SIZE+Ni*Nj
54     datatype=MPI_DOUBLE_PRECISION
55     dest=my_coupler_rank
56     tag=generate_tag(103,bibj,dataname)
57     comm=MPI_COMM_myglobal
58    
59     if (VERB) then
60     write(LogUnit,*)
61     & 'compsend_r8tiles: calling MPI_Send dest=',dest,' tile=',bibj
62     write(LogUnit,*) 'compsend_r8tiles: dataname=',dataname
63     call flush(LogUnit)
64     endif
65     call MPI_Send( r8buf, count, datatype, dest, tag, comm, ierr )
66     if (VERB) then
67     write(LogUnit,*) 'compsend_r8tiles: returned ierr=',ierr
68     call flush(LogUnit)
69     endif
70    
71     if (ierr.ne.0) then
72     write(LogUnit,*) 'compsend_r8tiles: rank(W,G,L)=',
73     & my_rank_in_world,my_rank_in_global,my_rank_in_local,
74     & ' ierr=',ierr
75     stop 'compsend_r8tiles: MPI_Send failed'
76     endif
77    
78     enddo ! bibj
79    
80     ! ------------------------------------------------------------------
81     return
82     end
83     !=======================================================================

  ViewVC Help
Powered by ViewVC 1.1.22