![]() |
Bsoft 2.1.4
Bernard's software package
|
Library routines used for model and component selection. More...
#include "rwmodel.h"
#include "model_poly.h"
#include "model_views.h"
#include "model_compare.h"
#include "model_util.h"
#include "rwimg.h"
#include "qsort_functions.h"
#include "linked_list.h"
#include "utilities.h"
Functions | |
long | model_selection_stats (Bmodel *model) |
Calculates selection statistics. More... | |
long | model_show_selection (Bmodel *model) |
Shows the distrubution of selections. More... | |
long | model_select (Bmodel *model, Bstring &selstr) |
Selects models, components and component types. More... | |
long | model_select (Bmodel *model, long number) |
Resets the selection to all models. More... | |
long | model_select_all (Bmodel *model) |
Resets the selection to all models. More... | |
long | model_select_unknowns (Bmodel *model) |
Resets the selection to all unknown models. More... | |
long | model_reset_selection (Bmodel *model) |
Resets the selection to all components. More... | |
long | model_unset_selection (Bmodel *model) |
Unsets the selection to all components but not models. More... | |
long | model_invert_selection (Bmodel *model) |
Inverts the selection of all components. More... | |
long | model_select_sets (Bmodel *model, int size, int flag) |
Selects sets of components, each set with the same size. More... | |
long | models_select_within_bounds (Bmodel *model, Vector3< double > &start, Vector3< double > &end) |
Selects components within bounds for a list of models. More... | |
long | model_select_number_of_components (Bmodel *model, int ncomp_min, int ncomp_max) |
Selects models within a range of the number of components. More... | |
int | model_select_random (Bmodel *model, long number) |
Selects a random number of components. More... | |
long | model_select_closed (Bmodel *model, int closure_rule, int val_order) |
Determines if a polyhedron is closed given a rule. More... | |
long | model_select_fullerene (Bmodel *model) |
Selects fullerene type models. More... | |
long | model_select_non_fullerene (Bmodel *model) |
Selects non-fullerene type models. More... | |
long | model_select_valence (Bmodel *model, int valence) |
Selects fullerene type models. More... | |
long | model_select_polygons (Bmodel *model, int order) |
Selects model polygons with a given order. More... | |
long | model_select_first (Bmodel *model, int first) |
Selects the first number of components. More... | |
long | model_select_within_shell (Bmodel *model, Vector3< double > center, double minrad, double maxrad) |
Selects components within a shell. More... | |
long | model_select_in_mask (Bmodel *model, Bimage *pmask) |
Selects components within a mask. More... | |
long | model_delete (Bmodel **model) |
Deletes tagged models, components and links. More... | |
long | model_delete_comp_type (Bmodel *model, Bstring &comptype) |
Deletes components based on their types. More... | |
long | model_delete_non_selected (Bmodel **model) |
Deletes selected components and associated links from the model. More... | |
long | model_type_from_selection (Bmodel *model, Bstring *comp_type, Bstring &filename) |
Converts selection numbers to different types. More... | |
long | model_fom_deselect (Bmodel *model, double fom_cutoff) |
Deselects components based on the FOM. More... | |
long | model_fom_max_fraction_deselect (Bmodel *model, double fom_fraction) |
Deselects components based on the FOM relative to the maximum FOM. More... | |
long | models_radius_deselect (Bmodel *model, double minrad, double maxrad) |
Deselects components based on the distance from the origin. More... | |
long | model_fom_histogram (Bmodel *model, double fom_step) |
Selects components based on the FOM. More... | |
long | model_fom_ranking (Bmodel *model, int nrank) |
Ranks components based on the FOM. More... | |
long | model_delete_overlapped_components (Bmodel **model, double distance) |
Deletes components that are too close to others. More... | |
long | model_average_overlapped_components (Bmodel *model, double distance) |
Averages components within a distance. More... | |
long | model_prune_simple (Bmodel *model, double mindist) |
Prunes overlapping components based on a distance criterion. More... | |
long | model_prune_fom_old (Bmodel *model, double distance) |
Deselects components based on overlap distance and the FOM. More... | |
long | model_prune_fom (Bmodel *model, double distance) |
long | model_prune_similar (Bmodel *model) |
Averages similar components and reduces redundancy. More... | |
Bimage * | model_composite_map (Bmodel *model, Vector3< int > size, Vector3< double > origin, Vector3< double > sampling) |
long | comp_select_one (Bcomponent *comp, long i, long n, double *ovm) |
double * | model_overlap_matrix (Bmodel *model, double distance) |
long | model_prune_fit (Bmodel *model, double distance) |
Deselects components based on overlap distance and the FOM. More... | |
long | model_prune_large (Bmodel *model, double sampling) |
Deselects components based on overlap distance and the FOM. More... | |
long | model_find_overlap (Bmodel *model, Bstring &reffile, double distance) |
Finds the components that overlap with the reference. More... | |
Variables | |
int | verbose |
Library routines used for model and component selection.
long comp_select_one | ( | Bcomponent * | comp, |
long | i, | ||
long | n, | ||
double * | ovm | ||
) |
long model_average_overlapped_components | ( | Bmodel * | model, |
double | distance | ||
) |
Averages components within a distance.
*model | model. |
distance | overlap distance. |
Components within the given distance are averaged and the first instance is retained while the others are deleted.
Bimage * model_composite_map | ( | Bmodel * | model, |
Vector3< int > | size, | ||
Vector3< double > | origin, | ||
Vector3< double > | sampling | ||
) |
long model_delete | ( | Bmodel ** | model | ) |
Deletes tagged models, components and links.
**model | pointer to model parameters. |
A model, component or link is tagged for deletion by a negative selection flag.
Deletes components based on their types.
*model | model parameters. |
&comptype | component type. |
long model_delete_non_selected | ( | Bmodel ** | model | ) |
Deletes selected components and associated links from the model.
**model | pointer to model parameters. |
long model_delete_overlapped_components | ( | Bmodel ** | model, |
double | distance | ||
) |
Deletes components that are too close to others.
**model | model handle. |
distance | distance to test for. |
If a component is too close to another, the component further down in the list is deleted.
Finds the components that overlap with the reference.
*model | model parameters. |
&reffile | reference molecule file name. |
distance | distance for overlap. |
All components that are within the distance criterion of the components in the reference model are deselected.
long model_fom_deselect | ( | Bmodel * | model, |
double | fom_cutoff | ||
) |
Deselects components based on the FOM.
*model | model parameters. |
fom_cutoff | FOM threshold. |
The FOM is assumed to be a value from 0 to 1 and components with FOM's below the given cutoff are deselected.
long model_fom_histogram | ( | Bmodel * | model, |
double | fom_step | ||
) |
Selects components based on the FOM.
*model | model parameters. |
fom_step | FOM step size. |
The FOM is assumed to represent energy and therefore values below the given cutoff are selected.
long model_fom_max_fraction_deselect | ( | Bmodel * | model, |
double | fom_fraction | ||
) |
Deselects components based on the FOM relative to the maximum FOM.
*model | model parameters. |
fom_fraction | FOM threshold. |
The FOM/FOMmax ratio is calculated and components with ratios below the given fraction are deselected.
long model_fom_ranking | ( | Bmodel * | model, |
int | nrank | ||
) |
Ranks components based on the FOM.
*model | model parameters. |
nrank | number of levels. |
Each selected component is ranked into one of a number of levels.
long model_invert_selection | ( | Bmodel * | model | ) |
Inverts the selection of all components.
*model | model parameters. |
double * model_overlap_matrix | ( | Bmodel * | model, |
double | distance | ||
) |
long model_prune_fit | ( | Bmodel * | model, |
double | distance | ||
) |
Deselects components based on overlap distance and the FOM.
*model | model parameters. |
distance | overlap distance. |
The selected component with the highest FOM is chosen as the first component in the solution set and all components within the distance criterion from it are deselected. The selected component with the next highest FOM is then chosen and again, all components within the distance criterion from it are deselected. This is repeated until there are no more selected components to choose. Only the first model is processed.
long model_prune_fom | ( | Bmodel * | model, |
double | distance | ||
) |
long model_prune_fom_old | ( | Bmodel * | model, |
double | distance | ||
) |
Deselects components based on overlap distance and the FOM.
*model | model parameters. |
distance | overlap distance. |
The selected component with the highest FOM is chosen as the first component in the solution set and all components within the distance criterion from it are deselected. The selected component with the next highest FOM is then chosen and again, all components within the distance criterion from it are deselected. This is repeated until there are no more selected components to choose. Only the first model is processed.
long model_prune_large | ( | Bmodel * | model, |
double | sampling | ||
) |
Deselects components based on overlap distance and the FOM.
*model | model parameters. |
sampling | grid sampling. |
Each component location is mapped to a grid point and the component FOM is assigned to the grid point if it is larger. Peaks are determined within each 3x3x3 kernel: a peak is defined as having a higher value than any of its 26 neigbors. The component associated with the value at a peak is selected. Only the first model is processed.
long model_prune_similar | ( | Bmodel * | model | ) |
Averages similar components and reduces redundancy.
*model | model parameters. |
A view distance matrix is calculated for all the components. All components within a small distance (~10% of model radius) are averaged into the first component and the rest are deselected. Only the first model is processed.
long model_prune_simple | ( | Bmodel * | model, |
double | mindist | ||
) |
Prunes overlapping components based on a distance criterion.
*model | model structure to be modified. |
mindist | distance criterion. |
The first component in any pair of overlapping components is kept.
long model_reset_selection | ( | Bmodel * | model | ) |
Resets the selection to all components.
*model | model parameters. |
Only one model is modified.
Selects models, components and component types.
*model | model parameters. |
&selstr | selection string. |
The selection string can have one of the following formats: #model@component #model%comp_type ^model_type@component ^model_type%comp_type Only elements originally selected is considered, except where the "." is used.
long model_select | ( | Bmodel * | model, |
long | number | ||
) |
Resets the selection to all models.
*model | model parameters. |
number | selection number to select. |
long model_select_all | ( | Bmodel * | model | ) |
Resets the selection to all models.
*model | model parameters. |
long model_select_closed | ( | Bmodel * | model, |
int | closure_rule, | ||
int | val_order | ||
) |
Determines if a polyhedron is closed given a rule.
*model | model parameters. |
closure_rule | 1=valency, 2=order. |
val_order | magnitude of valency or order. |
Polyhedron closure is arbitrarily decided by a rule: 1. Fixed valency. 3. Fixed polygon order. Models that fail the given rule are deselected. If the given valency or order is zero, it is set to the value for the first component or polygon.
long model_select_first | ( | Bmodel * | model, |
int | first | ||
) |
Selects the first number of components.
*model | model parameters. |
first | first number of components to select. |
long model_select_fullerene | ( | Bmodel * | model | ) |
Selects fullerene type models.
*model | model structure. |
Each model is tested for the presence of polygons that have orders only of five or six.
Selects components within a mask.
*model | model parameters. |
pmask | first number of components to select. |
long model_select_non_fullerene | ( | Bmodel * | model | ) |
Selects non-fullerene type models.
*model | model structure. |
Each model is tested for the presence of polygons that have orders different from five or six.
long model_select_number_of_components | ( | Bmodel * | model, |
int | ncomp_min, | ||
int | ncomp_max | ||
) |
Selects models within a range of the number of components.
*model | model parameters. |
ncomp_min | minimum number of components. |
ncomp_max | maximum number of components. |
long model_select_polygons | ( | Bmodel * | model, |
int | order | ||
) |
Selects model polygons with a given order.
*model | model to color. |
order | polygon order. |
int model_select_random | ( | Bmodel * | model, |
long | number | ||
) |
Selects a random number of components.
*model | model parameters. |
number | number of components. |
long model_select_sets | ( | Bmodel * | model, |
int | size, | ||
int | flag | ||
) |
Selects sets of components, each set with the same size.
*model | parameter structure with all parameters. |
size | number of components in each set. |
flag | flag to not count across model boundaries. |
Sets up sets of components, each set identified as a number in the selection array.
long model_select_unknowns | ( | Bmodel * | model | ) |
Resets the selection to all unknown models.
*model | model parameters. |
long model_select_valence | ( | Bmodel * | model, |
int | valence | ||
) |
Selects fullerene type models.
*model | model structure. |
valence | component valence |
Each model is tested for the presence of polygons that have orders only of five or six.
long model_select_within_shell | ( | Bmodel * | model, |
Vector3< double > | center, | ||
double | minrad, | ||
double | maxrad | ||
) |
Selects components within a shell.
*model | model parameters. |
center | center of the shell. |
minrad | minimum radius of the sphere. |
maxrad | radius of the sphere. |
long model_selection_stats | ( | Bmodel * | model | ) |
Calculates selection statistics.
*model | model parameters. |
The FOM is assumed to be a value from 0 to 1.
long model_show_selection | ( | Bmodel * | model | ) |
Shows the distrubution of selections.
*model | model parameters. |
Converts selection numbers to different types.
*model | model parameters. |
*comp_type | linked list of strings. |
&filename | component type file name. |
long model_unset_selection | ( | Bmodel * | model | ) |
Unsets the selection to all components but not models.
*model | model parameters. |
long models_radius_deselect | ( | Bmodel * | model, |
double | minrad, | ||
double | maxrad | ||
) |
Deselects components based on the distance from the origin.
*model | model parameters. |
minrad | minimum distance from origin. |
maxrad | maximum distance from the origin. |
All components inside a minmum radius and beyond a maximum radius are deselected.
long models_select_within_bounds | ( | Bmodel * | model, |
Vector3< double > & | start, | ||
Vector3< double > & | end | ||
) |
Selects components within bounds for a list of models.
*model | list of models. |
&start | minimum coordinates. |
&end | maximum coordinates. |
|
extern |