Bsoft 2.1.4
Bernard's software package
model_create.cpp File Reference

A tool to expand models. More...

#include "rwmodel.h"
#include "model_create.h"
#include "model_transform.h"
#include "model_symmetry.h"
#include "model_util.h"
#include "random_numbers.h"
#include "utilities.h"

Functions

Bcomponentmodel_add_component (Bmodel *model, Bstring &id, Bstring &type, Vector3< double > loc)
 Adds a component to a model. More...
 
Bcomponentmodel_add_components (Bmodel *model, Bstring &id, Bstring &type, vector< Vector3< double > > loc)
 Adds a set of components to a model. More...
 
Bmodelmodel_platonic (Bsymmetry &sym, double radius)
 Generates a platonic polyhedron. More...
 
Bmodelmodel_helix (double radius, double helix_rise, double helix_angle, long ncomp)
 Generates a helix model. More...
 
Bmodelmodel_random (long ncomp, double comp_radius, double max_radius)
 Generates components at random non-overlapping locations and with random views. More...
 
Bmodelmodel_random (long ncomp, double comp_radius, Vector3< double > min, Vector3< double > max)
 Generates components at random non-overlapping locations and with random views. More...
 
Bmodelmodel_random_gaussian (long ncomp, double std)
 Generates components at random gaussian-distributed locations and with random views. More...
 
Bmodelmodel_random_shell (long ncomp, double radius)
 Generates random components on a shell. More...
 
Bmodelmodel_random_shell (long ncomp, double radius, double separation)
 Generates random components on a shell. More...
 
Bmodelmodel_create_shell (long number, double radius, double separation)
 Creates a spherical shell point model. More...
 
Bmodelmodel_create_fibonacci_sphere (long number, double radius)
 Creates a Fibonacci sphere point model. More...
 
long model_subdivide (Bmodel *model, long division)
 
long model_spherize (Bmodel *model, double sphere_fraction)
 
Bmodelmodel_create_dodecahedron (double radius, long divisions, double sphere_fraction)
 Creates an dodecahedron point model. More...
 
Bmodelmodel_create_icosahedron (double radius, long divisions, double sphere_fraction)
 Creates an icosahedron point model. More...
 
double ellipsoid_slice_distance_error (Vector3< double > axes, double z1, double z2, double distance)
 Calculates how far the slice separation is from ideal. More...
 
double ellipsoid_next_z (Vector3< double > axes, double z1, double distance)
 Calculates the next slice's z component to evenly space slices. More...
 
double ellipse_next_x (Vector3< double > axes, double x1, double distance)
 Calculates the next slice's x component to evenly space points. More...
 
Bmodelmodel_create_circle (double radius, double z, double distance)
 Creates a circle with the given radius, height, and component separation. More...
 
Bmodelmodel_create_ellipse (Vector3< double > axes, double distance)
 Creates an ellipse with given semiaxes and component separation. More...
 
Bmodelmodel_create_ellipsoid (Vector3< double > axes, double distance)
 Creates an ellipsoid with given semiaxes and component separation. More...
 
Bmodelmodel_create_cylinder (Vector3< double > direction, double radius, double length, double distance)
 Creates a cylinder with the given radius, length, and component separation. More...
 
Bmodelmodel_create_plane (double length, double width, double z, double distance)
 Creates a plane with the given length, width, height, and component separation. More...
 
Bmodelmodel_create_spindle (Vector3< double > direction, double radius, double separation, double packing)
 Creates a spindle packed inside a sphere. More...
 
Bmodelmodel_create_cubic_lattice (Vector3< long > lattice, double separation)
 Creates a cubic lattice with the given length, width, height, and component separation. More...
 
Bmodelmodel_create_hexagonal_lattice (Vector3< long > lattice, double separation)
 Creates a hexagonal lattice with the given length, width, height, and component separation. More...
 

Variables

int verbose
 

Detailed Description

A tool to expand models.

Author
Bernard Heymann
Date
Created: 20090714
Modified: 20220210

Function Documentation

◆ ellipse_next_x()

double ellipse_next_x ( Vector3< double >  axes,
double  x1,
double  distance 
)

Calculates the next slice's x component to evenly space points.

Author
Jonathan Luo
Parameters
axesthe ellipse axis lengths and z-height
x1???
distanceseparation distance
Returns
double next x-coordinate
Solves the 4th order polynomial by bracketing.

◆ ellipsoid_next_z()

double ellipsoid_next_z ( Vector3< double >  axes,
double  z1,
double  distance 
)

Calculates the next slice's z component to evenly space slices.

Author
Jonathan Luo
Parameters
axesellipsoid axis lengths
z1previous slice's z
distanceseparation distance
Returns
double next slice's z-coordinate

