libcrystfel 0.11.0
Loading...
Searching...
No Matches
Data Structures | Macros | Typedefs | Enumerations | Functions | Variables
utils.h File Reference

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)   (((val) | (bit)) ^ (bit))
 
#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))
 
#define ph_en_to_lambda(a)   ((PLANCK*C_VACUO)/(a))
 
#define ph_lambda_to_en(a)   ((PLANCK*C_VACUO)/(a))
 
#define eV_to_J(a)   ((a)*ELECTRON_CHARGE)
 
#define J_to_eV(a)   ((a)/ELECTRON_CHARGE)
 
#define ph_lambda_to_eV(a)   J_to_eV(ph_lambda_to_en(a))
 
#define ph_eV_to_lambda(a)   ph_en_to_lambda(eV_to_J(a))
 
#define ph_eV_to_k(a)   ((a)*ELECTRON_CHARGE/PLANCK/C_VACUO)
 
#define UNUSED
 
#define likely(x)   (x)
 
#define unlikely(x)   (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_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_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)
 
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
 

Detailed Description

Miscellaneous utility functions

Enumeration Type Documentation

◆ AssplodeFlag

Controls the behaviour of assplode.

Enumerator
ASSPLODE_NONE 

Nothing

ASSPLODE_DUPS 

Don't merge deliminators

Function Documentation

◆ assplode()

int assplode ( const char *  a,
const char *  delims,
char ***  pbits,
AssplodeFlag  flags 
)

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 free() when finished with. The array of bits also needs to be freed with free() when finished with, unless n=0 in which case bits==NULL

◆ normalise_quaternion()

struct quaternion normalise_quaternion ( struct quaternion  q)
Parameters
qA quaternion

Rescales the quaternion such that its modulus is unity.

Returns
The normalised version of q

◆ quat_rot()

struct rvec quat_rot ( struct rvec  q,
struct quaternion  z 
)
Parameters
qA vector (in the form of an rvec struct)
zA quaternion

Rotates a vector according to a quaternion.

Returns
a rotated version of p.

◆ quaternion_modulus()

double quaternion_modulus ( struct quaternion  q)
Parameters
qA quaternion

If a quaternion represents a pure rotation, its modulus should be unity.

Returns
The modulus of the given quaternion.

◆ quaternion_valid()

int quaternion_valid ( struct quaternion  q)
Parameters
qA 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.

Returns
1 if the quaternion is normalised, 0 if not.

◆ random_quaternion()

struct quaternion random_quaternion ( gsl_rng *  rng)
Parameters
rngA GSL random number generator to use
Returns
A randomly generated, normalised, quaternion.

◆ show_matrix()

void show_matrix ( gsl_matrix *  M)
Parameters
MA matrix

Displays a matrix.

◆ show_matrix_eqn()

void show_matrix_eqn ( gsl_matrix *  M,
gsl_vector *  v 
)
Parameters
MA matrix
vA vector

Displays a matrix equation of the form M.a = v.

◆ solve_svd()

gsl_vector * solve_svd ( gsl_vector *  v,
gsl_matrix *  M,
int *  pn_filt,
int  verbose 
)
Parameters
va gsl_vector
Ma gsl_matrix
pn_filtpointer to store the number of filtered eigenvalues
verboseflag for verbosity on the terminal

Solves the matrix equation M.x = v, returning x. Performs rescaling and eigenvalue filtering.