GADGET-4
ewald Class Reference

Ewald correction functionality. More...

#include <ewald.h>

Inherits intposconvert, io_streamcount, and setcomm.

Classes

struct  ewald_header
 

Public Types

enum  interpolate_options { POINTMASS , MULTIPOLES , LONGRANGE_FORCETEST }
 

Public Member Functions

 ewald (void)
 
void ewald_init (void)
 This function initializes tables with the correction force and the correction potential due to the periodic images of a point mass located at the origin. More...
 
void ewald_corr (double dx, double dy, double dz, enum interpolate_options, ewald_data &fper)
 
void ewald_corr_exact (double dx, double dy, double dz, enum interpolate_options flag, ewald_data &fper)
 
void ewald_gridlookup (const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag, ewald_data &fper)
 
double ewald_D0 (double x, double y, double z)
 This function computes the potential correction term by means of Ewald summation. More...
 
vector< double > ewald_D1 (double x, double y, double z)
 This function computes the force correction term (difference between full force of infinite lattice and nearest image) by Ewald summation. More...
 
symtensor2< double > ewald_D2 (double x, double y, double z)
 
symtensor3< double > ewald_D3 (double x, double y, double z)
 
symtensor4< double > ewald_D4 (double x, double y, double z)
 
symtensor5< double > ewald_D5 (double x, double y, double z)
 
symtensor6< double > ewald_D6 (double x, double y, double z)
 
symtensor7< double > ewald_D7 (double x, double y, double z)
 
ewaldtensor6< double > ewald_P6 (void)
 
ewaldtensor8< double > ewald_P8 (void)
 
ewaldtensor10< double > ewald_P10 (void)
 
- Public Member Functions inherited from intposconvert
void constrain_intpos (MyIntPosType *pos)
 
template<typename T >
void diff_intpos_to_pos (MyIntPosType *a, MyIntPosType *b, T *posdiff, offset_tuple off=0)
 
MyIntPosType nearest_image_intpos_to_intpos_X (const MyIntPosType a, const MyIntPosType b)
 
MyIntPosType nearest_image_intpos_to_intpos_Y (const MyIntPosType a, const MyIntPosType b)
 
MyIntPosType nearest_image_intpos_to_intpos_Z (const MyIntPosType a, const MyIntPosType b)
 
template<typename T >
void nearest_image_intpos_to_pos (const MyIntPosType *const a, const MyIntPosType *const b, T *posdiff)
 
void nearest_image_intpos_to_absolute_intdist (const MyIntPosType *a, const MyIntPosType *b, MyIntPosType *delta)
 
template<typename T >
MySignedIntPosType pos_to_signedintpos (T posdiff)
 
template<typename T >
void pos_to_signedintpos (T *posdiff, MySignedIntPosType *intpos)
 
template<typename T >
void signedintpos_to_pos (MySignedIntPosType *intpos, T *pos)
 
double signedintpos_to_distanceorigin (MySignedIntPosType *intpos)
 
template<typename T >
void intpos_to_pos (MyIntPosType *intpos, T *posdiff)
 
void intpos_to_intpos (MyIntPosType *intpos, MyIntPosType *xyz)
 
template<typename T >
constrain_pos (T pos)
 
template<typename T >
void pos_to_intpos (T *posdiff, MyIntPosType *intpos)
 
- Public Member Functions inherited from io_streamcount
size_t my_fwrite (const void *ptr, size_t size, size_t nmemb, FILE *stream)
 A wrapper for the fwrite() function. More...
 
size_t my_fread (void *ptr, size_t size, size_t nmemb, FILE *stream)
 A wrapper for the fread() function. More...
 
void reset_io_byte_count (void)
 
long long get_io_byte_count (void)
 
- Public Member Functions inherited from setcomm
 setcomm (MPI_Comm Comm)
 
 setcomm (const char *str)
 
void initcomm (MPI_Comm Comm)
 
void mpi_printf (const char *fmt,...)
 
void determine_compute_nodes (void)
 

Additional Inherited Members

- Public Attributes inherited from intposconvert
MyReal FacIntToCoord
 
MyReal FacCoordToInt
 
MyReal RegionLen
 
MyReal RegionCorner [3]
 
MyReal RegionCenter [3]
 
- Public Attributes inherited from io_streamcount
long long byte_count = 0
 
- Public Attributes inherited from setcomm
MPI_Comm Communicator
 
int NTask
 
int ThisTask
 
int PTask
 
int ThisNode
 
int NumNodes = 0
 
int TasksInThisNode
 
