param.cpp 27.2 KB
Newer Older
Pilar Cossio's avatar
License  
Pilar Cossio committed
1
/* ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 3 4
   < BioEM software for Bayesian inference of Electron Microscopy images>
   Copyright (C) 2014 Pilar Cossio, David Rohr and Gerhard Hummer.
   Max Planck Institute of Biophysics, Frankfurt, Germany.
5

6
   See license statement for terms of distribution.
Pilar Cossio's avatar
License  
Pilar Cossio committed
7 8 9

   ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++*/

qon's avatar
qon committed
10 11 12 13 14 15 16 17
#include <stdio.h>
#include <stdlib.h>
#include <iostream>
#include <fstream>
#include <cstring>
#include <math.h>
#include <fftw3.h>

18 19 20 21
#ifdef WITH_OPENMP
#include <omp.h>
#endif

qon's avatar
qon committed
22 23 24 25 26 27 28
#include "param.h"
#include "map.h"

using namespace std;

bioem_param::bioem_param()
{
29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51
  //Number of Pixels
  param_device.NumberPixels = 0;
  param_device.NumberFFTPixels1D = 0;
  // Euler angle grid spacing
  angleGridPointsAlpha = 0;
  angleGridPointsBeta = 0;
  //Envelop function paramters
  numberGridPointsEnvelop = 0;
  //Contrast transfer function paramters
  numberGridPointsCTF_amp = 0;
  numberGridPointsCTF_phase = 0;

  // ****center displacement paramters Equal in both directions***
  param_device.maxDisplaceCenter = 0;
  numberGridPointsDisplaceCenter = 0;

  fft_plans_created = 0;

  refCTF = NULL;
  CtfParam = NULL;
  angles = NULL;

  printModel = false;
qon's avatar
qon committed
52 53
}

