00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018 #include "largefiles.h"
00019
00020 #include <stdio.h>
00021 #include <stdlib.h>
00022 #include <math.h>
00023 #include <string.h>
00024
00025 #include "molfile_plugin.h"
00026
00027 #ifndef M_PI_2
00028 #define M_PI_2 1.57079632679489661922
00029 #endif
00030
00031 #include "periodic_table.h"
00032
00033 static const char *xsf_symtab[] = {
00034 "(unknown keyword)", "#",
00035 "BEGIN_INFO", "END_INFO",
00036 "BEGIN_BLOCK_DATAGRID_2D", "END_BLOCK_DATAGRID_2D",
00037 "BEGIN_DATAGRID_2D", "END_DATAGRID_2D",
00038 "BEGIN_BLOCK_DATAGRID_3D", "END_BLOCK_DATAGRID_3D",
00039 "BEGIN_DATAGRID_3D", "END_DATAGRID_3D",
00040 "BEGIN_BLOCK_BANDGRID_3D", "END_BLOCK_BANDGRID_3D",
00041 "ATOMS", "ANIMSTEPS", "BAND",
00042 "MOLECULE", "POLYMER", "SLAB", "CRYSTAL",
00043 "PRIMVEC", "CONVVEC", "PRIMCOORD", "CONVCOORD"
00044 };
00045
00046 typedef enum {
00047 xsf_UNKNOWN = 0, xsf_COMMENT,
00048 xsf_BEGINFO, xsf_ENDINFO,
00049 xsf_BEG_2D_BLOCK, xsf_END_2D_BLOCK,
00050 xsf_BEG_2D_DATA, xsf_END_2D_DATA,
00051 xsf_BEG_3D_BLOCK, xsf_END_3D_BLOCK,
00052 xsf_BEG_3D_DATA, xsf_END_3D_DATA,
00053 xsf_BEG_3D_BAND, xsf_END_3D_BAND,
00054 xsf_ATOMS, xsf_ANIMSTEPS, xsf_BAND,
00055 xsf_MOLECULE, xsf_POLYMER, xsf_SLAB, xsf_CRYSTAL,
00056 xsf_PRIMVEC, xsf_CONVVEC, xsf_PRIMCOORD, xsf_CONVCOORD,
00057 xsf_NR_KEYWORDS
00058 } xsf_keyword_t;
00059
00060
00061
00062 static const struct {
00063 const char *name;
00064 xsf_keyword_t kw;
00065 } xsf_aliases[] = {
00066 { "DATAGRID_2D", xsf_BEG_2D_DATA },
00067 { "DATAGRID_3D", xsf_BEG_3D_DATA },
00068 { "BEGIN_BLOCK_DATAGRID2D", xsf_BEG_2D_BLOCK },
00069 { "BEGIN_BLOCK_DATAGRID3D", xsf_BEG_3D_BLOCK },
00070 { "END_BLOCK_DATAGRID2D", xsf_END_2D_BLOCK },
00071 { "END_BLOCK_DATAGRID3D", xsf_END_3D_BLOCK },
00072 { NULL, xsf_UNKNOWN }
00073 };
00074
00075 static xsf_keyword_t lookup_keyword(const char* word)
00076 {
00077 int i, j;
00078
00079 if (word == 0) return xsf_UNKNOWN;
00080
00081
00082 j=0;
00083 for (i=0; i < (int)strlen(word); ++i) {
00084 j=i;
00085 if (! isspace(word[i])) break;
00086 }
00087
00088 for (i=1; i < xsf_NR_KEYWORDS; ++i) {
00089 if (0 == strncmp(word + j, xsf_symtab[i], strlen(xsf_symtab[i])))
00090 return (xsf_keyword_t) i;
00091 }
00092
00093
00094 for (i=0; xsf_aliases[i].kw != xsf_UNKNOWN; ++i) {
00095 const char *name = xsf_aliases[i].name;
00096
00097 if (0 == strncmp(word + j, name, strlen(name)))
00098 return xsf_aliases[i].kw;
00099
00100 }
00101
00102 return xsf_UNKNOWN;
00103 }
00104
00105
00106 typedef struct {
00107 float A, B, C, alpha, beta, gamma, cell[3][3];
00108 } xsf_box;
00109
00110 typedef struct {
00111 FILE *fd;
00112 int nvolsets;
00113 int numatoms;
00114 int animsteps;
00115 int numsteps;
00116 bool coord;
00117 char *file_name;
00118 xsf_keyword_t pbctype;
00119 molfile_volumetric_t *vol;
00120 int numvolmeta;
00121 float origin[3];
00122 float rotmat[3][3];
00123 float invmat[3][3];
00124 xsf_box box;
00125 } xsf_t;
00126
00127
00128
00129
00130
00131 static int xsf_readbox(xsf_box *box, float *x, float *y, float *z) {
00132 float A, B, C;
00133 int i;
00134
00135 if (!box) {
00136 return 1;
00137 }
00138
00139
00140 box->A = 10.0;
00141 box->B = 10.0;
00142 box->C = 10.0;
00143 box->alpha = 90.0;
00144 box->beta = 90.0;
00145 box->gamma = 90.0;
00146
00147
00148 A = sqrt( x[0]*x[0] + x[1]*x[1] + x[2]*x[2] );
00149 B = sqrt( y[0]*y[0] + y[1]*y[1] + y[2]*y[2] );
00150 C = sqrt( z[0]*z[0] + z[1]*z[1] + z[2]*z[2] );
00151 if ((A<=0) || (B<=0) || (C<=0)) {
00152 return 1;
00153 }
00154 box->A = A;
00155 box->B = B;
00156 box->C = C;
00157
00158
00159
00160 box->gamma = acos( (x[0]*y[0]+x[1]*y[1]+x[2]*y[2])/(A*B) ) * 90.0/M_PI_2;
00161 box->beta = acos( (x[0]*z[0]+x[1]*z[1]+x[2]*z[2])/(A*C) ) * 90.0/M_PI_2;
00162 box->alpha = acos( (y[0]*z[0]+y[1]*z[1]+y[2]*z[2])/(B*C) ) * 90.0/M_PI_2;
00163
00164
00165 for (i=0; i<3; ++i) {
00166 box->cell[0][i] = x[i];
00167 box->cell[1][i] = y[i];
00168 box->cell[2][i] = z[i];
00169 }
00170
00171 return 0;
00172 }
00173
00174
00175 static void xsf_buildrotmat(xsf_t *xsf, float *a, float *b)
00176 {
00177
00178 const double len = sqrt(a[0]*a[0] + a[1]*a[1]);
00179 const double phi = atan2((double) a[2], (double) len);
00180 const double theta = atan2((double) a[1], (double) a[0]);
00181
00182 const double cph = cos(phi);
00183 const double cth = cos(theta);
00184 const double sph = sin(phi);
00185 const double sth = sin(theta);
00186
00187
00188 const double psi = atan2(-sph*cth*b[0] - sph*sth*b[1] + cph*b[2],-sth*b[0] + cth*b[1]);
00189 const double cps = cos(psi);
00190 const double sps = sin(psi);
00191
00192 const double r[3][3] = {
00193 { cph*cth, cph*sth, sph},
00194 {-sth*cps - sph*cth*sps, cth*cps - sph*sth*sps, cph*sps},
00195 { sth*sps - sph*cth*cps, -cth*sps - sph*sth*cps, cph*cps}
00196 };
00197
00198 for (int i=0; i<3; ++i) {
00199 for (int j=0; j<3; ++j) {
00200 xsf->rotmat[i][j] = r[i][j];
00201 }
00202 }
00203 }
00204
00205 static void xsf_buildinvmat(xsf_t *xsf, float *a, float *b, float *c)
00206 {
00207 float det, id;
00208
00209 det = a[0]*b[1]*c[2] + b[0]*c[1]*a[2] + c[0]*a[1]*b[2]
00210 - a[0]*c[1]*b[2] - b[0]*a[1]*c[2] - c[0]*b[1]*a[2];
00211
00212 id = 1.0 / det;
00213 xsf->invmat[0][0] = id * ( b[1]*c[2]-b[2]*c[1] );
00214 xsf->invmat[1][0] = id * ( a[2]*c[1]-a[1]*c[2] );
00215 xsf->invmat[2][0] = id * ( a[1]*b[2]-a[2]*b[1] );
00216 xsf->invmat[0][1] = id * ( b[2]*c[0]-b[0]*c[2] );
00217 xsf->invmat[1][1] = id * ( a[0]*c[2]-a[2]*c[0] );
00218 xsf->invmat[2][1] = id * ( a[2]*b[0]-a[0]*b[2] );
00219 xsf->invmat[0][2] = id * ( b[0]*c[1]-b[1]*c[0] );
00220 xsf->invmat[1][2] = id * ( a[1]*c[0]-a[0]*c[1] );
00221 xsf->invmat[2][2] = id * ( a[0]*b[1]-a[1]*b[0] );
00222 }
00223
00224
00225
00226 static void eatline(FILE *fd) {
00227 char readbuf[1025];
00228 fgets(readbuf, 1024, fd);
00229 }
00230
00231 static bool xsf_read_cell(FILE *fd, float *a, float *b, float *c)
00232 {
00233 return (9 == fscanf(fd, "%f%f%f%f%f%f%f%f%f",
00234 &a[0],&a[1],&a[2],
00235 &b[0],&b[1],&b[2],
00236 &c[0],&c[1],&c[2]));
00237 }
00238
00239 static void close_xsf_read(void *v);
00240
00241 static void *open_xsf_read(const char *filepath, const char *filetype,
00242 int *natoms) {
00243 FILE *fd;
00244 xsf_t *xsf;
00245 int i,j;
00246
00247 fd = fopen(filepath, "rb");
00248 if (!fd)
00249 return NULL;
00250
00251 xsf = new xsf_t;
00252 xsf->fd = fd;
00253 xsf->vol = NULL;
00254 xsf->numvolmeta = 0;
00255 xsf->coord = false;
00256 xsf->nvolsets = 0;
00257 xsf->numatoms = 0;
00258 xsf->numsteps = 0;
00259 xsf->file_name = strdup(filepath);
00260
00261 xsf->pbctype = xsf_MOLECULE;
00262
00263
00264 for (i=0; i<3; ++i) {
00265 for (j=0; j<3; ++j) {
00266 xsf->rotmat[i][j] = 0.0;
00267 }
00268 }
00269 for (i=0; i<3; ++i) {
00270 xsf->origin[i] = 0.0;
00271 xsf->rotmat[i][i] = 1.0;
00272 }
00273
00274
00275
00276
00277 char readbuf[256];
00278 xsf_keyword_t kw;
00279
00280
00281 do {
00282 if (NULL == fgets(readbuf, 256, xsf->fd)) break;
00283
00284 again:
00285 kw = lookup_keyword(readbuf);
00286 #ifdef TEST_PLUGIN
00287 fprintf(stderr, "keyword: %d / %s", kw, readbuf);
00288 #endif
00289
00290 switch (kw) {
00291 case xsf_ANIMSTEPS:
00292 #ifdef TEST_PLUGIN
00293 {
00294 int n;
00295 if (1 == sscanf(readbuf, "%*s%d", &n)) {
00296 fprintf(stderr, "ANIMSTEPS: found %d steps\n", n);
00297 }
00298 }
00299 #endif
00300 break;
00301
00302 case xsf_ATOMS:
00303
00304 ++ xsf->numsteps;
00305 if (xsf->numatoms == 0) {
00306 while (fgets(readbuf, 256, xsf->fd)) {
00307 float x,y,z;
00308
00309 if (3 == sscanf(readbuf, "%*s%f%f%f", &x, &y, &z)) {
00310 ++ xsf->numatoms;
00311 } else {
00312
00313
00314 goto again;
00315 break;
00316 }
00317 }
00318 #ifdef TEST_PLUGIN
00319 fprintf(stderr, "ATOMS: found %d atoms\n", xsf->numatoms);
00320 #endif
00321 } else {
00322 int n;
00323 for (n=0; n < xsf->numatoms; ++n) eatline(xsf->fd);
00324 }
00325 break;
00326
00327 case xsf_PRIMCOORD:
00328
00329 if(fgets(readbuf, 256, xsf->fd)) {
00330 int mol, mult;
00331
00332 if (xsf->numatoms == 0) {
00333 if (2 == sscanf(readbuf, "%d%d", &mol, &mult)) {
00334 xsf->numatoms = mol * mult;
00335 } else {
00336 xsf->numatoms = mol;
00337 }
00338 }
00339
00340 int n;
00341 for (n=0; n < xsf->numatoms; ++n) eatline(xsf->fd);
00342 ++ xsf->numsteps;
00343
00344 #ifdef TEST_PLUGIN
00345 fprintf(stderr, "PRIMCOORD: found %d atoms\n", xsf->numatoms);
00346 #endif
00347 }
00348 break;
00349
00350 case xsf_CONVCOORD:
00351
00352 if(fgets(readbuf, 256, xsf->fd)) {
00353 int mol, mult, num;
00354
00355 num = 0;
00356 if (2 == sscanf(readbuf, "%d%d", &mol, &mult)) {
00357 num = mol * mult;
00358 }
00359
00360
00361 int n;
00362 for (n=0; n < num; ++n) eatline(xsf->fd);
00363 #ifdef TEST_PLUGIN
00364 fprintf(stderr, "CONVCOORD: ignoring %d atoms\n", num);
00365 #endif
00366 }
00367 break;
00368
00369 case xsf_PRIMVEC:
00370 {
00371 float a[3], b[3], c[3];
00372
00373 if (xsf_read_cell(xsf->fd, a, b, c)) {
00374 xsf_buildrotmat(xsf, a, b);
00375 } else {
00376 fprintf(stderr, "xsfplugin) WARNING: error reading unit cell. ignoring unit cell info.\n");
00377 }
00378 }
00379 break;
00380
00381 case xsf_CONVVEC:
00382 {
00383 int n;
00384 for (n=0; n < 3; ++n) eatline(xsf->fd);
00385 }
00386 break;
00387
00388 case xsf_BEG_3D_BLOCK:
00389
00390
00391
00392
00393 if (xsf->vol == NULL) {
00394 xsf->numvolmeta = 32;
00395 xsf->vol = new molfile_volumetric_t[xsf->numvolmeta];
00396 }
00397
00398
00399 fgets(readbuf, 256, xsf->fd);
00400 printf("xsfplugin) found grid data block: %s", readbuf);
00401
00402 do {
00403
00404 if (NULL == fgets(readbuf, 256, xsf->fd)) break;
00405 switch (lookup_keyword(readbuf)) {
00406 case xsf_BEG_3D_DATA:
00407 {
00408 int n;
00409 molfile_volumetric_t *set;
00410 float a[3], b[3], c[3];
00411
00412 ++ xsf->nvolsets;
00413
00414
00415 if (xsf->nvolsets > xsf->numvolmeta) {
00416 molfile_volumetric_t *ptr = xsf->vol;
00417 xsf->vol = new molfile_volumetric_t[2 * xsf->numvolmeta];
00418 memcpy((void *)xsf->vol, (void *)ptr, xsf->numvolmeta*sizeof(molfile_volumetric_t));
00419 xsf->numvolmeta *= 2;
00420 delete[] ptr;
00421 }
00422
00423
00424 set = &(xsf->vol[xsf->nvolsets - 1]);
00425 set->has_color = 0;
00426
00427
00428
00429 strncpy(set->dataname, readbuf, 255);
00430
00431
00432
00433 fgets(readbuf, 256, xsf->fd);
00434 sscanf(readbuf, "%d%d%d", &(set->xsize), &(set->ysize), &(set->zsize));
00435 fgets(readbuf, 256, xsf->fd);
00436 sscanf(readbuf, "%f%f%f", &(set->origin[0]), &(set->origin[1]), &(set->origin[2]));
00437 fgets(readbuf, 256, xsf->fd);
00438 sscanf(readbuf, "%f%f%f", &a[0], &a[1], &a[2]);
00439 fgets(readbuf, 256, xsf->fd);
00440 sscanf(readbuf, "%f%f%f", &b[0], &b[1], &b[2]);
00441 fgets(readbuf, 256, xsf->fd);
00442 sscanf(readbuf, "%f%f%f", &c[0], &c[1], &c[2]);
00443
00444
00445
00446 -- set->xsize; -- set->ysize; -- set->zsize;
00447
00448
00449 for (n=0; n<3; ++n) {
00450 set->xaxis[n] = xsf->rotmat[n][0] * a[0]
00451 + xsf->rotmat[n][1] * a[1] + xsf->rotmat[n][2] * a[2];
00452
00453 set->yaxis[n] = xsf->rotmat[n][0] * b[0]
00454 + xsf->rotmat[n][1] * b[1] + xsf->rotmat[n][2] * b[2];
00455
00456 set->zaxis[n] = xsf->rotmat[n][0] * c[0]
00457 + xsf->rotmat[n][1] * c[1] + xsf->rotmat[n][2] * c[2];
00458 }
00459
00460 do {
00461 fgets(readbuf, 256, xsf->fd);
00462 } while (xsf_END_3D_DATA != lookup_keyword(readbuf));
00463
00464 #if 1
00465
00466
00467
00468
00469
00470
00471
00472
00473 set->origin[0] -= 0.5 * ( set->xaxis[0] / (double) set->xsize
00474 + set->yaxis[0] / (double) set->ysize
00475 + set->zaxis[0] / (double) set->zsize);
00476 set->origin[1] -= 0.5 * ( set->xaxis[1] / (double) set->xsize
00477 + set->yaxis[1] / (double) set->ysize
00478 + set->zaxis[1] / (double) set->zsize);
00479 set->origin[2] -= 0.5 * ( set->xaxis[2] / (double) set->xsize
00480 + set->yaxis[2] / (double) set->ysize
00481 + set->zaxis[2] / (double) set->zsize);
00482 #endif
00483 }
00484 break;
00485
00486 default:
00487 break;
00488 }
00489 } while (xsf_END_3D_BLOCK != lookup_keyword(readbuf));
00490
00491 #ifdef TEST_PLUGIN
00492 fprintf(stderr, "found %d volumetric data sets\n", xsf->nvolsets);
00493 #endif
00494 break;
00495
00496 case xsf_BEG_2D_BLOCK:
00497 do {
00498 fgets(readbuf, 256, xsf->fd);
00499 } while (xsf_END_2D_BLOCK != lookup_keyword(readbuf));
00500 break;
00501
00502
00503 case xsf_MOLECULE:
00504 case xsf_SLAB:
00505 case xsf_POLYMER:
00506 case xsf_CRYSTAL:
00507 xsf->pbctype = kw;
00508 break;
00509
00510 case xsf_COMMENT:
00511 case xsf_UNKNOWN:
00512 default:
00513 break;
00514
00515 }
00516 } while (! (feof(xsf->fd) || ferror(xsf->fd)));
00517 #ifdef TEST_PLUGIN
00518 fprintf(stderr, "total of %d coordinate sets\n", xsf->numsteps);
00519 #endif
00520 rewind(xsf->fd);
00521 *natoms = xsf->numatoms;
00522 return xsf;
00523 }
00524
00525
00526 static int read_xsf_structure(void *v, int *optflags, molfile_atom_t *atoms) {
00527 int i;
00528 xsf_t *xsf = (xsf_t *)v;
00529
00530
00531 if (xsf->numatoms < 1) return MOLFILE_SUCCESS;
00532
00533
00534
00535 rewind(xsf->fd);
00536
00537
00538
00539
00540
00541 char readbuf[1024];
00542 do {
00543 if (NULL == fgets(readbuf, 256, xsf->fd)) break;
00544
00545 switch (lookup_keyword(readbuf)) {
00546
00547 case xsf_PRIMCOORD:
00548 eatline(xsf->fd);
00549
00550 case xsf_ATOMS:
00551
00552
00553
00554 *optflags = MOLFILE_ATOMICNUMBER | MOLFILE_MASS | MOLFILE_RADIUS;
00555
00556 for(i=0; i<xsf->numatoms; ++i) {
00557 int j;
00558 char *k;
00559 float coord;
00560 molfile_atom_t *atom;
00561
00562 char buffer[1024];
00563 k = fgets(readbuf, 1024, xsf->fd);
00564 atom = atoms + i;
00565 j=sscanf(readbuf, "%s %f %f %f", buffer, &coord, &coord, &coord);
00566 if (k == NULL) {
00567 fprintf(stderr, "xsf structure) missing atom(s) in file '%s'\n",xsf->file_name);
00568 fprintf(stderr, "xsf structure) expecting '%d' atoms, found only '%d'\n",xsf->numatoms,i+1);
00569 return MOLFILE_ERROR;
00570 } else if (j < 4) {
00571 fprintf(stderr, "xsf structure) missing type or coordinate(s) in file '%s' for atom '%d'\n",
00572 xsf->file_name, i+1);
00573 return MOLFILE_ERROR;
00574 }
00575
00576
00577
00578 if (isdigit(buffer[0])) {
00579 int idx;
00580 idx = atoi(buffer);
00581 strncpy(atom->name, get_pte_label(idx), sizeof(atom->name));
00582 atom->atomicnumber = idx;
00583 atom->mass = get_pte_mass(idx);
00584 atom->radius = get_pte_vdw_radius(idx);
00585 } else {
00586 int idx;
00587 strncpy(atom->name, buffer, sizeof(atom->name));
00588 idx = get_pte_idx(buffer);
00589 atom->atomicnumber = idx;
00590 atom->mass = get_pte_mass(idx);
00591 atom->radius = get_pte_vdw_radius(idx);
00592 }
00593 strncpy(atom->type, atom->name, sizeof(atom->type));
00594 atom->resname[0] = '\0';
00595 atom->resid = 1;
00596 atom->chain[0] = '\0';
00597 atom->segid[0] = '\0';
00598 #ifdef TEST_PLUGIN
00599 fprintf(stderr,"xsf structure) atom %4d: %s / mass= %f\n", i, atom->name, atom->mass);
00600 #endif
00601 }
00602
00603
00604 rewind(xsf->fd);
00605 return MOLFILE_SUCCESS;
00606 break;
00607
00608
00609 case xsf_PRIMVEC:
00610 {
00611 float a[3], b[3], c[3];
00612
00613 if (xsf_read_cell(xsf->fd, a, b, c)) {
00614 xsf_readbox(&(xsf->box), a, b, c);
00615 xsf_buildrotmat(xsf, a, b);
00616
00617 if ((fabs((double) a[1]) + fabs((double) a[2]) + fabs((double) b[2]))
00618 > 0.001) {
00619 fprintf(stderr, "xsfplugin) WARNING: Coordinates will be rotated to comply \n"
00620 "xsfplugin) with VMD's conventions for periodic display...\n");
00621 }
00622 xsf_buildinvmat(xsf, a, b, c);
00623
00624 #if defined(TEST_PLUGIN)
00625 printf("cell vectors:\n");
00626 printf("<a>: %12.8f %12.8f %12.8f\n", a[0], a[1], a[2]);
00627 printf("<b>: %12.8f %12.8f %12.8f\n", b[0], b[1], b[2]);
00628 printf("<c>: %12.8f %12.8f %12.8f\n", c[0], c[1], c[2]);
00629 printf("cell dimensions:\n");
00630 printf("a= %12.8f b= %12.8f c= %12.8f\n", xsf->box.A, xsf->box.B, xsf->box.C);
00631 printf("alpha= %6.2f beta= %6.2f gamma= %6.2f\n", xsf->box.alpha, xsf->box.beta, xsf->box.gamma);
00632 printf("reciprocal cell vectors:\n");
00633 printf("i: %12.8f %12.8f %12.8f\n", xsf->invmat[0][0], xsf->invmat[0][1], xsf->invmat[0][2]);
00634 printf("k: %12.8f %12.8f %12.8f\n", xsf->invmat[1][0], xsf->invmat[1][1], xsf->invmat[1][2]);
00635 printf("l: %12.8f %12.8f %12.8f\n", xsf->invmat[2][0], xsf->invmat[2][1], xsf->invmat[2][2]);
00636 printf("cell rotation matrix:\n");
00637 printf("x: %12.8f %12.8f %12.8f\n", xsf->rotmat[0][0], xsf->rotmat[0][1], xsf->rotmat[0][2]);
00638 printf("y: %12.8f %12.8f %12.8f\n", xsf->rotmat[1][0], xsf->rotmat[1][1], xsf->rotmat[1][2]);
00639 printf("z: %12.8f %12.8f %12.8f\n", xsf->rotmat[2][0], xsf->rotmat[2][1], xsf->rotmat[2][2]);
00640 #endif
00641 }
00642 }
00643 break;
00644
00645 default:
00646 break;
00647 }
00648 } while (! (feof(xsf->fd) || ferror(xsf->fd)));
00649
00650
00651 return MOLFILE_ERROR;
00652
00653 }
00654
00655
00656 static int read_xsf_timestep(void *v, int natoms, molfile_timestep_t *ts) {
00657 int i;
00658
00659 xsf_t *xsf = (xsf_t *)v;
00660
00661
00662
00663
00664
00665
00666 char readbuf[1024];
00667 do {
00668 if (NULL == fgets(readbuf, 256, xsf->fd)) break;
00669
00670 switch (lookup_keyword(readbuf)) {
00671
00672 case xsf_PRIMCOORD:
00673 eatline(xsf->fd);
00674
00675 case xsf_ATOMS:
00676
00677
00678 for(i=0; i<natoms; ++i) {
00679 int j, n;
00680 char *k, buffer[1024];
00681 float x, y, z;
00682
00683 k = fgets(readbuf, 1024, xsf->fd);
00684 j=sscanf(readbuf, "%s %f %f %f", buffer, &x, &y, &z);
00685
00686 if (k == NULL) {
00687 return MOLFILE_ERROR;
00688 } else if (j < 4) {
00689 fprintf(stderr, "xsf structure) missing type or coordinate(s) in file '%s' for atom '%d'\n",
00690 xsf->file_name, i+1);
00691 return MOLFILE_ERROR;
00692 } else if (j>=3) {
00693 if (ts != NULL) {
00694
00695
00696 float xf, yf, zf;
00697
00698
00699 #ifdef TEST_PLUGIN
00700 printf("wrap PBC: before: %12.6f %12.6f %12.6f\n", x, y, z);
00701 #endif
00702 switch(xsf->pbctype) {
00703 case xsf_CRYSTAL:
00704 xf = xsf->invmat[0][0] * x + xsf->invmat[0][1] * y + xsf->invmat[0][2] * z;
00705 xf = xf - floor((double)xf);
00706 yf = xsf->invmat[1][0] * x + xsf->invmat[1][1] * y + xsf->invmat[1][2] * z;
00707 yf = yf - floor((double)yf);
00708 zf = xsf->invmat[2][0] * x + xsf->invmat[2][1] * y + xsf->invmat[2][2] * z;
00709 zf = zf - floor((double)zf);
00710 x = xsf->box.cell[0][0] * xf + xsf->box.cell[0][1] * yf + xsf->box.cell[0][2] * zf;
00711 y = xsf->box.cell[1][0] * xf + xsf->box.cell[1][1] * yf + xsf->box.cell[1][2] * zf;
00712 z = xsf->box.cell[2][0] * xf + xsf->box.cell[2][1] * yf + xsf->box.cell[2][2] * zf;
00713 break;
00714
00715 case xsf_SLAB:
00716 xf = xsf->invmat[0][0] * x + xsf->invmat[0][1] * y + xsf->invmat[0][2] * z;
00717 xf = xf - floor((double)xf);
00718 yf = xsf->invmat[1][0] * x + xsf->invmat[1][1] * y + xsf->invmat[1][2] * z;
00719 yf = yf - floor((double)yf);
00720 zf = xsf->invmat[2][0] * x + xsf->invmat[2][1] * y + xsf->invmat[2][2] * z;
00721 x = xsf->box.cell[0][0] * xf + xsf->box.cell[0][1] * yf + xsf->box.cell[0][2] * zf;
00722 y = xsf->box.cell[1][0] * xf + xsf->box.cell[1][1] * yf + xsf->box.cell[1][2] * zf;
00723 z = xsf->box.cell[2][0] * xf + xsf->box.cell[2][1] * yf + xsf->box.cell[2][2] * zf;
00724 break;
00725
00726 case xsf_POLYMER:
00727 xf = xsf->invmat[0][0] * x + xsf->invmat[0][1] * y + xsf->invmat[0][2] * z;
00728 xf = xf - floor((double)xf);
00729 yf = xsf->invmat[1][0] * x + xsf->invmat[1][1] * y + xsf->invmat[1][2] * z;
00730 zf = xsf->invmat[2][0] * x + xsf->invmat[2][1] * y + xsf->invmat[2][2] * z;
00731 x = xsf->box.cell[0][0] * xf + xsf->box.cell[0][1] * yf + xsf->box.cell[0][2] * zf;
00732 y = xsf->box.cell[1][0] * xf + xsf->box.cell[1][1] * yf + xsf->box.cell[1][2] * zf;
00733 z = xsf->box.cell[2][0] * xf + xsf->box.cell[2][1] * yf + xsf->box.cell[2][2] * zf;
00734 break;
00735
00736 case xsf_MOLECULE:
00737 xf = x;
00738 yf = y;
00739 zf = z;
00740 break;
00741
00742 default:
00743 break;
00744 }
00745
00746 #ifdef TEST_PLUGIN
00747 printf("wrap PBC: fract: %12.6f %12.6f %12.6f\n", xf, yf, zf);
00748 printf("wrap PBC: after: %12.6f %12.6f %12.6f\n", x, y, z);
00749 #endif
00750
00751
00752
00753
00754
00755
00756 x -= xsf->origin[0];
00757 y -= xsf->origin[1];
00758 z -= xsf->origin[2];
00759
00760 for (n=0; n<3; ++n) {
00761 ts->coords[3*i + n] = (xsf->origin[n] + xsf->rotmat[n][0] * x
00762 + xsf->rotmat[n][1] * y + xsf->rotmat[n][2] * z);
00763 }
00764 }
00765 } else {
00766 break;
00767 }
00768 }
00769
00770 if (ts != NULL) {
00771 ts->A = xsf->box.A;
00772 ts->B = xsf->box.B;
00773 ts->C = xsf->box.C;
00774 ts->alpha = xsf->box.alpha;
00775 ts->beta = xsf->box.beta;
00776 ts->gamma = xsf->box.gamma;
00777 }
00778
00779 return MOLFILE_SUCCESS;
00780 break;
00781
00782
00783 case xsf_PRIMVEC:
00784 {
00785 float a[3], b[3], c[3];
00786
00787 if (xsf_read_cell(xsf->fd, a, b, c)) {
00788 xsf_readbox(&(xsf->box), a,b,c);
00789 xsf_buildrotmat(xsf, a, b);
00790
00791 if ((fabs((double) a[1]) + fabs((double) a[2]) + fabs((double) b[2]))
00792 > 0.001) {
00793 fprintf(stderr, "xsfplugin) WARNING: Coordinates will be rotated to comply \n"
00794 "xsfplugin) with VMD's conventions for periodic display...\n");
00795 }
00796 xsf_buildinvmat(xsf, a, b, c);
00797
00798 #if defined(TEST_PLUGIN)
00799 printf("new cell vectors:\n");
00800 printf("<a>: %12.8f %12.8f %12.8f\n", a[0], a[1], a[2]);
00801 printf("<b>: %12.8f %12.8f %12.8f\n", b[0], b[1], b[2]);
00802 printf("<c>: %12.8f %12.8f %12.8f\n", c[0], c[1], c[2]);
00803 printf("new cell dimensions:\n");
00804 printf("a= %12.8f b= %12.8f c= %12.8f\n", xsf->box.A, xsf->box.B, xsf->box.C);
00805 printf("alpha= %6.2f beta= %6.2f gamma= %6.2f\n", xsf->box.alpha, xsf->box.beta, xsf->box.gamma);
00806 printf("new reciprocal cell vectors:\n");
00807 printf("i: %12.8f %12.8f %12.8f\n", xsf->invmat[0][0], xsf->invmat[0][1], xsf->invmat[0][2]);
00808 printf("k: %12.8f %12.8f %12.8f\n", xsf->invmat[1][0], xsf->invmat[1][1], xsf->invmat[1][2]);
00809 printf("l: %12.8f %12.8f %12.8f\n", xsf->invmat[2][0], xsf->invmat[2][1], xsf->invmat[2][2]);
00810 printf("new cell rotation matrix:\n");
00811 printf("x: %12.8f %12.8f %12.8f\n", xsf->rotmat[0][0], xsf->rotmat[0][1], xsf->rotmat[0][2]);
00812 printf("y: %12.8f %12.8f %12.8f\n", xsf->rotmat[1][0], xsf->rotmat[1][1], xsf->rotmat[1][2]);
00813 printf("z: %12.8f %12.8f %12.8f\n", xsf->rotmat[2][0], xsf->rotmat[2][1], xsf->rotmat[2][2]);
00814 #endif
00815 }
00816 }
00817 break;
00818
00819 default:
00820 break;
00821 }
00822 } while (! (feof(xsf->fd) || ferror(xsf->fd)));
00823
00824
00825 return MOLFILE_ERROR;
00826 }
00827
00828 static int read_xsf_metadata(void *v, int *nvolsets,
00829 molfile_volumetric_t **metadata) {
00830 xsf_t *xsf = (xsf_t *)v;
00831 *nvolsets = xsf->nvolsets;
00832 *metadata = xsf->vol;
00833
00834 return MOLFILE_SUCCESS;
00835 }
00836
00837 static int read_xsf_data(void *v, int set, float *datablock, float *colorblock) {
00838 xsf_t *xsf = (xsf_t *)v;
00839 const char *block = xsf->vol[set].dataname;
00840
00841 fprintf(stderr, "xsfplugin) trying to read xsf data set %d: %s\n", set, block);
00842
00843 int xsize = xsf->vol[set].xsize;
00844 int ysize = xsf->vol[set].ysize;
00845 int zsize = xsf->vol[set].zsize;
00846 int x, y, z;
00847 int n;
00848 char readbuf[1024];
00849 float dummy;
00850
00851
00852 rewind(xsf->fd);
00853 do {
00854 if (NULL == fgets(readbuf, 1024, xsf->fd)) return MOLFILE_ERROR;
00855 } while (strncmp(readbuf, block, 1024));
00856
00857 eatline(xsf->fd);
00858 eatline(xsf->fd);
00859 eatline(xsf->fd);
00860 eatline(xsf->fd);
00861 eatline(xsf->fd);
00862
00863
00864
00865
00866 n = 0;
00867 for (z=0; z<zsize+1; z++) {
00868 for (y=0; y<ysize+1; y++) {
00869 for (x=0; x<xsize+1; x++) {
00870 if((x>=xsize) || (y>=ysize) || (z>=zsize)) {
00871 if (fscanf(xsf->fd, "%f", &dummy) != 1) return MOLFILE_ERROR;
00872 } else {
00873 if (fscanf(xsf->fd, "%f", &datablock[n]) != 1) return MOLFILE_ERROR;
00874 ++ n;
00875 }
00876 }
00877 }
00878 }
00879 rewind(xsf->fd);
00880 return MOLFILE_SUCCESS;
00881 }
00882
00883 static void close_xsf_read(void *v) {
00884 xsf_t *xsf = (xsf_t *)v;
00885
00886 fclose(xsf->fd);
00887 if (xsf->vol) {
00888 delete[] xsf->vol;
00889 }
00890 free(xsf->file_name);
00891 delete xsf;
00892 }
00893
00894
00895
00896
00897 static molfile_plugin_t plugin = {
00898 vmdplugin_ABIVERSION,
00899 MOLFILE_PLUGIN_TYPE,
00900 "xsf",
00901 "XSF",
00902 "Axel Kohlmeyer, John E. Stone",
00903 0,
00904 5,
00905 VMDPLUGIN_THREADSAFE,
00906 "xsf",
00907 open_xsf_read,
00908 read_xsf_structure,
00909 0,
00910 read_xsf_timestep,
00911 close_xsf_read,
00912 0,
00913 0,
00914 0,
00915 0,
00916 read_xsf_metadata,
00917 read_xsf_data,
00918 0
00919 };
00920
00921 VMDPLUGIN_API int VMDPLUGIN_init(void) { return VMDPLUGIN_SUCCESS; }
00922 VMDPLUGIN_API int VMDPLUGIN_fini(void) { return VMDPLUGIN_SUCCESS; }
00923 VMDPLUGIN_API int VMDPLUGIN_register(void *v, vmdplugin_register_cb cb) {
00924 (*cb)(v, (vmdplugin_t *)&plugin);
00925 return VMDPLUGIN_SUCCESS;
00926 }
00927
00928
00929
00930 #ifdef TEST_PLUGIN
00931
00932 int main(int argc, char *argv[]) {
00933 int natoms, optflags;
00934 void *v;
00935 int i, nvolsets, set;
00936 molfile_volumetric_t * meta;
00937
00938
00939 printf("got index: %d\n", lookup_keyword(" ATOMS "));
00940
00941 while (--argc) {
00942 ++argv;
00943 v = open_xsf_read(*argv, "xsf", &natoms);
00944 if (!v) {
00945 fprintf(stderr, "open_xsf_read failed for file %s\n", *argv);
00946 return 1;
00947 }
00948 fprintf(stderr, "open_xsf_read succeeded for file %s\n", *argv);
00949 fprintf(stderr, "input contains %d atoms\n", natoms);
00950
00951 molfile_atom_t atoms[natoms];
00952
00953
00954 i = read_xsf_structure(v, &optflags, atoms);
00955 if (i) {
00956 fprintf(stderr, "read_xsf_structure failed for file %s\n", *argv);
00957 return 1;
00958 }
00959 fprintf(stderr, "read_xsf_structure succeeded for file %s\n", *argv);
00960
00961
00962 if (read_xsf_metadata(v, &nvolsets, &meta)) {
00963 return 1;
00964 }
00965 fprintf(stderr, "read_xsf_metadata succeeded for file %s\n", *argv);
00966
00967 for (set=0; set<nvolsets; set++) {
00968 printf("Loading volume set: %d\n", set);
00969
00970 int elements = meta[set].xsize * meta[set].ysize * meta[set].zsize;
00971 printf(" Grid Elements: %d\n", elements);
00972 printf(" Grid dimensions: X: %d Y: %d Z: %d\n",
00973 meta[set].xsize, meta[set].ysize, meta[set].zsize);
00974
00975 float * voldata = (float *) malloc(sizeof(float) * elements);
00976 float * coldata = NULL;
00977
00978 if (meta[set].has_color) {
00979 coldata = (float *) malloc(sizeof(float) * elements * 3);
00980 }
00981
00982
00983 if (read_xsf_data(v, set, voldata, coldata)) {
00984 return 1;
00985 }
00986
00987 printf("First 6 elements:\n ");
00988 for (i=0; i<6; i++) {
00989 printf("%g, ", voldata[i]);
00990 }
00991 printf("\n");
00992
00993 printf("Last 6 elements:\n ");
00994 for (i=elements - 6; i<elements; i++) {
00995 printf("%g, ", voldata[i]);
00996 }
00997 printf("\n");
00998 }
00999
01000 molfile_timestep_t ts;
01001 ts.coords = new float[3*natoms];
01002
01003 if (read_xsf_timestep(v, natoms, &ts)) {
01004 printf("read_xsf_timestep 0: failed\n");
01005 } else {
01006 printf("read_xsf_timestep 0: success\n");
01007 }
01008
01009 if (read_xsf_timestep(v, natoms, &ts)) {
01010 printf("read_xsf_timestep 1: failed\n");
01011 } else {
01012 printf("read_xsf_timestep 1: success\n");
01013 }
01014
01015 close_xsf_read(v);
01016 }
01017 return 0;
01018 }
01019
01020 #endif
01021