TABLE OF CONTENTS
- ABINIT/m_harmonics_terms
- m_harmonics_terms/harmonics_terms_evaluateElastic
- m_harmonics_terms/harmonics_terms_evaluateIFC
- m_harmonics_terms/harmonics_terms_init
- m_harmonics_terms/harmonics_terms_type
ABINIT/m_harmonics_terms [ Functions ]
NAME
m_harmonics_term
FUNCTION
Module with datatype and tools for the harmonics terms
COPYRIGHT
Copyright (C) 2010-2024 ABINIT group (AM) This file is distributed under the terms of the GNU General Public Licence, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
SOURCE
19 #if defined HAVE_CONFIG_H 20 #include "config.h" 21 #endif 22 23 #include "abi_common.h" 24 25 module m_harmonics_terms 26 27 use defs_basis 28 use m_errors 29 use m_abicore 30 use m_supercell,only: getPBCIndexes_supercell 31 use m_xmpi,only : xmpi_sum 32 use m_ifc 33 34 implicit none 35 36 public :: harmonics_terms_init 37 public :: harmonics_terms_free 38 public :: harmonics_terms_applySumRule 39 public :: harmonics_terms_evaluateIFC 40 public :: harmonics_terms_evaluateElastic 41 public :: harmonics_terms_setEffectiveCharges 42 public :: harmonics_terms_setDynmat 43 public :: harmonics_terms_setInternalStrain
m_harmonics_terms/harmonics_terms_evaluateElastic [ Functions ]
[ Top ] [ m_harmonics_terms ] [ Functions ]
NAME
harmonics_terms_evaluateElastic
FUNCTION
Compute the energy, forces and stresses related to the application of strain
INPUTS
elastic_constants(6,6) = elastic constants in Hartree disp(3,natom_sc) = atomics displacement between configuration and the reference natom = number of atoms in the supercell natom_uc = number of atoms in the unit cell ncell = total number of cell strain_coupling(6,3,natom) = internal strain coupling parameters strain(6) = strain between configuration and the reference
OUTPUT
energy = contribution to the energy fcart(3,natom) = contribution to the forces strten(6) = contribution to the stress tensor
SOURCE
642 subroutine harmonics_terms_evaluateElastic(elastic_constants,disp,energy,fcart,natom,natom_uc,ncell,& 643 & strain_coupling,strten,strain) 644 645 real(dp),intent(out):: energy 646 integer, intent(in) :: natom,natom_uc,ncell 647 ! array 648 real(dp),intent(in) :: elastic_constants(6,6),strain_coupling(6,3,natom) 649 real(dp),intent(out):: strten(6) 650 real(dp),intent(out):: fcart(3,natom) 651 real(dp),intent(in) :: disp(3,natom) 652 real(dp),intent(in) :: strain(6) 653 654 !Local variables------------------------------- 655 ! scalar 656 integer :: ia,ii,mu,alpha,beta 657 real(dp):: cij 658 ! array 659 ! ************************************************************************* 660 661 energy = zero 662 fcart = zero 663 strten = zero 664 665 ! write(*,*) "----- STRAIN -----" 666 ! write(*,*) strain 667 668 !1- Part due to elastic constants 669 do alpha=1,6 670 do beta=1,6 671 ! write(*,*) "--- cij --- alpha: ", alpha, " beta: ", beta 672 cij = ncell*elastic_constants(alpha,beta) 673 ! write(*,*) cij 674 energy = energy + half*cij*strain(alpha)*strain(beta) 675 strten(alpha) = strten(alpha) + cij*strain(beta) 676 end do 677 ! write(*,*) "strten(",alpha,"): ", strten(alpha) 678 end do 679 680 !2-Part due to the internal strain coupling parameters 681 ii = 1 682 do ia = 1,natom 683 do mu = 1,3 684 do alpha=1,6 685 cij = strain_coupling(alpha,mu,ii) 686 ! Accumulte for this atom 687 energy = energy + half*cij*strain(alpha)*disp(mu,ia) 688 fcart(mu,ia) = fcart(mu,ia) + half*cij*strain(alpha) 689 strten(alpha) = strten(alpha) + half*cij*disp(mu,ia) 690 end do 691 end do 692 ii = ii +1 693 ! Reset to 1 if the number of atoms is superior than in the initial cell 694 if(ii==natom_uc+1) ii = 1 695 end do 696 697 ! write(*,*) "--- STRTEN at the end --- " 698 ! write(*,*) strten(:) 699 700 end subroutine harmonics_terms_evaluateElastic
m_harmonics_terms/harmonics_terms_evaluateIFC [ Functions ]
[ Top ] [ m_harmonics_terms ] [ Functions ]
NAME
harmonics_terms_evaluateIFC
FUNCTION
This fonction compute the contribution of the ifc harmonic part of the energy and forces.
INPUTS
atmfrc(3,natom_uc,3,natom_uc,nrpt) = atomic force constants disp(3,natom_sc) = atomics displacement between configuration and the reference ncell = total number of cell to treat nrpt = total number of rpt to treat natom_sc = number of atoms in the supercell natom_uc = number of atoms in the unit cell nrpt = number of rpt atmrpt_index(nrpt,cell) = For each cell in the supercell and each rpt, give the index of the first atoms in the rpt cell rpt(nrpt) = index of rpt in atmfrc (6th dimension) index_cells(3,ncell) = indexes of the cells into supercell (-1 -1 -1 ,...,1 1 1) comm=MPI communicator
OUTPUT
energy = contribution of the ifc to the energy fcart(3,natom) = contribution of the ifc to the forces PARENT effective_potential_evaluate
SOURCE
536 subroutine harmonics_terms_evaluateIFC(atmfrc,disp,energy,fcart,natom_sc,natom_uc,& 537 & ncell,nrpt,atmrpt_index,index_cells,sc_size,rpt,comm) 538 539 implicit none 540 541 !Arguments ------------------------------- 542 ! scalars 543 real(dp),intent(out) :: energy 544 integer,intent(in) :: natom_uc,natom_sc,ncell,nrpt 545 integer,intent(in) :: comm 546 ! array 547 integer,intent(in) :: sc_size(3),atmrpt_index(nrpt,ncell) 548 integer,intent(in) :: index_cells(4,ncell),rpt(nrpt) 549 real(dp),intent(in) :: atmfrc(3,natom_uc,3,natom_uc,nrpt) 550 real(dp),intent(in) :: disp(3,natom_sc) 551 real(dp),intent(out) :: fcart(3,natom_sc) 552 553 !Local variables------------------------------- 554 ! scalar 555 integer :: i1,i2,i3,ia,ib,icell,ierr,irpt,irpt_tmp,ii,jj,kk,ll 556 integer :: mu,nu 557 real(dp):: disp1,disp2,ifc,tmp_etot1,tmp_etot2 558 !Variables for separation of short and dipdip ifc contribution 559 !real(dp):: short_ifc,ewald_ifc 560 !real(dp):: tmp_ewald1,tmp_ewald2,tmp_short1,tmp_short2 561 ! array 562 character(500) :: msg 563 564 ! ************************************************************************* 565 566 if (any(sc_size <= 0)) then 567 write(msg,'(a,a)')' sc_size can not be inferior or equal to zero' 568 ABI_ERROR(msg) 569 end if 570 571 ! Initialisation of variables 572 energy = zero 573 fcart(:,:) = zero 574 575 do icell = 1,ncell 576 i1 = index_cells(1,icell) 577 i2 = index_cells(2,icell) 578 i3 = index_cells(3,icell) 579 ! index of the first atom in the current cell 580 ii = index_cells(4,icell) 581 do irpt_tmp = 1,nrpt 582 irpt = rpt(irpt_tmp) 583 ! index of the first atom in the irpt cell 584 jj = atmrpt_index(irpt_tmp,icell) 585 ! Loop over the atom in the cell 586 do ib = 1, natom_uc 587 ll = jj + ib 588 do nu=1,3 589 disp2 = disp(nu,ll) 590 do ia = 1, natom_uc 591 kk = ii + ia 592 do mu=1,3 593 disp1 = disp(mu,kk) 594 ifc = atmfrc(mu,ia,nu,ib,irpt) 595 596 ! if(abs(ifc) > tol10)then 597 tmp_etot1 = disp2 * ifc 598 ! accumule energy 599 tmp_etot2 = disp1*tmp_etot1 600 energy = energy + tmp_etot2 601 ! accumule forces 602 fcart(mu,kk) = fcart(mu,kk) + tmp_etot1 603 ! end if 604 end do 605 end do 606 end do 607 end do 608 end do 609 end do 610 611 energy = half * energy 612 ! MPI_SUM 613 call xmpi_sum(energy, comm, ierr) 614 call xmpi_sum(fcart , comm, ierr) 615 616 end subroutine harmonics_terms_evaluateIFC
m_harmonics_terms/harmonics_terms_init [ Functions ]
[ Top ] [ m_harmonics_terms ] [ Functions ]
NAME
harmonics_terms_init
FUNCTION
Initialize harmonics_terms datatype
INPUTS
ifc<type(ifc_type)> = interatomic forces constants natom = number of atoms in primitive cell nrpt = number rpt (cell) in the ifc dynmat(2,3,natom,3,natom,3,nqpt) = optional, dynamical matricies for each q-point epsilon_inf(3,3) = optional, dielectric tensor elastic_constant(6,6) = optional, elastic constant strain_coupling(6,3,natom) = optional, internal strain coupling parameters nqpt = optional, number of q-points phfrq(3*natom,nqpt) = optional,phonons frequencies for each q points in Hartree/cm qpoints(3,nqpt) = list of qpoints wavevectors zeff(3,3,natom) = optional,effective charges
OUTPUT
harmonics_terms<type(harmonics_terms_type)> = harmonics_terms datatype to be initialized
SOURCE
123 subroutine harmonics_terms_init(harmonics_terms,ifcs,natom,nrpt,& 124 & dynmat,epsilon_inf,elastic_constants,strain_coupling,& 125 & nqpt,phfrq,qpoints,zeff) 126 127 implicit none 128 129 !Arguments ------------------------------------ 130 !scalars 131 integer, intent(in) :: natom,nrpt 132 !arrays 133 type(ifc_type),intent(in) :: ifcs 134 type(harmonics_terms_type), intent(out) :: harmonics_terms 135 integer, optional,intent(in) :: nqpt 136 real(dp),optional,intent(in) :: epsilon_inf(3,3),dynmat(:,:,:,:,:,:) 137 real(dp),optional,intent(in) :: elastic_constants(6,6) 138 real(dp),optional,intent(in) :: strain_coupling(6,3,natom),zeff(3,3,natom) 139 real(dp),optional,intent(in) :: phfrq(:,:),qpoints(:,:) 140 !Local variables------------------------------- 141 !scalar 142 !arrays 143 character(len=500) :: msg 144 145 ! ************************************************************************* 146 147 call harmonics_terms_free(harmonics_terms) 148 149 ! Do some Checks 150 if (natom < 1) then 151 write(msg, '(a,a,a,i10,a)' )& 152 & 'The cell must have at least one atom.',ch10,& 153 & 'The number of atom is ',natom,'.' 154 ABI_BUG(msg) 155 end if 156 157 if (nrpt < 1) then 158 write(msg, '(a,a,a,i10,a)' )& 159 & 'The cell must have at least one rpt point.',ch10,& 160 & 'The number of rpt points is ',nrpt,'.' 161 ABI_BUG(msg) 162 end if 163 164 if (nrpt /= ifcs%nrpt) then 165 write(msg, '(3a,i5,a,i5,a)' )& 166 & 'nrpt must have the same dimension as ifcs.',ch10,& 167 & 'The number of cell is ',nrpt,' instead of ',ifcs%nrpt,'.' 168 ABI_BUG(msg) 169 end if 170 171 if(present(nqpt).and.(.not.present(dynmat).or.& 172 & .not.present(qpoints) .or.& 173 & .not.present(phfrq)))then 174 write(msg, '(a)' )& 175 & 'nqpt is specified but dynamt,qpoints or phfrq are not.' 176 ABI_BUG(msg) 177 end if 178 179 if(.not.present(nqpt).and.(present(dynmat).or.& 180 & present(qpoints) .or.& 181 & present(phfrq)))then 182 write(msg, '(a)' )& 183 & ' dynamt,qpoints or phfrq are specified but nqpt is not.' 184 ABI_BUG(msg) 185 end if 186 187 !Set number of cell 188 harmonics_terms%ifcs%nrpt = nrpt 189 190 !Allocation of total ifc 191 ABI_MALLOC(harmonics_terms%ifcs%atmfrc,(3,natom,3,natom,nrpt)) 192 harmonics_terms%ifcs%atmfrc(:,:,:,:,:) = ifcs%atmfrc(:,:,:,:,:) 193 194 !Allocation of ewald part of ifc 195 ABI_MALLOC(harmonics_terms%ifcs%ewald_atmfrc,(3,natom,3,natom,nrpt)) 196 harmonics_terms%ifcs%ewald_atmfrc(:,:,:,:,:) = ifcs%ewald_atmfrc(:,:,:,:,:) 197 198 !Allocation of short range part of ifc 199 ABI_MALLOC(harmonics_terms%ifcs%short_atmfrc,(3,natom,3,natom,nrpt)) 200 harmonics_terms%ifcs%short_atmfrc(:,:,:,:,:) = ifcs%short_atmfrc(:,:,:,:,:) 201 202 !Allocation of cell of ifc 203 ABI_MALLOC(harmonics_terms%ifcs%cell,(3,nrpt)) 204 harmonics_terms%ifcs%cell(:,:) = ifcs%cell(:,:) 205 206 !Allocation of the dynamical matrix 207 harmonics_terms%nqpt = 0 208 if(present(nqpt).and.present(dynmat).and.present(qpoints).and.present(phfrq))then 209 call harmonics_terms_setDynmat(dynmat,harmonics_terms,natom,nqpt,& 210 & harmonics_terms%phfrq,harmonics_terms%qpoints) 211 end if 212 213 !Allocation of the elastic constants 214 harmonics_terms%elastic_constants = zero 215 if (present(elastic_constants)) then 216 harmonics_terms%elastic_constants = elastic_constants 217 end if 218 219 !Allication of the dielectric tensor 220 harmonics_terms%epsilon_inf = zero 221 if (present(epsilon_inf)) then 222 harmonics_terms%epsilon_inf = epsilon_inf 223 end if 224 225 !Allocation of Effective charges array 226 ABI_MALLOC(harmonics_terms%zeff,(3,3,natom)) 227 harmonics_terms%zeff = zero 228 if (present(zeff)) then 229 call harmonics_terms_setEffectiveCharges(harmonics_terms,natom,zeff) 230 harmonics_terms%zeff = zeff 231 end if 232 233 !Allocation of internal strain tensor 234 ABI_MALLOC(harmonics_terms%strain_coupling,(6,3,natom)) 235 harmonics_terms%strain_coupling = zero 236 if (present(strain_coupling)) then 237 call harmonics_terms_setInternalStrain(harmonics_terms,natom,strain_coupling) 238 end if 239 240 end subroutine harmonics_terms_init
m_harmonics_terms/harmonics_terms_type [ Types ]
[ Top ] [ m_harmonics_terms ] [ Types ]
NAME
harmonics_terms_type
FUNCTION
datatype for harmonic part of effective potential.
SOURCE
55 type, public :: harmonics_terms_type 56 57 integer :: nqpt 58 ! Number of qpoints 59 60 real(dp) :: epsilon_inf(3,3) 61 ! epsilon_inf(3,3) 62 ! Dielectric tensor 63 64 real(dp) :: elastic_constants(6,6) 65 ! elastic_constant(6,6) 66 ! Elastic tensor Hartree 67 68 real(dp), allocatable :: strain_coupling(:,:,:) 69 ! strain_coupling(6,3,natom) 70 ! internal strain tensor 71 72 real(dp), allocatable :: zeff(:,:,:) 73 ! zeff(3,3,natom) Effective charges 74 75 type(ifc_type) :: ifcs 76 ! type with ifcs constants (short + ewald) 77 ! also contains the number of cell and the indexes 78 79 real(dp), allocatable :: qpoints(:,:) 80 ! qph1l(3,nqpt) 81 ! List of qpoints wavevectors 82 83 real(dp), allocatable :: dynmat(:,:,:,:,:,:) 84 ! dynmat(2,3,natom,3,natom,nqpt) 85 ! dynamical matrix for each q points 86 87 real(dp), allocatable :: phfrq(:,:) 88 ! phfrq(3*natom,nqpt) 89 ! array with all phonons frequencies for each q points in Hartree/cm 90 91 end type harmonics_terms_type