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

#include <volume_local_sharpening.h>

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

Public Member Functions

void defineParams ()
 
void readParams ()
 
void produceSideInfo ()
 
void lowPassFilterFunction (const MultidimArray< std::complex< double > > &myfftV, double w, double wL, MultidimArray< double > &filteredVol, int count)
 
void bandPassFilterFunction (const MultidimArray< std::complex< double > > &myfftV, double w, double wL, MultidimArray< double > &filteredVol, int count)
 
void wideBandPassFilter (const MultidimArray< std::complex< double > > &myfftV, double wmin, double wmax, double wL, MultidimArray< double > &filteredVol)
 
void maxMinResolution (MultidimArray< double > &resVol, double &maxRes, double &minRes)
 
void computeAvgStdev_within_binary_mask (const MultidimArray< double > &resVol, const MultidimArray< double > &vol, double &stddev, bool outside=false)
 
void localfiltering (MultidimArray< std::complex< double > > &myfftV, MultidimArray< double > &localfilteredVol, double &minFreq, double &maxFreq, double &step)
 
void amplitudeMonogenicSignalBP (MultidimArray< std::complex< double > > &myfftV, double w1, double w1l, MultidimArray< double > &amplitude, int count)
 
void sameEnergy (MultidimArray< std::complex< double > > &myfftV, MultidimArray< double > &localfilteredVol, double &minFreq, double &maxFreq, double &step)
 
void run ()
 
- Public Member Functions inherited from XmippProgram
const char * getParam (const char *param, int arg=0)
 
const char * getParam (const char *param, const char *subparam, int arg=0)
 
int getIntParam (const char *param, int arg=0)
 
int getIntParam (const char *param, const char *subparam, int arg=0)
 
double getDoubleParam (const char *param, int arg=0)
 
double getDoubleParam (const char *param, const char *subparam, int arg=0)
 
float getFloatParam (const char *param, int arg=0)
 
float getFloatParam (const char *param, const char *subparam, int arg=0)
 
void getListParam (const char *param, StringVector &list)
 
int getCountParam (const char *param)
 
bool checkParam (const char *param)
 
bool existsParam (const char *param)
 
void addParamsLine (const String &line)
 
void addParamsLine (const char *line)
 
ParamDefgetParamDef (const char *param) const
 
virtual void quit (int exit_code=0) const
 
virtual int tryRun ()
 
void initProgress (size_t total, size_t stepBin=60)
 
void setProgress (size_t value=0)
 
void endProgress ()
 
void processDefaultComment (const char *param, const char *left)
 
void setDefaultComment (const char *param, const char *comment)
 
virtual void initComments ()
 
void setProgramName (const char *name)
 
void addUsageLine (const char *line, bool verbatim=false)
 
void clearUsage ()
 
void addExampleLine (const char *example, bool verbatim=true)
 
void addSeeAlsoLine (const char *seeAlso)
 
void addKeywords (const char *keywords)
 
const char * name () const
 
virtual void usage (int verb=0) const
 
virtual void usage (const String &param, int verb=2)
 
int version () const
 
virtual void show () const
 
virtual void read (int argc, const char **argv, bool reportErrors=true)
 
virtual void read (int argc, char **argv, bool reportErrors=true)
 
void read (const String &argumentsLine)
 
 XmippProgram ()
 
 XmippProgram (int argc, const char **argv)
 
virtual ~XmippProgram ()
 

Public Attributes

FileName fnOut
 
FileName fnVol
 
FileName fnRes
 
FileName fnMD
 
double sampling
 
double maxRes
 
double minRes
 
double lambda
 
double K
 
double maxFreq
 
double minFreq
 
double desv_Vorig
 
double desvOutside_Vorig
 
int Niter
 
int Nthread
 
MultidimArray< double > Vorig
 
MultidimArray< int > mask
 
MultidimArray< double > resVol
 
MultidimArray< double > iu
 
MultidimArray< double > sharpenedMap
 
MultidimArray< std::complex< double > > fftV
 
MultidimArray< std::complex< double > > fftVfilter
 
FourierTransformer transformer
 
FourierTransformer transformer_inv
 
FourierFilter FilterBand
 
- Public Attributes inherited from XmippProgram
bool doRun
 
bool runWithoutArgs
 
int verbose
 Verbosity level. More...
 
int debug
 

Additional Inherited Members

- Protected Member Functions inherited from XmippProgram
void defineCommons ()
 
- Protected Attributes inherited from XmippProgram
int errorCode
 