54 55
int bioem_param::readParameters(const char* fileinput,const char* fileangles)
{	// **************************************************************************************
David Rohr's avatar
David Rohr committed
56 57
	// ***************************** Reading Input Parameters ******************************
	// **************************************************************************************
58

David Rohr's avatar
David Rohr committed
59
	// Control for Parameters
60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78
  bool yesPixSi = false;
  bool yesNumPix = false;
  bool yesGPal = false;
  bool yesGPbe = false;
  bool yesGPEnv = false;
  bool yesGPamp = false;
  bool yesGPpha = false;
  bool yesSTEnv = false;
  bool yesSTamp = false;
  bool yesSTpha = false;
  bool yesGSPamp = false ;
  bool yesGSPEnv = false ;
  bool yesGSPpha = false ;
  bool yesMDC = false ;
  bool yesGCen = false ;

  //Default not flipped micrographs
  param_device.flipped=false;
  param_device.debugterm=false;
79
  param_device.writeCC=false;
Pilar Cossio's avatar
Pilar Cossio committed
80
  param_device.CCwithBayes=true;
81
  writeCTF=false;
82
  elecwavel=220;
83
  ignoreCCoff=false;
84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145
  //Assigning name of angles
    
  //Storing angle file name if existing. 
  inanglef=std::string(fileangles);
  NotUn_angles=0;

  ifstream input(fileinput);
  if (!input.good())
    {
      cout << "Failed to open file: " << fileinput << "\n";
      exit(0);
    }

  char line[512] = {0};
  char saveline[512];

  cout << "\n +++++++++++++++++++++++++++++++++++++++++ \n";
  cout << "\n   READING BioEM PARAMETERS             \n\n";
  cout << " +++++++++++++++++++++++++++++++++++++++++ \n";
	
  while (!input.eof())
    {
      input.getline(line, 512);
      strcpy(saveline, line);
      char *token = strtok(line, " ");

      if (token == NULL || line[0] == '#' || strlen(token) == 0)
	{
	  // comment or blank line
	}
      else if (strcmp(token, "PIXEL_SIZE") == 0)
	{
	  token = strtok(NULL, " ");
	  pixelSize = atof(token);
	  if (pixelSize < 0 ) { cout << "*** Error: Negative pixelSize "; exit(1);}
	  cout << "Pixel Sixe " << pixelSize << "\n";
	  yesPixSi= true;
	}
      else if (strcmp(token, "NUMBER_PIXELS") == 0)
	{
	  token = strtok(NULL, " ");
	  param_device.NumberPixels = int(atoi(token));
	  if (param_device.NumberPixels < 0 ) { cout << "*** Error: Negative Number of Pixels "; exit(1);}
	  cout << "Number of Pixels " << param_device.NumberPixels << "\n";
	  yesNumPix= true ;
	}
      else if (strcmp(token, "GRIDPOINTS_ALPHA") == 0)
	{
	  token = strtok(NULL, " ");
	  angleGridPointsAlpha = int(atoi(token));
	  if (angleGridPointsAlpha < 0 ) { cout << "*** Error: Negative GRIDPOINTS_ALPHA "; exit(1);}
	  cout << "Grid points alpha " << angleGridPointsAlpha << "\n";
	  yesGPal= true;
	}
      else if (strcmp(token, "GRIDPOINTS_BETA") == 0)
	{
	  token = strtok(NULL, " ");
	  angleGridPointsBeta = int(atoi(token));
	  if (angleGridPointsBeta < 0 ) { cout << "*** Error: Negative GRIDPOINTS_BETA "; exit(1);}
	  cout << "Grid points in Cosine ( beta ) " << angleGridPointsBeta << "\n";
	  yesGPbe= true;
	}
146

147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238
      else if (strcmp(token, "GRIDPOINTS_ENVELOPE") == 0)
	{
	  token = strtok(NULL, " ");
	  numberGridPointsEnvelop = int(atoi(token));
	  if (numberGridPointsDisplaceCenter < 0 ) { cout << "*** Error: Negative GRIDPOINTS_ENVELOPE "; exit(1);}
	  cout << "Grid points envelope " << numberGridPointsEnvelop << "\n";
	  yesGPEnv = true;
	}
      else if (strcmp(token, "START_ENVELOPE") == 0)
	{
	  token = strtok(NULL, " ");
	  startGridEnvelop = atof(token);
	  if (startGridEnvelop < 0 ) { cout << "*** Error: Negative START_ENVELOPE "; exit(1);}
	  cout << "Start Envelope " << startGridEnvelop << "\n";
	  yesSTEnv = true ;
	}
      else if (strcmp(token, "GRIDSPACE_ENVELOPE") == 0)
	{
	  token = strtok(NULL, " ");
	  gridEnvelop = atof(token);
	  if (gridEnvelop < 0 ) { cout << "*** Error: Negative GRIDSPACE_ENVELOPE "; exit(1);}
	  cout << "Grid spacing Envelope " << gridEnvelop << "\n";
	  yesGSPEnv = true ;
	}
      else if (strcmp(token,"GRIDPOINTS_PSF_PHASE")==0)
	{
	  token = strtok(NULL," ");
	  numberGridPointsCTF_phase=int(atoi(token));
	  cout << "Grid points PSF " << numberGridPointsCTF_phase << "\n";
	  yesGPpha = true;
	}
      else if (strcmp(token,"START_PSF_PHASE")==0)
	{
	  token = strtok(NULL," ");
	  startGridCTF_phase=atof(token);
	  cout << "Start PSF " << startGridCTF_phase << "\n";
	  yesSTpha = true ;
	}
      else if (strcmp(token,"GRIDSPACE_PSF_PHASE")==0)
	{
	  token = strtok(NULL," ");
	  gridCTF_phase=atof(token);
	  cout << "Grid Space PSF " << gridCTF_phase << "\n";
	  yesGSPpha = true ;
	}
      else if (strcmp(token,"GRIDPOINTS_PSF_AMP")==0)
	{
	  token = strtok(NULL," ");
	  numberGridPointsCTF_amp=int(atoi(token));
	  cout << "Grid points PSF " << numberGridPointsCTF_amp << "\n";
	  yesGPamp = true ;
	}
      else if (strcmp(token,"START_PSF_AMP")==0)
	{
	  token = strtok(NULL," ");
	  startGridCTF_amp=atof(token);
	  cout << "Start PSF " << startGridCTF_amp << "\n";
	  yesSTamp = true ;
	}
      else if (strcmp(token,"GRIDSPACE_PSF_AMP")==0)
	{
	  token = strtok(NULL," ");
	  gridCTF_amp=atof(token);
	  cout << "Grid Space PSF " << gridCTF_amp << "\n";
	  yesGSPamp = true ;
	}
      else if (strcmp(token, "MAX_D_CENTER") == 0)
	{
	  token = strtok(NULL, " ");
	  param_device.maxDisplaceCenter = int(atoi(token));
	  if (param_device.maxDisplaceCenter < 0 ) { cout << "*** Error: Negative MAX_D_CENTER "; exit(1);}
	  cout << "Maximum displacement Center " <<  param_device.maxDisplaceCenter << "\n";
	  yesMDC = true;
	}
      else if (strcmp(token, "PIXEL_GRID_CENTER") == 0)
	{
	  token = strtok(NULL, " ");
	  param_device.GridSpaceCenter = int(atoi(token));
	  if (param_device.GridSpaceCenter < 0 ) { cout << "*** Error: Negative PIXEL_GRID_CENTER "; exit(1);}
	  cout << "Grid space displacement center " <<   param_device.GridSpaceCenter << "\n";
	  yesGCen = true;
	}
      else if (strcmp(token, "WRITE_PROB_ANGLES") == 0) //Key word if writing down each angle probabilities
	{
	  param_device.writeAngles = true;
	  cout << "Writing Probabilies of each angle \n";
	}
      else if (strcmp(token, "WRITE_CROSSCOR") == 0)//Key word if writing down full micrograph cross correlation
	{
	  param_device.writeCC = true;
	  cout << "Writing CrossCorrelations \n";
	}
Pilar Cossio's avatar
Pilar Cossio committed
239
      else if (strcmp(token, "CROSSCOR_GRID_SPACE") == 0)
240 241 242 243 244 245
	{
	  token = strtok(NULL, " ");
	  param_device.CCdisplace=int(atoi(token));
	  if (param_device.CCdisplace < 0 ) { cout << "*** Error: Negative CROSSCOR_DISPLACE "; exit(1);}
	  cout << "Writing Cross Correlation Displacement " <<  param_device.CCdisplace << "  \n";
	}
Pilar Cossio's avatar
Pilar Cossio committed
246
        else if (strcmp(token, "CROSSCOR_NOTBAYESIAN") == 0)
247
        {
Pilar Cossio's avatar
Pilar Cossio committed
248
          param_device.CCwithBayes=false;
249 250
          cout << "Using Bayesian Analysis to write Cross Correlation  \n";
        }
251 252 253 254 255
      else if (strcmp(token, "FLIPPED") == 0) //Key word if images are flipped for cross-correlation
	{
	  param_device.flipped = true;
	  cout << "Micrograph Flipped Intensities \n";
	}
256 257 258 259 260
        else if (strcmp(token, "IGNORE_CROSSCORR_OFFSET") == 0) //Key word if images are flipped for cross-correlation
        {
          ignoreCCoff = true;
          cout << "Ignoring Cross-Correlation offset \n";
        }
Pilar Cossio's avatar
Pilar Cossio committed
261
      else if (strcmp(token, "PROJECT_RADIUS") == 0) //If projecting CA with amino-acid radius
262 263 264 265 266 267 268 269 270 271 272 273 274 275 276
	{
	  doaaradius = true;
	  cout << "Projecting corresponding radius \n";
	}
      else if (strcmp(token, "DEBUG_INDI_PROB_TERM") == 0)//writing out each term of the probability
	{
	  param_device.debugterm = true;
	  cout << "Debugging Individual Probability Terms \n";
	}
      else if (strcmp(token, "NOT_UNIFORM_TOTAL_ANGS") == 0)//Number of Euler angle tripplets in non uniform Euler angle sampling
	{
	  token = strtok(NULL, " ");
	  NotUn_angles=int(atoi(token));
	  cout << "Not uniform total angel: " <<  NotUn_angles<< "\n";
	}
277 278 279 280
	 else if (strcmp(token, "WRITE_CTF_PARAM") == 0)//Number of Euler angle tripplets in non uniform Euler angle sampling
        {
         writeCTF=true;
	 token = strtok(NULL," ");
281 282 283 284 285 286 287 288

         cout << "Writing CTF parameters from PSF param that maximize the posterior \n";
//         cout << "tt " << token;
       	if(token!=NULL){	
		elecwavel=atof(token);
          	if(elecwavel < 150 ){
			cout << "Wrong electron wave length " << elecwavel << "\n";
			 cout << "Has to be in Angstrom (A)\n"; 
289
 		exit(1);} 
290 291 292 293
	        cout << "Electron wave length in (A) is: " << elecwavel << "\n";
 	}else{
		 cout << "Using default electron wave length: 220 (A)\n";
 		};
294
        }
295

296

297 298 299 300 301 302
    }
  input.close();

  // Checks for ALL INPUT
  if( not ( yesPixSi ) ){ cout << "**** INPUT MISSING: Please provide PIXEL_SIZE\n" ; exit (1);};
  if( not ( yesNumPix ) ){ cout << "**** INPUT MISSING: Please provide NUMBER_PIXELS \n" ; exit (1);};
Pilar Cossio's avatar
Pilar Cossio committed
303 304 305 306
  if(!notuniformangles){
  	if( not ( yesGPal) ) { cout << "**** INPUT MISSING: Please provide GRIDPOINTS_ALPHA \n" ; exit (1);};
	if( not ( yesGPbe )) { cout << "**** INPUT MISSING: Please provide GRIDPOINTS_BETA \n" ; exit (1);};
  }
307 308 309 310 311 312 313 314 315 316 317 318
  if( not (  yesGPEnv ) ) { cout << "**** INPUT MISSING: Please provide GRIDPOINTS_ENVELOPE \n" ; exit (1);};
  if( not (  yesGPamp ) ) { cout << "**** INPUT MISSING: Please provide GRIDPOINTS_PSF_AMP \n" ; exit (1);};
  if( not (  yesGPpha ) ) { cout << "**** INPUT MISSING: Please provide GRIDPOINTS_PSF_PHASE \n" ; exit (1);};
  if( not (  yesSTEnv  ) ) { cout << "**** INPUT MISSING: Please provide START_ENVELOPE \n" ; exit (1);};
  if( not (  yesSTamp  ) ) { cout << "**** INPUT MISSING: Please provide START_PSF_AMP \n" ; exit (1);};
  if( not (  yesSTpha  ) ) { cout << "**** INPUT MISSING: Please provide START_PSF_PHASE \n" ; exit (1);};
  if( not (  yesGSPamp  ) ) { cout << "**** INPUT MISSING: Please provide GRIDSPACE_PSF_AMP \n" ; exit (1);};
  if( not ( yesGSPEnv  ) ) { cout << "**** INPUT MISSING: Please provide GRIDSPACE_ENVELOPE \n" ; exit (1);};
  if( not ( yesGSPpha  ) ) { cout << "**** INPUT MISSING: Please provide GRIDSPACE_PSF_PHASE \n" ; exit (1);};
  if( not (  yesMDC  ) ) { cout << "**** INPUT MISSING: Please provide MAX_D_CENTER \n" ; exit (1);};
  if( not (  yesGCen  ) ) { cout << "**** INPUT MISSING: Please provide PIXEL_GRID_CENTER \n" ; exit (1);};
  if( param_device.writeCC && param_device.CCdisplace < 1 ){ cout << "**** INPUT MISSING: Please provide CROSSCOR_DISPLACE \n" ; exit (1);};
319 320 321 322
  if( !param_device.writeCC && param_device.CCwithBayes ){  cout << "**** INPUT MISSING: WRITE_CROSSCOR keyword \n"; exit(1);}
  if( param_device.writeCC) {if(!param_device.CCwithBayes ){  cout << "Remark:: Not Using Bayesian method to store Cross-Correlation.\n Only Printing out Maximum\n";}
       if(param_device.flipped){ cout << "Remark:: Micrographs are Flipped = Particles are white\n";} else {  cout << "Remark:: Micrographs are NOT Flipped = Particles are dark\n";}
	}
323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359

  //if only one grid point for PSF kernel:
  if( (myfloat_t) numberGridPointsCTF_amp == 1 ) gridCTF_amp = startGridCTF_amp;
  if( (myfloat_t) numberGridPointsCTF_phase == 1 ) gridCTF_phase = startGridCTF_phase;
  if( (myfloat_t) numberGridPointsEnvelop == 1 ) gridEnvelop = startGridEnvelop;

  //More checks with input parameters
  // Envelope should not have a standard deviation greater than Npix/2
  if(sqrt(1./( (myfloat_t) numberGridPointsDisplaceCenter  * gridEnvelop + startGridEnvelop))> float(param_device.NumberPixels)/2.0) {
    cout << "MAX Standard deviation of envelope is larger than Allowed KERNEL Length\n";
    exit(1);
  }
  // Envelop param should be positive
  if( startGridCTF_amp < 0 || startGridCTF_amp > 1){
    cout << "PSF Amplitud should be between 0 and 1\n";
    exit(1);
  }

  if( (myfloat_t) numberGridPointsCTF_amp * gridCTF_amp + startGridCTF_amp < 0 || (myfloat_t) (numberGridPointsCTF_amp - 1) * gridCTF_amp + startGridCTF_amp > 1){
    cout << "Value should be between 0 and 1\n" ;
    exit(1);
  }

  param_device.NumberFFTPixels1D = param_device.NumberPixels / 2 + 1;
  FFTMapSize = param_device.NumberPixels * param_device.NumberFFTPixels1D;

  //Does currently not work with custom alignment on GPU
  /*if (FFTMapSize % Alignment)
    {
    FFTMapSize += Alignment - FFTMapSize % Alignment;
    }
    cout << "Using MAP Size " << FFTMapSize << " (Alignment " << Alignment << ", Unaligned Size " << param_device.NumberPixels * param_device.NumberFFTPixels1D << ")\n";*/

  cout << " +++++++++++++++++++++++++++++++++++++++++ \n";

  return(0);
}
360

