TABLE OF CONTENTS


ABINIT/gettag [ Functions ]

[ Top ] [ Functions ]

NAME

 gettag

FUNCTION

 Set the tag associated to each atom,

INPUTS

 prtallatoms = Logical for PRTint ALL ATOMS
 atlist      = ATom LIST
 index       = index for each atom
 natom       = Number of ATOMs

OUTPUT

  tag = The string to put for each atom

PARENTS

      prtxfase

CHILDREN

      gettag,wrtout

SOURCE

370 subroutine gettag(atlist,index,natom,prtallatoms,tag)
371 
372 
373 !This section has been created automatically by the script Abilint (TD).
374 !Do not modify the following lines by hand.
375 #undef ABI_FUNC
376 #define ABI_FUNC 'gettag'
377 !End of the abilint section
378 
379 implicit none
380 
381 !Arguments ------------------------------------
382 !scalars
383   logical,intent(in) :: prtallatoms
384   integer,intent(in) :: natom
385   logical,intent(in) :: atlist(natom)
386   integer,intent(in) :: index
387   character(len=7),intent(out)   :: tag
388 
389 !Local variables -------------------------
390 
391 ! *********************************************************************
392 !The numbering will be from (1) to (9999)
393 
394    if (prtallatoms)then
395      tag=''
396    elseif (atlist(index)) then
397      if (natom<10) then
398        write(tag, '(a,I1.1,a)') ' (',index,')'
399      elseif (natom<100) then
400        write(tag, '(a,I2.2,a)') ' (',index,')'
401      elseif (natom<1000) then
402        write(tag, '(a,I3.3,a)') ' (',index,')'
403      elseif (natom<10000) then
404        write(tag, '(a,I4.4,a)') ' (',index,')'
405      end if
406    end if
407 
408  end subroutine gettag

ABINIT/prtnatom [ Functions ]

[ Top ] [ Functions ]

NAME

 prtnatom

FUNCTION

 Print information for N atoms

INPUTS

 prtallatoms = Logical for PRTint ALL ATOMS
 atlist      = ATom LIST
 index       = index for each atom
 natom       = Number of ATOMs

OUTPUT

  tag = The string to put for aech atom

PARENTS

      prtxfase

CHILDREN

      gettag,wrtout

SOURCE

438 subroutine prtnatom(atlist,iout,message,natom,prtallatoms,thearray)
439 
440 
441 !This section has been created automatically by the script Abilint (TD).
442 !Do not modify the following lines by hand.
443 #undef ABI_FUNC
444 #define ABI_FUNC 'prtnatom'
445  use interfaces_14_hidewrite
446 !End of the abilint section
447 
448 implicit none
449 
450 !Arguments ------------------------------------
451 !scalars
452   logical,intent(in) :: prtallatoms
453   integer,intent(in) :: natom
454   logical,intent(in) :: atlist(natom)
455   integer,intent(in) :: iout
456   character(len=*),intent(inout) :: message
457 !arrays
458   real(dp) :: thearray(3,natom)
459 
460 !Local variables-------------------------------
461 !scalars
462   integer :: kk
463   character(len=7)   :: tag ! Maximal ' (9999)'
464   character(len=18)   :: fmt
465 
466 ! *********************************************************************
467 
468    fmt='(a,a,1p,3e22.14,a)'
469 
470    do kk=1,natom
471 
472      if (atlist(kk)) then
473        call gettag(atlist,kk,natom,prtallatoms,tag)
474        write(message,fmt)&
475 &       TRIM(message),ch10,&
476 &       thearray(:,kk),&
477 &       tag
478      end if
479 
480    end do
481  !MGNAG
482  ! Runtime Error: wrtout_cpp.f90, line 896: Buffer overflow on output
483    call wrtout(iout,message,'COLL')
484 
485  end subroutine prtnatom

ABINIT/prtxfase [ Functions ]

[ Top ] [ Functions ]

NAME

 prtxfase

FUNCTION

 Print the values of xcart (X), forces (F)
 acell (A), Stresses (S), and energy (E)
 All values come from the history hist
 Also compute and print max and rms forces.
 Also compute absolute and relative differences
 with previous calculation

