 |
ProSHADE
0.7.6.6 (JUL 2022)
Protein Shape Detection
|
This namespace contains all the functions required for map overlays.
More...
|
| void | getOptimalRotation (ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *eulA, proshade_double *eulB, proshade_double *eulG) |
| | This function finds the optimal rotation between two structures as described by the settings object. More...
|
| |
| void | getOptimalTranslation (ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ, proshade_double eulA, proshade_double eulB, proshade_double eulG) |
| | This function finds the optimal translation between two structures as described by the settings object given a rotation between the two objects. More...
|
| |
| void | computeTranslationsFromPeak (ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ) |
| | This function computes the translation in Angstroms that corresponds to the translation function peak. More...
|
| |
| void | computeBeforeAfterZeroCounts (proshade_unsign *addXPre, proshade_unsign *addYPre, proshade_unsign *addZPre, proshade_unsign *addXPost, proshade_unsign *addYPost, proshade_unsign *addZPost, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_unsign xDimIndices, proshade_unsign yDimIndices, proshade_unsign zDimIndices) |
| | This function finds the number of zeroes to be added after and before the structure along each dimension. More...
|
| |
| void | paddMapWithZeroes (proshade_double *oldMap, proshade_double *&newMap, proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim, proshade_unsign xDimIndices, proshade_unsign yDimIndices, proshade_unsign zDimIndices, proshade_unsign addXPre, proshade_unsign addYPre, proshade_unsign addZPre) |
| | This function adds zeroes before and after the central map and copies the central map values into a new map. More...
|
| |
| void | allocateTranslationFunctionMemory (fftw_complex *&tmpIn1, fftw_complex *&tmpOut1, fftw_complex *&tmpIn2, fftw_complex *&tmpOut2, fftw_complex *&resIn, fftw_complex *&resOut, fftw_plan &forwardFourierObj1, fftw_plan &forwardFourierObj2, fftw_plan &inverseFourierCombo, proshade_unsign xD, proshade_unsign yD, proshade_unsign zD) |
| | This function allocates the memory for the Fourier transforms required for translation function computation. More...
|
| |
| void | freeTranslationFunctionMemory (fftw_complex *&tmpIn1, fftw_complex *&tmpOut1, fftw_complex *&tmpIn2, fftw_complex *&tmpOut2, fftw_complex *&resOut, fftw_plan &forwardFourierObj1, fftw_plan &forwardFourierObj2, fftw_plan &inverseFourierCombo) |
| | This function releases the memory for the Fourier transforms required for translation function computation. More...
|
| |
| void | computeAngularThreshold (std::vector< proshade_double > *lonCO, std::vector< proshade_double > *latCO, proshade_unsign angRes) |
| | This function computes the angular thresholds for longitude and lattitude angles. More...
|
| |
| void | initialiseInverseSHComputation (proshade_unsign shBand, double *&sigR, double *&sigI, double *&rcoeffs, double *&icoeffs, double *&weights, double *&workspace, fftw_plan &idctPlan, fftw_plan &ifftPlan) |
| | This function initialises internal variables for inverse Spherical Harmonics computation. More...
|
| |
This namespace contains all the functions required for map overlays.
The map overlay task does have some specific functions that are betters groupped together for higher clarity of the code and also simple modifications later. This is where these live.
◆ allocateTranslationFunctionMemory()
| void ProSHADE_internal_overlay::allocateTranslationFunctionMemory |
( |
fftw_complex *& |
tmpIn1, |
|
|
fftw_complex *& |
tmpOut1, |
|
|
fftw_complex *& |
tmpIn2, |
|
|
fftw_complex *& |
tmpOut2, |
|
|
fftw_complex *& |
resIn, |
|
|
fftw_complex *& |
resOut, |
|
|
fftw_plan & |
forwardFourierObj1, |
|
|
fftw_plan & |
forwardFourierObj2, |
|
|
fftw_plan & |
inverseFourierCombo, |
|
|
proshade_unsign |
xD, |
|
|
proshade_unsign |
yD, |
|
|
proshade_unsign |
zD |
|
) |
| |
This function allocates the memory for the Fourier transforms required for translation function computation.
- Parameters
-
| [in] | tmpIn1 | Array to hold the static structure Fourier inputs. |
| [in] | tmpOut1 | Array to hold the static structure Fourier outputs. |
| [in] | tmpIn2 | Array to hold the moving structure Fourier inputs. |
| [in] | tmpOut2 | Array to hold the moving structure Fourier outputs. |
| [in] | resIn | Array to hold the final translation function values. |
| [in] | resOut | Array to hold the combined Fourier coefficients of both structures. |
| [in] | forwardFourierObj1 | FFTW plan for the forward Fourier of the static structure. |
| [in] | forwardFourierObj2 | FFTW plan for the forward Fourier of the moving structure. |
| [in] | inverseFourierCombo | FFTW plan for the backward Fourier of the combined Fourier factors of both structures. |
| [in] | xD | The dimension of the X axis of the structures (assumes both structures have the same sizes and sampling). |
| [in] | yD | The dimension of the Y axis of the structures (assumes both structures have the same sizes and sampling). |
| [in] | zD | The dimension of the Z axis of the structures (assumes both structures have the same sizes and sampling). |
Definition at line 367 of file ProSHADE_overlay.cpp.
370 tmpIn1 =
reinterpret_cast< fftw_complex*
> ( fftw_malloc (
sizeof ( fftw_complex ) * xD * yD * zD ) );
371 tmpOut1 =
reinterpret_cast< fftw_complex*
> ( fftw_malloc (
sizeof ( fftw_complex ) * xD * yD * zD ) );
372 tmpIn2 =
reinterpret_cast< fftw_complex*
> ( fftw_malloc (
sizeof ( fftw_complex ) * xD * yD * zD ) );
373 tmpOut2 =
reinterpret_cast< fftw_complex*
> ( fftw_malloc (
sizeof ( fftw_complex ) * xD * yD * zD ) );
374 resIn =
reinterpret_cast< fftw_complex*
> ( fftw_malloc (
sizeof ( fftw_complex ) * xD * yD * zD ) );
375 resOut =
reinterpret_cast< fftw_complex*
> ( fftw_malloc (
sizeof ( fftw_complex ) * xD * yD * zD ) );
386 forwardFourierObj1 = fftw_plan_dft_3d (
static_cast< int > ( xD ),
static_cast< int > ( yD ),
static_cast< int > ( zD ), tmpIn1, tmpOut1, FFTW_FORWARD , FFTW_ESTIMATE );
387 forwardFourierObj2 = fftw_plan_dft_3d (
static_cast< int > ( xD ),
static_cast< int > ( yD ),
static_cast< int > ( zD ), tmpIn2, tmpOut2, FFTW_FORWARD , FFTW_ESTIMATE );
388 inverseFourierCombo = fftw_plan_dft_3d (
static_cast< int > ( xD ),
static_cast< int > ( yD ),
static_cast< int > ( zD ), resOut, resIn , FFTW_BACKWARD, FFTW_ESTIMATE );
◆ computeAngularThreshold()
| void ProSHADE_internal_overlay::computeAngularThreshold |
( |
std::vector< proshade_double > * |
lonCO, |
|
|
std::vector< proshade_double > * |
latCO, |
|
|
proshade_unsign |
angRes |
|
) |
| |
This function computes the angular thresholds for longitude and lattitude angles.
- Parameters
-
| [in] | lonCO | Pointer to vector where longitude thresholds are to be stored at. |
| [in] | latCO | Pointer to vector where lattitude thresholds are to be stored at. |
| [in] | angRes | The angular resolution of the computation. |
Definition at line 1287 of file ProSHADE_overlay.cpp.
1290 for ( proshade_unsign iter = 0; iter <= angRes; iter++ )
1292 ProSHADE_internal_misc::addToDoubleVector ( lonCO,
static_cast<proshade_double
> ( iter ) * ( (
static_cast<proshade_double
> ( M_PI ) * 2.0 ) /
static_cast<proshade_double
> ( angRes ) ) - (
static_cast<double> ( M_PI ) ) );
1296 for ( proshade_unsign iter = 0; iter <= angRes; iter++ )
1298 ProSHADE_internal_misc::addToDoubleVector ( latCO,
static_cast<proshade_double
> ( iter ) * (
static_cast<proshade_double
> ( M_PI ) /
static_cast<proshade_double
> ( angRes ) ) - (
static_cast<proshade_double
> ( M_PI ) / 2.0 ) );
◆ computeBeforeAfterZeroCounts()
| void ProSHADE_internal_overlay::computeBeforeAfterZeroCounts |
( |
proshade_unsign * |
addXPre, |
|
|
proshade_unsign * |
addYPre, |
|
|
proshade_unsign * |
addZPre, |
|
|
proshade_unsign * |
addXPost, |
|
|
proshade_unsign * |
addYPost, |
|
|
proshade_unsign * |
addZPost, |
|
|
proshade_unsign |
xDim, |
|
|
proshade_unsign |
yDim, |
|
|
proshade_unsign |
zDim, |
|
|
proshade_unsign |
xDimIndices, |
|
|
proshade_unsign |
yDimIndices, |
|
|
proshade_unsign |
zDimIndices |
|
) |
| |
This function finds the number of zeroes to be added after and before the structure along each dimension.
- Parameters
-
| [in] | addXPre | Variable pointer where the number of zeroes to be added before the map along X axis is saved. |
| [in] | addYPre | Variable pointer where the number of zeroes to be added before the map along Y axis is saved. |
| [in] | addZPre | Variable pointer where the number of zeroes to be added before the map along Z axis is saved. |
| [in] | addXPost | Variable pointer where the number of zeroes to be added after the map along X axis is saved. |
| [in] | addYPost | Variable pointer where the number of zeroes to be added after the map along Y axis is saved. |
| [in] | addZPost | Variable pointer where the number of zeroes to be added after the map along Z axis is saved. |
| [in] | xDim | The X dimension size in indices of the new map. |
| [in] | yDim | The Y dimension size in indices of the new map. |
| [in] | zDim | The Z dimension size in indices of the new map. |
| [in] | xDimIndices | The X dimension size in indices of the old map. |
| [in] | yDimIndices | The Y dimension size in indices of the old map. |
| [in] | zDimIndices | The Z dimension size in indices of the old map. |
Definition at line 492 of file ProSHADE_overlay.cpp.
495 *addXPre = ( xDim - xDimIndices ) / 2;
496 *addYPre = ( yDim - yDimIndices ) / 2;
497 *addZPre = ( zDim - zDimIndices ) / 2;
498 *addXPost = *addXPre;
if ( ( xDim - xDimIndices ) % 2 == 1 ) { *addXPost += 1; }
499 *addYPost = *addYPre;
if ( ( yDim - yDimIndices ) % 2 == 1 ) { *addYPost += 1; }
500 *addZPost = *addZPre;
if ( ( zDim - zDimIndices ) % 2 == 1 ) { *addZPost += 1; }
◆ computeTranslationsFromPeak()
This function computes the translation in Angstroms that corresponds to the translation function peak.
- Parameters
-
| [in] | staticStructure | A pointer to the data class object of the static structure. |
| [in] | movingStructure | A pointer to the data class object of the moving structure. |
| [in] | trsX | The variable to which the best X-axis position value will be saved to. |
| [in] | trsY | The variable to which the best Y-axis position value will be saved to. |
| [in] | trsZ | The variable to which the best Z-axis position value will be saved to. |
Definition at line 220 of file ProSHADE_overlay.cpp.
223 proshade_single xSamplRate = (
static_cast< proshade_single
> ( staticStructure->
getXDimSize() ) /
static_cast< proshade_single
> ( staticStructure->
getXDim() ) );
224 proshade_single ySamplRate = (
static_cast< proshade_single
> ( staticStructure->
getYDimSize() ) /
static_cast< proshade_single
> ( staticStructure->
getYDim() ) );
225 proshade_single zSamplRate = (
static_cast< proshade_single
> ( staticStructure->
getZDimSize() ) /
static_cast< proshade_single
> ( staticStructure->
getZDim() ) );
227 proshade_single xCor =
static_cast< proshade_single
> (
static_cast< proshade_single
> ( staticStructure->
xFrom - movingStructure->
xFrom ) * xSamplRate );
228 proshade_single yCor =
static_cast< proshade_single
> (
static_cast< proshade_single
> ( staticStructure->
yFrom - movingStructure->
yFrom ) * ySamplRate );
229 proshade_single zCor =
static_cast< proshade_single
> (
static_cast< proshade_single
> ( staticStructure->
zFrom - movingStructure->
zFrom ) * zSamplRate );
231 proshade_single xMov =
static_cast< proshade_single
> ( 0.0 -
static_cast< proshade_double
> ( xCor ) - ( *trsX *
static_cast< proshade_double
> ( xSamplRate ) ) );
232 proshade_single yMov =
static_cast< proshade_single
> ( 0.0 -
static_cast< proshade_double
> ( yCor ) - ( *trsY *
static_cast< proshade_double
> ( ySamplRate ) ) );
233 proshade_single zMov =
static_cast< proshade_single
> ( 0.0 -
static_cast< proshade_double
> ( zCor ) - ( *trsZ *
static_cast< proshade_double
> ( zSamplRate ) ) );
236 *trsX =
static_cast< proshade_double
> ( -xMov );
237 *trsY =
static_cast< proshade_double
> ( -yMov );
238 *trsZ =
static_cast< proshade_double
> ( -zMov );
◆ freeTranslationFunctionMemory()
| void ProSHADE_internal_overlay::freeTranslationFunctionMemory |
( |
fftw_complex *& |
tmpIn1, |
|
|
fftw_complex *& |
tmpOut1, |
|
|
fftw_complex *& |
tmpIn2, |
|
|
fftw_complex *& |
tmpOut2, |
|
|
fftw_complex *& |
resOut, |
|
|
fftw_plan & |
forwardFourierObj1, |
|
|
fftw_plan & |
forwardFourierObj2, |
|
|
fftw_plan & |
inverseFourierCombo |
|
) |
| |
This function releases the memory for the Fourier transforms required for translation function computation.
- Parameters
-
| [in] | tmpIn1 | Array to hold the static structure Fourier inputs. |
| [in] | tmpOut1 | Array to hold the static structure Fourier outputs. |
| [in] | tmpIn2 | Array to hold the moving structure Fourier inputs. |
| [in] | tmpOut2 | Array to hold the moving structure Fourier outputs. |
| [in] | resOut | Array to hold the combined Fourier coefficients of both structures. |
| [in] | forwardFourierObj1 | FFTW plan for the forward Fourier of the static structure. |
| [in] | forwardFourierObj2 | FFTW plan for the forward Fourier of the moving structure. |
| [in] | inverseFourierCombo | FFTW plan for the backward Fourier of the combined Fourier factors of both structures. |
Definition at line 406 of file ProSHADE_overlay.cpp.
409 fftw_destroy_plan ( forwardFourierObj1 );
410 fftw_destroy_plan ( forwardFourierObj2 );
411 fftw_destroy_plan ( inverseFourierCombo );
412 fftw_free ( tmpIn1 );
413 fftw_free ( tmpIn2 );
414 fftw_free ( tmpOut1 );
415 fftw_free ( tmpOut2 );
416 fftw_free ( resOut );
◆ getOptimalRotation()
This function finds the optimal rotation between two structures as described by the settings object.
This function takes the settings and two structure classes. It then reads in and processes both structures so that the globally optimal rotation overlay Euler angles are detected. However, this is only the case for rotation along the centre of the map of the second structure; therefore, either use Patterson data (usePhase = false), or be aware that better rotation may exist for different centre of rotation.
- Parameters
-
| [in] | settings | A pointer to settings class containing all the information required for map symmetry detection. |
| [in] | obj1 | A pointer to the data class object of the other ( static ) structure. |
| [in] | obj2 | A pointer to the data class object of the first ( moving ) structure. |
| [in] | eulA | The variable to which the best Euler alpha angle value will be saved to. |
| [in] | eulB | The variable to which the best Euler beta angle value will be saved to. |
| [in] | eulG | The variable to which the best Euler gamma angle value will be saved to. |
Definition at line 74 of file ProSHADE_overlay.cpp.
98 eulA, eulB, eulG, settings );
◆ getOptimalTranslation()
This function finds the optimal translation between two structures as described by the settings object given a rotation between the two objects.
This function starts by loading and processing the structures according to the settings object (keeping phase is assumed, but callers responsibility). It then applies the required rotation to the second (moing) strucutre and then it follows with zero padding to make sure the structures have the same dimensions (again, it assumes map re-sampling was done, but setting to is callers responsibility). It then computes the translation function, finds the highest peak and returns the positions as well as height of this peak.
- Parameters
-
| [in] | settings | A pointer to settings class containing all the information required for map overlay computation. |
| [in] | staticStructure | A pointer to the data class object of the other ( static ) structure. |
| [in] | movingStructure | A pointer to the data class object of the first ( moving ) structure. |
| [in] | trsX | The variable to which the best X-axis position value will be saved to. |
| [in] | trsY | The variable to which the best Y-axis position value will be saved to. |
| [in] | trsZ | The variable to which the best Z-axis position value will be saved to. |
| [in] | eulA | The Euler alpha angle value, by which the moving structure is to be rotated by. |
| [in] | eulB | The Euler beta angle value, by which the moving structure is to be rotated by. |
| [in] | eulG | The Euler gamma angle value, by which the moving structure is to be rotated by. |
Definition at line 123 of file ProSHADE_overlay.cpp.
142 std::max ( staticStructure->
getYDim(), movingStructure->
getYDim() ),
143 std::max ( staticStructure->
getZDim(), movingStructure->
getZDim() ) );
145 std::max ( staticStructure->
getYDim(), movingStructure->
getYDim() ),
146 std::max ( staticStructure->
getZDim(), movingStructure->
getZDim() ) );
155 proshade_double mapPeak = 0.0;
156 proshade_unsign xDimS = staticStructure->
getXDim();
157 proshade_unsign yDimS = staticStructure->
getYDim();
158 proshade_unsign zDimS = staticStructure->
getZDim();
162 if ( *trsX > (
static_cast< proshade_double
> ( xDimS ) / 2.0 ) ) { *trsX = *trsX -
static_cast< proshade_double
> ( xDimS ); }
163 if ( *trsY > (
static_cast< proshade_double
> ( yDimS ) / 2.0 ) ) { *trsY = *trsY -
static_cast< proshade_double
> ( yDimS ); }
164 if ( *trsZ > (
static_cast< proshade_double
> ( zDimS ) / 2.0 ) ) { *trsZ = *trsZ -
static_cast< proshade_double
> ( zDimS ); }
173 proshade_single xMov =
static_cast< proshade_single
> ( -(*trsX) );
174 proshade_single yMov =
static_cast< proshade_single
> ( -(*trsY) );
175 proshade_single zMov =
static_cast< proshade_single
> ( -(*trsZ) );
178 std::stringstream hlpSS;
179 hlpSS <<
"Optimal map translation distances are " << *trsX <<
" ; " << *trsY <<
" ; " << *trsZ <<
" Angstroms with peak height " << mapPeak / (
static_cast< proshade_double
> ( xDimS ) *
static_cast< proshade_double
> ( yDimS ) *
static_cast< proshade_double
> ( zDimS ) );
195 static_cast< proshade_signed
> ( movingStructure->
getXDim() ),
static_cast< proshade_signed
> ( movingStructure->
getYDim() ),
196 static_cast< proshade_signed
> ( movingStructure->
getZDim() ) );
◆ initialiseInverseSHComputation()
| void ProSHADE_internal_overlay::initialiseInverseSHComputation |
( |
proshade_unsign |
shBand, |
|
|
double *& |
sigR, |
|
|
double *& |
sigI, |
|
|
double *& |
rcoeffs, |
|
|
double *& |
icoeffs, |
|
|
double *& |
weights, |
|
|
double *& |
workspace, |
|
|
fftw_plan & |
idctPlan, |
|
|
fftw_plan & |
ifftPlan |
|
) |
| |
This function initialises internal variables for inverse Spherical Harmonics computation.
- Parameters
-
| [in] | shBand | The bandwidth for this particular shell. |
| [in] | sigR | Pointer to be initialised for the real signal values. |
| [in] | sigI | Pointer to be initialised for the imaginary signal values. |
| [in] | rcoeffs | Pointer to be initialised for the real coefficient values. |
| [in] | icoeffs | Pointer to be initialised for the imaginary coefficient values. |
| [in] | weights | Pointer to be initialised for the transform weight values. |
| [in] | workspace | Pointer to be initialised for the computation screatch space. |
| [in] | idctPlan | Pointer reference to the cosine/sine transform plan to be created. |
| [in] | ifftPlan | Pointer reference to the discrete 3D Fourier transform plan to be created. |
Definition at line 1171 of file ProSHADE_overlay.cpp.
1174 proshade_unsign oneDim = shBand * 2;
1177 sigR =
new double [(oneDim * oneDim * 4)];
1178 sigI = sigR + (oneDim * oneDim * 2);
1179 rcoeffs =
new double [(oneDim * oneDim * 2)];
1180 icoeffs = rcoeffs + (oneDim * oneDim);
1181 weights =
new double [4 * shBand];
1182 workspace =
new double [( 20 * shBand * shBand ) + ( 24 * shBand )];
1191 idctPlan = fftw_plan_r2r_1d (
static_cast< int > ( oneDim ), weights, workspace, FFTW_REDFT01, FFTW_ESTIMATE );
1194 int rank, howmany_rank;
1195 fftw_iodim dims[1], howmany_dims[1];
1198 dims[0].n = 2 *
static_cast< int > ( shBand );
1199 dims[0].is = 2 *
static_cast< int > ( shBand );
1202 howmany_dims[0].n = 2 *
static_cast< int > ( shBand );
1203 howmany_dims[0].is = 1;
1204 howmany_dims[0].os = 2 *
static_cast< int > ( shBand );
1207 ifftPlan = fftw_plan_guru_split_dft ( rank, dims, howmany_rank, howmany_dims, sigR, sigI, rcoeffs, icoeffs, FFTW_ESTIMATE );
◆ paddMapWithZeroes()
| void ProSHADE_internal_overlay::paddMapWithZeroes |
( |
proshade_double * |
oldMap, |
|
|
proshade_double *& |
newMap, |
|
|
proshade_unsign |
xDim, |
|
|
proshade_unsign |
yDim, |
|
|
proshade_unsign |
zDim, |
|
|
proshade_unsign |
xDimIndices, |
|
|
proshade_unsign |
yDimIndices, |
|
|
proshade_unsign |
zDimIndices, |
|
|
proshade_unsign |
addXPre, |
|
|
proshade_unsign |
addYPre, |
|
|
proshade_unsign |
addZPre |
|
) |
| |
This function adds zeroes before and after the central map and copies the central map values into a new map.
- Parameters
-
| [in] | oldMap | The original, unpadded map. |
| [in] | newMap | The array to which the new map will be saved into. |
| [in] | xDim | The X dimension size of the new map. |
| [in] | yDim | The Y dimension size of the new map. |
| [in] | zDim | The Z dimension size of the new map. |
| [in] | xDimIndices | The X dimension size in indices of the old map. |
| [in] | yDimIndices | The Y dimension size in indices of the old map. |
| [in] | zDimIndices | The Z dimension size in indices of the old map. |
| [in] | addXPre | How many zeroes are to be added before the central map along the X axis? |
| [in] | addYPre | How many zeroes are to be added before the central map along the Y axis? |
| [in] | addZPre | How many zeroes are to be added before the central map along the Z axis? |
Definition at line 521 of file ProSHADE_overlay.cpp.
524 proshade_unsign newMapIndex = 0;
525 proshade_unsign oldMapIndex = 0;
528 for ( proshade_unsign xIt = 0; xIt < xDim; xIt++ )
530 for ( proshade_unsign yIt = 0; yIt < yDim; yIt++ )
532 for ( proshade_unsign zIt = 0; zIt < zDim; zIt++ )
535 newMapIndex = zIt + zDim * ( yIt + yDim * xIt );
538 if ( xIt < addXPre ) { newMap[newMapIndex] = 0.0;
continue; }
539 if ( yIt < addYPre ) { newMap[newMapIndex] = 0.0;
continue; }
540 if ( zIt < addZPre ) { newMap[newMapIndex] = 0.0;
continue; }
543 if ( xIt >= ( addXPre + xDimIndices ) ) { newMap[newMapIndex] = 0.0;
continue; }
544 if ( yIt >= ( addYPre + yDimIndices ) ) { newMap[newMapIndex] = 0.0;
continue; }
545 if ( zIt >= ( addZPre + zDimIndices ) ) { newMap[newMapIndex] = 0.0;
continue; }
548 oldMapIndex = (zIt-addZPre) + zDimIndices * ( (yIt-addYPre) + yDimIndices * (xIt-addXPre) );
549 newMap[newMapIndex] = oldMap[oldMapIndex];
void computeSphericalHarmonics(ProSHADE_settings *settings)
This function computes the spherical harmonics decomposition for the whole structure.
proshade_double mapMovFromsChangeX
When the map is translated, the xFrom and xTo values are changed. This variable holds how much they h...
void getOverlayRotationFunction(ProSHADE_settings *settings, ProSHADE_internal_data::ProSHADE_data *obj2)
This function computes the overlay rotation function (i.e. the correlation function in SO(3) space).
proshade_signed zFrom
This is the starting index along the z axis.
proshade_unsign xDimIndices
This is the size of the map cell x dimension in indices.
proshade_single getZDimSize(void)
This function allows access to the map size in angstroms along the Z axis.
void determineAllSHValues(proshade_unsign xDim, proshade_unsign yDim, proshade_single xDimAngs, proshade_single yDimAngs, proshade_single zDimAngs)
This function determines all the required values for spherical harmonics computation.
proshade_signed * getXAxisOrigin(void)
This function allows access to the map X axis origin value.
proshade_single getXDimSize(void)
This function allows access to the map size in angstroms along the X axis.
proshade_double *& getInternalMap(void)
This function allows access to the first map array value address.
proshade_single getYDimSize(void)
This function allows access to the map size in angstroms along the Y axis.
proshade_single zDimSize
This is the size of the map cell z dimension in Angstroms.
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.
proshade_signed * getXFromPtr(void)
This function allows access to the map start along the X axis.
proshade_signed verbose
Should the software report on the progress, or just be quiet? Value between -1 (nothing) and 4 (loud)
proshade_unsign yDimIndices
This is the size of the map cell y dimension in indices.
proshade_signed messageShift
This value allows shifting the messages to create more readable log for sub-processes.
void addToDoubleVector(std::vector< proshade_double > *vecToAddTo, proshade_double elementToAdd)
Adds the element to the vector.
proshade_unsign getXDim(void)
This function allows access to the map size in indices along the X axis.
proshade_double originalPdbTransX
The optimal translation vector as it relates to the original PDB positions (and not the ProSHADE inte...
void readInStructure(std::string fName, proshade_unsign inputO, ProSHADE_settings *settings, proshade_double *maskArr=nullptr, proshade_unsign maskXDim=0, proshade_unsign maskYDim=0, proshade_unsign maskZDim=0, proshade_double *weightsArr=nullptr, proshade_unsign weigXDim=0, proshade_unsign weigYDim=0, proshade_unsign weigZDim=0)
This function initialises the basic ProSHADE_data variables and reads in a single structure.
void mapToSpheres(ProSHADE_settings *settings)
This function converts the internal map onto a set of concentric spheres.
proshade_double originalPdbTransY
The optimal translation vector as it relates to the original PDB positions (and not the ProSHADE inte...
proshade_single xDimSize
This is the size of the map cell x dimension in Angstroms.
proshade_double mapMovFromsChangeY
When the map is translated, the yFrom and yTo values are changed. This variable holds how much they h...
void getBestPeakEulerAngsNaive(proshade_complex *map, proshade_unsign dim, proshade_double *eulA, proshade_double *eulB, proshade_double *eulG, ProSHADE_settings *settings)
This function finds the highest peaks optimised Euler angles using the "naive" approach.
void computeTranslationMap(ProSHADE_internal_data::ProSHADE_data *obj1)
This function does the computation of the translation map and saves results internally.
proshade_signed xFrom
This is the starting index along the x axis.
proshade_unsign getYDim(void)
This function allows access to the map size in indices along the Y axis.
void findHighestValueInMap(fftw_complex *resIn, proshade_unsign xD, proshade_unsign yD, proshade_unsign zD, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ, proshade_double *mapPeak)
This function simply finds the highest value in fftw_complex map and returns its position and value.
proshade_double originalPdbTransZ
The optimal translation vector as it relates to the original PDB positions (and not the ProSHADE inte...
void processInternalMap(ProSHADE_settings *settings)
This function simply clusters several other functions which should be called together.
proshade_complex * getTranslationFnPointer(void)
This function allows access to the translation function through a pointer.
std::vector< proshade_double > rotateMapRealSpaceInPlace(proshade_double eulA, proshade_double eulB, proshade_double eulG)
This function rotates a map based on the given Euler angles in place.
proshade_unsign getZDim(void)
This function allows access to the map size in indices along the Z axis.
proshade_double mapMovFromsChangeZ
When the map is translated, the zFrom and zTo values are changed. This variable holds how much they h...
proshade_signed * getYToPtr(void)
This function allows access to the map last position along the Y axis.
proshade_signed * getZAxisOrigin(void)
This function allows access to the map Z axis origin value.
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.
proshade_signed * getXToPtr(void)
This function allows access to the map last position along the X axis.
proshade_signed * getYFromPtr(void)
This function allows access to the map start along the Y axis.
proshade_complex * getInvSO3Coeffs(void)
This function allows access to the inverse SO(3) coefficients array.
proshade_unsign getEMatDim(void)
This function allows access to the maximum band for the E matrix.
proshade_signed yFrom
This is the starting index along the y axis.
std::vector< std::string > inputFiles
This vector contains the filenames of all input structure files.
proshade_signed * getYAxisOrigin(void)
This function allows access to the map Y axis origin value.
void computeTranslationsFromPeak(ProSHADE_internal_data::ProSHADE_data *staticStructure, ProSHADE_internal_data::ProSHADE_data *movingStructure, proshade_double *trsX, proshade_double *trsY, proshade_double *trsZ)
This function computes the translation in Angstroms that corresponds to the translation function peak...
void printProgressMessage(proshade_signed verbose, proshade_signed messageLevel, std::string message, proshade_signed messageShift=0)
General stdout message printing.
proshade_single yDimSize
This is the size of the map cell y dimension in Angstroms.
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.
void zeroPaddToDims(proshade_unsign xDim, proshade_unsign yDim, proshade_unsign zDim)
This function changes the size of a structure to fit the supplied new limits.
proshade_signed * getZFromPtr(void)
This function allows access to the map start along the Z axis.
proshade_signed * getZToPtr(void)
This function allows access to the map last position along the Z axis.