Hawaii Hybrid
|
#include <assert.h>
#include <errno.h>
#include <stdarg.h>
#include <stdbool.h>
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <regex.h>
#include <time.h>
#include <unistd.h>
#include <gsl/gsl_histogram.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_multifit_nlinear.h>
#include <gsl/gsl_spline.h>
#include <gsl/gsl_fft_real.h>
#include <gsl/gsl_fft_complex.h>
#include <gsl/gsl_fft_halfcomplex.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_randist.h>
#include <math.h>
#include <complex.h>
#include "mtwist.h"
#include <cvode/cvode.h>
#include <nvector/nvector_serial.h>
#include "array.h"
#include "constants.h"
#include "arena.h"
Go to the source code of this file.
Classes | |
struct | StoredHistogram |
struct | Monomer |
Monomer represents a monomer in a pair with associated dynamic variables. The order of variables in the qp array is specified by the following indices: More... | |
struct | MoleculeSystem |
Angular variables and conjugated momenta are stored within the MoleculeSystem struct in the same order as within the qp of Monomer struct. More... | |
struct | CalcParams |
struct | CFnc |
struct | CFncs |
struct | CFncArray |
struct | SFnc |
struct | Spectrum |
struct | Spectra |
struct | WingParams |
struct | WingData |
struct | Tracker |
struct | String_Builder |
Macros | |
#define | OMPI_SKIP_MPICXX |
#define | _GNU_SOURCE |
#define | MT_GENERATE_CODE_IN_HEADER 0 |
#define | IPHI 0 |
= 0 | |
#define | IPPHI 1 |
= 1 | |
#define | ITHETA 2 |
= 2 | |
#define | IPTHETA 3 |
= 3 | |
#define | IR 4 |
= 4 | |
#define | IPR 5 |
= 5 | |
#define | IPSI 4 |
= 4 | |
#define | IPPSI 5 |
= 5 | |
#define | UNUSED(x) |
#define | TODO(message) |
#define | UNREACHABLE(message) |
#define | return_defer(value) |
#define | DA_INIT_CAP 256 |
#define | da_append(da, item) |
#define | da_insert(da, i, item) |
#define | da_last(da) |
#define | INIT_WRANK |
#define | INIT_WSIZE |
#define | INFO(...) |
#define | WARNING(...) |
#define | ERROR(...) |
#define | PRINT0(...) |
#define | MONOMER_COUNT 6 |
#define | MODULO_BASE 100 |
#define | average_correlation_functions_ext(average, ...) |
Typedefs | |
typedef void(* | dipolePtr) (double *, double[3]) |
typedef void(* | dipoleFree) (void) |
typedef double(* | pesPtr) (double *) |
typedef void(* | dpesPtr) (double *, double *) |
Enumerations | |
enum | MonomerType { ATOM = 0 , LINEAR_MOLECULE = 4 , LINEAR_MOLECULE_REQ_INTEGER = MODULO_BASE + 4 , LINEAR_MOLECULE_REQ_HALFINTEGER = 2*MODULO_BASE + 4 , ROTOR = 6 , ROTOR_REQUANTIZED_ROTATION = MODULO_BASE + 6 } |
MonomerType enum is used to distinguish between systems of different types and stores the phase point size. More... | |
enum | PairState { PAIR_STATE_NONE , PAIR_STATE_FREE_AND_METASTABLE , PAIR_STATE_BOUND , PAIR_STATE_ALL , PAIR_STATE_COUNT } |
PairState is an enumeration representing the possible states of a pair. More... | |
enum | CalculationType { CALCULATION_NONE , CALCULATION_PR_MU , CALCULATION_CORRELATION_SINGLE , CALCULATION_CORRELATION_ARRAY , CALCULATION_PROCESSING , CALCULATION_PHASE_SPACE_M0 , CALCULATION_PHASE_SPACE_M2 , CALCULATION_TYPES_COUNT } |
enum | SamplingType { REJECT , MCMC } |
Functions | |
void | sincos (double, double *, double *) |
int | syncfs (int) |
MoleculeSystem | init_ms (double mu, MonomerType t1, MonomerType t2, double *I1, double *I2, size_t seed) |
Function init_ms prepares the MoleculeSystem based on the specified monomer types, sets up memory buffers and initializes random number generator. | |
MoleculeSystem | init_ms_from_monomers (double mu, Monomer *m1, Monomer *m2, size_t seed) |
Function init_ms_from_monomers prepares the MoleculeSystem based on the provided Monomer structs, sets up memory buffers and initializes random number generator. | |
void | free_ms (MoleculeSystem *ms) |
const char * | display_monomer_type (MonomerType t) |
void | put_qp_into_ms (MoleculeSystem *ms, Array qp) |
void | get_qp_from_ms (MoleculeSystem *ms, Array *qp) |
void | extract_q_and_write_into_ms (MoleculeSystem *ms) |
void | extract_dVdq_and_write_into_monomers (MoleculeSystem *ms) |
int | rhs (realtype t, N_Vector y, N_Vector ydot, void *data) |
rhs function is passed to CVode library to propagate the trajectory. First, the phase-point coordinates are stored into MoleculeSystem struct. A contiguous vector of coordinates is assembled via extract_q_and_write_into_ms. Next, by calling the external function dpes, the derivatives of potential energy are computed and stored into MoleculeSystem::dVdq field. The components of derivative vector are then copied into the field Monomer::dVdq of the corresponding Monomer via the call to extract_dVdq_and_write_into_monomers. The right-hand side of Hamilton's equations with respect to intermolecular degrees of freedom are readily obtained and filled into output ydot, while the derivatives with respect to monomer's coordinates are handled by rhsMonomer function. | |
Array | compute_numerical_derivatives (double(*f)(double *q), double *q, size_t len, size_t order) |
compute_numerical_derivatives | |
Array | compute_numerical_rhs (MoleculeSystem *ms, size_t order) |
compute_numerical_rhs Computes the right-hand side of Hamilton's equations of motion using finite-difference formulas. The returned Array contains derivatives whose order corresponds to variable ordering employed in the MoleculeSystem structure. | |
gsl_matrix * | compute_numerical_jac (void(*transform_angles)(double *qlab, double *qmol), double *qlab, size_t ninput_coordinates, size_t noutput_coordinates, size_t order) |
compute_numerical_jac | |
double | kinetic_energy (MoleculeSystem *ms) |
Function kinetic_energy computes the kinetic energy at the phase-point stored in MoleculeSystem. | |
double | Hamiltonian (MoleculeSystem *ms) |
Function Hamiltonian calls kinetic_energy function to compute kinetic energy, assembles a contiguous vector of coordinates via calling extract_q_and_write_into_ms and passes it to external pes function. The obtained values are added and returned. | |
double | generate_normal (double sigma) |
void | q_generator (MoleculeSystem *ms, CalcParams *params) |
Function q_generator generates ![]() ![]() ![]() ![]() ![]() ![]() | |
void | p_generator (MoleculeSystem *ms, double Temperature) |
Function p_generator samples momenta ![]() ![]() | |
bool | reject (MoleculeSystem *ms, double Temperature, double pesmin) |
reject applies the rejection step to the phase-point that is stored in the MoleculeSystem. It presupposes that the provided phase-point is sampled from ![]() ![]() ![]() | |
void | try_applying_requantization_for_monomer (Monomer *m, size_t step_counter) |
void | j_monomer (Monomer *m, double j[3]) |
j_monomer computes the magnitude of angular momentum of passed-in monomer. Currently, implemented only for linear molecules. | |
double | torque_monomer (Monomer *m) |
torque_monomer computes the magnitude of torque (time-derivative of angular momentum) of passed-in monomer. Requires the derivative of potential Monomer::dVdq to be set within the monomer. Currently, implemented only for linear molecules. | |
double | find_closest_integer (double j) |
find_closest_integer | |
double | find_closest_half_integer (double j) |
find_closest_half_integer | |
void | invert_momenta (MoleculeSystem *ms) |
invert_momenta | |
double | analytic_full_partition_function_by_V (MoleculeSystem *ms, double Temperature) |
double | integrate_composite_simpson (double *x, double *y, size_t len) |
bool | compute_Mn_from_cf_using_classical_detailed_balance (CFnc cf, size_t n, double *result) |
double | compute_Mn_from_sf_using_classical_detailed_balance (SFnc sf, size_t n) |
double | compute_Mn_from_sf_using_quantum_detailed_balance (SFnc sf, size_t n) |
void | calculate_M0 (MoleculeSystem *ms, CalcParams *params, double Temperature, double *m, double *q) |
calculate_M0 Running mean/variance formulas taken from GSL 1.15 https://github.com/ampl/gsl/blob/master/monte/plain.c | |
void | compute_dHdp (MoleculeSystem *ms, gsl_matrix *dHdp) |
compute_dHdp | |
void | calculate_M2 (MoleculeSystem *ms, CalcParams *params, double Temperature, double *m, double *q) |
calculate_M2 Running mean/variance formulas taken from GSL 1.15 https://github.com/ampl/gsl/blob/master/monte/plain.c | |
int | assert_float_is_equal_to (double estimate, double true_value, double abs_tolerance) |
double * | linspace (double start, double end, size_t n) |
double * | arena_linspace (Arena *a, double start, double end, size_t n) |
size_t * | arena_linspace_size_t (Arena *a, size_t start, size_t end, size_t n) |
void | send_histogram_and_reset (gsl_histogram *h) |
void | recv_histogram_and_append (Arena *a, int source, gsl_histogram **h) |
gsl_histogram * | gsl_histogram_extend_right (gsl_histogram *h, size_t add_bins) |
int | write_histogram_ext (FILE *fp, gsl_histogram *h, int count) |
void | sb_reserve (String_Builder *sb, size_t n) |
void | sb_append (String_Builder *sb, const char *line, size_t n) |
void | sb_append_null (String_Builder *sb) |
void | sb_append_cstring (String_Builder *sb, const char *line) |
void | sb_append_format (String_Builder *sb, const char *format,...) |
void | sb_append_seconds_as_datetime_string (String_Builder *sb, int seconds) |
void | sb_reset (String_Builder *sb) |
void | sb_free (String_Builder *sb) |
CFnc | copy_cfnc (CFnc cf) |
SFnc | copy_sfnc (SFnc sf) |
Spectrum | copy_spectrum (Spectrum sp) |
void | normalize_cfnc (CFnc *cf) |
void | free_cfnc (CFnc cf) |
void | free_cfnc_array (CFncArray ca) |
void | free_sfnc (SFnc cf) |
void | free_spectrum (Spectrum sp) |
bool | writetxt (const char *filename, double *x, double *y, size_t len, const char *header) |
bool | write_correlation_function (const char *filename, CFnc cf) |
int | write_correlation_function_ext (FILE *fp, CFnc cf) |
bool | write_spectral_function (const char *filename, SFnc sf) |
int | write_spectral_function_ext (FILE *fp, SFnc sf) |
bool | write_spectrum (const char *filename, Spectrum sp) |
int | write_spectrum_ext (FILE *fp, Spectrum sp) |
bool | parse_number_of_trajectories_from_header (const char *header, double *ntraj) |
bool | parse_temperature_from_header (const char *header, double *Temperature) |
bool | read_correlation_function (const char *filename, String_Builder *sb, CFnc *cf) |
bool | read_spectral_function (const char *filename, String_Builder *sb, SFnc *sf) |
bool | read_spectrum (const char *filename, String_Builder *sb, Spectrum *sp) |
bool | average_correlation_functions (CFnc *average, CFncs cfncs) |
int | average_correlation_functions__impl (CFnc *average, int arg_count,...) |
void | connes_apodization (Array a, double sampling_time) |
double * | dct (double *v, size_t len) |
double * | idct (double *v, size_t len) |
CFnc | dct_sf_to_cf (SFnc sf) |
SFnc | idct_cf_to_sf (CFnc cf) |
Spectrum | compute_alpha (SFnc sf) |
SFnc | desymmetrize_d1 (SFnc sf) |
SFnc | desymmetrize_d2_sf (SFnc sf) |
Spectrum | desymmetrize_d2_sp (Spectrum sp) |
SFnc | desymmetrize_schofield_sf (SFnc sf) |
Spectrum | desymmetrize_schofield_sp (Spectrum sp) |
SFnc | desymmetrize_egelstaff (SFnc sf) |
SFnc | desymmetrize_egelstaff_from_cf (CFnc cf) |
SFnc | desymmetrize_frommhold (SFnc sf) |
SFnc | desymmetrize_frommhold_from_cf (CFnc cf) |
CFnc | egelstaff_time_transform (CFnc cf, bool frommhold_renormalization) |
Spectrum | inv_desymmetrize_schofield (Spectrum sp) |
Spectrum | inv_desymmetrize_d2 (Spectrum sp) |
double * | pad_to_power_of_two (double *v, size_t len, size_t *padded_len) |
double | wingmodel (WingParams *wp, double t) |
double | wingmodel_image (WingParams *wp, double nu) |
int | wingmodel_f (const gsl_vector *x, void *data, gsl_vector *f) |
int | wingmodel_df (const gsl_vector *x, void *data, gsl_matrix *J) |
void | gsl_nonlinear_opt (size_t n, double *x, double *y, WingParams *wing_params) |
WingParams | fit_baseline (CFnc *cf, size_t EXT_RANGE_MIN) |
Variables | |
int | _wrank |
int | _wsize |
bool | _print0_suppress_info |
int | _print0_margin |
MonomerType | MONOMER_TYPES [MONOMER_COUNT] |
const char * | PAIR_STATES [PAIR_STATE_COUNT] |
const char * | CALCULATION_TYPES [CALCULATION_TYPES_COUNT] |
dipolePtr | dipole_1 |
dipolePtr | dipole_2 |
dipoleFree | free_dipole_1 |
dipoleFree | free_dipole_2 |
pesPtr | pes |
dpesPtr | dpes |
WingParams | INIT_WP |
#define _GNU_SOURCE |
#define average_correlation_functions_ext | ( | average, | |
... ) |
#define da_append | ( | da, | |
item ) |
#define DA_INIT_CAP 256 |
#define da_insert | ( | da, | |
i, | |||
item ) |
#define da_last | ( | da | ) |
#define ERROR | ( | ... | ) |
#define INFO | ( | ... | ) |
#define INIT_WRANK |
#define INIT_WSIZE |
#define IPHI 0 |
= 0
#define IPPHI 1 |
= 1
#define IPPSI 5 |
= 5
#define IPR 5 |
= 5
#define IPSI 4 |
= 4
#define IPTHETA 3 |
= 3
#define IR 4 |
= 4
#define ITHETA 2 |
= 2
#define MODULO_BASE 100 |
#define MONOMER_COUNT 6 |
#define MT_GENERATE_CODE_IN_HEADER 0 |
#define OMPI_SKIP_MPICXX |
#define PRINT0 | ( | ... | ) |
#define return_defer | ( | value | ) |
#define TODO | ( | message | ) |
#define UNREACHABLE | ( | message | ) |
#define UNUSED | ( | x | ) |
#define WARNING | ( | ... | ) |
typedef void(* dipoleFree) (void) |
typedef void(* dipolePtr) (double *, double[3]) |
typedef void(* dpesPtr) (double *, double *) |
dpes function accepts an array of Q_SIZE variables. MoleculeSystem describes the order of variables. The pointer to dpes is a global variable in hawaii.c file. The derivative of potential energy function is loaded from the dynamic library (see keyword SO_POTENTIAL in Input File Description). Intermolecular distance measured in Bohrs, while angles are measured in radians. The function is expected to return an array of derivatives of energy in units of Hartree/Bohr of length Q_SIZE.
typedef double(* pesPtr) (double *) |
pes function accepts an array of Q_SIZE variables. MoleculeSystem describes the order of variables. The pointer to pes is a global variable in hawaii.c file. The potential energy function is loaded from the dynamic library (see keyword SO_POTENTIAL in Input File Description). Intermolecular distance is measured in Bohrs, while angles are measured in radians. The function is expected to return energy in units of Hartree.
enum CalculationType |
enum MonomerType |
MonomerType enum is used to distinguish between systems of different types and stores the phase point size.
The size of phase point is calculated as size(phase_point) = MonomerType % MODULO_BASE, where MODULO_BASE is defined to 100 by default.
enum PairState |
PairState is an enumeration representing the possible states of a pair.
enum SamplingType |
double analytic_full_partition_function_by_V | ( | MoleculeSystem * | ms, |
double | Temperature ) |
double * arena_linspace | ( | Arena * | a, |
double | start, | ||
double | end, | ||
size_t | n ) |
size_t * arena_linspace_size_t | ( | Arena * | a, |
size_t | start, | ||
size_t | end, | ||
size_t | n ) |
int assert_float_is_equal_to | ( | double | estimate, |
double | true_value, | ||
double | abs_tolerance ) |
int average_correlation_functions__impl | ( | CFnc * | average, |
int | arg_count, | ||
... ) |
void calculate_M0 | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | Temperature, | ||
double * | m, | ||
double * | q ) |
calculate_M0 Running mean/variance formulas taken from GSL 1.15 https://github.com/ampl/gsl/blob/master/monte/plain.c
void calculate_M2 | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | Temperature, | ||
double * | m, | ||
double * | q ) |
calculate_M2 Running mean/variance formulas taken from GSL 1.15 https://github.com/ampl/gsl/blob/master/monte/plain.c
void compute_dHdp | ( | MoleculeSystem * | ms, |
gsl_matrix * | dHdp ) |
bool compute_Mn_from_cf_using_classical_detailed_balance | ( | CFnc | cf, |
size_t | n, | ||
double * | result ) |
double compute_Mn_from_sf_using_classical_detailed_balance | ( | SFnc | sf, |
size_t | n ) |
double compute_Mn_from_sf_using_quantum_detailed_balance | ( | SFnc | sf, |
size_t | n ) |
Array compute_numerical_derivatives | ( | double(* | f )(double *q), |
double * | q, | ||
size_t | len, | ||
size_t | order ) |
gsl_matrix * compute_numerical_jac | ( | void(* | transform_angles )(double *qlab, double *qmol), |
double * | qlab, | ||
size_t | ninput_coordinates, | ||
size_t | noutput_coordinates, | ||
size_t | order ) |
Array compute_numerical_rhs | ( | MoleculeSystem * | ms, |
size_t | order ) |
compute_numerical_rhs Computes the right-hand side of Hamilton's equations of motion using finite-difference formulas. The returned Array contains derivatives whose order corresponds to variable ordering employed in the MoleculeSystem structure.
ms | Pointer to the MoleculeSystem for which to compute the derivatives. |
order | Order of the finite-difference formula to employ for calculation. |
void connes_apodization | ( | Array | a, |
double | sampling_time ) |
double * dct | ( | double * | v, |
size_t | len ) |
const char * display_monomer_type | ( | MonomerType | t | ) |
void extract_dVdq_and_write_into_monomers | ( | MoleculeSystem * | ms | ) |
void extract_q_and_write_into_ms | ( | MoleculeSystem * | ms | ) |
double find_closest_half_integer | ( | double | j | ) |
requantization to: 0.0, 0.5, 1.5, 2.5, 3.5 ...
double find_closest_integer | ( | double | j | ) |
requantization to: 0.0, 1.0, 2.0, 3.0 ...
WingParams fit_baseline | ( | CFnc * | cf, |
size_t | EXT_RANGE_MIN ) |
void free_cfnc | ( | CFnc | cf | ) |
void free_cfnc_array | ( | CFncArray | ca | ) |
void free_ms | ( | MoleculeSystem * | ms | ) |
void free_sfnc | ( | SFnc | cf | ) |
void free_spectrum | ( | Spectrum | sp | ) |
double generate_normal | ( | double | sigma | ) |
void get_qp_from_ms | ( | MoleculeSystem * | ms, |
Array * | qp ) |
gsl_histogram * gsl_histogram_extend_right | ( | gsl_histogram * | h, |
size_t | add_bins ) |
void gsl_nonlinear_opt | ( | size_t | n, |
double * | x, | ||
double * | y, | ||
WingParams * | wing_params ) |
double Hamiltonian | ( | MoleculeSystem * | ms | ) |
Function Hamiltonian calls kinetic_energy function to compute kinetic energy, assembles a contiguous vector of coordinates via calling extract_q_and_write_into_ms and passes it to external pes function. The obtained values are added and returned.
double * idct | ( | double * | v, |
size_t | len ) |
MoleculeSystem init_ms | ( | double | mu, |
MonomerType | t1, | ||
MonomerType | t2, | ||
double * | II1, | ||
double * | II2, | ||
size_t | seed ) |
Function init_ms prepares the MoleculeSystem based on the specified monomer types, sets up memory buffers and initializes random number generator.
mu | Reduced mass of the molecule pair |
t1 | Type of the first monomer |
t2 | Type of the second monomer |
I1 | Inertia tensor values for the first monomer. If the monomer is MonomerType::ATOM, no values will be read from the pointer, so NULL can be passed. Two and three values are expected for linear molecule (MonomerType::LINEAR_MOLECULE, MonomerType::LINEAR_MOLECULE_REQ_INTEGER and MonomerType::LINEAR_MOLECULE_REQ_HALFINTEGER) and rotor (MonomerType::ROTOR and MonomerType::ROTOR_REQUANTIZED_ROTATION). |
I2 | Inertia tensor values for the second monomer. |
seed | Seed to set up random number generator (Mersenne Twister). A unique seed will be produced if 0 value is passed. |
MoleculeSystem init_ms_from_monomers | ( | double | mu, |
Monomer * | m1, | ||
Monomer * | m2, | ||
size_t | seed ) |
Function init_ms_from_monomers prepares the MoleculeSystem based on the provided Monomer structs, sets up memory buffers and initializes random number generator.
double integrate_composite_simpson | ( | double * | x, |
double * | y, | ||
size_t | len ) |
void invert_momenta | ( | MoleculeSystem * | ms | ) |
void j_monomer | ( | Monomer * | m, |
double | j[3] ) |
j_monomer computes the magnitude of angular momentum of passed-in monomer. Currently, implemented only for linear molecules.
double kinetic_energy | ( | MoleculeSystem * | ms | ) |
Function kinetic_energy computes the kinetic energy at the phase-point stored in MoleculeSystem.
The kinetic energy in the laboratory frame of reference is calculated in parts. The part corresponding to intermolecular degrees of freedom:
to linear molecule:
and to rotor:
double * linspace | ( | double | start, |
double | end, | ||
size_t | n ) |
void normalize_cfnc | ( | CFnc * | cf | ) |
void p_generator | ( | MoleculeSystem * | ms, |
double | Temperature ) |
Function p_generator samples momenta
double * pad_to_power_of_two | ( | double * | v, |
size_t | len, | ||
size_t * | padded_len ) |
bool parse_number_of_trajectories_from_header | ( | const char * | header, |
double * | ntraj ) |
bool parse_temperature_from_header | ( | const char * | header, |
double * | Temperature ) |
void put_qp_into_ms | ( | MoleculeSystem * | ms, |
Array | qp ) |
void q_generator | ( | MoleculeSystem * | ms, |
CalcParams * | params ) |
Function q_generator generates
bool read_correlation_function | ( | const char * | filename, |
String_Builder * | sb, | ||
CFnc * | cf ) |
bool read_spectral_function | ( | const char * | filename, |
String_Builder * | sb, | ||
SFnc * | sf ) |
bool read_spectrum | ( | const char * | filename, |
String_Builder * | sb, | ||
Spectrum * | sp ) |
void recv_histogram_and_append | ( | Arena * | a, |
int | source, | ||
gsl_histogram ** | h ) |
bool reject | ( | MoleculeSystem * | ms, |
double | Temperature, | ||
double | pesmin ) |
reject applies the rejection step to the phase-point that is stored in the MoleculeSystem. It presupposes that the provided phase-point is sampled from
int rhs | ( | realtype | t, |
N_Vector | y, | ||
N_Vector | ydot, | ||
void * | data ) |
rhs function is passed to CVode library to propagate the trajectory. First, the phase-point coordinates are stored into MoleculeSystem struct. A contiguous vector of coordinates is assembled via extract_q_and_write_into_ms. Next, by calling the external function dpes, the derivatives of potential energy are computed and stored into MoleculeSystem::dVdq field. The components of derivative vector are then copied into the field Monomer::dVdq of the corresponding Monomer via the call to extract_dVdq_and_write_into_monomers. The right-hand side of Hamilton's equations with respect to intermolecular degrees of freedom are readily obtained and filled into output ydot, while the derivatives with respect to monomer's coordinates are handled by rhsMonomer function.
void sb_append | ( | String_Builder * | sb, |
const char * | line, | ||
size_t | n ) |
void sb_append_cstring | ( | String_Builder * | sb, |
const char * | line ) |
void sb_append_format | ( | String_Builder * | sb, |
const char * | format, | ||
... ) |
void sb_append_null | ( | String_Builder * | sb | ) |
void sb_append_seconds_as_datetime_string | ( | String_Builder * | sb, |
int | seconds ) |
void sb_free | ( | String_Builder * | sb | ) |
void sb_reserve | ( | String_Builder * | sb, |
size_t | n ) |
void sb_reset | ( | String_Builder * | sb | ) |
void send_histogram_and_reset | ( | gsl_histogram * | h | ) |
void sincos | ( | double | , |
double * | , | ||
double * | ) |
int syncfs | ( | int | ) |
double torque_monomer | ( | Monomer * | m | ) |
torque_monomer computes the magnitude of torque (time-derivative of angular momentum) of passed-in monomer. Requires the derivative of potential Monomer::dVdq to be set within the monomer. Currently, implemented only for linear molecules.
void try_applying_requantization_for_monomer | ( | Monomer * | m, |
size_t | step_counter ) |
double wingmodel | ( | WingParams * | wp, |
double | t ) |
int wingmodel_df | ( | const gsl_vector * | x, |
void * | data, | ||
gsl_matrix * | J ) |
int wingmodel_f | ( | const gsl_vector * | x, |
void * | data, | ||
gsl_vector * | f ) |
double wingmodel_image | ( | WingParams * | wp, |
double | nu ) |
bool write_correlation_function | ( | const char * | filename, |
CFnc | cf ) |
int write_correlation_function_ext | ( | FILE * | fp, |
CFnc | cf ) |
int write_histogram_ext | ( | FILE * | fp, |
gsl_histogram * | h, | ||
int | count ) |
bool write_spectral_function | ( | const char * | filename, |
SFnc | sf ) |
int write_spectral_function_ext | ( | FILE * | fp, |
SFnc | sf ) |
bool write_spectrum | ( | const char * | filename, |
Spectrum | sp ) |
int write_spectrum_ext | ( | FILE * | fp, |
Spectrum | sp ) |
bool writetxt | ( | const char * | filename, |
double * | x, | ||
double * | y, | ||
size_t | len, | ||
const char * | header ) |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |
|
extern |