source: branches/jerome/src_flexwrf_v3.1/sendreal2d_mpi.f90 @ 16

Last change on this file since 16 was 16, checked in by jebri, 11 years ago

sources for flexwrf v3.1

File size: 3.8 KB
Line 
1!***********************************************************************
2!* Copyright 2012,2013                                                *
3!* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine,       *
4!* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,*
5!* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso,          *
6!* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann,   *
7!*                                                                     *
8!* This file is part of FLEXPART WRF                                   *
9!*                                                                     *
10!* FLEXPART is free software: you can redistribute it and/or modify    *
11!* it under the terms of the GNU General Public License as published by*
12!* the Free Software Foundation, either version 3 of the License, or   *
13!* (at your option) any later version.                                 *
14!*                                                                     *
15!* FLEXPART is distributed in the hope that it will be useful,         *
16!* but WITHOUT ANY WARRANTY; without even the implied warranty of      *
17!* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
18!* GNU General Public License for more details.                        *
19!*                                                                     *
20!* You should have received a copy of the GNU General Public License   *
21!* along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
22!***********************************************************************
23
24  !   routine used to send 2D real array by MPI
25  !    Author: J. Brioude                                                      *
26  !    March 2012                                                           *
27        subroutine sendreal2d_mpi(tag2,numpart2,nspec2,chunksize,direc)
28
29      use mpi_mod
30      use com_mod
31          implicit none
32      include 'mpif.h'
33
34!      character :: varname*20
35       integer :: chunksize,numpart2,ks,jj1
36!      real :: dummyr(numpart2,nspec2),
37!      real :: dummyr2(chunksize)
38       integer :: myid,ierr,ntasks,ii,jdeb,jfin,jj,tag,direc
39!      integer :: MPI_COMM_WORLD
40
41       integer :: tag2,nspec2
42       integer :: jj2,from,jj3   
43  integer, dimension(MPI_STATUS_SIZE) :: status
44
45      call MPI_COMM_RANK ( MPI_COMM_WORLD, myid, ierr )
46      call MPI_COMM_SIZE ( MPI_COMM_WORLD, ntasks, ierr )
47
48      do ks=1,nspec2
49      tag=100+ks
50
51     if (direc.eq.0) then ! the slaves get
52
53      if (myid.eq.0) then
54       do ii=1,ntasks-1
55       do jj2=1,chunksize
56        jj=(jj2-1)*ntasks+ii+1
57!       do jj=ii+1,numpart2+ii,ntasks
58!        jj2=(jj-ii-1)/ntasks+1
59!        dummyr2(jj2)=dummyr(jj,ks)
60         dummyr2(jj2)=xmass1(jj,ks)
61        enddo
62       call MPI_SEND(dummyr2, chunksize, MPI_REAL, ii,tag,MPI_COMM_WORLD, ierr)
63       enddo
64
65!     chunksize2=int((numpart2-1)/ntasks)+1
66!     chunksize2=chunksize
67      ii=0
68      do jj=1,numpart2,ntasks
69      ii=ii+1
70      jj2=jj
71      enddo
72      chunksize2=ii+numpart2-jj2
73     do jj=1,numpart2,ntasks
74      jj3=(jj-1)/ntasks+1
75     mpi_xmass1(jj3,ks)=xmass1(jj,ks)
76     enddo
77     mpi_xmass1(jj3:chunksize2,ks)=xmass1(jj2:numpart2,ks)
78
79   else ! slave receive
80    call MPI_RECV(mpi_xmass1(1,ks), chunksize, MPI_REAL, 0,tag, MPI_COMM_WORLD,status, ierr)
81   endif
82
83   else ! the master is going to get
84
85   if (myid.gt.0) then !slaves send
86
87   call MPI_SEND(mpi_xmass1(1,ks), chunksize, MPI_REAL,0,tag, MPI_COMM_WORLD, ierr)
88
89    else ! the master gets
90     do from =1,ntasks-1
91    call MPI_RECV(dummyr2, chunksize, MPI_REAL, from,tag, MPI_COMM_WORLD,status,ierr)
92          jj1=(from-1)*chunksize+1
93          jj2=from*chunksize
94        xmass1(jj1:jj2,ks)=dummyr2(1:chunksize)
95     enddo
96     xmass1(jj2+1:numpart2,ks)=mpi_xmass1(1:chunksize2,ks)
97
98    endif
99
100   endif
101
102
103    enddo
104       end subroutine sendreal2d_mpi
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG