Implementation of the Δ_n factor as a series; error estimates.
[qpms.git] / qpms / ewald.h
blob9485a2fe58e3071db4dcdff3c91cb3addc1aa6a7
1 /*! \file ewald.h
2 * \brief Lattice sums of spherical waves.
4 * Implementation of two-dimensional lattice sum in three dimensions
5 * according to:
6 * - [1] C.M. Linton, I. Thompson
7 * Journal of Computational Physics 228 (2009) 1815–1829
8 * - [2] C.M.Linton
9 * SIAM Review Vol 52, No. 4, pp. 630–674
11 * N.B.!!! currently, the long-range parts are calculated
12 * not according to [1,(4.5)], but rather
13 * according to the spherical-harmonic-normalisation-independent
14 * formulation in my notes notes/ewald.lyx.
15 * Both parts of lattice sums are then calculated with
16 * the \f$ P_n^{|m|} e^{im\phi} \f$
17 * (N.B. or \f$ P_n^{|m|} e^{imf} (-1)^m \f$ for negative m)
18 * substituted in place of \f$ Y_n^m \f$
19 * (this is quite a weird normalisation especially
20 * for negative \f$ |m| \f$, but it is consistent
21 * with the current implementation of the translation coefficients in
22 * @ref translations.c;
23 * in the long run, it might make more sense to replace it everywhere with normalised
24 * Legendre polynomials).
27 #ifndef EWALD_H
28 #define EWALD_H
29 #include <gsl/gsl_sf_result.h>
30 #include <stdlib.h>
31 #include <gsl/gsl_sf_legendre.h>
32 #include <gsl/gsl_errno.h>
33 #include <math.h> // for inlined lilgamma
34 #include <complex.h>
35 #include "qpms_types.h"
36 #include "lattices.h"
39 typedef enum {
40 QPMS_EWALD_LONG_RANGE = 1,
41 QPMS_EWALD_SHORT_RANGE = 2,
42 QPMS_EWALD_0TERM = 4,
43 QPMS_EWALD_FULL = QPMS_EWALD_LONG_RANGE | QPMS_EWALD_SHORT_RANGE | QPMS_EWALD_0TERM,
44 } qpms_ewald_part;
47 /// Use this handler to ignore underflows of incomplete gamma.
48 gsl_error_handler_t IgnoreUnderflowsGSLErrorHandler;
51 /// Object holding the Ewald sum constant factors.
52 /**
53 * Used internally by qpms_translation_calculator_t.
54 * Initialised by qpms_ewald3_constants_init() and freed by qpms_ewald3_constants_free().
56 typedef struct qpms_ewald3_constants_t {
57 qpms_l_t lMax;
58 qpms_y_t nelem_sc;
59 /// The values of maximum \a j's in the long-range part summation, `[(l-|m|/2)]`.
60 qpms_l_t *s1_jMaxes;
61 /// The constant factors for the long range part of a 2D Ewald sum.
62 complex double **s1_constfacs; // indices [y][j] where j is same as in [1, (4.5)]
63 /* These are the actual numbers now: (in the EWALD32_CONSTANTS_AGNOSTIC version)
64 * for m + n EVEN:
66 * s1_constfacs[y(m,n)][j] =
68 * -2 * I**(n+1) * sqrt(π) * ((n-m)/2)! * ((n+m)/2)! * (-1)**j
69 * -----------------------------------------------------------
70 * j! * ((n-m)/2 - j)! * ((n+m)/2 + j)!
72 * for m + n ODD:
74 * s1_constfacs[y(m,n)][j] = 0
76 complex double *s1_constfacs_base; ///< Internal pointer holding memory for the 2D Ewald sum constant factors.
77 // similarly for the 1D z-axis aligned case; now the indices are [n][j] (as m == 0)
78 /// The constant factors for the long range part of a 1D Ewald sum along the \a z axis.
79 /** If the summation points lie along a different direction, use the formula for
80 * 2D sum with additional factor of
81 * \f$ \sqrt{pi} \kappa \gamma(\abs{\vect{k}+\vect{K}}/\kappa) \f$.
83 complex double **s1_constfacs_1Dz;
84 /* These are the actual numbers now:
85 * s1_constfacs_1Dz[n][j] =
87 * -I**(n+1) (-1)**j * n!
88 * --------------------------
89 * j! * 2**(2*j) * (n - 2*j)!
91 complex double *s1_constfacs_1Dz_base; ///<Internal pointer holding memory for the 1D Ewald sum constant factors.
93 double *legendre0; /* now with GSL_SF_LEGENDRE_NONE normalisation, because this is what is
94 * what the multipliers from translations.c count with.
96 double *legendre_plus1; // needed? TODO; in any case, nonzero only for m=0
97 double *legendre_minus1; // needed? TODO; in any case, nonzero only for m=0
98 gsl_sf_legendre_t legendre_normconv;
99 int legendre_csphase; /* 1 or -1; csphase of the Legendre polynomials saved in legendre0 etc.
100 This is because I dont't actually consider this fixed in
101 translations.c */
103 } qpms_ewald3_constants_t;
105 /// Constructor for qpms_ewald3_constants_t.
106 qpms_ewald3_constants_t *qpms_ewald3_constants_init(qpms_l_t lMax, int csphase);
107 /// Destructor for qpms_ewald3_constants_t.
108 void qpms_ewald3_constants_free(qpms_ewald3_constants_t *);
111 /// Structure for holding complex-valued result of computation and an error estimate.
112 /** Similar to gsl_sf_result, but with complex val. */
113 typedef struct qpms_csf_result {
114 complex double val; ///< Calculation result.
115 double err; ///< Error estimate.
116 } qpms_csf_result;
119 // [1, (A.9)]
120 static inline complex double lilgamma(double t) {
121 t = fabs(t);
122 if (t >= 1)
123 return sqrt(t*t - 1);
124 else
125 return -I * sqrt(1 - t*t);
128 // [1, (A.8)], complex version of lilgamma()
129 static inline complex double clilgamma(complex double z) {
130 complex double a1 = z - 1, a2 = z + 1;
131 // ensure -pi/2 < arg(z + 1) < 3*pi/2
132 if (creal(a2) < 0 && cimag(a2) <= 0)
133 a2 = -csqrt(a2);
134 else
135 a2 = csqrt(a2);
136 // ensure -3*pi/2 < arg(z - 1) < pi/2
137 if (creal(a1) < 0 && cimag(a1) >= 0)
138 a1 = -csqrt(a1);
139 else
140 a1 = csqrt(a1);
141 return a1 * a2;
144 /// Incomplete Gamma function as a series.
145 /** DLMF 8.7.3 (latter expression) for complex second argument.
147 * The principal value is calculated. On the negative real axis
148 * (where the function has branch cut), the sign of the imaginary
149 * part is what matters (even if it is zero). Therefore one
150 * can have
151 * `cx_gamma_inc_series_e(a, z1) != cx_gamma_inc_series_e(a, z2)`
152 * even if `z1 == z2`, because `-0 == 0` according to IEEE 754.
153 * The side of the branch cut can be determined using `signbit(creal(z))`.
155 int cx_gamma_inc_series_e(double a, complex double z, qpms_csf_result * result);
157 /// Incomplete Gamma function as continued fractions.
158 /**
159 * The principal value is calculated. On the negative real axis
160 * (where the function has branch cut), the sign of the imaginary
161 * part is what matters (even if it is zero). Therefore one
162 * can have
163 * `cx_gamma_inc_CF_e(a, z1) != cx_gamma_inc_CF_e(a, z2)`
164 * even if `z1 == z2`, because `-0 == 0` according to IEEE 754.
165 * The side of the branch cut can be determined using `signbit(creal(z))`.
167 int cx_gamma_inc_CF_e(double a, complex double z, qpms_csf_result * result);
169 /// Incomplete gamma for complex second argument.
170 /**
171 * If x is (almost) real, it just uses gsl_sf_gamma_inc_e().
173 * On the negative real axis
174 * (where the function has branch cut), the sign of the imaginary
175 * part is what matters (even if it is zero). Therefore one
176 * can have
177 * `complex_gamma_inc_e(a, z1, m) != complex_gamma_inc_e(a, z2, m)`
178 * even if `z1 == z2`, because `-0 == 0` according to IEEE 754.
179 * The side of the branch cut can be determined using `signbit(creal(z))`.
181 * Another than principal branch can be selected using non-zero \a m
182 * argument.
184 int complex_gamma_inc_e(double a, complex double x,
185 /// Branch index.
186 /** If zero, the principal value is calculated.
187 * Other branches might be chosen using non-zero \a m.
188 * In such case, the returned value corresponds to \f[
189 * \Gamma(a,ze^{2\pi mi})=e^{2\pi mia} \Gamma(a,z)
190 * + (1-e^{2\pi mia}) \Gamma(a).
191 * \f]
193 * If \a a is non-positive integer, the limiting value should
194 * be used, but this is not yet implemented!
196 int m,
197 qpms_csf_result *result);
199 /// Exponential integral for complex second argument.
200 /** If x is (almost) positive real, it just uses gsl_sf_expint_En_e(). */
201 int complex_expint_n_e(int n, complex double x, qpms_csf_result *result);
204 /// Hypergeometric 2F2, used to calculate some errors.
205 int hyperg_2F2_series(double a, double b, double c, double d, double x,
206 gsl_sf_result *result);
208 #if 0
209 // The integral from (4.6); maybe should be static and not here.
210 int ewald32_sr_integral(double r, double k, double n, double eta, double *result, double *err, gsl_integration_workspace *workspace);
211 #endif
213 /// The Delta_n factor from [Kambe II], Appendix 3, used in 2D-in-3D long range sum.
214 /** \f[ \Delta_n = \int_n^\infty t^{-1/2 - n} \exp(-t + z^2/(4t))\ud t \f]
216 * \bug The current choice of method, based purely on the value of \a z, might be
217 * unsuitable especially for big values of \a maxn.
220 void ewald3_2_sigma_long_Delta(complex double *target, double *target_err, int maxn, complex double x, complex double z);
222 /// The Delta_n factor from [Kambe II], Appendix 3, used in 2D-in-3D long range sum.
223 /** This function always uses Kambe's (corrected) recurrent formula.
224 * For production, use ewald3_2_sigma_long_Delta() instead.
226 void ewald3_2_sigma_long_Delta_recurrent(complex double *target, double *target_err, int maxn, complex double x, complex double z);
228 /// The Delta_n factor from [Kambe II], Appendix 3, used in 2D-in-3D long range sum.
229 /** This function always uses Taylor expansion in \a z.
230 * For production, use ewald3_2_sigma_long_Delta() instead.
232 * \bug The error estimate seems to be wrong (too small) at least in some cases: try
233 * parameters maxn = 40, z = 0.5, x = -3. This might be related to the exponential growth
234 * of the error.
236 void ewald3_2_sigma_long_Delta_series(complex double *target, double *target_err, int maxn, complex double x, complex double z);
238 // General functions acc. to [2], sec. 4.6 – currently valid for 2D and 1D lattices in 3D space
240 /// The Ewald sum "self-interaction" term that appears in the lattice sums with zero (direct-space) Bravais lattice displacement.
241 int ewald3_sigma0(complex double *result, ///< Pointer to save the result (single complex double).
242 double *err, ///< Pointer to save the result error estimate (single double).
243 const qpms_ewald3_constants_t *c, ///< Constant factors structure initialised by qpms_ewald3_constants_init().
244 double eta, ///< Ewald parameter.
245 complex double wavenumber ///< Wavenumber of the background medium.
248 /// Short-range part of outgoing scalar spherical wavefunctions' lattice sum \f$ \sigma_{l,m}^\mathrm{S}(\vect k,\vect s)\f$.
249 int ewald3_sigma_short(
250 complex double *target_sigmasr_y, ///< Target array for \f$ \sigma_{l,m}^\mathrm{S} \f$, must be `c->nelem_sc` long.
251 double *target_sigmasr_y_err, ///< Target array for error estimates, must be `c->nelem_sc` long or `NULL`.
252 const qpms_ewald3_constants_t *c, ///< Constant factors structure initialised by qpms_ewald3_constants_init().
253 double eta, ///< Ewald parameter.
254 complex double wavenumber, ///< Wavenumber of the background medium.
255 /// Lattice dimensionality.
256 /** Ignored apart from asserts and possible optimisations, as the SR formula stays the same. */
257 LatticeDimensionality latdim,
258 /// Lattice point generator for the direct Bravais lattice.
259 /** There is a possibility that the whole PGen is not consumed
260 * (this might happen if the summand start to be consistently smaller
261 * than the (partial) sums * DBL_EPSILON.
262 * In such case, it is the responsibility of the caller to deallocate
263 * the generator.
265 PGen *pgen_R,
266 /// Indicates whether pgen_R already generates shifted points.
267 /** If false, the behaviour corresponds to the old ewald32_sigma_short_points_and_shift(),
268 * so the function assumes that the generated points correspond to the unshifted Bravais lattice,
269 * and adds particle_shift to the generated points before calculations.
270 * If true, it assumes that they are already shifted (if calculating interaction between
271 * different particles in the unit cell).
273 bool pgen_generates_shifted_points,
274 /// Wave vector \f$\vect k\f$.
275 cart3_t k,
276 /// Lattice offset \f$\vect s\f$ wrt. the Bravais lattice.
277 cart3_t particle_shift
280 /// Long-range part of outgoing scalar spherical wavefunctions' lattice sum \f$ \sigma_{l,m}^\mathrm{L}(\vect k,\vect s)\f$.
281 int ewald3_sigma_long( // calls ewald3_21_sigma_long or ewald3_3_sigma_long, depending on latdim
282 complex double *target_sigmalr_y, ///< Target array for \f$ \sigma_{l,m}^\mathrm{L} \f$, must be `c->nelem_sc` long.
283 double *target_sigmalr_y_err, ///< Target array for error estimates, must be `c->nelem_sc` long or `NULL`.
284 const qpms_ewald3_constants_t *c, ///< Constant factors structure initialised by qpms_ewald3_constants_init().
285 double eta, ///< Ewald parameter.
286 complex double wavenumber, ///< Wavenumber of the background medium.
287 double unitcell_volume, ///< Volume of the (direct lattice) unit cell (with dimension corresponding to the lattice dimensionality).
288 /// Lattice dimensionality.
289 LatticeDimensionality latdim,
290 /// Lattice point generator for the reciprocal lattice.
291 /** There is a possibility that the whole PGen is not consumed
292 * (this might happen if the summand start to be consistently smaller
293 * than the (partial) sums * DBL_EPSILON.
294 * In such case, it is the responsibility of the caller to deallocate
295 * the generator.
297 PGen *pgen_K,
298 /// Indicates whether pgen_K already generates shifted points.
299 /** If false, the behaviour corresponds to the old ewald32_sigma_long_points_and_shift(),
300 * so the function assumes that the generated points correspond to the unshifted reciprocal Bravais lattice,
301 * and adds beta to the generated points before calculations.
302 * If true, it assumes that they are already shifted.
304 bool pgen_generates_shifted_points,
305 /// Wave vector \f$\vect k\f$.
306 cart3_t k,
307 /// Lattice offset \f$\vect s\f$ wrt. the Bravais lattice.
308 cart3_t particle_shift
311 #endif //EWALD_H