Bsoft 2.1.4
Bernard's software package
model_util.h File Reference

Library routines used for model processing. More...

#include "rwmodel.h"
#include "rwmolecule.h"
#include "Bstring.h"

Functions

long models_process (Bmodel *model, long(modfunc)(Bmodel *))
 Process a list of models using the specified function. More...
 
long models_process (Bmodel *model, long i, long(modfunc)(Bmodel *, long))
 
long models_process (Bmodel *model, double d, long(modfunc)(Bmodel *, double))
 
long models_process (Bmodel *model, Bstring &str, long(modfunc)(Bmodel *, Bstring &str))
 
long models_process (Bmodel *model, string &str, long(modfunc)(Bmodel *, string &str))
 
long model_list (Bmodel *model)
 Lists models in table form. More...
 
long model_list_comp (Bmodel *model)
 Lists models with component counts in table form. More...
 
long model_replace_components (Bmodel *model, Bmodel *modref)
 Replaces all components in a model with reference components. More...
 
long model_merge (Bmodel *model)
 Merges components from all models into one. More...
 
long model_number_ids (Bmodel *model)
 Add a number to each model id. More...
 
long model_rename (Bmodel *model, char first_name)
 Rename models with alphabetical letters. More...
 
long model_rename_components (Bmodel *model)
 Rename components. More...
 
double model_mass (Bmodel *model)
 Calculates the mass of a model from component masses. More...
 
long model_mass_all (Bmodel *model)
 Calculates the masses of all the models in the list. More...
 
Vector3< double > model_center_of_mass (Bmodel *model)
 Calculates the center-of-mass of a model. More...
 
Vector3< double > models_center_of_coordinates (Bmodel *model)
 Calculates the center-of-mass of a list of models. More...
 
Vector3< double > model_geometric_median (Bmodel *model)
 Calculates the geometric median of model components. More...
 
double model_gyration_radius (Bmodel *model)
 Calculates the radius of gyration for a model. More...
 
Vector3< double > model_principal_axes (Bmodel *model, Vector3< double > *eigenvec)
 Calculates the principal axes of a model. More...
 
Vector3< double > model_principal_axes (Bmodel *model, Matrix &eigenvec)
 
long model_principal_axes (Bmodel *model)
 
long model_radial_distribution (Bmodel *model, double interval)
 Calculates the radial distribution function of a model. More...
 
Bmolgroupmodel_assemble (Bmodel *model, Bstring &paramfile, int separate)
 Concatenates selected molecules into one group. More...
 
Bmodelmodel_generate_com (Bmolgroup *molgroup)
 Calculates the centers-of-mass of molecule group components and generates a new model. More...
 
long model_update_centers_of_mass (Bmodel *model, Bmolgroup *molgroup)
 Updates the centers-of-mass of molecule group components. More...
 
long model_average_components (Bmodel *model, int number)
 Averages sequential components. More...
 
Bcomponent ** component_get_array (Bmodel *model, long &ncomp)
 Generates an array of pointers to component structures. More...
 
Vector3< double > component_plane (vector< Bcomponent * > &comparray, double &offset)
 Calculates a plane through an array of components. More...
 
vector< Vector3< double > > models_calculate_bounds (Bmodel *model)
 Calculates the bounds of a list of models. More...
 
vector< vector< Bcomponent * > > model_component_grid (Bmodel *model, Vector3< long > &size, Vector3< double > &origin, Vector3< double > &sampling)
 Generates lists of atoms based on a grid. More...
 
vector< Bcomponent * > models_get_component_array (Bmodel *model)
 Generates an array of pointers to model components. More...
 
vector< vector< Bcomponent * > > model_split_into_slices (Bmodel *model, double bottom, double top, double thickness)
 Splits the components of models into slices. More...
 

Detailed Description

Library routines used for model processing.

Author
Bernard Heymann
Date
Created: 20060908
Modified: 20230120

Function Documentation

◆ component_get_array()

