Removing scripts and introducing wrapper binaries instead.
[gromacs.git] / include / fftgrid.h
blobf0892c71857b879a6569957b304cc7ba9bd854ba
1 /*
2 * $Id$
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gromacs Runs On Most of All Computer Systems
37 #ifndef _fftgrid_h
38 #define _fftgrid_h
40 #ifdef HAVE_CONFIG_H
41 #include <config.h>
42 #endif
44 #include <stdio.h>
45 #include "typedefs.h"
46 #include "fftw_wrapper.h"
47 #include "gmxcomplex.h"
48 #include "network.h"
50 /* Use FFTW */
52 #ifndef WITHOUT_FFTW
53 typedef t_complex t_fft_c;
54 typedef real t_fft_r;
55 #else /* NO FFTW PRESENT! */
56 typedef struct {
57 real re, im;
58 } t_fft_c;
59 typedef real t_fft_r;
60 typedef enum {
61 FFTW_FORWARD = -1, FFTW_BACKWARD = 1
62 } fftw_direction;
63 #endif
65 #define INDEX(i,j,k) ((i)*la12+(j)*la2+(k))
67 typedef struct {
68 int local_nx,local_x_start,local_ny_after_transpose;
69 int local_y_start_after_transpose,total_local_size;
70 } t_parfft;
72 typedef struct {
73 t_fft_r *ptr;
74 t_fft_r *localptr;
75 t_fft_r *workspace;
76 int nx,ny,nz,la2r,la2c,la12r,la12c;
77 int nptr,nxyz;
78 #ifndef WITHOUT_FFTW
79 rfftwnd_plan plan_fw;
80 rfftwnd_plan plan_bw;
81 #ifdef USE_MPI
82 rfftwnd_mpi_plan plan_mpi_fw;
83 rfftwnd_mpi_plan plan_mpi_bw;
84 t_parfft pfft;
85 #endif
86 #endif
87 } t_fftgrid;
89 extern t_fftgrid *mk_fftgrid(FILE *fp,bool bParallel,int nx,int ny,
90 int nz,bool bOptFFT);
91 /* Create an FFT grid (1 Dimensional), to be indexed by the INDEX macro
92 * Setup FFTW plans and extract local sizes for the grid.
93 * If the file pointer is given, information is printed to it.
96 extern void done_fftgrid(t_fftgrid *grid);
97 /* And throw it away again */
99 extern void gmxfft3D(t_fftgrid *grid,int dir,t_commrec *cr);
100 /* Do the FFT, direction may be either
101 * FFTW_FORWARD (sign -1) for real -> complex transform
102 * FFTW_BACKWARD (sign 1) for complex -> real transform
105 extern void clear_fftgrid(t_fftgrid *grid);
106 /* Set it to zero */
108 extern void unpack_fftgrid(t_fftgrid *grid,int *nx,int *ny,int *nz,
109 int *nx2,int *ny2,int *nz2,
110 int *la2, int *la12,bool bReal, t_fft_r **ptr);
112 /* Get the values for the constants into local copies */
117 /************************************************************************
119 * For backward compatibility (for testing the ewald code vs. PPPM etc)
120 * some old grid routines are retained here.
122 ************************************************************************/
124 extern real ***mk_rgrid(int nx,int ny,int nz);
126 extern void free_rgrid(real ***grid,int nx,int ny);
128 extern real print_rgrid(FILE *fp,char *title,int nx,int ny,int nz,
129 real ***grid);
131 extern void print_rgrid_pdb(char *fn,int nx,int ny,int nz,real ***grid);
133 extern t_complex ***mk_cgrid(int nx,int ny,int nz);
135 extern void free_cgrid(t_complex ***grid,int nx,int ny);
137 extern t_complex print_cgrid(FILE *fp,char *title,int nx,int ny,int nz,
138 t_complex ***grid);
140 extern void clear_cgrid(int nx,int ny,int nz,t_complex ***grid);
142 extern void clear_rgrid(int nx,int ny,int nz,real ***grid);
144 #endif