Main Page   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Related Pages  

VolumetricData Class Reference

Volumetric data class for potential maps, electron density maps, etc. More...

#include <VolumetricData.h>

List of all members.

Public Methods

 VolumetricData (const char *name, const float *origin, const float *xaxis, const float *yaxis, const float *zaxis, int xs, int ys, int zs, float *dataptr)
 constructors for both single- and double-precision axes. More...

 VolumetricData (const char *name, const double *origin, const double *xaxis, const double *yaxis, const double *zaxis, int xs, int ys, int zs, float *dataptr)
 constructor for double precision origin/axes info. More...

 ~VolumetricData ()
 destructor. More...

ptrdiff_t gridsize () const
 return total number of gridpoints. More...

void datarange (float &min, float &max)
 return min/max data values. More...

void set_name (const char *name)
 Sets data name to an internal copy of the provided string. More...

float mean ()
 return the mean value of the density map, implemented via O(1) access to a cached value. More...

float sigma ()
 return the standard deviation (sigma) of the density map, implemented via O(1) access to a cached value. More...

double integral ()
 return integral of map as product of mean() times total volume. More...

void cell_lengths (float *xl, float *yl, float *zl) const
 return cell side lengths. More...

void cell_axes (float *xax, float *yax, float *zax) const
 return cell axes. More...

void cell_axes (double *xax, double *yax, double *zax) const
 return cell axes. More...

void cell_dirs (float *xax, float *yax, float *zax) const
 return cell axes directions. More...

double cell_volume () const
 return cell volume. More...

void set_volume_axes (float *xax, float *yax, float *zax)
 set volume/cell lengths. More...

void set_volume_axes (double *xax, double *yax, double *zax)
void scale_volume (double scalex, double scaley, double scalez)
 scale volume/cell lengths. More...

void set_volume_origin (float *org)
 set volume origin. More...

void set_volume_origin (double *org)
void voxel_coord_from_cartesian_coord (const float *carcoord, float *voxcoord, int shiftflag) const
 return volumetric coordinate from cartesian coordinate. More...

ptrdiff_t voxel_index_from_coord (float xpos, float ypos, float zpos) const
 return index of the voxel nearest to a cartesian coordinate. More...

float voxel_value (int x, int y, int z) const
 return voxel at requested index, no safety checks. More...

float voxel_value_safe (int x, int y, int z) const
 return voxel, after safely clamping index to valid range. More...

float voxel_value_interpolate (float xv, float yv, float zv) const
 return interpolated value from 8 nearest neighbor voxels. More...

float voxel_value_from_coord (float xpos, float ypos, float zpos) const
 return voxel value based on cartesian coordinates. More...

float voxel_value_interpolate_from_coord (float xpos, float ypos, float zpos) const
 return interpolated value of voxel, based on cartesian coords. More...

float voxel_value_from_coord_safe (float xpos, float ypos, float zpos) const
 return voxel value based on cartesian coordinates. safe versions return zero if coordinates are outside the map. More...

float voxel_value_interpolate_from_coord_safe (float xpos, float ypos, float zpos) const
 return interpolated value of voxel, based on cartesian coords. this safe version returns zero if the coordinates are outside the map. More...

const float * access_volume_gradient ()
 get read-only access to the gradient. More...

void set_volume_gradient (float *gradient)
 provide the volume gradient. More...

void compute_volume_gradient (void)
 (re)compute the volume gradient. More...

void voxel_gradient_fast (int x, int y, int z, float *grad) const
 return gradient at requested index, no safety checks. More...

void voxel_gradient_safe (int x, int y, int z, float *grad) const
 return gradient, after safely clamping index to valid range. More...

void voxel_gradient_interpolate (const float *voxcoord, float *gradient) const
 interpolate the gradient between the eight neighboring voxels. More...

void voxel_gradient_from_coord (const float *coord, float *gradient) const
 return voxel gradient based on cartesian coordinates. More...

void voxel_gradient_interpolate_from_coord (const float *coord, float *gradient) const
 return interpolated voxel gradient for cartesian coordinate. More...

