TABLE OF CONTENTS
- ABINIT/m_fit_data
- m_fit_data/fit_data_compute
- m_fit_data/fit_data_free
- m_fit_data/fit_data_init
- m_fit_data/fit_data_type
- m_fit_data/training_set_free
- m_fit_data/training_set_init
- m_fit_data/training_set_type
ABINIT/m_fit_data [ Modules ]
NAME
m_fit_data
FUNCTION
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
18 #if defined HAVE_CONFIG_H 19 #include "config.h" 20 #endif 21 22 #include "abi_common.h" 23 24 module m_fit_data 25 26 use defs_basis 27 use m_errors 28 use m_abicore 29 30 use m_geometry, only : metric 31 32 implicit none
m_fit_data/fit_data_compute [ Functions ]
[ Top ] [ m_fit_data ] [ Functions ]
NAME
fit_data_compute
FUNCTION
Conpute the strain of each configuration. Compute the displacmeent of each configuration. Compute the variation of the displacement due to strain of each configuration. Compute fixed forces and stresse and get the standard deviation. Compute Sheppard and al Factors \Omega^{2} see J.Chem Phys 136, 074103 (2012) [[cite:Sheppard2012]]
INPUTS
eff_pot<type(effective_potential)> = effective potential hist<type(abihist)> = The history of the MD (or snapshot of DFT comm = MPI communicator verbose = optional, flag for the verbose mode
OUTPUT
fit_data<fit_data_type> = fit_data is now filled
SOURCE
279 subroutine fit_data_compute(fit_data,eff_pot,hist,comm,verbose) 280 281 use m_strain,only : strain_type,strain_get 282 use m_effective_potential,only : effective_potential_type,effective_potential_evaluate 283 use m_effective_potential,only : effective_potential_getDisp 284 use m_abihist, only : abihist 285 use m_strain,only : strain_type,strain_get 286 implicit none 287 288 !Arguments ------------------------------------ 289 !scalars 290 integer,intent(in) :: comm 291 logical,optional,intent(in) :: verbose 292 !arrays 293 type(fit_data_type),intent(inout) :: fit_data 294 type(effective_potential_type),intent(in) :: eff_pot 295 type(abihist),intent(in) :: hist 296 !Local variables------------------------------- 297 !scalar 298 integer :: ii,itime,natom,ntime 299 real(dp):: energy 300 logical :: need_verbose 301 !arrays 302 character(len=500) :: message 303 real(dp) :: gmet(3,3),gprimd(3,3),rmet(3,3) 304 real(dp),allocatable :: energy_diff(:) 305 real(dp),allocatable :: du_delta(:,:,:,:),displacement(:,:,:),strain(:,:) 306 real(dp),allocatable :: fcart_diff(:,:,:),gred_fixed(:,:,:),fcart_fixed(:,:,:) 307 real(dp),allocatable :: strten_diff(:,:),strten_fixed(:,:),sqomega(:),ucvol(:) 308 type(strain_type) :: strain_t 309 type(training_set_type) :: ts 310 ! ************************************************************************* 311 312 !Initialisation of optional arguments 313 need_verbose = .TRUE. 314 if(present(verbose)) need_verbose = verbose 315 316 !Checks 317 natom = eff_pot%supercell%natom 318 ntime = hist%mxhist 319 if(natom /= size(hist%xred,2))then 320 write(message, '(5a)' )& 321 & 'The number of atoms in the hist file does not correspond to the supercell.',ch10,& 322 & 'You should call the routine fit_polynomial_coeff_mapHistToRef before',ch10,& 323 & 'Action: Contact abinit group' 324 ABI_BUG(message) 325 end if 326 327 !Allocation of temporary arrays 328 ABI_MALLOC(displacement,(3,natom,ntime)) 329 ABI_MALLOC(du_delta,(6,3,natom,ntime)) 330 ABI_MALLOC(energy_diff,(ntime)) 331 ABI_MALLOC(fcart_fixed,(3,natom,ntime)) 332 ABI_MALLOC(fcart_diff,(3,natom,ntime)) 333 ABI_MALLOC(gred_fixed,(3,natom,ntime)) 334 ABI_MALLOC(strain,(6,ntime)) 335 ABI_MALLOC(strten_fixed,(6,ntime)) 336 ABI_MALLOC(strten_diff,(6,ntime)) 337 ABI_MALLOC(sqomega,(ntime)) 338 ABI_MALLOC(ucvol,(ntime)) 339 340 displacement = zero 341 du_delta = zero 342 strain = zero; 343 fcart_fixed = zero 344 gred_fixed = zero 345 strain = zero 346 strten_fixed = zero 347 strten_diff = zero 348 sqomega = zero 349 ucvol = zero 350 351 do itime=1,ntime 352 ! Get strain 353 call strain_get(strain_t,rprim=eff_pot%supercell%rprimd,& 354 & rprim_def=hist%rprimd(:,:,itime),symmetrized=.FALSE.) 355 if (strain_t%name /= "reference") then 356 do ii=1,3 357 strain(ii,itime) = strain_t%strain(ii,ii) 358 end do 359 strain(4,itime) = (strain_t%strain(2,3) + strain_t%strain(3,2)) 360 strain(5,itime) = (strain_t%strain(3,1) + strain_t%strain(1,3)) 361 strain(6,itime) = (strain_t%strain(2,1) + strain_t%strain(1,2)) 362 else 363 strain(:,itime) = zero 364 end if 365 366 ! Get displacement and du_delta 367 call effective_potential_getDisp(displacement(:,:,itime),du_delta(:,:,:,itime),natom,& 368 & hist%rprimd(:,:,itime),eff_pot%supercell%rprimd,comm,& 369 & xred_hist=hist%xred(:,:,itime),xcart_ref=eff_pot%supercell%xcart,& 370 & compute_displacement=.TRUE.,compute_duDelta=.TRUE.) 371 372 ! Get forces and stresses from harmonic part (fixed part) 373 call effective_potential_evaluate(eff_pot,energy,fcart_fixed(:,:,itime),gred_fixed(:,:,itime),& 374 & strten_fixed(:,itime),natom,hist%rprimd(:,:,itime),& 375 & displacement=displacement(:,:,itime),& 376 & du_delta=du_delta(:,:,:,itime),strain=strain(:,itime),& 377 & compute_anharmonic=.true.,verbose=.FALSE.) 378 379 ! Compute \Omega^{2} and ucvol for each time 380 call metric(gmet,gprimd,-1,rmet,hist%rprimd(:,:,itime),ucvol(itime)) 381 ! Formula: sqomega(itime) = (((ucvol(itime)**(-2.))* ((natom)**(0.5)))**(-1.0/3.0))**2 382 ! Compact form: 383 sqomega(itime) = ((ucvol(itime)**(4.0/3.0)) / ((natom)**(1/3.0))) 384 385 ! Compute the difference between History and model (fixed part) 386 fcart_diff(:,:,itime) = hist%fcart(:,:,itime) - fcart_fixed(:,:,itime) 387 energy_diff(itime) = hist%etot(itime) - energy 388 strten_diff(:,itime) = hist%strten(:,itime) - strten_fixed(:,itime) 389 end do ! End Loop itime 390 391 !Set the training set 392 call training_set_init(ts,displacement,du_delta,natom,ntime,strain,sqomega,ucvol) 393 !Set the fit_data 394 call fit_data_init(fit_data,energy_diff,fcart_diff,natom,ntime,strten_diff,ts) 395 396 !Free space 397 call training_set_free(ts) 398 ABI_FREE(displacement) 399 ABI_FREE(du_delta) 400 ABI_FREE(energy_diff) 401 ABI_FREE(fcart_fixed) 402 ABI_FREE(fcart_diff) 403 ABI_FREE(gred_fixed) 404 ABI_FREE(strain) 405 ABI_FREE(strten_fixed) 406 ABI_FREE(strten_diff) 407 ABI_FREE(sqomega) 408 ABI_FREE(ucvol) 409 410 end subroutine fit_data_compute
m_fit_data/fit_data_free [ Functions ]
[ Top ] [ m_fit_data ] [ Functions ]
NAME
fit_data_free
FUNCTION
Free the fit_data datatype
INPUTS
fit_data<fit_data_type> = fit_data to be free
OUTPUT
SOURCE
224 subroutine fit_data_free(fit_data) 225 226 implicit none 227 228 !Arguments ------------------------------------ 229 !scalars 230 !arrays 231 type(fit_data_type),intent(inout) :: fit_data 232 !Local variables------------------------------- 233 !scalar 234 !arrays 235 ! ************************************************************************* 236 237 ! Reset integer values 238 fit_data%ntime = 0 239 fit_data%natom = 0 240 241 ! Deallocate arrays 242 if(allocated(fit_data%energy_diff)) then 243 ABI_FREE(fit_data%energy_diff) 244 end if 245 if(allocated(fit_data%fcart_diff)) then 246 ABI_FREE(fit_data%fcart_diff) 247 end if 248 if(allocated(fit_data%strten_diff)) then 249 ABI_FREE(fit_data%strten_diff) 250 end if 251 call training_set_free(fit_data%training_set) 252 253 end subroutine fit_data_free
m_fit_data/fit_data_init [ Functions ]
[ Top ] [ m_fit_data ] [ Functions ]
NAME
fit_data_init
FUNCTION
Initialize fit_data datatype
INPUTS
energy_diff(3,natom,ntime) = Difference of energy between DFT calculation and fixed part of the model (more often harmonic part) fcart_diff(3,natom,ntime) = Difference of cartesian forces between DFT calculation and fixed part of the model (more often harmonic part) natom = Number of atoms ntime = Number of time (number of snapshot, number of md step...) strten_diff(6,natom) = Difference of stress tensor between DFT calculation and fixed part of the model (more often harmonic part) sqomega(ntime) = Sheppard and al Factors \Omega^{2} see J.Chem Phys 136, 074103 (2012) [[cite:Sheppard2012]] ucvol(ntime) = Volume of the system for each time ts<training_set_type> = datatype with the information about the training set
OUTPUT
fit_data<fit_data_type> = fit_data datatype to be initialized
SOURCE
159 subroutine fit_data_init(fit_data,energy_diff,fcart_diff,natom,ntime,strten_diff,ts) 160 161 !Arguments ------------------------------------ 162 !scalars 163 integer,intent(in) :: natom,ntime 164 !arrays 165 real(dp),intent(in) :: energy_diff(ntime),fcart_diff(3,natom,ntime) 166 real(dp),intent(in) :: strten_diff(6,ntime) 167 type(training_set_type),intent(in) :: ts 168 type(fit_data_type),intent(inout) :: fit_data 169 !Local variables------------------------------- 170 !scalar 171 !arrays 172 character(len=500) :: message 173 ! ************************************************************************* 174 175 if(natom /= ts%natom)then 176 write(message, '(a)')& 177 & ' The number of atoms does not correspond to the training set' 178 ABI_BUG(message) 179 end if 180 181 if(ntime /= ts%ntime)then 182 write(message, '(a)')& 183 & ' The number of time does not correspond to the training set' 184 ABI_BUG(message) 185 end if 186 187 !Free the output 188 call fit_data_free(fit_data) 189 190 !Set integer values 191 fit_data%ntime = ntime 192 fit_data%natom = natom 193 194 !allocate arrays 195 ABI_MALLOC(fit_data%fcart_diff,(3,natom,ntime)) 196 fit_data%fcart_diff(:,:,:) = fcart_diff(:,:,:) 197 198 ABI_MALLOC(fit_data%strten_diff,(6,ntime)) 199 fit_data%strten_diff(:,:) = strten_diff(:,:) 200 201 ABI_MALLOC(fit_data%energy_diff,(ntime)) 202 fit_data%energy_diff(:) = energy_diff 203 204 call training_set_init(fit_data%training_set,ts%displacement,ts%du_delta,& 205 & natom,ntime,ts%strain,ts%sqomega,ts%ucvol) 206 207 end subroutine fit_data_init
m_fit_data/fit_data_type [ Types ]
[ Top ] [ m_fit_data ] [ Types ]
NAME
fit_data_type
FUNCTION
SOURCE
97 type, public :: fit_data_type 98 99 integer :: ntime 100 ! Number of time in the training set 101 102 integer :: natom 103 ! Number of atoms in the training set 104 105 real(dp),allocatable :: energy_diff(:) 106 ! energy(ntime) 107 ! Array with the diffence between energy from training set and energy from initial model. 108 ! The model constains only harmonic part 109 110 real(dp),allocatable :: fcart_diff(:,:,:) 111 ! fcart_diff(3,natom,ntime) 112 ! Array with the diffence between cartesian forces from training set and forces from initial model. 113 ! The model constains only harmonic part 114 115 real(dp),allocatable :: strten_diff(:,:) 116 ! strten_diff(6,ntime) 117 ! Array with the diffence between strain from training set and strain from initial model. 118 ! The model constains only harmonic part 119 120 type(training_set_type) :: training_set 121 ! datatype with the information of the training set 122 123 end type fit_data_type 124 125 !routine for fit_data 126 public :: fit_data_compute 127 public :: fit_data_init 128 public :: fit_data_free
m_fit_data/training_set_free [ Functions ]
[ Top ] [ m_fit_data ] [ Functions ]
NAME
training_set_free
FUNCTION
Free the training_set datatype
INPUTS
training_set<training_set_type> = training_set to be free
OUTPUT
SOURCE
492 subroutine training_set_free(ts) 493 494 implicit none 495 496 !Arguments ------------------------------------ 497 !scalars 498 !arrays 499 type(training_set_type),intent(inout) :: ts 500 !Local variables------------------------------- 501 !scalar 502 !arrays 503 ! ************************************************************************* 504 505 ! Reset integer values 506 ts%ntime = 0 507 ts%natom = 0 508 509 ! Deallocate arrays 510 if(allocated(ts%displacement)) then 511 ABI_FREE(ts%displacement) 512 end if 513 if(allocated(ts%du_delta)) then 514 ABI_FREE(ts%du_delta) 515 end if 516 if(allocated(ts%strain)) then 517 ABI_FREE(ts%strain) 518 end if 519 if(allocated(ts%sqomega))then 520 ABI_FREE(ts%sqomega) 521 end if 522 if(allocated(ts%ucvol)) then 523 ABI_FREE(ts%ucvol) 524 end if 525 526 end subroutine training_set_free
m_fit_data/training_set_init [ Functions ]
[ Top ] [ m_fit_data ] [ Functions ]
NAME
training_set_init
FUNCTION
Initialize training_set datatype
INPUTS
du_delta(6,3,natom,ntime) = Variation to displacements wrt to the strain (Bohr) displacement(3,natom,ntime)= Atomic displacement wrt to the reference (Bohr) natom = number of atoms ntime = number of time step strain(6,ntime) = Strain sqomega = Sheppard and al Factors \Omega^{2} see J.Chem Phys 136, 074103 (2012) [[cite:Sheppard2012]] ucvol(ntime) = Volume of the supercell for each time (Bohr^3)
OUTPUT
ts<training_set_type> = training set to be initialized
SOURCE
438 subroutine training_set_init(ts,displacement,du_delta,natom,ntime,strain,sqomega,ucvol) 439 440 !Arguments ------------------------------------ 441 !scalars 442 integer,intent(in) :: natom,ntime 443 !arrays 444 real(dp),intent(in) :: displacement(3,natom,ntime),du_delta(6,3,natom,ntime) 445 real(dp),intent(in) :: strain(6,ntime),sqomega(ntime),ucvol(ntime) 446 type(training_set_type),intent(inout) :: ts 447 !Local variables------------------------------- 448 !scalar 449 !arrays 450 ! ************************************************************************* 451 452 !Free the output 453 call training_set_free(ts) 454 455 !Set integer values 456 ts%ntime = ntime 457 ts%natom = natom 458 459 !allocate arrays 460 ABI_MALLOC(ts%displacement,(3,natom,ntime)) 461 ts%displacement(:,:,:) = displacement(:,:,:) 462 463 ABI_MALLOC(ts%du_delta,(6,3,natom,ntime)) 464 ts%du_delta(:,:,:,:) = du_delta(:,:,:,:) 465 466 ABI_MALLOC(ts%strain,(6,ntime)) 467 ts%strain(:,:) = strain(:,:) 468 469 ABI_MALLOC(ts%sqomega,(ntime)) 470 ts%sqomega(:) = sqomega(:) 471 472 ABI_MALLOC(ts%ucvol,(ntime)) 473 ts%ucvol(:) = ucvol(:) 474 475 end subroutine training_set_init
m_fit_data/training_set_type [ Types ]
[ Top ] [ m_fit_data ] [ Types ]
NAME
training_set_type
FUNCTION
datatype for with all the information for the fit process This structure contains all the information of the training set: - ntime - displacement - du_delta - strain - sqomega - ucvol
SOURCE
51 type, public :: training_set_type 52 53 integer :: ntime 54 ! Number of time in the training set 55 56 integer :: natom 57 ! Number of atoms in the training set 58 59 real(dp), allocatable :: displacement(:,:,:) 60 ! displacement(3,natom,ntime) 61 ! displacement array, difference of the position in cartisian coordinates 62 ! between training set and reference 63 64 real(dp), allocatable :: du_delta(:,:,:,:) 65 ! du_delta(6,3,natom,ntime) 66 ! du_delta array, variation of displacement wrt to the strain 67 68 real(dp), allocatable :: strain(:,:) 69 ! strain(6,ntime) 70 ! strain array, strain in the training set 71 72 real(dp), allocatable :: sqomega(:) 73 ! sqomega(ntime) 74 ! sqomega(itime) = (((ucvol(itime)**(-2.))* ((natom_sc)**(0.5)))**(-1.0/3.0))**2 75 76 real(dp), allocatable :: ucvol(:) 77 ! ucvol(ntime) 78 ! ucvol array, volume of the cell in the training set 79 80 end type training_set_type 81 82 !routine for training_set 83 public :: training_set_init 84 public :: training_set_free