| 
    IMP Reference Guide
    2.23.0
    
   The Integrative Modeling Platform 
   | 
 
#include <IMP/saxs/Profile.h>
 Inheritance diagram for IMP::saxs::Profile:Basic profile class, can be initialized from the input file (experimental or theoretical) or computed from a set of Model Particles (theoretical)
Public Member Functions | |
| Profile (const std::string &file_name, bool fit_file=false, double max_q=0.0, int units=1) | |
| init from file  More... | |
| Profile (double qmin=0.0, double qmax=0.5, double delta=0.005) | |
| init for theoretical profile  More... | |
| void | add (const Profile *other_profile, double weight=1.0) | 
| add another profile - useful for rigid bodies  More... | |
| void | add (const Vector< Profile * > &profiles, const Vector< double > &weights=Vector< double >()) | 
| add other profiles - useful for weighted ensembles  More... | |
| void | add_errors () | 
| add simulated error  More... | |
| void | add_noise (double percentage=0.03) | 
| add simulated noise  More... | |
| void | add_partial_profiles (const Profile *other_profile, double weight=1.0) | 
| add partial profiles  More... | |
| void | add_partial_profiles (const Vector< Profile * > &profiles, const Vector< double > &weights=Vector< double >()) | 
| add other partial profiles  More... | |
| void | background_adjust (double start_q) | 
| background adjustment option  More... | |
| double | calculate_I0 (const Particles &particles, FormFactorType ff_type=HEAVY_ATOMS) | 
| calculate Intensity at zero (= squared number of electrons)  More... | |
| void | calculate_profile (const Particles &particles, FormFactorType ff_type=HEAVY_ATOMS, bool reciprocal=false) | 
| computes theoretical profile  More... | |
| void | calculate_profile (const Particles &particles1, const Particles &particles2, FormFactorType ff_type=HEAVY_ATOMS) | 
| void | calculate_profile_constant_form_factor (const Particles &particles, double form_factor=1.0) | 
| calculate profile for any type of Particles that have coordinates  More... | |
| void | calculate_profile_partial (const Particles &particles, const Vector< double > &surface=Vector< double >(), FormFactorType ff_type=HEAVY_ATOMS) | 
| compute profile for fitting with hydration layer and excluded volume  More... | |
| void | calculate_profile_partial (const Particles &particles1, const Particles &particles2, const Vector< double > &surface1=Vector< double >(), const Vector< double > &surface2=Vector< double >(), FormFactorType ff_type=HEAVY_ATOMS) | 
| compute profile for fitting with hydration layer and excluded volume  More... | |
| void | calculate_profile_reciprocal_partial (const Particles &particles, const Vector< double > &surface=Vector< double >(), FormFactorType ff_type=HEAVY_ATOMS) | 
| void | calculate_profile_symmetric (const Particles &particles, unsigned int n, FormFactorType ff_type=HEAVY_ATOMS) | 
| void | copy_errors (const Profile *exp_profile) | 
| void | distribution_2_profile (const RadialDistributionFunction &r_dist) | 
| convert to reciprocal space I(q) = Sum(P(r)*sin(qr)/qr)  More... | |
| void | downsample (Profile *downsampled_profile, unsigned int point_number) const | 
| downsample the profile to a given number of points  More... | |
| double | get_average_radius () const | 
| double | get_delta_q () const | 
| return sampling resolution  More... | |
| double | get_error (unsigned int i) const | 
| const Eigen::VectorXd & | get_errors () const | 
| unsigned int | get_id () const | 
| const Eigen::VectorXd & | get_intensities () const | 
| double | get_intensity (unsigned int i) const | 
| double | get_max_q () const | 
| return maximal sampling point  More... | |
| double | get_min_q () const | 
| return minimal sampling point  More... | |
| std::string | get_name () const | 
| double | get_q (unsigned int i) const | 
| const Eigen::VectorXd & | get_qs () const | 
| virtual std::string | get_type_name () const override | 
| virtual ::IMP::VersionInfo | get_version_info () const override | 
| Get information about the module and version of the object.  More... | |
| double | get_weight (unsigned int i) const | 
| bool | is_partial_profile () const | 
| bool | is_uniform_sampling () const | 
| checks the sampling of experimental profile  More... | |
| double | mean_intensity () const | 
| calculate mean intensity  More... | |
| void | offset (double c) | 
| offset profile by c, I(q) = I(q) - c  More... | |
| void | profile_2_distribution (RadialDistributionFunction &rd, double max_distance) const | 
| convert to real space P(r) function P(r) = 1/2PI^2 Sum(I(q)*qr*sin(qr))  More... | |
| double | radius_of_gyration (double end_q_rg=1.3) const | 
| compute radius of gyration with Guinier approximation  More... | |
| void | read_partial_profiles (const std::string &file_name) | 
| read a partial profile from file (7 columns)  More... | |
| void | read_SAXS_file (const std::string &file_name, bool fit_file=false, double max_q=0.0, int units=1) | 
| reads SAXS profile from file  More... | |
| void | resample (const Profile *exp_profile, Profile *resampled_profile) const | 
| return a profile that is sampled on the q values of the exp_profile  More... | |
| void | scale (double c) | 
| scale  More... | |
| void | set_average_radius (double r) | 
| void | set_average_volume (double v) | 
| void | set_beam_profile (std::string beam_profile_file) | 
| void | set_errors (const Eigen::VectorXd &e) | 
| void | set_ff_table (FormFactorTable *ff_table) | 
| required for reciprocal space calculation  More... | |
| void | set_id (unsigned int id) | 
| void | set_intensities (const Eigen::VectorXd &i) | 
| void | set_intensity (unsigned int i, double iq) | 
| void | set_name (std::string name) | 
| void | set_qs (const Eigen::VectorXd &q) | 
| unsigned int | size () const | 
| return number of entries in SAXS profile  More... | |
| void | sum_partial_profiles (double c1, double c2, bool check_cashed=true) | 
| computes full profile for given fitting parameters  More... | |
| void | write_partial_profiles (const std::string &file_name) const | 
| write a partial profile to file  More... | |
| void | write_SAXS_file (const std::string &file_name, double max_q=0.0) const | 
| print to file  More... | |
  Public Member Functions inherited from IMP::Object | |
