TABLE OF CONTENTS
ABINIT/prtimg [ Functions ]
NAME
prtimg
FUNCTION
Print out results obtained by as ground-state calculation of an image of the cell. The printing format is condensed in order to facilitate the reading.
COPYRIGHT
Copyright (C) 2011-2018 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 . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
dynimage(nimagetot)=flags defining static/dynamic state of images imagealgo_str=name of the algorithm (with images) used imgmov=index of algorithm (with images) used iout=unit number for output mpi_enreg=MPI-parallelisation information nimage=number of images stored on current proc nimage_tot=total number of images (should be dtset%nimage) prt_all_images=true if all images have to be printed out (ignoring dynimage) prtvolimg=printing volume for each image <0 : nothing 0 : only a title 1 : energy, residuals, forces, stresses, velocities, atomic positions 2 : energy, residuals resimg(nimage) <type(results_img_type)>=results of the ground-state computations for all images treated by current proc
OUTPUT
(data written to unit iout)
PARENTS
gstateimg
CHILDREN
destroy_results_img,gather_results_img,metric,prtxvf,wrtout,xred2xcart
SOURCE
48 #if defined HAVE_CONFIG_H 49 #include "config.h" 50 #endif 51 52 #include "abi_common.h" 53 54 subroutine prtimg(dynimage,imagealgo_str,imgmov,iout,mpi_enreg,nimage,nimage_tot,& 55 & prt_all_images,prtvolimg,resimg) 56 57 use defs_basis 58 use defs_abitypes 59 use m_results_img 60 use m_errors 61 use m_profiling_abi 62 63 use m_geometry, only : xred2xcart, metric 64 65 !This section has been created automatically by the script Abilint (TD). 66 !Do not modify the following lines by hand. 67 #undef ABI_FUNC 68 #define ABI_FUNC 'prtimg' 69 use interfaces_14_hidewrite 70 use interfaces_67_common, except_this_one => prtimg 71 !End of the abilint section 72 73 implicit none 74 75 !Arguments ------------------------------------ 76 !scalars 77 integer,intent(in) :: nimage_tot,dynimage(nimage_tot),imgmov,iout,nimage,prtvolimg !vz_d 78 logical,intent(in) :: prt_all_images 79 character(len=60),intent(in) :: imagealgo_str 80 type(MPI_type),intent(in) :: mpi_enreg 81 !arrays 82 type(results_img_type),target,intent(inout) :: resimg(nimage) 83 84 !Local variables------------------------------- 85 !scalars 86 integer :: ii,prtvel 87 logical :: test_img 88 real(dp) :: ucvol_img 89 character(len=500) :: msg 90 !arrays 91 integer,allocatable :: iatfix_img(:,:) 92 real(dp),allocatable :: gmet_img(:,:),gprimd_img(:,:),rmet_img(:,:),xcart_img(:,:) 93 type(results_img_type),pointer :: resimg_all(:) 94 95 ! **************************************************************** 96 97 DBG_ENTER('COLL') 98 99 if (prtvolimg<=0) return 100 if (mpi_enreg%me_cell/=0) return 101 102 !Gather data 103 if (prtvolimg==1.or.prtvolimg==2) then 104 test_img=(nimage_tot/=1.and.mpi_enreg%paral_img==1) 105 if (test_img) then 106 if (mpi_enreg%me==0) then 107 ABI_DATATYPE_ALLOCATE(resimg_all,(nimage_tot)) 108 end if 109 call gather_results_img(mpi_enreg,resimg,resimg_all,master=0,& 110 & allgather=.false.,only_one_per_img=.true.) 111 else 112 resimg_all => resimg 113 end if 114 end if 115 116 !===== First option for the printing volume === 117 if (prtvolimg==1.and.mpi_enreg%me==0) then 118 119 prtvel=0;if (imgmov==0.or.imgmov==9.or.imgmov==10.or.imgmov==13) prtvel=1 120 121 do ii=1,nimage_tot 122 if (dynimage(ii)==1.or.prt_all_images) then 123 124 ! Title 125 write(msg,'(6a,i4,a,i4,2a)') ch10,& 126 & '----------------------------------------------------------------------',ch10,& 127 & ' ',trim(imagealgo_str),' - CELL # ',ii,'/',nimage_tot,ch10,& 128 & '----------------------------------------------------------------------' 129 call wrtout(iout,msg,'COLL') 130 131 ! Total energy 132 write(msg,'(2a,es20.12)') ch10,' Total energy for the cell [Ha]: ',resimg_all(ii)%results_gs%etotal 133 call wrtout(iout,msg,'COLL') 134 135 ! Residuals of the SCF cycle 136 write(msg,'(3a,4(a,es16.8,a))') ch10,& 137 & ' Residuals from SCF cycle: ',ch10,& 138 & ' Total energy difference =',resimg_all(ii)%results_gs%deltae,ch10,& 139 & ' Maximal forces difference =',resimg_all(ii)%results_gs%diffor,ch10,& 140 & ' Max. residual of wave-functions=',resimg_all(ii)%results_gs%residm,ch10,& 141 & ' Density/potential residual (^2)=',resimg_all(ii)%results_gs%res2,ch10 142 call wrtout(iout,msg,'COLL') 143 144 ! Cell parameters 145 ABI_ALLOCATE(rmet_img,(3,3)) 146 ABI_ALLOCATE(gmet_img,(3,3)) 147 ABI_ALLOCATE(gprimd_img,(3,3)) 148 call metric(gmet_img,gprimd_img,iout,rmet_img,resimg_all(ii)%rprim,ucvol_img) 149 ABI_DEALLOCATE(rmet_img) 150 ABI_DEALLOCATE(gmet_img) 151 ABI_DEALLOCATE(gprimd_img) 152 153 ! Positions, forces and velocities 154 ABI_ALLOCATE(iatfix_img,(3,resimg_all(ii)%natom)) 155 ABI_ALLOCATE(xcart_img,(3,resimg_all(ii)%natom)) 156 iatfix_img=0 157 call xred2xcart(resimg_all(ii)%natom,resimg_all(ii)%rprim,xcart_img,resimg_all(ii)%xred) 158 call prtxvf(resimg_all(ii)%results_gs%fcart,resimg_all(ii)%results_gs%fred,& 159 & iatfix_img,iout,resimg_all(ii)%natom,prtvel,& 160 & resimg_all(ii)%vel,xcart_img,resimg_all(ii)%xred) 161 ABI_DEALLOCATE(iatfix_img) 162 ABI_DEALLOCATE(xcart_img) 163 164 ! Stress tensor 165 write(msg, '(a,es12.4,a)' ) & 166 & '-Cartesian components of stress tensor (GPa) [Pressure=',& 167 & -(resimg_all(ii)%results_gs%strten(1)+resimg_all(ii)%results_gs%strten(2) & 168 & +resimg_all(ii)%results_gs%strten(3))*HaBohr3_GPa/three,' GPa]' 169 call wrtout(iout,msg,'COLL') 170 write(msg, '(2(a,1p,e16.8))' ) '- sigma(1 1)=',resimg_all(ii)%results_gs%strten(1)*HaBohr3_GPa,& 171 & ' sigma(3 2)=',resimg_all(ii)%results_gs%strten(4)*HaBohr3_GPa 172 call wrtout(iout,msg,'COLL') 173 write(msg, '(2(a,1p,e16.8))' ) '- sigma(2 2)=',resimg_all(ii)%results_gs%strten(2)*HaBohr3_GPa,& 174 & ' sigma(3 1)=',resimg_all(ii)%results_gs%strten(5)*HaBohr3_GPa 175 call wrtout(iout,msg,'COLL') 176 write(msg, '(2(a,1p,e16.8))' ) '- sigma(3 3)=',resimg_all(ii)%results_gs%strten(3)*HaBohr3_GPa,& 177 & ' sigma(2 1)=',resimg_all(ii)%results_gs%strten(6)*HaBohr3_GPa 178 call wrtout(iout,msg,'COLL') 179 end if 180 end do 181 end if 182 183 184 !===== 2nd option for the printing volume === 185 if (prtvolimg==2.and.mpi_enreg%me==0) then 186 write(msg,'(a,1x,a)') ch10,& 187 & 'Cell Total_energy[Ha] deltae diffor residm res2' 188 call wrtout(iout,msg,'COLL') 189 do ii=1,nimage_tot 190 if (dynimage(ii)==1.or.prt_all_images) then 191 write(msg,'(1x,i4,2x,es16.8,4(1x,es13.5))') & 192 & ii,resimg_all(ii)%results_gs%etotal,resimg_all(ii)%results_gs%deltae,& 193 & resimg_all(ii)%results_gs%diffor,resimg_all(ii)%results_gs%residm,& 194 & resimg_all(ii)%results_gs%res2 195 call wrtout(iout,msg,'COLL') 196 end if 197 end do 198 end if 199 200 !===== 201 if (prtvolimg==1.or.prtvolimg==2) then 202 if (test_img.and.mpi_enreg%me==0) then 203 call destroy_results_img(resimg_all) 204 ABI_DATATYPE_DEALLOCATE(resimg_all) 205 end if 206 nullify(resimg_all) 207 end if 208 209 DBG_EXIT('COLL') 210 211 end subroutine prtimg