param.cpp 12.5 KB
Newer Older
qon's avatar
qon committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
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
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
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
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
360
361
362
363
364
365
366
367
368
369
370
371
372
373
#include <stdio.h>
#include <stdlib.h>
#include <iostream>
#include <fstream>
#include <cstring>
#include <math.h>
#include <fftw3.h>

#include "param.h"
#include "map.h"

using namespace std;

bioem_param::bioem_param()
{
    //Number of Pixels
    param_device.NumberPixels=0;
    // Pixel Size
    pixelSize = NULL;
    // Euler angle grid spacing
    angleGridPointsAlpha = 0;
    angleGridPointsBeta = 0;
    //Envelop function paramters
    startGridEnvelop = NULL;
    numberGridPointsEnvelop = 0;
    gridEnvelop = NULL;
    //Contrast transfer function paramters
    startGridCTF_amp = NULL;
    numberGridPointsCTF_amp = 0;
    gridCTF_amp = NULL;
    startGridCTF_phase = NULL;
    numberGridPointsCTF_phase = 0;
    gridCTF_phase = NULL;

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

int bioem_param::readParameters()
{
    /**************************************************************************************/
    /***************************** Reading Input Parameters ******************************/
    /**************************************************************************************/

    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 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 = new myfloat_t;
            *pixelSize=atof(token);
            cout << "Pixel Sixe " << *pixelSize << "\n";

        }
        else if (strcmp(token,"NUMBER_PIXELS")==0)
        {
            token = strtok(NULL," ");
            param_device.NumberPixels=int(atoi(token));
            cout << "Number of Pixels " << param_device.NumberPixels << "\n";

        }
        else if (strcmp(token,"GRIDPOINTS_ALPHA")==0)
        {
            token = strtok(NULL," ");
            angleGridPointsAlpha=int(atoi(token));
            cout << "Grid points alpha " << angleGridPointsAlpha << "\n";

        }
        else if (strcmp(token,"GRIDPOINTS_BETA")==0)
        {
            token = strtok(NULL," ");
            angleGridPointsBeta=int(atoi(token));
            cout << "Grid points beta " << angleGridPointsBeta << "\n";

        }

        else if (strcmp(token,"GRIDPOINTS_ENVELOPE")==0)
        {
            token = strtok(NULL," ");
            numberGridPointsEnvelop=int(atoi(token));
            cout << "Grid points envelope " << numberGridPointsEnvelop << "\n";

        }
        else if (strcmp(token,"START_ENVELOPE")==0)
        {
            token = strtok(NULL," ");
            startGridEnvelop = new myfloat_t;
            *startGridEnvelop=atof(token);
            cout << "Start Envelope " << *startGridEnvelop << "\n";

        }
        else if (strcmp(token,"GRIDSPACE_ENVELOPE")==0)
        {
            token = strtok(NULL," ");
            gridEnvelop= new myfloat_t;
            *gridEnvelop=atof(token);
            cout << "Grid spacing Envelope " << *gridEnvelop << "\n";

        }
        else if (strcmp(token,"GRIDPOINTS_CTF_PHASE")==0)
        {
            token = strtok(NULL," ");
            numberGridPointsCTF_phase=int(atoi(token));
            cout << "Grid points CTF " << numberGridPointsCTF_phase << "\n";

        }
        else if (strcmp(token,"START_CTF_PHASE")==0)
        {
            token = strtok(NULL," ");
            startGridCTF_phase= new myfloat_t;
            *startGridCTF_phase=atof(token);
            cout << "Start CTF " << *startGridCTF_phase << "\n";

        }
        else if (strcmp(token,"GRIDSPACE_CTF_PHASE")==0)
        {
            token = strtok(NULL," ");
            gridCTF_phase=new myfloat_t;
            *gridCTF_phase=atof(token);
            cout << "Grid Space CTF " << *gridCTF_phase << "\n";

        } else if (strcmp(token,"GRIDPOINTS_CTF_AMP")==0)
        {
            token = strtok(NULL," ");
            numberGridPointsCTF_amp=int(atoi(token));
            cout << "Grid points CTF " << numberGridPointsCTF_amp << "\n";

        }
        else if (strcmp(token,"START_CTF_AMP")==0)
        {
            token = strtok(NULL," ");
            startGridCTF_amp= new myfloat_t;
            *startGridCTF_amp=atof(token);
            cout << "Start CTF " << *startGridCTF_amp << "\n";

        }
        else if (strcmp(token,"GRIDSPACE_CTF_AMP")==0)
        {
            token = strtok(NULL," ");
            gridCTF_amp=new myfloat_t;
            *gridCTF_amp=atof(token);
            cout << "Grid Space CTF " << *gridCTF_amp << "\n";

        }
        else if (strcmp(token,"MAX_D_CENTER")==0)
        {
            token = strtok(NULL," ");
            param_device.maxDisplaceCenter=int(atoi(token));
            cout << "Maximum displacement Center " <<  param_device.maxDisplaceCenter << "\n";

        }
        else if (strcmp(token,"PIXEL_GRID_CENTER")==0)
        {
            token = strtok(NULL," ");
            param_device.GridSpaceCenter=int(atoi(token));
            cout << "Grid space displacement center " <<   param_device.GridSpaceCenter << "\n";

        }
        else if (strcmp(token,"WRITE_PROB_ANGLES")==0)
        {
            writeAngles=true;
            cout << "Writing Probabilies of each angle \n";

        }

    }
    input.close();
    cout << " +++++++++++++++++++++++++++++++++++++++++ \n";
    return(0);
}