void voxel_coord (int x, int y, int z, float &gx, float &gy, float &gz) const
 get the cartesian coordinate of a voxel given its x,y,z indices. More...

void voxel_coord (ptrdiff_t i, float &x, float &y, float &z) const
 get the cartesian coordinate of a voxel given its 1-D index. More...

void pad (int padxm, int padxp, int padym, int padyp, int padzm, int padzp)
 add or remove voxels in the given axis directions. More...

void crop (double crop_minx, double crop_miny, double crop_minz, double crop_maxx, double crop_maxy, double crop_maxz)
 crop a volumetric data to a given set of cartesian coordinates. More...

void clamp (float min_value, float max_value)
 clamp out of range voxel values. More...

void scale_by (float ff)
 scales voxel data by given amount. More...

void scalar_add (float ff)
 add scalar value to to all voxels. More...

void rescale_voxel_value_range (float min_value, float max_value)
 rescale voxel data to a given range. More...

void downsample ()
 decimate/dowmnsample voxels by 2 in each dimension (x8 total reduction). More...

void supersample ()
 refine/supersample voxels by 2 in each dimension (x8 total increase). More...

void sigma_scale ()
 Transform map to a sigma scale, so that isovalues in VMD correspond to number of sigmas above the mean. More...

void binmask (float threshold=0.0f)
 Make a binary mask out of a map, i.e. map values > threshold are set to 1, and all others are set to 0. More...

void gaussian_blur (double sigma)
 Guassian blur by sigma. More...

void mdff_potential (float threshold)
 Create a potential for use with MDFF. More...


Public Attributes

double origin [3]
 origin of volume (x=0, y=0, z=0 corner). More...

double xaxis [3]
 direction and length for X axis (non-unit). More...

double yaxis [3]
 direction and length for Y axis (non-unit). More...

double zaxis [3]
 direction and length for Z axis (non-unit). More...

int xsize
int ysize
int zsize
 number of samples along each axis. More...

char * name
 human-readable volume dataset identifier. More...

float * data
 raw data, total of xsize*ysize*zsize voxels stored x varying fastest, then y, then z. More...


Detailed Description

Volumetric data class for potential maps, electron density maps, etc.

Definition at line 35 of file VolumetricData.h.


Constructor & Destructor Documentation

VolumetricData::VolumetricData const char *    name,
const float *    origin,
const float *    xaxis,
const float *    yaxis,
const float *    zaxis,
int    xs,
int    ys,
int    zs,
float *    dataptr
 

constructors for both single- and double-precision axes.

Definition at line 41 of file VolumetricData.C.

References data, name, NULL, origin, stringdup, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

VolumetricData::VolumetricData const char *    name,
const double *    origin,
const double *    xaxis,
const double *    yaxis,
const double *    zaxis,
int    xs,
int    ys,
int    zs,
float *    dataptr
 

constructor for double precision origin/axes info.

Definition at line 68 of file VolumetricData.C.

References data, name, NULL, origin, stringdup, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

VolumetricData::~VolumetricData  
 

destructor.

Definition at line 95 of file VolumetricData.C.

References data, and name.


Member Function Documentation

const float * VolumetricData::access_volume_gradient  
 

get read-only access to the gradient.

Definition at line 415 of file VolumetricData.C.

References compute_volume_gradient, and NULL.

Referenced by IsoSurface::DoCellGeneral, and IsoSurface::DoGridPosNorms.

void VolumetricData::binmask float    threshold = 0.0f
 

Make a binary mask out of a map, i.e. map values > threshold are set to 1, and all others are set to 0.

Definition at line 1030 of file VolumetricData.C.

References data, gridsize, and threshold.

Referenced by density_binmask.

void VolumetricData::cell_axes double *    xax,
double *    yax,
double *    zax
const
 

return cell axes.

Definition at line 172 of file VolumetricData.C.

References xaxis, xsize, yaxis, ysize, zaxis, and zsize.

void VolumetricData::cell_axes float *    xax,
float *    yax,
float *    zax
const
 

return cell axes.

Definition at line 138 of file VolumetricData.C.

References xaxis, xsize, yaxis, ysize, zaxis, and zsize.

