ProSHADE  0.7.6.6 (JUL 2022)
Protein Shape Detection
ProSHADE_io.cpp
Go to the documentation of this file.
1 
22 //==================================================== ProSHADE
23 #include "ProSHADE_io.hpp"
24 
25 //==================================================== Forward declarations
27 {
28  void getNonZeroBounds ( proshade_double* map, proshade_signed xDim, proshade_signed yDim, proshade_signed zDim, proshade_signed*& ret );
29 }
30 
38 bool ProSHADE_internal_io::isFilePDB ( std::string fName )
39 {
40  //================================================ Try reading the file using Gemmi
41  try
42  {
43  gemmi::Structure structure = gemmi::read_structure ( gemmi::MaybeGzipped ( fName ) );
44  }
45  catch ( std::runtime_error& e )
46  {
47  //============================================ Supress MSVC C4101 Unferenced variable warning / clang and gcc -Wunused-exception-parameter
48  (void)e;
49 
50  //============================================ Read failed. Done
51  return ( false );
52  }
53 
54  //================================================ Read successfull. Done
55  return ( true );
56 
57 }
58 
66 bool ProSHADE_internal_io::isFileMAP ( std::string fName )
67 {
68  gemmi::Ccp4<float> map;
69  try
70  {
71  map.read_ccp4 ( gemmi::MaybeGzipped (fName.c_str() ) );
72  }
73  catch ( std::runtime_error& e )
74  {
75  //============================================ Supress MSVC C4101 Unferenced variable warning / clang and gcc -Wunused-exception-parameter
76  (void)e;
77 
78  //============================================ Failed to read the map
79  return ( false );
80  }
81 
82  //================================================ Done
83  return ( true );
84 
85 }
86 
115 void ProSHADE_internal_io::readInMapHeader ( gemmi::Ccp4<float> *map, proshade_unsign *xDimInds, proshade_unsign *yDimInds, proshade_unsign *zDimInds, proshade_single *xDim, proshade_single *yDim, proshade_single *zDim, proshade_single *aAng, proshade_single *bAng, proshade_single *cAng, proshade_signed *xFrom, proshade_signed *yFrom, proshade_signed *zFrom, proshade_signed *xAxOrigin, proshade_signed *yAxOrigin, proshade_signed *zAxOrigin, proshade_unsign *xAxOrder, proshade_unsign *yAxOrder, proshade_unsign *zAxOrder, proshade_unsign *xGridInds, proshade_unsign *yGridInds, proshade_unsign *zGridInds )
116 {
117  //================================================ Read in the map file header
118  *xDimInds = static_cast<proshade_unsign> ( map->header_i32 ( 1 ) );
119  *yDimInds = static_cast<proshade_unsign> ( map->header_i32 ( 2 ) );
120  *zDimInds = static_cast<proshade_unsign> ( map->header_i32 ( 3 ) );
121 
122  *xFrom = static_cast<proshade_signed> ( map->header_i32 ( 5 ) );
123  *yFrom = static_cast<proshade_signed> ( map->header_i32 ( 6 ) );
124  *zFrom = static_cast<proshade_signed> ( map->header_i32 ( 7 ) );
125 
126  *xDim = static_cast<proshade_single> ( map->header_float ( 11 ) );
127  *yDim = static_cast<proshade_single> ( map->header_float ( 12 ) );
128  *zDim = static_cast<proshade_single> ( map->header_float ( 13 ) );
129 
130  *aAng = static_cast<proshade_single> ( map->header_float ( 14 ) );
131  *bAng = static_cast<proshade_single> ( map->header_float ( 15 ) );
132  *cAng = static_cast<proshade_single> ( map->header_float ( 16 ) );
133 
134  *xAxOrigin = static_cast<proshade_signed> ( map->header_i32 ( 50 ) ) + (*xFrom);
135  *yAxOrigin = static_cast<proshade_signed> ( map->header_i32 ( 51 ) ) + (*yFrom);
136  *zAxOrigin = static_cast<proshade_signed> ( map->header_i32 ( 52 ) ) + (*zFrom);
137 
138  *xAxOrder = static_cast<proshade_unsign> ( map->header_i32 ( 17 ) );
139  *yAxOrder = static_cast<proshade_unsign> ( map->header_i32 ( 18 ) );
140  *zAxOrder = static_cast<proshade_unsign> ( map->header_i32 ( 19 ) );
141 
142  *xGridInds = static_cast<proshade_unsign> ( map->header_i32 ( 8 ) );
143  *yGridInds = static_cast<proshade_unsign> ( map->header_i32 ( 9 ) );
144  *zGridInds = static_cast<proshade_unsign> ( map->header_i32 ( 10 ) );
145 
146  //================================================ Deal with sampling being different from cell size
147  if ( *xGridInds != *xDimInds )
148  {
149  *xDim = *xDim * ( static_cast<proshade_single> ( *xDimInds ) / static_cast<proshade_single> ( *xGridInds ) );
150  *xGridInds = *xDimInds;
151  }
152 
153  if ( *yGridInds != *yDimInds )
154  {
155  *yDim = *yDim * ( static_cast<proshade_single> ( *yDimInds ) / static_cast<proshade_single> ( *yGridInds ) );
156  *yGridInds = *yDimInds;
157  }
158 
159  if ( *zGridInds != *zDimInds )
160  {
161  *zDim = *zDim * ( static_cast<proshade_single> ( *zDimInds ) / static_cast<proshade_single> ( *zGridInds ) );
162  *zGridInds = *zDimInds;
163  }
164 
165  //================================================ Done
166  return ;
167 
168 }
169 
184 void ProSHADE_internal_io::readInMapData ( gemmi::Ccp4<float> *gemmiMap, proshade_double*& map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_unsign xAxOrder, proshade_unsign yAxOrder, proshade_unsign zAxOrder )
185 {
186  //================================================ Allocate internal variables
187  proshade_unsign *axOrdArr = new proshade_unsign[3];
188  proshade_unsign *axDimArr = new proshade_unsign[3];
189  proshade_unsign arrPos = 0;
190 
191  //================================================ Check memory allocation and fill in values
192  ProSHADE_internal_misc::checkMemoryAllocation ( axOrdArr, __FILE__, __LINE__, __func__ );
193  ProSHADE_internal_misc::checkMemoryAllocation ( axDimArr, __FILE__, __LINE__, __func__ );
194  axDimArr[0] = xDimInds;
195  axDimArr[1] = yDimInds;
196  axDimArr[2] = zDimInds;
197 
198  //================================================ Allocate the ProSHADE internal map variable memory
199  map = new proshade_double [xDimInds * yDimInds * zDimInds];
200  ProSHADE_internal_misc::checkMemoryAllocation ( map, __FILE__, __LINE__, __func__ );
201 
202  //================================================ Copy read in data to internal map variable
203  for ( axOrdArr[0] = 0; axOrdArr[0] < axDimArr[xAxOrder-1]; axOrdArr[0]++ )
204  {
205  for ( axOrdArr[1] = 0; axOrdArr[1] < axDimArr[yAxOrder-1]; axOrdArr[1]++ )
206  {
207  for ( axOrdArr[2] = 0; axOrdArr[2] < axDimArr[zAxOrder-1]; axOrdArr[2]++ )
208  {
209  arrPos = axOrdArr[2] + axDimArr[zAxOrder-1] * ( axOrdArr[1] + axDimArr[yAxOrder-1] * axOrdArr[0] );
210  map[arrPos] = static_cast< proshade_double > ( gemmiMap->grid.get_value_q( static_cast< int > ( axOrdArr[xAxOrder-1] ),
211  static_cast< int > ( axOrdArr[yAxOrder-1] ),
212  static_cast< int > ( axOrdArr[zAxOrder-1] ) ) );
213  }
214  }
215  }
216 
217  //================================================ Release internal variables memory
218  delete[] axDimArr;
219  delete[] axOrdArr;
220 
221  //================================================ Done
222  return ;
223 
224 }
225 
240 void ProSHADE_internal_io::readInMapData ( gemmi::Ccp4<int8_t> *gemmiMap, proshade_double*& map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_unsign xAxOrder, proshade_unsign yAxOrder, proshade_unsign zAxOrder )
241 {
242  //================================================ Allocate internal variables
243  proshade_unsign *axOrdArr = new proshade_unsign[3];
244  proshade_unsign *axDimArr = new proshade_unsign[3];
245  proshade_unsign arrPos = 0;
246 
247  //================================================ Check memory allocation and fill in values
248  ProSHADE_internal_misc::checkMemoryAllocation ( axOrdArr, __FILE__, __LINE__, __func__ );
249  ProSHADE_internal_misc::checkMemoryAllocation ( axDimArr, __FILE__, __LINE__, __func__ );
250  axDimArr[0] = xDimInds;
251  axDimArr[1] = yDimInds;
252  axDimArr[2] = zDimInds;
253 
254  //================================================ Allocate the ProSHADE internal map variable memory
255  map = new proshade_double [xDimInds * yDimInds * zDimInds];
256  ProSHADE_internal_misc::checkMemoryAllocation ( map, __FILE__, __LINE__, __func__ );
257 
258  //================================================ Copy read in data to internal map variable
259  for ( axOrdArr[0] = 0; axOrdArr[0] < axDimArr[xAxOrder-1]; axOrdArr[0]++ )
260  {
261  for ( axOrdArr[1] = 0; axOrdArr[1] < axDimArr[yAxOrder-1]; axOrdArr[1]++ )
262  {
263  for ( axOrdArr[2] = 0; axOrdArr[2] < axDimArr[zAxOrder-1]; axOrdArr[2]++ )
264  {
265  arrPos = axOrdArr[2] + axDimArr[zAxOrder-1] * ( axOrdArr[1] + axDimArr[yAxOrder-1] * axOrdArr[0] );
266  map[arrPos] = static_cast< proshade_double > ( gemmiMap->grid.get_value_q( static_cast< int > ( axOrdArr[xAxOrder-1] ),
267  static_cast< int > ( axOrdArr[yAxOrder-1] ),
268  static_cast< int > ( axOrdArr[zAxOrder-1] ) ) );
269  }
270  }
271  }
272 
273  //================================================ Release internal variables memory
274  delete[] axDimArr;
275  delete[] axOrdArr;
276 
277  //================================================ Done
278  return ;
279 
280 }
281 
302 void ProSHADE_internal_io::applyMask ( proshade_double*& map, std::string maskFile, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_signed verbose, proshade_signed messageShift, std::vector< proshade_double >* calcBounds, proshade_double* maskArray, proshade_unsign maXInds, proshade_unsign maYInds, proshade_unsign maZInds )
303 {
304  //================================================ Report progress
305  std::stringstream hlpSS;
306  hlpSS << "Reading mask " << maskFile;
307  ProSHADE_internal_messages::printProgressMessage ( verbose, 2, hlpSS.str(), messageShift );
308 
309  //================================================ Are we reading from array or from file?
310  if ( ( maskArray != nullptr ) && ( maXInds != 0 ) && ( maYInds != 0 ) && ( maZInds != 0 ) )
311  {
312  //============================================ Array it is!
313  ProSHADE_internal_io::applyMaskFromArray ( map, xDimInds, yDimInds, zDimInds, maskArray, maXInds, maYInds, maZInds, calcBounds, verbose, messageShift );
314  }
315  else
316  {
317  //============================================ Check if filename was given
318  if ( maskFile == "" ) { ProSHADE_internal_messages::printProgressMessage ( verbose, 3, "No mask supplied." ); return; }
319 
320  //============================================ File it is! Open the mask
321  gemmi::Ccp4<float> mask;
322  mask.read_ccp4 ( gemmi::MaybeGzipped ( maskFile.c_str() ) );
323 
324  //============================================ Convert to XYZ and create complete mask, if need be
325  mask.setup ( 0.0f, gemmi::MapSetup::ReorderOnly );
326 
327  //============================================ Read in the rest of the mask file header
328  proshade_unsign xDI, yDI, zDI, xAOR, yAOR, zAOR, xGI, yGI, zGI;
329  proshade_single xDS, yDS, zDS, aA, bA, cA;
330  proshade_signed xF, yF, zF, xAO, yAO, zAO;
332  &xDI, &yDI, &zDI,
333  &xDS, &yDS, &zDS,
334  &aA, &bA, &cA,
335  &xF, &yF, &zF,
336  &xAO, &yAO, &zAO,
337  &xAOR, &yAOR, &zAOR,
338  &xGI, &yGI, &zGI );
339 
340  //============================================ Save the mask values to ProSHADE variable
341  proshade_double* internalMask = nullptr;
342  ProSHADE_internal_io::readInMapData ( &mask, internalMask, xDI, yDI, zDI, xAOR, yAOR, zAOR );
343 
344  //============================================ Apply mask from array
345  ProSHADE_internal_io::applyMaskFromArray ( map, xDimInds, yDimInds, zDimInds, internalMask, xDI, yDI, zDI, calcBounds, verbose, messageShift );
346 
347  //============================================ Release the memory
348  delete[] internalMask;
349  }
350 
351  //================================================ Done
352  return ;
353 
354 }
355 
375 void ProSHADE_internal_io::applyMaskFromArray ( proshade_double*& map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_double*& mask, proshade_unsign xDimIndsMsk, proshade_unsign yDimIndsMsk, proshade_unsign zDimIndsMsk, std::vector< proshade_double >* calcBounds, proshade_signed verbose, proshade_signed messageShift )
376 {
377  //================================================ Initialise local variables
378  size_t origVolume = xDimInds * yDimInds * zDimInds;
379  size_t newVolume = xDimIndsMsk * yDimIndsMsk * zDimIndsMsk;
380  proshade_double* maskFinal;
381 
382  //================================================ If mask has different number of indices than map, then re-sample mask
383  if ( ( xDimIndsMsk != xDimInds ) || ( yDimIndsMsk != yDimInds ) || ( zDimIndsMsk != zDimInds ) )
384  {
385  //============================================ Initialise variables
386  fftw_complex* origCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * newVolume ) );
387  fftw_complex* origCoeffsHKL = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * newVolume ) );
388  fftw_complex* modifCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
389  fftw_complex* modifCoeffsHKL = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
390  fftw_complex* inMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * newVolume ) );
391  fftw_complex* outMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
392 
393  //============================================ Check memory allocation
394  ProSHADE_internal_misc::checkMemoryAllocation ( inMap, __FILE__, __LINE__, __func__ );
395  ProSHADE_internal_misc::checkMemoryAllocation ( outMap, __FILE__, __LINE__, __func__ );
396  ProSHADE_internal_misc::checkMemoryAllocation ( origCoeffs, __FILE__, __LINE__, __func__ );
397  ProSHADE_internal_misc::checkMemoryAllocation ( modifCoeffs, __FILE__, __LINE__, __func__ );
398  ProSHADE_internal_misc::checkMemoryAllocation ( origCoeffsHKL, __FILE__, __LINE__, __func__ );
399  ProSHADE_internal_misc::checkMemoryAllocation ( modifCoeffsHKL, __FILE__, __LINE__, __func__ );
400 
401  //============================================ Set array to zeroes
402  for ( size_t iter = 0; iter < origVolume; iter++ ) { modifCoeffsHKL[iter][0] = 0.0; modifCoeffsHKL[iter][1] = 0.0; }
403 
404  //============================================ Cope mask to Fourier input array
405  for ( size_t iter = 0; iter < newVolume; iter++ ) { inMap[iter][0] = mask[iter]; inMap[iter][1] = 0.0; }
406 
407  //============================================ Prepare Fourier transform plans
408  fftw_plan planForwardFourier = fftw_plan_dft_3d ( static_cast< int > ( xDimIndsMsk ), static_cast< int > ( yDimIndsMsk ), static_cast< int > ( zDimIndsMsk ), inMap, origCoeffs, FFTW_FORWARD, FFTW_ESTIMATE );
409  fftw_plan inverseFoourier = fftw_plan_dft_3d ( static_cast< int > ( xDimInds ), static_cast< int > ( yDimInds ), static_cast< int > ( zDimInds ), modifCoeffs, outMap, FFTW_BACKWARD, FFTW_ESTIMATE );
410 
411  //============================================ Compute pre and post changes
412  proshade_signed xPre, yPre, zPre;
413  xPre = std::abs ( ( static_cast< proshade_signed > ( xDimInds ) - static_cast< proshade_signed > ( xDimIndsMsk ) ) / 2 );
414  yPre = std::abs ( ( static_cast< proshade_signed > ( yDimInds ) - static_cast< proshade_signed > ( yDimIndsMsk ) ) / 2 );
415  zPre = std::abs ( ( static_cast< proshade_signed > ( zDimInds ) - static_cast< proshade_signed > ( zDimIndsMsk ) ) / 2 );
416 
417  if ( ( ( static_cast< proshade_signed > ( xDimInds ) - static_cast< proshade_signed > ( xDimIndsMsk ) ) % 2 ) == 1 ) { xPre -= 1; }
418  if ( ( ( static_cast< proshade_signed > ( yDimInds ) - static_cast< proshade_signed > ( yDimIndsMsk ) ) % 2 ) == 1 ) { yPre -= 1; }
419  if ( ( ( static_cast< proshade_signed > ( zDimInds ) - static_cast< proshade_signed > ( zDimIndsMsk ) ) % 2 ) == 1 ) { zPre -= 1; }
420 
421  //============================================ Run forward Fourier
422  fftw_execute ( planForwardFourier );
423 
424  //============================================ Initialise local variables
425  proshade_signed maskMapIndex = 0;
426  proshade_signed densMapIndex = 0;
427  proshade_signed xMaskPos, yMaskPos, zMaskPos, xDensPos, yDensPos, zDensPos;
428  proshade_signed maskH, maskK, maskL;
429 
430  //============================================ Convert mask to HKL for re-boxing
431  for ( proshade_signed xIt = 0; xIt < static_cast< proshade_signed > ( xDimIndsMsk ); xIt++ )
432  {
433  for ( proshade_signed yIt = 0; yIt < static_cast< proshade_signed > ( yDimIndsMsk ); yIt++ )
434  {
435  for ( proshade_signed zIt = 0; zIt < static_cast< proshade_signed > ( zDimIndsMsk ); zIt++ )
436  {
437  //================================ Convert to HKL
438  maskH = xIt + static_cast< proshade_signed > ( (xDimIndsMsk+1) / 2 ); if ( maskH >= static_cast< proshade_signed > ( xDimIndsMsk ) ) { maskH -= xDimIndsMsk; }
439  maskK = yIt + static_cast< proshade_signed > ( (yDimIndsMsk+1) / 2 ); if ( maskK >= static_cast< proshade_signed > ( yDimIndsMsk ) ) { maskK -= yDimIndsMsk; }
440  maskL = zIt + static_cast< proshade_signed > ( (zDimIndsMsk+1) / 2 ); if ( maskL >= static_cast< proshade_signed > ( zDimIndsMsk ) ) { maskL -= zDimIndsMsk; }
441 
442  //================================ Find the positions
443  maskMapIndex = zIt + static_cast< proshade_signed > ( zDimIndsMsk ) * ( yIt + static_cast< proshade_signed > ( yDimIndsMsk ) * xIt );
444  densMapIndex = maskL + static_cast< proshade_signed > ( zDimIndsMsk ) * ( maskK + static_cast< proshade_signed > ( yDimIndsMsk ) * maskH );
445 
446  //================================ Save the values
447  origCoeffsHKL[densMapIndex][0] = origCoeffs[maskMapIndex][0];
448  origCoeffsHKL[densMapIndex][1] = origCoeffs[maskMapIndex][1];
449  }
450  }
451  }
452 
453  //============================================ Rebox
454  for ( proshade_signed xIt = 0; xIt < static_cast< proshade_signed > ( xDimInds ); xIt++ )
455  {
456  for ( proshade_signed yIt = 0; yIt < static_cast< proshade_signed > ( yDimInds ); yIt++ )
457  {
458  for ( proshade_signed zIt = 0; zIt < static_cast< proshade_signed > ( zDimInds ); zIt++ )
459  {
460  //================================ Deal with X
461  if ( xDimIndsMsk >= xDimInds ) { xMaskPos = xIt + xPre; }
462  else { xMaskPos = xIt - xPre; }
463  xDensPos = xIt;
464 
465  //================================ Deal with Y
466  if ( yDimIndsMsk >= yDimInds ) { yMaskPos = yIt + yPre; }
467  else { yMaskPos = yIt - yPre; }
468  yDensPos = yIt;
469 
470  //================================ Deal with Z
471  if ( zDimIndsMsk >= zDimInds ) { zMaskPos = zIt + zPre; }
472  else { zMaskPos = zIt - zPre; }
473  zDensPos = zIt;
474 
475  //================================ Skip if mask value not available (because the modifCoeffsHKL array is zeroed, we do not need to do anything here)
476  if ( ( xMaskPos < 0 ) || ( xMaskPos >= static_cast< proshade_signed > ( xDimIndsMsk ) ) ) { continue; }
477  if ( ( yMaskPos < 0 ) || ( yMaskPos >= static_cast< proshade_signed > ( yDimIndsMsk ) ) ) { continue; }
478  if ( ( zMaskPos < 0 ) || ( zMaskPos >= static_cast< proshade_signed > ( zDimIndsMsk ) ) ) { continue; }
479 
480  //================================ Find the positions
481  maskMapIndex = zMaskPos + static_cast< proshade_signed > ( zDimIndsMsk ) * ( yMaskPos + static_cast< proshade_signed > ( yDimIndsMsk ) * xMaskPos );
482  densMapIndex = zDensPos + static_cast< proshade_signed > ( zDimInds ) * ( yDensPos + static_cast< proshade_signed > ( yDimInds ) * xDensPos );
483 
484  //================================ Copy values
485  modifCoeffsHKL[densMapIndex][0] = origCoeffsHKL[maskMapIndex][0];
486  modifCoeffsHKL[densMapIndex][1] = origCoeffsHKL[maskMapIndex][1];
487  }
488  }
489  }
490 
491  //============================================ Convert mask back to FFTW order
492  for ( proshade_signed xIt = 0; xIt < static_cast< proshade_signed > ( xDimInds ); xIt++ )
493  {
494  for ( proshade_signed yIt = 0; yIt < static_cast< proshade_signed > ( yDimInds ); yIt++ )
495  {
496  for ( proshade_signed zIt = 0; zIt < static_cast< proshade_signed > ( zDimInds ); zIt++ )
497  {
498  //================================ Convert to HKL
499  maskH = xIt + static_cast< proshade_signed > ( xDimInds / 2 ); if ( maskH >= static_cast< proshade_signed > ( xDimInds ) ) { maskH -= xDimInds; }
500  maskK = yIt + static_cast< proshade_signed > ( yDimInds / 2 ); if ( maskK >= static_cast< proshade_signed > ( yDimInds ) ) { maskK -= yDimInds; }
501  maskL = zIt + static_cast< proshade_signed > ( zDimInds / 2 ); if ( maskL >= static_cast< proshade_signed > ( zDimInds ) ) { maskL -= zDimInds; }
502 
503  //================================ Find the positions
504  maskMapIndex = zIt + static_cast< proshade_signed > ( zDimInds ) * ( yIt + static_cast< proshade_signed > ( yDimInds ) * xIt );
505  densMapIndex = maskL + static_cast< proshade_signed > ( zDimInds ) * ( maskK + static_cast< proshade_signed > ( yDimInds ) * maskH );
506 
507  //================================ Save the values
508  modifCoeffs[densMapIndex][0] = modifCoeffsHKL[maskMapIndex][0];
509  modifCoeffs[densMapIndex][1] = modifCoeffsHKL[maskMapIndex][1];
510  }
511  }
512  }
513 
514  //============================================ Run inverse Fourier on the modified coefficients
515  fftw_execute ( inverseFoourier );
516 
517  //============================================ Delete old mask and allocate memory for the new, re-sampled mask
518  maskFinal = new proshade_double [origVolume];
519  ProSHADE_internal_misc::checkMemoryAllocation ( maskFinal, __FILE__, __LINE__, __func__ );
520 
521  //============================================ Copy results into a new, properly sampled mask
522  for ( size_t iter = 0; iter < origVolume; iter++ ) { maskFinal[iter] = outMap[iter][0]; }
523 
524  //============================================ Release remaining memory
525  fftw_destroy_plan ( planForwardFourier );
526  fftw_destroy_plan ( inverseFoourier );
527  fftw_free ( origCoeffs );
528  fftw_free ( modifCoeffs );
529  fftw_free ( origCoeffsHKL );
530  fftw_free ( modifCoeffsHKL );
531  fftw_free ( inMap );
532  fftw_free ( outMap );
533  }
534  else
535  {
536  maskFinal = new proshade_double [origVolume];
537  ProSHADE_internal_misc::checkMemoryAllocation ( maskFinal, __FILE__, __LINE__, __func__ );
538  for ( size_t iter = 0; iter < origVolume; iter++ ) { maskFinal[iter] = mask[iter]; }
539  }
540 
541  //================================================ Apply the mask to the map
542  for ( size_t iter = 0; iter < static_cast< size_t > ( xDimInds * yDimInds * zDimInds ); iter++ ) { map[iter] *= maskFinal[iter]; }
543 
544  //================================================ Report progress
545  ProSHADE_internal_messages::printProgressMessage ( verbose, 3, "Mask read in and applied successfully.", messageShift );
546 
547  //================================================ Determine bounds for computation
548  proshade_signed* bnds = new proshade_signed[6];
549  ProSHADE_internal_misc::checkMemoryAllocation ( bnds, __FILE__, __LINE__, __func__ );
550  ProSHADE_internal_mapManip::getNonZeroBounds ( maskFinal, static_cast< proshade_signed > ( xDimInds ), static_cast< proshade_signed > ( yDimInds ), static_cast< proshade_signed > ( zDimInds ), bnds );
551 
552  //================================================ Copy to settings
553  calcBounds->at(0) = static_cast< proshade_double > ( bnds[1] - bnds[0] );
554  calcBounds->at(1) = static_cast< proshade_double > ( bnds[3] - bnds[2] );
555  calcBounds->at(2) = static_cast< proshade_double > ( bnds[5] - bnds[4] );
556 
557  //================================================ Report progress
558  ProSHADE_internal_messages::printProgressMessage ( verbose, 4, "Computation boundaries determined.", messageShift );
559 
560  //================================================ Release memory
561  delete[] bnds;
562  delete[] maskFinal;
563 
564  //================================================ Done
565  return ;
566 
567 }
568 
588 void ProSHADE_internal_io::applyWeights ( proshade_double*& map, std::string weightsFile, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_signed verbose, proshade_signed messageShift, proshade_double* weightsArray, proshade_unsign waXInds, proshade_unsign waYInds, proshade_unsign waZInds )
589 {
590  //================================================ Report progress
591  std::stringstream hlpSS;
592  hlpSS << "Reading weights " << weightsFile;
593  ProSHADE_internal_messages::printProgressMessage ( verbose, 2, hlpSS.str(), messageShift );
594 
595  //================================================ Are we reading from file, or from array?
596  if ( ( weightsArray != nullptr ) && ( waXInds != 0 ) && ( waYInds != 0 ) && ( waZInds != 0 ) )
597  {
598  //============================================ From array it is!
599  ProSHADE_internal_io::applyWeightsFromArray ( map, xDimInds, yDimInds, zDimInds, weightsArray, waXInds, waYInds, waZInds, verbose, messageShift );
600  }
601  else
602  {
603  //============================================ Check if weights file was given
604  if ( weightsFile == "" ) { ProSHADE_internal_messages::printProgressMessage ( verbose, 3, "No weights supplied. Assuming all weights to be 1.0.", messageShift ); return; }
605 
606  //============================================ From file it is! Open the weights file
607  gemmi::Ccp4<float> weights;
608  weights.read_ccp4 ( gemmi::MaybeGzipped ( weightsFile.c_str() ) );
609 
610  //============================================ Convert to XYZ and create complete weights, if need be
611  weights.setup ( 0.0f, gemmi::MapSetup::ReorderOnly );
612 
613  //============================================ Read in the rest of the weights file header
614  proshade_unsign xDI, yDI, zDI, xAOR, yAOR, zAOR, xGI, yGI, zGI;
615  proshade_single xDS, yDS, zDS, aA, bA, cA;
616  proshade_signed xF, yF, zF, xAO, yAO, zAO;
618  &xDI, &yDI, &zDI,
619  &xDS, &yDS, &zDS,
620  &aA, &bA, &cA,
621  &xF, &yF, &zF,
622  &xAO, &yAO, &zAO,
623  &xAOR, &yAOR, &zAOR,
624  &xGI, &yGI, &zGI );
625 
626  //============================================ Save the weights values to local variable
627  proshade_double* internalWeights = nullptr;
628  ProSHADE_internal_io::readInMapData ( &weights, internalWeights, xDI, yDI, zDI, xAOR, yAOR, zAOR );
629 
630  //============================================ Apply weights from array
631  ProSHADE_internal_io::applyWeightsFromArray ( map, xDimInds, yDimInds, zDimInds, internalWeights, xDI, yDI, zDI, verbose, messageShift );
632 
633  //============================================ Release the memory
634  delete[] internalWeights;
635  }
636 
637  //================================================ Done
638  return ;
639 
640 }
641 
660 void ProSHADE_internal_io::applyWeightsFromArray ( proshade_double*& map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_double*& weights, proshade_unsign xDimIndsWgh, proshade_unsign yDimIndsWgh, proshade_unsign zDimIndsWgh, proshade_signed verbose, proshade_signed messageShift )
661 {
662  //================================================ Initialise local variables
663  proshade_double* weightsFinal;
664  size_t origVolume = xDimInds * yDimInds * zDimInds;
665  size_t newVolume = xDimIndsWgh * yDimIndsWgh * zDimIndsWgh;
666 
667  //================================================ If weights have different number of indices than map, then re-sample weights in supplied space
668  if ( ( xDimIndsWgh != xDimInds ) || ( yDimIndsWgh != yDimInds ) || ( zDimIndsWgh != zDimInds ) )
669  {
670  //============================================ Initialise variables
671  fftw_complex* origCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * newVolume ) );
672  fftw_complex* origCoeffsHKL = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * newVolume ) );
673  fftw_complex* modifCoeffs = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
674  fftw_complex* modifCoeffsHKL = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
675  fftw_complex* inMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * newVolume ) );
676  fftw_complex* outMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
677 
678  //============================================ Check memory allocation
679  ProSHADE_internal_misc::checkMemoryAllocation ( inMap, __FILE__, __LINE__, __func__ );
680  ProSHADE_internal_misc::checkMemoryAllocation ( outMap, __FILE__, __LINE__, __func__ );
681  ProSHADE_internal_misc::checkMemoryAllocation ( origCoeffs, __FILE__, __LINE__, __func__ );
682  ProSHADE_internal_misc::checkMemoryAllocation ( modifCoeffs, __FILE__, __LINE__, __func__ );
683  ProSHADE_internal_misc::checkMemoryAllocation ( origCoeffsHKL, __FILE__, __LINE__, __func__ );
684  ProSHADE_internal_misc::checkMemoryAllocation ( modifCoeffsHKL, __FILE__, __LINE__, __func__ );
685 
686  //============================================ Set array to zeroes
687  for ( size_t iter = 0; iter < origVolume; iter++ ) { modifCoeffsHKL[iter][0] = 0.0; modifCoeffsHKL[iter][1] = 0.0; }
688 
689  //============================================ Copy weights to Fourier input array
690  for ( size_t iter = 0; iter < newVolume; iter++ ) { inMap[iter][0] = weights[iter]; inMap[iter][1] = 0.0; }
691 
692  //============================================ Prepare Fourier transform plans
693  fftw_plan planForwardFourier = fftw_plan_dft_3d ( static_cast< int > ( xDimIndsWgh ), static_cast< int > ( yDimIndsWgh ), static_cast< int > ( zDimIndsWgh ), inMap, origCoeffs, FFTW_FORWARD, FFTW_ESTIMATE );
694  fftw_plan inverseFoourier = fftw_plan_dft_3d ( static_cast< int > ( xDimInds ), static_cast< int > ( yDimInds ), static_cast< int > ( zDimInds ), modifCoeffs, outMap, FFTW_BACKWARD, FFTW_ESTIMATE );
695 
696  //============================================ Compute pre and post changes
697  proshade_signed xPre, yPre, zPre;
698  xPre = std::abs ( ( static_cast< proshade_signed > ( xDimInds ) - static_cast< proshade_signed > ( xDimIndsWgh ) ) / 2 );
699  yPre = std::abs ( ( static_cast< proshade_signed > ( yDimInds ) - static_cast< proshade_signed > ( yDimIndsWgh ) ) / 2 );
700  zPre = std::abs ( ( static_cast< proshade_signed > ( zDimInds ) - static_cast< proshade_signed > ( zDimIndsWgh ) ) / 2 );
701 
702  if ( ( ( static_cast< proshade_signed > ( xDimInds ) - static_cast< proshade_signed > ( xDimIndsWgh ) ) % 2 ) == 1 ) { xPre -= 1; }
703  if ( ( ( static_cast< proshade_signed > ( yDimInds ) - static_cast< proshade_signed > ( yDimIndsWgh ) ) % 2 ) == 1 ) { yPre -= 1; }
704  if ( ( ( static_cast< proshade_signed > ( zDimInds ) - static_cast< proshade_signed > ( zDimIndsWgh ) ) % 2 ) == 1 ) { zPre -= 1; }
705 
706  //============================================ Run forward Fourier
707  fftw_execute ( planForwardFourier );
708 
709  //============================================ Initialise local variables
710  proshade_signed maskMapIndex = 0;
711  proshade_signed densMapIndex = 0;
712  proshade_signed xMaskPos, yMaskPos, zMaskPos, xDensPos, yDensPos, zDensPos;
713  proshade_signed maskH, maskK, maskL;
714 
715  //============================================ Convert weights to HKL for re-boxing
716  for ( proshade_signed xIt = 0; xIt < static_cast< proshade_signed > ( xDimIndsWgh ); xIt++ )
717  {
718  for ( proshade_signed yIt = 0; yIt < static_cast< proshade_signed > ( yDimIndsWgh ); yIt++ )
719  {
720  for ( proshade_signed zIt = 0; zIt < static_cast< proshade_signed > ( zDimIndsWgh ); zIt++ )
721  {
722  //================================ Convert to HKL
723  maskH = xIt + static_cast< proshade_signed > ( (xDimIndsWgh+1) / 2 ); if ( maskH >= static_cast< proshade_signed > ( xDimIndsWgh ) ) { maskH -= xDimIndsWgh; }
724  maskK = yIt + static_cast< proshade_signed > ( (yDimIndsWgh+1) / 2 ); if ( maskK >= static_cast< proshade_signed > ( yDimIndsWgh ) ) { maskK -= yDimIndsWgh; }
725  maskL = zIt + static_cast< proshade_signed > ( (zDimIndsWgh+1) / 2 ); if ( maskL >= static_cast< proshade_signed > ( zDimIndsWgh ) ) { maskL -= zDimIndsWgh; }
726 
727  //================================ Find the positions
728  maskMapIndex = zIt + static_cast< proshade_signed > ( zDimIndsWgh ) * ( yIt + static_cast< proshade_signed > ( yDimIndsWgh ) * xIt );
729  densMapIndex = maskL + static_cast< proshade_signed > ( zDimIndsWgh ) * ( maskK + static_cast< proshade_signed > ( yDimIndsWgh ) * maskH );
730 
731  //================================ Save the values
732  origCoeffsHKL[densMapIndex][0] = origCoeffs[maskMapIndex][0];
733  origCoeffsHKL[densMapIndex][1] = origCoeffs[maskMapIndex][1];
734  }
735  }
736  }
737 
738  //============================================ Rebox
739  for ( proshade_signed xIt = 0; xIt < static_cast< proshade_signed > ( xDimInds ); xIt++ )
740  {
741  for ( proshade_signed yIt = 0; yIt < static_cast< proshade_signed > ( yDimInds ); yIt++ )
742  {
743  for ( proshade_signed zIt = 0; zIt < static_cast< proshade_signed > ( zDimInds ); zIt++ )
744  {
745  //================================ Deal with X
746  if ( xDimIndsWgh >= xDimInds ) { xMaskPos = xIt + xPre; }
747  else { xMaskPos = xIt - xPre; }
748  xDensPos = xIt;
749 
750  //================================ Deal with Y
751  if ( yDimIndsWgh >= yDimInds ) { yMaskPos = yIt + yPre; }
752  else { yMaskPos = yIt - yPre; }
753  yDensPos = yIt;
754 
755  //================================ Deal with Z
756  if ( zDimIndsWgh >= zDimInds ) { zMaskPos = zIt + zPre; }
757  else { zMaskPos = zIt - zPre; }
758  zDensPos = zIt;
759 
760  //================================ Skip if weights value not available (because the modifCoeffsHKL array is zeroed, we do not need to do anything here)
761  if ( ( xMaskPos < 0 ) || ( xMaskPos >= static_cast< proshade_signed > ( xDimIndsWgh ) ) ) { continue; }
762  if ( ( yMaskPos < 0 ) || ( yMaskPos >= static_cast< proshade_signed > ( yDimIndsWgh ) ) ) { continue; }
763  if ( ( zMaskPos < 0 ) || ( zMaskPos >= static_cast< proshade_signed > ( zDimIndsWgh ) ) ) { continue; }
764 
765  //================================ Find the positions
766  maskMapIndex = zMaskPos + static_cast< proshade_signed > ( zDimIndsWgh ) * ( yMaskPos + static_cast< proshade_signed > ( yDimIndsWgh ) * xMaskPos );
767  densMapIndex = zDensPos + static_cast< proshade_signed > ( zDimInds ) * ( yDensPos + static_cast< proshade_signed > ( yDimInds ) * xDensPos );
768 
769  //================================ Copy values
770  modifCoeffsHKL[densMapIndex][0] = origCoeffsHKL[maskMapIndex][0];
771  modifCoeffsHKL[densMapIndex][1] = origCoeffsHKL[maskMapIndex][1];
772  }
773  }
774  }
775 
776  //============================================ Convert weights back to FFTW order
777  for ( proshade_signed xIt = 0; xIt < static_cast< proshade_signed > ( xDimInds ); xIt++ )
778  {
779  for ( proshade_signed yIt = 0; yIt < static_cast< proshade_signed > ( yDimInds ); yIt++ )
780  {
781  for ( proshade_signed zIt = 0; zIt < static_cast< proshade_signed > ( zDimInds ); zIt++ )
782  {
783  //================================ Convert to HKL
784  maskH = xIt + static_cast< proshade_signed > ( xDimInds / 2 ); if ( maskH >= static_cast< proshade_signed > ( xDimInds ) ) { maskH -= xDimInds; }
785  maskK = yIt + static_cast< proshade_signed > ( yDimInds / 2 ); if ( maskK >= static_cast< proshade_signed > ( yDimInds ) ) { maskK -= yDimInds; }
786  maskL = zIt + static_cast< proshade_signed > ( zDimInds / 2 ); if ( maskL >= static_cast< proshade_signed > ( zDimInds ) ) { maskL -= zDimInds; }
787 
788  //================================ Find the positions
789  maskMapIndex = zIt + static_cast< proshade_signed > ( zDimInds ) * ( yIt + static_cast< proshade_signed > ( yDimInds ) * xIt );
790  densMapIndex = maskL + static_cast< proshade_signed > ( zDimInds ) * ( maskK + static_cast< proshade_signed > ( yDimInds ) * maskH );
791 
792  //================================ Save the values
793  modifCoeffs[densMapIndex][0] = modifCoeffsHKL[maskMapIndex][0];
794  modifCoeffs[densMapIndex][1] = modifCoeffsHKL[maskMapIndex][1];
795  }
796  }
797  }
798 
799  //============================================ Run inverse Fourier on the modified coefficients
800  fftw_execute ( inverseFoourier );
801 
802  //============================================ Delete old weights and allocate memory for the new, re-sampled weights
803  weightsFinal = new proshade_double [origVolume];
804  ProSHADE_internal_misc::checkMemoryAllocation ( weightsFinal, __FILE__, __LINE__, __func__ );
805 
806  //============================================ Copy results into a new, properly sampled weights
807  for ( size_t iter = 0; iter < origVolume; iter++ ) { weightsFinal[iter] = outMap[iter][0]; }
808 
809  //============================================ Release remaining memory
810  fftw_destroy_plan ( planForwardFourier );
811  fftw_destroy_plan ( inverseFoourier );
812  fftw_free ( origCoeffs );
813  fftw_free ( modifCoeffs );
814  fftw_free ( origCoeffsHKL );
815  fftw_free ( modifCoeffsHKL );
816  fftw_free ( inMap );
817  fftw_free ( outMap );
818  }
819  else
820  {
821  weightsFinal = new proshade_double [origVolume];
822  ProSHADE_internal_misc::checkMemoryAllocation ( weightsFinal, __FILE__, __LINE__, __func__ );
823  for ( size_t iter = 0; iter < origVolume; iter++ ) { weightsFinal[iter] = weights[iter]; }
824  }
825 
826  //================================================ Allocate memory for map Fourier transform
827  fftw_complex* inMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
828  fftw_complex* outMap = reinterpret_cast< fftw_complex* > ( fftw_malloc ( sizeof ( fftw_complex ) * origVolume ) );
829  ProSHADE_internal_misc::checkMemoryAllocation ( inMap, __FILE__, __LINE__, __func__ );
830  ProSHADE_internal_misc::checkMemoryAllocation ( outMap, __FILE__, __LINE__, __func__ );
831  fftw_plan planForwardFourier = fftw_plan_dft_3d ( static_cast< int > ( xDimInds ), static_cast< int > ( yDimInds ), static_cast< int > ( zDimInds ), inMap, outMap, FFTW_FORWARD, FFTW_ESTIMATE );
832  fftw_plan inverseFoourier = fftw_plan_dft_3d ( static_cast< int > ( xDimInds ), static_cast< int > ( yDimInds ), static_cast< int > ( zDimInds ), outMap, inMap, FFTW_BACKWARD, FFTW_ESTIMATE );
833 
834  //================================================ Set data
835  for ( size_t iter = 0; iter < static_cast< size_t > ( origVolume ); iter++ ) { inMap[iter][0] = map[iter]; inMap[iter][1] = 0.0; }
836 
837  //================================================ Convert map to Fourier space
838  fftw_execute ( planForwardFourier );
839 
840  //================================================ Apply the weights to the map in Fourier space
841  proshade_double normFactor = static_cast<proshade_double> ( origVolume );
842  for ( size_t iter = 0; iter < static_cast< size_t > ( origVolume ); iter++ ) { outMap[iter][0] *= weightsFinal[iter] / normFactor; outMap[iter][1] *= weightsFinal[iter] / normFactor; }
843 
844  //================================================ Convert weighted map from Fourier space
845  fftw_execute ( inverseFoourier );
846 
847  //================================================ Copy results to map
848  for ( size_t iter = 0; iter < static_cast< size_t > ( xDimInds * yDimInds * zDimInds ); iter++ ) { map[iter] = inMap[iter][0]; }
849 
850  //================================================ Release memory
851  delete[] weightsFinal;
852  fftw_free ( inMap );
853  fftw_free ( outMap );
854  fftw_destroy_plan ( planForwardFourier );
855  fftw_destroy_plan ( inverseFoourier );
856 
857  //================================================ Report progress
858  ProSHADE_internal_messages::printProgressMessage ( verbose, 3, "Mask read in and applied successfully.", messageShift );
859 
860  //================================================ Done
861  return ;
862 
863 }
864 
895 void ProSHADE_internal_io::writeOutMapHeader ( gemmi::Ccp4<float> *map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_single xDim, proshade_single yDim, proshade_single zDim, proshade_single aAng, proshade_single bAng, proshade_single cAng, proshade_signed xFrom, proshade_signed yFrom, proshade_signed zFrom, proshade_signed xAxOrigin, proshade_signed yAxOrigin, proshade_signed zAxOrigin, proshade_unsign xAxOrder, proshade_unsign yAxOrder, proshade_unsign zAxOrder, proshade_unsign xGridInds, proshade_unsign yGridInds, proshade_unsign zGridInds, std::string title, int mode )
896 {
897  //================================================ Fill in the map file header
898  map->set_header_i32 ( 1 , static_cast<int32_t> ( xDimInds ) ); // Number of columns in 3D data array (fast axis)
899  map->set_header_i32 ( 2 , static_cast<int32_t> ( yDimInds ) ); // Number of columns in 3D data array (medium axis)
900  map->set_header_i32 ( 3 , static_cast<int32_t> ( zDimInds ) ); // Number of columns in 3D data array (slow axis)
901  map->set_header_i32 ( 4 , static_cast<int32_t> ( mode ) ); // Map mode
902  map->set_header_i32 ( 5 , static_cast<int32_t> ( xFrom ) ); // Starting index (fast axis)
903  map->set_header_i32 ( 6 , static_cast<int32_t> ( yFrom ) ); // Starting index (medium axis)
904  map->set_header_i32 ( 7 , static_cast<int32_t> ( zFrom ) ); // Starting index (slow axis)
905  map->set_header_i32 ( 8 , static_cast<int32_t> ( xGridInds ) ); // Grid sampling (fast axis)
906  map->set_header_i32 ( 9 , static_cast<int32_t> ( yGridInds ) ); // Grid sampling (medium axis)
907  map->set_header_i32 ( 10, static_cast<int32_t> ( zGridInds ) ); // Grid sampling (slow axis)
908  map->set_header_float ( 11, static_cast<float> ( xDim ) ); // Grid dimension in Angstrom (fast axis)
909  map->set_header_float ( 12, static_cast<float> ( yDim ) ); // Grid dimension in Angstrom (medium axis)
910  map->set_header_float ( 13, static_cast<float> ( zDim ) ); // Grid dimension in Angstrom (slow axis)
911  map->set_header_float ( 14, static_cast<float> ( aAng ) ); // Alpha angle in degrees
912  map->set_header_float ( 15, static_cast<float> ( bAng ) ); // Beta angle in degrees
913  map->set_header_float ( 16, static_cast<float> ( cAng ) ); // Gamma angle in degrees
914  map->set_header_i32 ( 17, static_cast<int32_t> ( xAxOrder ) ); // MAPC
915  map->set_header_i32 ( 18, static_cast<int32_t> ( yAxOrder ) ); // MAPR
916  map->set_header_i32 ( 19, static_cast<int32_t> ( zAxOrder ) ); // MAPS
917  if ( map->grid.spacegroup ) { map->set_header_i32 ( 23, static_cast<int32_t> ( map->grid.spacegroup->ccp4 ) ); } // Space group
918  else { map->set_header_i32 ( 23, static_cast<int32_t> ( 1 ) ); }
919  map->set_header_i32 ( 24, static_cast<int32_t> ( map->grid.spacegroup->operations().order() * 80 ) ); // NSYMBT - size of extended header (which follows main header) in bytes
920  map->set_header_str ( 27, "CCP4" ); // Code for the type of extended header
921  map->set_header_i32 ( 28, static_cast<int32_t> ( 20140 ) ); // Version
922  map->set_header_i32 ( 50, static_cast<int32_t> ( xAxOrigin ) ); // Origin of the map (fast axis)
923  map->set_header_i32 ( 51, static_cast<int32_t> ( yAxOrigin ) ); // Origin of the map (medium axis)
924  map->set_header_i32 ( 52, static_cast<int32_t> ( zAxOrigin ) ); // Origin of the map (slow axis)
925  map->set_header_str ( 53, "MAP" ); // File format
926  if ( gemmi::is_little_endian() ) { map->set_header_i32 ( 54, static_cast<int32_t> ( 0x00004144 ) ); } // Machine stamp encoding byte ordering of data
927  else { map->set_header_i32 ( 54, static_cast<int32_t> ( 0x11110000 ) ); }
928  map->set_header_i32 ( 56, static_cast<int32_t> ( 1 ) ); // Number of labels used
929  std::memset ( reinterpret_cast<void*> ( &(map->ccp4_header.at( 56 )) ), ' ', static_cast< size_t > ( 800 + map->grid.spacegroup->operations().order() * 80 ) ); // 56 is used because the vector is indexed from 0
930  map->set_header_str ( 57, title ); // Title
931 
932  //================================================ Done
933  return ;
934 
935 }
936 
945 ProSHADE_internal_io::InputType ProSHADE_internal_io::figureDataType ( std::string fName )
946 {
947  //================================================ Try readin as PDB
948  if ( isFilePDB ( fName ) )
949  {
950  return ( PDB );
951  }
952 
953  //================================================ If not, try readin as MAP
954  if ( isFileMAP ( fName ) )
955  {
956  return ( MAP );
957  }
958 
959  //================================================ No luck? UNKNOWN it is ...
960  return ( UNKNOWN );
961 
962  //================================================ Done
963 
964 }
965 
987 void ProSHADE_internal_io::writeRotationTranslationJSON ( proshade_double trsX1, proshade_double trsY1, proshade_double trsZ1, proshade_double eulA, proshade_double eulB, proshade_double eulG, proshade_double trsX2, proshade_double trsY2, proshade_double trsZ2, std::string fileName )
988 {
989  //================================================ Open file for writing
990  std::ofstream jsonFile;
991  jsonFile.open ( fileName );
992 
993  //================================================ Check file opening success
994  if ( !jsonFile.is_open( ) )
995  {
996  throw ProSHADE_exception ( "Failed to open JSON output file.", "E000056", __FILE__, __LINE__, __func__, "Failed to open json file to which the rotation and\n : translation would be written into. Most likely cause is\n : lack of rights to write in the current folder." );
997  }
998 
999  //================================================ Get rotation matrix from Euler angles
1000  proshade_double* rotMat = new proshade_double[9];
1001  ProSHADE_internal_misc::checkMemoryAllocation ( rotMat, __FILE__, __LINE__, __func__ );
1003 
1004  //================================================ Write the info
1005  jsonFile << "{\n";
1006  jsonFile << " \"translationToOrigin\" : [ " << trsX1 << ", " << trsY1 << ", " << trsZ1 << " ], \n";
1007 
1008  jsonFile << " \"rotationMatrix:\" : [ " << rotMat[0] << ", " << rotMat[1] << ", " << rotMat[2] << ", \n";
1009  jsonFile << " " << rotMat[3] << ", " << rotMat[4] << ", " << rotMat[5] << ", \n";
1010  jsonFile << " " << rotMat[6] << ", " << rotMat[7] << ", " << rotMat[8] << "], \n";
1011 
1012  jsonFile << " \"translationFromRotCenToOverlay\" : [ " << trsX2 << ", " << trsY2 << ", " << trsZ2 << " ] \n";
1013  jsonFile << "}\n";
1014 
1015  //================================================ Close file
1016  jsonFile.close ( );
1017 
1018  //================================================ Release memory
1019  delete[] rotMat;
1020 
1021  //================================================ Done
1022  return ;
1023 
1024 }
ProSHADE_internal_io::isFilePDB
bool isFilePDB(std::string fName)
Function determining if the input data type is PDB.
Definition: ProSHADE_io.cpp:38
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_io::applyWeightsFromArray
void applyWeightsFromArray(proshade_double *&map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_double *&weights, proshade_unsign xDimIndsWgh, proshade_unsign yDimIndsWgh, proshade_unsign zDimIndsWgh, proshade_signed verbose, proshade_signed messageShift)
This function applies the weights to the map.
Definition: ProSHADE_io.cpp:660
ProSHADE_internal_io::isFileMAP
bool isFileMAP(std::string fName)
Function determining if the input data type is MAP.
Definition: ProSHADE_io.cpp:66
ProSHADE_internal_io::figureDataType
InputType figureDataType(std::string fName)
Function determining input data type.
Definition: ProSHADE_io.cpp:945
ProSHADE_io.hpp
This header file declares all the functions required for low level file format access.
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_io::readInMapData
void readInMapData(gemmi::Ccp4< float > *gemmiMap, proshade_double *&map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_unsign xAxOrder, proshade_unsign yAxOrder, proshade_unsign zAxOrder)
This function converts the gemmi Ccp4 object data to ProSHADE internal map representation.
Definition: ProSHADE_io.cpp:184
ProSHADE_internal_mapManip
This namespace contains the internal functions for manipulating maps already present in the internal ...
Definition: ProSHADE_io.cpp:27
ProSHADE_internal_io::applyMaskFromArray
void applyMaskFromArray(proshade_double *&map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_double *&mask, proshade_unsign xDimIndsMsk, proshade_unsign yDimIndsMsk, proshade_unsign zDimIndsMsk, std::vector< proshade_double > *calcBounds, proshade_signed verbose, proshade_signed messageShift)
This function applies the mask to the map.
Definition: ProSHADE_io.cpp:375
ProSHADE_internal_io::writeRotationTranslationJSON
void writeRotationTranslationJSON(proshade_double trsX1, proshade_double trsY1, proshade_double trsZ1, proshade_double eulA, proshade_double eulB, proshade_double eulG, proshade_double trsX2, proshade_double trsY2, proshade_double trsZ2, std::string fileName)
Function for writing out the optimal rotation and translation into a JSON file.
Definition: ProSHADE_io.cpp:987
ProSHADE_internal_io::writeOutMapHeader
void writeOutMapHeader(gemmi::Ccp4< float > *map, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_single xDim, proshade_single yDim, proshade_single zDim, proshade_single aAng, proshade_single bAng, proshade_single cAng, proshade_signed xFrom, proshade_signed yFrom, proshade_signed zFrom, proshade_signed xAxOrigin, proshade_signed yAxOrigin, proshade_signed zAxOrigin, proshade_unsign xAxOrder, proshade_unsign yAxOrder, proshade_unsign zAxOrder, proshade_unsign xGridInds, proshade_unsign yGridInds, proshade_unsign zGridInds, std::string title, int mode)
This function parses the CCP4 MAP file header as read in by gemmi.
Definition: ProSHADE_io.cpp:895
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_io::applyMask
void applyMask(proshade_double *&map, std::string maskFile, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_signed verbose, proshade_signed messageShift, std::vector< proshade_double > *calcBounds, proshade_double *maskArray=nullptr, proshade_unsign maXInds=0, proshade_unsign maYInds=0, proshade_unsign maZInds=0)
This function reads and applies the mask to the map.
Definition: ProSHADE_io.cpp:302
ProSHADE_internal_io::applyWeights
void applyWeights(proshade_double *&map, std::string weightsFile, proshade_unsign xDimInds, proshade_unsign yDimInds, proshade_unsign zDimInds, proshade_signed verbose, proshade_signed messageShift, proshade_double *weightsArray=nullptr, proshade_unsign waXInds=0, proshade_unsign waYInds=0, proshade_unsign waZInds=0)
This function reads and applies the Fourier weights to the map.
Definition: ProSHADE_io.cpp:588
ProSHADE_internal_messages::printProgressMessage
void printProgressMessage(proshade_signed verbose, proshade_signed messageLevel, std::string message, proshade_signed messageShift=0)
General stdout message printing.
Definition: ProSHADE_messages.cpp:71
ProSHADE_internal_io::readInMapHeader
void readInMapHeader(gemmi::Ccp4< float > *map, proshade_unsign *xDimInds, proshade_unsign *yDimInds, proshade_unsign *zDimInds, proshade_single *xDim, proshade_single *yDim, proshade_single *zDim, proshade_single *aAng, proshade_single *bAng, proshade_single *cAng, proshade_signed *xFrom, proshade_signed *yFrom, proshade_signed *zFrom, proshade_signed *xAxOrigin, proshade_signed *yAxOrigin, proshade_signed *zAxOrigin, proshade_unsign *xAxOrder, proshade_unsign *yAxOrder, proshade_unsign *zAxOrder, proshade_unsign *xGridInds, proshade_unsign *yGridInds, proshade_unsign *zGridInds)
This function parses the CCP4 MAP file header as read in by gemmi.
Definition: ProSHADE_io.cpp:115