361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391
int bioem_param::forprintBest(const char* fileinput)
{
  // **************************************************************************************
  // **********Alternative parameter routine for only printing out a map ******************

  ifstream input(fileinput);
  withnoise=false;

  if (!input.good())
    {
      cout << "Failed to open Best Parameter file: " << fileinput << "\n";
      exit(0);
    }

  delete[] angles;
  angles = new myfloat3_t[ 1 ] ; //Only best orientation

  char line[512] = {0};
  char saveline[512];

  cout << "\n +++++++++++++++++++++++++++++++++++++++++ \n";
  cout << "\n     ONLY READING BEST PARAMETERS \n";
  cout << "\n     FOR PRINTING MAXIMIZED MAP \n";
  cout << " +++++++++++++++++++++++++++++++++++++++++ \n";
  while (!input.eof())
    {
      input.getline(line, 512);
      strcpy(saveline, line);
      char *token = strtok(line, " ");

      if (token == NULL || line[0] == '#' || strlen(token) == 0)
392
	{
393
	  // comment or blank line
394
	}
395
      else if (strcmp(token, "PIXEL_SIZE") == 0)
396
	{
397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419
	  token = strtok(NULL, " ");
	  pixelSize = atof(token);
	  if (pixelSize < 0 ) { cout << "*** Error: Negative pixelSize "; exit(1);}
	  cout << "Pixel Sixe " << pixelSize << "\n";
	}
      else if (strcmp(token, "NUMBER_PIXELS") == 0)
	{
	  token = strtok(NULL, " ");
	  param_device.NumberPixels = int(atoi(token));
	  if (param_device.NumberPixels < 0 ) { cout << "*** Error: Negative Number of Pixels "; exit(1);}
	  cout << "Number of Pixels " << param_device.NumberPixels << "\n";
	}
      else if (strcmp(token, "BEST_ALPHA") == 0)
	{
	  token = strtok(NULL, " ");
	  angles[0].pos[0] = atof(token);
	  cout << "Best Alpha " <<  angles[0].pos[0] << "\n";
	}
      else if (strcmp(token, "BEST_BETA") == 0)
	{
	  token = strtok(NULL, " ");
	  angles[0].pos[1] = atof(token);
	  cout << "Best Alpha " <<  angles[0].pos[1] << "\n";
420
	}
421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500
      else if (strcmp(token, "BEST_GAMMA") == 0)
	{
	  token = strtok(NULL, " ");
	  angles[0].pos[2] = atof(token);
	  cout << "Best Gamma " <<  angles[0].pos[2] << "\n";
	}
               
      else if (strcmp(token, "BEST_ENVELOPE") == 0)
	{
	  token = strtok(NULL, " ");
	  startGridEnvelop = atof(token);
	  if (startGridEnvelop < 0 ) { cout << "*** Error: Negative START_ENVELOPE "; exit(1);}
	  cout << "Best Envelope PSF " << startGridEnvelop << "\n";
	}
      else if (strcmp(token,"BEST_PHASE")==0)
	{
	  token = strtok(NULL," ");
	  startGridCTF_phase=atof(token);
	  cout << "Best Phase PSF " << startGridCTF_phase << "\n";
	}
      else if (strcmp(token,"BEST_AMP")==0)
	{
	  token = strtok(NULL," ");
	  startGridCTF_amp=atof(token);
	  cout << "Best Amplitud PSF " << startGridCTF_amp << "\n";
	}
      else if (strcmp(token, "BEST_DX") == 0)
	{
	  token = strtok(NULL, " ");
	  ddx = atoi(token);
	  cout << "Best dx " << ddx << "\n";
	}
      else if (strcmp(token, "BEST_DY") == 0)
	{
	  token = strtok(NULL, " ");
	  ddy = atoi(token);
	  cout << "Best dy " << ddy << "\n";
	}
      else if (strcmp(token, "BEST_NORM") == 0)
	{
	  token = strtok(NULL, " ");
	  bestnorm= atof(token);
	  cout << "Best norm " << bestnorm << "\n";
	}
      else if (strcmp(token, "BEST_OFFSET") == 0)
	{
	  token = strtok(NULL, " ");
	  bestoff = atof(token);
	  cout << "Best offset " << bestoff << "\n";
	}
      else if (strcmp(token, "WITHNOISE") == 0)
	{
	  token = strtok(NULL, " ");
	  stnoise = atof(token);
	  withnoise=true;
	  cout << "Including noise with standard deviation " << stnoise << "\n";
	}
      else if (strcmp(token, "DO_AA_RADIUS") == 0) //If projecting CA with amino-acid radius
	{
	  doaaradius = true;
	  cout << "Projecting corresponding radius \n";
	}


    }

  input.close();

  //Automatic definitions
  numberGridPointsCTF_amp = 1 ;
  gridCTF_amp = startGridCTF_amp;
  numberGridPointsCTF_phase = 1;
  gridCTF_phase = startGridCTF_phase;
  numberGridPointsEnvelop = 1 ;
  gridEnvelop = startGridEnvelop;

  param_device.NumberFFTPixels1D = param_device.NumberPixels / 2 + 1;
  FFTMapSize = param_device.NumberPixels * param_device.NumberFFTPixels1D;

  return 0;
501

qon's avatar
qon committed
502 503
}