Referenced by cell_dirs, cell_volume, IsoSurface::compute, VolMapCreateCoulombPotentialMSM::compute_frame, VolMapCreateCoulombPotential::compute_frame, VolMapCreateOccupancy::compute_frame, VolMapCreateDensity::compute_frame, VolMapCreateMask::compute_frame, density_info, volmap_write_dx_file, voxel_coord, and voxel_coord.

void VolumetricData::cell_dirs float *    xax,
float *    yax,
float *    zax
const
 

return cell axes directions.

Definition at line 206 of file VolumetricData.C.

References cell_axes, and cell_lengths.

Referenced by IsoSurface::compute.

void VolumetricData::cell_lengths float *    xl,
float *    yl,
float *    zl
const
 

return cell side lengths.

Definition at line 112 of file VolumetricData.C.

References dot_prod, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

Referenced by cell_dirs, and compute_volume_gradient.

double VolumetricData::cell_volume   const
 

return cell volume.

Definition at line 227 of file VolumetricData.C.

References cell_axes.

Referenced by density_info, and integral.

void VolumetricData::clamp float    min_value,
float    max_value
 

clamp out of range voxel values.

Definition at line 693 of file VolumetricData.C.

References data, and gridsize.

Referenced by density_clamp.

void VolumetricData::compute_volume_gradient void   
 

(re)compute the volume gradient.

Definition at line 442 of file VolumetricData.C.

References cell_lengths, clamp_int, data, gridsize, voxel_value, xsize, ysize, and zsize.

Referenced by access_volume_gradient, BaseMolecule::add_volume_data, and vmd_measure_volinterior.

void VolumetricData::crop double    crop_minx,
double    crop_miny,
double    crop_minz,
double    crop_maxx,
double    crop_maxy,
double    crop_maxz
 

crop a volumetric data to a given set of cartesian coordinates.

Definition at line 674 of file VolumetricData.C.

References origin, pad, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

Referenced by density_crop.

void VolumetricData::datarange float &    min,
float &    max
 

return min/max data values.

Definition at line 1116 of file VolumetricData.C.

Referenced by BaseMolecule::add_volume_data, density_info, AtomColor::find, VolumeTexture::generateContourLineTexture, histogram, and GraphicsFltkMenu::update_molchooser.

void VolumetricData::downsample  
 

decimate/dowmnsample voxels by 2 in each dimension (x8 total reduction).

Definition at line 792 of file VolumetricData.C.

References data, n, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

Referenced by density_downsample.

void VolumetricData::gaussian_blur double    sigma
 

Guassian blur by sigma.

Definition at line 1051 of file VolumetricData.C.

References GaussianBlur::blur, data, GaussianBlur::get_image, sigma, xsize, ysize, and zsize.

Referenced by density_smooth.

ptrdiff_t VolumetricData::gridsize   const [inline]
 

return total number of gridpoints.

Definition at line 59 of file VolumetricData.h.

References xsize, ysize, and zsize.

Referenced by add, BaseMolecule::add_volume_data, average, binmask, clamp, VolMapCreateILS::compute, compute_volume_gradient, countIsoGrids, histogram, init_from_identity, init_from_intersection, init_from_union, markIsoGrid, mdff_potential, multiply, normalize_pmap, process_pmap, rescale_voxel_value_range, scalar_add, scale_by, segment_volume, sigma_scale, subtract, vmd_volmap_compare, vol_com, and vol_probability.

double VolumetricData::integral  
 

return integral of map as product of mean() times total volume.

Definition at line 1190 of file VolumetricData.C.

References cell_volume, mean, xsize, ysize, and zsize.

Referenced by density_info.

void VolumetricData::mdff_potential float    threshold
 

Create a potential for use with MDFF.

Definition at line 1073 of file VolumetricData.C.

References data, gridsize, and threshold.

Referenced by density_mdff_potential.

float VolumetricData::mean  
 

return the mean value of the density map, implemented via O(1) access to a cached value.

Definition at line 1159 of file VolumetricData.C.

Referenced by density_info, integral, and sigma_scale.

