The Elk Code
cfmtconj.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2016 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 pure subroutine cfmtconj(nr,nri,np,cfmt)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: nr,nri,np
11 complex(4), intent(inout) :: cfmt(np)
12 ! local variables
13 integer i
14 ! automatic arrays
15 complex(4) cfmt1(np)
16 cfmt1(1:np)=cfmt(1:np)
17 call cflmnconj(lmaxi,nri,lmmaxi,cfmt1,cfmt)
18 i=lmmaxi*nri+1
19 call cflmnconj(lmaxo,nr-nri,lmmaxo,cfmt1(i),cfmt(i))
20 
21 contains
22 
23 !BOP
24 ! !ROUTINE: cflmnconj
25 ! !INTERFACE:
26 pure subroutine cflmnconj(lmax,n,ld,cflm1,cflm2)
27 ! !INPUT/OUTPUT PARAMETERS:
28 ! lmax : maximum angular momentum (in,integer)
29 ! n : number of functions to conjugate (in,integer)
30 ! ld : leading dimension (in,integer)
31 ! cflm1 : coefficients of input complex spherical harmonic expansion
32 ! (in,complex((lmax+1)**2)))
33 ! cflm2 : coefficients of output complex spherical harmonic expansion
34 ! (out,complex((lmax+1)**2)))
35 ! !DESCRIPTION:
36 ! Returns the complex conjugate of a function expanded in spherical harmonics.
37 ! In other words, given the input function coefficients $c_{lm}$, the routine
38 ! returns $c'_{lm}=(-1)^m c^*_{l-m}$ so that
39 ! $$ \sum_{lm}c'_{lm}Y_{lm}(\theta,\phi)=\left(\sum_{lm}c_{lm}Y_{lm}
40 ! (\theta,\phi)\right)^* $$
41 ! for all $(\theta,\phi)$.
42 !
43 ! !REVISION HISTORY:
44 ! Created April 2004 (JKD)
45 !EOP
46 !BOC
47 implicit none
48 ! arguments
49 integer, intent(in) :: lmax,n,ld
50 complex(4), intent(in) :: cflm1(ld,n)
51 complex(4), intent(out) :: cflm2(ld,n)
52 ! local variables
53 integer l,m,lm1,lm2
54 do l=0,lmax
55  lm1=l**2
56  lm2=(l+1)**2+1
57  do m=-l,-1
58  lm1=lm1+1
59  lm2=lm2-1
60  if (mod(m,2) == 0) then
61  cflm2(lm1,1:n)=conjg(cflm1(lm2,1:n))
62  cflm2(lm2,1:n)=conjg(cflm1(lm1,1:n))
63  else
64  cflm2(lm1,1:n)=-conjg(cflm1(lm2,1:n))
65  cflm2(lm2,1:n)=-conjg(cflm1(lm1,1:n))
66  end if
67  end do
68 ! m = 0 case
69  lm1=lm1+1
70  cflm2(lm1,1:n)=conjg(cflm1(lm1,1:n))
71 end do
72 end subroutine
73 !EOC
74 
75 end subroutine
76 
integer lmmaxo
Definition: modmain.f90:203
integer lmaxo
Definition: modmain.f90:201
pure subroutine cfmtconj(nr, nri, np, cfmt)
Definition: cfmtconj.f90:7
integer lmmaxi
Definition: modmain.f90:207
pure subroutine cflmnconj(lmax, n, ld, cflm1, cflm2)
Definition: cfmtconj.f90:27
integer lmaxi
Definition: modmain.f90:205