◆ ellipsoid_slice_distance_error()

double ellipsoid_slice_distance_error ( Vector3< double >  axes,
double  z1,
double  z2,
double  distance 
)

Calculates how far the slice separation is from ideal.

Author
Jonathan Luo
Parameters
axesellipsoid axis lengths
z1previous slice z value
z2next slice z value
distancetarget separation distance
Returns
double the distance error (2.0*distance is target sum for the 2 distances)
For the generalized ellipsoid, the distance between points of 2 slices varies
if the x and y axes are not equal. The best slice will have a component some
distance+error separated at the shorter axis and distance-error at the longer axis,
distributing the error. Solve the ellipsoid equation for those two points and distances. 

◆ model_add_component()

Bcomponent * model_add_component ( Bmodel model,
Bstring id,
Bstring type,
Vector3< double >  loc 
)

Adds a component to a model.

Parameters
*modelmodel list.
&idmodel identifier.
&typecomponent type.
loccomponent location.
Returns
Bcomponent* new component.
The component is added to the model indicated by the ID, the type and
the location given.

◆ model_add_components()

Bcomponent * model_add_components ( Bmodel model,
Bstring id,
Bstring type,
vector< Vector3< double > >  loc 
)

Adds a set of components to a model.

Parameters
*modelmodel list.
&idmodel identifier.
&typecomponent type.
loccomponent locations.
Returns
Bcomponent* last new component.
The component is added to the model indicated by the ID, the type and
the location given.

◆ model_create_circle()

Bmodel * model_create_circle ( double  radius,
double  z,
double  distance 
)

Creates a circle with the given radius, height, and component separation.

Author
Jonathan Luo
Parameters
radiusradius length.
zheight.
distanceseparation distance.
Returns
Bmodel* new circle model.
Find angle associated with separation distance using Law of Cosines.

◆ model_create_cubic_lattice()

Bmodel * model_create_cubic_lattice ( Vector3< long >  lattice,
double  separation 
)

Creates a cubic lattice with the given length, width, height, and component separation.

Parameters
latticelattice dimensions.
separationdistances between components.
Returns
Bmodel* new lattice model.
The components are placed in a cubic arrangement.

◆ model_create_cylinder()

Bmodel * model_create_cylinder ( Vector3< double >  direction,
double  radius,
double  length,
double  distance 
)

Creates a cylinder with the given radius, length, and component separation.

Author
Jonathan Luo
Parameters
directioncylinder axis direction
radiusradius length.
lengthcylinder length.
distanceseparation distance.
Returns
Bmodel* new cylinder model.
Find angle associated with separation distance using Law of Cosines.

◆ model_create_dodecahedron()

Bmodel * model_create_dodecahedron ( double  radius,
long  divisions,
double  sphere_fraction 
)

Creates an dodecahedron point model.

Parameters
radiussphere radius.
divisionsnumber of divisions from a base dodecahedron.
sphere_fractionspherical fraction: 0=dodecahedral, 1=spherical.
Returns
Bmodel* new sphere model.

◆ model_create_ellipse()

Bmodel * model_create_ellipse ( Vector3< double >  axes,
double  distance 
)

Creates an ellipse with given semiaxes and component separation.

Author
Jonathan Luo
Parameters
axesellipse axis x and y lengths and z height
distanceseparation distance
Returns
Bmodel* new ellipse model.

◆ model_create_ellipsoid()

Bmodel * model_create_ellipsoid ( Vector3< double >  axes,
double  distance 
)

Creates an ellipsoid with given semiaxes and component separation.

Author
Jonathan Luo
Parameters
axesellipsoid axis lengths
distanceseparation distance
Returns
Bmodel* new ellipsoid model.

◆ model_create_fibonacci_sphere()

Bmodel * model_create_fibonacci_sphere ( long  number,
double  radius 
)

Creates a Fibonacci sphere point model.

Parameters
numbernumber of points.
radiussphere radius.
Returns
Bmodel* new sphere model.

◆ model_create_hexagonal_lattice()

Bmodel * model_create_hexagonal_lattice ( Vector3< long >  lattice,
double  separation 
)

Creates a hexagonal lattice with the given length, width, height, and component separation.

Parameters
latticelattice dimensions.
separationdistances between components.
Returns
Bmodel* new lattice model.
The components are placed in a hexagonal arrangement to simulate the
most commonly expected close packing.

◆ model_create_icosahedron()

Bmodel * model_create_icosahedron ( double  radius,
long  divisions,
double  sphere_fraction 
)

Creates an icosahedron point model.

