changed reading hint
[gromacs/adressmacs.git] / include / fftgrid.h
blobe26842d6ec4943972f466676d55a40c4288203bd
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 2.0
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
16 * Please refer to:
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
23 * or e-mail to:
24 * gromacs@chem.rug.nl
26 * And Hey:
27 * Good ROcking Metal Altar for Chronical Sinners
30 #ifndef _fftgrid_h
31 #define _fftgrid_h
33 static char *SRCID_fftgrid_h = "$Id$";
35 #include <stdio.h>
36 #include "typedefs.h"
37 #include "fftw.h"
38 #include "rfftw.h"
39 #ifdef USE_MPI
40 #include "mpi.h"
41 #include "rfftw_mpi.h"
42 #endif
43 #include "complex.h"
44 #include "network.h"
46 /* Use FFTW */
48 typedef t_complex t_fft_c;
49 typedef real t_fft_r;
51 #define INDEX(i,j,k) ((i)*la12+(j)*la2+(k))
54 typedef struct {
55 int local_nx,local_x_start,local_ny_after_transpose;
56 int local_y_start_after_transpose,total_local_size;
57 } t_parfft;
59 typedef struct {
60 t_fft_r *ptr;
61 t_fft_r *localptr;
62 t_fft_r *workspace;
63 int nx,ny,nz,la2r,la2c,la12r,la12c;
64 int nptr,nxyz;
65 rfftwnd_plan plan_fw;
66 rfftwnd_plan plan_bw;
67 #ifdef USE_MPI
68 rfftwnd_mpi_plan plan_mpi_fw;
69 rfftwnd_mpi_plan plan_mpi_bw;
70 t_parfft pfft;
71 #endif
72 } t_fftgrid;
74 extern t_fftgrid *mk_fftgrid(FILE *fp,bool bParallel,int nx,int ny,
75 int nz,bool bOptFFT);
76 /* Create an FFT grid (1 Dimensional), to be indexed by the INDEX macro
77 * Setup FFTW plans and extract local sizes for the grid.
78 * If the file pointer is given, information is printed to it.
81 extern void done_fftgrid(t_fftgrid *grid);
82 /* And throw it away again */
84 extern void gmxfft3D(t_fftgrid *grid,int dir,t_commrec *cr);
85 /* Do the FFT, direction may be either
86 * FFTW_FORWARD (sign -1) for real -> complex transform
87 * FFTW_BACKWARD (sign 1) for complex -> real transform
90 extern void clear_fftgrid(t_fftgrid *grid);
91 /* Set it to zero */
93 extern void unpack_fftgrid(t_fftgrid *grid,int *nx,int *ny,int *nz,
94 int *la2, int *la12,bool bReal, t_fft_r **ptr);
96 /* Get the values for the constants into local copies */
101 /************************************************************************
103 * For backward compatibility (for testing the ewald code vs. PPPM etc)
104 * some old grid routines are retained here.
106 ************************************************************************/
108 extern real ***mk_rgrid(int nx,int ny,int nz);
110 extern void free_rgrid(real ***grid,int nx,int ny);
112 extern real print_rgrid(FILE *fp,char *title,int nx,int ny,int nz,
113 real ***grid);
115 extern void print_rgrid_pdb(char *fn,int nx,int ny,int nz,real ***grid);
117 extern t_complex ***mk_cgrid(int nx,int ny,int nz);
119 extern void free_cgrid(t_complex ***grid,int nx,int ny);
121 extern t_complex print_cgrid(FILE *fp,char *title,int nx,int ny,int nz,
122 t_complex ***grid);
124 extern void clear_cgrid(int nx,int ny,int nz,t_complex ***grid);
126 extern void clear_rgrid(int nx,int ny,int nz,real ***grid);
128 #endif