-
Notifications
You must be signed in to change notification settings - Fork 8
/
Copy pathcalcul_new_velocity.f90
245 lines (143 loc) · 5.4 KB
/
calcul_new_velocity.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
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
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
subroutine calcul_new_velocity()
use variables
implicit none
!-------------------------------------------------!
! Calculation of velocity field at t = dt*n+1 !
!-------------------------------------------------!
!!! In x direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
u1(i,j,k) = u_star(i,j,k) - dt*(p(i+1,j,k)-p(i,j,k)) / Dxs(i)
end do;end do
!$OMP END PARALLEL DO
end do
!!! In y direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
v1(i,j,k) = v_star(i,j,k) - dt*(p(i,j+1,k)-p(i,j,k)) / Dys(j)
end do;end do
!$OMP END PARALLEL DO
end do
!!! In w direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
w1(i,j,k) = w_star(i,j,k) - dt*(p(i,j,k+1)-p(i,j,k)) / Dzs(k)
end do;end do
!$OMP END PARALLEL DO
end do
!-------------------------------------------------!
! Calculation of velocity field for DFIB !
!-------------------------------------------------!
!!! In x direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
u2(i,j,k) = ETA(i,j,k) *u_solid + (1-ETA(i,j,k)) * u1(i,j,k)
FX(i,j,k) = (u2(i,j,k) - u1(i,j,k)) / dt
end do;end do
!$OMP END PARALLEL DO
end do
!!! In y direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
v2(i,j,k) = ETA(i,j,k) *v_solid + (1-ETA(i,j,k)) * v1(i,j,k)
FY(i,j,k) = (v2(i,j,k) - v1(i,j,k)) / dt
end do;end do
!$OMP END PARALLEL DO
end do
!!! In w direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
w2(i,j,k) = ETA(i,j,k) *w_solid + (1-ETA(i,j,k)) * w1(i,j,k)
FZ(i,j,k) = (w2(i,j,k) - w1(i,j,k)) / dt
end do;end do
!$OMP END PARALLEL DO
end do
!write(*,*) myid, 'gcount = ', gcount(myid), 'igcount = ', igcount
!----------data collect among nodes----------!
icount = igcount*(nx+4)*(ny+4)
!Send my results back to the master
if(myid>master)then
itag = 350
call MPI_SEND( FX(-1,-1,istart), icount, MPI_REAL8, master, itag, MPI_COMM_WORLD, ierr )
itag = 360
call MPI_SEND( FY(-1,-1,istart), icount, MPI_REAL8, master, itag, MPI_COMM_WORLD, ierr )
itag = 370
call MPI_SEND( FZ(-1,-1,istart), icount, MPI_REAL8, master, itag, MPI_COMM_WORLD, ierr )
end if
!Wait to receive results from each task
if(myid==master)then
do i = 1, (nproc-1)
itag = 350
call MPI_RECV( FX(-1,-1,gstart(i)), icount, MPI_REAL8, i, itag, MPI_COMM_WORLD, status, ierr )
itag = 360
call MPI_RECV( FY(-1,-1,gstart(i)), icount, MPI_REAL8, i, itag, MPI_COMM_WORLD, status, ierr )
itag = 370
call MPI_RECV( FZ(-1,-1,gstart(i)), icount, MPI_REAL8, i, itag, MPI_COMM_WORLD, status, ierr )
end do
end if
!icount = (nx+4)*(ny+4)*(nz+4)
!call MPI_BARRIER(MPI_COMM_WORLD, ierr)
!call MPI_BCAST( FX(-1,-1,-1), icount, MPI_REAL8, master, MPI_COMM_WORLD, ierr )
!call MPI_BARRIER(MPI_COMM_WORLD, ierr)
!call MPI_BCAST( FY(-1,-1,-1), icount, MPI_REAL8, master, MPI_COMM_WORLD, ierr )
!call MPI_BARRIER(MPI_COMM_WORLD, ierr)
!call MPI_BCAST( FZ(-1,-1,-1), icount, MPI_REAL8, master, MPI_COMM_WORLD, ierr )
!call MPI_BARRIER(MPI_COMM_WORLD, ierr)
!----------data collect among nodes----------!
!----------------------------------------------------------------------------------------!
! Centering variables to get matrix (1:nx,1:ny,1:nz) for Tecplot !
!----------------------------------------------------------------------------------------!
!----------data transformation among nodes----------!
icount = (nx+4)*(ny+4)
itag = 380
call MPI_SENDRECV( w2(-1,-1,iend), icount, MPI_REAL8, r_nbr, itag, &
w2(-1,-1,istart-1), icount, MPI_REAL8, l_nbr, itag, MPI_COMM_WORLD, status, ierr )
call MPI_BARRIER(MPI_COMM_WORLD, ierr)
!----------data transformation among nodes----------!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
uc(i,j,k) = 0.5*(u2(i,j,k)+u2(i-1,j,k))
vc(i,j,k) = 0.5*(v2(i,j,k)+v2(i,j-1,k))
wc(i,j,k) = 0.5*(w2(i,j,k)+w2(i,j,k-1))
end do;end do
!$OMP END PARALLEL DO
end do
end subroutine calcul_new_velocity
subroutine Updating_velocity()
use variables
implicit none
!---------------------------------------------------------!
! loops to update the main u and v velocity fields !
!---------------------------------------------------------!
!!! In x direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
u(i,j,k) = u2(i,j,k)
end do;end do
!$OMP END PARALLEL DO
end do
!!! In y direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
v(i,j,k) = v2(i,j,k)
end do;end do
!$OMP END PARALLEL DO
end do
!!! In w direction !!!
do k=istart,iend
!$OMP PARALLEL DO PRIVATE(i)
do j=1,ny;do i=1,nx
w(i,j,k) = w2(i,j,k)
end do;end do
!$OMP END PARALLEL DO
end do
end subroutine Updating_velocity