ProgramDefprogDef
 Program definition and arguments parser. More...
 
std::map< String, CommentListdefaultComments
 
int argc
 Original command line arguments. More...
 
const char ** argv
 

Detailed Description

SSNR parameters.

Definition at line 49 of file volume_local_sharpening.h.

Member Function Documentation

◆ amplitudeMonogenicSignalBP()

void ProgLocSharpening::amplitudeMonogenicSignalBP ( MultidimArray< std::complex< double > > &  myfftV,
double  w1,
double  w1l,
MultidimArray< double > &  amplitude,
int  count 
)

◆ bandPassFilterFunction()

void ProgLocSharpening::bandPassFilterFunction ( const MultidimArray< std::complex< double > > &  myfftV,
double  w,
double  wL,
MultidimArray< double > &  filteredVol,
int  count 
)

Definition at line 181 of file volume_local_sharpening.cpp.

183 {
184  fftVfilter.initZeros(myfftV);
185  size_t xdim, ydim, zdim, ndim;
186  Vorig.getDimensions(xdim, ydim, zdim, ndim);
187 
188 
189  double delta = wL-w;
190  double w_inf = w-delta;
191  // Filter the input volume and add it to amplitude
192  long n=0;
193  double ideltal=PI/delta;
195  {
196  double un=DIRECT_MULTIDIM_ELEM(iu,n);
197  if (un>=w && un<=wL)
198  {
200  DIRECT_MULTIDIM_ELEM(fftVfilter, n) *= 0.5*(1+cos((un-w)*ideltal));//H;
201  } else if (un<=w && un>=w_inf)
202  {
204  DIRECT_MULTIDIM_ELEM(fftVfilter, n) *= 0.5*(1+cos((un-w)*ideltal));//H;
205  }
206  }
207 
208  filteredVol.resizeNoCopy(Vorig);
209 
211 
212 // #ifdef DEBUG_FILTER
213 // Image<double> filteredvolume;
214 // filteredvolume() = filteredVol;
215 // filteredvolume.write(formatString("Volumen_filtrado_%i.vol", count));
216 // #endif
217 }
void inverseFourierTransform()
Definition: xmipp_fftw.cpp:329
void resizeNoCopy(const MultidimArray< T1 > &v)
doublereal * w
MultidimArray< std::complex< double > > fftVfilter
#define FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(v)
#define DIRECT_MULTIDIM_ELEM(v, n)
FourierTransformer transformer_inv
MultidimArray< double > Vorig
MultidimArray< double > iu
void initZeros(const MultidimArray< T1 > &op)
#define PI
Definition: tools.h:43
int * n
double * delta
void getDimensions(size_t &Xdim, size_t &Ydim, size_t &Zdim, size_t &Ndim) const

◆ computeAvgStdev_within_binary_mask()

void ProgLocSharpening::computeAvgStdev_within_binary_mask ( const MultidimArray< double > &  resVol,
const MultidimArray< double > &  vol,
double &  stddev,
bool  outside = false 
)

Definition at line 151 of file volume_local_sharpening.cpp.

153 {
154 
156  double sum1 = 0;
157  double sum2 = 0;
158  int N = 0;
159  double avg=0;
160 
162  {
163  if ((not outside && A3D_ELEM(resVol, k, i, j) > 2*sampling) ||
164  (outside && A3D_ELEM(resVol, k, i, j) < 2*sampling))
165  {
166  ++N;
167  double aux=A3D_ELEM(vol, k, i, j);
168  sum1 += aux;
169  sum2 += aux*aux;
170  }
171  }
172 
173  // average and standard deviation
174  avg = sum1 / (double) N;
175  if (N > 1)
176  stddev = sqrt(fabs(sum2 / N - avg * avg) * N / (N - 1));
177  else
178  stddev = 0;
179 }
void sqrt(Image< double > &op)
#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 A3D_ELEM(V, k, i, j)
#define j
#define SPEED_UP_tempsInt
Definition: xmipp_macros.h:408
#define FOR_ALL_ELEMENTS_IN_COMMON_IN_ARRAY3D(V1, V2)

◆ defineParams()

void ProgLocSharpening::defineParams ( )
virtual

Function in which the param of each Program are defined.

Reimplemented from XmippProgram.

Definition at line 44 of file volume_local_sharpening.cpp.

