The Elk Code
symrf.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: symrf
8 ! !INTERFACE:
9 subroutine symrf(nrmt_,nrmti_,npmt_,ngridg_,ngtot_,ngvec_,nfgrz_,igfft_,igrzf_,&
10  ld,rfmt,rfir)
11 ! !USES:
12 use modmain
13 use modomp
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! nrmt_ : number of radial points for each species (in,integer(nspecies))
16 ! nrmti_ : number of radial points on the inner part (in,integer(nspecies))
17 ! npmt_ : total number of points in each muffin-tin (in,integer(nspecies))
18 ! ngridg_ : G-vector grid sizes (in,integer(3))
19 ! ngtot_ : total number of G-vectors (in,integer)
20 ! ngvec_ : number of G-vectors within cut-off (in,integer)
21 ! nfgrz_ : number of FFT elements for real-complex transforms (in,integer)
22 ! igfft_ : map from G-vector index to FFT array (in,integer(ngvec_))
23 ! igrzf_ : map from real-complex FFT index to G-point index
24 ! (in,integer(nfgrz_)
25 ! ld : leading dimension (in,integer)
26 ! rfmt : real muffin-tin function (inout,real(ld,natmtot))
27 ! rfir : real intersitial function (inout,real(ngtot_))
28 ! !DESCRIPTION:
29 ! Symmetrises a real scalar function defined over the entire unit cell using
30 ! the full set of crystal symmetries. In the muffin-tin of a particular atom
31 ! the spherical harmonic coefficients of every equivlent atom are rotated and
32 ! averaged. The interstitial part of the function is first Fourier transformed
33 ! to $G$-space, and then averaged over each symmetry by rotating the Fourier
34 ! coefficients and multiplying them by a phase factor corresponding to the
35 ! symmetry translation.
36 !
37 ! !REVISION HISTORY:
38 ! Created May 2007 (JKD)
39 !EOP
40 !BOC
41 implicit none
42 ! arguments
43 integer, intent(in) :: nrmt_(nspecies),nrmti_(nspecies),npmt_(nspecies)
44 integer, intent(in) :: ngridg_(3),ngtot_,ngvec_,nfgrz_
45 integer, intent(in) :: igfft_(ngvec_),igrzf_(nfgrz_),ld
46 real(8), intent(inout) :: rfmt(ld,natmtot),rfir(ngtot_)
47 ! local variables
48 integer nthd
49 call holdthd(2,nthd)
50 !$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
51 !$OMP NUM_THREADS(nthd)
52 !$OMP SECTION
53 call symrfmt(nrmt_,nrmti_,npmt_,ld,rfmt)
54 !$OMP SECTION
55 call symrfir(ngridg_,ngtot_,ngvec_,nfgrz_,igfft_,igrzf_,rfir)
56 !$OMP END PARALLEL SECTIONS
57 call freethd(nthd)
58 end subroutine
59 !EOC
60 
subroutine symrfmt(nrmt_, nrmti_, npmt_, ld, rfmt)
Definition: symrfmt.f90:7
subroutine symrfir(ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, rfir)
Definition: symrfir.f90:10
Definition: modomp.f90:6
subroutine symrf(nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rfmt, rfir)
Definition: symrf.f90:11
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78