-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathget_lambda.f90
193 lines (179 loc) · 5.9 KB
/
get_lambda.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
module getlambda
use get_neighbour_lists
implicit none
contains
!!------------------------------------------------------------------------
!! Get the value of lambda from a linear combination of two kernels
!!------------------------------------------------------------------------
subroutine get_lambda(lambda,ntot)
use dimen_mhd, only:ndim,ndimV
use debug, only:trace
use loguns, only:iprint
use kernels, only:interpolate_kernels2,radkern2
use linklist
use part
use setup_params
use timestep, only:time
!
!--define local variables
!
implicit none
integer, intent(in) :: ntot
real, dimension(:), intent(out) :: lambda
integer :: i,j,n
integer :: icell,iprev,nneigh
integer, dimension(ntot) :: listneigh ! neighbour list
integer :: idone
!
! (particle properties - local copies)
!
real :: rij,rij2
real :: hi,hi1,hj,hj1,h2,hi2,hj2
real :: hfacwabi,hfacwabj
real :: rho1i,rho1j,rho21i, rho21j
real :: rhoalt1i,rhoalt1j,rhoalt21i,rhoalt21j
real, dimension(ndim) :: dx
real, dimension(ndimv) :: dr
!
! (kernel quantities)
!
real :: q2,q2i,q2j
real :: wabi,wabj
real :: wabalti,wabaltj
real :: grkerni,grkernj,grkernalti,grkernaltj
real :: term,termalt,rhochecki,rhocheckalti
real, dimension(ndimV) :: gradunityi,gradunityalti
!
!--allow for tracing flow
!
if (trace) write(iprint,*) ' Entering subroutine get_lambda, ntot=',ntot
!
!--initialise quantities
!
!!!print*,'ntot = ',ntot, 'size= ',size(divBonrho)
do i=1,ntot
lambda(i) = 0.
enddo
listneigh = 0
!
!--Loop over all the link-list cells
!
loop_over_cells: do icell=1,ncellsloop ! step through all cells
!
!--get the list of neighbours for this cell
! (common to all particles in the cell)
!
call get_neighbour_list_partial(icell,listneigh,nneigh)
!
!--now loop over all particles in the current cell
!
i = ifirstincell(icell)
idone = -1 ! note density summation includes current particle
if (i.NE.-1) iprev = i
loop_over_cell_particles: do while (i.NE.-1) ! loop over home cell particles
! PRINT*,'Doing particle ',i,nneigh,' neighbours',hh(i)
idone = idone + 1
hi = hh(i)
hi1 = 1./hi
hi2 = hi*hi
hfacwabi = hi1**ndim
rho1i = 1./rho(i)
rho21i = rho1i*rho1i
rhoalt1i = 1./rhoalt(i)
rhoalt21i = rhoalt1i*rhoalt1i
rhochecki = 0.
rhocheckalti = 0.
gradunityi(:) = 0.
gradunityalti(:) = 0.
if (i.gt.npart) then
i = ll(i)
cycle loop_over_cell_particles
endif
!
!--for each particle in the current cell, loop over its neighbours
!
loop_over_neighbours: DO n = 1,nneigh
j = listneigh(n)
if (j.lt.0 .or. j.gt.ntotal) then
print*,i,j
print*,listneigh(1:nneigh)
print*,'n=',n,' listneigh(n) = ',listneigh(n)
endif
dx(:) = x(:,i) - x(:,j)
hj = hh(j)
hj1 = 1./hj
hj2 = hj*hj
hfacwabj = hj1**ndim
rho1j = 1./rho(j)
rho21j = rho1j*rho1j
rhoalt1j = 1./rhoalt(j)
rhoalt21j = rhoalt1j*rhoalt1j
rij2 = DOT_PRODUCT(dx,dx)
q2 = rij2/h2
q2i = rij2/hi2
q2j = rij2/hj2
!
!--do interaction if r/h < compact support size
! don't calculate interactions between ghost particles
!
if (((q2i.lt.radkern2).OR.(q2j.lt.radkern2)) &
.and. .not.(i.gt.npart.and.j.gt.npart)) then
if (rij2.gt.0.) then
rij = SQRT(rij2)
dr(1:ndim) = dx(1:ndim)/rij ! unit vector
else
dr(:) = 0.
endif
if (ndimV.gt.ndim) dr(ndim+1:ndimV) = 0.
!
!--interpolate from kernel table
! (use either average h or average kernel gradient)
!
! (using hi)
call interpolate_kernels2(q2i,wabi,wabalti,grkerni,grkernalti)
wabi = wabi*hfacwabi
wabalti = wabalti*hfacwabi
grkerni = grkerni*hfacwabi*hi1
grkernalti = grkernalti*hfacwabi*hi1
! (using hj)
call interpolate_kernels2(q2j,wabj,wabaltj,grkernj,grkernaltj)
wabj = wabj*hfacwabj
wabaltj = wabaltj*hfacwabj
grkernj = grkernj*hfacwabj*hj1
grkernaltj = grkernaltj*hfacwabj*hj1
!
!--calculate term
!
rhochecki = rhochecki + pmass(j)*wabi
rhocheckalti = rhocheckalti + pmass(j)*wabalti
!term = pmass(j)*(rho21i*grkerni + rho21j*grkernj)
!termalt = pmass(j)*(rhoalt21i*grkernalti + rhoalt21j*grkernaltj)
term = -pmass(j)*rho1i*dx(1)*dr(1)*grkerni
termalt = -pmass(j)*rhoalt1i*dx(1)*dr(1)*grkernalti
gradunityi(:) = gradunityi(:) + term !*dr(:)
gradunityalti(:) = gradunityalti(:) + termalt !*dr(:)
endif
enddo loop_over_neighbours
!
!--compute lambda for this particle
!
lambda(i) = 0.8 ! max(min((1. - gradunityalti(1))/(gradunityi(1) - gradunityalti(1)),1.),0.)
if (i.eq.50) then
print*,i,'h = ',hh(i),' lambda = ',lambda(i)
print*,i,'terms = ',gradunityi(1),gradunityalti(1),' rho = ',rho(i),rhoalt(i)
print*,i,'corrected term = ',lambda(i)*gradunityi(1) + (1. - lambda(i))*gradunityalti(1)
!print*,i,'corrected rho = ',lambda(i)*rho(i) + (1. - lambda(i))*rhoalt(i)
if (abs(rhoalt(i)-rhocheckalti).gt.1.e-10) then
print*,' error, rhoalt = ',rhoalt(i),' should be ',rhocheckalti
stop
endif
!print*,' rhoalt is ',rhoalt(i),' should be ',rhocheckalti
write(50,*) time,lambda(i)*rho(i) + (1. - lambda(i))*rhoalt(i),rho(i),rhoalt(i)
endif
iprev = i
if (iprev.NE.-1) i = ll(i) ! possibly should be only if (iprev.NE.-1)
enddo loop_over_cell_particles
enddo loop_over_cells
return
end subroutine get_lambda
end module getlambda