The Elk Code
symrvf.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 General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: symrvf
8 ! !INTERFACE:
9 subroutine symrvf(tspin,tnc,nrmt_,nrmti_,npmt_,ngridg_,ngtot_,ngvec_,nfgrz_, &
10  igfft_,igrzf_,ld1,rvfmt,ld2,rvfir)
11 ! !USES:
12 use modmain
13 use modomp
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! tspin : .true. if spin rotations should be used (in,logical)
16 ! tnc : .true. if the vector field is non-collinear, otherwise it is
17 ! collinear along the z-axis (in,logical)
18 ! nrmt_ : number of radial points for each species (in,integer(nspecies))
19 ! nrmti_ : number of radial points on the inner part (in,integer(nspecies))
20 ! npmt_ : total number of points in each muffin-tin (in,integer(nspecies))
21 ! ngridg_ : G-vector grid sizes (in,integer(3))
22 ! ngtot_ : total number of G-vectors (in,integer)
23 ! ngvec_ : number of G-vectors within cut-off (in,integer)
24 ! nfgrz_ : number of FFT elements for real-complex transforms (in,integer)
25 ! igfft_ : map from G-vector index to FFT array (in,integer(ngvec_))
26 ! igrzf_ : map from real-complex FFT index to G-point index
27 ! (in,integer(nfgrz_)
28 ! ld1 : leading dimension (in,integer)
29 ! rvfmt : real muffin-tin vector field (in,real(ld1,natmtot,*))
30 ! ld2 : leading dimension (in,integer)
31 ! rvfir : real interstitial vector field (in,real(ld2,*))
32 ! !DESCRIPTION:
33 ! Symmetrises a vector field defined over the entire unit cell using the full
34 ! set of crystal symmetries. If a particular symmetry involves rotating atom
35 ! 1 into atom 2, then the spatial and spin rotations of that symmetry are
36 ! applied to the vector field in atom 2 (expressed in spherical harmonic
37 ! coefficients), which is then added to the field in atom 1. This is repeated
38 ! for all symmetry operations. The fully symmetrised field in atom 1 is then
39 ! rotated and copied to atom 2. Symmetrisation of the interstitial part of the
40 ! field is performed by {\tt symrvfir}. See also {\tt symrfmt} and
41 ! {\tt findsym}.
42 !
43 ! !REVISION HISTORY:
44 ! Created May 2007 (JKD)
45 ! Fixed problem with improper rotations, February 2008 (L. Nordstrom,
46 ! F. Bultmark and F. Cricchio)
47 !EOP
48 !BOC
49 implicit none
50 ! arguments
51 logical, intent(in) :: tspin,tnc
52 integer, intent(in) :: nrmt_(nspecies),nrmti_(nspecies),npmt_(nspecies)
53 integer, intent(in) :: ngridg_(3),ngtot_,ngvec_,nfgrz_
54 integer, intent(in) :: igfft_(ngvec_),igrzf_(nfgrz_)
55 integer, intent(in) :: ld1
56 real(8), intent(inout) :: rvfmt(ld1,natmtot,*)
57 integer, intent(in) :: ld2
58 real(8), intent(inout) :: rvfir(ld2,*)
59 ! local variables
60 integer nthd
61 call holdthd(2,nthd)
62 !$OMP PARALLEL SECTIONS DEFAULT(SHARED) &
63 !$OMP NUM_THREADS(nthd)
64 !$OMP SECTION
65 call symrvfmt(tspin,tnc,nrmt_,nrmti_,npmt_,ld1,rvfmt)
66 !$OMP SECTION
67 call symrvfir(tspin,tnc,ngridg_,ngtot_,ngvec_,nfgrz_,igfft_,igrzf_,ld2,rvfir)
68 !$OMP END PARALLEL SECTIONS
69 call freethd(nthd)
70 end subroutine
71 !EOC
72 
subroutine symrvfir(tspin, tnc, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rvfir)
Definition: symrvfir.f90:11
subroutine symrvfmt(tspin, tnc, nrmt_, nrmti_, npmt_, ld, rvfmt)
Definition: symrvfmt.f90:7
Definition: modomp.f90:6
subroutine symrvf(tspin, tnc, nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld1, rvfmt, ld2, rvfir)
Definition: symrvf.f90:11
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78