int bioem_param::CalculateGridsParam() //TO DO FOR QUATERNIONS
{
    /**************************************************************************************/
    /**************** Routine that pre-calculates Euler angle grids **********************/
    /************************************************************************************/
    myfloat_t grid_alpha,cos_grid_beta;
    int n=0;

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

    //cosine beta is uniform in -1,1
    cos_grid_beta=2.0/float(angleGridPointsBeta);

    // Euler Angle Array
    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]= float(ialpha)*grid_alpha-M_PI+grid_alpha*0.5; //ALPHA centered in the middle
                angles[n].pos[1]= acos(float(ibeta)*cos_grid_beta-1+cos_grid_beta*0.5); //BETA centered in the middle
                angles[n].pos[2]= float(igamma)*grid_alpha-M_PI+grid_alpha*0.5; //GAMMA centered in the middle
                n++;
            }
        }
    }
    nTotGridAngles=n;

    /********** Calculating normalized volumen element *********/

    param_device.volu=grid_alpha*grid_alpha*cos_grid_beta*float(param_device.GridSpaceCenter)*(*pixelSize)*float(param_device.GridSpaceCenter)*(*pixelSize)
         *(*gridCTF_phase)*(*gridCTF_amp)*(*gridEnvelop)/(2*M_PI)/(2*M_PI)/2/(2*float(param_device.maxDisplaceCenter))/(2*float(param_device.maxDisplaceCenter))/(float(numberGridPointsCTF_amp)*(*gridCTF_amp)+*startGridCTF_amp)
         /(float(numberGridPointsCTF_phase)*(*gridCTF_phase)+*startGridCTF_phase)/(float(numberGridPointsEnvelop)*(*gridEnvelop)+*startGridEnvelop);

    /*** Number of total pixels***/

    param_device.Ntotpi=float(param_device.NumberPixels*param_device.NumberPixels);
    return(0);

}


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

    myfloat_t amp,env,phase,ctf,radsq;
    fftw_plan plan;
    mycomplex_t* localCTF;
    int nctfmax=int(float(param_device.NumberPixels)/2.0);
    int n=0;

    localCTF= (mycomplex_t *) fftw_malloc(sizeof(mycomplex_t) *4*param_device.NumberPixels*param_device.NumberPixels);
    refCTF = new bioem_map_forFFT[MAX_REF_CTF];


    for (int iamp = 0; iamp <  numberGridPointsCTF_amp ; iamp++) //Loop over amplitud
    {
        amp=float(iamp)*(*gridCTF_amp) + *startGridCTF_amp;

        for (int iphase = 0; iphase <  numberGridPointsCTF_phase ; iphase++)//Loop over phase
        {
            phase=float(iphase)*(*gridCTF_phase) + *startGridCTF_phase;

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

                memset(localCTF[0],0,4*param_device.NumberPixels*param_device.NumberPixels*sizeof(double));
                memset(localCTF[1],0,4*param_device.NumberPixels*param_device.NumberPixels*sizeof(double));

                //Assigning CTF values
                for(int i=0; i< nctfmax; i++)
                {
                    for(int j=0; j< nctfmax; j++)
                    {
                        radsq=float(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));

                        localCTF[i*2*param_device.NumberPixels+j][0]=double(ctf);
                        localCTF[i*2*param_device.NumberPixels+2*param_device.NumberPixels-j-1][0]=double(ctf);
                        localCTF[(2*param_device.NumberPixels-i-1)*2*param_device.NumberPixels+j][0]=double(ctf);
                        localCTF[(2*param_device.NumberPixels-i-1)*2*param_device.NumberPixels+2*param_device.NumberPixels-j-1][0]=double(ctf);

                    }
                }
                //Creatting FFT_Forward of Kernel to store
                mycomplex_t* localout;
                localout= (mycomplex_t *) fftw_malloc(sizeof(mycomplex_t) *4*param_device.NumberPixels*param_device.NumberPixels);
                //Calling FFT_Forward
                plan = fftw_plan_dft_2d(2*param_device.NumberPixels,2*param_device.NumberPixels,localCTF,localout,FFTW_FORWARD,FFTW_ESTIMATE);
                fftw_execute(plan);
                fftw_destroy_plan(plan);

                // Normalizing and saving the Reference CTFs
                for(int i=0; i < 2*param_device.NumberPixels ; i++ )
                {
                    for(int j=0; j < 2*param_device.NumberPixels ; j++ )
                    {
                        refCTF[n].cpoints[i*2*param_device.NumberPixels+j][0]=localout[i*2*param_device.NumberPixels+j][0]/double(4*param_device.NumberPixels*param_device.NumberPixels);
                        refCTF[n].cpoints[i*2*param_device.NumberPixels+j][1]=localout[i*2*param_device.NumberPixels+j][1]/double(4*param_device.NumberPixels*param_device.NumberPixels);
                    }
                }
                CtfParam[n].pos[0]=amp;
                CtfParam[n].pos[1]=phase;
                CtfParam[n].pos[2]=env;
                n++;
                fftw_free(localout);
                if(n>MAX_REF_CTF)
                {   cout << n << "PROBLEM WITH CTF KERNEL PARAMETERS AND MAX NUMBER ALLOWED\n";
                    exit(1);
                }
            }
        }
    }

    fftw_free(localCTF);
    nTotCTFs=n;

    return(0);
}


bioem_param::~bioem_param()
{
    // delete all parameters
    if (pixelSize)
    {
        delete[] pixelSize;
        pixelSize = NULL;
    }

    if (startGridEnvelop)
    {
        delete[] startGridEnvelop;
        startGridEnvelop = NULL;
    }
    if (gridEnvelop)
    {
        delete[] gridEnvelop;
        gridEnvelop = NULL;
    }
    if (startGridCTF_amp)
    {
        delete[] startGridCTF_amp;
        startGridCTF_amp = NULL;
    }
    if (gridCTF_amp)
    {
        delete[] gridCTF_amp;
        gridCTF_amp = NULL;
    }

    if (startGridCTF_phase)
    {
        delete[] startGridCTF_phase;
        startGridCTF_phase = NULL;
    }
    if (gridCTF_amp)
    {
        delete[] gridCTF_phase;
        gridCTF_phase = NULL;
    }

    param_device.NumberPixels=0;
    angleGridPointsAlpha = 0;
    angleGridPointsBeta = 0;
    numberGridPointsEnvelop = 0;
    numberGridPointsCTF_amp = 0;
    numberGridPointsCTF_phase = 0;
    param_device.maxDisplaceCenter = 0;
    numberGridPointsDisplaceCenter = 0;

}