COPYRIGHT

 Copyright (C) 1998-2018 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

 ab_mover<type abimover>=Subset of dtset only related with
          |                 movement of ions and acell, contains:
          | dtion:  Time step
          ! natom:  Number of atoms
          | vis:    viscosity
          | iatfix: Index of atoms and directions fixed
          | amass:  Mass of ions
 hist<type abihist>=Historical record of positions, forces,
                               stresses, cell and energies,
 itime= time step
 iout=unit number for printing

OUTPUT

  (only writing)

PARENTS

      mover

CHILDREN

      gettag,wrtout

SOURCE

 47 #if defined HAVE_CONFIG_H
 48 #include "config.h"
 49 #endif
 50 
 51 #include "abi_common.h"
 52 
 53 subroutine prtxfase(ab_mover,hist,itime,iout,pos)
 54 
 55  use defs_basis
 56  use m_profiling_abi
 57  use m_abimover
 58  use m_abihist
 59 
 60  use m_geometry,  only : fcart2fred, xred2xcart
 61 
 62 !This section has been created automatically by the script Abilint (TD).
 63 !Do not modify the following lines by hand.
 64 #undef ABI_FUNC
 65 #define ABI_FUNC 'prtxfase'
 66  use interfaces_14_hidewrite
 67 !End of the abilint section
 68 
 69 implicit none
 70 
 71 !Arguments ------------------------------------
 72 !scalars
 73  type(abimover),intent(in) :: ab_mover
 74  type(abihist),intent(in),target :: hist
 75  integer,intent(in) :: itime,iout
 76  integer,intent(in) :: pos
 77 !arrays
 78 
 79 !Local variables-------------------------------
 80 !scalars
 81  integer :: jj,kk,unfixd,iprt
 82  real(dp) :: val_max,val_rms,ucvol ! Values maximal and RMS, Volume of Unitary cell
 83  real(dp) :: dEabs,dErel ! Diff of energy absolute and relative
 84  real(dp) :: ekin
 85  real(dp) :: angle(3),rmet(3,3)
 86 !character(len=80*(max(ab_mover%natom,3)+1)) :: message
 87 !MGNAG: This is not very safe. One should use line-based output istead of appending chars
 88 ! and then outputting everything! For the time being I use this temporary hack to solve the problem with NAG
 89  character(len=max(80*(max(ab_mover%natom,3)+1),50000)) :: message
 90  character(len=18)   :: fmt1
 91  logical :: prtallatoms
 92 !arrays
 93  logical :: atlist(ab_mover%natom)
 94  real(dp),allocatable :: fred(:,:),xcart(:,:)
 95  real(dp),pointer :: acell(:),fcart(:,:),rprimd(:,:),strten(:),vel(:,:),xred(:,:)
 96 
 97 ! ***********************************************************
 98 
 99  fmt1='(a,a,1p,3e22.14)'