int RankInThisNode
 
int MinTasksPerNode
 
int MaxTasksPerNode
 
long long MemoryOnNode
 
long long SharedMemoryOnNode
 

Detailed Description

Ewald correction functionality.

This class collects all the functionality provided for Ewald table lookups.

Definition at line 53 of file ewald.h.

Member Enumeration Documentation

◆ interpolate_options

Enumerator
POINTMASS 
MULTIPOLES 
LONGRANGE_FORCETEST 

Definition at line 65 of file ewald.h.

Constructor & Destructor Documentation

◆ ewald()

ewald ( void  )
inline

Definition at line 56 of file ewald.h.

Member Function Documentation

◆ ewald_corr()

void ewald_corr ( double  dx,
double  dy,
double  dz,
enum  interpolate_options,
ewald_data fper 
)

◆ ewald_corr_exact()

void ewald_corr_exact ( double  dx,
double  dy,
double  dz,
enum interpolate_options  flag,
ewald_data fper 
)

◆ ewald_D0()

double ewald_D0 ( double  x,
double  y,
double  z 
)

This function computes the potential correction term by means of Ewald summation.

Parameters
x,y,zcontains the distance vector for which the correction term should be computed
Returns
the correction term

Definition at line 705 of file ewald.cc.

◆ ewald_D1()

vector< double > ewald_D1 ( double  x,
double  y,
double  z 
)

This function computes the force correction term (difference between full force of infinite lattice and nearest image) by Ewald summation.

Parameters
x,y,zcontains the distance vector for which the correction force should be computed
forcearray will containing the correction force

Definition at line 898 of file ewald.cc.

◆ ewald_D2()

symtensor2< double > ewald_D2 ( double  x,
double  y,
double  z 
)

Definition at line 1090 of file ewald.cc.

◆ ewald_D3()

symtensor3< double > ewald_D3 ( double  x,
double  y,
double  z 
)

Definition at line 1312 of file ewald.cc.

◆ ewald_D4()

symtensor4< double > ewald_D4 ( double  x,
double  y,
double  z 
)

Definition at line 1550 of file ewald.cc.

◆ ewald_D5()

symtensor5< double > ewald_D5 ( double  x,
double  y,
double  z 
)

Definition at line 1694 of file ewald.cc.

◆ ewald_D6()

symtensor6< double > ewald_D6 ( double  x,
double  y,
double  z 
)

Definition at line 1845 of file ewald.cc.

◆ ewald_D7()

symtensor7< double > ewald_D7 ( double  x,
double  y,
double  z 
)

Definition at line 2011 of file ewald.cc.

◆ ewald_gridlookup()

void ewald_gridlookup ( const MyIntPosType p_intpos,
const MyIntPosType target_intpos,
enum interpolate_options  flag,
ewald_data fper 
)

Definition at line 355 of file ewald.cc.

◆ ewald_init()

void ewald_init ( void  )

This function initializes tables with the correction force and the correction potential due to the periodic images of a point mass located at the origin.

This file contains the computation of the Ewald correction table, and the corresponding lookup functions.

in D0phi we store the correction potential:

phi = 1/x + pi/alpha^2 - sum_q (erfc(alpha |x-q|)/|x-q|) - 4pi/V sum_k exp(-k^2/(4alpha^2))/k^2 cos(k*x)

in D1phi we store the first derivative of correction potential

dphi/dx_i

in D2phi we store the correction tensor (second derivatives of correction potential)

d2phi/(dx_i dx_j)

in D3phi we store the third order correction tensor (third derivatives of correction potential)

d3phi/(dx_i dx_j dx_k)

and so on also for D4phi and D5phi

These corrections are obtained by Ewald summation. (See for example Hernquist, Bouchet, Suto, ApJS, 1991, 75, 231) The correction fields are used to obtain the full periodic force if periodic boundaries combined with the pure tree algorithm are used. For the TreePM/FMM-PM algorithms, the Ewald correction is not used.

The correction fields are stored on disk once they are computed. If a corresponding file is found, they are loaded from disk to speed up the initialization. The Ewald summation issrc/gravtree_forcetest.c done in parallel, i.e. the processors share the work to compute the tables if needed.

Definition at line 69 of file ewald.cc.

◆ ewald_P10()

ewaldtensor10< double > ewald_P10 ( void  )

◆ ewald_P6()

ewaldtensor6< double > ewald_P6 ( void  )

◆ ewald_P8()

ewaldtensor8< double > ewald_P8 ( void  )

The documentation for this class was generated from the following files: