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 double precision vectors by MPI |
---|
25 | ! Author: J. Brioude * |
---|
26 | ! March 2012 * |
---|
27 | ! * |
---|
28 | |
---|
29 | subroutine senddouble_mpi(tag,numpart2,chunksize,direc) |
---|
30 | |
---|
31 | use mpi_mod |
---|
32 | use com_mod |
---|
33 | implicit none |
---|
34 | include 'mpif.h' |
---|
35 | |
---|
36 | ! character :: varname*20 |
---|
37 | integer :: chunksize,numpart2,jj1 |
---|
38 | ! real(kind=dp) :: dummyr(numpart2), |
---|
39 | ! real(kind=dp) :: dummyr22(chunksize) |
---|
40 | integer :: myid,ierr,ntasks,ii,jdeb,jfin,jj,direc,tag |
---|
41 | ! integer :: MPI_COMM_WORLD |
---|
42 | |
---|
43 | integer :: jj2,from,jj3 |
---|
44 | integer, dimension(MPI_STATUS_SIZE) :: status |
---|
45 | |
---|
46 | call MPI_COMM_RANK ( MPI_COMM_WORLD, myid, ierr ) |
---|
47 | call MPI_COMM_SIZE ( MPI_COMM_WORLD, ntasks, ierr ) |
---|
48 | if (direc.eq.0) then ! the slaves get |
---|
49 | |
---|
50 | if (myid.eq.0) then |
---|
51 | do ii=1,ntasks-1 |
---|
52 | do jj2=1,chunksize |
---|
53 | jj=(jj2-1)*ntasks+ii+1 |
---|
54 | ! do jj=ii+1,numpart2+ii,ntasks |
---|
55 | ! jj2=(jj-ii-1)/ntasks+1 |
---|
56 | ! dummyr2(jj2)=dummyr(jj) |
---|
57 | if (tag.eq.11) dummyr22(jj2)=xtra1(jj) |
---|
58 | if (tag.eq.12) dummyr22(jj2)=ytra1(jj) |
---|
59 | enddo |
---|
60 | call MPI_SEND(dummyr22, chunksize, MPI_DOUBLE_PRECISION, ii,tag, MPI_COMM_WORLD, ierr) |
---|
61 | enddo |
---|
62 | ! chunksize2=int((numpart2-1)/ntasks)+1 |
---|
63 | ! chunksize2=chunksize |
---|
64 | ii=0 |
---|
65 | do jj=1,numpart2,ntasks |
---|
66 | ii=ii+1 |
---|
67 | jj2=jj |
---|
68 | enddo |
---|
69 | chunksize2=ii+numpart2-jj2 |
---|
70 | |
---|
71 | if (tag.eq.11) then |
---|
72 | do jj=1,numpart2,ntasks |
---|
73 | jj3=(jj-1)/ntasks+1 |
---|
74 | mpi_xtra1(jj3)=xtra1(jj) |
---|
75 | enddo |
---|
76 | mpi_xtra1(jj3:chunksize2)=xtra1(jj2:numpart2) |
---|
77 | elseif (tag.eq.12) then |
---|
78 | do jj=1,numpart2,ntasks |
---|
79 | jj3=(jj-1)/ntasks+1 |
---|
80 | mpi_ytra1(jj3)=ytra1(jj) |
---|
81 | enddo |
---|
82 | mpi_ytra1(jj3:chunksize2)=ytra1(jj2:numpart2) |
---|
83 | endif |
---|
84 | |
---|
85 | else ! the slaves receive |
---|
86 | if (tag.eq.11) call MPI_RECV(mpi_xtra1, chunksize, MPI_DOUBLE_PRECISION, 0,11, MPI_COMM_WORLD,status, ierr) |
---|
87 | if (tag.eq.12) call MPI_RECV(mpi_ytra1, chunksize, MPI_DOUBLE_PRECISION, 0,12, MPI_COMM_WORLD,status, ierr) |
---|
88 | endif |
---|
89 | |
---|
90 | else ! the master is going to get |
---|
91 | |
---|
92 | if (myid.gt.0) then !slaves send |
---|
93 | if (tag.eq.11) call MPI_SEND(mpi_xtra1, chunksize, MPI_DOUBLE_PRECISION, 0,11, MPI_COMM_WORLD, ierr) |
---|
94 | if (tag.eq.12) call MPI_SEND(mpi_ytra1, chunksize, MPI_DOUBLE_PRECISION, 0,12, MPI_COMM_WORLD, ierr) |
---|
95 | |
---|
96 | else ! the master gets |
---|
97 | |
---|
98 | do from =1,ntasks-1 |
---|
99 | call MPI_RECV(dummyr22, chunksize, MPI_DOUBLE_PRECISION, from,tag, MPI_COMM_WORLD,status,ierr) |
---|
100 | jj1=(from-1)*chunksize+1 |
---|
101 | jj2=from*chunksize |
---|
102 | if (tag.eq.11) xtra1(jj1:jj2)=dummyr22(1:chunksize) |
---|
103 | if (tag.eq.12) ytra1(jj1:jj2)=dummyr22(1:chunksize) |
---|
104 | enddo |
---|
105 | |
---|
106 | if (tag.eq.11) xtra1(jj2+1:numpart2)=mpi_xtra1(1:chunksize2) |
---|
107 | if (tag.eq.12) ytra1(jj2+1:numpart2)=mpi_ytra1(1:chunksize2) |
---|
108 | |
---|
109 | |
---|
110 | endif |
---|
111 | |
---|
112 | endif |
---|
113 | |
---|
114 | end subroutine senddouble_mpi |
---|