TABLE OF CONTENTS


ABINIT/pimd_nosehoover_npt [ Functions ]

[ Top ] [ Functions ]

NAME

 pimd_nosehoover_npt

FUNCTION

 Predicts new positions in Path Integral Molecular Dynamics using Nose-Hoover in the NPT ensemble.
 Given the positions at time t and t-dtion, an estimation of the velocities at time t,
 the forces and an estimation of the stress at time t, and an estimation of the cell at time t,
 computes in the Path Integral Molecular Dynamics framework the new positions at time t+dtion,
 computes self-consistently the velocities, the stress and the cell at time t and produces
 an estimation of the velocities, stress and new cell at time t+dtion
 No change of acell and rprim at present.

COPYRIGHT

 Copyright (C) 2011-2018 ABINIT group (GG,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

  etotal(trotter)=electronic total energy for all images
  itimimage=number of the current time for image propagation (itimimage+1 is to be predicted here)
  natom=dimension of vel_timimage and xred_timimage
  pimd_param=datastructure that contains all the parameters necessary to Path-Integral MD
  prtvolimg=printing volume
  rprimd(3,3)=dimensionless unit cell vectors (common to all images) at time t (present time step)
  rprimd_prev(3,3)=dimensionless unit cell vectors (common to all images) at time t-dt (previous time step)
  stressin(3,3,trotter)=electronic stress tensor for each image
  trotter=Trotter number (total number of images)
  volume=voume of unit cell (common to all images)
  xred(3,natom,trotter)=reduced coordinates of atoms for all images at time t (present time step)
  xred_prev(3,natom,trotter)=reduced coordinates of atoms for all images at time t-dt (previous time step)

OUTPUT

  rprimd_next(3,3)=dimensionless unit cell vectors (common to all images) at time t+dt (next time step)
  xred_next(3,natom,trotter)=reduced coordinates of atoms for all images at time t+dt (next time step)

SIDE EFFECTS

  forces(3,natom,trotter)=forces over atoms for all images
    at input,  electronic forces
    at output, electronic forces + Langevin contribution
  vel(3,natom,trotter)=velocies of atoms for all images
    at input,  values at time t
    at output, values at time t+dt
  vel_cell(3,3,trotter)=time derivative of cell parameters
    at input,  values at time t
    at output, values at time t+dt

PARENTS

      predict_pimd

CHILDREN

      pimd_energies,pimd_initvel,pimd_print,pimd_stresses,xcart2xred
      xred2xcart

SOURCE

 60 #if defined HAVE_CONFIG_H
 61 #include "config.h"
 62 #endif
 63 
 64 #include "abi_common.h"
 65 
 66 subroutine pimd_nosehoover_npt(etotal,forces,itimimage,natom,pimd_param,prtvolimg,&
 67 &          rprimd,rprimd_next,rprimd_prev,stressin,trotter,vel,vel_cell,&
 68 &          volume,xred,xred_next,xred_prev)
 69 
 70  use defs_basis
 71  use m_pimd
 72  use m_profiling_abi
 73 
 74  use m_geometry,  only : xcart2xred, xred2xcart
 75 
 76 !This section has been created automatically by the script Abilint (TD).
 77 !Do not modify the following lines by hand.
 78 #undef ABI_FUNC
 79 #define ABI_FUNC 'pimd_nosehoover_npt'
 80 !End of the abilint section
 81 
 82  implicit none
 83 
 84 !Arguments ------------------------------------
 85 !scalars
 86  integer,intent(in) :: itimimage,natom,prtvolimg,trotter
 87  real(dp),intent(in) :: volume
 88  type(pimd_type),intent(in) :: pimd_param
 89 !arrays
 90  real(dp),intent(in) :: etotal(trotter),rprimd(3,3),rprimd_prev(3,3),stressin(3,3,trotter)
 91  real(dp),intent(in),target :: xred(3,natom,trotter),xred_prev(3,natom,trotter)
 92  real(dp),intent(out) :: rprimd_next(3,3),xred_next(3,natom,trotter)
 93  real(dp),intent(inout) :: forces(3,natom,trotter),vel(3,natom,trotter),vel_cell(3,3,trotter)
 94 
 95 !Local variables-------------------------------
 96 !Options
 97  real(dp),parameter :: tolerance=tol7 ! SCF tolerance
 98 !scalars
 99  integer :: idum=-5
100  integer :: constraint,iimage,irestart,ndof,nnos,pitransform,prtstress
101  real(dp) :: dtion,eharm,eharm2,epot,initemp,kt,temperature1,temperature2,thermtemp
102 !arrays
103  real(dp) :: constraint_output(2),ddh(3,3),stress_pimd(3,3,3)
104  real(dp),allocatable :: dzeta(:,:,:,:),forces_orig(:,:,:),forces_pimd(:,:,:)
105  real(dp),allocatable :: inertmass(:),masseff(:,:),qmass(:),quantummass(:),springeff(:,:)
106  real(dp),allocatable :: xcart(:,:,:),xcart_next(:,:,:),xcart_prev(:,:,:),zeta(:)
107 
108 ! *************************************************************************
109 
110 !############# Initializations ###########################
111 
112 !Allocation of local arrays
113  ABI_ALLOCATE(xcart,(3,natom,trotter))
114  ABI_ALLOCATE(xcart_prev,(3,natom,trotter))
115  ABI_ALLOCATE(xcart_next,(3,natom,trotter))
116  ABI_ALLOCATE(forces_orig,(3,natom,trotter))
117  ABI_ALLOCATE(forces_pimd,(3,natom,trotter))
118  ABI_ALLOCATE(inertmass,(natom))
119  ABI_ALLOCATE(quantummass,(natom))
120 
121 !Fill in the local variables
122  ndof=3*natom*trotter
123  quantummass(1:natom)=pimd_param%amu   (pimd_param%typat(1:natom))*amu_emass
124  inertmass  (1:natom)=pimd_param%pimass(pimd_param%typat(1:natom))*amu_emass
125  initemp=pimd_param%mdtemp(1);thermtemp=pimd_param%mdtemp(2)
126  dtion=pimd_param%dtion;pitransform=pimd_param%pitransform
127  kt=thermtemp*kb_HaK
128  forces_orig=forces
129 
130 !Allocation/initialization of local variables used for Nose-Hoover chains
131 !Associated variables:
132 !nnos = number of thermostats
133 !dzeta(3,natom,trotter,nnos) = variables of thermostats, in (atomic time unit)^(-1)
134 !qmass(nnos) = masses of thermostats
135 !specific to PIMD: pitransform = coordinate transformation (0:no; 1:normal mode; 2:staging)
136  nnos=pimd_param%nnos
137  ABI_ALLOCATE(qmass,(nnos))
138  ABI_ALLOCATE(zeta,(nnos))
139  ABI_ALLOCATE(dzeta,(3,natom,trotter,nnos))
140  qmass(1:nnos)=pimd_param%qmass(1:nnos)
141  zeta=zero;dzeta=zero
142 
143 !Compute cartesian coordinates
144  do iimage=1,trotter
145    call xred2xcart(natom,rprimd,xcart     (:,:,iimage),xred(:,:,iimage))
146    call xred2xcart(natom,rprimd,xcart_prev(:,:,iimage),xred_prev(:,:,iimage))
147  end do
148 
149 !Determine if it is a restart or not
150 !If this is a calculation from scratch,generate random distribution of velocities
151  irestart=1;if (itimimage==1) irestart=pimd_is_restart(masseff,vel,vel_cell)
152 
153 !Initialize derivatives
154  if (mod(irestart,10)==0) then
155    call pimd_initvel(idum,masseff,natom,initemp,trotter,vel,pimd_param%constraint,pimd_param%wtatcon)
156  end if
157 !vel_cell does not depend on Trotter...
158  ddh=vel_cell(:,:,1);if (irestart<10) ddh=zero
159 
160 !Compute temperature at t
161  temperature1=pimd_temperature(masseff,vel)
162 
163 !################## Images evolution #####################
164 
165 !This is temporary
166  xcart_next=zero
167  rprimd_next=rprimd_prev
168  temperature2=pimd_temperature(masseff,vel)
169 
170 !Compute contributions to energy
171  call pimd_energies(eharm,eharm2,epot,etotal,forces_orig,natom,springeff,trotter,xcart)
172 
173 !Compute stress tensor at t from virial theorem
174  call pimd_stresses(masseff,natom,quantummass,stress_pimd,stressin,thermtemp,thermtemp,trotter,vel,volume,xcart)
175 
176 
177 !############# Final operations ############################
178 
179 !Print messages
180  prtstress=1
181  call pimd_print(constraint,constraint_output,&
182 & eharm,eharm2,epot,forces_pimd,inertmass,irestart,&
183 & itimimage,kt,natom,pimd_param%optcell,prtstress,prtvolimg,rprimd,&
184 & stress_pimd,temperature1,temperature2,&
185 & pimd_param%traj_unit,trotter,vel,ddh,xcart,xred)
186 
187 !Come back to reduced coordinates
188  do iimage=1,trotter
189    call xcart2xred(natom,rprimd,xcart_next(:,:,iimage),xred_next(:,:,iimage))
190  end do
191 
192 !Return cell velocities (does not depend on Trotter)
193  do iimage=1,trotter
194    vel_cell(:,:,iimage)=ddh(:,:)
195  end do
196 
197 !Free memory
198  ABI_DEALLOCATE(xcart)
199  ABI_DEALLOCATE(xcart_prev)
200  ABI_DEALLOCATE(xcart_next)
201  ABI_DEALLOCATE(forces_orig)
202  ABI_DEALLOCATE(forces_pimd)
203  ABI_DEALLOCATE(inertmass)
204  ABI_DEALLOCATE(quantummass)
205  ABI_DEALLOCATE(masseff)
206  ABI_DEALLOCATE(springeff)
207  ABI_DEALLOCATE(qmass)
208  ABI_DEALLOCATE(dzeta)
209  ABI_DEALLOCATE(zeta)
210 
211 end subroutine pimd_nosehoover_npt