Skip to content

Commit

Permalink
Adding explanation
Browse files Browse the repository at this point in the history
  • Loading branch information
cpinte committed Nov 21, 2023
1 parent 31d5cde commit 066beeb
Show file tree
Hide file tree
Showing 2 changed files with 29 additions and 1 deletion.
14 changes: 14 additions & 0 deletions src/cylindrical_grid.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1201,6 +1201,20 @@ real(dp) function distance_to_closest_wall_cyl(id,icell,x,y,z) result(s)

if (l3D) then
! phi walls
!
! Distance between a line with equation a*x+by+c=0, ie y=-a/b*x-c/b and point (x0,y0) is:
!
! d = |a*x0+b*y0+c|/sqrt(a**2+b**2)
!
! In case of phi walls, we have y = x * tan(phi) so:
!
! a = -sin(phi)
! b = cos(phi)
! c = 0
!
! and therefore
!
! d = |x0*sin(phi) - y0*cos(theta)| / sqrt[cos^2 + sin^2]
s5 = abs(x*sin_phi_lim(k0) - y*cos_phi_lim(k0))
s6 = abs(x*sin_phi_lim(k0-1) - y*cos_phi_lim(k0-1))
s = min(s1,s2,s3,s4,s5,s6)
Expand Down
16 changes: 15 additions & 1 deletion src/spherical_grid.f90
Original file line number Diff line number Diff line change
Expand Up @@ -466,13 +466,27 @@ real(dp) function distance_to_closest_wall_sph(id,icell,x,y,z) result(s)
s1 = r_lim(ri0) - r
s2 = r - r_lim(ri0-1)

! theta walls
! theta walls (same method as phi)
z0 = abs(z)
s3 = abs(rcyl*w_lim(thetaj0) - z0*cos_phi_lim(thetaj0))
s4 = abs(rcyl*w_lim(thetaj0-1) - z0*cos_phi_lim(thetaj0-1))

if (l3D) then
! phi walls
!
! Distance between a line with equation a*x+by+c=0, ie y=-a/b*x-c/b and point (x0,y0) is:
!
! d = |a*x0+b*y0+c|/sqrt(a**2+b**2)
!
! In case of phi walls, we have y = x * tan(phi) so:
!
! a = -sin(phi)
! b = cos(phi)
! c = 0
!
! and therefore
!
! d = |x0*sin(phi) - y0*cos(theta)| / sqrt[cos^2 + sin^2]
s5 = abs(x*sin_phi_lim(k0) - y*cos_phi_lim(k0))
s6 = abs(x*sin_phi_lim(k0-1) - y*cos_phi_lim(k0-1))
s = min(s1,s2,s3,s4,s5,s6)
Expand Down

0 comments on commit 066beeb

Please sign in to comment.