TPCCLIB
Loading...
Searching...
No Matches
mrp.c
Go to the documentation of this file.
1
9/*****************************************************************************/
10#include "libtpcrec.h"
11#ifdef HAVE_OMP_H
12#include <omp.h>
13#endif
14/*****************************************************************************/
15
16/*****************************************************************************/
24 IMG *scn,
26 IMG *img,
28 int imgDim,
30 float zoom,
32 float shiftX,
34 float shiftY,
36 float rotation,
38 int maxIterNr,
40 int skipPriorNr,
42 float beta,
44 int maskDim,
46 int osSetNr,
48 int verbose
49) {
50 if(verbose>1)
51 printf("imgMRP(scn, img, %d, %g, %g, %g, %g, %d, %d, %g, %d, %d)\n",
52 imgDim, zoom, shiftX, shiftY, rotation, maxIterNr, skipPriorNr, beta,
53 maskDim, osSetNr);
54
55
56 /* Check the arguments */
57 if(scn->status!=IMG_STATUS_OCCUPIED) return(1);
58 if(imgDim<2 || imgDim>4096 || imgDim%2) return(1);
59 if(zoom<0.01 || zoom>1000.) return(1);
60 if(scn->dimx<=1 || scn->dimx>16384) return(1);
61 if(beta<0.0) return(1);
62 if(maskDim!=3 && maskDim!=5) return(1);
63
64
65 /*
66 * Allocate output image
67 */
68 if(verbose>1) printf("allocating memory for the image\n");
69 imgEmpty(img);
70 if(imgAllocate(img, scn->dimz, imgDim, imgDim, scn->dimt)!=0) return(3);
71
72 /* Set image "header" information */
73 if(verbose>1) printf("setting image header\n");
75 img->unit=CUNIT_CPS; /* (cnts/sec) */
76 img->scanStart=scn->scanStart;
77 img->axialFOV=scn->axialFOV;
79 img->sizez=scn->sizez;
80 strcpy(img->studyNr, scn->studyNr);
82 if(scn->sampleDistance<=0.0) {
83 if(scn->dimx==281) { // GE Advance
84 img->sizez=4.25;
85 scn->sampleDistance=1.95730;
86 scn->axialFOV=150.; scn->transaxialFOV=550.;
87 } else { // ECAT 931
88 img->sizez=6.75;
89 scn->sampleDistance=3.12932;
90 scn->axialFOV=108.; scn->transaxialFOV=600.829;
91 }
92 }
93 float pixSize; /* note that pixSize is in cm in the ECAT image */
94 pixSize=scn->sampleDistance*(float)scn->dimx/((float)imgDim*zoom);
95 img->sizex=img->sizey=pixSize;
96 int plane, frame;
97 for(plane=0; plane<scn->dimz; plane++)
98 img->planeNumber[plane]=scn->planeNumber[plane];
102 for(frame=0; frame<scn->dimt; frame++) {
103 img->start[frame]=scn->start[frame]; img->end[frame]=scn->end[frame];
104 img->mid[frame]=0.5*(img->start[frame]+img->end[frame]);
105 }
106 img->isWeight=0;
107
108 /*
109 * Preparations for reconstruction
110 */
111
112 /* Pre-compute the sine tables for back-projection (note the rotation!) */
113 if(verbose>1) printf("computing sine tables for back-projection\n");
114 float sinB[3*scn->dimy/2];
115 float sinBrot[3*scn->dimy/2];
116 recSinTables(scn->dimy, sinB, sinBrot, rotation);
117
118 /* Set the backprojection zoom and inverse (globals) */
119 float bpZoom, bpZoomInv;
120 bpZoom=zoom*(float)imgDim/(float)scn->dimx;
121 bpZoomInv=1.0/bpZoom;
122 if(verbose>2) printf("bpZoom=%g bpZoomInv=%g\n", bpZoom, bpZoomInv);
123
124 /* Initialize variables used by back-projection */
125 if(verbose>1) printf("initialize variables for back-projection\n");
126 float offsX, offsY;
127 int halfDim;
128 halfDim=imgDim/2;
129 offsX=shiftX/pixSize; offsY=shiftY/pixSize;
130 for(int i=0; i<3*(scn->dimy)/2; i++) {
131 sinB[i]*=bpZoomInv;
132 sinBrot[i]*=bpZoomInv;
133 }
134 if(verbose>2) {
135 printf("halfDim=%d offsX=%g offsY=%g\n", halfDim, offsX, offsY);
136 }
137
138
139 /*
140 * Reconstruct one matrix at a time
141 */
142 if(verbose>1) printf("reconstruct one matrix at a time...\n");
143 int failed=0;
144#pragma omp parallel for private(frame)
145 for(plane=0; plane<scn->dimz; plane++) {
146 for(frame=0; frame<scn->dimt; frame++) {
147 if(failed) break;
148
149 if(verbose>3) {
150 printf("reconstructing plane %d frame %d\n",
151 scn->planeNumber[plane], frame+1);
152 fflush(stdout);
153 }
154 int i, j, k, ret;
155
156 /* Copy scan data into the array */
157 float scnData[scn->dimx*scn->dimy];
158 float imgData[img->dimx*img->dimy];
159 //float *imgOrigin=imgData+imgDim*(halfDim-1)+halfDim;
160
161 for(i=0, k=0; i<scn->dimy; i++)
162 for(j=0; j<scn->dimx; j++)
163 scnData[k++]=scn->m[plane][i][j][frame];
164 /* Initiate image buffer */
165 for(i=0; i<imgDim*imgDim; i++) imgData[i]=0.0;
166
167 /* Reconstruct */
168 ret=mrp(scnData, scn->dimx, scn->dimy, maxIterNr, osSetNr, maskDim, zoom,
169 beta, skipPriorNr, imgDim, imgData, verbose-3);
170 if(ret!=0) {
171 if(verbose>0) fprintf(stderr, "mrp() return value %d\n", ret);
172 failed=ret; break;
173 }
174
175 /* Copy the image array to matrix data */
176 /* At the same time, correct for the frame length */
177 float f;
178 f=img->end[frame]-img->start[frame];
179 if(f<1.0) f=1.0;
180 f=1.0/f;
181 for(i=0, k=0; i<img->dimy; i++)
182 for(j=0; j<img->dimx; j++)
183 img->m[plane][i][j][frame] = imgData[k++]*f;
184 if(verbose==1) {fprintf(stdout, "."); fflush(stdout);}
185 } /* next frame */
186 } /* next plane */
187 if(verbose==1) {fprintf(stdout, "\n"); fflush(stdout);}
188 if(failed) return(8);
189
190 if(verbose>1) printf("imgMRP() done.\n");
191 return(0);
192}
193/*****************************************************************************/
194
195/*****************************************************************************/
203 float *coef,
205 float *img,
207 float *oimg,
209 int n
210) {
211 for(int i=0; i<n; i++) {
212 oimg[i]=coef[i]*img[i];
213 if(oimg[i]<0.0) oimg[i]=0.0;
214 }
215}
216/*****************************************************************************/
217
218/*****************************************************************************/
231 float *measured,
233 float *proj,
235 float *correct,
237 int os_sets,
239 int rays,
241 int views
242) {
243 for(int i=0; i<rays*views/os_sets; i++) {
244 if(fabs(proj[i])>1.0E-40) {
245 correct[i]=(float)os_sets*measured[i]/proj[i];
246 if(correct[i]>100.) correct[i]=100.;
247 else if(correct[i]<0.0) correct[i]=0.0;
248 } else {
249 correct[i]=0.0;
250 }
251 }
252}
253/*****************************************************************************/
254
255/*****************************************************************************/
262int mrp(
265 float *sinogram,
267 int rays,
269 int views,
271 int iter,
273 int os_sets,
275 int maskdim,
277 float zoom,
279 float beta,
281 int skip_prior,
285 int dim,
287 float *image,
289 int verbose
290) {
291 if(verbose>0)
292 printf("mrp(%d, %d, %d, %d, %d, %g, %g, %d, %d)\n",
293 rays, views, iter, os_sets, maskdim, zoom, beta, skip_prior, dim);
294 if(sinogram==NULL || image==NULL) return(1);
295 if(rays<2 || views<2 || iter<1 || os_sets<1 || dim<2) return(1);
296 if(dim%2 || zoom<0.05) return(2);
297 if(maskdim!=3 && maskdim!=5) return(2);
298
299 int i, k;
300 int imgsize=dim*dim;
301 int halfdim=dim/2;
302 int views_in_set=views/os_sets;
303
304 /* Set scale */
305 int recrays;
306 float scale, bp_zoom;
307 recrays=(int)((float)dim*zoom);
308 bp_zoom=zoom*(float)dim/(float)recrays;
309 scale=((float)recrays/(float)rays)*bp_zoom*bp_zoom/(float)views;
310 if(verbose>1) {
311 printf(" recrays := %d\n", recrays);
312 printf(" bp_zoom := %g\n", bp_zoom);
313 printf(" scale := %g\n", scale);
314 }
315
316 /* Make the Ordered Subset process order (bit-reversed sequence) */
317 int seq[os_sets]; set_os_set(os_sets, seq);
318 if(verbose>2) {
319 printf("os_sets :=");
320 for(i=0; i<os_sets; i++) printf(" %d", seq[i]);
321 printf("\n"); fflush(stdout);
322 }
323
324 /* Arrange the sinogram and interpolate so that
325 bin width equals pixel width */
326 float sino[recrays*views];
327 {
328 float scnset[rays*views];
329 /* arrange */
330 for(int s=0; s<os_sets; s++) {
331 for(int j=0; j<views/os_sets; j++) {
332 memcpy((char*)(scnset + s*rays*views/os_sets + j*rays),
333 (char*)(sinogram + j*rays*os_sets + s*rays), rays*sizeof(float));
334 }
335 }
336 /* interpolate */
337 recInterpolateSinogram(scnset, sino, rays, recrays, views);
338 }
339
340 /* Get some statistics from the sinogram matrix */
341 int scnsize=recrays*views;
342 int nonzeroNr;
343 float counts;
344 nonzeroNr=recGetStatistics(sino, scnsize, &counts, NULL, NULL, 1);
345 if(verbose>1) {
346 printf(" total_counts := %g\n", counts);
347 printf(" non-zeroes := %d/%d\n", nonzeroNr, scnsize);
348 }
349
350
351 /* Make an initial image: an uniform disk enclosed by rays with a value
352 matching the total count */
353 if(verbose>1) printf("creating initial %dx%d image\n", dim, dim);
354 float current_img[imgsize];
355 for(i=0; i<imgsize; i++) image[i]=current_img[i]=0.0;
356 float init;
357 k=0;
358 init=counts/(M_PI*(float)((halfdim-1)*(halfdim-1)));
359 for(int j=halfdim-1; j>=-halfdim; j--) {
360 for(i=-halfdim; i<halfdim; i++) {
361 if((int)hypot((double)i, (double)j) < halfdim-1)
362 current_img[k]=init;
363 k++;
364 }
365 }
366
367
368 /* Pre-compute the sine tables for back-projection */
369 if(verbose>1) printf("computing sine table\n");
370 float sinB[3*views/2];
371 recSinTables(views, sinB, NULL, 0.0);
372 for(i=0; i<3*views/2; i++) sinB[i]/=bp_zoom;
373
374 /* Iterations */
375 if(verbose>1) {printf("iterations\n"); fflush(stdout);}
376 float current_proj[recrays*views_in_set+1];
377 float correction[recrays*views/os_sets];
378 float medcoefs[imgsize], mlcoefs[imgsize];
379 int s, itercount=1, view;
380 do {
381 if(verbose>3) {printf(" iteration %d\n", itercount); fflush(stdout);}
382
383 for(s=0; s<os_sets; s++) {
384 if(verbose>4) {
385 printf(" os_set %d; seq=%d\n", 1+s, seq[s]);
386 fflush(stdout);
387 }
388 /* Image reprojection */
389 for(i=0; i<recrays*views_in_set; i++) current_proj[i]=0.0;
390 for(i=0; i<views_in_set; i++) {
391 view=seq[s]+i*os_sets;
392 viewReprojection(current_img, current_proj+i*recrays, view, dim,
393 views, recrays, sinB, sinB, 0.0, 0.0, bp_zoom);
394 }
395 /* Calculate correction = measured / re-projected */
396 mrpProjectionCorrection(sino+seq[s]*recrays*views_in_set, current_proj,
397 correction, os_sets, recrays, views);
398 /* Make the ML coefficients */
399 for(i=0; i<imgsize; i++) mlcoefs[i]=0.0;
400 for(i=0; i<views_in_set; i++)
401 viewBackprojection(correction+i*recrays, mlcoefs, dim,
402 seq[s]+i*os_sets, views, recrays, sinB, sinB, 0.0, 0.0, bp_zoom);
403 /* Make the prior coefficients */
404 if(skip_prior<=0 && beta>0.0) {
405 if(verbose>3) {printf(" applying prior\n"); fflush(stdout);}
406 float maxv, maxm;
407 fMinMaxFin(current_img, imgsize, NULL, &maxv);
408 do_prior(current_img, beta, medcoefs, dim,
409 1.0E-08*maxv, maskdim, &maxm);
410 if(verbose>3) {
411 printf(" max value in current image := %g\n", maxv);
412 printf(" max median coefficient := %g\n", maxm);
413 }
414 /* Adjust ML coefficients */
415 mrpUpdate(medcoefs, mlcoefs, mlcoefs, imgsize);
416 }
417
418 /* Calculate next image */
419 mrpUpdate(mlcoefs, current_img, current_img, imgsize);
420
421 } // next set
422 itercount++; skip_prior--;
423 if(verbose>3) {printf(" -> next iteration maybe\n"); fflush(stdout);}
424 } while(itercount<iter);
425 if(verbose>2) {printf(" iterations done.\n"); fflush(stdout);}
426
427 for(i=0; i<imgsize; i++) image[i]=current_img[i]*=scale;
428
429 if(verbose>1) {printf("mrp() done.\n"); fflush(stdout);}
430 return(0);
431} // mrp()
432/*****************************************************************************/
433
434/*****************************************************************************/
int imgAllocate(IMG *image, int planes, int rows, int columns, int frames)
Definition img.c:194
void imgEmpty(IMG *image)
Definition img.c:121
void fMinMaxFin(float *data, long long int n, float *fmin, float *fmax)
Definition imgminmax.c:649
#define IMG_STATUS_OCCUPIED
#define IMG_DC_NONCORRECTED
#define IMG_TYPE_IMAGE
Header file for libtpcrec.
void recSinTables(int views, float *sinB, float *sinBrot, float rotation)
Definition recutil.c:12
void viewBackprojection(float *prj, float *idata, int dim, int view, int viewNr, int rayNr, float *sinB, float *sinBrot, float offsX, float offsY, float bpZoom)
void viewReprojection(float *idata, float *sdata, int view, int dim, int viewNr, int rayNr, float *sinB, float *sinBrot, float offsX, float offsY, float bpZoom)
void recInterpolateSinogram(float *srcsino, float *newsino, int srcrays, int newrays, int views)
Definition recutil.c:37
void set_os_set(int os_sets, int *set_seq)
Definition recutil.c:113
int recGetStatistics(float *buf, int n, float *osum, float *omin, float *omax, int skip_zero_mins)
Definition recutil.c:154
void do_prior(float *img, float beta, float *med_coef, int dim, float small, int maskdim, float *maxm)
Definition mrprior.c:77
int mrp(float *sinogram, int rays, int views, int iter, int os_sets, int maskdim, float zoom, float beta, int skip_prior, int dim, float *image, int verbose)
Definition mrp.c:262
int imgMRP(IMG *scn, IMG *img, int imgDim, float zoom, float shiftX, float shiftY, float rotation, int maxIterNr, int skipPriorNr, float beta, int maskDim, int osSetNr, int verbose)
Definition mrp.c:21
void mrpProjectionCorrection(float *measured, float *proj, float *correct, int os_sets, int rays, int views)
Definition mrp.c:229
void mrpUpdate(float *coef, float *img, float *oimg, int n)
Definition mrp.c:201
float sizex
unsigned short int dimx
char type
float sampleDistance
float **** m
char decayCorrection
float transaxialFOV
char unit
char status
time_t scanStart
unsigned short int dimt
int * planeNumber
float sizey
float * start
unsigned short int dimz
unsigned short int dimy
float * end
char radiopharmaceutical[32]
float isotopeHalflife
char studyNr[MAX_STUDYNR_LEN+1]
float axialFOV
float * mid
char isWeight
float sizez