-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathphi_im.f90
executable file
·409 lines (340 loc) · 10.6 KB
/
phi_im.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
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
module mod_phi_im
use mod_var
use mod_global
use mod_bou
use mod_cons
use mod_chem
use decomp_2d
implicit none
private init_bigphi, solver2ddd
public updatephi_im
contains
!
subroutine updatephi_im
implicit none
integer i,j,k
!real :: factor1,factor2
real::advectivex,advectivey, advectivez,advective,source,phi_lap
!real :: t11,t22,t33,t44,t55
!
real :: middle, SSSS,aaaa,xishu
real, target, dimension(-2:i1+2,-2:j1+2,-2:k1+2) :: Q,dPFM
real, pointer, dimension(:,:,:) :: qqq
real, dimension(-2:i1+2,-2:j1+2,-2:k1+2) :: g_prime
real :: xst(itot), yst(jtot)
real, dimension(imax,jmax) :: xyst
!real :: bb
real(mytype) :: si(itot+15), sj(jtot+15)
real, dimension(kmax) :: a0,b0,c0
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
middle=3.0*sqrt(2.)*PFM_Sigma_ref*PFM_l !this is the gamma1 in the paper
SSSS=1.0*PFM_l**2*sqrt(4.0*1.0/mobility/middle/dt) !this is the S in the paper
aaaa=-SSSS/2.0/PFM_l**2*( 1+sqrt(1-4.0*1.0/middle/mobility/dt*PFM_l**4/SSSS**2) ) !this is the alpha in the paper
xst=0
yst=0
xyst=0
a0=0.
b0=0.
c0=0.
si=0
sj=0
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! here start my own code
xishu= -(aaaa+SSSS/PFM_l**2) !it is negative
call init_bigphi(xishu,xst,yst,xyst,si,sj,a0,b0,c0)
qqq => Q(1:imax,1:jmax,1:kmax)
do k=0,k1
do j=0,j1
do i=0,i1
!do k=-2,k1+2
! do j=-2,j1+2
! do i=-2,i1+2
g_prime(i,j,k)=(2.0*PFM_phi(i,j,k)*(PFM_phi(i,j,k)-1)*(2.0*PFM_phi(i,j,k)-1)-SSSS*PFM_phi(i,j,k))/PFM_l**2 !this has an additional term
enddo
enddo
enddo
do k=1,kmax
do j=1,jmax
do i=1,imax
source = PFM_phi(i,j,k)*( (unew(i,j,k)-unew(i-1,j,k))/dx + &
(vnew(i,j,k)-vnew(i,j-1,k))/dy + (wnew(i,j,k)-wnew(i,j,k-1))/dz )
advectivex=( unew(i,j,k)*(PFM_phi(i+1,j,k)+PFM_phi(i,j,k))/2.0-unew(i-1,j,k)*(PFM_phi(i,j,k)+PFM_phi(i-1,j,k))/2.0 )/dx !conservative
advectivey=( vnew(i,j,k)*(PFM_phi(i,j+1,k)+PFM_phi(i,j,k))/2.0-vnew(i,j-1,k)*(PFM_phi(i,j,k)+PFM_phi(i,j-1,k))/2.0 )/dy
advectivez=( wnew(i,j,k)*(PFM_phi(i,j,k+1)+PFM_phi(i,j,k))/2.0-wnew(i,j,k-1)*(PFM_phi(i,j,k)+PFM_phi(i,j,k-1))/2.0 )/dz
advective=advectivex+advectivey+advectivez
phi_lap= (g_prime(i+1,j,k)+g_prime(i-1,j,k)-2.0*g_prime(i,j,k))/dx2 + &
(g_prime(i,j+1,k)+g_prime(i,j-1,k)-2.0*g_prime(i,j,k))/dy2 + &
(g_prime(i,j,k+1)+g_prime(i,j,k-1)-2.0*g_prime(i,j,k))/dz2
Q(i,j,k)=( PFM_phi(i,j,k)/dt-advective+source )/middle/mobility + phi_lap
!debug1(i,j,k)=Q(i,j,k)
enddo
enddo
enddo
call solver2ddd(qqq,xyst,si,sj,a0,b0,c0) !then here we get the value of bigphi
!call boundChem(Q) !here Q has the value of bigphi, this is not essential
!do k=1,kmax
! do j=1,jmax
! do i=1,imax
! debug2(i,j,k)=(Q(i-1,j,k)+Q(i+1,j,k)-2.0*Q(i,j,k))/dx2+&
! (Q(i,j-1,k)+Q(i,j+1,k)-2.0*Q(i,j,k))/dy2+&
! (Q(i,j,k-1)+Q(i,j,k+1)-2.0*Q(i,j,k))/dz2+&
! xishu*Q(i,j,k)-debug1(i,j,k)
! debug3(i,j,k)=Q(i,j,k)
! enddo
! enddo
!enddo
xishu= aaaa !it is negative
call init_bigphi(xishu,xst,yst,xyst,si,sj,a0,b0,c0)
call solver2ddd(qqq,xyst,si,sj,a0,b0,c0) !then here we get the value of bigphi
!call boundChem(Q) !here Q has the value of bigphi, this is not essential
!unless different boundary condition is used, not the Neumann scheme
!do k=1,kmax
! do j=1,jmax
! do i=1,imax
! debug4(i,j,k)=(Q(i-1,j,k)+Q(i+1,j,k)-2.0*Q(i,j,k))/dx2+&
! (Q(i,j-1,k)+Q(i,j+1,k)-2.0*Q(i,j,k))/dy2+&
! (Q(i,j,k-1)+Q(i,j,k+1)-2.0*Q(i,j,k))/dz2+&
! xishu*Q(i,j,k)-debug3(i,j,k)
! debug5(i,j,k)=Q(i,j,k)
! enddo
! enddo
!enddo
do k=1,kmax
do j=1,jmax
do i=1,imax
PFM_phi(i,j,k)=Q(i,j,k)
!debug1(i,j,k)=PFM_phi(i,j,k)
enddo
enddo
enddo
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!calculate the overall residual
!call boundPFM(PFM_phi,dPFM_bound_old,dPFM_bound_old_top,unew,vnew,wnew)
!call ChemicalPotential_im
!do k=1,kmax
! do j=1,jmax
! do i=1,imax
! dPFM(i,j,k)=mobility*middle*SSSS/PFM_l**2*(PFM_phi(i,j,k)-PFM_phi_old(i,j,k))
!! t11 = 2.0*PFM_phi_old(i,j,k)*(PFM_phi_old(i,j,k)-1)*(2.0*PFM_phi_old(i,j,k)-1)
!! t22 = (PFM_phi(i+1,j,k)-2.*PFM_phi(i,j,k)+PFM_phi(i-1,j,k))/dx2
!! t33 = (PFM_phi(i,j+1,k)-2.*PFM_phi(i,j,k)+PFM_phi(i,j-1,k))/dy2
!! t44 = (PFM_phi(i,j,k+1)-2.*PFM_phi(i,j,k)+PFM_phi(i,j,k-1))/dz2
!! t55= t22+t33+t44
!! Chem_Pot(i,j,k) =3.0*sqrt(2.0)*(PFM_Sigma_ref/PFM_l*t11-PFM_Sigma_ref*PFM_l*t55)
! enddo
! enddo
!enddo
!call boundChem(dPFM)
!call boundChem(chem_pot)
!do k=1,kmax
! do j=1,jmax
! do i=1,imax
! source = PFM_phi_old(i,j,k)*( (unew(i,j,k)-unew(i-1,j,k))/dx + &
! (vnew(i,j,k)-vnew(i,j-1,k))/dy + (wnew(i,j,k)-wnew(i,j,k-1))/dz )
!
! advectivex=( unew(i,j,k)*(PFM_phi_old(i+1,j,k)+PFM_phi_old(i,j,k))/2.0-&
! unew(i-1,j,k)*(PFM_phi_old(i,j,k)+PFM_phi_old(i-1,j,k))/2.0 )/dx !conservative
! advectivey=( vnew(i,j,k)*(PFM_phi_old(i,j+1,k)+PFM_phi_old(i,j,k))/2.0-&
! vnew(i,j-1,k)*(PFM_phi_old(i,j,k)+PFM_phi_old(i,j-1,k))/2.0 )/dy
! advectivez=( wnew(i,j,k)*(PFM_phi_old(i,j,k+1)+PFM_phi_old(i,j,k))/2.0-&
! wnew(i,j,k-1)*(PFM_phi_old(i,j,k)+PFM_phi_old(i,j,k-1))/2.0 )/dz
!
! advective=advectivex+advectivey+advectivez
!
! debug6(i,j,k)=(PFM_phi(i,j,k)-PFM_phi_old(i,j,k))/dt+&
! advective-source- mobility*&
! ((Chem_pot(i-1,j,k)+Q(i+1,j,k)-2.0*Chem_pot(i,j,k))/dx2+&
! (Chem_pot(i,j-1,k)+Q(i,j+1,k)-2.0*Chem_pot(i,j,k))/dy2+&
! (Chem_pot(i,j,k-1)+Q(i,j,k+1)-2.0*Chem_pot(i,j,k))/dz2 )
! !((dPFM(i-1,j,k)+dPFM(i+1,j,k)-2.0*dPFM(i,j,k))/dx2+&
! !(dPFM(i,j-1,k)+dPFM(i,j+1,k)-2.0*dPFM(i,j,k))/dy2+&
! !(dPFM(i,j,k-1)+dPFM(i,j,k+1)-2.0*dPFM(i,j,k))/dz2 )
! enddo
! enddo
!enddo
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!update the boundary to n+1
!do j=1,jmax
! do i=1,imax
! Phi_wall_n1=phi_wall_n(i,j)+factor1*dPFM_bound_n(i,j)+factor2*dPFM_bound_old(i,j)
! Phi_wall_n1_top=phi_wall_n_top(i,j)+factor1*dPFM_bound_n_top(i,j)+factor2*dPFM_bound_old_top(i,j)
! if (abs(PFM_mu_f).gt.1e-8) then
! PFM_phi(i,j,0) = 2*Phi_wall_n1-PFM_phi(i,j,1)
! dPFM_bound_old(i,j)=dPFM_bound_n(i,j)
! endif
!
! enddo
!enddo
!
return
end subroutine updatephi_im
!!&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
!!&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
!
subroutine init_bigphi(xishu,xst,yst,xyst,si,sj,a0,b0,c0)
implicit none
integer :: i,j,k,iv,jv
real, dimension(itot),intent(inout) :: xst
real, dimension(jtot),intent(inout) :: yst
real, dimension(imax,jmax),intent(inout) :: xyst
real(mytype),dimension(itot+15),intent(inout) :: si
real(mytype),dimension(jtot+15),intent(inout) :: sj
real, dimension(kmax),intent(inout) :: a0,b0,c0
real :: xishu
!
! Generate tridiagonal matrix
!
!print*, xishu
do k=1,kmax
a0(k) = dzi*dzi
c0(k) = dzi*dzi
b0(k) = -(a0(k) + c0(k))
enddo
! Neumann boundary condition for correction big_phi in z-direction;
! For big_phi, I only need Neumann, but Q need to be changed as well, because this is no a
! simple Neumann
!
b0(1) = b0(1) + a0(1)
b0(kmax) = b0(kmax) + c0(kmax)
a0(1) = 0.
c0(kmax) = 0.
!
! set lookup tables.
!
call vrffti(itot,si)
call vrffti(jtot,sj)
!
! generate eigenvalues ( xrt and yrt ).
!
!
! x direction
!
do i=3,itot,2
xst(i-1) = -4.*dxi*dxi*(sin(float((i-1))*pi/(2.*itot)))**2
xst(i) = xst(i-1)
enddo
xst(1 ) = 0.
xst(itot) = -4.*dxi*dxi
!
! y direction
!
do j=3,jtot,2
yst(j-1) = -4.*dyi*dyi*(sin(float((j-1))*pi/(2.*jtot)))**2
yst(j) = yst(j-1)
enddo
yst(1 ) = 0.
yst(jtot) = -4.*dyi*dyi
!!Neumann-Neumann
!do j=1,jtot
! yst(j) = -4.*dyi*dyi*(sin(float((j-1))*pi/(2.*jtot)))**2
!enddo
do j=1,jmax
jv = j + zstart(2) - 1
do i=1,imax
iv = i + zstart(1) - 1
xyst(i,j) = xst(iv)+yst(jv)+xishu !this is the difference between Poisson and Helmholtz
enddo
enddo
!
return
end subroutine init_bigphi
!
subroutine solver2ddd(phiz,xyst,si,sj,a0,b0,c0)
use decomp_2d
implicit none
real, intent(inout), dimension(1:,1:,1:) :: phiz
real, dimension(itot/dims(1),jtot,ktot/dims(2)) :: phiy
real, dimension(itot,jtot/dims(1),ktot/dims(2)) :: phix
real :: bb
real :: z,d(imax,jmax,kmax)
real :: di(itot),dj(jtot)
integer :: i,j,k
!real, intent(in) :: xst(itot), yst(jtot)
real, dimension(imax,jmax),intent(in) :: xyst
real(mytype),intent(in) :: si(itot+15), sj(jtot+15)
real, dimension(kmax),intent(in) :: a0,b0,c0
!
call transpose_z_to_y(phiz,phiy)
call transpose_y_to_x(phiy,phix)
!
do k=1,xsize(3)
do j=1,xsize(2)
call vrfftf(1,itot,phix(1:itot,j,k),di,1,si)
enddo
enddo
!
!!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
!do i=1,itot
! phix(1:itot,j,k)=phix(1:itot,j,k)*0.5*(1+cos(pi*float(i-1)/itot))
!enddo
!!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
call transpose_x_to_y(phix,phiy)
!
do k=1,ysize(3)
do i=1,ysize(1)
call vrfftf(1,jtot,phiy(i,1:jtot,k),dj,1,sj)
enddo
enddo
!
!!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
!do i=1,itot
! phiy(i,1:jtot,k)=phiy(i,1:jtot,k)*0.5*(1+cos(pi*float(i-1)/jtot))
!enddo
!!^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
call transpose_y_to_z(phiy,phiz)
!
do j=1,jmax
do i=1,imax
z = 1./(b0(1)+xyst(i,j))
d(i,j,1) = c0(1)*z
phiz(i,j,1) = phiz(i,j,1)*z
enddo
enddo
do k=2,kmax-1
do j=1,jmax
do i=1,imax
bb = b0(k)+xyst(i,j)
z = 1./(bb-a0(k)*d(i,j,k-1))
d(i,j,k) = c0(k)*z
phiz(i,j,k) = (phiz(i,j,k)-a0(k)*phiz(i,j,k-1))*z
enddo
enddo
enddo
do j=1,jmax
do i=1,imax
bb = b0(kmax)+xyst(i,j)
z = bb-a0(kmax)*d(i,j,kmax-1)
if(z.ne.0.) then
phiz(i,j,kmax) = (phiz(i,j,kmax)-a0(kmax)*phiz(i,j,kmax-1))/z
else
phiz(i,j,kmax) =0.
endif
enddo
enddo
do k=kmax-1,1,-1
do j=1,jmax
do i=1,imax
phiz(i,j,k) = phiz(i,j,k)-d(i,j,k)*phiz(i,j,k+1)
enddo
enddo
enddo
call transpose_z_to_y(phiz,phiy)
do k=1,ysize(3)
do i=1,ysize(1)
call vrfftb(1,jtot,phiy(i,1:jtot,k),dj,1,sj)
enddo
enddo
!
call transpose_y_to_x(phiy,phix)
!
do k=1,xsize(3)
do j=1,xsize(2)
call vrfftb(1,itot,phix(1:itot,j,k),di,1,si)
enddo
enddo
!
call transpose_x_to_y(phix,phiy)
call transpose_y_to_z(phiy,phiz)
!
return
end subroutine solver2ddd
end module mod_phi_im