504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543
void bioem_param::PrepareFFTs()
{
  if (mpi_rank == 0) cout << "Preparing FFTs\n";
  releaseFFTPlans();
  mycomplex_t *tmp_map, *tmp_map2;
  tmp_map = (mycomplex_t *) myfftw_malloc(sizeof(mycomplex_t) * param_device.NumberPixels * param_device.NumberPixels);
  tmp_map2 = (mycomplex_t *) myfftw_malloc(sizeof(mycomplex_t) * param_device.NumberPixels * param_device.NumberPixels);
  Alignment = 64;

  fft_plan_c2c_forward = myfftw_plan_dft_2d(param_device.NumberPixels, param_device.NumberPixels, tmp_map, tmp_map2, FFTW_FORWARD, FFTW_MEASURE | FFTW_DESTROY_INPUT);
  fft_plan_c2c_backward = myfftw_plan_dft_2d(param_device.NumberPixels, param_device.NumberPixels, tmp_map, tmp_map2, FFTW_BACKWARD, FFTW_MEASURE | FFTW_DESTROY_INPUT);
  fft_plan_r2c_forward = myfftw_plan_dft_r2c_2d(param_device.NumberPixels, param_device.NumberPixels, (myfloat_t*) tmp_map, tmp_map2, FFTW_MEASURE | FFTW_DESTROY_INPUT);
  fft_plan_c2r_backward = myfftw_plan_dft_c2r_2d(param_device.NumberPixels, param_device.NumberPixels, tmp_map, (myfloat_t*) tmp_map2, FFTW_MEASURE | FFTW_DESTROY_INPUT);

  if (fft_plan_c2c_forward == 0 || fft_plan_c2c_backward == 0 || fft_plan_r2c_forward == 0 || fft_plan_c2r_backward == 0)
    {
      cout << "Error planing FFTs\n";
      exit(1);
    }

  myfftw_free(tmp_map);
  myfftw_free(tmp_map2);

  const int count = omp_get_max_threads();
  fft_scratch_complex = new mycomplex_t*[count];
  fft_scratch_real = new myfloat_t*[count];
#pragma omp parallel
  {
#pragma omp critical
    {
      const int i = omp_get_thread_num();
      fft_scratch_complex[i] = (mycomplex_t *) myfftw_malloc(sizeof(mycomplex_t) * param_device.NumberPixels * param_device.NumberFFTPixels1D);
      fft_scratch_real[i] = (myfloat_t *) myfftw_malloc(sizeof(myfloat_t) * param_device.NumberPixels * param_device.NumberPixels);
    }
  }

  fft_plans_created = 1;
}


