ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/openbabel/obutil.cpp
Revision: 1081
Committed: Thu Oct 19 20:49:05 2006 UTC (18 years, 6 months ago) by gezelter
File size: 32501 byte(s)
Log Message:
updated OpenBabel to version 2.0.2

File Contents

# Content
1 /**********************************************************************
2 obutil.cpp - Various utility methods.
3
4 Copyright (C) 1998-2001 by OpenEye Scientific Software, Inc.
5 Some portions Copyright (C) 2001-2005 by Geoffrey R. Hutchison
6
7 This file is part of the Open Babel project.
8 For more information, see <http://openbabel.sourceforge.net/>
9
10 This program is free software; you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation version 2 of the License.
13
14 This program is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU General Public License for more details.
18 ***********************************************************************/
19
20 #include "config.h"
21 #include "matrix3x3.hpp"
22 #include "vector3.hpp"
23 #include "mol.hpp"
24 #include "obutil.hpp"
25
26 #if HAVE_CONIO_H
27 #include <conio.h>
28 #endif
29
30 using namespace std;
31 namespace OpenBabel
32 {
33
34 /*! \class OBStopwatch
35 \brief Stopwatch class used for timing length of execution
36
37 The OBStopwatch class makes timing the execution of blocks of
38 code to microsecond accuracy very simple. The class effectively
39 has two functions, Start() and Elapsed(). The usage of the
40 OBStopwatch class is demonstrated by the following code:
41 \code
42 OBStopwatch sw;
43 sw.Start();
44 //insert code here
45 cout << "Elapsed time = " << sw.Elapsed() << endl;
46 \endcode
47 */
48
49 //! Deprecated: use the OBMessageHandler class instead
50 //! \deprecated Throw an error through the OpenBabel::OBMessageHandler class
51 OBAPI void ThrowError(char *str)
52 {
53 obErrorLog.ThrowError("", str, obInfo);
54 }
55
56 //! Deprecated: use the OBMessageHandler class instead
57 //! \deprecated Throw an error through the OpenBabel::OBMessageHandler class
58 OBAPI void ThrowError(std::string &str)
59 {
60 obErrorLog.ThrowError("", str, obInfo);
61 }
62
63 // Comparison function (for sorting ints) returns a < b
64 OBAPI bool OBCompareInt(const int &a,const int &b)
65 {
66 return(a<b);
67 }
68
69 // Comparison function (for sorting unsigned ints) returns a < b
70 OBAPI bool OBCompareUnsigned(const unsigned int &a,const unsigned int &b)
71 {
72 return(a<b);
73 }
74
75 // Comparison for doubles: returns a < (b + epsilon)
76 OBAPI bool IsNear(const double &a, const double &b, const double epsilon)
77 {
78 return (fabs(a - b) < epsilon);
79 }
80
81 // Comparison for doubles: returns a < (0.0 + epsilon)
82 OBAPI bool IsNearZero(const double &a, const double epsilon)
83 {
84 return (fabs(a) < epsilon);
85 }
86
87 //! Utility function: replace the last extension in string &src with new extension char *ext.
88 OBAPI string NewExtension(string &src,char *ext)
89 {
90 unsigned int pos = (unsigned int)src.find_last_of(".");
91 string dst;
92
93 if (pos != string::npos)
94 dst = src.substr(0,pos+1);
95 else
96 {
97 dst = src;
98 dst += ".";
99 }
100
101 dst += ext;
102 return(dst);
103 }
104
105 //! Return the geometric centroid to an array of coordinates in double* format
106 //! and center the coordinates to the origin. Operates on the first "size"
107 //! coordinates in the array.
108 OBAPI vector3 center_coords(double *c, unsigned int size)
109 {
110 if (size == 0)
111 {
112 vector3 v(0.0f, 0.0f, 0.0f);
113 return(v);
114 }
115 unsigned int i;
116 double x=0,y=0,z=0;
117 for (i = 0;i < size;i++)
118 {
119 x += c[i*3];
120 y += c[i*3+1];
121 z += c[i*3+2];
122 }
123 x /= (double) size;
124 y /= (double) size;
125 z /= (double) size;
126 for (i = 0;i < size;i++)
127 {
128 c[i*3] -= x;
129 c[i*3+1] -= y;
130 c[i*3+2] -= z;
131 }
132 vector3 v(x,y,z);
133 return(v);
134 }
135
136 //! Rotates the coordinate set *c by the transformation matrix m[3][3]
137 //! Operates on the first "size" coordinates in the array.
138 OBAPI void rotate_coords(double *c,double m[3][3],unsigned int size)
139 {
140 double x,y,z;
141 for (unsigned int i = 0;i < size;i++)
142 {
143 x = c[i*3]*m[0][0] + c[i*3+1]*m[0][1] + c[i*3+2]*m[0][2];
144 y = c[i*3]*m[1][0] + c[i*3+1]*m[1][1] + c[i*3+2]*m[1][2];
145 z = c[i*3]*m[2][0] + c[i*3+1]*m[2][1] + c[i*3+2]*m[2][2];
146 c[i*3] = x;
147 c[i*3+1] = y;
148 c[i*3+2] = z;
149 }
150 }
151
152 //! Calculate the RMS deviation between the first N coordinates of *r and *f
153 OBAPI double calc_rms(double *r,double *f, unsigned int N)
154 {
155 if (N == 0)
156 return 0.0f; // no RMS deviation between two empty sets
157
158 double d2=0.0;
159 for (unsigned int i = 0;i < N;i++)
160 {
161 d2 += SQUARE(r[i*3] - f[i*3]) +
162 SQUARE(r[i*3+1] - f[i*3+1]) +
163 SQUARE(r[i*3+2] - f[i*3+2]);
164 }
165
166 d2 /= (double) N;
167 return(sqrt(d2));
168 }
169
170 //! Rotate the coordinates of 'atoms'
171 //! such that tor == ang - atoms in 'tor' should be ordered such
172 //! that the 3rd atom is the pivot around which atoms rotate
173 OBAPI void SetRotorToAngle(double *c,vector<int> &tor,double ang,vector<int> &atoms)
174 {
175 double v1x,v1y,v1z,v2x,v2y,v2z,v3x,v3y,v3z;
176 double c1x,c1y,c1z,c2x,c2y,c2z,c3x,c3y,c3z;
177 double c1mag,c2mag,radang,costheta,m[9];
178 double x,y,z,mag,rotang,sn,cs,t,tx,ty,tz;
179
180 //
181 //calculate the torsion angle
182 //
183 v1x = c[tor[0]] - c[tor[1]];
184 v2x = c[tor[1]] - c[tor[2]];
185 v1y = c[tor[0]+1] - c[tor[1]+1];
186 v2y = c[tor[1]+1] - c[tor[2]+1];
187 v1z = c[tor[0]+2] - c[tor[1]+2];
188 v2z = c[tor[1]+2] - c[tor[2]+2];
189 v3x = c[tor[2]] - c[tor[3]];
190 v3y = c[tor[2]+1] - c[tor[3]+1];
191 v3z = c[tor[2]+2] - c[tor[3]+2];
192
193 c1x = v1y*v2z - v1z*v2y;
194 c2x = v2y*v3z - v2z*v3y;
195 c1y = -v1x*v2z + v1z*v2x;
196 c2y = -v2x*v3z + v2z*v3x;
197 c1z = v1x*v2y - v1y*v2x;
198 c2z = v2x*v3y - v2y*v3x;
199 c3x = c1y*c2z - c1z*c2y;
200 c3y = -c1x*c2z + c1z*c2x;
201 c3z = c1x*c2y - c1y*c2x;
202
203 c1mag = SQUARE(c1x)+SQUARE(c1y)+SQUARE(c1z);
204 c2mag = SQUARE(c2x)+SQUARE(c2y)+SQUARE(c2z);
205 if (c1mag*c2mag < 0.01)
206 costheta = 1.0; //avoid div by zero error
207 else
208 costheta = (c1x*c2x + c1y*c2y + c1z*c2z)/(sqrt(c1mag*c2mag));
209
210 if (costheta < -0.999999)
211 costheta = -0.999999;
212 if (costheta > 0.999999)
213 costheta = 0.999999;
214
215 if ((v2x*c3x + v2y*c3y + v2z*c3z) > 0.0)
216 radang = -acos(costheta);
217 else
218 radang = acos(costheta);
219
220 //
221 // now we have the torsion angle (radang) - set up the rot matrix
222 //
223
224 //find the difference between current and requested
225 rotang = ang - radang;
226
227 sn = sin(rotang);
228 cs = cos(rotang);
229 t = 1 - cs;
230 //normalize the rotation vector
231 mag = sqrt(SQUARE(v2x)+SQUARE(v2y)+SQUARE(v2z));
232 x = v2x/mag;
233 y = v2y/mag;
234 z = v2z/mag;
235
236 //set up the rotation matrix
237 m[0]= t*x*x + cs;
238 m[1] = t*x*y + sn*z;
239 m[2] = t*x*z - sn*y;
240 m[3] = t*x*y - sn*z;
241 m[4] = t*y*y + cs;
242 m[5] = t*y*z + sn*x;
243 m[6] = t*x*z + sn*y;
244 m[7] = t*y*z - sn*x;
245 m[8] = t*z*z + cs;
246
247 //
248 //now the matrix is set - time to rotate the atoms
249 //
250 tx = c[tor[1]];
251 ty = c[tor[1]+1];
252 tz = c[tor[1]+2];
253 vector<int>::iterator i;
254 int j;
255 for (i = atoms.begin();i != atoms.end();i++)
256 {
257 j = *i;
258 c[j] -= tx;
259 c[j+1] -= ty;
260 c[j+2]-= tz;
261 x = c[j]*m[0] + c[j+1]*m[1] + c[j+2]*m[2];
262 y = c[j]*m[3] + c[j+1]*m[4] + c[j+2]*m[5];
263 z = c[j]*m[6] + c[j+1]*m[7] + c[j+2]*m[8];
264 c[j] = x;
265 c[j+1] = y;
266 c[j+2] = z;
267 c[j] += tx;
268 c[j+1] += ty;
269 c[j+2] += tz;
270 }
271 }
272
273 //! Safely open the supplied filename and return an ifstream, throwing an error
274 //! to the default OBMessageHandler error log if it fails.
275 OBAPI bool SafeOpen(ifstream &fs,char *filename)
276 {
277 #ifdef WIN32
278 string s = filename;
279 if (s.find(".bin") != string::npos)
280 fs.open(filename,ios::binary);
281 else
282 #endif
283
284 fs.open(filename);
285
286 if (!fs)
287 {
288 string error = "Unable to open file \'";
289 error += filename;
290 error += "\' in read mode";
291 obErrorLog.ThrowError(__func__, error, obError);
292 return(false);
293 }
294
295 return(true);
296 }
297
298
299 //! Safely open the supplied filename and return an ofstream, throwing an error
300 //! to the default OBMessageHandler error log if it fails.
301 OBAPI bool SafeOpen(ofstream &fs,char *filename)
302 {
303 #ifdef WIN32
304 string s = filename;
305 if (s.find(".bin") != string::npos)
306 fs.open(filename,ios::binary);
307 else
308 #endif
309
310 fs.open(filename);
311
312 if (!fs)
313 {
314 string error = "Unable to open file \'";
315 error += filename;
316 error += "\' in write mode";
317 obErrorLog.ThrowError(__func__, error, obError);
318 return(false);
319 }
320
321 return(true);
322 }
323
324 //! Safely open the supplied filename and return an ifstream, throwing an error
325 //! to the default OBMessageHandler error log if it fails.
326 OBAPI bool SafeOpen(ifstream &fs,string &filename)
327 {
328 return(SafeOpen(fs,(char*)filename.c_str()));
329 }
330
331 //! Safely open the supplied filename and return an ofstream, throwing an error
332 //! to the default OBMessageHandler error log if it fails.
333 OBAPI bool SafeOpen(ofstream &fs,string &filename)
334 {
335 return(SafeOpen(fs,(char*)filename.c_str()));
336 }
337
338 //! Shift the supplied string to uppercase
339 OBAPI void ToUpper(std::string &s)
340 {
341 if (s.empty())
342 return;
343 unsigned int i;
344 for (i = 0;i < s.size();i++)
345 if (isalpha(s[i]) && !isdigit(s[i]))
346 s[i] = toupper(s[i]);
347 }
348
349 //! Shift the supplied char* to uppercase
350 OBAPI void ToUpper(char *cptr)
351 {
352 char *c;
353 for (c = cptr;*c != '\0';c++)
354 if (isalpha(*c) && !isdigit(*c))
355 *c = toupper(*c);
356 }
357
358 //! Shift the supplied string to lowercase
359 OBAPI void ToLower(std::string &s)
360 {
361 if (s.empty())
362 return;
363 unsigned int i;
364 for (i = 0;i < s.size();i++)
365 if (isalpha(s[i]) && !isdigit(s[i]))
366 s[i] = tolower(s[i]);
367 }
368
369 //! Shift the supplied char* to lowercase
370 OBAPI void ToLower(char *cptr)
371 {
372 char *c;
373 for (c = cptr;*c != '\0';c++)
374 if (isalpha(*c) && !isdigit(*c))
375 *c = tolower(*c);
376 }
377
378 //! "Clean" the supplied atom type, shifting the first character to uppercase,
379 //! the second character (if it's a letter) to lowercase, and terminating with a NULL
380 //! to strip off any trailing characters
381 OBAPI void CleanAtomType(char *id)
382 {
383 id[0] = toupper(id[0]);
384 if (isalpha(id[1]) == 0)
385 id[1] = '\0';
386 else
387 {
388 id[1] = tolower(id[1]);
389 id[2] = '\0';
390 }
391 }
392
393 //! Transform the supplied vector<OBInternalCoord*> into cartesian and update
394 //! the OBMol accordingly.
395 //! Implements <a href="http://qsar.sourceforge.net/dicts/blue-obelisk/index.xhtml#zmatrixCoordinatesIntoCartesianCoordinates">blue-obelisk:zmatrixCoordinatesIntoCartesianCoordinates</a>
396 OBAPI void InternalToCartesian(std::vector<OBInternalCoord*> &vic,OBMol &mol)
397 {
398 vector3 n,nn,v1,v2,v3,avec,bvec,cvec;
399 double dst = 0.0, ang = 0.0, tor = 0.0;
400 OBAtom *atom;
401 vector<OBNodeBase*>::iterator i;
402 int index;
403
404 if (vic.empty())
405 return;
406
407 obErrorLog.ThrowError(__func__,
408 "Ran OpenBabel::InternalToCartesian", obAuditMsg);
409
410 for (atom = mol.BeginAtom(i);atom;atom = mol.NextAtom(i))
411 {
412 index = atom->GetIdx();
413
414 // make sure we always have valid pointers
415 if (index >= vic.size() || !vic[index])
416 return;
417
418 if (vic[index]->_a) // make sure we have a valid ptr
419 {
420 avec = vic[index]->_a->GetVector();
421 dst = vic[index]->_dst;
422 }
423 else
424 {
425 // atom 1
426 atom->SetVector(0.0, 0.0, 0.0);
427 continue;
428 }
429
430 if (vic[index]->_b)
431 {
432 bvec = vic[index]->_b->GetVector();
433 ang = vic[index]->_ang * DEG_TO_RAD;
434 }
435 else
436 {
437 // atom 2
438 atom->SetVector(dst, 0.0, 0.0);
439 continue;
440 }
441
442 if (vic[index]->_c)
443 {
444 cvec = vic[index]->_c->GetVector();
445 tor = vic[index]->_tor * DEG_TO_RAD;
446 }
447 else
448 {
449 // atom 3
450 cvec = VY;
451 tor = 90.0 * DEG_TO_RAD;
452 }
453
454 v1 = avec - bvec;
455 v2 = avec - cvec;
456 n = cross(v1,v2);
457 nn = cross(v1,n);
458 n.normalize();
459 nn.normalize();
460
461 n *= -sin(tor);
462 nn *= cos(tor);
463 v3 = n + nn;
464 v3.normalize();
465 v3 *= dst * sin(ang);
466 v1.normalize();
467 v1 *= dst * cos(ang);
468 v2 = avec + v3 - v1;
469
470 atom->SetVector(v2);
471 }
472
473 // Delete dummy atoms
474 for (atom = mol.BeginAtom(i);atom;atom = mol.NextAtom(i))
475 if (atom->GetAtomicNum() == 0)
476 mol.DeleteAtom(atom);
477 }
478
479 //! Use the supplied OBMol and its Cartesian coordinates to generate
480 //! a set of internal (z-matrix) coordinates as supplied in the
481 //! vector<OBInternalCoord*> argument.
482 //! Implements <a href="http://qsar.sourceforge.net/dicts/blue-obelisk/index.xhtml#cartesianCoordinatesIntoZmatrixCoordinates">blue-obelisk:cartesianCoordinatesIntoZmatrixCoordinates</a>.
483 OBAPI void CartesianToInternal(std::vector<OBInternalCoord*> &vic,OBMol &mol)
484 {
485 double r,sum;
486 OBAtom *atom,*nbr,*ref;
487 vector<OBNodeBase*>::iterator i,j,m;
488
489 obErrorLog.ThrowError(__func__,
490 "Ran OpenBabel::CartesianToInternal", obAuditMsg);
491
492 //set reference atoms
493 for (atom = mol.BeginAtom(i);atom;atom = mol.NextAtom(i))
494 {
495 if (atom->GetIdx() == 1)
496 continue;
497 else if (atom->GetIdx() == 2)
498 {
499 vic[atom->GetIdx()]->_a = mol.GetAtom(1);
500 continue;
501 }
502 else if (atom->GetIdx() == 3)
503 {
504 if( (atom->GetVector()-mol.GetAtom(2)->GetVector()).length_2()
505 <(atom->GetVector()-mol.GetAtom(1)->GetVector()).length_2())
506 {
507 vic[atom->GetIdx()]->_a = mol.GetAtom(2);
508 vic[atom->GetIdx()]->_b = mol.GetAtom(1);
509 }
510 else
511 {
512 vic[atom->GetIdx()]->_a = mol.GetAtom(1);
513 vic[atom->GetIdx()]->_b = mol.GetAtom(2);
514 }
515 continue;
516 }
517 sum=1.0E10;
518 ref = mol.GetAtom(1);
519 for(nbr = mol.BeginAtom(j);nbr && (i != j);nbr = mol.NextAtom(j))
520 {
521 r = (atom->GetVector()-nbr->GetVector()).length_2();
522 if((r < sum) && (vic[nbr->GetIdx()]->_a != nbr) &&
523 (vic[nbr->GetIdx()]->_b != nbr))
524 {
525 sum = r;
526 ref = nbr;
527 }
528 }
529
530 vic[atom->GetIdx()]->_a = ref;
531 if (ref->GetIdx() >= 3)
532 {
533 vic[atom->GetIdx()]->_b = vic[ref->GetIdx()]->_a;
534 vic[atom->GetIdx()]->_c = vic[ref->GetIdx()]->_b;
535 }
536 else
537 {
538 if(ref->GetIdx()== 1)
539 {
540 vic[atom->GetIdx()]->_b = mol.GetAtom(2);
541 vic[atom->GetIdx()]->_c = mol.GetAtom(3);
542 }
543 else
544 {//ref->GetIdx()== 2
545 vic[atom->GetIdx()]->_b = mol.GetAtom(1);
546 vic[atom->GetIdx()]->_c = mol.GetAtom(3);
547 }
548 }
549 }
550
551 //fill in geometries
552 unsigned int k;
553 vector3 v1,v2;
554 OBAtom *a,*b,*c;
555 for (k = 2;k <= mol.NumAtoms();k++)
556 {
557 atom = mol.GetAtom(k);
558 a = vic[k]->_a;
559 b = vic[k]->_b;
560 c = vic[k]->_c;
561 if (k == 2)
562 {
563 vic[k]->_dst = (atom->GetVector() - a->GetVector()).length();
564 continue;
565 }
566
567 v1 = atom->GetVector() - a->GetVector();
568 v2 = b->GetVector() - a->GetVector();
569 vic[k]->_dst = v1.length();
570 vic[k]->_ang = vectorAngle(v1,v2);
571
572 if (k == 3)
573 continue;
574 vic[k]->_tor = CalcTorsionAngle(atom->GetVector(),
575 a->GetVector(),
576 b->GetVector(),
577 c->GetVector());
578 }
579
580 //check for linear geometries and try to correct if possible
581 bool done;
582 double ang;
583 for (k = 2;k <= mol.NumAtoms();k++)
584 {
585 ang = fabs(vic[k]->_ang);
586 if (ang > 5.0 && ang < 175.0)
587 continue;
588 atom = mol.GetAtom(k);
589 done = false;
590 for (a = mol.BeginAtom(i);a && a->GetIdx() < k && !done;a = mol.NextAtom(i))
591 for (b=mol.BeginAtom(j);b && b->GetIdx()<a->GetIdx() && !done;b = mol.NextAtom(j))
592 {
593 v1 = atom->GetVector() - a->GetVector();
594 v2 = b->GetVector() - a->GetVector();
595 ang = fabs(vectorAngle(v1,v2));
596 if (ang < 5.0 || ang > 175.0)
597 continue;
598
599 for (c = mol.BeginAtom(m);c && c->GetIdx() < atom->GetIdx();c = mol.NextAtom(m))
600 if (c != atom && c != a && c != b)
601 break;
602 if (!c)
603 continue;
604
605 vic[k]->_a = a;
606 vic[k]->_b = b;
607 vic[k]->_c = c;
608 vic[k]->_dst = v1.length();
609 vic[k]->_ang = vectorAngle(v1,v2);
610 vic[k]->_tor = CalcTorsionAngle(atom->GetVector(),
611 a->GetVector(),
612 b->GetVector(),
613 c->GetVector());
614 done = true;
615 }
616 }
617 }
618
619 OBAPI void qtrfit (double *r,double *f,int size, double u[3][3])
620 {
621 register int i;
622 double xxyx, xxyy, xxyz;
623 double xyyx, xyyy, xyyz;
624 double xzyx, xzyy, xzyz;
625 double d[4],q[4];
626 double c[16],v[16];
627 double rx,ry,rz,fx,fy,fz;
628
629 /* generate the upper triangle of the quadratic form matrix */
630
631 xxyx = 0.0;
632 xxyy = 0.0;
633 xxyz = 0.0;
634 xyyx = 0.0;
635 xyyy = 0.0;
636 xyyz = 0.0;
637 xzyx = 0.0;
638 xzyy = 0.0;
639 xzyz = 0.0;
640
641 for (i = 0; i < size; i++)
642 {
643 rx = r[i*3];
644 ry = r[i*3+1];
645 rz = r[i*3+2];
646 fx = f[i*3];
647 fy = f[i*3+1];
648 fz = f[i*3+2];
649
650 xxyx += fx * rx;
651 xxyy += fx * ry;
652 xxyz += fx * rz;
653 xyyx += fy * rx;
654 xyyy += fy * ry;
655 xyyz += fy * rz;
656 xzyx += fz * rx;
657 xzyy += fz * ry;
658 xzyz += fz * rz;
659 }
660
661 c[4*0+0] = xxyx + xyyy + xzyz;
662
663 c[4*0+1] = xzyy - xyyz;
664 c[4*1+1] = xxyx - xyyy - xzyz;
665
666 c[4*0+2] = xxyz - xzyx;
667 c[4*1+2] = xxyy + xyyx;
668 c[4*2+2] = xyyy - xzyz - xxyx;
669
670 c[4*0+3] = xyyx - xxyy;
671 c[4*1+3] = xzyx + xxyz;
672 c[4*2+3] = xyyz + xzyy;
673 c[4*3+3] = xzyz - xxyx - xyyy;
674
675 /* diagonalize c */
676
677 matrix3x3::jacobi(4, c, d, v);
678
679 /* extract the desired quaternion */
680
681 q[0] = v[4*0+3];
682 q[1] = v[4*1+3];
683 q[2] = v[4*2+3];
684 q[3] = v[4*3+3];
685
686 /* generate the rotation matrix */
687
688 u[0][0] = q[0]*q[0] + q[1]*q[1] - q[2]*q[2] - q[3]*q[3];
689 u[1][0] = 2.0 * (q[1] * q[2] - q[0] * q[3]);
690 u[2][0] = 2.0 * (q[1] * q[3] + q[0] * q[2]);
691
692 u[0][1] = 2.0 * (q[2] * q[1] + q[0] * q[3]);
693 u[1][1] = q[0]*q[0] - q[1]*q[1] + q[2]*q[2] - q[3]*q[3];
694 u[2][1] = 2.0 * (q[2] * q[3] - q[0] * q[1]);
695
696 u[0][2] = 2.0 * (q[3] * q[1] - q[0] * q[2]);
697 u[1][2] = 2.0 * (q[3] * q[2] + q[0] * q[1]);
698 u[2][2] = q[0]*q[0] - q[1]*q[1] - q[2]*q[2] + q[3]*q[3];
699 }
700
701
702
703 static double Roots[4];
704
705 #define ApproxZero 1E-7
706 #define IsZero(x) ((double)fabs(x)<ApproxZero)
707 #ifndef PI
708 #define PI 3.14159265358979323846226433
709 #endif
710 #define OneThird (1.0/3.0)
711 #define FourThirdsPI (4.0*PI/3.0)
712 #define TwoThirdsPI (2.0*PI/3.0)
713
714 #ifdef OLD_RMAT
715
716 /*FUNCTION */
717 /* recieves: the co-efficients for a general
718 * equation of degree one.
719 * Ax + B = 0 !!
720 */
721 OBAPI static int SolveLinear(double A,double B)
722 {
723 if( IsZero(A) )
724 return( 0 );
725 Roots[0] = -B/A;
726 return( 1 );
727 }
728
729 /*FUNCTION */
730 /* recieves: the co-efficients for a general
731 * linear equation of degree two.
732 * Ax^2 + Bx + C = 0 !!
733 */
734 OBAPI static int SolveQuadratic(double A,double B,double C)
735 {
736 register double Descr, Temp, TwoA;
737
738 if( IsZero(A) )
739 return( SolveLinear(B,C) );
740
741 TwoA = A+A;
742 Temp = TwoA*C;
743 Descr = B*B - (Temp+Temp);
744 if( Descr<0.0 )
745 return( 0 );
746
747 if( Descr>0.0 )
748 {
749 Descr = sqrt(Descr);
750 #ifdef ORIG
751
752 Roots[0] = (-B-Descr)/TwoA;
753 Roots[1] = (-B+Descr)/TwoA;
754 #else
755 /* W. Press, B. Flannery, S. Teukolsky and W. Vetterling,
756 * "Quadratic and Cubic Equations", Numerical Recipes in C,
757 * Chapter 5, pp. 156-157, 1989.
758 */
759 Temp = (B<0.0)? -0.5*(B-Descr) : -0.5*(B+Descr);
760 Roots[0] = Temp/A;
761 Roots[1] = C/Temp;
762 #endif
763
764 return( 2 );
765 }
766 Roots[0] = -B/TwoA;
767 return( 1 );
768 }
769
770 /*FUNCTION */
771 /* task: to return the cube root of the
772 * given value taking into account
773 * that it may be negative.
774 */
775 OBAPI static double CubeRoot(double X)
776 {
777 if( X>=0.0 )
778 {
779 return pow( X, OneThird );
780 }
781 else
782 return -pow( -X, OneThird );
783 }
784
785 OBAPI static int SolveCubic(double A,double B,double C,double D)
786 {
787 register double TwoA, ThreeA, BOver3A;
788 register double Temp, POver3, QOver2;
789 register double Desc, Rho, Psi;
790
791
792 if( IsZero(A) )
793 {
794 return( SolveQuadratic(B,C,D) );
795 }
796
797 TwoA = A+A;
798 ThreeA = TwoA+A;
799 BOver3A = B/ThreeA;
800 QOver2 = ((TwoA*BOver3A*BOver3A-C)*BOver3A+D)/TwoA;
801 POver3 = (C-B*BOver3A)/ThreeA;
802
803
804 Rho = POver3*POver3*POver3;
805 Desc = QOver2*QOver2 + Rho;
806
807 if( Desc<=0.0 )
808 {
809 Rho = sqrt( -Rho );
810 Psi = OneThird*acos(-QOver2/Rho);
811 Temp = CubeRoot( Rho );
812 Temp = Temp+Temp;
813
814 Roots[0] = Temp*cos( Psi )-BOver3A;
815 Roots[1] = Temp*cos( Psi+TwoThirdsPI )-BOver3A;
816 Roots[2] = Temp*cos( Psi+FourThirdsPI )-BOver3A;
817 return( 3 );
818 }
819
820 if( Desc> 0.0 )
821 {
822 Temp = CubeRoot( -QOver2 );
823 Roots[0] = Temp+Temp-BOver3A;
824 Roots[1] = -Temp-BOver3A;
825 return( 2 );
826 }
827
828 Desc = sqrt( Desc );
829 Roots[0] = CubeRoot(Desc-QOver2)-CubeRoot(Desc+QOver2) - BOver3A;
830
831 return( 1 );
832 }
833 #endif
834
835
836 #define MAX_SWEEPS 50
837
838 OBAPI void ob_make_rmat(double a[3][3],double rmat[9])
839 {
840 double onorm, dnorm;
841 double b, dma, q, t, c, s,d[3];
842 double atemp, vtemp, dtemp,v[3][3];
843 double r1[3],r2[3],v1[3],v2[3],v3[3];
844 int i, j, k, l;
845
846 memset((char*)d,'\0',sizeof(double)*3);
847
848 for (j = 0; j < 3; j++)
849 {
850 for (i = 0; i < 3; i++)
851 v[i][j] = 0.0;
852
853 v[j][j] = 1.0;
854 d[j] = a[j][j];
855 }
856
857 for (l = 1; l <= MAX_SWEEPS; l++)
858 {
859 dnorm = 0.0;
860 onorm = 0.0;
861 for (j = 0; j < 3; j++)
862 {
863 dnorm = dnorm + (double)fabs(d[j]);
864 for (i = 0; i <= j - 1; i++)
865 {
866 onorm = onorm + (double)fabs(a[i][j]);
867 }
868 }
869
870 if((onorm/dnorm) <= 1.0e-12)
871 goto Exit_now;
872 for (j = 1; j < 3; j++)
873 {
874 for (i = 0; i <= j - 1; i++)
875 {
876 b = a[i][j];
877 if(fabs(b) > 0.0)
878 {
879 dma = d[j] - d[i];
880 if((fabs(dma) + fabs(b)) <= fabs(dma))
881 t = b / dma;
882 else
883 {
884 q = 0.5 * dma / b;
885 t = 1.0/((double)fabs(q) + (double)sqrt(1.0+q*q));
886 if(q < 0.0)
887 t = -t;
888 }
889 c = 1.0/(double)sqrt(t * t + 1.0);
890 s = t * c;
891 a[i][j] = 0.0;
892 for (k = 0; k <= i-1; k++)
893 {
894 atemp = c * a[k][i] - s * a[k][j];
895 a[k][j] = s * a[k][i] + c * a[k][j];
896 a[k][i] = atemp;
897 }
898 for (k = i+1; k <= j-1; k++)
899 {
900 atemp = c * a[i][k] - s * a[k][j];
901 a[k][j] = s * a[i][k] + c * a[k][j];
902 a[i][k] = atemp;
903 }
904 for (k = j+1; k < 3; k++)
905 {
906 atemp = c * a[i][k] - s * a[j][k];
907 a[j][k] = s * a[i][k] + c * a[j][k];
908 a[i][k] = atemp;
909 }
910 for (k = 0; k < 3; k++)
911 {
912 vtemp = c * v[k][i] - s * v[k][j];
913 v[k][j] = s * v[k][i] + c * v[k][j];
914 v[k][i] = vtemp;
915 }
916 dtemp = c*c*d[i] + s*s*d[j] - 2.0*c*s*b;
917 d[j] = s*s*d[i] + c*c*d[j] + 2.0*c*s*b;
918 d[i] = dtemp;
919 } /* end if */
920 } /* end for i */
921 } /* end for j */
922 } /* end for l */
923
924 Exit_now:
925
926 /* max_sweeps = l;*/
927
928 for (j = 0; j < 3-1; j++)
929 {
930 k = j;
931 dtemp = d[k];
932 for (i = j+1; i < 3; i++)
933 if(d[i] < dtemp)
934 {
935 k = i;
936 dtemp = d[k];
937 }
938
939 if(k > j)
940 {
941 d[k] = d[j];
942 d[j] = dtemp;
943 for (i = 0; i < 3 ; i++)
944 {
945 dtemp = v[i][k];
946 v[i][k] = v[i][j];
947 v[i][j] = dtemp;
948 }
949 }
950 }
951
952 r1[0] = v[0][0];
953 r1[1] = v[1][0];
954 r1[2] = v[2][0];
955 r2[0] = v[0][1];
956 r2[1] = v[1][1];
957 r2[2] = v[2][1];
958
959 v3[0] = r1[1]*r2[2] - r1[2]*r2[1];
960 v3[1] = -r1[0]*r2[2] + r1[2]*r2[0];
961 v3[2] = r1[0]*r2[1] - r1[1]*r2[0];
962 s = (double)sqrt(v3[0]*v3[0] + v3[1]*v3[1] + v3[2]*v3[2]);
963 v3[0] /= s;
964 v3[0] /= s;
965 v3[0] /= s;
966
967 v2[0] = v3[1]*r1[2] - v3[2]*r1[1];
968 v2[1] = -v3[0]*r1[2] + v3[2]*r1[0];
969 v2[2] = v3[0]*r1[1] - v3[1]*r1[0];
970 s = (double)sqrt(v2[0]*v2[0] + v2[1]*v2[1] + v2[2]*v2[2]);
971 v2[0] /= s;
972 v2[0] /= s;
973 v2[0] /= s;
974
975 v1[0] = v2[1]*v3[2] - v2[2]*v3[1];
976 v1[1] = -v2[0]*v3[2] + v2[2]*v3[0];
977 v1[2] = v2[0]*v3[1] - v2[1]*v3[0];
978 s = (double)sqrt(v1[0]*v1[0] + v1[1]*v1[1] + v1[2]*v1[2]);
979 v1[0] /= s;
980 v1[0] /= s;
981 v1[0] /= s;
982
983 rmat[0] = v1[0];
984 rmat[1] = v1[1];
985 rmat[2] = v1[2];
986 rmat[3] = v2[0];
987 rmat[4] = v2[1];
988 rmat[5] = v2[2];
989 rmat[6] = v3[0];
990 rmat[7] = v3[1];
991 rmat[8] = v3[2];
992 }
993
994 static int get_roots_3_3(double mat[3][3], double roots[3])
995 {
996 double rmat[9];
997
998 ob_make_rmat(mat,rmat);
999
1000 mat[0][0]=rmat[0];
1001 mat[0][1]=rmat[3];
1002 mat[0][2]=rmat[6];
1003 mat[1][0]=rmat[1];
1004 mat[1][1]=rmat[4];
1005 mat[1][2]=rmat[7];
1006 mat[2][0]=rmat[2];
1007 mat[2][1]=rmat[5];
1008 mat[2][2]=rmat[8];
1009
1010 roots[0]=(double)Roots[0];
1011 roots[1]=(double)Roots[1];
1012 roots[2]=(double)Roots[2];
1013
1014 return 1;
1015 }
1016
1017 OBAPI double superimpose(double *r,double *f,int size)
1018 {
1019 int i,j;
1020 double x,y,z,d2;
1021 double mat[3][3],rmat[3][3],mat2[3][3],roots[3];
1022
1023 /* make inertial cross tensor */
1024 for(i=0;i<3;i++)
1025 for(j=0;j<3;j++)
1026 mat[i][j]=0.0;
1027
1028 for(i=0;i < size;i++)
1029 {
1030 mat[0][0]+=r[3*i] *f[3*i];
1031 mat[1][0]+=r[3*i+1]*f[3*i];
1032 mat[2][0]+=r[3*i+2]*f[3*i];
1033 mat[0][1]+=r[3*i] *f[3*i+1];
1034 mat[1][1]+=r[3*i+1]*f[3*i+1];
1035 mat[2][1]+=r[3*i+2]*f[3*i+1];
1036 mat[0][2]+=r[3*i] *f[3*i+2];
1037 mat[1][2]+=r[3*i+1]*f[3*i+2];
1038 mat[2][2]+=r[3*i+2]*f[3*i+2];
1039 }
1040
1041 d2=mat[0][0]*(mat[1][1]*mat[2][2]-mat[1][2]*mat[2][1])
1042 -mat[0][1]*(mat[1][0]*mat[2][2]-mat[1][2]*mat[2][0])
1043 +mat[0][2]*(mat[1][0]*mat[2][1]-mat[1][1]*mat[2][0]);
1044
1045
1046 /* square matrix= ((mat transpose) * mat) */
1047 for(i=0;i<3;i++)
1048 for(j=0;j<3;j++)
1049 {
1050 x=mat[0][i]*mat[0][j]+mat[1][i]*mat[1][j]+mat[2][i]*mat[2][j];
1051 mat2[i][j]=mat[i][j];
1052 rmat[i][j]=x;
1053 }
1054 get_roots_3_3(rmat,roots);
1055
1056 roots[0]=(roots[0]<0.0001) ? 0.0: (roots[0]);
1057 roots[1]=(roots[1]<0.0001) ? 0.0: (roots[1]);
1058 roots[2]=(roots[2]<0.0001) ? 0.0: (roots[2]);
1059
1060 /* make sqrt of rmat, store in mat*/
1061
1062 roots[0]=roots[0]<0.0001? 0.0: 1.0/(double)sqrt(roots[0]);
1063 roots[1]=roots[1]<0.0001? 0.0: 1.0/(double)sqrt(roots[1]);
1064 roots[2]=roots[2]<0.0001? 0.0: 1.0/(double)sqrt(roots[2]);
1065
1066 if(d2<0.0)
1067 {
1068 if( (roots[0]>=roots[1]) && (roots[0]>=roots[2]) )
1069 roots[0]*=-1.0;
1070 if( (roots[1]>roots[0]) && (roots[1]>=roots[2]) )
1071 roots[1]*=-1.0;
1072 if( (roots[2]>roots[1]) && (roots[2]>roots[0]) )
1073 roots[2]*=-1.0;
1074 }
1075
1076 for(i=0;i<3;i++)
1077 for(j=0;j<3;j++)
1078 mat[i][j]=roots[0]*rmat[i][0]*rmat[j][0]+
1079 roots[1]*rmat[i][1]*rmat[j][1]+
1080 roots[2]*rmat[i][2]*rmat[j][2];
1081
1082 /* and multiply into original inertial cross matrix, mat2 */
1083 for(i=0;i<3;i++)
1084 for(j=0;j<3;j++)
1085 rmat[i][j]=mat[0][j]*mat2[i][0]+
1086 mat[1][j]*mat2[i][1]+
1087 mat[2][j]*mat2[i][2];
1088
1089 /* rotate all coordinates */
1090 d2 = 0.0;
1091 for(i=0;i<size;i++)
1092 {
1093 x=f[3*i]*rmat[0][0]+f[3*i+1]*rmat[0][1]+f[3*i+2]*rmat[0][2];
1094 y=f[3*i]*rmat[1][0]+f[3*i+1]*rmat[1][1]+f[3*i+2]*rmat[1][2];
1095 z=f[3*i]*rmat[2][0]+f[3*i+1]*rmat[2][1]+f[3*i+2]*rmat[2][2];
1096 f[3*i ]=x;
1097 f[3*i+1]=y;
1098 f[3*i+2]=z;
1099
1100 x = r[i*3] - f[i*3];
1101 y = r[i*3+1] - f[i*3+1];
1102 z = r[i*3+2] - f[i*3+2];
1103 d2 += x*x+y*y+z*z;
1104 }
1105
1106 d2 /= (double) size;
1107
1108 return((double)sqrt(d2));
1109 }
1110
1111 OBAPI void get_rmat(double *rvec,double *r,double *f,int size)
1112 {
1113 int i,j;
1114 double x,d2;
1115 double mat[3][3],rmat[3][3],mat2[3][3],roots[3];
1116
1117 /* make inertial cross tensor */
1118 for(i=0;i<3;i++)
1119 for(j=0;j<3;j++)
1120 mat[i][j]=0.0;
1121
1122 for(i=0;i < size;i++)
1123 {
1124 mat[0][0]+=r[3*i] *f[3*i];
1125 mat[1][0]+=r[3*i+1]*f[3*i];
1126 mat[2][0]+=r[3*i+2]*f[3*i];
1127 mat[0][1]+=r[3*i] *f[3*i+1];
1128 mat[1][1]+=r[3*i+1]*f[3*i+1];
1129 mat[2][1]+=r[3*i+2]*f[3*i+1];
1130 mat[0][2]+=r[3*i] *f[3*i+2];
1131 mat[1][2]+=r[3*i+1]*f[3*i+2];
1132 mat[2][2]+=r[3*i+2]*f[3*i+2];
1133 }
1134
1135 d2=mat[0][0]*(mat[1][1]*mat[2][2]-mat[1][2]*mat[2][1])
1136 -mat[0][1]*(mat[1][0]*mat[2][2]-mat[1][2]*mat[2][0])
1137 +mat[0][2]*(mat[1][0]*mat[2][1]-mat[1][1]*mat[2][0]);
1138
1139 /* square matrix= ((mat transpose) * mat) */
1140 for(i=0;i<3;i++)
1141 for(j=0;j<3;j++)
1142 {
1143 x=mat[0][i]*mat[0][j]+mat[1][i]*mat[1][j]+mat[2][i]*mat[2][j];
1144 mat2[i][j]=mat[i][j];
1145 rmat[i][j]=x;
1146 }
1147 get_roots_3_3(rmat,roots);
1148
1149 roots[0]=(roots[0]<0.0001) ? 0.0: (roots[0]);
1150 roots[1]=(roots[1]<0.0001) ? 0.0: (roots[1]);
1151 roots[2]=(roots[2]<0.0001) ? 0.0: (roots[2]);
1152
1153 /* make sqrt of rmat, store in mat*/
1154
1155 roots[0]=(roots[0]<0.0001) ? 0.0: 1.0/(double)sqrt(roots[0]);
1156 roots[1]=(roots[1]<0.0001) ? 0.0: 1.0/(double)sqrt(roots[1]);
1157 roots[2]=(roots[2]<0.0001) ? 0.0: 1.0/(double)sqrt(roots[2]);
1158
1159 if(d2<0.0)
1160 {
1161 if( (roots[0]>=roots[1]) && (roots[0]>=roots[2]) )
1162 roots[0]*=-1.0;
1163 if( (roots[1]>roots[0]) && (roots[1]>=roots[2]) )
1164 roots[1]*=-1.0;
1165 if( (roots[2]>roots[1]) && (roots[2]>roots[0]) )
1166 roots[2]*=-1.0;
1167 }
1168
1169 for(i=0;i<3;i++)
1170 for(j=0;j<3;j++)
1171 mat[i][j]=roots[0]*rmat[i][0]*rmat[j][0]+
1172 roots[1]*rmat[i][1]*rmat[j][1]+
1173 roots[2]*rmat[i][2]*rmat[j][2];
1174
1175 /* and multiply into original inertial cross matrix, mat2 */
1176 for(i=0;i<3;i++)
1177 for(j=0;j<3;j++)
1178 rmat[i][j]=mat[0][j]*mat2[i][0]+
1179 mat[1][j]*mat2[i][1]+
1180 mat[2][j]*mat2[i][2];
1181
1182 rvec[0] = rmat[0][0];
1183 rvec[1] = rmat[0][1];
1184 rvec[2] = rmat[0][2];
1185 rvec[3] = rmat[1][0];
1186 rvec[4] = rmat[1][1];
1187 rvec[5] = rmat[1][2];
1188 rvec[6] = rmat[2][0];
1189 rvec[7] = rmat[2][1];
1190 rvec[8] = rmat[2][2];
1191 }
1192
1193 } // end namespace OpenBabel
1194
1195 //! \file obutil.cpp
1196 //! \brief Various utility methods.