Main Page   Alphabetical List   Compound List   File List   Compound Members   File Members   Related Pages  

cubeplugin.C

Go to the documentation of this file.
00001 /***************************************************************************
00002  *cr
00003  *cr            (C) Copyright 1995-2006 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: cubeplugin.C,v $
00013  *      $Author: johns $       $Locker:  $             $State: Exp $
00014  *      $Revision: 1.23 $       $Date: 2006/01/05 00:05:52 $
00015  *
00016  ***************************************************************************/
00017 
00018 //
00019 // Plugin reader for Gaussian "cube" files
00020 //
00021 // Gaussian "cube" file format described here:
00022 //   http://www.gaussian.com/00000430.htm
00023 //
00024 
00025 #include <stdio.h>
00026 #include <stdlib.h>
00027 #include <math.h>
00028 #include <string.h>
00029 
00030 #include "molfile_plugin.h"
00031 
00032 #ifndef M_PI_2
00033 #define M_PI_2 1.57079632679489661922
00034 #endif
00035 
00036 #include "periodic_table.h"
00037 
00038 static const float bohr = 0.529177249;
00039 
00040 // A format-independent structure to hold unit cell data
00041 typedef struct {
00042   float A, B, C, alpha, beta, gamma;
00043 } cube_box;
00044 
00045 typedef struct {
00046   FILE *fd;              // file descriptor
00047   int nsets;             // number of volume datasets
00048   int numatoms;          // number of atoms
00049   bool coord;            // has coordinate data
00050   long crdpos, datapos;  // seek offsets for coords and data
00051   char *file_name;       // original filename 
00052   float *datacache;      // temporary cache of orbital data prior to conversion
00053   molfile_volumetric_t *vol; // volume data
00054   float origin[3];       // origin, stored for periodic display hack 
00055   float rotmat[3][3];    // rotation matrix, stored for periodic display hack
00056   cube_box box;          // unit cell dimensions
00057 } cube_t;
00058 
00059 // Converts box basis vectors to A, B, C, alpha, beta, and gamma.  
00060 // Stores values in cube_box struct, which should be allocated before calling
00061 // this function.
00062 static int cube_readbox(cube_box *box, float *x, float *y, float *z) {
00063   float A, B, C;
00064 
00065   if (!box) {
00066     return 1;
00067   }
00068 
00069   // provide defaults
00070   box->A = 10.0;
00071   box->B = 10.0;
00072   box->C = 10.0;
00073   box->alpha = 90.0;
00074   box->beta  = 90.0;
00075   box->gamma = 90.0;
00076 
00077   // A, B, C are the lengths of the x, y, z vectors, respectively
00078   A = sqrt( x[0]*x[0] + x[1]*x[1] + x[2]*x[2] );
00079   B = sqrt( y[0]*y[0] + y[1]*y[1] + y[2]*y[2] );
00080   C = sqrt( z[0]*z[0] + z[1]*z[1] + z[2]*z[2] );
00081   if ((A<=0) || (B<=0) || (C<=0)) {
00082     return 1;
00083   }
00084   box->A = A;
00085   box->B = B;
00086   box->C = C;
00087 
00088   // gamma, beta, alpha are the angles between the x & y, x & z, y & z
00089   // vectors, respectively
00090   box->gamma = acos( (x[0]*y[0]+x[1]*y[1]+x[2]*y[2])/(A*B) ) * 90.0/M_PI_2;
00091   box->beta = acos( (x[0]*z[0]+x[1]*z[1]+x[2]*z[2])/(A*C) ) * 90.0/M_PI_2;
00092   box->alpha = acos( (y[0]*z[0]+y[1]*z[1]+y[2]*z[2])/(B*C) ) * 90.0/M_PI_2; 
00093 
00094   return 0;
00095 }
00096 
00097 // calculate and store origin and rotation matrix to realign everything later.
00098 static void cube_buildrotmat(cube_t *cube, float *o, float *a, float *b)
00099 {
00100   // we rotate first around y and z to align a along the x-axis...
00101   const double len   = sqrt(a[0]*a[0] + a[1]*a[1]);
00102   const double phi   = atan2((double) a[2], (double) len);
00103   const double theta = atan2((double) a[1], (double) a[0]);
00104 
00105   const double cph = cos(phi);
00106   const double cth = cos(theta);
00107   const double sph = sin(phi);
00108   const double sth = sin(theta);
00109 
00110   // ...then we rotate around x to put b into the xy-plane.
00111   const double psi = atan2(-sph*cth*b[0] - sph*sth*b[1] + cph*b[2],-sth*b[0] + cth*b[1]);
00112   const double cps = cos(psi);
00113   const double sps = sin(psi);
00114 
00115   const double r[3][3] = { 
00116     {               cph*cth,                    cph*sth,      sph},
00117     {-sth*cps - sph*cth*sps,      cth*cps - sph*sth*sps,  cph*sps},
00118     { sth*sps - sph*cth*cps,     -cth*sps - sph*sth*cps,  cph*cps}
00119   };
00120 
00121   for (int i=0; i<3; ++i) {
00122     cube->origin[i] = o[i];
00123     for (int j=0; j<3; ++j) {
00124       cube->rotmat[i][j] = r[i][j];
00125     }
00126   }
00127 }
00128 
00129 
00130 static void eatline(FILE * fd) {
00131   char readbuf[1025];
00132   fgets(readbuf, 1024, fd);    // go on to next line
00133 }  
00134 
00135 // prototype.
00136 static void close_cube_read(void *v);
00137 
00138 static void *open_cube_read(const char *filepath, const char *filetype,
00139     int *natoms) {
00140   FILE *fd;
00141   cube_t *cube;
00142   int xsize, ysize, zsize;
00143   float a[3], b[3], c[3];
00144   int i;
00145  
00146   fd = fopen(filepath, "rb");
00147   if (!fd) 
00148     return NULL;
00149 
00150   cube = new cube_t;
00151   cube->fd = fd;
00152   cube->vol = NULL;
00153   cube->coord = false;
00154   cube->file_name = strdup(filepath);
00155   cube->datacache = NULL;
00156 
00157   // initialize origin and rotmat to sensible defaults.
00158   for (i=0; i<3; ++i) {
00159       for (int j=0; j<3; ++j) {
00160           cube->rotmat[i][j] = 0.0;
00161       }
00162   }
00163   for (i=0; i<3; ++i) {
00164       cube->origin[i] = 0.0;
00165       cube->rotmat[i][i] = 1.0;
00166   }
00167 
00168   molfile_volumetric_t voltmpl; // base information for all data sets.
00169 
00170   // read in cube file header information
00171   char readbuf[256]; 
00172   fgets(readbuf, 256, cube->fd);    // go on to next line
00173 
00174   // identify this file, and read title string into dataset info
00175   strcpy(voltmpl.dataname, "Gaussian Cube: ");
00176   strncat(voltmpl.dataname, readbuf, 240);      // 240 is max space left after
00177                                                 //   "Gaussian Cube: "
00178   eatline(cube->fd);          // skip second header line 
00179 
00180   // read in number of atoms
00181   if (fscanf(cube->fd, "%d", &cube->numatoms) != 1) {
00182     close_cube_read(cube);
00183     return NULL;
00184   }
00185 
00186   if (cube->numatoms > 0) {   // density cube file
00187     cube->nsets = 1;          // this cube file contains only one data set
00188   } else {
00189     // cube file with orbitals => multiple densities.
00190     cube->numatoms = - cube->numatoms;
00191     cube->nsets = 0;          // we don't know yet how many data sets.
00192   }
00193   *natoms = cube->numatoms;
00194 
00195   // read in cube origin
00196   if ( fscanf(cube->fd, "%f %f %f", 
00197          &voltmpl.origin[0], 
00198          &voltmpl.origin[1], 
00199          &voltmpl.origin[2]) != 3 ) {
00200     close_cube_read(cube);
00201     return NULL;
00202   }
00203 
00204   // read in cube axes and sizes
00205   if ((fscanf(cube->fd, "%d", &xsize) != 1) ||
00206       (fscanf(cube->fd, "%f %f %f", &a[0], &a[1], &a[2]) != 3)) {
00207     close_cube_read(cube);
00208     return NULL;
00209   }
00210 
00211   if ((fscanf(cube->fd, "%d", &ysize) != 1) ||
00212       (fscanf(cube->fd, "%f %f %f", &b[0], &b[1], &b[2]) != 3)) {
00213     close_cube_read(cube);
00214     return NULL;
00215   }
00216 
00217   if ((fscanf(cube->fd, "%d", &zsize) != 1) ||
00218       (fscanf(cube->fd, "%f %f %f", &c[0], &c[1], &c[2]) != 3)) {
00219     close_cube_read(cube);
00220     return NULL;
00221   }
00222 
00223   // calculate number of samples in each dimension
00224   voltmpl.xsize = xsize;
00225   voltmpl.ysize = ysize;
00226   voltmpl.zsize = zsize;
00227   voltmpl.has_color = 0;
00228 
00229   eatline(cube->fd);     // skip remaining end of line characters
00230 
00231   // to make the periodic display work, we need to rotate
00232   // the cellvectors (and the coordinates) in such a way,
00233   // that the a-vector is collinear with the x-axis and
00234   // the b-vector is in the xy-plane. 
00235   cube_buildrotmat(cube, voltmpl.origin, a, b);
00236   // print warning, if the rotation will be significant:
00237   if ((fabs((double) a[1]) + fabs((double) a[2]) + fabs((double) b[2]))
00238       > 0.001) {
00239     fprintf(stderr, "cubeplugin) WARNING: Coordinates will be rotated to comply \n"
00240             "cubeplugin) with VMD's conventions for periodic display...\n");
00241   }
00242   
00243 
00244   // all dimensional units are always in Bohr
00245   // so scale axes and origin correctly.
00246   // NOTE: the angstroms are only allowed in input.
00247   voltmpl.origin[0] *= bohr; 
00248   voltmpl.origin[1] *= bohr;
00249   voltmpl.origin[2] *= bohr;
00250 
00251   // store aligned axes.
00252   for (i=0; i<3; ++i) {
00253     voltmpl.xaxis[i] = cube->rotmat[i][0] * a[0] 
00254       + cube->rotmat[i][1] * a[1] + cube->rotmat[i][2] * a[2];
00255 
00256     voltmpl.yaxis[i] = cube->rotmat[i][0] * b[0] 
00257       + cube->rotmat[i][1] * b[1] + cube->rotmat[i][2] * b[2];
00258     
00259     voltmpl.zaxis[i] = cube->rotmat[i][0] * c[0] 
00260       + cube->rotmat[i][1] * c[1] + cube->rotmat[i][2] * c[2];
00261   }
00262 
00263   voltmpl.xaxis[0] *= bohr * xsize;
00264   voltmpl.xaxis[1] *= bohr * xsize; 
00265   voltmpl.xaxis[2] *= bohr * xsize; 
00266 
00267   voltmpl.yaxis[0] *= bohr * ysize; 
00268   voltmpl.yaxis[1] *= bohr * ysize; 
00269   voltmpl.yaxis[2] *= bohr * ysize; 
00270 
00271   voltmpl.zaxis[0] *= bohr * zsize; 
00272   voltmpl.zaxis[1] *= bohr * zsize; 
00273   voltmpl.zaxis[2] *= bohr * zsize; 
00274 
00275 #if 1
00276   /*   As of VMD version 1.8.3, volumetric data points are 
00277    *   expected to represent the center of a grid box. cube format 
00278    *   volumetric data represents the value at the edges of the 
00279    *   grid boxes, so we need to shift the internal origin by half
00280    *   a grid box diagonal to have the data at the correct position.
00281    *   This will need to be changed again when the plugin interface 
00282    *   is updated to explicitly allow point/face-centered data sets.
00283    */
00284   voltmpl.origin[0] -= 0.5 * ( voltmpl.xaxis[0] / (double) xsize
00285                         + voltmpl.yaxis[0] / (double) ysize
00286                         + voltmpl.zaxis[0] / (double) zsize);
00287   voltmpl.origin[1] -= 0.5 * ( voltmpl.xaxis[1] / (double) xsize
00288                         + voltmpl.yaxis[1] / (double) ysize
00289                         + voltmpl.zaxis[1] / (double) zsize);
00290   voltmpl.origin[2] -= 0.5 * ( voltmpl.xaxis[2] / (double) xsize
00291                         + voltmpl.yaxis[2] / (double) ysize
00292                         + voltmpl.zaxis[2] / (double) zsize);
00293 #endif
00294 
00295 #if defined(TEST_PLUGIN)
00296   printf("cell before rotation:\n");
00297   printf("a: %12.8f %12.8f %12.8f\n", a[0]*xsize*bohr, a[1]*ysize*bohr, a[2]*zsize*bohr);
00298   printf("b: %12.8f %12.8f %12.8f\n", b[0]*xsize*bohr, b[1]*ysize*bohr, b[2]*zsize*bohr);
00299   printf("c: %12.8f %12.8f %12.8f\n", c[0]*xsize*bohr, c[1]*ysize*bohr, c[2]*zsize*bohr);
00300 
00301   printf("cell after rotation:\n");
00302   printf("x: %12.8f %12.8f %12.8f\n", voltmpl.xaxis[0], voltmpl.xaxis[1], voltmpl.xaxis[2]);
00303   printf("y: %12.8f %12.8f %12.8f\n", voltmpl.yaxis[0], voltmpl.yaxis[1], voltmpl.yaxis[2]);
00304   printf("z: %12.8f %12.8f %12.8f\n", voltmpl.zaxis[0], voltmpl.zaxis[1], voltmpl.zaxis[2]);
00305 #endif
00306 
00307   // store the unit cell information for later perusal.
00308   if (cube_readbox(&(cube->box), voltmpl.xaxis, voltmpl.yaxis, voltmpl.zaxis)) {
00309           fprintf(stderr, "cubeplugin) Calculation of unit cell size failed. Trying to continue...\n");
00310   }
00311 
00312   cube->crdpos = ftell(cube->fd); // and record beginning of coordinates
00313   // XXX fseek()/ftell() are incompatible with 64-bit LFS I/O implementations, 
00314   // hope we don't read any files >= 2GB...
00315 
00316   if (cube->nsets >0) { 
00317     int i;
00318 
00319     // density cube file, copy voltmpl into the cube struct.
00320     cube->vol = new molfile_volumetric_t[1];
00321     memcpy(cube->vol, &voltmpl, sizeof(voltmpl));
00322 
00323     // skip over coordinates to find the start of volumetric data
00324     for (i=0; i < cube->numatoms; i++) {
00325       eatline(cube->fd);
00326     }
00327 
00328     cube->datapos = ftell(cube->fd); // and record beginning of data
00329     // XXX fseek()/ftell() are incompatible with 64-bit LFS I/O, 
00330     // hope we don't read any files >= 2GB...
00331   } else {              
00332     int i;
00333 
00334     // orbital cube file. we now have to read the orbitals section
00335     // skip over coordinates
00336     for (i=0; i < cube->numatoms; i++) {
00337       eatline(cube->fd);
00338     }
00339       
00340     fscanf(cube->fd, "%d", &cube->nsets);
00341     fprintf(stderr, "\ncubeplugin) found %d orbitals\n", cube->nsets);
00342     cube->vol = new molfile_volumetric_t[cube->nsets];
00343     for (i=0; i < cube->nsets; ++i) {
00344       int orb;
00345       fscanf(cube->fd, "%d", &orb);
00346       memcpy(&cube->vol[i], &voltmpl, sizeof(voltmpl));
00347       sprintf(cube->vol[i].dataname, "Gaussian Cube: Orbital %d", orb);
00348     }
00349       
00350     eatline(cube->fd);        // gobble up rest of line
00351     cube->datapos = ftell(cube->fd); // and record beginning of data
00352     // XXX fseek()/ftell() are incompatible with 64-bit LFS I/O, 
00353     // hope we don't read any files >= 2GB...
00354   }
00355 
00356   return cube;
00357 }
00358 
00359   
00360 static int read_cube_structure(void *v, int *optflags, molfile_atom_t *atoms) {
00361   int i, j;
00362   char *k;
00363   molfile_atom_t *atom;
00364 
00365   cube_t *cube = (cube_t *)v;
00366 
00367   // go to coordinates
00368   fseek(cube->fd, cube->crdpos, SEEK_SET);
00369   // XXX fseek()/ftell() are incompatible with 64-bit LFS I/O implementations, 
00370   // hope we don't read any files >= 2GB...
00371 
00372   /* set mass and radius from PTE. */
00373   *optflags = MOLFILE_ATOMICNUMBER | MOLFILE_MASS | MOLFILE_RADIUS; 
00374 
00375   for(i=0;i<cube->numatoms;i++) {
00376     int idx;
00377     float chrg, coord;
00378     char fbuffer[1024];
00379 
00380     atom = atoms + i;
00381 
00382     k = fgets(fbuffer, 1024, cube->fd);
00383     j=sscanf(fbuffer, "%d %f %f %f %f", &idx, &chrg, &coord, &coord, &coord);
00384     if (k == NULL) {
00385       fprintf(stderr, "cube structure) missing atom(s) in "
00386               "file '%s'\n",cube->file_name);
00387       fprintf(stderr, "cube structure) expecting '%d' atoms,"
00388               " found only '%d'\n",cube->numatoms,i+1);
00389       return MOLFILE_ERROR;
00390     } else if (j < 5) {
00391       fprintf(stderr, "cube structure) missing type or coordinate(s)"
00392               " in file '%s' for atom '%d'\n",cube->file_name,i+1);
00393       return MOLFILE_ERROR;
00394     }
00395 
00396     // assign atom symbol to number. flag unknown or dummy atoms with X.
00397     atom->atomicnumber = idx;
00398     strncpy(atom->name, get_pte_label(idx), sizeof(atom->name));
00399     strncpy(atom->type, atom->name, sizeof(atom->type));
00400     atom->mass = get_pte_mass(idx);
00401     atom->radius = get_pte_vdw_radius(idx);
00402     atom->resname[0] = '\0';
00403     atom->resid = 1;
00404     atom->chain[0] = '\0';
00405     atom->segid[0] = '\0';
00406     /* skip to the end of line */
00407   }
00408 
00409   return MOLFILE_SUCCESS;
00410 }
00411 
00412 
00413 static int read_cube_timestep(void *v, int natoms, molfile_timestep_t *ts) {
00414   int i, j, n;
00415   char fbuffer[1024];
00416   float x, y, z;
00417   char *k;
00418   
00419   cube_t *cube = (cube_t *)v;
00420 
00421   // there is only one set of coordinates
00422   if (cube->coord) return MOLFILE_EOF;
00423   cube->coord = true;
00424 
00425   // jump to coordinate position
00426   fseek(cube->fd, cube->crdpos, SEEK_SET);
00427   // XXX fseek()/ftell() are incompatible with 64-bit LFS I/O implementations, 
00428   // hope we don't read any files >= 2GB...
00429  
00430   /* read the coordinates */
00431   for (i=0; i<cube->numatoms; i++) {
00432     k = fgets(fbuffer, 1024, cube->fd);
00433     j = sscanf(fbuffer, "%*d %*f %f %f %f", &x, &y, &z);
00434     
00435     if (k == NULL) {
00436       return MOLFILE_ERROR;
00437     } else if (j < 3) {
00438       fprintf(stderr, "cube timestep) missing type or coordinate(s) in file '%s' for atom '%d'\n",cube->file_name,i+1);
00439       return MOLFILE_ERROR;
00440     } else if (j>=3) {
00441       if (ts != NULL) { 
00442         // Only save coords if we're given a timestep pointer, 
00443         // otherwise assume that VMD wants us to skip past it.
00444         
00445         // In order to make the periodic display work, we need to
00446         // rotate the coordinates around the origin by the stored
00447         // rotation matrix. All coordinates are in Bohr, so they
00448         // must be converted to Angstrom, too.
00449         x -= cube->origin[0];
00450         y -= cube->origin[1];
00451         z -= cube->origin[2];
00452         
00453         for (n=0; n<3; ++n) {
00454           ts->coords[3*i + n] = bohr*(cube->origin[n] 
00455                                       + cube->rotmat[n][0] * x
00456                                       + cube->rotmat[n][1] * y
00457                                       + cube->rotmat[n][2] * z);
00458         }
00459       }
00460     } else {
00461       break;
00462     }
00463   }
00464   // set unitcell dimensions from cached data.
00465   if (ts != NULL) { 
00466       ts->A = cube->box.A;
00467       ts->B = cube->box.B;
00468       ts->C = cube->box.C;
00469       ts->alpha = cube->box.alpha;
00470       ts->beta  = cube->box.beta;
00471       ts->gamma = cube->box.gamma;
00472   }
00473   
00474   return MOLFILE_SUCCESS;
00475 }
00476 
00477 static int read_cube_metadata(void *v, int *nsets, 
00478   molfile_volumetric_t **metadata) {
00479   cube_t *cube = (cube_t *)v;
00480   *nsets = cube->nsets; 
00481   *metadata = cube->vol;  
00482 
00483   return MOLFILE_SUCCESS;
00484 }
00485 
00486 static int read_cube_data(void *v, int set, float *datablock, float *colorblock) {
00487   cube_t *cube = (cube_t *)v;
00488 
00489   fprintf(stderr, "cubeplugin) trying to read cube data set %d\n", set);
00490 
00491   int xsize = cube->vol[set].xsize; 
00492   int ysize = cube->vol[set].ysize;
00493   int zsize = cube->vol[set].zsize;
00494   int xysize = xsize*ysize;
00495   int nsize = cube->nsets;
00496   int nzsize = nsize*zsize;
00497   int nyzsize = nsize*zsize*ysize;
00498   int x, y, z;
00499 
00500   // go to data
00501   fseek(cube->fd, cube->datapos, SEEK_SET);
00502   // XXX fseek()/ftell() are incompatible with 64-bit LFS I/O implementations, 
00503   // hope we don't read any files >= 2GB...
00504 
00505   // read the data values in 
00506   if (cube->nsets == 1) { // density cube file
00507     for (x=0; x<xsize; x++) {
00508       for (y=0; y<ysize; y++) {
00509         for (z=0; z<zsize; z++) {
00510           if (fscanf(cube->fd, "%f", &datablock[z*xysize + y*xsize + x]) != 1) {
00511               return MOLFILE_ERROR;
00512           }
00513         }
00514       } 
00515     }
00516   } else {
00517     // XXX we may wish to examine this strategy for alternatives that provide
00518     // the same performance but without the extra copy, but it makes sense 
00519     // for now.  
00520 
00521     // Since the orbital cube file stores the data orb1(a1,b1,c1), orb2(a1,b1,c1), 
00522     // ... orbn(a1,b1,c1), orb1(a1,b1,c2), orb2(a1,a1,c2), ..., orbn(ai,bj,ck)
00523     // we have to cache the whole data set of have any kind of reasonable performance.
00524     // otherwise we would have to read (and parse!) the whole file over and over again.
00525     // this way we have to do it only once at the temporary expense of some memory.
00526     if (cube->datacache == NULL) {
00527       int points = xsize*ysize*zsize * nsize; // number of grid cells * nsets
00528       int i;
00529 
00530       // let people know what is going on.
00531       fprintf(stderr, "\ncubeplugin) creating %d MByte cube orbital cache.\n", 
00532               (int) (points*sizeof(float)) / 1048576);
00533       cube->datacache = new float[points];
00534             
00535       for (i=0; i < points; ++i) {
00536         if (fscanf(cube->fd, "%f", &cube->datacache[i]) != 1) {
00537           return MOLFILE_ERROR;
00538         }
00539 
00540         // print an ascii progress bar so impatient people do not get scared.
00541         if ((i % (1048576/sizeof(float))) == 0) {  // one dot per MB.
00542           fprintf(stderr, "."); 
00543         }
00544       }
00545     }
00546       
00547     for (x=0; x<xsize; x++) {
00548       for (y=0; y<ysize; y++) {
00549         for (z=0; z<zsize; z++) {
00550           datablock[z*xysize + y*xsize + x] = cube->datacache[x*nyzsize + y*nzsize + z*nsize + set];
00551         }
00552       }
00553     }
00554   }
00555 
00556   return MOLFILE_SUCCESS;
00557 }
00558 
00559 static void close_cube_read(void *v) {
00560   cube_t *cube = (cube_t *)v;
00561 
00562   fclose(cube->fd);
00563   if (cube->vol) {
00564     delete[] cube->vol; 
00565   }
00566   free(cube->file_name);
00567   if (cube->datacache) { 
00568     fprintf(stderr, "\ncubeplugin) freeing cube orbital cache.\n");
00569     delete[] cube->datacache; 
00570   }  
00571 
00572   delete cube;
00573 }
00574 
00575 /*
00576  * Initialization stuff here
00577  */
00578 static molfile_plugin_t plugin = {
00579   vmdplugin_ABIVERSION,   // ABI version
00580   MOLFILE_PLUGIN_TYPE,    // plugin type
00581   "cube",                 // short file format description
00582   "Gaussian Cube",        // pretty file format description
00583   "Axel Kohlmeyer, John E. Stone", // author(s)
00584   0,                      // major version
00585   8,                      // minor version
00586   VMDPLUGIN_THREADSAFE,   // is reentrant
00587   "cube",                 // filename extension
00588   open_cube_read,               
00589   read_cube_structure,
00590   0,                      // read_bonds
00591   read_cube_timestep,
00592   close_cube_read,
00593   0,                      // open_file_write
00594   0,                      // write_structure
00595   0,                      // write_timestep
00596   0,                      // close_file_write
00597   read_cube_metadata,
00598   read_cube_data,
00599   0                       // read_rawgraphics
00600 };
00601 
00602 int VMDPLUGIN_init(void) { return VMDPLUGIN_SUCCESS; }
00603 int VMDPLUGIN_fini(void) { return VMDPLUGIN_SUCCESS; }
00604 int VMDPLUGIN_register(void *v, vmdplugin_register_cb cb) {
00605   (*cb)(v, (vmdplugin_t *)&plugin);
00606   return VMDPLUGIN_SUCCESS;
00607 }
00608 
00609 
00610 
00611 #ifdef TEST_PLUGIN
00612 
00613 int main(int argc, char *argv[]) {
00614   int natoms;
00615   void *v;
00616   int i, nsets, set;
00617   molfile_volumetric_t * meta;
00618 
00619   while (--argc) {
00620     ++argv;
00621     v = open_cube_read(*argv, "cube", &natoms);
00622     if (!v) {
00623       fprintf(stderr, "open_cube_read failed for file %s\n", *argv);
00624       return 1;
00625     }
00626     fprintf(stderr, "open_cube_read succeeded for file %s\n", *argv);
00627 
00628     // try loading the EDM metadata now
00629     if (read_cube_metadata(v, &nsets, &meta)) {
00630       return 1; // failed to load cube file
00631     }
00632     fprintf(stderr, "read_cube_metadata succeeded for file %s\n", *argv);
00633 
00634     for (set=0; set<nsets; set++) {
00635       printf("Loading volume set: %d\n", set);   
00636       
00637       int elements = meta[set].xsize * meta[set].ysize * meta[set].zsize;
00638       printf("   Grid Elements: %d\n", elements);
00639       printf(" Grid dimensions: X: %d Y: %d Z: %d\n", 
00640              meta[set].xsize, meta[set].ysize, meta[set].zsize);
00641 
00642       float * voldata = (float *) malloc(sizeof(float) * elements);
00643       float * coldata = NULL;
00644 
00645       if (meta[set].has_color) {
00646         coldata = (float *) malloc(sizeof(float) * elements * 3);
00647       }
00648 
00649       // try loading the data sets now
00650       if (read_cube_data(v, set, voldata, coldata)) {
00651         return 1; // failed to load cube file
00652       }
00653 
00654       printf("First 6 elements:\n   ");
00655       for (i=0; i<6; i++) {
00656         printf("%g, ", voldata[i]);
00657       }
00658       printf("\n"); 
00659 
00660       printf("Last 6 elements:\n   ");
00661       for (i=elements - 6; i<elements; i++) {
00662         printf("%g, ", voldata[i]);
00663       }
00664       printf("\n"); 
00665     }
00666 
00667     close_cube_read(v);
00668   }
00669   return 0;
00670 }
00671 
00672 #endif

Generated on Wed Mar 22 13:15:28 2006 for VMD Plugins (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002