b2cc6609907d4f4b1671d4e9cec688a1a74ce3c6
[alexxy/gromacs.git] / src / contrib / timefft.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  * GROwing Monsters And Cloning Shrimps
34  */
35 #include <math.h>
36 #include <stdio.h>
37 #include <time.h>
38 #include "typedefs.h"
39 #include "macros.h"
40 #include "smalloc.h"
41 #include "xvgr.h"
42 #include "complex.h"
43 #include "copyrite.h"
44 #include "gmx_fft.h"
45 #include "mdrun.h"
46 #include "main.h"
47 #include "statutil.h"
48
49 #ifdef GMX_MPI
50 #include "gmx_parallel_3dfft.h"
51 #endif
52
53 #include "fftgrid.h"
54
55
56 int main(int argc,char *argv[])
57 {
58   int       mmm[] = { 8, 10, 12, 15, 16, 18, 20, 24, 25, 27, 30, 32, 36, 40,
59                       45, 48, 50, 54, 60, 64, 72, 75, 80, 81, 90, 100 };
60   int       nnn[] = { 24, 32, 48, 60, 72, 84, 96 };
61 #define NNN asize(nnn)
62   FILE      *fp,*fplog;
63   int       *niter;
64   int       i,j,n,nit,ntot,n3,rsize;
65   double    t,nflop,start;
66   double    *rt,*ct;
67   t_fftgrid *g;
68   t_commrec *cr;
69   static gmx_bool bReproducible = FALSE;
70   static int  nnode    = 1;
71   static int  nitfac  = 1;
72   t_pargs pa[] = {
73     { "-reproducible",   FALSE, etBOOL, {&bReproducible}, 
74       "Request binary reproducible results" },
75     { "-np",    FALSE, etINT, {&nnode},
76       "Number of NODEs" },
77     { "-itfac", FALSE, etINT, {&nitfac},
78       "Multiply number of iterations by this" }
79   };
80   static t_filenm fnm[] = {
81     { efLOG, "-g", "fft",      ffWRITE },
82     { efXVG, "-o", "fft",      ffWRITE }
83   };
84 #define NFILE asize(fnm)
85   
86   cr = init_par(&argc,&argv);
87   if (MASTER(cr))
88     CopyRight(stdout,argv[0]);
89   parse_common_args(&argc,argv,
90                     PCA_CAN_SET_DEFFNM | (MASTER(cr) ? 0 : PCA_QUIET),
91                     NFILE,fnm,asize(pa),pa,0,NULL,0,NULL);
92   gmx_log_open(ftp2fn(efLOG,NFILE,fnm),cr,1,0,&fplog);
93
94   snew(niter,NNN);
95   snew(ct,NNN);
96   snew(rt,NNN);
97   rsize = sizeof(real);
98   for(i=0; (i<NNN); i++) {
99     n  = nnn[i];
100     if (n < 16)
101       niter[i] = 50;
102     else if (n < 26)
103       niter[i] = 20;
104     else if (n < 51)
105       niter[i] = 10;
106     else
107       niter[i] = 5;
108     niter[i] *= nitfac;
109     nit = niter[i];
110     
111     if (MASTER(cr))
112       fprintf(stderr,"\r3D FFT (%s precision) %3d^3, niter %3d     ",
113               (rsize == 8) ? "Double" : "Single",n,nit);
114     
115     g  = mk_fftgrid(n,n,n,NULL,NULL,cr,bReproducible);
116
117     if (PAR(cr))
118       start = time(NULL);
119     else
120       start_time();
121     for(j=0; (j<nit); j++) {
122       gmxfft3D(g,GMX_FFT_REAL_TO_COMPLEX,cr);
123       gmxfft3D(g,GMX_FFT_COMPLEX_TO_REAL,cr);
124     }
125     if (PAR(cr)) 
126       rt[i] = time(NULL)-start;
127     else {
128       update_time();
129       rt[i] = node_time();
130     }
131     done_fftgrid(g);
132     sfree(g);
133   }
134   if (MASTER(cr)) {
135     fprintf(stderr,"\n");
136     fp=xvgropen(ftp2fn(efXVG,NFILE,fnm),
137                 "FFT timings","n^3","t (s)");
138     for(i=0; (i<NNN); i++) {
139       n3 = 2*niter[i]*nnn[i]*nnn[i]*nnn[i];
140       fprintf(fp,"%10d  %10g\n",nnn[i],rt[i]/(2*niter[i]));
141     }
142     gmx_fio_fclose(fp);
143   }
144   return 0;
145 }