TPCCLIB
Loading...
Searching...
No Matches
trmrp.c File Reference

Transmission image reconstruction using MRP. More...

#include "libtpcrec.h"

Go to the source code of this file.

Functions

int trmrp (float *bla, float *tra, int dim, float *image, int iter, int os_sets, int rays, int views, int maskdim, float zoom, float beta, float axial_fov, float sample_distance, int skip_prior, int osl, float shiftX, float shiftY, float rotation, int verbose)

Detailed Description

Transmission image reconstruction using MRP.

Remarks
Based on the program fbprec (Feb 1998) written by Sakari Alenius for Sun UNIX workstations.
Author
Vesa Oikonen

Definition in file trmrp.c.

Function Documentation

◆ trmrp()

int trmrp ( float * bla,
float * tra,
int dim,
float * image,
int iter,
int os_sets,
int rays,
int views,
int maskdim,
float zoom,
float beta,
float axial_fov,
float sample_distance,
int skip_prior,
int osl,
float shiftX,
float shiftY,
float rotation,
int verbose )

Median Root Prior (MRP) reconstruction of one 2D data matrix given as an array of floats.

See also
fbp, reprojection
Returns
Returns 0 if ok.
Parameters
blaFloat array containing rays*views blank sinogram values. Data must be normalization-corrected.
traFloat array containing rays*views transmission sinogram values. Data must be normalization-corrected.
dimImage x and y dimensions.
imagePointer to pre-allocated image data; size must be at least dim*dim; log-transformed attenuation correction factors will be written in here.
iterNr of iterations.
os_setsLength of ordered subset process order array; 1, 2, 4, ... 128.
raysNr of rays (bins or columns) in sinogram data.
viewsNr of views (rows) in sinogram data.
maskdimMask dimension; 3 or 5 (9 or 21 pixels).
zoomReconstruction zoom.
betaBeta.
axial_fovAxial field-of-view in mm (found in transmission mainheader in cm).
sample_distanceSample distance in mm (found in transmission subheader in cm).
skip_priorNumber of iteration before prior; usually 1.
oslUse OSL-type (0 or 1).
shiftXPossible shifting in x dimension (mm).
shiftYPossible shifting in y dimension (mm).
rotationPossible image rotation, -180 - +180 (in degrees).
verboseVerbose level; if zero, then nothing is printed to stderr or stdout.

Definition at line 21 of file trmrp.c.

