CPMD VMD
CPMD VMD
CPMD VMD
Theoretische Chemie, Ruhr-Universitt Bochum,Germany http://www.theochem.ruhr-uni-bochum.de/go/cpmd-vmd.html Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
1. Introduction
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD using Density-Functional Perturbation Theory, BOMD, and for determination of energy barriers and reaction paths using Nudged Elastic Band (NEB) and Fourier String Method Dynamics (SMD). Visual Molecular Dynamics (VMD) is a very powerful and flexible toolkit for visualization and analysis of molecular dynamics simulation data. Although it was initially developed with focus on large biomolecular systems, newer versions contain many features that make it a very attractive choice for visualizing results from CPMD and other electronic structure calculations. Further information can be found on the VMD homepage at http://www.ks.uiuc.edu/Research/vmd/. There also is a VMD Script Library with more examples of VMD scripting.
1.3. Notes
Most of the examples presented here assume, were tested with and may only run with VMD version 1.8.3 or later on a UNIX(TM)/Linux machine. Contact axel.kohlmeyer@theochem.ruhr-uni-bochum.de if there are compatibility problems with newer versions of VMD or on other platforms. Any other form of feedback, e.g. corrections, enhancements, further examples, or new visualization problems, is also highly welcome.
2. Table of Contents
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
1. Introduction 1.1. About this Tutorial 1.2. About the Programs 1.3. Notes 1.4. Recent Changes 2. Table of Contents 3. Preparation and Installation Issues 3.1. Customizing the VMD Setup 3.2. Predefining Additional Items 3.3. Extending the Script and Plugin Search Path 4. Loading and Displaying Configurations and Trajectories 4.1. Loading a Single Geometry 4.2. Creating a Visualization 4.3. Saving and Re-using a Visualization 4.4. Viewing a Trajectory 4.5. How to Show Breaking and Formation of Bonds 4.6. Adding Graphics to a Visualization 5. Adding Dynamic Graphics to a Trajectory 5.1. Adding a Progress Bar for the Elapsed Time 5.2. Display the Total Dipole Moment of the Simulation Cell 5.3. Visualizing Changing Atom Properties with Color 5.4. Modify an Atom Property Dynamically from an External File 6. Dynamic Atom Selection 6.1. Display a Changing Number of Molecules 6.2. Keeping Atoms or a Molecule in the Center and Aligned 6.3. Modify a Selection During a Trajectory 6.4. Using the User Field for Computed Selections 6.5. Tracing a Dynamic Property 7. Visualizing Volumetric Data from Cube-Files 7.1. Electron Density and Electrostatic Potential 7.2. Canonical and Localized Orbitals 7.3. Electron Localization Function (ELF) 7.4. Manipulation of Cube Files / Response to an External Potential 7.5. Bulk Systems 7.6. Animations with Isosurfaces 7.7. Volumetric data from Gaussian 8. Using Data Processing to Tailor Data for VMD 8.1. Visualizing Path-Integral Trajectories 8.2. Extracting the Geometry Information from old CPMD Output Files 8.3. Removing Unneeded Parts From a Cube File 8.4. Extract Some Coordinates with Bounding Box Information 8.5. Creating 3d-Ramachandran Histograms 9. Misc Tips and Tricks 9.1. Collected 'draw' Extensions 9.2. Using a Different Default Visualization 9.3. Changing the Default vdW Radii 9.4. Reloading the Current Trajectory
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD 9.5. Visualize a Trajectory With a Changing Number of Atoms or Bonds 9.6. Set the Unit Cell Information for a Whole Trajectory 9.7. Directly Print the Current Visualization 9.8. Transferring a Visualization from a Molecule to Others 9.9. Adding a TCL-Plugin to the Extensions Menu 9.10. Turn Off Output in Analysis Scripts 9.11. Automatically Turn on TCL mode in (X)Emacs for .vmd Files 10. Credits 11. Script distribution policy
Many of the examples presented here utilize the excellent scripting capabilities of VMD and consist in part of reusable subroutines written in the VMD/Tcl dialect. To adapt the VMD scripting examples presented here, you should therefore first have a look at the excellent VMD User's Guide and also get a little bit acquainted with programming in Tcl/Tk. There is a Python scripting interface in VMD as well, but it is not (yet) covered here. Via the autoloading feature of the Tcl/Tk script engine embedded in VMD users can have a repository of Tcl subroutines which will be loaded on demand. This way you don't have to copy the files into every single VMD script and you can update your scripts by just replacing the files in the repository (of course only as long as the calling sequence stays intact).
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
There are several ways of how you may want to look at configurations generated by Single Point and MD calculations. The following sections will give you some examples on how to visualize single geometries, whole trajectories, Wannier centers, and (optionally) copies of the simulation unit cell, which can be very helpful for systems with periodic boundary conditions.
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD The easiest way to load molecules from CPMD calculations into VMD would be to load the GEOMETRY.xyz file that recent CPMD versions (3.5.x and later) generate by default. Beginning with VMD version 1.8.2 CPMD style xyz files are supported. If you have an older version of VMD, you have to first convert the file to a supported format, e.g. PDB, with Babel or Open Babel (VMD has internal Babel support, if installed correctly). For PWscf calculations you can use the scripts pwi2xsf.sh and/or pwo2xsf.sh to convert the respective input and output files to the xsf format, which is supported by VMD since version 1.8.3. Alternatively, you can use the simple pwscf2xyz.awk awk script to convert (some) PWscf output files to xyz format (this script may need some adjustments to support older/newer versions of PWscf). Please note, that the xsf format also usually contains the information of the supercell, which is not available for xyz files. For the configuration on the left and the following examples, use the file al-h2o_6.xyz. You can load this file into VMD by just typing vmd al-h2o_6.xyz and you should see something very similar to the picture on the left. In the next sections we are trying to improve it a little bit.
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD used for a different way to show the octahedron.
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD individual groups of atoms. If your trajectory includes Wannier centers, the distance based bond selection will produce strange results. You could use the selections "not name X" for the atom representations and "name X" for the Wannier centers (looks pretty cool, if you use the Glass material and a VDW representation for the 'Wannier-atoms'. See the next chapter for some examples).
# put a small yellow ball in the center of the box draw color yellow draw sphere {0 0 0} radius 0.25 resolution 20 # add simulation cell object draw_cubic_unitcell 9.956
Sometimes you want to draw additional graphics the screen that show how some properties change during a trajectory. VMD offers you access to graphic primitves like lines, spheres, cylinders, cones, trangles and text with which you can add graphical elements to your molecule. To make these change with the evolution of the trajectory you can use the powerful Tcl variable tracing mechanism. In this case you want to trace the variable vmd_frame(). The following sections will give a few examples. The strategy is always the same: 1. write a small function that changes the graphics depending on the current frame number and an additional property 2. read a list of values for the property into an array indexed by the corresponding frame number. 3. trace vmd_frame() with that function
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD Although VMD shows the current frame number in the main window, it is sometimes nicer to have an indicator of the elapsed time in the visualization window. One big problem is, that the progress indicator should stay fixed, even if the molecule itself is rotated. As described above, usually you 'draw' in the coordinate system of the molecule and so the added graphics would rotate with the visualize molecule. This problem is avoided by creating an (empty) stub molecule (by calling mol new without a file name), fixing the position of that 'molecule' in the global coordinate system, and adding the time bar graphics to the fixed 'molecule'. This is done by the script function display_time_bar. After you have loaded the trajectory from the files zundel.xyz and zundel.dcd and created a nice representation you run the following code: # function to draw the time bar with the required options proc do_time {args} { # see source for explanations of the options display_time_bar 1.0 50.0 "fs" 0 } # connect do_time to vmd_frame trace variable vmd_frame(0) w do_time # go back to the start of the trajectory animate goto start Now you can animate the trajectory as you like and the red line in the time bar will always the display the current time of the simulation. You can also move the mouse and the time display should stay fixed. The file zundel-time.vmd contains all the commands. You need to have the display_time_bar.tcl file installed properly to make this work. Alternatively you can copy the file to the current directory and add the command source display_time_bar.tcl to the .vmd script (or execute it manually) before you define do_time.
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD if {[info exists dipgraph]} then { foreach g $dipgraph { graphics $molid delete $g } } graphics $molid color yellow lassign $dipdata($frame) cnt vec # scaling res set dipgraph [vmd_draw_vector $molid $cnt $vec 100000.0 10 } } # read in the dipole data from file set dip [open "zundel.dip" r] # define selection to calculate the center of mass set sel [atomselect 0 "not name X"] set n [molinfo 0 get numframes] for {set i 0} {$i < $n} {incr i} { # advance selection to the current frame and update. $sel frame $n $sel update set dipdata($i) [list [center_of_mass $sel] [gets $dip]] } close $dip # connect to vmd_frame trace variable vmd_frame(0) w do_dipdraw animate goto 0 Again the whole vmd script file is also available under the name zundel-dip.vmd. Of course, you can combine the time bar with the dipole vector display. A vmd script file to do that is available under the name zundel-all.vmd NOTE: for this simple example to work, the trajectory and the dipole file have to be written with the same frequency.
radius 0.08]
10
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD code: set mol 0 set left [atomselect $mol {index < 60}] set right [atomselect $mol {index>= 60}] set n [ molinfo $mol get numframes ] for {set i 0} {$i < $n} {incr i} { # compute for first molecule $left frame $i set com [measure center $left] set dlist "" # loop over each {x y z} vector in $c foreach c [$left get {x y z}] { lappend dlist [veclength [vecsub $c $com]] } # set distance from COM for all atoms in selection. $left set user $dlist # now repeat for second molecule $right frame $i set com [measure center $right] set dlist "" foreach c [$right get {x y z}] { lappend dlist [veclength [vecsub $c $com]] } $right set user $dlist } The complete script of this visualization script is also available for download under the name bucky.vmd.
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD set fp [open "ar3plus-charge.dat" r] for {set i 0} {$i < $n} {incr i} { set chrg($i) [gets $fp] } close $fp # procedure to change the charge field from the data in $chrg proc do_charge {args} { global chrg molid set a [molinfo $molid get numatoms] set f [molinfo $molid get frame] for {set i 0} {$i < $a} {incr i} { set s [atomselect $molid "index $i"] $s set charge [lindex $chrg($f) $i] } } # turn on update of the charge info and coloring in each frame. trace variable vmd_frame($molid) w do_charge mol colupdate 0 $molid on # set color mapping parameters mol scaleminmax 0 $molid 0.0 0.7 color scale method RGB color scale midpoint 0.25 color scale min 0.0 color scale max 0.7 # go back to the beginning and activate the additional feartures. animate goto start do_charge The complete script of this visualization is also available for download under the name ar3plus-charge.vmd.
12
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD shared between both ions. Now how was this done? First we need to load the trajectory from the files 800h2o-16nacl.xyz and 800h2o-16nacl.dcd and delete the default representation. Next we define a set of atomselect macros that make the selection of the atoms much easier. Here is the block for the sodium hydration shells (the cutoff distances 3.3 and 5.5 were taken from the gNaO(r) function for this simulation). # selection macros. this is the first part of the magic atomselect macro hyd {name H} atomselect macro oxy {name O} # first and second solvation shell of Na+ and oxygens in transition atomselect macro naox1 {oxy within 3.3 of index 2405} atomselect macro naox2 {(oxy within 5.5 of index 2405) and not naox1} atomselect macro naotr {(oxy within 3.35 of index 2405) and not (oxy within 3.25 of index 2405)} atomselect macro nahy1 {hyd within 1.31 of naox1} atomselect macro nahy2 {hyd within 1.31 of naox2} Note how the hydrogens are selected with a distance criterium around the oxygen selections. This way you make sure that you select only complete water molecules for each shell. The selection for the chloride follows the same pattern. The next step is to create the representations utilizing the selection macros. Finally you have to turn on recalculation of the selection for each representation (and we add a little trajectory smoothing to have a less jumpy animation). This is done with the following code: # let all selections be recalculated for each frame # and smooth the trajectory a little bit for all representations # that's part two of the magic. set molid 0 set n [molinfo $molid get numreps] for {set i 0} {$i < $n} {incr i} { mol selupdate $i $molid on mol smoothrep $molid $i 2 } # go back to the start of the trajectory. animate goto start That is it. How to display all hydrogen bonds for all visible water molecule is left to you as an exercise. The full VMD script file is also available for download under the name nacl-shell.vmd. There also is a somewhat larger (640x480 pixel, 4.7 MByte ) rendered MPEG-1 Movie of the same animation.
13
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD center and direction of the ion pair directly would give a very jumpy picture. So we average both parameters over several frames (here 20) and used the smoothed realignment instead. See the movie on the left as an example. The code in nacl-follow.vmd differs from the previous example primarily by the following, additional code: proc do_realign {args} { global chist dhist hcount hoffs molid # this is the axis to align to set avec [vecnorm {1.0 0.0 0.0}] # this is the sliding window size set asize 20 # initialize the cache counters if ([info exists hcount]) { } else { set hcount 0 set hoffs -1 } # calculate center and direction of the ion pair set sel [atomselect $molid "index 2405 or index 2431"] lassign [$sel get {x y z}] na cl set cent [vecscale [vecadd $na $cl] 0.5] set dir [vecsub $na $cl] # store data in cache for sliding window averaring if {$hcount < $asize} then { incr hcount } incr hoffs if {$hoffs>= $hcount} then { set hoffs 0 } set chist($hoffs) $cent set dhist($hoffs) $dir # calculate averages set csum [veczero] set dsum [veczero] for {set i 0} {$i < $hcount} {incr i} { set csum [vecadd $csum $chist($i)] set dsum [vecadd $dsum $dhist($i)] } set csum [vecscale [expr 1.0/[expr $hcount * 1.0]] $csum] set dsum [vecnorm $dsum] # get rotation axis. set rvec [vecnorm [veccross $dsum $avec]] # set origin and rotation molinfo $molid set { center_matrix rotate_matrix} \ [list [trans origin $csum] \ [trans axis $rvec [expr acos([vecdot $dsum $avec])] rad]]
14
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD # clean up selections $sel delete } trace variable vmd_frame($molid) w do_realign # zoom to fit the display when vmd is started with '-size 800 600' scale to 0.15 # the smoothing works better with this (no discontinuity) animate style rock
15
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD set highlight($i) [gets $fp] } close $fp # hook highlight function into animation trace variable vmd_frame($molid) w do_highlight animate goto start For method b) we put the code to compute the selection within the highlight function, but first check the array whether we have already computed the selection and then use the cached value instead. This example is available under the name h3oplus-cache.vmd: # subroutine to update the selection of # representation $selid and cache the results proc do_highlight {args} { global highlight molid selid # get the current frame number set frame [molinfo $molid get frame] if {[info exists highlight($frame)]} then { mol modselect $selid $molid "$highlight($frame)" } else { set hl {none} # loop over all oxygens to count hydrogens within bondlength set osel [atomselect $molid "name O"] foreach ox [$osel list] { set sel [atomselect $molid "name H and within 1.32 of index $ox"] $sel frame $frame $sel update if {[$sel num] != 2 } { set hl "$hl or index $ox" } $sel delete } $osel delete # apply new selection and cache it mol modselect $selid $molid "$hl" set highlight($frame) "$hl" } } set molid [molinfo top] set selid 2 trace variable vmd_frame($molid) w do_highlight animate goto start
16
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD be accessed in selection expressions. Using this method on the previous example one would first need to compute the number of bonded Hydrogens and store it in the user field. Using selections for that avoids the time consuming computations in Tcl and thus combines the advantages of both methods presented above. # prep set num [molinfo $mol get numframes] set ox [atomselect $mol {name O}] set all [atomselect $mol {all}] # create a selection for each oxygen atom to compute foreach i [$ox get index] { set sel($i) [atomselect $mol "exwithin 1.30 of index $i"] } # loop over all frames, and for each frame loop over # all oxygens and store the number of hydrogens in user for {set n 0} {$n < $num} {incr n} { set bc {} foreach i [$ox get index] { $sel($i) frame $n $sel($i) update $all frame $n $all set user 0 lappend bc [$sel($i) num] } $ox frame $n $ox set user $bc unset bc } # clean up selections foreach i [$ox get index] { $sel($i) delete } $ox delete $all delete unset ox all sel i n Now to recognize all Oxygen atoms with more than two bonded Hydrogens, we simply need to define a selection with user > 2 and make sure that the selection is evaluated for every frame during an animation (-> Trajectory tab in the Graphical Representations menu. mol representation VDW 0.60000 20.000000 mol color ColorID 4 mol selection {name O and user > 2} mol material Transparent mol addrep top mol selupdate 3 $mol on The complete example is available for download: 32h2o_h3op_user.vmd.
17
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
# create set hyd1 set hyd2 set hyd3 set hyd4 set hyd5
selection for tracing hydrogens and the h3o-plus [atomselect $mol {index 85}] [atomselect $mol {index 99}] [atomselect $mol {index 86}] [atomselect $mol {index 97}] [atomselect $mol {index 98}]
18
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD set h3op [atomselect $mol {name O and user > 2}] # draw trajectory paths trajectory_path $hyd1 blue 0 4 trajectory_path $hyd2 green 0 4 trajectory_path $hyd3 orange 0 4 trajectory_path $hyd4 purple 0 4 trajectory_path $hyd5 ochre 0 4 trajectory_path $h3op yellow 1 4 The image demonstrates nicely and without any animation, that the excess proton 'defect' is much more movable than the individual Hydrogen atoms due to the Grotthuss structural diffusion mechanism.
Apart from calculating trajectories various volumetric properties can be calculated with the CPMD program: e.g. electron densities, spin densities, electrostatic potentials, electron localization functions (ELF), localized and canonical (occupied) orbitals. This data can be visualized with VMD via the Gaussian cube file format. Beginning with version 1.8.2 VMD fully supports reading atom coordinates and volumetric data in the cube file format. Please note that some of the examples need a very recent CPMD version (3.9.1 or newer) to work properly, as i found and fixed a few cubefile related bugs while creating this part of the tutorial (the same goes for the cpmd2cube program, a patch relative to the latest released version from cpmd.org is available on request). VMD version 1.8.3 and CPMD version 3.9.2 as well as the cpmd2cube version relased with it contain additional improvements, mainly for handling non-orthogonal supercells. For some of the aforementioned properties cube file can be directly generated from cpmd, the rest is written in a native format which can be converted into a cube file with the help of the cpmd2cube program (available at the download area of http://www.cpmd.org). Since the CPMD manual is not very detailed about these tasks, the following section will give a few CPMD input file examples alongside the suggestions for visualizations.
19
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD CUTOFF 120.0 ANGSTROM CELL 6.0 1.0 1.0 0.0 0.0 0.0 &END &ATOMS *O_MT_PBE KLEINMAN-BYLANDER LMAX=P 1 2.904516 3.000000 2.926732 *H_MT_PBE KLEINMAN-BYLANDER LMAX=S 2 2.900437 3.000000 3.897528 3.841176 3.000000 2.671532 &END
By adding the keywords RHOOUT and ELECTROSTATIC POTENTIAL the electron density and the electrostatic potential will be written to files named DENSITY and ELPOT, respectively, at the end of the wavefunction optimization. Note the rather high cutoff of 120ryd which is not really required, but helps to get smooth surfaces. Since we want to re-use this restart file as base for further calculations, please rename it from RESTART.1 to RESTART. With RESTART WAVEFUNCTION COORDINATES (Note: no LATEST) all further calculation will always read in this 'high quality' restart and it will not be overwritten. Finally we need to convert the two volumetric files to cube files by: cpmd2cube.x -o h2o-dens -dens DENSITY cpmd2cube.x -o h2o-pot -dens ELPOT As a result you will get the files h2o-dens.cube and h2o-pot.cube. To visualize volumetric data in VMD you currently have two options: Isosurfaces and Volume Slices. The two images in the top left of this section give an example for both styles. The upper image shows in addition to a CPK model of the water molecule an isosurface of the electron density (for = 0.05), the lower image a volume slice through the electrostatic potential (for y = 0.5). The CPMD input (h2o-dens-pot.in), the pseudopotential files (O_MT_PBE, H_MT_PBE), the resulting cube files (h2o-dens.cube, h2o-pot.cube) and a combined VMD script file (h2o-dens-pot.vmd) are available for download, so that you can experiment with it. VMD version 1.8.3 introduces the Volume coloring method (for OpenGL implementations that support it), which can be used to colorcode the surface by the value of the electrostatic potential at the position of the surface, i.e. color-map the electrostatic potential to the surface. The value range of the colormap can be adjusted via the Color Scale Data Range fields in the Trajectory tab of the Isosurface representation. Note: that for a 'properly symmetric' behavior the electrostatic potential from CPMD calculations needs to be corrected, e.g. via the trimcube utility. Downloadable example files: h2o-dens.cube, h2o-pot-norm.cube, and h2o-pot-map.vmd.
20
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD To calculate the canonical (i.e. projected to an atomic basis set) and localized orbitals for our water molecule we need to do a properties calculation starting from the previously generated wavefunction. Both sets of orbitals can be calculated individually as well as simultaneously (as done here). Note that when calculating orbital cube files, especially for localized orbitals, the default cube files contain a lot of unneeded data, so that using the trimcube utility, or the -trim option to cpmd2cube is highly recommended to save disk space and reduce the memory requirements of VMD. For the localized orbitals the keywords LOCALIZE and WANNIER WFNOUT are required. This generates a series of files with the names WANNIER_*.* which have to be converted to cube format with cpmd2cube.x -o h2o-local WANNIER_1.1. Note that you must not use the -dens option here. For the projected orbitals the CUBEFILE ORBITALS keyword (plus the two additional lines specifying how many and which orbitals shall be written) is required. This will generate a series of cube files with the names PSI.*.cube, that can be read into VMD directly. Alltogether the first part of the CPMD input now contains: &CPMD PROPERTIES RESTART WAVEFUNCTION COORDINATES WANNIER WFNOUT ALL &END &PROP LOCALIZE CUBEFILE ORBITALS HALFMESH 4 1 2 3 4 &END
The picture on the left demonstrates the localized orbitals. Here all cube files have been read in on top of each other and visualized with an isosurface (blue for the lone pairs, green for the OH-bonds). The picture on the right shows the HOMO of the water molecule. Here the two different phases are visualized by creating two repesentations from the same data set and just using isovalues with opposite sign for each of them. The input and cube files used in this section are: h2o-orbs.in, h2o-local-1.cube, h2o-local-2.cube, h2o-local-3.cube, h2o-local-4.cube, h2o-orbs-local.vmd, h2o-homo.cube, and h2o-homo.vmd.
21
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD atomic shells and other elements of chemical bonding. For more details please consult the ELF homepage at http://www.cpfs.mpg.de/ELF/. To calculate ELF with CPMD we need to add the keyword ELF PARAMETER to the &CPMD section of our input file. The resulting file named ELF is a density file similar to DENSITY and has to be converted to cube format with cpmd2cube. The picture on the left shows the ELF of our water example with an isovalue of 0.85. When increasing this value you can see, that the two lone pairs of the water are only partially localized. The relevant downloadable files are: h2o-elf.in, h2o-elf.cube, and h2o-elf.vmd.
# cpmd-vmd/files> cubman Action [Add, Copy, Difference, Properties, SUbtract, SCale]? su First input? h2o-nopot-dens.cube Is it formatted [no,yes,old]? yes Opened special file h2o-nopot-dens.cube.
22
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD Second input? h2o-extpot-dens.cube Is it formatted [no,yes,old]? yes Opened special file h2o-extpot-dens.cube. Output file? h2o-extpot-resp.cube Should it be formatted [no,yes,old]? yes Opened special file h2o-extpot-resp.cube.
The image on the left shows the result. The arrow illustrates the direction of the electrostatic field represented by the external potential and the two colored isosurfaced show areas of reduced (green) and increased (yellow) electron density, the transparent isosurface represents the total electron density. The formation of an induced dipole moment in the water molecule is clearly visible. The figure was created with the VMD script h2o-extpot-resp.vmd.
(click here or on the image to view a larger version) Creating and visualizing volumetric data is not restricted to isolated molecules, but can also be used for bulk systems, for example a bulk water system. Occasionally you have to make sure, that the cube file is properly centered and that all atoms are inside the density/simulation box. The command lines to get the cube files used here were: cpmd2cube.x -o h2o-dens -inbox -center -dens DENSITY cpmd2cube.x -o h2o-elpot -inbox -center -dens ELPOT After visualizing the water molecules you can add one or more isosurface representations in order to display different isosurfaces from the same data set: blue (= negative electrostatic potential), red (= same value only positive), pink (=even more positive isovalue), and green (= electron density). Again the full saved state and the data file are available for download with the filenames h2o-cube.vmd, h2o-dens.cube, and h2o-elpot.cube.
23
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD (click here or on the image to view a larger version) With the help of a little bit of VMD scripting we can also do an animation with volumetric data, where the the isosurface is changed in each frame. In the present example, we have a simulation, in which a hydrogen molecule is shot at a double layer of gold atoms with very high velocity. We then need a cube file for each frame, that should be animated, therefor we need to instruct the CPMD program to write a (different) restart for each of these frames. This is done by using the STORE and RESTFILE keywords in the CPMD input file (for this example STORE was set to 20 and RESTFILE to 500 with a TIMESTEP of 2 a.u.). Now during the simulation (which should not exceed 10000 steps), a new restart is written every 20 steps. The next step is to convert all of these restarts into cube files. This can done by doing one step of wavefunction optimization and using the keyword RHOOUT to write a DENSITY file, which then can be converted into a cube file with cpmd2cube.x. The first part of the CPMD input thus is: &CPMD OPTIMIZE WAVEFUNCTION RESTART WAVEFUNCTION COORDINATES MAXSTEP 1 RHOOUT &END The conversion can be automated with: for s in `seq 1 25` do \ rm -f RESTART ln -s RESTART.$s RESTART && ./cpmd.x au-dens.in && ./cpmd2cube.x -rho -o au-dens-$s DENSITY done We now load all those cube files in the correct order into VMD and create a visualization including an isosurface. To switch the volumetric data set during the animation, we write a small tcl procedure, that updates the representation and hook it into the animation loop by tracing vmd_frame. To make this as transparent as possible, we record the molecule id and the (unique) representation name of the isosurface in question in two global variables. set updmol [mol new {au-dens-0.cube} type cube waitfor all] ... set updrep [mol repname top 3] ... proc update_iso {args} { global updmol global updrep set repid [mol repindex $updmol $updrep] if {$repid < 0} { return } set frame [molinfo $updmol get frame] lassign [molinfo $updmol get "{rep $repid}"] rep
24
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD mol representation [lreplace $rep 2 2 $frame] mol modrep $repid $updmol } trace variable vmd_frame(0) w update_iso There is one drawback: you have to know, which index the isosurface visualization has. But that is easily done (if you can count and remember, that the first representation id is 0), since you start from a script anyways. Any subsequent changes to the representations should be transparent to the script. You can download the full VMD script au-iso.vmd, and an archive with the cubefiles au-dens-cube.tar.gz (25MB).
25
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD electrons of the quartett state are rather delocalized.
In some cases the existing data is not (well) suited to be read in for a VMD visualization, so it needs to be augmented or converted using external programs. Here are some examples.
8.2. Extracting the Geometry Information from old CPMD Output Files
When running geometry optimizations, having a look at the 'trajectory' of the optmization is often desirable. Newer CPMD versions support the XYZ flag to the OPTIMIZE GEOMETRY keyword to create an xyz-file that can be easily visualized, but for outputs from older CPMD versions (or in cases where setting the XYZ keyword was overlook or forgotten), one may want to extract the coordinates from the output file. Doing this with a text editor can be quite tiresome. The perl script out2xyz.pl will try to do this for you. Of course this script also works (or at least it should) for other jobs that contain geometry information.
26
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD usage: trimcube [options] <input cube> [<output cube>] available options: -h[elp] -t[hresh] <value> -n[orm]
print this info set trim threshold to <value>, (default 0.005) normalize so that the integral over the cube is zero.
use '-' to read from standard input program writes to standard output if not output file is given
27
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD The ramaplot tool in VMD is very useful for tracking changes of the protein backbone angles from a trajectory. But to get an impression of the statistical distribution of some configurations, the 'show all timesteps' mode in ramaplot can be easily misleading (see the picture on the right). In this example we have a short trajectory of a 12-ALA alpha-helix which is in fact stable during the time of the simulation. A histogram of the data would put an equal weight to each data point and with a so-called 'rubbersheet' representation (see picture on the left) one could could identify the regions of statistical relevance much better and see, that the alpha-helix does indeed not fall apart significantly. To achieve this, we have to first create a histogram of the alpha-carbon angles and then visualize it. Since VMD does not support this internally, we create an empty dummy molecule and add the surface to it using VMDs graphics primitives. With the following script code, we can create the raw histogram.
set sel [atomselect top {protein and name CA}] set res 36 set w [expr ($res - 1.0)/360.0] set n [molinfo [$sel molid] get numframes] for {set i 0 } { $i < $n } { incr i } { $sel frame $i $sel update foreach a [$sel get set phi [lindex set psi [lindex incr data([expr } } {phi psi}] { $a 0] $a 1] int(($phi + 180.0) * $w)],[expr int(($psi + 180.0) * $w)])
We then have to normalize the histogram and create the surface by drawing four triangles between the corners of each square of data points and their mid-point. We also use the z-value of the mid-point to colorize the triangles according to a BGR color scale. so that the bottom is blue and the peaks will be green with red tips: color scale method BGR color scale max 0.9 color scale midpoint 0.3 for {set i 0} {$i < $res} {incr i} { for {set j 0} {$j < $res} {incr j} { set i2 [expr $i + 1] set j2 [expr $j + 1] set x1 [expr ($i - (0.5 * $res)) * $len] set x2 [expr ($i2 - (0.5 * $res)) * $len] set xm [expr 0.5 * ($x1 + $x2)]
28
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD set y1 [expr ($j - (0.5 * $res)) * $len] set y2 [expr ($j2 - (0.5 * $res)) * $len] set ym [expr 0.5 * ($y1 + $y2)] set zm [expr ($data($i,$j) + $data($i2,$j2) \ + $data($i2,$j) + $data($i,$j2)) / 4.0] graphics $mol color [expr 17 + int (200 * $zm)] graphics $mol triangle "$x1 $y1 $data($i,$j)" "$xm $ym $zm" "$x2 $y1 $data($i2,$j)" graphics $mol triangle "$x1 $y1 $data($i,$j)" " $x1 $y2 $data($i,$j2)" "$xm $ym $zm" graphics $mol triangle "$x2 $y2 $data($i2,$j2)" "$x2 $y1 $data($i2,$j)" "$xm $ym $zm" graphics $mol triangle "$x2 $y2 $data($i2,$j2)" "$xm $ym $zm" "$x1 $y2 $data($i,$j2)" } } \ \ \ \
Completed with a nice border and some labels the full code can be put into a separate subroutine, so that the code to create the picture on the left becomes: mol new {12-ala.pdb} type pdb waitfor all mol addfile {12-ala.dcd} type dcd waitfor all set sel [atomselect top {resid > 1 and resid < 12 and name CA}] mk3drama $sel
The VMD script code is available for download under mk3drama.tcl. This subroutines also turns off, but does not delete the originally loaded molecule(s), so you can in fact create multiple histograms from different selections and, e.g., view them in turn by clicking on the 'D' symbols in the main VMD window. The visualization and the data files are available under 12-ala-rama.vmd, 12-ala.pdb, and 12-ala.dcd.
This chapter is a collection of (hopefully) useful odds and ends that came up while working on these pages and with VMD in general.
29
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD The following is a collection of extensions to the draw command. All of them return a list of graphics ids for the graphics primitives they were build of. When stored in a variable, they can be selectively deleted like in the following example (a more elaborate version is availiable for download as test-draw-ext.vmd):
30
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD # draw one blue prism and a red vector on top of it set gidlist [draw prism {-0.7 -0.5 0.0} {0.7 -0.5 0.0} {0.0 1.0 0.0}] draw color red append gidlist " " [draw vector2 {0.0 0.0 0.0} {0.0 0.0 2.0} 1.0 20]
31
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD # draw a green sphere draw color green set sphere [draw sphere {0.0 0.0 -2.0} radius 2.0 resolution 30] # delete only the sphere.
32
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD draw delete $sphere The individual extensions are: delete [<gid> [<gid>] ...] (vmd_draw_delete.tcl) This replaces the regular draw delete command with a version that can handle multiple ids and lists as arguments. arrow {x1 y1 z1} {x2 y2 z2} [<scale> <resolution> {<radius>}] (vmd_draw_arrow.tcl) This is an updated version of the draw arrow example from the VMD User's Guide, which returns the graphics ids of the components. The optional third argument allows to scale the size of the vector. vector {x1 y1 z1} {x2 y2 z2} [<scale> <resolution> {<radius>}] (vmd_draw_vector.tcl) This is similar to draw arrow, but the two vectors specify the center of the vector and the direction. vector2 {x1 y1 z1} {x2 y2 z2} [<scale> <resolution> {<radius>}] (vmd_draw_vector.tcl) This is another version of an arrow, but now the two vectors specify the basepoint and the direction. prism {x1 y1 z1} {x2 y2 z2} {x3 y3 z3} [<thickness>] (vmd_draw_prism.tcl) This is basically like draw triangle but the optional parameter gives the "thickness" of the triangle. vecfield <list of lists> [<scale> <resolution> {<radius>}] (vmd_draw_vector.tcl) This is a wrapper around draw vector that processes a list of pairs of coordinate triples ({{{x1_1 y1_1 z1_1} {x2_1 y2_1 z2_1}} {{x1_2 y1_2 z1_2} {x2_2 y2_2 z2_2}} ...}) for creating a large number of vectors conveniently. arrowfield <list of lists> [<scale> <resolution> {<radius>}] (vmd_draw_arrow.tcl) This is is like draw vecfield only that the arguments follow the draw arrow conventions. unitcell <options> (vmd_draw_unitcell.tcl) This will draw add a unitcell graph to the top molecule. Available options: cell (vmd|auto|{a b c alpha beta gamma}), default:"vmd" "vmd" will use the internal values "auto" will build a orthogonal unitcell based on 'measure minmax {all}' else a list of a,b,c,alpha,beta,gamma will be assumed. origin ({x y z}|auto) default: {0.0 0.0 0.0} or "auto" with "cell auto" style (lines|dashed|rod) default: line width 'width' of the lines/rods, default: 1.0 resolution resolution of cylinders/spheres for 'style rod',default: 8
33
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD mol representation CPK 1.000000 0.300000 19.000000 16.000000 mol addrep $molid } # add 'default visualization' hack for the first three molecules. trace variable vmd_read_trajectory(0) w my_def_viz trace variable vmd_read_trajectory(1) w my_def_viz trace variable vmd_read_trajectory(2) w my_def_viz Note that this 'hack' does not work for molecules loaded from the commandline as they are loaded, before .vmdrc is processed. Also, in case of vmd_init_structure the default visualization (Lines) will still be added to the molecule, as it will be added after my_def_viz has been processed so that there is no chance to delete there. In the case of vmd_read_trajectory the new visualization will be created after the whole file has been read, which may not be what you want, if you want to see the new visualization alread during the loading of the trajectory. WARNING: Visualizations other then Line can put a severe strain on your available graphics and memory resources when loading very large molecules or structures.
hing radii 1.5 } {name N } { 1.4 } {name O } { 1.3 } {name F } { 1.2 } {name P } { 2.0 }
{na
$selstr"]
34
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
35
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD # the molmovie procedure creates an animation by # sequentially turning the molecules on and off # default is to wait a 500 miliseconds between 'frames' proc molmovie {{loops 10} {delay 500}} { global molmovie_last set nmols [molinfo num] if {![info exists molmovie_last]} { set molmovie_last [expr $nmols - 1] } for {set i 0} {$i < $loops} {incr i} { for {set n 0} {$n < $nmols} {incr n} { display update display update ui mol on $n mol off $molmovie_last set molmovie_last $n after $delay } } } # example use of molmovie command # load frames as separate molecules for {set i 0} {$i < 15} {incr i} { mol new [format "zundel-frame-%02d.pdb" $i] } # now disable displaying all molecules for {set i 0} {$i < [molinfo num]} {incr i} { mol off $i } # run animation for two loops molmovie 2 The upper example are just some frames of the zundel ion movie from the beginning of the tutorial (zundel-frames.tar.gz), the lower example is more elaborate and from a simulation, where the crossing of the periodic boundaries of long molecules requires new connectivity information in each step after wrapping all atom positions back into the unitcell (psf-molmovie-example.vmd). To get a uniform visualization for all frames, it utilizes the clone_reps command. The data for this example is also available: psf-molmovie-example.tar.gz (1 MByte).
36
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD following set_unitcell procedure (set_unitcell.tcl) will make this more convenient for simulations with fixed cell dimensions. proc set_unitcell {a b c {molid top} {alpha 90.0} {beta 90.0} {gamma 90.0}} { if {![string compare $molid top]} { set molid [molinfo top] } set n [molinfo $molid get numframes] for {set i 0} {$i < $n} {incr i} { molinfo $molid set frame $i molinfo $molid set {a b c alpha beta gamma} \ [list $a $b $c $alpha $beta $gamma] } }
37
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
} Now we need to provide a procedure ::CloneRep::clonerepgui that will create the GUI and a few callback functions for the active elements of the GUI (e.g. the buttons). Since this should not become a general Tk/TCL tutorial, we will concentrate on one small example (see the comments in the script file clonerepgui.tcl for some more info): the callback for the Clone button. We import the two global variables with the molecule ids and then call the clone_reps procedure explained above. proc ::CloneRep::CloneRepDoClone { args } { variable fromid variable toid clone_reps $fromid $toid }
38
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD Since we want to register the GUI in the Extensions menu, we also have to provide a callback function, that will create the GUI and return a handle to the main widget. proc clonerepgui_tk_cb {} { ::CloneRep::clonerepgui return $::CloneRep::w } To finally register the plugin with the main VMD GUI, we have to put the following code into our .vmdrc file. if { [catch "package require clonerepgui" msg] } { puts "CloneRep plugin could not be found:\n$msg" } elseif { [catch {menu tk register clonerepgui clonerepgui_tk_cb} msg] } { puts "CloneRep plugin could not be registered:\n$msg" } else { puts "CloneRep activated" } As usual the files clone_reps.tcl and clonerepgui.tcl are available for download.
39
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD proc source_silent {args} { foreach f $args { source $f } } # execute the original script source_silent {get-his.vmd}
. tcl-mode)) auto-mode-alist))
10. Credits The following is a list of tools and persons (in no specific order) that played an important role in making this page reality. This page was compiled by Axel Kohlmeyer, Lehrstuhl fr Theoretische Chemie, Ruhr-Universitt Bochum, axel.kohlmeyer@theochem.ruhr-uni-bochum.de Additional simulation data was contributed by Holger Langer, Volker Kleinschmidt, Marcel Baer Lehrstuhl fr Theoretische Chemie, Ruhr-Universitt Bochum Heshe Peshkin Helpful suggesions and encouragement John E. Stone, Dominik Marx, Roger Rousseau, Volker Kleinschmidt, Holger Langer, Eduard Schreiner, Amalendu Chandra, Marcel Baer VMD Humphrey, W., Dalke, A. and Schulten, K., VMD - Visual Molecular Dynamics" J. Molec. Graphics 1996, 14.1, 33-38. http://www.ks.uiuc.edu/Research/vmd/ CPMD Copyright IBM Corp 1990-2005, Copyright MPI fr Festkrperforschung Stuttgart 1997-2001. http://www.cpmd.org/ ESPRESSO / PWscf
40
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD Stefano Baroni, Stefano de Gironcoli, Andrea Dal Corso (SISSA, Trieste), Paolo Giannozzi (Scuola Normale, Pisa) and others. http://www.pwscf.org/ Gaussian 98 (Revision A.11) M. J. Frisch, G. W. Trucks, H. B. Schlegel, G. E. Scuseria, M. A. Robb, J. R. Cheeseman, V. G. Zakrzewski, J. A. Montgomery, R. E. Stratmann, J. C. Burant, S. Dapprich, J. M. Millam, A. D. Daniels, K. N. Kudin, M. C. Strain, O. Farkas, J. Tomasi, V. Barone, M. Cossi, R. Cammi, B. Mennucci, C. Pomelli, C. Adamo, S. Clifford, J. Ochterski, G. A. Petersson, P. Y. Ayala, Q. Cui, K. Morokuma, D. K. Malick, A. D. Rabuck, K. Raghavachari, J. B. Foresman, J. Cioslowski, J. V. Ortiz, B. B. Stefanov, G. Liu, A. Liashenko, P. Piskorz, I. Komaromi, R. Gomperts, R. L. Martin, D. J. Fox, T. Keith, M. A. Al-Laham, C. Y. Peng, A. Nanayakkara, C. Gonzalez, M. Challacombe, P. M. W. Gill, B. G. Johnson, W. Chen, M. W. Wong, J. L. Andres, M. Head-Gordon, E. S. Replogle and J. A. Pople, Gaussian, Inc., Pittsburgh PA, 1998. Copyright 1988, 1990, 1992-1999, Gaussian, Inc., and copyright 1983, Carnegie Mellon University http://www.gaussian.com/ Tachyon Copyright 1994-2005 John E. Stone http://jedi.ks.uiuc.edu/~johns/raytracer/ Raster3D Merritt & Bacon (1997) Meth. Enzymol. 277, 505-524. http://www.bmsc.washington.edu/raster3d/ ImageMagick Copyright 1999-2004 ImageMagick Studio LLC. http://www.imagemagick.org/ Gifsicle Copyright 1997-2003 by Eddie Kohler. http://www.lcdf.org/gifsicle/ htmldoc Copyright 1997-2004 by Easy Software Products. http://www.easysw.com/htmldoc/ Mjpegtools by Rainer Johanni, Gernot Ziegler, Andrew Stevens, Bernhard Praschinger, Ronald Bultje, Xavier Biquard, Matthew Marjanovic, pHilipp Zabel, Kawamata/Hitoshi, Stefan Fendt, Scott Moser, Shawn Sulma, Mike Bernson, James Klicman http://mjpeg.sourceforge.net/ WML -- Website META Language Copyright 1996-2001 Ralf S. Engelschall Copyright 1999-2001 Denis Barbier http://thewml.org/
41
Visualization and Analysis of Quantum Chemical and Molecular Dynamics Data with VMD
If you want to cite or bookmark this page please use the URL http://www.theochem.ruhr-uni-bochum.de/go/cpmd-vmd.html as the underlying link might change in the future.
The scripts linked to on this page are made available free of charge for personal use. They still are copyrighted by the author(s) mentioned in the individual files. They may not be used in any commercial software without prior agreement of the author(s). These scripts are distributed in the hope that they will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. All trademarks and registered trademarks are the properties of their respective holders. All data files used for the examples presented here are available for download. You can look a all the files, if you follow the link to the files section.
42