544 545
void bioem_param::releaseFFTPlans()
{
546 547 548 549
  if (fft_plans_created)
    {
      const int count = omp_get_max_threads();
      for (int i = 0;i < count;i++)
550
	{
551 552
	  myfftw_free(fft_scratch_complex[i]);
	  myfftw_free(fft_scratch_real[i]);
553
	}
554 555 556 557 558 559 560 561 562 563
      delete[] fft_scratch_complex;
      delete[] fft_scratch_real;

      myfftw_destroy_plan(fft_plan_c2c_forward);
      myfftw_destroy_plan(fft_plan_c2c_backward);
      myfftw_destroy_plan(fft_plan_r2c_forward);
      myfftw_destroy_plan(fft_plan_c2r_backward);
      myfftw_cleanup();
    }
  fft_plans_created = 0;
564 565
}

qon's avatar
qon committed
566 567
int bioem_param::CalculateGridsParam() //TO DO FOR QUATERNIONS
{
568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651
  // **************************************************************************************
  // **************** Routine that pre-calculates Euler angle grids **********************
  // ************************************************************************************

  myfloat_t voluang;

  if(!notuniformangles){
    myfloat_t grid_alpha, cos_grid_beta;
    int n = 0;

    //alpha and gamma are uniform in -PI,PI
    grid_alpha = 2.f * M_PI / (myfloat_t) angleGridPointsAlpha;

    //cosine beta is uniform in -1,1
    cos_grid_beta = 2.f / (myfloat_t) angleGridPointsBeta;

    // Euler Angle Array
    delete[] angles;
    angles = new myfloat3_t[angleGridPointsAlpha * angleGridPointsBeta * angleGridPointsAlpha];
    for (int ialpha = 0; ialpha < angleGridPointsAlpha; ialpha ++)
      {
	for (int ibeta = 0; ibeta < angleGridPointsBeta; ibeta ++)
	  {
	    for (int igamma = 0; igamma < angleGridPointsAlpha; igamma ++)
	      {
		angles[n].pos[0] = (myfloat_t) ialpha * grid_alpha - M_PI + grid_alpha * 0.5f; //ALPHA centered in the middle
		angles[n].pos[1] = acos((myfloat_t) ibeta * cos_grid_beta - 1 + cos_grid_beta * 0.5f); //BETA centered in the middle
		angles[n].pos[2] = (myfloat_t) igamma * grid_alpha - M_PI + grid_alpha * 0.5f; //GAMMA centered in the middle
		n++;
	      }
	  }
      }
    nTotGridAngles = n;
    voluang= grid_alpha * grid_alpha * cos_grid_beta / (2.f * M_PI) / (2.f * M_PI) / 2.f ;

  } else{

    ifstream input(inanglef.c_str());

    if (!input.good())
      {
	cout << "Euler Angle File Failed to open file " <<  inanglef.c_str() << " " << endl ;
	exit(1);
      }

    char line[512] = {0};
    int n=0;

    if(NotUn_angles<1) {
      cout << "\nNot defined number of Euler angles in INPUT file:" << endl ;
      cout << "Use key word: NOT_UNIFORM_TOTAL_ANGS\n";
      exit(1);
    }

    delete[] angles;
    angles = new myfloat3_t[ NotUn_angles] ;

    while (!input.eof())
      {
	input.getline(line, 511);
	if(n< NotUn_angles){
	  float a,b,g;

	  char tmpVals[36]  = {0};

	  strncpy (tmpVals, line, 12);
	  sscanf (tmpVals, "%f", &a);

	  strncpy (tmpVals, line + 12, 12);
	  sscanf (tmpVals, "%f", &b);

	  strncpy (tmpVals, line + 24, 12);
	  sscanf (tmpVals, "%f", &g);

	  angles[n].pos[0] = a;
	  angles[n].pos[1] = b;
	  angles[n].pos[2] = g;
	  //cout << "angs " << angles[n].pos[0] << " " << angles[n].pos[1] << " " << angles[n].pos[2] << "\n";
	}
	n++;
	if(NotUn_angles+1 < n) {
	  cout << "Not properly defined total Euler angles " << n << " instead of " << NotUn_angles << "\n";
	  exit(1);
	}
652

653 654 655 656 657
	
      }
    nTotGridAngles = NotUn_angles;
    voluang= 1./ (myfloat_t) NotUn_angles;
  }
658 659


660
  // ********** Calculating normalized volumen element *********
661

662 663 664
  param_device.volu = voluang * (myfloat_t) param_device.GridSpaceCenter * pixelSize * (myfloat_t) param_device.GridSpaceCenter * pixelSize
    * gridCTF_phase * gridCTF_amp * gridEnvelop / ((2.f * (myfloat_t) param_device.maxDisplaceCenter+1.)) / (2.f * (myfloat_t) (param_device.maxDisplaceCenter + 1.)) / ((myfloat_t) numberGridPointsCTF_amp * gridCTF_amp )
    / ((myfloat_t) numberGridPointsCTF_phase * gridCTF_phase) / ((myfloat_t) numberGridPointsEnvelop * gridEnvelop );
665

666 667
  //  cout << "VOLU " << param_device.volu  << " " << gridCTF_amp << "\n";
  // *** Number of total pixels***
668

669 670
  param_device.Ntotpi = (myfloat_t) (param_device.NumberPixels * param_device.NumberPixels);
  param_device.NtotDist = (2 * (int) (param_device.maxDisplaceCenter / param_device.GridSpaceCenter) + 1 ) * (2 * (int) (param_device.maxDisplaceCenter / param_device.GridSpaceCenter) + 1);
671

672
  nTotCC = (int) ((myfloat_t) param_device.NumberPixels / (myfloat_t) param_device.CCdisplace + 1) * (int) ((myfloat_t) param_device.NumberPixels / (myfloat_t) param_device.CCdisplace + 1);  
673
  return(0);
qon's avatar
qon committed
674 675 676 677 678

}

