/[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.3 - (hide annotations) (download)
Mon Sep 14 16:22:02 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/compsend_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 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 jmc 1.3 call mitcplr_char2dbl( dataname, r8buf(9) )
42 jmc 1.1
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