void VolumetricData::pad int    padxm,
int    padxp,
int    padym,
int    padyp,
int    padzm,
int    padzp
 

add or remove voxels in the given axis directions.

Definition at line 601 of file VolumetricData.C.

References data, MAX, MIN, origin, vec_scaled_add, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

Referenced by crop, and density_trim.

void VolumetricData::rescale_voxel_value_range float    min_value,
float    max_value
 

rescale voxel data to a given range.

Definition at line 767 of file VolumetricData.C.

References data, and gridsize.

Referenced by density_range.

void VolumetricData::scalar_add float    ff
 

add scalar value to to all voxels.

Definition at line 749 of file VolumetricData.C.

References data, and gridsize.

Referenced by density_sadd.

void VolumetricData::scale_by float    ff
 

scales voxel data by given amount.

Definition at line 724 of file VolumetricData.C.

References data, and gridsize.

Referenced by density_smult.

void VolumetricData::scale_volume double    scalex,
double    scaley,
double    scalez
[inline]
 

scale volume/cell lengths.

Definition at line 124 of file VolumetricData.h.

References xaxis, yaxis, and zaxis.

Referenced by text_cmd_mol.

void VolumetricData::set_name const char *    name
 

Sets data name to an internal copy of the provided string.

Definition at line 104 of file VolumetricData.C.

References name.

Referenced by VolMapCreateCoulombPotentialMSM::compute_init, VolMapCreateCoulombPotential::compute_init, VolMapCreateDistance::compute_init, VolMapCreateOccupancy::compute_init, VolMapCreateInterp::compute_init, VolMapCreateDensity::compute_init, VolMapCreateMask::compute_init, and VolMapCreateILS::write_map.

void VolumetricData::set_volume_axes double *    xax,
double *    yax,
double *    zax
[inline]
 

Definition at line 108 of file VolumetricData.h.

References xaxis, yaxis, and zaxis.

void VolumetricData::set_volume_axes float *    xax,
float *    yax,
float *    zax
[inline]
 

set volume/cell lengths.

Definition at line 93 of file VolumetricData.h.

References xaxis, yaxis, and zaxis.

Referenced by text_cmd_mol.

void VolumetricData::set_volume_gradient float *    gradient
 

provide the volume gradient.

Definition at line 424 of file VolumetricData.C.

References NULL, xsize, ysize, and zsize.

Referenced by BaseMolecule::add_volume_data.

void VolumetricData::set_volume_origin double *    org [inline]
 

Definition at line 147 of file VolumetricData.h.

References origin.

void VolumetricData::set_volume_origin float *    org [inline]
 

set volume origin.

Definition at line 140 of file VolumetricData.h.

References origin.

Referenced by text_cmd_mol.

float VolumetricData::sigma  
 

return the standard deviation (sigma) of the density map, implemented via O(1) access to a cached value.

Definition at line 1176 of file VolumetricData.C.

Referenced by density_info, gaussian_blur, and sigma_scale.

void VolumetricData::sigma_scale  
 

Transform map to a sigma scale, so that isovalues in VMD correspond to number of sigmas above the mean.

Definition at line 1009 of file VolumetricData.C.

References data, gridsize, mean, and sigma.

Referenced by density_sigma.

void VolumetricData::supersample  
 

refine/supersample voxels by 2 in each dimension (x8 total increase).

Definition at line 845 of file VolumetricData.C.

References cubic_interp, data, xsize, ysize, and zsize.

Referenced by density_supersample.

void VolumetricData::voxel_coord ptrdiff_t    i,
float &    x,
float &    y,
float &    z
const [inline]
 

get the cartesian coordinate of a voxel given its 1-D index.

Definition at line 221 of file VolumetricData.h.

References cell_axes, origin, xsize, ysize, and z.

void VolumetricData::voxel_coord int    x,
int    y,
int    z,
float &    gx,
float &    gy,
float &    gz
const [inline]
 

get the cartesian coordinate of a voxel given its x,y,z indices.

Definition at line 210 of file VolumetricData.h.

References cell_axes, origin, and z.