int bioem_param::CalculateRefCTF()
{
679 680 681
  // **************************************************************************************
  // ********** Routine that pre-calculates Kernels for Convolution **********************
  // ************************************************************************************
682

683 684 685 686 687
  myfloat_t amp, env, phase, ctf, radsq;
  myfloat_t* localCTF;
  mycomplex_t* localout;
  int nctfmax = param_device.NumberPixels / 2;
  int n = 0;
688

689 690
  localCTF = (myfloat_t *) myfftw_malloc(sizeof(myfloat_t) * param_device.NumberPixels * param_device.NumberPixels);
  localout = (mycomplex_t *) myfftw_malloc(sizeof(mycomplex_t) * param_device.NumberPixels * param_device.NumberFFTPixels1D);
691

692 693 694 695 696
  nTotCTFs = numberGridPointsCTF_amp * numberGridPointsCTF_phase * numberGridPointsEnvelop;
  delete[] refCTF;
  refCTF = new mycomplex_t[getRefCtfCount()];
  delete[] CtfParam;
  CtfParam = new myfloat3_t[getCtfParamCount()];
David Rohr's avatar
David Rohr committed
697

698
  myfloat_t normctf;
Pilar Cossio's avatar
Pilar Cossio committed
699

700 701 702 703 704
  for (int iamp = 0; iamp <  numberGridPointsCTF_amp ; iamp++) //Loop over amplitud
    {
      amp = (myfloat_t) iamp * gridCTF_amp + startGridCTF_amp;

      for (int iphase = 0; iphase <  numberGridPointsCTF_phase ; iphase++)//Loop over phase
705
	{
706 707 708 709 710 711 712
	  phase = (myfloat_t) iphase * gridCTF_phase + startGridCTF_phase;

	  for (int ienv = 0; ienv <  numberGridPointsEnvelop ; ienv++)//Loop over envelope
	    {
	      env = (myfloat_t) ienv * gridEnvelop + startGridEnvelop;

	      memset(localCTF, 0, param_device.NumberPixels * param_device.NumberPixels * sizeof(myfloat_t));
713

714 715 716 717
	      normctf = 0.0;

	      //                                cout << "CTF param " << amp << " " << phase << " " << env  << " \n" ;
	      for(int i = 0; i < nctfmax; i++)
718
		{
719 720 721 722 723 724
		  for(int j = 0; j < nctfmax; j++)
		    {
		      radsq = (myfloat_t) (i * i + j * j) * pixelSize * pixelSize;
		      ctf = exp(-radsq * env / 2.0) * (amp * cos(radsq * phase / 2.0) - sqrt((1 - amp * amp)) * sin(radsq * phase / 2.0)) ;
		      normctf += 4.0 * (myfloat_t) ctf ;
		    }
725 726
		}

727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756

	      //Assigning CTF values
	      for(int i = 0; i < nctfmax; i++)
		{
		  for(int j = 0; j < nctfmax; j++)
		    {
		      radsq = (myfloat_t) (i * i + j * j) * pixelSize * pixelSize;
		      ctf = exp(-radsq * env / 2.0) * (amp * cos(radsq * phase / 2.0) - sqrt((1 - amp * amp)) * sin(radsq * phase / 2.0)) / normctf ;

		      localCTF[i * param_device.NumberPixels + j] = (myfloat_t) ctf;
		      localCTF[i * param_device.NumberPixels + param_device.NumberPixels - j - 1] = (myfloat_t) ctf;
		      localCTF[(param_device.NumberPixels - i - 1)*param_device.NumberPixels + j] = (myfloat_t) ctf;
		      localCTF[(param_device.NumberPixels - i - 1)*param_device.NumberPixels + param_device.NumberPixels - j - 1] = (myfloat_t) ctf;
		    }
		}
	      //Calling FFT_Forward
	      myfftw_execute_dft_r2c(fft_plan_r2c_forward, localCTF, localout);

	      // Normalizing and saving the Reference CTFs
	      mycomplex_t* curRef = &refCTF[n * FFTMapSize];
	      for(int i = 0; i < param_device.NumberPixels * param_device.NumberFFTPixels1D; i++ )
		{
		  curRef[i][0] = localout[i][0];
		  curRef[i][1] = localout[i][1];
		}
	      CtfParam[n].pos[0] = amp;
	      CtfParam[n].pos[1] = phase;
	      CtfParam[n].pos[2] = env;
	      n++;
	    }
757
	}
758 759 760 761 762 763 764 765 766 767
    }


  myfftw_free(localCTF);
  myfftw_free(localout);
  if (nTotCTFs != n)
    {
      cout << "Internal error during CTF preparation\n";
      exit(1);
    }
768

769
  return(0);
qon's avatar
qon committed
770 771 772 773
}

bioem_param::~bioem_param()
{
774 775 776 777 778 779 780 781 782 783 784 785
  releaseFFTPlans();
  param_device.NumberPixels = 0;
  angleGridPointsAlpha = 0;
  angleGridPointsBeta = 0;
  numberGridPointsEnvelop = 0;
  numberGridPointsCTF_amp = 0;
  numberGridPointsCTF_phase = 0;
  param_device.maxDisplaceCenter = 0;
  numberGridPointsDisplaceCenter = 0;
  if (refCTF) delete[] refCTF;
  if (CtfParam) delete[] CtfParam;
  if (angles) delete[] angles;
qon's avatar
qon committed
786
}