Random cosmetic fixes
[alexxy/gromacs.git] / src / external / vmd_molfile / molfile_plugin.h
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  * This file is part of Gromacs        Copyright (c) 1991-2008
5  * David van der Spoel, Erik Lindahl, Berk Hess, University of Groningen.
6  *
7  * This program is free software; you can redistribute it and/or
8  * modify it under the terms of the GNU General Public License
9  * as published by the Free Software Foundation; either version 2
10  * of the License, or (at your option) any later version.
11  *
12  * To help us fund GROMACS development, we humbly ask that you cite
13  * the research papers on the package. Check out http://www.gromacs.org
14  * 
15  * And Hey:
16  * Gnomes, ROck Monsters And Chili Sauce
17  */
18
19 /***************************************************************************
20  *cr
21  *cr            (C) Copyright 1995-2006 The Board of Trustees of the
22  *cr                        University of Illinois
23  *cr                         All Rights Reserved
24  *cr
25
26 Developed by:           Theoretical and Computational Biophysics Group
27                         University of Illinois at Urbana-Champaign
28                         http://www.ks.uiuc.edu/
29
30 Permission is hereby granted, free of charge, to any person obtaining a copy of
31 this software and associated documentation files (the Software), to deal with
32 the Software without restriction, including without limitation the rights to
33 use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies
34 of the Software, and to permit persons to whom the Software is furnished to
35 do so, subject to the following conditions:
36
37 Redistributions of source code must retain the above copyright notice,
38 this list of conditions and the following disclaimers.
39
40 Redistributions in binary form must reproduce the above copyright notice,
41 this list of conditions and the following disclaimers in the documentation
42 and/or other materials provided with the distribution.
43
44 Neither the names of Theoretical and Computational Biophysics Group,
45 University of Illinois at Urbana-Champaign, nor the names of its contributors
46 may be used to endorse or promote products derived from this Software without
47 specific prior written permission.
48
49 THE SOFTWARE IS PROVIDED AS IS, WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
50 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
51 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.  IN NO EVENT SHALL
52 THE CONTRIBUTORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR
53 OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE,
54 ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR
55 OTHER DEALINGS WITH THE SOFTWARE.
56  ***************************************************************************/
57
58 /***************************************************************************
59  * RCS INFORMATION:
60  *
61  *      $RCSfile: molfile_plugin.h,v $
62  *      $Author: saam $       $Locker:  $             $State: Exp $
63  *      $Revision: 1.91 $       $Date: 2009/07/28 21:54:30 $
64  *
65  ***************************************************************************/
66
67 /** @file 
68  * API for C extensions to define a way to load structure, coordinate,
69  * trajectory, and volumetric data files  
70  */ 
71
72 #ifndef MOL_FILE_PLUGIN_H
73 #define MOL_FILE_PLUGIN_H
74
75 #include "vmdplugin.h"
76
77 /**
78  * Define a common plugin type to be used when registering the plugin.
79  */
80 #define MOLFILE_PLUGIN_TYPE "mol file reader"
81
82 /**
83  * File converter plugins use the same API  but register under a different
84  * type so that regular file readers can have priority.
85  */
86 #define MOLFILE_CONVERTER_PLUGIN_TYPE "mol file converter"
87
88 /* File plugin symbolic constants for better code readability */
89 #define MOLFILE_SUCCESS           0   /**< succeeded in reading file      */
90 #define MOLFILE_EOF              (-1)   /**< end of file                    */
91 #define MOLFILE_ERROR            (-1)   /**< error reading/opening a file   */
92 #define MOLFILE_NOSTRUCTUREDATA  (-2)   /**< no structure data in this file */
93
94 #define MOLFILE_NUMATOMS_UNKNOWN (-1)   /**< unknown number of atoms       */
95 #define MOLFILE_NUMATOMS_NONE     0   /**< no atoms in this file type    */
96
97 /**
98  * Maximum string size macro
99  */
100 #define MOLFILE_BUFSIZ           81   /**< maximum chars in string data  */
101 #define MOLFILE_BIGBUFSIZ        4096 /**< maximum chars in long strings */
102
103 #define MOLFILE_MAXWAVEPERTS     25   /**< maximum number of wavefunctions
104                                        *   per timestep */
105
106 #define MOLFILE_QM_STATUS_UNKNOWN  (-1)
107 #define MOLFILE_QM_OPT_CONVERGED    0
108 #define MOLFILE_QM_SCF_NOT_CONV     1
109 #define MOLFILE_QM_OPT_NOT_CONV     2
110 #define MOLFILE_QM_FILE_TRUNCATED   3
111
112
113 /**
114  * File level comments, origin information, and annotations.
115  */
116 typedef struct {
117   char database[81];   /**< database of origin, if any        */
118   char accession[81];  /**< database accession code, if any   */
119   char date[81];       /**< date/time stamp for this data     */
120   char title[81];      /**< brief title for this data         */
121   int remarklen;       /**< length of remarks string          */
122   char *remarks;       /**< free-form remarks about data      */
123 } molfile_metadata_t;
124
125
126 /* 
127  * Struct for specifying atoms in a molecular structure.  The first 
128  * six components are required, the rest are optional and their presence is 
129  * indicating by setting the corresponding bit in optsflag.  When omitted,
130  * the application (for read_structure) or plugin (for write_structure) 
131  * must be able to supply default values if the missing parameters are 
132  * part of its internal data structure.
133  * Note that it is not possible to specify coordinates with this structure.
134  * This is intentional; all coordinate I/O is done with the read_timestep and 
135  * write_timestep functions. 
136  */
137
138 /**
139  * Per-atom attributes and information.
140  */
141 typedef struct {
142   /* these fields absolutely must be set or initialized to empty */
143   char name[16];      /**< required atom name string             */
144   char type[16];      /**< required atom type string             */
145   char resname[8];    /**< required residue name string          */
146   int resid;          /**< required integer residue ID           */
147   char segid[8];      /**< required segment name string, or ""   */
148   char chain[2];      /**< required chain name, or ""            */
149
150   /* rest are optional; use optflags to specify what's present   */
151   char altloc[2];     /**< optional PDB alternate location code  */
152   char insertion[2];  /**< optional PDB insertion code           */
153   float occupancy;    /**< optional occupancy value              */
154   float bfactor;      /**< optional B-factor value               */
155   float mass;         /**< optional mass value                   */
156   float charge;       /**< optional charge value                 */
157   float radius;       /**< optional radius value                 */
158   int atomicnumber;   /**< optional element atomic number        */
159 } molfile_atom_t;
160
161 /*@{*/
162 /** Plugin optional data field availability flag */
163 #define MOLFILE_NOOPTIONS     0x0000 /**< no optional data                 */
164 #define MOLFILE_INSERTION     0x0001 /**< insertion codes provided         */
165 #define MOLFILE_OCCUPANCY     0x0002 /**< occupancy data provided          */
166 #define MOLFILE_BFACTOR       0x0004 /**< B-factor data provided           */
167 #define MOLFILE_MASS          0x0008 /**< Atomic mass provided             */
168 #define MOLFILE_CHARGE        0x0010 /**< Atomic charge provided           */
169 #define MOLFILE_RADIUS        0x0020 /**< Atomic VDW radius provided       */
170 #define MOLFILE_ALTLOC        0x0040 /**< Multiple conformations present   */
171 #define MOLFILE_ATOMICNUMBER  0x0080 /**< Atomic element number provided   */
172 #define MOLFILE_BONDSSPECIAL  0x0100 /**< Only non-standard bonds provided */
173 #define MOLFILE_BADOPTIONS    0xFFFFFFFF /**< Detect badly behaved plugins */
174                               
175 /*@}*/
176
177 /*@{*/
178 /** Plugin optional data field availability flag */
179 #define MOLFILE_QMTS_NOOPTIONS     0x0000 /**< no optional data               */
180 #define MOLFILE_QMTS_GRADIENT      0x0001 /**< energy gradients provided      */
181 #define MOLFILE_QMTS_SCFITER       0x0002
182 /*@}*/
183
184 #if vmdplugin_ABIVERSION > 10
185 typedef struct molfile_timestep_metadata {
186   unsigned int count;                  /**< total # timesteps; -1 if unknown */
187   unsigned int avg_bytes_per_timestep; /** bytes per timestep                */
188   int has_velocities;                  /**< if timesteps have velocities     */
189 } molfile_timestep_metadata_t;
190 #endif
191
192 #if vmdplugin_ABIVERSION > 11
193 typedef struct molfile_qm_timestep_metadata {
194   unsigned int count;                  /**< total # timesteps; -1 if unknown */
195   unsigned int avg_bytes_per_timestep; /** bytes per timestep                */
196   int has_gradient;                    /**< if timestep contains gradient    */
197   int num_scfiter;                     /**< # scf iterations for this ts     */
198   int num_orbitals_per_wavef[MOLFILE_MAXWAVEPERTS]; /**< # orbitals for each wavefunction */
199   int has_orben_per_wavef[MOLFILE_MAXWAVEPERTS]; /**< orbital energy flags */
200   int has_occup_per_wavef[MOLFILE_MAXWAVEPERTS]; /**< orbital occupancy flags */
201   int num_wavef ;                      /**< # wavefunctions in this ts     */
202   int wavef_size;                      /**< size of one wavefunction 
203                                         *   (# of gaussian basis fctns)    */
204   int num_charge_sets;            /**< # of charge values per atom */
205 } molfile_qm_timestep_metadata_t;
206 #endif
207
208 /*
209  * Per-timestep atom coordinates and periodic cell information
210  */ 
211 typedef struct {
212   float *coords;        /**< coordinates of all atoms, arranged xyzxyzxyz   */
213 #if vmdplugin_ABIVERSION > 10
214   float *velocities;    /**< space for velocities of all atoms; same layout */
215                         /**< NULL unless has_velocities is set              */
216 #endif
217
218   /*@{*/   
219   /**
220    * Unit cell specification of the form A, B, C, alpha, beta, gamma.
221    * notes: A, B, C are side lengths of the unit cell
222    * alpha = angle between b and c
223    *  beta = angle between a and c
224    * gamma = angle between a and b
225    */ 
226   float A, B, C, alpha, beta, gamma; 
227   /*@}*/   
228
229 #if vmdplugin_ABIVERSION > 10
230   double physical_time; /**< physical time point associated with this frame */
231 #endif
232 } molfile_timestep_t;
233
234
235 /**
236  * Metadata for volumetric datasets, read initially and used for subsequent
237  * memory allocations and file loading.  
238  */
239 typedef struct {
240   char dataname[256];   /**< name of volumetric data set                    */
241   float origin[3];      /**< origin: origin of volume (x=0, y=0, z=0 corner */
242
243   /*
244    * x/y/z axis:
245    * These the three cell sides, providing both direction and length
246    * (not unit vectors) for the x, y, and z axes.  In the simplest
247    * case, these would be <size,0,0> <0,size,0> and <0,0,size) for 
248    * an orthogonal cubic volume set.  For other cell shapes these
249    * axes can be oriented non-orthogonally, and the parallelpiped
250    * may have different side lengths, not just a cube/rhombus.
251    */
252   float xaxis[3];       /**< direction (and length) for X axis              */ 
253   float yaxis[3];       /**< direction (and length) for Y axis              */
254   float zaxis[3];       /**< direction (and length) for Z axis              */
255
256   /*
257    * x/y/z size: 
258    * Number of grid cells along each axis.  This is _not_ the
259    * physical size of the box, this is the number of voxels in each
260    * direction, independent of the shape of the volume set. 
261    */
262   int xsize;            /**< number of grid cells along the X axis          */
263   int ysize;            /**< number of grid cells along the Y axis          */
264   int zsize;            /**< number of grid cells along the Z axis          */
265
266   int has_color;        /**< flag indicating presence of voxel color data   */
267 } molfile_volumetric_t;
268
269
270 #if vmdplugin_ABIVERSION > 9
271
272 /**
273  * Sizes of various QM-related data arrays which must be allocated by
274  * the caller (VMD) so that the plugin can fill in the arrays with data.
275  */
276 typedef struct {
277   /* hessian data */
278   int nimag;                    /**< number of imaginary modes */
279   int nintcoords;               /**< number internal coordinates */
280   int ncart;                    /**< number cartesian coordinates */
281
282   /* orbital/basisset data */
283   int num_basis_funcs;          /**< number of uncontracted basis functions in basis array */
284   int num_basis_atoms;          /**< number of atoms in basis set */
285   int num_shells;               /**< total number of atomic shells */
286   int wavef_size;               /**< size of the wavefunction
287                                  *   i.e. size of secular eq. or
288                                  *   # of cartesian contracted
289                                  *   gaussian basis functions */
290
291   /* everything else */
292   int have_sysinfo;
293   int have_carthessian;         /**< hessian in cartesian coords available  */
294   int have_inthessian;          /**< hessian in internal coords available  */
295   int have_normalmodes;         /**< normal modes available  */
296 } molfile_qm_metadata_t;
297
298
299 /**
300  * struct holding the data of hessian/normal mode runs
301  * needed to calculate bond/angle/dihedral force constants
302  * XXX: do we really need doubles here??
303  */
304 typedef struct {
305   double *carthessian;  /**< hessian matrix in cartesian coordinates (ncart)*(ncart)
306                          *   as a single array of doubles (row(1), ...,row(natoms)) */
307   int    *imag_modes;   /**< list(nimag) of imaginary modes */
308   double *inthessian;   /**< hessian matrix in internal coordinates
309                          *   (nintcoords*nintcoords) as a single array of
310                          *   doubles (row(1), ...,row(nintcoords)) */
311   float *wavenumbers;   /**< array(ncart) of wavenumbers of normal modes */
312   float *intensities;   /**< array(ncart) of intensities of normal modes */
313   float *normalmodes;   /**< matrix(ncart*ncart) of normal modes  */
314 } molfile_qm_hessian_t;
315
316
317 /**
318  * struct holding the data for wavefunction/orbitals
319  * needed to generate the volumetric orbital data
320  */
321 typedef struct {
322   int *num_shells_per_atom; /**< number of shells per atom */
323   int *num_prim_per_shell;  /**< number of shell primitives shell */
324
325   float *basis;             /**< contraction coeffients and exponents for
326                              *   the basis functions in the form
327                              *   { exp(1), c-coeff(1), exp(2), c-coeff(2), ....};
328                              *   size=2*num_basis_funcs */
329   int *atomic_number;       /**< atomic numbers (chem. element) of atoms in basis set */
330   int *angular_momentum;    /**< 3 ints per wave function coefficient do describe the 
331                              *   cartesian components of the angular momentum.
332                              *   E.g. S={0 0 0}, Px={1 0 0}, Dxy={1 1 0}, or Fyyz={0 2 1}. 
333                              */
334   int *shell_symmetry;      /**< symmetry type per shell in basis */
335 } molfile_qm_basis_t;
336
337
338 /**
339  * QM run info. Parameters that stay unchanged during a single file.
340  */ 
341 typedef struct {
342   int nproc;             /**< number of processors used. XXX:? */
343   int memory;            /**< amount of memory used in Mbyte. XXX:? */ 
344   int runtype;           /**< flag indicating the calculation method. */
345   int scftype;           /**< SCF type: RHF, UHF, ROHF, GVB or MCSCF wfn. */
346   int status;            /**< indicates whether SCF and geometry optimization
347                           *   have converged properly. */
348   int num_electrons;     /**< number of electrons.    XXX: can be fractional in some DFT codes */
349   int totalcharge;       /**< total charge of system. XXX: can be fractional in some DFT codes */
350   int num_occupied_A;    /**< number of occupied alpha orbitals */
351   int num_occupied_B;    /**< number of occupied beta orbitals */
352
353   double *nuc_charge;    /**< array(natom) containing the nuclear charge of atom i */
354
355   char basis_string[MOLFILE_BUFSIZ];    /**< basis name as "nice" string. */
356   char runtitle[MOLFILE_BIGBUFSIZ];     /**< title of run.                */
357   char geometry[MOLFILE_BUFSIZ];        /**< type of provided geometry,   XXX: remove?
358                                          * e.g. UNIQUE, ZMT, CART, ...    */
359   char version_string[MOLFILE_BUFSIZ];  /**< QM code version information. */
360 } molfile_qm_sysinfo_t;
361
362
363 typedef struct {
364   int   type;               /**< CANONICAL, LOCALIZED, OTHER */
365   int   spin;               /**< 1 for alpha, -1 for beta */
366   int   excitation;         /**< 0 for ground state, 1,2,3,... for excited states */
367   int   multiplicity;       /**< spin multiplicity of the state, zero if unknown */
368   char info[MOLFILE_BUFSIZ]; /**< string for additional type info */
369
370   double energy;            /**< energy of the electronic state.
371                              *   i.e. HF-SCF energy, CI state energy,
372                              *   MCSCF energy, etc. */
373
374   float *wave_coeffs;       /**< expansion coefficients for wavefunction in the
375                              *   form {orbital1(c1),orbital1(c2),.....,orbitalM(cN)} */
376   float *orbital_energies;  /**< list of orbital energies for wavefunction */
377   float *occupancies;       /**< orbital occupancies */
378   int   *orbital_ids;       /**< orbital ID numbers; If NULL then VMD will
379                              *   assume 1,2,3,...num_orbs.     */
380 } molfile_qm_wavefunction_t;
381
382
383 /**
384  * QM per trajectory timestep info
385  */
386 typedef struct {
387   molfile_qm_wavefunction_t *wave; /**< array of wavefunction objects */
388   float  *gradient;         /**< force on each atom (=gradient of energy) */
389
390   double *scfenergies;      /**< scfenergy per trajectory point. */
391   double *charges;          /**< per-atom charges */
392   int    *charge_types;     /**< type of each charge set */
393 } molfile_qm_timestep_t;
394
395
396 /**
397  * QM wavefunctions, and related information 
398  */
399 typedef struct {
400   molfile_qm_hessian_t hess;            /* hessian info */
401   molfile_qm_basis_t   basis;           /* basis set info */
402   molfile_qm_sysinfo_t run;             /* system info  */
403 } molfile_qm_t;
404
405
406 /**
407  *  Enumeration of all of the wavefunction types that can be read
408  *  from QM file reader plugins.
409  *
410  * CANON    = canonical (i.e diagonalized) wavefunction
411  * GEMINAL  = GVB-ROHF geminal pairs
412  * MCSCFNAT = Multi-Configuration SCF natural orbitals
413  * MCSCFOPT = Multi-Configuration SCF optimized orbitals
414  * CINATUR  = Configuration-Interaction natural orbitals
415  * BOYS     = Boys localization
416  * RUEDEN   = Ruedenberg localization
417  * PIPEK    = Pipek-Mezey population localization
418  *
419  * NBO related localizations:
420  * --------------------------
421  * NAO      = Natural Atomic Orbitals
422  * PNAO     = pre-orthogonal NAOs
423  * NBO      = Natural Bond Orbitals
424  * PNBO     = pre-orthogonal NBOs
425  * NHO      = Natural Hybrid Orbitals
426  * PNHO     = pre-orthogonal NHOs
427  * NLMO     = Natural Localized Molecular Orbitals
428  * PNLMO    = pre-orthogonal NLMOs
429  *
430  * UNKNOWN  = Use this for any type not listed here
431  *            You can use the string field for description
432  */
433 enum molfile_qm_wavefunc_type {
434   MOLFILE_WAVE_CANON,    MOLFILE_WAVE_GEMINAL,
435   MOLFILE_WAVE_MCSCFNAT, MOLFILE_WAVE_MCSCFOPT,
436   MOLFILE_WAVE_CINATUR,
437   MOLFILE_WAVE_PIPEK,  MOLFILE_WAVE_BOYS, MOLFILE_WAVE_RUEDEN,
438   MOLFILE_WAVE_NAO,    MOLFILE_WAVE_PNAO, MOLFILE_WAVE_NHO, 
439   MOLFILE_WAVE_PNHO,   MOLFILE_WAVE_NBO,  MOLFILE_WAVE_PNBO, 
440   MOLFILE_WAVE_PNLMO,  MOLFILE_WAVE_NLMO, MOLFILE_WAVE_MOAO, 
441   MOLFILE_WAVE_NATO,   MOLFILE_WAVE_UNKNOWN
442 };
443
444 enum molfile_qm_charge_type {
445   MOLFILE_QMCHARGE_MULLIKEN, MOLFILE_QMCHARGE_LOWDIN,
446   MOLFILE_QMCHARGE_ESP, MOLFILE_QMCHARGE_NPA,
447   MOLFILE_QMCHARGE_UNKNOWN
448 };
449 #endif
450
451
452 /**
453  *  Enumeration of all of the supported graphics objects that can be read
454  *  from graphics file reader plugins.
455  */
456 enum molfile_graphics_type {
457   MOLFILE_POINT,  MOLFILE_TRIANGLE, MOLFILE_TRINORM, MOLFILE_NORMS, 
458   MOLFILE_LINE,   MOLFILE_CYLINDER, MOLFILE_CAPCYL,  MOLFILE_CONE,    
459   MOLFILE_SPHERE, MOLFILE_TEXT,     MOLFILE_COLOR,   MOLFILE_TRICOLOR
460 };
461
462 /**
463  *  Individual graphics object/element data
464  */ 
465 typedef struct {
466   int type;             /* One of molfile_graphics_type */
467   int style;            /* A general style parameter    */
468   float size;           /* A general size parameter     */
469   float data[9];        /* All data for the element     */
470 } molfile_graphics_t;
471
472
473 /*
474  * Types for raw graphics elements stored in files.  Data for each type
475  * should be stored by the plugin as follows:
476
477 type        data                                     style       size
478 ----        ----                                     -----       ----
479 point       x, y, z                                              pixel size
480 triangle    x1,y1,z1,x2,y2,z2,x3,y3,z3                 
481 trinorm     x1,y1,z1,x2,y2,z2,x3,y3,z3                 
482             the next array element must be NORMS
483 tricolor    x1,y1,z1,x2,y2,z2,x3,y3,z3                 
484             the next array elements must be NORMS
485             the following element must be COLOR, with three RGB triples
486 norms       x1,y1,z1,x2,y2,z2,x3,y3,z3                 
487 line        x1,y1,z1,x2,y2,z2                        0=solid     pixel width
488                                                      1=stippled
489 cylinder    x1,y1,z1,x2,y2,z2                        resolution  radius
490 capcyl      x1,y1,z1,x2,y2,z2                        resolution  radius
491 sphere      x1,y1,z1                                 resolution  radius
492 text        x, y, z, up to 24 bytes of text                      pixel size
493 color       r, g, b
494 */
495
496
497 /**
498  * Main file reader API.  Any function in this struct may be NULL
499  * if not implemented by the plugin; the application checks this to determine
500  * what functionality is present in the plugin. 
501  */ 
502 typedef struct {
503   /**
504    * Required header 
505    */
506   vmdplugin_HEAD;
507
508   /**
509    * Filename extension for this file type.  May be NULL if no filename 
510    * extension exists and/or is known.  For file types that match several
511    * common extensions, list them in a comma separated list such as:
512    *  "pdb,ent,foo,bar,baz,ban"
513    * The comma separated list will be expanded when filename extension matching
514    * is performed.  If multiple plugins solicit the same filename extensions,
515    * the one that lists the extension earliest in its list is selected. In the 
516    * case of a "tie", the first one tried/checked "wins".
517    */
518   const char *filename_extension;
519
520   /**
521    * Try to open the file for reading.  Return an opaque handle, or NULL on
522    * failure. Set the number of atoms; if the number of atoms cannot be 
523    * determined, set natoms to MOLFILE_NUMATOMS_UNKNOWN. 
524    * Filetype should be the name under which this plugin was registered;
525    * this is provided so that plugins can provide the same function pointer
526    * to handle multiple file types.
527    */
528   void *(* open_file_read)(const char *filepath, const char *filetype, 
529       int *natoms);
530   
531   /**
532    * Read molecular structure from the given file handle.  atoms is allocated
533    * by the caller and points to space for natoms.
534    * On success, place atom information in the passed-in pointer.  
535    * optflags specifies which optional fields in the atoms will be set by
536    * the plugin.
537    */
538   int (*read_structure)(void *, int *optflags, molfile_atom_t *atoms);
539
540   /**
541    * Read bond information for the molecule.  On success the arrays from
542    * and to should point to the (one-based) indices of bonded atoms.
543    * Each unique bond should be specified only once, so file formats that list
544    * bonds twice will need post-processing before the results are returned to
545    * the caller.
546    * If the plugin provides bond information, but the file loaded doesn't 
547    * actually contain any bond info, the nbonds parameter should be
548    * set to 0 and from/to should be set to NULL to indicate that no bond
549    * information was actually present, and automatic bond search should be
550    * performed.  
551    *
552    * If the plugin provides bond order information, the bondorder array
553    * will contain the bond order for each from/to pair.  If not, the bondorder
554    * pointer should be set to NULL, in which case the caller will provide a 
555    * default bond order value of 1.0.
556    *
557    * If the plugin provides bond type information, the bondtype array
558    * will contain the bond type index for each from/to pair. These numbers
559    * are consecutive integers starting from 0.
560    * the bondtypenames list, contains the corresponding names, if available,
561    * as a NULL string terminated list. nbondtypes is provided for convenience
562    * and consistency checking.
563    *
564    * These arrays must be freed by the plugin in the close_file_read function.
565    * This function can be called only after read_structure().  
566    * Return MOLFILE_SUCCESS if no errors occur. 
567    */
568 #if vmdplugin_ABIVERSION > 14
569   int (*read_bonds)(void *, int *nbonds, int **from, int **to, float **bondorder, 
570                     int **bondtype, int *nbondtypes, char ***bondtypename);
571 #else
572   int (*read_bonds)(void *, int *nbonds, int **from, int **to, float **bondorder);
573 #endif
574
575   /**
576    * XXX this function will be augmented and possibly superceded by a 
577    *     new QM-capable version named read_timestep(), when finished.
578    *
579    * Read the next timestep from the file.  Return MOLFILE_SUCCESS, or 
580    * MOLFILE_EOF on EOF.  If the molfile_timestep_t argument is NULL, then 
581    * the frame should be skipped.  Otherwise, the application must prepare 
582    * molfile_timestep_t by allocating space in coords for the corresponding 
583    * number of coordinates.  
584    * The natoms parameter exists because some coordinate file formats 
585    * (like CRD) cannot determine for themselves how many atoms are in a 
586    * timestep; the app must therefore obtain this information elsewhere
587    * and provide it to the plugin.
588    */
589   int (* read_next_timestep)(void *, int natoms, molfile_timestep_t *);
590
591   /**
592    * Close the file and release all data.  The handle cannot be reused.
593    */
594   void (* close_file_read)(void *);
595    
596   /**
597    * Open a coordinate file for writing using the given header information.
598    * Return an opaque handle, or NULL on failure.  The application must
599    * specify the number of atoms to be written. 
600    * filetype should be the name under which this plugin was registered.
601    */
602   void *(* open_file_write)(const char *filepath, const char *filetype, 
603       int natoms);
604   
605   /**
606    * Write structure information.  Return success.
607    */
608   int (* write_structure)(void *, int optflags, const molfile_atom_t *atoms);
609
610   /**
611    * Write a timestep to the coordinate file.  Return MOLFILE_SUCCESS if no
612    * errors occur.  If the file contains structure information in each 
613    * timestep (like a multi-entry PDB), it will have to cache the information 
614    * from the initial calls from write_structure.
615    */
616   int (* write_timestep)(void *, const molfile_timestep_t *);
617   
618   /**
619    * Close the file and release all data.  The handle cannot be reused.
620    */
621   void (* close_file_write)(void *);
622
623   /**
624    * Retrieve metadata pertaining to volumetric datasets in this file.
625    * Set nsets to the number of volumetric data sets, and set *metadata
626    * to point to an array of molfile_volumetric_t.  The array is owned by
627    * the plugin and should be freed by close_file_read().  The application
628    * may call this function any number of times.
629    */
630   int (* read_volumetric_metadata)(void *, int *nsets, 
631         molfile_volumetric_t **metadata);
632
633   /** 
634    * Read the specified volumetric data set into the space pointed to by 
635    * datablock.  The set is specified with a zero-based index.  The space 
636    * allocated for the datablock must be equal to
637    * xsize * ysize * zsize.  No space will be allocated for colorblock 
638    * unless has_color is nonzero; in that case, colorblock should be
639    * filled in with three RGB floats per datapoint.
640    */
641   int (* read_volumetric_data)(void *, int set, float *datablock, 
642         float *colorblock);
643
644   /**
645    * Read raw graphics data stored in this file.   Return the number of data
646    * elements and the data itself as an array of molfile_graphics_t in the 
647    * pointer provided by the application.  The plugin is responsible for 
648    * freeing the data when the file is closed.
649    */
650   int (* read_rawgraphics)(void *, int *nelem, const molfile_graphics_t **data);
651
652   /**
653    * Read molecule metadata such as what database (if any) this file/data
654    * came from, what the accession code for the database is, textual remarks
655    * and other notes pertaining to the contained structure/trajectory/volume
656    * and anything else that's informative at the whole file level.
657    */ 
658   int (* read_molecule_metadata)(void *, molfile_metadata_t **metadata);
659   
660   /**
661    * Write bond information for the molecule.  The arrays from
662    * and to point to the (one-based) indices of bonded atoms.
663    * Each unique bond will be specified only once by the caller. 
664    * File formats that list bonds twice will need to emit both the 
665    * from/to and to/from versions of each.
666    * This function must be called before write_structure().  
667    *
668    * Like the read_bonds() routine, the bondorder pointer is set to NULL
669    * if the caller doesn't have such information, in which case the 
670    * plugin should assume a bond order of 1.0 if the file format requires
671    * bond order information.
672    *
673    * Support for bond types follows the bondorder rules. bondtype is
674    * an integer array of the size nbonds that contains the bond type
675    * index (consecutive integers starting from 0) and bondtypenames
676    * contain the corresponding strings, in case the naming/numbering
677    * scheme is different from the index numbers.
678    * if the pointers are set to NULL, then this information is not available.
679    * bondtypenames can only be used of bondtypes is also given.
680    * Return MOLFILE_SUCCESS if no errors occur. 
681    */
682 #if vmdplugin_ABIVERSION > 14
683   int (* write_bonds)(void *, int nbonds, int *from, int *to, float *bondorder, 
684                      int *bondtype, int nbondtypes, char **bondtypename);
685 #else
686   int (* write_bonds)(void *, int nbonds, int *from, int *to, float *bondorder);
687 #endif
688
689 #if vmdplugin_ABIVERSION > 9
690   /**
691    * Write the specified volumetric data set into the space pointed to by 
692    * datablock.  The * allocated for the datablock must be equal to
693    * xsize * ysize * zsize.  No space will be allocated for colorblock 
694    * unless has_color is nonzero; in that case, colorblock should be
695    * filled in with three RGB floats per datapoint.
696    */
697   int (* write_volumetric_data)(void *, molfile_volumetric_t *metadata,
698                                 float *datablock, float *colorblock);
699
700 #if vmdplugin_ABIVERSION > 15
701   /** 
702    * Read in Angles, Dihedrals, Impropers, and Cross Terms and optionally types.
703    * (Cross terms pertain to the CHARMM/NAMD CMAP feature) 
704    */
705   int (* read_angles)(void *handle, int *numangles, int **angles, int **angletypes,
706                       int *numangletypes, char ***angletypenames, int *numdihedrals,
707                       int **dihedrals, int **dihedraltypes, int *numdihedraltypes,
708                       char ***dihedraltypenames, int *numimpropers, int **impropers,        
709                       int **impropertypes, int *numimpropertypes, char ***impropertypenames,
710                       int *numcterms, int **cterms, int *ctermcols, int *ctermrows);
711
712   /** 
713    * Write out Angles, Dihedrals, Impropers, and Cross Terms
714    * (Cross terms pertain to the CHARMM/NAMD CMAP feature) 
715    */
716   int (* write_angles)(void *handle, int numangles, const int *angles, const int *angletypes,
717                        int numangletypes, const char **angletypenames, int numdihedrals,
718                        const int *dihedrals, const int *dihedraltypes, int numdihedraltypes,
719                        const char **dihedraltypenames, int numimpropers, 
720                        const int *impropers, const int *impropertypes, int numimpropertypes,
721                        const char **impropertypenames, int numcterms,  const int *cterms, 
722                        int ctermcols, int ctermrows);
723 #else
724   /** 
725    * Read in Angles, Dihedrals, Impropers, and Cross Terms
726    * Forces are in Kcal/mol
727    * (Cross terms pertain to the CHARMM/NAMD CMAP feature, forces are given
728    *  as a 2-D matrix)
729    */
730   int (* read_angles)(void *,
731                 int *numangles,    int **angles,    double **angleforces,
732                 int *numdihedrals, int **dihedrals, double **dihedralforces,
733                 int *numimpropers, int **impropers, double **improperforces,
734                 int *numcterms,    int **cterms, 
735                 int *ctermcols,    int *ctermrows,  double **ctermforces);
736
737   /** 
738    * Write out Angles, Dihedrals, Impropers, and Cross Terms
739    * Forces are in Kcal/mol
740    * (Cross terms pertain to the CHARMM/NAMD CMAP feature, forces are given
741    *  as a 2-D matrix)
742    */
743   int (* write_angles)(void *,
744         int numangles,    const int *angles,    const double *angleforces,
745         int numdihedrals, const int *dihedrals, const double *dihedralforces,
746         int numimpropers, const int *impropers, const double *improperforces,
747         int numcterms,   const int *cterms,    
748         int ctermcols, int ctermrows, const double *ctermforces);
749 #endif
750
751   /**
752    * Retrieve metadata pertaining to QM datasets in this file.
753    */
754   int (* read_qm_metadata)(void *, molfile_qm_metadata_t *metadata);
755
756   /**
757    * Read QM data
758    */
759   int (* read_qm_rundata)(void *, molfile_qm_t *qmdata);
760
761   /**
762    * Read the next timestep from the file.  Return MOLFILE_SUCCESS, or 
763    * MOLFILE_EOF on EOF.  If the molfile_timestep_t or molfile_qm_metadata_t
764    * arguments are NULL, then the coordinate or qm data should be skipped.  
765    * Otherwise, the application must prepare molfile_timestep_t and 
766    * molfile_qm_timestep_t by allocating space for the corresponding 
767    * number of coordinates, orbital wavefunction coefficients, etc.
768    * Since it is common for users to want to load only the final timestep 
769    * data from a QM run, the application may provide any combination of
770    * valid, or NULL pointers for the molfile_timestep_t and 
771    * molfile_qm_timestep_t parameters, depending on what information the
772    * user is interested in.
773    * The natoms and qm metadata parameters exist because some file formats 
774    * cannot determine for themselves how many atoms etc are in a 
775    * timestep; the app must therefore obtain this information elsewhere
776    * and provide it to the plugin.
777    */
778   int (* read_timestep)(void *, int natoms, molfile_timestep_t *,
779                         molfile_qm_metadata_t *, molfile_qm_timestep_t *);
780 #endif
781
782 #if vmdplugin_ABIVERSION > 10
783   int (* read_timestep_metadata)(void *, molfile_timestep_metadata_t *);
784 #endif
785 #if vmdplugin_ABIVERSION > 11
786   int (* read_qm_timestep_metadata)(void *, molfile_qm_timestep_metadata_t *);
787 #endif
788
789 #if vmdplugin_ABIVERSION > 13
790   /**
791    *  Console output, READ-ONLY function pointer.
792    *  Function pointer that plugins can use for printing to the host
793    *  application's text console.  This provides a clean way for plugins
794    *  to send message strings back to the calling application, giving the
795    *  caller the ability to prioritize, buffer, and redirect console messages
796    *  to an appropriate output channel, window, etc.  This enables the use of
797    *  graphical consoles like TkCon without losing console output from plugins.
798    *  If the function pointer is NULL, no console output service is provided
799    *  by the calling application, and the output should default to stdout
800    *  stream.  If the function pointer is non-NULL, all output will be
801    *  subsequently dealt with by the calling application.
802    *
803    *  XXX this should really be put into a separate block of
804    *      application-provided read-only function pointers for any
805    *      application-provided services
806    */
807   int (* cons_fputs)(const int, const char*);
808 #endif
809
810 } molfile_plugin_t;
811
812 #endif
813