100 
101 !##########################################################
102 !### 1. Organize list of atoms to print
103 
104  prtallatoms=.TRUE.
105  do kk=1,ab_mover%natom
106    if (ab_mover%prtatlist(kk)/=kk) prtallatoms=.FALSE.
107  end do
108 
109  atlist(:)=.FALSE.
110  do iprt=1,ab_mover%natom
111    if (ab_mover%prtatlist(iprt)>0.and.ab_mover%prtatlist(iprt)<=ab_mover%natom) atlist(ab_mover%prtatlist(iprt))=.TRUE.
112  end do
113 
114  acell  => hist%acell(:,hist%ihist)
115  rprimd => hist%rprimd(:,:,hist%ihist)
116  xred   => hist%xred(:,:,hist%ihist)
117  fcart  => hist%fcart(:,:,hist%ihist)
118  strten => hist%strten(:,hist%ihist)
119  vel    => hist%vel(:,:,hist%ihist)
120 
121 !###########################################################
122 !### 1. Positions
123 
124  ABI_ALLOCATE(xcart,(3,ab_mover%natom))
125  call xred2xcart(ab_mover%natom,rprimd,xcart,xred)
126 
127  write(message, '(a,a)' )&
128 & ch10,' Cartesian coordinates (xcart) [bohr]'
129  call prtnatom(atlist,iout,message,ab_mover%natom,prtallatoms,xcart)
130 
131  write(message, '(a)' )&
132 & ' Reduced coordinates (xred)'
133  call prtnatom(atlist,iout,message,ab_mover%natom,prtallatoms,xred)
134 
135  ABI_DEALLOCATE(xcart)
136 
137 !###########################################################
138 !### 2. Forces
139 
140  if(pos==mover_AFTER)then
141 
142    ABI_ALLOCATE(fred,(3,ab_mover%natom))
143    call fcart2fred(fcart,fred,rprimd,ab_mover%natom)
144 
145 !  Compute max |f| and rms f,
146 !  EXCLUDING the components determined by iatfix
147    val_max=0.0_dp
148    val_rms=0.0_dp
149    unfixd=0
150    do kk=1,ab_mover%natom
151      do jj=1,3
152        if (ab_mover%iatfix(jj,kk) /= 1) then
153          unfixd=unfixd+1
154          val_rms=val_rms+fcart(jj,kk)**2
155          val_max=max(val_max,abs(fcart(jj,kk)**2))
156        end if
157      end do
158    end do
159    if ( unfixd /= 0 ) val_rms=sqrt(val_rms/dble(unfixd))
160 
161    write(message, '(a,1p,2e12.5,a)' ) &
162 &   ' Cartesian forces (fcart) [Ha/bohr]; max,rms=',&
163 &   sqrt(val_max),val_rms,' (free atoms)'
164    call prtnatom(atlist,iout,message,ab_mover%natom,prtallatoms,fcart)
165 
166    write(message, '(a)' )&
167 &   ' Reduced forces (fred)'
168    call prtnatom(atlist,iout,message,ab_mover%natom,prtallatoms,fred)
169 
170    ABI_DEALLOCATE(fred)
171 
172  end if
173 
174 !###########################################################
175 !### 3. Velocities
176 
177 !Only if the velocities are being used
178  if (hist%isVused)then
179 !  Only if velocities are recorded in a history
180    if (allocated(hist%vel))then
181 !    Compute max |v| and rms v,
182 !    EXCLUDING the components determined by iatfix
183      val_max=0.0_dp
184      val_rms=0.0_dp
185      unfixd=0
186      do kk=1,ab_mover%natom
187        do jj=1,3
188          if (ab_mover%iatfix(jj,kk) /= 1) then
189            unfixd=unfixd+1
190            val_rms=val_rms+vel(jj,kk)**2
191            val_max=max(val_max,abs(vel(jj,kk)**2))
192          end if
193        end do
194      end do
195      if ( unfixd /= 0 ) val_rms=sqrt(val_rms/dble(unfixd))
196 
197      write(message, '(a,1p,2e12.5,a)' ) &
198 &     ' Cartesian velocities (vel) [bohr*Ha/hbar]; max,rms=',&
199 &     sqrt(val_max),val_rms,' (free atoms)'
200      call prtnatom(atlist,iout,message,ab_mover%natom,prtallatoms,vel)
201 
202 !    Compute the ionic kinetic energy (no cell shape kinetic energy yet)
203      ekin=0.0_dp
204      do kk=1,ab_mover%natom
205        do jj=1,3
206 !        Warning : the fixing of atoms is implemented in reduced
207 !        coordinates, so that this expression is wrong
208          if (ab_mover%iatfix(jj,kk) == 0) then
209            ekin=ekin+0.5_dp*ab_mover%amass(kk)*vel(jj,kk)**2
210          end if
211        end do
212      end do
213      write(message, '(a,1p,e22.14,a)' )&
214 &     ' Kinetic energy of ions (ekin) [Ha]=',ekin
215      call wrtout(iout,message,'COLL')
216 
217    end if
218  end if
219 
220 !###########################################################
221 !### 3. ACELL
222 
223 !Only if the acell is being used
224  if (hist%isARused)then
225 !  Only if acell is recorded in a history
226    if (allocated(hist%acell))then
227 
228      write(message, '(a)' ) &
229 &     ' Scale of Primitive Cell (acell) [bohr]'
230      write(message,fmt1)&
231 &     TRIM(message),ch10,acell(:)
232      call wrtout(iout,message,'COLL')
233    end if
234  end if
235 
236 !###########################################################
237 !### 4. RPRIMD
238 
239 !Only if the acell is being used
240  if (hist%isARused)then
241 !  Only if rprimd is recorded in a history
242    if (allocated(hist%rprimd))then
243      write(message, '(a)' ) &
244 &     ' Real space primitive translations (rprimd) [bohr]'
245      do kk=1,3
246        write(message,fmt1)&
247 &       TRIM(message),ch10,&
248 &       rprimd(:,kk)
249      end do
250      call wrtout(iout,message,'COLL')
251    end if
252  end if
253 
254 !###########################################################
255 !### 5. Unitary cell volume
256 
257  if (ab_mover%optcell/=0)then
258 
259    ucvol=&
260 &   rprimd(1,1)*(rprimd(2,2)*rprimd(3,3)-rprimd(3,2)*rprimd(2,3))+&
261 &   rprimd(2,1)*(rprimd(3,2)*rprimd(1,3)-rprimd(1,2)*rprimd(3,3))+&
262 &   rprimd(3,1)*(rprimd(1,2)*rprimd(2,3)-rprimd(2,2)*rprimd(1,3))
263 
264    write(message, '(a,1p,e22.14)' )&
265 &   ' Unitary Cell Volume (ucvol) [Bohr^3]=',&
266 &   ucvol
267    call wrtout(iout,message,'COLL')
268 
269 !  ###########################################################
270 !  ### 5. Angles and lengths
271 
272 !  Compute real space metric.
273    rmet = MATMUL(TRANSPOSE(rprimd),rprimd)
274 
275    angle(1)=acos(rmet(2,3)/sqrt(rmet(2,2)*rmet(3,3)))/two_pi*360.0d0
276    angle(2)=acos(rmet(1,3)/sqrt(rmet(1,1)*rmet(3,3)))/two_pi*360.0d0
277    angle(3)=acos(rmet(1,2)/sqrt(rmet(1,1)*rmet(2,2)))/two_pi*360.0d0
278 
279    write(message, '(a,a)' ) &
280 &   ' Angles (23,13,12)= [degrees]'
281    write(message,fmt1)&
282 &   TRIM(message),ch10,&
283 &   angle(:)
284    call wrtout(iout,message,'COLL')
285 
286    write(message, '(a,a)' ) &
287 &   ' Lengths [Bohr]'
288    write(message,fmt1)&
289 &   TRIM(message),ch10,&
290 &   sqrt(rmet(1,1)),sqrt(rmet(2,2)),sqrt(rmet(3,3))
291    call wrtout(iout,message,'COLL')
292 
293 
294 !  ###########################################################
295 !  ### 5. Stress Tensor
296 
297    if(pos==mover_AFTER)then
298 !    Only if strten is recorded in a history
299      if (allocated(hist%strten))then
300 
301        write(message, '(a)' ) &
302 &       ' Stress tensor in cartesian coordinates (strten) [Ha/bohr^3]'
303 
304        write(message,fmt1)&
305 &       TRIM(message),ch10,&
306 &       strten(1),strten(6),strten(5)
307        write(message,fmt1)&
308 &       TRIM(message),ch10,&
309 &       strten(6),strten(2),strten(4)
310        write(message,fmt1)&
311 &       TRIM(message),ch10,&
312 &       strten(5),strten(4),strten(3)
313        call wrtout(iout,message,'COLL')
314      end if
315    end if
316  end if
317 
318 !###########################################################
319 !### 6. Energy
320 
321  if(pos==mover_AFTER)then
322    write(message, '(a,1p,e22.14)' )&
323 &   ' Total energy (etotal) [Ha]=',&
324 &   hist%etot(hist%ihist)
325 
326    if (itime>1)then
327      jj = abihist_findIndex(hist,-1)
328      dEabs=hist%etot(hist%ihist)-hist%etot(jj)
329      dErel=2*dEabs/(abs(hist%etot(hist%ihist))+&
330 &     abs(hist%etot(jj)))
331      write(message, '(a,a,a,a)' )&
332 &     TRIM(message),ch10,ch10,&
333 &     ' Difference of energy with previous step (new-old):'
334      write(message, '(a,a,10a,a,1p,e12.5,a,10a,a,1p,e12.5)')&
335 &     TRIM(message),ch10,&
336 &     (' ',jj=1,10),' Absolute (Ha)=',dEabs,ch10,&
337 &     (' ',jj=1,10),' Relative     =',dErel
338    end if
339    call wrtout(iout,message,'COLL')
340  end if
341 
342  contains