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

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

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


Revision 1.3 - (show annotations) (download)
Wed Nov 4 17:04:21 2015 UTC (8 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, HEAD
Changes since 1.2: +9 -7 lines
fix error message

1 C $Header: /u/gcmpack/MITgcm/pkg/compon_communic/coupsend_r4tiles.F,v 1.2 2007/10/08 23:58:20 jmc Exp $
2 C $Name: $
3
4 !=======================================================================
5 subroutine coupsend_r4tiles( component, dataname, Nx, Ny, arr )
6 implicit none
7 ! Predefined constants/arrays
8 #include "CPLR_SIG.h"
9 ! MPI variables
10 #include "mpif.h"
11 ! Arguments
12 character*(*) component
13 character*(*) dataname
14 integer Nx,Ny
15 real*4 arr(Nx,Ny)
16 ! Functions
17 integer mitcplr_match_comp
18 integer generate_tag
19 external mitcplr_match_comp
20 external generate_tag
21 ! Local
22 integer count,dtype,dest,tag,comm,ierr
23 integer compind,numprocs
24 integer i,j,ij,n,bibj
25 integer Ni,Io,Nj,Jo
26 ! ------------------------------------------------------------------
27
28 ! Establish who I am communicating with
29 compind=mitcplr_match_comp( component )
30 if (compind.le.0) stop 'coupsend_r4tiles: Bad component id'
31 comm=MPI_COMM_compcplr( compind )
32 numprocs=num_component_procs(compind)
33 if (numprocs.lt.1) then
34 write(LogUnit,*) 'coupsend_r4tiles: compind = ',compind
35 stop 'coupsend_r4tiles: numprocs < 1'
36 endif
37 if (VERB)
38 & write(LogUnit,*) 'coupsend_r4tiles: ',component_Name(compind)
39 if (VERB)
40 & write(LogUnit,*) 'coupsend_r4tiles: dataname=',dataname
41
42 ! Foreach component process
43 do n=1,numprocs
44
45 ! Foreach tile on that process
46 do bibj=1,component_num_tiles(n,compind)
47
48 ! Create header
49 Io=component_tile_i0(bibj,n,compind)
50 Jo=component_tile_j0(bibj,n,compind)
51 Ni=component_tile_nx(bibj,n,compind)
52 Nj=component_tile_ny(bibj,n,compind)
53 r4buf(1)=float( Io )
54 r4buf(2)=float( Jo )
55 r4buf(3)=float( Ni )
56 r4buf(4)=float( Nj )
57 call mitcplr_char2real( dataname, r4buf(9) )
58
59 ! Pack data
60 do j=1,Nj
61 do i=1,Ni
62 ij=HEADER_SIZE+i+Ni*(j-1)
63 r4buf(ij)=arr(Io+i-1,Jo+j-1)
64 enddo
65 enddo
66
67 ! Send message
68 count=HEADER_SIZE+Ni*Nj
69 dtype=MPI_REAL
70 tag=generate_tag(123,bibj,dataname)
71 dest=rank_component_procs(n,compind)
72
73 if (VERB) then
74 write(LogUnit,*) 'coupsend_r4tiles: calling MPI_Send dest=',
75 & dest,' proc=',n,'/',numprocs,' tile=',bibj
76 call flush(LogUnit)
77 endif
78 call MPI_Send( r4buf, count, dtype, dest, tag, comm, ierr )
79 if (VERB) then
80 write(LogUnit,*) 'coupsend_r4tiles: returned ierr=',ierr
81 call flush(LogUnit)
82 endif
83
84 if (ierr.ne.0) then
85 write(LogUnit,*) 'coupsend_r4tiles: rank(W,G)=',
86 & my_rank_in_world,my_rank_in_global,
87 & ' ierr=',ierr
88 stop 'coupsend_r4tiles: MPI_Send failed'
89 endif
90
91 enddo ! bibj
92
93 enddo ! n
94
95 ! ------------------------------------------------------------------
96 return
97 end
98 !=======================================================================

  ViewVC Help
Powered by ViewVC 1.1.22