45 {
46  addUsageLine("This function performs local sharpening");
47  addParamsLine(" --vol <vol_file=\"\"> : Input volume");
48  addParamsLine(" --resolution_map <vol_file=\"\">: Resolution map");
49  addParamsLine(" --sampling <s=1>: sampling");
50  addParamsLine(" -o <output=\"Sharpening.vol\">: sharpening volume");
51  addParamsLine(" [--md <output=\"params.xmd\">]: sharpening params");
52  addParamsLine(" [-l <lambda=1>]: regularization param");
53  addParamsLine(" [-k <K=0.025>]: K param");
54  addParamsLine(" [-i <Niter=50>]: iteration");
55  addParamsLine(" [-n <Nthread=1>]: threads number");
56 }
void addUsageLine(const char *line, bool verbatim=false)
void addParamsLine(const String &line)

◆ localfiltering()

void ProgLocSharpening::localfiltering ( MultidimArray< std::complex< double > > &  myfftV,
MultidimArray< double > &  localfilteredVol,
double &  minFreq,
double &  maxFreq,
double &  step 
)

Definition at line 219 of file volume_local_sharpening.cpp.

222 {
223  MultidimArray<double> filteredVol, lastweight, weight;
224  localfilteredVol.initZeros(Vorig);
225  weight.initZeros(Vorig);
226  lastweight.initZeros(Vorig);
227 
228  double freq, lastResolution=1e38;
229  int idx, lastidx = -1;
230 
231  for(double res = minRes; res<maxRes; res+=step)
232  {
233  freq = sampling/res;
234 
235  DIGFREQ2FFT_IDX(freq, ZSIZE(myfftV), idx);
236 
237  if (idx == lastidx)
238  {
239  continue;
240  }
241 
242  double wL = sampling/(res - step);
243 
244  bandPassFilterFunction(myfftV, freq, wL, filteredVol, idx);
245 
247  {
248 
250  {
251  DIRECT_MULTIDIM_ELEM(filteredVol, n)=0;
252  }
253  else
254  {
255  double res_map = DIRECT_MULTIDIM_ELEM(resVol, n);//+1e-38;
256  DIRECT_MULTIDIM_ELEM(weight, n) = (exp(-K*(res-res_map)*(res-res_map)));
257  DIRECT_MULTIDIM_ELEM(filteredVol, n) *= DIRECT_MULTIDIM_ELEM(weight, n);
258  }
259 
260  }
261 
262  localfilteredVol += filteredVol;
263  lastweight += weight;
264  lastResolution = res;
265  lastidx = idx;
266  }
267 // double sigmaBefore=0;
268 // computeAvgStdev_within_binary_mask(resVol, localfilteredVol, sigmaBefore);
269 
271  {
272  if (DIRECT_MULTIDIM_ELEM(lastweight, n)>0)
273  DIRECT_MULTIDIM_ELEM(localfilteredVol, n) /=DIRECT_MULTIDIM_ELEM(lastweight, n);
274  }
275 // double sigmaAfter=0;
276 // computeAvgStdev_within_binary_mask(resVol, localfilteredVol, sigmaAfter);
277 // std::cout << "sigmaBefore div=" << sigmaBefore << " sigmaAfter div=" << sigmaAfter << std::endl;
278 // FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(localfilteredVol)
279 // {
280 // if (DIRECT_MULTIDIM_ELEM(lastweight, n)>0)
281 // DIRECT_MULTIDIM_ELEM(localfilteredVol, n) *=0.01*sigmaBefore/sigmaAfter;
282 // }
283 }
void bandPassFilterFunction(const MultidimArray< std::complex< double > > &myfftV, double w, double wL, MultidimArray< double > &filteredVol, int count)
#define DIGFREQ2FFT_IDX(freq, size, idx)
Definition: xmipp_fft.h:61
#define FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(v)
#define ZSIZE(v)
#define DIRECT_MULTIDIM_ELEM(v, n)
MultidimArray< double > Vorig
void initZeros(const MultidimArray< T1 > &op)
int * n
MultidimArray< double > resVol

◆ lowPassFilterFunction()

void ProgLocSharpening::lowPassFilterFunction ( const MultidimArray< std::complex< double > > &  myfftV,
double  w,
double  wL,
MultidimArray< double > &  filteredVol,
int  count 
)

◆ maxMinResolution()

void ProgLocSharpening::maxMinResolution ( MultidimArray< double > &  resVol,
double &  maxRes,
double &  minRes 
)

Definition at line 132 of file volume_local_sharpening.cpp.

