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

Measure.C

Go to the documentation of this file.
00001 /***************************************************************************
00002  *cr
00003  *cr            (C) Copyright 1995-2011 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: Measure.C,v $
00013  *      $Author: johns $        $Locker:  $             $State: Exp $
00014  *      $Revision: 1.134 $       $Date: 2011/10/20 22:46:19 $
00015  *
00016  ***************************************************************************
00017  * DESCRIPTION:
00018  *   Code to measure atom distances, angles, dihedrals, etc.
00019  ***************************************************************************/
00020 
00021 #include <stdio.h>
00022 #include <stdlib.h>
00023 #include <math.h>
00024 #include "Measure.h"
00025 #include "AtomSel.h"
00026 #include "Matrix4.h"
00027 #include "utilities.h"
00028 #include "fitrms.h"
00029 #include "ResizeArray.h"
00030 #include "SpatialSearch.h"  // for find_within()
00031 #include "MoleculeList.h"
00032 #include "Inform.h"
00033 #include "Timestep.h"
00034 #include "VMDApp.h"
00035 #include "WKFThreads.h"
00036 #include "WKFUtils.h"
00037 
00038 // Standard functions available to everyone
00039 static const char *measure_error_messages[] = {
00040   "no atom selection",                                  // -1
00041   "no atoms in selection",                              // -2
00042   "incorrect number of weights for selection",          // -3
00043   "internal error: NULL pointer given",                 // -4
00044   "bad weight sum, would cause divide by zero",         // -5
00045   "molecule was deleted(?)",                            // -6
00046   "cannot understand weight parameter",                 // -7
00047   "non-number given as parameter",                      // -8
00048   "two selections don't have the same number of atoms", // -9
00049   "internal error: out of range",                       // -10
00050   "no coordinates in selection",                        // -11
00051   "couldn't compute eigenvalue/vectors",                // -12
00052   "unknown Tcl problem",                                // -13
00053   "no atom radii",                                      // -14
00054   "order parameter contains out-of-range atom index",   // -15
00055   "order parameter not supported in old RMS fit mode",  // -16
00056   "specified frames are out of range",                  // -17
00057   "both selections must reference the same molecule",   // -18
00058   "one atom appears twice in list",                     // -19
00059   "molecule contains no frames",                        // -20
00060   "invalid atom id",                                    // -21
00061   "cutoff must be smaller than cell dimension",         // -22
00062   "Zero volmap gridsize"                                // -23
00063 };
00064   
00065 const char *measure_error(int errnum) {
00066   if (errnum >= 0 || errnum < -23) 
00067     return "bad error number";
00068   return measure_error_messages[-errnum - 1];
00069 }
00070 
00071 int measure_move(const AtomSel *sel, float *framepos, const Matrix4 &mat) {
00072   if (!sel)       return MEASURE_ERR_NOSEL;
00073   if (!framepos)  return MEASURE_ERR_NOFRAMEPOS;
00074 
00075   // and apply it to the coordinates
00076   int i;
00077   float *pos = framepos;
00078   if (mat.mat[3]==0 && mat.mat[7]==0 && mat.mat[11]==0 && mat.mat[15]==1) {
00079     // then this is just a rotation followed by a translation.  Do it
00080     // the fast way.
00081     const float dx = mat.mat[12];
00082     const float dy = mat.mat[13];
00083     const float dz = mat.mat[14];
00084     const float *m = mat.mat;
00085     if (sel->selected == sel->num_atoms) {
00086       for (i=0; i<sel->num_atoms; i++) {
00087         float x = pos[0]*m[0] + pos[1]*m[4] + pos[2]*m[8] + dx;
00088         float y = pos[0]*m[1] + pos[1]*m[5] + pos[2]*m[9] + dy;
00089         float z = pos[0]*m[2] + pos[1]*m[6] + pos[2]*m[10] + dz;
00090         pos[0] = x;
00091         pos[1] = y;
00092         pos[2] = z;
00093         pos += 3;
00094       }
00095     } else {
00096       pos += sel->firstsel*3;
00097       for (i=sel->firstsel; i<=sel->lastsel; i++) {
00098         if (sel->on[i]) {
00099           float x = pos[0]*m[0] + pos[1]*m[4] + pos[2]*m[8] + dx;
00100           float y = pos[0]*m[1] + pos[1]*m[5] + pos[2]*m[9] + dy;
00101           float z = pos[0]*m[2] + pos[1]*m[6] + pos[2]*m[10] + dz;
00102           pos[0] = x;
00103           pos[1] = y;
00104           pos[2] = z;
00105         }
00106         pos += 3;
00107       }
00108     }
00109   } else {
00110     pos += sel->firstsel*3;
00111     for (i=sel->firstsel; i<=sel->lastsel; i++) {
00112       if (sel->on[i]) {
00113         mat.multpoint3d(pos, pos);
00114       }
00115       pos += 3;
00116     }
00117   }
00118   return MEASURE_NOERR;
00119 }
00120 
00121 // compute the sum of a set of weights which correspond either to
00122 // all atoms, or to selected atoms
00123 int measure_sumweights(const AtomSel *sel, int numweights, 
00124                        const float *weights, float *weightsum) {
00125   if (!sel)       return MEASURE_ERR_NOSEL;
00126   if (!weights)   return MEASURE_ERR_NOWEIGHT;
00127   if (!weightsum) return MEASURE_ERR_GENERAL;
00128 
00129   double sum = 0;
00130 
00131   if (numweights == sel->num_atoms) {
00132     int i;
00133     for (i=sel->firstsel; i<=sel->lastsel; i++) {
00134       if (sel->on[i]) {
00135         sum += weights[i];
00136       }
00137     }
00138   } else if (numweights == sel->selected) {
00139     int i, j;
00140     for (j=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00141       if (sel->on[i]) {
00142         sum += weights[j];
00143         j++;
00144       }
00145     }
00146   } else {
00147     return MEASURE_ERR_BADWEIGHTPARM;  
00148   }
00149 
00150   *weightsum = (float) sum;
00151   return MEASURE_NOERR;
00152 }
00153 
00154 
00155 // compute the center of mass
00156 // return -5 if total weight == 0, otherwise 0 for success.
00157 int measure_center(const AtomSel *sel, const float *framepos,
00158                    const float *weight, float *com) {
00159   if (!sel)      return MEASURE_ERR_NOSEL;
00160   if (!framepos) return MEASURE_ERR_NOFRAMEPOS;
00161   if (!weight)   return MEASURE_ERR_NOWEIGHT;
00162   if (!com)      return MEASURE_ERR_NOCOM;   
00163 
00164   int i, j;
00165   float x, y, z, w;
00166   j=0;
00167   w=x=y=z=0;
00168   for (i=sel->firstsel; i<=sel->lastsel; i++) {
00169     if (sel->on[i]) {
00170       float tw = weight[j];
00171       w += tw;
00172       x += tw * framepos[3*i  ];
00173       y += tw * framepos[3*i+1];
00174       z += tw * framepos[3*i+2];
00175       j++;
00176     }
00177   }
00178 
00179   if (w == 0) {
00180     return MEASURE_ERR_BADWEIGHTSUM;
00181   }
00182 
00183   com[0] = x / w;
00184   com[1] = y / w;
00185   com[2] = z / w;
00186 
00187   return MEASURE_NOERR;
00188 }
00189 
00190 
00191 // compute the axis aligned aligned bounding box for the selected atoms
00192 int measure_minmax(int num, const int *on, const float *framepos, 
00193                    const float *radii, float *min_coord, float *max_coord) {
00194   int i;
00195   float minx, miny, minz;
00196   float maxx, maxy, maxz;
00197 
00198   if (!on)                      return MEASURE_ERR_NOSEL;
00199   if (num == 0)                 return MEASURE_ERR_NOATOMS;
00200   if (!min_coord || !max_coord) return MEASURE_ERR_NOMINMAXCOORDS;
00201 
00202   vec_zero(min_coord);
00203   vec_zero(max_coord);
00204 
00205   // find first selected atom
00206   for (i=0; i<num; i++) if (on[i]) break;
00207 
00208   // if no atoms are selected return { 0 0 0 } for min and max
00209   if (i==num) return MEASURE_NOERR;
00210 
00211   // initialize minmax coords to the first selected atom
00212   if (radii == NULL) {
00213     // calculate bounding box of atom centers
00214     minx = maxx = framepos[i*3  ];
00215     miny = maxy = framepos[i*3+1];
00216     minz = maxz = framepos[i*3+2];
00217   } else {
00218     // calculate bounding box for atoms /w given radii 
00219     minx = framepos[i*3  ] - radii[i];
00220     maxx = framepos[i*3  ] + radii[i];
00221     miny = framepos[i*3+1] - radii[i];
00222     maxy = framepos[i*3+1] + radii[i];
00223     minz = framepos[i*3+2] - radii[i];
00224     maxz = framepos[i*3+2] + radii[i];
00225   }
00226 
00227   // continue looping from there until we finish
00228   if (radii == NULL) {
00229     // calculate bounding box of atom centers
00230     for (i++; i<num; i++) {
00231       if (on[i]) {
00232         int ind = i * 3;
00233         float tmpx = framepos[ind  ];
00234         if (tmpx < minx) minx = tmpx; 
00235         if (tmpx > maxx) maxx = tmpx;
00236   
00237         float tmpy = framepos[ind+1];
00238         if (tmpy < miny) miny = tmpy; 
00239         if (tmpy > maxy) maxy = tmpy;
00240   
00241         float tmpz = framepos[ind+2];
00242         if (tmpz < minz) minz = tmpz; 
00243         if (tmpz > maxz) maxz = tmpz;
00244       }
00245     }
00246   } else {
00247     // calculate bounding box for atoms /w given radii 
00248     for (i++; i<num; i++) {
00249       if (on[i]) {
00250         int ind = i * 3;
00251         float mintmpx = framepos[ind  ] - radii[i];
00252         float maxtmpx = framepos[ind  ] + radii[i];
00253         if (mintmpx < minx) minx = mintmpx;
00254         if (maxtmpx > maxx) maxx = maxtmpx;
00255   
00256         float mintmpy = framepos[ind+1] - radii[i];
00257         float maxtmpy = framepos[ind+1] + radii[i];
00258         if (mintmpy < miny) miny = mintmpy; 
00259         if (maxtmpy > maxy) maxy = maxtmpy;
00260   
00261         float mintmpz = framepos[ind+2] - radii[i];
00262         float maxtmpz = framepos[ind+2] + radii[i];
00263         if (mintmpz < minz) minz = mintmpz; 
00264         if (maxtmpz > maxz) maxz = maxtmpz;
00265       }
00266     }
00267   }
00268 
00269   // set the final min/max output values
00270   min_coord[0]=minx;
00271   min_coord[1]=miny;
00272   min_coord[2]=minz;
00273   max_coord[0]=maxx;
00274   max_coord[1]=maxy;
00275   max_coord[2]=maxz;
00276 
00277   return MEASURE_NOERR;
00278 }
00279 
00280 
00281 // Calculate average position of selected atoms over selected frames
00282 extern int measure_avpos(const AtomSel *sel, MoleculeList *mlist, 
00283                          int start, int end, int step, float *avpos) {
00284   if (!sel)                     return MEASURE_ERR_NOSEL;
00285   if (sel->num_atoms == 0)      return MEASURE_ERR_NOATOMS;
00286 
00287   Molecule *mymol = mlist->mol_from_id(sel->molid());
00288   int maxframes = mymol->numframes();
00289   
00290   // accept value of -1 meaning "all" frames
00291   if (end == -1)
00292     end = maxframes-1;
00293 
00294   if (maxframes == 0 || start < 0 || start > end || 
00295       end >= maxframes || step <= 0)
00296     return MEASURE_ERR_BADFRAMERANGE;
00297 
00298   int i;
00299   for (i=0; i<(3*sel->selected); i++)
00300     avpos[i] = 0.0f;
00301 
00302   int frame, avcount, j;
00303   for (avcount=0,frame=start; frame<=end; avcount++,frame+=step) {
00304     const float *framepos = (mymol->get_frame(frame))->pos;
00305     for (j=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00306       if (sel->on[i]) {
00307         avpos[j*3    ] += framepos[i*3    ];
00308         avpos[j*3 + 1] += framepos[i*3 + 1];
00309         avpos[j*3 + 2] += framepos[i*3 + 2];
00310         j++;
00311       }
00312     } 
00313   }
00314 
00315   float avinv = 1.0f / (float) avcount;
00316   for (j=0; j<(3*sel->selected); j++) {
00317     avpos[j] *= avinv;
00318   } 
00319 
00320   return MEASURE_NOERR;
00321 }
00322 
00323 
00324 // Calculate dipole moment for selected atoms
00325 extern int measure_dipole(const AtomSel *sel, MoleculeList *mlist, 
00326                           float *dipole, int unitsdebye, int usecenter) {
00327   if (!sel)                     return MEASURE_ERR_NOSEL;
00328   if (sel->num_atoms == 0)      return MEASURE_ERR_NOATOMS;
00329 
00330   Molecule *mymol = mlist->mol_from_id(sel->molid());
00331   double  rvec[3] = {0, 0, 0};
00332   double qrvec[3] = {0, 0, 0};
00333   double mrvec[3] = {0, 0, 0};
00334   double totalq = 0.0;
00335   double totalm = 0.0;
00336   int i;
00337 
00338   // get atom coordinates
00339   const float *framepos = sel->coordinates(mlist);
00340 
00341   // get atom charges
00342   const float *q = mymol->charge();
00343   const float *m = mymol->mass();
00344 
00345   for (i=sel->firstsel; i<=sel->lastsel; i++) {
00346     if (sel->on[i]) {
00347       int ind = i * 3;
00348       rvec[0] += framepos[ind    ];
00349       rvec[1] += framepos[ind + 1];
00350       rvec[2] += framepos[ind + 2];
00351 
00352       qrvec[0] += framepos[ind    ] * q[i];
00353       qrvec[1] += framepos[ind + 1] * q[i];
00354       qrvec[2] += framepos[ind + 2] * q[i];
00355 
00356       mrvec[0] += framepos[ind    ] * m[i];
00357       mrvec[1] += framepos[ind + 1] * m[i];
00358       mrvec[2] += framepos[ind + 2] * m[i];
00359 
00360       totalq += q[i];
00361       totalm += m[i];
00362     }
00363   }
00364 
00365   // fall back to geometrical center when bad or no masses
00366   if (totalm < 0.0001)
00367     usecenter=1;
00368   
00369   switch (usecenter) {
00370     case 1:
00371     {
00372         double rscale = totalq / sel->selected; 
00373         dipole[0] = (float) (qrvec[0] - (rvec[0] * rscale)); 
00374         dipole[1] = (float) (qrvec[1] - (rvec[1] * rscale)); 
00375         dipole[2] = (float) (qrvec[2] - (rvec[2] * rscale)); 
00376         break;
00377     }
00378 
00379     case -1: 
00380     {
00381         double mscale = totalq / totalm; 
00382         dipole[0] = (float) (qrvec[0] - (mrvec[0] * mscale)); 
00383         dipole[1] = (float) (qrvec[1] - (mrvec[1] * mscale)); 
00384         dipole[2] = (float) (qrvec[2] - (mrvec[2] * mscale)); 
00385         break;
00386     }
00387     
00388     case 0: // fallthrough
00389     default: 
00390     {
00391         dipole[0] = (float) qrvec[0];
00392         dipole[1] = (float) qrvec[1];
00393         dipole[2] = (float) qrvec[2];
00394         break;
00395     }
00396   }
00397 
00398   // According to the values in
00399   // http://www.physics.nist.gov/cuu/Constants/index.html
00400   // 1 e*A = 4.80320440079 D
00401   // 1 D = 1E-18 Fr*cm = 0.208194346224 e*A
00402   // 1 e*A = 1.60217653 E-29 C*m
00403   // 1 C*m = 6.24150947961 E+28 e*A
00404   // 1 e*A = 1.88972613458 e*a0
00405   // 1 e*a0 = 0.529177208115 e*A
00406 
00407   if (unitsdebye) {
00408     // 1 e*A = 4.80320440079 D 
00409     // latest CODATA (2006) gives:
00410     //         4.80320425132073913031
00411     dipole[0] *= 4.80320425132f;
00412     dipole[1] *= 4.80320425132f;
00413     dipole[2] *= 4.80320425132f;
00414   }
00415  
00416   return MEASURE_NOERR;
00417 }
00418 
00419 
00420 // Calculate RMS fluctuation of selected atoms over selected frames
00421 extern int measure_rmsf(const AtomSel *sel, MoleculeList *mlist, 
00422                         int start, int end, int step, float *rmsf) {
00423   if (!sel)                     return MEASURE_ERR_NOSEL;
00424   if (sel->num_atoms == 0)      return MEASURE_ERR_NOATOMS;
00425 
00426   Molecule *mymol = mlist->mol_from_id(sel->molid());
00427   int maxframes = mymol->numframes();
00428 
00429   // accept value of -1 meaning "all" frames
00430   if (end == -1)
00431     end = maxframes-1;
00432 
00433   if (maxframes == 0 || start < 0 || start > end ||
00434       end >= maxframes || step <= 0)
00435     return MEASURE_ERR_BADFRAMERANGE;
00436 
00437   int i;
00438   for (i=0; i<sel->selected; i++)
00439     rmsf[i] = 0.0f;
00440 
00441   int rc; 
00442   float *avpos = new float[3*sel->selected];
00443   rc = measure_avpos(sel, mlist, start, end, step, avpos);
00444 
00445   if (rc != MEASURE_NOERR) {
00446     delete [] avpos;
00447     return rc;
00448   }
00449 
00450   // calculate per-atom variance here
00451   int frame, avcount, j;
00452   for (avcount=0,frame=start; frame<=end; avcount++,frame+=step) {
00453     const float *framepos = (mymol->get_frame(frame))->pos;
00454     for (j=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00455       if (sel->on[i]) {
00456         rmsf[j] += distance2(&avpos[3*j], &framepos[3*i]);
00457         j++;
00458       }
00459     }
00460   }
00461 
00462   float avinv = 1.0f / (float) avcount;
00463   for (j=0; j<sel->selected; j++) {
00464     rmsf[j] = sqrtf(rmsf[j] * avinv);
00465   }
00466 
00467   delete [] avpos;
00468   return MEASURE_NOERR;
00469 }
00470 
00471 
00473 //  rgyr := sqrt(sum (mass(n) ( r(n) - r(com) )^2)/sum(mass(n)))
00474 //  The return value, a float, is put in 'float *rgyr'
00475 //  The function return value is 0 if ok, <0 if not
00476 int measure_rgyr(const AtomSel *sel, MoleculeList *mlist, const float *weight, 
00477                  float *rgyr) {
00478   if (!rgyr) return MEASURE_ERR_NORGYR;
00479   if (!sel) return MEASURE_ERR_NOSEL;
00480   if (!mlist) return MEASURE_ERR_GENERAL;
00481 
00482   const float *framepos = sel->coordinates(mlist);
00483 
00484   // compute the center of mass with the current weights
00485   float com[3];
00486   int ret_val = measure_center(sel, framepos, weight, com);
00487   if (ret_val < 0) 
00488     return ret_val;
00489   
00490   // measure center of gyration
00491   int i, j;
00492   float total_w, sum;
00493 
00494   total_w=sum=0;
00495   for (j=0,i=sel->firstsel; i<=sel->lastsel; i++) {
00496     if (sel->on[i]) {
00497       float w = weight[j];
00498       total_w += w;
00499       sum += w * distance2(framepos + 3*i, com);
00500       j++;
00501     } 
00502   }
00503 
00504   if (total_w == 0) {
00505     return MEASURE_ERR_BADWEIGHTSUM;
00506   }
00507 
00508   // and finalize the computation
00509   *rgyr = sqrtf(sum/total_w);
00510 
00511   return MEASURE_NOERR;
00512 }
00513 
00514 
00516 //  1) if num == sel.selected ; assumes there is one weight per 
00517 //           selected atom
00518 //  2) if num == sel.num_atoms; assumes weight[i] is for atom[i]
00519 //  returns 0 and value in rmsd if good
00520 //   return < 0 if invalid
00521 //  Function is::=  rmsd = 
00522 //    sqrt(sum(weight(n) * sqr(r1(i(n))-r2(i(n))))/sum(weight(n)) / N
00523 int measure_rmsd(const AtomSel *sel1, const AtomSel *sel2,
00524                  int num, const float *framepos1, const float *framepos2,
00525                  float *weight, float *rmsd) {
00526   if (!sel1 || !sel2)   return MEASURE_ERR_NOSEL;
00527   if (sel1->selected < 1 || sel2->selected < 1) return MEASURE_ERR_NOSEL;
00528   if (!weight || !rmsd) return MEASURE_ERR_NOWEIGHT;
00529 
00530   // the number of selected atoms must be the same
00531   if (sel1->selected != sel2->selected) return MEASURE_ERR_MISMATCHEDCNT;
00532 
00533   // need to know how to traverse the list of weights
00534   // there could be 1 weight per atom (sel_flg == 1) or 
00535   // 1 weight per selected atom (sel_flg == 0)
00536   int sel_flg;
00537   if (num == sel1->num_atoms) {
00538     sel_flg = 1; // using all elements
00539   } else {
00540     sel_flg = 0; // using elements from selection
00541   }
00542 
00543   // temporary variables
00544   float tmp_w;
00545   int w_index = 0;              // the term in the weight field to use
00546   int sel1ind = sel1->firstsel; // start from the first selected atom
00547   int sel2ind = sel2->firstsel; // start from the first selected atom
00548   float wsum = 0;
00549   float rmsdsum = 0;
00550 
00551   *rmsd = 10000000; // if we bail out, return a huge value
00552 
00553   // compute the rmsd
00554   int count = sel1->selected;
00555   while (count-- > 0) {
00556     // find next 'on' atom in sel1 and sel2
00557     // loop is safe since we already stop the on count > 0 above
00558     while (!sel1->on[sel1ind])
00559       sel1ind++;
00560     while (!sel2->on[sel2ind])
00561       sel2ind++;
00562 
00563     // the weight offset to use depends on how many terms there are
00564     if (sel_flg == 0) {
00565       tmp_w = weight[w_index++];
00566     } else {
00567       tmp_w = weight[sel1ind]; // use the first selection for the weights
00568     }
00569 
00570     // sum the calculated rmsd and weight values
00571     rmsdsum += tmp_w * distance2(framepos1 + 3*sel1ind, framepos2 + 3*sel2ind);
00572     wsum += tmp_w;
00573 
00574     // and advance to the next atom pair
00575     sel1ind++;
00576     sel2ind++;
00577   }
00578 
00579   // check weight sum
00580   if (wsum == 0) {
00581     return MEASURE_ERR_BADWEIGHTSUM;
00582   }
00583 
00584   // finish the rmsd calcs
00585   *rmsd = sqrtf(rmsdsum / wsum);
00586 
00587   return MEASURE_NOERR; // and say rmsd is OK
00588 }
00589 
00590 /* jacobi.C, taken from Numerical Recipes and specialized to 3x3 case */
00591 
00592 #define ROTATE(a,i,j,k,l) g=a[i][j];h=a[k][l];a[i][j]=g-s*(h+g*tau);\
00593         a[k][l]=h+s*(g-h*tau);
00594 
00595 static int jacobi(float a[4][4], float d[3], float v[3][3])
00596 {
00597   int n=3;
00598   int j,iq,ip,i;
00599   float tresh,theta,tau,t,sm,s,h,g,c,*b,*z;
00600   
00601   b=new float[n];
00602   z=new float[n];
00603   for (ip=0;ip<n;ip++) {
00604     for (iq=0;iq<n;iq++) v[ip][iq]=0.0;
00605     v[ip][ip]=1.0;
00606   }
00607   for (ip=0;ip<n;ip++) {
00608     b[ip]=d[ip]=a[ip][ip];
00609     z[ip]=0.0;
00610   }
00611   for (i=1;i<=50;i++) {
00612     sm=0.0;
00613     for (ip=0;ip<n-1;ip++) {
00614       for (iq=ip+1;iq<n;iq++)
00615         sm += (float) fabs(a[ip][iq]);
00616     }
00617     if (sm == 0.0) {
00618       delete [] z;
00619       delete [] b;
00620       return 0; // Exit normally
00621     }
00622     if (i < 4)
00623       tresh=0.2f*sm/(n*n);
00624     else
00625       tresh=0.0f;
00626     for (ip=0;ip<n-1;ip++) {
00627       for (iq=ip+1;iq<n;iq++) {
00628         g=100.0f * fabsf(a[ip][iq]);
00629         if (i > 4 && (float)(fabs(d[ip])+g) == (float)fabs(d[ip])
00630             && (float)(fabs(d[iq])+g) == (float)fabs(d[iq]))
00631           a[ip][iq]=0.0;
00632         else if (fabs(a[ip][iq]) > tresh) {
00633           h=d[iq]-d[ip];
00634           if ((float)(fabs(h)+g) == (float)fabs(h))
00635             t=(a[ip][iq])/h;
00636           else {
00637             theta=0.5f*h/(a[ip][iq]);
00638             t=1.0f/(fabsf(theta)+sqrtf(1.0f+theta*theta));
00639             if (theta < 0.0f) t = -t;
00640           }
00641           c=1.0f/sqrtf(1+t*t);
00642           s=t*c;
00643           tau=s/(1.0f+c);
00644           h=t*a[ip][iq];
00645           z[ip] -= h;
00646           z[iq] += h;
00647           d[ip] -= h;
00648           d[iq] += h;
00649           a[ip][iq]=0.0;
00650           for (j=0;j<=ip-1;j++) {
00651             ROTATE(a,j,ip,j,iq)
00652               }
00653           for (j=ip+1;j<=iq-1;j++) {
00654             ROTATE(a,ip,j,j,iq)
00655               }
00656           for (j=iq+1;j<n;j++) {
00657             ROTATE(a,ip,j,iq,j)
00658               }
00659           for (j=0;j<n;j++) {
00660             ROTATE(v,j,ip,j,iq)
00661               }
00662         }
00663       }
00664     }
00665     for (ip=0;ip<n;ip++) {
00666       b[ip] += z[ip];
00667       d[ip]=b[ip];
00668       z[ip]=0.0;
00669     }
00670   }
00671   delete [] b;
00672   delete [] z;
00673   return 1; // Failed to converge
00674 }
00675 #undef ROTATE
00676 
00677 static int transvecinv(const double *v, Matrix4 &res) {
00678   double x, y, z;
00679   x=v[0];
00680   y=v[1];
00681   z=v[2];
00682   if (x == 0.0 && y == 0.0) {
00683     if (z == 0.0) {
00684       return -1;
00685     }
00686     if (z > 0)
00687       res.rot(90, 'y');
00688     else
00689       res.rot(-90, 'y');
00690     return 0;
00691   }
00692   double theta = atan2(x,y);
00693   double length = sqrt(x*x + y*y);
00694   double phi = atan2(z,length);
00695   res.rot((float) RADTODEG(theta), 'z');
00696   res.rot((float) RADTODEG(-phi),  'y');
00697   return 0;
00698 }
00699  
00700 static int transvec(const double *v, Matrix4 &res) {
00701   double x, y, z;
00702   x=v[0];
00703   y=v[1];
00704   z=v[2];
00705   if (x == 0.0 && y == 0.0) {
00706     if (z == 0.0) {
00707       return -1;
00708     }
00709     if (z > 0)
00710       res.rot(-90, 'y');
00711     else
00712       res.rot(90, 'y');
00713     return 0;
00714   }
00715   double theta = atan2(x,y);
00716   double length = sqrt(x*x + y*y);
00717   double phi = atan2(z,length);
00718   res.rot((float) RADTODEG(phi),  'y');
00719   res.rot((float) RADTODEG(-theta), 'z');
00720   return 0;
00721 }
00722 
00723 static Matrix4 myfit2(const float *x, const float *y, 
00724                const float *comx, const float *comy) {
00725 
00726   Matrix4 res;
00727   double dx[3], dy[3];
00728   dx[0] = x[0] - comx[0];
00729   dx[1] = x[1] - comx[1];
00730   dx[2] = x[2] - comx[2];
00731   dy[0] = y[0] - comy[0];
00732   dy[1] = y[1] - comy[1];
00733   dy[2] = y[2] - comy[2];
00734 
00735   res.translate(-comx[0], -comx[1], -comx[2]);
00736   transvecinv(dx, res); 
00737   transvec(dy, res);
00738   res.translate(comy[0], comy[1], comy[2]); 
00739   return res;
00740 }
00741 
00742 static Matrix4 myfit3(const float *x1, const float *x2, 
00743                       const float *y1, const float *y2,
00744                       const float *comx, const float *comy) {
00745    
00746   Matrix4 mx, my, rx1, ry1;
00747   double dx1[3], dy1[3], angle;
00748   float dx2[3], dy2[3], x2t[3], y2t[3];
00749 
00750   for (int i=0; i<3; i++) {
00751     dx1[i] = x1[i] - comx[i];
00752     dx2[i] = x2[i] - comx[i];
00753     dy1[i] = y1[i] - comy[i];
00754     dy2[i] = y2[i] - comy[i];
00755   }
00756 
00757   transvecinv(dx1, rx1);
00758   rx1.multpoint3d(dx2, x2t);
00759   angle = atan2(x2t[2], x2t[1]);
00760   mx.translate(-comx[0], -comx[1], -comy[2]);
00761   mx.multmatrix(rx1);
00762   mx.rot((float) RADTODEG(angle), 'x'); 
00763 
00764   transvecinv(dy1, ry1);
00765   ry1.multpoint3d(dy2, y2t);
00766   angle = atan2(y2t[2], y2t[1]);
00767   my.translate(-comy[0], -comy[1], -comy[2]);
00768   my.multmatrix(ry1);
00769   my.rot((float) RADTODEG(angle), 'x');
00770   
00771   my.inverse();
00772   mx.multmatrix(my);
00773   return mx;
00774 }
00775 
00776 // find the best fit alignment to take the first structure into the second
00777 // Put the result in the matrix 'mat'
00778 //  This algorithm comes from Kabsch, Acta Cryst. (1978) A34, 827-828.
00779 // Need the 2nd weight for the com calculation
00780 int measure_fit(const AtomSel *sel1, const AtomSel *sel2, const float *x, 
00781                 const float *y, const float *weight, 
00782                 const int *order, Matrix4 *mat) {
00783   float comx[3];
00784   float comy[3];
00785   int num = sel1->selected;
00786   int ret_val;
00787   ret_val = measure_center(sel1, x, weight, comx);
00788   if (ret_val < 0) {
00789     return ret_val;
00790   }
00791   ret_val = measure_center(sel2, y, weight, comy);
00792   if (ret_val < 0) {
00793     return ret_val;
00794   }
00795 
00796   // the Kabsch method won't work of the number of atoms is less than 4
00797   // (and won't work in some cases of n > 4; I think it works so long as
00798   // three or more planes are needed to intersect all the data points
00799   switch (sel1->selected) {
00800   case 1: { // simple center of mass alignment
00801     Matrix4 tmp;
00802     tmp.translate(-comx[0], -comx[1], -comx[2]);
00803     tmp.translate(comy[0], comy[1], comy[2]);
00804     memcpy(mat->mat, tmp.mat, 16*sizeof(float));
00805     return MEASURE_NOERR;
00806   }
00807   case 3:
00808   case 2: { 
00809     // find the first (n-1) points (from each molecule)
00810     int pts[6], count = 0;
00811     int n;
00812     for (n=sel1->firstsel; n<=sel1->lastsel; n++) {
00813       if (sel1->on[n]) {
00814         pts[count++] = n;
00815         if (sel1->selected == 2) {
00816           count++;                   // will put y data in pts[3]
00817           break;
00818         }
00819         if (count == 2) break;
00820       }
00821     }
00822     for (n=sel2->firstsel; n<=sel2->lastsel; n++) {
00823       if (sel2->on[n]) {
00824         pts[count++] = n;
00825         if (sel1->selected == 2) {
00826           count++;
00827           break;
00828         }
00829         if (count == 4) break;
00830       }
00831     }
00832     if (count != 4) {
00833       return MEASURE_ERR_MISMATCHEDCNT;
00834     }
00835     
00836     // reorder the sel2 atoms according to the order parameter
00837     if (order != NULL) {
00838       int i; 
00839       int tmp[6];
00840       memcpy(tmp, pts, sizeof(pts));
00841       for (i=0; i<num; i++) {
00842         pts[i + num] = tmp[num + order[i]]; // order indices are 0-based
00843       } 
00844     }
00845 
00846     if (sel1->selected == 2) {
00847       *mat = myfit2(x+3*pts[0],y+3*pts[2], comx, comy);
00848       ret_val = 0;
00849     } else {
00850       *mat = myfit3(x+3*pts[0],x+3*pts[1],y+3*pts[2],y+3*pts[3],comx,comy);
00851       ret_val = 0;
00852     }  
00853     if (ret_val != 0) {
00854       return MEASURE_ERR_GENERAL;
00855     }
00856 
00857     return 0;
00858   }
00859   default:
00860     break;
00861   }
00862   // at this point I know all the data values are good
00863  
00864 
00865   // use the new RMS fit implementation by default unless told otherwise 
00866   char *opt = getenv("VMDRMSFITMETHOD");
00867   if (!opt || strcmp(opt, "oldvmd")) {
00868     int i, k;
00869     float *v1, *v2;
00870     v1 = new float[3*num];
00871     v2 = new float[3*num];
00872     for (k=0, i=sel1->firstsel; i<=sel1->lastsel; i++) {
00873       if (sel1->on[i]) {
00874         int ind = 3 * i;
00875         v1[k++] = x[ind    ];
00876         v1[k++] = x[ind + 1];
00877         v1[k++] = x[ind + 2];
00878       }
00879     }
00880     for (k=0, i=sel2->firstsel; i<=sel2->lastsel; i++) {
00881       if (sel2->on[i]) {
00882         int ind = 3 * i;
00883         v2[k++] = y[ind    ];
00884         v2[k++] = y[ind + 1];
00885         v2[k++] = y[ind + 2];
00886       }
00887     }
00888 
00889     // reorder the sel2 atoms according to the order parameter
00890     if (order != NULL) {
00891       int i; 
00892       float *tmp = new float[3*num];
00893       memcpy(tmp, v2, 3*num*sizeof(float));
00894       for (i=0; i<num; i++) {
00895         int ind = 3 * i;
00896         int idx = 3 * order[i]; // order indices are 0-based
00897         v2[ind    ] = tmp[idx    ];
00898         v2[ind + 1] = tmp[idx + 1];
00899         v2[ind + 2] = tmp[idx + 2];
00900       } 
00901       delete[] tmp;
00902     }
00903 
00904     float tmp[16];
00905     // MatrixFitRMS returns RMS distance of fitted molecule.  Would be nice
00906     // to return this information to the user since it would make computing
00907     // the fitted RMSD much faster (no need to get the matrix, apply the
00908     // transformation, recompute RMS).
00909     MatrixFitRMS(num, v1, v2, weight, tmp);
00910 
00911     delete [] v1;
00912     delete [] v2;
00913     //fprintf(stderr, "got err %f\n", err);
00914     // output from FitRMS is a 3x3 matrix, plus a pre-translation stored in
00915     // row 3, and a post-translation stored in column 3.
00916     float pre[3], post[3];
00917     for (i=0; i<3; i++) {
00918       post[i] = tmp[4*i+3];
00919       tmp[4*i+3] = 0;
00920     }
00921     for (i=0; i<3; i++) {
00922       pre[i] = tmp[12+i];
00923       tmp[12+i] = 0;
00924     }
00925     Matrix4 result;
00926     result.translate(pre);
00927     result.multmatrix(Matrix4(tmp));
00928     result.translate(post);
00929     memcpy(mat->mat, result.mat, 16*sizeof(float));
00930     return 0;
00931   }
00932 
00933   // the old RMS fit code doesn't support reordering of sel2 currently
00934   if (order != NULL) {
00935     return MEASURE_ERR_NOTSUP;
00936   }
00937 
00938   // a) compute R = r(i,j) = sum( w(n) * (y(n,i)-comy(i)) * (x(n,j)-comx(j)))
00939   Matrix4 R;
00940   int i,j;
00941   float scale = (float) num * num;
00942   for (i=0; i<3; i++) {
00943     for (j=0; j<3; j++) {
00944       float tmp = 0;
00945       int nx = 0, ny = 0, k = 0; 
00946       while (nx < sel1->num_atoms && ny < sel2->num_atoms) { 
00947         if (!sel1->on[nx]) {
00948           nx++;
00949           continue;
00950         }
00951         if (!sel2->on[ny]) {
00952           ny++;
00953           continue;
00954         }
00955 
00956         // found both, so get data
00957         
00958         tmp += weight[k] * (y[3*ny+i] - comy[i]) * (x[3*nx+j] - comx[j]) /
00959           scale;
00960         nx++;
00961         ny++;
00962         k++;
00963       }
00964       R.mat[4*i+j] = tmp;
00965     }
00966   }
00967 
00968   // b) 1) form R~R
00969   Matrix4 Rt;
00970   for (i=0; i<3; i++) {
00971     for (j=0; j<3; j++) {
00972       Rt.mat[4*i+j] = R.mat[4*j+i];
00973     }
00974   }
00975   Matrix4 RtR(R);
00976   RtR.multmatrix(Rt);
00977 
00978   // b) 2) find the eigenvalues and eigenvectors
00979   float evalue[3];
00980   float evector[3][3];
00981   float tmpmat[4][4];
00982   for (i=0; i<4; i++)
00983     for (j=0; j<4; j++)
00984       tmpmat[i][j]=RtR.mat[4*i+j];
00985 
00986   if(jacobi(tmpmat,evalue,evector) != 0) return MEASURE_ERR_NONZEROJACOBI;
00987   // transposition the evector matrix to put the vectors in rows
00988   float vectmp;
00989   vectmp=evector[0][1]; evector[0][1]=evector[1][0]; evector[1][0]=vectmp;
00990   vectmp=evector[0][2]; evector[0][2]=evector[2][0]; evector[2][0]=vectmp;
00991   vectmp=evector[2][1]; evector[2][1]=evector[1][2]; evector[1][2]=vectmp;
00992 
00993 
00994   // b) 4) sort so that the eigenvalues are from largest to smallest
00995   //      (or rather so a[0] is eigenvector with largest eigenvalue, ...)
00996   float *a[3];
00997   a[0] = evector[0];
00998   a[1] = evector[1];
00999   a[2] = evector[2];
01000 #define SWAP(qq,ww) {                                           \
01001     float v; float *v1;                                         \
01002     v = evalue[qq]; evalue[qq] = evalue[ww]; evalue[ww] = v;    \
01003     v1 = a[qq]; a[qq] = a[ww]; a[ww] = v1;                      \
01004 }
01005   if (evalue[0] < evalue[1]) {
01006     SWAP(0, 1);
01007   }
01008   if (evalue[0] < evalue[2]) {
01009     SWAP(0, 2);
01010   }
01011   if (evalue[1] < evalue[2]) {
01012     SWAP(1, 2);
01013   }
01014 
01015   
01016   // c) determine b(i) = R*a(i)
01017   float b[3][3];
01018 
01019   Rt.multpoint3d(a[0], b[0]);
01020   vec_normalize(b[0]);
01021 
01022   Rt.multpoint3d(a[1], b[1]);
01023   vec_normalize(b[1]);
01024 
01025   Rt.multpoint3d(a[2], b[2]);
01026   vec_normalize(b[2]);
01027 
01028   // d) compute U = u(i,j) = sum(b(k,i) * a(k,j))
01029   Matrix4 U;
01030   for (i=0; i<3; i++) {
01031     for (j=0; j<3; j++) {
01032       float *tmp = &(U.mat[4*j+i]);
01033       *tmp = 0;
01034       for (int k=0; k<3; k++) {
01035         *tmp += b[k][i] * a[k][j];
01036       }
01037     }
01038   }
01039 
01040   // Check the determinant of U.  If it's negative, we need to
01041   // flip the sign of the last row.
01042   float *um = U.mat;
01043   float det = 
01044     um[0]*(um[4+1]*um[8+2] - um[4+2]*um[8+1]) -
01045     um[1]*(um[4+0]*um[8+2] - um[4+2]*um[8+0]) +
01046     um[2]*(um[4+0]*um[8+1] - um[4+1]*um[8+0]);
01047   if (det < 0) {
01048     for (int q=0; q<3; q++) um[8+q] = -um[8+q];
01049   }
01050 
01051   // e) apply the offset for com
01052   Matrix4 tx;
01053   tx.translate(-comx[0], -comx[1], -comx[2]);
01054   Matrix4 ty;
01055   ty.translate(comy[0], comy[1], comy[2]);
01056   //  U.multmatrix(com);
01057   ty.multmatrix(U);
01058   ty.multmatrix(tx);
01059   memcpy(mat->mat, ty.mat, 16*sizeof(float));
01060   return MEASURE_NOERR;
01061 }
01062 
01063 // For different values of the random seed, the computed SASA's of brH.pdb 
01064 // converge to within 1% of each other when the number of points is about
01065 // 500.  We therefore use 500 as the default number.
01066 #define NPTS 500 
01067 extern int measure_sasa(const AtomSel *sel, const float *framepos,
01068     const float *radius, float srad, float *sasa, 
01069     ResizeArray<float> *sasapts, const AtomSel *restrictsel,
01070     const int *nsamples) {
01071 
01072   // check arguments
01073   if (!sel) return MEASURE_ERR_NOSEL;
01074   if (!sel->selected) {
01075     *sasa = 0;
01076     return MEASURE_NOERR;
01077   }
01078   if (!framepos) return MEASURE_ERR_NOFRAMEPOS;
01079   if (!radius)   return MEASURE_ERR_NORADII;
01080   if (restrictsel && restrictsel->num_atoms != sel->num_atoms)
01081     return MEASURE_ERR_MISMATCHEDCNT;
01082 
01083   int i;
01084   int npts = nsamples ? *nsamples : NPTS;
01085   float maxrad = -1;
01086 
01087 
01088   // find biggest atom radius 
01089   for (i=0; i<sel->num_atoms; i++) {
01090     float rad = radius[i];
01091     if (maxrad < rad) maxrad = rad;
01092   }
01093 
01094   // Find atoms within maxrad of each other.  
01095   // build a list of pairs for each atom
01096   ResizeArray<int> *pairlist = new ResizeArray<int>[sel->num_atoms];
01097   {
01098     GridSearchPair *pairs;
01099     pairs = vmd_gridsearch1(framepos, sel->num_atoms, sel->on, 
01100                             2.0f * (maxrad + srad), 0, sel->num_atoms * 1000);
01101 
01102     GridSearchPair *p, *tmp; 
01103     for (p = pairs; p != NULL; p = tmp) {
01104       int ind1=p->ind1;
01105       int ind2=p->ind2;
01106       pairlist[ind1].append(ind2);
01107       pairlist[ind2].append(ind1);
01108       tmp = p->next;
01109       free(p);
01110     }
01111   }
01112 
01113   static const float RAND_MAX_INV = 1.0f/VMD_RAND_MAX;
01114   // Seed the random number generator before each calculation.  This gives
01115   // reproducible results and still allows a more accurate answer to be
01116   // obtained by increasing the samples size.  I don't know if this is a
01117   // "good" seed value or not, I just picked something random-looking.
01118   vmd_srandom(38572111);
01119 
01120   // All the spheres use the same random points.  
01121   float *spherepts = new float[3*npts];
01122   for (i=0; i<npts; i++) {
01123     float u1 = (float) vmd_random();
01124     float u2 = (float) vmd_random();
01125     float z = 2.0f*u1*RAND_MAX_INV -1.0f;
01126     float phi = (float) (2.0f*VMD_PI*u2*RAND_MAX_INV);
01127     float R = sqrtf(1.0f-z*z);
01128     spherepts[3*i  ] = R*cosf(phi);
01129     spherepts[3*i+1] = R*sinf(phi);
01130     spherepts[3*i+2] = z;
01131   }
01132 
01133   const float prefac = (float) (4 * VMD_PI / npts);
01134   float totarea = 0.0f;
01135   // compute area for each atom based on its pairlist
01136   for (i=sel->firstsel; i<=sel->lastsel; i++) {
01137     if (sel->on[i]) {
01138       // only atoms in restrictsel contribute
01139       if (restrictsel && !restrictsel->on[i]) continue;
01140       const float *loc = framepos+3*i;
01141       float rad = radius[i]+srad;
01142       float surfpos[3];
01143       int surfpts = npts;
01144       const ResizeArray<int> &nbrs = pairlist[i];
01145       for (int j=0; j<npts; j++) {
01146         surfpos[0] = loc[0] + rad*spherepts[3*j  ];
01147         surfpos[1] = loc[1] + rad*spherepts[3*j+1];
01148         surfpos[2] = loc[2] + rad*spherepts[3*j+2];
01149         int on = 1;
01150         for (int k=0; k<nbrs.num(); k++) {
01151           int ind = nbrs[k];
01152           const float *nbrloc = framepos+3*ind;
01153           float radsq = radius[ind]+srad; radsq *= radsq;
01154           float dx = surfpos[0]-nbrloc[0];
01155           float dy = surfpos[1]-nbrloc[1];
01156           float dz = surfpos[2]-nbrloc[2];
01157           if (dx*dx + dy*dy + dz*dz <= radsq) {
01158             on = 0;
01159             break;
01160           }
01161         }
01162         if (on) {
01163           if (sasapts) {
01164             sasapts->append(surfpos[0]);
01165             sasapts->append(surfpos[1]);
01166             sasapts->append(surfpos[2]);
01167           }
01168         } else {
01169           surfpts--;
01170         }
01171       }
01172       float atomarea = prefac * rad * rad * surfpts;
01173       totarea += atomarea;
01174     }
01175   }
01176 
01177   delete [] pairlist;
01178   delete [] spherepts;
01179   *sasa = totarea;
01180   return 0;
01181 }
01182 
01183 
01184 //
01185 // Calculate g(r)
01186 //
01187 
01188 // find the minimum image distance for one coordinate component 
01189 // and square the result (orthogonal cells only).
01190 static float fix_pbc_n_sqr(float delta, const float boxby2) {
01191   while (delta >= boxby2) { delta -= 2.0f * boxby2; }
01192   while (delta < -boxby2) { delta += 2.0f * boxby2; }
01193   return delta * delta;
01194 }
01195 
01196 // calculate the minimum distance between two positions with respect 
01197 // to the periodic cell (orthogonal cells only).
01198 static float min_dist_with_pbc(const float *a, const float *b, 
01199                                 const float *boxby2) {
01200   float distsqr;
01201   distsqr  = fix_pbc_n_sqr(a[0] - b[0], boxby2[0]);
01202   distsqr += fix_pbc_n_sqr(a[1] - b[1], boxby2[1]);
01203   distsqr += fix_pbc_n_sqr(a[2] - b[2], boxby2[2]);
01204   return sqrtf(distsqr);
01205 }
01206 
01211 static inline double spherical_cap(const double &radius, const double &boxby2) {
01212   return (VMD_PI / 3.0 * (radius - boxby2) * (radius - boxby2)
01213           * ( 2.0 * radius + boxby2));
01214 }
01215 
01216 
01217 typedef struct {
01218   int threadid;
01219   int threadcount;
01220   int count_o_start;
01221   int count_o_end;
01222   const float *olist;
01223   int count_i;
01224   const float *ilist;
01225   int count_h;
01226   int *hlist;
01227   float delta;
01228   const float *boxby2;
01229   wkfmsgtimer *msgtp;
01230   int curframe;
01231   int maxframe;
01232 } gofrparms_t;
01233     
01234 // calculate the non-normalized pair-distribution function
01235 // for two lists of atom coordinates and store the resulting
01236 // histogram in the hlist array. orthogonal cell version.
01237 //
01238 // NOTE: this is just the workhorse function. special issues related
01239 // to atoms present in both lists have to be dealt with in the uplevel
01240 // functions, that then also can do various post-processing steps.
01241 extern "C" void * measure_gofr_orth(void *voidparms) {
01242   // handle per-thread parameters
01243   gofrparms_t *parms = (gofrparms_t *) voidparms;
01244   const int count_o_start = parms->count_o_start;
01245   const int count_o_end   = parms->count_o_end;
01246   const int count_i       = parms->count_i;
01247   const int count_h       = parms->count_h;
01248   const float *olist      = parms->olist;
01249   const float *ilist      = parms->ilist;
01250   const float *boxby2     = parms->boxby2;
01251   wkfmsgtimer *msgtp         = parms->msgtp;
01252   const int curframe      = parms->curframe;
01253   const int maxframe      = parms->maxframe;
01254   int *hlist              = parms->hlist;
01255   // other local variables
01256   int i, j, idx;
01257   float dist;
01258   const float deltascale = 1.0f / parms->delta;
01259   int msgcnt=0;
01260 
01261   // loop over the chunk of pairs that was associated to this thread.
01262   for (i=count_o_start; i<count_o_end; ++i) {
01263 
01264     // print progress messages only for thread(s) that have
01265     // a timer defined (usually only tid==0).
01266     if (msgtp && wkf_msg_timer_timeout(msgtp)) {
01267       char tmpbuf[1024];
01268       if (msgcnt==0) 
01269         sprintf(tmpbuf, "timestep %d of %d", curframe, maxframe);
01270       else
01271         sprintf(tmpbuf, "timestep %d of %d: (%6.2f %% complete)", curframe, maxframe, 
01272                 (100.0f * i) / (float) (count_o_end - count_o_start + 1));
01273       msgInfo << "vmd_measure_gofr_orth: " << tmpbuf << sendmsg;
01274       msgcnt++;
01275       // XXX we should update the display here...
01276     }
01277 
01278     for (j=0; j<count_i; ++j) {
01279       // calculate distance and add to histogram
01280       dist = min_dist_with_pbc(&olist[i*3], &ilist[j*3], boxby2);
01281       idx = (int) (dist * deltascale);
01282       if ((idx >= 0) && (idx < count_h)) 
01283         ++hlist[idx];
01284     }
01285   }
01286 
01287   return MEASURE_NOERR;
01288 }
01289 
01290 // main entry point for 'measure gofr'.
01291 // tasks:
01292 // - sanity check on arguments.
01293 // - select proper algorithm for PBC treatment.
01294 int measure_gofr(AtomSel *sel1, AtomSel *sel2, MoleculeList *mlist,
01295                  const int count_h, double *gofr, double *numint, double *histog,
01296                  const float delta, int first, int last, int step, int *framecntr,
01297                  int usepbc, int selupdate) {
01298   int i, j, frame;
01299   float a, b, c, alpha, beta, gamma;
01300   int isortho=0;    // orthogonal unit cell not assumed by default.
01301   int duplicates=0; // counter for duplicate atoms in both selections.
01302 
01303   // initialize a/b/c/alpha/beta/gamma to arbitrary defaults to please the compiler.
01304   a=b=c=9999999.0;
01305   alpha=beta=gamma=90.0;
01306 
01307   // reset counter for total, skipped, and _orth processed frames.
01308   framecntr[0]=framecntr[1]=framecntr[2]=0;
01309 
01310   // First round of sanity checks.
01311   // neither list can be undefined
01312   if (!sel1 || !sel2 ) {
01313     return MEASURE_ERR_NOSEL;
01314   }
01315 
01316   // make sure that both selections are from the same molecule
01317   // so that we know that PBC unit cell info is the same for both
01318   if (sel2->molid() != sel1->molid()) {
01319     return MEASURE_ERR_MISMATCHEDMOLS;
01320   }
01321 
01322   Molecule *mymol = mlist->mol_from_id(sel1->molid());
01323   int maxframe = mymol->numframes() - 1;
01324   int nframes = 0;
01325 
01326   if (last == -1)
01327     last = maxframe;
01328 
01329   if ((last < first) || (last < 0) || (step <=0) || (first < -1)
01330       || (maxframe < 0) || (last > maxframe)) {
01331       msgErr << "measure gofr: bad frame range given." 
01332              << " max. allowed frame#: " << maxframe << sendmsg;
01333     return MEASURE_ERR_BADFRAMERANGE;
01334   }
01335 
01336   // test for non-orthogonal PBC cells, zero volume cells, etc.
01337   if (usepbc) { 
01338     for (isortho=1, nframes=0, frame=first; frame <=last; ++nframes, frame += step) {
01339       const Timestep *ts;
01340 
01341       if (first == -1) {
01342         // use current frame only. don't loop.
01343         ts = sel1->timestep(mlist);
01344         frame=last;
01345       } else {
01346         ts = mymol->get_frame(frame);
01347       }
01348       // get periodic cell information for current frame
01349       a = ts->a_length;
01350       b = ts->b_length;
01351       c = ts->c_length;
01352       alpha = ts->alpha;
01353       beta = ts->beta;
01354       gamma = ts->gamma;
01355 
01356       // check validity of PBC cell side lengths
01357       if (fabsf(a*b*c) < 0.0001) {
01358         msgErr << "measure gofr: unit cell volume is zero." << sendmsg;
01359         return MEASURE_ERR_GENERAL;
01360       }
01361 
01362       // check PBC unit cell shape to select proper low level algorithm.
01363       if ((alpha != 90.0) || (beta != 90.0) || (gamma != 90.0)) {
01364         isortho=0;
01365       }
01366     }
01367   } else {
01368     // initialize a/b/c/alpha/beta/gamma to arbitrary defaults
01369     isortho=1;
01370     a=b=c=9999999.0;
01371     alpha=beta=gamma=90.0;
01372   }
01373 
01374   // until we can handle non-orthogonal periodic cells, this is fatal
01375   if (!isortho) {
01376     msgErr << "measure gofr: only orthorhombic cells are supported (for now)." << sendmsg;
01377     return MEASURE_ERR_GENERAL;
01378   }
01379 
01380   // clear the result arrays
01381   for (i=0; i<count_h; ++i) {
01382     gofr[i] = numint[i] = histog[i] = 0.0;
01383   }
01384 
01385   // pre-allocate coordinate buffers of the max size we'll
01386   // ever need, so we don't have to reallocate if/when atom
01387   // selections are updated on-the-fly
01388   float *sel1coords = new float[3*sel1->num_atoms];
01389   float *sel2coords = new float[3*sel2->num_atoms];
01390 
01391   // setup status message timer
01392   wkfmsgtimer *msgt = wkf_msg_timer_create(5);
01393 
01394   // threading setup.
01395   wkf_thread_t *threads;
01396   gofrparms_t  *parms;
01397 #if defined(VMDTHREADS)
01398   int numprocs = wkf_thread_numprocessors();
01399 #else
01400   int numprocs = 1;
01401 #endif
01402 
01403   threads = new wkf_thread_t[numprocs];
01404   memset(threads, 0, numprocs * sizeof(wkf_thread_t));
01405 
01406   // allocate and (partially) initialize array of per-thread parameters
01407   parms = new gofrparms_t[numprocs];
01408   for (i=0; i<numprocs; ++i) {
01409     parms[i].threadid = i;
01410     parms[i].threadcount = numprocs;
01411     parms[i].delta = (float) delta;
01412     parms[i].msgtp = NULL;
01413     parms[i].count_h = count_h;
01414     parms[i].hlist = new int[count_h];
01415   }
01416 
01417   msgInfo << "measure gofr: using multi-threaded implementation with " 
01418           << numprocs << ((numprocs > 1) ? " CPUs" : " CPU") << sendmsg;
01419 
01420   for (nframes=0,frame=first; frame <=last; ++nframes, frame += step) {
01421     const Timestep *ts1, *ts2;
01422 
01423     if (frame  == -1) {
01424       // use current frame only. don't loop.
01425       ts1 = sel1->timestep(mlist);
01426       ts2 = sel2->timestep(mlist);
01427       frame=last;
01428     } else {
01429       sel1->which_frame = frame;
01430       sel2->which_frame = frame;
01431       ts1 = ts2 = mymol->get_frame(frame); // requires sels from same mol
01432     }
01433 
01434     if (usepbc) {
01435       // get periodic cell information for current frame
01436       a     = ts1->a_length;
01437       b     = ts1->b_length;
01438       c     = ts1->c_length;
01439       alpha = ts1->alpha;
01440       beta  = ts1->beta;
01441       gamma = ts1->gamma;
01442     }
01443 
01444     // compute half periodic cell size
01445     float boxby2[3];
01446     boxby2[0] = 0.5f * a;
01447     boxby2[1] = 0.5f * b;
01448     boxby2[2] = 0.5f * c;
01449 
01450     // update the selections if the user desires it
01451     if (selupdate) {
01452       if (sel1->change(NULL, mymol) != AtomSel::PARSE_SUCCESS)
01453         msgErr << "measure gofr: failed to evaluate atom selection update";
01454       if (sel2->change(NULL, mymol) != AtomSel::PARSE_SUCCESS)
01455         msgErr << "measure gofr: failed to evaluate atom selection update";
01456     }
01457 
01458     // check for duplicate atoms in the two lists, as these will have
01459     // to be subtracted back out of the first histogram slot
01460     if (sel2->molid() == sel1->molid()) {
01461       int i;
01462       for (i=0, duplicates=0; i<sel1->num_atoms; ++i) {
01463         if (sel1->on[i] && sel2->on[i])
01464           ++duplicates;
01465       }
01466     }
01467 
01468     // copy selected atoms to the two coordinate lists
01469     // requires that selections come from the same molecule
01470     const float *framepos = ts1->pos;
01471     for (i=0, j=0; i<sel1->num_atoms; ++i) {
01472       if (sel1->on[i]) {
01473         int a = i*3;
01474         sel1coords[j    ] = framepos[a    ];
01475         sel1coords[j + 1] = framepos[a + 1];
01476         sel1coords[j + 2] = framepos[a + 2];
01477         j+=3;
01478       }
01479     }
01480     framepos = ts2->pos;
01481     for (i=0, j=0; i<sel2->num_atoms; ++i) {
01482       if (sel2->on[i]) {
01483         int a = i*3;
01484         sel2coords[j    ] = framepos[a    ];
01485         sel2coords[j + 1] = framepos[a + 1];
01486         sel2coords[j + 2] = framepos[a + 2];
01487         j+=3;
01488       }
01489     }
01490 
01491     // clear the histogram for this frame
01492     // and set up argument structure for threaded execution.
01493     int maxframe = (int) ((last - first + 1) / ((float) step));
01494     for (i=0; i<numprocs; ++i) {
01495       memset(parms[i].hlist, 0, count_h * sizeof(int));
01496       parms[i].boxby2 = boxby2;
01497       parms[i].curframe = frame;
01498       parms[i].maxframe = maxframe;
01499     }
01500     parms[0].msgtp = msgt;
01501     
01502     if (isortho && sel1->selected && sel2->selected) {
01503       int count_o = sel1->selected;
01504       int count_i = sel2->selected;
01505       const float *olist = sel1coords;
01506       const float *ilist = sel2coords;
01507       // make sure the outer loop is the longer one to have 
01508       // better threading efficiency and cache utilization.
01509       if (count_o < count_i) {
01510         count_o = sel2->selected;
01511         count_i = sel1->selected;
01512         olist = sel2coords;
01513         ilist = sel1coords;
01514       }
01515 
01516       // distribute outer loop across threads in fixed size chunks.
01517       // this should work very well for small numbers of threads.
01518       // thrdelta is the chunk size and we need it to be at least 1 
01519       // _and_ numprocs*thrdelta >= count_o.
01520       int thrdelta = (count_o + (numprocs-1)) / numprocs;
01521       int o_min = 0;
01522       int o_max = thrdelta;
01523       for (i=0; i<numprocs; ++i, o_min += thrdelta, o_max += thrdelta) {
01524         if (o_max >  count_o)  o_max = count_o; // curb loop to max
01525         if (o_min >= count_o)  o_max = - 1;     // no work for this thread. too little data.
01526         parms[i].count_o_start = o_min;
01527         parms[i].count_o_end   = o_max;
01528         parms[i].count_i       = count_i;
01529         parms[i].olist         = olist;
01530         parms[i].ilist         = ilist;
01531       }
01532       
01533       // do the gofr calculation for orthogonal boxes.
01534       // XXX. non-orthogonal box not supported yet. detected and handled above.
01535 #if defined(VMDTHREADS)
01536       for (i=0; i<numprocs; ++i) {
01537         wkf_thread_create(&threads[i], measure_gofr_orth, &parms[i]);
01538       }
01539       for (i=0; i<numprocs; ++i) {
01540         wkf_thread_join(threads[i], NULL);
01541       } 
01542 #else
01543       measure_gofr_orth((void *) &parms[0]);
01544 #endif
01545       ++framecntr[2]; // frame processed with _orth algorithm
01546     } else {
01547       ++framecntr[1]; // frame skipped
01548     }
01549     ++framecntr[0];   // total frames.
01550 
01551     // correct the first histogram slot for the number of atoms that are 
01552     // present in both lists. they'll end up in the first histogram bin. 
01553     // we subtract only from the first thread histogram which is always defined.
01554     parms[0].hlist[0] -= duplicates;
01555 
01556     // in case of going 'into the edges', we should cut 
01557     // off the part that is not properly normalized to 
01558     // not confuse people that don't know about this.
01559     int h_max=count_h;
01560     float smallside=a;
01561     if (isortho && usepbc) {
01562       if(b < smallside) {
01563         smallside=b;
01564       }
01565       if(c < smallside) {
01566         smallside=c;
01567       }
01568       h_max=(int) (sqrtf(0.5f)*smallside/delta) +1;
01569       if (h_max > count_h) {
01570         h_max=count_h;
01571       }
01572     }
01573     // compute normalization function.
01574     double all=0.0;
01575     double pair_dens = 0.0;
01576     if (sel1->selected && sel2->selected) {
01577       if (usepbc) {
01578         pair_dens = a * b * c / ((double)sel1->selected * (double)sel2->selected - (double)duplicates);
01579       } else { // assume a particle volume of 30 \AA^3 (~ 1 water).
01580         pair_dens = 30.0 * (double)sel1->selected / 
01581           ((double)sel1->selected * (double)sel2->selected - (double)duplicates);
01582       }
01583     }
01584     
01585     // XXX for orthogonal boxes, we can reduce this to rmax < sqrt(0.5)*smallest side
01586     for (i=0; i<h_max; ++i) {
01587       // radius of inner and outer sphere that form the spherical slice
01588       double r_in  = delta * (double)i;
01589       double r_out = delta * (double)(i+1);
01590       double slice_vol = 4.0 / 3.0 * VMD_PI
01591         * ((r_out * r_out * r_out) - (r_in * r_in * r_in));
01592 
01593       if (isortho && usepbc) {
01594         // add correction for 0.5*box < r <= sqrt(0.5)*box
01595         if (r_out > boxby2[0]) {
01596           slice_vol -= 2.0 * spherical_cap(r_out, boxby2[0]);
01597         }
01598         if (r_out > boxby2[1]) {
01599           slice_vol -= 2.0 * spherical_cap(r_out, boxby2[1]);
01600         }
01601         if (r_out > boxby2[2]) {
01602           slice_vol -= 2.0 * spherical_cap(r_out, boxby2[2]);
01603         }
01604         if (r_in > boxby2[0]) {
01605           slice_vol += 2.0 * spherical_cap(r_in, boxby2[0]);
01606         }
01607         if (r_in > boxby2[1]) {
01608           slice_vol += 2.0 * spherical_cap(r_in, boxby2[1]);
01609         }
01610         if (r_in > boxby2[2]) {
01611           slice_vol += 2.0 * spherical_cap(r_in, boxby2[2]);
01612         }
01613       }
01614 
01615       double normf = pair_dens / slice_vol;
01616       double histv = 0.0;
01617       for (j=0; j<numprocs; ++j) {
01618         histv += (double) parms[j].hlist[i];
01619       }
01620       gofr[i] += normf * histv;
01621       all     += histv;
01622       if (sel1->selected) {
01623         numint[i] += all / (double)(sel1->selected);
01624       }
01625       histog[i] += histv;
01626     }
01627   }
01628   delete [] sel1coords;
01629   delete [] sel2coords;
01630 
01631   double norm = 1.0 / (double) nframes;
01632   for (i=0; i<count_h; ++i) {
01633     gofr[i]   *= norm;
01634     numint[i] *= norm;
01635     histog[i] *= norm;
01636   }
01637   wkf_msg_timer_destroy(msgt);
01638 
01639   // release thread-related storage
01640   for (i=0; i<numprocs; ++i) {
01641     delete [] parms[i].hlist;
01642   }
01643   delete [] threads;
01644   delete [] parms;
01645 
01646   return MEASURE_NOERR;
01647 }
01648 
01649 
01650 int measure_geom(MoleculeList *mlist, int *molid, int *atmid, ResizeArray<float> *gValues,
01651                  int frame, int first, int last, int defmolid, int geomtype) {
01652   int i, ret_val;
01653   for(i=0; i < geomtype; i++) {
01654     // make sure an atom is not repeated in this list
01655     if(i > 0 && molid[i-1]==molid[i] && atmid[i-1]==atmid[i]) {
01656       printf("measure_geom: %i/%i %i/%i\n", molid[i-1],atmid[i-1],molid[i],atmid[i]);
01657       return MEASURE_ERR_REPEATEDATOM;
01658     }
01659   }
01660 
01661   float value;
01662   int max_ts, orig_ts;
01663 
01664   // use the default molecule to determine which frames to cycle through
01665   Molecule *mol = mlist->mol_from_id(defmolid);
01666   if( !mol )
01667     return MEASURE_ERR_NOMOLECULE;
01668   
01669   // get current frame number and make sure there are frames
01670   if((orig_ts = mol->frame()) < 0)
01671     return MEASURE_ERR_NOFRAMES;
01672   
01673   // get the max frame number and determine frame range
01674   max_ts = mol->numframes()-1;
01675   if (frame<0) {
01676     if (first<0 && last<0) first = last = orig_ts; 
01677     if (last<0 || last>max_ts) last = max_ts;
01678     if (first<0) first = 0;
01679   } else {
01680     if (frame>max_ts) frame = max_ts;
01681     first = last = frame; 
01682   }
01683   
01684   // go through all the frames, calculating values
01685   for(i=first; i <= last; i++) {
01686     mol->override_current_frame(i);
01687     switch (geomtype) {
01688     case MEASURE_BOND:
01689       if ((ret_val=calculate_bond(mlist, molid, atmid, &value))<0)
01690         return ret_val;
01691       gValues->append(value);
01692       break;
01693     case MEASURE_ANGLE:
01694       if ((ret_val=calculate_angle(mlist, molid, atmid, &value))<0)
01695         return ret_val;
01696       gValues->append(value);
01697       break;
01698     case MEASURE_DIHED:
01699       if ((ret_val=calculate_dihed(mlist, molid, atmid, &value))<0)
01700         return ret_val;
01701       gValues->append(value);
01702       break;
01703     }
01704   }
01705   
01706   // reset the current frame
01707   mol->override_current_frame(orig_ts);
01708   
01709   return MEASURE_NOERR;
01710 }
01711   
01712   
01713 // calculate the value of this geometry, and return it
01714 int calculate_bond(MoleculeList *mlist, int *molid, int *atmid, float *value) {
01715 
01716   // get coords to calculate distance 
01717   int ret_val;
01718   float pos1[3], pos2[3];
01719   if ((ret_val=normal_atom_coord(mlist->mol_from_id(molid[0]), atmid[0], pos1))<0)
01720     return ret_val;
01721   if ((ret_val=normal_atom_coord(mlist->mol_from_id(molid[1]), atmid[1], pos2))<0)
01722     return ret_val;
01723   
01724   vec_sub(pos2, pos2, pos1);
01725   *value = norm(pos2);
01726 
01727   return MEASURE_NOERR;
01728 }
01729 
01730 // calculate the value of this geometry, and return it
01731 int calculate_angle(MoleculeList *mlist, int *molid, int *atmid, float *value) {
01732 
01733   // get coords to calculate distance 
01734   int ret_val;
01735   float pos1[3], pos2[3], pos3[3], r1[3], r2[3];
01736   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[0]), atmid[0], pos1))<0)
01737     return ret_val;
01738   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[1]), atmid[1], pos2))<0)
01739     return ret_val;
01740   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[2]), atmid[2], pos3))<0)
01741     return ret_val;
01742 
01743   vec_sub(r1, pos1, pos2);
01744   vec_sub(r2, pos3, pos2);
01745   *value = angle(r1, r2);
01746 
01747   return MEASURE_NOERR;
01748 }
01749 
01750 // calculate the value of this geometry, and return it
01751 int calculate_dihed(MoleculeList *mlist, int *molid, int *atmid, float *value) {
01752 
01753   // get coords to calculate distance 
01754   int ret_val;
01755   float pos1[3], pos2[3], pos3[3], pos4[3]; 
01756   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[0]), atmid[0], pos1))<0)
01757     return ret_val;
01758   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[1]), atmid[1], pos2))<0)
01759     return ret_val;
01760   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[2]), atmid[2], pos3))<0)
01761     return ret_val;
01762   if((ret_val=normal_atom_coord(mlist->mol_from_id(molid[3]), atmid[3], pos4))<0)
01763     return ret_val;
01764 
01765   *value = dihedral(pos1, pos2, pos3, pos4);
01766 
01767   return MEASURE_NOERR;
01768 }
01769 
01770 
01771 // for the given Molecule, find the UNTRANSFORMED coords for the given atom
01772 // return Molecule pointer if successful, NULL otherwise.
01773 int normal_atom_coord(Molecule *mol, int a, float *pos) {
01774   Timestep *now;
01775 
01776   int cell[3];
01777   memset(cell, 0, 3*sizeof(int));
01778 
01779   // get the molecule pointer, and get the coords for the current timestep
01780   int ret_val = check_mol(mol, a);
01781   if (ret_val<0) 
01782     return ret_val;
01783 
01784   if ((now = mol->current())) {
01785     memcpy((void *)pos, (void *)(now->pos + 3*a), 3*sizeof(float));
01786     
01787     // Apply periodic image transformation before returning
01788     Matrix4 mat;
01789     now->get_transform_from_cell(cell, mat);
01790     mat.multpoint3d(pos, pos);
01791     
01792     return MEASURE_NOERR;
01793   }
01794   
01795   // if here, error (i.e. molecule contains no frames)
01796   return MEASURE_ERR_NOFRAMES;
01797 }
01798 
01799 
01800 // check whether the given molecule & atom index is OK
01801 // if OK, return Molecule pointer; otherwise, return NULL
01802 int check_mol(Molecule *mol, int a) {
01803 
01804   if (!mol)
01805     return MEASURE_ERR_NOMOLECULE;
01806   if (a < 0 || a >= mol->nAtoms)
01807     return MEASURE_ERR_BADATOMID;
01808   
01809   return MEASURE_NOERR;
01810 }
01811 
01812 
01813 int measure_energy(MoleculeList *mlist, int *molid, int *atmid, int natoms, ResizeArray<float> *gValues,
01814                  int frame, int first, int last, int defmolid, double *params, int geomtype) {
01815   int i, ret_val;
01816   for(i=0; i < natoms; i++) {
01817     // make sure an atom is not repeated in this list
01818     if(i > 0 && molid[i-1]==molid[i] && atmid[i-1]==atmid[i]) {
01819       printf("measure_energy: %i/%i %i/%i\n", molid[i-1],atmid[i-1],molid[i],atmid[i]);
01820       return MEASURE_ERR_REPEATEDATOM;
01821     }
01822   }
01823 
01824   float value;
01825   int max_ts, orig_ts;
01826 
01827   // use the default molecule to determine which frames to cycle through
01828   Molecule *mol = mlist->mol_from_id(defmolid);
01829   if( !mol )
01830     return MEASURE_ERR_NOMOLECULE;
01831   
01832   // get current frame number and make sure there are frames
01833   if((orig_ts = mol->frame()) < 0)
01834     return MEASURE_ERR_NOFRAMES;
01835   
01836   // get the max frame number and determine frame range
01837   max_ts = mol->numframes()-1;
01838   if (frame==-1) {
01839     if (first<0 && last<0) first = last = orig_ts; 
01840     if (last<0 || last>max_ts) last = max_ts;
01841     if (first<0) first = 0;
01842   } else {
01843     if (frame>max_ts || frame==-2) frame = max_ts;
01844     first = last = frame; 
01845   }
01846   
01847   // go through all the frames, calculating values
01848   for(i=first; i <= last; i++) {
01849     mol->override_current_frame(i);
01850     switch (geomtype) {
01851     case MEASURE_BOND:
01852       if ((ret_val=compute_bond_energy(mlist, molid, atmid, &value, (float) params[0], (float) params[1]))<0)
01853         return ret_val;
01854       gValues->append(value);
01855       break;
01856     case MEASURE_ANGLE:
01857       if ((ret_val=compute_angle_energy(mlist, molid, atmid, &value, (float) params[0], (float) params[1], (float) params[2], (float) params[3]))<0)
01858         return ret_val;
01859       gValues->append(value);
01860       break;
01861     case MEASURE_DIHED:
01862       if ((ret_val=compute_dihed_energy(mlist, molid, atmid, &value, (float) params[0], int(params[1]), (float) params[2]))<0)
01863         return ret_val;
01864       gValues->append(value);
01865       break;
01866     case MEASURE_IMPRP:
01867       if ((ret_val=compute_imprp_energy(mlist, molid, atmid, &value, (float) params[0], (float) params[1]))<0)
01868         return ret_val;
01869       gValues->append(value);
01870       break;
01871     case MEASURE_VDW:
01872       if ((ret_val=compute_vdw_energy(mlist, molid, atmid, &value, (float) params[0], (float) params[1], (float) params[2], (float) params[3], (float) params[4], (float) params[5]))<0)
01873         return ret_val;
01874       gValues->append(value);
01875       break;
01876     case MEASURE_ELECT:
01877       if ((ret_val=compute_elect_energy(mlist, molid, atmid, &value, (float) params[0], (float) params[1], (bool) params[2], (bool) params[3], (float) params[4]))<0)
01878         return ret_val;
01879       gValues->append(value);
01880       break;
01881     }
01882   }
01883   
01884   // reset the current frame
01885   mol->override_current_frame(orig_ts);
01886   
01887   return MEASURE_NOERR;
01888 }
01889   
01890 // calculate the energy of this geometry
01891 int compute_bond_energy(MoleculeList *mlist, int *molid, int *atmid, float *energy, float k, float x0) {
01892   int ret_val;
01893   float dist;
01894 
01895   // Get the coordinates
01896   if ((ret_val=calculate_bond(mlist, molid, atmid, &dist))<0)
01897         return ret_val;
01898   float x = dist-x0;
01899   *energy = k*x*x;
01900 
01901   return MEASURE_NOERR;
01902 }
01903 
01904 // calculate the energy of this geometry
01905 int compute_angle_energy(MoleculeList *mlist, int *molid, int *atmid, float *energy,
01906                          float k, float x0, float kub, float s0) {
01907   int ret_val;
01908   float value;
01909 
01910   // Get the coordinates
01911   if ((ret_val=calculate_angle(mlist, molid, atmid, &value))<0)
01912         return ret_val;
01913   float x = (float) DEGTORAD((value-x0));
01914   float s = 0.0f;
01915 
01916   if (kub>0.0f) {
01917     int twoatoms[2];
01918     twoatoms[0] = atmid[0];
01919     twoatoms[1] = atmid[2];
01920     if ((ret_val=calculate_bond(mlist, molid, twoatoms, &value))<0)
01921       return ret_val;
01922     s = value-s0;
01923   }
01924 
01925   *energy = k*x*x + kub*s*s;
01926 
01927   return MEASURE_NOERR;
01928 }
01929 
01930 // calculate the energy of this geometry
01931 int compute_dihed_energy(MoleculeList *mlist, int *molid, int *atmid, float *energy,
01932                          float k, int n, float delta) {
01933   int ret_val;
01934   float value;
01935 
01936   // Get the coordinates
01937   if ((ret_val=calculate_dihed(mlist, molid, atmid, &value))<0)
01938         return ret_val;
01939   *energy = k*(1+cosf((float) (DEGTORAD((n*value-delta)))));
01940 
01941   return MEASURE_NOERR;
01942 }
01943 
01944 // calculate the energy of this geometry
01945 int compute_imprp_energy(MoleculeList *mlist, int *molid, int *atmid, float *energy,
01946                          float k, float x0) {
01947   int ret_val;
01948   float value;
01949 
01950   // Get the coordinates
01951   if ((ret_val=calculate_dihed(mlist, molid, atmid, &value))<0)
01952         return ret_val;
01953   float x = (float) (DEGTORAD((value-x0)));
01954   *energy = k*x*x;
01955 
01956   return MEASURE_NOERR;
01957 }
01958 
01959 // Calculate the VDW energy for specified pair of atoms
01960 // VDW energy:                                               
01961 // Evdw = eps * ((Rmin/dist)**12 - 2*(Rmin/dist)**6)         
01962 // eps = sqrt(eps1*eps2),  Rmin = Rmin1+Rmin2                
01963 int compute_vdw_energy(MoleculeList *mlist, int *molid, int *atmid, float *energy, float eps1, float rmin1,
01964                          float eps2, float rmin2, float cutoff, float switchdist) {
01965   int ret_val;
01966   float dist;
01967 
01968   // Get the coordinates
01969   if ((ret_val=calculate_bond(mlist, molid, atmid, &dist))<0)
01970     return ret_val;
01971 
01972   float sw=1.0;
01973   if (switchdist>0.0 && cutoff>0.0) {
01974     if (dist>=cutoff) {
01975       sw = 0.0;
01976     } else if (dist>=switchdist) {
01977       // This is the CHARMM switching function
01978       float dist2 = dist*dist;
01979       float cut2 = cutoff*cutoff;
01980       float switch2 = switchdist*switchdist;
01981       float s = cut2-dist2;
01982       float range = cut2-switch2;
01983       sw = s*s*(cut2+2*dist2-3*switch2)/(range*range*range);
01984     }
01985   }
01986 
01987   float term6 = (float) powf((rmin1+rmin2)/dist,6);
01988   *energy = sqrtf(eps1*eps2)*(term6*term6 - 2.0f*term6)*sw;
01989 
01990   return MEASURE_NOERR;
01991 }
01992 
01993 int compute_elect_energy(MoleculeList *mlist, int *molid, int *atmid, float *energy, float q1, float q2,
01994                          bool flag1, bool flag2, float cutoff) {
01995   int ret_val;
01996   float dist;
01997 
01998   // Get the coordinates
01999   if ((ret_val=calculate_bond(mlist, molid, atmid, &dist))<0)
02000     return ret_val;
02001 
02002   // Get atom charges
02003   if (!flag1) q1 = mlist->mol_from_id(molid[0])->charge()[atmid[0]];
02004   if (!flag2) q2 = mlist->mol_from_id(molid[0])->charge()[atmid[1]];
02005 
02006   if (cutoff>0.0) {
02007     if (dist<cutoff) {
02008       float efac = 1.0f-dist*dist/(cutoff*cutoff);
02009       *energy = 332.0636f*q1*q2/dist*efac*efac;
02010     } else {
02011       *energy = 0.0f;
02012     }
02013   } else {
02014     *energy = 332.0636f*q1*q2/dist;
02015   }
02016 
02017   return MEASURE_NOERR;
02018 }
02019  
02020 
02021 // Compute the center of mass for a given selection.
02022 // The result is put in rcom which has to have a size of at least 3.
02023 static void center_of_mass(AtomSel *sel, MoleculeList *mlist, float *rcom) {
02024   int i;
02025   float m = 0, mtot = 0;
02026   Molecule *mol = mlist->mol_from_id(sel->molid());
02027 
02028   // get atom masses
02029   const float *mass = mol->mass();
02030 
02031   // get atom coordinates
02032   const float *pos = sel->coordinates(mlist);
02033 
02034   memset(rcom, 0, 3*sizeof(float));
02035 
02036   // center of mass
02037   for (i=sel->firstsel; i<=sel->lastsel; i++) {
02038     if (sel->on[i]) {
02039       int ind = i * 3;
02040 
02041       m = mass[i];
02042 
02043       rcom[0] += m*pos[ind    ];
02044       rcom[1] += m*pos[ind + 1];
02045       rcom[2] += m*pos[ind + 2];
02046 
02047       // total mass
02048       mtot += m;
02049     }
02050   }
02051 
02052   rcom[0] /= mtot;
02053   rcom[1] /= mtot;
02054   rcom[2] /= mtot;
02055 }
02056 
02057 
02058 // Calculate principle axes and moments of inertia for selected atoms.
02059 // The corresponding eigenvalues are also returned, they can be used
02060 // to see if two axes are equivalent. The center of mass will be put
02061 // in parameter rcom.
02062 // The user can provide his own set of coordinates in coor. If this
02063 // parameter is NULL then the coordinates from the selection are used.
02064 extern int measure_inertia(AtomSel *sel, MoleculeList *mlist, const float *coor, float rcom[3],
02065                            float priaxes[3][3], float itensor[4][4], float evalue[3]) {
02066   if (!sel)                     return MEASURE_ERR_NOSEL;
02067   if (sel->num_atoms == 0)      return MEASURE_ERR_NOATOMS;
02068 
02069   Molecule *mol = mlist->mol_from_id(sel->molid());
02070 
02071   float x, y, z, m;
02072   float Ixx=0, Iyy=0, Izz=0, Ixy=0, Ixz=0, Iyz=0;
02073   int i,j=0;
02074 
02075   // need to put 3x3 inertia tensor into 4x4 matrix for jacobi eigensolver
02076   // itensor = {{0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 0}, {0, 0, 0, 1}};
02077   memset(itensor, 0, 16*sizeof(float));
02078   itensor[3][3] = 1.0;
02079 
02080   // compute center of mass
02081   center_of_mass(sel, mlist, rcom);
02082 
02083   // get atom coordinates
02084   const float *pos = sel->coordinates(mlist);
02085 
02086   // get atom masses
02087   const float *mass = mol->mass();
02088 
02089 
02090   // moments of inertia tensor
02091   for (i=sel->firstsel; i<=sel->lastsel; i++) {
02092     if (sel->on[i]) {
02093       // position relative to COM
02094       if (coor) {
02095         // use user provided coordinates
02096         x = coor[j*3    ] - rcom[0];
02097         y = coor[j*3 + 1] - rcom[1];
02098         z = coor[j*3 + 2] - rcom[2];
02099         j++;
02100       } else {
02101         // use coordinates from selection
02102         x = pos[i*3    ] - rcom[0];
02103         y = pos[i*3 + 1] - rcom[1];
02104         z = pos[i*3 + 2] - rcom[2];
02105       }
02106 
02107       m = mass[i];
02108 
02109       Ixx += m*(y*y+z*z);
02110       Iyy += m*(x*x+z*z);
02111       Izz += m*(x*x+y*y);
02112       Ixy -= m*x*y;
02113       Ixz -= m*x*z;
02114       Iyz -= m*y*z;
02115     }
02116   }
02117 
02118   itensor[0][0] = Ixx;
02119   itensor[1][1] = Iyy;
02120   itensor[2][2] = Izz;
02121   itensor[0][1] = Ixy;
02122   itensor[1][0] = Ixy;
02123   itensor[0][2] = Ixz;
02124   itensor[2][0] = Ixz;
02125   itensor[1][2] = Iyz;
02126   itensor[2][1] = Iyz;
02127 
02128   // Find the eigenvalues and eigenvectors of moments of inertia tensor.
02129   // The eigenvectors correspond to the principle axes of inertia.
02130   float evector[3][3];
02131   if (jacobi(itensor,evalue,evector) != 0) return MEASURE_ERR_NONZEROJACOBI;
02132 
02133   // transpose the evector matrix to put the vectors in rows
02134   float vectmp;
02135   vectmp=evector[0][1]; evector[0][1]=evector[1][0]; evector[1][0]=vectmp;
02136   vectmp=evector[0][2]; evector[0][2]=evector[2][0]; evector[2][0]=vectmp;
02137   vectmp=evector[2][1]; evector[2][1]=evector[1][2]; evector[1][2]=vectmp;
02138 
02139 
02140   // sort so that the eigenvalues are from largest to smallest
02141   // (or rather so a[0] is eigenvector with largest eigenvalue, ...)
02142   float *a[3];
02143   a[0] = evector[0];
02144   a[1] = evector[1];
02145   a[2] = evector[2];
02146   // The code for SWAP is copied from measure_fit().
02147   // It swaps rows in the eigenvector matrix.
02148 #define SWAP(qq,ww) {                                           \
02149     float v; float *v1;                                         \
02150     v = evalue[qq]; evalue[qq] = evalue[ww]; evalue[ww] = v;    \
02151     v1 = a[qq]; a[qq] = a[ww]; a[ww] = v1;                      \
02152 }
02153   if (evalue[0] < evalue[1]) {
02154     SWAP(0, 1);
02155   }
02156   if (evalue[0] < evalue[2]) {
02157     SWAP(0, 2);
02158   }
02159   if (evalue[1] < evalue[2]) {
02160     SWAP(1, 2);
02161   }
02162 
02163 #if 0
02164   // If the 2nd and 3rd eigenvalues are identical and not close to zero
02165   // then the corresponding axes are not unique. 
02166   if (evalue[1]/evalue[0]>0.1 && fabs(evalue[1]-evalue[2])/evalue[0]<0.05
02167       && fabs(evalue[0]-evalue[1])/evalue[0]>0.05) {
02168     msgInfo << "Principal axes of inertia 2 and 3 are not unique!" << sendmsg;
02169   }
02170 #endif
02171 
02172   for (i=0; i<3; i++) {
02173     for (j=0; j<3; j++) 
02174       priaxes[i][j] = a[i][j];
02175   }
02176 
02177   return MEASURE_NOERR;
02178 }
02179 

Generated on Mon May 21 01:50:57 2012 for VMD (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002