| virtual void | clear_caches () | 
| CheckLevel | get_check_level () const | 
| LogLevel | get_log_level () const | 
| void | set_check_level (CheckLevel l) | 
| void | set_log_level (LogLevel l) | 
| Set the logging level used in this object.  More... | |
| void | set_was_used (bool tf) const | 
| void | show (std::ostream &out=std::cout) const | 
| const std::string & | get_name () const | 
| void | set_name (std::string name) | 
Static Public Attributes | |
| static const double | modulation_function_parameter_ | 
Protected Member Functions | |
| void | init (unsigned int size=0, unsigned int partial_profiles_size=0) | 
  Protected Member Functions inherited from IMP::Object | |
| Object (std::string name) | |
| Construct an object with the given name.  More... | |
| virtual void | do_destroy () | 
Protected Attributes | |
| double | average_radius_ | 
| double | average_volume_ | 
| Pointer< Profile > | beam_profile_ | 
| double | c1_ | 
| double | c2_ | 
| double | delta_q_ | 
| Eigen::VectorXd | error_ | 
| bool | experimental_ | 
| FormFactorTable * | ff_table_ | 
| unsigned int | id_ | 
| Eigen::VectorXd | intensity_ | 
| double | max_q_ | 
| double | min_q_ | 
| std::string | name_ | 
| std::vector< Eigen::VectorXd > | partial_profiles_ | 
| Eigen::VectorXd | q_ | 
| std::map< double, unsigned int > | q_mapping_ | 
| IMP::saxs::Profile::Profile | ( | const std::string & | file_name, | 
| bool | fit_file = false,  | 
        ||
| double | max_q = 0.0,  | 
        ||
| int | units = 1  | 
        ||
| ) | 
init from file
| [in] | file_name | profile file name | 
| [in] | fit_file | if true, intensities are read from column 3 | 
| [in] | max_q | read till maximal q value = max_q, or all if max_q<=0 | 
| [in] | units | gets 1, 2, or 3 for unknown q units, 1/A, or 1/nm | 
| IMP::saxs::Profile::Profile | ( | double | qmin = 0.0,  | 
        
| double | qmax = 0.5,  | 
        ||
| double | delta = 0.005  | 
        ||
| ) | 
init for theoretical profile
| void IMP::saxs::Profile::add | ( | const Profile * | other_profile, | 
| double | weight = 1.0  | 
        ||
| ) | 
add another profile - useful for rigid bodies
| void IMP::saxs::Profile::add | ( | const Vector< Profile * > & | profiles, | 
| const Vector< double > & | weights = Vector< double >()  | 
        ||
| ) | 
add other profiles - useful for weighted ensembles
| void IMP::saxs::Profile::add_errors | ( | ) | 
add simulated error
| void IMP::saxs::Profile::add_noise | ( | double | percentage = 0.03 | ) | 
add simulated noise
| void IMP::saxs::Profile::add_partial_profiles | ( | const Profile * | other_profile, | 
| double | weight = 1.0  | 
        ||
| ) | 
add partial profiles
| void IMP::saxs::Profile::add_partial_profiles | ( | const Vector< Profile * > & | profiles, | 
| const Vector< double > & | weights = Vector< double >()  | 
        ||
| ) | 
add other partial profiles
| void IMP::saxs::Profile::background_adjust | ( | double | start_q | ) | 
background adjustment option
| double IMP::saxs::Profile::calculate_I0 | ( | const Particles & | particles, | 
| FormFactorType | ff_type = HEAVY_ATOMS  | 
        ||
| ) | 
calculate Intensity at zero (= squared number of electrons)
| void IMP::saxs::Profile::calculate_profile | ( | const Particles & | particles, | 
| FormFactorType | ff_type = HEAVY_ATOMS,  | 
        ||
| bool | reciprocal = false  | 
        ||
| ) | 
| void IMP::saxs::Profile::calculate_profile | ( | const Particles & | particles1, | 
| const Particles & | particles2, | ||
| FormFactorType | ff_type = HEAVY_ATOMS  | 
        ||
| ) | 
| void IMP::saxs::Profile::calculate_profile_constant_form_factor | ( | const Particles & | particles, | 
| double | form_factor = 1.0  | 
        ||
| ) | 
calculate profile for any type of Particles that have coordinates
| void IMP::saxs::Profile::calculate_profile_partial | ( | const Particles & | particles, | 
| const Vector< double > & | surface = Vector< double >(),  | 
        ||
| FormFactorType | ff_type = HEAVY_ATOMS  | 
        ||
| ) | 
compute profile for fitting with hydration layer and excluded volume
A partial profile is a pre-computed profile, where intensities are split into 6 parts that can be summed up into a regular profile given a pair of c1/c2 values by sum_partial_profiles function. see FoXS paper for details.
| void IMP::saxs::Profile::calculate_profile_partial | ( | const Particles & | particles1, | 
| const Particles & | particles2, | ||
| const Vector< double > & | surface1 = Vector< double >(),  | 
        ||
| const Vector< double > & | surface2 = Vector< double >(),  | 
        ||
| FormFactorType | ff_type = HEAVY_ATOMS  | 
        ||
| ) | 
compute profile for fitting with hydration layer and excluded volume
| void IMP::saxs::Profile::distribution_2_profile | ( | const RadialDistributionFunction & | r_dist | ) | 
convert to reciprocal space I(q) = Sum(P(r)*sin(qr)/qr)
| void IMP::saxs::Profile::downsample | ( | Profile * | downsampled_profile, | 
| unsigned int | point_number | ||
| ) | const | 
downsample the profile to a given number of points
| double IMP::saxs::Profile::get_delta_q | ( | ) | const | 
| double IMP::saxs::Profile::get_max_q | ( | ) | const | 
| double IMP::saxs::Profile::get_min_q | ( | ) | const | 
      
  | 
  overridevirtual | 
