Merge "Improve master-specific CMake behavior."
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_confrms.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "filenm.h"
40 #include "smalloc.h"
41 #include "macros.h"
42 #include <math.h>
43 #include "typedefs.h"
44 #include "xvgr.h"
45 #include "copyrite.h"
46 #include "statutil.h"
47 #include "tpxio.h"
48 #include "string2.h"
49 #include "vec.h"
50 #include "index.h"
51 #include "pbc.h"
52 #include "gmx_fatal.h"
53 #include "futil.h"
54 #include "confio.h"
55 #include "pdbio.h"
56 #include "txtdump.h"
57 #include "do_fit.h"
58 #include "viewit.h"
59 #include "rmpbc.h"
60 #include "gmx_ana.h"
61
62
63 void calc_rm_cm(int isize, atom_id index[], t_atoms *atoms, rvec x[], rvec xcm)
64 {
65     int  i, d;
66     real tm, m;
67
68     /* calculate and remove center of mass of reference structure */
69     tm = 0;
70     clear_rvec(xcm);
71     for (i = 0; i < isize; i++)
72     {
73         m = atoms->atom[index[i]].m;
74         for (d = 0; d < DIM; d++)
75         {
76             xcm[d] += m*x[index[i]][d];
77         }
78         tm += m;
79     }
80     svmul(1/tm, xcm, xcm);
81     for (i = 0; i < atoms->nr; i++)
82     {
83         rvec_dec(x[i], xcm);
84     }
85 }
86
87 int build_res_index(int isize, atom_id index[], t_atom atom[], int rindex[])
88 {
89     int i, r;
90
91     r         = 0;
92     rindex[r] = atom[index[0]].resind;
93     r++;
94     for (i = 1; i < isize; i++)
95     {
96         if (atom[index[i]].resind != rindex[r-1])
97         {
98             rindex[r] = atom[index[i]].resind;
99             r++;
100         }
101     }
102
103     return r;
104 }
105
106 int find_res_end(int i, int isize, atom_id index[], t_atoms *atoms)
107 {
108     int rnr;
109
110     rnr = atoms->atom[index[i]].resind;
111     while (i < isize && atoms->atom[index[i]].resind == rnr)
112     {
113         i++;
114     }
115     return i;
116 }
117
118 int debug_strcmp(char s1[], char s2[])
119 {
120     if (debug)
121     {
122         fprintf(debug, " %s-%s", s1, s2);
123     }
124     return strcmp(s1, s2);
125 }
126
127 int find_next_match_atoms_in_res(int *i1, int isize1, atom_id index1[],
128                                  int m1, char **atnms1[],
129                                  int *i2, int isize2, atom_id index2[],
130                                  int m2, char **atnms2[])
131 {
132     int      dx, dy, dmax, cmp;
133     gmx_bool bFW = FALSE;
134
135     dx   = dy = 0;
136     cmp  = NOTSET;
137     dmax = max(m1-*i1, m2-*i2);
138     for (dx = 0; dx < dmax && cmp != 0; dx++)
139     {
140         for (dy = dx; dy < dmax && cmp != 0; dy++)
141         {
142             if (dx || dy)
143             {
144                 if (debug)
145                 {
146                     fprintf(debug, ".");
147                 }
148                 cmp = NOTSET;
149                 if (*i1+dx < m1 && *i2+dy < m2)
150                 {
151                     bFW = TRUE;
152                     cmp = debug_strcmp(*atnms1[index1[*i1+dx]], *atnms2[index2[*i2+dy]]);
153                     if (debug)
154                     {
155                         fprintf(debug, "(%d %d)", *i1+dx, *i2+dy);
156                     }
157                 }
158                 if (cmp != 0 && *i1+dy < m1 && *i2+dx < m2)
159                 {
160                     bFW = FALSE;
161                     cmp = debug_strcmp(*atnms1[index1[*i1+dy]], *atnms2[index2[*i2+dx]]);
162                     if (debug)
163                     {
164                         fprintf(debug, "(%d %d)", *i1+dy, *i2+dx);
165                     }
166                 }
167             }
168         }
169     }
170     /* apparently, dx and dy are incremented one more time
171        as the loops terminate: we correct this here */
172     dx--;
173     dy--;
174     if (cmp == 0)
175     {
176         if (debug)
177         {
178             fprintf(debug, "{%d %d}", *i1+bFW ? dx : dy, *i2+bFW ? dy : dx );
179         }
180         if (bFW)
181         {
182             *i1 += dx;
183             *i2 += dy;
184         }
185         else
186         {
187             *i1 += dy;
188             *i2 += dx;
189         }
190     }
191
192     return cmp;
193 }
194
195 static int find_next_match_res(int *rnr1, int isize1,
196                                int index1[], t_resinfo *resinfo1,
197                                int *rnr2, int isize2,
198                                int index2[], t_resinfo *resinfo2)
199 {
200     int      dx, dy, dmax, cmp, rr1, rr2;
201     gmx_bool bFW = FALSE, bFF = FALSE;
202
203     dx  = dy = 0;
204     rr1 = 0;
205     while (rr1 < isize1 && *rnr1 != index1[rr1])
206     {
207         rr1++;
208     }
209     rr2 = 0;
210     while (rr2 < isize2 && *rnr2 != index2[rr2])
211     {
212         rr2++;
213     }
214
215     cmp  = NOTSET;
216     dmax = max(isize1-rr1, isize2-rr2);
217     if (debug)
218     {
219         fprintf(debug, " R:%d-%d:%d-%d:%d ",
220                 rr1, isize1, rr2, isize2, dmax);
221     }
222     for (dx = 0; dx < dmax && cmp != 0; dx++)
223     {
224         for (dy = 0; dy <= dx && cmp != 0; dy++)
225         {
226             if (dx != dy)
227             {
228                 cmp = NOTSET;
229                 if (rr1+dx < isize1 && rr2+dy < isize2)
230                 {
231                     bFW = TRUE;
232                     cmp = debug_strcmp(*resinfo1[index1[rr1+dx]].name,
233                                        *resinfo2[index2[rr2+dy]].name);
234                     if (debug)
235                     {
236                         fprintf(debug, "(%d %d)", rr1+dx, rr2+dy);
237                     }
238                 }
239                 if (cmp != 0 && rr1+dy < isize1 && rr2+dx < isize2)
240                 {
241                     bFW = FALSE;
242                     cmp = debug_strcmp(*resinfo1[index1[rr1+dy]].name,
243                                        *resinfo2[index2[rr2+dx]].name);
244                     if (debug)
245                     {
246                         fprintf(debug, "(%d %d)", rr1+dy, rr2+dx);
247                     }
248                 }
249                 if (dx != 0 && cmp != 0 && rr1+dx < isize1 && rr2+dx < isize2)
250                 {
251                     bFF = TRUE;
252                     cmp = debug_strcmp(*resinfo1[index1[rr1+dx]].name,
253                                        *resinfo2[index2[rr2+dx]].name);
254                     if (debug)
255                     {
256                         fprintf(debug, "(%d %d)", rr1+dx, rr2+dx);
257                     }
258                 }
259                 else
260                 {
261                     bFF = FALSE;
262                 }
263             }
264         }
265     }
266     /* apparently, dx and dy are incremented one more time
267        as the loops terminate: we correct this here */
268     dx--;
269     dy--;
270     /* if we skipped equal on both sides, only skip one residue
271        to allow for single mutations of residues... */
272     if (bFF)
273     {
274         if (debug)
275         {
276             fprintf(debug, "%d.%d.%dX%sX%s", dx, rr1, rr2,
277                     *resinfo1[index1[rr1+1]].name,
278                     *resinfo2[index2[rr2+1]].name);
279         }
280         dx = 1;
281     }
282     if (cmp == 0)
283     {
284         if (debug)
285         {
286             fprintf(debug, "!");
287         }
288         if (bFF)
289         {
290             rr1 += dx;
291             rr2 += dx;
292         }
293         else
294         if (bFW)
295         {
296             rr1 += dx;
297             rr2 += dy;
298         }
299         else
300         {
301             rr1 += dy;
302             rr2 += dx;
303         }
304         *rnr1 = index1[rr1];
305         *rnr2 = index2[rr2];
306     }
307
308     return cmp;
309 }
310
311 int find_first_atom_in_res(int rnr, int isize, atom_id index[], t_atom atom[])
312 {
313     int i;
314
315     i = 0;
316     while (i < isize && atom[index[i]].resind != rnr)
317     {
318         i++;
319     }
320
321     if (atom[index[i]].resind == rnr)
322     {
323         return i;
324     }
325     else
326     {
327         return NOTSET;
328     }
329 }
330
331 void find_matching_names(int *isize1, atom_id index1[], t_atoms *atoms1,
332                          int *isize2, atom_id index2[], t_atoms *atoms2)
333 {
334     int        i, i1, i2, ii1, ii2, m1, m2;
335     int        atcmp, rescmp;
336     int        r, rnr1, rnr2, prnr1, prnr2;
337     int        rsize1, rsize2;
338     int       *rindex1, *rindex2;
339     char      *resnm1, *resnm2, *atnm1, *atnm2;
340     char    ***atnms1, ***atnms2;
341     t_resinfo *resinfo1, *resinfo2;
342
343     /* set some handy shortcuts */
344     resinfo1 = atoms1->resinfo;
345     atnms1   = atoms1->atomname;
346     resinfo2 = atoms2->resinfo;
347     atnms2   = atoms2->atomname;
348
349     /* build indexes of selected residues */
350     snew(rindex1, atoms1->nres);
351     rsize1 = build_res_index(*isize1, index1, atoms1->atom, rindex1);
352     snew(rindex2, atoms2->nres);
353     rsize2 = build_res_index(*isize2, index2, atoms2->atom, rindex2);
354
355     i1    = i2 = 0;
356     ii1   = ii2 = 0;
357     atcmp = rescmp = 0;
358     prnr1 = prnr2 = NOTSET;
359     if (debug)
360     {
361         fprintf(debug, "Find matching names: %d, %d\n", *isize1, *isize2);
362     }
363     while (atcmp == 0 && i1 < *isize1 && i2 < *isize2)
364     {
365         /* prologue */
366         rnr1 = atoms1->atom[index1[i1]].resind;
367         rnr2 = atoms2->atom[index2[i2]].resind;
368         if (rnr1 != prnr1 || rnr2 != prnr2)
369         {
370             if (debug)
371             {
372                 fprintf(debug, "R: %s%d %s%d\n",
373                         *resinfo1[rnr1].name, rnr1, *resinfo2[rnr2].name, rnr2);
374             }
375             rescmp = strcmp(*resinfo1[rnr1].name, *resinfo2[rnr2].name);
376         }
377         if (debug)
378         {
379             fprintf(debug, "comparing %d %d", i1, i2);
380         }
381         atcmp = debug_strcmp(*atnms1[index1[i1]], *atnms2[index2[i2]]);
382
383         /* the works */
384         if (atcmp != 0) /* no match -> find match within residues */
385         {
386             m1 = find_res_end(i1, *isize1, index1, atoms1);
387             m2 = find_res_end(i2, *isize2, index2, atoms2);
388             if (debug)
389             {
390                 fprintf(debug, " [%d<%d %d<%d]", i1, m1, i2, m2);
391             }
392             atcmp = find_next_match_atoms_in_res(&i1, *isize1, index1, m1, atnms1,
393                                                  &i2, *isize2, index2, m2, atnms2);
394             if (debug)
395             {
396                 fprintf(debug, " -> %d %d %s-%s", i1, i2,
397                         *atnms1[index1[i1]], *atnms2[index2[i2]]);
398             }
399
400         }
401         if (atcmp != 0) /* still no match -> next residue(s) */
402         {
403             prnr1  = rnr1;
404             prnr2  = rnr2;
405             rescmp = find_next_match_res(&rnr1, rsize1, rindex1, resinfo1,
406                                          &rnr2, rsize2, rindex2, resinfo2);
407             if (rnr1 != prnr1)
408             {
409                 i1 = find_first_atom_in_res(rnr1, *isize1, index1, atoms1->atom);
410             }
411             if (rnr2 != prnr2)
412             {
413                 i2 = find_first_atom_in_res(rnr2, *isize2, index2, atoms2->atom);
414             }
415             if (debug)
416             {
417                 fprintf(debug, " -> %s%d-%s%d %s%d-%s%d",
418                         *resinfo1[rnr1].name, rnr1,
419                         *atnms1[index1[i1]], index1[i1],
420                         *resinfo2[rnr2].name, rnr2,
421                         *atnms2[index2[i2]], index2[i2]);
422             }
423             m1 = find_res_end(i1, *isize1, index1, atoms1);
424             m2 = find_res_end(i2, *isize2, index2, atoms2);
425             if (debug)
426             {
427                 fprintf(debug, " [%d<%d %d<%d]", i1, m1, i2, m2);
428             }
429             atcmp = find_next_match_atoms_in_res(&i1, *isize1, index1, m1, atnms1,
430                                                  &i2, *isize2, index2, m2, atnms2);
431             if (debug)
432             {
433                 fprintf(debug, " -> %d %d %s-%s", i1, i2,
434                         *atnms1[index1[i1]], *atnms2[index2[i2]]);
435             }
436         }
437         if (debug)
438         {
439             fprintf(debug, "(%d %d): %d %d\n", ii1, ii2, atcmp, rescmp);
440         }
441         if (atcmp == 0) /* if match -> store indices */
442         {
443             index1[ii1++] = index1[i1];
444             index2[ii2++] = index2[i2];
445         }
446         i1++;
447         i2++;
448
449         /* epilogue */
450         prnr1 = rnr1;
451         prnr2 = rnr2;
452     }
453
454     if (ii1 != ii2)
455     {
456         gmx_fatal(FARGS, "DEATH HORROR: non-equal number of matching atoms!\n");
457     }
458     if (ii1 == i1 && ii2 == i2)
459     {
460         printf("All atoms in index groups 1 and 2 match\n");
461     }
462     else
463     {
464         if (i1 == i2 && ii1 == ii2)
465         {
466             printf("Both index groups modified from %d to %d atoms\n", i1, ii1);
467         }
468         else
469         {
470             if (ii1 != i1)
471             {
472                 printf("Index group 1 modified from %d to %d atoms\n", i1, ii1);
473             }
474             if (ii2 != i2)
475             {
476                 printf("Index group 2 modified from %d to %d atoms\n", i2, ii2);
477             }
478         }
479         *isize1 = ii1;
480         *isize2 = ii2;
481     }
482 }
483 /* 1 */
484
485 int gmx_confrms(int argc, char *argv[])
486 {
487     const char     *desc[] = {
488         "[TT]g_confrms[tt] computes the root mean square deviation (RMSD) of two",
489         "structures after least-squares fitting the second structure on the first one.",
490         "The two structures do NOT need to have the same number of atoms,",
491         "only the two index groups used for the fit need to be identical.",
492         "With [TT]-name[tt] only matching atom names from the selected groups",
493         "will be used for the fit and RMSD calculation. This can be useful ",
494         "when comparing mutants of a protein.",
495         "[PAR]",
496         "The superimposed structures are written to file. In a [TT].pdb[tt] file",
497         "the two structures will be written as separate models",
498         "(use [TT]rasmol -nmrpdb[tt]). Also in a [TT].pdb[tt] file, B-factors",
499         "calculated from the atomic MSD values can be written with [TT]-bfac[tt].",
500     };
501     static gmx_bool bOne  = FALSE, bRmpbc = FALSE, bMW = TRUE, bName = FALSE,
502                     bBfac = FALSE, bFit = TRUE, bLabel = FALSE;
503
504     t_pargs  pa[] = {
505         { "-one", FALSE, etBOOL, {&bOne},   "Only write the fitted structure to file" },
506         { "-mw",  FALSE, etBOOL, {&bMW},    "Mass-weighted fitting and RMSD" },
507         { "-pbc", FALSE, etBOOL, {&bRmpbc}, "Try to make molecules whole again" },
508         { "-fit", FALSE, etBOOL, {&bFit},
509           "Do least squares superposition of the target structure to the reference" },
510         { "-name", FALSE, etBOOL, {&bName},
511           "Only compare matching atom names" },
512         { "-label", FALSE, etBOOL, {&bLabel},
513           "Added chain labels A for first and B for second structure"},
514         { "-bfac", FALSE, etBOOL, {&bBfac},
515           "Output B-factors from atomic MSD values" }
516     };
517     t_filenm fnm[] = {
518         { efTPS, "-f1",  "conf1.gro", ffREAD  },
519         { efSTX, "-f2",  "conf2",     ffREAD  },
520         { efSTO, "-o",   "fit.pdb",   ffWRITE },
521         { efNDX, "-n1", "fit1.ndx",  ffOPTRD },
522         { efNDX, "-n2", "fit2.ndx",  ffOPTRD },
523         { efNDX, "-no", "match.ndx", ffOPTWR }
524     };
525 #define NFILE asize(fnm)
526
527     /* the two structure files */
528     const char  *conf1file, *conf2file, *matchndxfile, *outfile;
529     FILE        *fp;
530     char         title1[STRLEN], title2[STRLEN], *name1, *name2;
531     t_topology  *top1, *top2;
532     int          ePBC1, ePBC2;
533     t_atoms     *atoms1, *atoms2;
534     int          warn = 0;
535     atom_id      at;
536     real        *w_rls, mass, totmass;
537     rvec        *x1, *v1, *x2, *v2, *fit_x;
538     matrix       box1, box2;
539
540     output_env_t oenv;
541
542     /* counters */
543     int     i, j, m;
544
545     /* center of mass calculation */
546     real    tmas1, tmas2;
547     rvec    xcm1, xcm2;
548
549     /* variables for fit */
550     char    *groupnames1, *groupnames2;
551     int      isize1, isize2;
552     atom_id *index1, *index2;
553     real     rms, msd, minmsd, maxmsd;
554     real    *msds;
555
556
557     parse_common_args(&argc, argv, PCA_BE_NICE | PCA_CAN_VIEW,
558                       NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
559     matchndxfile = opt2fn_null("-no", NFILE, fnm);
560     conf1file    = ftp2fn(efTPS, NFILE, fnm);
561     conf2file    = ftp2fn(efSTX, NFILE, fnm);
562
563     /* reading reference structure from first structure file */
564     fprintf(stderr, "\nReading first structure file\n");
565     snew(top1, 1);
566     read_tps_conf(conf1file, title1, top1, &ePBC1, &x1, &v1, box1, TRUE);
567     atoms1 = &(top1->atoms);
568     fprintf(stderr, "%s\nContaining %d atoms in %d residues\n",
569             title1, atoms1->nr, atoms1->nres);
570
571     if (bRmpbc)
572     {
573         rm_gropbc(atoms1, x1, box1);
574     }
575
576     fprintf(stderr, "Select group from first structure\n");
577     get_index(atoms1, opt2fn_null("-n1", NFILE, fnm),
578               1, &isize1, &index1, &groupnames1);
579     printf("\n");
580
581     if (bFit && (isize1 < 3))
582     {
583         gmx_fatal(FARGS, "Need >= 3 points to fit!\n");
584     }
585
586     /* reading second structure file */
587     fprintf(stderr, "\nReading second structure file\n");
588     snew(top2, 1);
589     read_tps_conf(conf2file, title2, top2, &ePBC2, &x2, &v2, box2, TRUE);
590     atoms2 = &(top2->atoms);
591     fprintf(stderr, "%s\nContaining %d atoms in %d residues\n",
592             title2, atoms2->nr, atoms2->nres);
593
594     if (bRmpbc)
595     {
596         rm_gropbc(atoms2, x2, box2);
597     }
598
599     fprintf(stderr, "Select group from second structure\n");
600     get_index(atoms2, opt2fn_null("-n2", NFILE, fnm),
601               1, &isize2, &index2, &groupnames2);
602
603     if (bName)
604     {
605         find_matching_names(&isize1, index1, atoms1, &isize2, index2, atoms2);
606         if (matchndxfile)
607         {
608             fp = ffopen(matchndxfile, "w");
609             fprintf(fp, "; Matching atoms between %s from %s and %s from %s\n",
610                     groupnames1, conf1file, groupnames2, conf2file);
611             fprintf(fp, "[ Match_%s_%s ]\n", conf1file, groupnames1);
612             for (i = 0; i < isize1; i++)
613             {
614                 fprintf(fp, "%4u%s", index1[i]+1, (i%15 == 14 || i == isize1-1) ? "\n" : " ");
615             }
616             fprintf(fp, "[ Match_%s_%s ]\n", conf2file, groupnames2);
617             for (i = 0; i < isize2; i++)
618             {
619                 fprintf(fp, "%4u%s", index2[i]+1, (i%15 == 14 || i == isize2-1) ? "\n" : " ");
620             }
621         }
622     }
623
624     /* check isizes, must be equal */
625     if (isize2 != isize1)
626     {
627         gmx_fatal(FARGS, "You selected groups with differen number of atoms.\n");
628     }
629
630     for (i = 0; i < isize1; i++)
631     {
632         name1 = *atoms1->atomname[index1[i]];
633         name2 = *atoms2->atomname[index2[i]];
634         if (strcmp(name1, name2))
635         {
636             if (warn < 20)
637             {
638                 fprintf(stderr,
639                         "Warning: atomnames at index %d don't match: %u %s, %u %s\n",
640                         i+1, index1[i]+1, name1, index2[i]+1, name2);
641             }
642             warn++;
643         }
644         if (!bMW)
645         {
646             atoms1->atom[index1[i]].m = 1;
647             atoms2->atom[index2[i]].m = 1;
648         }
649     }
650     if (warn)
651     {
652         fprintf(stderr, "%d atomname%s did not match\n", warn, (warn == 1) ? "" : "s");
653     }
654
655     if (bFit)
656     {
657         /* calculate and remove center of mass of structures */
658         calc_rm_cm(isize1, index1, atoms1, x1, xcm1);
659         calc_rm_cm(isize2, index2, atoms2, x2, xcm2);
660
661         snew(w_rls, atoms2->nr);
662         snew(fit_x, atoms2->nr);
663         for (at = 0; (at < isize1); at++)
664         {
665             w_rls[index2[at]] = atoms1->atom[index1[at]].m;
666             copy_rvec(x1[index1[at]], fit_x[index2[at]]);
667         }
668
669         /* do the least squares fit to the reference structure */
670         do_fit(atoms2->nr, w_rls, fit_x, x2);
671
672         sfree(fit_x);
673         sfree(w_rls);
674         w_rls = NULL;
675     }
676     else
677     {
678         clear_rvec(xcm1);
679         clear_rvec(xcm2);
680         w_rls = NULL;
681     }
682
683     /* calculate the rms deviation */
684     rms     = 0;
685     totmass = 0;
686     maxmsd  = -1e18;
687     minmsd  =  1e18;
688     snew(msds, isize1);
689     for (at = 0; at < isize1; at++)
690     {
691         mass = atoms1->atom[index1[at]].m;
692         for (m = 0; m < DIM; m++)
693         {
694             msd       = sqr(x1[index1[at]][m] - x2[index2[at]][m]);
695             rms      += msd*mass;
696             msds[at] += msd;
697         }
698         maxmsd   = max(maxmsd, msds[at]);
699         minmsd   = min(minmsd, msds[at]);
700         totmass += mass;
701     }
702     rms = sqrt(rms/totmass);
703
704     printf("Root mean square deviation after lsq fit = %g nm\n", rms);
705     if (bBfac)
706     {
707         printf("Atomic MSD's range from %g to %g nm^2\n", minmsd, maxmsd);
708     }
709
710     if (bFit)
711     {
712         /* reset coordinates of reference and fitted structure */
713         for (i = 0; i < atoms1->nr; i++)
714         {
715             for (m = 0; m < DIM; m++)
716             {
717                 x1[i][m] += xcm1[m];
718             }
719         }
720         for (i = 0; i < atoms2->nr; i++)
721         {
722             for (m = 0; m < DIM; m++)
723             {
724                 x2[i][m] += xcm1[m];
725             }
726         }
727     }
728
729     outfile = ftp2fn(efSTO, NFILE, fnm);
730     switch (fn2ftp(outfile))
731     {
732         case efPDB:
733         case efBRK:
734         case efENT:
735             if (bBfac || bLabel)
736             {
737                 srenew(atoms1->pdbinfo, atoms1->nr);
738                 srenew(atoms1->atom, atoms1->nr); /* Why renew atom? */
739
740                 /* Avoid segfaults when writing the pdb-file */
741                 for (i = 0; i < atoms1->nr; i++)
742                 {
743                     atoms1->pdbinfo[i].type         = eptAtom;
744                     atoms1->pdbinfo[i].occup        = 1.00;
745                     atoms1->pdbinfo[i].bAnisotropic = FALSE;
746                     if (bBfac)
747                     {
748                         atoms1->pdbinfo[i].bfac = 0;
749                     }
750                     if (bLabel)
751                     {
752                         atoms1->resinfo[atoms1->atom[i].resind].chainid = 'A';
753                     }
754                 }
755
756                 for (i = 0; i < isize1; i++)
757                 {
758                     /* atoms1->pdbinfo[index1[i]].type = eptAtom; */
759 /*  atoms1->pdbinfo[index1[i]].bAnisotropic = FALSE; */
760                     if (bBfac)
761                     {
762                         atoms1->pdbinfo[index1[i]].bfac = (800*M_PI*M_PI/3.0)*msds[i];
763                     }
764 /*  if (bLabel) */
765 /*    atoms1->resinfo[atoms1->atom[index1[i]].resind].chain = 'A'; */
766                 }
767                 srenew(atoms2->pdbinfo, atoms2->nr);
768                 srenew(atoms2->atom, atoms2->nr); /* Why renew atom? */
769
770                 for (i = 0; i < atoms2->nr; i++)
771                 {
772                     atoms2->pdbinfo[i].type         = eptAtom;
773                     atoms2->pdbinfo[i].occup        = 1.00;
774                     atoms2->pdbinfo[i].bAnisotropic = FALSE;
775                     if (bBfac)
776                     {
777                         atoms2->pdbinfo[i].bfac = 0;
778                     }
779                     if (bLabel)
780                     {
781                         atoms2->resinfo[atoms1->atom[i].resind].chainid = 'B';
782                     }
783                 }
784
785                 for (i = 0; i < isize2; i++)
786                 {
787                     /* atoms2->pdbinfo[index2[i]].type = eptAtom; */
788 /*  atoms2->pdbinfo[index2[i]].bAnisotropic = FALSE; */
789                     if (bBfac)
790                     {
791                         atoms2->pdbinfo[index2[i]].bfac = (800*M_PI*M_PI/3.0)*msds[i];
792                     }
793 /*  if (bLabel) */
794 /*    atoms2->resinfo[atoms2->atom[index2[i]].resind].chain = 'B'; */
795                 }
796             }
797             fp = ffopen(outfile, "w");
798             if (!bOne)
799             {
800                 write_pdbfile(fp, title1, atoms1, x1, ePBC1, box1, ' ', 1, NULL, TRUE);
801             }
802             write_pdbfile(fp, title2, atoms2, x2, ePBC2, box2, ' ', bOne ? -1 : 2, NULL, TRUE);
803             ffclose(fp);
804             break;
805         case efGRO:
806             if (bBfac)
807             {
808                 fprintf(stderr, "WARNING: cannot write B-factor values to gro file\n");
809             }
810             fp = ffopen(outfile, "w");
811             if (!bOne)
812             {
813                 write_hconf_p(fp, title1, atoms1, 3, x1, v1, box1);
814             }
815             write_hconf_p(fp, title2, atoms2, 3, x2, v2, box2);
816             ffclose(fp);
817             break;
818         default:
819             if (bBfac)
820             {
821                 fprintf(stderr, "WARNING: cannot write B-factor values to %s file\n",
822                         ftp2ext(fn2ftp(outfile)));
823             }
824             if (!bOne)
825             {
826                 fprintf(stderr,
827                         "WARNING: cannot write the reference structure to %s file\n",
828                         ftp2ext(fn2ftp(outfile)));
829             }
830             write_sto_conf(outfile, title2, atoms2, x2, v2, ePBC2, box2);
831             break;
832     }
833
834     view_all(oenv, NFILE, fnm);
835
836     thanx(stderr);
837
838     return 0;
839 }