Initial qpms version metadata in library and output files
[qpms.git] / oldtests / ewaldshift3g_vargeom.c
blobf6ed8ffb874d821da15ea557144dd76d5338f0e3
1 // c99 -o ewaldshift3g_vargeom -ggdb -Wall -I ../ ewaldshift3g_vargeom.c ../qpms/ewald.c ../qpms/ewaldsf.c ../qpms/lattices2d.c ../qpms/latticegens.c -lgsl -lm -lblas
3 // implementation of the [LT(4.16)] test
4 #include <math.h>
5 #define M_SQRTPI 1.7724538509055160272981674833411452
6 #define M_SQRT3 1.7320508075688772935274463415058724
7 #include <qpms/ewald.h>
8 #include <qpms/tiny_inlines.h>
9 #include <qpms/indexing.h>
10 #include <stdlib.h>
11 #include <stdio.h>
12 #include <float.h>
13 #include <gsl/gsl_sf_legendre.h>
14 typedef struct ewaldtest2d_params {
15 qpms_l_t lMax;
16 point2d b1, b2;
17 point2d beta;
18 point2d particle_shift;
19 double k;
20 //double a;
21 double eta;
22 double maxR;
23 double maxK;
24 double csphase;
25 } ewaldtest2d_params;
27 typedef struct ewaldtest2d_results {
28 ewaldtest2d_params p;
29 complex double *sigmas_short,
30 *sigmas_long,
31 sigma0,
32 *sigmas_total;
33 double *err_sigmas_short,
34 *err_sigmas_long,
35 err_sigma0,
36 *err_sigmas_total;
37 complex double *regsigmas_416;
38 } ewaldtest2d_results;
42 const double a = 582e-9;
43 const double inv_a = 4*M_PI/a/M_SQRT3;
44 const double Klen = 4*M_PI/a/3;
47 #define AA (582.e-9)
48 #define INV_A (4*M_PI/AA/M_SQRT3)
49 #define KLEN (4*M_PI/AA/3)
51 ewaldtest2d_params paramslist[] = {
52 // lMax, b1, b2, beta, shift, k, eta, maxR, maxK, csphase
53 { 2, {1,0}, {0,1}, {2.7, 1}, {0.5,0.1325}, 2.3, 0.5, 20, 160, 1.},
54 { 2, {1,0}, {0,1}, {2.7, 1}, {0.5,0.1325}, 2.3, 1.5, 20, 160, 1.},
55 { 2, {1,0}, {0,1}, {2.7, 1}, {0.5,0.1325}, 2.3, 2.5, 20, 160, 1.},
56 { 2, {1,0}, {0,1}, {2.7, 1}, {0.5,0.1325}, 2.3, 3.5, 20, 160, 1.},
58 { 2, {1,0}, {0,1}, {1.7, 1}, {0.5,0.1325}, 2.3, 0.5, 20, 160, 1.},
59 { 2, {1,0}, {0,1}, {1.7, 1}, {0.5,0.1325}, 2.3, 1.5, 20, 160, 1.},
60 { 2, {1,0}, {0,1}, {1.7, 1}, {0.5,0.1325}, 2.3, 2.5, 20, 160, 1.},
61 { 2, {1,0}, {0,1}, {1.7, 1}, {0.5,0.1325}, 2.3, 3.5, 20, 160, 1.},
63 //{ 3, (AA,0}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 0.5 / AA, 20*AA, 160/AA, 1.},
64 { 3, {AA,0}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 2.5 / AA, 20*AA, 160/AA, 1.},
65 { 3, {AA,0}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 4.5 / AA, 20*AA, 160/AA, 1.},
66 { 3, {AA,0}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 6.5 / AA, 20*AA, 160/AA, 1.},
67 { 3, {AA,0}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 8.5 / AA, 20*AA, 160/AA, 1.},
69 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 2.5 / AA, 20*AA, 160/AA, 1.},
70 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 4.5 / AA, 20*AA, 160/AA, 1.},
71 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 6.5 / AA, 20*AA, 160/AA, 1.},
72 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, 8.5 / AA, 20*AA, 160/AA, 1.},
74 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {.1*AA,.03*AA}, 2.62 * 4 * M_PI/3/AA, 2.5 / AA, 20*AA, 160/AA, 1.},
75 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {.1*AA,.03*AA}, 2.62 * 4 * M_PI/3/AA, 4.5 / AA, 20*AA, 160/AA, 1.},
76 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {.1*AA,.03*AA}, 2.62 * 4 * M_PI/3/AA, 6.5 / AA, 20*AA, 160/AA, 1.},
77 { 3, {AA,.1*AA}, {0,AA}, {0.3*KLEN,KLEN}, {.1*AA,.03*AA}, 2.62 * 4 * M_PI/3/AA, 8.5 / AA, 20*AA, 160/AA, 1.},
82 #if 0
84 { 3, {0,4198609.6394310603}, {0,0}, 11255786.828366444, 9.9766126515967311e-07, 29088820.866572164, 20*9.9766126515967311e-07, 20*7272205.21664304, 1., TRIANGULAR_VERTICAL},
85 { 3, {0,4198609.6394310603}, {0,0}, 11255786.828366444, 9.9766126515967311e-07, 0.18*29088820.866572164, 20*9.9766126515967311e-07, 20*7272205.21664304, 1., TRIANGULAR_VERTICAL},
86 { 3, {0,4198609.6394310603}, {0,0}, 11255786.828366444, 9.9766126515967311e-07, 0.13*29088820.866572164, 20*9.9766126515967311e-07, 20*7272205.21664304, 1., TRIANGULAR_VERTICAL},
87 { 3, {0,4198609.6394310603}, {0,0}, 11255786.828366444, 9.9766126515967311e-07, 0.07*29088820.866572164, 20*9.9766126515967311e-07, 20*7272205.21664304, 1., TRIANGULAR_VERTICAL},
88 { 3, {0,4198609.6394310603}, {0,0}, 11255786.828366444, 9.9766126515967311e-07, 0.03*29088820.866572164, 20*9.9766126515967311e-07, 20*7272205.21664304, 1., TRIANGULAR_VERTICAL},
90 // { 3, {0,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, AA, 0.5 / AA, 20*AA, 160/AA, 1., TRIANGULAR_VERTICAL},
91 { 3, {0,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, AA, 2.5 / AA, 20*AA, 160/AA, 1., TRIANGULAR_VERTICAL},
92 { 3, {0,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, AA, 4.5 / AA, 20*AA, 160/AA, 1., TRIANGULAR_VERTICAL},
93 { 3, {0,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, AA, 6.5 / AA, 20*AA, 160/AA, 1., TRIANGULAR_VERTICAL},
94 { 3, {0,KLEN}, {0,0}, 2.62 * 4 * M_PI/3/AA, AA, 8.5 / AA, 20*AA, 160/AA, 1., TRIANGULAR_VERTICAL},
96 { 3, {0,Klen}, {0,0}, 2.62 * 4 * M_PI/3/a, a, 0.5 / a, 20*a, 2*M_PI*160/a, 1., TRIANGULAR_VERTICAL},
97 { 3, {0,Klen}, {0,0}, 2.62 * 4 * M_PI/3/a, a, 2.5 / a, 20*a, 2*M_PI*160/a, 1., TRIANGULAR_VERTICAL},
98 { 3, {0,Klen}, {0,0}, 2.62 * 4 * M_PI/3/a, a, 4.5 / a, 20*a, 2*M_PI*160/a, 1., TRIANGULAR_VERTICAL},
99 { 3, {0,Klen}, {0,0}, 2.62 * 4 * M_PI/3/a, a, 6.5 / a, 20*a, 2*M_PI*160/a, 1., TRIANGULAR_VERTICAL},
100 { 3, {0,Klen}, {0,0}, 2.62 * 4 * M_PI/3/a, a, 8.5 / a, 20*a, 2*M_PI*160/a, 1., TRIANGULAR_VERTICAL},
104 { 2, {1.1, 2.1}, {0.5,0.1325}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
105 { 2, {1.1, 2.1}, {0.5,0.1325}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
106 { 2, {1.1, 2.1}, {0.5,0.1325}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
107 { 2, {1.1, 2.1}, {0.5,0.1325}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
109 { 2, {0, 3.1}, {0.5,0}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
110 { 2, {0, 3.1}, {0.5,0}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
111 { 2, {0, 3.1}, {0.5,0}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
112 { 2, {0, 3.1}, {0.5,0}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
114 { 2, {0, 1.1}, {0.5,0}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
115 { 2, {0, 1.1}, {0.5,0}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
116 { 2, {0, 1.1}, {0.5,0}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
117 { 2, {0, 1.1}, {0.5,0}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
119 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
120 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
121 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
122 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
124 { 2, {1.1,0}, {0,0.5}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
125 { 2, {1.1,0}, {0,0.5}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
126 { 2, {1.1,0}, {0,0.5}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
127 { 2, {1.1,0}, {0,0.5}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
129 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
130 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
131 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
132 { 2, {3.1,0}, {0,0.5}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
134 { 2, {3.1*0.5,-3.1*0.8}, {0.8,0.5}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
135 { 2, {3.1*0.5,-3.1*0.8}, {0.8,0.5}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
136 { 2, {3.1*0.5,-3.1*0.8}, {0.8,0.5}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
137 { 2, {3.1*0.5,-3.1*0.8}, {0.8,0.5}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
139 { 2, {1.1*0.5,-1.1*0.8}, {0.8,0.5}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
140 { 2, {1.1*0.5,-1.1*0.8}, {0.8,0.5}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
141 { 2, {1.1*0.5,-1.1*0.8}, {0.8,0.5}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
142 { 2, {1.1*0.5,-1.1*0.8}, {0.8,0.5}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
144 // Poloviční posun oproti přodchozímu
145 { 2, {3.1*0.5,-3.1*0.8}, {0.4,0.25}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
146 { 2, {3.1*0.5,-3.1*0.8}, {0.4,0.25}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
147 { 2, {3.1*0.5,-3.1*0.8}, {0.4,0.25}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
148 { 2, {3.1*0.5,-3.1*0.8}, {0.4,0.25}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
150 { 2, {1.1*0.5,-1.1*0.8}, {0.4,0.25}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
151 { 2, {1.1*0.5,-1.1*0.8}, {0.4,0.25}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
152 { 2, {1.1*0.5,-1.1*0.8}, {0.4,0.25}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
153 { 2, {1.1*0.5,-1.1*0.8}, {0.4,0.25}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
155 // miniposun
156 { 3, {3.1*0.5,-3.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
157 { 3, {3.1*0.5,-3.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
158 { 3, {3.1*0.5,-3.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
159 { 3, {3.1*0.5,-3.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
161 { 3, {1.1*0.5,-1.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
162 { 3, {1.1*0.5,-1.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
163 { 3, {1.1*0.5,-1.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
164 { 3, {1.1*0.5,-1.1*0.8}, {0.004,0.0025}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
166 { 3, {3.1*0.5,-3.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
167 { 3, {3.1*0.5,-3.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
168 { 3, {3.1*0.5,-3.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
169 { 3, {3.1*0.5,-3.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
171 { 3, {1.1*0.5,-1.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
172 { 3, {1.1*0.5,-1.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
173 { 3, {1.1*0.5,-1.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
174 { 3, {1.1*0.5,-1.1*0.8}, {-0.004,-0.0025}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
180 { 2, {0, 3.1}, {0,0.5}, 2.3, 0.97, 0.5, 20, 160, 1., TRIANGULAR_VERTICAL},
181 { 2, {0, 3.1}, {0,0.5}, 2.3, 0.97, 1.5, 20, 160, 1., TRIANGULAR_VERTICAL},
182 { 2, {0, 3.1}, {0,0.5}, 2.3, 0.97, 2.5, 20, 160, 1., TRIANGULAR_VERTICAL},
183 { 2, {0, 3.1}, {0,0.5}, 2.3, 0.97, 3.5, 20, 160, 1., TRIANGULAR_VERTICAL},
187 #endif
189 // end:
190 // { 0, {0, 0}, 0, 0, 0, 0, 0, 0, 0}
193 void ewaldtest2d_results_free(ewaldtest2d_results *r) {
194 free(r->sigmas_short);
195 free(r->sigmas_long);
196 free(r->sigmas_total);
197 free(r->err_sigmas_long);
198 free(r->err_sigmas_total);
199 free(r->err_sigmas_short);
200 free(r->regsigmas_416);
201 free(r);
204 static inline double san(double x) {
205 return fabs(x) < 1e-13 ? 0 : x;
208 ewaldtest2d_results *ewaldtest2d(const ewaldtest2d_params p);
210 int main() {
211 gsl_set_error_handler(IgnoreUnderflowsGSLErrorHandler);
212 for (size_t i = 0; i < sizeof(paramslist)/sizeof(ewaldtest2d_params); ++i) {
213 ewaldtest2d_params p = paramslist[i];
214 ewaldtest2d_results *r = ewaldtest2d(p);
215 // TODO print per-test header here
216 printf("===============================\n");
217 printf("b1 = (%g, %g), b2 = (%g, %g)," /* "K1 = (%g, %g), K2 = (%g, %g),"*/ " Kmax = %g, Rmax = %g, lMax = %d, eta = %g, k = %g, beta = (%g,%g), ps = (%g,%g), csphase = %g\n",
218 p.b1.x, p.b1.y, p.b2.x, p.b2.y,/*TODO K1, K2*/ p.maxK, p.maxR, p.lMax, p.eta, p.k, p.beta.x, p.beta.y, p.particle_shift.x, p.particle_shift.y, p.csphase);
219 printf("sigma0: %.16g%+.16gj\n", creal(r->sigma0), cimag(r->sigma0));
220 for (qpms_l_t n = 0; n <= p.lMax; ++n) {
221 for (qpms_m_t m = -n; m <= n; ++m){
222 if ((m+n)%2) continue;
223 qpms_y_t y = qpms_mn2y_sc(m,n);
224 qpms_y_t y_conj = qpms_mn2y_sc(-m,n);
225 // y n m sigma_total (err), regsigmas_416 regsigmas_415_recon
226 printf("%zd %d %d: T:%.16g%+.16gj(%.3g) L:%.16g%+.16gj(%.3g) S:%.16g%+.16gj(%.3g) \n"
227 //"| predict %.16g%+.16gj \n| actual %.16g%+.16gj\n"
229 y, n, m, creal(san(r->sigmas_total[y])), san(cimag(r->sigmas_total[y])),
230 r->err_sigmas_total[y],
231 san(creal(r->sigmas_long[y])), san(cimag(r->sigmas_long[y])),
232 r->err_sigmas_long[y],
233 san(creal(r->sigmas_short[y])), san(cimag(r->sigmas_short[y])),
234 r->err_sigmas_short[y]
235 //san(creal(r->regsigmas_416[y])), san(cimag(r->regsigmas_416[y])),
236 //san(creal(r->sigmas_total[y]) + creal(r->sigmas_total[y_conj])),
237 //san(cimag(r->sigmas_total[y]) - cimag(r->sigmas_total[y_conj]))
241 ewaldtest2d_results_free(r);
243 return 0;
247 int ewaldtest_counter = 0;
250 ewaldtest2d_results *ewaldtest2d(const ewaldtest2d_params p) {
251 cart3_t beta3 = cart22cart3xy(p.beta);
252 cart3_t particle_shift3 = cart22cart3xy(p.particle_shift);
254 cart2_t b1 = p.b1, b2 = p.b2, rb1, rb2;
255 if (QPMS_SUCCESS != l2d_reciprocalBasis2pi(b1, b2, &rb1, &rb2))
256 abort();
258 const double A = l2d_unitcell_area(b1, b2); // sqrt(3) * a * a / 2.; // unit cell size
259 const double K_len = cart2norm(rb1)+cart2norm(rb2); //4*M_PI/a/sqrt(3); // reciprocal vector length
261 ewaldtest2d_results *results = malloc(sizeof(ewaldtest2d_results));
262 results->p = p;
264 // skip zeroth point if it coincides with origin
265 bool include_origin = !(fabs(p.particle_shift.x) == 0
266 && fabs(p.particle_shift.y) == 0);
268 PGen Rlgen = PGen_xyWeb_new(b1, b2, BASIS_RTOL, CART2_ZERO, 0, include_origin, p.maxR, false);
269 //PGen Rlgen_plus_shift = PGen_xyWeb_new(b1, b2, BASIS_RTOL, cart2_scale(-1 /* CHECKSIGN */, particle_shift2), 0, include_origin, p.maxR + a, false);
270 PGen Klgen = PGen_xyWeb_new(rb1, rb2, BASIS_RTOL, CART2_ZERO, 0, true, p.maxK + K_len, false);
271 //PGen Klgen_plus_beta = PGen_xyWeb_new(rb1, rb2, BASIS_RTOL, beta2, 0, true, p.maxK + K_len, false);
273 qpms_y_t nelem_sc = qpms_lMax2nelem_sc(p.lMax);
275 results->sigmas_short = malloc(sizeof(complex double)*nelem_sc);
276 results->sigmas_long = malloc(sizeof(complex double)*nelem_sc);
277 results->sigmas_total = malloc(sizeof(complex double)*nelem_sc);
278 results->err_sigmas_short = malloc(sizeof(double)*nelem_sc);
279 results->err_sigmas_long = malloc(sizeof(double)*nelem_sc);
280 results->err_sigmas_total = malloc(sizeof(double)*nelem_sc);
282 qpms_ewald32_constants_t *c = qpms_ewald32_constants_init(p.lMax, p.csphase);
284 if (0!=ewald3_sigma_long(results->sigmas_long,
285 results->err_sigmas_long, c, p.eta, p.k, A,
286 LAT_2D_IN_3D_XYONLY, &Klgen, false, beta3, particle_shift3))
287 abort();
288 if (0!=ewald3_sigma_short(
289 results->sigmas_short, results->err_sigmas_short, c,
290 p.eta, p.k, LAT_2D_IN_3D_XYONLY, &Rlgen, false, beta3, particle_shift3))
291 abort();
292 if (0!=ewald3_sigma0(&(results->sigma0), &(results->err_sigma0), c, p.eta, p.k))
293 abort();
294 for(qpms_y_t y = 0; y < nelem_sc; ++y) {
295 results->sigmas_total[y] = results->sigmas_short[y] + results->sigmas_long[y];
296 results->err_sigmas_total[y] = results->err_sigmas_short[y] + results->err_sigmas_long[y];
298 results->sigmas_total[0] += results->sigma0;
299 results->err_sigmas_total[0] += results->err_sigma0;
301 // Now calculate the reference values [LT(4.16)]
302 results->regsigmas_416 = calloc(nelem_sc, sizeof(complex double));
303 results->regsigmas_416[0] = -2 * c->legendre0[gsl_sf_legendre_array_index(0,0)];
305 #if 0 // not yet implemented for the new API
307 double legendres[gsl_sf_legendre_array_n(p.lMax)];
308 points2d_rordered_t sel =
309 points2d_rordered_annulus(Kpoints_plus_beta, 0, true, p.k, false);
310 if (0 != sel.nrs)
312 point2d *beta_pq_lessthan_k = sel.base + sel.r_offsets[0];
313 size_t beta_pq_lessthan_k_count = sel.r_offsets[sel.nrs] - sel.r_offsets[0];
314 for(size_t i = 0; i < beta_pq_lessthan_k_count; ++i) {
315 point2d beta_pq = beta_pq_lessthan_k[i];
316 double rbeta_pq = cart2norm(beta_pq);
317 double arg_pq = atan2(beta_pq.y, beta_pq.x);
318 double denom = sqrt(p.k*p.k - rbeta_pq*rbeta_pq);
319 if( gsl_sf_legendre_array_e(GSL_SF_LEGENDRE_NONE,
320 p.lMax, denom/p.k, p.csphase, legendres) != 0)
321 abort();
322 for (qpms_y_t y = 0; y < nelem_sc; ++y) {
323 qpms_l_t n; qpms_m_t m;
324 qpms_y2mn_sc_p(y, &m, &n);
325 if ((m+n)%2 != 0)
326 continue;
327 complex double eimf = cexp(I*m*arg_pq);
328 results->regsigmas_416[y] +=
329 4*M_PI*ipow(n)/p.k/A
330 * eimf * legendres[gsl_sf_legendre_array_index(n,abs(m))] * min1pow_m_neg(m)
331 / denom;
336 #else
337 for(qpms_y_t y = 0; y < nelem_sc; ++y) {
338 qpms_l_t n; qpms_m_t m;
339 qpms_y2mn_sc_p(y, &m, &n);
340 if ((m+n)%2 != 0)
341 continue;
342 results->regsigmas_416[y] = NAN;
344 #endif
347 qpms_ewald32_constants_free(c);
348 ++ewaldtest_counter;
349 return results;