TABLE OF CONTENTS
- ABINIT/m_results_gs
- m_results_gs/copy_results_gs
- m_results_gs/destroy_results_gs
- m_results_gs/destroy_results_gs_array
- m_results_gs/init_results_gs
- m_results_gs/init_results_gs_array
- m_results_gs/results_gs_ncwrite
- m_results_gs/results_gs_type
- m_results_gs/results_gs_yaml_write
ABINIT/m_results_gs [ Modules ]
NAME
m_results_gs
FUNCTION
This module provides the definition of the results_gs_type used to store results from GS calculations.
COPYRIGHT
Copyright (C) 2011-2024 ABINIT group (MT) 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
17 #if defined HAVE_CONFIG_H 18 #include "config.h" 19 #endif 20 21 #include "abi_common.h" 22 23 MODULE m_results_gs 24 25 use defs_basis 26 use m_abicore 27 use m_xmpi 28 use m_energies 29 use m_errors 30 use m_yaml 31 use m_crystal 32 use m_stream_string 33 use m_pair_list 34 use m_nctk 35 #ifdef HAVE_NETCDF 36 use netcdf 37 #endif 38 39 use m_io_tools, only : file_exists 40 use m_fstrings, only : sjoin 41 use m_numeric_tools, only : get_trace 42 use m_geometry, only : stress_voigt_to_mat 43 44 implicit none 45 46 private
m_results_gs/copy_results_gs [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
copy_results_gs
FUNCTION
Copy a results_gs datastructure into another
INPUTS
results_gs_in=<type(results_gs_type)>=input results_gs datastructure
OUTPUT
results_gs_out=<type(results_gs_type)>=output results_gs datastructure
SOURCE
558 subroutine copy_results_gs(results_gs_in,results_gs_out) 559 560 !Arguments ------------------------------------ 561 !arrays 562 class(results_gs_type),intent(in) :: results_gs_in 563 type(results_gs_type),intent(inout) :: results_gs_out !vz_i 564 565 !Local variables------------------------------- 566 !scalars 567 integer :: natom_in,natom_out,ngrvdw_in,nspden_in,nspden_out,nsppol_in,nsppol_out 568 569 !************************************************************************ 570 571 !@results_gs_type 572 573 natom_in =results_gs_in%natom 574 natom_out =results_gs_out%natom 575 ngrvdw_in =results_gs_in%ngrvdw 576 nspden_in =results_gs_in%nspden 577 nspden_out=results_gs_out%nspden 578 nsppol_in =results_gs_in%nsppol 579 nsppol_out=results_gs_out%nsppol 580 581 if (natom_in>natom_out) then 582 ABI_SFREE(results_gs_out%fcart) 583 ABI_SFREE(results_gs_out%gred) 584 ABI_SFREE(results_gs_out%grchempottn) 585 ABI_SFREE(results_gs_out%grcondft) 586 ABI_SFREE(results_gs_out%gresid) 587 ABI_SFREE(results_gs_out%grewtn) 588 ABI_SFREE(results_gs_out%grvdw) 589 ABI_SFREE(results_gs_out%grxc) 590 ABI_SFREE(results_gs_out%intgres) 591 ABI_SFREE(results_gs_out%synlgr) 592 593 if (allocated(results_gs_in%fcart)) then 594 ABI_MALLOC(results_gs_out%fcart,(3,natom_in)) 595 end if 596 if (allocated(results_gs_in%gred)) then 597 ABI_MALLOC(results_gs_out%gred,(3,natom_in)) 598 end if 599 if (allocated(results_gs_in%gresid)) then 600 ABI_MALLOC(results_gs_out%gresid,(3,natom_in)) 601 end if 602 if (allocated(results_gs_in%grchempottn)) then 603 ABI_MALLOC(results_gs_out%grchempottn,(3,natom_in)) 604 end if 605 if (allocated(results_gs_in%grcondft)) then 606 ABI_MALLOC(results_gs_out%grcondft,(3,natom_in)) 607 end if 608 if (allocated(results_gs_in%grewtn)) then 609 ABI_MALLOC(results_gs_out%grewtn,(3,natom_in)) 610 end if 611 if (allocated(results_gs_in%grvdw)) then 612 ABI_MALLOC(results_gs_out%grvdw,(3,ngrvdw_in)) 613 end if 614 if (allocated(results_gs_in%grxc)) then 615 ABI_MALLOC(results_gs_out%grxc,(3,natom_in)) 616 end if 617 if (allocated(results_gs_in%synlgr)) then 618 ABI_MALLOC(results_gs_out%synlgr,(3,natom_in)) 619 end if 620 end if 621 622 if (nsppol_in>nsppol_out) then 623 ABI_SFREE(results_gs_out%gaps) 624 if (allocated(results_gs_in%gaps)) then 625 ABI_MALLOC(results_gs_out%gaps,(3,nsppol_in)) 626 end if 627 endif 628 629 if (nspden_in>nspden_out .or. natom_in>natom_out) then 630 ABI_SFREE(results_gs_out%intgres) 631 if (allocated(results_gs_in%intgres)) then 632 ABI_MALLOC(results_gs_out%intgres,(max(nspden_in,nspden_out),max(natom_in,natom_out))) 633 end if 634 endif 635 636 637 results_gs_out%natom =results_gs_in%natom 638 results_gs_out%ngrvdw =results_gs_in%ngrvdw 639 results_gs_out%nspden =results_gs_in%nspden 640 results_gs_out%nsppol =results_gs_in%nsppol 641 results_gs_out%berryopt=results_gs_in%berryopt 642 results_gs_out%deltae =results_gs_in%deltae 643 results_gs_out%diffor =results_gs_in%diffor 644 results_gs_out%entropy=results_gs_in%entropy 645 results_gs_out%entropy_extfpmd=results_gs_in%entropy_extfpmd 646 results_gs_out%etotal =results_gs_in%etotal 647 results_gs_out%fermie =results_gs_in%fermie 648 results_gs_out%fermih =results_gs_in%fermih 649 results_gs_out%nelect_extfpmd=results_gs_in%nelect_extfpmd 650 results_gs_out%residm =results_gs_in%residm 651 results_gs_out%res2 =results_gs_in%res2 652 results_gs_out%shiftfactor_extfpmd=results_gs_in%shiftfactor_extfpmd 653 results_gs_out%vxcavg =results_gs_in%vxcavg 654 655 call energies_copy(results_gs_in%energies,results_gs_out%energies) 656 657 results_gs_out%pel(:)=results_gs_in%pel(:) 658 results_gs_out%pion(:)=results_gs_in%pion(:) 659 results_gs_out%strten(:)=results_gs_in%strten(:) 660 661 if (allocated(results_gs_in%fcart)) results_gs_out%fcart(:,1:natom_in) =results_gs_in%fcart(:,1:natom_in) 662 if (allocated(results_gs_in%gred)) results_gs_out%gred(:,1:natom_in) =results_gs_in%gred(:,1:natom_in) 663 if (allocated(results_gs_in%gaps)) results_gs_out%gaps(:,1:nsppol_in) =results_gs_in%gaps(:,1:nsppol_in) 664 if (allocated(results_gs_in%grchempottn))& 665 & results_gs_out%grchempottn(:,1:natom_in)=results_gs_in%grchempottn(:,1:natom_in) 666 if (allocated(results_gs_in%grcondft)) results_gs_out%grcondft(:,1:natom_in)=results_gs_in%grcondft(:,1:natom_in) 667 if (allocated(results_gs_in%gresid)) results_gs_out%gresid(:,1:natom_in)=results_gs_in%gresid(:,1:natom_in) 668 if (allocated(results_gs_in%grewtn)) results_gs_out%grewtn(:,1:natom_in)=results_gs_in%grewtn(:,1:natom_in) 669 if (allocated(results_gs_in%grxc)) results_gs_out%grxc(:,1:natom_in) =results_gs_in%grxc(:,1:natom_in) 670 if (allocated(results_gs_in%intgres))results_gs_out%intgres(1:nspden_in,1:natom_in) =results_gs_in%intgres(1:nspden_in,1:natom_in) 671 if (allocated(results_gs_in%synlgr)) results_gs_out%synlgr(:,1:natom_in)=results_gs_in%synlgr(:,1:natom_in) 672 if (allocated(results_gs_in%grvdw).and.ngrvdw_in>0) then 673 results_gs_out%grvdw(:,1:ngrvdw_in)=results_gs_in%grvdw(:,1:ngrvdw_in) 674 end if 675 676 end subroutine copy_results_gs
m_results_gs/destroy_results_gs [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
destroy_results_gs
FUNCTION
Clean and destroy a results_gs datastructure
INPUTS
OUTPUT
SIDE EFFECTS
results_gs(:)=<type(results_gs_type)>=results_gs datastructure
SOURCE
449 subroutine destroy_results_gs(results_gs) 450 451 !Arguments ------------------------------------ 452 !arrays 453 type(results_gs_type),intent(inout) :: results_gs 454 455 !************************************************************************ 456 457 !@results_gs_type 458 459 results_gs%natom =0 460 results_gs%ngrvdw=0 461 results_gs%nspden=0 462 results_gs%nsppol=0 463 results_gs%berryopt=0 464 465 ABI_SFREE(results_gs%fcart) 466 ABI_SFREE(results_gs%gred) 467 ABI_SFREE(results_gs%gaps) 468 ABI_SFREE(results_gs%grcondft) 469 ABI_SFREE(results_gs%gresid) 470 ABI_SFREE(results_gs%grewtn) 471 ABI_SFREE(results_gs%grchempottn) 472 ABI_SFREE(results_gs%grvdw) 473 ABI_SFREE(results_gs%grxc) 474 ABI_SFREE(results_gs%intgres) 475 ABI_SFREE(results_gs%synlgr) 476 477 end subroutine destroy_results_gs
m_results_gs/destroy_results_gs_array [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
destroy_results_gs_array
FUNCTION
Clean and destroy a 2D-array of results_gs datastructures
INPUTS
OUTPUT
SIDE EFFECTS
results_gs(:)=<type(results_gs_type)>=results_gs datastructure 2D-array
SOURCE
498 subroutine destroy_results_gs_array(results_gs) 499 500 !Arguments ------------------------------------ 501 !arrays 502 type(results_gs_type),intent(inout) :: results_gs(:,:) 503 !Local variables------------------------------- 504 !scalars 505 integer :: ii,jj,results_gs_size1,results_gs_size2 506 507 !************************************************************************ 508 509 !@results_gs_type 510 511 results_gs_size1=size(results_gs,1) 512 results_gs_size2=size(results_gs,2) 513 if (results_gs_size1>0.and.results_gs_size2>0) then 514 515 do ii=1,results_gs_size2 516 do jj=1,results_gs_size1 517 results_gs(jj,ii)%natom =0 518 results_gs(jj,ii)%ngrvdw=0 519 results_gs(jj,ii)%nspden=0 520 results_gs(jj,ii)%nsppol=0 521 results_gs(jj,ii)%berryopt=0 522 523 ABI_SFREE(results_gs(jj,ii)%fcart) 524 ABI_SFREE(results_gs(jj,ii)%gred) 525 ABI_SFREE(results_gs(jj,ii)%gaps) 526 ABI_SFREE(results_gs(jj,ii)%grchempottn) 527 ABI_SFREE(results_gs(jj,ii)%grcondft) 528 ABI_SFREE(results_gs(jj,ii)%gresid) 529 ABI_SFREE(results_gs(jj,ii)%grewtn) 530 ABI_SFREE(results_gs(jj,ii)%grvdw) 531 ABI_SFREE(results_gs(jj,ii)%grxc) 532 ABI_SFREE(results_gs(jj,ii)%intgres) 533 ABI_SFREE(results_gs(jj,ii)%synlgr) 534 end do 535 end do 536 537 end if 538 539 end subroutine destroy_results_gs_array
m_results_gs/init_results_gs [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
init_results_gs
FUNCTION
Init all (or part of) scalars and allocatables in a results_gs datastructure
INPUTS
natom=number of atoms in cell nsppol=number of spin channels for this dataset only_part= --optional, default=false-- if this flag is activated only the following parts of results_gs are initalized: all scalars, fcart,gred,strten
OUTPUT
SIDE EFFECTS
results_gs=<type(results_gs_type)>=results_gs datastructure
SOURCE
252 subroutine init_results_gs(natom,nspden,nsppol,results_gs,only_part) 253 254 !Arguments ------------------------------------ 255 !scalars 256 integer,intent(in) :: natom,nspden,nsppol 257 logical,optional,intent(in) :: only_part 258 !arrays 259 type(results_gs_type),intent(inout) :: results_gs 260 !Local variables------------------------------- 261 !scalars 262 logical :: full_init 263 264 !************************************************************************ 265 266 !@results_gs_type 267 268 full_init=.true.;if (present(only_part)) full_init=(.not.only_part) 269 270 results_gs%berryopt=0 271 results_gs%natom =natom 272 results_gs%ngrvdw =0 273 results_gs%nspden =nspden 274 results_gs%nsppol =nsppol 275 276 results_gs%deltae =zero 277 results_gs%diffor =zero 278 results_gs%entropy=zero 279 results_gs%entropy_extfpmd=zero 280 results_gs%etotal =zero 281 results_gs%fermie =zero 282 results_gs%fermih =zero 283 results_gs%nelect_extfpmd=zero 284 results_gs%residm =zero 285 results_gs%res2 =zero 286 results_gs%shiftfactor_extfpmd=zero 287 results_gs%vxcavg =zero 288 289 call energies_init(results_gs%energies) 290 291 results_gs%strten=zero 292 ABI_MALLOC(results_gs%fcart,(3,natom)) 293 results_gs%fcart=zero 294 ABI_MALLOC(results_gs%gred,(3,natom)) 295 results_gs%gred =zero 296 ABI_MALLOC(results_gs%gaps,(3,nsppol)) 297 results_gs%gaps =zero 298 ABI_MALLOC(results_gs%intgres,(nspden,natom)) 299 results_gs%intgres=zero 300 301 if (full_init) then 302 results_gs%pel=zero 303 results_gs%pion=zero 304 305 ABI_MALLOC(results_gs%grchempottn,(3,natom)) 306 results_gs%grchempottn=zero 307 ABI_MALLOC(results_gs%grcondft,(3,natom)) 308 results_gs%grcondft=zero 309 ABI_MALLOC(results_gs%gresid,(3,natom)) 310 results_gs%gresid=zero 311 ABI_MALLOC(results_gs%grewtn,(3,natom)) 312 results_gs%grewtn=zero 313 ABI_MALLOC(results_gs%grvdw,(3,natom)) 314 results_gs%grvdw=zero 315 ABI_MALLOC(results_gs%grxc,(3,natom)) 316 results_gs%grxc =zero 317 ABI_MALLOC(results_gs%synlgr,(3,natom)) 318 results_gs%synlgr=zero 319 end if 320 321 end subroutine init_results_gs
m_results_gs/init_results_gs_array [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
init_results_gs_array
FUNCTION
Init all (or part of) scalars and allocatables in a 2D-array of results_gs datastructures
INPUTS
natom=number of atoms in cell nsppol=number of spin channels for this dataset only_part= --optional, default=false-- if this flag is activated only the following parts of results_gs are initalized: all scalars, fcart,gred,strten
OUTPUT
SIDE EFFECTS
results_gs(:)=<type(results_gs_type)>=results_gs datastructure 2Darray
SOURCE
347 subroutine init_results_gs_array(natom,nspden,nsppol,results_gs,only_part) 348 349 !Arguments ------------------------------------ 350 !scalars 351 integer,intent(in) :: natom,nspden,nsppol 352 logical,optional,intent(in) :: only_part 353 !arrays 354 type(results_gs_type),intent(inout) :: results_gs(:,:) 355 !Local variables------------------------------- 356 !scalars 357 integer :: ii,jj,results_gs_size1,results_gs_size2 358 logical :: full_init 359 !arrays 360 361 !************************************************************************ 362 363 !@results_gs_type 364 365 results_gs_size1=size(results_gs,1) 366 results_gs_size2=size(results_gs,2) 367 full_init=.true.;if (present(only_part)) full_init=(.not.only_part) 368 369 if (results_gs_size1>0.and.results_gs_size2>0) then 370 371 do ii=1,results_gs_size2 372 do jj=1,results_gs_size1 373 374 results_gs(jj,ii)%berryopt=0 375 results_gs(jj,ii)%natom =natom 376 results_gs(jj,ii)%ngrvdw =0 377 results_gs(jj,ii)%nspden =nspden 378 results_gs(jj,ii)%nsppol =nsppol 379 380 results_gs(jj,ii)%deltae =zero 381 results_gs(jj,ii)%diffor =zero 382 results_gs(jj,ii)%entropy=zero 383 results_gs(jj,ii)%entropy_extfpmd=zero 384 results_gs(jj,ii)%etotal =zero 385 results_gs(jj,ii)%fermie =zero 386 results_gs(jj,ii)%fermih =zero 387 results_gs(jj,ii)%nelect_extfpmd=zero 388 results_gs(jj,ii)%residm =zero 389 results_gs(jj,ii)%res2 =zero 390 results_gs(jj,ii)%shiftfactor_extfpmd=zero 391 results_gs(jj,ii)%vxcavg =zero 392 393 call energies_init(results_gs(jj,ii)%energies) 394 395 results_gs(jj,ii)%strten=zero 396 ABI_MALLOC(results_gs(jj,ii)%fcart,(3,natom)) 397 results_gs(jj,ii)%fcart=zero 398 ABI_MALLOC(results_gs(jj,ii)%gred,(3,natom)) 399 results_gs(jj,ii)%gred =zero 400 ABI_MALLOC(results_gs(jj,ii)%gaps,(3,nsppol)) 401 results_gs(jj,ii)%gaps =zero 402 ABI_MALLOC(results_gs(jj,ii)%intgres,(nspden,natom)) 403 results_gs(jj,ii)%intgres =zero 404 405 if (full_init) then 406 results_gs(jj,ii)%pel=zero 407 results_gs(jj,ii)%pion=zero 408 ABI_MALLOC(results_gs(jj,ii)%grchempottn,(3,natom)) 409 results_gs(jj,ii)%grchempottn=zero 410 ABI_MALLOC(results_gs(jj,ii)%grcondft,(3,natom)) 411 results_gs(jj,ii)%grcondft=zero 412 ABI_MALLOC(results_gs(jj,ii)%gresid,(3,natom)) 413 results_gs(jj,ii)%gresid=zero 414 ABI_MALLOC(results_gs(jj,ii)%grewtn,(3,natom)) 415 results_gs(jj,ii)%grewtn=zero 416 ABI_MALLOC(results_gs(jj,ii)%grxc,(3,natom)) 417 results_gs(jj,ii)%grxc =zero 418 ABI_MALLOC(results_gs(jj,ii)%grvdw,(3,results_gs(jj,ii)%ngrvdw)) 419 results_gs(jj,ii)%grvdw =zero 420 ABI_MALLOC(results_gs(jj,ii)%synlgr,(3,natom)) 421 results_gs(jj,ii)%synlgr=zero 422 end if 423 424 end do 425 end do 426 end if 427 428 end subroutine init_results_gs_array
m_results_gs/results_gs_ncwrite [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
results_gs_ncwrite
FUNCTION
INPUTS
ncid=NC file handle ecut, pawecutdg= Input cutoff energies in Ha.
OUTPUT
SOURCE
695 integer function results_gs_ncwrite(res, ncid, ecut, pawecutdg) result(ncerr) 696 697 !Arguments ------------------------------------ 698 !scalars 699 integer,intent(in) :: ncid 700 real(dp),intent(in) :: ecut,pawecutdg 701 class(results_gs_type),intent(in) :: res 702 703 !Local variables------------------------------- 704 !scalars 705 #ifdef HAVE_NETCDF 706 707 ! ************************************************************************* 708 709 ! ============================================== 710 ! === Write the dimensions specified by ETSF === 711 ! ============================================== 712 !FIXME: do not handle k_dependent = 1 713 ncerr = nctk_def_dims(ncid, [nctkdim_t("six", 6), nctkdim_t("number_of_cartesian_dimensions", 3), & 714 nctkdim_t("number_of_atoms", res%natom), nctkdim_t("number_of_spins", res%nsppol)], defmode=.True.) 715 NCF_CHECK(ncerr) 716 717 ! Define variables. 718 ! scalars passed in input (not belonging to results_gs) as well as scalars defined in results_gs 719 ncerr = nctk_def_dpscalars(ncid, [character(len=nctk_slen) :: & 720 "ecut", "pawecutdg", "deltae", "diffor", "entropy", "entropy_extfpmd", "etotal", "fermie", "fermih",& 721 & "nelect_extfpmd", "residm", "res2", "shiftfactor_extfpmd"]) 722 NCF_CHECK(ncerr) 723 724 ! arrays 725 ! 726 ! Note: unlike gred, this array has been corrected by enforcing 727 ! the translational symmetry, namely that the sum of force on all atoms is zero. 728 729 ncerr = nctk_def_arrays(ncid, [& 730 nctkarr_t('cartesian_forces', "dp", "number_of_cartesian_dimensions, number_of_atoms"),& 731 nctkarr_t('cartesian_stress_tensor', "dp", 'six')]) 732 NCF_CHECK(ncerr) 733 734 ! In case of a Berry phase calculation output the polarization 735 if (res%berryopt/=0) then 736 ncerr = nctk_def_arrays(ncid, [& 737 nctkarr_t('reduced_electronic_polarization', "dp", "number_of_cartesian_dimensions"),& 738 nctkarr_t('reduced_ionic_polarization', "dp", "number_of_cartesian_dimensions")]) 739 NCF_CHECK(ncerr) 740 end if 741 742 ! Write data 743 ! Write variables 744 ncerr = nctk_write_dpscalars(ncid, [character(len=nctk_slen) :: & 745 & 'ecut', 'pawecutdg', 'deltae', 'diffor', 'entropy', 'entropy_extfpmd', 'etotal', 'fermie', 'fermih',& 746 & 'nelect_extfpmd', 'residm', 'res2', 'shiftfactor_extfpmd'],& 747 & [ecut, pawecutdg, res%deltae, res%diffor, res%entropy, res%entropy_extfpmd, res%etotal, res%fermie, res%fermih,& 748 & res%nelect_extfpmd, res%residm, res%res2, res%shiftfactor_extfpmd],& 749 & datamode=.True.) 750 NCF_CHECK(ncerr) 751 752 NCF_CHECK(nctk_set_datamode(ncid)) 753 NCF_CHECK(nf90_put_var(ncid, vid("cartesian_forces"), res%fcart)) 754 NCF_CHECK(nf90_put_var(ncid, vid("cartesian_stress_tensor"), res%strten)) 755 756 if (res%berryopt/=0) then 757 NCF_CHECK(nf90_put_var(ncid, vid("reduced_electronic_polarization"), res%pel)) 758 NCF_CHECK(nf90_put_var(ncid, vid("reduced_ionic_polarization"), res%pion)) 759 end if 760 761 ! Add energies 762 call energies_ncwrite(res%energies, ncid) 763 764 #else 765 ABI_ERROR("netcdf support is not activated.") 766 #endif 767 768 contains 769 integer function vid(vname) 770 771 character(len=*),intent(in) :: vname 772 vid = nctk_idname(ncid, vname) 773 end function vid 774 775 end function results_gs_ncwrite
m_results_gs/results_gs_type [ Types ]
[ Top ] [ m_results_gs ] [ Types ]
NAME
results_gs_type
FUNCTION
This structured datatype contains the results of a GS calculation : energy and its decomposition, forces and their decompositions, stresses and their decompositions
SOURCE
60 type, public :: results_gs_type 61 62 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines, 63 ! declared in another part of ABINIT, that might need to take into account your modification. 64 65 ! Integer scalar 66 67 integer :: natom 68 ! The number of atoms for this dataset 69 70 integer :: nspden 71 ! The number of spin-density components for this dataset 72 73 integer :: nsppol 74 ! The number of spin channels for this dataset 75 76 integer :: ngrvdw 77 ! Size of grvdw array 78 ! Can be 0 (not allocated) or natom 79 80 integer :: berryopt 81 ! Store the Berry phase option to know whether to use pel and pion (0 means no, otherwise yes) 82 83 ! Real (real(dp)) scalars 84 85 real(dp) :: deltae 86 ! change in energy (Hartree) 87 88 real(dp) :: diffor 89 ! maximal absolute value of changes in the components of force 90 91 real(dp) :: nelect_extfpmd 92 ! Contribution of the Extended FPMD model to the number of electrons for high temperature simulations 93 94 ! All the energies are in Hartree, obtained "per unit cell". 95 type(energies_type) :: energies 96 !!! real(dp) :: eei ! local pseudopotential energy (Hartree) 97 !!! real(dp) :: eeig ! sum of eigenvalue energy (Hartree) 98 !!! real(dp) :: eew ! Ewald energy (Hartree) 99 !!! real(dp) :: ehart ! Hartree part of total energy (Hartree) 100 !!! real(dp) :: eii ! pseudopotential core-core energy 101 !!! real(dp) :: ek ! kinetic energy (Hartree) 102 !!! real(dp) :: enefield ! the term of the energy functional that depends 103 !!! ! explicitely on the electric field 104 !!! ! enefield = -ucvol*E*P 105 !!! real(dp) :: enl ! nonlocal pseudopotential energy (Hartree) 106 real(dp) :: entropy ! entropy (Hartree) 107 !!! real(dp) :: enxc ! exchange-correlation energy (Hartree) 108 !!! real(dp) :: enxcdc ! exchange-correlation double-counting energy (Hartree) 109 !!! real(dp) :: epaw ! PAW spherical energy (Hartree) 110 !!! real(dp) :: epawdc ! PAW spherical double-counting energy (Hartree) 111 real(dp) :: entropy_extfpmd ! Entropy contribution of the Extended FPMD model 112 ! for high temperature simulations 113 real(dp) :: etotal ! total energy (Hartree) 114 ! for fixed occupation numbers (occopt==0,1,or 2): 115 ! etotal=ek+ehart+enxc+eei+eew+eii+enl+PAW_spherical_part 116 ! for varying occupation numbers (occopt>=3): 117 ! etotal=ek+ehart+enxc+eei+eew+eii+enl - tsmear*entropy +PAW_spherical_part 118 real(dp) :: fermie ! Fermi energy (Hartree) 119 real(dp) :: fermih ! Fermi energy (Hartree) for excited holes in case occopt 9 120 real(dp) :: residm ! maximum value for the residual over all bands, all k points, 121 ! and all spins (Hartree or Hartree**2, to be checked !) 122 real(dp) :: res2 ! density/potential residual (squared) 123 real(dp) :: vxcavg ! Average of the exchange-correlation energy. The average 124 ! of the local psp pot and the Hartree pot is set to zero (due 125 ! to the usual problem at G=0 for Coulombic system, so vxcavg 126 ! is also the average of the local part of the Hamiltonian 127 128 ! Real (real(dp)) arrays 129 130 real(dp), allocatable :: fcart(:,:) 131 ! fcart(3,natom) 132 ! Cartesian forces (Hartree/Bohr) 133 ! Note: unlike gred, this array has been corrected by enforcing 134 ! the translational symmetry, namely that the sum of force 135 ! on all atoms is zero. 136 137 real(dp), allocatable :: gred(:,:) 138 ! gred(3,natom) 139 ! Forces in reduced coordinates (Hartree) 140 ! Actually, gradient of the total energy with respect 141 ! to change of reduced coordinates 142 143 real(dp),allocatable :: gaps(:,:) 144 ! gaps(3,nsppol) 145 ! gaps(1,:) : fundamental gap 146 ! gaps(2,:) : optical gap 147 ! gaps(3,:) : "status" for each channel: 148 ! 0.0dp if the gap was not computed (because there are only valence bands); 149 ! -1.0dp if the system (or spin-channel) is metallic; 150 ! 1.0dp if the gap was computed 151 152 real(dp), allocatable :: grchempottn(:,:) 153 ! grchempottn(3,natom) 154 ! Part of the gradient of the total energy (Hartree) with respect 155 ! to change of reduced coordinates, that comes from the spatially-varying chemical potential 156 157 real(dp), allocatable :: grcondft(:,:) 158 ! grcondft(nspden,natom) 159 ! Part of the gradient of the total energy (Hartree) with respect 160 ! to change of reduced coordinates, that comes from the constrained DFT contribution 161 162 real(dp), allocatable :: gresid(:,:) 163 ! gresid(3,natom) 164 ! Part of the gradient of the total energy (Hartree) with respect 165 ! to change of reduced coordinates, that comes from the residual 166 ! of the potential 167 168 real(dp), allocatable :: grewtn(:,:) 169 ! grewtn(3,natom) 170 ! Part of the gradient of the total energy (Hartree) with respect 171 ! to change of reduced coordinates, that comes from the Ewald energy 172 173 real(dp), allocatable :: grvdw(:,:) 174 ! grvdw(3,natom) 175 ! Part of the gradient of the total energy (Hartree) with respect 176 ! to change of reduced coordinates, that comes from 177 ! Van der Waals DFT-D dispersion (hartree) 178 ! ngrvdw can be 0 or natom 179 180 real(dp), allocatable :: grxc(:,:) 181 ! grxc(3,natom) 182 ! Part of the gradient of the total energy (Hartree) with respect 183 ! to change of reduced coordinates, that comes from the XC energy 184 185 real(dp), allocatable :: intgres(:,:) 186 ! intgres(nspden,natom) 187 ! Derivative of the total energy with respect to changes of constraints, in constrained DFT. 188 189 real(dp) :: pel(3) 190 ! ucvol times the electronic polarization in reduced coordinates 191 192 real(dp) :: pion(3) 193 ! ucvol times the ionic polarization in reduced coordinates 194 195 real(dp) :: shiftfactor_extfpmd 196 ! Energy shift factor of the Extended FPMD model for high temperature simulations 197 198 real(dp) :: strten(6) 199 ! Stress tensor in cartesian coordinates (Hartree/Bohr^3) 200 ! 6 unique components of this symmetric 3x3 tensor: 201 ! Given in order (1,1), (2,2), (3,3), (3,2), (3,1), (2,1). 202 203 real(dp), allocatable :: synlgr(:,:) 204 ! synlgr(3,natom) 205 ! Part of the gradient of the total energy (Hartree) with respect 206 ! to change of reduced coordinates, that comes from the non-local energy 207 ! The "sy" prefix refer to the fact that this gradient has been 208 ! symmetrized. 209 210 contains 211 212 procedure :: yaml_write => results_gs_yaml_write 213 ! Write the most important results in Yaml format. 214 215 end type results_gs_type 216 217 !public procedures. 218 public :: init_results_gs 219 public :: init_results_gs_array 220 public :: destroy_results_gs 221 public :: destroy_results_gs_array 222 public :: copy_results_gs 223 public :: results_gs_ncwrite
m_results_gs/results_gs_yaml_write [ Functions ]
[ Top ] [ m_results_gs ] [ Functions ]
NAME
results_gs_yaml_write
FUNCTION
Write results_gs in yaml format to unit.
INPUTS
results <type(results_gs_type)>=miscellaneous information about the system after ground state computation unit= unit of output file [cryst]: optional Crystal structure [info]: optional info for the final document [occopt]: optional Input variable occopt [with_conv]: optional True if the convergence dictionary with residuals and diffs should be written.
SOURCE
796 subroutine results_gs_yaml_write(results, unit, cryst, info, occopt, with_conv) 797 798 class(results_gs_type),intent(in) :: results 799 integer,intent(in) :: unit 800 type(crystal_t),intent(in),optional :: cryst 801 integer, intent(in),optional :: occopt 802 logical,intent(in),optional :: with_conv 803 character(len=*),intent(in),optional :: info 804 805 !Local variables------------------------------- 806 integer,parameter :: width=10 807 integer :: ii 808 type(yamldoc_t) :: ydoc 809 !arrays 810 real(dp) :: strten(3,3), abc(3), fnorms(results%natom) 811 character(len=2) :: species(results%natom) 812 813 !************************************************************************ 814 815 if (unit == dev_null) return 816 817 if (present(info)) then 818 ydoc = yamldoc_open('ResultsGS', info=info, width=width) 819 else 820 ydoc = yamldoc_open('ResultsGS', width=width) 821 end if 822 823 ! Write lattice parameters 824 if (present(cryst)) then 825 call ydoc%add_real2d('lattice_vectors', cryst%rprimd, real_fmt="(f11.7)") 826 !ori abc = [(sqrt(sum(cryst%rprimd(:, ii) ** 2)), ii=1,3)] 827 !replace the implicit loop by an explicit one 828 !workaround works with both ifort and ifx on oneapi 2024 829 do ii=1,3 830 abc(ii) = sqrt(sum(cryst%rprimd(:, ii) ** 2)) 831 end do 832 call ydoc%add_real1d('lattice_lengths', abc, real_fmt="(f10.5)") 833 call ydoc%add_real1d('lattice_angles', cryst%angdeg, real_fmt="(f7.3)", comment="degrees, (23, 13, 12)") 834 call ydoc%add_real('lattice_volume', cryst%ucvol + tol10, real_fmt="(es15.7)") 835 endif 836 837 ! Write convergence degree. 838 ! It seems there's a portability problem for residm computed with nstep = 0 and iscf -3 839 ! because one may get very small value e.g. 7.91-323. residm with nstep > 0 are OK though 840 ! so print zero if residm < tol30 or allow the caller not to write the convergence dict. 841 if(present(with_conv))then 842 if (with_conv) then 843 call ydoc%add_reals( & 844 "deltae, res2, residm, diffor", & 845 [results%deltae, results%res2, merge(results%residm, zero, results%residm > tol30), results%diffor], & 846 real_fmt="(es10.3)", dict_key="convergence") 847 else 848 call ydoc%set_keys_to_string("deltae, res2, residm, diffor", "null", dict_key="convergence") 849 end if 850 else 851 call ydoc%set_keys_to_string("deltae, res2, residm, diffor", "null", dict_key="convergence") 852 endif 853 854 ! Write energies. 855 if (present(occopt))then 856 if (occopt == 9) then 857 call ydoc%add_reals("etotal, entropy, fermie, fermih", [results%etotal, results%entropy, results%fermie, results%fermih]) 858 else 859 call ydoc%add_reals("etotal, entropy, fermie", [results%etotal, results%entropy, results%fermie]) 860 endif 861 else 862 call ydoc%add_reals("etotal, entropy, fermie", [results%etotal, results%entropy, results%fermie]) 863 endif 864 865 ! Cartesian stress tensor and forces. 866 call stress_voigt_to_mat(results%strten, strten) 867 if (strten(1,1) /= MAGIC_UNDEF) then 868 call ydoc%add_real2d('cartesian_stress_tensor', strten, comment="hartree/bohr^3") 869 call ydoc%add_real('pressure_GPa', - get_trace(strten) * HaBohr3_GPa / three, real_fmt="(es12.4)") 870 else 871 call ydoc%set_keys_to_string("cartesian_stress_tensor, pressure_GPa", "null") 872 end if 873 874 if(present(cryst))then 875 species = [(cryst%symbol_iatom(ii), ii=1,cryst%natom)] 876 call ydoc%add_real2d('xred', cryst%xred, slist=species, real_fmt="(es12.4)") 877 endif 878 879 if (results%fcart(1,1) /= MAGIC_UNDEF) then 880 call ydoc%add_real2d('cartesian_forces', results%fcart, comment="hartree/bohr") 881 fnorms = [(sqrt(sum(results%fcart(:, ii) ** 2)), ii=1,results%natom)] 882 ! Write force statistics 883 call ydoc%add_reals('min, max, mean', & 884 values=[minval(fnorms), maxval(fnorms), sum(fnorms) / results%natom], dict_key="force_length_stats") 885 else 886 ! Set entries to null (python None) to facilitate life to the parsing routines! 887 call ydoc%add_string('cartesian_forces', "null") 888 call ydoc%set_keys_to_string("min, max, mean", "null", dict_key="force_length_stats") 889 end if 890 891 call ydoc%write_and_free(unit) 892 893 end subroutine results_gs_yaml_write