Xmipp  v3.23.11-Nereus
blobs.h
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: Carlos Oscar S. Sorzano (coss@cnb.csic.es)
4  *
5  * Unidad de Bioinformatica of Centro Nacional de Biotecnologia , CSIC
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software
19  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA
20  * 02111-1307 USA
21  *
22  * All comments concerning this program package may be sent to the
23  * e-mail address 'xmipp@cnb.csic.es'
24  ***************************************************************************/
25 
26 #ifndef _CORE_BLOBS_HH
27 #define _CORE_BLOBS_HH
28 
29 #include "grids.h"
30 
31 class ImageOver;
32 
33 // Thread declaration
34 void * blobs2voxels_SimpleGrid( void * data );
35 
36 // Mutex used to make load balancing thread-safe
37 extern pthread_mutex_t blobs_conv_mutex;
38 
39 // There exist one integer per slice. If value == 0, this slice
40 // can be processed. If value == -1, this slice has already been processed and
41 // if value > 0, this slice is influenced by a neighbour slice being processed and,
42 // at this moment, can not be scheduled for processing. In the latest case, the
43 // threads will have to look for a value == 0 slice.
44 extern int * slices_status;
45 
46 // Number of already processed slices.
47 extern int slices_processed;
48 
49 // This structure is needed to pass parameters to threads.
50 typedef struct{
52  const SimpleGrid *grid;
53  const struct blobtype *blob;
56  int istep;
59  bool FORW;
60  int eq_mode;
61  int thread_id;
65 
66 /* ========================================================================= */
67 /* BLOBS */
68 /* ========================================================================= */
72 // Blob structure ----------------------------------------------------------
112 struct blobtype
113 {
115  double radius;
116 
118  int order;
119 
121  double alpha;
122 };
123 
124 // Blob value --------------------------------------------------------------
139 #define blob_val(r, blob) kaiser_value(r, blob.radius, blob.alpha, blob.order)
140 
142 double kaiser_value(double r, double a, double alpha, int m);
143 
144 // Blob projection ---------------------------------------------------------
161 #define blob_proj(r, blob) kaiser_proj(r, blob.radius, blob.alpha, blob.order)
162 
164 double kaiser_proj(double r, double a, double alpha, int m);
165 
174 #define blob_Fourier_val(w, blob) \
175  kaiser_Fourier_value(w, blob.radius, blob.alpha, blob.order)
176 
178 double kaiser_Fourier_value(double w, double a, double alpha, int m);
179 
181 #define blob_mass(blob) \
182  basvolume(blob.radius, blob.alpha, blob.order,3)
183 
185 double basvolume(double a, double alpha, int m, int n);
186 
188 double in_zeroarg(int n);
189 
191 double inph_zeroarg(int n);
192 
194 double i_nph(int n, double x);
195 
198 double i_n(int n, double x);
199 
202 double blob_freq_zero(struct blobtype b);
203 
208 double blob_att(double w, struct blobtype b);
209 
213 double blob_ops(double w, struct blobtype b);
214 
219 
224 
229 
230 /* ========================================================================= */
231 /* TOOLS */
232 /* ========================================================================= */
240 blobtype best_blob(double alpha_0, double alpha_F, double inc_alpha,
241  double a_0, double a_F, double inc_a, double w, double *target_att,
242  int target_length = 1);
243 
271 void footprint_blob(ImageOver &blobprint, const struct blobtype &blob,
272  int istep = 50, int normalise = 0);
273 
299 double sum_blob_Grid(const struct blobtype &blob, const Grid &grid,
300  const Matrix2D<double> *D = nullptr);
301 
302 
307 void voxel_volume_shape(const GridVolume &vol_blobs,
308  const struct blobtype &blob, const Matrix2D<double> *D,
309  Matrix1D<int> &corner1, Matrix1D<int> &size);
310 
343 void blobs2voxels(const GridVolume &vol_blobs,
344  const struct blobtype &blob, MultidimArray<double> *vol_voxels,
345  const Matrix2D<double> *D = nullptr, int threads = 1, int Zdim = 0, int Ydim = 0, int Xdim = 0);
346 
350 void blobs2space_coefficients(const GridVolume &vol_blobs, const struct blobtype &blob,
351  MultidimArray<double> *vol_coefs);
352 
353 constexpr int SHOW_CONVERSION = 1;
371 void voxels2blobs(const MultidimArray<double> *vol_voxels,
372  const struct blobtype &blob, GridVolume &vol_blobs,
373  int grid_type, double grid_relative_size,
374  double lambda, const MultidimArray<double> *vol_mask = nullptr,
375  const Matrix2D<double> *D = nullptr,
376  double final_error_change = 0.01, int tell = 0, double R = -1, int threads = 1);
377 
378 constexpr int VARTK = 1;
379 constexpr int VMAXARTK = 2;
380 
404  GridVolume &vol_in, GridVolume *vol_out, const struct blobtype &blob,
405  const Matrix2D<double> *D, double lambda, MultidimArray<double> *theo_vol,
407  const MultidimArray<double> *mask_vol,
408  double &mean_error, double &max_error, int eq_mode = VARTK, int threads=1);
410 #endif
double kaiser_value(double r, double a, double alpha, int m)
Definition: blobs.cpp:37
void voxel_volume_shape(const GridVolume &vol_blobs, const struct blobtype &blob, const Matrix2D< double > *D, Matrix1D< int > &corner1, Matrix1D< int > &size)
Definition: blobs.cpp:851
double alpha
Smoothness parameter.
Definition: blobs.h:121
constexpr int SHOW_CONVERSION
Definition: blobs.h:353
double optimal_FCC_grid_relative_size(struct blobtype b)
Definition: blobs.cpp:356
void blobs2voxels(const GridVolume &vol_blobs, const struct blobtype &blob, MultidimArray< double > *vol_voxels, const Matrix2D< double > *D=nullptr, int threads=1, int Zdim=0, int Ydim=0, int Xdim=0)
Definition: blobs.cpp:926
MultidimArray< double > * vol_corr
Definition: blobs.h:57
void ART_voxels2blobs_single_step(GridVolume &vol_in, GridVolume *vol_out, const struct blobtype &blob, const Matrix2D< double > *D, double lambda, MultidimArray< double > *theo_vol, const MultidimArray< double > *read_vol, MultidimArray< double > *corr_vol, const MultidimArray< double > *mask_vol, double &mean_error, double &max_error, int eq_mode=VARTK, int threads=1)
Definition: blobs.cpp:1094
constexpr int VARTK
Definition: blobs.h:378
void blobs2space_coefficients(const GridVolume &vol_blobs, const struct blobtype &blob, MultidimArray< double > *vol_coefs)
Definition: blobs.cpp:1029
void read_vol(char *vol_file, double *width, double *origx, double *origy, double *origz, unsigned *extx, unsigned *exty, unsigned *extz, double **phi)
Definition: lib_vio.cpp:25
doublereal * w
double blob_freq_zero(struct blobtype b)
Definition: blobs.cpp:330
const SimpleGrid * grid
Definition: blobs.h:52
int * slices_status
Definition: blobs.cpp:33
Definition: grids.h:479
double i_nph(int n, double x)
Definition: blobs.cpp:270
doublereal * x
MultidimArray< double > data
Definition: xmipp_image.h:55
doublereal * b
double sum_blob_Grid(const struct blobtype &blob, const Grid &grid, const Matrix2D< double > *D=nullptr)
Definition: blobs.cpp:320
double * lambda
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=nullptr, const Matrix2D< double > *D=nullptr, double final_error_change=0.01, int tell=0, double R=-1, int threads=1)
Definition: blobs.cpp:1316
double blob_ops(double w, struct blobtype b)
Definition: blobs.cpp:342
const struct blobtype * blob
Definition: blobs.h:53
int slices_processed
Definition: blobs.cpp:34
blobtype best_blob(double alpha_0, double alpha_F, double inc_alpha, double a_0, double a_F, double inc_a, double w, double *target_att, int target_length=1)
Definition: blobs.cpp:363
int min_separation
Definition: blobs.h:63
double blob_att(double w, struct blobtype b)
Definition: blobs.cpp:336
pthread_mutex_t blobs_conv_mutex
Definition: blobs.cpp:31
void footprint_blob(ImageOver &blobprint, const struct blobtype &blob, int istep=50, int normalise=0)
Definition: blobs.cpp:417
double kaiser_Fourier_value(double w, double a, double alpha, int m)
Definition: blobs.cpp:144
int m
double i_n(int n, double x)
Definition: blobs.cpp:244
void * blobs2voxels_SimpleGrid(void *data)
Definition: blobs.cpp:452
const MultidimArray< double > * vol_blobs
Definition: blobs.h:51
const Matrix2D< double > * D
Definition: blobs.h:55
const MultidimArray< double > * vol_mask
Definition: blobs.h:58
int order
Derivation order and Bessel function order.
Definition: blobs.h:118
double in_zeroarg(int n)
Definition: blobs.cpp:294
double kaiser_proj(double r, double a, double alpha, int m)
Definition: blobs.cpp:94
double inph_zeroarg(int n)
Definition: blobs.cpp:307
double optimal_CC_grid_relative_size(struct blobtype b)
Definition: blobs.cpp:348
double optimal_BCC_grid_relative_size(struct blobtype b)
Definition: blobs.cpp:352
double radius
Spatial radius in Universal System units.
Definition: blobs.h:115
int * n
doublereal * a
constexpr int VMAXARTK
Definition: blobs.h:379
MultidimArray< double > * vol_voxels
Definition: blobs.h:54
double basvolume(double a, double alpha, int m, int n)
Definition: blobs.cpp:215