PAWpySeed
Parallel C/Python package for numerical analysis of PAW DFT wavefunctions
utils.h
Go to the documentation of this file.
1 
17 #ifndef UTILS_H
18 #define UTILS_H
19 #include <complex.h>
20 #include <math.h>
21 
30 typedef struct funcset {
31  int l;
32  double* proj;
33  double** proj_spline;
34  double* aewave;
35  double** aewave_spline;
36  double* pswave;
37  double** pswave_spline;
38  double* diffwave;
39  double** diffwave_spline;
40  double* kwave;
41  double** kwave_spline;
42  double* smooth_diffwave;
44  double* dense_kwave;
46 } funcset_t;
47 
48 typedef struct ppot {
49  int num_projs;
51  int lmax;
53  double rmax; //< maximum radius of the projector functions
54  double wave_rmax;
61  double* wave_grid;
62  double* kwave_grid;
63  double* proj_grid;
64  double* smooth_grid;
65  double* dense_kgrid;
66 } ppot_t;
67 
68 typedef struct projection {
69  int num_projs;
71  int* ns;
72  int* ls;
73  int* ms;
74  double complex* overlaps;
75 } projection_t;
76 
82 typedef struct band {
83  int n;
84  int num_waves;
85  double occ;
86  double N;
87  double complex energy;
88  float complex* Cs;
89  double complex* CRs;
92 } band_t;
93 
94 typedef struct rayleigh_set {
95  int l;
96  double complex* terms;
98 
99 typedef struct kpoint {
100  short int up;
101  int num_waves;
102  int* Gs;
103  double* k;
104  double weight;
105  int num_bands;
108 } kpoint_t;
109 
110 typedef struct pswf {
114  int* G_bounds;
116  int nspin;
117  int nband;
118  int nwk;
119  double* lattice;
120  double* reclattice;
121  int* fftg;
123  double* dcoords;
124  double complex** overlaps;
125 } pswf_t;
126 
127 typedef struct projgrid {
128  double complex* values;
129 } projgrid_t;
130 
131 typedef struct real_proj {
132  int l;
133  int m;
134  int func_num;
135  double* paths;
136  double complex* values;
137 } real_proj_t;
138 
139 typedef struct real_proj_site {
140  int index;
141  int elem;
145  int gridsize;
146  double rmax;
147  double* coord;
148  int* indices;
151 
155 int min(int a, int b);
156 
160 int max(int a, int b);
161 
167 void vcross(double* res, double* top, double* bottom);
168 
173 double dot(double* x1, double* x2);
174 
178 double mag(double* x1);
179 
184 double determinant(double* m);
185 
194 double dist_from_frac(double* coords1, double* coords2, double* lattice);
195 
196 void frac_to_cartesian(double* coord, double* lattice);
197 
198 void cartesian_to_frac(double* coord, double* reclattice);
199 
200 void min_cart_path(double* coord, double* center, double* lattice, double* path, double* r);
201 
202 double complex trilinear_interpolate(double complex* c, double* frac, int* fftg);
203 
204 void free_kpoint(kpoint_t* kpt, int num_elems, ppot_t* pps);
205 
206 void free_ppot(ppot_t* pp);
207 
209 
211 
212 void free_pswf(pswf_t* wf);
213 
214 void free_ptr(void* ptr);
215 
216 void free_real_proj_site_list(real_proj_site_t* sites, int length);
217 
218 void free_ppot_list(ppot_t* pps, int length);
219 
227 double* get_occs(pswf_t* wf);
228 
230 int get_nband(pswf_t* wf);
231 
233 int get_nwk(pswf_t* wf);
234 
236 int get_nspin(pswf_t* wf);
237 
239 void set_num_sites(pswf_t* wf, int nsites);
240 
242 double legendre(int l, int m, double x);
243 
245 double fac(int n);
246 
251 double complex Ylm(int l, int m, double theta, double phi);
252 
256 double complex Ylm2(int l, int m, double costheta, double phi);
257 
264 double proj_interpolate(double r, double rmax, int size, double* x, double* proj, double** proj_spline);
265 
271 double wave_interpolate(double r, int size, double* x, double* f, double** wave_spline);
272 
276 double complex proj_value_helper(double r, double rmax, int size,
277  double* pos, double* x, double* f, double** s, int l, int m);
278 
283 double complex proj_value(funcset_t funcs, double* x, int m, double rmax,
284  int size, double* ion_pos, double* pos, double* lattice);
285 
290 double complex smooth_wave_value(funcset_t funcs, double* x, int m, double rmax,
291  int size, double* ion_pos, double* pos, double* lattice);
292 
297 double complex wave_value(funcset_t funcs, int size, double* x, int m,
298  double* ion_pos, double* pos, double* lattice);
299 
304 double complex wave_value2(double* x, double* wave, double** spline, int size,
305  int l, int m, double* pos);
306 
310 void setup_site(real_proj_site_t* sites, ppot_t* pps, int num_sites, int* site_nums,
311  int* labels, double* coords, double* lattice, int* fftg, int pr0_pw1);
312 
319 double** spline_coeff(double* x, double* y, int N);
320 
327 double spline_integral(double* x, double* a, double** s, int size);
328 
329 void frac_from_index(int index, double* coord, int* fftg);
330 
331 double sph_bessel(double k, double r, int l);
332 
336 double sbf(double x, int l);
337 
343 double complex rayexp(double* kpt, int* Gs, float complex* Cs, int l, int m,
344  int num_waves, double complex* sum_terms, double* ionp);
345 
350 double complex* rayexp_terms(double* kpt, int* Gs, int num_waves,
351  int l, int wave_gridsize, double* grid,
352  double* wave, double** spline, double* reclattice);
353 
359 void generate_rayleigh_expansion_terms(pswf_t* wf, ppot_t* pps, int num_elems);
360 
364 void copy_rayleigh_expansion_terms(pswf_t* wf, ppot_t* pps, int num_elems, pswf_t* wf_R);
365 
370 void CHECK_ALLOCATION(void* ptr);
371 
375 void ALLOCATION_FAILED();
376 
377 #endif
int * G_bounds
Definition: utils.h:114
int get_nband(pswf_t *wf)
int total_projs
number of projector functions
Definition: utils.h:70
double rmax
Definition: utils.h:146
band_t ** bands
bands with this k-point
Definition: utils.h:106
double * aewave
all electron partial wave
Definition: utils.h:34
real_proj_t * projs
Definition: utils.h:149
struct projection projection_t
double dist_from_frac(double *coords1, double *coords2, double *lattice)
short int up
spin
Definition: utils.h:100
void ALLOCATION_FAILED()
grid
Definition: rayleigh.py:30
double ** spline_coeff(double *x, double *y, int N)
double * lattice
Definition: utils.h:119
double complex Ylm2(int l, int m, double costheta, double phi)
double * smooth_grid
Definition: utils.h:64
double occ
occupancy of the band
Definition: utils.h:85
double sbf(double x, int l)
int nspin
Definition: utils.h:116
double * proj_grid
real radial grid for projector functions
Definition: utils.h:63
double ** smooth_diffwave_spline
spline coefficients for smooth_diffwave
Definition: utils.h:43
double complex * overlaps
list of <p_i|psi>
Definition: utils.h:74
double * pspw_overlap_matrix
overlap matrix for pseudo partial waves
Definition: utils.h:55
double complex wave_value(funcset_t funcs, int size, double *x, int m, double *ion_pos, double *pos, double *lattice)
int num_waves
number of plane waves
Definition: utils.h:84
double fac(int n)
double complex Ylm(int l, int m, double theta, double phi)
double * smooth_diffwave
diffwave on linear grid with high-frequency components removed
Definition: utils.h:42
int index
Definition: utils.h:140
double * dense_kwave
Definition: utils.h:44
double dot(double *x1, double *x2)
double complex proj_value(funcset_t funcs, double *x, int m, double rmax, int size, double *ion_pos, double *pos, double *lattice)
double weight
k-point weight
Definition: utils.h:104
int * Gs
plane wave coefficients, in sets of three
Definition: utils.h:102
Definition: utils.h:99
double * coord
Definition: utils.h:147
int nwk
Definition: utils.h:118
int n
band number
Definition: utils.h:83
Definition: utils.h:68
Definition: utils.h:48
y
Definition: rayleigh.py:27
int num_projs
Definition: utils.h:142
int num_projs
number of radial projector functions
Definition: utils.h:49
int total_projs
Definition: utils.h:143
double * diffwave
aewave-pswave
Definition: utils.h:38
double complex proj_value_helper(double r, double rmax, int size, double *pos, double *x, double *f, double **s, int l, int m)
void free_ppot(ppot_t *pp)
double * get_occs(pswf_t *wf)
projection_t * projections
length==number of sites in structure
Definition: utils.h:90
int m
Definition: utils.h:133
double N
Definition: utils.h:86
f
Definition: gaunt.py:28
double determinant(double *m)
int num_projs
number of radial projector functions
Definition: utils.h:69
void setup_site(real_proj_site_t *sites, ppot_t *pps, int num_sites, int *site_nums, int *labels, double *coords, double *lattice, int *fftg, int pr0_pw1)
double ** diffwave_spline
spline coefficients for diffwave
Definition: utils.h:39
double mag(double *x1)
double * proj
projector function
Definition: utils.h:32
Definition: utils.h:131
void generate_rayleigh_expansion_terms(pswf_t *wf, ppot_t *pps, int num_elems)
double * k
k-point vector
Definition: utils.h:103
double complex trilinear_interpolate(double complex *c, double *frac, int *fftg)
void free_real_proj_site(real_proj_site_t *site)
double complex rayexp(double *kpt, int *Gs, float complex *Cs, int l, int m, int num_waves, double complex *sum_terms, double *ionp)
int * ns
radial projector index
Definition: utils.h:71
double complex energy
energy of the band
Definition: utils.h:87
void free_ptr(void *ptr)
void vcross(double *res, double *top, double *bottom)
float complex * Cs
plane wave coefficients (normalized to 1)
Definition: utils.h:88
int l
angular momentum quantum number
Definition: utils.h:95
struct real_proj real_proj_t
list x
Definition: quadrature.py:9
double complex * rayexp_terms(double *kpt, int *Gs, int num_waves, int l, int wave_gridsize, double *grid, double *wave, double **spline, double *reclattice)
int nband
Definition: utils.h:117
struct band band_t
Definition: utils.h:127
int num_indices
Definition: utils.h:144
void free_kpoint(kpoint_t *kpt, int num_elems, ppot_t *pps)
void set_num_sites(pswf_t *wf, int nsites)
int num_bands
number of bands
Definition: utils.h:105
rayleigh_set_t ** expansion
Definition: utils.h:107
struct real_proj_site real_proj_site_t
double complex ** overlaps
Definition: utils.h:124
double complex * values
Definition: utils.h:136
double legendre(int l, int m, double x)
int gridsize
Definition: utils.h:145
void free_pswf(pswf_t *wf)
double ** pswave_spline
ps partial wave spline coefficients
Definition: utils.h:37
int elem
Definition: utils.h:141
struct funcset funcset_t
int num_cart_gridpts
number of real space grid points that can fit in the projector sphere
Definition: utils.h:60
double complex * values
Definition: utils.h:128
void free_real_proj_site_list(real_proj_site_t *sites, int length)
double ** aewave_spline
ae partial wave spline coefficients
Definition: utils.h:35
int num_aug_overlap_sites
Definition: utils.h:122
double * paths
Definition: utils.h:135
double sph_bessel(double k, double r, int l)
double ** kwave_spline
spline coefficients for kwave
Definition: utils.h:41
struct projgrid projgrid_t
Definition: utils.h:30
double rmax
Definition: utils.h:53
void cartesian_to_frac(double *coord, double *reclattice)
int max(int a, int b)
int lmax
maximum l-value of any projector
Definition: utils.h:51
Definition: utils.h:110
int get_nwk(pswf_t *wf)
double * kwave_grid
reciprocal radial grid for partial waves
Definition: utils.h:62
void copy_rayleigh_expansion_terms(pswf_t *wf, ppot_t *pps, int num_elems, pswf_t *wf_R)
double complex wave_value2(double *x, double *wave, double **spline, int size, int l, int m, double *pos)
int * fftg
Definition: utils.h:121
void free_ppot_list(ppot_t *pps, int length)
void frac_from_index(int index, double *coord, int *fftg)
struct kpoint kpoint_t
int proj_gridsize
number of points on projector radial grid
Definition: utils.h:58
struct rayleigh_set rayleigh_set_t
int func_num
Definition: utils.h:134
int total_projs
number of projector functions
Definition: utils.h:50
int * ls
l values of projectors
Definition: utils.h:72
struct pswf pswf_t
struct ppot ppot_t
double proj_interpolate(double r, double rmax, int size, double *x, double *proj, double **proj_spline)
int l
Definition: utils.h:132
double spline_integral(double *x, double *a, double **s, int size)
double * dense_kgrid
Definition: utils.h:65
int * indices
Definition: utils.h:148
Definition: utils.h:82
int wave_gridsize
number of points on partial wave radial grid
Definition: utils.h:59
int k
Definition: rayleigh.py:4
void free_real_proj(real_proj_t *proj)
int num_elems
Definition: utils.h:111
r
Definition: rayleigh.py:38
double * reclattice
Definition: utils.h:120
int l
l quantum number
Definition: utils.h:31
double * aepw_overlap_matrix
overlap matrix for all electron partial waves
Definition: utils.h:56
double * dcoords
Definition: utils.h:123
Definition: utils.h:94
kpoint_t ** kpts
Definition: utils.h:115
double * kwave
Expansion of diffwave in spherical Bessel functions.
Definition: utils.h:40
int * ms
m values of projectors
Definition: utils.h:73
void CHECK_ALLOCATION(void *ptr)
double ** dense_kwave_spline
Definition: utils.h:45
void min_cart_path(double *coord, double *center, double *lattice, double *path, double *r)
double wave_interpolate(double r, int size, double *x, double *f, double **wave_spline)
ppot_t * pps
Definition: utils.h:113
funcset_t * funcs
funcset for each projector, see funcset
Definition: utils.h:52
double complex * CRs
wavefunction in real space
Definition: utils.h:89
int min(int a, int b)
int num_sites
Definition: utils.h:112
double * wave_grid
real radial grid for partial waves
Definition: utils.h:61
Definition: utils.h:139
double * pswave
pseudo partial wave
Definition: utils.h:36
double complex smooth_wave_value(funcset_t funcs, double *x, int m, double rmax, int size, double *ion_pos, double *pos, double *lattice)
void frac_to_cartesian(double *coord, double *lattice)
double wave_rmax
maximum radius of the partial waves
Definition: utils.h:54
double complex * terms
rayleigh epansion terms
Definition: utils.h:96
double * diff_overlap_matrix
overlap matrix of difference between all electron and partial waves
Definition: utils.h:57
projection_t * wave_projections
length==number of sites in structure
Definition: utils.h:91
int get_nspin(pswf_t *wf)
double ** proj_spline
projector function spline
Definition: utils.h:33
int num_waves
number of plane waves in a band
Definition: utils.h:101