/[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.1 - (show annotations) (download)
Mon Dec 15 02:35:29 2003 UTC (20 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint57t_post, checkpoint57o_post, checkpoint52l_pre, checkpoint52e_pre, hrcube4, checkpoint58e_post, checkpoint57v_post, checkpoint52n_post, checkpoint52j_post, checkpoint53d_post, checkpoint58u_post, checkpoint58w_post, checkpoint54a_pre, checkpoint57m_post, checkpoint55c_post, checkpoint54e_post, checkpoint52e_post, checkpoint57s_post, checkpoint54a_post, checkpoint53c_post, checkpoint57k_post, checkpoint55d_pre, checkpoint57d_post, checkpoint57g_post, checkpoint57b_post, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint56b_post, checkpoint57i_post, checkpoint57y_post, checkpoint57e_post, checkpoint52l_post, checkpoint55h_post, checkpoint58n_post, checkpoint58x_post, checkpoint52k_post, checkpoint57g_pre, checkpoint54b_post, checkpoint53b_pre, checkpoint55b_post, checkpoint58t_post, checkpoint58h_post, checkpoint54d_post, checkpoint56c_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55, checkpoint53a_post, checkpoint57f_pre, checkpoint57a_post, checkpoint54, checkpoint58q_post, checkpoint54f_post, checkpoint53b_post, checkpoint55g_post, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint55f_post, checkpoint59c, checkpoint59b, checkpoint59h, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, checkpoint53, checkpoint52d_post, eckpoint57e_pre, checkpoint57h_done, checkpoint58f_post, checkpoint53g_post, checkpoint52f_post, checkpoint57x_post, checkpoint57n_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint57f_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, hrcube5, checkpoint58o_post, checkpoint57z_post, checkpoint57c_post, checkpoint58y_post, checkpoint55e_post, checkpoint58k_post, checkpoint52i_post, checkpoint52j_pre, checkpoint58v_post, checkpoint53f_post, checkpoint55a_post, checkpoint53d_pre, checkpoint54c_post, checkpoint58s_post, checkpoint58p_post, checkpoint57j_post, checkpoint58b_post, checkpoint57h_pre, checkpoint58m_post, checkpoint57l_post, checkpoint52i_pre, checkpoint52h_pre, checkpoint52f_pre, checkpoint57h_post, hrcube_2, hrcube_3, checkpoint56a_post, checkpoint55d_post
coupler - multi-components interface: low-level communication S/R package.

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

  ViewVC Help
Powered by ViewVC 1.1.22