Updated documentation for electric fields.
authorDavid van der Spoel <spoel@xray.bmc.uu.se>
Sun, 15 Feb 2015 20:04:48 +0000 (21:04 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 2 Mar 2015 15:47:24 +0000 (16:47 +0100)
Simulations with alternating electric fields have been possible for
8 years without this being documented. Now it is documented.
A please cite call has been added as well.

Change-Id: Iec60972622af3b2f9c7f1009e54e24bbd5600cd2

docs/user-guide/mdp-options.rst
src/gromacs/gmxlib/copyrite.cpp
src/gromacs/mdlib/sim_util.cpp

index 5129754f26ec07eddb69fdcee5771ba24baa3052..86fdfef6105a8ce74ef67ae216e13002d17a4dc8 100644 (file)
@@ -2544,7 +2544,25 @@ Electric fields
 
 .. mdp:: E-xt; E-yt; E-zt:
 
-   not implemented yet
+   Here you can specify a pulsed alternating electric field. The field
+   has the form of a gaussian laser pulse:
+
+   E(t) = E0 exp ( -(t-t0)^2/(2 sigma^2) ) cos(omega (t-t0))
+
+   For example, the four parameters for direction x are set in the
+   three fields of :mdp:`E-x` and :mdp:`E-xt` like
+
+   E-x  = 1 E0 0
+
+   E-xt = omega t0 sigma
+
+   In the special case that sigma = 0, the exponential term is omitted
+   and only the cosine term is used.
+
+   More details in Carl Caleman and David van der Spoel: Picosecond
+   Melting of Ice by an Infrared Laser Pulse - A Simulation Study
+   Angew. Chem. Intl. Ed. 47 pp. 14 17-1420 (2008)
+
 
 
 Mixed quantum/classical molecular dynamics
index ab622947ca6c38d168f1c639875e0a1680577900..25225b882e8b06b98d46cd1f8b439b470b9257aa 100644 (file)
@@ -384,6 +384,11 @@ void please_cite(FILE *fp, const char *key)
           "Thermodynamics of liquids: standard molar entropies and heat capacities of common solvents from 2PT molecular dynamics",
           "Phys. Chem. Chem. Phys.",
           13, 2011, "169-181" },
+        { "Caleman2008a",
+          "C. Caleman and D. van der Spoel",
+          "Picosecond Melting of Ice by an Infrared Laser Pulse: A Simulation Study",
+          "Angew. Chem. Int. Ed",
+          47, 2008, "1417-1420" },
         { "Caleman2011b",
           "C. Caleman and P. J. van Maaren and M. Hong and J. S. Hub and L. T. da Costa and D. van der Spoel",
           "Force Field Benchmark of Organic Liquids: Density, Enthalpy of Vaporization, Heat Capacities, Surface Tension, Isothermal Compressibility, Volumetric Expansion Coefficient, and Dielectric Constant",
index 9e2fc4de1b61d1feec81c05982a9392a0c22e269..2df9b2a73e9a76a033dfc92ae4fcf1ff756de455 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -211,11 +211,9 @@ static void sum_forces(int start, int end, rvec f[], rvec flr[])
  * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
  *
  * Et[] contains the parameters for the time dependent
- * part of the field (not yet used).
+ * part of the field.
  * Ex[] contains the parameters for
- * the spatial dependent part of the field. You can have cool periodic
- * fields in principle, but only a constant field is supported
- * now.
+ * the spatial dependent part of the field.
  * The function should return the energy due to the electric field
  * (if any) but for now returns 0.
  *
@@ -226,7 +224,7 @@ static void sum_forces(int start, int end, rvec f[], rvec flr[])
  * For neutral systems with many charged molecules the error is small.
  * But for systems with a net charge or a few charged molecules
  * the error can be significant when the field is high.
- * Solution: implement a self-consitent electric field into PME.
+ * Solution: implement a self-consistent electric field into PME.
  */
 static void calc_f_el(FILE *fp, int  start, int homenr,
                       real charge[], rvec f[],
@@ -2785,7 +2783,10 @@ void init_md(FILE *fplog,
             please_cite(fplog, "Goga2012");
         }
     }
-
+    if ((ir->et[XX].n > 0) || (ir->et[YY].n > 0) || (ir->et[ZZ].n > 0))
+    {
+        please_cite(fplog, "Caleman2008a");
+    }
     init_nrnb(nrnb);
 
     if (nfile != -1)