Remove repeat sfree(state->nosehoover_xi)
[alexxy/gromacs.git] / src / gromacs / gmxlib / typedefs.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
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.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 /* This file is completely threadsafe - keep it that way! */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include "smalloc.h"
42 #include "symtab.h"
43 #include "vec.h"
44 #include "pbc.h"
45 #include "macros.h"
46 #include <string.h>
47
48 #ifdef GMX_THREAD_MPI
49 #include "thread_mpi.h"
50 #endif
51
52 /* The source code in this file should be thread-safe.
53       Please keep it that way. */
54
55
56
57 static gmx_bool            bOverAllocDD = FALSE;
58 #ifdef GMX_THREAD_MPI
59 static tMPI_Thread_mutex_t over_alloc_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
60 #endif
61
62
63 void set_over_alloc_dd(gmx_bool set)
64 {
65 #ifdef GMX_THREAD_MPI
66     tMPI_Thread_mutex_lock(&over_alloc_mutex);
67     /* we just make sure that we don't set this at the same time.
68        We don't worry too much about reading this rarely-set variable */
69 #endif
70     bOverAllocDD = set;
71 #ifdef GMX_THREAD_MPI
72     tMPI_Thread_mutex_unlock(&over_alloc_mutex);
73 #endif
74 }
75
76 int over_alloc_dd(int n)
77 {
78     if (bOverAllocDD)
79     {
80         return OVER_ALLOC_FAC*n + 100;
81     }
82     else
83     {
84         return n;
85     }
86 }
87
88 int gmx_large_int_to_int(gmx_large_int_t step, const char *warn)
89 {
90     int i;
91
92     i = (int)step;
93
94     if (warn != NULL && (step < INT_MIN || step > INT_MAX))
95     {
96         fprintf(stderr, "\nWARNING during %s:\n", warn);
97         fprintf(stderr, "step value ");
98         fprintf(stderr, gmx_large_int_pfmt, step);
99         fprintf(stderr, " does not fit in int, converted to %d\n\n", i);
100     }
101
102     return i;
103 }
104
105 char *gmx_step_str(gmx_large_int_t i, char *buf)
106 {
107     sprintf(buf, gmx_large_int_pfmt, i);
108
109     return buf;
110 }
111
112 void init_block(t_block *block)
113 {
114     int i;
115
116     block->nr           = 0;
117     block->nalloc_index = 1;
118     snew(block->index, block->nalloc_index);
119     block->index[0]     = 0;
120 }
121
122 void init_blocka(t_blocka *block)
123 {
124     int i;
125
126     block->nr           = 0;
127     block->nra          = 0;
128     block->nalloc_index = 1;
129     snew(block->index, block->nalloc_index);
130     block->index[0]     = 0;
131     block->nalloc_a     = 0;
132     block->a            = NULL;
133 }
134
135 void init_atom(t_atoms *at)
136 {
137     int i;
138
139     at->nr        = 0;
140     at->nres      = 0;
141     at->atom      = NULL;
142     at->resinfo   = NULL;
143     at->atomname  = NULL;
144     at->atomtype  = NULL;
145     at->atomtypeB = NULL;
146     at->pdbinfo   = NULL;
147 }
148
149 void init_atomtypes(t_atomtypes *at)
150 {
151     at->nr         = 0;
152     at->radius     = NULL;
153     at->vol        = NULL;
154     at->atomnumber = NULL;
155     at->gb_radius  = NULL;
156     at->S_hct      = NULL;
157 }
158
159 void init_groups(gmx_groups_t *groups)
160 {
161     int g;
162
163     groups->ngrpname = 0;
164     groups->grpname  = NULL;
165     for (g = 0; (g < egcNR); g++)
166     {
167         groups->grps[g].nm_ind = NULL;
168         groups->ngrpnr[g]      = 0;
169         groups->grpnr[g]       = NULL;
170     }
171
172 }
173
174 void init_mtop(gmx_mtop_t *mtop)
175 {
176     mtop->name         = NULL;
177     mtop->nmoltype     = 0;
178     mtop->moltype      = NULL;
179     mtop->nmolblock    = 0;
180     mtop->molblock     = NULL;
181     mtop->maxres_renum = 0;
182     mtop->maxresnr     = -1;
183     init_groups(&mtop->groups);
184     init_block(&mtop->mols);
185     open_symtab(&mtop->symtab);
186 }
187
188 void init_top (t_topology *top)
189 {
190     int i;
191
192     top->name = NULL;
193     init_atom (&(top->atoms));
194     init_atomtypes(&(top->atomtypes));
195     init_block(&top->cgs);
196     init_block(&top->mols);
197     init_blocka(&top->excls);
198     open_symtab(&top->symtab);
199 }
200
201 void init_inputrec(t_inputrec *ir)
202 {
203     memset(ir, 0, (size_t)sizeof(*ir));
204     snew(ir->fepvals, 1);
205     snew(ir->expandedvals, 1);
206     snew(ir->simtempvals, 1);
207 }
208
209 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup)
210 {
211     int i;
212
213     if (bOneIndexGroup)
214     {
215         grp->nalloc_index = 2;
216         snew(grp->index, grp->nalloc_index);
217         grp->index[0] = 0;
218         grp->index[1] = natom;
219         grp->nr       = 1;
220     }
221     else
222     {
223         grp->nalloc_index = natom+1;
224         snew(grp->index, grp->nalloc_index);
225         snew(grp->index, natom+1);
226         for (i = 0; (i <= natom); i++)
227         {
228             grp->index[i] = i;
229         }
230         grp->nr = natom;
231     }
232 }
233
234 void stupid_fill_blocka(t_blocka *grp, int natom)
235 {
236     int i;
237
238     grp->nalloc_a = natom;
239     snew(grp->a, grp->nalloc_a);
240     for (i = 0; (i < natom); i++)
241     {
242         grp->a[i] = i;
243     }
244     grp->nra = natom;
245
246     grp->nalloc_index = natom + 1;
247     snew(grp->index, grp->nalloc_index);
248     for (i = 0; (i <= natom); i++)
249     {
250         grp->index[i] = i;
251     }
252     grp->nr = natom;
253 }
254
255 void copy_blocka(const t_blocka *src, t_blocka *dest)
256 {
257     int i;
258
259     dest->nr           = src->nr;
260     dest->nalloc_index = dest->nr + 1;
261     snew(dest->index, dest->nalloc_index);
262     for (i = 0; i < dest->nr+1; i++)
263     {
264         dest->index[i] = src->index[i];
265     }
266     dest->nra      = src->nra;
267     dest->nalloc_a = dest->nra + 1;
268     snew(dest->a, dest->nalloc_a);
269     for (i = 0; i < dest->nra+1; i++)
270     {
271         dest->a[i] = src->a[i];
272     }
273 }
274
275 void done_block(t_block *block)
276 {
277     block->nr    = 0;
278     sfree(block->index);
279     block->nalloc_index = 0;
280 }
281
282 void done_blocka(t_blocka *block)
283 {
284     block->nr    = 0;
285     block->nra   = 0;
286     sfree(block->index);
287     sfree(block->a);
288     block->index        = NULL;
289     block->a            = NULL;
290     block->nalloc_index = 0;
291     block->nalloc_a     = 0;
292 }
293
294 void done_atom (t_atoms *at)
295 {
296     at->nr       = 0;
297     at->nres     = 0;
298     sfree(at->atom);
299     sfree(at->resinfo);
300     sfree(at->atomname);
301     sfree(at->atomtype);
302     sfree(at->atomtypeB);
303     if (at->pdbinfo)
304     {
305         sfree(at->pdbinfo);
306     }
307 }
308
309 void done_atomtypes(t_atomtypes *atype)
310 {
311     atype->nr = 0;
312     sfree(atype->radius);
313     sfree(atype->vol);
314     sfree(atype->surftens);
315     sfree(atype->atomnumber);
316     sfree(atype->gb_radius);
317     sfree(atype->S_hct);
318 }
319
320 void done_moltype(gmx_moltype_t *molt)
321 {
322     int f;
323
324     done_atom(&molt->atoms);
325     done_block(&molt->cgs);
326     done_blocka(&molt->excls);
327
328     for (f = 0; f < F_NRE; f++)
329     {
330         sfree(molt->ilist[f].iatoms);
331         molt->ilist[f].nalloc = 0;
332     }
333 }
334
335 void done_molblock(gmx_molblock_t *molb)
336 {
337     if (molb->nposres_xA > 0)
338     {
339         molb->nposres_xA = 0;
340         free(molb->posres_xA);
341     }
342     if (molb->nposres_xB > 0)
343     {
344         molb->nposres_xB = 0;
345         free(molb->posres_xB);
346     }
347 }
348
349 void done_mtop(gmx_mtop_t *mtop, gmx_bool bDoneSymtab)
350 {
351     int i;
352
353     if (bDoneSymtab)
354     {
355         done_symtab(&mtop->symtab);
356     }
357
358     sfree(mtop->ffparams.functype);
359     sfree(mtop->ffparams.iparams);
360
361     for (i = 0; i < mtop->nmoltype; i++)
362     {
363         done_moltype(&mtop->moltype[i]);
364     }
365     sfree(mtop->moltype);
366     for (i = 0; i < mtop->nmolblock; i++)
367     {
368         done_molblock(&mtop->molblock[i]);
369     }
370     sfree(mtop->molblock);
371     done_block(&mtop->mols);
372 }
373
374 void done_top(t_topology *top)
375 {
376     int f;
377
378     sfree(top->idef.functype);
379     sfree(top->idef.iparams);
380     for (f = 0; f < F_NRE; ++f)
381     {
382         sfree(top->idef.il[f].iatoms);
383         top->idef.il[f].iatoms = NULL;
384         top->idef.il[f].nalloc = 0;
385     }
386
387     done_atom (&(top->atoms));
388
389     /* For GB */
390     done_atomtypes(&(top->atomtypes));
391
392     done_symtab(&(top->symtab));
393     done_block(&(top->cgs));
394     done_block(&(top->mols));
395     done_blocka(&(top->excls));
396 }
397
398 static void done_pullgrp(t_pullgrp *pgrp)
399 {
400     sfree(pgrp->ind);
401     sfree(pgrp->ind_loc);
402     sfree(pgrp->weight);
403     sfree(pgrp->weight_loc);
404 }
405
406 static void done_pull(t_pull *pull)
407 {
408     int i;
409
410     for (i = 0; i < pull->ngrp+1; i++)
411     {
412         done_pullgrp(pull->grp);
413         done_pullgrp(pull->dyna);
414     }
415 }
416
417 void done_inputrec(t_inputrec *ir)
418 {
419     int m;
420
421     for (m = 0; (m < DIM); m++)
422     {
423         if (ir->ex[m].a)
424         {
425             sfree(ir->ex[m].a);
426         }
427         if (ir->ex[m].phi)
428         {
429             sfree(ir->ex[m].phi);
430         }
431         if (ir->et[m].a)
432         {
433             sfree(ir->et[m].a);
434         }
435         if (ir->et[m].phi)
436         {
437             sfree(ir->et[m].phi);
438         }
439     }
440
441     sfree(ir->opts.nrdf);
442     sfree(ir->opts.ref_t);
443     sfree(ir->opts.annealing);
444     sfree(ir->opts.anneal_npoints);
445     sfree(ir->opts.anneal_time);
446     sfree(ir->opts.anneal_temp);
447     sfree(ir->opts.tau_t);
448     sfree(ir->opts.acc);
449     sfree(ir->opts.nFreeze);
450     sfree(ir->opts.QMmethod);
451     sfree(ir->opts.QMbasis);
452     sfree(ir->opts.QMcharge);
453     sfree(ir->opts.QMmult);
454     sfree(ir->opts.bSH);
455     sfree(ir->opts.CASorbitals);
456     sfree(ir->opts.CASelectrons);
457     sfree(ir->opts.SAon);
458     sfree(ir->opts.SAoff);
459     sfree(ir->opts.SAsteps);
460     sfree(ir->opts.bOPT);
461     sfree(ir->opts.bTS);
462
463     if (ir->pull)
464     {
465         done_pull(ir->pull);
466         sfree(ir->pull);
467     }
468 }
469
470 static void zero_ekinstate(ekinstate_t *eks)
471 {
472     eks->ekin_n         = 0;
473     eks->ekinh          = NULL;
474     eks->ekinf          = NULL;
475     eks->ekinh_old      = NULL;
476     eks->ekinscalef_nhc = NULL;
477     eks->ekinscaleh_nhc = NULL;
478     eks->vscale_nhc     = NULL;
479     eks->dekindl        = 0;
480     eks->mvcos          = 0;
481 }
482
483 void init_energyhistory(energyhistory_t * enerhist)
484 {
485     enerhist->nener = 0;
486
487     enerhist->ener_ave     = NULL;
488     enerhist->ener_sum     = NULL;
489     enerhist->ener_sum_sim = NULL;
490     enerhist->dht          = NULL;
491
492     enerhist->nsteps     = 0;
493     enerhist->nsum       = 0;
494     enerhist->nsteps_sim = 0;
495     enerhist->nsum_sim   = 0;
496
497     enerhist->dht = NULL;
498 }
499
500 static void done_delta_h_history(delta_h_history_t *dht)
501 {
502     int i;
503
504     for (i = 0; i < dht->nndh; i++)
505     {
506         sfree(dht->dh[i]);
507     }
508     sfree(dht->dh);
509     sfree(dht->ndh);
510 }
511
512 void done_energyhistory(energyhistory_t * enerhist)
513 {
514     sfree(enerhist->ener_ave);
515     sfree(enerhist->ener_sum);
516     sfree(enerhist->ener_sum_sim);
517
518     if (enerhist->dht != NULL)
519     {
520         done_delta_h_history(enerhist->dht);
521         sfree(enerhist->dht);
522     }
523 }
524
525 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
526 {
527     int i, j;
528
529     state->ngtc          = ngtc;
530     state->nnhpres       = nnhpres;
531     state->nhchainlength = nhchainlength;
532     if (state->ngtc > 0)
533     {
534         snew(state->nosehoover_xi, state->nhchainlength*state->ngtc);
535         snew(state->nosehoover_vxi, state->nhchainlength*state->ngtc);
536         snew(state->therm_integral, state->ngtc);
537         for (i = 0; i < state->ngtc; i++)
538         {
539             for (j = 0; j < state->nhchainlength; j++)
540             {
541                 state->nosehoover_xi[i*state->nhchainlength + j]   = 0.0;
542                 state->nosehoover_vxi[i*state->nhchainlength + j]  = 0.0;
543             }
544         }
545         for (i = 0; i < state->ngtc; i++)
546         {
547             state->therm_integral[i]  = 0.0;
548         }
549     }
550     else
551     {
552         state->nosehoover_xi  = NULL;
553         state->nosehoover_vxi = NULL;
554         state->therm_integral = NULL;
555     }
556
557     if (state->nnhpres > 0)
558     {
559         snew(state->nhpres_xi, state->nhchainlength*nnhpres);
560         snew(state->nhpres_vxi, state->nhchainlength*nnhpres);
561         for (i = 0; i < nnhpres; i++)
562         {
563             for (j = 0; j < state->nhchainlength; j++)
564             {
565                 state->nhpres_xi[i*nhchainlength + j]   = 0.0;
566                 state->nhpres_vxi[i*nhchainlength + j]  = 0.0;
567             }
568         }
569     }
570     else
571     {
572         state->nhpres_xi  = NULL;
573         state->nhpres_vxi = NULL;
574     }
575 }
576
577
578 void init_state(t_state *state, int natoms, int ngtc, int nnhpres, int nhchainlength, int nlambda)
579 {
580     int i;
581
582     state->natoms = natoms;
583     state->nrng   = 0;
584     state->flags  = 0;
585     state->lambda = 0;
586     snew(state->lambda, efptNR);
587     for (i = 0; i < efptNR; i++)
588     {
589         state->lambda[i] = 0;
590     }
591     state->veta   = 0;
592     clear_mat(state->box);
593     clear_mat(state->box_rel);
594     clear_mat(state->boxv);
595     clear_mat(state->pres_prev);
596     clear_mat(state->svir_prev);
597     clear_mat(state->fvir_prev);
598     init_gtc_state(state, ngtc, nnhpres, nhchainlength);
599     state->nalloc = state->natoms;
600     if (state->nalloc > 0)
601     {
602         snew(state->x, state->nalloc);
603         snew(state->v, state->nalloc);
604     }
605     else
606     {
607         state->x = NULL;
608         state->v = NULL;
609     }
610     state->sd_X = NULL;
611     state->cg_p = NULL;
612
613     zero_ekinstate(&state->ekinstate);
614
615     init_energyhistory(&state->enerhist);
616
617     init_df_history(&state->dfhist, nlambda, 0);
618
619     state->ddp_count       = 0;
620     state->ddp_count_cg_gl = 0;
621     state->cg_gl           = NULL;
622     state->cg_gl_nalloc    = 0;
623 }
624
625 void done_state(t_state *state)
626 {
627     if (state->x)
628     {
629         sfree(state->x);
630     }
631     if (state->v)
632     {
633         sfree(state->v);
634     }
635     if (state->sd_X)
636     {
637         sfree(state->sd_X);
638     }
639     if (state->cg_p)
640     {
641         sfree(state->cg_p);
642     }
643     state->nalloc = 0;
644     if (state->cg_gl)
645     {
646         sfree(state->cg_gl);
647     }
648     state->cg_gl_nalloc = 0;
649     if (state->lambda)
650     {
651         sfree(state->lambda);
652     }
653     if (state->ngtc > 0)
654     {
655         sfree(state->nosehoover_xi);
656         sfree(state->nosehoover_vxi);
657         sfree(state->therm_integral);
658     }
659 }
660
661 static void do_box_rel(t_inputrec *ir, matrix box_rel, matrix b, gmx_bool bInit)
662 {
663     int d, d2;
664
665     for (d = YY; d <= ZZ; d++)
666     {
667         for (d2 = XX; d2 <= (ir->epct == epctSEMIISOTROPIC ? YY : ZZ); d2++)
668         {
669             /* We need to check if this box component is deformed
670              * or if deformation of another component might cause
671              * changes in this component due to box corrections.
672              */
673             if (ir->deform[d][d2] == 0 &&
674                 !(d == ZZ && d2 == XX && ir->deform[d][YY] != 0 &&
675                   (b[YY][d2] != 0 || ir->deform[YY][d2] != 0)))
676             {
677                 if (bInit)
678                 {
679                     box_rel[d][d2] = b[d][d2]/b[XX][XX];
680                 }
681                 else
682                 {
683                     b[d][d2] = b[XX][XX]*box_rel[d][d2];
684                 }
685             }
686         }
687     }
688 }
689
690 void set_box_rel(t_inputrec *ir, t_state *state)
691 {
692     /* Make sure the box obeys the restrictions before we fix the ratios */
693     correct_box(NULL, 0, state->box, NULL);
694
695     clear_mat(state->box_rel);
696
697     if (PRESERVE_SHAPE(*ir))
698     {
699         do_box_rel(ir, state->box_rel, state->box, TRUE);
700     }
701 }
702
703 void preserve_box_shape(t_inputrec *ir, matrix box_rel, matrix b)
704 {
705     if (PRESERVE_SHAPE(*ir))
706     {
707         do_box_rel(ir, box_rel, b, FALSE);
708     }
709 }
710
711 void add_t_atoms(t_atoms *atoms, int natom_extra, int nres_extra)
712 {
713     int i;
714
715     if (natom_extra > 0)
716     {
717         srenew(atoms->atomname, atoms->nr+natom_extra);
718         srenew(atoms->atom, atoms->nr+natom_extra);
719         if (NULL != atoms->pdbinfo)
720         {
721             srenew(atoms->pdbinfo, atoms->nr+natom_extra);
722         }
723         if (NULL != atoms->atomtype)
724         {
725             srenew(atoms->atomtype, atoms->nr+natom_extra);
726         }
727         if (NULL != atoms->atomtypeB)
728         {
729             srenew(atoms->atomtypeB, atoms->nr+natom_extra);
730         }
731         for (i = atoms->nr; (i < atoms->nr+natom_extra); i++)
732         {
733             atoms->atomname[i] = NULL;
734             memset(&atoms->atom[i], 0, sizeof(atoms->atom[i]));
735             if (NULL != atoms->pdbinfo)
736             {
737                 memset(&atoms->pdbinfo[i], 0, sizeof(atoms->pdbinfo[i]));
738             }
739             if (NULL != atoms->atomtype)
740             {
741                 atoms->atomtype[i] = NULL;
742             }
743             if (NULL != atoms->atomtypeB)
744             {
745                 atoms->atomtypeB[i] = NULL;
746             }
747         }
748         atoms->nr += natom_extra;
749     }
750     if (nres_extra > 0)
751     {
752         srenew(atoms->resinfo, atoms->nres+nres_extra);
753         for (i = atoms->nres; (i < atoms->nres+nres_extra); i++)
754         {
755             memset(&atoms->resinfo[i], 0, sizeof(atoms->resinfo[i]));
756         }
757         atoms->nres += nres_extra;
758     }
759 }
760
761 void init_t_atoms(t_atoms *atoms, int natoms, gmx_bool bPdbinfo)
762 {
763     atoms->nr   = natoms;
764     atoms->nres = 0;
765     snew(atoms->atomname, natoms);
766     atoms->atomtype  = NULL;
767     atoms->atomtypeB = NULL;
768     snew(atoms->resinfo, natoms);
769     snew(atoms->atom, natoms);
770     if (bPdbinfo)
771     {
772         snew(atoms->pdbinfo, natoms);
773     }
774     else
775     {
776         atoms->pdbinfo = NULL;
777     }
778 }
779
780 t_atoms *copy_t_atoms(t_atoms *src)
781 {
782     t_atoms *dst;
783     int      i;
784
785     snew(dst, 1);
786     init_t_atoms(dst, src->nr, (NULL != src->pdbinfo));
787     dst->nr = src->nr;
788     if (NULL != src->atomname)
789     {
790         snew(dst->atomname, src->nr);
791     }
792     if (NULL != src->atomtype)
793     {
794         snew(dst->atomtype, src->nr);
795     }
796     if (NULL != src->atomtypeB)
797     {
798         snew(dst->atomtypeB, src->nr);
799     }
800     for (i = 0; (i < src->nr); i++)
801     {
802         dst->atom[i] = src->atom[i];
803         if (NULL != src->pdbinfo)
804         {
805             dst->pdbinfo[i] = src->pdbinfo[i];
806         }
807         if (NULL != src->atomname)
808         {
809             dst->atomname[i]  = src->atomname[i];
810         }
811         if (NULL != src->atomtype)
812         {
813             dst->atomtype[i] = src->atomtype[i];
814         }
815         if (NULL != src->atomtypeB)
816         {
817             dst->atomtypeB[i] = src->atomtypeB[i];
818         }
819     }
820     dst->nres = src->nres;
821     for (i = 0; (i < src->nres); i++)
822     {
823         dst->resinfo[i] = src->resinfo[i];
824     }
825     return dst;
826 }
827
828 void t_atoms_set_resinfo(t_atoms *atoms, int atom_ind, t_symtab *symtab,
829                          const char *resname, int resnr, unsigned char ic,
830                          int chainnum, char chainid)
831 {
832     t_resinfo *ri;
833
834     ri           = &atoms->resinfo[atoms->atom[atom_ind].resind];
835     ri->name     = put_symtab(symtab, resname);
836     ri->rtp      = NULL;
837     ri->nr       = resnr;
838     ri->ic       = ic;
839     ri->chainnum = chainnum;
840     ri->chainid  = chainid;
841 }
842
843 void free_t_atoms(t_atoms *atoms, gmx_bool bFreeNames)
844 {
845     int i;
846
847     if (bFreeNames && atoms->atomname != NULL)
848     {
849         for (i = 0; i < atoms->nr; i++)
850         {
851             if (atoms->atomname[i] != NULL)
852             {
853                 sfree(*atoms->atomname[i]);
854                 *atoms->atomname[i] = NULL;
855             }
856         }
857     }
858     if (bFreeNames && atoms->resinfo != NULL)
859     {
860         for (i = 0; i < atoms->nres; i++)
861         {
862             if (atoms->resinfo[i].name != NULL)
863             {
864                 sfree(*atoms->resinfo[i].name);
865                 *atoms->resinfo[i].name = NULL;
866             }
867         }
868     }
869     if (bFreeNames && atoms->atomtype != NULL)
870     {
871         for (i = 0; i < atoms->nr; i++)
872         {
873             if (atoms->atomtype[i] != NULL)
874             {
875                 sfree(*atoms->atomtype[i]);
876                 *atoms->atomtype[i] = NULL;
877             }
878         }
879     }
880     if (bFreeNames && atoms->atomtypeB != NULL)
881     {
882         for (i = 0; i < atoms->nr; i++)
883         {
884             if (atoms->atomtypeB[i] != NULL)
885             {
886                 sfree(*atoms->atomtypeB[i]);
887                 *atoms->atomtypeB[i] = NULL;
888             }
889         }
890     }
891     sfree(atoms->atomname);
892     sfree(atoms->atomtype);
893     sfree(atoms->atomtypeB);
894     sfree(atoms->resinfo);
895     sfree(atoms->atom);
896     sfree(atoms->pdbinfo);
897     atoms->nr        = 0;
898     atoms->nres      = 0;
899     atoms->atomname  = NULL;
900     atoms->atomtype  = NULL;
901     atoms->atomtypeB = NULL;
902     atoms->resinfo   = NULL;
903     atoms->atom      = NULL;
904     atoms->pdbinfo   = NULL;
905 }
906
907 real max_cutoff(real cutoff1, real cutoff2)
908 {
909     if (cutoff1 == 0 || cutoff2 == 0)
910     {
911         return 0;
912     }
913     else
914     {
915         return max(cutoff1, cutoff2);
916     }
917 }
918
919 extern void init_df_history(df_history_t *dfhist, int nlambda, real wl_delta)
920 {
921     int i;
922
923     dfhist->bEquil   = 0;
924     dfhist->nlambda  = nlambda;
925     dfhist->wl_delta = wl_delta;
926     snew(dfhist->sum_weights, dfhist->nlambda);
927     snew(dfhist->sum_dg, dfhist->nlambda);
928     snew(dfhist->sum_minvar, dfhist->nlambda);
929     snew(dfhist->sum_variance, dfhist->nlambda);
930     snew(dfhist->n_at_lam, dfhist->nlambda);
931     snew(dfhist->wl_histo, dfhist->nlambda);
932
933     /* allocate transition matrices here */
934     snew(dfhist->Tij, dfhist->nlambda);
935     snew(dfhist->Tij_empirical, dfhist->nlambda);
936
937     for (i = 0; i < dfhist->nlambda; i++)
938     {
939         snew(dfhist->Tij[i], dfhist->nlambda);
940         snew(dfhist->Tij_empirical[i], dfhist->nlambda);
941     }
942
943     snew(dfhist->accum_p, dfhist->nlambda);
944     snew(dfhist->accum_m, dfhist->nlambda);
945     snew(dfhist->accum_p2, dfhist->nlambda);
946     snew(dfhist->accum_m2, dfhist->nlambda);
947
948     for (i = 0; i < dfhist->nlambda; i++)
949     {
950         snew((dfhist->accum_p)[i], dfhist->nlambda);
951         snew((dfhist->accum_m)[i], dfhist->nlambda);
952         snew((dfhist->accum_p2)[i], dfhist->nlambda);
953         snew((dfhist->accum_m2)[i], dfhist->nlambda);
954     }
955 }
956
957 extern void copy_df_history(df_history_t *df_dest, df_history_t *df_source)
958 {
959     int i, j;
960
961     init_df_history(df_dest, df_source->nlambda, df_source->wl_delta);
962     df_dest->nlambda = df_source->nlambda;
963     df_dest->bEquil  = df_source->bEquil;
964     for (i = 0; i < df_dest->nlambda; i++)
965     {
966         df_dest->sum_weights[i]  = df_source->sum_weights[i];
967         df_dest->sum_dg[i]       = df_source->sum_dg[i];
968         df_dest->sum_minvar[i]   = df_source->sum_minvar[i];
969         df_dest->sum_variance[i] = df_source->sum_variance[i];
970         df_dest->n_at_lam[i]     = df_source->n_at_lam[i];
971         df_dest->wl_histo[i]     = df_source->wl_histo[i];
972         df_dest->accum_p[i]      = df_source->accum_p[i];
973         df_dest->accum_m[i]      = df_source->accum_m[i];
974         df_dest->accum_p2[i]     = df_source->accum_p2[i];
975         df_dest->accum_m2[i]     = df_source->accum_m2[i];
976     }
977
978     for (i = 0; i < df_dest->nlambda; i++)
979     {
980         for (j = 0; j < df_dest->nlambda; j++)
981         {
982             df_dest->Tij[i][j]            = df_source->Tij[i][j];
983             df_dest->Tij_empirical[i][j]  = df_source->Tij_empirical[i][j];
984         }
985     }
986 }