Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / gmxlib / copyrite.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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #ifdef GMX_THREAD_MPI
43 #include <thread_mpi.h>
44 #endif
45
46 #ifdef HAVE_LIBMKL
47 #include <mkl.h>
48 #endif
49 #ifdef GMX_GPU
50 #include <cuda.h>
51 #include <cuda_runtime_api.h>
52 #endif
53 #ifdef GMX_FFT_FFTW3
54 #include <fftw3.h>
55 #endif
56
57 /* This file is completely threadsafe - keep it that way! */
58
59 #include <string.h>
60 #include <ctype.h>
61 #include "sysstuff.h"
62 #include "smalloc.h"
63 #include "string2.h"
64 #include "macros.h"
65 #include <time.h>
66 #include "random.h"
67 #include "statutil.h"
68 #include "copyrite.h"
69 #include "strdb.h"
70 #include "futil.h"
71 #include "vec.h"
72 #include "buildinfo.h"
73 #include "gmx_cpuid.h"
74
75 static void pr_two(FILE *out,int c,int i)
76 {
77   if (i < 10)
78     fprintf(out,"%c0%1d",c,i);
79   else
80     fprintf(out,"%c%2d",c,i);
81 }
82
83 void pr_difftime(FILE *out,double dt)
84 {
85   int    ndays,nhours,nmins,nsecs;
86   gmx_bool   bPrint,bPrinted;
87
88   ndays = dt/(24*3600);
89   dt    = dt-24*3600*ndays;
90   nhours= dt/3600;
91   dt    = dt-3600*nhours;
92   nmins = dt/60;
93   dt    = dt-nmins*60;
94   nsecs = dt;
95   bPrint= (ndays > 0);
96   bPrinted=bPrint;
97   if (bPrint) 
98     fprintf(out,"%d",ndays);
99   bPrint=bPrint || (nhours > 0);
100   if (bPrint) {
101     if (bPrinted)
102       pr_two(out,'d',nhours);
103     else 
104       fprintf(out,"%d",nhours);
105   }
106   bPrinted=bPrinted || bPrint;
107   bPrint=bPrint || (nmins > 0);
108   if (bPrint) {
109     if (bPrinted)
110       pr_two(out,'h',nmins);
111     else 
112       fprintf(out,"%d",nmins);
113   }
114   bPrinted=bPrinted || bPrint;
115   if (bPrinted)
116     pr_two(out,':',nsecs);
117   else
118     fprintf(out,"%ds",nsecs);
119   fprintf(out,"\n");
120 }
121
122
123 gmx_bool be_cool(void)
124 {
125   /* Yes, it is bad to check the environment variable every call,
126    * but we dont call this routine often, and it avoids using 
127    * a mutex for locking the variable...
128    */
129 #ifdef GMX_FAHCORE
130   /*be uncool*/
131   return FALSE;
132 #else
133   return (getenv("GMX_NO_QUOTES") == NULL);
134 #endif
135 }
136
137 void space(FILE *out, int n)
138 {
139   fprintf(out,"%*s",n,"");
140 }
141
142 void f(char *a)
143 {
144     int i;
145     int len=strlen(a);
146     
147     for(i=0;i<len;i++)
148         a[i]=~a[i]; 
149 }
150
151 static void sp_print(FILE *out,const char *s)
152 {
153   int slen;
154   
155   slen=strlen(s);
156   space(out,(80-slen)/2);
157   fprintf(out,"%s\n",s);
158 }
159
160 static void ster_print(FILE *out,const char *s)
161 {
162   int  slen;
163   char buf[128];
164   
165   snprintf(buf,128,":-)  %s  (-:",s);
166   slen=strlen(buf);
167   space(out,(80-slen)/2);
168   fprintf(out,"%s\n",buf);
169 }
170
171
172 static void pukeit(const char *db,const char *defstring, char *retstring, 
173                    int retsize, int *cqnum)
174 {
175   FILE *fp;
176   char **help;
177   int  i,nhlp;
178   int  seed;
179  
180   if (be_cool() && ((fp = low_libopen(db,FALSE)) != NULL)) {
181     nhlp=fget_lines(fp,&help);
182     /* for libraries we can use the low-level close routines */
183     ffclose(fp);
184     seed=time(NULL);
185     *cqnum=nhlp*rando(&seed);
186     if (strlen(help[*cqnum]) >= STRLEN)
187       help[*cqnum][STRLEN-1] = '\0';
188     strncpy(retstring,help[*cqnum],retsize);
189     f(retstring);
190     for(i=0; (i<nhlp); i++)
191       sfree(help[i]);
192     sfree(help);
193   }
194   else 
195     strncpy(retstring,defstring,retsize);
196 }
197
198 void bromacs(char *retstring, int retsize)
199 {
200   int dum;
201
202   pukeit("bromacs.dat",
203          "Groningen Machine for Chemical Simulation",
204          retstring,retsize,&dum);
205 }
206
207 void cool_quote(char *retstring, int retsize, int *cqnum)
208 {
209   char *tmpstr;
210   char *s,*ptr;
211   int tmpcq,*p;
212   
213   if (cqnum!=NULL)
214     p = cqnum;
215   else
216     p = &tmpcq;
217   
218   /* protect audience from explicit lyrics */
219   snew(tmpstr,retsize+1);
220   pukeit("gurgle.dat","Thanx for Using GROMACS - Have a Nice Day",
221          tmpstr,retsize-2,p);
222
223   if ((ptr = strchr(tmpstr,'_')) != NULL) {
224     *ptr='\0';
225     ptr++;
226     sprintf(retstring,"\"%s\" %s",tmpstr,ptr);
227   }
228   else {
229     strcpy(retstring,tmpstr);
230   }
231   sfree(tmpstr);
232 }
233
234 void CopyRight(FILE *out,const char *szProgram)
235 {
236   /* Dont change szProgram arbitrarily - it must be argv[0], i.e. the 
237    * name of a file. Otherwise, we won't be able to find the library dir.
238    */
239 #define NCR (int)asize(CopyrightText)
240 #ifdef GMX_FAHCORE
241 #define NGPL 0 /*FAH has an exception permission from GPL to allow digital signatures in Gromacs*/
242 #else
243 #define NGPL (int)asize(GPLText)
244 #endif
245
246   char buf[256],tmpstr[1024];
247   int i;
248
249 #ifdef GMX_FAHCORE
250   set_program_name("Gromacs");
251 #else
252   set_program_name(szProgram);
253 #endif
254
255   ster_print(out,"G  R  O  M  A  C  S");
256   fprintf(out,"\n");
257   
258   bromacs(tmpstr,1023);
259   sp_print(out,tmpstr); 
260   fprintf(out,"\n");
261
262   ster_print(out,GromacsVersion());
263   fprintf(out,"\n");
264
265   /* fprintf(out,"\n");*/
266
267   /* sp_print(out,"PLEASE NOTE: THIS IS A BETA VERSION\n");
268   
269   fprintf(out,"\n"); */
270
271   for(i=0; (i<NCR); i++) 
272     sp_print(out,CopyrightText[i]);
273   for(i=0; (i<NGPL); i++)
274     sp_print(out,GPLText[i]);
275
276   fprintf(out,"\n");
277
278   snprintf(buf,256,"%s",Program());
279 #ifdef GMX_DOUBLE
280   strcat(buf," (double precision)");
281 #endif
282   ster_print(out,buf);
283   fprintf(out,"\n");
284 }
285
286
287 void thanx(FILE *fp)
288 {
289   char cq[1024];
290   int  cqnum;
291
292   /* protect the audience from suggestive discussions */
293   cool_quote(cq,1023,&cqnum);
294   
295   if (be_cool()) 
296     fprintf(fp,"\ngcq#%d: %s\n\n",cqnum,cq);
297   else
298     fprintf(fp,"\n%s\n\n",cq);
299 }
300
301 typedef struct {
302   const char *key;
303   const char *author;
304   const char *title;
305   const char *journal;
306   int volume,year;
307   const char *pages;
308 } t_citerec;
309
310 void please_cite(FILE *fp,const char *key)
311 {
312   static const t_citerec citedb[] = {
313     { "Allen1987a",
314       "M. P. Allen and D. J. Tildesley",
315       "Computer simulation of liquids",
316       "Oxford Science Publications",
317       1, 1987, "1" },
318     { "Berendsen95a",
319       "H. J. C. Berendsen, D. van der Spoel and R. van Drunen",
320       "GROMACS: A message-passing parallel molecular dynamics implementation",
321       "Comp. Phys. Comm.",
322       91, 1995, "43-56" },
323     { "Berendsen84a",
324       "H. J. C. Berendsen, J. P. M. Postma, A. DiNola and J. R. Haak",
325       "Molecular dynamics with coupling to an external bath",
326       "J. Chem. Phys.",
327       81, 1984, "3684-3690" },
328     { "Ryckaert77a",
329       "J. P. Ryckaert and G. Ciccotti and H. J. C. Berendsen",
330       "Numerical Integration of the Cartesian Equations of Motion of a System with Constraints; Molecular Dynamics of n-Alkanes",
331       "J. Comp. Phys.",
332       23, 1977, "327-341" },
333     { "Miyamoto92a",
334       "S. Miyamoto and P. A. Kollman",
335       "SETTLE: An Analytical Version of the SHAKE and RATTLE Algorithms for Rigid Water Models",
336       "J. Comp. Chem.",
337       13, 1992, "952-962" },
338     { "Cromer1968a",
339       "D. T. Cromer & J. B. Mann",
340       "X-ray scattering factors computed from numerical Hartree-Fock wave functions",
341       "Acta Cryst. A",
342       24, 1968, "321" },
343     { "Barth95a",
344       "E. Barth and K. Kuczera and B. Leimkuhler and R. D. Skeel",
345       "Algorithms for Constrained Molecular Dynamics",
346       "J. Comp. Chem.",
347       16, 1995, "1192-1209" },
348     { "Essmann95a",
349       "U. Essmann, L. Perera, M. L. Berkowitz, T. Darden, H. Lee and L. G. Pedersen ",
350       "A smooth particle mesh Ewald method",
351       "J. Chem. Phys.",
352       103, 1995, "8577-8592" },
353     { "Torda89a",
354       "A. E. Torda and R. M. Scheek and W. F. van Gunsteren",
355       "Time-dependent distance restraints in molecular dynamics simulations",
356       "Chem. Phys. Lett.",
357       157, 1989, "289-294" },
358     { "Tironi95a",
359       "I. G. Tironi and R. Sperb and P. E. Smith and W. F. van Gunsteren",
360       "Generalized reaction field method for molecular dynamics simulations",
361       "J. Chem. Phys",
362       102, 1995, "5451-5459" },
363     { "Hess97a",
364       "B. Hess and H. Bekker and H. J. C. Berendsen and J. G. E. M. Fraaije",
365       "LINCS: A Linear Constraint Solver for molecular simulations",
366       "J. Comp. Chem.",
367       18, 1997, "1463-1472" },
368     { "Hess2008a",
369       "B. Hess",
370       "P-LINCS: A Parallel Linear Constraint Solver for molecular simulation",
371       "J. Chem. Theory Comput.",
372       4, 2008, "116-122" },
373     { "Hess2008b",
374       "B. Hess and C. Kutzner and D. van der Spoel and E. Lindahl",
375       "GROMACS 4: Algorithms for highly efficient, load-balanced, and scalable molecular simulation",
376       "J. Chem. Theory Comput.",
377       4, 2008, "435-447" },
378     { "Hub2010",
379       "J. S. Hub, B. L. de Groot and D. van der Spoel",
380       "g_wham - A free weighted histogram analysis implementation including robust error and autocorrelation estimates",
381       "J. Chem. Theory Comput.",
382       6, 2010, "3713-3720"}, 
383     { "In-Chul99a",
384       "Y. In-Chul and M. L. Berkowitz",
385       "Ewald summation for systems with slab geometry",
386       "J. Chem. Phys.",
387       111, 1999, "3155-3162" },
388     { "DeGroot97a",
389       "B. L. de Groot and D. M. F. van Aalten and R. M. Scheek and A. Amadei and G. Vriend and H. J. C. Berendsen",
390       "Prediction of Protein Conformational Freedom From Distance Constrains",
391       "Proteins",
392       29, 1997, "240-251" },
393     { "Spoel98a",
394       "D. van der Spoel and P. J. van Maaren and H. J. C. Berendsen",
395       "A systematic study of water models for molecular simulation. Derivation of models optimized for use with a reaction-field.",
396       "J. Chem. Phys.",
397       108, 1998, "10220-10230" },
398     { "Wishart98a",
399       "D. S. Wishart and A. M. Nip",
400       "Protein Chemical Shift Analysis: A Practical Guide",
401       "Biochem. Cell Biol.",
402       76, 1998, "153-163" },
403     { "Maiorov95",
404       "V. N. Maiorov and G. M. Crippen",
405       "Size-Independent Comparison of Protein Three-Dimensional Structures",
406       "PROTEINS: Struct. Funct. Gen.",
407       22, 1995, "273-283" },
408     { "Feenstra99",
409       "K. A. Feenstra and B. Hess and H. J. C. Berendsen",
410       "Improving Efficiency of Large Time-scale Molecular Dynamics Simulations of Hydrogen-rich Systems",
411       "J. Comput. Chem.",
412       20, 1999, "786-798" },
413     { "Timneanu2004a",
414       "N. Timneanu and C. Caleman and J. Hajdu and D. van der Spoel",
415       "Auger Electron Cascades in Water and Ice",
416       "Chem. Phys.",
417       299, 2004, "277-283" },
418     { "Pascal2011a",
419       "T. A. Pascal and S. T. Lin and W. A. Goddard III",
420       "Thermodynamics of liquids: standard molar entropies and heat capacities of common solvents from 2PT molecular dynamics",
421       "Phys. Chem. Chem. Phys.",
422       13, 2011, "169-181" },
423     { "Caleman2011b",
424       "C. Caleman and P. J. van Maaren and M. Hong and J. S. Hub and L. T. da Costa and D. van der Spoel",
425       "Force Field Benchmark of Organic Liquids: Density, Enthalpy of Vaporization, Heat Capacities, Surface Tension, Isothermal Compressibility, Volumetric Expansion Coefficient, and Dielectric Constant",
426       "J. Chem. Theo. Comp.",
427       8, 2012, "61" },
428     { "Lindahl2001a",
429       "E. Lindahl and B. Hess and D. van der Spoel",
430       "GROMACS 3.0: A package for molecular simulation and trajectory analysis",
431       "J. Mol. Mod.",
432       7, 2001, "306-317" },
433     { "Wang2001a",
434       "J. Wang and W. Wang and S. Huo and M. Lee and P. A. Kollman",
435       "Solvation model based on weighted solvent accessible surface area",
436       "J. Phys. Chem. B",
437       105, 2001, "5055-5067" },
438     { "Eisenberg86a",
439       "D. Eisenberg and A. D. McLachlan",
440       "Solvation energy in protein folding and binding",
441       "Nature",
442       319, 1986, "199-203" },
443     { "Eisenhaber95",
444       "Frank Eisenhaber and Philip Lijnzaad and Patrick Argos and Chris Sander and Michael Scharf",
445       "The Double Cube Lattice Method: Efficient Approaches to Numerical Integration of Surface Area and Volume and to Dot Surface Contouring of Molecular Assemblies",
446       "J. Comp. Chem.",
447       16, 1995, "273-284" },
448     { "Hess2002",
449       "B. Hess, H. Saint-Martin and H.J.C. Berendsen",
450       "Flexible constraints: an adiabatic treatment of quantum degrees of freedom, with application to the flexible and polarizable MCDHO model for water",
451       "J. Chem. Phys.",
452       116, 2002, "9602-9610" },
453     { "Hetenyi2002b",
454       "Csaba Hetenyi and David van der Spoel",
455       "Efficient docking of peptides to proteins without prior knowledge of the binding site.",
456       "Prot. Sci.",
457       11, 2002, "1729-1737" },
458     { "Hess2003",
459       "B. Hess and R.M. Scheek",
460       "Orientation restraints in molecular dynamics simulations using time and ensemble averaging",
461       "J. Magn. Res.",
462       164, 2003, "19-27" },
463     { "Rappe1991a",
464       "A. K. Rappe and W. A. Goddard III",
465       "Charge Equillibration for Molecular Dynamics Simulations",
466       "J. Phys. Chem.",
467       95, 1991, "3358-3363" },
468     { "Mu2005a",
469       "Y. Mu, P. H. Nguyen and G. Stock",
470       "Energy landscape of a small peptide revelaed by dihedral angle principal component analysis",
471       "Prot. Struct. Funct. Bioinf.",
472       58, 2005, "45-52" },
473     { "Okabe2001a",
474       "T. Okabe and M. Kawata and Y. Okamoto and M. Mikami",
475       "Replica-exchange {M}onte {C}arlo method for the isobaric-isothermal ensemble",
476       "Chem. Phys. Lett.",
477       335, 2001, "435-439" },
478     { "Hukushima96a",
479       "K. Hukushima and K. Nemoto",
480       "Exchange Monte Carlo Method and Application to Spin Glass Simulations",
481       "J. Phys. Soc. Jpn.",
482       65, 1996, "1604-1608" },
483     { "Tropp80a",
484       "J. Tropp",
485       "Dipolar Relaxation and Nuclear Overhauser effects in nonrigid molecules: The effect of fluctuating internuclear distances",
486       "J. Chem. Phys.",
487       72, 1980, "6035-6043" },
488     { "Bultinck2002a",
489        "P. Bultinck and W. Langenaeker and P. Lahorte and F. De Proft and P. Geerlings and M. Waroquier and J. P. Tollenaere",
490       "The electronegativity equalization method I: Parametrization and validation for atomic charge calculations",
491       "J. Phys. Chem. A",
492       106, 2002, "7887-7894" },
493     { "Yang2006b",
494       "Q. Y. Yang and K. A. Sharp",
495       "Atomic charge parameters for the finite difference Poisson-Boltzmann method using electronegativity neutralization",
496       "J. Chem. Theory Comput.",
497       2, 2006, "1152-1167" },
498     { "Spoel2005a",
499       "D. van der Spoel, E. Lindahl, B. Hess, G. Groenhof, A. E. Mark and H. J. C. Berendsen",
500       "GROMACS: Fast, Flexible and Free",
501       "J. Comp. Chem.",
502       26, 2005, "1701-1719" },
503     { "Spoel2006b",
504       "D. van der Spoel, P. J. van Maaren, P. Larsson and N. Timneanu",
505       "Thermodynamics of hydrogen bonding in hydrophilic and hydrophobic media",
506       "J. Phys. Chem. B",
507       110, 2006, "4393-4398" },
508     { "Spoel2006d",
509       "D. van der Spoel and M. M. Seibert",
510       "Protein folding kinetics and thermodynamics from atomistic simulations",
511       "Phys. Rev. Letters",
512       96, 2006, "238102" },
513     { "Palmer94a",
514       "B. J. Palmer",
515       "Transverse-current autocorrelation-function calculations of the shear viscosity for molecular liquids",
516       "Phys. Rev. E",
517       49, 1994, "359-366" },
518     { "Bussi2007a",
519       "G. Bussi, D. Donadio and M. Parrinello",
520       "Canonical sampling through velocity rescaling",
521       "J. Chem. Phys.",
522       126, 2007, "014101" },
523     { "Hub2006",
524       "J. S. Hub and B. L. de Groot",
525       "Does CO2 permeate through Aquaporin-1?",
526       "Biophys. J.",
527       91, 2006, "842-848" },
528     { "Hub2008",
529       "J. S. Hub and B. L. de Groot",
530       "Mechanism of selectivity in aquaporins and aquaglyceroporins",
531       "PNAS",
532       105, 2008, "1198-1203" },
533     { "Friedrich2009",
534       "M. S. Friedrichs, P. Eastman, V. Vaidyanathan, M. Houston, S. LeGrand, A. L. Beberg, D. L. Ensign, C. M. Bruns, and V. S. Pande",
535       "Accelerating Molecular Dynamic Simulation on Graphics Processing Units",
536       "J. Comp. Chem.",
537       30, 2009, "864-872" },
538     { "Engin2010",
539       "O. Engin, A. Villa, M. Sayar and B. Hess",
540       "Driving Forces for Adsorption of Amphiphilic Peptides to Air-Water Interface",
541       "J. Phys. Chem. B",
542       114, 2010, "11093" },
543     { "Fritsch12",
544       "S. Fritsch, C. Junghans and K. Kremer",
545       "Adaptive molecular simulation study on structure formation of toluene around C60 using Gromacs",
546       "J. Chem. Theo. Comp.",
547       8, 2012, "398" },
548     { "Junghans10",
549       "C. Junghans and S. Poblete",
550       "A reference implementation of the adaptive resolution scheme in ESPResSo",
551       "Comp. Phys. Comm.",
552       181, 2010, "1449" },
553     { "Wang2010",
554       "H. Wang, F. Dommert, C.Holm",
555       "Optimizing working parameters of the smooth particle mesh Ewald algorithm in terms of accuracy and efficiency",
556       "J. Chem. Phys. B",
557       133, 2010, "034117" },
558     { "Sugita1999a",
559       "Y. Sugita, Y. Okamoto",
560       "Replica-exchange molecular dynamics method for protein folding",
561       "Chem. Phys. Lett.",
562       314, 1999, "141-151" },
563     { "Kutzner2011",
564       "C. Kutzner and J. Czub and H. Grubmuller",
565       "Keep it Flexible: Driving Macromolecular Rotary Motions in Atomistic Simulations with GROMACS",
566       "J. Chem. Theory Comput.",
567       7, 2011, "1381-1393" },
568     { "Hoefling2011",
569       "M. Hoefling, N. Lima, D. Haenni, C.A.M. Seidel, B. Schuler, H. Grubmuller",
570       "Structural Heterogeneity and Quantitative FRET Efficiency Distributions of Polyprolines through a Hybrid Atomistic Simulation and Monte Carlo Approach",
571       "PLoS ONE",
572       6, 2011, "e19791" },
573     { "Hockney1988",
574       "R. W. Hockney and J. W. Eastwood",
575       "Computer simulation using particles",
576       "IOP, Bristol",
577       1, 1988, "1" },
578     { "Ballenegger2012",
579       "V. Ballenegger, J.J. Cerda, and C. Holm",
580       "How to Convert SPME to P3M: Influence Functions and Error Estimates",
581       "J. Chem. Theory Comput.",
582       8, 2012, "936-947" },
583     { "Garmay2012",
584       "Garmay Yu, Shvetsov A, Karelov D, Lebedev D, Radulescu A, Petukhov M, Isaev-Ivanov V",
585       "Correlated motion of protein subdomains and large-scale conformational flexibility of RecA protein filament",
586       "Journal of Physics: Conference Series",
587       340, 2012, "012094" }
588   };
589 #define NSTR (int)asize(citedb)
590   
591   int  j,index;
592   char *author;
593   char *title;
594 #define LINE_WIDTH 79
595   
596   if (fp == NULL)
597     return;
598
599   for(index=0; (index<NSTR) && (strcmp(citedb[index].key,key) != 0); index++)
600     ;
601   
602   fprintf(fp,"\n++++ PLEASE READ AND CITE THE FOLLOWING REFERENCE ++++\n");
603   if (index < NSTR) {
604     /* Insert newlines */
605     author = wrap_lines(citedb[index].author,LINE_WIDTH,0,FALSE);
606     title  = wrap_lines(citedb[index].title,LINE_WIDTH,0,FALSE);
607     fprintf(fp,"%s\n%s\n%s %d (%d) pp. %s\n",
608             author,title,citedb[index].journal,
609             citedb[index].volume,citedb[index].year,
610             citedb[index].pages);
611     sfree(author);
612     sfree(title);
613   }
614   else {
615     fprintf(fp,"Entry %s not found in citation database\n",key);
616   }
617   fprintf(fp,"-------- -------- --- Thank You --- -------- --------\n\n");
618   fflush(fp);
619 }
620
621 #ifdef USE_VERSION_H
622 /* Version information generated at compile time. */
623 #include "version.h"
624 #else
625 /* Fall back to statically defined version. */
626 static const char _gmx_ver_string[]="VERSION " VERSION;
627 #endif
628
629 const char *GromacsVersion()
630 {
631   return _gmx_ver_string;
632 }
633
634 void gmx_print_version_info(FILE *fp)
635 {
636 #ifdef GMX_GPU
637     int cuda_driver,cuda_runtime;
638 #endif
639
640     fprintf(fp, "Gromacs version:    %s\n", _gmx_ver_string);
641 #ifdef USE_VERSION_H
642     fprintf(fp, "GIT SHA1 hash:      %s\n", _gmx_full_git_hash);
643     /* Only print out the branch information if present.
644      * The generating script checks whether the branch point actually
645      * coincides with the hash reported above, and produces an empty string
646      * in such cases. */
647     if (_gmx_central_base_hash[0] != 0)
648     {
649         fprintf(fp, "Branched from:      %s\n", _gmx_central_base_hash);
650     }
651 #endif
652
653 #ifdef GMX_DOUBLE
654     fprintf(fp, "Precision:          double\n");
655 #else
656     fprintf(fp, "Precision:          single\n");
657 #endif
658
659 #ifdef GMX_THREAD_MPI
660     fprintf(fp, "MPI library:        thread_mpi\n");
661 #elif defined(GMX_MPI)
662     fprintf(fp, "MPI library:        MPI\n");
663 #else
664     fprintf(fp, "MPI library:        none\n");
665 #endif
666 #ifdef GMX_OPENMP
667     fprintf(fp, "OpenMP support:     enabled\n");
668 #else
669     fprintf(fp, "OpenMP support:     disabled\n");
670 #endif
671 #ifdef GMX_GPU
672     fprintf(fp, "GPU support:        enabled\n");
673 #else
674     fprintf(fp, "GPU support:        disabled\n");
675 #endif
676     /* A preprocessor trick to avoid duplicating logic from vec.h */
677 #define gmx_stringify2(x) #x
678 #define gmx_stringify(x) gmx_stringify2(x)
679     fprintf(fp, "invsqrt routine:    %s\n", gmx_stringify(gmx_invsqrt(x)));
680     fprintf(fp, "CPU acceleration:   %s\n", GMX_CPU_ACCELERATION_STRING);
681
682     /* TODO: Would be nicer to wrap this in a gmx_fft_version() call, but
683      * since that is currently in mdlib, can wait for master. */
684 #ifdef GMX_FFT_FFTPACK
685     fprintf(fp, "FFT library:        fftpack (built-in)\n");
686 #elif defined(GMX_FFT_FFTW3) && defined(GMX_NATIVE_WINDOWS)
687     fprintf(fp, "FFT library:        %s\n", "fftw3");
688 #elif defined(GMX_FFT_FFTW3) && defined(GMX_DOUBLE)
689     fprintf(fp, "FFT library:        %s\n", fftw_version);
690 #elif defined(GMX_FFT_FFTW3)
691     fprintf(fp, "FFT library:        %s\n", fftwf_version);
692 #elif defined(GMX_FFT_MKL)
693     fprintf(fp, "FFT library:        MKL\n");
694 #else
695     fprintf(fp, "FFT library:        unknown\n");
696 #endif
697 #ifdef GMX_LARGEFILES
698     fprintf(fp, "Large file support: enabled\n");
699 #else
700     fprintf(fp, "Large file support: disabled\n");
701 #endif
702 #ifdef HAVE_RDTSCP
703     fprintf(fp, "RDTSCP usage:       enabled\n");
704 #else
705     fprintf(fp, "RDTSCP usage:       disabled\n");
706 #endif
707
708     fprintf(fp, "Built on:           %s\n", BUILD_TIME);
709     fprintf(fp, "Built by:           %s\n", BUILD_USER);
710     fprintf(fp, "Build OS/arch:      %s\n", BUILD_HOST);
711     fprintf(fp, "Build CPU vendor:   %s\n", BUILD_CPU_VENDOR);
712     fprintf(fp, "Build CPU brand:    %s\n", BUILD_CPU_BRAND);
713     fprintf(fp, "Build CPU family:   %d   Model: %d   Stepping: %d\n",
714             BUILD_CPU_FAMILY, BUILD_CPU_MODEL, BUILD_CPU_STEPPING);
715     /* TODO: The below strings can be quite long, so it would be nice to wrap
716      * them. Can wait for later, as the master branch has ready code to do all
717      * that. */
718     fprintf(fp, "Build CPU features: %s\n", BUILD_CPU_FEATURES);
719     fprintf(fp, "C compiler:         %s\n", BUILD_C_COMPILER);
720     fprintf(fp, "C compiler flags:   %s\n", BUILD_CFLAGS);
721     if (BUILD_CXX_COMPILER[0] != '\0')
722     {
723         fprintf(fp, "C++ compiler:       %s\n", BUILD_CXX_COMPILER);
724         fprintf(fp, "C++ compiler flags: %s\n", BUILD_CXXFLAGS);
725     }
726 #ifdef HAVE_LIBMKL
727     /* MKL might be used for LAPACK/BLAS even if FFTs use FFTW, so keep it separate */
728     fprintf(fp, "Linked with Intel MKL version %s.%s.%s.\n",
729             __INTEL_MKL__,__INTEL_MKL_MINOR__,__INTEL_MKL_UPDATE__);
730 #endif
731 #ifdef GMX_GPU
732     fprintf(fp, "CUDA compiler:      %s\n",CUDA_NVCC_COMPILER_INFO);
733     cuda_driver = 0;
734     cudaDriverGetVersion(&cuda_driver);
735     cuda_runtime = 0;
736     cudaRuntimeGetVersion(&cuda_runtime);
737     fprintf(fp, "CUDA driver:        %d.%d\n",cuda_driver/1000, cuda_driver%100);
738     fprintf(fp, "CUDA runtime:       %d.%d\n",cuda_runtime/1000, cuda_runtime%100);
739 #endif
740 }