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

#include <basis.h>

Collaboration diagram for Basis:
Collaboration graph
[legend]

Public Types

enum  tBasisFunction { blobs, voxels, splines }
 Type of basis function. More...
 

Public Member Functions

 Basis ()
 Empty constructor. By default, blobs. More...
 
void setDefault ()
 Default values. More...
 
String basisName () const
 Basis name. More...
 
void read (int argc, char **argv)
 
void read (const FileName &fn)
 
void readParams (XmippProgram *program)
 
void produceSideInfo (const Grid &grid)
 
void setSamplingRate (double _Tm)
 
void setD (Matrix2D< double > *_D)
 
double maxLength () const
 
void changeToVoxels (GridVolume &vol_basis, MultidimArray< double > *vol_voxels, int Zdim, int Ydim, int Xdim, int threads=1) const
 
void changeFromVoxels (const MultidimArray< double > &vol_voxels, GridVolume &vol_basis, int grid_type, double grid_relative_size, const MultidimArray< double > *vol_mask, const Matrix2D< double > *D, double R, int threads=1) const
 
double valueAt (const Matrix1D< double > &r) const
 
double projectionAt (const Matrix1D< double > &u, const Matrix1D< double > &r) const
 

Static Public Member Functions

static void defineParams (XmippProgram *program, const char *prefix=NULL, const char *comment=NULL)
 

Public Attributes

tBasisFunction type
 Basis function to use. More...
 
MultidimArray< double > * VolPSF
 Footprint is convolved with a volume PSF // At this moment only used with blobs. More...
 
double Tm
 Sampling rate. More...
 
struct blobtype blob
 Blob parameters. More...
 
double grid_relative_size
 Relative size for the grid. More...
 
Matrix2D< double > * D
 
ImageOver blobprint
 Blob footprint. More...
 
ImageOver blobprint2
 Square of the footprint. More...
 
double sum_on_grid
 Sum of the basis on the grid points. More...
 
Matrix1D< double > aux
 Auxiliary vector for projections. More...
 

Friends

std::ostream & operator<< (std::ostream &out, const Basis &basis)
 Show. More...
 

Detailed Description

Basis class.

Definition at line 45 of file basis.h.

Member Enumeration Documentation

◆ tBasisFunction

Type of basis function.

Enumerator
blobs 
voxels 
splines 

Definition at line 49 of file basis.h.

tBasisFunction
Type of basis function.
Definition: basis.h:49

Constructor & Destructor Documentation

◆ Basis()

Basis::Basis ( )

Empty constructor. By default, blobs.

Definition at line 30 of file basis.cpp.

31 {
32  setDefault();
33 }
void setDefault()
Default values.
Definition: basis.cpp:36

Member Function Documentation

◆ basisName()

String Basis::basisName ( ) const

Basis name.

Definition at line 50 of file basis.cpp.

51 {
52  switch (type)
53  {
54  case Basis::blobs:
55  return "blobs";
56  break;
57  case Basis::voxels:
58  return "voxels";
59  break;
60  case Basis::splines:
61  return "splines";
62  break;
63  default:
64  return "Unknown";
65  break;
66  }
67 }
tBasisFunction type
Basis function to use.
Definition: basis.h:52

◆ changeFromVoxels()

void Basis::changeFromVoxels ( const MultidimArray< double > &  vol_voxels,
GridVolume vol_basis,
int  grid_type,
double  grid_relative_size,
const MultidimArray< double > *  vol_mask,
const Matrix2D< double > *  D,
double  R,
int  threads = 1 
) const

Change basis from voxels. A voxel volume is provided, then the output vol_basis will be shaped to represent this voxel volume. If the volume is already in voxels, then the mask and the radius mask are applied.

Definition at line 292 of file basis.cpp.

