git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@16053 f3b2605a-c512-4ea7-a41b...
[lammps.git] / doc / src / fix_pimd.txt
blob38022e4c7dd3d6d2e2bb1b0a56177fae2918604c
1 "LAMMPS WWW Site"_lws - "LAMMPS Documentation"_ld - "LAMMPS Commands"_lc :c
3 :link(lws,http://lammps.sandia.gov)
4 :link(ld,Manual.html)
5 :link(lc,Section_commands.html#comm)
7 :line
9 fix pimd command :h3
11 [Syntax:]
13 fix ID group-ID pimd keyword value ... :pre
15 ID, group-ID are documented in "fix"_fix.html command :ulb,l
16 pimd = style name of this fix command :l
17 zero or more keyword/value pairs may be appended :l
18 keyword = {method} or {fmass} or {sp} or {temp} or {nhc} :l
19   {method} value = {pimd} or {nmpimd} or {cmd}
20   {fmass} value = scaling factor on mass
21   {sp} value = scaling factor on Planck constant
22   {temp} value = temperature (temperarate units)
23   {nhc} value = Nc = number of chains in Nose-Hoover thermostat :pre
24 :ule
26 [Examples:]
28 fix 1 all pimd method nmpimd fmass 1.0 sp 2.0 temp 300.0 nhc 4 :pre
30 [Description:]
32 This command performs quantum molecular dynamics simulations based on
33 the Feynman path integral to include effects of tunneling and
34 zero-point motion.  In this formalism, the isomorphism of a quantum
35 partition function for the original system to a classical partition
36 function for a ring-polymer system is exploited, to efficiently sample
37 configurations from the canonical ensemble "(Feynman)"_#Feynman.
38 The classical partition function and its components are given
39 by the following equations:
41 :c,image(Eqs/fix_pimd.jpg)
43 The interested user is referred to any of the numerous references on
44 this methodology, but briefly, each quantum particle in a path
45 integral simulation is represented by a ring-polymer of P quasi-beads,
46 labeled from 1 to P.  During the simulation, each quasi-bead interacts
47 with beads on the other ring-polymers with the same imaginary time
48 index (the second term in the effective potential above).  The
49 quasi-beads also interact with the two neighboring quasi-beads through
50 the spring potential in imaginary-time space (first term in effective
51 potential).  To sample the canonical ensemble, a Nose-Hoover massive
52 chain thermostat is applied "(Tuckerman)"_#pimd-Tuckerman.  With the
53 massive chain algorithm, a chain of NH thermostats is coupled to each
54 degree of freedom for each quasi-bead.  The keyword {temp} sets the
55 target temperature for the system and the keyword {nhc} sets the
56 number {Nc} of thermostats in each chain.  For example, for a
57 simulation of N particles with P beads in each ring-polymer, the total
58 number of NH thermostats would be 3 x N x P x Nc.
60 NOTE: This fix implements a complete velocity-verlet integrator
61 combined with NH massive chain thermostat, so no other time
62 integration fix should be used.
64 The {method} keyword determines what style of PIMD is performed.  A
65 value of {pimd} is standard PIMD.  A value of {nmpimd} is for
66 normal-mode PIMD.  A value of {cmd} is for centroid molecular dynamics
67 (CMD).  The difference between the styles is as follows.
69 In standard PIMD, the value used for a bead's fictitious mass is
70 arbitrary.  A common choice is to use Mi = m/P, which results in the
71 mass of the entire ring-polymer being equal to the real quantum
72 particle.  But it can be difficult to efficiently integrate the
73 equations of motion for the stiff harmonic interactions in the ring
74 polymers.
76 A useful way to resolve this issue is to integrate the equations of
77 motion in a normal mode representation, using Normal Mode
78 Path-Integral Molecular Dynamics (NMPIMD) "(Cao1)"_#Cao1.  In NMPIMD,
79 the NH chains are attached to each normal mode of the ring-polymer and
80 the fictitious mass of each mode is chosen as Mk = the eigenvalue of
81 the Kth normal mode for k > 0. The k = 0 mode, referred to as the
82 zero-frequency mode or centroid, corresponds to overall translation of
83 the ring-polymer and is assigned the mass of the real particle.
85 Motion of the centroid can be effectively uncoupled from the other
86 normal modes by scaling the fictitious masses to achieve a partial
87 adiabatic separation.  This is called a Centroid Molecular Dynamics
88 (CMD) approximation "(Cao2)"_#Cao2.  The time-evolution (and resulting
89 dynamics) of the quantum particles can be used to obtain centroid time
90 correlation functions, which can be further used to obtain the true
91 quantum correlation function for the original system.  The CMD method
92 also uses normal modes to evolve the system, except only the k > 0
93 modes are thermostatted, not the centroid degrees of freedom.
95 The keyword {fmass} sets a further scaling factor for the fictitious
96 masses of beads, which can be used for the Partial Adiabatic CMD
97 "(Hone)"_#Hone, or to be set as P, which results in the fictitious
98 masses to be equal to the real particle masses.
100 The keyword {sp} is a scaling factor on Planck's constant, which can
101 be useful for debugging or other purposes.  The default value of 1.0
102 is appropriate for most situations.
104 The PIMD algorithm in LAMMPS is implemented as a hyper-parallel scheme
105 as described in "(Calhoun)"_#Calhoun.  In LAMMPS this is done by using
106 "multi-replica feature"_Section_howto.html#howto_5 in LAMMPS, where
107 each quasi-particle system is stored and simulated on a separate
108 partition of processors.  The following diagram illustrates this
109 approach.  The original system with 2 ring polymers is shown in red.
110 Since each ring has 4 quasi-beads (imaginary time slices), there are 4
111 replicas of the system, each running on one of the 4 partitions of
112 processors.  Each replica (shown in green) owns one quasi-bead in each
113 ring.
115 :c,image(JPG/pimd.jpg)
117 To run a PIMD simulation with M quasi-beads in each ring polymer using
118 N MPI tasks for each partition's domain-decomposition, you would use P
119 = MxN processors (cores) and run the simulation as follows:
121 mpirun -np P lmp_mpi -partition MxN -in script :pre
123 Note that in the LAMMPS input script for a multi-partition simulation,
124 it is often very useful to define a "uloop-style
125 variable"_variable.html such as
127 variable ibead uloop M pad :pre
129 where M is the number of quasi-beads (partitions) used in the
130 calculation.  The uloop variable can then be used to manage I/O
131 related tasks for each of the partitions, e.g.
133 dump dcd all dcd 10 system_$\{ibead\}.dcd
134 restart 1000 system_$\{ibead\}.restart1 system_$\{ibead\}.restart2
135 read_restart system_$\{ibead\}.restart2 :pre
137 [Restrictions:]
139 This fix is part of the USER-MISC package.  It is only enabled if
140 LAMMPS was built with that package.  See the "Making
141 LAMMPS"_Section_start.html#start_3 section for more info.
143 A PIMD simulation can be initialized with a single data file read via
144 the "read_data"_read_data.html command.  However, this means all
145 quasi-beads in a ring polymer will have identical positions and
146 velocities, resulting in identical trajectories for all quasi-beads.
147 To avoid this, users can simply initialize velocities with different
148 random number seeds assigned to each partition, as defined by the
149 uloop variable, e.g.
151 velocity all create 300.0 1234$\{ibead\} rot yes dist gaussian :pre
153 [Default:]
155 The keyword defaults are method = pimd, fmass = 1.0, sp = 1.0, temp = 300.0,
156 and nhc = 2.
158 :line
160 :link(Feynman)
161 [(Feynman)] R. Feynman and A. Hibbs, Chapter 7, Quantum Mechanics and
162 Path Integrals, McGraw-Hill, New York (1965).
164 :link(pimd-Tuckerman)
165 [(Tuckerman)] M. Tuckerman and B. Berne, J Chem Phys, 99, 2796 (1993).
167 :link(Cao1)
168 [(Cao1)] J. Cao and B. Berne, J Chem Phys, 99, 2902 (1993).
170 :link(Cao2)
171 [(Cao2)] J. Cao and G. Voth, J Chem Phys, 100, 5093 (1994).
173 :link(Hone)
174 [(Hone)] T. Hone, P. Rossky, G. Voth, J Chem Phys, 124,
175 154103 (2006).
177 :link(Calhoun)
178 [(Calhoun)] A. Calhoun, M. Pavese, G. Voth, Chem Phys Letters, 262,
179 415 (1996).