diff --git a/doc/index.html b/doc/index.html
index e4d27c4..01b067f 100644
--- a/doc/index.html
+++ b/doc/index.html
@@ -53,9 +53,11 @@
5 Utility Functions
@@ -100,10 +102,7 @@ Scope and usage
larger systems by replicating a given unitcell.
Version
-This documentation describes version 1.10 of the topotools plugin.
-This is the final version of this plugin unless a new maintainer is found. If you
-are interested in maintaining TopoTools, please contact the original author,
-Axel Kohlmeyer via email.
+This documentation describes the version of the topotools plugin that is displayed in the title.
Main command interface
This is the middleware part of the package that provides abstract operations on top of the low-level API. This is modeled after the example of the internal mol or molinfo commands, or the pbc command from the PBCTools package to provide a somewhat consistent interface to the functionality. All command lines start with the topo keyword and then take a subcommand name to determine the actual functionality that is requested.
topo <command> [args...] <flags>
@@ -208,7 +207,11 @@ readlammpsdata <fil
Read in atom coordinates, properties, bond, angle, dihedral and other related topology info from a LAMMPS data file, i.e. a file suitable for the read_data command. This can be used to check a data file for its validity, for manipulations from within VMD, or to generate a .psf file to be used for visualization of .dcd or .xtc format trajectory files in VMD. The 'atom style' is the value given to the atom_style command in the LAMMPS input file (by default TopoTools will try to infer the atom style from information embedded in the data file as comments; if no such hints are present, it will use 'full'). This subcommand creates a new molecule in VMD and returns its molecule id or -1 in case of failure. The -sel parameter is currently ignored.
writelammpsdata <file name> [typelabels] [<atom style>]
Write out atom coordinates, properties, bond, angle, dihedral and other related topology info stored inside VMD to a LAMMPS data file, i.e. a file suitable for the read_data command. Using the optional 'typelabels' flag will trigger writing a data file with typelabel support requiring LAMMPS version 15Sep2022 or later. By default a traditional data file with numerical types will be written. This this way VMD can be used to build LAMMPS input with Tcl scripting and convert existing inputs from other MD codes to be used in LAMMPS. For some examples, please see the TopoTools tutorials. The 'atom style' is the value you want to give to the atom_style command in the LAMMPS input file (default is 'full'). Only data that is present will be written and non-zero box sizes are required.
-readvarxyz <file name>
+readlammpsmol <file name>
+Read in atom coordinates, properties, bond, angle, dihedral and other related topology info from a LAMMPS molecule template file, i.e. a file suitable for the molecule command. This can be used to check a molecule template for its validity (e.g. for fixes that utilize molecule templates), or for manipulations from within VMD. This subcommand creates a new molecule in VMD and returns its molecule id or -1 in case of failure. The -sel parameter is currently ignored.
+writelammpsmol <file name> [typelabels]
+Write out atom coordinates, properties, bond, angle, dihedral and other related topology info stored inside VMD to a LAMMPS molecule template file, i.e. a file suitable for the molecule command. Using the optional 'typelabels' flag will trigger writing a molecule template file with typelabel support requiring LAMMPS version 15Sep2022 or later. By default a traditional molecule template with numerical types will be written. This way VMD can be used to build LAMMPS input with Tcl scripting and convert existing inputs from other MD codes to be used in LAMMPS. Only data that is present will be written.
+readvarxyz <file name> [first|last|step <frame]
Read in an xyz-format trajectory file (in xmol style) with a varying number of atoms per frame. This format is normally not supported in VMD and the script circumvents the restriction by automatically adding a sufficient number of dummy particles. Whether an atom is actually present in a given frame or not is flagged by the value of the corresponding user field, which is set to either 1.0 or -1.0, respectively. For efficiency reasons the atoms are sorted by type, thus atom order and bonding is not preserved. This subcommand creates a new molecule and returns its molecule id or -1, in case of failure.
writevarxyz <file name> [selmod <sel>] [first|last|step <frame]
Write out an xyz-format trajectory file with a varying number of atoms per frame. This is the counterpart to the readvarxyz subcommand. The optional selection string defines how atoms for each frame have to be selected. If not given, as selection string of "user > 0" is assumed.
diff --git a/topolammps.tcl b/topolammps.tcl
index ab353e8..27f3278 100644
--- a/topolammps.tcl
+++ b/topolammps.tcl
@@ -1463,3 +1463,465 @@ proc ::TopoTools::writelammpscoeffhint {fp sel typelabels type} {
puts $fp ""
return
}
+
+
+# read lammps molecule file
+proc ::TopoTools::readlammpsmol {filename {flags none}} {
+ global M_PI
+ if {[catch {open $filename r} fp]} {
+ vmdcon -err "readlammpsmol: problem opening data file: $fp\n"
+ return -1
+ }
+
+ # parse lammps header section.
+ array set lammps [readlammpsmolheader $fp]
+ if {$lammps(atoms) < 1} {
+ vmdcon -err "readlammpsdata: failed to parse lammps data header. abort."
+ return -1
+ }
+
+ # create an empty molecule and timestep
+ set mol -1
+ if {[catch {mol new atoms $lammps(atoms)} mol]} {
+ vmdcon -err "readlammpsdata: problem creating empty molecule: $mol"
+ return -1
+ } else {
+ animate dup $mol
+ }
+ mol rename $mol [file tail $filename]
+ set sel [atomselect $mol all]
+ set boxdim {}
+ set atomidmap {}
+ set atommasses {}
+ set atomlabels {}
+ set bondlabels {}
+ set anglelabels {}
+ set dihedrallabels {}
+ set improperlabels {}
+
+ # now loop through the file until we find a known section header.
+ # then call a subroutine that parses this section. those subroutines
+ # all take the current line number as last argument and return the
+ # new value, so we can keep track of the current line and print more
+ # useful error messages or warnings.
+ set lineno $lammps(lineno)
+ while {[gets $fp line] >= 0} {
+ incr lineno
+ if {[regexp {^\s*Types} $line ]} {
+ set lineno [readlammpsmoltypes $fp $sel $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Types section."
+ return -1
+ }
+ set atomidmap {}
+ foreach id [$sel get user] {
+ lappend atomidmap [expr {int($id + 0.5)}]
+ }
+ } elseif {[regexp {^\s*Coords} $line ]} {
+ set lineno [readlammpsmolcoords $fp $sel $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Coords section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Charges} $line ]} {
+ set lineno [readlammpsmolcharges $fp $sel $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Charges section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Molecules} $line ]} {
+ set lineno [readlammpsmolmolecules $fp $sel $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Molecules section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Masses} $line ]} {
+ set lineno [readlammpsmolmasses $fp $sel $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Masses section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Bonds} $line ]} {
+ if {[llength $atomidmap] < 1} {
+ vmdcon -err "readlammpsdata: Atoms section must come before Bonds in data file"
+ return -1
+ }
+ set lineno [readlammpsbonds $fp $sel $lammps(bonds) $atomidmap bondlabels $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Bonds section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Angles} $line ]} {
+ if {[llength $atomidmap] < 1} {
+ vmdcon -err "readlammpsdata: Atoms section must come before Angles in data file"
+ return -1
+ }
+ set lineno [readlammpsangles $fp $sel $lammps(angles) $atomidmap anglelabels $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Angles section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Dihedrals} $line ]} {
+ if {[llength $atomidmap] < 1} {
+ vmdcon -err "readlammpsdata: Atoms section must come before Dihedrals in data file"
+ return -1
+ }
+ set lineno [readlammpsdihedrals $fp $sel $lammps(dihedrals) $atomidmap dihedrallabels $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Dihedrals section."
+ return -1
+ }
+ } elseif {[regexp {^\s*Impropers} $line ]} {
+ if {[llength $atomidmap] < 1} {
+ vmdcon -err "readlammpsdata: Atoms section must come before Impropers in data file"
+ return -1
+ }
+ set lineno [readlammpsimpropers $fp $sel $lammps(impropers) $atomidmap improperlabels $lineno]
+ if {$lineno < 0} {
+ vmdcon -err "readlammpsdata: error reading Impropers section."
+ return -1
+ }
+ } elseif { [regexp {^\s*(\#.*|)$} $line ] } {
+ # skip empty lines silently
+ } elseif {[regexp {^\s*Shake Flags|^\s*Shake Atoms|^\s*Diameters \
+ |^\s*Special Bond Counts|^\s*Special Bond Counts|^\s*Shake Bond Types} $line ]} {
+ vmdcon -err "readlammpsmol: section $line not yet supported"
+ return -1
+ } else {
+ vmdcon -err "readlammpsdata: unkown header line: $lineno : $line "
+ vmdcon -err "readlammpsdata: cannot continue. "
+ return -1
+ }
+ set lammps(lineno) $lineno
+ }
+ close $fp
+
+ # apply masses. Atoms section sets a default of 1.0.
+ # since the Masses section can appear before the Atoms section
+ # we have to set it here after the parsing.
+ if {[llength atommasses] > 0} {
+ foreach {t m} $atommasses {
+ set msel [atomselect $mol "type '$t'"]
+ $msel set mass $m
+ $msel delete
+ }
+ }
+ mol reanalyze $mol
+ variable newaddsrep
+ if {$newaddsrep} {
+ adddefaultrep $mol
+ }
+ $sel delete
+ return $mol
+}
+
+
+# read lammps molecule header section from opened file
+# and return as an array.
+proc ::TopoTools::readlammpsmolheader {fp} {
+ array set lammps {
+ atoms 0 atomtypes 0 bonds 0 bondtypes 0 angles 0 angletypes 0
+ dihedrals 0 dihedraltypes 0 impropers 0 impropertypes 0 xtrabond 0
+ xlo -0.5 xhi 0.5 ylo -0.5 yhi 0.5 zlo -0.5 zhi 0.5 xy {} xz {} yz {}
+ lineno 0 cgcmm 0 triclinic 0 style unknown typemass 1
+ }
+ set x {}
+
+ vmdcon -info "parsing LAMMPS header."
+
+ # first header line is skipped by LAMMPS. so we put a flag to
+ # detect CMM style CG data files with additional information.
+ gets $fp line
+ set lineno 1
+ set offs [tell $fp]
+ set lammps(lineno) $lineno
+
+ #
+ while {[gets $fp line] >= 0} {
+ incr lineno
+ if { [ regexp {^\s*(\d+)\s+atoms} $line x lammps(atoms) ] } {
+ } elseif { [regexp {^\s*(\d+)\s+bonds} $line x lammps(bonds) ] } {
+ } elseif { [regexp {^\s*(\d+)\s+angles} $line x lammps(angles)] } {
+ } elseif { [regexp {^\s*(\d+)\s+dihedrals} $line x lammps(dihedrals)] } {
+ } elseif { [regexp {^\s*(\d+)\s+impropers} $line x lammps(impropers)] } {
+ } elseif { [regexp {^\s*(\#.*|)$} $line ] } {
+ } elseif {[regexp {^\s*(Types|Coords|Charges|Molecules|Masses|Shapes|Dipoles|Bonds|Angles|Dihedrals|Impropers)} $line ]} {
+ seek $fp $offs start
+ break
+ } else {
+ vmdcon -warn "readlammpsheader: skipping unkown header line: $lineno : $line "
+ }
+ set offs [tell $fp]
+ set lammps(lineno) $lineno
+ }
+
+ return [array get lammps]
+}
+
+
+# parse mol coords section
+proc ::TopoTools::readlammpsmolcoords {fp sel lineno} {
+ set numatoms [$sel num]
+
+ vmdcon -info "parsing LAMMPS Coords section."
+
+ set atomdata {}
+ set x 0.0
+ set y 0.0
+ set z 0.0
+
+ set curatoms 0
+ while {[gets $fp line] >= 0} {
+ incr lineno
+
+ set atomid 0
+
+ if { [regexp {^\s*(\#.*|)$} $line ] } {
+ # skip empty lines.
+ } else {
+ incr curatoms
+ lassign $line atomid x y z
+ lappend atomdata [list $atomid $x $y $z]
+ }
+ if {$curatoms >= $numatoms} break
+ }
+ $sel set {user x y z} \
+ [lsort -integer -index 0 $atomdata]
+ return $lineno
+}
+
+
+# parse mol types section
+proc ::TopoTools::readlammpsmoltypes {fp sel lineno} {
+ set numatoms [$sel num]
+
+ vmdcon -info "parsing LAMMPS Types section."
+
+ set atomdata {}
+
+ set curatoms 0
+ while {[gets $fp line] >= 0} {
+ incr lineno
+
+ set atomid 0
+ set atomtype 0
+
+ if { [regexp {^\s*(\#.*|)$} $line ] } {
+ # skip empty lines.
+ } else {
+ incr curatoms
+ lassign $line atomid atomtype
+ lappend atomdata [list $atomid $atomtype]
+ }
+ if {$curatoms >= $numatoms} break
+ }
+ vmdcon -info "applying atoms data. sorted by atom id."
+ $sel set {user type} \
+ [lsort -integer -index 0 $atomdata]
+ return $lineno
+}
+
+
+# parse mol charges section
+proc ::TopoTools::readlammpsmolcharges {fp sel lineno} {
+ set numatoms [$sel num]
+
+ vmdcon -info "parsing LAMMPS Charges section."
+
+ set atomdata {}
+ set charge 0
+
+ set curatoms 0
+ while {[gets $fp line] >= 0} {
+ incr lineno
+
+ set atomid 0
+
+ if { [regexp {^\s*(\#.*|)$} $line ] } {
+ # skip empty lines.
+ } else {
+ incr curatoms
+ lassign $line atomid charge
+ lappend atomdata [list $atomid $charge]
+ }
+ if {$curatoms >= $numatoms} break
+ }
+ $sel set {user charge} \
+ [lsort -integer -index 0 $atomdata]
+ return $lineno
+}
+
+# parse mol Molecules section. put LAMMPS molecule ID into VMD resid
+proc ::TopoTools::readlammpsmolmolecules {fp sel lineno} {
+ set numatoms [$sel num]
+
+ vmdcon -info "parsing LAMMPS Molecules section."
+
+ set atomdata {}
+ set resid 0
+
+ set curatoms 0
+ while {[gets $fp line] >= 0} {
+ incr lineno
+
+ set atomid 0
+
+ if { [regexp {^\s*(\#.*|)$} $line ] } {
+ # skip empty lines.
+ } else {
+ incr curatoms
+ lassign $line atomid resid
+ lappend atomdata [list $atomid $resid]
+ }
+ if {$curatoms >= $numatoms} break
+ }
+ $sel set {user resid} \
+ [lsort -integer -index 0 $atomdata]
+ return $lineno
+}
+
+# parse mol masses section
+proc ::TopoTools::readlammpsmolmasses {fp sel lineno} {
+ set numatoms [$sel num]
+
+ vmdcon -info "parsing LAMMPS Masses section."
+
+ set mass 0
+
+ set curatoms 0
+ while {[gets $fp line] >= 0} {
+ incr lineno
+
+ set atomid 0
+
+ if { [regexp {^\s*(\#.*|)$} $line ] } {
+ # skip empty lines.
+ } else {
+ incr curatoms
+ lassign $line atomid mass
+ }
+ if {$curatoms >= $numatoms} break
+ }
+ return $lineno
+}
+
+# export internal structure data to a LAMMPS molecule file.
+# Arguments:
+# mol = molecule id with matching coordinate data
+# filename = name of data file
+# sel = selection function for the subset to be written out.
+# flags = more flags. (currently not used)
+proc ::TopoTools::writelammpsmol {mol filename typelabels sel {flags none}} {
+ if {[catch {open $filename w} fp]} {
+ vmdcon -err "writelammpsmol: problem opening molecule file: $fp\n"
+ return -1
+ }
+
+ # initialize system default settings
+ array set lammps {
+ atoms 0 bonds 0 angles 0 dihedrals 0 impropers 0 typelabels 0
+ }
+
+ # gather available system information
+ set lammps(atoms) [$sel num]
+ set lammps(bonds) [bondinfo numbonds $sel]
+ set lammps(angles) [angleinfo numangles $sel]
+ set lammps(dihedrals) [dihedralinfo numdihedrals $sel]
+ set lammps(impropers) [improperinfo numimpropers $sel]
+ set lammps(atomtypes) [llength [lsort -ascii -unique [$sel get {type}]]]
+ set lammps(bondtypes) [bondinfo numbondtypes $sel]
+ set lammps(angletypes) [angleinfo numangletypes $sel]
+ set lammps(dihedraltypes) [dihedralinfo numdihedraltypes $sel]
+ set lammps(impropertypes) [improperinfo numimpropertypes $sel]
+ set lammps(typelabels) $typelabels
+
+ # write out supported data file sections
+ writelammpsmolheader $fp [array get lammps]
+
+ set atomidmap [$sel list]
+
+ writelammpsmolcoords $fp $sel
+
+ writelammpsmoltypes $fp $sel
+
+ writelammpsmolcharges $fp $sel
+
+ if {$lammps(bonds) > 0} {
+ writelammpsbonds $fp $sel $atomidmap $typelabels
+ }
+ if {$lammps(angles) > 0} {
+ writelammpsangles $fp $sel $atomidmap $typelabels
+ }
+ if {$lammps(dihedrals) > 0} {
+ writelammpsdihedrals $fp $sel $atomidmap $typelabels
+ }
+ if {$lammps(impropers) > 0} {
+ writelammpsimpropers $fp $sel $atomidmap $typelabels
+ }
+ close $fp
+ return 0
+}
+
+# write lammps mol header section to open file
+proc ::TopoTools::writelammpsmolheader {fp flags} {
+ variable version
+ array set lammps $flags
+ # first header line is skipped.
+ puts $fp "LAMMPS molecule file generated by VMD/TopoTools v$version on [clock format [clock seconds]]"
+
+ foreach key {atoms bonds angles dihedrals impropers} {
+ puts $fp [format " %d %s" $lammps($key) $key]
+ }
+
+ puts $fp ""
+ return
+}
+
+# write mol types section
+proc ::TopoTools::writelammpsmoltypes {fp sel} {
+ global M_PI
+
+ puts $fp " Types\n"
+ set typemap [lsort -unique -ascii [$sel get type]]
+ set atomid 0
+ foreach adat [$sel get {type}] {
+ lassign $adat type
+ set atomtype [lsearch -sorted -ascii $typemap $type]
+ incr atomid
+ incr atomtype
+ # why $atomtype not correct? --jrg
+ puts $fp [format "%d %s" $atomid $type]
+ }
+ puts $fp ""
+ return
+}
+
+# write mol charges section
+proc ::TopoTools::writelammpsmolcharges {fp sel} {
+ global M_PI
+
+ puts $fp " Charges\n"
+ set atomid 0
+ foreach adat [$sel get {charge}] {
+ lassign $adat charge
+ incr atomid
+ puts $fp [format "%d %.6f" $atomid $charge]
+ }
+ puts $fp ""
+ return
+}
+
+# write mol coords section
+proc ::TopoTools::writelammpsmolcoords {fp sel} {
+ global M_PI
+
+ puts $fp " Coords\n"
+ set atomid 0
+ foreach adat [$sel get {x y z}] {
+ lassign $adat x y z
+ incr atomid
+ puts $fp [format "%d %.6f %.6f %.6f" $atomid $x $y $z]
+ }
+ puts $fp ""
+ return
+}
diff --git a/topotools.tcl b/topotools.tcl
index a2b2a12..91f74ab 100644
--- a/topotools.tcl
+++ b/topotools.tcl
@@ -183,7 +183,7 @@ proc ::TopoTools::usage {} {
vmdcon -info " section is the value given to the 'atom_style' LAMMPS keyword. Default value is 'full'."
vmdcon -info " Only data that is present is written. "
vmdcon -info ""
- vmdcon -info " readvarxyz "
+ vmdcon -info " readvarxyz \[first|last|step \]"
vmdcon -info " read an xmol/xyz format trajectory with a varying numer of particles."
vmdcon -info " This is normally not supported by VMD and the script circumvents this"
vmdcon -info " restriction by automatically adding dummy particles and then indicating"
@@ -314,19 +314,20 @@ proc ::TopoTools::topo { args } {
}
# check whether we have a valid command.
- set validcmd {readvarxyz writevarxyz readlammpsdata writelammpsdata
- writegmxtop help numatoms numatomtypes atomtypenames guessatom
- getbondlist bondtypenames numbondtypes numbonds setbondlist
- retypebonds clearbonds guessbonds addbond delbond getanglelist
- angletypenames numangletypes numangles setanglelist retypeangles
- clearangles guessangles addangle delangle sortangles getdihedrallist
+ set validcmd {readlammpsmol writelammpsmol readvarxyz
+ writevarxyz readlammpsdata writelammpsdata writegmxtop help
+ numatoms numatomtypes atomtypenames guessatom getbondlist
+ bondtypenames numbondtypes numbonds setbondlist retypebonds
+ clearbonds guessbonds addbond delbond getanglelist angletypenames
+ numangletypes numangles setanglelist retypeangles clearangles
+ guessangles addangle delangle sortangles getdihedrallist
dihedraltypenames numdihedraltypes numdihedrals setdihedrallist
retypedihedrals cleardihedrals guessdihedrals adddihedral
deldihedral sortdihedrals getimproperlist impropertypenames
numimpropertypes numimpropers setimproperlist retypeimpropers
- clearimpropers guessimpropers addimproper delimproper sortimpropers
- getcrosstermlist numcrossterms setcrosstermlist clearcrossterms
- addcrossterm delcrossterm}
+ clearimpropers guessimpropers addimproper delimproper
+ sortimpropers getcrosstermlist numcrossterms setcrosstermlist
+ clearcrossterms addcrossterm delcrossterm}
if {[lsearch -exact $validcmd $cmd] < 0} {
vmdcon -err "Unknown topotools command '$cmd'"
usage
@@ -355,9 +356,23 @@ proc ::TopoTools::topo { args } {
return $retval
}
+ if {[string equal $cmd readlammpsmol]} {
+ if {[llength $newargs] < 1} {
+ vmdcon -err "Not enough arguments for 'topo readlammpsdata'"
+ usage
+ return
+ }
+ set fname [lindex $newargs 0]
+
+ set retval [readlammpsmol $fname]
+ citation_reminder
+ return $retval
+ }
+
+
if {[string equal $cmd readvarxyz]} {
set fname [lindex $newargs 0]
- set retval [readvarxyz $fname]
+ set retval [readvarxyz $fname [lrange $newargs 1 end]]
citation_reminder
return $retval
}
@@ -711,6 +726,22 @@ proc ::TopoTools::topo { args } {
set retval [writelammpsdata $molid $fname $typelabels $style $sel]
}
+ writelammpsmol { ;# NOTE: readlammpsdata is handled above to bypass check for sel/molid.
+ if {[llength $newargs] < 1} {
+ vmdcon -err "Not enough arguments for 'topo writelammpsmol'"
+ usage
+ return
+ }
+ set typelabels 0
+ set fname [lindex $newargs 0]
+ if {[llength $newargs] > 1} {
+ if {[lindex $newargs 1] == "typelabels"} {
+ set typelabels 1
+ }
+ }
+ set retval [writelammpsmol $molid $fname $typelabels $sel]
+ }
+
writevarxyz { ;# NOTE: readvarxyz is handled above to bypass check for sel/molid.
if {[llength $newargs] < 1} {
vmdcon -err "Not enough arguments for 'topo writevarxyz'"
diff --git a/topovarxyz.tcl b/topovarxyz.tcl
index 22d6be2..65f41d2 100644
--- a/topovarxyz.tcl
+++ b/topovarxyz.tcl
@@ -19,19 +19,36 @@
#
# Arguments:
# filename = name of data file
-# flags = more flags. (currently not used)
-proc ::TopoTools::readvarxyz {filename {flags none}} {
+# flags = more flags
+proc ::TopoTools::readvarxyz {filename {flags {}}} {
if {[catch {open $filename r} fp]} {
vmdcon -err "readvarxyz: problem opening xyz file: $fp\n"
return -1
}
# initialize local variables
- set nframes 0 ; # total number of frames
+ set nframes 0 ; # number of frames read in
+ set tframes 0 ; # total number of frames in file
+ set step 1 ; # default step size when reading trajectory
+ set first 1 ; # default first frame to read in
+ set last -1 ; # default last frame to read in
set typemap {} ; # atom type map
set typecount {} ; # atom type count for one frame
set maxcount {} ; # max. atom type count for all frames
array set traj {} ; # temporary trajectory storage
+
+ # parse optional flags
+ foreach {key value} $flags {
+ switch -- $key {
+ step {set step $value}
+ first {set first $value}
+ last {set last $value}
+ default {
+ vmdcon -err "writevarxyz: unknown flag: $key"
+ return -1
+ }
+ }
+ }
# to be able to determine the number of dummy atoms we first
# have to parse and store away the whole trajectory and while
@@ -51,6 +68,18 @@ proc ::TopoTools::readvarxyz {filename {flags none}} {
break
}
+ incr tframes
+ if {$tframes == $last+1} break
+ if {fmod($tframes-$first+1,$step) != 0 || $tframes < $first} {
+ for {set i 0} {$i < $numlines} {incr i} {
+ if {[catch {gets $fp line} msg]} {
+ vmdcon -err "readvarxyz: error reading frame $nframes of xyz file: $msg. "
+ break
+ }
+ }
+ continue
+ }
+
# collect data for this frame.
set frame {}
for {set i 0} {$i < $numlines} {incr i} {
@@ -171,7 +200,7 @@ proc ::TopoTools::readvarxyz {filename {flags none}} {
# Arguments:
# filename = name of data file
-# flags = more flags. (currently not used)
+# flags = more flags
proc ::TopoTools::writevarxyz {filename mol sel {flags {}}} {
if {[catch {open $filename w} fp]} {
vmdcon -err "writevarxyz: problem opening xyz file: $fp\n"