296 {
297  Grid grid;
298  Matrix1D<double> corner1(3);
299  Matrix1D<double> corner2(3);
300  double R2 = R * R;
301  switch (type)
302  {
303  case blobs:
304  voxels2blobs(&vol_voxels, blob, vol_basis,
305  grid_type, grid_relative_size, 0.05, vol_mask,
306  D, 0.01, 0, R, threads);
307  break;
308  case voxels:
309  VECTOR_R3(corner1, STARTINGX(vol_voxels),
310  STARTINGY(vol_voxels), STARTINGZ(vol_voxels));
311  VECTOR_R3(corner2, FINISHINGX(vol_voxels),
312  FINISHINGY(vol_voxels), FINISHINGZ(vol_voxels));
313  grid = Create_CC_grid(1, corner1, corner2);
314  vol_basis.adapt_to_grid(grid);
315  vol_basis(0)() = vol_voxels;
316  if (R != -1)
317  FOR_ALL_ELEMENTS_IN_ARRAY3D(vol_voxels)
318  if (k*k + i*i + j*j > R2)
319  vol_basis(0)()(k, i, j) = 0;
320  if (vol_mask != nullptr)
321  FOR_ALL_ELEMENTS_IN_ARRAY3D(vol_voxels)
322  if ((*vol_mask)(k, i, j) == 0)
323  vol_basis(0)()(k, i, j) = 0;
324  break;
325  case splines:
327  /* TODO */
328  break;
329  }
330 }
Case or algorithm not implemented yet.
Definition: xmipp_error.h:177
#define FINISHINGX(v)
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
SimpleGrid Create_CC_grid(double relative_size, const Matrix1D< double > &corner1, const Matrix1D< double > &corner2, const Matrix1D< double > &origin)
Definition: grids.cpp:196
double grid_relative_size
Relative size for the grid.
Definition: basis.h:64
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
Definition: grids.h:479
#define FINISHINGZ(v)
#define STARTINGX(v)
#define i
ql0001_ & k(htemp+1),(cvec+1),(atemp+1),(bj+1),(bl+1),(bu+1),(x+1),(clamda+1), &iout, infoqp, &zero,(w+1), &lenw,(iw+1), &leniw, &glob_grd.epsmac
#define STARTINGY(v)
void adapt_to_grid(const Grid &_grid)
Definition: grids.h:838
#define FOR_ALL_ELEMENTS_IN_ARRAY3D(V)
void voxels2blobs(const MultidimArray< double > *vol_voxels, const struct blobtype &blob, GridVolume &vol_blobs, int grid_type, double grid_relative_size, double lambda, const MultidimArray< double > *vol_mask, const Matrix2D< double > *D, double final_error_change, int tell, double R, int threads)
Definition: blobs.cpp:1316
#define j
#define FINISHINGY(v)
#define VECTOR_R3(v, x, y, z)
Definition: matrix1d.h:124
#define STARTINGZ(v)

◆ changeToVoxels()

void Basis::changeToVoxels ( GridVolume vol_basis,
MultidimArray< double > *  vol_voxels,
int  Zdim,
int  Ydim,
int  Xdim,
int  threads = 1 
) const

Change basis to voxels. This function takes a grid volume in the basis indicated in this object and translates it into voxels of the given size. If the volume is already in voxels a padding is done so that the output is of the given size and the basis volume is in the center.

Definition at line 261 of file basis.cpp.

263 {
264  int xdiff;
265  int ydiff;
266  int zdiff;
267  switch (type)
268  {
269  case blobs:
270  blobs2voxels(vol_basis, blob, vol_voxels, D, threads , Zdim, Ydim, Xdim );
271  break;
272  case voxels:
273  *vol_voxels = vol_basis(0)();
274  xdiff = (Xdim - XSIZE(*vol_voxels)) / 2;
275  ydiff = (Ydim - YSIZE(*vol_voxels)) / 2;
276  zdiff = (Zdim - ZSIZE(*vol_voxels)) / 2;
277  vol_voxels->selfWindow(
278  STARTINGZ(*vol_voxels) - zdiff,
279  STARTINGY(*vol_voxels) - ydiff,
280  STARTINGX(*vol_voxels) - xdiff,
281  (int)(FINISHINGZ(*vol_voxels) + Zdim - ZSIZE(*vol_voxels) - zdiff),
282  (int)(FINISHINGY(*vol_voxels) + Ydim - YSIZE(*vol_voxels) - ydiff),
283  (int)(FINISHINGX(*vol_voxels) + Xdim - XSIZE(*vol_voxels) - xdiff));
284  break;
285  case splines:
286  spatial_Bspline032voxels(vol_basis, vol_voxels, Zdim, Ydim, Xdim);
287  break;
288  }
289 }
void selfWindow(int n0, int z0, int y0, int x0, int nF, int zF, int yF, int xF, T init_value=0)
#define YSIZE(v)
#define FINISHINGX(v)
void blobs2voxels(const GridVolume &vol_blobs, const struct blobtype &blob, MultidimArray< double > *vol_voxels, const Matrix2D< double > *D, int threads, int Zdim, int Ydim, int Xdim)
Definition: blobs.cpp:926
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
#define FINISHINGZ(v)
Matrix2D< double > * D
Definition: basis.h:68
#define STARTINGX(v)
#define STARTINGY(v)
#define XSIZE(v)
#define ZSIZE(v)
#define FINISHINGY(v)
void spatial_Bspline032voxels(const GridVolume &vol_splines, MultidimArray< double > *vol_voxels, int Zdim, int Ydim, int Xdim)
Definition: splines.cpp:386
#define STARTINGZ(v)