Parameters
radiussphere radius.
divisionsnumber of divisions from a base icosahedron.
sphere_fractionspherical fraction: 0=dodecahedral, 1=spherical.
Returns
Bmodel* new sphere model.

◆ model_create_plane()

Bmodel * model_create_plane ( double  length,
double  width,
double  z,
double  distance 
)

Creates a plane with the given length, width, height, and component separation.

Parameters
lengthdimension in x.
widthdimension in y.
zheight.
distanceseparation distance.
Returns
Bmodel* new plane model.
The components are placed in a hexagonal arrangement to simulate the
most commonly expected close packingx.

◆ model_create_shell()

Bmodel * model_create_shell ( long  number,
double  radius,
double  separation 
)

Creates a spherical shell point model.

Parameters
numbernumber of points.
radiusshell radius.
separationdistance between points.
Returns
Bmodel* new shell model.
Two of the three arguments need to be given (the other zero).
If all three arguments are given, only the radius and distance is used.

◆ model_create_spindle()

Bmodel * model_create_spindle ( Vector3< double >  direction,
double  radius,
double  separation,
double  packing 
)

Creates a spindle packed inside a sphere.

Parameters
directionspindle axis direction
radiussphere radius.
separationseparation distance between successive components.
packingdistance between spindle strands.
Returns
Bmodel* new spindle model.

The spindle starts at the minimum Z-pole of the sphere, winding up along the spherical wall to the other pole, and back down inside the first shell. The process is repeated until the whole sphere is packed.

◆ model_helix()

Bmodel * model_helix ( double  radius,
double  helix_rise,
double  helix_angle,
long  ncomp 
)

Generates a helix model.

Parameters
radiusdistance of each vertex from helical axis.
helix_risehelical rise.
helix_anglehelical rotation angle.
ncompnumber of components.
Returns
Bmodel* new model.

◆ model_platonic()

Bmodel * model_platonic ( Bsymmetry sym,
double  radius 
)

Generates a platonic polyhedron.

Parameters
*symsymmetry.
radiusdistance of each vertex from origin.
Returns
Bmodel* new model.
Polyhedral components are generated on symmetry axes.
Symmetry designations supported:
    T-3     tetrahedron.
    O-2     truncated octahedron.
    O-3     cube.
    O-4     octahedron.
    I-2     truncated icosahedron.
    I-3     dodecahedron.
    I-5     icosahedron.

◆ model_random() [1/2]

Bmodel * model_random ( long  ncomp,
double  comp_radius,
double  max_radius 
)

Generates components at random non-overlapping locations and with random views.

Parameters
ncompnumber of components.
comp_radiuscomponent radius.
max_radiusmaximum radius of components.
Returns
Bmodel* new model.
If a new component overlaps within an existing component, as defined
by the component radius, new random coordinates are generated for it.

◆ model_random() [2/2]

Bmodel * model_random ( long  ncomp,
double  comp_radius,
Vector3< double >  min,
Vector3< double >  max 
)

Generates components at random non-overlapping locations and with random views.

Parameters
ncompnumber of components.
comp_radiuscomponent radius.
minminimum bounds.
maxmaximum bounds.
Returns
Bmodel* new model.
If a new component overlaps within an existing component, as defined
by the component radius, new random coordinates are generated for it.

◆ model_random_gaussian()

Bmodel * model_random_gaussian ( long  ncomp,
double  std 
)

Generates components at random gaussian-distributed locations and with random views.

Parameters
ncompnumber of components.
stdstandard deviation of gaussian distribution.
Returns
Bmodel* new model.
Overlapping components are generated.

◆ model_random_shell() [1/2]

Bmodel * model_random_shell ( long  ncomp,
double  radius 
)

Generates random components on a shell.

Parameters
ncompnumber of components.
radiusshell radius.
Returns
Bmodel* new model.
Overlapping components are generated.

◆ model_random_shell() [2/2]

Bmodel * model_random_shell ( long  ncomp,
double  radius,
double  separation 
)

Generates random components on a shell.

Parameters
ncompnumber of components.
radiusshell radius.
separationminimum vertex separation distance.
Returns
Bmodel* new model.
The vertices are created with a desired separation, which means that
if too many vertices are requested with too large separation, the 
function may never return. The following must therefore be true:
    n*(d/r)^2 < 8
where
    n: number of vertices
    d: minimum separation distance
    r: polyhedron radius
If this not true, the number of vertices is decreased.

◆ model_spherize()

long model_spherize ( Bmodel model,
double  sphere_fraction 
)

◆ model_subdivide()

long model_subdivide ( Bmodel model,
long  division 
)

Variable Documentation

◆ verbose

int verbose
extern