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

Contents of /MITgcm/pkg/compon_communic/compsend_r8tiles.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 compsend_r8tiles( dataname,Ni,Oi,Nj,Oj,Nk,Ti,Tj, arr )
3 implicit none
4 ! Arguments
5 character*(*) dataname
6 integer Ni,Oi,Nj,Oj,Nk,Ti,Tj
7 real*8 arr(1-Oi:Ni+Oi,1-Oj:Nj+Oj,Nk,Ti,Tj)
8 ! Predefined constants/arrays
9 #include "CPLR_SIG.h"
10 ! MPI variables
11 #include "mpif.h"
12 integer count,datatype,dest,tag,comm,ierr
13 ! Functions
14 integer generate_tag
15 ! Local
16 integer i,j,ij,bi,bj,k,bibj
17 ! ------------------------------------------------------------------
18
19 if (HEADER_SIZE+Ni*Nj.gt.MAX_R8_BUFLEN)
20 & stop 'compsend_r8tiles: Nx*Ny too big'
21
22 ! Foreach tile which is non-blank
23 do bibj=1,my_num_tiles
24
25 bi=my_tile_bi(bibj)
26 bj=my_tile_bj(bibj)
27
28 if (Ni.ne.my_tile_nx(bibj))
29 & stop 'compsend_r8tiles: ni != my_tile_nx'
30 if (Nj.ne.my_tile_ny(bibj))
31 & stop 'compsend_r8tiles: nj != my_tile_ny'
32
33 ! Set up buffer
34 r8buf(1)=float( my_tile_i0(bibj) )
35 r8buf(2)=float( my_tile_nx(bibj) )
36 r8buf(3)=float( my_tile_j0(bibj) )
37 r8buf(4)=float( my_tile_ny(bibj) )
38 call mitcplr_char2real( dataname, r8buf(9) )
39
40 ! Copy interior of tile to buffer
41 k=1
42 do j=1,Nj
43 do i=1,Ni
44 ij=HEADER_SIZE+i+Ni*(j-1)
45 r8buf(ij)=arr(i,j,k,bi,bj)
46 enddo
47 enddo
48
49 ! Send message
50 count=HEADER_SIZE+Ni*Nj
51 datatype=MPI_DOUBLE_PRECISION
52 dest=my_coupler_rank
53 tag=generate_tag(103,bibj,dataname)
54 comm=MPI_COMM_myglobal
55
56 if (VERB) then
57 write(LogUnit,*)
58 & 'compsend_r8tiles: calling MPI_Send dest=',dest,' tile=',bibj
59 write(LogUnit,*) 'compsend_r8tiles: dataname=',dataname
60 call flush(LogUnit)
61 endif
62 call MPI_Send( r8buf, count, datatype, dest, tag, comm, ierr )
63 if (VERB) then
64 write(LogUnit,*) 'compsend_r8tiles: returned ierr=',ierr
65 call flush(LogUnit)
66 endif
67
68 if (ierr.ne.0) then
69 write(LogUnit,*) 'compsend_r8tiles: rank(W,G,L)=',
70 & my_rank_in_world,my_rank_in_global,my_rank_in_local,
71 & ' ierr=',ierr
72 stop 'compsend_r8tiles: MPI_Send failed'
73 endif
74
75 enddo ! bibj
76
77 ! ------------------------------------------------------------------
78 return
79 end
80 !=======================================================================

  ViewVC Help
Powered by ViewVC 1.1.22