◆ defineParams()

void Basis::defineParams ( XmippProgram program,
const char *  prefix = NULL,
const char *  comment = NULL 
)
static

Definition of paramaters

Definition at line 69 of file basis.cpp.

70 {
71  char tempLine[256];
72  char lineOut[512];
73 
74  if(prefix == nullptr)
75  sprintf(tempLine, " [--basis <basis_type=blobs>] ");
76  else
77  sprintf(tempLine,"%s --basis <basis_type=blobs> ", prefix);
78 
79  //std::cerr << "DEBUG_JM: tempLine: " << tempLine << std::endl;
80 
81 
82  if (comment != nullptr)
83  sprintf(lineOut, "%s : %s", tempLine, comment);
84  else
85  sprintf(lineOut, "%s : Basis function to use for the reconstruction", tempLine);
86 
87  program->addParamsLine(lineOut);
88 
89  //std::cerr << "DEBUG_JM: tempLine: " << lineOut << std::endl;
90 
91  program->addParamsLine(" where <basis_type>");
92  program->addParamsLine(" blobs <radius=2> <Bessel_order=2> <alpha_param=10.4> : Default blob parameters and grid relative size adjusted to use small blobs");
93  program->addParamsLine(" voxels");
94  program->addParamsLine(" splines");
95  program->addParamsLine(" or --big_blobs : blob parameters and grid relative size adjusted to use big blobs (r=2, m=2, alpha=2.26, g=2.26)");
96  program->addParamsLine(" or --visual_blobs : blobs optimal for direct visualization (r=2.4, m=2, alpha=13.3633, g=1.41)");
97  program->addParamsLine(" or --tomo_blobs : blobs with high overlapping (r=3.469269, m=2, alpha=13.7385, g=1/sqrt(2))");
98 }
void addParamsLine(const String &line)

◆ maxLength()

double Basis::maxLength ( ) const

Max length of the basis. This is the maximum distance between the center of the basis and its further point.

Definition at line 242 of file basis.cpp.

243 {
244  switch (type)
245  {
246  case blobs:
247  return blob.radius;
248  break;
249  case voxels:
250  return sqrt(3.0) * 0.5;
251  break;
252  case splines:
253  return sqrt(3.0) * 2.0;
254  break;
255  }
256  return 0.0;
257 }
void sqrt(Image< double > &op)
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
double radius
Spatial radius in Universal System units.
Definition: blobs.h:115

◆ produceSideInfo()

void Basis::produceSideInfo ( const Grid grid)

Produce side information. You must provide the grid in which this basis function will live

Definition at line 162 of file basis.cpp.

