ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
communication.F90
Go to the documentation of this file.
2 !
3 ! Contains the variable definition of the processor number,
4 ! myID and the number of processors, nProc, which belong to the
5 ! group defined by the communicator ADflow_comm_world. The range
6 ! of processor numbers is <0..Nproc-1>, i.e. the numbering
7 ! starts at 0. This is done for compatibility with MPI.
8 ! Furthermore this module contains the communication pattern for
9 ! all the multigrid levels.
10 !
11  use constants, only: inttype, realtype
12  implicit none
13  save
14 !
15 ! The definition of the derived data type commListType, which
16 ! stores the i,j and k indices as well as the block id of the
17 ! data to be communicated. Send lists may contain interpolants
18 ! since the indices may refer to a stencil, while the receive
19 ! list does not. All interpolations should be done on the send
20 ! side to keep message sizes to a minimum.
21 !
22 #ifndef USE_TAPENADE
24 
25  ! block(..): Local block id to which the cell/node belongs.
26  ! The dimension is equal to the number of entities
27  ! to be communicated with the particular processor.
28  ! indices(..,3): I, j and k indices of the data to be communicated.
29  ! For the first dimension, see block.
30  ! interp(..,..): Interpolants for indices that represent a cell
31  ! stencil (allocated only when needed, e.g. send
32  ! list for overset communication).
33  ! For the first dimension, see block.
34 
35  integer(kind=intType), pointer, dimension(:) :: block
36  integer(kind=intType), pointer, dimension(:, :) :: indices
37  real(kind=realtype), pointer, dimension(:, :) :: interp, interpd
38  real(kind=realtype), pointer, dimension(:, :) :: xcen
39  end type sendcommlisttype
40 
42 
43  ! block(..): Local block id to which the cell/node belongs.
44  ! The dimension is equal to the number of entities
45  ! to be communicated with the particular processor.
46  ! indices(..,3): I, j and k indices of the data to be communicated.
47  ! For the first dimension, see block.
48 
49  integer(kind=intType), pointer, dimension(:) :: block
50  integer(kind=intType), pointer, dimension(:, :) :: indices
51 
52  end type recvcommlisttype
53 !
54 ! The definition of the derived data type periodicDataType,
55 ! which stores the rotation matrix, the rotation center and the
56 ! translation vector of the periodic transformation, as well as
57 ! the halos to which this transformation must be applied.
58 !
60 
61  ! rotMatrix(3,3): Rotation matrix.
62  ! rotCenter(3): Coordinates of center of rotation.
63  ! translation(3): Translation vector.
64 
65  real(kind=realtype), dimension(3, 3) :: rotmatrix
66  real(kind=realtype), dimension(3) :: rotcenter, translation
67 
68  ! nHalos: # of halos to which this periodic
69  ! transformation must be applied.
70  ! block(nHalos): Local block id to which the halos belong.
71  ! indices(nHalos,3): I, j and k indices of the halos in
72  ! this block.
73 
74  integer(kind=intType) :: nhalos
75 
76  integer(kind=intType), pointer, dimension(:) :: block
77  integer(kind=intType), pointer, dimension(:, :) :: indices
78 
79  end type periodicdatatype
80 !
81 ! The definition of the derived data type commType, which
82 ! stores the communication pattern for a certain halo type for a
83 ! certain grid level.
84 !
85  type commtype
86 
87  ! nProcSend: # of procs, to whom this proc will send
88  ! nProcRecv: # of procs, from whom this proc will
89  ! receive.
90  ! sendProc(nProcSend): Send processor numbers.
91  ! recvProc(nProcRecv): Receive processor numbers.
92 
93  integer(kind=intType) :: nprocsend, nprocrecv
94  integer(kind=intType), pointer, dimension(:) :: sendproc, recvproc
95 
96  ! nsend(nProcSend): # of entities to send to other processors.
97  ! nrecv(nProcRecv): # of entities to receive from other processors.
98 
99  integer(kind=intType), pointer, dimension(:) :: nsend, nrecv
100 
101  ! nsendCum(0:NprocSend): cumulative version of nsend.
102  ! nrecvCum(0:NprocRecv): cumulative version of nrecv.
103 
104  integer(kind=intType), pointer, dimension(:) :: nsendcum
105  integer(kind=intType), pointer, dimension(:) :: nrecvcum
106 
107  ! indexSendProc(0:Nproc-1): index of the processors in sendProc.
108  ! If nothing is to be sent to a
109  ! processor this value is 0.
110  ! indexRecvProc(0:Nproc-1): index of the processors in recvProc.
111  ! If nothing is to be received from a
112  ! processor this value is 0.
113 
114  integer(kind=intType), pointer, dimension(:) :: indexsendproc
115  integer(kind=intType), pointer, dimension(:) :: indexrecvproc
116 
117  ! sendList(nProcSend): Indices and block ids to send to these
118  ! processors.
119  ! recvList(nProcRecv): Indices and block ids to receive from
120  ! these processors.
121 
122  type(sendcommlisttype), pointer, dimension(:) :: sendlist
123  type(recvcommlisttype), pointer, dimension(:) :: recvlist
124 
125  ! nPeriodic: # of periodic data arrays.
126  ! periodicData(nPeriodic): Periodic data and entities to which
127  ! the transformation must be applied.
128 
129  integer(kind=intType) :: nperiodic
130  type(periodicdatatype), pointer, dimension(:) :: periodicdata
131 
132  end type commtype
133 !
134 ! The definition of the derived data type internalCommType,
135 ! which stores the memory to memory copy on this processor for a
136 ! certain halo type.
137 !
139 
140  ! ncopy: # of entities to copy internally.
141 
142  ! donorBlock(ncopy): Local block id of the donor cell.
143  ! donorIndices(ncopy,3): The indices of the donor cell.
144  ! donorInterp(ncopy,3): Interpolants of the donor stencil.
145  ! Only allocated when needed, e.g.
146  ! overset communication).
147 
148  ! haloBlock(ncopy): Local block id of the halo cell.
149  ! haloIndices(ncopy,3): The indices of the halo cell.
150 
151  integer(kind=intType) :: ncopy
152 
153  integer(kind=intType), pointer, dimension(:) :: donorblock
154  integer(kind=intType), pointer, dimension(:, :) :: donorindices
155  real(kind=realtype), pointer, dimension(:, :) :: donorinterp, donorinterpd
156  real(kind=realtype), pointer, dimension(:, :) :: xcen
157 
158  integer(kind=intType), pointer, dimension(:) :: haloblock
159  integer(kind=intType), pointer, dimension(:, :) :: haloindices
160 
161  ! nPeriodic: # of periodic data arrays.
162  ! periodicData(nPeriodic): Periodic data and entities to which
163  ! the transformation must be applied.
164 
165  integer(kind=intType) :: nperiodic
166  type(periodicdatatype), pointer, dimension(:) :: periodicdata
167 
168  end type internalcommtype
169 !
170 ! Variables stored in this module.
171 !
172  ! ADflow_comm_world: The communicator of this processor group.
173  ! myID: My processor number in ADflow_comm_world.
174  ! nProc: The number of processors in ADflow_comm_world.
176 
177  integer :: myid, nproc
178 
179  ! commPatternCell_1st(nLevel): The communication pattern for 1st
180  ! level cell halo's on the multiple
181  ! grids.
182  ! commPatternCell_2nd(nLevel): The communication pattern for 2nd
183  ! level cell halo's (including the
184  ! 1st level) on the multiple grids.
185  ! commPatternNode_1st(nLevel): The communication pattern for 1st
186  ! level node halo's on the multiple
187  ! grids.
188 
189  type(commtype), allocatable, dimension(:) :: commpatterncell_1st
190  type(commtype), allocatable, dimension(:) :: commpatterncell_2nd
191  type(commtype), allocatable, dimension(:) :: commpatternnode_1st
192  type(commtype), allocatable, target, dimension(:, :) :: commpatternoverset
193 
194  ! internalCell_1st(nLevel): Memory to memory copies for 1st level
195  ! cell halo's on the multiple grids.
196  ! internalCell_2nd(nLevel): Memory to memory copies for 2nd level
197  ! cell halo's on the multiple grids.
198  ! internalNode_1st(nLevel): Memory to memory copies for 1st level
199  ! node halo's on the multiple grids.
200 
201  type(internalcommtype), allocatable, dimension(:) :: internalcell_1st
202  type(internalcommtype), allocatable, dimension(:) :: internalcell_2nd
203  type(internalcommtype), allocatable, dimension(:) :: internalnode_1st
204  type(internalcommtype), allocatable, target, dimension(:, :) :: internaloverset
205 
206  ! sendBufferSize_1to1: Size of the send buffer needed to perform
207  ! all 1 to 1 communication.
208  ! recvBufferSize_1to1: Idem for the receive buffer.
209  ! sendBufferSizeOver: Size of the send buffer needed to perform
210  ! all overset communication.
211  ! recvBufferSizeOver: Idem for the receive buffer.
212  ! sendBufferSize: Size of the send buffer to perform all
213  ! possible communication.
214  ! recvBufferSize: Idem for the receive buffer.
215 
216  integer(kind=intType) :: sendbuffersize_1to1, sendbuffersize
217  integer(kind=intType) :: recvbuffersize_1to1, recvbuffersize
218  integer(kind=intType) :: sendbuffersizeover, recvbuffersizeover
219 
220  ! sendBuffer: Buffer used to store the info to be send during
221  ! a nonblocking communication.
222  ! recvBuffer: Buffer used to store the info to be received
223  ! during a nonblocking communication.
224  ! sendRequests: Array of requests for the nonblocking sends.
225  ! recvRequests: Array of requests for the nonblocking receives.
226 
227  real(kind=realtype), allocatable, dimension(:) :: sendbuffer
228  real(kind=realtype), allocatable, dimension(:) :: recvbuffer
229 
230  integer, allocatable, dimension(:) :: sendrequests, recvrequests
231 #endif
232 end module communication
Definition: block.F90:1
integer(kind=inttype) recvbuffersizeover
integer(kind=inttype) sendbuffersize
real(kind=realtype), dimension(:), allocatable recvbuffer
integer, dimension(:), allocatable recvrequests
real(kind=realtype), dimension(:), allocatable sendbuffer
type(internalcommtype), dimension(:, :), allocatable, target internaloverset
type(internalcommtype), dimension(:), allocatable internalcell_1st
integer(kind=inttype) sendbuffersizeover
type(internalcommtype), dimension(:), allocatable internalcell_2nd
integer, dimension(:), allocatable sendrequests
type(commtype), dimension(:, :), allocatable, target commpatternoverset
type(commtype), dimension(:), allocatable commpatterncell_1st
integer(kind=inttype) recvbuffersize
integer adflow_comm_self
type(commtype), dimension(:), allocatable commpatterncell_2nd
type(commtype), dimension(:), allocatable commpatternnode_1st
integer(kind=inttype) sendbuffersize_1to1
type(internalcommtype), dimension(:), allocatable internalnode_1st
integer(kind=inttype) recvbuffersize_1to1
integer adflow_comm_world