Hawaii Hybrid
|
#include "hawaii.h"
#include "arena.h"
#include "trajectory.h"
#include "hep_hawaii.h"
#include "angles_handler.hpp"
Macros | |
#define | ARENA_REGION_DEFAULT_CAPACITY (8*1024) |
#define | HISTOGRAM_MAX_TPS 50 |
#define | R_HISTOGRAM_BINS 100 |
#define | R_HISTOGRAM_MAX 10000000.0 |
#define | DEFAULT_JINI_HISTOGRAM_FILENAME1 "jini.dat.1" |
#define | DEFAULT_JINI_HISTOGRAM_FILENAME2 "jini.dat.2" |
#define | DEFAULT_JINI_HISTOGRAM_BINS 352 |
#define | DEFAULT_JINI_HISTOGRAM_MAX 35.0 |
#define | DEFAULT_JFIN_HISTOGRAM_FILENAME1 "jfin.dat.1" |
#define | DEFAULT_JFIN_HISTOGRAM_FILENAME2 "jfin.dat.2" |
#define | DEFAULT_JFIN_HISTOGRAM_BINS 352 |
#define | DEFAULT_JFIN_HISTOGRAM_MAX 35.0 |
#define | DEFAULT_NSWITCH_HISTOGRAM_FILENAME1 "nswitch.dat.1" |
#define | DEFAULT_NSWITCH_HISTOGRAM_FILENAME2 "nswitch.dat.2" |
#define | DEFAULT_NSWITCH_HISTOGRAM_BINS 20 |
#define | DEFAULT_NSWITCH_HISTOGRAM_MAX 20.0 |
#define | ARENA_IMPLEMENTATION |
#define | nparams 3 |
#define | REAL(z, i) |
#define | IMAG(z, i) |
Functions | |
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. | |
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. | |
void | free_stored_histogram (StoredHistogram *sh) |
void | free_ms (MoleculeSystem *ms) |
const char * | display_monomer_type (MonomerType t) |
double | find_closest_integer (double j) |
find_closest_integer | |
double | find_closest_half_integer (double j) |
find_closest_half_integer | |
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. | |
void | rhsMonomer (Monomer *m, double *deriv) |
rhsMonomer dispatches between computing the right-hand side of Hamilton's equations of motion depending on the type of monomer. In addition to differentiating the kinetic energy, the derivatives of potential energy, which are taken from Monomer::dVdq, are also added to compute the right-hand side. rhsMonomer does not account for requantization which should happen after this function has been called. | |
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. | |
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 | |
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. | |
double | kinetic_energy (MoleculeSystem *ms) |
Function kinetic_energy computes the kinetic energy at the phase-point stored in MoleculeSystem. | |
void | put_qp_into_ms (MoleculeSystem *ms, Array qp) |
void | get_qp_from_ms (MoleculeSystem *ms, Array *qp) |
void | extract_q (double *qp, double *q, size_t QP_SIZE) |
void | extract_q_and_write_into_ms (MoleculeSystem *ms) |
void | invert_momenta (MoleculeSystem *ms) |
invert_momenta | |
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 | 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 | |
void | mpi_calculate_M0 (MoleculeSystem *ms, CalcParams *params, double Temperature, double *m, double *q) |
mpi_calculate_M0 | |
void | mpi_calculate_M2 (MoleculeSystem *ms, CalcParams *params, double Temperature, double *m, double *q) |
mpi_calculate_M2 | |
void | track_turning_points (Tracker *tr, double R) |
int | correlation_eval_zimmerman_trick (MoleculeSystem *ms, Trajectory *traj, CalcParams *params, double *crln, size_t *tps) |
correlation_eval_zimmerman_trick | |
void | try_applying_requantization_for_monomer (Monomer *m, size_t step_counter) |
int | correlation_eval (MoleculeSystem *ms, Trajectory *traj, CalcParams *params, double *crln, size_t *tps) |
correlation_eval | |
void | gsl_fft_square (double *farr, size_t N) |
CFncArray | calculate_correlation_array_and_save (MoleculeSystem *ms, CalcParams *params, double base_temperature) |
calculate_correlation_array_and_save | |
void | setup_nswitch_histogram_for_monomer (Monomer *m) |
void | setup_jini_histogram_for_monomer (Monomer *m) |
void | setup_jfin_histogram_for_monomer (Monomer *m) |
void | send_histogram_and_reset (gsl_histogram *h) |
void | normalize_cfnc (CFnc *cf) |
CFnc | calculate_correlation_and_save (MoleculeSystem *ms, CalcParams *params, double Temperature) |
calculate_correlation_and_save | |
void | recv_histogram_and_append (Arena *a, int source, gsl_histogram **h) |
SFnc | calculate_spectral_function_using_prmu_representation_and_save (MoleculeSystem *ms, CalcParams *params, double Temperature) |
calculate_spectral_function_using_prmu_representation_and_save | |
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) |
bool | write_spectrum (const char *filename, Spectrum sp) |
int | write_spectrum_ext (FILE *fp, Spectrum sp) |
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) |
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_reset (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 s) |
void | sb_free (String_Builder *sb) |
bool | read_correlation_function (const char *filename, String_Builder *sb, CFnc *cf) |
bool | average_correlation_functions (CFnc *average, CFncs cfncs) |
int | average_correlation_functions__impl (CFnc *average, int arg_count,...) |
bool | read_spectral_function (const char *filename, String_Builder *sb, SFnc *sf) |
bool | parse_temperature_from_header (const char *header, double *Temperature) |
bool | parse_number_of_trajectories_from_header (const char *header, double *ntraj) |
bool | read_spectrum (const char *filename, String_Builder *sb, Spectrum *sp) |
bool | writetxt (const char *filename, double *x, double *y, size_t len, const char *header) |
double | analytic_full_partition_function_by_V (MoleculeSystem *ms, double Temperature) |
gsl_histogram * | gsl_histogram_extend_right (gsl_histogram *h, size_t add_bins) |
void | free_cfnc (CFnc cf) |
void | free_cfnc_array (CFncArray ca) |
void | free_sfnc (SFnc sf) |
void | free_spectrum (Spectrum sp) |
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_multifit_callback (const size_t iter, void *params, const gsl_multifit_nlinear_workspace *w) |
void | gsl_nonlinear_opt (size_t n, double *x, double *y, WingParams *wing_params) |
WingParams | fit_baseline (CFnc *cf, size_t EXT_RANGE_MIN) |
void | connes_apodization (Array a, double sampling_time) |
double * | dct (double *v, size_t len) |
double * | idct (double *v, size_t len) |
double * | pad_to_power_of_two (double *v, size_t len, size_t *padded_len) |
SFnc | idct_cf_to_sf (CFnc cf) |
CFnc | dct_sf_to_cf (SFnc sf) |
SFnc | desymmetrize_schofield_sf (SFnc sf) |
Spectrum | desymmetrize_schofield_sp (Spectrum sp) |
Spectrum | inv_desymmetrize_schofield (Spectrum sp) |
SFnc | desymmetrize_d2_sf (SFnc sf) |
Spectrum | desymmetrize_d2_sp (Spectrum sp) |
Spectrum | inv_desymmetrize_d2 (Spectrum sp) |
SFnc | desymmetrize_d1 (SFnc sf) |
CFnc | egelstaff_time_transform (CFnc cf, bool frommhold_renormalization) |
SFnc | desymmetrize_frommhold (SFnc sf) |
SFnc | desymmetrize_egelstaff_from_cf (CFnc cf) |
SFnc | desymmetrize_frommhold_from_cf (CFnc cf) |
SFnc | desymmetrize_egelstaff (SFnc sf) |
Spectrum | compute_alpha (SFnc sf) |
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) |
CFnc | copy_cfnc (CFnc cf) |
SFnc | copy_sfnc (SFnc sf) |
Spectrum | copy_spectrum (Spectrum sp) |
Variables | |
dipolePtr | dipole_1 = NULL |
dipoleFree | free_dipole_1 = NULL |
dipolePtr | dipole_2 = NULL |
dipoleFree | free_dipole_2 = NULL |
pesPtr | pes = NULL |
dpesPtr | dpes = NULL |
int | _wrank = 0 |
int | _wsize = 1 |
bool | _print0_suppress_info = false |
int | _print0_margin = 0 |
const char * | PAIR_STATES [PAIR_STATE_COUNT] |
const char * | CALCULATION_TYPES [CALCULATION_TYPES_COUNT] |
MonomerType | MONOMER_TYPES [MONOMER_COUNT] |
WingParams | INIT_WP |
#define ARENA_IMPLEMENTATION |
#define ARENA_REGION_DEFAULT_CAPACITY (8*1024) |
#define DEFAULT_JFIN_HISTOGRAM_BINS 352 |
#define DEFAULT_JFIN_HISTOGRAM_FILENAME1 "jfin.dat.1" |
#define DEFAULT_JFIN_HISTOGRAM_FILENAME2 "jfin.dat.2" |
#define DEFAULT_JFIN_HISTOGRAM_MAX 35.0 |
#define DEFAULT_JINI_HISTOGRAM_BINS 352 |
#define DEFAULT_JINI_HISTOGRAM_FILENAME1 "jini.dat.1" |
#define DEFAULT_JINI_HISTOGRAM_FILENAME2 "jini.dat.2" |
#define DEFAULT_JINI_HISTOGRAM_MAX 35.0 |
#define DEFAULT_NSWITCH_HISTOGRAM_BINS 20 |
#define DEFAULT_NSWITCH_HISTOGRAM_FILENAME1 "nswitch.dat.1" |
#define DEFAULT_NSWITCH_HISTOGRAM_FILENAME2 "nswitch.dat.2" |
#define DEFAULT_NSWITCH_HISTOGRAM_MAX 20.0 |
#define HISTOGRAM_MAX_TPS 50 |
#define IMAG | ( | z, | |
i ) |
#define nparams 3 |
#define R_HISTOGRAM_BINS 100 |
#define R_HISTOGRAM_MAX 10000000.0 |
#define REAL | ( | z, | |
i ) |
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, | ||
... ) |
CFnc calculate_correlation_and_save | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | Temperature ) |
calculate_correlation_and_save
TODO: check 'sampler_Rmin': we want to catch the situation when it's too low for given Temperature
CFncArray calculate_correlation_array_and_save | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | base_temperature ) |
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
SFnc calculate_spectral_function_using_prmu_representation_and_save | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | Temperature ) |
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 ) |
int correlation_eval | ( | MoleculeSystem * | ms, |
Trajectory * | traj, | ||
CalcParams * | params, | ||
double * | crln, | ||
size_t * | tps ) |
int correlation_eval_zimmerman_trick | ( | MoleculeSystem * | ms, |
Trajectory * | traj, | ||
CalcParams * | params, | ||
double * | crln, | ||
size_t * | tps ) |
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 | ( | double * | qp, |
double * | q, | ||
size_t | QP_SIZE ) |
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 | sf | ) |
void free_spectrum | ( | Spectrum | sp | ) |
void free_stored_histogram | ( | StoredHistogram * | sh | ) |
double generate_normal | ( | double | sigma | ) |
void get_qp_from_ms | ( | MoleculeSystem * | ms, |
Array * | qp ) |
void gsl_fft_square | ( | double * | farr, |
size_t | N ) |
gsl_histogram * gsl_histogram_extend_right | ( | gsl_histogram * | h, |
size_t | add_bins ) |
void gsl_multifit_callback | ( | const size_t | iter, |
void * | params, | ||
const gsl_multifit_nlinear_workspace * | w ) |
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 mpi_calculate_M0 | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | Temperature, | ||
double * | m, | ||
double * | q ) |
void mpi_calculate_M2 | ( | MoleculeSystem * | ms, |
CalcParams * | params, | ||
double | Temperature, | ||
double * | m, | ||
double * | q ) |
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 rhsMonomer | ( | Monomer * | m, |
double * | deriv ) |
rhsMonomer dispatches between computing the right-hand side of Hamilton's equations of motion depending on the type of monomer. In addition to differentiating the kinetic energy, the derivatives of potential energy, which are taken from Monomer::dVdq, are also added to compute the right-hand side.
rhsMonomer does not account for requantization which should happen after this function has been called.
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 | s ) |
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 setup_jfin_histogram_for_monomer | ( | Monomer * | m | ) |
void setup_jini_histogram_for_monomer | ( | Monomer * | m | ) |
void setup_nswitch_histogram_for_monomer | ( | Monomer * | m | ) |
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 track_turning_points | ( | Tracker * | tr, |
double | R ) |
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 ) |
int _print0_margin = 0 |
bool _print0_suppress_info = false |
int _wrank = 0 |
int _wsize = 1 |
const char* CALCULATION_TYPES[CALCULATION_TYPES_COUNT] |
dipolePtr dipole_1 = NULL |
dipolePtr dipole_2 = NULL |
dpesPtr dpes = NULL |
dipoleFree free_dipole_1 = NULL |
dipoleFree free_dipole_2 = NULL |
WingParams INIT_WP |
MonomerType MONOMER_TYPES[MONOMER_COUNT] |
const char* PAIR_STATES[PAIR_STATE_COUNT] |
pesPtr pes = NULL |