163 {
164  switch (type)
165  {
166  case blobs:
167  {
168  int subsampling = (VolPSF == nullptr )? BLOB_SUBSAMPLING : PIXEL_SUBSAMPLING;
169 
170  footprint_blob(blobprint, blob, subsampling);
171  sum_on_grid = sum_blob_Grid(blob, grid, D);
172  blobprint() /= sum_on_grid;
173 
174  if (VolPSF != nullptr)
175  {
176  // let adjust to the same resolution and size both blobprint and VolPSF
177  // selfScaleToSize(LINEAR, *VolPSF, XSIZE(*VolPSF)*BLOB_SUBSAMPLING,
178  // YSIZE(*VolPSF)*BLOB_SUBSAMPLING, ZSIZE(*VolPSF));
179 
180  blobprint().setXmippOrigin();
181 
182  if (XSIZE(blobprint()) < XSIZE(*VolPSF))
183  blobprint().selfWindow(STARTINGY(*VolPSF),
184  STARTINGX(*VolPSF),
186  STARTINGX(*VolPSF)+XSIZE(*VolPSF)-1);
187  else if (XSIZE(blobprint()) > XSIZE(*VolPSF))
189  STARTINGX(blobprint()),
192 
193  // creation of the 3D blobprint from convolution of blobprint and PSF
194  ImageOver footprintT;
195  footprintT.init(*VolPSF, subsampling);
196 
197  convolutionFFTStack(*VolPSF, blobprint(), footprintT());
198 
199  blobprint = footprintT;
200  }
201 
202  blobprint2() = blobprint();
203  blobprint2() *= blobprint();
204  break;
205  }
206  case voxels: sum_on_grid = 1;
207  break;
209  break;
210  }
211 
212 #ifdef DEBUG
213  std::cout << "Sum of a basis on the grid=" << sum_on_grid << std::endl;
214  std::cout << "D\n" << D << std::endl;
215  ImageXmipp save;
216  save() = blobprint();
217  save.write("footprint.xmp");
218 #endif
219 }
void selfWindow(int n0, int z0, int y0, int x0, int nF, int zF, int yF, int xF, T init_value=0)
#define YSIZE(v)
void init(int _umin, int _umax, int _uistep, int _vmin, int _vmax, int _vistep, int _wmin=0, int _wmax=0, int _wistep=1)
double sum_on_grid
Sum of the basis on the grid points.
Definition: basis.h:77
struct blobtype blob
Blob parameters.
Definition: basis.h:61
void convolutionFFTStack(const MultidimArray< double > &img, const MultidimArray< double > &kernel, MultidimArray< double > &result)
Definition: xmipp_fftw.cpp:429
tBasisFunction type
Basis function to use.
Definition: basis.h:52
ImageOver blobprint2
Square of the footprint.
Definition: basis.h:74
Matrix2D< double > * D
Definition: basis.h:68
#define STARTINGX(v)
#define STARTINGY(v)
double sum_blob_Grid(const struct blobtype &blob, const Grid &grid, const Matrix2D< double > *D)
Definition: blobs.cpp:320
MultidimArray< double > * VolPSF
Footprint is convolved with a volume PSF // At this moment only used with blobs.
Definition: basis.h:55
#define XSIZE(v)
void footprint_blob(ImageOver &blobprint, const struct blobtype &blob, int istep, int normalise)
Definition: blobs.cpp:417
double sum_spatial_Bspline03_Grid(const Grid &grid)
Definition: splines.cpp:72
const int BLOB_SUBSAMPLING
Definition: basis.h:34
const int PIXEL_SUBSAMPLING
Definition: basis.h:35
ImageOver blobprint
Blob footprint.
Definition: basis.h:71

◆ projectionAt()

double Basis::projectionAt ( const Matrix1D< double > &  u,
const Matrix1D< double > &  r 
) const

Projection at a given direction (u) with a given point (r).

Definition at line 367 of file basis.cpp.

368 {
369  const double p0 = 1.0 / (2 * PIXEL_SUBSAMPLING) - 0.5;
370  const double pStep = 1.0 / PIXEL_SUBSAMPLING;
371  const double pAvg = 1.0 / (PIXEL_SUBSAMPLING * PIXEL_SUBSAMPLING);
372  double module_r;
373  double px;
374  double py;
375  int i;
376  int j;
377  switch (type)
378  {
379  case blobs:
380  module_r = sqrt(XX(r) * XX(r) + YY(r) * YY(r) + ZZ(r) * ZZ(r));
381  return blob_proj(module_r, blob);
382  break;
383  case voxels:
384  {
385  double retval = 0;
386  ZZ(aux) = ZZ(r);
387  for (i = 0, px = p0; i < PIXEL_SUBSAMPLING; i++, px += pStep)
388  {
389  XX(aux) = XX(r) + px;
390  for (j = 0, py = p0; j < PIXEL_SUBSAMPLING; j++, py += pStep)
391  {
392  YY(aux) = YY(r) + py;
393  retval += intersection_unit_cube(u, aux);
394  }
395  }
396  return retval*pAvg;
397  break;
398  }
399  case splines:
400  return spatial_Bspline03_proj(r, u);
401  break;
402  }
403  return 0.0;
404 }
void sqrt(Image< double > &op)
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
#define i
#define XX(v)
Definition: matrix1d.h:85
double intersection_unit_cube(const Matrix1D< double > &u, const Matrix1D< double > &r)
Definition: geometry.cpp:1190
Matrix1D< double > aux
Auxiliary vector for projections.
Definition: basis.h:80
#define j
#define YY(v)
Definition: matrix1d.h:93
#define blob_proj(r, blob)
Definition: blobs.h:161
double spatial_Bspline03_proj(const Matrix1D< double > &r, const Matrix1D< double > &u)
Definition: splines.cpp:103
const int PIXEL_SUBSAMPLING
Definition: basis.h:35
#define ZZ(v)
Definition: matrix1d.h:101

