1 #ifndef PMESOLVERUTIL_H 2 #define PMESOLVERUTIL_H 33 const int jblock,
const int kblock,
34 int& i0,
int& i1,
int& j0,
int& j1,
int& k0,
int& k1) {
36 int isize, jsize, ksize;
48 isize = pmeGrid.
K1/2+1;
57 ksize = pmeGrid.
K1/2+1;
63 jsize = pmeGrid.
K1/2+1;
69 NAMD_bug(
"getPencilDim, invalid permutation");
73 if (jblock < 0 || jblock >= jblocks || kblock < 0 || kblock >= kblocks)
74 NAMD_bug(
"getPencilDim, invalid block indices");
79 j0 = jsize*jblock/jblocks;
80 j1 = jsize*(jblock+1)/jblocks - 1;
82 k0 = ksize*kblock/kblocks;
83 k1 = ksize*(kblock+1)/kblocks - 1;
90 const int iblock,
const int jblock,
const int kblock,
91 int& i0,
int& i1,
int& j0,
int& j1,
int& k0,
int& k1) {
93 getPencilDim(pmeGrid, permutation, jblock, kblock, i0, i1, j0, j1, k0, k1);
111 NAMD_bug(
"getBlockDim, invalid permutation");
115 if (iblock < 0 || iblock >= iblocks)
116 NAMD_bug(
"getBlockDim, invalid block index");
120 i0 = isize*iblock/iblocks;
121 i1 = isize*(iblock+1)/iblocks - 1;
139 float* dataSrc_in,
int dataSrcSize_in,
140 float* dataDst_in,
int dataDstSize_in,
141 int permutation,
PmeGrid pmeGrid,
144 if (dataSrc_in != NULL && dataSrc_in == dataDst_in)
145 NAMD_bug(
"FFTCompute::init, only out-of-place FFTs supported");
147 int permutationDst = permutation;
152 if (dataSrc_in == NULL) {
158 if (dataSrcSize_in < getDataSizeRequired(permutation, pmeGrid,
jblock,
kblock))
159 NAMD_bug(
"FFTCompute::init, invalid dataSrcSize_in");
165 if (dataDst_in == NULL) {
171 if (dataDstSize_in < getDataSizeRequired(permutationDst, pmeGrid,
jblock,
kblock))
172 NAMD_bug(
"FFTCompute::init, invalid dataDstSize_in");
182 NAMD_bug(
"FFTCompute::init, error setting up data buffers");
185 if (pmePencilType == 3) {
188 NAMD_bug(
"FFTCompute::init, 3D FFT requires a single pencil");
189 int n[3] = {pmeGrid.
K1, pmeGrid.
K2, pmeGrid.
K3};
191 }
else if (pmePencilType == 1) {
192 int i0, i1, j0, j1, k0, k1;
194 int n[1] = {i1-i0+1};
195 int howmany = (j1-j0+1)*(k1-k0+1);
197 plan1DX(n, howmany, flags);
199 plan1DY(n, howmany, flags);
201 plan1DZ(n, howmany, flags);
203 NAMD_bug(
"FFTCompute::init, invalid permutation");
205 }
else if (pmePencilType == 2) {
207 int i0, i1, j0, j1, k0, k1;
209 int n[2] = {pmeGrid.
K1, pmeGrid.
K2};
210 int howmany = k1-k0+1;
211 plan2D(n, howmany, flags);
213 NAMD_bug(
"FFTCompute::init, invalid pmePencilType");
235 int getDataSizeRequired(
int permutation,
PmeGrid pmeGrid,
241 int i0, i1, j0, j1, k0, k1;
243 i0, i1, j0, j1, k0, k1);
245 int size = (i1-i0+1)*(j1-j0+1)*(k1-k0+1);
257 virtual float* allocateData(
const int dataSizeRequired)=0;
258 virtual void plan3D(
int *n,
int flags)=0;
259 virtual void plan2D(
int *n,
int howmany,
int flags)=0;
260 virtual void plan1DX(
int *n,
int howmany,
int flags)=0;
261 virtual void plan1DY(
int *n,
int howmany,
int flags)=0;
262 virtual void plan1DZ(
int *n,
int howmany,
int flags)=0;
306 virtual void solve(
Lattice &lattice,
const bool doEnergy,
const bool doVirial,
float* data)=0;
308 virtual void getVirial(
double *virial)=0;
367 static inline double calcGridCoord(
const double x,
const double recip11,
const int nfftx) {
370 return (
double)((double)nfftx*(w - (floor(w + 0.5) - 0.5)));
373 static inline void calcGridCoord(
const double x,
const double y,
const double z,
374 const double recip11,
const double recip22,
const double recip33,
375 const int nfftx,
const int nffty,
const int nfftz,
376 double &frx,
double &fry,
double &frz) {
379 frx = (double)((
double)nfftx*(w - (floor(w + 0.5) - 0.5)));
381 fry = (double)((
double)nffty*(w - (floor(w + 0.5) - 0.5)));
383 frz = (double)((
double)nfftz*(w - (floor(w + 0.5) - 0.5)));
386 static inline void calcGridCoord(
const float x,
const float y,
const float z,
387 const float recip11,
const float recip22,
const float recip33,
388 const int nfftx,
const int nffty,
const int nfftz,
389 float &frx,
float &fry,
float &frz) {
391 w = x*recip11 + 2.0f;
392 frx = (float)(nfftx*(w - (floorf(w + 0.5f) - 0.5f)));
393 w = y*recip22 + 2.0f;
394 fry = (float)(nffty*(w - (floorf(w + 0.5f) - 0.5f)));
395 w = z*recip33 + 2.0f;
396 frz = (float)(nfftz*(w - (floorf(w + 0.5f) - 0.5f)));
399 static inline void calcGridCoord(
const float x,
const float y,
const float z,
400 const int nfftx,
const int nffty,
const int nfftz,
float &frx,
float &fry,
float &frz) {
401 frx = (float)(nfftx)*x;
402 fry = (float)(nffty)*y;
403 frz = (float)(nfftz)*z;
406 static inline void calcGridCoord(
const double x,
const double y,
const double z,
407 const int nfftx,
const int nffty,
const int nfftz,
double &frx,
double &fry,
double &frz) {
408 frx = (double)(nfftx)*x;
409 fry = (double)(nffty)*y;
410 frz = (double)(nfftz)*z;
417 #include <hip/hip_vector_types.h> 437 int i0, i1, j0, j1, k0, k1;
442 dataSize = (i1-i0+1)*(j1-j0+1)*(k1-k0+1);
455 NAMD_bug(
"PmeTranspose::PmeTranspose, invalid permutation");
463 for (
int iblock=0;iblock <
nblock;iblock++) {
467 int x0, y0dummy, y1dummy, z0dummy, z1dummy;
468 getBlockDim(
pmeGrid,
permutation, iblock,
jblock,
kblock, x0, x1, y0dummy, y1dummy, z0dummy, z1dummy);
480 #endif // PMESOLVERUTIL_H virtual void gatherForce(Lattice &lattice, CudaForce *force)=0
void compute_b_moduli(double *bm, int K, int order)
virtual void solve(Lattice &lattice, const bool doEnergy, const bool doVirial, float *data)=0
PmeRealSpaceCompute(PmeGrid pmeGrid, const int jblock, const int kblock, unsigned int grid=0)
virtual void backward()=0
static int getPencilIndexY(const PmeGrid &pmeGrid, const int y)
static void getPencilDim(const PmeGrid &pmeGrid, const int permutation, const int jblock, const int kblock, int &i0, int &i1, int &j0, int &j1, int &k0, int &k1)
static void calcGridCoord(const float x, const float y, const float z, const float recip11, const float recip22, const float recip33, const int nfftx, const int nffty, const int nfftz, float &frx, float &fry, float &frz)
const unsigned int NUM_GRID_MAX
virtual void transposeXYZtoYZX(const float2 *data)=0
static void calcGridCoord(const float x, const float y, const float z, const int nfftx, const int nffty, const int nfftz, float &frx, float &fry, float &frz)
virtual void transposeXYZtoZXY(const float2 *data)=0
PmeTranspose(PmeGrid pmeGrid, const int permutation, const int jblock, const int kblock)
static double calcGridCoord(const double x, const double recip11, const int nfftx)
void NAMD_bug(const char *err_msg)
virtual void getVirial(double *virial)=0
void setGrid(unsigned int i)
static void getBlockDim(const PmeGrid &pmeGrid, const int permutation, const int iblock, const int jblock, const int kblock, int &i0, int &i1, int &j0, int &j1, int &k0, int &k1)
static void calcGridCoord(const double x, const double y, const double z, const double recip11, const double recip22, const double recip33, const int nfftx, const int nffty, const int nfftz, double &frx, double &fry, double &frz)
unsigned int multipleGridIndex
void init(float *dataSrc_in, int dataSrcSize_in, float *dataDst_in, int dataDstSize_in, int permutation, PmeGrid pmeGrid, int pmePencilType, int jblock, int kblock, int flags)
static void calcGridCoord(const double x, const double y, const double z, const int nfftx, const int nffty, const int nfftz, double &frx, double &fry, double &frz)
virtual void spreadCharge(Lattice &lattice)=0
virtual void setGrid(unsigned int iGrid)
PmeKSpaceCompute(PmeGrid pmeGrid, const int permutation, const int jblock, const int kblock, double kappa, unsigned int multipleGridIndex=0)
virtual ~PmeKSpaceCompute()
PatchInfo(int x, int y, int z, int w)
virtual ~PmeRealSpaceCompute()
static int getPencilIndexZ(const PmeGrid &pmeGrid, const int z)
virtual double getEnergy()=0
virtual void copyAtoms(const int numAtoms, const CudaAtom *atoms)=0
unsigned int multipleGridIndex