TABLE OF CONTENTS


ABINIT/clnup1 [ Functions ]

[ Top ] [ Functions ]

NAME

 clnup1

FUNCTION

 Perform "cleanup" at end of execution of gstate routine.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  acell(3)=length scales of primitive translations (bohr)
  dosdeltae=DOS delta of Energy
  dtset <type(dataset_type)>=all input variables in this dataset
  eigen(mband*nkpt*nsppol)=eigenvalues (hartree) for all bands
                           at each k point
  enunit=choice for units of output eigenvalues: 0=>hartree,
   1=> eV, 2=> hartree and eV
  fermie=fermi energy (Hartree)
  fnameabo_dos=filename of output DOS file
  fnameabo_eig=filename of output EIG file
  fred(3,natom)=d(E)/d(xred) (hartree)
  iatfix(3,natom)=0 if not fixed along specified direction,
                  1 if fixed
  iscf=parameter controlling scf or non-scf choice
  kptopt=option for the generation of k points
  kptns(3,nkpt)=k points in terms of recip primitive translations
  mband=maximum number of bands
  mpi_enreg=information about MPI parallelization
  natom=number of atoms in unit cell
  nband(nkpt*nsppol)=number of bands
  nfft=(effective) number of FFT grid points (for this processor)
  ngfft(18)=contain all needed information about 3D FFT,
            see ~abinit/doc/input_variables/vargs.htm#ngfft
  nkpt=number of k points
  nspden=number of spin-density components
  nsppol=1 for unpolarized, 2 for spin-polarized
  nstep=desired number of electron iteration steps
  occ(maxval(nband(:))*nkpt*nsppol)=occupancies for each band and k point
  occopt=option for occupancies
  prtdos= if == 1, will print the density of states
  prtfor= if >0, will print the forces
  prtstm= input variable prtstm
  prtvol=control print volume and debugging
  resid(mband*nkpt*nsppol)=squared residuals for each band and k point where
                     resid(n,k)=|<C(n,k)|(H-e(n,k))|C(n,k)>|^2
  rhor(nfft,nspden)=electron density (electrons/bohr^3)
  rprimd(3,3)=dimensional real space primitive translations (bohr)
  tphysel="physical" electronic temperature with FD occupations
  tsmear=smearing energy or temperature (if metal)
  vxcavg=average of vxc potential
  wtk(nkpt)=real(dp) array of k-point weights
  xred(3,natom)=reduced atomic coordinates

OUTPUT

  (only print and write to disk)

PARENTS

      gstate

CHILDREN

      getnel,metric,prteigrs,prtrhomxmn,prtxf,write_eig,wrtout

SOURCE

 72 #if defined HAVE_CONFIG_H
 73 #include "config.h"
 74 #endif
 75 
 76 #include "abi_common.h"
 77 
 78 subroutine clnup1(acell,dtset,eigen,fermie,&
 79   & fnameabo_dos,fnameabo_eig,fred,&
 80   & mpi_enreg,nfft,ngfft,occ,prtfor,&
 81   & resid,rhor,rprimd,vxcavg,xred)
 82 
 83  use defs_basis
 84  use defs_abitypes
 85  use defs_wvltypes
 86  use m_profiling_abi
 87  use m_xmpi
 88  use m_errors
 89 
 90  use m_io_tools,     only : open_file
 91 
 92 !This section has been created automatically by the script Abilint (TD).
 93 !Do not modify the following lines by hand.
 94 #undef ABI_FUNC
 95 #define ABI_FUNC 'clnup1'
 96  use interfaces_14_hidewrite
 97  use interfaces_41_geometry
 98  use interfaces_59_ionetcdf
 99  use interfaces_61_occeig