134 {
135  // Count number of voxels with resolution
136  size_t n=0;
137  double lastMinRes=1e38, lastMaxRes=1e-38, value;
139  {
140  value = DIRECT_MULTIDIM_ELEM(resVol, n);
141  if (value>lastMaxRes)
142  lastMaxRes = value;
143  if (value<lastMinRes && value>0)
144  lastMinRes = value;
145  }
146 
147  maxRes = lastMaxRes;
148  minRes = lastMinRes;
149 }
#define FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(v)
#define DIRECT_MULTIDIM_ELEM(v, n)
int * n

◆ produceSideInfo()

void ProgLocSharpening::produceSideInfo ( )

Definition at line 58 of file volume_local_sharpening.cpp.

59 {
60  std::cout << "Starting..." << std::endl;
61  Image<double> V;
62  V.read(fnVol);
63  V().setXmippOrigin();
64 
65 
66  if (Nthread>1)
67  {
68  std::cout << "used procesors = " << Nthread << std::endl;
71  }
72 
73  MultidimArray<double> &inputVol = V();
74  Vorig = inputVol;
75 
77 
78  iu.initZeros(fftV);
79  double uz, uy, ux, uz2, u2, uz2y2;
80  long n=0;
81  for(size_t k=0; k<ZSIZE(fftV); ++k)
82  {
83  FFT_IDX2DIGFREQ(k,ZSIZE(inputVol),uz);
84  uz2=uz*uz;
85 
86  for(size_t i=0; i<YSIZE(fftV); ++i)
87  {
88  FFT_IDX2DIGFREQ(i,YSIZE(inputVol),uy);
89  uz2y2=uz2+uy*uy;
90 
91  for(size_t j=0; j<XSIZE(fftV); ++j)
92  {
93  FFT_IDX2DIGFREQ(j,XSIZE(inputVol),ux);
94  u2=uz2y2+ux*ux;
95  if ((k != 0) || (i != 0) || (j != 0))
96  DIRECT_MULTIDIM_ELEM(iu,n) = sqrt(u2);
97  else
98  DIRECT_MULTIDIM_ELEM(iu,n) = 1e-38;
99  ++n;
100  }
101  }
102  }
103 
104  inputVol.clear();
105 
106  Image<double> resolutionVolume;
107  resolutionVolume.read(fnRes);
108 
109  resVol = resolutionVolume();
110  resolutionVolume().clear();
111 
113 // std::cout << "maxRes = " << maxRes << " minRes = " << minRes << std::endl;
114 
116  {
117 
119  {
121  }
122 
123  }
124 
126 
128  //std::cout << "desvOutside_Vorig = " << desvOutside_Vorig << std::endl;
129 
130 }
#define YSIZE(v)
void setThreadsNumber(int tNumber)
Definition: xmipp_fftw.h:119
void sqrt(Image< double > &op)
MultidimArray< std::complex< double > > fftV
#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 maxMinResolution(MultidimArray< double > &resVol, double &maxRes, double &minRes)
#define XSIZE(v)
#define FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(v)
#define ZSIZE(v)
#define DIRECT_MULTIDIM_ELEM(v, n)
void FourierTransform(T &v, T1 &V, bool getCopy=true)
Definition: xmipp_fftw.h:166
__device__ float FFT_IDX2DIGFREQ(int idx, int size)
#define j
FourierTransformer transformer_inv
void computeAvgStdev_within_binary_mask(const MultidimArray< double > &resVol, const MultidimArray< double > &vol, double &stddev, bool outside=false)
MultidimArray< double > Vorig
MultidimArray< double > iu
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
void initZeros(const MultidimArray< T1 > &op)
FourierTransformer transformer
int * n
MultidimArray< double > resVol
void clear()
Definition: xmipp_image.h:144

◆ readParams()

void ProgLocSharpening::readParams ( )
virtual

Function in which each program will read parameters that it need. If some error occurs the usage will be printed out.

Reimplemented from XmippProgram.

Definition at line 31 of file volume_local_sharpening.cpp.

32 {
33  fnVol = getParam("--vol");
34  fnRes = getParam("--resolution_map");
35  sampling = getDoubleParam("--sampling");
36  lambda = getDoubleParam("-l");
37  K= getDoubleParam("-k");
38  Niter = getIntParam("-i");
39  Nthread = getIntParam("-n");
40  fnOut = getParam("-o");
41  fnMD = getParam("--md");
42 }
double getDoubleParam(const char *param, int arg=0)
const char * getParam(const char *param, int arg=0)
int getIntParam(const char *param, int arg=0)

