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