-
Notifications
You must be signed in to change notification settings - Fork 1
/
comms5.F90
157 lines (137 loc) · 6.72 KB
/
comms5.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
module comms5
use comms_common
use params
use mpi
implicit none
integer :: halo_5_xup_send(4:4), halo_5_xdn_send(4:4)
integer :: halo_5_xup_recv(4:4), halo_5_xdn_recv(4:4)
integer :: halo_5_yup_send(4:4), halo_5_ydn_send(4:4)
integer :: halo_5_yup_recv(4:4), halo_5_ydn_recv(4:4)
integer :: halo_5_tup_send(4:4), halo_5_tdn_send(4:4)
integer :: halo_5_tup_recv(4:4), halo_5_tdn_recv(4:4)
contains
! Initialise a single halo type.
! direction: direction to be communicated; x=0, y=1, t=2
! position: coordinate in direction of the region
! i5, i6: size in extra directions
! datatype: datatype to transfer
! typetarget: array element to put type into;
! should be halo_{4,4_real,5,6}_{x,y,t}{up,dn}_{send,recv}(size{4,5}[, size6])
subroutine init_single_halo_type_5(direction, position, size5, typetarget)
integer, intent(in) :: direction, position, size5
integer, intent(out) :: typetarget
integer, dimension(5) :: sizes, subsizes, starts
integer :: ierr
sizes = (/kthird, ksizex_l + 2, ksizey_l + 2, ksizet_l + 2, size5/)
subsizes = (/kthird, ksizex_l, ksizey_l, ksizet_l, size5/)
subsizes(direction + 2) = 1
starts = (/0, 1, 1, 1, 0/)
starts(direction + 2) = position
call MPI_Type_Create_Subarray(5, sizes, subsizes, starts, MPI_Order_Fortran, &
& MPI_Double_Complex, typetarget, ierr)
call MPI_Type_Commit(typetarget, ierr)
return
end subroutine init_single_halo_type_5
subroutine init_halo_types_5()
integer, parameter :: nsize5 = 1
integer :: size5(nsize5) = (/4/)
integer :: i5
do i5 = 1, nsize5
call init_single_halo_type_5(0, 1, size5(i5), halo_5_xdn_send(size5(i5)))
call init_single_halo_type_5(0, 0, size5(i5), halo_5_xdn_recv(size5(i5)))
call init_single_halo_type_5(0, ksizex_l, size5(i5), halo_5_xup_send(size5(i5)))
call init_single_halo_type_5(0, ksizex_l + 1, size5(i5), halo_5_xup_recv(size5(i5)))
call init_single_halo_type_5(1, 1, size5(i5), halo_5_ydn_send(size5(i5)))
call init_single_halo_type_5(1, 0, size5(i5), halo_5_ydn_recv(size5(i5)))
call init_single_halo_type_5(1, ksizey_l, size5(i5), halo_5_yup_send(size5(i5)))
call init_single_halo_type_5(1, ksizey_l + 1, size5(i5), halo_5_yup_recv(size5(i5)))
call init_single_halo_type_5(2, 1, size5(i5), halo_5_tdn_send(size5(i5)))
call init_single_halo_type_5(2, 0, size5(i5), halo_5_tdn_recv(size5(i5)))
call init_single_halo_type_5(2, ksizet_l, size5(i5), halo_5_tup_send(size5(i5)))
call init_single_halo_type_5(2, ksizet_l + 1, size5(i5), halo_5_tup_recv(size5(i5)))
end do
end subroutine init_halo_types_5
subroutine init_halo_update_5(size5, Array, tag, reqs)
!
integer, intent(in) :: size5, tag
complex(dp), intent(inout) :: Array(kthird, 0:ksizex_l + 1, 0:ksizey_l + 1, &
& 0:ksizet_l + 1, size5)
integer, intent(out) :: reqs(12)
integer :: ip_xup, ip_xdn, ip_yup, ip_ydn, ip_tup, ip_tdn
integer :: tag_offset
integer :: ierr
tag_offset = 6*tag
! Start send and receive in x direction
call MPI_Cart_Shift(comm, 0, 1, ip_xdn, ip_xup, ierr)
call MPI_Send_init(Array, 1, halo_5_xup_send(size5), ip_xup, 0 + tag_offset, comm, &
& reqs(1), ierr)
call MPI_Recv_init(Array, 1, halo_5_xdn_recv(size5), ip_xdn, 0 + tag_offset, comm, &
& reqs(2), ierr)
call MPI_Send_init(Array, 1, halo_5_xdn_send(size5), ip_xdn, 1 + tag_offset, comm, &
& reqs(3), ierr)
call MPI_Recv_init(Array, 1, halo_5_xup_recv(size5), ip_xup, 1 + tag_offset, comm, &
& reqs(4), ierr)
! Start send and receive in y direction
call MPI_Cart_Shift(comm, 1, 1, ip_ydn, ip_yup, ierr)
call MPI_Send_init(Array, 1, halo_5_yup_send(size5), ip_yup, 2 + tag_offset, comm, &
& reqs(5), ierr)
call MPI_Recv_init(Array, 1, halo_5_ydn_recv(size5), ip_ydn, 2 + tag_offset, comm, &
& reqs(6), ierr)
call MPI_Send_init(Array, 1, halo_5_ydn_send(size5), ip_ydn, 3 + tag_offset, comm, &
& reqs(7), ierr)
call MPI_Recv_init(Array, 1, halo_5_yup_recv(size5), ip_yup, 3 + tag_offset, comm, &
& reqs(8), ierr)
! Start send and receive in t direction
call MPI_Cart_Shift(comm, 2, 1, ip_tdn, ip_tup, ierr)
call MPI_Send_init(Array, 1, halo_5_tup_send(size5), ip_tup, 4 + tag_offset, comm, &
& reqs(9), ierr)
call MPI_Recv_init(Array, 1, halo_5_tdn_recv(size5), ip_tdn, 4 + tag_offset, comm, &
& reqs(10), ierr)
call MPI_Send_init(Array, 1, halo_5_tdn_send(size5), ip_tdn, 5 + tag_offset, comm, &
& reqs(11), ierr)
call MPI_Recv_init(Array, 1, halo_5_tup_recv(size5), ip_tup, 5 + tag_offset, comm, &
& reqs(12), ierr)
!
end subroutine init_halo_update_5
subroutine start_halo_update_5(size5, Array, tag, reqs)
!
integer, intent(in) :: size5, tag
complex(dp), intent(inout) :: Array(kthird, 0:ksizex_l + 1, 0:ksizey_l + 1, &
& 0:ksizet_l + 1, size5)
integer, intent(out) :: reqs(12)
integer :: ip_xup, ip_xdn, ip_yup, ip_ydn, ip_tup, ip_tdn
integer :: tag_offset
integer :: ierr
tag_offset = 6*tag
! Start send and receive in x direction
call MPI_Cart_Shift(comm, 0, 1, ip_xdn, ip_xup, ierr)
call MPI_Isend(Array, 1, halo_5_xup_send(size5), ip_xup, 0 + tag_offset, comm, &
& reqs(1), ierr)
call MPI_Irecv(Array, 1, halo_5_xdn_recv(size5), ip_xdn, 0 + tag_offset, comm, &
& reqs(2), ierr)
call MPI_Isend(Array, 1, halo_5_xdn_send(size5), ip_xdn, 1 + tag_offset, comm, &
& reqs(3), ierr)
call MPI_Irecv(Array, 1, halo_5_xup_recv(size5), ip_xup, 1 + tag_offset, comm, &
& reqs(4), ierr)
! Start send and receive in y direction
call MPI_Cart_Shift(comm, 1, 1, ip_ydn, ip_yup, ierr)
call MPI_Isend(Array, 1, halo_5_yup_send(size5), ip_yup, 2 + tag_offset, comm, &
& reqs(5), ierr)
call MPI_Irecv(Array, 1, halo_5_ydn_recv(size5), ip_ydn, 2 + tag_offset, comm, &
& reqs(6), ierr)
call MPI_Isend(Array, 1, halo_5_ydn_send(size5), ip_ydn, 3 + tag_offset, comm, &
& reqs(7), ierr)
call MPI_Irecv(Array, 1, halo_5_yup_recv(size5), ip_yup, 3 + tag_offset, comm, &
& reqs(8), ierr)
! Start send and receive in t direction
call MPI_Cart_Shift(comm, 2, 1, ip_tdn, ip_tup, ierr)
call MPI_Isend(Array, 1, halo_5_tup_send(size5), ip_tup, 4 + tag_offset, comm, &
& reqs(9), ierr)
call MPI_Irecv(Array, 1, halo_5_tdn_recv(size5), ip_tdn, 4 + tag_offset, comm, &
& reqs(10), ierr)
call MPI_Isend(Array, 1, halo_5_tdn_send(size5), ip_tdn, 5 + tag_offset, comm, &
& reqs(11), ierr)
call MPI_Irecv(Array, 1, halo_5_tup_recv(size5), ip_tup, 5 + tag_offset, comm, &
& reqs(12), ierr)
end subroutine start_halo_update_5
end module comms5