◆ read() [1/2]

void Basis::read ( int  argc,
char **  argv 
)

Read parameters from a command line. This function reads the parameters from a command line defined by argc and argv. An exception might be thrown by any of the internal conversions, this would mean that there is an error in the command line and you might show a usage message.

◆ read() [2/2]

void Basis::read ( const FileName fn)

Read parameters from a file. An exception is thrown if the file cannot be open

◆ readParams()

void Basis::readParams ( XmippProgram program)

Read the parameters from the command line

Definition at line 100 of file basis.cpp.

101 {
102  String basisType = program->getParam("--basis");
103 
104  if (basisType == "blobs")
105  {
106  blob.radius = program->getDoubleParam("--basis", "blobs", 0);
107  blob.order = program->getIntParam("--basis", "blobs", 1);
108  blob.alpha = program->getDoubleParam("--basis", "blobs", 2);
109  }
110  if (!program->checkParam("--basis")) // Default is for small blobs
111  grid_relative_size = 1.41;
112 
113  if (program->checkParam("--big_blobs"))
114  {
115  blob.radius = 2;
116  blob.order = 2;
117  blob.alpha = 3.6;
118 
119  grid_relative_size = 2.26;
120  }
121  else if (program->checkParam("--visual_blobs"))
122  {
123  blob.radius = 2.4;
124  blob.order = 2;
125  blob.alpha = 13.3633;
126 
127  grid_relative_size = 1.41;
128  }
129  else if (program->checkParam("--tomo_blobs"))
130  {
131  blob.radius = 3.469269392209083;
132  blob.order = 2;
133  blob.alpha = 13.7385068372842;
134 
135  grid_relative_size = 1.0/sqrt(2);
136  }
137  else
138  {
139  if (basisType == "voxels")
140  {
141  type = voxels;
142  grid_relative_size = 1;
143  }
144  else if (basisType == "splines")
145  {
146  type = splines;
147  grid_relative_size = 1;
148  }
149  }
150 }
double alpha
Smoothness parameter.
Definition: blobs.h:121
double getDoubleParam(const char *param, int arg=0)
double grid_relative_size
Relative size for the grid.
Definition: basis.h:64
void sqrt(Image< double > &op)
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
const char * getParam(const char *param, int arg=0)
std::string String
Definition: xmipp_strings.h:34
int order
Derivation order and Bessel function order.
Definition: blobs.h:118
bool checkParam(const char *param)
int getIntParam(const char *param, int arg=0)
double radius
Spatial radius in Universal System units.
Definition: blobs.h:115

◆ setD()

void Basis::setD ( Matrix2D< double > *  _D)
inline

Set D. D is the deformation matrix used for crystals.

Definition at line 123 of file basis.h.

124  {
125  D = _D;
126  }
Matrix2D< double > * D
Definition: basis.h:68

◆ setDefault()

void Basis::setDefault ( )

Default values.

Definition at line 36 of file basis.cpp.

37 {
38  type = blobs;
39  blob.radius = 2;
40  blob.order = 2;
41  blob.alpha = 10.4;
42  VolPSF = nullptr;
43  D = nullptr;
44  blobprint.clear();
45  blobprint2.clear();
46  aux.resizeNoCopy(3);
47 }
double alpha
Smoothness parameter.
Definition: blobs.h:121
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
ImageOver blobprint2
Square of the footprint.
Definition: basis.h:74
Matrix2D< double > * D
Definition: basis.h:68
MultidimArray< double > * VolPSF
Footprint is convolved with a volume PSF // At this moment only used with blobs.
Definition: basis.h:55
Matrix1D< double > aux
Auxiliary vector for projections.
Definition: basis.h:80
int order
Derivation order and Bessel function order.
Definition: blobs.h:118
void resizeNoCopy(int Xdim)
Definition: matrix1d.h:458
double radius
Spatial radius in Universal System units.
Definition: blobs.h:115
ImageOver blobprint
Blob footprint.
Definition: basis.h:71

◆ setSamplingRate()