Referenced by add, average, multiply, subtract, and vol_com.

void VolumetricData::voxel_coord_from_cartesian_coord const float *    carcoord,
float *    voxcoord,
int    shiftflag
const
 

return volumetric coordinate from cartesian coordinate.

Definition at line 242 of file VolumetricData.C.

References Matrix4::inverse, Matrix4::multpoint3d, origin, xaxis, xsize, yaxis, ysize, zaxis, and zsize.

Referenced by colvarproxy_vmd::compute_voldata, voxel_gradient_from_coord, voxel_gradient_interpolate_from_coord, voxel_index_from_coord, voxel_value_interpolate_from_coord, and voxel_value_interpolate_from_coord_safe.

void VolumetricData::voxel_gradient_fast int    x,
int    y,
int    z,
float *    grad
const [inline]
 

return gradient at requested index, no safety checks.

Definition at line 192 of file VolumetricData.h.

References xsize, ysize, and z.

void VolumetricData::voxel_gradient_from_coord const float *    coord,
float *    gradient
const
 

return voxel gradient based on cartesian coordinates.

Definition at line 570 of file VolumetricData.C.

References voxel_coord_from_cartesian_coord, and voxel_gradient_safe.

void VolumetricData::voxel_gradient_interpolate const float *    voxcoord,
float *    gradient
const
 

interpolate the gradient between the eight neighboring voxels.

Definition at line 531 of file VolumetricData.C.

References vec_lerp, voxel_gradient_safe, and z.

Referenced by colvarproxy_vmd::compute_voldata, and voxel_gradient_interpolate_from_coord.

void VolumetricData::voxel_gradient_interpolate_from_coord const float *    coord,
float *    gradient
const
 

return interpolated voxel gradient for cartesian coordinate.

Definition at line 582 of file VolumetricData.C.

References NAN, voxel_coord_from_cartesian_coord, voxel_gradient_interpolate, xsize, ysize, and zsize.

void VolumetricData::voxel_gradient_safe int    x,
int    y,
int    z,
float *    grad
const
 

return gradient, after safely clamping index to valid range.

Definition at line 518 of file VolumetricData.C.

References xsize, ysize, z, and zsize.

Referenced by voxel_gradient_from_coord, and voxel_gradient_interpolate.

ptrdiff_t VolumetricData::voxel_index_from_coord float    xpos,
float    ypos,
float    zpos
const
 

return index of the voxel nearest to a cartesian coordinate.

Definition at line 280 of file VolumetricData.C.

References voxel_coord_from_cartesian_coord, xsize, ysize, and zsize.

Referenced by atomsel_gridindex_array, voxel_value_from_coord, and voxel_value_from_coord_safe.

float VolumetricData::voxel_value int    x,
int    y,
int    z
const [inline]
 

return voxel at requested index, no safety checks.

Definition at line 162 of file VolumetricData.h.

References data, xsize, ysize, and z.

Referenced by compute_volume_gradient, and volmap_write_dx_file.

float VolumetricData::voxel_value_from_coord float    xpos,
float    ypos,
float    zpos
const
 

return voxel value based on cartesian coordinates.

Definition at line 343 of file VolumetricData.C.

References data, NAN, and voxel_index_from_coord.

Referenced by atomsel_volume_array, and multiply.

float VolumetricData::voxel_value_from_coord_safe float    xpos,
float    ypos,
float    zpos
const
 

return voxel value based on cartesian coordinates. safe versions return zero if coordinates are outside the map.

Definition at line 375 of file VolumetricData.C.

References data, and voxel_index_from_coord.

Referenced by add, average, and subtract.

float VolumetricData::voxel_value_interpolate float    xv,
float    yv,
float    zv
const
 

return interpolated value from 8 nearest neighbor voxels.

Definition at line 312 of file VolumetricData.C.

References voxel_value_safe, and z.

Referenced by colvarproxy_vmd::compute_voldata, VolumeTexture::generateContourLineTexture, voxel_value_interpolate_from_coord, and voxel_value_interpolate_from_coord_safe.

float VolumetricData::voxel_value_interpolate_from_coord float    xpos,
float    ypos,
float    zpos
const
 

return interpolated value of voxel, based on cartesian coords.

Definition at line 353 of file VolumetricData.C.

References NAN, voxel_coord_from_cartesian_coord, voxel_value_interpolate, xsize, ysize, and zsize.

Referenced by atomsel_interp_volume_array, and multiply.

float VolumetricData::voxel_value_interpolate_from_coord_safe float    xpos,
float    ypos,
float    zpos
const
 

return interpolated value of voxel, based on cartesian coords. this safe version returns zero if the coordinates are outside the map.

Definition at line 386 of file VolumetricData.C.

References voxel_coord_from_cartesian_coord, voxel_value_interpolate, xsize, ysize, and zsize.

Referenced by add, average, and subtract.

float VolumetricData::voxel_value_safe int    x,
int    y,
int    z
const
 

return voxel, after safely clamping index to valid range.

Definition at line 301 of file VolumetricData.C.

References data, xsize, ysize, z, and zsize.

Referenced by voxel_value_interpolate.


Member Data Documentation

float* VolumetricData::data
 

raw data, total of xsize*ysize*zsize voxels stored x varying fastest, then y, then z.

Definition at line 43 of file VolumetricData.h.

Referenced by add, VolMapCreateILS::add_map_to_molecule, average, binmask, cc_threaded, clamp, VolMapCreate::combo_begin, VolMapCreate::combo_export, VolMapCreateILS::compute, VolMapCreate::compute_all, compute_volume_gradient, countIsoGrids, IsoContour::DoCell, IsoSurface::DoCellGeneral, IsoSurface::DoGridPosNorms, downsample, gaussian_blur, VolumeTexture::generateChargeTexture, VolumeTexture::generateColorScaleTexture, VolumeTexture::generateHSVTexture, histogram, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, markIsoGrid, mdff_cc, mdff_potential, multiply, normalize_pmap, pad, process_pmap, py_mol_get_volumetric, RaycastGrid, rescale_voxel_value_range, scalar_add, scale_by, segment_volume, Segmentation::Segmentation, sigma_scale, subtract, supersample, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, vol_com, vol_probability, VolIn_CleanGrid, volin_threaded, volin_threaded_prob, volmap_write_dx_file, VolumetricData, voxel_value, voxel_value_from_coord, voxel_value_from_coord_safe, voxel_value_safe, MolFilePlugin::write_volumetric, and ~VolumetricData.

char* VolumetricData::name
 

human-readable volume dataset identifier.

Definition at line 42 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, BaseMolecule::add_volume_data, mdff_cc, SaveTrajectoryFltkMenu::molchooser_activate_selection, set_name, GraphicsFltkMenu::update_molchooser, vmd_measure_volinterior, vmd_volmap_new_fromtype, GraphicsFltkMenu::volindex_update, volmap_write_dx_file, VolumetricData, MolFilePlugin::write_volumetric, and ~VolumetricData.

double VolumetricData::origin[3]
 

origin of volume (x=0, y=0, z=0 corner).

Definition at line 37 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, VolMapCreateCoulombPotentialMSM::compute_frame, VolMapCreateCoulombPotential::compute_frame, VolMapCreateDistance::compute_frame, VolMapCreateOccupancy::compute_frame, VolMapCreateInterp::compute_frame, VolMapCreateDensity::compute_frame, VolMapCreateMask::compute_frame, VolMapCreate::compute_init, DrawMolecule::cov, CreateEmptyGrid, CreateProbGrid, crop, density_info, density_rotate, IsoContour::DoCell, IsoSurface::DoCellGeneral, IsoSurface::DoGridPosNorms, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, mdff_cc, normalize_pmap, pad, prepare_texture_coordinates, process_pmap, py_mol_get_volumetric, segment_volume, IsoSurface::set_color_voltex_rgb3fv, set_volume_origin, text_cmd_mol, vmd_cuda_calc_density, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, vol_move, vol_moveto, volmap_write_dx_file, VolumetricData, voxel_coord, voxel_coord, voxel_coord_from_cartesian_coord, and MolFilePlugin::write_volumetric.

double VolumetricData::xaxis[3]
 

direction and length for X axis (non-unit).

Definition at line 38 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, cell_axes, cell_lengths, IsoContour::compute, VolMapCreate::compute_init, DrawMolecule::cov, CreateEmptyGrid, CreateProbGrid, crop, density_rotate, downsample, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, mdff_cc, normalize_pmap, pad, prepare_texture_coordinates, process_pmap, py_mol_get_volumetric, DrawMolecule::scale_factor, scale_volume, segment_volume, set_volume_axes, text_cmd_mol, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, vol_move, VolumetricData, voxel_coord_from_cartesian_coord, and MolFilePlugin::write_volumetric.

int VolumetricData::xsize
 

Definition at line 41 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, cc_threaded, cell_axes, cell_lengths, VolMapCreate::combo_addframe, VolMapCreate::combo_begin, VolMapCreate::combo_export, VolMapCreateILS::compute, IsoSurface::compute, IsoContour::compute, VolMapCreate::compute_all, VolMapCreateCoulombPotentialMSM::compute_frame, VolMapCreateCoulombPotential::compute_frame, VolMapCreateDistance::compute_frame, VolMapCreateOccupancy::compute_frame, VolMapCreateInterp::compute_frame, VolMapCreateDensity::compute_frame, VolMapCreateMask::compute_frame, VolMapCreate::compute_init, colvarproxy_vmd::compute_voldata, compute_volume_gradient, CreateEmptyGrid, CreateProbGrid, crop, density_info, IsoContour::DoCell, IsoSurface::DoCellGeneral, IsoSurface::DoGridPosNorms, downsample, gaussian_blur, VolumeTexture::generateChargeTexture, VolumeTexture::generateColorScaleTexture, VolumeTexture::generateContourLineTexture, VolumeTexture::generateHSVTexture, gridsize, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, integral, mdff_cc, normalize_pmap, pad, process_pmap, py_mol_get_volumetric, RaycastGrid, segment_volume, Segmentation::Segmentation, IsoSurface::set_color_voltex_rgb3fv, set_volume_gradient, supersample, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, VolIn_CleanGrid, volin_threaded, volin_threaded_prob, volmap_write_dx_file, VolumetricData, voxel_coord, voxel_coord, voxel_coord_from_cartesian_coord, voxel_gradient_fast, voxel_gradient_interpolate_from_coord, voxel_gradient_safe, voxel_index_from_coord, voxel_value, voxel_value_interpolate_from_coord, voxel_value_interpolate_from_coord_safe, voxel_value_safe, and MolFilePlugin::write_volumetric.

double VolumetricData::yaxis[3]
 

direction and length for Y axis (non-unit).

Definition at line 39 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, cell_axes, cell_lengths, IsoContour::compute, VolMapCreate::compute_init, DrawMolecule::cov, CreateEmptyGrid, CreateProbGrid, crop, density_rotate, downsample, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, mdff_cc, normalize_pmap, pad, prepare_texture_coordinates, process_pmap, py_mol_get_volumetric, DrawMolecule::scale_factor, scale_volume, segment_volume, set_volume_axes, text_cmd_mol, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, vol_move, VolumetricData, voxel_coord_from_cartesian_coord, and MolFilePlugin::write_volumetric.

int VolumetricData::ysize
 

