TABLE OF CONTENTS
ABINIT/generate_training_set [ Functions ]
NAME
generate_training_set
FUNCTION
INPUTS
acell(3)=length scales by which rprim is to be multiplied amplitudes(namplitude) = list of the amplitudes of the unstable phonons amplitudes(1:3,iamplitude) = qpt amplitudes(4,iamplitude) = mode amplitudes(5,iamplitude) = amplitude natom= Number of atoms in the supercell nconfig = Number of configuration to generate hist<type(abihist)> = The history of the MD namplitude = number of amplitude provided by the user ngqpt(3)= Grid of qpoint in the DDB nqshft=number of shift vectors in the repeated cell option = option to deal with negative frequency -> Bose factor explodes (eg acoustic at Gamma) several philosophies to be implemented for the unstable modes: option == 1 => ignore option == 2 => populate them according to a default amplitude option == 3 => populate according to their modulus squared option == 4 => USER defined value(s), require namplitude and amplitude qshft(3,nqshft)=vectors that will be used to determine rlatt(3,3)= size of the supercell rprimd(3,3)=dimensional primitive translations (bohr) temperature_K = temperature in Kelvin xred(3,natom) = reduced dimensionless atomic coordinates in the supercell comm=MPI communicator
OUTPUT
SIDE EFFECTS
hist<type abihist>=Historical record of positions
PARENTS
multibinit
CHILDREN
SOURCE
85 subroutine generate_training_set(acell,add_strain,amplitudes,filename,hist,natom,namplitude,nconfig,& 86 & ngqpt,nqshift,option,qshift,rlatt,rprimd,temperature_k,xred,comm,DEBUG) 87 88 use defs_basis 89 use m_errors 90 use m_abicore 91 use m_xmpi 92 use m_strain 93 use m_abihist, only : abihist,var2hist,abihist_findIndex 94 use m_ifc, only : ifc_type,ifc_init_fromFile,ifc_free 95 use m_crystal, only : crystal_t,crystal_free 96 use m_supercell, only : supercell_type 97 use m_geometry, only : xcart2xred 98 use m_phonons ,only :thermal_supercell_make,thermal_supercell_free 99 100 !This section has been created automatically by the script Abilint (TD). 101 !Do not modify the following lines by hand. 102 #undef ABI_FUNC 103 #define ABI_FUNC 'generate_training_set' 104 !End of the abilint section 105 106 implicit none 107 108 !Arguments ------------------------------------ 109 !scalars 110 integer,intent(in) :: natom,nconfig,nqshift,option,comm 111 logical,intent(in) :: DEBUG 112 real(dp),intent(in):: temperature_k 113 integer,intent(in) :: namplitude 114 logical,intent(in) :: add_strain 115 !arrays 116 integer,intent(in) :: ngqpt(3) 117 integer,intent(in) :: rlatt(3,3) 118 real(dp),intent(in) :: amplitudes(5,namplitude) 119 real(dp),intent(in) :: qshift(3,nqshift) 120 real(dp),intent(in) :: acell(3) 121 real(dp), intent(in) :: rprimd(3,3) 122 real(dp),intent(in) :: xred(3,natom) 123 character(len=*),intent(in) :: filename 124 type(abihist),intent(inout) :: hist 125 !Local variables------------------------------- 126 !scalar 127 real(dp):: delta,rand 128 integer :: ii,iconfig,natom_uc,direction 129 character(len=500) :: message 130 131 !arrays 132 real(dp) :: dielt(3,3) 133 real(dp),allocatable :: zeff(:,:,:) 134 real(dp) :: acell_next(3),xred_next(3,natom),rprimd_next(3,3) 135 type(ifc_type) :: ifc 136 type(crystal_t) :: crystal 137 type(supercell_type),allocatable :: thm_scells(:) 138 type(strain_type) :: strain 139 140 ! ************************************************************************* 141 142 write(message,'(a,(80a),a)') ch10,('=',ii=1,80),ch10 143 call wrtout(ab_out,message,'COLL') 144 call wrtout(std_out,message,'COLL') 145 146 write(message, '(a,a,a)' )' Generation of all the configuration for the training set',ch10 147 call wrtout(std_out,message,'COLL') 148 call wrtout(ab_out,message,'COLL') 149 150 call ifc_init_fromFile(dielt,trim(filename),ifc,natom_uc,ngqpt,nqshift,qshift,crystal,zeff,comm) 151 152 write(message, '(a,I0,a,f10.2,02a)' )' Generation of ',nconfig,' at the temperature ',& 153 & temperature_K,' K from the phonons',ch10 154 155 call wrtout(std_out,message,'COLL') 156 call wrtout(ab_out,message,'COLL') 157 158 ABI_DATATYPE_ALLOCATE(thm_scells,(nconfig)) 159 160 call thermal_supercell_make(amplitudes,crystal, Ifc,namplitude, nconfig, option,int(rlatt),& 161 & temperature_K, thm_scells) 162 163 do iconfig = 1,nconfig 164 ! Get the atomic position for the new configuration 165 call xcart2xred(thm_scells(iconfig)%natom,thm_scells(iconfig)%rprimd,& 166 & thm_scells(iconfig)%xcart,xred_next) 167 ! Just fill acell with the reference value, we apply strain on rprimd 168 acell_next(:) = acell(:) 169 170 ! Get the rprim for the new configuration 171 if(.not.add_strain)then 172 rprimd_next(:,:) = thm_scells(iconfig)%rprimd(:,:) 173 else 174 call random_number(rand) 175 direction = int(rand*6+1) 176 call random_number(rand) 177 delta = rand/500 178 call random_number(rand) 179 delta = delta*(two*rand-1) 180 181 call strain_init(strain,delta=delta,direction=direction) 182 call strain_apply(thm_scells(iconfig)%rprimd,rprimd_next,strain) 183 end if 184 185 ! Fill history with the values of xred, acell and rprimd 186 call var2hist(acell_next,hist,natom,rprimd_next,xred_next,DEBUG) 187 hist%ihist = abihist_findIndex(hist,+1) 188 189 end do 190 191 ! Restart ihist before to leave 192 hist%ihist = 1 193 194 call ifc_free(ifc) 195 call crystal_free(crystal) 196 call thermal_supercell_free(nconfig,thm_scells) 197 ABI_DATATYPE_DEALLOCATE(thm_scells) 198 ABI_DEALLOCATE(zeff) 199 200 end subroutine generate_training_set
ABINIT/m_generate_training_set [ Modules ]
NAME
m_generate_training_set
FUNCTION
COPYRIGHT
Copyright (C) 2008-2018 ABINIT group () This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
PARENTS
CHILDREN
SOURCE
21 #if defined HAVE_CONFIG_H 22 #include "config.h" 23 #endif 24 25 #include "abi_common.h" 26 27 module m_generate_training_set 28 29 implicit none 30 31 private