Newer
Older
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
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
! Field metadata for heat equation solver
module heat
use iso_fortran_env, only : REAL64
implicit none
integer, parameter :: dp = REAL64
real(dp), parameter :: DX = 0.01, DY = 0.01 ! Fixed grid spacing
type :: field
integer :: nx ! local dimension of the field
integer :: ny
integer :: nx_full ! global dimension of the field
integer :: ny_full
real(dp) :: dx
real(dp) :: dy
real(dp), dimension(:,:), allocatable :: data
end type field
type :: parallel_data
integer :: size
integer :: rank
integer :: nup, ndown, nleft, nright ! Ranks of neighbouring MPI tasks
integer :: comm
integer :: requests(8) ! Non-blocking communication handles
integer :: rowtype ! MPI Datatype for communication of
! rows
integer :: columntype ! MPI Datatype for communication of
! columns
integer :: subarraytype ! MPI Datatype for text I/O
integer :: restarttype ! MPI Datatype for restart I/O
integer :: filetype ! MPI Datatype for file view in
! restart I/O
end type parallel_data
contains
! Initialize the field type metadata
! Arguments:
! field0 (type(field)): input field
! nx, ny, dx, dy: field dimensions and spatial step size
subroutine set_field_dimensions(field0, nx, ny, parallel)
implicit none
type(field), intent(out) :: field0
integer, intent(in) :: nx, ny
type(parallel_data), intent(in) :: parallel
integer :: nx_local, ny_local
integer, dimension(2) :: dims, coords
logical :: periods(2)
integer :: ierr
call mpi_cart_get(parallel%comm, 2, dims, periods, coords, ierr)
nx_local = nx / dims(1)
ny_local = ny / dims(2)
field0%dx = DX
field0%dy = DY
field0%nx = nx_local
field0%ny = ny_local
field0%nx_full = nx
field0%ny_full = ny
end subroutine set_field_dimensions
subroutine parallel_setup(parallel, nx, ny)
use mpi
implicit none
type(parallel_data), intent(out) :: parallel
integer, intent(in), optional :: nx, ny
integer :: nx_local, ny_local
integer :: world_size
integer :: dims(2) = (/0, 0/)
logical :: periods(2) = (/.false., .false./)
integer :: coords(2)
integer, dimension(2) :: sizes, subsizes, offsets
integer :: ierr
call mpi_comm_size(MPI_COMM_WORLD, world_size, ierr)
! Set grid dimensions
call mpi_dims_create(world_size, 2, dims, ierr)
nx_local = nx / dims(1)
ny_local = ny / dims(2)
! Ensure that the grid is divisible to the MPI tasks
if (nx_local * dims(1) /= nx) then
write(*,*) 'Cannot divide grid evenly to processors in x-direction', &
& nx_local, dims(1), nx
call mpi_abort(MPI_COMM_WORLD, -2, ierr)
end if
if (ny_local * dims(2) /= ny) then
write(*,*) 'Cannot divide grid evenly to processors in y-direction', &
& ny_local, dims(2), ny
call mpi_abort(MPI_COMM_WORLD, -2, ierr)
end if
! Create cartesian communicator
call mpi_cart_create(MPI_COMM_WORLD, 2, dims, periods, .true., parallel%comm, ierr)
call mpi_cart_shift(parallel%comm, 0, 1, parallel%nup, parallel%ndown, ierr)
call mpi_cart_shift(parallel%comm, 1, 1, parallel%nleft, &
& parallel%nright, ierr)
call mpi_comm_size(parallel%comm, parallel%size, ierr)
call mpi_comm_rank(parallel%comm, parallel%rank, ierr)
if (parallel%rank == 0) then
write(*,'(A, I3, A3, I3)') 'Using domain decomposition', dims(1), 'x', &
& dims(2)
write(*,'(A, I5, A3, I5)') 'Local domain size', nx_local, 'x', ny_local
end if
! Create datatypes for halo exchange
call mpi_type_vector(ny_local + 2, 1, nx_local + 2, MPI_DOUBLE_PRECISION, &
& parallel%rowtype, ierr)
call mpi_type_contiguous(nx_local + 2, MPI_DOUBLE_PRECISION, &
& parallel%columntype, ierr)
call mpi_type_commit(parallel%rowtype, ierr)
call mpi_type_commit(parallel%columntype, ierr)
! Create datatype for subblock needed in I/O
! Rank 0 uses datatype for receiving data into full array while
! other ranks use datatype for sending the inner part of array
sizes(1) = nx_local + 2
sizes(2) = ny_local + 2
subsizes(1) = nx_local
subsizes(2) = ny_local
offsets(1) = 1
offsets(2) = 1
if (parallel%rank == 0) then
sizes(1) = nx
sizes(2) = ny
offsets(1) = 0
offsets(2) = 0
end if
call mpi_type_create_subarray(2, sizes, subsizes, offsets, &
& MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, parallel%subarraytype, &
& ierr)
call mpi_type_commit(parallel%subarraytype, ierr)
! Create datatypes for restart I/O
! For boundary ranks also the ghost layer (boundary condition)
! is written
call mpi_cart_coords(parallel%comm, parallel%rank, 2, coords, ierr);
sizes(1) = nx + 2
sizes(2) = ny + 2
offsets(1) = 1 + coords(1) * nx_local
offsets(2) = 1 + coords(2) * ny_local
if (coords(1) == 0) then
offsets(1) = offsets(1) - 1
subsizes(1) = subsizes(1) + 1
end if
if (coords(1) == dims(1) - 1) then
subsizes(1) = subsizes(1) + 1;
end if
if (coords(2) == 0) then
offsets(2) = offsets(2) - 1
subsizes(2) = subsizes(2) + 1
end if
if (coords(2) == dims(2) - 1) then
subsizes(2) = subsizes(2) + 1
end if
call mpi_type_create_subarray(2, sizes, subsizes, offsets, &
& MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, parallel%filetype, ierr)
call mpi_type_commit(parallel%filetype, ierr)
sizes(1) = nx_local + 2
sizes(2) = ny_local + 2
offsets(1) = 1
offsets(2) = 1
if (coords(1) == 0) then
offsets(1) = 0
end if
if (coords(2) == 0) then
offsets(2) = 0
end if
call mpi_type_create_subarray(2, sizes, subsizes, offsets, &
& MPI_ORDER_FORTRAN, MPI_DOUBLE_PRECISION, parallel%restarttype, ierr)
call mpi_type_commit(parallel%restarttype, ierr)
end subroutine parallel_setup
end module heat