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

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

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


Revision 1.3 - (show annotations) (download)
Mon Sep 14 16:22:03 2009 UTC (14 years, 9 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, checkpoint65p, 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
Changes since 1.2: +2 -2 lines
change conversion S/R calls to use the right type (real*8).

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

  ViewVC Help
Powered by ViewVC 1.1.22