Xmipp  v3.23.11-Nereus
Public Member Functions | Public Attributes | Friends | List of all members
Ellipsoid Class Reference

#include <phantom.h>

Inheritance diagram for Ellipsoid:
Inheritance graph
[legend]
Collaboration diagram for Ellipsoid:
Collaboration graph
[legend]

Public Member Functions

void prepare ()
 
void assign (const Ellipsoid &F)
 
int point_inside (const Matrix1D< double > &r, Matrix1D< double > &aux) const
 
double density_inside (const Matrix1D< double > &r, Matrix1D< double > &aux) const
 
Featurescale (double factor) const
 
 __attribute__ ((deprecated)) void scale(double factor
 
double intersection (const Matrix1D< double > &direction, const Matrix1D< double > &passing_point, Matrix1D< double > &r, Matrix1D< double > &u) const
 
double volume () const
 
void read_specific (char *line)
 
void read_specific (const std::vector< double > &vect)
 
void feat_printf (FILE *fh) const
 
void feat_printm (MetaData &MD, size_t id)
 
void init_rnd (double minXradius, double maxXradius, double minYradius, double maxYradius, double minZradius, double maxZradius, double minden=1.0, double maxden=1.0, double minx0=0, double maxx0=0, double miny0=0, double maxy0=0, double minz0=0, double maxz0=0, double minrot=0, double maxrot=360, double mintilt=0, double maxtilt=180, double minpsi=0, double maxpsi=360)
 
- Public Member Functions inherited from Oriented_Feature
void prepare_Euler ()
 
void assign (const Oriented_Feature &OF)
 
virtual void rotate (const Matrix2D< double > &E)
 
- Public Member Functions inherited from Feature
void assign (const Feature &F)
 
virtual void rotate_center (const Matrix2D< double > &E)
 
int point_inside (const Matrix1D< double > &r) const
 
int voxel_inside (const Matrix1D< double > &r, Matrix1D< double > &aux1, Matrix1D< double > &aux2) const
 
int voxel_inside (const Matrix1D< double > &r) const
 
double voxel_inside_by_normalized_density (const Matrix1D< double > &r, Matrix1D< double > &aux1, Matrix1D< double > &aux2) const
 
int intersects_sphere (const Matrix1D< double > &r, double radius, Matrix1D< double > &aux1, Matrix1D< double > &aux2, Matrix1D< double > &aux3) const
 
int intersects_sphere (const Matrix1D< double > &r, double radius) const
 
Featureencircle (double radius=0) const
 
Featurebackground (int back_mode, double back_param) const
 
double intersection (const Matrix1D< double > &direction, const Matrix1D< double > &passing_point) const
 
void mean_variance_in_plane (Image< double > *V, double z, double &mean, double &var)
 
void project_to (Projection &P, const Matrix2D< double > &VP, const Matrix2D< double > &PV) const
 
void corners (const MultidimArray< double > &V, Matrix1D< double > &corner1, Matrix1D< double > &corner2)
 
void draw_in (MultidimArray< double > &V, int color_mode=INTERNAL, double colour=-1)
 
void sketch_in (MultidimArray< double > &V, double colour=2)
 
void shift (double shiftX, double shiftY, double shiftZ)
 
void selfApplyGeometry (const Matrix2D< double > &A)
 
void readCommon (char *line)
 
void readCommon (MDRow &row)
 
void read (MDRow &row)
 

Public Attributes

double xradius
 X radius before rotating. More...
 
double yradius
 Y radius before rotating. More...
 
double zradius
 Z radius before rotating. More...
 
Feature **_f const
 
- Public Attributes inherited from Oriented_Feature
double rot
 First Euler angle. More...
 
double tilt
 Second Euler angle. More...
 
double psi
 Third Euler angle. More...
 
Matrix2D< double > euler
 Euler matrix. More...
 
Matrix2D< double > eulert
 Inverse Euler matrix. More...
 
- Public Attributes inherited from Feature
std::string type
 
char add_assign
 
double density
 
Matrix1D< double > center
 
double max_distance
 

Friends

std::ostream & operator<< (std::ostream &o, const Ellipsoid &f)
 

Detailed Description

Ellipsoid. An elliposid is defined by its center, and the radius in the 3 dimensions (X, Y and Z, before rotating), label="ell". The center of the ellipsoid is at the geometrical center, and the dimensions define the elliposid from "-Xradius/2" to "Xradius/2", ... around the center. Then the ellipsoid is rotated after the three Euler angles. The corresponding Euler matrix is stored in the field "euler" while its inverse is in "eulert".

A point (r) in the universal coordinate system, where elliposids are defined in general, can be expressed (rp) with respect to a system where the ellipsoid is centered at the origin, all its radii are unity and its axes are aligned with XYZ by

V3_MINUS_V3(rp,r,ell.Center);
M3x3_BY_V3x1(rp,ell.euler,rp);
XX(rp) /= ell.xradius;
YY(rp) /= ell.yradius;
ZZ(rp) /= ell.zradius;

Directions (free vectors) are transformed in the same fashion except that for the first vector subtraction (referring the origin).

Definition at line 1127 of file phantom.h.

Member Function Documentation

◆ __attribute__()

Ellipsoid::__attribute__ ( (deprecated)  )

Another function for return a scaled elliposoid.

◆ assign()

void Ellipsoid::assign ( const Ellipsoid F)

Another function for assignment.

Definition at line 138 of file phantom.cpp.

139 {
140  *this = F;
141 }

◆ density_inside()

double Ellipsoid::density_inside ( const Matrix1D< double > &  r,
Matrix1D< double > &  aux 
) const
inlinevirtual

Density inside a Ellipsoid. This function tells you the density of the ellipsoid at point r for constant valued features is trivial

Implements Feature.

Definition at line 1156 of file phantom.h.

1157  {
1158  return 1.;
1159  }

◆ feat_printf()

void Ellipsoid::feat_printf ( FILE *  fh) const
virtual

Print ellipsoid in the standard feature format. \ See {Feature::feat_printf}, Phantoms

Implements Feature.

Definition at line 523 of file phantom.cpp.

524 {
525  fprintf(fh, "ell %c %1.4f % 7.2f % 7.2f % 7.2f % 7.2f "
526  "% 7.2f % 7.2f % 7.2f % 7.2f % 7.2f\n",
529  rot, tilt, psi);
530 }
double tilt
Second Euler angle.
Definition: phantom.h:430
double xradius
X radius before rotating.
Definition: phantom.h:1131
Matrix1D< double > center
Definition: phantom.h:119
double density
Definition: phantom.h:114
char add_assign
Definition: phantom.h:108
#define XX(v)
Definition: matrix1d.h:85
#define YY(v)
Definition: matrix1d.h:93
double rot
First Euler angle.
Definition: phantom.h:427
double psi
Third Euler angle.
Definition: phantom.h:433
fprintf(glob_prnt.io, "\)
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134
#define ZZ(v)
Definition: matrix1d.h:101

◆ feat_printm()

void Ellipsoid::feat_printm ( MetaData MD,
size_t  id 
)
virtual

Implements Feature.

Definition at line 533 of file phantom.cpp.

534 {
535  std::vector<double> FSVect;
536  FSVect.push_back(xradius);
537  FSVect.push_back(yradius);
538  FSVect.push_back(zradius);
539  FSVect.push_back(rot);
540  FSVect.push_back(tilt);
541  FSVect.push_back(psi);
543 }
double tilt
Second Euler angle.
Definition: phantom.h:430
double xradius
X radius before rotating.
Definition: phantom.h:1131
Specific parameters for a feature (vector double)
bool setValue(const MDLabel label, const T &valueIn, size_t id)
double rot
First Euler angle.
Definition: phantom.h:427
double psi
Third Euler angle.
Definition: phantom.h:433
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134

◆ init_rnd()

void Ellipsoid::init_rnd ( double  minXradius,
double  maxXradius,
double  minYradius,
double  maxYradius,
double  minZradius,
double  maxZradius,
double  minden = 1.0,
double  maxden = 1.0,
double  minx0 = 0,
double  maxx0 = 0,
double  miny0 = 0,
double  maxy0 = 0,
double  minz0 = 0,
double  maxz0 = 0,
double  minrot = 0,
double  maxrot = 360,
double  mintilt = 0,
double  maxtilt = 180,
double  minpsi = 0,
double  maxpsi = 360 
)

Generate a random ellipsoid. An ellipsoid is generated randomly within the range of the parameters given. Notice that most of them are set by default. The exact parameters are picked uniformly from the range. The maximum distance is adjusted properly according to the randomly generated feature. 'den' stands for density and (x0,y0,z0) is the center of the feature. The behaviour of the new sphere is always Add.

Definition at line 1854 of file phantom.cpp.

1865 {
1867  center.resize(3);
1868  type = "ell";
1869  add_assign = '+';
1870  density = rnd_unif(minden, maxden);
1871  XX(center) = rnd_unif(minx0, maxx0);
1872  YY(center) = rnd_unif(miny0, maxy0);
1873  ZZ(center) = rnd_unif(minz0, maxz0);
1874 
1875  if (minYradius == 0)
1876  minYradius = minXradius;
1877  if (minZradius == 0)
1878  minZradius = minXradius;
1879  if (maxYradius == 0)
1880  maxYradius = maxXradius;
1881  if (maxZradius == 0)
1882  maxZradius = maxXradius;
1883  xradius = rnd_unif(minXradius, maxXradius);
1884  yradius = rnd_unif(minYradius, maxYradius);
1885  zradius = rnd_unif(minZradius, maxZradius);
1886  rot = rnd_unif(minrot, maxrot);
1887  tilt = rnd_unif(mintilt, maxtilt);
1888  psi = rnd_unif(minpsi, maxpsi);
1890  eulert = euler.transpose();
1891 
1893 }
#define XMIPP_MAX(x, y)
Definition: xmipp_macros.h:193
Matrix2D< double > eulert
Inverse Euler matrix.
Definition: phantom.h:439
void Euler_angles2matrix(T alpha, T beta, T gamma, Matrix2D< T > &A, bool homogeneous)
Definition: geometry.cpp:624
double tilt
Second Euler angle.
Definition: phantom.h:430
double xradius
X radius before rotating.
Definition: phantom.h:1131
Matrix2D< double > euler
Euler matrix.
Definition: phantom.h:436
Matrix1D< double > center
Definition: phantom.h:119
double density
Definition: phantom.h:114
Matrix2D< T > transpose() const
Definition: matrix2d.cpp:1314
char add_assign
Definition: phantom.h:108
double rnd_unif()
#define XX(v)
Definition: matrix1d.h:85
void resize(size_t Xdim, bool copy=true)
Definition: matrix1d.h:410
double max_distance
Definition: phantom.h:126
#define YY(v)
Definition: matrix1d.h:93
std::string type
Definition: phantom.h:101
double rot
First Euler angle.
Definition: phantom.h:427
double psi
Third Euler angle.
Definition: phantom.h:433
unsigned int randomize_random_generator()
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134
#define ZZ(v)
Definition: matrix1d.h:101

◆ intersection()

double Ellipsoid::intersection ( const Matrix1D< double > &  direction,
const Matrix1D< double > &  passing_point,
Matrix1D< double > &  r,
Matrix1D< double > &  u 
) const
virtual

Intersection of a ray with an ellipsoid. See Feature::intersection to know more about the parameters meaning.

Implements Feature.

Definition at line 1250 of file phantom.cpp.

1255 {
1256  double norm = direction.module();
1258 
1259  // Set the passing point in the ellipsoid coordinate system
1260  // and normalise to a unit sphere
1261  V3_MINUS_V3(r, passing_point, center);
1262  M3x3_BY_V3x1(r, euler, r);
1263  XX(r) /= xradius;
1264  YY(r) /= yradius;
1265  ZZ(r) /= zradius;
1266 
1267  // Express also the direction in the ellipsoid coordinate system
1268  // and normalise to a unit sphere
1269  M3x3_BY_V3x1(u, euler, direction);
1270  XX(u) /= xradius;
1271  YY(u) /= yradius;
1272  ZZ(u) /= zradius;
1273 
1274  return intersection_unit_sphere(u, r) / norm;
1275 }
double module() const
Definition: matrix1d.h:983
double xradius
X radius before rotating.
Definition: phantom.h:1131
Matrix2D< double > euler
Euler matrix.
Definition: phantom.h:436
#define V3_MINUS_V3(a, b, c)
Definition: matrix1d.h:202
Matrix1D< double > center
Definition: phantom.h:119
T norm(const std::vector< T > &v)
Definition: vector_ops.h:399
#define XX(v)
Definition: matrix1d.h:85
double intersection_unit_sphere(const Matrix1D< double > &u, const Matrix1D< double > &r)
Definition: geometry.cpp:1122
#define M3x3_BY_V3x1(a, M, b)
Definition: matrix2d.h:170
#define YY(v)
Definition: matrix1d.h:93
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134
#define ZZ(v)
Definition: matrix1d.h:101
#define SPEED_UP_temps012
Definition: xmipp_macros.h:403

◆ point_inside()

int Ellipsoid::point_inside ( const Matrix1D< double > &  r,
Matrix1D< double > &  aux 
) const
virtual

Speeded up point inside an ellipsoid. This function tells you if a point is inside the elliposoid or not. See Feature::point_inside

Implements Feature.

Definition at line 794 of file phantom.cpp.

795 {
797  double tx;
798  double ty;
799  double tz;
800 
801  // Express r in the feature coord. system
802  V3_MINUS_V3(aux, r, center);
803  M3x3_BY_V3x1(aux, euler, aux);
804 
805  // Check if inside
806  tx = XX(aux) / xradius;
807  ty = YY(aux) / yradius;
808  tz = ZZ(aux) / zradius;
809  if (tx*tx + ty*ty + tz*tz <= 1.0)
810  return 1;
811  return 0;
812 }
double xradius
X radius before rotating.
Definition: phantom.h:1131
Matrix2D< double > euler
Euler matrix.
Definition: phantom.h:436
#define V3_MINUS_V3(a, b, c)
Definition: matrix1d.h:202
Matrix1D< double > center
Definition: phantom.h:119
#define XX(v)
Definition: matrix1d.h:85
#define M3x3_BY_V3x1(a, M, b)
Definition: matrix2d.h:170
#define YY(v)
Definition: matrix1d.h:93
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134
#define ZZ(v)
Definition: matrix1d.h:101
#define SPEED_UP_temps012
Definition: xmipp_macros.h:403

◆ prepare()

void Ellipsoid::prepare ( )
virtual

Prepare ellipsoid for work. Computes the maximum distance, in this case, equal to "MAX(MAX(xradius,yradius),zradius)", and computes the Euler and inverse Euler matrices as a function of the Euler angles

Implements Feature.

Definition at line 77 of file phantom.cpp.

78 {
79  prepare_Euler();
81 }
#define XMIPP_MAX(x, y)
Definition: xmipp_macros.h:193
double xradius
X radius before rotating.
Definition: phantom.h:1131
double max_distance
Definition: phantom.h:126
void prepare_Euler()
Definition: phantom.cpp:102
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134

◆ read_specific() [1/2]

void Ellipsoid::read_specific ( char *  line)

Read specific description for an ellipsoid. An exception is thrown if the line doesn't conform the standard specification. See Feature::read_specific

Definition at line 357 of file phantom.cpp.

358 {
359  int stat;
360  stat = sscanf(line, "%*s %*c %*f %*f %*f %*f %lf %lf %lf %lf %lf %lf",
361  &xradius, &yradius, &zradius, &rot, &tilt, &psi);
362  if (stat != 6)
363  REPORT_ERROR(ERR_IO_NOREAD, (std::string)"Error when reading an ellipsoid" + line);
364  prepare();
365 }
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
double tilt
Second Euler angle.
Definition: phantom.h:430
double xradius
X radius before rotating.
Definition: phantom.h:1131
void prepare()
Definition: phantom.cpp:77
Couldn&#39;t read from file.
Definition: xmipp_error.h:139
double rot
First Euler angle.
Definition: phantom.h:427
double psi
Third Euler angle.
Definition: phantom.h:433
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134

◆ read_specific() [2/2]

void Ellipsoid::read_specific ( const std::vector< double > &  vector)
virtual

Read a feature from a file, VIRTUAL!!!. The format must be the one given in Phantoms, and each subclass must implement its own I/O routines. These routines must fill only the non common part of the feature description, but they receive the whole line with the description.

Implements Feature.

Definition at line 368 of file phantom.cpp.

369 {
370  if (vect.size() != 6)
371  REPORT_ERROR(ERR_IO_NOREAD, (std::string)"Error when reading a ellipsoid");
372  xradius = vect[0];
373  yradius = vect[1];
374  zradius = vect[2];
375  rot = vect[3];
376  tilt = vect[4];
377  psi = vect[5];
378  prepare();
379 }
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
double tilt
Second Euler angle.
Definition: phantom.h:430
double xradius
X radius before rotating.
Definition: phantom.h:1131
void prepare()
Definition: phantom.cpp:77
Couldn&#39;t read from file.
Definition: xmipp_error.h:139
double rot
First Euler angle.
Definition: phantom.h:427
double psi
Third Euler angle.
Definition: phantom.h:433
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134

◆ scale()

void Ellipsoid::scale ( double  factor) const
virtual

Return a scaled elliposoid. The center, density, angles and behaviour of the new ellipsoid is exactly the same as the actual one. The dimensions are multiplied by the scale factor and the maximum distance is recalculated for the new ellipsoid.

Implements Feature.

Definition at line 1599 of file phantom.cpp.

1600 {
1601  Ellipsoid *f;
1602  f = new Ellipsoid;
1604  COPY_ANGLES;
1605 
1606  f->xradius = factor * xradius;
1607  f->yradius = factor * yradius;
1608  f->zradius = factor * zradius;
1609  f->prepare();
1610  return (Feature *)f;
1611 }
double xradius
X radius before rotating.
Definition: phantom.h:1131
void prepare()
Definition: phantom.cpp:77
#define COPY_COMMON_PART
Definition: phantom.cpp:1473
double * f
#define COPY_ANGLES
Definition: phantom.cpp:1479
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134

◆ volume()

double Ellipsoid::volume ( ) const
inlinevirtual

Volume of an ellipsoid. This function returns 4/3*PI*xradius*yradius*zradius. See Feature::volume

Implements Feature.

Definition at line 1180 of file phantom.h.

1181  {
1182  return 4 / 3*PI*xradius*yradius*zradius;
1183  }
double xradius
X radius before rotating.
Definition: phantom.h:1131
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134
#define PI
Definition: tools.h:43

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  o,
const Ellipsoid f 
)
friend

Show feature not in the standard format but more informatively. This function only shows the non common part of the feature. Use the << operator of Feature to show the whole feature.

Definition at line 658 of file phantom.cpp.

659 {
660  o << " Xradius: " << f.xradius << std::endl;
661  o << " Yradius: " << f.yradius << std::endl;
662  o << " Zradius: " << f.zradius << std::endl;
663  o << " Rot: " << f.rot << std::endl;
664  o << " Tilt: " << f.tilt << std::endl;
665  o << " Psi: " << f.psi << std::endl;
666  return o;
667 }
double tilt
Second Euler angle.
Definition: phantom.h:430
double xradius
X radius before rotating.
Definition: phantom.h:1131
double rot
First Euler angle.
Definition: phantom.h:427
double psi
Third Euler angle.
Definition: phantom.h:433
double zradius
Z radius before rotating.
Definition: phantom.h:1137
double yradius
Y radius before rotating.
Definition: phantom.h:1134

Member Data Documentation

◆ const

Feature** _f Ellipsoid::const

Definition at line 1168 of file phantom.h.

◆ xradius

double Ellipsoid::xradius

X radius before rotating.

Definition at line 1131 of file phantom.h.

◆ yradius

double Ellipsoid::yradius

Y radius before rotating.

Definition at line 1134 of file phantom.h.

◆ zradius

double Ellipsoid::zradius

Z radius before rotating.

Definition at line 1137 of file phantom.h.


The documentation for this class was generated from the following files: