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

#include <volume_apply_deform_sph.h>

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

Public Member Functions

void defineParams () override
 
void readParams () override
 
void show () const override
 
void run () override
 
std::string readNthLine (int N) const
 
std::vector< double > string2vector (std::string const &s) const
 
void fillVectorTerms ()
 
- 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 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 ()
 

Additional Inherited Members

- Public Attributes inherited from XmippProgram
bool doRun
 
bool runWithoutArgs
 
int verbose
 Verbosity level. More...
 
int debug
 
- 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

Definition at line 31 of file volume_apply_deform_sph.h.

Member Function Documentation

◆ defineParams()

void ProgApplyVolDeformSph::defineParams ( )
overridevirtual

Params definitions

Reimplemented from XmippProgram.

Definition at line 32 of file volume_apply_deform_sph.cpp.

33 {
34  addUsageLine("Deform a PDB according to a list of SPH deformation coefficients");
35  addParamsLine("-i <volume> : Volume to deform");
36  addParamsLine("--clnm <metadata_file> : List of deformation coefficients");
37  addParamsLine("-o <volume> : Deformed volume");
38  addExampleLine("xmipp_apply_deform_sph -i input.vol -o volume_deformed.vol --clnm coefficients.txt");
39 }
void addExampleLine(const char *example, bool verbatim=true)
void addUsageLine(const char *line, bool verbatim=false)
void addParamsLine(const String &line)

◆ fillVectorTerms()

void ProgApplyVolDeformSph::fillVectorTerms ( )

Fill degree and order vectors

Definition at line 145 of file volume_apply_deform_sph.cpp.

146 {
147  int idx = 0;
148  auto vecSize = (int)((clnm.size())/3);
149  vL1.initZeros(vecSize);
150  vN.initZeros(vecSize);
151  vL2.initZeros(vecSize);
152  vM.initZeros(vecSize);
153  for (int h=0; h<=basisParams[1]; h++)
154  {
155  int totalSPH = 2*h+1;
156  auto aux = (int)(std::floor(totalSPH/2));
157  for (int l=h; l<=basisParams[0]; l+=2)
158  {
159  for (int m=0; m<totalSPH; m++)
160  {
161  VEC_ELEM(vL1,idx) = l;
162  VEC_ELEM(vN,idx) = h;
163  VEC_ELEM(vL2,idx) = h;
164  VEC_ELEM(vM,idx) = m-aux;
165  idx++;
166  }
167  }
168  }
169 }
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
__host__ __device__ float2 floor(const float2 v)
void initZeros()
Definition: matrix1d.h:592
int m

◆ readNthLine()

std::string ProgApplyVolDeformSph::readNthLine ( int  N) const

Read Nth line of file

Definition at line 122 of file volume_apply_deform_sph.cpp.

123 {
124  std::ifstream in(fn_sph.getString());
125  std::string s;
126 
127  //skip N lines
128  for(int i = 0; i < N; ++i)
129  std::getline(in, s);
130 
131  std::getline(in,s);
132  return s;
133 }
#define i
int in
String getString() const

◆ readParams()

void ProgApplyVolDeformSph::readParams ( )
overridevirtual

Read from a command line.

Reimplemented from XmippProgram.

Definition at line 41 of file volume_apply_deform_sph.cpp.

42 {
43  fn_vol=getParam("-i");
44  fn_sph=getParam("--clnm");
45  fn_out=getParam("-o");
46 }
const char * getParam(const char *param, int arg=0)

◆ run()

void ProgApplyVolDeformSph::run ( )
overridevirtual

Run.

Reimplemented from XmippProgram.

Definition at line 59 of file volume_apply_deform_sph.cpp.

60 {
61  Image<double> VI;
62  Image<double> VO;
63  VI.read(fn_vol);
64  VI().setXmippOrigin();
65  VO().initZeros(VI());
66  VO().setXmippOrigin();
67  std::string line;
68  line = readNthLine(0);
69  basisParams = string2vector(line);
70  line = readNthLine(1);
71  clnm = string2vector(line);
73  size_t idxY0=(clnm.size())/3;
74  size_t idxZ0=2*idxY0;
75  const MultidimArray<double> &mVI=VI();
76  double voxelI=0.0;
77  double Rmax=basisParams[2];
78  double Rmax2=Rmax*Rmax;
79  double iRmax=1.0/Rmax;
80  double zsph=0.0;
81  for (int k=STARTINGZ(mVI); k<=FINISHINGZ(mVI); k++)
82  {
83  for (int i=STARTINGY(mVI); i<=FINISHINGY(mVI); i++)
84  {
85  for (int j=STARTINGX(mVI); j<=FINISHINGX(mVI); j++)
86  {
87  double gx=0.0;
88  double gy=0.0;
89  double gz=0.0;
90  double k2=k*k;
91  double kr=k*iRmax;
92  double k2i2=k2+i*i;
93  double ir=i*iRmax;
94  double r2=k2i2+j*j;
95  double jr=j*iRmax;
96  double rr=std::sqrt(r2)*iRmax;
97  if (r2<Rmax2)
98  {
99  for (size_t idx=0; idx<idxY0; idx++)
100  {
101  int l1 = VEC_ELEM(vL1,idx);
102  int n = VEC_ELEM(vN,idx);
103  int l2 = VEC_ELEM(vL2,idx);
104  int m = VEC_ELEM(vM,idx);
105  zsph=ZernikeSphericalHarmonics(l1,n,l2,m,jr,ir,kr,rr);
106  if (rr>0 && l2!=0)
107  {
108  gx += clnm[idx] *zsph;
109  gy += clnm[idx+idxY0] *zsph;
110  gz += clnm[idx+idxZ0] *zsph;
111  }
112  }
113  }
114  voxelI=mVI.interpolatedElement3D(j+gx,i+gy,k+gz);
115  A3D_ELEM(VO(), k, i, j)=voxelI;
116  }
117  }
118  }
119  VO.write(fn_out);
120 }
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
#define FINISHINGX(v)
void sqrt(Image< double > &op)
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)
#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)
#define A3D_ELEM(V, k, i, j)
std::string readNthLine(int N) const
T interpolatedElement3D(double x, double y, double z, T outside_value=(T) 0) const
#define j
int m
#define FINISHINGY(v)
double ZernikeSphericalHarmonics(int n, int m, double xr, double yr, double zr, double r)
int read(const FileName &name, DataMode datamode=DATA, size_t select_img=ALL_IMAGES, bool mapData=false, int mode=WRITE_READONLY)
float r2
std::vector< double > string2vector(std::string const &s) const
#define STARTINGZ(v)
int * n
int ir

◆ show()

void ProgApplyVolDeformSph::show ( ) const
overridevirtual

Show parameters.

Reimplemented from XmippProgram.

Definition at line 48 of file volume_apply_deform_sph.cpp.

49 {
50  if (verbose==0)
51  return;
52  std::cout
53  << "Volume: " << fn_vol << std::endl
54  << "Coefficient list: " << fn_sph << std::endl
55  << "Output: " << fn_out << std::endl
56  ;
57 }
int verbose
Verbosity level.

◆ string2vector()

std::vector< double > ProgApplyVolDeformSph::string2vector ( std::string const &  s) const

Convert String to Vector

Definition at line 135 of file volume_apply_deform_sph.cpp.

136 {
137  std::stringstream iss(s);
138  double number;
139  std::vector<double> v;
140  while (iss >> number)
141  v.push_back(number);
142  return v;
143 }

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