1 | module columnslab |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | contains |
---|
6 | |
---|
7 | subroutine slabs_to_columns(slab,col) |
---|
8 | ! In the cfd code, data is partitioned across processors in equal |
---|
9 | ! sub-blocks in dimension 3, i.e. in horizontal slabs. It is sometimes |
---|
10 | ! convenient to exchange data so that it is partitioned along dimension |
---|
11 | ! 1, i.e. in vertical columns. This routine redistributes the data |
---|
12 | ! array "slab" into column orientation "col". |
---|
13 | |
---|
14 | use mpi_parameters |
---|
15 | use boundary_information |
---|
16 | include '../input/problem_size.h' |
---|
17 | include 'mpif.h' |
---|
18 | |
---|
19 | integer :: tag,ierr,status_array(MPI_STATUS_SIZE,4*numprocs+1) |
---|
20 | integer :: dest,source,istart,iend,kstart,kend,reqs(4*numprocs) |
---|
21 | integer :: sizeofreal,nreqs |
---|
22 | integer :: iproc,numwords |
---|
23 | integer :: m1,m2,m3 |
---|
24 | real,dimension(:,:,:),intent(in) :: slab |
---|
25 | real,dimension(:,:,:),intent(out) :: col |
---|
26 | |
---|
27 | if( boundary_type2(1,1) .ne. 'periodic' .or. boundary_type4(1,1) .ne. 'periodic' ) then |
---|
28 | write(0,*) 'myid,boundary_type2(1,1) = ',myid,boundary_type2(1,1) |
---|
29 | write(0,*) ' USING PERIODIC ENDPOINT LOGIC IN SLABS_TO_COLUMNS W/O PERIODIC BCS !!!' |
---|
30 | stop |
---|
31 | endif |
---|
32 | |
---|
33 | m1=size(col,1) |
---|
34 | m2=size(slab,2) |
---|
35 | m3=size(slab,3) |
---|
36 | |
---|
37 | !!******************************************************************************* |
---|
38 | !!************* DO THE DATA EXCHANGES **************************************** |
---|
39 | !!******************************************************************************* |
---|
40 | !! Ship portions of my local data "slab" slab needed by other nodes |
---|
41 | !! and recieve remote data required to construct my local data "column" col . |
---|
42 | |
---|
43 | !! All processors initiatze nonblocking receives, establishing buffers for the |
---|
44 | !! sends |
---|
45 | nreqs=0 |
---|
46 | dest=myid |
---|
47 | numwords=m1*m2*m3 |
---|
48 | tag=1 |
---|
49 | |
---|
50 | do iproc=0,numprocs-1 |
---|
51 | source=iproc |
---|
52 | kstart=(iproc*m3)+1 |
---|
53 | if ( dest .ne. source ) then !! bypass mpi when source=dest=myid |
---|
54 | nreqs=nreqs+1 |
---|
55 | call MPI_IRECV(col(1,1,kstart),numwords,MPI_REAL,source,tag,comm,reqs(nreqs),ierr) |
---|
56 | if ( ierr .ne. 0 ) then |
---|
57 | stop 'MPI IRECV ERROR IN SLABS TO COLS' |
---|
58 | endif |
---|
59 | endif |
---|
60 | enddo |
---|
61 | |
---|
62 | call MPI_BARRIER(comm,ierr) !! make sure all processors receives are posted |
---|
63 | !! to avoid buffering |
---|
64 | |
---|
65 | !! All processors now initiate nonblocking ready sends, writing to the |
---|
66 | !! pre-established dest locations |
---|
67 | source=myid |
---|
68 | numwords=1 |
---|
69 | tag = 1 |
---|
70 | |
---|
71 | do iproc=0,numprocs-1 |
---|
72 | dest=iproc |
---|
73 | istart=(iproc*m1)+1 |
---|
74 | if ( dest .ne. source ) then !! bypass mpi when source=dest=myid |
---|
75 | nreqs=nreqs+1 |
---|
76 | call MPI_IRSEND(slab(istart,1,1),numwords,subslice,dest,tag,comm,reqs(nreqs),ierr) |
---|
77 | if ( ierr .ne. 0 ) then |
---|
78 | stop 'MPI IRSEND ERROR IN SLABS TO COLS' |
---|
79 | endif |
---|
80 | endif |
---|
81 | enddo |
---|
82 | |
---|
83 | !! case where source=dest=myid |
---|
84 | kstart=myid*m3+1 |
---|
85 | kend=kstart+m3-1 |
---|
86 | istart=myid*m1+1 |
---|
87 | iend=istart+m1-1 |
---|
88 | col(:,:,kstart:kend)=slab(istart:iend,:,:) |
---|
89 | |
---|
90 | call MPI_WAITALL(nreqs,reqs(1),status_array,ierr) |
---|
91 | |
---|
92 | end subroutine slabs_to_columns |
---|
93 | |
---|
94 | |
---|
95 | subroutine columns_to_slabs(col,slab) |
---|
96 | ! In the cfd code, data is partitioned across processors in equal |
---|
97 | ! sub-blocks in dimension 3, i.e. in horizontal slabs. It is sometimes |
---|
98 | ! convenient to exchange data so that it is partitioned along dimension |
---|
99 | ! 1, i.e. in vertical columns. This routine redistributes the data |
---|
100 | ! array "col" into slab orientation "slab". |
---|
101 | |
---|
102 | use mpi_parameters |
---|
103 | use boundary_information |
---|
104 | include '../input/problem_size.h' |
---|
105 | include 'mpif.h' |
---|
106 | |
---|
107 | integer :: tag,ierr,status_array(MPI_STATUS_SIZE,2*numprocs) |
---|
108 | integer :: dest,source,istart,iend,kstart,kend,reqs(4*numprocs+1) |
---|
109 | integer :: sizeofreal,nreqs,j |
---|
110 | integer :: iproc,numwords |
---|
111 | integer :: m1,m2,m3 |
---|
112 | real,dimension(:,:,:),intent(out) :: slab |
---|
113 | real,dimension(:,:,:),intent(in) :: col |
---|
114 | |
---|
115 | if( boundary_type2(1,1) .ne. 'periodic' .or. boundary_type4(1,1) .ne. 'periodic' ) then |
---|
116 | write(0,*) ' USING PERIODIC ENDPOINT LOGIC IN SLABS_TO_COLUMNS W/O PERIODIC BCS !!!' |
---|
117 | stop |
---|
118 | endif |
---|
119 | |
---|
120 | m1=size(col,1) |
---|
121 | m2=size(slab,2) |
---|
122 | m3=size(slab,3) |
---|
123 | |
---|
124 | !!****************************************************************************** |
---|
125 | !!************* DO THE DATA EXCHANGES *************************************** |
---|
126 | !!****************************************************************************** |
---|
127 | !!Ship portions of my local data "column" x%ColVals needed by other nodes |
---|
128 | !!and recieve remote data required to construct my local data "slab" x%SlabVals. |
---|
129 | |
---|
130 | !!Post the receives first....immediate mode, ie return before data is actually transferred |
---|
131 | nreqs=0 |
---|
132 | dest=myid |
---|
133 | numwords=1 |
---|
134 | tag = 2 |
---|
135 | |
---|
136 | do iproc=0,numprocs-1 |
---|
137 | source=iproc |
---|
138 | istart=iproc*m1+1 !! start index for receiving "subslab" data x |
---|
139 | if( dest .ne. source ) then !! bypass mpi when source=dest=myid |
---|
140 | nreqs=nreqs+1 |
---|
141 | call MPI_IRECV(slab(istart,1,1),numwords,subslice,source,tag,comm,reqs(nreqs),ierr) |
---|
142 | if( ierr .ne. 0 ) then |
---|
143 | stop 'MPI IRECV ERROR IN COLS TO SLABS' |
---|
144 | endif |
---|
145 | endif |
---|
146 | enddo |
---|
147 | |
---|
148 | call MPI_BARRIER(comm,ierr) !! make sure all processors receives are posted |
---|
149 | !! to avoid buffering |
---|
150 | |
---|
151 | !!All processors have issued all their receives, now the sends can be done w/o buffering. |
---|
152 | !!Use immediate-ready-send meaning, return before transfer is finished and guarantee that |
---|
153 | !!the target locations on the destination processors are ready for writing. Faster. |
---|
154 | source=myid |
---|
155 | numwords=m1*m2*m3 |
---|
156 | tag = 2 |
---|
157 | |
---|
158 | do iproc=0,numprocs-1 |
---|
159 | dest=iproc |
---|
160 | kstart=iproc*m3+1 !! start index for outgoing column data in x%ColVals |
---|
161 | if( dest .ne. source ) then !! bypass mpi when source=dest=myid |
---|
162 | nreqs=nreqs+1 |
---|
163 | call MPI_IRSEND(col(1,1,kstart),numwords,MPI_REAL,dest,tag,comm,reqs(nreqs),ierr) |
---|
164 | if( ierr .ne. 0 ) then |
---|
165 | stop 'MPI IRSEND ERROR IN COLS TO SLABS' |
---|
166 | endif |
---|
167 | endif |
---|
168 | enddo |
---|
169 | |
---|
170 | !! case where source=dest=myid |
---|
171 | istart = myid*m1 + 1 |
---|
172 | iend = istart + m1 - 1 |
---|
173 | kstart = myid*m3 + 1 |
---|
174 | kend = kstart + m3 - 1 |
---|
175 | slab(istart:iend,:,:) = col(:,:,kstart:kend) |
---|
176 | |
---|
177 | !!make sure all communications are done before proceeding |
---|
178 | call MPI_WAITALL(nreqs,reqs,status_array,ierr) |
---|
179 | !****************************************************************************** |
---|
180 | !****************************************************************************** |
---|
181 | !****************************************************************************** |
---|
182 | |
---|
183 | end subroutine columns_to_slabs |
---|
184 | |
---|
185 | |
---|
186 | end module |
---|
187 | |
---|
188 | |
---|
189 | |
---|