63 {
64 if(verbose>0) printf("trmrp()\n");
65
66 if(bla==NULL || tra==NULL || image==NULL) return(1);
67 if(rays<2 || views<2 || dim<2 || iter<1 || os_sets<1) return(1);
68 if(maskdim!=3 && maskdim!=5) return(1);
69 if(zoom<0.05) return(1);
70 //if(dim%2) return(2);
71
72
73 /* Set scale */
74 int recrays;
75 float scale, bp_zoom;
76 recrays=(int)((float)dim*zoom);
77 bp_zoom=zoom*(float)dim/(float)recrays;
78 scale=((float)recrays/(float)rays)*bp_zoom*bp_zoom/(float)views;
79 if(verbose>1) {
80 printf(" recrays := %d\n", recrays);
81 printf(" bp_zoom := %g\n", bp_zoom);
82 printf(" scale := %g\n", scale);
83 }
84
85 int views_in_set;
86 views_in_set=views/os_sets;
87 if(verbose>1) printf(" views_in_set := %d\n", views_in_set);
88
89 /* Make the Ordered Subset process order (bit-reversed sequence) */
90 int seq[os_sets]; set_os_set(os_sets, seq);
91 if(verbose>2) {
92 printf("os_sets :=");
93 for(int i=0; i<os_sets; i++) printf(" %d", seq[i]);
94 printf("\n");
95 }
96 /* Arrange transmission and blank sinograms, and interpolate so that
97 bin width equals pixel width */
98 float recbla[recrays*views], rectra[recrays*views];
99 {
100 float blaset[rays*views], traset[rays*views];
101 /* arrange */
102 for(int s=0; s<os_sets; s++) {
103 for(int j=0; j<views_in_set; j++) {
104 memcpy((char*)(blaset + s*rays*views_in_set + j*rays),
105 (char*)(bla + j*rays*os_sets + s*rays), rays*sizeof(float));
106 memcpy((char*)(traset + s*rays*views_in_set + j*rays),
107 (char*)(tra + j*rays*os_sets + s*rays), rays*sizeof(float));
108 }
109 }
110 /* interpolate */
111 recInterpolateSinogram(blaset, recbla, rays, recrays, views);
112 recInterpolateSinogram(traset, rectra, rays, recrays, views);
113 }
114
115
116 /* Sum of blank and transmission */
117 float bsum=0.0, tsum=0.0;
118 for(int i=0; i<recrays*views_in_set; i++) bsum+=recbla[i];
119 for(int i=0; i<recrays*views_in_set; i++) tsum+=rectra[i];
120 if(verbose>2) {
121 printf(" blank_sum := %g\n", bsum);
122 printf(" transmission_sum := %g\n", tsum);
123 }
124
125 /* Make an initial image: an uniform disk enclosed by rays with a value
126 matching the total count */
127 if(verbose>1) printf("creating initial %dx%d image\n", dim, dim);
128 int imgsize=dim*dim;
129 float current_img[imgsize];
130 for(int i=0; i<imgsize; i++) image[i]=0.0;
131 float init;
132 init=-logf(tsum/bsum)*sample_distance/axial_fov;
133 if(verbose>2) printf(" init := %g\n", init);
134
135 for(int k=0; k<imgsize; k++) current_img[k]=0.0;
136 {
137 int j, k=0;
138 for(j=dim/2-1; j>=-dim/2; j--) {
139 for(int i=-dim/2; i<dim/2; i++) {
140 if((int)hypot((double)i, (double)j) < dim/2-1) current_img[k]=init;
141 k++;
142 }
143 }
144 }
145
146 /* Pre-compute the sine tables for back-projection */
147 if(verbose>1) printf("computing sine table\n");
148 float sinB[3*views/2], sinBrot[3*views/2];
149 recSinTables(views, sinB, sinBrot, rotation);
150 for(int i=0; i<3*views/2; i++) sinB[i]/=bp_zoom;
151 for(int i=0; i<3*views/2; i++) sinBrot[i]/=bp_zoom;
152
153 /* Calculate pixel size */
154 float pixsize;
155 pixsize=sample_distance*(float)rays/(zoom*(float)dim);
156 if(verbose>2) printf(" pixsize := %g\n", pixsize);
157 /* ... and convert shifts from mm to pixels */
158 shiftX/=pixsize;
159 shiftY/=pixsize;
160
161
162 /* Iterations */
163 if(verbose>1) {printf("iterations\n"); fflush(stdout);}
164 float muproj[recrays*views_in_set+1];
165 float atnproj[recrays*views_in_set];
166 float projdiff[recrays*views_in_set];
167 float numerator[imgsize], denominator[imgsize];
168 float med_img[imgsize], oslcoefs[imgsize];
169 int itercount=1;
170 do {
171 if(verbose>3) {printf(" iteration %d\n", itercount); fflush(stdout);}
172
173 if(verbose>1) {
174 float mi, ma;
175 fMinMaxFin(current_img, imgsize, &mi, &ma);
176 printf(" min=%g max=%g iter=%i\n", mi, ma, itercount);
177 for(int i=0; i<imgsize; i++)
178 if(!isfinite(current_img[i])) {
179 printf(" inf in current image! index=%d, iter=%d\n", i, itercount);
180 break;
181 }
182 }
183
184 for(int s=0; s<os_sets; s++) {
185 if(verbose>4) {printf(" os_set %d; seq=%d\n", 1+s, seq[s]); fflush(stdout);}
186 /* Image reprojection */
187 for(int i=0; i<recrays*views_in_set; i++) muproj[i]=0.0;
188#pragma omp parallel for
189 for(int i=0; i<views_in_set; i++) {
190 int view=seq[s]+i*os_sets;
191 if(verbose>8 && (i==0 || i==views_in_set-1)) printf(" reprojecting view %d\n", view);
192 viewReprojection(current_img, muproj+i*recrays, view, dim,
193 views, recrays, sinB, sinBrot, shiftX, shiftY, bp_zoom);
194 //re_proj(current_img, muproj+i*recrays, seq[s]+i*sets, dim, views);
195 }
196
197 /* Calculate correction = measured / re-projected */
198#pragma omp parallel for
199 for(int i=0; i<recrays*views_in_set; i++) {
200 float mup=muproj[i]*scale;
201 float est_count=expf(-mup)*recbla[i];
202 atnproj[i]=est_count*mup;
203 projdiff[i]=est_count - rectra[i];
204 }
205 /* Make numerator and denominator images */
206 for(int i=0; i<imgsize; i++) numerator[i]=0.0;
207 for(int i=0; i<imgsize; i++) denominator[i]=0.0;
208#pragma omp parallel for
209 for(int i=0; i<views_in_set; i++) {
210 int view=seq[s]+i*os_sets;
211 viewBackprojection(projdiff+i*recrays, numerator, dim,
212 view, views, recrays, sinB, sinBrot, shiftX, shiftY, bp_zoom);
213 viewBackprojection(atnproj+i*recrays, denominator, dim,
214 view, views, recrays, sinB, sinBrot, shiftX, shiftY, bp_zoom);
215 }
216 if(verbose>4) {
217 float mi, ma;
218 fMinMaxFin(numerator, imgsize, &mi, &ma);
219 printf(" numerator_range := %g - %g\n", mi, ma);
220 fMinMaxFin(denominator, imgsize, &mi, &ma);
221 printf(" denominator_range := %g - %g\n", mi, ma);
222 }
223
224 /* Apply the prior */
225 if(skip_prior<=0 && beta>0.0) {
226 if(verbose>3) {printf(" applying prior\n"); fflush(stdout);}
227 if(osl) {
228 float maxv, maxm;
229 fMinMaxFin(current_img, imgsize, NULL, &maxv);
230 do_prior(current_img, beta, oslcoefs, dim, 1.0E-08*maxv, maskdim, &maxm);
231 if(verbose>6) {
232 printf(" max value in current image := %g\n", maxv);
233 printf(" max median coefficient := %g\n", maxm);
234 }
235 } else {
236 float *iptr, *mptr;
237 if(maskdim==3) {
238 iptr=current_img+dim+1;
239 mptr=med_img+dim+1;
240 for(int i=2*(dim+1); i<imgsize; i++) *mptr++=med9(iptr++, dim);
241 } else if(maskdim==5) {
242 iptr=current_img+2*dim+2;
243 mptr=med_img+2*dim+2;
244 for(int i=2*(2*dim+2); i<imgsize; i++) *mptr++=med21(iptr++, dim);
245 }
246 for(int i=0; i<imgsize; i++) numerator[i]*=med_img[i];
247 for(int i=0; i<imgsize; i++) numerator[i]-=beta*(current_img[i]-med_img[i]);
248 for(int i=0; i<imgsize; i++) denominator[i]*=med_img[i];
249 for(int i=0; i<imgsize; i++) denominator[i]+=beta*current_img[i];
250 }
251 }
252 /* Calculate the next image */
253 if(verbose>3) {printf(" calculating next image\n"); fflush(stdout);}
254 if(osl && skip_prior<=0 && beta>0.0) { /* convex-OSL */
255 if(verbose>4) printf(" convex-OSL\n");
256 for(int i=0; i<imgsize; i++) {
257 float f=fmaf(numerator[i], 1.0/denominator[i], 1.0);
258 f*=oslcoefs[i];
259 f*=current_img[i];
260 if(isfinite(f)) current_img[i]=f;
261 }
262 } else { /* convex */
263 if(verbose>4) printf(" convex\n");
264 for(int i=0; i<imgsize; i++) {
265 float f=fmaf(numerator[i], 1.0/denominator[i], 1.0);
266 f*=current_img[i];
267 if(isfinite(f)) current_img[i]=f;
268 }
269 }
270 if(verbose>4) {printf(" -> next os_set maybe\n"); fflush(stdout);}
271 } // next set
272 itercount++; skip_prior--;
273 if(verbose>3) {printf(" -> next iteration maybe\n"); fflush(stdout);}
274 } while(itercount<iter);
275 if(verbose>2) {printf(" iterations done.\n"); fflush(stdout);}
276
277 for(int i=0; i<imgsize; i++) image[i]=current_img[i]*=scale;
278
279 if(verbose>1) {printf("trmrp() done.\n"); fflush(stdout);}
280 return(0);
281}
void fMinMaxFin(float *data, long long int n, float *fmin, float *fmax)
Definition imgminmax.c:649
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
float med21(float *inp, int dim)
Definition mrprior.c:49
void set_os_set(int os_sets, int *set_seq)
Definition recutil.c:113
float med9(float *inp, int dim)
Definition mrprior.c:15
void do_prior(float *img, float beta, float *med_coef, int dim, float small, int maskdim, float *maxm)
Definition mrprior.c:77

Referenced by atnMake().