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
SOURCE
77 subroutine generate_training_set(acell,add_strain,amplitudes,filename,hist,natom,namplitude,nconfig,& 78 & ngqpt,nqshift,option,qshift,rlatt,rprimd,temperature_k,xred,comm,DEBUG) 79 80 use defs_basis 81 use m_errors 82 use m_abicore 83 use m_xmpi 84 use m_strain 85 use m_abihist, only : abihist,var2hist,abihist_findIndex 86 use m_ifc, only : ifc_type 87 use m_crystal, only : crystal_t 88 use m_supercell, only : supercell_type 89 use m_geometry, only : xcart2xred 90 use m_phonons ,only :thermal_supercell_make,thermal_supercell_free 91 use m_xmpi 92 93 !Arguments ------------------------------------ 94 !scalars 95 integer,intent(in) :: natom,nconfig,nqshift,option,comm 96 logical,intent(in) :: DEBUG 97 real(dp),intent(in):: temperature_k 98 integer,intent(in) :: namplitude 99 logical,intent(in) :: add_strain 100 !arrays 101 integer,intent(in) :: ngqpt(3) 102 integer,intent(in) :: rlatt(3,3) 103 real(dp),intent(in) :: amplitudes(5,namplitude) 104 real(dp),intent(in) :: qshift(3,nqshift) 105 real(dp),intent(in) :: acell(3) 106 real(dp), intent(in) :: rprimd(3,3) 107 real(dp),intent(in) :: xred(3,natom) 108 character(len=*),intent(in) :: filename 109 type(abihist),intent(inout) :: hist 110 !Local variables------------------------------- 111 !scalar 112 real(dp):: delta,rand 113 integer :: ii,iconfig,natom_uc,direction 114 character(len=500) :: message 115 INTEGER :: n 116 INTEGER :: i, iseed 117 INTEGER, DIMENSION(:), ALLOCATABLE :: seed 118 119 !arrays 120 real(dp) :: dielt(3,3) 121 real(dp),allocatable :: zeff(:,:,:),qdrp_cart(:,:,:,:) 122 real(dp) :: acell_next(3),xred_next(3,natom),rprimd_next(3,3) 123 type(ifc_type) :: ifc 124 type(crystal_t) :: crystal 125 type(supercell_type),allocatable :: thm_scells(:) 126 type(strain_type) :: strain 127 128 ! ************************************************************************* 129 130 ABI_UNUSED((/rprimd(1,1), xred(1,1)/)) 131 132 if ( .not. DEBUG ) then 133 CALL RANDOM_SEED(size = n) 134 ABI_MALLOC(seed,(n)) 135 seed = iseed + (/ (i - 1, i = 1, n) /) 136 137 CALL RANDOM_SEED(PUT = seed+xmpi_comm_rank(xmpi_world)) 138 end if 139 140 write(message,'(a,(80a),a)') ch10,('=',ii=1,80),ch10 141 call wrtout(ab_out,message,'COLL') 142 call wrtout(std_out,message,'COLL') 143 144 write(message, '(a,a,a)' )' Generation of all the configuration for the training set',ch10 145 call wrtout(std_out,message,'COLL') 146 call wrtout(ab_out,message,'COLL') 147 148 call ifc%from_file(dielt,trim(filename),natom_uc,ngqpt,nqshift,qshift,crystal,zeff,qdrp_cart,comm) 149 150 write(message, '(a,I0,a,f10.2,02a)' )' Generation of ',nconfig,' at the temperature ',& 151 & temperature_K,' K from the phonons',ch10 152 153 call wrtout(std_out,message,'COLL') 154 call wrtout(ab_out,message,'COLL') 155 156 ABI_MALLOC(thm_scells,(nconfig)) 157 158 call thermal_supercell_make(amplitudes,crystal, Ifc,namplitude, nconfig, option,int(rlatt),& 159 & temperature_K, thm_scells) 160 161 do iconfig = 1,nconfig 162 ! Get the atomic position for the new configuration 163 call xcart2xred(thm_scells(iconfig)%natom,thm_scells(iconfig)%rprimd,& 164 & thm_scells(iconfig)%xcart,xred_next) 165 ! Just fill acell with the reference value, we apply strain on rprimd 166 acell_next(:) = acell(:) 167 168 ! Get the rprim for the new configuration 169 if(.not.add_strain)then 170 rprimd_next(:,:) = thm_scells(iconfig)%rprimd(:,:) 171 else 172 call random_number(rand) 173 direction = int(rand*6+1) 174 call random_number(rand) 175 delta = rand/500 176 call random_number(rand) 177 delta = delta*(two*rand-1) 178 179 call strain_init(strain,delta=delta,direction=direction) 180 call strain_apply(thm_scells(iconfig)%rprimd,rprimd_next,strain) 181 end if 182 183 ! Fill history with the values of xred, acell and rprimd 184 call var2hist(acell_next,hist,natom,rprimd_next,xred_next,DEBUG) 185 hist%ihist = abihist_findIndex(hist,+1) 186 187 end do 188 189 ! Restart ihist before to leave 190 hist%ihist = 1 191 192 call ifc%free() 193 call crystal%free() 194 call thermal_supercell_free(nconfig,thm_scells) 195 ABI_FREE(thm_scells) 196 ABI_FREE(zeff) 197 ABI_FREE(qdrp_cart) 198 199 end subroutine generate_training_set
ABINIT/m_generate_training_set [ Modules ]
NAME
m_generate_training_set
FUNCTION
XG20200322 : Was initiated by Alexandre Martin, but not tested, not documented, perhaps not completed ?!? So, should be removed at some point ...
COPYRIGHT
Copyright (C) 2008-2024 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 .
SOURCE
18 #if defined HAVE_CONFIG_H 19 #include "config.h" 20 #endif 21 22 #include "abi_common.h" 23 24 module m_generate_training_set 25 26 implicit none 27 28 private