◆ run()

void ProgLocSharpening::run ( )
virtual

This function will be start running the program. it also should be implemented by derived classes.

Reimplemented from XmippProgram.

Definition at line 285 of file volume_local_sharpening.cpp.

286 {
287  produceSideInfo();
288 
289  MultidimArray<double> auxVol;
290  MultidimArray<double> operatedfiltered, Vk, filteredVol;
291  double lastnorm=0, lastporc=1;
292  double freq;
293  double step = 0.2;
294  int idx, bool1=1, bool2=1;
295  int lastidx = -1;
296 
297  minRes = 2*sampling;
298  maxRes=maxRes+2;
299 
300  //std::cout << "Resolutions between " << minRes << " and " << maxRes << std::endl;
301 
302  filteredVol = Vorig;
303 
305  double normOrig=0;
306 
307  MetaDataVec mditer;
308  size_t objId;
309  objId = mditer.addObject();
310 
311  for (size_t i = 1; i<=Niter; ++i)
312  {
313  //std::cout << "----------------Iteration " << i << "----------------" << std::endl;
314  mditer.setValue(MDL_ITER, (int) i, objId);
315  auxVol = filteredVol;
317 
318  localfiltering(fftV, operatedfiltered, minRes, maxRes, step);
319 
320  filteredVol = Vorig;
321 
322  filteredVol -= operatedfiltered;
323 
324  //calculate norm for Vorig
325  if (i==1)
326  {
328  {
330  }
331  normOrig = sqrt(normOrig);
332  //std::cout << "norma del original " << normOrig << std::endl;
333  }
334 
335 
336  //calculate norm for operatedfiltered
337  double norm=0;
339  {
340  norm +=(DIRECT_MULTIDIM_ELEM(operatedfiltered,n)*DIRECT_MULTIDIM_ELEM(operatedfiltered,n));
341  }
342  norm=sqrt(norm);
343 
344 
345  double porc=lastnorm*100/norm;
346  //std::cout << "norm " << norm << " percetage " << porc << std::endl;
347 
348  double subst=porc-lastporc;
349 
350  if ((subst<1)&&(bool1==1)&&(i>2))
351  {
352  bool1=2;
353  //std::cout << "-----iteration completed-----" << std::endl;
354 
355  }
356 
357  lastnorm=norm;
358  lastporc=porc;
359 
360  if (i==1 && lambda==1)
361  {
362  lambda=(normOrig/norm)/12;
363  std::cout << " lambda " << lambda << std::endl;
364  }
365 
367  transformer.FourierTransform(filteredVol, fftV);
368  localfiltering(fftV, filteredVol, minRes, maxRes, step);
369 
370  if (i == 1)
371  Vk = Vorig;
372  else
373  Vk = sharpenedMap;
374 
375  //sharpenedMap=Vk+lambda*(filteredVol);
377  {
379  lambda*DIRECT_MULTIDIM_ELEM(filteredVol,n);
380  //-0.01*DIRECT_MULTIDIM_ELEM(Vk,n)*SGN(DIRECT_MULTIDIM_ELEM(Vk,n));
383  }
384 
385 // double desv_sharp=0;
386 // computeAvgStdev_within_binary_mask(resVol, sharpenedMap, desv_sharp);
387 // std::cout << "desv_sharp = " << desv_sharp << std::endl;
388 
389  filteredVol = sharpenedMap;
390 
391  if (bool1==2)
392  {
393  Image<double> filteredvolume;
394  filteredvolume() = sharpenedMap;
395  filteredvolume.write(fnOut);
396  break;
397  }
398  }
399 
400  mditer.setValue(MDL_COST, lambda, objId);
401  mditer.write(fnMD);
402 
403  Image<double> filteredvolume;
404  filteredvolume() = sharpenedMap;
405  filteredvolume.write(fnOut);
406 
407 }
void resizeNoCopy(const MultidimArray< T1 > &v)
void sqrt(Image< double > &op)
MultidimArray< std::complex< double > > fftV
void write(const FileName &name="", size_t select_img=ALL_IMAGES, bool isStack=false, int mode=WRITE_OVERWRITE, CastWriteMode castMode=CW_CAST, int _swapWrite=0)
void write(const FileName &outFile, WriteModeMetaData mode=MD_OVERWRITE) const
T norm(const std::vector< T > &v)
Definition: vector_ops.h:399
#define i
bool setValue(const MDObject &mdValueIn, size_t id)
size_t addObject() override
Cost for the image (double)
void localfiltering(MultidimArray< std::complex< double > > &myfftV, MultidimArray< double > &localfilteredVol, double &minFreq, double &maxFreq, double &step)
#define FOR_ALL_DIRECT_ELEMENTS_IN_MULTIDIMARRAY(v)
#define DIRECT_MULTIDIM_ELEM(v, n)
void FourierTransform(T &v, T1 &V, bool getCopy=true)
Definition: xmipp_fftw.h:166
MultidimArray< double > Vorig
FourierTransformer transformer
MultidimArray< double > sharpenedMap
Current iteration number (int)
int * n