100  use interfaces_67_common, except_this_one => clnup1
101 !End of the abilint section
102 
103  implicit none
104 
105 !Arguments ------------------------------------
106 !scalars
107  integer,intent(in) :: nfft
108  integer,intent(in) :: prtfor
109  real(dp),intent(in) :: fermie
110  real(dp),intent(in) :: vxcavg
111  character(len=*),intent(in) :: fnameabo_dos,fnameabo_eig
112  type(dataset_type),intent(in) :: dtset
113  type(MPI_type),intent(in) :: mpi_enreg
114 !arrays
115  integer,intent(in)  :: ngfft(18)
116  real(dp),intent(in) :: acell(3)
117  real(dp),intent(in) :: eigen(dtset%mband*dtset%nkpt*dtset%nsppol)
118  real(dp),intent(in) :: fred(3,dtset%natom)
119  real(dp),intent(in) :: resid(dtset%mband*dtset%nkpt*dtset%nsppol)
120  real(dp),intent(in) :: rhor(nfft,dtset%nspden)
121  real(dp),intent(in) :: rprimd(3,3)
122  real(dp),intent(in) :: xred(3,dtset%natom)
123  real(dp),intent(inout) :: occ(dtset%mband*dtset%nkpt*dtset%nsppol)
124 
125 !Local variables-------------------------------
126 !scalars
127  integer,parameter :: master=0
128  integer :: comm,iatom,ii,iscf_dum,iwfrc,me,nnonsc,option,unitdos
129  real(dp) :: entropy,grmax,grsum,maxocc,nelect,tolwf,ucvol
130  real(dp) :: gmet(3,3),gprimd(3,3),rmet(3,3)
131  character(len=500) :: message
132  character(len=fnlen) filename
133 !arrays
134  real(dp),allocatable :: doccde(:)
135 
136 ! ****************************************************************
137 
138  comm=mpi_enreg%comm_cell; me=xmpi_comm_rank(comm)
139 
140  if(dtset%prtstm==0)then ! Write reduced coordinates xred
141    write(message, '(a,i5,a)' )' reduced coordinates (array xred) for',dtset%natom,' atoms'
142    call wrtout(ab_out,message,'COLL')
143    do iatom=1,dtset%natom
144      write(message, '(1x,3f20.12)' ) xred(:,iatom)
145      call wrtout(ab_out,message,'COLL')
146    end do
147  end if
148 
149 !Write reduced gradients if iscf > 0 and dtset%nstep>0 and
150  if (dtset%iscf>=0.and.dtset%nstep>0.and.dtset%prtstm==0) then
151 
152 !  Compute absolute maximum and root mean square value of gradients
153    grmax=0.0_dp
154    grsum=0.0_dp
155    do iatom=1,dtset%natom
156      do ii=1,3
157 !      To be activated in v5.5
158 !      grmax=max(grmax,abs(fred(ii,iatom)))
159        grmax=max(grmax,fred(ii,iatom))
160        grsum=grsum+fred(ii,iatom)**2
161      end do
162    end do
163    grsum=sqrt(grsum/dble(3*dtset%natom))
164 
165    write(message, '(1x,a,1p,e12.4,a,e12.4,a)' )'rms dE/dt=',grsum,'; max dE/dt=',grmax,'; dE/dt below (all hartree)'
166    call wrtout(ab_out,message,'COLL')
167    do iatom=1,dtset%natom
168      write(message, '(i5,1x,3f20.12)' ) iatom,fred(1:3,iatom)
169      call wrtout(ab_out,message,'COLL')
170    end do
171 
172  end if
173 
174  if(dtset%prtstm==0)then
175 
176 !  Compute and write out dimensional cartesian coords and forces:
177    call wrtout(ab_out,' ','COLL')
178 
179 !  (only write forces if iscf > 0 and dtset%nstep>0)
180    if (dtset%iscf<0.or.dtset%nstep<=0.or.prtfor==0) then
181      iwfrc=0
182    else
183      iwfrc=1
184    end if
185 
186    call prtxf(fred,dtset%iatfix,ab_out,iwfrc,dtset%natom,rprimd,xred)
187 
188 !  Write length scales
189    write(message, '(1x,a,3f16.12,a)' )'length scales=',acell,' bohr'
190    call wrtout(ab_out,message,'COLL')
191    write(message, '(14x,a,3f16.12,a)' )'=',Bohr_Ang*acell(1:3),' angstroms'
192    call wrtout(ab_out,message,'COLL')
193 
194  end if
195 
196  option=1; nnonsc=0; tolwf=0.0_dp
197 
198  if(dtset%iscf<0 .and. dtset%iscf/=-3)option=3
199  iscf_dum=dtset%iscf
200  if(dtset%nstep==0)iscf_dum=-1
201 
202  if(dtset%tfkinfunc==0)then
203    call prteigrs(eigen,dtset%enunit,fermie,fnameabo_eig,ab_out,&
204 &   iscf_dum,dtset%kptns,dtset%kptopt,dtset%mband,&
205 &   dtset%nband,dtset%nkpt,nnonsc,dtset%nsppol,occ,&
206 &   dtset%occopt,option,dtset%prteig,dtset%prtvol,resid,tolwf,&
207 &   vxcavg,dtset%wtk)
208    call prteigrs(eigen,dtset%enunit,fermie,fnameabo_eig,std_out,&
209 &   iscf_dum,dtset%kptns,dtset%kptopt,dtset%mband,&
210 &   dtset%nband,dtset%nkpt,nnonsc,dtset%nsppol,occ,&
211 &   dtset%occopt,option,dtset%prteig,dtset%prtvol,resid,tolwf,&
212 &   vxcavg,dtset%wtk)
213 
214 #if defined HAVE_NETCDF
215    if (dtset%prteig==1 .and. me == master) then
216      filename=trim(fnameabo_eig)//'.nc'
217      call write_eig(eigen,filename,dtset%kptns,dtset%mband,dtset%nband,dtset%nkpt,dtset%nsppol)
218    end if
219 #endif
220 
221  end if
222 
223 !Compute and print location of maximal and minimal density
224  call metric(gmet,gprimd,-1,rmet,rprimd,ucvol)
225  call prtrhomxmn(std_out,mpi_enreg,nfft,ngfft,dtset%nspden,2,rhor,ucvol=ucvol)
226  if( dtset%prtvol>1)then
227    call prtrhomxmn(ab_out,mpi_enreg,nfft,ngfft,dtset%nspden,2,rhor,ucvol=ucvol)
228  end if
229 
230 !If needed, print DOS (unitdos is closed in getnel, occ is not changed if option == 2
231  if (dtset%prtdos==1 .and. me == master) then
232    if (open_file(fnameabo_dos,message, newunit=unitdos, status='unknown', action="write", form='formatted') /= 0) then
233      MSG_ERROR(message)
234    end if
235    rewind(unitdos)
236    maxocc=two/(dtset%nspinor*dtset%nsppol)  ! Will not work in the fixed moment case
237    option=2
238    ABI_ALLOCATE(doccde,(dtset%mband*dtset%nkpt*dtset%nsppol))
239    call getnel(doccde,dtset%dosdeltae,eigen,entropy,fermie,&
240 &   maxocc,dtset%mband,dtset%nband,nelect,dtset%nkpt,&
241 &   dtset%nsppol,occ,dtset%occopt,option,dtset%tphysel,&
242 &   dtset%tsmear,unitdos,dtset%wtk)
243    ABI_DEALLOCATE(doccde)
244  end if
245 
246 end subroutine clnup1