The Elk Code
sphcrd.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: sphcrd
8 ! !INTERFACE:
9 pure subroutine sphcrd(v,r,tp)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! v : input vector (in,real(3))
12 ! r : length of v (out,real)
13 ! tp : (theta, phi) coordinates (out,real(2))
14 ! !DESCRIPTION:
15 ! Returns the spherical coordinates $(r,\theta,\phi)$ of a vector
16 ! $$ {\bf v}=(r\sin(\theta)\cos(\phi), r\sin(\theta)\sin(\phi),
17 ! r\cos(\theta)). $$
18 !
19 ! !REVISION HISTORY:
20 ! Created October 2002 (JKD)
21 !EOP
22 !BOC
23 implicit none
24 ! arguments
25 real(8), intent(in) :: v(3)
26 real(8), intent(out) :: r,tp(2)
27 ! local variables
28 real(8), parameter :: eps=1.d-14
29 real(8) t1
30 r=sqrt(v(1)**2+v(2)**2+v(3)**2)
31 if (r > eps) then
32  t1=v(3)/r
33  if (t1 >= 1.d0) then
34  tp(1)=0.d0
35  else if (t1 <= -1.d0) then
36  tp(1)=3.1415926535897932385d0
37  else
38  tp(1)=acos(t1)
39  end if
40  if ((abs(v(1)) > eps).or.(abs(v(2)) > eps)) then
41  tp(2)=atan2(v(2),v(1))
42  else
43  tp(2)=0.d0
44  end if
45 else
46  tp(1)=0.d0
47  tp(2)=0.d0
48 end if
49 end subroutine
50 !EOC
51 
pure subroutine sphcrd(v, r, tp)
Definition: sphcrd.f90:10