1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
|
/*
* utils.h
*
* Utility stuff
*
* Copyright © 2012-2014 Deutsches Elektronen-Synchrotron DESY,
* a research centre of the Helmholtz Association.
*
* Authors:
* 2009-2014 Thomas White <taw@physics.org>
*
* This file is part of CrystFEL.
*
* CrystFEL is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* CrystFEL is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with CrystFEL. If not, see <http://www.gnu.org/licenses/>.
*
*/
#ifndef UTILS_H
#define UTILS_H
#ifdef HAVE_CONFIG_H
#include <config.h>
#endif
#include <math.h>
#include <complex.h>
#include <string.h>
#include <stdlib.h>
#include <pthread.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_rng.h>
#include "thread-pool.h"
/* -------------------------- Fundamental constants ------------------------ */
/* Electron charge in C */
#define ELECTRON_CHARGE (1.6021773e-19)
/* Planck's constant (Js) */
#define PLANCK (6.62606896e-34)
/* Speed of light in vacuo (m/s) */
#define C_VACUO (299792458)
/* Thomson scattering length (m) */
#define THOMSON_LENGTH (2.81794e-15)
/* ------------------------------ Quaternions ------------------------------- */
/**
* quaternion:
*
* <programlisting>
* struct quaternion
* {
* double w
* double x
* double y
* double z
* };
* </programlisting>
*
* A structure representing a quaternion.
*
**/
struct quaternion;
struct quaternion {
double w;
double x;
double y;
double z;
};
#ifdef __cplusplus
extern "C" {
#endif
extern struct quaternion normalise_quaternion(struct quaternion q);
extern double quaternion_modulus(struct quaternion q);
extern struct quaternion random_quaternion(gsl_rng *rng);
extern int quaternion_valid(struct quaternion q);
extern struct rvec quat_rot(struct rvec q, struct quaternion z);
/* --------------------------- Useful functions ----------------------------- */
extern void show_matrix_eqn(gsl_matrix *M, gsl_vector *v);
extern void show_matrix(gsl_matrix *M);
extern size_t notrail(char *s);
extern void chomp(char *s);
/**
* AssplodeFlag:
* @ASSPLODE_NONE: Nothing
* @ASSPLODE_DUPS: Don't merge deliminators
**/
typedef enum {
ASSPLODE_NONE = 0,
ASSPLODE_DUPS = 1<<0
} AssplodeFlag;
extern int assplode(const char *a, const char *delims, char ***pbits,
AssplodeFlag flags);
extern void progress_bar(int val, int total, const char *text);
extern double random_flat(gsl_rng *rng, double max);
extern double flat_noise(gsl_rng *rng, double expected, double width);
extern double gaussian_noise(gsl_rng *rng, double expected, double stddev);
extern int poisson_noise(gsl_rng *rng, double expected);
/* Keep these ones inline, to avoid function call overhead */
static inline struct quaternion invalid_quaternion(void)
{
struct quaternion quat;
quat.w = 0.0;
quat.x = 0.0;
quat.y = 0.0;
quat.z = 0.0;
return quat;
}
static inline double distance(double x1, double y1, double x2, double y2)
{
return sqrt((x2-x1)*(x2-x1) + (y2-y1)*(y2-y1));
}
static inline double modulus(double x, double y, double z)
{
return sqrt(x*x + y*y + z*z);
}
static inline double modulus2d(double x, double y)
{
return sqrt(x*x + y*y);
}
static inline double modulus_squared(double x, double y, double z) {
return x*x + y*y + z*z;
}
static inline double distance3d(double x1, double y1, double z1,
double x2, double y2, double z2)
{
double d = modulus(x1-x2, y1-y2, z1-z2);
return d;
}
/* Answer in radians */
static inline double angle_between(double x1, double y1, double z1,
double x2, double y2, double z2)
{
double mod1 = modulus(x1, y1, z1);
double mod2 = modulus(x2, y2, z2);
double cosine = (x1*x2 + y1*y2 + z1*z2) / (mod1*mod2);
/* Fix domain if outside due to rounding */
if ( cosine > 1.0 ) cosine = 1.0;
if ( cosine < -1.0 ) cosine = -1.0;
return acos(cosine);
}
/* Answer in radians */
static inline double angle_between_2d(double x1, double y1,
double x2, double y2)
{
double mod1 = modulus2d(x1, y1);
double mod2 = modulus2d(x2, y2);
double cosine = (x1*x2 + y1*y2) / (mod1*mod2);
/* Fix domain if outside due to rounding */
if ( cosine > 1.0 ) cosine = 1.0;
if ( cosine < -1.0 ) cosine = -1.0;
return acos(cosine);
}
static inline int within_tolerance(double a, double b, double percent)
{
double tol = fabs(a) * (percent/100.0);
if ( fabs(b-a) < tol ) return 1;
return 0;
}
/* ----------------------------- Useful macros ------------------------------ */
#define rad2deg(a) ((a)*180/M_PI)
#define deg2rad(a) ((a)*M_PI/180)
#define is_odd(a) ((a)%2==1)
#define biggest(a,b) ((a>b) ? (a) : (b))
#define smallest(a,b) ((a<b) ? (a) : (b))
/* Photon energy (J) to wavelength (m) */
#define ph_en_to_lambda(a) ((PLANCK*C_VACUO)/(a))
/* Photon wavelength (m) to energy (J) */
#define ph_lambda_to_en(a) ((PLANCK*C_VACUO)/(a))
/* eV to Joules */
#define eV_to_J(a) ((a)*ELECTRON_CHARGE)
/* Joules to eV */
#define J_to_eV(a) ((a)/ELECTRON_CHARGE)
/* Photon wavelength (m) to energy (eV) */
#define ph_lambda_to_eV(a) J_to_eV(ph_lambda_to_en(a))
/* Photon energy (eV) to wavelength (m) */
#define ph_eV_to_lambda(a) ph_en_to_lambda(eV_to_J(a))
#define UNUSED __attribute__((unused))
/* ------------------------------ Message macros ---------------------------- */
extern pthread_mutex_t stderr_lock;
#define ERROR(...) { \
int error_print_val = get_status_label(); \
pthread_mutex_lock(&stderr_lock); \
if ( error_print_val >= 0 ) { \
fprintf(stderr, "%3i: ", error_print_val); \
} \
fprintf(stderr, __VA_ARGS__); \
pthread_mutex_unlock(&stderr_lock); \
}
#define STATUS(...) { \
int status_print_val = get_status_label(); \
pthread_mutex_lock(&stderr_lock); \
if ( status_print_val >= 0 ) { \
fprintf(stderr, "%3i: ", status_print_val); \
} \
fprintf(stderr, __VA_ARGS__); \
pthread_mutex_unlock(&stderr_lock); \
}
/* ------------------------------ File handling ----------------------------- */
extern char *check_prefix(char *prefix);
extern char *safe_basename(const char *in);
#ifdef __cplusplus
}
#endif
#endif /* UTILS_H */
|