Main Page   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages  

mol2plugin.C

Go to the documentation of this file.
00001 /***************************************************************************
00002  *cr
00003  *cr            (C) Copyright 1995-2016 The Board of Trustees of the
00004  *cr                        University of Illinois
00005  *cr                         All Rights Reserved
00006  *cr
00007  ***************************************************************************/
00008 
00009 /***************************************************************************
00010  * RCS INFORMATION:
00011  *
00012  *      $RCSfile: mol2plugin.C,v $
00013  *      $Author: johns $       $Locker:  $             $State: Exp $
00014  *      $Revision: 1.38 $       $Date: 2016/11/28 05:01:54 $
00015  *
00016  ***************************************************************************/
00017   
00018 /*
00019  * mol2 file reader
00020  * More information on this format is available at
00021  *   http://www.tripos.com/data/support/mol2.pdf
00022  *   http://www.tripos.com/mol2/
00023  *
00024  *   DOCK mol2 page: 
00025  *     http://www.csb.yale.edu/userguides/datamanip/dock/DOCK_4.0.1/html/Manual.41.html
00026  *
00027  * This plugin currently reads the following record types:
00028  *  MOLECULE
00029  *  ATOM
00030  *  BOND
00031  *
00032  */
00033 
00034 #include "molfile_plugin.h"
00035 
00036 #include <stdlib.h>
00037 #include <stdio.h>
00038 #include <string.h>
00039 
00040 #if defined(_AIX)
00041 #include <strings.h>
00042 #endif
00043 
00044 #define LINESIZE 256
00045 
00046 typedef struct {
00047   FILE *file;
00048   molfile_atom_t *atomlist;
00049   int natoms, nbonds, optflags, coords_read;
00050   int *from, *to;
00051   float *bondorder;
00052 } mol2data;
00053 
00054 // Open the file and create the mol2 struct used to pass data to the other
00055 // functions.
00056 static void *open_mol2_read(const char *path, const char *filetype, 
00057     int *natoms) {
00058   FILE *fd;
00059   mol2data *mol2;
00060   char line[LINESIZE]; 
00061   int match, nbonds, optflags;
00062 
00063   fd = fopen(path, "r");
00064   if (!fd)
00065     return NULL;
00066   
00067   // Find and read the MOLECULE record
00068   do {
00069     fgets(line, LINESIZE, fd);
00070     if ( ferror(fd) || feof(fd) ) {
00071       fprintf(stderr, "mol2plugin) No molecule record found in file.\n");
00072       return NULL;
00073     }
00074   } while ( strncmp(line, "@<TRIPOS>MOLECULE", 17) );
00075 
00076   fgets(line, LINESIZE, fd);  // Read and ignore the mol_name
00077   fgets(line, LINESIZE, fd);  // Read the molecule info
00078   match = sscanf(line, " %d %d", natoms, &nbonds);
00079   if (match == 1) {
00080     nbonds = 0;
00081   }
00082   else if (match != 2) {
00083     fprintf(stderr, "mol2plugin) Cannot determine the number of atoms.\n");
00084     return NULL;
00085   }
00086   fgets(line, LINESIZE, fd);  // Read and ignore the mol_type
00087   fgets(line, LINESIZE, fd);  // Read the charge_type
00088   if ( strncmp(line, "NO_CHARGES", 10) == 0 ) {
00089     optflags = MOLFILE_NOOPTIONS;
00090   }
00091   else {
00092     optflags = MOLFILE_CHARGE;
00093   }
00094 
00095   // Allocate and initialize the mol2 structure
00096   mol2 =  (mol2data *) malloc(sizeof(mol2data));
00097   memset(mol2, 0, sizeof(mol2data));
00098   mol2->file = fd;
00099   mol2->natoms = *natoms;
00100   mol2->nbonds = nbonds;
00101   mol2->optflags = optflags;
00102   mol2->coords_read = 0;
00103   mol2->from = NULL;
00104   mol2->to = NULL;
00105   mol2->bondorder = NULL;
00106 
00107   return mol2;
00108 }
00109 
00110 // Read atom information, but not coordinates.
00111 static int read_mol2(void *v, int *optflags, molfile_atom_t *atoms) {
00112   mol2data *mol2 = (mol2data *)v;
00113   char line[LINESIZE]; 
00114   int i, match;
00115   molfile_atom_t *atom;
00116 
00117   *optflags = mol2->optflags;
00118 
00119   // Find and read the ATOM record
00120   rewind(mol2->file);
00121   do {
00122     fgets(line, LINESIZE, mol2->file);
00123     if ( ferror(mol2->file) || feof(mol2->file) ) {
00124       fprintf(stderr, "mol2plugin) No atom record found in file.\n");
00125       return MOLFILE_ERROR;
00126     }
00127   } while ( strncmp(line, "@<TRIPOS>ATOM", 13) );
00128 
00129   // Read the atoms
00130   for (i = 0; i < mol2->natoms; i++) {
00131     atom = atoms+i;
00132 
00133     fgets(line, LINESIZE, mol2->file);
00134     if ( ferror(mol2->file) || feof(mol2->file) ) {
00135       fprintf(stderr, "mol2plugin) Error occurred reading atom record.\n");
00136       return MOLFILE_ERROR;
00137     }
00138 
00139     match = sscanf(line, " %*d %s %*f %*f %*f %s %d %s %f", 
00140       atom->name, atom->type, &atom->resid, atom->resname, &atom->charge);
00141 
00142     // The last three records are optional for mol2 files, supply values if
00143     // any are missing. Note that these cases are meant to fall through.
00144     switch (match) {
00145       case 0: 
00146         fprintf(stderr, "mol2plugin) Improperly formatted atom record.\n");
00147         return MOLFILE_ERROR;
00148 
00149       case 1:
00150         atom->resid = 0;
00151 
00152       case 2:
00153         sprintf(atom->resname, "%d", atom->resid);
00154 
00155       case 3:
00156         atom->charge = 0.0;
00157 
00158       default:
00159         break;
00160     }
00161 
00162     // Leave these blank when not provided by the file.
00163     atom->chain[0] = '\0';
00164     atom->segid[0] = '\0';
00165   }
00166 
00167   rewind(mol2->file);
00168 
00169   return MOLFILE_SUCCESS;
00170 }
00171 
00172 
00173 
00174 // Create arrays of one-based bond indicies.
00175 static int read_mol2_bonds_aux(void *v, int *nbonds, int **fromptr, int **toptr, 
00176                                float **bondorderptr) {
00177   mol2data *mol2 = (mol2data *)v;
00178   char line[LINESIZE], bond_type[16]; 
00179   int i, match, bond_from, bond_to, bond_index, current_nbonds;
00180   float curr_order;
00181 
00182   if (mol2->nbonds == 0) {
00183     *nbonds = 0;
00184     *fromptr = NULL;
00185     *toptr = NULL;
00186     return MOLFILE_SUCCESS;
00187   }
00188 
00189   current_nbonds = mol2->nbonds;
00190 
00191   // Find and read the BOND record
00192   rewind(mol2->file);
00193   do {
00194     fgets(line, LINESIZE, mol2->file);
00195     if ( ferror(mol2->file) || feof(mol2->file) ) {
00196       fprintf(stderr, "mol2plugin) No bond record found in file.\n");
00197       return MOLFILE_ERROR;
00198     }
00199   } while ( strncmp(line, "@<TRIPOS>BOND", 13) );
00200 
00201   // Read the bonds
00202   bond_index = 0;
00203   for (i = 0; i < mol2->nbonds; i++) {
00204     fgets(line, LINESIZE, mol2->file);
00205     if ( ferror(mol2->file) || feof(mol2->file) ) {
00206       fprintf(stderr, "mol2plugin) Error occurred reading bond record.\n");
00207       return MOLFILE_ERROR;
00208     }
00209 
00210     // Move on if the next line is a header
00211     if (strncmp(line, "@", 1) == 0) {
00212       //Then the bonds are over
00213       break;
00214     }
00215 
00216     match = sscanf(line, " %*d %d %d %s", &bond_from, &bond_to, bond_type);
00217     if (match < 3) {
00218       fprintf(stderr, "mol2plugin) Improperly formatted bond record.\n");
00219       continue;
00220     }
00221     if ( strncmp(bond_type, "nc", 2) == 0 ) {
00222       // Not an actual bond, don't add it to the list
00223       current_nbonds--;
00224     } else if ( strncmp(bond_type, "ar", 2) == 0) {
00225        // Aromatic; consider it order 1.5
00226        curr_order = 1.5;
00227        mol2->from[bond_index] = bond_from;
00228        mol2->to[bond_index] = bond_to;
00229        mol2->bondorder[bond_index]=curr_order;
00230        bond_index++;
00231      } else {
00232       // Add the bond to the list
00233       curr_order=strtod(bond_type,NULL);
00234       if (curr_order<1.0 || curr_order>4.0) curr_order=1;
00235       //fprintf(stdout,"mol2plugin: Bond from %d to %d of order %f\n", bond_from, bond_to, curr_order);
00236       fflush(stdout);
00237       mol2->from[bond_index] = bond_from;
00238       mol2->to[bond_index] = bond_to;
00239       mol2->bondorder[bond_index]=curr_order;
00240       bond_index++;
00241     }
00242   }
00243   if (bond_index > 0) {
00244     *nbonds = current_nbonds;
00245     *fromptr = mol2->from;
00246     *toptr = mol2->to;
00247     *bondorderptr = mol2->bondorder; 
00248   } else {
00249     printf("mol2plugin) WARNING: no bonds defined in mol2 file\n");
00250     *nbonds = 0;
00251     *fromptr = NULL;
00252     *toptr = NULL;
00253     *bondorderptr = NULL; 
00254   }
00255     
00256   rewind(mol2->file);
00257   return MOLFILE_SUCCESS;
00258 }
00259 
00260 
00261 // Read atom coordinates
00262 static int read_mol2_timestep(void *v, int natoms, molfile_timestep_t *ts) {
00263   mol2data *mol2 = (mol2data *)v;
00264   char line[LINESIZE];
00265   int i, match;
00266   float x, y, z;
00267 
00268   // Find and read the ATOM record
00269   do {
00270     fgets(line, LINESIZE, mol2->file);
00271     if ( ferror(mol2->file) || feof(mol2->file) ) {
00272 
00273       //This is a problem, unless we've read at least one timestep
00274       if (mol2->coords_read == 0) {
00275         fprintf(stderr, "mol2plugin) No atom record found in file.\n");
00276       }
00277 
00278       return MOLFILE_ERROR;
00279     }
00280   } while ( strncmp(line, "@<TRIPOS>ATOM", 13) );
00281 
00282   // Read the atoms
00283   for (i = 0; i < mol2->natoms; i++) {
00284     fgets(line, LINESIZE, mol2->file);
00285     if ( ferror(mol2->file) || feof(mol2->file) ) {
00286       fprintf(stderr, "mol2plugin) Error occurred reading atom coordinates.\n");
00287       return MOLFILE_ERROR;
00288     }
00289 
00290 
00291     match = sscanf(line, " %*d %*s %f %f %f", &x, &y, &z);
00292     if (match < 3) {
00293       fprintf(stderr, "mol2plugin) Improperly formatted atom coordinates.\n");
00294       return MOLFILE_ERROR;
00295     }
00296 
00297     if (ts) {
00298       ts->coords[3*i  ] = x;
00299       ts->coords[3*i+1] = y;
00300       ts->coords[3*i+2] = z;
00301     }
00302   }
00303 
00304   mol2->coords_read = 1;
00305   return MOLFILE_SUCCESS;
00306 }
00307 
00308 
00309 static void *open_mol2_write(const char *filename, const char *filetype, 
00310                            int natoms) {
00311   FILE *fd;
00312   mol2data *data;
00313 
00314   fd = fopen(filename, "w");
00315   if (!fd) { 
00316     fprintf(stderr, "mol2plugin) Error: unable to open mol2 file %s for writing\n",
00317             filename);
00318     return NULL;
00319   }
00320   
00321   data = (mol2data *) malloc(sizeof(mol2data));
00322   memset(data, 0, sizeof(mol2data));
00323   data->natoms = natoms;
00324   data->file = fd;
00325   data->nbonds = 0;
00326 //  data->file_name = strdup(filename);
00327   return data;
00328 }
00329 
00330 
00331 static int write_mol2_structure(void *mydata, int optflags, 
00332                                const molfile_atom_t *atoms) {
00333   mol2data *data = (mol2data *)mydata;
00334   data->atomlist = (molfile_atom_t *)malloc(data->natoms*sizeof(molfile_atom_t));
00335   memcpy(data->atomlist, atoms, data->natoms*sizeof(molfile_atom_t));
00336   return MOLFILE_SUCCESS;
00337 }
00338 
00339 /*
00340 void getmol2ff(char* outputtype, const char* psftype) {
00341   //fprintf(stdout,"Doing ff typing on %s\n",psftype);
00342   if (strncmp(psftype,"H",1)==0) {
00343     //It's a hydrogen
00344     strncpy(outputtype, "H   ",4);
00345     return;
00346   } else if (strncmp(psftype,"C",1)==0) {
00347     //It's a carbon... probably
00348     if (strncmp(psftype,"C ",2)==0  || strncmp(psftype,"CA ",3)==0 || 
00349         strncmp(psftype,"CPH",3)==0 || strncmp(psftype,"CPT",3)==0 || 
00350         strncmp(psftype,"CC ",3)==0 || strncmp(psftype,"CD ",3)==0 || 
00351         strncmp(psftype,"CN1",3)==0 || strncmp(psftype,"CN2",3)==0 || 
00352         strncmp(psftype,"CN3",3)==0 || strncmp(psftype,"CN4",3)==0 || 
00353         strncmp(psftype,"CN5",3)==0 || strncmp(psftype,"CNA",3)==0) {
00354       strncpy(outputtype, "C.2 ",4);
00355       return;
00356     } else {
00357       strncpy(outputtype, "C.3 ",4);
00358       return;
00359     }  
00360   } else if (strncmp(psftype,"N",1)==0) {
00361      //It"s probably nitrogen
00362      if (strncmp(psftype,"NR",2)==0 || strncmp(psftype,"NH1",3)==0 || 
00363          strncmp(psftype,"NH2",3)==0 || strncmp(psftype,"NC2",3)==0 || 
00364          strncmp(psftype,"NY",2)==0 || 
00365          (strncmp(psftype,"NN",2)==0 && strncmp(psftype,"NN6",3)!=0)) {
00366        strncpy(outputtype, "N.am",4);
00367        return;
00368        } else {
00369        strncpy(outputtype, "N.3 ",4);
00370        return;
00371        }
00372   } else if (strncmp(psftype,"O",1)==0) {
00373      //Probably an oxygen
00374      if (strncmp(psftype,"OH1",3)==0 || strncmp(psftype,"OS",2)==0 || 
00375          strncmp(psftype,"OT ",3)==0 || strncmp(psftype,"ON4",3)==0 || 
00376          strncmp(psftype,"ON5",3)==0 || strncmp(psftype,"ON6",3)==0) {
00377        strncpy(outputtype, "O.3 ",4);
00378        return;
00379      } else {
00380        strncpy(outputtype, "O.2 ",4);
00381        return;
00382      } 
00383   } else if (strncmp(psftype,"S",1)==0) {
00384      strncpy(outputtype, "S.3 ",4);
00385      return;
00386   } else if (strncmp(psftype,"P",1)==0) {
00387      strncpy(outputtype, "P.3 ",4);
00388      return;
00389   } else {
00390      strncpy(outputtype, "X.  ",4);
00391      return;
00392   }
00393 }
00394 */
00395 
00396 
00397 
00398 
00399 
00400 
00401 static int write_mol2_timestep(void *mydata, const molfile_timestep_t *ts) {
00402   mol2data *data = (mol2data *)mydata; 
00403   const molfile_atom_t *atom;
00404   const float *pos;
00405   float chrgsq;
00406   int i;
00407 
00408   // try to guess whether we have charge information.
00409   chrgsq=0.0;
00410   atom = data->atomlist;
00411   for (i = 0; i < data->natoms; i++) {
00412       chrgsq += atom->charge*atom->charge;
00413       ++atom;
00414   }
00415 
00416   //print header block
00417   fprintf(data->file, "@<TRIPOS>MOLECULE\n");
00418   fprintf(data->file, "generated by VMD\n");
00419   fprintf(data->file, " %4d %4d    1    0    0\n", data->natoms, data->nbonds);
00420   fprintf(data->file, "SMALL\n");
00421   // educated guess
00422   if (chrgsq > 0.0001) {
00423       fprintf(data->file, "USER_CHARGES\n");
00424   } else {
00425       fprintf(data->file, "NO_CHARGES\n");
00426   }
00427   fprintf(data->file, "****\n");
00428   fprintf(data->file, "Energy = 0\n\n");
00429   
00430   //print atoms block
00431   fprintf(data->file, "@<TRIPOS>ATOM\n");
00432   atom = data->atomlist;
00433   pos = ts->coords;
00434   //char mol2fftype[5];
00435   //mol2fftype[4] = '\0';
00436   for (i = 0; i < data->natoms; i++) {
00437 #if 0
00438     getmol2ff(mol2fftype, atom->type);
00439     fprintf(data->file, "%7d %-4s      %8.4f  %8.4f  %8.4f %4s %4d  %3s        %8.6f\n", i+1, atom->name, pos[0], pos[1], pos[2], mol2fftype, atom->resid, atom->resname, atom->charge);
00440 #else
00441     fprintf(data->file, "%7d %-4s      %8.4f  %8.4f  %8.4f %4s %4d  %3s        %8.6f\n", i+1, atom->name, pos[0], pos[1], pos[2], atom->type, atom->resid, atom->resname, atom->charge);
00442 #endif
00443     ++atom; 
00444     pos += 3;
00445   }
00446 
00447   //print bond info
00448 
00449   int l=1; //number of bond record
00450   printf("mol2plugin) numbonds: %d\n", data->nbonds);
00451   if (data->nbonds>0) fprintf(data->file, "@<TRIPOS>BOND\n");
00452   for (i=0; i<data->nbonds; i++) {
00453     // For mol2, only write bonds for fromptr[i]<toptr[i]
00454     // bondorder is either 1, 2, 3 or a textual representation: am,ar,du,un,nc
00455     // we don't have the info for the text, so we truncate to integer.
00456     if (data->bondorder != NULL) {
00457       fprintf(data->file, "%5d %5d %5d %2d\n", l ,data->from[i], data->to[i],
00458               (int)data->bondorder[i]);
00459     } else {
00460       fprintf(data->file, "%5d %5d %5d %2d\n", l ,data->from[i], data->to[i], 1);
00461     }
00462 
00463     l++;
00464   } 
00465 
00466   // Print out substructure info to keep some programs sane
00467   fprintf(data->file,"\n@<TRIPOS>SUBSTRUCTURE\n");
00468   fprintf(data->file,"1 ****        1 TEMP                        ");
00469   fprintf(data->file,"0 ****  **** 0 ROOT\n");
00470 
00471   return MOLFILE_SUCCESS;
00472 }
00473 
00474 static int write_mol2_bonds(void *v, int nbonds, int *fromptr, int *toptr, 
00475                             float *bondorderptr,  int *bondtype, 
00476                             int nbondtypes, char **bondtypename) {
00477   printf("*** RUNNING WRITE_MOL2_BONDS\n");
00478   mol2data *data = (mol2data *)v;
00479   data->nbonds = nbonds;
00480   data->from = (int *) malloc(nbonds * sizeof(int));
00481   data->to = (int *) malloc(nbonds * sizeof(int));
00482   //set the pointers for use later
00483   for (int i=0;i<nbonds;i++) {
00484     data->from[i]=fromptr[i];
00485     data->to[i]=toptr[i];
00486   }
00487   printf("*** I THINK nbonds is %i\n", nbonds);
00488   data->nbonds = nbonds;
00489   if (bondorderptr != NULL) {
00490     data->bondorder = (float *) malloc(nbonds * sizeof(float));
00491     for (int i=0;i<nbonds;i++) {
00492       data->bondorder[i]=bondorderptr[i];
00493     }
00494   }
00495 
00496   return MOLFILE_SUCCESS;
00497 }
00498 
00499 
00500 static void close_mol2_write(void *mydata) {
00501   mol2data *data = (mol2data *)mydata;
00502   if (data) {
00503     if (data->file) fclose(data->file);
00504     if (data->from != NULL) free(data->from);
00505     data->from = NULL;
00506     if (data->to != NULL)   free(data->to);
00507     data->to = NULL;
00508     if (data->bondorder != NULL)   free(data->bondorder);
00509     data->bondorder = NULL;
00510     if (data->atomlist != NULL) free(data->atomlist);
00511     data->atomlist = NULL;
00512     free(data);
00513   }
00514 }
00515 
00516 //
00517 // Free the memory used by the mol2 structure
00518 static void close_mol2_read(void *v) {
00519   mol2data *mol2 = (mol2data *)v;
00520   if (mol2) {
00521     if (mol2->file) fclose(mol2->file);
00522     if (mol2->from != NULL) free(mol2->from);
00523     mol2->from = NULL;
00524     if (mol2->to != NULL)   free(mol2->to);
00525     mol2->to = NULL;
00526     if (mol2->bondorder != NULL)   free(mol2->bondorder);
00527     mol2->bondorder = NULL;
00528     free(mol2);
00529   }
00530 }
00531 
00532 
00533 static int read_mol2_bonds(void *v, int *nbonds, int **fromptr, int **toptr, 
00534                            float **bondorderptr, int **bondtype, 
00535                       int *nbondtypes, char ***bondtypename) {
00536   mol2data *mol2 = (mol2data *)v;
00537 
00538   /* now read bond data */
00539 //  *nbonds = start_psf_bonds(psf->fp);
00540 
00541   if (mol2->nbonds > 0) {
00542     mol2->from = (int *) malloc(mol2->nbonds*sizeof(int));
00543     mol2->to = (int *) malloc(mol2->nbonds*sizeof(int));
00544     mol2->bondorder = (float *) malloc(mol2->nbonds*sizeof(float));
00545     if (mol2->from == NULL || mol2->to == NULL || mol2->bondorder == NULL) {
00546       fprintf(stderr, "mol2plugin) ERROR: Failed to allocate memory for bonds\n");
00547       fclose(mol2->file);
00548       mol2->file = NULL;
00549       return MOLFILE_ERROR;
00550     }
00551     if ((read_mol2_bonds_aux(mol2, nbonds, &(mol2->from), &(mol2->to), &(mol2->bondorder))) != MOLFILE_SUCCESS) {
00552       fclose(mol2->file);
00553       mol2->file = NULL;
00554       return MOLFILE_ERROR;
00555     }
00556     *fromptr = mol2->from;
00557     *toptr = mol2->to;
00558     *bondorderptr = mol2->bondorder; 
00559     *bondtype = NULL;
00560     *nbondtypes = 0;
00561     *bondtypename = NULL;
00562   } else {
00563     printf("mol2plugin) WARNING: zero bonds defined in mol2 file.\n");
00564     *nbonds = 0;
00565     *fromptr=NULL;
00566     *toptr=NULL;
00567     *bondorderptr=NULL;
00568     *bondtype = NULL;
00569     *nbondtypes = 0;
00570     *bondtypename = NULL;
00571   }
00572   return MOLFILE_SUCCESS;
00573 }
00574 
00575 static molfile_plugin_t plugin;
00576 
00577 VMDPLUGIN_API int VMDPLUGIN_init() {
00578   memset(&plugin, 0, sizeof(molfile_plugin_t));
00579   plugin.abiversion = vmdplugin_ABIVERSION;
00580   plugin.type = MOLFILE_PLUGIN_TYPE;
00581   plugin.name = "mol2";
00582   plugin.prettyname = "MDL mol2";
00583   plugin.author = "Peter Freddolino, Eamon Caddigan";
00584   plugin.majorv = 0;
00585   plugin.minorv = 17;
00586   plugin.is_reentrant = VMDPLUGIN_THREADSAFE;
00587   plugin.filename_extension = "mol2";
00588   plugin.open_file_read = open_mol2_read;
00589   plugin.read_structure = read_mol2;
00590   plugin.read_bonds = read_mol2_bonds;
00591   plugin.read_next_timestep = read_mol2_timestep;
00592   plugin.close_file_read = close_mol2_read;
00593   plugin.open_file_write = open_mol2_write;
00594   plugin.write_structure = write_mol2_structure;
00595   plugin.write_timestep = write_mol2_timestep;
00596   plugin.close_file_write = close_mol2_write;
00597   plugin.write_bonds = write_mol2_bonds;
00598   return VMDPLUGIN_SUCCESS;
00599 }
00600 
00601 VMDPLUGIN_API int VMDPLUGIN_register(void *v, vmdplugin_register_cb cb) {
00602   (*cb)(v, (vmdplugin_t *)&plugin);
00603   return VMDPLUGIN_SUCCESS;
00604 }
00605 
00606 VMDPLUGIN_API int VMDPLUGIN_fini() {
00607   return VMDPLUGIN_SUCCESS;
00608 }
00609 

Generated on Fri Mar 29 03:10:12 2024 for VMD Plugins (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002