OpenAtom
Version1.5a
|
Add type declarations for simulationConstants class (readonly vars) and once class for each type of object we currently have (GStateXYslab, RealStateXZslab, RealRhoSlab, GRhoSlab, etc.) More...
Go to the source code of this file.
Classes | |
class | GStateSlab |
class | RealStateSlab |
class | RhoRealSlab |
class | RhoGSlab |
struct | essl_work |
= Holder classes for the plans : Allows many fft libaries to be used More... | |
struct | fftplanholder |
struct | rfftplanholder |
class | FFTcache |
class | GSlabInfo |
Typedefs | |
typedef struct essl_work | ESSL_WORK |
= Holder classes for the plans : Allows many fft libaries to be used | |
typedef struct fftplanholder | FFTplanHolder |
typedef struct rfftplanholder | RFFTplanHolder |
Functions | |
void | initGStateSlab (GStateSlab *gs, int sizeX, int sizeY, int sizeZ, int gSpaceUnits, int realSpaceUnits, int s_grain, int iplane_ind, int istate_ind, int len_nhc_cp, int num_nhc_cp, int nck_nhc_cp) |
= slab initialization helpers More... | |
void | initRealStateSlab (RealStateSlab *rs, int ngrid_a, int ngrid_b, int ngrid_c, int gSpaceUnits, int realSpaceUnits, int stateIndex, int thisPlane) |
void | initRhoRealSlab (RhoRealSlab *rho_rs, int xdim, int ydim, int zdim, int xdimA, int ydimA, int myIndexX, int myIndexY, int rhoRsubplanes) |
void | fft_split (FFTplanHolder *fftplanholder, int howmany, fftw_complex *in, int istride, int idist, fftw_complex *out, int ostride, int odist, int split, int index) |
= Eric's really cool BG/L progress callers More... | |
void | rfftwnd_complex_to_real_split (RFFTplanHolder *rfftplanholder, int howmany, fftw_complex *in, int istride, int idist, fftw_real *out, int ostride, int odist, int split, int index) |
split up an fft call into multiple invocations with cmiprogress calls between them. | |
void | rfftwnd_real_to_complex_split (RFFTplanHolder *rfftplanholder, int howmany, fftw_real *in, int istride, int idist, fftw_complex *out, int ostride, int odist, int split, int index) |
split up an fft call into multiple invocations with cmiprogress calls between them. | |
void | initFFTholder (FFTplanHolder *, int *, int *, int *, double *, int *, int *, int *, int *, int, int *, int *) |
void | initRCFFTholder (RFFTplanHolder *, int *, int *, int *, double *, int *, int *, int *, int *, int, int *, int *) |
void | initCRFFTholder (RFFTplanHolder *, int *, int *, int *, double *, int *, int *, int *, int *, int, int *, int *) |
void | make_essl_work_map (int, int *, int *, int *, int *, int *, int *, int) |
Add type declarations for simulationConstants class (readonly vars) and once class for each type of object we currently have (GStateXYslab, RealStateXZslab, RealRhoSlab, GRhoSlab, etc.)
Definition in file fftCacheSlab.h.
void fft_split | ( | FFTplanHolder * | fftplanholder, |
int | howmany, | ||
fftw_complex * | in, | ||
int | istride_in, | ||
int | idist_in, | ||
fftw_complex * | out_in, | ||
int | ostride_in, | ||
int | odist_in, | ||
int | split, | ||
int | index | ||
) |
= Eric's really cool BG/L progress callers
= Eric's really cool BG/L progress callers
Definition at line 1645 of file fftCache.C.
Referenced by FFTcache::doEextFFTGtoR_Gchare(), FFTcache::doEextFFTGtoR_Rchare(), FFTcache::doEextFFTGyToRy_Rchare(), FFTcache::doEextFFTRtoG_Gchare(), FFTcache::doEextFFTRtoG_Rchare(), FFTcache::doEextFFTRyToGy_Rchare(), FFTcache::doHartFFTGtoR_Gchare(), FFTcache::doNlFFTGtoR_Gchare(), FFTcache::doNlFFTGtoR_Rchare(), FFTcache::doNlFFTRtoG_Gchare(), FFTcache::doNlFFTRtoG_Rchare(), FFTcache::doRhoFFTGtoR_Gchare(), FFTcache::doRhoFFTGtoR_Rchare(), FFTcache::doRhoFFTGyToRy_Rchare(), FFTcache::doRhoFFTRtoG_Gchare(), FFTcache::doRhoFFTRtoG_Rchare(), FFTcache::doRhoFFTRyToGy_Rchare(), FFTcache::doStpFFTGtoR_Gchare(), FFTcache::doStpFFTGtoR_Rchare(), FFTcache::doStpFFTRtoG_Gchare(), and FFTcache::doStpFFTRtoG_Rchare().
void initGStateSlab | ( | GStateSlab * | gs, |
int | sizeX, | ||
int | sizeY, | ||
int | sizeZ, | ||
int | gSpaceUnits, | ||
int | realSpaceUnits, | ||
int | s_grain, | ||
int | iplane_ind, | ||
int | istate_ind, | ||
int | len_nhc_cp, | ||
int | num_nhc_cp, | ||
int | nck_nhc_cp | ||
) |
= slab initialization helpers
= Explanation of organization of data: A point is psi(kx,ky,kz)
Each gStateSlab is a collection of pts in lines of constant kx,ky. The maximum number of pts in any line is size[1] = nfftz Some lines are longer/shorter to spherical trunction (|k|<k_cut). Also, kx>=0 when at the Gamma point or doublePack=1, the only thing that has been recently tested. Collections of lines are parallelized. Each collection does not correspond to a unique plane of kx. The parameter gSpaceUnits is obsolete and must be unity. In order to create the state in real space, psi(x,y,z), an fft of all kz lines is performed followed by a transpose. The result of the transpose is parallelized by planes of z. This is different than it used to be. The number of points AFTER the fft is nlines*nfftz will be different for each chare array
Definition at line 35 of file stateSlab.C.
Referenced by CP_State_GSpacePlane::CP_State_GSpacePlane().
void initRealStateSlab | ( | RealStateSlab * | rs, |
int | ngrid_a, | ||
int | ngrid_b, | ||
int | ngrid_c, | ||
int | gSpaceUnits, | ||
int | realSpaceUnits, | ||
int | stateIndex, | ||
int | thisPlane | ||
) |
= Explanation of the organization: Each CP_State_RealSpacePlane instance is actually a planes of psi_I(x,y,z). We index into planeArr using the z-coordinate and each chare is an x-y plane. Data arrives from the gspaceplane of total size : nlines_tot*nfftz
Definition at line 392 of file stateSlab.C.
Referenced by CP_State_RealSpacePlane::CP_State_RealSpacePlane().