void Basis::setSamplingRate ( double  _Tm)

Set sampling rate.

Definition at line 153 of file basis.cpp.

154 {
155  if (_Tm == 0)
156  REPORT_ERROR(ERR_VALUE_INCORRECT, "Basis::set_sampling_rate: Sampling rate cannot be 0");
157  Tm = _Tm;
158  blob.radius /= Tm;
159 }
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
struct blobtype blob
Blob parameters.
Definition: basis.h:61
double Tm
Sampling rate.
Definition: basis.h:58
Incorrect value received.
Definition: xmipp_error.h:195
double radius
Spatial radius in Universal System units.
Definition: blobs.h:115

◆ valueAt()

double Basis::valueAt ( const Matrix1D< double > &  r) const

Basis value at a given point.

Definition at line 332 of file basis.cpp.

333 {
334  double module_r;
335  switch (type)
336  {
337  case blobs:
338  {
339  module_r = sqrt(XX(r) * XX(r) + YY(r) * YY(r) + ZZ(r) * ZZ(r));
340  return blob_val(module_r, blob);
341  break;
342  }
343  case voxels:
344  {
345  if (-0.5 <= XX(r) && XX(r) < 0.5 &&
346  -0.5 <= YY(r) && YY(r) < 0.5 &&
347  -0.5 <= ZZ(r) && ZZ(r) < 0.5)
348  return 1.0;
349  else
350  return 0.0;
351  break;
352  }
353  case splines:
354  {
355  if (-2 <= XX(r) && XX(r) < 2 &&
356  -2 <= YY(r) && YY(r) < 2 &&
357  -2 <= ZZ(r) && ZZ(r) < 2)
358  return spatial_Bspline03LUT(r);
359  else
360  return 0.0;
361  break;
362  }
363  }
364  return 0.0;
365 }
void sqrt(Image< double > &op)
double spatial_Bspline03LUT(const Matrix1D< double > &r)
Definition: splines.h:60
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
#define XX(v)
Definition: matrix1d.h:85
#define YY(v)
Definition: matrix1d.h:93
#define blob_val(r, blob)
Definition: blobs.h:139
#define ZZ(v)
Definition: matrix1d.h:101

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  out,
const Basis basis 
)
friend

Show.

Definition at line 223 of file basis.cpp.

224 {
225  switch (basis.type)
226  {
227  case Basis::blobs:
228  out << " Blobs: radius=" << basis.blob.radius << " pixels"
229  << " alpha=" << basis.blob.alpha
230  << " order=" << basis.blob.order << std::endl;
231  break;
232  case Basis::voxels:
233  out << "Voxels\n";
234  break;
235  case Basis::splines:
236  out << "Splines\n";
237  break;
238  }
239  return out;
240 }
double alpha
Smoothness parameter.
Definition: blobs.h:121
struct blobtype blob
Blob parameters.
Definition: basis.h:61
tBasisFunction type
Basis function to use.
Definition: basis.h:52
int order
Derivation order and Bessel function order.
Definition: blobs.h:118
double radius
Spatial radius in Universal System units.
Definition: blobs.h:115

Member Data Documentation

◆ aux

Matrix1D<double> Basis::aux

Auxiliary vector for projections.

Definition at line 80 of file basis.h.

◆ blob

struct blobtype Basis::blob

Blob parameters.

Definition at line 61 of file basis.h.

◆ blobprint

ImageOver Basis::blobprint

Blob footprint.

Definition at line 71 of file basis.h.

◆ blobprint2

ImageOver Basis::blobprint2

Square of the footprint.

Definition at line 74 of file basis.h.

◆ D

Matrix2D<double>* Basis::D

Volume deformation matrix. See the documentation of BasicARTParameters for further explanation.

Definition at line 68 of file basis.h.

◆ grid_relative_size

double Basis::grid_relative_size

Relative size for the grid.

Definition at line 64 of file basis.h.

◆ sum_on_grid

double Basis::sum_on_grid

Sum of the basis on the grid points.

Definition at line 77 of file basis.h.

◆ Tm

double Basis::Tm

Sampling rate.

Definition at line 58 of file basis.h.

◆ type

tBasisFunction Basis::type

Basis function to use.

Definition at line 52 of file basis.h.

◆ VolPSF

MultidimArray<double>* Basis::VolPSF

Footprint is convolved with a volume PSF // At this moment only used with blobs.

Definition at line 55 of file basis.h.


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