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

vaspposcarplugin.c

Go to the documentation of this file.
00001 /***************************************************************************
00002  * RCS INFORMATION:
00003  *
00004  *      $RCSfile: vaspposcarplugin.c,v $
00005  *      $Author: johns $       $Locker:  $             $State: Exp $
00006  *      $Revision: 1.10 $       $Date: 2017/10/27 18:02:32 $
00007  *
00008  ***************************************************************************/
00009 
00010 /*
00011  *  VASP plugins for VMD
00012  *  Sung Sakong, Dept. of Phys., Univsity Duisburg-Essen
00013  *  
00014  *  VASP manual   
00015  *  http://cms.mpi.univie.ac.at/vasp/
00016  * 
00017  *  LINUX
00018  *  g++ -O2 -Wall -I. -I$VMDBASEDIR/plugins/include -c vaspposcarplugin.c
00019  *  ld -shared -o vaspposcarplugin.so vaspposcarplugin.o
00020  *
00021  *  MACOSX
00022  *  c++ -O2 -Wall -I. -I$VMDBASEDIR/plugins/include -c vaspposcarplugin.c
00023  *  c++ -bundle -o vaspposcarplugin.so vaspposcarplugin.o
00024  *
00025  *  Install
00026  *  copy vaspposcarplugin.so $VMDBASEDIR/plugins/$ARCH/molfile
00027  */
00028 
00029 #include <stdio.h>
00030 #include <stdlib.h>
00031 #include <string.h>
00032 
00033 #include "molfile_plugin.h"
00034 #include "vaspplugin.h"
00035 #include "periodic_table.h"
00036 
00037 
00038 static void *open_vaspposcar_read(const char *filename, const char *filetype, int *natoms)
00039 {
00040   vasp_plugindata_t *data;
00041   char lineptr[LINESIZE];
00042   int i;
00043 
00044   /* Verify that input is OK */
00045   if (!filename || !natoms) return NULL;
00046 
00047   /* Start with undefined value; set it after successful read */
00048   *natoms = MOLFILE_NUMATOMS_UNKNOWN;
00049 
00050   data = vasp_plugindata_malloc();
00051   if (!data) return NULL;
00052 
00053   /* VASP4 is assumed in default */
00054   data->version = 4;
00055   data->file = fopen(filename, "rb");
00056   if (!data->file) {
00057     vasp_plugindata_free(data);
00058     return NULL;
00059   }
00060 
00061   data->filename = strdup(filename);
00062 
00063   /* Read title line */
00064   fgets(lineptr, LINESIZE, data->file);
00065   data->titleline = strdup(lineptr);
00066 
00067   /* Ignore rest of header up to the line with atom numbers */
00068   for (i = 0; i < 5; ++i) fgets(lineptr, LINESIZE, data->file);
00069 
00070   /* Read the number of atoms per atom type */
00071   data->numatoms = 0;
00072   for (i = 0; i < MAXATOMTYPES; ++i) {
00073     char const *tmplineptr = lineptr;
00074     char const *token = (i == 0 ? strtok(lineptr, " ") : strtok(NULL, " "));
00075     int const n = (token ? atoi(token) : -1);
00076 
00077     /* if fails to read number of atoms, then assume VASP5 */
00078     if (i == 0 && n <= 0) {
00079       data->version = 5;
00080       free(data->titleline);
00081       data->titleline =  strdup(tmplineptr);
00082       fgets(lineptr, LINESIZE, data->file);
00083       break;
00084     } else if (n <= 0) 
00085       break;
00086 
00087     data->eachatom[i] = n;
00088     data->numatoms += n;
00089   }
00090 
00091   if (data->version == 5) {
00092     data->numatoms = 0;
00093     for (i = 0; i < MAXATOMTYPES; ++i) {
00094       char const *token = (i == 0 ? strtok(lineptr, " ") : strtok(NULL, " "));
00095       int const n = (token ? atoi(token) : -1);
00096       
00097       if (n <= 0) break;
00098       
00099       data->eachatom[i] = n;
00100       data->numatoms += n;
00101     }
00102   }
00103 
00104   if (data->numatoms == 0) {
00105     vasp_plugindata_free(data);
00106     fprintf(stderr, "\n\nVASP POSCAR read) ERROR: file '%s' does not have list of atom numbers.\n", filename);
00107     return NULL;
00108   }
00109 
00110   *natoms = data->numatoms;
00111   rewind(data->file);
00112 
00113   return data;
00114 }
00115 
00116 
00117 static int read_vaspposcar_structure(void *mydata, int *optflags, molfile_atom_t *atoms)
00118 {
00119   vasp_plugindata_t *data = (vasp_plugindata_t *)mydata;
00120   FILE *potcar = NULL;
00121   int atomcount, i;
00122   char lineptr[LINESIZE], potcarfile[1000], *cp;
00123  
00124   if (!data || !optflags || !atoms) return MOLFILE_ERROR;
00125 
00126   *optflags = MOLFILE_MASS; /* we set atom mass from the PTE. */
00127   *optflags |= MOLFILE_ATOMICNUMBER | MOLFILE_RADIUS; 
00128 
00129   /* This plugin can read both POSCAR and CONTCAR files */
00130   strcpy(potcarfile, data->filename);
00131   cp = strstr(potcarfile, "POSCAR");
00132   if (!cp) cp = strstr(potcarfile, "CONTCAR");
00133 
00134   if (cp) {
00135     strcpy(cp, "POTCAR");
00136     potcar = fopen(potcarfile, "r");
00137   }
00138 
00139   /* Read POTCAR file to determine atom types.
00140    * Each atom type section in POTCAR starts with a line
00141    * that contains the name of the element (H, He, C etc.).
00142    * Otherwise try the title line instead.
00143    */
00144   for (atomcount = i = 0; atomcount < data->numatoms; ++i) {
00145     int idx, j;
00146     char const *label;
00147     float mass, radius;
00148 
00149     if (potcar) {
00150        char atomtype[5] = "X";
00151        /* Obtain atom types from POTCAR file */
00152        if (fgets(lineptr, LINESIZE, potcar)) sscanf(lineptr, "%*s %4[^_. 0-9]", atomtype);
00153        idx = get_pte_idx(atomtype);
00154        /* Skip lines in potcar file until next element */
00155        while (fgets(lineptr, LINESIZE, potcar)) if (strstr(lineptr, "End of Dataset")) break;
00156     } else {
00157        /* Try to obtain atom types from title line */
00158        char const *token = (i == 0 ? strtok(data->titleline, " ") : strtok(NULL, " "));
00159        idx = get_pte_idx(token);
00160     }
00161 
00162     label = get_pte_label(idx);
00163     mass = get_pte_mass(idx);
00164     radius = get_pte_vdw_radius(idx);
00165     for (j = 0; j < data->eachatom[i]; ++j, ++atomcount) {
00166       molfile_atom_t *const atom = &(atoms[atomcount]);
00167 
00168       /* Required settings */
00169       strncpy(atom->name, label, sizeof(atom->name));
00170       strncpy(atom->type, atom->name, sizeof(atom->type));
00171       atom->resname[0] = '\0';
00172       atom->resid = 1;
00173       atom->segid[0]='\0';
00174       atom->chain[0]='\0';
00175 
00176 
00177       /* Optional flags (as defined in *optflags) */
00178       atom->mass = mass;
00179       atom->radius = radius;
00180       atom->atomicnumber = idx;
00181     }
00182   }
00183   if (potcar) fclose(potcar);
00184 
00185   if (atomcount != data->numatoms) {
00186     fprintf(stderr, "\n\nVASP POSCAR read) ERROR: file '%s' doesn't seem to have list of atoms.\n", data->filename);
00187     return MOLFILE_ERROR;
00188   }
00189 
00190  /* Ignore header until X,Y,Z-coordinates */
00191  for (i = 0; i < data->version + 3; ++i) fgets(lineptr, LINESIZE, data->file);
00192 
00193  /* Ignore selective tag-line, starting with either 's' or 'S'. */
00194  if (tolower(lineptr[0]) == 's') fgets(lineptr, LINESIZE, data->file);
00195 
00196  /* Check whether all coordinates are present in the file */
00197  for (i = 0; i < data->numatoms; ++i) {
00198    float coord;
00199    fgets(lineptr, LINESIZE, data->file);
00200    if (3 != sscanf(lineptr, "%f %f %f", &coord, &coord, &coord)) {
00201      fprintf(stderr, "\n\nVASP POSCAR read) ERROR: structure is missing type or coordinate(s) in file '%s' for atom '%d'\n", data->filename, i+1);
00202      return MOLFILE_ERROR;
00203    }
00204  }
00205 
00206  rewind(data->file);
00207 
00208  return MOLFILE_SUCCESS;
00209 }
00210 
00211 
00212 static int read_vaspposcar_timestep(void *mydata, int natoms, molfile_timestep_t *ts)
00213 {
00214   int i, direct;
00215   char lineptr[LINESIZE];
00216   float lc;
00217   
00218   vasp_plugindata_t *data = (vasp_plugindata_t *)mydata;
00219 
00220   /* Save coords only if we're given a timestep pointer,
00221    * otherwise assume that VMD wants us to skip past it.
00222    */
00223   if (!ts || !data) return MOLFILE_EOF;
00224 
00225   /* VMD keeps calling for a next timestep, until we reach End-Of-File here */
00226   if (fgets(lineptr, LINESIZE, data->file) == NULL) return MOLFILE_EOF;
00227 
00228   fgets(lineptr, LINESIZE, data->file);
00229   sscanf(lineptr, "%f", &lc);
00230 
00231   for (i = 0; i < 3; ++i) {
00232     float x, y, z;
00233     fgets(lineptr, LINESIZE, data->file);
00234     sscanf(lineptr, "%f %f %f", &x, &y, &z);
00235     data->cell[i][0] = x*lc;
00236     data->cell[i][1] = y*lc;
00237     data->cell[i][2] = z*lc;
00238   }
00239   vasp_buildrotmat(data);
00240 
00241   /* Skip numbers of atom types */
00242   for (i = 0; i < data->version - 2; ++i) fgets(lineptr, LINESIZE, data->file);
00243 
00244   /* Skip selective tag-line, starting with 's' or 'S'. */
00245   if (tolower(lineptr[0]) == 's') fgets(lineptr, LINESIZE, data->file);
00246 
00247   /* Detect direct coordinates tag, starting with 'd' or 'D'. */
00248   direct = (tolower(lineptr[0]) == 'd' ? 1 : 0);
00249 
00250   for (i = 0; i < data->numatoms; ++i) {
00251     float x, y, z, rotx, roty, rotz;
00252     fgets(lineptr, LINESIZE, data->file);
00253     if (3 != sscanf(lineptr, "%f %f %f", &x, &y, &z)) {
00254       fprintf(stderr, "VASP POSCAR read) missing type or coordinate(s) in file '%s' for atom '%d'\n", data->filename, i+1);
00255       return MOLFILE_EOF;
00256     }
00257 
00258     if (direct) {
00259        rotx = x*data->cell[0][0]+y*data->cell[1][0]+z*data->cell[2][0];
00260        roty = x*data->cell[0][1]+y*data->cell[1][1]+z*data->cell[2][1];
00261        rotz = x*data->cell[0][2]+y*data->cell[1][2]+z*data->cell[2][2];
00262     } else {
00263        rotx = x*lc;
00264        roty = y*lc;
00265        rotz = z*lc;
00266     }
00267     ts->coords[3*i  ] = data->rotmat[0][0]*rotx+data->rotmat[0][1]*roty+data->rotmat[0][2]*rotz;
00268     ts->coords[3*i+1] = data->rotmat[1][0]*rotx+data->rotmat[1][1]*roty+data->rotmat[1][2]*rotz;
00269     ts->coords[3*i+2] = data->rotmat[2][0]*rotx+data->rotmat[2][1]*roty+data->rotmat[2][2]*rotz;
00270   }
00271 
00272   vasp_timestep_unitcell(ts, data);
00273 
00274   /* POSCAR type files have only one single timestep.
00275    * Therefore proceed till end of file after reading the coordinates.
00276    */
00277   fseek(data->file, 0, SEEK_END);
00278 
00279   return MOLFILE_SUCCESS;
00280 }
00281 
00282 
00283 static void close_vaspposcar_read(void *mydata)
00284 {
00285   vasp_plugindata_t *data = (vasp_plugindata_t *)mydata;
00286 
00287   vasp_plugindata_free(data);
00288 }
00289 
00290 
00291 static void *open_vaspposcar_write(const char *filename, const char *filetype, int natoms)
00292 {
00293   vasp_plugindata_t *data;
00294 
00295   data = vasp_plugindata_malloc();
00296   if (!data) return NULL;
00297 
00298   data->file = fopen(filename, "w");
00299   if (!data->file) {
00300     vasp_plugindata_free(data);
00301     fprintf(stderr, "VASP POSCAR write) ERROR: Unable to open vaspposcar file '%s' for writing\n", filename);
00302     return NULL;
00303   }
00304 
00305   data->filename = strdup(filename);
00306   data->numatoms = natoms;
00307 
00308   return data;
00309 }
00310 
00311 
00312 static int write_vaspposcar_structure(void *mydata, int optflags, const molfile_atom_t *atoms)
00313 {
00314   vasp_plugindata_t *data = (vasp_plugindata_t *)mydata;
00315 
00316   if (!data || !atoms) return MOLFILE_ERROR;
00317 
00318   data->atomlist = (molfile_atom_t *)malloc(data->numatoms*sizeof(molfile_atom_t));
00319   if (!data->atomlist) return MOLFILE_ERROR;
00320 
00321   memcpy(data->atomlist, atoms, data->numatoms*sizeof(molfile_atom_t));
00322 
00323   return MOLFILE_SUCCESS;
00324 }
00325 
00326 
00327 static int write_vaspposcar_timestep(void *mydata, const molfile_timestep_t *ts)
00328 {
00329   vasp_plugindata_t *data = (vasp_plugindata_t *)mydata; 
00330   molfile_atom_t const *atom;
00331   int i, maxtype, eachatom[MAXATOMTYPES];
00332   float x1, x2, y2, x3, y3, z3;
00333   char tmptype[LINESIZE] = "";
00334 
00335   if (!data || !ts) {
00336     fprintf(stderr, "VASP POSCAR write) ERROR: Wrong input for writing POSCAR file\n");
00337     return MOLFILE_ERROR;
00338   }
00339 
00340   /* restore unit cell information in Cartesian */
00341   x1 = ts->A;
00342   x2 = ts->B*cos(ts->gamma*M_PI/180.0);
00343   y2 = ts->B*sin(ts->gamma*M_PI/180.0);
00344   x3 = ts->C*cos(ts->beta*M_PI/180.0);
00345   y3 = (ts->B*ts->C*cos(ts->alpha*M_PI/180.0)-x2*x3)/y2;
00346   z3 = sqrt(ts->C*ts->C - x3*x3 - y3*y3);
00347 
00348   maxtype = -1;
00349   atom = data->atomlist;
00350   for (i = 0; i < data->numatoms && maxtype < MAXATOMTYPES - 1; ++i, ++atom) {
00351     if (strcmp(tmptype, atom->type) != 0) {
00352       fprintf(data->file, "%-2s  ", atom->type);
00353       eachatom[++maxtype] = 0;
00354     }
00355     eachatom[maxtype]++;
00356     strncpy(tmptype, atom->type, sizeof(atom->type));
00357   }
00358 
00359   fprintf(data->file, "\n%20.12f\n", 1.0);  
00360   fprintf(data->file, "%20.12f  %20.12f  %20.12f\n", x1, 0.0, 0.0);
00361   fprintf(data->file, "%20.12f  %20.12f  %20.12f\n", x2, y2,  0.0);
00362   fprintf(data->file, "%20.12f  %20.12f  %20.12f\n", x3, y3,  z3);
00363 
00364   for (i = 0; i <= maxtype; ++i) fprintf(data->file, " %d ", eachatom[i]);
00365 
00366   fprintf(data->file, "\nDirect\n");
00367 
00368   for (i = 0; i < data->numatoms; ++i) {
00369     float const *pos = ts->coords + 3*i;
00370     /* in direct coordinate */
00371     fprintf(data->file, "%20.12f %20.12f %20.12f \n", pos[0]/x1, -x2*pos[0]/(x1*y2) + pos[1]/y2, (-y2*x3+x2*y3)*pos[0]/(x1*y2*z3) - y3*pos[1]/(y2*z3) + pos[2]/z3);
00372   }
00373 
00374   return MOLFILE_SUCCESS;
00375 }
00376 
00377 
00378 static void close_vaspposcar_write(void *mydata)
00379 {
00380   vasp_plugindata_t *data = (vasp_plugindata_t *)mydata;
00381   vasp_plugindata_free(data);
00382 }
00383 
00384 
00385 /* registration stuff */
00386 static molfile_plugin_t plugin;
00387 
00388 int VMDPLUGIN_init() {
00389   memset(&plugin, 0, sizeof(molfile_plugin_t));
00390   plugin.abiversion = vmdplugin_ABIVERSION;
00391   plugin.type = MOLFILE_PLUGIN_TYPE;
00392   plugin.name = "POSCAR";
00393   plugin.prettyname = "VASP_POSCAR";
00394   plugin.author = "Sung Sakong";
00395   plugin.majorv = 0;
00396   plugin.minorv = 8;
00397   plugin.is_reentrant = VMDPLUGIN_THREADUNSAFE;
00398   plugin.filename_extension = "POSCAR";
00399   plugin.open_file_read = open_vaspposcar_read;
00400   plugin.read_structure = read_vaspposcar_structure;
00401   plugin.read_next_timestep = read_vaspposcar_timestep;
00402   plugin.close_file_read = close_vaspposcar_read;
00403   plugin.open_file_write = open_vaspposcar_write;
00404   plugin.write_structure = write_vaspposcar_structure;
00405   plugin.write_timestep = write_vaspposcar_timestep;
00406   plugin.close_file_write = close_vaspposcar_write;
00407   return VMDPLUGIN_SUCCESS;
00408 }
00409 
00410 int VMDPLUGIN_register(void *v, vmdplugin_register_cb cb) {
00411   (*cb)(v, (vmdplugin_t *)&plugin);
00412   return VMDPLUGIN_SUCCESS;
00413 }
00414 
00415 int VMDPLUGIN_fini() {
00416   return VMDPLUGIN_SUCCESS;
00417 }

Generated on Thu Apr 25 03:07:44 2024 for VMD Plugins (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002