![]() |
libcrystfel 0.11.1-244-g4019144e+
|
Data Structures | |
| struct | quaternion |
Macros | |
| #define | ELECTRON_CHARGE (1.6021773e-19) |
| #define | ELECTRON_MASS (9.1093837015e-31) |
| #define | PLANCK (6.62606896e-34) |
| #define | C_VACUO (299792458) |
| #define | THOMSON_LENGTH (2.81794e-15) |
| #define | CLEAR_BIT(val, bit) |
| #define | rad2deg(a) |
| #define | deg2rad(a) |
| #define | is_odd(a) |
| #define | biggest(a, b) |
| #define | smallest(a, b) |
| #define | ph_en_to_lambda(a) |
| #define | ph_lambda_to_en(a) |
| #define | eV_to_J(a) |
| #define | J_to_eV(a) |
| #define | ph_lambda_to_eV(a) |
| #define | ph_eV_to_lambda(a) |
| #define | ph_eV_to_k(a) |
| #define | UNUSED |
| #define | likely(x) |
| #define | unlikely(x) |
Typedefs | |
| typedef void(* | LogMsgFunc) (enum log_msg_type type, const char *msg, void *vp) |
Enumerations | |
| enum | AssplodeFlag { ASSPLODE_NONE = 0 , ASSPLODE_DUPS = 1<<0 } |
| enum | log_msg_type { LOG_MSG_STATUS , LOG_MSG_ERROR } |
Functions | |
| void | show_matrix_eqn (gsl_matrix *M, gsl_vector *v) |
| void | show_matrix (gsl_matrix *M) |
| void | show_vector (gsl_vector *M) |
| gsl_vector * | solve_inv (gsl_vector *v, gsl_matrix *M) |
| gsl_vector * | solve_svd (gsl_vector *v, gsl_matrix *M, int *n_filt, int verbose) |
| gsl_matrix * | matrix_mult2 (gsl_matrix *A, gsl_matrix *B) |
| gsl_matrix * | matrix_mult3 (gsl_matrix *A, gsl_matrix *B, gsl_matrix *C) |
| gsl_matrix * | matrix_invert (gsl_matrix *m) |
| size_t | notrail (char *s) |
| int | convert_int (const char *str, int *pval) |
| int | convert_long_int (const char *str, long long int *pval) |
| int | convert_float (const char *str, double *pval) |
| void | chomp (char *s) |
| void * | srealloc (void *arr, size_t new_size) |
| int | assplode (const char *a, const char *delims, char ***pbits, AssplodeFlag flags) |
| void | progress_bar (int val, int total, const char *text) |
| double | random_flat (gsl_rng *rng, double max) |
| double | flat_noise (gsl_rng *rng, double expected, double width) |
| double | gaussian_noise (gsl_rng *rng, double expected, double stddev) |
| int | poisson_noise (gsl_rng *rng, double expected) |
| void | rotate2d (double *x, double *y, double cx, double cy, double ang) |
| void | STATUS (const char *format,...) |
| void | ERROR (const char *format,...) |
| void | set_log_message_func (LogMsgFunc new_log_msg_func, void *vp) |
| void * | cfmalloc (size_t size) |
| void | cffree (void *ptr) |
| void * | cfcalloc (size_t nmemb, size_t size) |
| void * | cfrealloc (void *ptr, size_t size) |
| char * | cfstrdup (const char *s) |
| char * | cfstrndup (const char *s, size_t n) |
| int | set_mm_funcs (void *(*cfmalloc)(size_t size), void(*cffree)(void *ptr), void *(*cfcalloc)(size_t nmemb, size_t size), void *(*cfrealloc)(void *ptr, size_t size)) |
| void | set_last_task (const char *task) |
| int | notify_alive (void) |
| int | set_debug_funcs (void(*slt)(const char *, void *vp), int(*ping)(void *vp), void *debug_data) |
| char * | check_prefix (char *prefix) |
| char * | safe_basename (const char *in) |
| char * | safe_strdup (const char *in) |
| void | strip_extension (char *bfn) |
| char * | load_entire_file (const char *filename) |
| int | file_exists (const char *filename) |
| int | is_dir (const char *filename) |
| const char * | filename_extension (const char *fn, const char **ext2) |
| struct quaternion | normalise_quaternion (struct quaternion q) |
| double | quaternion_modulus (struct quaternion q) |
| struct quaternion | random_quaternion (gsl_rng *rng) |
| int | quaternion_valid (struct quaternion q) |
| struct rvec | quat_rot (struct rvec q, struct quaternion z) |
| int | compare_double (const void *av, const void *bv) |
| int | crystfel_has_peakfinder9 (void) |
Variables | |
| pthread_mutex_t | stderr_lock |
Miscellaneous utility functions
| #define biggest | ( | a, | |
| b ) |
| #define CLEAR_BIT | ( | val, | |
| bit ) |
| #define deg2rad | ( | a | ) |
| #define eV_to_J | ( | a | ) |
| #define is_odd | ( | a | ) |
| #define J_to_eV | ( | a | ) |
| #define likely | ( | x | ) |
| #define ph_en_to_lambda | ( | a | ) |
| #define ph_eV_to_k | ( | a | ) |
| #define ph_eV_to_lambda | ( | a | ) |
| #define ph_lambda_to_en | ( | a | ) |
| #define ph_lambda_to_eV | ( | a | ) |
| #define rad2deg | ( | a | ) |
| #define smallest | ( | a, | |
| b ) |
| #define unlikely | ( | x | ) |
| enum AssplodeFlag |
|
extern |
Split the string 'a' using 'delims' as a zero-terminated list of deliminators. Store each segment in bits[0...n] where n is the number of segments and is the return value. pbits = &bits Each segment needs to be freed with cffree() when finished with. The array of bits also needs to be freed with cffree() when finished with, unless n=0 in which case bits==NULL
|
extern |
| q | A quaternion |
Rescales the quaternion such that its modulus is unity.
q
|
extern |
| q | A vector (in the form of an rvec struct) |
| z | A quaternion |
Rotates a vector according to a quaternion.
p.
|
extern |
| q | A quaternion |
If a quaternion represents a pure rotation, its modulus should be unity.
|
extern |
| q | A quaternion |
Checks if the given quaternion is normalised.
This function performs a nasty floating point comparison of the form (modulus > 0.999) && (modulus < 1.001), and so should not be relied upon to spot anything other than the most obvious input error.
|
extern |
| rng | A GSL random number generator to use |
|
extern |
| M | A matrix |
Displays a matrix.
|
extern |
| M | A matrix |
| v | A vector |
Displays a matrix equation of the form M.a = v.
|
extern |
| v | a gsl_vector |
| M | a gsl_matrix |
| pn_filt | pointer to store the number of filtered eigenvalues |
| verbose | flag for verbosity on the terminal |
Solves the matrix equation M.x = v, returning x. Performs rescaling and eigenvalue filtering.