Xmipp  v3.23.11-Nereus
multidim_array_generic.h
Go to the documentation of this file.
1 /***************************************************************************
2  * Authors: Joaquin Oton (joton@cnb.csic.es)
3  *
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_MULTIDIM_ARRAY_GENERIC_H_
27 #define CORE_MULTIDIM_ARRAY_GENERIC_H_
28 
29 #include "xmipp_datatype.h"
30 #include "multidim_array.h"
31 #include "utils/half.hpp"
32 
33 /* Switch among different datatypes.
34  *
35  * This macro replicates the code for the different data type options.
36  *
37  *@code
38  *
39  *#define MYFUNC(type) getSlice(k, *(MultidimArray<type>*)image, axis, n)
40  *
41  *SWITCHDATATYPE(datatype,MYFUNC)
42  *
43  *@endcode
44  */
45 #define SWITCHDATATYPE(datatype,OP) \
46  switch (datatype)\
47  {\
48  case DT_Double:\
49  {OP(double)};\
50  break;\
51  case DT_Float:\
52  {OP(float)};\
53  break;\
54  case DT_UInt:\
55  {OP(unsigned int)};\
56  break;\
57  case DT_Int:\
58  {OP(int)};\
59  break;\
60  case DT_Short:\
61  {OP(short)};\
62  break;\
63  case DT_UShort:\
64  {OP(unsigned short)};\
65  break;\
66  case DT_SChar:\
67  {OP(char)};\
68  break;\
69  case DT_UHalfByte:\
70  case DT_UChar:\
71  {OP(unsigned char)};\
72  break;\
73  case DT_Long:\
74  {OP(long)};\
75  break;\
76  case DT_ULong:\
77  {OP(unsigned long)};\
78  break;\
79  case DT_HalfFloat:\
80  {OP(half_float::half)};\
81  break;\
82  default:\
83  REPORT_ERROR(ERR_ARG_INCORRECT,"Do not know how to handle this type at this point");\
84  }
85 
87 
89 
96 #ifndef MULTIDIM_ARRAY_BASE
97 #define MULTIDIM_ARRAY_BASE(v) (*((v).data->im))
98 #endif
99 
100 #ifndef MULTIDIM_ARRAY_GENERIC
101 #define MULTIDIM_ARRAY_GENERIC(v) (*((v).data))
102 #endif
103 
104 #ifndef MULTIDIM_ARRAY_TYPE
105 #define MULTIDIM_ARRAY_TYPE(v, type) (*((MultidimArray<type>*)(v.im)))
106 #endif
107 
112 {
113 public:
116 
117 protected:
118  // Flag to allow destroy MultidimArrayBase im data.
120 
121 public:
122 
123  /* Empty constructor */
125  {
126  init();
127  }
128 
134 
135  /* Get an aliasSlice of the selected slice from the multidimarray
136  */
137  MultidimArrayGeneric(MultidimArrayGeneric &mdim, int select_slice);
138 
143 
144  /* Initialize
145  */
146  void init();
147 
148  /* Clear the MultidimArrayBase and others
149  */
150  void clear();
151 
152  /* Set the datatype of the multidimarray object
153  */
154  void setDatatype(DataType imgType);
155 
156  /* Return mdim pointing to a specific slice inside this mda generic
157  */
158  void aliasSlice(MultidimArrayGeneric &mdim, int select_slice);
159 
163  void link(MultidimArrayBase* array);
164 
168  void resize(size_t Ndim, int Zdim, int Ydim, int Xdim, bool copy=true)
169  {
170  im->resize(Ndim,Zdim,Ydim,Xdim,copy);
171  }
172 
175  void setDimensions(int Xdim, int Ydim, int Zdim, size_t Ndim)
176  {
177  im->setDimensions(Xdim,Ydim,Zdim,Ndim);
178  }
179 
180  void resize(ArrayDim &adim, bool copy=true)
181  {
182  im->resize(adim, copy);
183  }
184 
185  void resize(MultidimArrayGeneric &mdim, bool copy=true)
186  {
187  ArrayDim adim;
188  mdim.getDimensions(adim);
189  im->resize(adim, copy);
190  }
191 
194  {
195  im->selfReverseX();
196  }
199  {
200  im->selfReverseY();
201  }
204  {
205  im->selfReverseZ();
206  }
207 
211  template <typename T>
212  void getArrayPointer(T* &M) const
213  {
214 #define GETMULTIDIMARRAY(type) M = (T*) (((MultidimArray<type>*) im)->data);
216 #undef GETMULTIDIMARRAY
217 
218  }
219 
220  void * getArrayPointer() const
221  {
222  void *res;
223 #define GETMULTIDIMARRAY(type) res = (void *) (((MultidimArray<type>*) im)->data);
225 #undef GETMULTIDIMARRAY
226  return res;
227  }
228 
232  template <typename T>
234  {
235 #define GETMULTIDIMARRAY(type) M = (MultidimArray<T>*) ((MultidimArray<type>*) im);
237 #undef GETMULTIDIMARRAY
238 
239  }
240 
242  void window(MultidimArrayGeneric &result, int z0, int y0, int x0,
243  int zF, int yF, int xF,
244  double init_value = 0.) const
245  {
246 #define WINDOW(type) _window(*((MultidimArray<type>*)(result.im)), z0,y0,x0,zF,yF,xF,init_value);
248 #undef WINDOW
249 
250  }
251 
252  template <class T>
253  void _window(MultidimArray<T> &result, int z0, int y0, int x0,
254  int zF, int yF, int xF,
255  double init_value = 0.) const
256  {
257 #define WINDOW(type) ((MultidimArray<type>*)im)->window(result, z0,y0,x0,zF,yF,xF,(T)init_value);
258  SWITCHDATATYPE(datatype,WINDOW)
259 #undef WINDOW
260 
261  }
262 
263  void selfWindow(int z0, int y0, int x0,
264  int zF, int yF, int xF,
265  double init_value = 0.)
266  {
267  if (im->mmapOn)
268  REPORT_ERROR(ERR_MMAP, "Cannot resize the image when it is mapped to file.");
269 
270 #define WINDOW(type) ((MultidimArray<type>*)im)->selfWindow(z0,y0,x0,zF,yF,xF,(type)init_value);
271 
272  SWITCHDATATYPE(datatype,WINDOW)
273 #undef WINDOW
274 
275  }
276 
277  void patch(MultidimArrayGeneric &patchArray, int x, int y)
278  {
279 #define PATCH(type) ((MultidimArray<type>*)im)->patch(*((MultidimArray<type>*)(patchArray.im)), x, y);
280  SWITCHDATATYPE(datatype, PATCH)
281 #undef PATCH
282 
283  }
284 
288  template <typename T>
289  void getSlice(int k, MultidimArray<T> &M, char axis = 'Z', bool reverse = false, size_t n = 0) const
290  {
291 #define GETSLICE(type) ((MultidimArray<type>*) im)->getSlice(k, M, axis, reverse, n);
292  SWITCHDATATYPE(datatype,GETSLICE)
293 #undef GETSLICE
294 
295  }
296 
300  void getSlice(int k, MultidimArrayGeneric* M, char axis = 'Z', bool reverse = false, size_t n = 0) const
301  {
302 #define GETSLICE(type) getSlice(k, *(MultidimArray<type>*)M->im, axis, reverse, n);
304 #undef GETSLICE
305 
306  }
307 
311  template <typename T1>
312  void setSlice(int k, const MultidimArray <T1>& v, size_t n = 0)
313  {
314 #define SETSLICE(type) ((MultidimArray<type>*) im)->setSlice(k, v, n);
315 
316  SWITCHDATATYPE(datatype,SETSLICE)
317 
318 #undef SETSLICE
319 
320  }
321 
325  void setSlice(int k, const MultidimArrayGeneric* v, size_t n = 0)
326  {
327 #define SETSLICE(type) setSlice(k,*(MultidimArray<type>*) v->im, n);
329 #undef SETSLICE
330 
331  }
332 
336  void getDimensions(size_t& Xdim, size_t& Ydim, size_t& Zdim, size_t &Ndim) const
337  {
338  im->getDimensions(Xdim,Ydim,Zdim,Ndim);
339  }
340 
341  void getDimensions(size_t& Xdim, size_t& Ydim, size_t& Zdim) const
342  {
343  size_t Ndim;
344  im->getDimensions(Xdim,Ydim,Zdim,Ndim);
345  }
346 
347  void getDimensions(size_t& Xdim, size_t& Ydim) const
348  {
349  size_t Zdim, Ndim;
350  im->getDimensions(Xdim,Ydim,Zdim,Ndim);
351  }
352 
354  {
355  im->getDimensions(adim);
356  }
357 
361  inline void setXmippOrigin()
362  {
363  im->setXmippOrigin();
364  }
365 
366  void maxIndex(ArrayCoord &pos)
367  {
368  im->maxIndex(pos);
369  }
370 
372  double computeAvg() const
373  {
374  return im->computeAvg();
375  }
376 
382  void computeStats(double& avg, double& stddev, double& minval, double& maxval) const
383  {
384 #define COMPUTESTATS(type) type Tminval(0); \
385  type Tmaxval(0); \
386  ((MultidimArray<type>*)(im))->computeStats(avg, stddev, Tminval, Tmaxval);\
387  minval = Tminval;\
388  maxval = Tmaxval;
389 
390  SWITCHDATATYPE(datatype, COMPUTESTATS)
391 #undef COMPUTESTATS
392 
393  }
394 
397  void computeDoubleMinMax(double& minval, double& maxval) const
398  {
399 #define COMPUTESDOUBLEMINMAX(type) ((MultidimArray<type>*)(im))->computeDoubleMinMax(minval, maxval);
401 #undef COMPUTESDOUBLEMINMAX
402 
403  }
404 
408  void rangeAdjust(const MultidimArrayGeneric &example, const MultidimArray<int> *mask=NULL)
409  {
410 #define RANGEADJUST(type) ((MultidimArray<type>*)(im))->rangeAdjust(*(MultidimArray<type>*)(example.im),mask);
411  SWITCHDATATYPE(datatype,RANGEADJUST)
412 #undef RANGEADJUST
413  }
414 
417  {
418  if (&input != this && input.datatype != DT_Unknown)
419  {
420  setDatatype(input.datatype);
421  *im = *input.im;
422  }
423  return *this;
424  }
425 
427  void operator*=(double op1)
428  {
429 #define OPERATORTIMESEQUAL(type) *((MultidimArray<type>*)(im)) *= (type)op1;
431 #undef OPERATORTIMESEQUAL
432  }
433 
436  {
437 #define OPERATORPLUSEQUAL(type) *((MultidimArray<type>*)(im)) += *((MultidimArray<type>*)(op1.im));
439 #undef OPERATORPLUSEQUAL
440  }
441 
443  double operator()(size_t n, int k, int i, int j) const
444  {
445  double ret;
446 #define GETVALUE(type) ret = NZYX_ELEM(*(MultidimArray<type>*)im,n,k,i,j);
447  SWITCHDATATYPE(datatype,GETVALUE)
448 #undef GETVALUE
449  return ret;
450  }
454  bool operator==(const MultidimArrayGeneric &mdA) const;
455 
461  bool equal(const MultidimArrayGeneric &op,
462  double accuracy = XMIPP_EQUAL_ACCURACY) const;
463 
465  double operator()(int i, int j) const
466  {
467  double ret;
468 #define GETVALUE(type) ret = A2D_ELEM(*(MultidimArray<type>*)im,i,j);
469  SWITCHDATATYPE(datatype,GETVALUE)
470 #undef GETVALUE
471  return ret;
472  }
473 
476  {
477  return *im;
478  }
479 
482  {
483  return *im;
484  }
485 
486 
489  template <typename T>
490  void getImage(MultidimArray<T> &M) const
491  {
492 #define TYPECAST(type) typeCast(*(MultidimArray<type>*)(im), M);
493  SWITCHDATATYPE(datatype, TYPECAST)
494 #undef TYPECAST
495 
496  }
497 
500  template <typename T>
502  {
503 #define TYPECAST(type) typeCast(M, *(MultidimArray<type>*)(im));
504  SWITCHDATATYPE(datatype, TYPECAST)
505 #undef TYPECAST
506 
507  }
508 
511  template <typename T1>
512  void window(MultidimArray<T1> &result, int n0,int z0, int y0, int x0,
513  int nF,int zF, int yF, int xF,
514  T1 init_value = 0) const
515  {
516 #define WINDOW(type) ((MultidimArray<type>*)(im))->window(result, n0, z0, y0, x0, nF, zF, yF, xF, init_value);
517  SWITCHDATATYPE(datatype, WINDOW)
518 #undef WINDOW
519 
520  }
521 }
522 ;
524 
525 #endif /* MULTIDIM_ARRAY_GENERIC_H_ */
void window(MultidimArrayGeneric &result, int z0, int y0, int x0, int zF, int yF, int xF, double init_value=0.) const
void window(MultidimArray< T1 > &result, int n0, int z0, int y0, int x0, int nF, int zF, int yF, int xF, T1 init_value=0) const
#define WINDOW(type)
MultidimArrayGeneric & operator=(const MultidimArrayGeneric &input)
void aliasSlice(MultidimArrayGeneric &mdim, int select_slice)
#define COMPUTESTATS(type)
#define REPORT_ERROR(nerr, ErrormMsg)
Definition: xmipp_error.h:211
void computeStats(double &avg, double &stddev, double &minval, double &maxval) const
void selfWindow(int z0, int y0, int x0, int zF, int yF, int xF, double init_value=0.)
void computeDoubleMinMax(double &minval, double &maxval) const
Global mmap error.
Definition: xmipp_error.h:170
void setSlice(int k, const MultidimArray< T1 > &v, size_t n=0)
virtual double computeAvg() const =0
void getMultidimArrayPointer(MultidimArray< T > *&M) const
static double * y
#define PATCH(type)
virtual void selfReverseX()=0
void resize(ArrayDim &adim, bool copy=true)
#define zF
#define z0
doublereal * x
#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
void setDatatype(DataType imgType)
#define GETSLICE(type)
#define COMPUTESDOUBLEMINMAX(type)
void getSlice(int k, MultidimArray< T > &M, char axis='Z', bool reverse=false, size_t n=0) const
void getImage(MultidimArray< T > &M) const
void _window(MultidimArray< T > &result, int z0, int y0, int x0, int zF, int yF, int xF, double init_value=0.) const
#define y0
char axis
#define x0
#define OPERATORPLUSEQUAL(type)
virtual void selfReverseY()=0
void getSlice(int k, MultidimArrayGeneric *M, char axis='Z', bool reverse=false, size_t n=0) const
void getDimensions(size_t &Xdim, size_t &Ydim) const
void resize(size_t Ndim, int Zdim, int Ydim, int Xdim, bool copy=true)
#define XMIPP_EQUAL_ACCURACY
Definition: xmipp_macros.h:119
double operator()(size_t n, int k, int i, int j) const
bool equal(const MultidimArrayGeneric &op, double accuracy=XMIPP_EQUAL_ACCURACY) const
void setImage(MultidimArray< T > &M)
void resize(MultidimArrayGeneric &mdim, bool copy=true)
void patch(MultidimArrayGeneric &patchArray, int x, int y)
DataType
void link(MultidimArrayBase *array)
void maxIndex(ArrayCoord &pos)
#define xF
bool operator==(const MultidimArrayGeneric &mdA) const
MultidimArrayBase & operator()()
void getArrayPointer(T *&M) const
void setDimensions(int Xdim, int Ydim, int Zdim, size_t Ndim)
#define j
virtual void resize(size_t Ndim, size_t Zdim, size_t Ydim, size_t Xdim, bool copy=true)=0
#define TYPECAST(type)
double operator()(int i, int j) const
void rangeAdjust(const MultidimArrayGeneric &example, const MultidimArray< int > *mask=NULL)
void getDimensions(size_t &Xdim, size_t &Ydim, size_t &Zdim) const
virtual void selfReverseZ()=0
#define OPERATORTIMESEQUAL(type)
void setDimensions(int Xdim, int Ydim, int Zdim, size_t Ndim)
#define GETMULTIDIMARRAY(type)
void getDimensions(size_t &Xdim, size_t &Ydim, size_t &Zdim, size_t &Ndim) const
void operator+=(const MultidimArrayGeneric &op1)
virtual void maxIndex(size_t &lmax, int &kmax, int &imax, int &jmax) const =0
void setSlice(int k, const MultidimArrayGeneric *v, size_t n=0)
#define GETVALUE(type)
#define yF
int * n
#define SETSLICE(type)
#define SWITCHDATATYPE(datatype, OP)
#define RANGEADJUST(type)
void getDimensions(ArrayDim &adim)
void getDimensions(size_t &Xdim, size_t &Ydim, size_t &Zdim, size_t &Ndim) const
const MultidimArrayBase & operator()() const