◆ sameEnergy()

void ProgLocSharpening::sameEnergy ( MultidimArray< std::complex< double > > &  myfftV,
MultidimArray< double > &  localfilteredVol,
double &  minFreq,
double &  maxFreq,
double &  step 
)

◆ wideBandPassFilter()

void ProgLocSharpening::wideBandPassFilter ( const MultidimArray< std::complex< double > > &  myfftV,
double  wmin,
double  wmax,
double  wL,
MultidimArray< double > &  filteredVol 
)

Member Data Documentation

◆ desv_Vorig

double ProgLocSharpening::desv_Vorig

Definition at line 56 of file volume_local_sharpening.h.

◆ desvOutside_Vorig

double ProgLocSharpening::desvOutside_Vorig

Definition at line 56 of file volume_local_sharpening.h.

◆ fftV

MultidimArray< std::complex<double> > ProgLocSharpening::fftV

Definition at line 99 of file volume_local_sharpening.h.

◆ fftVfilter

MultidimArray< std::complex<double> > ProgLocSharpening::fftVfilter

Definition at line 99 of file volume_local_sharpening.h.

◆ FilterBand

FourierFilter ProgLocSharpening::FilterBand

Definition at line 101 of file volume_local_sharpening.h.

◆ fnMD

FileName ProgLocSharpening::fnMD

Definition at line 53 of file volume_local_sharpening.h.

◆ fnOut

FileName ProgLocSharpening::fnOut

Filenames

Definition at line 53 of file volume_local_sharpening.h.

◆ fnRes

FileName ProgLocSharpening::fnRes

Definition at line 53 of file volume_local_sharpening.h.

◆ fnVol

FileName ProgLocSharpening::fnVol

Definition at line 53 of file volume_local_sharpening.h.

◆ iu

MultidimArray<double> ProgLocSharpening::iu

Definition at line 98 of file volume_local_sharpening.h.

◆ K

double ProgLocSharpening::K

Definition at line 56 of file volume_local_sharpening.h.

◆ lambda

double ProgLocSharpening::lambda

Definition at line 56 of file volume_local_sharpening.h.

◆ mask

MultidimArray<int> ProgLocSharpening::mask

Definition at line 96 of file volume_local_sharpening.h.

◆ maxFreq

double ProgLocSharpening::maxFreq

Definition at line 56 of file volume_local_sharpening.h.

◆ maxRes

double ProgLocSharpening::maxRes

Definition at line 56 of file volume_local_sharpening.h.

◆ minFreq

double ProgLocSharpening::minFreq

Definition at line 56 of file volume_local_sharpening.h.

◆ minRes

double ProgLocSharpening::minRes

Definition at line 56 of file volume_local_sharpening.h.

◆ Niter

int ProgLocSharpening::Niter

Definition at line 57 of file volume_local_sharpening.h.

◆ Nthread

int ProgLocSharpening::Nthread

Definition at line 57 of file volume_local_sharpening.h.

◆ resVol

MultidimArray<double> ProgLocSharpening::resVol

Definition at line 97 of file volume_local_sharpening.h.

◆ sampling

double ProgLocSharpening::sampling

sampling rate, minimum resolution, and maximum resolution

Definition at line 56 of file volume_local_sharpening.h.

◆ sharpenedMap

MultidimArray<double> ProgLocSharpening::sharpenedMap

Definition at line 98 of file volume_local_sharpening.h.

◆ transformer

FourierTransformer ProgLocSharpening::transformer

Definition at line 100 of file volume_local_sharpening.h.

◆ transformer_inv

FourierTransformer ProgLocSharpening::transformer_inv

Definition at line 100 of file volume_local_sharpening.h.

◆ Vorig

MultidimArray<double> ProgLocSharpening::Vorig

Definition at line 95 of file volume_local_sharpening.h.


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