Remove unnecessary config.h includes
[alexxy/gromacs.git] / src / gromacs / gmxlib / typedefs.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "gromacs/legacyheaders/typedefs.h"
41
42 #include <string.h>
43
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/pbcutil/pbc.h"
47 #include "gromacs/random/random.h"
48 #include "gromacs/utility/smalloc.h"
49
50 /* The source code in this file should be thread-safe.
51       Please keep it that way. */
52
53 int gmx_int64_to_int(gmx_int64_t step, const char *warn)
54 {
55     int i;
56
57     i = (int)step;
58
59     if (warn != NULL && (step < INT_MIN || step > INT_MAX))
60     {
61         fprintf(stderr, "\nWARNING during %s:\n", warn);
62         fprintf(stderr, "step value ");
63         fprintf(stderr, "%"GMX_PRId64, step);
64         fprintf(stderr, " does not fit in int, converted to %d\n\n", i);
65     }
66
67     return i;
68 }
69
70 void init_inputrec(t_inputrec *ir)
71 {
72     memset(ir, 0, (size_t)sizeof(*ir));
73     snew(ir->fepvals, 1);
74     snew(ir->expandedvals, 1);
75     snew(ir->simtempvals, 1);
76 }
77
78 static void done_pull_group(t_pull_group *pgrp)
79 {
80     if (pgrp->nat > 0)
81     {
82         sfree(pgrp->ind);
83         sfree(pgrp->ind_loc);
84         sfree(pgrp->weight);
85         sfree(pgrp->weight_loc);
86     }
87 }
88
89 static void done_pull(t_pull *pull)
90 {
91     int i;
92
93     for (i = 0; i < pull->ngroup+1; i++)
94     {
95         done_pull_group(pull->group);
96         done_pull_group(pull->dyna);
97     }
98 }
99
100 void done_inputrec(t_inputrec *ir)
101 {
102     int m;
103
104     for (m = 0; (m < DIM); m++)
105     {
106         if (ir->ex[m].a)
107         {
108             sfree(ir->ex[m].a);
109         }
110         if (ir->ex[m].phi)
111         {
112             sfree(ir->ex[m].phi);
113         }
114         if (ir->et[m].a)
115         {
116             sfree(ir->et[m].a);
117         }
118         if (ir->et[m].phi)
119         {
120             sfree(ir->et[m].phi);
121         }
122     }
123
124     sfree(ir->opts.nrdf);
125     sfree(ir->opts.ref_t);
126     sfree(ir->opts.annealing);
127     sfree(ir->opts.anneal_npoints);
128     sfree(ir->opts.anneal_time);
129     sfree(ir->opts.anneal_temp);
130     sfree(ir->opts.tau_t);
131     sfree(ir->opts.acc);
132     sfree(ir->opts.nFreeze);
133     sfree(ir->opts.QMmethod);
134     sfree(ir->opts.QMbasis);
135     sfree(ir->opts.QMcharge);
136     sfree(ir->opts.QMmult);
137     sfree(ir->opts.bSH);
138     sfree(ir->opts.CASorbitals);
139     sfree(ir->opts.CASelectrons);
140     sfree(ir->opts.SAon);
141     sfree(ir->opts.SAoff);
142     sfree(ir->opts.SAsteps);
143     sfree(ir->opts.bOPT);
144     sfree(ir->opts.bTS);
145
146     if (ir->pull)
147     {
148         done_pull(ir->pull);
149         sfree(ir->pull);
150     }
151 }
152
153 static void zero_history(history_t *hist)
154 {
155     hist->disre_initf  = 0;
156     hist->ndisrepairs  = 0;
157     hist->disre_rm3tav = NULL;
158     hist->orire_initf  = 0;
159     hist->norire_Dtav  = 0;
160     hist->orire_Dtav   = NULL;
161 }
162
163 static void zero_ekinstate(ekinstate_t *eks)
164 {
165     eks->ekin_n         = 0;
166     eks->ekinh          = NULL;
167     eks->ekinf          = NULL;
168     eks->ekinh_old      = NULL;
169     eks->ekinscalef_nhc = NULL;
170     eks->ekinscaleh_nhc = NULL;
171     eks->vscale_nhc     = NULL;
172     eks->dekindl        = 0;
173     eks->mvcos          = 0;
174 }
175
176 static void init_swapstate(swapstate_t *swapstate)
177 {
178     int ii, ic;
179
180     swapstate->eSwapCoords = 0;
181     swapstate->nAverage    = 0;
182
183     /* Ion/water position swapping */
184     for (ic = 0; ic < eCompNR; ic++)
185     {
186         for (ii = 0; ii < eIonNR; ii++)
187         {
188             swapstate->nat_req[ic][ii]        = 0;
189             swapstate->nat_req_p[ic][ii]      = NULL;
190             swapstate->inflow_netto[ic][ii]   = 0;
191             swapstate->inflow_netto_p[ic][ii] = NULL;
192             swapstate->nat_past[ic][ii]       = NULL;
193             swapstate->nat_past_p[ic][ii]     = NULL;
194             swapstate->fluxfromAtoB[ic][ii]   = 0;
195             swapstate->fluxfromAtoB_p[ic][ii] = NULL;
196         }
197     }
198     swapstate->fluxleak               = NULL;
199     swapstate->nions                  = 0;
200     swapstate->comp_from              = NULL;
201     swapstate->channel_label          = NULL;
202     swapstate->bFromCpt               = 0;
203     swapstate->nat[eChan0]            = 0;
204     swapstate->nat[eChan1]            = 0;
205     swapstate->xc_old_whole[eChan0]   = NULL;
206     swapstate->xc_old_whole[eChan1]   = NULL;
207     swapstate->xc_old_whole_p[eChan0] = NULL;
208     swapstate->xc_old_whole_p[eChan1] = NULL;
209 }
210
211 void init_energyhistory(energyhistory_t * enerhist)
212 {
213     enerhist->nener = 0;
214
215     enerhist->ener_ave     = NULL;
216     enerhist->ener_sum     = NULL;
217     enerhist->ener_sum_sim = NULL;
218     enerhist->dht          = NULL;
219
220     enerhist->nsteps     = 0;
221     enerhist->nsum       = 0;
222     enerhist->nsteps_sim = 0;
223     enerhist->nsum_sim   = 0;
224
225     enerhist->dht = NULL;
226 }
227
228 static void done_delta_h_history(delta_h_history_t *dht)
229 {
230     int i;
231
232     for (i = 0; i < dht->nndh; i++)
233     {
234         sfree(dht->dh[i]);
235     }
236     sfree(dht->dh);
237     sfree(dht->ndh);
238 }
239
240 void done_energyhistory(energyhistory_t * enerhist)
241 {
242     sfree(enerhist->ener_ave);
243     sfree(enerhist->ener_sum);
244     sfree(enerhist->ener_sum_sim);
245
246     if (enerhist->dht != NULL)
247     {
248         done_delta_h_history(enerhist->dht);
249         sfree(enerhist->dht);
250     }
251 }
252
253 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
254 {
255     int i, j;
256
257     state->ngtc          = ngtc;
258     state->nnhpres       = nnhpres;
259     state->nhchainlength = nhchainlength;
260     if (state->ngtc > 0)
261     {
262         snew(state->nosehoover_xi, state->nhchainlength*state->ngtc);
263         snew(state->nosehoover_vxi, state->nhchainlength*state->ngtc);
264         snew(state->therm_integral, state->ngtc);
265         for (i = 0; i < state->ngtc; i++)
266         {
267             for (j = 0; j < state->nhchainlength; j++)
268             {
269                 state->nosehoover_xi[i*state->nhchainlength + j]   = 0.0;
270                 state->nosehoover_vxi[i*state->nhchainlength + j]  = 0.0;
271             }
272         }
273         for (i = 0; i < state->ngtc; i++)
274         {
275             state->therm_integral[i]  = 0.0;
276         }
277     }
278     else
279     {
280         state->nosehoover_xi  = NULL;
281         state->nosehoover_vxi = NULL;
282         state->therm_integral = NULL;
283     }
284
285     if (state->nnhpres > 0)
286     {
287         snew(state->nhpres_xi, state->nhchainlength*nnhpres);
288         snew(state->nhpres_vxi, state->nhchainlength*nnhpres);
289         for (i = 0; i < nnhpres; i++)
290         {
291             for (j = 0; j < state->nhchainlength; j++)
292             {
293                 state->nhpres_xi[i*nhchainlength + j]   = 0.0;
294                 state->nhpres_vxi[i*nhchainlength + j]  = 0.0;
295             }
296         }
297     }
298     else
299     {
300         state->nhpres_xi  = NULL;
301         state->nhpres_vxi = NULL;
302     }
303 }
304
305
306 void init_state(t_state *state, int natoms, int ngtc, int nnhpres, int nhchainlength, int nlambda)
307 {
308     int i;
309
310     state->natoms = natoms;
311     state->flags  = 0;
312     state->lambda = 0;
313     snew(state->lambda, efptNR);
314     for (i = 0; i < efptNR; i++)
315     {
316         state->lambda[i] = 0;
317     }
318     state->veta   = 0;
319     clear_mat(state->box);
320     clear_mat(state->box_rel);
321     clear_mat(state->boxv);
322     clear_mat(state->pres_prev);
323     clear_mat(state->svir_prev);
324     clear_mat(state->fvir_prev);
325     init_gtc_state(state, ngtc, nnhpres, nhchainlength);
326     state->nalloc = state->natoms;
327     if (state->nalloc > 0)
328     {
329         snew(state->x, state->nalloc);
330         snew(state->v, state->nalloc);
331     }
332     else
333     {
334         state->x = NULL;
335         state->v = NULL;
336     }
337     state->sd_X = NULL;
338     state->cg_p = NULL;
339     zero_history(&state->hist);
340     zero_ekinstate(&state->ekinstate);
341     init_energyhistory(&state->enerhist);
342     init_df_history(&state->dfhist, nlambda);
343     init_swapstate(&state->swapstate);
344     state->ddp_count       = 0;
345     state->ddp_count_cg_gl = 0;
346     state->cg_gl           = NULL;
347     state->cg_gl_nalloc    = 0;
348 }
349
350 void done_state(t_state *state)
351 {
352     if (state->x)
353     {
354         sfree(state->x);
355     }
356     if (state->v)
357     {
358         sfree(state->v);
359     }
360     if (state->sd_X)
361     {
362         sfree(state->sd_X);
363     }
364     if (state->cg_p)
365     {
366         sfree(state->cg_p);
367     }
368     state->nalloc = 0;
369     if (state->cg_gl)
370     {
371         sfree(state->cg_gl);
372     }
373     state->cg_gl_nalloc = 0;
374     if (state->lambda)
375     {
376         sfree(state->lambda);
377     }
378     if (state->ngtc > 0)
379     {
380         sfree(state->nosehoover_xi);
381         sfree(state->nosehoover_vxi);
382         sfree(state->therm_integral);
383     }
384 }
385
386 t_state *serial_init_local_state(t_state *state_global)
387 {
388     int      i;
389     t_state *state_local;
390
391     snew(state_local, 1);
392
393     /* Copy all the contents */
394     *state_local = *state_global;
395     snew(state_local->lambda, efptNR);
396     /* local storage for lambda */
397     for (i = 0; i < efptNR; i++)
398     {
399         state_local->lambda[i] = state_global->lambda[i];
400     }
401
402     return state_local;
403 }
404
405 static void do_box_rel(t_inputrec *ir, matrix box_rel, matrix b, gmx_bool bInit)
406 {
407     int d, d2;
408
409     for (d = YY; d <= ZZ; d++)
410     {
411         for (d2 = XX; d2 <= (ir->epct == epctSEMIISOTROPIC ? YY : ZZ); d2++)
412         {
413             /* We need to check if this box component is deformed
414              * or if deformation of another component might cause
415              * changes in this component due to box corrections.
416              */
417             if (ir->deform[d][d2] == 0 &&
418                 !(d == ZZ && d2 == XX && ir->deform[d][YY] != 0 &&
419                   (b[YY][d2] != 0 || ir->deform[YY][d2] != 0)))
420             {
421                 if (bInit)
422                 {
423                     box_rel[d][d2] = b[d][d2]/b[XX][XX];
424                 }
425                 else
426                 {
427                     b[d][d2] = b[XX][XX]*box_rel[d][d2];
428                 }
429             }
430         }
431     }
432 }
433
434 void set_box_rel(t_inputrec *ir, t_state *state)
435 {
436     /* Make sure the box obeys the restrictions before we fix the ratios */
437     correct_box(NULL, 0, state->box, NULL);
438
439     clear_mat(state->box_rel);
440
441     if (PRESERVE_SHAPE(*ir))
442     {
443         do_box_rel(ir, state->box_rel, state->box, TRUE);
444     }
445 }
446
447 void preserve_box_shape(t_inputrec *ir, matrix box_rel, matrix b)
448 {
449     if (PRESERVE_SHAPE(*ir))
450     {
451         do_box_rel(ir, box_rel, b, FALSE);
452     }
453 }
454
455 real max_cutoff(real cutoff1, real cutoff2)
456 {
457     if (cutoff1 == 0 || cutoff2 == 0)
458     {
459         return 0;
460     }
461     else
462     {
463         return max(cutoff1, cutoff2);
464     }
465 }
466
467 void init_df_history(df_history_t *dfhist, int nlambda)
468 {
469     int i;
470
471     dfhist->nlambda  = nlambda;
472     dfhist->bEquil   = 0;
473     dfhist->wl_delta = 0;
474
475     if (nlambda > 0)
476     {
477         snew(dfhist->sum_weights, dfhist->nlambda);
478         snew(dfhist->sum_dg, dfhist->nlambda);
479         snew(dfhist->sum_minvar, dfhist->nlambda);
480         snew(dfhist->sum_variance, dfhist->nlambda);
481         snew(dfhist->n_at_lam, dfhist->nlambda);
482         snew(dfhist->wl_histo, dfhist->nlambda);
483
484         /* allocate transition matrices here */
485         snew(dfhist->Tij, dfhist->nlambda);
486         snew(dfhist->Tij_empirical, dfhist->nlambda);
487
488         /* allocate accumulators for various transition matrix
489            free energy methods here */
490         snew(dfhist->accum_p, dfhist->nlambda);
491         snew(dfhist->accum_m, dfhist->nlambda);
492         snew(dfhist->accum_p2, dfhist->nlambda);
493         snew(dfhist->accum_m2, dfhist->nlambda);
494
495         for (i = 0; i < dfhist->nlambda; i++)
496         {
497             snew(dfhist->Tij[i], dfhist->nlambda);
498             snew(dfhist->Tij_empirical[i], dfhist->nlambda);
499             snew((dfhist->accum_p)[i], dfhist->nlambda);
500             snew((dfhist->accum_m)[i], dfhist->nlambda);
501             snew((dfhist->accum_p2)[i], dfhist->nlambda);
502             snew((dfhist->accum_m2)[i], dfhist->nlambda);
503         }
504     }
505 }
506
507 extern void copy_df_history(df_history_t *df_dest, df_history_t *df_source)
508 {
509     int i, j;
510
511     /* Currently, there should not be any difference in nlambda between the two,
512        but this is included for completeness for potential later functionality */
513     df_dest->nlambda  = df_source->nlambda;
514     df_dest->bEquil   = df_source->bEquil;
515     df_dest->wl_delta = df_source->wl_delta;
516
517     for (i = 0; i < df_dest->nlambda; i++)
518     {
519         df_dest->sum_weights[i]  = df_source->sum_weights[i];
520         df_dest->sum_dg[i]       = df_source->sum_dg[i];
521         df_dest->sum_minvar[i]   = df_source->sum_minvar[i];
522         df_dest->sum_variance[i] = df_source->sum_variance[i];
523         df_dest->n_at_lam[i]     = df_source->n_at_lam[i];
524         df_dest->wl_histo[i]     = df_source->wl_histo[i];
525     }
526
527     for (i = 0; i < df_dest->nlambda; i++)
528     {
529         for (j = 0; j < df_dest->nlambda; j++)
530         {
531             df_dest->accum_p[i][j]        = df_source->accum_p[i][j];
532             df_dest->accum_m[i][j]        = df_source->accum_m[i][j];
533             df_dest->accum_p2[i][j]       = df_source->accum_p2[i][j];
534             df_dest->accum_m2[i][j]       = df_source->accum_m2[i][j];
535             df_dest->Tij[i][j]            = df_source->Tij[i][j];
536             df_dest->Tij_empirical[i][j]  = df_source->Tij_empirical[i][j];
537         }
538     }
539 }
540
541 void done_df_history(df_history_t *dfhist)
542 {
543     int i;
544
545     if (dfhist->nlambda > 0)
546     {
547         sfree(dfhist->n_at_lam);
548         sfree(dfhist->wl_histo);
549         sfree(dfhist->sum_weights);
550         sfree(dfhist->sum_dg);
551         sfree(dfhist->sum_minvar);
552         sfree(dfhist->sum_variance);
553
554         for (i = 0; i < dfhist->nlambda; i++)
555         {
556             sfree(dfhist->Tij[i]);
557             sfree(dfhist->Tij_empirical[i]);
558             sfree(dfhist->accum_p[i]);
559             sfree(dfhist->accum_m[i]);
560             sfree(dfhist->accum_p2[i]);
561             sfree(dfhist->accum_m2[i]);
562         }
563     }
564     dfhist->bEquil   = 0;
565     dfhist->nlambda  = 0;
566     dfhist->wl_delta = 0;
567 }