Bcomponent ** component_get_array ( Bmodel model,
long &  n 
)

Generates an array of pointers to component structures.

Parameters
*modelmodel structure.
&npointer to number of comps found.
Returns
Bcomponent** array of pointers to components.

◆ component_plane()

Vector3< double > component_plane ( vector< Bcomponent * > &  comparray,
double &  offset 
)

Calculates a plane through an array of components.

Parameters
**comparrayarray of components.
&offsetoffset from plane.
Returns
Vector3<double> plane normal.
A plane is fit through the polygon vertices and the normal calculated from:
    n•p = d
where n is the normal vector, p is a point in the plane, and d is the offset.
The polygon planarity is defined as the root-mean-square-deviation from 
the fitted plane.

◆ model_assemble()

Bmolgroup * model_assemble ( Bmodel model,
Bstring paramfile,
int  separate 
)

Concatenates selected molecules into one group.

Parameters
*modelmodel parameters.
&paramfileatomic parameter file.
separateflag to generate separate molecule groups.
Returns
Bmolgroup* list of molecule groups.
Only the first model in the linked list is processed.

◆ model_average_components()

long model_average_components ( Bmodel model,
int  number 
)

Averages sequential components.

Parameters
*modelmodel structure to be modified.
numbernumber of components to average.
Returns
long number of remaining components.
Only the first component in each set with modified coordinates is kept.

◆ model_center_of_mass()

Vector3< double > model_center_of_mass ( Bmodel model)

Calculates the center-of-mass of a model.

Parameters
*modelmodel parameters.
Returns
Vector3<double> center-of-mass.
Only the first model in the list is processed.

◆ model_component_grid()

vector< vector< Bcomponent * > > model_component_grid ( Bmodel model,
Vector3< long > &  size,
Vector3< double > &  origin,
Vector3< double > &  sampling 
)

Generates lists of atoms based on a grid.

Parameters
*modelmodel list.
sizesize of grid.
originorigin of grid.
samplingspacing in each dimension.
Returns
vector<vector<Bcomponent*>> array of component arrays.
The goal is to fit all the components within the grid boundaries.
Components located outside the grid will be added to the edges.

◆ model_generate_com()

Bmodel * model_generate_com ( Bmolgroup molgroup)

Calculates the centers-of-mass of molecule group components and generates a new model.

Parameters
*molgrouplist of molecule groups.
Returns
Bmodel* new model.
Each molecule is assumed to be a component.

◆ model_geometric_median()

Vector3< double > model_geometric_median ( Bmodel model)

Calculates the geometric median of model components.

Parameters
*modelmodel parameters.
Returns
Vector3<double> geometric median.
Only the first model in the list is processed.
Based on Weiszfeld’s method.

◆ model_gyration_radius()

double model_gyration_radius ( Bmodel model)

Calculates the radius of gyration for a model.

Parameters
*modelmodel parameters.
Returns
double radius of gyration.
Only the first model in the list is processed.

◆ model_list()

long model_list ( Bmodel model)

Lists models in table form.

Parameters
*modelmodel parameters.
Returns
long number of models.

◆ model_list_comp()

long model_list_comp ( Bmodel model)

Lists models with component counts in table form.

Parameters
*modelmodel parameters.
Returns
long number of models.

◆ model_mass()

double model_mass ( Bmodel model)

Calculates the mass of a model from component masses.

Parameters
*modelmodel parameters.
Returns
double model mass.
The component type masses must be provided.
Only the first model in the list is processed.

◆ model_mass_all()

long model_mass_all ( Bmodel model)

Calculates the masses of all the models in the list.

Parameters
*modellinked list of model parameters.
Returns
long number of selected models.
The component type masses must be provided.

◆ model_merge()

long model_merge ( Bmodel model)

Merges components from all models into one.

Parameters
*modelmodel parameters.
Returns
long number of components.

◆ model_number_ids()

long model_number_ids ( Bmodel model)

Add a number to each model id.

