ProSHADE  0.7.6.6 (JUL 2022)
Protein Shape Detection
ProSHADE_mapManip.cpp
Go to the documentation of this file.
1 
22 //==================================================== ProSHADE
23 #include "ProSHADE_mapManip.hpp"
24 
25 //==================================================== Define round for C++98
31 proshade_signed ProSHADE_internal_mapManip::myRound ( proshade_double x )
32 {
33 #if __cplusplus >= 201103L
34  return ( static_cast< proshade_signed > ( std::round ( x ) ) );
35 #else
36  return ( static_cast< proshade_signed > ( round ( x ) ) );
37 #endif
38 }
39 
40 //==================================================== Define round for C++98
46 proshade_signed ProSHADE_internal_mapManip::myRound ( proshade_single x )
47 {
48 #if __cplusplus >= 201103L
49  return ( static_cast< proshade_signed > ( std::round ( x ) ) );
50 #else
51  return ( static_cast< proshade_signed > ( round ( x ) ) );
52 #endif
53 }
54 
72 void ProSHADE_internal_mapManip::determinePDBRanges ( gemmi::Structure pdbFile, proshade_single* xFrom, proshade_single* xTo, proshade_single* yFrom, proshade_single* yTo, proshade_single* zFrom, proshade_single* zTo, bool firstModel )
73 {
74  //================================================ Initialise structure crawl
75  bool firstAtom = true;
76 
77  //================================================ Use the first model, if it exists
78  if ( pdbFile.models.size() > 0 )
79  {
80  //============================================ For each model
81  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile.models.size() ); sIt++ )
82  {
83  //======================================== Check if multiple models are allowed
84  if ( firstModel && ( sIt != 0 ) ) { break; }
85 
86  //======================================== Get model
87  gemmi::Model model = pdbFile.models.at(sIt);
88 
89  //======================================== For each chain
90  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model.chains.size() ); mIt++ )
91  {
92  //==================================== Get chain
93  gemmi::Chain chain = model.chains.at(mIt);
94 
95  //==================================== For each residue
96  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( chain.residues.size() ); rIt++ )
97  {
98  //================================ Get residue
99  gemmi::Residue residue = chain.residues.at(rIt);
100 
101  //================================ For each atom
102  for ( proshade_unsign aIt = 0; aIt < static_cast<proshade_unsign> ( residue.atoms.size() ); aIt++ )
103  {
104  //============================ Get atom
105  gemmi::Atom atom = residue.atoms.at(aIt);
106 
107  //============================ Ignore hydrogens, map computations ignore them anyway and inclusion here causes map - co-ordinate mismatches.
108  if ( atom.is_hydrogen() ) { continue; }
109 
110  //============================ Find the coordinate ranges
111  if ( firstAtom )
112  {
113  *xTo = static_cast<proshade_single> ( atom.pos.x );
114  *xFrom = static_cast<proshade_single> ( atom.pos.x );
115  *yTo = static_cast<proshade_single> ( atom.pos.y );
116  *yFrom = static_cast<proshade_single> ( atom.pos.y );
117  *zTo = static_cast<proshade_single> ( atom.pos.z );
118  *zFrom = static_cast<proshade_single> ( atom.pos.z );
119  firstAtom = false;
120  }
121  else
122  {
123  if ( static_cast<proshade_single> ( atom.pos.x ) > *xTo ) { *xTo = static_cast<proshade_single> ( atom.pos.x ); }
124  if ( static_cast<proshade_single> ( atom.pos.x ) < *xFrom ) { *xFrom = static_cast<proshade_single> ( atom.pos.x ); }
125  if ( static_cast<proshade_single> ( atom.pos.y ) > *yTo ) { *yTo = static_cast<proshade_single> ( atom.pos.y ); }
126  if ( static_cast<proshade_single> ( atom.pos.y ) < *yFrom ) { *yFrom = static_cast<proshade_single> ( atom.pos.y ); }
127  if ( static_cast<proshade_single> ( atom.pos.z ) > *zTo ) { *zTo = static_cast<proshade_single> ( atom.pos.z ); }
128  if ( static_cast<proshade_single> ( atom.pos.z ) < *zFrom ) { *zFrom = static_cast<proshade_single> ( atom.pos.z ); }
129  }
130  }
131  }
132  }
133  }
134  }
135  else
136  {
137  std::stringstream hlpSS;
138  hlpSS << "Found 0 models in input file " << pdbFile.name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
139  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
140  }
141 
142  //================================================ Done
143  return ;
144 
145 }
146 
157 void ProSHADE_internal_mapManip::findPDBCOMValues ( gemmi::Structure* pdbFile, proshade_double *xCom, proshade_double *yCom, proshade_double *zCom, bool firstModel )
158 {
159  //================================================ Initialise structure crawl
160  proshade_double totAtoms = 0.0;
161  *xCom = 0.0;
162  *yCom = 0.0;
163  *zCom = 0.0;
164 
165  //================================================ Use the first model, if it exists
166  if ( pdbFile->models.size() > 0 )
167  {
168  //============================================ For each model
169  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile->models.size() ); sIt++ )
170  {
171  //======================================== Get model
172  gemmi::Model model = pdbFile->models.at(sIt);
173 
174  //======================================== Check if multiple models are allowed
175  if ( firstModel && ( sIt != 0 ) ) { break; }
176 
177  //======================================== For each chain
178  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model.chains.size() ); mIt++ )
179  {
180  //==================================== Get chain
181  gemmi::Chain chain = model.chains.at(mIt);
182 
183  gemmi::ResidueSpan span = chain.whole ( );
184 
185  //==================================== For each residue
186  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( span.size() ); rIt++ )
187  {
188  //================================ Get residue
189  gemmi::Residue residue = span.at(rIt);
190 
191  //================================ For each atom
192  for ( proshade_unsign aIt = 0; aIt < static_cast<proshade_unsign> ( residue.atoms.size() ); aIt++ )
193  {
194  //============================ Get atom
195  gemmi::Atom atom = residue.atoms.at(aIt);
196 
197  //============================ Save the COM sums
198  *xCom += atom.pos.x * atom.element.weight();
199  *yCom += atom.pos.y * atom.element.weight();
200  *zCom += atom.pos.z * atom.element.weight();
201  totAtoms += atom.element.weight();
202  }
203  }
204  }
205  }
206  }
207  else
208  {
209  std::stringstream hlpSS;
210  hlpSS << "Found 0 models in input file " << pdbFile->name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
211  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
212  }
213 
214  //================================================ Normalise sums to COM
215  *xCom /= totAtoms;
216  *yCom /= totAtoms;
217  *zCom /= totAtoms;
218 
219  //================================================ Done
220  return ;
221 
222 }
223 
243 void ProSHADE_internal_mapManip::findMAPCOMValues ( proshade_double* map, proshade_double *xCom, proshade_double *yCom, proshade_double *zCom, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xFrom, proshade_signed xTo, proshade_signed yFrom, proshade_signed yTo, proshade_signed zFrom, proshade_signed zTo, bool removeNegDens )
244 {
245  //================================================ Initialise computation
246  proshade_double totDensity = 0.0;
247  *xCom = 0.0;
248  *yCom = 0.0;
249  *zCom = 0.0;
250  proshade_signed arrPos = 0;
251  proshade_single xSampRate = xAngs / static_cast< proshade_single > ( xTo - xFrom + 1 );
252  proshade_single ySampRate = yAngs / static_cast< proshade_single > ( yTo - yFrom + 1 );
253  proshade_single zSampRate = zAngs / static_cast< proshade_single > ( zTo - zFrom + 1 );
254 
255  //================================================ For each map point
256  for ( proshade_signed xIt = 0; xIt <= ( xTo - xFrom ); xIt++ )
257  {
258  for ( proshade_signed yIt = 0; yIt <= ( yTo - yFrom ); yIt++ )
259  {
260  for ( proshade_signed zIt = 0; zIt <= ( zTo - zFrom ); zIt++ )
261  {
262  //==================================== Find array position
263  arrPos = zIt + ( zTo - zFrom + 1 ) * ( yIt + ( yTo - yFrom + 1 ) * xIt );
264 
265  //==================================== Ignore NaN's and NA's
266  const FloatingPoint< proshade_double > lhs ( map[arrPos] );
267  if ( !lhs.AlmostEquals ( lhs ) ) { map[arrPos] = 0.0; continue; }
268 
269  //==================================== If negative density is to be ignored, do so
270  if ( ( map[arrPos] < 0.0 ) && removeNegDens ) { continue; }
271 
272  //==================================== Compute position times weight
273  totDensity += map[arrPos];
274  *xCom += static_cast<proshade_double> ( xIt * map[arrPos] );
275  *yCom += static_cast<proshade_double> ( yIt * map[arrPos] );
276  *zCom += static_cast<proshade_double> ( zIt * map[arrPos] );
277  }
278  }
279  }
280 
281  //================================================ Divide by total weight
282  *xCom /= totDensity;
283  *yCom /= totDensity;
284  *zCom /= totDensity;
285 
286  //================================================ Convert from indices to Angstroms and into "real world" co-ordinate system.
287  *xCom = ( static_cast< proshade_double > ( xFrom ) + (*xCom) ) * static_cast< proshade_double > ( xSampRate );
288  *yCom = ( static_cast< proshade_double > ( yFrom ) + (*yCom) ) * static_cast< proshade_double > ( ySampRate );
289  *zCom = ( static_cast< proshade_double > ( zFrom ) + (*zCom) ) * static_cast< proshade_double > ( zSampRate );
290 
291  //================================================ Done
292  return ;
293 
294 }
295 
310 void ProSHADE_internal_mapManip::rotatePDBCoordinates ( gemmi::Structure *pdbFile, proshade_double euA, proshade_double euB, proshade_double euG, proshade_double xCom,
311 proshade_double yCom, proshade_double zCom, bool firstModel )
312 {
313  //================================================ Convert Euler angles to rotation matrix
314  proshade_double *rotMat = new proshade_double[9];
315  ProSHADE_internal_misc::checkMemoryAllocation ( rotMat, __FILE__, __LINE__, __func__ );
317 
318  //================================================ Initialise internal variables
319  proshade_double xTmp, yTmp, zTmp;
320 
321  //================================================ Use the first model, if it exists
322  if ( pdbFile->models.size() > 0 )
323  {
324  //============================================ For each model
325  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile->models.size() ); sIt++ )
326  {
327  //======================================== Get model
328  gemmi::Model *model = &pdbFile->models.at(sIt);
329 
330  //======================================== Check if multiple models are allowed
331  if ( firstModel && ( sIt != 0 ) ) { break; }
332 
333  //======================================== For each chain
334  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model->chains.size() ); mIt++ )
335  {
336  //==================================== Get chain
337  gemmi::Chain *chain = &model->chains.at(mIt);
338 
339  //==================================== For each residue
340  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( chain->residues.size() ); rIt++ )
341  {
342  //================================ Get residue
343  gemmi::Residue *residue = &chain->residues.at(rIt);
344 
345  //================================ For each atom
346  for ( proshade_unsign aIt = 0; aIt < static_cast<proshade_unsign> ( residue->atoms.size() ); aIt++ )
347  {
348  //============================ Get atom
349  gemmi::Atom *atom = &residue->atoms.at(aIt);
350 
351  //============================ Move to mid-point
352  xTmp = static_cast< proshade_double > ( atom->pos.x - xCom );
353  yTmp = static_cast< proshade_double > ( atom->pos.y - yCom );
354  zTmp = static_cast< proshade_double > ( atom->pos.z - zCom );
355 
356  //============================ Rotate the atom position
357  atom->pos.x = ( xTmp * rotMat[0] ) + ( yTmp * rotMat[1] ) + ( zTmp * rotMat[2] );
358  atom->pos.y = ( xTmp * rotMat[3] ) + ( yTmp * rotMat[4] ) + ( zTmp * rotMat[5] );
359  atom->pos.z = ( xTmp * rotMat[6] ) + ( yTmp * rotMat[7] ) + ( zTmp * rotMat[8] );
360 
361  //============================ Move back
362  atom->pos.x = atom->pos.x + xCom;
363  atom->pos.y = atom->pos.y + yCom;
364  atom->pos.z = atom->pos.z + zCom;
365  }
366  }
367  }
368  }
369  }
370  else
371  {
372  std::stringstream hlpSS;
373  hlpSS << "Found 0 models in input file " << pdbFile->name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
374  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
375  }
376 
377  //================================================ Release memory
378  delete[] rotMat;
379 
380  //================================================ Done
381  return ;
382 
383 }
384 
395 void ProSHADE_internal_mapManip::translatePDBCoordinates ( gemmi::Structure *pdbFile, proshade_double transX, proshade_double transY, proshade_double transZ, bool firstModel )
396 {
397  //================================================ Use the first model, if it exists
398  if ( pdbFile->models.size() > 0 )
399  {
400  //============================================ For each model
401  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile->models.size() ); sIt++ )
402  {
403  //======================================== Check if multiple models are allowed
404  if ( firstModel && ( sIt != 0 ) ) { break; }
405 
406  //======================================== Get model
407  gemmi::Model *model = &pdbFile->models.at(sIt);
408 
409  //======================================== For each chain
410  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model->chains.size() ); mIt++ )
411  {
412  //==================================== Get chain
413  gemmi::Chain *chain = &model->chains.at(mIt);
414 
415  //==================================== For each residue
416  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( chain->residues.size() ); rIt++ )
417  {
418  //================================ Get residue
419  gemmi::Residue *residue = &chain->residues.at(rIt);
420 
421  //================================ For each atom
422  for ( proshade_unsign aIt = 0; aIt < static_cast<proshade_unsign> ( residue->atoms.size() ); aIt++ )
423  {
424  //============================ Get atom
425  gemmi::Atom *atom = &residue->atoms.at(aIt);
426 
427  //============================ Translate
428  atom->pos.x += transX;
429  atom->pos.y += transY;
430  atom->pos.z += transZ;
431  }
432  }
433  }
434  }
435  }
436  else
437  {
438  std::stringstream hlpSS;
439  hlpSS << "Found 0 models in input file " << pdbFile->name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
440  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
441  }
442 
443  //================================================ Done
444  return ;
445 
446 }
447 
458 void ProSHADE_internal_mapManip::changePDBBFactors ( gemmi::Structure *pdbFile, proshade_double newBFactorValue, bool firstModel )
459 {
460  //================================================ Use the first model, if it exists
461  if ( pdbFile->models.size() > 0 )
462  {
463  //============================================ For each model
464  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile->models.size() ); sIt++ )
465  {
466  //======================================== Check if multiple models are allowed
467  if ( firstModel && ( sIt != 0 ) ) { break; }
468 
469  //======================================== Get model
470  gemmi::Model *model = &pdbFile->models.at(sIt);
471 
472  //======================================== For each chain
473  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model->chains.size() ); mIt++ )
474  {
475  //==================================== Get chain
476  gemmi::Chain *chain = &model->chains.at(mIt);
477 
478  //==================================== For each residue
479  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( chain->residues.size() ); rIt++ )
480  {
481  //================================ Get residue
482  gemmi::Residue *residue = &chain->residues.at(rIt);
483 
484  //================================ For each atom
485  for ( proshade_unsign aIt = 0; aIt < static_cast<proshade_unsign> ( residue->atoms.size() ); aIt++ )
486  {
487  //============================ Get atom
488  gemmi::Atom *atom = &residue->atoms.at(aIt);
489 
490  //============================ Change the B-factors
491  atom->b_iso = static_cast< float > ( newBFactorValue );
492  }
493  }
494  }
495  }
496  }
497  else
498  {
499  std::stringstream hlpSS;
500  hlpSS << "Found 0 models in input file " << pdbFile->name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
501  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
502  }
503 
504  //================================================ Done
505  return ;
506 
507 }
508 
518 void ProSHADE_internal_mapManip::removeWaters ( gemmi::Structure *pdbFile, bool firstModel )
519 {
520  //================================================ Use the first model, if it exists
521  if ( pdbFile->models.size() > 0 )
522  {
523  //============================================ For each model
524  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile->models.size() ); sIt++ )
525  {
526  //======================================== Check if multiple models are allowed
527  if ( firstModel && ( sIt != 0 ) ) { break; }
528 
529  //======================================== Get model
530  gemmi::Model *model = &pdbFile->models.at(sIt);
531 
532  //======================================== For each chain
533  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model->chains.size() ); mIt++ )
534  {
535  //==================================== Get chain
536  gemmi::Chain *chain = &model->chains.at(mIt);
537 
538  //==================================== Initialise del vector
539  std::vector< proshade_unsign > delVec;
540 
541  //==================================== For each residue
542  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( chain->residues.size() ); rIt++ )
543  {
544  //================================ Get residue
545  gemmi::Residue *residue = &chain->residues.at(rIt);
546 
547  //================================ If residue is water
548  if ( residue->is_water() )
549  {
551  }
552  }
553 
554  //==================================== Delete from end to avoid indexing issues
555  std::sort ( delVec.begin(), delVec.end(), std::greater<int>() );
556  for ( proshade_unsign vecIt = 0; vecIt < static_cast<proshade_unsign> ( delVec.size() ); vecIt++ )
557  {
558  chain->residues.erase ( chain->residues.begin() + static_cast< long int > ( delVec.at(vecIt) ) );
559  }
560  }
561  }
562  }
563  else
564  {
565  std::stringstream hlpSS;
566  hlpSS << "Found 0 models in input file " << pdbFile->name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
567  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
568  }
569 
570  //================================================ Done
571  return ;
572 
573 }
574 
587 void ProSHADE_internal_mapManip::movePDBForMapCalc ( gemmi::Structure *pdbFile, proshade_single xMov, proshade_single yMov, proshade_single zMov, bool firstModel )
588 {
589  //================================================ Use the first model, if it exists
590  if ( pdbFile->models.size() > 0 )
591  {
592  //============================================ For each model
593  for ( proshade_unsign sIt = 0; sIt < static_cast<proshade_unsign> ( pdbFile->models.size() ); sIt++ )
594  {
595  //======================================== Check if multiple models are allowed
596  if ( firstModel && ( sIt != 0 ) ) { break; }
597 
598  //======================================== Get model
599  gemmi::Model *model = &pdbFile->models.at(sIt);
600 
601  //======================================== For each chain
602  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( model->chains.size() ); mIt++ )
603  {
604  //==================================== Get chain
605  gemmi::Chain *chain = &model->chains.at(mIt);
606 
607  //==================================== For each residue
608  for ( proshade_unsign rIt = 0; rIt < static_cast<proshade_unsign> ( chain->residues.size() ); rIt++ )
609  {
610  //================================ Get residue
611  gemmi::Residue *residue = &chain->residues.at(rIt);
612 
613  //================================ For each atom
614  for ( proshade_unsign aIt = 0; aIt < static_cast<proshade_unsign> ( residue->atoms.size() ); aIt++ )
615  {
616  //============================ Get atom
617  gemmi::Atom *atom = &residue->atoms.at(aIt);
618 
619  //============================ Move the atoms
620  atom->pos = gemmi::Position ( atom->pos.x + static_cast< proshade_double > ( xMov ), atom->pos.y + static_cast< proshade_double > ( yMov ), atom->pos.z + static_cast< proshade_double > ( zMov ) );
621  }
622  }
623  }
624 
625  }
626  }
627  else
628  {
629  std::stringstream hlpSS;
630  hlpSS << "Found 0 models in input file " << pdbFile->name << ".\n : This suggests that the input co-ordinate file is\n : corrupted or mis-formatted.";
631  throw ProSHADE_exception ( "Found no model in co-ordinate file.", "EP00050", __FILE__, __LINE__, __func__, hlpSS.str() );
632  }
633 
634  //================================================ Done
635  return ;
636 }
637 
658 void ProSHADE_internal_mapManip::generateMapFromPDB ( gemmi::Structure pdbFile, proshade_double*& map, proshade_single requestedResolution, proshade_single xCell, proshade_single yCell, proshade_single zCell, proshade_signed* xTo, proshade_signed* yTo, proshade_signed* zTo, bool forceP1, bool firstModel )
659 {
660  //================================================ Set cell dimensions from the increased ranges (we need to add some space) and re-calculate cell properties
661  if ( forceP1 ) { pdbFile.cell = gemmi::UnitCell(); }
662  pdbFile.cell.a = static_cast< proshade_double > ( xCell );
663  pdbFile.cell.b = static_cast< proshade_double > ( yCell );
664  pdbFile.cell.c = static_cast< proshade_double > ( zCell );
665  pdbFile.cell.calculate_properties ( );
666 
667  //================================================ Get elements in Gemmi format
668  std::string totElString;
669  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( pdbFile.models.size() ); mIt++ )
670  {
671  //============================================ Check if multiple models are allowed
672  if ( firstModel && ( mIt != 0 ) )
673  {
674  std::stringstream hlpSS;
675  hlpSS << "!!! ProSHADE WARNING !!! Found multiple models (" << pdbFile.models.size() << ") in input file " << pdbFile.name << ", while the settings state that only the first PDB file model should be used. If all models should be used, please supply ProSHADE with the \"-x\" option.";
676  ProSHADE_internal_messages::printWarningMessage ( 0, hlpSS.str(), "WP00055" );
677  break;
678  }
679 
680  std::string hlpStr = pdbFile.models[mIt].present_elements ( ).to_string<char,std::char_traits<char>,std::allocator<char> >();
681  totElString = totElString + hlpStr;
682  }
683  std::bitset< static_cast< size_t > ( gemmi::El::END )> present_elems ( totElString );
684 
685  //================================================ Sanity checks
686  if ( present_elems[static_cast<int> ( gemmi::El::X )] )
687  {
688  throw ProSHADE_exception ( "Found unknown element in input file.", "EP00051", __FILE__, __LINE__, __func__, "Gemmi library does not recognise some of the elements in\n : the co-ordinate file. Please check the file for not being\n : corrupted and containing standard elements." );
689  }
690 
691  for ( proshade_unsign elIt = 0; elIt < static_cast<proshade_unsign> ( present_elems.size() ); elIt++ )
692  {
693  if ( present_elems[elIt] && !gemmi::IT92<double>::has ( static_cast<gemmi::El> ( elIt ) ) )
694  {
695  std::stringstream hlpSS;
696  hlpSS << "Missing form factor for element " << element_name ( static_cast<gemmi::El> ( elIt ) );
697  throw ProSHADE_exception ( hlpSS.str().c_str(), "EP00052", __FILE__, __LINE__, __func__, "Gemmi library does not have a form factor value for this\n : reported element. Please report this to the author." );
698  }
699  }
700 
701  //================================================ Compute the f's
702  double wavelength = 10.0;
703  double energy = gemmi::hc() / wavelength;
704 
705  //================================================ Create the density calculator object and fill it in
706  gemmi::DensityCalculator<gemmi::IT92<double>, float> dencalc;
707 
708  dencalc.d_min = static_cast< double > ( requestedResolution );
709  for ( size_t elIt = 0; elIt < present_elems.size(); elIt++ ) { if ( present_elems[elIt] ) { dencalc.addends.set ( static_cast< gemmi::El > ( elIt ), static_cast< float > ( gemmi::cromer_liberman ( static_cast< int > ( elIt ), energy, nullptr ) ) ); } }
710  dencalc.set_grid_cell_and_spacegroup ( pdbFile );
711 
712  //================================================ Force P1 spacegroup
713  if ( forceP1 ) { dencalc.grid.spacegroup = &gemmi::get_spacegroup_p1(); }
714 
715  //================================================ Compute the theoretical map for each model
716  dencalc.grid.data.clear ( );
717  dencalc.grid.set_size_from_spacing ( dencalc.d_min / ( 2.0 * dencalc.rate), true );
718  for ( proshade_unsign mIt = 0; mIt < static_cast<proshade_unsign> ( pdbFile.models.size() ); mIt++ )
719  {
720  if ( firstModel && ( mIt != 0 ) ) { break; }
721  dencalc.add_model_density_to_grid ( pdbFile.models[mIt] );
722  dencalc.grid.symmetrize ( [](float a, float b) { return a + b; } );
723  }
724 
725  //================================================ Get the map
726  const gemmi::Grid<float>& grid = dencalc.grid;
727 
728  //================================================ Save the map dimensions
729  *xTo = grid.nu;
730  *yTo = grid.nv;
731  *zTo = grid.nw;
732 
733  //================================================ Copy the gemmi::Grid to my map format
734  map = new proshade_double [(*xTo) * (*yTo) * (*zTo)];
735  ProSHADE_internal_misc::checkMemoryAllocation ( map, __FILE__, __LINE__, __func__ );
736 
737  proshade_signed arrPos = 0;
738  for ( proshade_signed uIt = 0; uIt < (*xTo); uIt++ )
739  {
740  for ( proshade_signed vIt = 0; vIt < (*yTo); vIt++ )
741  {
742  for ( proshade_signed wIt = 0; wIt < (*zTo); wIt++ )
743  {
744  arrPos = wIt + (*zTo) * ( vIt + (*yTo) * uIt );
745  map[arrPos] = static_cast< proshade_double > ( grid.get_value_q( static_cast< int > ( uIt ), static_cast< int > ( vIt ), static_cast< int > ( wIt ) ) );
746  }
747  }
748  }
749 
750  //================================================ Done
751  return ;
752 
753 }
754 
777 void ProSHADE_internal_mapManip::moveMapByIndices ( proshade_single* xMov, proshade_single* yMov, proshade_single* zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed* xFrom, proshade_signed* xTo, proshade_signed* yFrom, proshade_signed* yTo, proshade_signed* zFrom, proshade_signed* zTo, proshade_signed* xOrigin, proshade_signed* yOrigin, proshade_signed* zOrigin )
778 {
779  //================================================ Compute movement in indices
780  proshade_single xIndMove = std::floor ( -(*xMov) / ( xAngs / ( static_cast< proshade_single > ( *xTo ) - static_cast< proshade_single > ( *xFrom ) + 1.0f ) ) );
781  proshade_single yIndMove = std::floor ( -(*yMov) / ( yAngs / ( static_cast< proshade_single > ( *yTo ) - static_cast< proshade_single > ( *yFrom ) + 1.0f ) ) );
782  proshade_single zIndMove = std::floor ( -(*zMov) / ( zAngs / ( static_cast< proshade_single > ( *zTo ) - static_cast< proshade_single > ( *zFrom ) + 1.0f ) ) );
783 
784  //================================================ Set the movs to the remainder
785  *xMov = -( *xMov ) - ( xIndMove * ( xAngs / ( static_cast< proshade_single > ( *xTo ) - static_cast< proshade_single > ( *xFrom ) + 1.0f ) ) );
786  *yMov = -( *yMov ) - ( yIndMove * ( yAngs / ( static_cast< proshade_single > ( *yTo ) - static_cast< proshade_single > ( *yFrom ) + 1.0f ) ) );
787  *zMov = -( *zMov ) - ( zIndMove * ( zAngs / ( static_cast< proshade_single > ( *zTo ) - static_cast< proshade_single > ( *zFrom ) + 1.0f ) ) );
788 
789  //================================================ Move indices by as much
790  *xFrom += static_cast< proshade_signed > ( xIndMove );
791  *xTo += static_cast< proshade_signed > ( xIndMove );
792  *yFrom += static_cast< proshade_signed > ( yIndMove );
793  *yTo += static_cast< proshade_signed > ( yIndMove );
794  *zFrom += static_cast< proshade_signed > ( zIndMove );
795  *zTo += static_cast< proshade_signed > ( zIndMove );
796 
797  //================================================ And set origin to reflect the changes
798  *xOrigin = *xFrom;
799  *yOrigin = *yFrom;
800  *zOrigin = *zFrom;
801 
802  //================================================ Done
803  return ;
804 
805 }
806 
825 void ProSHADE_internal_mapManip::moveMapByFourier ( proshade_double*& map, proshade_single xMov, proshade_single yMov, proshade_single zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim )
826 {
827  //================================================ Local variables initialisation
828  proshade_unsign arrayPos = 0;
829 
830  //================================================ Create Fourier map variable
831  fftw_complex *fCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * static_cast< proshade_unsign > ( xDim * yDim * zDim ) ) );
832  fftw_complex *translatedMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * static_cast< proshade_unsign > ( xDim * yDim * zDim ) ) );
833 
834  //================================================ Check memory allocation
835  ProSHADE_internal_misc::checkMemoryAllocation ( fCoeffs, __FILE__, __LINE__, __func__ );
836  ProSHADE_internal_misc::checkMemoryAllocation ( translatedMap, __FILE__, __LINE__, __func__ );
837 
838  //================================================ Create plans
839  fftw_plan planForwardFourier = fftw_plan_dft_3d ( static_cast< int > ( xDim ), static_cast< int > ( yDim ), static_cast< int > ( zDim ), translatedMap, fCoeffs, FFTW_FORWARD, FFTW_ESTIMATE );
840  fftw_plan planBackwardFourier = fftw_plan_dft_3d ( static_cast< int > ( xDim ), static_cast< int > ( yDim ), static_cast< int > ( zDim ), fCoeffs, translatedMap, FFTW_BACKWARD, FFTW_ESTIMATE );
841 
842  //================================================ Copy map to complex format
843  for ( proshade_unsign uIt = 0; uIt < static_cast< proshade_unsign > ( xDim ); uIt++ )
844  {
845  for ( proshade_unsign vIt = 0; vIt < static_cast< proshade_unsign > ( yDim ); vIt++ )
846  {
847  for ( proshade_unsign wIt = 0; wIt < static_cast< proshade_unsign > ( zDim ); wIt++ )
848  {
849  arrayPos = wIt + static_cast< proshade_unsign > ( zDim ) * ( vIt + static_cast< proshade_unsign > ( yDim ) * uIt );
850 
851  const FloatingPoint< proshade_double > lhs ( map[arrayPos] ), rhs ( map[arrayPos] );
852  if ( lhs.AlmostEquals ( rhs ) ) { translatedMap[arrayPos][0] = map[arrayPos]; }
853  else { translatedMap[arrayPos][0] = 0.0; }
854  translatedMap[arrayPos][1] = 0.0;
855  }
856  }
857  }
858 
859  //================================================ Compute Forward Fourier
860  fftw_execute ( planForwardFourier );
861 
862  //================================================ Shift the Fourier coefficients
863  proshade_double *weight = nullptr;
864  moveMapByFourierInReci ( fCoeffs, weight, xMov, yMov, zMov, xAngs, yAngs, zAngs, xDim, yDim, zDim );
865 
866  //================================================ Compute inverse Fourier
867  fftw_execute ( planBackwardFourier );
868 
869  //================================================ Copy back to map
870  for ( proshade_unsign uIt = 0; uIt < static_cast< proshade_unsign > ( xDim ); uIt++ )
871  {
872  for ( proshade_unsign vIt = 0; vIt < static_cast< proshade_unsign > ( yDim ); vIt++ )
873  {
874  for ( proshade_unsign wIt = 0; wIt < static_cast< proshade_unsign > ( zDim ); wIt++ )
875  {
876  arrayPos = wIt + static_cast< proshade_unsign > ( zDim ) * ( vIt + static_cast< proshade_unsign > ( yDim ) * uIt );
877  map[arrayPos] = translatedMap[arrayPos][0];
878  }
879  }
880  }
881 
882  //================================================ Release memory
883  fftw_destroy_plan ( planForwardFourier );
884  fftw_destroy_plan ( planBackwardFourier );
885  fftw_free ( fCoeffs );
886  fftw_free ( translatedMap );
887 
888  //================================================ Done
889  return ;
890 
891 }
892 
910 void ProSHADE_internal_mapManip::moveMapByFourierInReci ( proshade_complex*& coeffs, proshade_double*& weights, proshade_single xMov, proshade_single yMov, proshade_single zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim )
911 {
912  //================================================ Local variables initialisation
913  proshade_unsign arrayPos = 0;
914  proshade_signed h, k, l;
915  proshade_double real = 0.0;
916  proshade_double imag = 0.0;
917  proshade_double trCoeffReal, trCoeffImag;
918  proshade_double normFactor = static_cast< proshade_double > ( xDim * yDim * zDim );
919  proshade_double exponent = 0.0;
920  proshade_double hlpArrReal;
921  proshade_double hlpArrImag;
922 
923  //================================================ If no weights are given, set them to 1.0, otherwise use supplied
924  proshade_double* wght = new proshade_double[xDim * yDim * zDim];
925  ProSHADE_internal_misc::checkMemoryAllocation ( wght, __FILE__, __LINE__, __func__ );
926  if ( weights == nullptr ) { for ( size_t iter = 0; iter < static_cast< size_t > ( xDim * yDim * zDim ); iter++ ) { wght[iter] = 1.0; } }
927  else { for ( size_t iter = 0; iter < static_cast< size_t > ( xDim * yDim * zDim ); iter++ ) { wght[iter] = weights[iter]; } }
928 
929  //================================================ Add the required shift
930  for ( proshade_unsign uIt = 0; uIt < static_cast<proshade_unsign> ( xDim ); uIt++ )
931  {
932  for ( proshade_unsign vIt = 0; vIt < static_cast<proshade_unsign> ( yDim ); vIt++ )
933  {
934  for ( proshade_unsign wIt = 0; wIt < static_cast<proshade_unsign> ( zDim ); wIt++ )
935  {
936  //==================================== Var init
937  arrayPos = wIt + static_cast< proshade_unsign > ( zDim ) * ( vIt + static_cast< proshade_unsign > ( yDim ) * uIt );
938  real = coeffs[arrayPos][0];
939  imag = coeffs[arrayPos][1];
940 
941  //==================================== Convert 0-max indices to HKL
942  if ( uIt > static_cast< proshade_unsign > ( (xDim+1) / 2) ) { h = static_cast < proshade_signed > ( uIt ) - xDim; } else { h = static_cast < proshade_signed > ( uIt ); }
943  if ( vIt > static_cast< proshade_unsign > ( (yDim+1) / 2) ) { k = static_cast < proshade_signed > ( vIt ) - yDim; } else { k = static_cast < proshade_signed > ( vIt ); }
944  if ( wIt > static_cast< proshade_unsign > ( (zDim+1) / 2) ) { l = static_cast < proshade_signed > ( wIt ) - zDim; } else { l = static_cast < proshade_signed > ( wIt ); }
945 
946  //==================================== Get translation coefficient change
947  exponent = ( ( ( static_cast <proshade_double> ( h ) / static_cast <proshade_double> ( xAngs ) ) * static_cast< proshade_double > ( -xMov ) ) +
948  ( ( static_cast <proshade_double> ( k ) / static_cast <proshade_double> ( yAngs ) ) * static_cast< proshade_double > ( -yMov ) ) +
949  ( ( static_cast <proshade_double> ( l ) / static_cast <proshade_double> ( zAngs ) ) * static_cast< proshade_double > ( -zMov ) ) ) * 2.0 * M_PI;
950 
951  trCoeffReal = cos ( exponent );
952  trCoeffImag = sin ( exponent );
953  ProSHADE_internal_maths::complexMultiplication ( &real, &imag, &trCoeffReal, &trCoeffImag, &hlpArrReal, &hlpArrImag );
954 
955  //==================================== Save the translated coefficient value and apply weights
956  coeffs[arrayPos][0] = ( hlpArrReal / normFactor ) * wght[arrayPos];
957  coeffs[arrayPos][1] = ( hlpArrImag / normFactor ) * wght[arrayPos];
958  }
959  }
960  }
961 
962  //================================================ Release weights
963  delete[] wght;
964 
965  //================================================ Done
966  return ;
967 
968 }
969 
986 void ProSHADE_internal_mapManip::blurSharpenMap ( proshade_double*& map, proshade_double*& blurredMap, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single blurringFactor )
987 {
988  //================================================ Set local variables
989  proshade_signed xDim = static_cast< proshade_signed > ( xDimS );
990  proshade_signed yDim = static_cast< proshade_signed > ( yDimS );
991  proshade_signed zDim = static_cast< proshade_signed > ( zDimS );
992  proshade_double real, imag, S, mag, phase;
993  proshade_signed h, k, l;
994  proshade_unsign arrayPos = 0;
995  proshade_double normFactor = static_cast<proshade_double> ( xDim * yDim * zDim );
996 
997  //================================================ Copy map for processing
998  fftw_complex* mapCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * static_cast< proshade_unsign > ( xDim * yDim * zDim ) ) );
999  fftw_complex* mapMask = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * static_cast< proshade_unsign > ( xDim * yDim * zDim ) ) );
1000 
1001  //================================================ Check memory allocation
1002  ProSHADE_internal_misc::checkMemoryAllocation ( mapCoeffs, __FILE__, __LINE__, __func__ );
1003  ProSHADE_internal_misc::checkMemoryAllocation ( mapMask, __FILE__, __LINE__, __func__ );
1004 
1005  //================================================ Copy data to mask
1006  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> (xDim * yDim * zDim); iter++ )
1007  {
1008  mapMask[iter][0] = map[iter];
1009  mapMask[iter][1] = 0.0;
1010  }
1011 
1012  //================================================ Prepare FFTW plans
1013  fftw_plan forward = fftw_plan_dft_3d ( static_cast< int > ( xDim ), static_cast< int > ( yDim ), static_cast< int > ( zDim ), mapMask, mapCoeffs, FFTW_FORWARD, FFTW_ESTIMATE );
1014  fftw_plan inverse = fftw_plan_dft_3d ( static_cast< int > ( xDim ), static_cast< int > ( yDim ), static_cast< int > ( zDim ), mapCoeffs, mapMask, FFTW_BACKWARD, FFTW_ESTIMATE );
1015 
1016  //================================================ Run forward Fourier
1017  fftw_execute ( forward );
1018 
1019  //================================================ Blur the coeffs
1020  for ( proshade_unsign uIt = 0; uIt < static_cast<proshade_unsign> ( xDim ); uIt++ )
1021  {
1022  for ( proshade_unsign vIt = 0; vIt < static_cast<proshade_unsign> ( yDim ); vIt++ )
1023  {
1024  for ( proshade_unsign wIt = 0; wIt < static_cast<proshade_unsign> ( zDim ); wIt++ )
1025  {
1026  //==================================== Var init
1027  arrayPos = wIt + static_cast< proshade_unsign > ( zDim ) * ( vIt + static_cast< proshade_unsign > ( yDim ) * uIt );
1028  real = mapCoeffs[arrayPos][0];
1029  imag = mapCoeffs[arrayPos][1];
1030 
1031  //==================================== Convert to HKL
1032  if ( uIt > static_cast< proshade_unsign > ( (xDim+1) / 2) ) { h = static_cast < proshade_signed > ( uIt ) - xDim; } else { h = static_cast < proshade_signed > ( uIt ); }
1033  if ( vIt > static_cast< proshade_unsign > ( (yDim+1) / 2) ) { k = static_cast < proshade_signed > ( vIt ) - yDim; } else { k = static_cast < proshade_signed > ( vIt ); }
1034  if ( wIt > static_cast< proshade_unsign > ( (zDim+1) / 2) ) { l = static_cast < proshade_signed > ( wIt ) - zDim; } else { l = static_cast < proshade_signed > ( wIt ); }
1035 
1036  //====================================Get magnitude and phase with mask parameters
1037  S = ( pow( static_cast< proshade_double > ( h ) / static_cast< proshade_double > ( xAngs ), 2.0 ) +
1038  pow( static_cast< proshade_double > ( k ) / static_cast< proshade_double > ( yAngs ), 2.0 ) +
1039  pow( static_cast< proshade_double > ( l ) / static_cast< proshade_double > ( zAngs ), 2.0 ) );
1040  mag = std::sqrt ( (real*real) + (imag*imag) ) * std::exp ( - ( ( static_cast< proshade_double > ( blurringFactor ) * S ) / 4.0 ) );
1041  phase = std::atan2 ( imag, real );
1042 
1043  //==================================== Save the mask data
1044  mapCoeffs[arrayPos][0] = ( mag * cos(phase) ) / normFactor;
1045  mapCoeffs[arrayPos][1] = ( mag * sin(phase) ) / normFactor;
1046  }
1047  }
1048  }
1049 
1050  //================================================ Run inverse Fourier
1051  fftw_execute ( inverse );
1052 
1053  //================================================ Save the results
1054  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> (xDim * yDim * zDim); iter++ )
1055  {
1056  blurredMap[iter] = mapMask[iter][0];
1057  }
1058 
1059  //================================================ Release memory
1060  fftw_free ( mapMask );
1061  fftw_free ( mapCoeffs );
1062 
1063  //================================================ Delete FFTW plans
1064  fftw_destroy_plan ( forward );
1065  fftw_destroy_plan ( inverse );
1066 
1067  //================================================ Done
1068  return ;
1069 
1070 }
1071 
1086 void ProSHADE_internal_mapManip::getMaskFromBlurr ( proshade_double*& blurMap, proshade_double*& outMap, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_single noIQRs )
1087 {
1088  //================================================ Initialise vector for map data
1089  std::vector<proshade_double> mapVals ( xDim * yDim * zDim, 0.0 );
1090 
1091  //================================================ Save map values in vector
1092  for ( proshade_unsign iter = 0; iter < ( xDim * yDim * zDim ); iter++ )
1093  {
1094  mapVals.at(iter) = blurMap[iter];
1095  }
1096 
1097  //================================================ Find median and IQRs
1098  proshade_double* medAndIQR = new proshade_double[2];
1099  ProSHADE_internal_maths::vectorMedianAndIQR ( &mapVals, medAndIQR );
1100 
1101  //================================================ Find the threshold
1102  proshade_double maskThreshold = medAndIQR[0] + ( medAndIQR[1] * static_cast<proshade_double> ( noIQRs ) );
1103 
1104  //================================================ Apply threshold
1105  for ( proshade_unsign iter = 0; iter < ( xDim * yDim * zDim ); iter++ )
1106  {
1107  if ( blurMap[iter] < maskThreshold )
1108  {
1109  outMap[iter] = 0.0;
1110  blurMap[iter] = 0.0;
1111  }
1112  }
1113 
1114  //================================================ Release vector values
1115  mapVals.clear ( );
1116 
1117  //================================================ Release memory
1118  delete[] medAndIQR;
1119 
1120  //================================================ Done
1121  return ;
1122 
1123 }
1124 
1136 void ProSHADE_internal_mapManip::getNonZeroBounds ( proshade_double* map, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, proshade_signed*& ret )
1137 {
1138  //================================================ Initialise local variables
1139  proshade_signed arrayPos = 0;
1140 
1141  //================================================ Initialise result variable
1142  ret[0] = xDim;
1143  ret[1] = 0;
1144  ret[2] = yDim;
1145  ret[3] = 0;
1146  ret[4] = zDim;
1147  ret[5] = 0;
1148 
1149  //================================================ Iterate through map and check bounds
1150  for ( proshade_signed xIt = 0; xIt < xDim; xIt++ )
1151  {
1152  for ( proshade_signed yIt = 0; yIt < yDim; yIt++ )
1153  {
1154  for ( proshade_signed zIt = 0; zIt < zDim; zIt++ )
1155  {
1156  //==================================== Var init
1157  arrayPos = zIt + zDim * ( yIt + yDim * xIt );
1158 
1159  //==================================== Check bounds
1160  if ( map[arrayPos] > 0.001 )
1161  {
1162  if ( xIt < ret[0] ) { ret[0] = xIt; }
1163  if ( xIt > ret[1] ) { ret[1] = xIt; }
1164  if ( yIt < ret[2] ) { ret[2] = yIt; }
1165  if ( yIt > ret[3] ) { ret[3] = yIt; }
1166  if ( zIt < ret[4] ) { ret[4] = zIt; }
1167  if ( zIt > ret[5] ) { ret[5] = zIt; }
1168  }
1169  }
1170  }
1171  }
1172 
1173  //================================================ Done
1174  return ;
1175 
1176 }
1177 
1193 void ProSHADE_internal_mapManip::addExtraBoundSpace ( proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed*& bounds, proshade_single extraSpace )
1194 {
1195  //================================================ Convert angstrom distance to indices
1196  proshade_signed xExtraInds = ProSHADE_internal_mapManip::myRound ( extraSpace / ( xAngs / static_cast<proshade_single> ( xDim ) ) );
1197  proshade_signed yExtraInds = ProSHADE_internal_mapManip::myRound ( extraSpace / ( yAngs / static_cast<proshade_single> ( yDim ) ) );
1198  proshade_signed zExtraInds = ProSHADE_internal_mapManip::myRound ( extraSpace / ( zAngs / static_cast<proshade_single> ( zDim ) ) );
1199 
1200  //=============================================== Changed bounds even if exceeding physical map - this will be deal with in the map creation part
1201  bounds[0] = bounds[0] - xExtraInds;
1202  bounds[1] = bounds[1] + xExtraInds;
1203  bounds[2] = bounds[2] - yExtraInds;
1204  bounds[3] = bounds[3] + yExtraInds;
1205  bounds[4] = bounds[4] - zExtraInds;
1206  bounds[5] = bounds[5] + zExtraInds;
1207 
1208  //================================================ Done
1209  return ;
1210 
1211 }
1212 
1229 void ProSHADE_internal_mapManip::reSampleMapToResolutionTrilinear ( proshade_double*& map, proshade_single resolution, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single*& corrs )
1230 {
1231  //================================================ Sanity check - the resolution needs to be set
1232  if ( resolution <= 0.0f )
1233  {
1234  throw ProSHADE_exception ( "Requested resolution not set for map re-sampling.", "EM00015", __FILE__, __LINE__, __func__, "There is no resolution value set, but map re-sampling to\n : this unset resolution value is required. This error\n : occurs when a task with no resolution requirement is\n : requested on a map data and the map resolution change is\n : set to \'on\'. Either supply a resolution value, or do not\n : re-sample the map." );
1235  }
1236 
1237  //================================================ Initialise local variables
1238  proshade_signed xDim = static_cast<proshade_signed> ( xDimS );
1239  proshade_signed yDim = static_cast<proshade_signed> ( yDimS );
1240  proshade_signed zDim = static_cast<proshade_signed> ( zDimS );
1241  proshade_single oldXSample = ( xAngs / static_cast<proshade_single> ( xDim ) );
1242  proshade_single oldYSample = ( yAngs / static_cast<proshade_single> ( yDim ) );
1243  proshade_single oldZSample = ( zAngs / static_cast<proshade_single> ( zDim ) );
1244  proshade_single newXSample = static_cast< proshade_single > ( resolution / 2.0f );
1245  proshade_single newYSample = static_cast< proshade_single > ( resolution / 2.0f );
1246  proshade_single newZSample = static_cast< proshade_single > ( resolution / 2.0f );
1247 
1248  //================================================ Compute required grid size
1249  proshade_signed newXDim = static_cast<proshade_signed> ( std::ceil ( xAngs / newXSample ) );
1250  proshade_signed newYDim = static_cast<proshade_signed> ( std::ceil ( yAngs / newYSample ) );
1251  proshade_signed newZDim = static_cast<proshade_signed> ( std::ceil ( zAngs / newZSample ) );
1252 
1253  //================================================ Create a new map variable
1254  proshade_double* newMap = new proshade_double [newXDim * newYDim * newZDim];
1255 
1256  //================================================ For each new map point
1257  proshade_signed xBottom = 0, xTop, yBottom = 0, yTop, zBottom = 0, zTop, oldMapIndex, newMapIndex;
1258  std::vector<proshade_double> c000 = std::vector<proshade_double> ( 4, 0.0 );
1259  std::vector<proshade_double> c001 = std::vector<proshade_double> ( 4, 0.0 );
1260  std::vector<proshade_double> c010 = std::vector<proshade_double> ( 4, 0.0 );
1261  std::vector<proshade_double> c011 = std::vector<proshade_double> ( 4, 0.0 );
1262  std::vector<proshade_double> c100 = std::vector<proshade_double> ( 4, 0.0 );
1263  std::vector<proshade_double> c101 = std::vector<proshade_double> ( 4, 0.0 );
1264  std::vector<proshade_double> c110 = std::vector<proshade_double> ( 4, 0.0 );
1265  std::vector<proshade_double> c111 = std::vector<proshade_double> ( 4, 0.0 );
1266  std::vector<proshade_double> c00 = std::vector<proshade_double> ( 4, 0.0 );
1267  std::vector<proshade_double> c01 = std::vector<proshade_double> ( 4, 0.0 );
1268  std::vector<proshade_double> c10 = std::vector<proshade_double> ( 4, 0.0 );
1269  std::vector<proshade_double> c11 = std::vector<proshade_double> ( 4, 0.0 );
1270  std::vector<proshade_double> c0 = std::vector<proshade_double> ( 4, 0.0 );
1271  std::vector<proshade_double> c1 = std::vector<proshade_double> ( 4, 0.0 );
1272  proshade_double xRelative, yRelative, zRelative;
1273 
1274  for ( proshade_signed xIt = 0; xIt < newXDim; xIt++ )
1275  {
1276  for ( proshade_signed yIt = 0; yIt < newYDim; yIt++ )
1277  {
1278  for ( proshade_signed zIt = 0; zIt < newZDim; zIt++ )
1279  {
1280  //==================================== Get this point's index
1281  newMapIndex = zIt + newZDim * ( yIt + newYDim * xIt );
1282 
1283  //==================================== Find this points bottom and top positions in the old map (including periodicity)
1284  for ( proshade_signed ox = 0; ox < ( static_cast< proshade_signed > ( xDimS ) - 1 ); ox++ ) { if ( ( ( static_cast< proshade_single > ( xIt ) * newXSample ) >= ( static_cast< proshade_single > ( ox ) * oldXSample ) ) && ( ( static_cast< proshade_single > ( xIt ) * newXSample ) <= ( ( static_cast< proshade_single > ( ox ) + 1 ) * oldXSample ) ) ) { xBottom = ox; break; } }
1285  for ( proshade_signed oy = 0; oy < ( static_cast< proshade_signed > ( yDimS ) - 1 ); oy++ ) { if ( ( ( static_cast< proshade_single > ( yIt ) * newYSample ) >= ( static_cast< proshade_single > ( oy ) * oldYSample ) ) && ( ( static_cast< proshade_single > ( yIt ) * newYSample ) <= ( ( static_cast< proshade_single > ( oy ) + 1 ) * oldYSample ) ) ) { yBottom = oy; break; } }
1286  for ( proshade_signed oz = 0; oz < ( static_cast< proshade_signed > ( zDimS ) - 1 ); oz++ ) { if ( ( ( static_cast< proshade_single > ( zIt ) * newZSample ) >= ( static_cast< proshade_single > ( oz ) * oldZSample ) ) && ( ( static_cast< proshade_single > ( zIt ) * newZSample ) <= ( ( static_cast< proshade_single > ( oz ) + 1 ) * oldZSample ) ) ) { zBottom = oz; break; } }
1287  xTop = xBottom + 1;
1288  yTop = yBottom + 1;
1289  zTop = zBottom + 1;
1290 
1291  //==================================== Find the surrounding point's values from the original map
1292  oldMapIndex = zBottom + static_cast< proshade_signed > ( zDimS ) * ( yBottom + static_cast< proshade_signed > ( yDimS ) * xBottom );
1293  c000.at(0) = static_cast<proshade_double> ( xBottom ) * static_cast<proshade_double> ( oldXSample );
1294  c000.at(1) = static_cast<proshade_double> ( yBottom ) * static_cast<proshade_double> ( oldYSample );
1295  c000.at(2) = static_cast<proshade_double> ( zBottom ) * static_cast<proshade_double> ( oldZSample );
1296  c000.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1297 
1298  oldMapIndex = zTop + static_cast< proshade_signed > ( zDimS ) * ( yBottom + static_cast< proshade_signed > ( yDimS ) * xBottom );
1299  c001.at(0) = static_cast<proshade_double> ( xBottom ) * static_cast<proshade_double> ( oldXSample );
1300  c001.at(1) = static_cast<proshade_double> ( yBottom ) * static_cast<proshade_double> ( oldYSample );
1301  c001.at(2) = static_cast<proshade_double> ( zTop ) * static_cast<proshade_double> ( oldZSample );
1302  c001.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1303 
1304  oldMapIndex = zBottom + static_cast< proshade_signed > ( zDimS ) * ( yTop + static_cast< proshade_signed > ( yDimS ) * xBottom );
1305  c010.at(0) = static_cast<proshade_double> ( xBottom ) * static_cast<proshade_double> ( oldXSample );
1306  c010.at(1) = static_cast<proshade_double> ( yTop ) * static_cast<proshade_double> ( oldYSample );
1307  c010.at(2) = static_cast<proshade_double> ( zBottom ) * static_cast<proshade_double> ( oldZSample );
1308  c010.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1309 
1310  oldMapIndex = zTop + static_cast< proshade_signed > ( zDimS ) * ( yTop + static_cast< proshade_signed > ( yDimS ) * xBottom );
1311  c011.at(0) = static_cast<proshade_double> ( xBottom ) * static_cast<proshade_double> ( oldXSample );
1312  c011.at(1) = static_cast<proshade_double> ( yTop ) * static_cast<proshade_double> ( oldYSample );
1313  c011.at(2) = static_cast<proshade_double> ( zTop ) * static_cast<proshade_double> ( oldZSample );
1314  c011.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1315 
1316  oldMapIndex = zBottom + static_cast< proshade_signed > ( zDimS ) * ( yBottom + static_cast< proshade_signed > ( yDimS ) * xTop );
1317  c100.at(0) = static_cast<proshade_double> ( xTop ) * static_cast<proshade_double> ( oldXSample );
1318  c100.at(1) = static_cast<proshade_double> ( yBottom ) * static_cast<proshade_double> ( oldYSample );
1319  c100.at(2) = static_cast<proshade_double> ( zBottom ) * static_cast<proshade_double> ( oldZSample );
1320  c100.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1321 
1322  oldMapIndex = zTop + static_cast< proshade_signed > ( zDimS ) * ( yBottom + static_cast< proshade_signed > ( yDimS ) * xTop );
1323  c101.at(0) = static_cast<proshade_double> ( xTop ) * static_cast<proshade_double> ( oldXSample );
1324  c101.at(1) = static_cast<proshade_double> ( yBottom ) * static_cast<proshade_double> ( oldYSample );
1325  c101.at(2) = static_cast<proshade_double> ( zTop ) * static_cast<proshade_double> ( oldZSample );
1326  c101.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1327 
1328  oldMapIndex = zBottom + static_cast< proshade_signed > ( zDimS ) * ( yTop + static_cast< proshade_signed > ( yDimS ) * xTop );
1329  c110.at(0) = static_cast<proshade_double> ( xTop ) * static_cast<proshade_double> ( oldXSample );
1330  c110.at(1) = static_cast<proshade_double> ( yTop ) * static_cast<proshade_double> ( oldYSample );
1331  c110.at(2) = static_cast<proshade_double> ( zBottom ) * static_cast<proshade_double> ( oldZSample );
1332  c110.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1333 
1334  oldMapIndex = zTop + static_cast< proshade_signed > ( zDimS ) * ( yTop + static_cast< proshade_signed > ( yDimS ) * xTop );
1335  c111.at(0) = static_cast<proshade_double> ( xTop ) * static_cast<proshade_double> ( oldXSample );
1336  c111.at(1) = static_cast<proshade_double> ( yTop ) * static_cast<proshade_double> ( oldYSample );
1337  c111.at(2) = static_cast<proshade_double> ( zTop ) * static_cast<proshade_double> ( oldZSample );
1338  c111.at(3) = static_cast<proshade_double> ( map[oldMapIndex] );
1339 
1340  //==================================== Interpolate to the new grid along X
1341  xRelative = ( ( static_cast<proshade_double> ( xIt ) * static_cast<proshade_double> ( newXSample ) ) - ( static_cast<proshade_double> ( xBottom ) * static_cast<proshade_double> ( oldXSample ) ) ) / ( ( static_cast<proshade_double> ( xTop ) * static_cast<proshade_double> ( oldXSample ) ) - ( static_cast<proshade_double> ( xBottom ) * static_cast<proshade_double> ( oldXSample ) ) );
1342 
1343  //==================================== Interpolate for the less less point
1344  c00.at(0) = ( static_cast< proshade_double > ( newXSample ) * xRelative ) + c000.at(0);
1345  c00.at(1) = c000.at(1);
1346  c00.at(2) = c000.at(2);
1347  c00.at(3) = ( c000.at(3) * ( 1.0 - xRelative ) ) + ( c100.at(3) * xRelative );
1348 
1349  //==================================== Interpolate for the less more point
1350  c01.at(0) = ( static_cast< proshade_double > ( newXSample ) * xRelative ) + c001.at(0);
1351  c01.at(1) = c001.at(1);
1352  c01.at(2) = c001.at(2);
1353  c01.at(3) = ( c001.at(3) * ( 1.0 - xRelative ) ) + ( c101.at(3) * xRelative );
1354 
1355  //==================================== Interpolate for the more less point
1356  c10.at(0) = ( static_cast< proshade_double > ( newXSample ) * xRelative ) + c010.at(0);
1357  c10.at(1) = c010.at(1);
1358  c10.at(2) = c010.at(2);
1359  c10.at(3) = ( c010.at(3) * ( 1.0 - xRelative ) ) + ( c110.at(3) * xRelative );
1360 
1361  //==================================== Interpolate for the more more point
1362  c11.at(0) = ( static_cast< proshade_double > ( newXSample ) * xRelative ) + c011.at(0);
1363  c11.at(1) = c011.at(1);
1364  c11.at(2) = c011.at(2);
1365  c11.at(3) = ( c011.at(3) * ( 1.0 - xRelative ) ) + ( c111.at(3) * xRelative );
1366 
1367  //==================================== Interpolate to the new grid along Y
1368  yRelative = ( ( static_cast<proshade_double> ( yIt ) * static_cast<proshade_double> ( newYSample ) ) - ( static_cast<proshade_double> ( yBottom ) * static_cast<proshade_double> ( oldYSample ) ) ) / ( ( static_cast<proshade_double> ( yTop ) * static_cast<proshade_double> ( oldYSample ) ) - ( static_cast<proshade_double> ( yBottom ) * static_cast<proshade_double> ( oldYSample ) ) );
1369 
1370  //==================================== Interpolate for the less point
1371  c0.at(0) = c00.at(0);
1372  c0.at(1) = ( static_cast< proshade_double > ( newYSample ) * yRelative ) + c00.at(1);
1373  c0.at(2) = c00.at(2);
1374  c0.at(3) = ( c00.at(3) * ( 1.0 - yRelative ) ) + ( c10.at(3) * yRelative );
1375 
1376  //==================================== Interpolate for the more point
1377  c1.at(0) = c01.at(0);
1378  c1.at(1) = ( static_cast< proshade_double > ( newYSample ) * yRelative ) + c01.at(1);
1379  c1.at(2) = c01.at(2);
1380  c1.at(3) = ( c01.at(3) * ( 1.0 - yRelative ) ) + ( c11.at(3) * yRelative );
1381 
1382  //==================================== Interpolate to the new grid along Z
1383  zRelative = ( ( static_cast<proshade_double> ( zIt ) * static_cast< proshade_double > ( newZSample ) ) - ( static_cast<proshade_double> ( zBottom ) * static_cast<proshade_double> ( oldZSample ) ) ) / static_cast< proshade_double > ( ( static_cast<proshade_double> ( zTop ) * static_cast<proshade_double> ( oldZSample ) ) - ( static_cast<proshade_double> ( zBottom ) * static_cast<proshade_double> ( oldZSample ) ) );
1384  newMap[newMapIndex] = ( c0.at(3) * ( 1.0 - zRelative ) ) + ( c1.at(3) * zRelative );
1385  }
1386  }
1387  }
1388 
1389  //================================================ Delete old map and allocate new memory
1390  delete[] map;
1391  map = new proshade_double [newXDim * newYDim * newZDim];
1392 
1393  //================================================ Copy map
1394  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( newXDim * newYDim * newZDim ); iter++ )
1395  {
1396  map[iter] = newMap[iter];
1397  }
1398 
1399  //================================================ Release memory
1400  delete[] newMap;
1401 
1402  //================================================ Define change in indices and return it
1403  corrs[0] = static_cast< proshade_single > ( newXDim - xDim );
1404  corrs[1] = static_cast< proshade_single > ( newYDim - yDim );
1405  corrs[2] = static_cast< proshade_single > ( newZDim - zDim );
1406  corrs[3] = static_cast< proshade_single > ( newXDim ) * static_cast< proshade_single > ( newXSample );
1407  corrs[4] = static_cast< proshade_single > ( newYDim ) * static_cast< proshade_single > ( newYSample );
1408  corrs[5] = static_cast< proshade_single > ( newZDim ) * static_cast< proshade_single > ( newZSample );
1409 
1410  //======================================== Done
1411  return ;
1412 
1413 }
1414 
1439 void ProSHADE_internal_mapManip::reSampleMapToResolutionFourier ( proshade_double*& map, proshade_single resolution, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single*& corrs )
1440 {
1441  //================================================ Sanity check - the resolution needs to be set
1442  if ( resolution <= 0.0f )
1443  {
1444  throw ProSHADE_exception ( "Requested resolution not set for map re-sampling.", "EM00015", __FILE__, __LINE__, __func__, "There is no resolution value set, but map re-sampling to\n : this unset resolution value is required. This error\n : occurs when a task with no resolution requirement is\n : requested on a map data and the map resolution change is\n : set to \'on\'. Either supply a resolution value, or do not\n : re-sample the map." );
1445  }
1446 
1447  //================================================ Initialise variables
1448  proshade_unsign newXDim = static_cast<proshade_unsign> ( ProSHADE_internal_mapManip::myRound ( xAngs / ( resolution / 2.0f ) ) );
1449  proshade_unsign newYDim = static_cast<proshade_unsign> ( ProSHADE_internal_mapManip::myRound ( yAngs / ( resolution / 2.0f ) ) );
1450  proshade_unsign newZDim = static_cast<proshade_unsign> ( ProSHADE_internal_mapManip::myRound ( zAngs / ( resolution / 2.0f ) ) );
1451 
1452  if ( newXDim % 2 != 0 ) { newXDim += 1; }
1453  if ( newYDim % 2 != 0 ) { newYDim += 1; }
1454  if ( newZDim % 2 != 0 ) { newZDim += 1; }
1455 
1456  proshade_signed preXChange, preYChange, preZChange;
1457  preXChange = static_cast< proshade_signed > ( std::floor ( ( static_cast<proshade_single> ( xDimS ) - static_cast<proshade_single> ( newXDim ) ) / 2.0f ) );
1458  preYChange = static_cast< proshade_signed > ( std::floor ( ( static_cast<proshade_single> ( yDimS ) - static_cast<proshade_single> ( newYDim ) ) / 2.0f ) );
1459  preZChange = static_cast< proshade_signed > ( std::floor ( ( static_cast<proshade_single> ( zDimS ) - static_cast<proshade_single> ( newZDim ) ) / 2.0f ) );
1460 
1461  proshade_signed postXChange = static_cast<proshade_signed> ( xDimS ) - ( preXChange + static_cast<proshade_signed> ( newXDim ) );
1462  proshade_signed postYChange = static_cast<proshade_signed> ( yDimS ) - ( preYChange + static_cast<proshade_signed> ( newYDim ) );
1463  proshade_signed postZChange = static_cast<proshade_signed> ( zDimS ) - ( preZChange + static_cast<proshade_signed> ( newZDim ) );
1464 
1465  proshade_unsign origSizeArr = 0, newSizeArr = 0;
1466  proshade_double normFactor = static_cast<proshade_double> ( xDimS * yDimS * zDimS );
1467 
1468  //================================================ Manage memory
1469  fftw_complex *origMap, *fCoeffs, *newFCoeffs, *newMap;
1470  fftw_plan planForwardFourier, planBackwardRescaledFourier;
1471  allocateResolutionFourierMemory ( origMap, fCoeffs, newFCoeffs, newMap, planForwardFourier, planBackwardRescaledFourier,
1472  xDimS, yDimS, zDimS, newXDim, newYDim, newZDim );
1473 
1474  //================================================ Fill maps with data and zeroes
1475  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( normFactor ); iter++ ) { origMap[iter][0] = map[iter]; origMap[iter][1] = 0.0; }
1476  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( newXDim * newYDim * newZDim ); iter++ ) { newFCoeffs[iter][0] = 0.0; newFCoeffs[iter][1] = 0.0; }
1477 
1478  //================================================ Get the Fourier coeffs
1479  fftw_execute ( planForwardFourier );
1480 
1481  //================================================ Change the order of Fourier coefficients
1482  changeFourierOrder ( fCoeffs, static_cast< proshade_signed > ( xDimS ), static_cast< proshade_signed > ( yDimS ), static_cast< proshade_signed > ( zDimS ), true );
1483 
1484  //================================================ Re-sample the coefficients by removing high frequencies or adding these with 0 values
1485  for ( proshade_unsign xIt = 0; xIt < newXDim; xIt++ )
1486  {
1487  for ( proshade_unsign yIt = 0; yIt < newYDim; yIt++ )
1488  {
1489  for ( proshade_unsign zIt = 0; zIt < newZDim; zIt++ )
1490  {
1491  //==================================== Find the array positions
1492  origSizeArr = ( ( zIt + static_cast< proshade_unsign > ( preZChange ) ) + zDimS *
1493  ( ( yIt + static_cast< proshade_unsign > ( preYChange ) ) + yDimS *
1494  ( xIt + static_cast< proshade_unsign > ( preXChange ) ) ) );
1495  newSizeArr = zIt + newZDim * ( yIt + newYDim * xIt );
1496 
1497  //==================================== If original coefficient for this new coefficient position exists, copy
1498  if ( ( ( -1 < static_cast< proshade_signed > ( xIt ) + preXChange ) && ( -1 < static_cast<proshade_signed> ( yIt ) + preYChange ) && ( -1 < static_cast<proshade_signed> ( zIt ) + preZChange ) ) &&
1499  ( ( xIt < newXDim + static_cast<proshade_unsign> ( postXChange ) ) && ( yIt < newYDim + static_cast<proshade_unsign> ( postYChange ) ) && ( zIt < newZDim + static_cast<proshade_unsign> ( postZChange ) ) ) )
1500  {
1501  //================================ Copy the Fourier coeff
1502  newFCoeffs[newSizeArr][0] = fCoeffs[origSizeArr][0] / normFactor;
1503  newFCoeffs[newSizeArr][1] = fCoeffs[origSizeArr][1] / normFactor;
1504  }
1505  }
1506  }
1507  }
1508 
1509  //================================================ Change the order of the re-sampled Fourier coefficients
1510  changeFourierOrder ( newFCoeffs, static_cast< proshade_signed > ( newXDim ), static_cast< proshade_signed > ( newYDim ), static_cast< proshade_signed > ( newZDim ), false );
1511 
1512  //================================================ Get the new map from the re-sized Fourier coefficients
1513  fftw_execute ( planBackwardRescaledFourier );
1514 
1515  //================================================ Delete the old map and create a new, re-sized one. Then copy the new map values into this new map memory.
1516  delete map;
1517  map = new proshade_double [newXDim * newYDim * newZDim];
1518  ProSHADE_internal_misc::checkMemoryAllocation ( map, __FILE__, __LINE__, __func__ );
1519  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( newXDim * newYDim * newZDim ); iter++ ) { map[iter] = newMap[iter][0]; }
1520 
1521  //================================================ Release memory
1522  releaseResolutionFourierMemory ( origMap, fCoeffs, newFCoeffs, newMap, planForwardFourier, planBackwardRescaledFourier );
1523 
1524  //================================================ Define change in indices and return it
1525  corrs[0] = static_cast< proshade_single > ( newXDim ) - static_cast< proshade_single > ( xDimS );
1526  corrs[1] = static_cast< proshade_single > ( newYDim ) - static_cast< proshade_single > ( yDimS );
1527  corrs[2] = static_cast< proshade_single > ( newZDim ) - static_cast< proshade_single > ( zDimS );
1528  corrs[3] = static_cast< proshade_single > ( newXDim ) * static_cast< proshade_single > ( resolution / 2.0f );
1529  corrs[4] = static_cast< proshade_single > ( newYDim ) * static_cast< proshade_single > ( resolution / 2.0f );
1530  corrs[5] = static_cast< proshade_single > ( newZDim ) * static_cast< proshade_single > ( resolution / 2.0f );
1531 
1532  //======================================== Done
1533  return ;
1534 
1535 }
1536 
1555 void ProSHADE_internal_mapManip::allocateResolutionFourierMemory ( fftw_complex*& origMap, fftw_complex*& fCoeffs, fftw_complex*& newFCoeffs, fftw_complex*& newMap, fftw_plan& planForwardFourier, fftw_plan& planBackwardRescaledFourier, proshade_unsign xDimOld, proshade_unsign yDimOld, proshade_unsign zDimOld, proshade_unsign xDimNew, proshade_unsign yDimNew, proshade_unsign zDimNew )
1556 {
1557  //================================================ Initialise memory
1558  origMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * xDimOld * yDimOld * zDimOld ) );
1559  fCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * xDimOld * yDimOld * zDimOld ) );
1560  newFCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * xDimNew * yDimNew * zDimNew ) );
1561  newMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * xDimNew * yDimNew * zDimNew ) );
1562 
1563  //================================================ Check memory allocation
1564  ProSHADE_internal_misc::checkMemoryAllocation ( origMap, __FILE__, __LINE__, __func__ );
1565  ProSHADE_internal_misc::checkMemoryAllocation ( fCoeffs, __FILE__, __LINE__, __func__ );
1566  ProSHADE_internal_misc::checkMemoryAllocation ( newFCoeffs, __FILE__, __LINE__, __func__ );
1567  ProSHADE_internal_misc::checkMemoryAllocation ( newMap, __FILE__, __LINE__, __func__ );
1568 
1569  //================================================ Create plans
1570  planForwardFourier = fftw_plan_dft_3d ( static_cast< int > ( xDimOld ), static_cast< int > ( yDimOld ), static_cast< int > ( zDimOld ), origMap, fCoeffs, FFTW_FORWARD, FFTW_ESTIMATE );
1571  planBackwardRescaledFourier = fftw_plan_dft_3d ( static_cast< int > ( xDimNew ), static_cast< int > ( yDimNew ), static_cast< int > ( zDimNew ), newFCoeffs, newMap, FFTW_BACKWARD, FFTW_ESTIMATE );
1572 
1573  //================================================ Done
1574  return ;
1575 
1576 }
1577 
1589 void ProSHADE_internal_mapManip::releaseResolutionFourierMemory ( fftw_complex*& origMap, fftw_complex*& fCoeffs, fftw_complex*& newFCoeffs, fftw_complex*& newMap, fftw_plan& planForwardFourier, fftw_plan& planBackwardRescaledFourier )
1590 {
1591  //================================================ Delete the FFTW plans
1592  fftw_destroy_plan ( planForwardFourier );
1593  fftw_destroy_plan ( planBackwardRescaledFourier );
1594 
1595  //================================================ Delete the complex arrays
1596  fftw_free ( origMap );
1597  fftw_free ( fCoeffs );
1598  fftw_free ( newFCoeffs );
1599  fftw_free ( newMap );
1600 
1601  //================================================ Done
1602  return ;
1603 
1604 }
1605 
1618 void ProSHADE_internal_mapManip::changeFourierOrder ( fftw_complex*& fCoeffs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, bool negativeFirst )
1619 {
1620  //================================================ Initialise local variables
1621  proshade_signed h = 0, k = 0, l = 0, origSizeArr = 0, newSizeArr = 0;
1622  proshade_signed xSeq1FreqStart, ySeq1FreqStart, zSeq1FreqStart, xSeq2FreqStart, ySeq2FreqStart, zSeq2FreqStart;
1623 
1624  //================================================ Find the positive and negative indices cut-offs
1625  if ( negativeFirst )
1626  {
1627  if ( ( xDim % 2 ) == 0 ) { xSeq1FreqStart = xDim / 2; xSeq2FreqStart = xDim / 2; } else { xSeq1FreqStart = (xDim / 2) + 1; xSeq2FreqStart = xDim / 2; }
1628  if ( ( yDim % 2 ) == 0 ) { ySeq1FreqStart = yDim / 2; ySeq2FreqStart = yDim / 2; } else { ySeq1FreqStart = (yDim / 2) + 1; ySeq2FreqStart = yDim / 2; }
1629  if ( ( zDim % 2 ) == 0 ) { zSeq1FreqStart = zDim / 2; zSeq2FreqStart = zDim / 2; } else { zSeq1FreqStart = (zDim / 2) + 1; zSeq2FreqStart = zDim / 2; }
1630  }
1631  else
1632  {
1633  if ( ( xDim % 2 ) == 0 ) { xSeq1FreqStart = xDim / 2; xSeq2FreqStart = xDim / 2; } else { xSeq1FreqStart = (xDim / 2); xSeq2FreqStart = xDim / 2 + 1; }
1634  if ( ( yDim % 2 ) == 0 ) { ySeq1FreqStart = yDim / 2; ySeq2FreqStart = yDim / 2; } else { ySeq1FreqStart = (yDim / 2); ySeq2FreqStart = yDim / 2 + 1; }
1635  if ( ( zDim % 2 ) == 0 ) { zSeq1FreqStart = zDim / 2; zSeq2FreqStart = zDim / 2; } else { zSeq1FreqStart = (zDim / 2); zSeq2FreqStart = zDim / 2 + 1; }
1636  }
1637 
1638  //================================================ Allocate helper array memory
1639  fftw_complex *hlpFCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * static_cast< proshade_unsign > ( xDim * yDim * zDim ) ) );
1640  ProSHADE_internal_misc::checkMemoryAllocation ( hlpFCoeffs, __FILE__, __LINE__, __func__ );
1641 
1642  //================================================ Change the coefficients order
1643  for ( proshade_signed xIt = 0; xIt < xDim; xIt++ )
1644  {
1645  //============================================ Find x frequency
1646  if ( xIt < xSeq1FreqStart ) { h = xIt + xSeq2FreqStart; } else { h = xIt - xSeq1FreqStart; }
1647 
1648  for ( proshade_signed yIt = 0; yIt < yDim; yIt++ )
1649  {
1650  //======================================== Find y frequency
1651  if ( yIt < ySeq1FreqStart ) { k = yIt + ySeq2FreqStart; } else { k = yIt - ySeq1FreqStart; }
1652 
1653  for ( proshade_signed zIt = 0; zIt < zDim; zIt++ )
1654  {
1655  //==================================== Find z frequency
1656  if ( zIt < zSeq1FreqStart ) { l = zIt + zSeq2FreqStart; } else { l = zIt - zSeq1FreqStart; }
1657 
1658  //==================================== Find array positions
1659  newSizeArr = l + zDim * ( k + yDim * h );
1660  origSizeArr = zIt + zDim * ( yIt + yDim * xIt );
1661 
1662  //==================================== Copy vals
1663  hlpFCoeffs[newSizeArr][0] = fCoeffs[origSizeArr][0];
1664  hlpFCoeffs[newSizeArr][1] = fCoeffs[origSizeArr][1];
1665  }
1666  }
1667  }
1668 
1669  //================================================ Copy the helper array to the input Fourier coefficients array
1670  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( xDim * yDim * zDim ); iter++ ) { fCoeffs[iter][0] = hlpFCoeffs[iter][0]; fCoeffs[iter][1] = hlpFCoeffs[iter][1]; }
1671 
1672  //================================================ Release helper array memory
1673  fftw_free ( hlpFCoeffs );
1674 
1675  //================================================ Done
1676  return ;
1677 
1678 }
1679 
1691 void ProSHADE_internal_mapManip::removeMapPhase ( fftw_complex*& mapCoeffs, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim )
1692 {
1693  //================================================ Set local variables
1694  proshade_double real, imag, mag, phase;
1695  proshade_unsign arrayPos = 0;
1696  proshade_double normFactor = static_cast<proshade_double> ( xDim * yDim * zDim );
1697 
1698  //================================================ Iterate through the map
1699  for ( proshade_unsign uIt = 0; uIt < xDim; uIt++ )
1700  {
1701  for ( proshade_unsign vIt = 0; vIt < yDim; vIt++ )
1702  {
1703  for ( proshade_unsign wIt = 0; wIt < zDim; wIt++ )
1704  {
1705  //==================================== Var init
1706  arrayPos = wIt + zDim * ( vIt + yDim * uIt );
1707  real = mapCoeffs[arrayPos][0];
1708  imag = mapCoeffs[arrayPos][1];
1709 
1710  //==================================== Get magnitude and phase with mask parameters
1711  mag = std::sqrt ( (real*real) + (imag*imag) );;
1712  phase = 0.0; // This would be std::atan2 ( imag, real ); - but here we remove the phase.
1713 
1714  //==================================== Save the phaseless data
1715  mapCoeffs[arrayPos][0] = ( mag * cos(phase) ) / normFactor;
1716  mapCoeffs[arrayPos][1] = ( mag * sin(phase) ) / normFactor;
1717  }
1718  }
1719  }
1720 
1721  //================================================ Done
1722  return ;
1723 
1724 }
1725 
1740 void ProSHADE_internal_mapManip::getFakeHalfMap ( proshade_double*& map, proshade_double*& fakeHalfMap, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_signed fakeMapKernel )
1741 {
1742  //================================================ Set local variables
1743  proshade_signed xDim = static_cast< proshade_signed > ( xDimS );
1744  proshade_signed yDim = static_cast< proshade_signed > ( yDimS );
1745  proshade_signed zDim = static_cast< proshade_signed > ( zDimS );
1746  proshade_signed currentPos, neighArrPos, neighXPos, neighYPos, neighZPos;
1747  proshade_double neighSum;
1748  proshade_double neighCount = pow ( ( ( fakeMapKernel * 2 ) + 1 ), 3.0 ) - 1.0;
1749 
1750  //================================================ Blur the coeffs
1751  for ( proshade_signed uIt = 0; uIt < xDim; uIt++ )
1752  {
1753  for ( proshade_signed vIt = 0; vIt < yDim; vIt++ )
1754  {
1755  for ( proshade_signed wIt = 0; wIt < zDim; wIt++ )
1756  {
1757  //==================================== Var init
1758  currentPos = wIt + zDim * ( vIt + yDim * uIt );
1759  neighSum = 0.0;
1760 
1761  //==================================== Average neighbours
1762  for ( proshade_signed xCh = -fakeMapKernel; xCh <= +fakeMapKernel; xCh++ )
1763  {
1764  for ( proshade_signed yCh = -fakeMapKernel; yCh <= +fakeMapKernel; yCh++ )
1765  {
1766  for ( proshade_signed zCh = -fakeMapKernel; zCh <= +fakeMapKernel; zCh++ )
1767  {
1768  if ( ( xCh == 0 ) && ( yCh == 0 ) && ( zCh == 0 ) ) { continue; }
1769 
1770  //======================== Find the nieghbour peak indices (with periodicity)
1771  neighXPos = uIt + xCh; if ( neighXPos >= xDim ) { neighXPos -= xDim; }; if ( neighXPos < 0 ) { neighXPos += xDim; }
1772  neighYPos = vIt + yCh; if ( neighYPos >= yDim ) { neighYPos -= yDim; }; if ( neighYPos < 0 ) { neighYPos += yDim; }
1773  neighZPos = wIt + zCh; if ( neighZPos >= zDim ) { neighZPos -= zDim; }; if ( neighZPos < 0 ) { neighZPos += zDim; }
1774  neighArrPos = neighZPos + zDim * ( neighYPos + yDim * neighXPos );
1775 
1776  //======================== Add to average
1777  neighSum += map[neighArrPos];
1778  }
1779  }
1780  }
1781 
1782  //==================================== Save the average to "fake" half-map
1783  fakeHalfMap[currentPos] = neighSum / neighCount;
1784  }
1785  }
1786  }
1787 
1788  //================================================ Done
1789  return ;
1790 
1791 }
1792 
1805 void ProSHADE_internal_mapManip::getCorrelationMapMask ( proshade_double*& map, proshade_double*& fakeHalfMap, proshade_double*& correlationMask, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_signed corrMaskKernel )
1806 {
1807  //================================================ Set local variables
1808  proshade_signed xDim = static_cast< proshade_signed > ( xDimS ), yDim = static_cast< proshade_signed > ( yDimS ), zDim = static_cast< proshade_signed > ( zDimS ), currentPos, neighArrPos, neighXPos, neighYPos, neighZPos, corrIter;
1809  proshade_unsign noCorrVals = static_cast<proshade_unsign> ( pow ( ( ( corrMaskKernel * 2 ) + 1 ), 3 ) );
1810 
1811  //================================================ Alocate memory
1812  proshade_double *origMap = new proshade_double [noCorrVals];
1813  proshade_double *fakeHM = new proshade_double [noCorrVals];
1814 
1815  //================================================ Check memory allocation
1816  ProSHADE_internal_misc::checkMemoryAllocation ( origMap, __FILE__, __LINE__, __func__ );
1817  ProSHADE_internal_misc::checkMemoryAllocation ( fakeHM, __FILE__, __LINE__, __func__ );
1818 
1819  //================================================ Blur the coeffs
1820  for ( proshade_signed uIt = 0; uIt < xDim; uIt++ )
1821  {
1822  for ( proshade_signed vIt = 0; vIt < yDim; vIt++ )
1823  {
1824  for ( proshade_signed wIt = 0; wIt < zDim; wIt++ )
1825  {
1826  //==================================== Var init
1827  currentPos = wIt + zDim * ( vIt + yDim * uIt );
1828  corrIter = 0;
1829 
1830  //==================================== Average neighbours
1831  for ( proshade_signed xCh = -corrMaskKernel; xCh <= +corrMaskKernel; xCh++ )
1832  {
1833  for ( proshade_signed yCh = -corrMaskKernel; yCh <= +corrMaskKernel; yCh++ )
1834  {
1835  for ( proshade_signed zCh = -corrMaskKernel; zCh <= +corrMaskKernel; zCh++ )
1836  {
1837  //======================== Find the nieghbour peak indices (with periodicity)
1838  neighXPos = uIt + xCh; if ( neighXPos >= xDim ) { neighXPos -= xDim; }; if ( neighXPos < 0 ) { neighXPos += xDim; }
1839  neighYPos = vIt + yCh; if ( neighYPos >= yDim ) { neighYPos -= yDim; }; if ( neighYPos < 0 ) { neighYPos += yDim; }
1840  neighZPos = wIt + zCh; if ( neighZPos >= zDim ) { neighZPos -= zDim; }; if ( neighZPos < 0 ) { neighZPos += zDim; }
1841  neighArrPos = neighZPos + zDim * ( neighYPos + yDim * neighXPos );
1842 
1843  //======================== Add to correct arrays
1844  origMap[corrIter] = map[neighArrPos];
1845  fakeHM[corrIter] = fakeHalfMap[neighArrPos];
1846  corrIter += 1;
1847  }
1848  }
1849  }
1850 
1851  //==================================== Save the correlation comparison result
1852  correlationMask[currentPos] = ProSHADE_internal_maths::pearsonCorrCoeff ( origMap, fakeHM, noCorrVals );
1853  }
1854  }
1855  }
1856 
1857  //================================================ Done
1858  return ;
1859 
1860 }
1861 
1876 proshade_single ProSHADE_internal_mapManip::getIndicesFromAngstroms ( proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single dist )
1877 {
1878  //================================================ Compute
1879  proshade_single ret = static_cast< proshade_single > ( ProSHADE_internal_mapManip::myRound ( std::max ( dist / ( xAngs / static_cast<proshade_single> ( xDim ) ),
1880  std::max ( dist / ( yAngs / static_cast<proshade_single> ( yDim ) ),
1881  dist / ( zAngs / static_cast<proshade_single> ( zDim ) ) ) ) ) );
1882 
1883  //================================================ Done
1884  return ( ret );
1885 
1886 }
1887 
1899 void ProSHADE_internal_mapManip::connectMaskBlobs ( proshade_double*& mask, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single maskThres )
1900 {
1901  //================================================ Initialise variables
1902  proshade_double* hlpMap = new proshade_double[xDim * yDim * zDim];
1903  proshade_signed addSurroundingPoints = static_cast< proshade_signed > ( std::max ( 3L, static_cast<proshade_signed> ( std::ceil ( getIndicesFromAngstroms( static_cast< proshade_unsign > ( xDim ), static_cast< proshade_unsign > ( yDim ), static_cast< proshade_unsign > ( zDim ), xAngs, yAngs, zAngs, static_cast< proshade_single > ( std::max( xAngs, std::max( yAngs, zAngs ) ) * 0.1f ) ) ) ) ) );
1904  proshade_signed currPos, neighXPos, neighYPos, neighZPos, neighArrPos;
1905 
1906  //================================================ Check memory allocation
1907  ProSHADE_internal_misc::checkMemoryAllocation ( hlpMap, __FILE__, __LINE__, __func__ );
1908 
1909  //================================================ Copy the mask map
1910  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( xDim * yDim * zDim ); iter++ ) { hlpMap[iter] = mask[iter]; }
1911 
1912  //================================================ Repeat as many times as needed
1913  for ( proshade_signed it = 0; it < addSurroundingPoints; it++ )
1914  {
1915  //============================================ For each x, y and z
1916  for ( proshade_signed xIt = 0; xIt < xDim; xIt++ )
1917  {
1918  for ( proshade_signed yIt = 0; yIt < yDim; yIt++ )
1919  {
1920  for ( proshade_signed zIt = 0; zIt < zDim; zIt++ )
1921  {
1922  //================================ Current position
1923  currPos = zIt + zDim * ( yIt + yDim * xIt );
1924 
1925  //================================ If zero, ignore
1926  if ( hlpMap[currPos] < static_cast< proshade_double > ( maskThres ) ) { continue; }
1927 
1928  //================================ Check neighbours
1929  for ( proshade_signed xCh = -1; xCh <= +1; xCh++ )
1930  {
1931  for ( proshade_signed yCh = -1; yCh <= +1; yCh++ )
1932  {
1933  for ( proshade_signed zCh = -1; zCh <= +1; zCh++ )
1934  {
1935  if ( ( xCh == 0 ) && ( yCh == 0 ) && ( zCh == 0 ) ) { continue; }
1936 
1937  //==================== Find the nieghbour indices (without periodicity)
1938  neighXPos = xIt + xCh; if ( neighXPos < 0 ) { continue; } if ( neighXPos >= xDim ) { continue; }
1939  neighYPos = yIt + yCh; if ( neighYPos < 0 ) { continue; } if ( neighYPos >= yDim ) { continue; }
1940  neighZPos = zIt + zCh; if ( neighZPos < 0 ) { continue; } if ( neighZPos >= zDim ) { continue; }
1941  neighArrPos = neighZPos + zDim * ( neighYPos + yDim * neighXPos );
1942 
1943  //==================== Add to mask if this point is below it (as it is a neighbour to a point which is part of the mask)
1944  if ( hlpMap[neighArrPos] < static_cast< proshade_double > ( maskThres ) ) { mask[neighArrPos] = static_cast< proshade_double > ( maskThres ); }
1945  }
1946  }
1947  }
1948  }
1949  }
1950  }
1951 
1952  //============================================ Now copy the updated mask to the temporary helper, unless last iteration
1953  for ( proshade_unsign iter = 0; iter < static_cast<proshade_unsign> ( xDim * yDim * zDim ); iter++ ) { hlpMap[iter] = mask[iter]; }
1954  }
1955 
1956  //================================================ Release memory
1957  delete[] hlpMap;
1958 
1959  //================================================ Done
1960  return ;
1961 
1962 }
1963 
1976 void ProSHADE_internal_mapManip::beautifyBoundaries ( proshade_signed*& bounds, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_signed boundsDiffThres )
1977 {
1978  //================================================ If extra bounds space pushed the bounds over the physical map, freely add up to 10 indices over the current value
1979  while ( bounds[1] >= static_cast<proshade_signed> ( xDim ) ) { xDim += 10; }
1980  while ( bounds[3] >= static_cast<proshade_signed> ( yDim ) ) { yDim += 10; }
1981  while ( bounds[5] >= static_cast<proshade_signed> ( zDim ) ) { zDim += 10; }
1982 
1983  //================================================ Find if better bouds exist in terms of prime numbers
1984  proshade_signed addToX = betterClosePrimeFactors ( bounds[1] - bounds[0] + 1, static_cast< proshade_signed > ( xDim ) );
1985  proshade_signed addToY = betterClosePrimeFactors ( bounds[3] - bounds[2] + 1, static_cast< proshade_signed > ( yDim ) );
1986  proshade_signed addToZ = betterClosePrimeFactors ( bounds[5] - bounds[4] + 1, static_cast< proshade_signed > ( zDim ) );
1987 
1988  //================================================ Figure if similar sizes should not be forced to be the same
1989  proshade_signed XtoY = std::abs ( addToX - addToY );
1990  proshade_signed XtoZ = std::abs ( addToX - addToZ );
1991  proshade_signed YtoZ = std::abs ( addToY - addToZ );
1992 
1993  if ( ( ( XtoY < boundsDiffThres ) && ( XtoZ < boundsDiffThres ) ) ||
1994  ( ( XtoY < boundsDiffThres ) && ( YtoZ < boundsDiffThres ) ) ||
1995  ( ( XtoZ < boundsDiffThres ) && ( YtoZ < boundsDiffThres ) ) )
1996  {
1997  //============================================ Dealing with chain ( a <= b <= c type ) where at least two dista are smaller than threshold
1998  proshade_signed maxSize = std::max ( addToX, std::max ( addToY, addToZ ) );
1999  addToX = maxSize;
2000  addToY = maxSize;
2001  addToZ = maxSize;
2002  }
2003  else
2004  {
2005  //============================================ Only two may be similar enough
2006  if ( XtoY <= boundsDiffThres )
2007  {
2008  proshade_signed maxSize = std::max ( addToX, addToY );
2009  addToX = maxSize;
2010  addToY = maxSize;
2011  }
2012  if ( XtoZ <= boundsDiffThres )
2013  {
2014  proshade_signed maxSize = std::max ( addToX, addToZ );
2015  addToX = maxSize;
2016  addToZ = maxSize;
2017  }
2018  if ( YtoZ <= boundsDiffThres )
2019  {
2020  proshade_signed maxSize = std::max ( addToY, addToZ );
2021  addToY = maxSize;
2022  addToZ = maxSize;
2023  }
2024  }
2025 
2026  //================================================ Add the extra space appropriately
2027  distributeSpaceToBoundaries ( bounds[0], bounds[1], ( bounds[1] - bounds[0] + 1 ), addToX );
2028  distributeSpaceToBoundaries ( bounds[2], bounds[3], ( bounds[3] - bounds[2] + 1 ), addToY );
2029  distributeSpaceToBoundaries ( bounds[4], bounds[5], ( bounds[5] - bounds[4] + 1 ), addToZ );
2030 
2031  //================================================ Done
2032  return ;
2033 
2034 }
2035 
2046 proshade_signed ProSHADE_internal_mapManip::betterClosePrimeFactors ( proshade_signed fromRange, proshade_signed toRange )
2047 {
2048  //================================================ Initialise variables
2049  proshade_signed ret = fromRange;
2050  std::vector < proshade_signed > posibles, hlp;
2051  proshade_signed sum;
2052 
2053  //================================================ Find the sums of prime factors for each number in the whole range
2054  for ( proshade_signed iter = fromRange; iter < toRange; iter++ )
2055  {
2056  sum = 0;
2058  for ( proshade_unsign i = 0; i < static_cast<proshade_unsign> ( hlp.size() ); i++ ) { sum += hlp.at(i); }
2059  hlp.clear ( );
2060  ProSHADE_internal_misc::addToSignedVector ( &posibles, sum );
2061  }
2062 
2063  //================================================ Find a better number
2064  for ( proshade_signed iter = fromRange; iter < toRange; iter++ )
2065  {
2066  //============================================ Ignore odd numbers
2067  if ( iter %2 != 0 ) { continue; }
2068 
2069  //============================================ Better number needs to have lower prime factor sum, but also needs not to be too far
2070  if ( posibles.at( static_cast< size_t > ( iter - fromRange ) ) < ( posibles.at( static_cast< size_t > ( ret - fromRange ) ) - ( iter - ret ) ) ) { ret = iter; }
2071  }
2072 
2073  //================================================ In the case of large prime number, just add one for even number
2074  if ( ( ret % 2 != 0 ) && ( ret < ( toRange - 1 ) ) ) { ret += 1; }
2075 
2076  //================================================ Done
2077  return ( ret );
2078 
2079 }
2080 
2092 void ProSHADE_internal_mapManip::distributeSpaceToBoundaries ( proshade_signed& minBound, proshade_signed& maxBound, proshade_signed oldBoundRange, proshade_signed newBoundRange )
2093 {
2094  //================================================ Only bother when sometings is to be added
2095  if ( newBoundRange > oldBoundRange )
2096  {
2097  //============================================ How much are we adding?
2098  proshade_signed distributeThis = newBoundRange - oldBoundRange;
2099 
2100  while ( distributeThis != 0 )
2101  {
2102  minBound -= 1;
2103  distributeThis -= 1;
2104 
2105  if ( distributeThis != 0 )
2106  {
2107  maxBound += 1;
2108  distributeThis -= 1;
2109  }
2110  }
2111  }
2112 
2113  //================================================ Done
2114  return ;
2115 
2116 }
2117 
2142 void ProSHADE_internal_mapManip::copyMapByBounds ( proshade_signed xFrom, proshade_signed xTo, proshade_signed yFrom, proshade_signed yTo, proshade_signed zFrom, proshade_signed zTo, proshade_signed origXFrom, proshade_signed origYFrom, proshade_signed origZFrom, proshade_unsign yDimIndices, proshade_unsign zDimIndices, proshade_unsign origXDimIndices, proshade_unsign origYDimIndices, proshade_unsign origZDimIndices, proshade_double*& newMap, proshade_double* origMap )
2143 {
2144  //================================================ Initialise variables
2145  proshade_signed newMapIndex, oldMapIndex, oldX, oldY, oldZ, newX, newY, newZ;
2146 
2147  //=============================================== For all values in the new map
2148  for ( proshade_signed xIt = xFrom; xIt <= xTo; xIt++ )
2149  {
2150  //============================================ Index position init
2151  newX = ( xIt - xFrom );
2152  oldX = ( newX + ( xFrom - origXFrom ) );
2153 
2154  for ( proshade_signed yIt = yFrom; yIt <= yTo; yIt++ )
2155  {
2156  //======================================== Index position init
2157  newY = ( yIt - yFrom );
2158  oldY = ( newY + ( yFrom - origYFrom ) );
2159 
2160  for ( proshade_signed zIt = zFrom; zIt <= zTo; zIt++ )
2161  {
2162  //==================================== Index position init
2163  newZ = ( zIt - zFrom );
2164  oldZ = ( newZ + ( zFrom - origZFrom ) );
2165 
2166  //==================================== Index arrays
2167  newMapIndex = newZ + static_cast< proshade_signed > ( zDimIndices ) * ( newY + static_cast< proshade_signed > ( yDimIndices ) * newX );
2168  oldMapIndex = oldZ + static_cast< proshade_signed > ( origZDimIndices ) * ( oldY + static_cast< proshade_signed > ( origYDimIndices ) * oldX );
2169 
2170  //============================ Check if we are adding new point
2171  if ( ( ( oldX < 0 ) || ( oldX >= static_cast< proshade_signed > ( origXDimIndices ) ) ) ||
2172  ( ( oldY < 0 ) || ( oldY >= static_cast< proshade_signed > ( origYDimIndices ) ) ) ||
2173  ( ( oldZ < 0 ) || ( oldZ >= static_cast< proshade_signed > ( origZDimIndices ) ) ) )
2174  {
2175  //================================ Yes, this is a new point, no known value for it
2176  newMap[newMapIndex] = 0.0;
2177  continue;
2178  }
2179 
2180  //==================================== If we got here, this should be within the old map, so just copy value
2181  newMap[newMapIndex] = origMap[oldMapIndex];
2182  }
2183  }
2184  }
2185 
2186  //================================================ Done
2187  return ;
2188 
2189 }
ProSHADE_internal_mapManip::addExtraBoundSpace
void addExtraBoundSpace(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed *&bounds, proshade_single extraSpace)
This function takes a set of bounds and adds a given number of Angstroms to them.
Definition: ProSHADE_mapManip.cpp:1193
ProSHADE_internal_mapManip::determinePDBRanges
void determinePDBRanges(gemmi::Structure pdbFile, proshade_single *xFrom, proshade_single *xTo, proshade_single *yFrom, proshade_single *yTo, proshade_single *zFrom, proshade_single *zTo, bool firstModel)
Function for finding the PDB file ranges.
Definition: ProSHADE_mapManip.cpp:72
ProSHADE_internal_mapManip::changePDBBFactors
void changePDBBFactors(gemmi::Structure *pdbFile, proshade_double newBFactorValue, bool firstModel)
Function for changing the PDB B-factor values to a specific single value.
Definition: ProSHADE_mapManip.cpp:458
ProSHADE_exception
This class is the representation of ProSHADE exception.
Definition: ProSHADE_exceptions.hpp:37
ProSHADE_internal_maths::getRotationMatrixFromEulerZYZAngles
void getRotationMatrixFromEulerZYZAngles(proshade_double eulerAlpha, proshade_double eulerBeta, proshade_double eulerGamma, proshade_double *matrix)
Function to find the rotation matrix from Euler angles (ZYZ convention).
Definition: ProSHADE_maths.cpp:1019
ProSHADE_internal_mapManip::distributeSpaceToBoundaries
void distributeSpaceToBoundaries(proshade_signed &minBound, proshade_signed &maxBound, proshade_signed oldBoundRange, proshade_signed newBoundRange)
Function for adding space to boundaries within the map confines.
Definition: ProSHADE_mapManip.cpp:2092
ProSHADE_internal_mapManip::getMaskFromBlurr
void getMaskFromBlurr(proshade_double *&blurMap, proshade_double *&outMap, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single noIQRs)
Function for computing mask from blurred map.
Definition: ProSHADE_mapManip.cpp:1086
ProSHADE_mapManip.hpp
This header file declares the ProSHADE_internal_mapManip namespace, which groups functions for intern...
ProSHADE_internal_mapManip::blurSharpenMap
void blurSharpenMap(proshade_double *&map, proshade_double *&maskedMap, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single blurringFactor)
Function for blurring/sharpening maps.
Definition: ProSHADE_mapManip.cpp:986
ProSHADE_internal_mapManip::moveMapByFourier
void moveMapByFourier(proshade_double *&map, proshade_single xMov, proshade_single yMov, proshade_single zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim)
Function for moving map back to original PDB location by using Fourier transformation.
Definition: ProSHADE_mapManip.cpp:825
ProSHADE_internal_mapManip::movePDBForMapCalc
void movePDBForMapCalc(gemmi::Structure *pdbFile, proshade_single xMov, proshade_single yMov, proshade_single zMov, bool firstModel)
Function for moving co-ordinate atoms to better suit theoretical map computation.
Definition: ProSHADE_mapManip.cpp:587
ProSHADE_internal_maths::complexMultiplication
void complexMultiplication(proshade_double *r1, proshade_double *i1, proshade_double *r2, proshade_double *i2, proshade_double *retReal, proshade_double *retImag)
Function to multiply two complex numbers.
Definition: ProSHADE_maths.cpp:38
ProSHADE_internal_messages::printWarningMessage
void printWarningMessage(proshade_signed verbose, std::string message, std::string warnCode)
General stderr message printing (used for warnings).
Definition: ProSHADE_messages.cpp:102
ProSHADE_internal_mapManip::reSampleMapToResolutionTrilinear
void reSampleMapToResolutionTrilinear(proshade_double *&map, proshade_single resolution, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single *&corrs)
This function re-samples a map to conform to given resolution using tri-linear interpolation.
Definition: ProSHADE_mapManip.cpp:1229
ProSHADE_internal_mapManip::getCorrelationMapMask
void getCorrelationMapMask(proshade_double *&map, proshade_double *&fakeHalfMap, proshade_double *&correlationMask, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_signed corrMaskKernel)
Function for creating the correlation mask.
Definition: ProSHADE_mapManip.cpp:1805
ProSHADE_internal_mapManip::changeFourierOrder
void changeFourierOrder(fftw_complex *&fCoeffs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, bool negativeFirst)
This function changes the order of Fourier coefficients in a 3D array between positive first (default...
Definition: ProSHADE_mapManip.cpp:1618
ProSHADE_internal_mapManip::releaseResolutionFourierMemory
void releaseResolutionFourierMemory(fftw_complex *&origMap, fftw_complex *&fCoeffs, fftw_complex *&newFCoeffs, fftw_complex *&newMap, fftw_plan &planForwardFourier, fftw_plan &planBackwardRescaledFourier)
This function releases the memory required by the Fourier resampling.
Definition: ProSHADE_mapManip.cpp:1589
ProSHADE_internal_mapManip::removeWaters
void removeWaters(gemmi::Structure *pdbFile, bool firstModel)
This function removed all waters from PDB input file.
Definition: ProSHADE_mapManip.cpp:518
ProSHADE_internal_maths::pearsonCorrCoeff
proshade_double pearsonCorrCoeff(proshade_double *valSet1, proshade_double *valSet2, proshade_unsign length)
Function for computing the Pearson's correlation coefficient.
Definition: ProSHADE_maths.cpp:246
ProSHADE_internal_misc::addToSignedVector
void addToSignedVector(std::vector< proshade_signed > *vecToAddTo, proshade_signed elementToAdd)
Adds the element to the vector.
Definition: ProSHADE_misc.cpp:121
ProSHADE_internal_mapManip::getIndicesFromAngstroms
proshade_single getIndicesFromAngstroms(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single dist)
This function converts distance in Angstroms to distance in map indices.
Definition: ProSHADE_mapManip.cpp:1876
ProSHADE_internal_mapManip::getNonZeroBounds
void getNonZeroBounds(proshade_double *map, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, proshade_signed *&ret)
Function for finding the map boundaries enclosing positive only values.
Definition: ProSHADE_mapManip.cpp:1136
ProSHADE_internal_mapManip::findMAPCOMValues
void findMAPCOMValues(proshade_double *map, proshade_double *xCom, proshade_double *yCom, proshade_double *zCom, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xFrom, proshade_signed xTo, proshade_signed yFrom, proshade_signed yTo, proshade_signed zFrom, proshade_signed zTo, bool removeNegDens)
This function finds the Centre of Mass for a map.
Definition: ProSHADE_mapManip.cpp:243
ProSHADE_internal_mapManip::removeMapPhase
void removeMapPhase(fftw_complex *&mapCoeffs, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim)
This function removes the phase from reciprocal (frequency) map.
Definition: ProSHADE_mapManip.cpp:1691
ProSHADE_internal_mapManip::allocateResolutionFourierMemory
void allocateResolutionFourierMemory(fftw_complex *&origMap, fftw_complex *&fCoeffs, fftw_complex *&newFCoeffs, fftw_complex *&newMap, fftw_plan &planForwardFourier, fftw_plan &planBackwardRescaledFourier, proshade_unsign xDimOld, proshade_unsign yDimOld, proshade_unsign zDimOld, proshade_unsign xDimNew, proshade_unsign yDimNew, proshade_unsign zDimNew)
This function allocates and checks the allocatio of the memory required by the Fourier resampling.
Definition: ProSHADE_mapManip.cpp:1555
ProSHADE_internal_maths::primeFactorsDecomp
std::vector< proshade_signed > primeFactorsDecomp(proshade_signed number)
Function to find prime factors of an integer.
Definition: ProSHADE_maths.cpp:1719
ProSHADE_internal_mapManip::betterClosePrimeFactors
proshade_signed betterClosePrimeFactors(proshade_signed fromRange, proshade_signed toRange)
Function for finding close numbers with better prime factors.
Definition: ProSHADE_mapManip.cpp:2046
ProSHADE_internal_mapManip::myRound
proshade_signed myRound(proshade_double x)
Calls the appropriate version of round function depending on compiler version.
Definition: ProSHADE_mapManip.cpp:31
ProSHADE_internal_mapManip::findPDBCOMValues
void findPDBCOMValues(gemmi::Structure *pdbFile, proshade_double *xCom, proshade_double *yCom, proshade_double *zCom, bool firstModel)
This function finds the Centre of Mass for the co-ordinate file.
Definition: ProSHADE_mapManip.cpp:157
ProSHADE_internal_mapManip::beautifyBoundaries
void beautifyBoundaries(proshade_signed *&bounds, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_signed boundsDiffThres)
Function for modifying boundaries to a mathematically more pleasant values.
Definition: ProSHADE_mapManip.cpp:1976
ProSHADE_internal_mapManip::rotatePDBCoordinates
void rotatePDBCoordinates(gemmi::Structure *pdbFile, proshade_double euA, proshade_double euB, proshade_double euG, proshade_double xCom, proshade_double yCom, proshade_double zCom, bool firstModel)
Function for rotating the PDB file co-ordinates by Euler angles.
Definition: ProSHADE_mapManip.cpp:310
ProSHADE_internal_mapManip::generateMapFromPDB
void generateMapFromPDB(gemmi::Structure pdbFile, proshade_double *&map, proshade_single requestedResolution, proshade_single xCell, proshade_single yCell, proshade_single zCell, proshade_signed *xTo, proshade_signed *yTo, proshade_signed *zTo, bool forceP1, bool firstModel)
This function generates a theoretical map from co-ordinate input files.
Definition: ProSHADE_mapManip.cpp:658
ProSHADE_internal_misc::checkMemoryAllocation
void checkMemoryAllocation(chVar checkVar, std::string fileP, unsigned int lineP, std::string funcP, std::string infoP="This error may occurs when ProSHADE requests memory to be\n : allocated to it and this operation fails. This could\n : happen when not enough memory is available, either due to\n : other processes using a lot of memory, or when the machine\n : does not have sufficient memory available. Re-run to see\n : if this problem persists.")
Checks if memory was allocated properly.
Definition: ProSHADE_misc.hpp:73
ProSHADE_internal_mapManip::translatePDBCoordinates
void translatePDBCoordinates(gemmi::Structure *pdbFile, proshade_double transX, proshade_double transY, proshade_double transZ, bool firstModel)
Function for translating the PDB file co-ordinates by given distances in Angstroms.
Definition: ProSHADE_mapManip.cpp:395
ProSHADE_internal_maths::vectorMedianAndIQR
void vectorMedianAndIQR(std::vector< proshade_double > *vec, proshade_double *&ret)
Function to get vector median and inter-quartile range.
Definition: ProSHADE_maths.cpp:149
ProSHADE_internal_mapManip::getFakeHalfMap
void getFakeHalfMap(proshade_double *&map, proshade_double *&fakeHalfMap, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_signed fakeMapKernel)
Function for creating "fake" half-maps.
Definition: ProSHADE_mapManip.cpp:1740
ProSHADE_internal_misc::addToUnsignVector
void addToUnsignVector(std::vector< proshade_unsign > *vecToAddTo, proshade_unsign elementToAdd)
Adds the element to the vector.
Definition: ProSHADE_misc.cpp:99
ProSHADE_internal_mapManip::reSampleMapToResolutionFourier
void reSampleMapToResolutionFourier(proshade_double *&map, proshade_single resolution, proshade_unsign xDimS, proshade_unsign yDimS, proshade_unsign zDimS, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single *&corrs)
This function re-samples a map to conform to given resolution using Fourier.
Definition: ProSHADE_mapManip.cpp:1439
ProSHADE_internal_mapManip::copyMapByBounds
void copyMapByBounds(proshade_signed xFrom, proshade_signed xTo, proshade_signed yFrom, proshade_signed yTo, proshade_signed zFrom, proshade_signed zTo, proshade_signed origXFrom, proshade_signed origYFrom, proshade_signed origZFrom, proshade_unsign yDimIndices, proshade_unsign zDimIndices, proshade_unsign origXDimIndices, proshade_unsign origYDimIndices, proshade_unsign origZDimIndices, proshade_double *&newMap, proshade_double *origMap)
This function copies an old map to a new map with different boundaries.
Definition: ProSHADE_mapManip.cpp:2142
ProSHADE_internal_mapManip::connectMaskBlobs
void connectMaskBlobs(proshade_double *&mask, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_single maskThres)
This function connects blobs in mask.
Definition: ProSHADE_mapManip.cpp:1899
ProSHADE_internal_mapManip::moveMapByFourierInReci
void moveMapByFourierInReci(proshade_complex *&coeffs, proshade_double *&weights, proshade_single xMov, proshade_single yMov, proshade_single zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim)
Function for moving map using Fourier coefficients in the reciprocal space only.
Definition: ProSHADE_mapManip.cpp:910
ProSHADE_internal_mapManip::moveMapByIndices
void moveMapByIndices(proshade_single *xMov, proshade_single *yMov, proshade_single *zMov, proshade_single xAngs, proshade_single yAngs, proshade_single zAngs, proshade_signed *xFrom, proshade_signed *xTo, proshade_signed *yFrom, proshade_signed *yTo, proshade_signed *zFrom, proshade_signed *zTo, proshade_signed *xOrigin, proshade_signed *yOrigin, proshade_signed *zOrigin)
Function for moving map back to original PDB location by changing the indices.
Definition: ProSHADE_mapManip.cpp:777