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

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

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


Revision 1.4 - (hide annotations) (download)
Wed Nov 4 17:04:21 2015 UTC (8 years, 7 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.3: +9 -7 lines
fix error message

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

  ViewVC Help
Powered by ViewVC 1.1.22