Parameters
*modelmodel parameters.
Returns
long number of models.
The intention is to give unique id's to models.

◆ model_principal_axes() [1/3]

long model_principal_axes ( Bmodel model)

◆ model_principal_axes() [2/3]

Vector3< double > model_principal_axes ( Bmodel model,
Matrix eigenvec 
)

◆ model_principal_axes() [3/3]

Vector3< double > model_principal_axes ( Bmodel model,
Vector3< double > *  eigenvec 
)

Calculates the principal axes of a model.

Parameters
*modelmodel structure.
*eigenveceigen vectors (can be NULL).
Returns
Vector3<double> 3-valued vector of principal axes.
Only the first model in the list is processed.

◆ model_radial_distribution()

long model_radial_distribution ( Bmodel model,
double  interval 
)

Calculates the radial distribution function of a model.

Parameters
*modelmodel structure.
intervalinterval between bins.
Returns
long 0.
Only the first model in the list is processed.

◆ model_rename()

long model_rename ( Bmodel model,
char  first_name 
)

Rename models with alphabetical letters.

Parameters
*modelmodel parameters.
first_nameletter of first model.
Returns
long number of models.

◆ model_rename_components()

long model_rename_components ( Bmodel model)

Rename components.

Parameters
*modelmodel parameters.
Returns
long number of components.
The number of links to a component determines its new name.
Only the first model is processed.

◆ model_replace_components()

long model_replace_components ( Bmodel model,
Bmodel modref 
)

Replaces all components in a model with reference components.

Parameters
*modelmodel structure to be modified.
*modrefreference model.
Returns
long number of components used in reference.
Only the first model in the reference is used.
All models will have identical sets of components.

◆ model_split_into_slices()

vector< vector< Bcomponent * > > model_split_into_slices ( Bmodel model,
double  bottom,
double  top,
double  thickness 
)

Splits the components of models into slices.

Parameters
*modelmodel.
bottomminimum coordinate in z.
topmaximum coordinate in z.
thicknessslice thickness.
Returns
vector<vector<Bcomponent*>> sets of arrays of pointers to components.

◆ model_update_centers_of_mass()

long model_update_centers_of_mass ( Bmodel model,
Bmolgroup molgroup 
)

Updates the centers-of-mass of molecule group components.

Parameters
*modelmodel parameters.
*molgrouplist of molecule groups.
Returns
long number of selected components.
The identifiers of the molecule groups must correspond to the component identifiers.

◆ models_calculate_bounds()

vector< Vector3< double > > models_calculate_bounds ( Bmodel model)

Calculates the bounds of a list of models.

Parameters
*modelmodel list.
Returns
vector<Vector3<double>> minimum and maximum bounds.

◆ models_center_of_coordinates()

Vector3< double > models_center_of_coordinates ( Bmodel model)

Calculates the center-of-mass of a list of models.

Parameters
*modelmodel parameters.
Returns
Vector3<double> center-of-mass.

◆ models_get_component_array()

vector< Bcomponent * > models_get_component_array ( Bmodel model)

Generates an array of pointers to model components.

Parameters
*modelmodel.
Returns
vector<Bcomponent*> array of pointers to components.

◆ models_process() [1/5]

long models_process ( Bmodel model,
Bstring str,
long(modfunc)(Bmodel *, Bstring &str)   
)

◆ models_process() [2/5]

long models_process ( Bmodel model,
double  d,
long(modfunc)(Bmodel *, double)   
)

◆ models_process() [3/5]

long models_process ( Bmodel model,
long  i,
long(modfunc)(Bmodel *, long)   
)

◆ models_process() [4/5]

long models_process ( Bmodel model,
long(modfunc)(Bmodel *)   
)

Process a list of models using the specified function.

Parameters
*modellist of models.
modfuncfunction to be called.
Returns
long aggregate number returned by function.

◆ models_process() [5/5]

long models_process ( Bmodel model,
string &  str,
long(modfunc)(Bmodel *, string &str)   
)