LCOV - code coverage report
Current view: top level - cdn - MagMoments.F90 (source / functions) Hit Total Coverage
Test: FLEUR test coverage Lines: 76 88 86.4 %
Date: 2024-04-27 04:44:07 Functions: 3 4 75.0 %

          Line data    Source code
       1             : !--------------------------------------------------------------------------------
       2             : ! Copyright (c) 2022 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
       3             : ! This file is part of FLEUR and available as free software under the conditions
       4             : ! of the MIT license as expressed in the LICENSE file in more detail.
       5             : !--------------------------------------------------------------------------------
       6             : 
       7             : module m_magmoments
       8             :     implicit none
       9             :     private
      10             :     public:: spinMoments,orbMoments
      11             :     !! Module to encapsulate routines to calculate and output the spin- and orbital Moments
      12             :     contains
      13         196 :     SUBROUTINE spinMoments(input,atoms,noco,nococonv,moments,den)
      14             :         !! Calculate and output the spin moments
      15             :         !! Either the moments have to be given as calculated in cdnval (moments-argument)
      16             :         !! or the moments are calculated from the density 
      17             :         USE m_types
      18             :         USE m_constants
      19             :         USE m_xmlOutput
      20             :         USE m_intgr, ONLY : intgr3
      21             :         TYPE(t_input), INTENT(IN)           :: input
      22             :         TYPE(t_atoms), INTENT(IN)           :: atoms
      23             :         TYPE(t_noco), INTENT(IN)            :: noco
      24             :         TYPE(t_nococonv), INTENT(IN)        :: nococonv
      25             :         TYPE(t_moments),INTENT(IN),OPTIONAL :: moments
      26             :         TYPE(t_potden),INTENT(IN),OPTIONAL  :: den
      27             :      
      28             :         INTEGER                       :: iType
      29             :         REAL                          :: up,down,local_m(0:3),global_m(0:3),off1,off2
      30             :         COMPLEX                       :: off_diag
      31             :         LOGICAL                       :: l_offdiag
      32             :      
      33         196 :         if (input%jspins==1) return ! No spin-pol calculation!
      34         196 :         if (.not.present(moments).and..not.present(den)) return !No data provided
      35         196 :         if (present(moments)) call priv_print_spin_density_at_nucleus(input,atoms,moments)
      36             :         
      37             : 
      38             :         !Header of output
      39         196 :         write(ounit,*)
      40         196 :         write(oUnit,'(a,19x,a,19x,a,19x,a)') "Spin Magn. mom.  |","Global Frame","  | ","Local Frame"
      41         196 :         write(oUnit,*) "------------------------------------------------------------------------------------------------------------------------"
      42         196 :         IF(.NOT.(noco%l_noco.or.noco%l_soc)) THEN
      43          84 :            write(oUnit,'(a,5x,a,2(" | ",5(a,5x)))') "Atom "," m    ","mx   ","my   ","mz   ","alpha","beta ","mx   ","my   ","mz   ","alpha","beta "
      44             :         ELSE
      45         112 :            write(oUnit,'(a,5x,a,2(" | ",5(a,5x)))') "Atom ","|m|   ","mx   ","my   ","mz   ","alpha","beta ","mx   ","my   ","mz   ","alpha","beta "
      46             :         END IF   
      47         588 :         CALL openXMLElement('magneticMomentsInMTSpheres',(/'units'/),(/'muBohr'/))
      48             : 
      49             : 
      50         524 :         DO iType = 1, atoms%ntype
      51             :            !Collect the data
      52         328 :            if (present(moments)) THEN
      53           0 :               up=moments%chmom(iType,1)
      54           0 :               down=moments%chmom(iType,2)
      55         328 :            elseif(present(den)) THEN
      56         328 :               CALL intgr3(den%mt(:,0,iType,1),atoms%rmsh(:,iType),atoms%dx(iType),atoms%jri(iType),up)
      57         328 :               CALL intgr3(den%mt(:,0,iType,2),atoms%rmsh(:,iType),atoms%dx(iType),atoms%jri(iType),down)
      58         328 :               up=up*sfp_const;down=down*sfp_const
      59             :            endif   
      60         328 :            off_diag=0.0
      61         328 :            l_offdiag=.false.
      62         328 :            if (present(moments).and.noco%l_mperp) THEN
      63           0 :                off_diag=moments%qa21(itype)
      64           0 :                l_offdiag=.true.
      65         328 :            elseif (present(den)) THEN
      66         328 :               if (size(den%mt,4)==4) THEN
      67          37 :                  CALL intgr3(den%mt(:,0,iType,3),atoms%rmsh(:,iType),atoms%dx(iType),atoms%jri(iType),off1)
      68          37 :                  CALL intgr3(den%mt(:,0,iType,4),atoms%rmsh(:,iType),atoms%dx(iType),atoms%jri(iType),off2)
      69          37 :                  off_diag=cmplx(off1,off2)*sfp_const
      70          37 :                  l_offdiag=.true.
      71             :               endif
      72             :            endif   
      73             :            !Determine local vector and angles
      74         328 :            local_m=nococonv%denmat_to_mag(up,down,off_diag)
      75             :            !Rotate to global frame
      76         328 :            if (noco%l_noco) then
      77         150 :                 call nococonv%rotdenmat(itype, up, down, off_diag, toGlobal=.true.)
      78             :            else
      79         178 :                 call nococonv%rotdenmat(nococonv%phi,nococonv%theta,up, down, off_diag, toGlobal=.true.)
      80             :            endif
      81         328 :            global_m=nococonv%denmat_to_mag(up,down,off_diag)
      82             : 
      83         524 :            call priv_print_mt_moment(itype,noco%l_noco,l_offdiag,noco%l_soc,local_m(1:3),global_m(1:3),"smm")
      84             :         enddo   
      85             :      
      86         196 :         CALL closeXMLElement('magneticMomentsInMTSpheres')
      87         196 :         write(oUnit,*) "------------------------------------------------------------------------------------------------------------------------"
      88         196 :         write(oUnit,*)
      89         196 :         if (noco%l_soc) write(ounit,*) "SOC calculation: Global frame aligned to lattice"
      90             :      
      91             :     END SUBROUTINE 
      92             : 
      93         196 :     SUBROUTINE orbMoments(input,atoms,noco,nococonv,moments)
      94             :         !! Calculate and output the orbital moments
      95             :         !! Either the moments have to be given as calculated in cdnval (moments-argument)
      96             :         !! or the moments are calculated from the density 
      97             :         USE m_types
      98             :         USE m_constants
      99             :         USE m_xmlOutput
     100             :         
     101             :         USE m_intgr, ONLY : intgr3
     102             :         TYPE(t_input), INTENT(IN)           :: input
     103             :         TYPE(t_atoms), INTENT(IN)           :: atoms
     104             :         TYPE(t_noco), INTENT(IN)            :: noco
     105             :         TYPE(t_nococonv), INTENT(IN)        :: nococonv
     106             :         TYPE(t_moments),INTENT(IN)          :: moments
     107             :      
     108             :         INTEGER                       :: iType
     109             :         REAL                          :: local_m(0:3),global_m(3),theta,phi
     110         196 :         if (input%jspins==1) return ! No spin-pol calculation!
     111             :         
     112             : 
     113             :         !Header of output
     114         196 :         write(ounit,*)
     115         196 :         write(oUnit,'(a,19x,a,19x,a,19x,a)') "Orbital  moments |","Global Frame","  | ","Local Frame"
     116         196 :         write(oUnit,*) "------------------------------------------------------------------------------------------------------------------------"
     117         196 :         write(oUnit,'(a,5x,a,2(" | ",5(a,5x)))') "Atom ","|m|   ","mx   ","my   ","mz   ","alpha","beta ","mx   ","my   ","mz   ","alpha","beta "
     118             :            
     119         588 :         CALL openXMLElement('orbitalMomentsInMTSpheres',(/'units'/),(/'muBohr'/))
     120             : 
     121             : 
     122         524 :         DO iType = 1, atoms%ntype
     123             :            !Collect the data
     124        1312 :             global_m(:)=moments%clmom(:,iType,1)+moments%clmom(:,iType,2)
     125         328 :             if (noco%l_noco) THEN
     126         150 :                 theta = nococonv%beta(iType)
     127         150 :                 phi   = nococonv%alph(iType)
     128             :             else
     129         178 :                 if (noco%l_soc) THEN    
     130          29 :                     theta = nococonv%theta
     131          29 :                     phi   = nococonv%phi
     132             :                 else
     133         149 :                     theta=0.0
     134         149 :                     phi=0.0
     135             :                 endif        
     136             :             endif
     137             :             
     138             :            !Rotate to local frame
     139        1312 :            local_m(1:3)=global_m
     140         328 :            local_m(0)=0 
     141         328 :            call nococonv%rot_magvec(phi,theta,local_m, toGlobal=.false.)
     142             :            
     143         524 :            call priv_print_mt_moment(itype,.true.,.true.,.true.,local_m(1:3),global_m,"omm")
     144             :         enddo   
     145             :      
     146         196 :         CALL closeXMLElement('orbitalMomentsInMTSpheres')
     147         196 :         write(oUnit,*) "------------------------------------------------------------------------------------------------------------------------"
     148         196 :         write(oUnit,*)
     149             :         
     150             :     END SUBROUTINE 
     151             :      
     152         656 :     subroutine priv_print_mt_moment(itype,l_noco,l_offdiag,l_soc,magmomL,magmom,grepstring)
     153             :         USE m_xmlOutput
     154             :         USE m_constants
     155             :         USE m_types_noco
     156             :         USE m_polangle
     157             :         integer,intent(in):: itype 
     158             :         logical,intent(in):: l_noco,l_offdiag,l_soc
     159             :         real,intent(in)   :: magmom(3)
     160             :         real,intent(in)   :: magmomL(3)
     161             :         character(len=*),INTENT(IN) :: grepstring
     162             :         
     163             :         
     164             :      
     165             :         character(len=15):: label 
     166             :         character(len=30):: attributes(2)      
     167             :         real :: alphal,betal,alpha,beta
     168             : 
     169             :         !Calculate angles
     170         656 :         CALL pol_angle(magmomL(1),magmomL(2),magmomL(3),betal,alphal,.true.)
     171         656 :         CALL pol_angle(magmom(1),magmom(2),magmom(3),beta,alpha,.true.)
     172             :                 
     173         656 :         if (l_offdiag) THEN
     174        1460 :            write(oUnit,'(i6,1x,f9.6,2(" | ",5(f9.6,1x)),"  mm ")') itype,sqrt(dot_product(magmom(1:3),magmom(1:3))),magmom,alpha,beta,magmomL,alphaL,betaL
     175         291 :         elseif (l_noco.or.l_soc) THEN
     176             :            write(oUnit,'(i6,1x,f9.6," | ",5(f9.6,1x)," |    ---       ---    ",f9.6,"    ---          ---           ",a)') &
     177         568 :             itype,sqrt(dot_product(magmom(1:3),magmom(1:3))),magmom,alpha,beta,magmomL(3),grepstring
     178             :         else
     179             :            write(oUnit,'(i6,1x,f9.6," |    ---       ---       ---       ---         ---   |    ---       ---       ---       ---          ---           ",a)') &
     180         149 :            itype,magmom(3),grepstring
     181             :         endif
     182             : 
     183         656 :         WRITE(attributes(1),'(i0)') iType
     184         656 :         WRITE(attributes(2),'(3(f9.5,1x))') magmom(1),magmom(2),magmom(3)
     185         656 :         label="globalMagMoment"
     186             :         CALL writeXMLElementFormPoly(label,(/'atomType','vec     '/),&
     187        1968 :                         attributes,reshape((/8,3,6,30/),(/2,2/)))
     188         656 :         WRITE(attributes(2),'(3(f9.5,1x))') magmomL(1),magmomL(2),magmomL(3)
     189         656 :         label="localMagMoment "
     190             :         CALL writeXMLElementFormPoly(label,(/'atomType','vec     '/),&
     191        1968 :                         attributes,reshape((/8,3,6,30/),(/2,2/)))
     192             :                                
     193             :      
     194         656 :     end subroutine
     195           0 :     subroutine priv_print_spin_density_at_nucleus(input,atoms,moments)
     196             :         USE m_types
     197             :         USE m_constants
     198             :         TYPE(t_input), INTENT(IN)           :: input
     199             :         TYPE(t_atoms), INTENT(IN)           :: atoms
     200             :         TYPE(t_moments),INTENT(IN)          :: moments
     201             :         INTEGER                       :: iType
     202             :         REAL                          :: sval,stot,scor
     203           0 :         WRITE (oUnit,FMT=8000)
     204           0 :         DO iType = 1,atoms%ntype
     205           0 :          sval = moments%svdn(iType,1) - moments%svdn(iType,input%jspins)
     206           0 :          stot = moments%stdn(iType,1) - moments%stdn(iType,input%jspins)
     207           0 :          scor = stot - sval
     208           0 :          WRITE (oUnit,FMT=8010) iType,stot,sval,scor,moments%svdn(iType,1),moments%stdn(iType,1)
     209             :         END DO
     210             :         8000 FORMAT (/,/,10x,'spin density at the nucleus:',/,10x,'type',t25,&
     211             :                  'total',t42,'valence',t65,'core',t90,&
     212             :                  'majority valence and total density',/)
     213             :         8010 FORMAT (i13,2x,3e20.8,5x,2e20.8)
     214           0 :     end subroutine
     215             : end module    

Generated by: LCOV version 1.14