-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathcim.f90
executable file
·273 lines (214 loc) · 9.39 KB
/
cim.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
module mod_cim
use mod_cons
use mod_var
use mod_global
use mod_bou
use mod_conserve
use mod_thomas
use mod_triperiodic
implicit none
private
public updatecim
contains
!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!1
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!I have a question: how to consider the boundary condition regariding the icing near wall
! c range now is 0 to 1, 1 means full water and 0 means full ice
subroutine updatecim
!implicit none
integer i,j,k
real:: advective,lap_c,c_star,temp,source,pterm
real:: t_step,num_s
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!for 2D pencil transfer
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
real, dimension(itot/dims(1),jtot/dims(2),ktot) :: AA_z,BB_z,CC_z,RHS_z !I think maybe it is not essential to be pointers
real, dimension(itot/dims(1),jtot,ktot/dims(2)) :: AA_y,BB_y,CC_y,RHS_y
!real, dimension(itot,jtot/dims(1),ktot/dims(2)) :: cox1,cox2,cox3,RHS_x
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
t_step=1.0 !we use two sub steps
num_s=1e-30
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
!first, let us calculate the z direction
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
!calculate RHS
do k=1,kmax
do j=1,jmax
do i=1,imax
advective = (m_x_phi_c(i,j,k)-m_x_phi_c(i-1,j,k))/dx + &
(m_y_phi_c(i,j,k)-m_y_phi_c(i,j-1,k))/dy + &
(m_z_phi_c(i,j,k)-m_z_phi_c(i,j,k-1))/dz
!advective =0.
!if (k==kmax) then
! advective = (m_x_phi_c(i,j,k)-m_x_phi_c(i-1,j,k))/dx + &
! (m_y_phi_c(i,j,k)-m_y_phi_c(i,j-1,k))/dy + &
! (m_z_phi_c(i,j,k-1)-m_z_phi_c(i,j,k-2))/dz
!endif
!use phi_old, because phi has already been updated
lap_c= -lamda_c *( ( (PFM_c(i+1,j,k)-PFM_c(i,j,k))/dx*(PFM_phi(i+1,j,k)+PFM_phi(i,j,k))/2.0 - &
(PFM_c(i,j,k)-PFM_c(i-1,j,k))/dx*(PFM_phi(i,j,k)+PFM_phi(i-1,j,k))/2.0 )/dx +&
( (PFM_c(i,j+1,k)-PFM_c(i,j,k))/dy*(PFM_phi(i,j+1,k)+PFM_phi(i,j,k))/2.0 - &
(PFM_c(i,j,k)-PFM_c(i,j-1,k))/dy*(PFM_phi(i,j,k)+PFM_phi(i,j-1,k))/2.0 )/dy &
)
c_star=2.0*PFM_c(i,j,k)**2*(3-4.0*PFM_c(i,j,k))*lamda_c/PFM_l_c**2*PFM_phi(i,j,k)
if (Tnew(i,j,k)>t_melt .and. PFM_c(i,j,k)==0.) then !here the temperature is in degree, not Kelvin
pterm=1.0
!pterm=0.0
else if (Tnew(i,j,k)<t_melt .and. PFM_c(i,j,k)==1 ) then
pterm=1.0
!pterm=0.0
else
pterm=-30.0*PFM_c(i,j,k)**2*(3.0*PFM_c(i,j,k)**2-4.0*PFM_c(i,j,k)+1.0)
!pterm=30.0*PFM_c(i,j,k)**2*(3.0*PFM_c(i,j,k)**2-4.0*PFM_c(i,j,k)+1.0)
endif
temp=rho_2*latent/t_melt*PFM_phi(i,j,k)*(t_melt-Tnew(i,j,k)) * pterm
source = Phi_c_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 )
RHS_z(i,j,k)=-advective-mobility_c*(lap_c+c_star+temp)+source+Phi_c_old(i,j,k)/dt*t_step
debug1(i,j,k)=RHS_z(i,j,k)
enddo
enddo
enddo
!calculate the coefficient of the matrix
!aa is coefficient for k-1, bb for k, cc for k+1
do k=1,kmax
do j=1,jmax
do i=1,imax !aa is coefficient for k-1, bb for k, cc for k+1
AA_z(i,j,k)=-mobility_c*lamda_c*(PFM_phi(i,j,k)+PFM_phi(i,j,k-1)+2.0*num_s)/2.0/dz2
BB_z(i,j,k)=(PFM_phi(i,j,k)+num_s)/dt*t_step &
+mobility_c*lamda_c*(PFM_phi(i,j,k+1)+2.0*PFM_phi(i,j,k)+PFM_phi(i,j,k-1)+4.0*num_s)/2.0/dz2 &
+mobility_c*lamda_c/PFM_l_c**2*(PFM_phi(i,j,k)+num_s)*(12.0*PFM_c(i,j,k)**2-12.0*PFM_c(i,j,k)+2.0) &
+mobility_c*rho_2*latent/t_melt*(PFM_phi(i,j,k)+num_s)*(t_melt-Tnew(i,j,k)) &
*60.0*PFM_c(i,j,k)*(2.0*PFM_c(i,j,k)-1.0)*(PFM_c(i,j,k)-1.0)
CC_z(i,j,k)=-mobility_c*lamda_c*(PFM_phi(i,j,k+1)+PFM_phi(i,j,k)+2.0*num_s)/2.0/dz2
debug2(i,j,k)=AA_z(i,j,k)
debug3(i,j,k)=BB_z(i,j,k)
debug4(i,j,k)=CC_z(i,j,k)
enddo
enddo
enddo
!revise the coefficient of the matrix according to the boundary condition
!aa is coefficient for k-1, bb for k, cc for k+1
do j=1,jmax
do i=1,imax
!for Neumann at top, just correct bb, as aa and cc are not used in thomas
BB_z(i,j,kmax)=BB_z(i,j,kmax)+CC_z(i,j,kmax)
!for Neumann at bottom, the RHS need to add the bottom value
BB_z(i,j,1) =BB_z(i,j,1) +AA_z(i,j,1)
enddo
enddo
!then calculate the C value
do j=1,jmax
do i=1,imax
call thomas(RHS_z(i,j,1:ktot),AA_z(i,j,1:ktot),BB_z(i,j,1:ktot),CC_z(i,j,1:ktot),ktot)
enddo
enddo
! then transfer the value from RHS to C
do k=1,kmax
do j=1,jmax
do i=1,imax
PFM_c(i,j,k)=RHS_z(i,j,k)
debug5(i,j,k)=PFM_c(i,j,k)
enddo
enddo
enddo
call boundc(PFM_c)
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
!then, let us calculate the y direction
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
!calculate RHS
do k=1,kmax
do j=1,jmax
do i=1,imax
! in fact, this advective does not need to computed again, should save it
advective = (m_x_phi_c(i,j,k)-m_x_phi_c(i-1,j,k))/dx + &
(m_y_phi_c(i,j,k)-m_y_phi_c(i,j-1,k))/dy + &
(m_z_phi_c(i,j,k)-m_z_phi_c(i,j,k-1))/dz
!use phi_old, because phi has already been updated
lap_c= -lamda_c *( ( (PFM_c(i+1,j,k)-PFM_c(i,j,k))/dx*(PFM_phi(i+1,j,k)+PFM_phi(i,j,k))/2.0 - &
(PFM_c(i,j,k)-PFM_c(i-1,j,k))/dx*(PFM_phi(i,j,k)+PFM_phi(i-1,j,k))/2.0 )/dx +&
( (PFM_c(i,j,k+1)-PFM_c(i,j,k))/dz*(PFM_phi(i,j,k+1)+PFM_phi(i,j,k))/2.0 - &
(PFM_c(i,j,k)-PFM_c(i,j,k-1))/dz*(PFM_phi(i,j,k)+PFM_phi(i,j,k-1))/2.0 )/dz &
)
c_star=2.0*PFM_c(i,j,k)**2*(3-4.0*PFM_c(i,j,k))*lamda_c/PFM_l_c**2*PFM_phi(i,j,k)
if (Tnew(i,j,k)>t_melt .and. PFM_c(i,j,k)==0.) then !here the temperature is in degree, not Kelvin
pterm=1.0
else if (Tnew(i,j,k)<t_melt .and. PFM_c(i,j,k)==1 ) then
pterm=1.0
else
pterm=-30.0*PFM_c(i,j,k)**2*(3.0*PFM_c(i,j,k)**2-4.0*PFM_c(i,j,k)+1.0)
endif
temp=rho_2*latent/t_melt*PFM_phi(i,j,k)*(t_melt-Tnew(i,j,k)) * pterm
source = Phi_c_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 )
RHS_z(i,j,k)=-advective-mobility_c*(lap_c+c_star+temp)+source+ Phi_c_old(i,j,k)/dt*t_step
enddo
enddo
enddo
!calculate the coefficient of the matrix
!aa is coefficient for j-1, bb for j, cc for j+1
do k=1,kmax
do j=1,jtot/dims(2)
do i=1,imax
AA_z(i,j,k)=-mobility_c*lamda_c*(PFM_phi(i,j,k)+PFM_phi(i,j-1,k)+2.0*num_s)/2.0/dy2
BB_z(i,j,k)=(PFM_phi(i,j,k)+num_s)/dt*t_step &
+mobility_c*lamda_c*(PFM_phi(i,j+1,k)+2.0*PFM_phi(i,j,k)+PFM_phi(i,j-1,k)+4.0*num_s)/2.0/dy2 &
+mobility_c*lamda_c/PFM_l_c**2*(PFM_phi(i,j,k)+num_s)*(12.0*PFM_c(i,j,k)**2-12.0*PFM_c(i,j,k)+2.0) &
+mobility_c*rho_2*latent/t_melt*(PFM_phi(i,j,k)+num_s)*(t_melt-Tnew(i,j,k)) &
!*(-60.0)*PFM_c(i,j,k)*(2.0*PFM_c(i,j,k)-1.0)*(PFM_c(i,j,k)-1.0)
*60.0*PFM_c(i,j,k)*(2.0*PFM_c(i,j,k)-1.0)*(PFM_c(i,j,k)-1.0)
CC_z(i,j,k)=-mobility_c*lamda_c*(PFM_phi(i,j+1,k)+PFM_phi(i,j,k)+2.0*num_s)/2.0/dy2
enddo
enddo
enddo
!for periodic boundary condition, no additional revision is needed for the coefficient, I will do it with triperiodic
! then I need to do 2D pencil
! then transfer the y and z direction
AA_y=0.
BB_y=0.
CC_y=0.
RHS_y=0.
call transpose_z_to_y(AA_z, AA_y)
call transpose_z_to_y(BB_z, BB_y)
call transpose_z_to_y(CC_z, CC_y)
call transpose_z_to_y(RHS_z, RHS_y)
! then calculate along y direction as implicit, Z as explicit
do i=1,imax !here the sequnce can be changed
do k=1,ktot/dims(2)
call triperiodic(RHS_y(i,1:jtot,k),AA_y(i,1:jtot,k),BB_y(i,1:jtot,k),CC_y(i,1:jtot,k),jtot)
enddo
enddo
!then transpose back
call transpose_y_to_z(RHS_y, RHS_z)
! then transfer the value from RHS to C
do k=1,kmax
do j=1,jmax
do i=1,imax
PFM_c(i,j,k)=RHS_z(i,j,k)
debug3(i,j,k)=RHS_z(i,j,k)
Phi_c(i,j,k)=PFM_c(i,j,k)*(PFM_phi(i,j,k)+1e-30) !this is right, or temperature will NAN
debug4(i,j,k)=Phi_c(i,j,k)
enddo
enddo
enddo
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
! at last, let us do the mapping
!$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$
do k=1,kmax
do j=1,jmax
do i=1,imax !here use new phi value
if ( (abs(PFM_phi(i,j,k)-Phi_c(i,j,k))<1e-12) .or. (PFM_phi(i,j,k)<1e-15) ) then
!if ( (abs(PFM_phi(i,j,k)-Phi_c(i,j,k))<1e-15) .or. (PFM_phi(i,j,k)<1e-15) ) then
PFM_c(i,j,k)=1.0
else
PFM_c(i,j,k)=MIN(MAX(PFM_c(i,j,k), 0.0),1.0)
endif
debug6(i,j,k)=PFM_c(i,j,k)
enddo
enddo
enddo
return
end subroutine updatecim
end module mod_cim