Get information about the module and version of the object.
Reimplemented from IMP::Object.
| bool IMP::saxs::Profile::is_uniform_sampling | ( | ) | const | 
checks the sampling of experimental profile
| double IMP::saxs::Profile::mean_intensity | ( | ) | const | 
calculate mean intensity
| void IMP::saxs::Profile::offset | ( | double | c | ) | 
offset profile by c, I(q) = I(q) - c
| void IMP::saxs::Profile::profile_2_distribution | ( | RadialDistributionFunction & | rd, | 
| double | max_distance | ||
| ) | const | 
convert to real space P(r) function P(r) = 1/2PI^2 Sum(I(q)*qr*sin(qr))
| double IMP::saxs::Profile::radius_of_gyration | ( | double | end_q_rg = 1.3 | ) | const | 
compute radius of gyration with Guinier approximation
ln[I(q)]=ln[I(0)] - (q^2*rg^2)/3
| [in] | end_q_rg | determines the range of profile used for approximation: i.e. q*rg < end_q_rg. Use 1.3 for globular proteins, 0.8 for elongated | 
| void IMP::saxs::Profile::read_partial_profiles | ( | const std::string & | file_name | ) | 
read a partial profile from file (7 columns)
| void IMP::saxs::Profile::read_SAXS_file | ( | const std::string & | file_name, | 
| bool | fit_file = false,  | 
        ||
| double | max_q = 0.0,  | 
        ||
| int | units = 1  | 
        ||
| ) | 
reads SAXS profile from file
| [in] | file_name | profile file name | 
| [in] | fit_file | if true, intensities are read from column 3 | 
| [in] | max_q | read till maximal q value = max_q, or all if max_q<=0 | 
| [in] | units | gets 1, 2, or 3 for unknown q units, 1/A, or 1/nm | 
| void IMP::saxs::Profile::resample | ( | const Profile * | exp_profile, | 
| Profile * | resampled_profile | ||
| ) | const | 
return a profile that is sampled on the q values of the exp_profile
| void IMP::saxs::Profile::scale | ( | double | c | ) | 
scale
| void IMP::saxs::Profile::set_ff_table | ( | FormFactorTable * | ff_table | ) | 
| unsigned int IMP::saxs::Profile::size | ( | ) | const | 
| void IMP::saxs::Profile::sum_partial_profiles | ( | double | c1, | 
| double | c2, | ||
| bool | check_cashed = true  | 
        ||
| ) | 
computes full profile for given fitting parameters
| void IMP::saxs::Profile::write_partial_profiles | ( | const std::string & | file_name | ) | const | 
write a partial profile to file
| void IMP::saxs::Profile::write_SAXS_file | ( | const std::string & | file_name, | 
| double | max_q = 0.0  | 
        ||
| ) | const | 
print to file
| [in] | file_name | output file name | 
| [in] | max_q | output till maximal q value = max_q, or all if max_q<=0 |