Definition at line 41 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, cc_threaded, cell_axes, cell_lengths, VolMapCreate::combo_addframe, VolMapCreate::combo_begin, VolMapCreate::combo_export, VolMapCreateILS::compute, IsoSurface::compute, IsoContour::compute, VolMapCreate::compute_all, VolMapCreateCoulombPotentialMSM::compute_frame, VolMapCreateCoulombPotential::compute_frame, VolMapCreateDistance::compute_frame, VolMapCreateOccupancy::compute_frame, VolMapCreateInterp::compute_frame, VolMapCreateDensity::compute_frame, VolMapCreateMask::compute_frame, VolMapCreate::compute_init, colvarproxy_vmd::compute_voldata, compute_volume_gradient, CreateEmptyGrid, CreateProbGrid, crop, density_info, IsoContour::DoCell, IsoSurface::DoCellGeneral, IsoSurface::DoGridPosNorms, downsample, gaussian_blur, VolumeTexture::generateChargeTexture, VolumeTexture::generateColorScaleTexture, VolumeTexture::generateContourLineTexture, VolumeTexture::generateHSVTexture, gridsize, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, integral, mdff_cc, normalize_pmap, pad, process_pmap, py_mol_get_volumetric, RaycastGrid, segment_volume, Segmentation::Segmentation, IsoSurface::set_color_voltex_rgb3fv, set_volume_gradient, supersample, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, VolIn_CleanGrid, volin_threaded, volin_threaded_prob, volmap_write_dx_file, VolumetricData, voxel_coord, voxel_coord, voxel_coord_from_cartesian_coord, voxel_gradient_fast, voxel_gradient_interpolate_from_coord, voxel_gradient_safe, voxel_index_from_coord, voxel_value, voxel_value_interpolate_from_coord, voxel_value_interpolate_from_coord_safe, voxel_value_safe, and MolFilePlugin::write_volumetric.

double VolumetricData::zaxis[3]
 

direction and length for Z axis (non-unit).

Definition at line 40 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, cell_axes, cell_lengths, IsoContour::compute, VolMapCreate::compute_init, DrawMolecule::cov, CreateEmptyGrid, CreateProbGrid, crop, density_rotate, downsample, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, mdff_cc, normalize_pmap, pad, prepare_texture_coordinates, process_pmap, py_mol_get_volumetric, DrawMolecule::scale_factor, scale_volume, segment_volume, set_volume_axes, text_cmd_mol, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, vol_move, VolumetricData, voxel_coord_from_cartesian_coord, and MolFilePlugin::write_volumetric.

int VolumetricData::zsize
 

number of samples along each axis.

Definition at line 41 of file VolumetricData.h.

Referenced by VolMapCreateILS::add_map_to_molecule, calc_density_bounds_overlapping_map, cc_threaded, cell_axes, cell_lengths, VolMapCreate::combo_addframe, VolMapCreate::combo_begin, VolMapCreate::combo_export, IsoSurface::compute, IsoContour::compute, VolMapCreate::compute_all, VolMapCreateCoulombPotentialMSM::compute_frame, VolMapCreateCoulombPotential::compute_frame, VolMapCreateDistance::compute_frame, VolMapCreateOccupancy::compute_frame, VolMapCreateInterp::compute_frame, VolMapCreateDensity::compute_frame, VolMapCreateMask::compute_frame, VolMapCreate::compute_init, colvarproxy_vmd::compute_voldata, compute_volume_gradient, CreateEmptyGrid, CreateProbGrid, crop, density_info, IsoSurface::DoGridPosNorms, downsample, gaussian_blur, VolumeTexture::generateChargeTexture, VolumeTexture::generateColorScaleTexture, VolumeTexture::generateContourLineTexture, VolumeTexture::generateHSVTexture, gridsize, init_from_identity, init_from_intersection, init_from_union, init_new_volume_molecule, integral, mdff_cc, normalize_pmap, pad, process_pmap, py_mol_get_volumetric, RaycastGrid, segment_volume, Segmentation::Segmentation, IsoSurface::set_color_voltex_rgb3fv, set_volume_gradient, supersample, vmd_cuda_calc_density, vmd_cuda_compare_sel_refmap, vmd_measure_volinterior, vmd_volmap_compare, vmd_volmap_new_fromtype, VolIn_CleanGrid, volin_threaded, volin_threaded_prob, volmap_write_dx_file, VolumetricData, voxel_coord_from_cartesian_coord, voxel_gradient_interpolate_from_coord, voxel_gradient_safe, voxel_index_from_coord, voxel_value_interpolate_from_coord, voxel_value_interpolate_from_coord_safe, voxel_value_safe, and MolFilePlugin::write_volumetric.


The documentation for this class was generated from the following files:
Generated on Fri Apr 19 02:47:18 2024 for VMD (current) by doxygen1.2.14 written by Dimitri van Heesch, © 1997-2002