Xmipp  v3.23.11-Nereus
pdb_sph_deform.cpp
Go to the documentation of this file.
1 /***************************************************************************
2  *
3  * Authors: David Herreros Calero (dherreros@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 #include "pdb_sph_deform.h"
27 #include <data/numerical_tools.h>
28 #include <fstream>
29 
31 {
32  addUsageLine("Deform a PDB according to a list of SPH deformation coefficients");
33  addParamsLine("--pdb <file> : PDB to deform");
34  addParamsLine("--clnm <metadata_file> : List of deformation coefficients");
35  addParamsLine("-o <file> : Deformed PDB");
36  addParamsLine("[--center_mass] : Center the PDB according with its center of mass");
37  addParamsLine("[--boxsize <b=0> ] : Box size of the volume where the PDB was fitted");
38  addParamsLine("[--sr <s=1> ] : Sampling rate of the volume where the PDB was fitted");
39  addExampleLine("xmipp_pdb_sph_deform --pdb 2tbv.pdb -o 2tbv_deformed.pdb --clnm coefficients.txt");
40 }
41 
43 {
44  fn_pdb = getParam("--pdb");
45  fn_sph = getParam("--clnm");
46  fn_out = getParam("-o");
47  center = checkParam("--center_mass");
48  boxSize= 0.5 * getDoubleParam("--boxsize") * getDoubleParam("--sr");
49 }
50 
52 {
53  if (verbose==0)
54  return;
55  std::cout
56  << "PDB: " << fn_pdb << std::endl
57  << "Coefficient list: " << fn_sph << std::endl
58  << "Output: " << fn_out << std::endl
59  ;
60 }
61 
63 {
64  PDBRichPhantom pdb;
65  pdb.read(fn_pdb);
66  std::string line;
67  line = readNthLine(0);
68  basisParams = string2vector(line);
69  line = readNthLine(1);
70  clnm = string2vector(line);
72  size_t idxY0=clnm.size()/3;
73  size_t idxZ0=2*idxY0;
75  cm.initZeros(3);
76  centerOfMass(pdb, cm);
77  for (size_t a=0; a<pdb.getNumberOfAtoms(); a++)
78  {
79  double gx=0.0;
80  double gy=0.0;
81  double gz=0.0;
82  RichAtom& atom_i=pdb.atomList[a];
83  auto k = atom_i.z;
84  auto i = atom_i.y;
85  auto j = atom_i.x;
86  if (center)
87  {
88  k -= VEC_ELEM(cm, 2);
89  i -= VEC_ELEM(cm, 1);
90  j -= VEC_ELEM(cm, 0);
91  }
92  else
93  {
94  k -= boxSize;
95  i -= boxSize;
96  j -= boxSize;
97  }
98  for (size_t idx=0; idx<idxY0; idx++)
99  {
100  double Rmax=basisParams[2];
101  double Rmax2=Rmax*Rmax;
102  double iRmax=1.0/Rmax;
103  double k2=k*k;
104  double kr=k*iRmax;
105  double k2i2=k2+i*i;
106  double ir=i*iRmax;
107  double r2=k2i2+j*j;
108  double jr=j*iRmax;
109  double rr=std::sqrt(r2)*iRmax;
110  double zsph=0.0;
111  if (r2<Rmax2)
112  {
113  int l1 = VEC_ELEM(vL1,idx);
114  int n = VEC_ELEM(vN,idx);
115  int l2 = VEC_ELEM(vL2,idx);
116  int m = VEC_ELEM(vM,idx);
117  zsph=ZernikeSphericalHarmonics(l1,n,l2,m,jr,ir,kr,rr);
118  if (rr>0 || l2==0)
119  {
120  gx += clnm[idx] *zsph;
121  gy += clnm[idx+idxY0] *zsph;
122  gz += clnm[idx+idxZ0] *zsph;
123  }
124  }
125  }
126  atom_i.x += gx;
127  atom_i.y += gy;
128  atom_i.z += gz;
129  }
130  pdb.write(fn_out);
131 }
132 
133 std::string ProgPdbSphDeform::readNthLine(int N) const
134 {
135  std::ifstream in(fn_sph.getString());
136  std::string s;
137 
138  //skip N lines
139  for(int i = 0; i < N; ++i)
140  std::getline(in, s);
141 
142  std::getline(in,s);
143  return s;
144 }
145 
146 std::vector<double> ProgPdbSphDeform::string2vector(std::string const &s) const
147 {
148  std::stringstream iss(s);
149  double number;
150  std::vector<double> v;
151  while (iss >> number)
152  v.push_back(number);
153  return v;
154 }
155 
157 {
158  int idx = 0;
159  auto vecSize = (int)(clnm.size()/3);
160  vL1.initZeros(vecSize);
161  vN.initZeros(vecSize);
162  vL2.initZeros(vecSize);
163  vM.initZeros(vecSize);
164  for (int h=0; h<=basisParams[1]; h++)
165  {
166  int totalSPH = 2*h+1;
167  auto aux = (int)(std::floor(totalSPH/2));
168  for (int l=h; l<=basisParams[0]; l+=2)
169  {
170  for (int m=0; m<totalSPH; m++)
171  {
172  VEC_ELEM(vL1,idx) = l;
173  VEC_ELEM(vN,idx) = h;
174  VEC_ELEM(vL2,idx) = h;
175  VEC_ELEM(vM,idx) = m-aux;
176  idx++;
177  }
178  }
179  }
180 }
181 
183 {
184  size_t numAtoms = pdb.getNumberOfAtoms();
185  for (size_t a=0; a<numAtoms; a++)
186  {
187  RichAtom& atom_i=pdb.atomList[a];
188  auto z = atom_i.z;
189  auto y = atom_i.y;
190  auto x = atom_i.x;
191  VEC_ELEM(cm, 0) += x;
192  VEC_ELEM(cm, 1) += y;
193  VEC_ELEM(cm, 2) += z;
194 
195  }
196  VEC_ELEM(cm, 0) /= numAtoms;
197  VEC_ELEM(cm, 1) /= numAtoms;
198  VEC_ELEM(cm, 2) /= numAtoms;
199 }
void defineParams() override
void read(const FileName &fnPDB, const bool pseudoatoms=false, const double threshold=0.0)
Read rich phantom from either a PDB of CIF file.
Definition: pdb.cpp:704
#define VEC_ELEM(v, i)
Definition: matrix1d.h:245
double getDoubleParam(const char *param, int arg=0)
__host__ __device__ float2 floor(const float2 v)
void write(const FileName &fnPDB, const bool renumber=false)
Write rich phantom to PDB or CIF file.
Definition: pdb.cpp:872
void sqrt(Image< double > &op)
static double * y
Definition: pdb.h:174
std::string readNthLine(int N) const
void centerOfMass(PDBRichPhantom pdb, Matrix1D< double > &cm)
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 show() const override
const char * getParam(const char *param, int arg=0)
double y
Position Y.
Definition: pdb.h:187
int in
std::vector< double > string2vector(std::string const &s) const
double z
void addExampleLine(const char *example, bool verbatim=true)
void run() override
int verbose
Verbosity level.
void initZeros()
Definition: matrix1d.h:592
void readParams() override
double x
Position X.
Definition: pdb.h:184
#define j
int m
std::vector< RichAtom > atomList
List of atoms.
Definition: pdb.h:257
size_t getNumberOfAtoms() const
Get number of atoms.
Definition: pdb.h:270
double ZernikeSphericalHarmonics(int n, int m, double xr, double yr, double zr, double r)
String getString() const
bool checkParam(const char *param)
float r2
void addUsageLine(const char *line, bool verbatim=false)
int * n
doublereal * a
void addParamsLine(const String &line)
int ir
double z
Position Z.
Definition: pdb.h:190