Partial commit of the project to remove all static variables.
[gromacs.git] / src / tools / smooth.c
bloba1305d481c28a4411cf1e7fc0df7d051c40255bc
1 /*
2 * $Id$
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.1
11 * Copyright (c) 1991-2001, University of Groningen, The Netherlands
12 * This program is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU General Public License
14 * as published by the Free Software Foundation; either version 2
15 * of the License, or (at your option) any later version.
17 * If you want to redistribute modifications, please consider that
18 * scientific software is very special. Version control is crucial -
19 * bugs must be traceable. We will be happy to consider code for
20 * inclusion in the official distribution, but derived work must not
21 * be called official GROMACS. Details are found in the README & COPYING
22 * files - if they are missing, get the official version at www.gromacs.org.
24 * To help us fund GROMACS development, we humbly ask that you cite
25 * the papers on the package - you can find them in the top README file.
27 * For more info, check our website at http://www.gromacs.org
29 * And Hey:
30 * Green Red Orange Magenta Azure Cyan Skyblue
33 #include "cdist.h"
35 #define NMRLEN 99.9
36 #define OVERLAP(a_lb,a_ub,b_lb,b_ub) !((a_ub)<(b_lb)) || ((b_ub)<(a_lb))
37 /* Bound smoothing for cdist*/
38 /* Adam Kirrander 990119 */
40 real _mysqrt(real x,char *fn,int line)
42 if (x < 0) {
43 if (debug)
44 fprintf(debug,"MYSQRT: negative argument %g called from %s, line %d\n",
45 x,fn,line);
46 return 0;
48 else
49 return sqrt(x);
52 #define mysqrt(x) _mysqrt(x,__FILE__,__LINE__)
54 void check_bounds(int ai,int aj,real lb,real ub)
56 if (lb > ub)
57 fatal_error(0,"Strange bounds for %d-%d: lb=%f, ub=%f",
58 ai+1,aj+1,lb,ub);
61 void check_len_(int ai,int aj,real lb,real ub,real len,char *file,int line)
63 if (len != NMRLEN)
64 if (((len > ub) || (len < lb)) && len > 0.0)
65 fatal_error(0,"%s, line %d: Ideal len. outside bounds for %d-%d:"
66 " lb=%f, ub=%f, len=%f",
67 file,line,ai+1,aj+1,lb,ub,len);
70 #define check_len(ai,aj,lb,ub,len) check_len_(ai,aj,lb,ub,len,__FILE__,__LINE__)
72 void check_tri(t_dist *d,int natoms,int i,int j,int k)
74 real a_lb,a_ub,a_len,b_lb,b_ub,b_len,c_lb,c_ub,c_len;
76 a_lb = d_lb(d,natoms,i,j);
77 a_ub = d_ub(d,natoms,i,j);
78 a_len = d_len(d,natoms,i,j);
80 b_lb = d_lb(d,natoms,i,k);
81 b_ub = d_ub(d,natoms,i,k);
82 b_len = d_len(d,natoms,i,k);
84 c_lb = d_lb(d,natoms,j,k);
85 c_ub = d_ub(d,natoms,j,k);
86 c_len = d_len(d,natoms,j,k);
88 check_len(i,j,a_lb,a_ub,a_len);
89 check_len(i,k,b_lb,b_ub,b_len);
90 check_len(j,k,c_lb,c_ub,c_len);
92 check_bounds(i,j,a_lb,a_ub);
93 check_bounds(i,k,b_lb,b_ub);
94 check_bounds(j,k,c_lb,c_ub);
96 /* Check if triangle inequality is fulfilled */
97 if ((a_lb > (b_ub+c_ub)) ||
98 (b_lb > (a_ub+c_ub)) ||
99 (c_lb > (a_ub+b_ub))) {
100 fprintf(stderr,"Impossible triangle (%d,%d,%d):\n"\
101 "a_lb=%10.5f, a_ub=%10.5f, a_len=%10.5f\n"\
102 "b_lb=%10.5f, b_ub=%10.5f, b_len=%10.5f\n"
103 "c_lb=%10.5f, c_ub=%10.5f, c_len=%10.5f\n",\
104 i,j,k,a_lb,a_ub,a_len,b_lb,b_ub,b_len,c_lb,c_ub,c_len);
105 exit(1);
109 double triangle_upper_bound (t_dist *d,int natoms,real tol)
111 double possible,ntotal=0;
112 int i,j,k,nsmooth,innerloop,count=0,nloops=0;
113 real a_lb,a_ub,a_len,b_lb,b_ub,b_len,c_lb,c_ub,c_len;
115 if (debug)
116 fprintf(debug,"Just entered triangle_upper_bound!\n");
118 /* Loop over all triangles until no smoothing occurs */
119 do {
120 nsmooth=0;
121 nloops=nloops+1;
122 for (i=0; (i < natoms); i++) {
123 for (j=i+1; (j < natoms); j++) {
124 if (!dist_set(d,natoms,i,j)) continue;
125 for (k=j+1; (k < natoms); k++) {
127 count=count+1;
129 /* Should make use of distances between zero-weighted if that can
130 help us to define better distances for e.g. loops. But make sure
131 we don't waste time smoothing zero-weighted tripples...*/
133 /* Are all distances defined? */
134 if (!dist_set(d,natoms,i,k) || !dist_set(d,natoms,j,k)) continue;
136 /* smooth the triangle */
137 do {
138 /* Reset inner counter */
139 innerloop=0;
141 check_tri(d,natoms,i,j,k);
143 /* Read bounds */
144 a_lb = d_lb(d,natoms,i,j);
145 a_ub = d_ub(d,natoms,i,j);
146 a_len = d_len(d,natoms,i,j);
148 b_lb = d_lb(d,natoms,i,k);
149 b_ub = d_ub(d,natoms,i,k);
150 b_len = d_len(d,natoms,i,k);
152 c_lb = d_lb(d,natoms,j,k);
153 c_ub = d_ub(d,natoms,j,k);
154 c_len = d_len(d,natoms,j,k);
156 /* Check if triangle inequality is fulfilled */
157 if (a_ub > (b_ub+c_ub+tol)) {
158 set_dist(d,natoms,i,j,a_lb,b_ub+c_ub,a_len);
159 a_ub = d_ub(d,natoms,i,j);
160 innerloop++;
162 if (b_ub > (a_ub+c_ub+tol)) {
163 set_dist(d,natoms,i,k,b_lb,a_ub+c_ub,b_len);
164 b_ub = d_ub(d,natoms,i,k);
165 innerloop++;
167 if (c_ub > (a_ub+b_ub+tol)) {
168 set_dist(d,natoms,j,k,c_lb,a_ub+b_ub,c_len);
169 c_ub = d_ub(d,natoms,j,k);
170 innerloop++;
173 nsmooth += innerloop;
175 while (innerloop > 0);
179 ntotal += nsmooth;
181 while (nsmooth>0);
183 possible = natoms*(natoms-1.0)*(natoms-2.0)/6;
184 fprintf(stderr,"Checked %d triples (of %.0f) in %d rounds of ub"
185 " triangle smoothing.\n",
186 count/nloops,possible,nloops);
187 fprintf(stderr,"Smoothed %g upper bounds with triagonal ineq.\n",ntotal);
189 return ntotal;
193 double triangle_lower_bound (t_dist *d,int natoms,real tol)
195 double ntotal = 0;
196 int i,j,k,nsmooth,innerloop;
197 real a_lb,a_ub,a_len,b_lb,b_ub,b_len,c_lb,c_ub,c_len,new_lb;
199 /*fprintf(stderr,"Just entered triangle_lower_bound!\n");*/
201 /* Loop over all triangles until no smoothing occurs */
202 do {
203 nsmooth=0;
204 for (i=0; (i < natoms); i++) {
205 for (j=i+1; (j < natoms); j++) {
206 if (!dist_set(d,natoms,i,j)) continue;
207 for (k=j+1; (k < natoms); k++) {
209 /* Are all distances defined? If not, continue */
210 if (!dist_set(d,natoms,i,k) || !dist_set(d,natoms,j,k)) continue;
213 /* smooth the triangle */
214 do {
215 /* Reset inner counter */
216 innerloop=0;
218 /* Read bounds */
219 a_lb = d_lb(d,natoms,i,j);
220 a_ub = d_ub(d,natoms,i,j);
221 a_len = d_len(d,natoms,i,j);
223 b_lb = d_lb(d,natoms,i,k);
224 b_ub = d_ub(d,natoms,i,k);
225 b_len = d_len(d,natoms,i,k);
227 c_lb = d_lb(d,natoms,j,k);
228 c_ub = d_ub(d,natoms,j,k);
229 c_len = d_len(d,natoms,j,k);
231 /* Smooth lower bound */
232 if (!OVERLAP(b_lb,b_ub,c_lb,c_ub)) {
233 new_lb = min(fabs(b_lb-c_ub),fabs(c_lb-b_ub));
234 if ( tol < new_lb - a_lb) {
235 set_dist(d,natoms,i,j,new_lb,a_ub,a_len);
236 nsmooth++;
237 innerloop++;
238 /*fprintf(stderr,"Smoothed %d-%d",i,j);*/
241 if (!OVERLAP(a_lb,a_ub,c_lb,c_ub)) {
242 new_lb = min(fabs(a_lb-c_ub),fabs(c_lb-a_ub));
243 if ( tol < new_lb - b_lb) {
244 set_dist(d,natoms,i,k,new_lb,b_ub,b_len);
245 nsmooth++;
246 innerloop++;
247 /*fprintf(stderr,"Smoothed %d-%d",i,k);*/
250 if (!OVERLAP(a_lb,a_ub,b_lb,b_ub)) {
251 new_lb = min(fabs(a_lb-b_ub),fabs(b_lb-a_ub));
252 /*Programming error? */
253 /* if ( (a_ub+b_ub) < c_ub ) { */
254 /* New if-statement 990609 Adam */
255 if ( tol < new_lb - c_lb) {
256 set_dist(d,natoms,j,k,new_lb,c_ub,c_len);
257 nsmooth++;
258 innerloop++;
259 /*fprintf(stderr,"Smoothed %d-%d",j,k);*/
263 /* Check that bounds make sense */
264 check_bounds(i,j,a_lb,a_ub);
265 check_bounds(i,k,b_lb,b_ub);
266 check_bounds(j,k,c_lb,c_ub);
268 while (innerloop > 0);
272 ntotal += nsmooth;
274 while (nsmooth>0);
275 fprintf(stderr,"Smoothed %g lower bounds with triagonal ineq.\n",ntotal);
276 return ntotal;
279 double do_triangle (t_dist *d,t_atoms *atoms,real tol)
281 double ntot=0;
282 /* Do triangular smoothing only */
283 int natoms = atoms->nr;
285 fprintf(stderr,"Starting triangle smoothing\n");
286 ntot += triangle_upper_bound (d,natoms,tol);
287 ntot += triangle_lower_bound (d,natoms,tol);
289 return ntot;
292 /* Tetrangle-bound smoothing for cdist.
293 Adam Kirrander 990209 */
295 /* According Havel,Kuntz and Crippen 1983 Bull.Math.Biol */
300 /* FUNCTIONS ASSOCIATED WITH TRIANGLE_LIMITED */
301 /* Routines to check if triangle inequalty limits are attainable between
302 points 1 and 4, given that the 5 other distances are defined.
303 Sets logical TUL=true if the upper limit of distance 1-4 is defined by
304 triangle inequalty, and TLL=true if corresponding lower limit is set by
305 the triangle ineq. Written to be used in conjunction with tetragonal
306 bound smoothing in cdist.
308 The main routine is trianglelimited, and it is the one that should be
309 called. The rest are supporting routines.
310 Ref.: Bull.Math.Biol.Vol45,Nr5,665-720,Havel,Crippen,Kuntz,1983 */
312 /* Adam Kirrander 990209 */
314 /*#define OVERLAP(a_lb,a_ub,b_lb,b_ub) !((a_ub<b_lb) || (b_ub<a_lb))*/
315 /* ((a_ub >= b_lb) && (b_ub >= a_lb)) */
316 #define apex(AC,AD,BC,BD) mysqrt((BC*(AD*AD-BD*BD)+BD*(AC*AC-BC*BC))/(BC+BD))
319 bool _triangle_bound_attainable(real BC,real BD,real AClb,real ACub,
320 real ADlb,real ADub,real ABlb,real ABub,
321 char *fn,int line)
323 /* Assume triangle with A,C and D at corners and B at C-D.
324 Given the distances AC,AD,BC,BD it calculates the distance AB
325 at its min and max, which corresponds to AC and AD being min and
326 max respectively.
327 If this intervall overlaps with the set bounds for AB, it follows that
328 the points C,B and D can be made colinear and thus are restricted by the
329 triangle inequality. */
331 real ABmin,ABmax;
333 if (debug)
334 fprintf(debug,"%s, line %d: ABlb: %g, ABub: %g, AClb: %g, ACub: %g, ADlb: %g, ADub: %g, BC: %g, BD: %g\n",
335 fn,line,ABlb,ABub,AClb,ACub,ADlb,ADub,BC,BD);
336 ABmax = apex(ACub,ADub,BC,BD);
337 /* Triangle can't be formed, i.e. there is a risk all 4 points can
338 become colinear */
339 if ( (AClb+ADlb) < (BC+BD) )
340 return ( ABlb <= ABmax );
342 /* I am not sure the above truly is the optimal solution, but it should
343 * be the safe solution
346 ABmin = apex(AClb,ADlb,BC,BD);
348 return (OVERLAP(ABlb,ABub,ABmin,ABmax));
351 #define triangle_bound_attainable(BC,BD,AClb,ACub,ADlb,ADub,ABlb,ABub) \
352 _triangle_bound_attainable(BC,BD,AClb,ACub,ADlb,ADub,ABlb,ABub, \
353 __FILE__,__LINE__)
356 void triangle_limited(int *ATMS,int natoms,t_dist *d,bool *TUL, bool *TLL)
358 /* Given 4 atoms it checks if the triangle inequality lower or upper bounds
359 for the distance 1-4 are attainable. */
360 /* Situation looks something like the following: */
361 /* Can: */
362 /* 1 4 */
363 /* \\\ //// */
364 /* 2------------3 */
365 /* become: */
366 /* 1---------------3------4 */
367 /* \\\\\\\ /////// */
368 /* 2 */
371 int N1=ATMS[0],N2=ATMS[1],N3=ATMS[2],N4=ATMS[3];
372 real d12[2],d13[2],d14[2],d23[2],d24[2],d34[2];
375 /* Initiate and read bounds */
376 *TUL = FALSE;
377 *TLL = FALSE;
379 d12[0] = d_lb(d,natoms,N1,N2);
380 d12[1] = d_ub(d,natoms,N1,N2);
382 d13[0] = d_lb(d,natoms,N1,N3);
383 d13[1] = d_ub(d,natoms,N1,N3);
385 d14[0] = d_lb(d,natoms,N1,N4);
386 d14[1] = d_ub(d,natoms,N1,N4);
388 d23[0] = d_lb(d,natoms,N2,N3);
389 d23[1] = d_ub(d,natoms,N2,N3);
391 d24[0] = d_lb(d,natoms,N2,N4);
392 d24[1] = d_ub(d,natoms,N2,N4);
394 d34[0] = d_lb(d,natoms,N3,N4);
395 d34[1] = d_ub(d,natoms,N3,N4);
397 /* Check if UPPER triangle inequality limit is attainable. */
398 if ( d12[1] + d24[1] < d13[1] + d34[1] ) {
399 /* Corresponds to N1,N2 and N4 colinear.
400 N1=D,N2=B,N3=A,N4=C ;in notation of subroutines */
401 *TUL = triangle_bound_attainable(d24[1],d12[1],d34[0],d34[1],
402 d13[0],d13[1],d23[0],d23[1]);
404 else if ( d12[1] + d24[1] > d13[1] + d34[1] ) {
405 /* Corresponds to N1,N3 and N4 colinear.
406 N1=D,N2=A,N3=B,N4=C ;in notation of subroutines */
407 *TUL = triangle_bound_attainable(d34[1],d13[1],d24[0],d24[1],
408 d12[0],d12[1],d23[0],d23[1]);
410 /* if N2 and N3 can superimpose TUL is true by necessity (AB=0) */
411 else if (d23[0] == 0) *TUL = TRUE;
413 /* Check if LOWER triangle inequality limit is attainable */
414 if ( (d13[0] - d34[1] <= 0) && (d34[0] - d13[1] <= 0) &&
415 (d24[0] - d12[1] <= 0) && (d12[0] - d24[1] <= 0) ) {
416 /* if all inverse triangle ineq. limits are zero then */
417 *TLL = TRUE;
419 else if (d12[0] - d24[1] > 0) {
420 /* Corresponds to N1,N4 and N2 colinear.
421 A=N3,B=N4,C=N1,D=N2 */
422 *TLL = triangle_bound_attainable(d24[1],d12[0]-d24[1],d23[0],d23[1],
423 d13[0],d13[1],d34[0],d34[1]);
425 else if (d24[0] - d12[1] > 0) {
426 /* Corresponds to N2,N1 and N4 colinear.
427 A=N3,B=N1,C=N2,D=N4 */
428 *TLL = triangle_bound_attainable(d12[1],d24[0]-d12[1],d23[0],d23[1],
429 d34[0],d34[1],d13[0],d13[1]);
431 else if (d13[0] - d34[1] > 0) {
432 /* Corresponds to N1,N4 and N3 colinear.
433 A=N2,B=N4,C=N3,D=N1 */
434 *TLL = triangle_bound_attainable(d34[1],d13[0]-d34[1],d23[0],d23[1],
435 d12[0],d12[1],d24[0],d24[1]);
437 else if (d34[0] - d13[1] > 0) {
438 /* Corresponds to N3,N1 and N4 colinear.
439 A=N2,B=N1,C=N3,D=N4 */
440 *TLL = triangle_bound_attainable(d13[1],d34[0]-d13[1],d23[0],d23[1],
441 d24[0],d24[1],d12[0],d12[1]);
445 /* END OF FUNCTIONS >TRIANGLE_LIMITED< */
449 /* FUNCTIONS ASSOCIATED WITH TETRAGONAL-SMOOTHING */
451 void minmax(real array[],int SIZE,int *min,int *max)
453 /* Finds min and max (indices thereof) in array.
454 Early version that only accepts real.
455 Fix it!
456 Adam Kirrander 990211 */
458 int i;
459 *min = 0;
460 *max = 0;
462 for (i=1; (i<SIZE); i++) {
463 if (array[i] < array[*min])
464 *min = i;
465 else if (array[i] > array[*max])
466 *max = i;
470 void swap (int *x,int *y)
472 /* Swaps x and y.
473 Early version that only can swap integers.
474 Fix it!
475 Adam Kirrander 990211 */
477 int temp;
478 temp = *x;
479 *x = *y;
480 *y = temp;
484 bool alldistset (int i,int j,int k,int l,int natoms,t_dist *d)
486 /* Returns FALSE if all distances are not set */
487 if (!dist_set(d,natoms,i,j)) return FALSE;
488 if (!dist_set(d,natoms,i,k)) return FALSE;
489 if (!dist_set(d,natoms,i,l)) return FALSE;
490 if (!dist_set(d,natoms,j,k)) return FALSE;
491 if (!dist_set(d,natoms,j,l)) return FALSE;
492 if (!dist_set(d,natoms,k,l)) return FALSE;
494 return TRUE;
497 int nr_nb (int i,int j,int k,int l,int natoms,t_dist *d)
499 /* Counts the number of nb's */
500 int nnb=0;
502 if (d_len(d,natoms,i,j) == 0.0) nnb++;
503 if (d_len(d,natoms,i,k) == 0.0) nnb++;
504 if (d_len(d,natoms,i,l) == 0.0) nnb++;
505 if (d_len(d,natoms,j,k) == 0.0) nnb++;
506 if (d_len(d,natoms,j,l) == 0.0) nnb++;
507 if (d_len(d,natoms,k,l) == 0.0) nnb++;
509 return nnb;
512 void sort_atoms(int natoms,t_dist *d,int *ATMS)
514 /* Sorts atoms. The nb-atoms end up in ATMS[0] and ATMS[3]. */
516 if (d_len(d,natoms,ATMS[0],ATMS[1]) == 0.0)
517 swap(&ATMS[1],&ATMS[3]);
518 else if (d_len(d,natoms,ATMS[0],ATMS[2]) == 0.0)
519 swap(&ATMS[2],&ATMS[3]);
520 else if (d_len(d,natoms,ATMS[1],ATMS[2]) == 0.0) {
521 swap(&ATMS[0],&ATMS[1]);
522 swap(&ATMS[2],&ATMS[3]);
524 else if (d_len(d,natoms,ATMS[1],ATMS[3]) == 0.0)
525 swap(&ATMS[1],&ATMS[0]);
526 else if (d_len(d,natoms,ATMS[2],ATMS[3]) == 0.0)
527 swap(&ATMS[2],&ATMS[0]);
529 /* Put the two middle ones in order (ambivalent procedure, necessary?) */
530 if (d_len(d,natoms,ATMS[0],ATMS[1]) > d_len(d,natoms,ATMS[0],ATMS[2]))
531 swap(&ATMS[1],&ATMS[2]);
534 real solve_tetrangle(real dAB,real dAC,real dBC,real dBD,real dCD,int cosphi)
536 /* Solves tetrangle eq. for distance AD.
537 cosphi=cos(phi), where phi is the dihedral angle
538 cosphi=1 corresponds to min, cosphi=-1 to max
539 eq. solved and optimized by Maple, watch out for bugs */
541 real t1,t2,t3,t4,t5,t7,t8,t9,t11,t12,t13,t14,t15,t17,t20;
543 t1 = dAB*dAB;
544 t2 = dBC*dBC;
545 t3 = t1*t2;
546 t4 = t1*t1;
547 t5 = dAC*dAC;
548 t7 = t2*t2;
549 t8 = t5*t2;
550 t9 = t5*t5;
551 t11 = dCD*dCD;
552 t12 = t2*t11;
553 t13 = dBD*dBD;
554 t14 = t2*t13;
555 t15 = t11*t11;
556 t17 = t13*t13;
557 t20 = mysqrt((-2.0*t3+t4-2.0*t1*t5+t7-2.0*t8+t9)*(-2.0*t12+t7-2.0*t14+t15
558 -2.0*t11*t13+t17));
559 return mysqrt(-(cosphi*t20-t8-t14+t7-t3-t1*t11+t1*t13-t12+t5*t11-t5*t13)/(2*t2));
560 /* t20 = mysqrt((-2.0*(t3+t1*t5+t8) + t4 + t7 + t9)*
561 (-2.0*(t12+t14+t11*t13) + t7 + t15 + t17));
562 return mysqrt(-0.5*(cosphi*t20-t8-t14+t7-t3-t1*t11+t1*t13-t12+t5*t11-t5*t13))/dBC;*/
566 int tetrangle_bound(int *ATMS,int cosphi,int natoms,t_dist *d,real tol)
568 int sol=0,D1,D2,D3,D4,D5,nsmooth=0,dmin,dmax;
569 real dAD[32],dAB,dAC,dBC,dBD,dCD,present_lb,present_ub,present_len;
571 /*fprintf(stderr,"entering tetrangle_bound\n");*/
573 /* Try all 32 combinations of bounds */
574 for (D1=0; (D1<2); D1++) {
575 for (D2=0; (D2<2); D2++) {
576 for (D3=0; (D3<2); D3++) {
577 for (D4=0; (D4<2); D4++) {
578 for (D5=0; (D5<2); D5++) {
579 dAB = (D1) ? d_lb(d,natoms,ATMS[0],ATMS[1]):
580 d_ub(d,natoms,ATMS[0],ATMS[1]);
581 dAC = (D2) ? d_lb(d,natoms,ATMS[0],ATMS[2]):
582 d_ub(d,natoms,ATMS[0],ATMS[2]);
583 dBC = (D3) ? d_lb(d,natoms,ATMS[1],ATMS[2]):
584 d_ub(d,natoms,ATMS[1],ATMS[2]);
585 dBD = (D4) ? d_lb(d,natoms,ATMS[1],ATMS[3]):
586 d_ub(d,natoms,ATMS[1],ATMS[3]);
587 dCD = (D5) ? d_lb(d,natoms,ATMS[2],ATMS[3]):
588 d_ub(d,natoms,ATMS[2],ATMS[3]);
590 dAD[sol] = solve_tetrangle(dAB,dAC,dBC,dBD,dCD,cosphi);
591 if (debug)
592 fprintf(debug,"dAD[%d]=%10.5f\n",sol,dAD[sol]);
593 sol++;
600 /* we need these to re-set one of the bounds for the distance */
601 present_len = d_len(d,natoms,ATMS[0],ATMS[3]);
602 present_lb = d_lb(d,natoms,ATMS[0],ATMS[3]);
603 present_ub = d_ub(d,natoms,ATMS[0],ATMS[3]);
605 /* Set the new bound(s) */
606 minmax(dAD,32,&dmin,&dmax);
607 if (debug)
608 fprintf(debug,"Max=dAD[%d] (%12g), Min=dAD[%d] (%12g)\n",
609 dmax,dAD[dmax],dmin,dAD[dmin]);
610 /* Set new upper limit */
611 if ((cosphi == -1) && ((present_ub - dAD[dmax]) > tol) ) {
612 set_dist(d,natoms,ATMS[0],ATMS[3],present_lb,dAD[dmax],present_len);
613 if (debug)
614 fprintf(debug,"Corrected %d-%d-%d-%d ub to %10.5f\n",ATMS[0]+1,
615 ATMS[1]+1,ATMS[2]+1,ATMS[3]+1,dAD[dmax]);
616 nsmooth++;
618 /* Set new lower limit */
619 else if ((cosphi == 1) && ( tol < (dAD[dmin] - present_lb))) {
620 set_dist(d,natoms,ATMS[0],ATMS[3],dAD[dmin],present_ub,present_len);
621 if (debug)
622 fprintf(debug,"Corrected %d-%d-%d-%d lb to %10.5f\n",ATMS[0]+1,
623 ATMS[1]+1,ATMS[2]+1,ATMS[3]+1,dAD[dmin]);
624 nsmooth++;
627 /* Check the new bounds */
628 if (d_ub(d,natoms,ATMS[0],ATMS[3]) < d_lb(d,natoms,ATMS[0],ATMS[3]))
629 fatal_error(0,"Fatal error for tetr. smooth of (%d,%d) ub=%f, lb=%f",
630 ATMS[0],ATMS[3],d_ub(d,natoms,ATMS[0],ATMS[3]),
631 d_lb(d,natoms,ATMS[0],ATMS[3]));
633 return nsmooth;
638 int tetrangle (t_dist *d,int natoms,real tol)
641 int i,j,k,l,AT[4],nsmooth=0,nnonbonds1=0,nnonbonds2=0;
642 bool TLL=FALSE,TUL=FALSE;
644 fprintf(stderr,"Starting tetrangle smoothing\n");
645 /* Loops over all sets of four atoms */
646 for (i=0;i<(natoms-3);i++) {
647 for (j=i+1;j<(natoms-2);j++) {
648 if (!dist_set(d,natoms,i,j)) continue;
649 for (k=j+1;k<(natoms-1);k++) {
650 if (!dist_set(d,natoms,i,k)) continue;
651 if (!dist_set(d,natoms,j,k)) continue;
652 nnonbonds1 = 0;
653 if (d_len(d,natoms,i,j) == 0.0) nnonbonds1++;
654 if (d_len(d,natoms,i,k) == 0.0) nnonbonds1++;
655 if (d_len(d,natoms,j,k) == 0.0) nnonbonds1++;
656 if (nnonbonds1 > 1) continue;
657 for (l=k+1;l<natoms;l++) {
658 if (!dist_set(d,natoms,i,l)) continue;
659 if (!dist_set(d,natoms,j,l)) continue;
660 if (!dist_set(d,natoms,k,l)) continue;
661 nnonbonds2 = 0;
662 if (d_len(d,natoms,i,l) == 0.0) nnonbonds2++;
663 if (d_len(d,natoms,j,l) == 0.0) nnonbonds2++;
664 if (d_len(d,natoms,k,l) == 0.0) nnonbonds2++;
665 if ( (nnonbonds1+nnonbonds2) != 1) continue;
667 /*fprintf(stderr,"Trying %d-%d-%d-%d\n",i+1,j+1,k+1,l+1);*/
669 /* The two following subroutines have been nested into the loops
670 above. Doesn't look as good, but should be substantially
671 faster. */
672 /* if (!alldistset(i,j,k,l,natoms,d)) continue; */
673 /* if (nr_nb(i,j,k,l,natoms,d) != 1) continue; */
675 /* Copy indices to array AT for easier handling & sort them */
676 AT[0] = i;
677 AT[1] = j;
678 AT[2] = k;
679 AT[3] = l;
680 /*fprintf(stderr,"OK'd atoms (%d,%d,%d,%d)\n",i+1,j+1,k+1,l+1);*/
681 sort_atoms(natoms,d,AT);
682 /*fprintf(stderr,"After sorting (%d,%d,%d,%d)\n",AT[0]+1,AT[1]+1,AT[2]+1,AT[3]+1);*/
684 /* Are the distances limited by the triangle-ineq.? */
685 triangle_limited(AT,natoms,d,&TUL,&TLL);
686 /*if (TUL) fprintf(stderr,"TUL=true\n");
687 else fprintf(stderr,"TUL=false\n");
688 if (TLL) fprintf(stderr,"TLL=true\n");
689 else fprintf(stderr,"TLL=false\n"); */
691 /* If not, correct bounds according to tetrangle-ineq. */
692 /* upper limit */
693 if (!TUL) nsmooth += tetrangle_bound(AT,-1,natoms,d,tol);
694 /* lower limit */
695 if (!TLL) nsmooth += tetrangle_bound(AT,1,natoms,d,tol);
697 /* Comment: If we were to calculate upper and lower bound
698 simultaneously we could make the code a bit faster */
703 return nsmooth;
706 /* END OF FUNCTIONS >TETRANGLE-SMOOTH< */
710 void do_smooth (t_dist *d,t_atoms *atoms,real tol)
712 /* Triangular and tetragonal bound smoothing.
713 Do we need to nestle triangular and tetragonal
714 smoothing or can we do them separately? THINK! */
716 /* In this case I think we can get away with this, because
717 we are not using distances that we correct with tetragonal
718 ineq. to smooth other distances. The triangle smoothing
719 should then spread the tighter bounds to other distances,
720 not corrected by the teragonal ineq.
721 (Of course, if you correct one distance by tetragonal ineq.
722 it will become shorter, and other distances that are involved
723 in triagonal ineq. with it will be affected.)
726 /* It should only do 2 loops, and in the second no tetragonal
727 correction will be made. This is because (in the current version)
728 we only correct nb's that are connected by bonded distances, and
729 these are unlikely to change during triagonal smoothing. */
731 /* If this turns out to be true => SKIP DO-LOOP AND JUST HAVE
732 TRIAGONAL-TETRAGONAL-TRIAGONAL */
734 int ntri;
735 int nloops = 0,nsmooth,natoms = atoms->nr;
737 if (debug)
738 fprintf(debug,"Starting tetragonal smoothing routine\n");
740 /* ntri = do_triangle(d,atoms,tol);
741 do {
742 nloops += 1;
744 nsmooth = tetrangle (d,natoms,tol);
745 fprintf(stderr,"Smoothed %d distances with tetr. ineq.\n",nsmooth);
746 if (nsmooth > 0)
747 ntri = do_triangle(d,atoms,tol);
748 else
749 ntri = 0;
751 fprintf(stderr,"Finished %d loop(s) of smoothing.\n",nloops);
752 } while (ntri > 0);
754 ntri = do_triangle(d,atoms,tol);
755 nsmooth = tetrangle (d,natoms,tol);
756 fprintf(stderr,"Smoothed %d distances with tetr. ineq.\n",nsmooth);
757 ntri = do_triangle(d,atoms,tol);
759 fprintf(stderr,"Finished smoothing.\n");