TPCCLIB
Loading...
Searching...
No Matches
imgbfh2o.c
Go to the documentation of this file.
1
8/*****************************************************************************/
9#include "tpcclibConfig.h"
10/*****************************************************************************/
11#include <stdio.h>
12#include <stdlib.h>
13#include <unistd.h>
14#include <string.h>
15#include <math.h>
16#include <time.h>
17/*****************************************************************************/
18#include "libtpccurveio.h"
19#include "libtpcmodext.h"
20#include "libtpcmisc.h"
21#include "libtpcmodel.h"
22#include "libtpcimgio.h"
23#include "libtpcimgp.h"
24/*****************************************************************************/
25#define MAX_N 2
26/*****************************************************************************/
27
28/*****************************************************************************/
29static char *info[] = {
30 "Estimation of rate constants K1, k2 and Va from dynamic PET image in",
31 "ECAT 6.3, ECAT 7.x, NIfTI-1, or Analyze 7.5 file format using",
32 "one-tissue compartment model (1), solved using the basis function",
33 "approach (2, 3).",
34 " ",
35 "When applied to dynamic [O-15]H2O studies, the resulting K1 image",
36 "equals perfusion (blood flow) image. K1 image can be divided by tissue",
37 "density (g/mL) (option -density) and multiplied by 100 (option -dL)",
38 "to achieve the blood flow image in units (mL blood)/((100 g tissue) * min).",
39 " ",
40 "When applied to dynamic [O-15]O2 brain studies, the resulting K1 image",
41 "can be converted to oxygen consumption image by multiplying it by",
42 "arterial oxygen concentration (4) (mL O2 / dL blood) to get the",
43 "parametric image in units mL O2 / ((100 ml tissue) * min).",
44 "The model assumptions hold only when oxygen consumption is 1-6.7",
45 "mL O2/(100g * min) and fit time is set to 300 s or less (4).",
46 " ",
47 "Arterial blood TAC must be corrected for decay and delay, with sample times",
48 "in seconds. Dynamic PET image must be corrected for decay. Fit time must",
49 "be given in seconds.",
50 " ",
51 "Usage: @P [Options] btacfile imgfile fittime flowfile",
52 " ",
53 "Options:",
54 " -mL or -dL",
55 " Units in flow and Va images will be given per mL or per dL,",
56 " respectively. By default, units are per mL.",
57 " -density[=<value>]",
58 " With option -density the flow is calculated per gram or 100g tissue.",
59 " Tissue density can be changed from the default 1.04 g/mL.",
60 " -Vd=<filename>",
61 " Parametric K1/k2 (Vd, apparent p) image is saved.",
62 " -k2=<filename>",
63 " Parametric k2 image is saved; in some situations perfusion calculation",
64 " from k2 can be more accurate than the default assumption of f=K1.",
65 " Perfusion can be calculated from k2 using equation f=k2*pH2O, where",
66 " pH2O is the physiological partition coefficient of water in tissue.",
67 " -Va=<filename>",
68 " Parametric Va image is saved.",
69 " Set -Va=0, if Va=0 is assumed (pre-corrected); otherwise Va is fitted.",
70 " -wss=<filename>",
71 " Weighted sum-of-squares are written in specified image file.",
72 " -thr=<threshold%>",
73 " Pixels with AUC less than (threshold/100 x input AUC) are set to zero;",
74 " default is 0%",
75 " -k2min=<Min k2> and -k2max=<Max k2>",
76 " Enter the minimum and maximum k2 in units 1/min, applying to decay",
77 " corrected data.",
78 " -fmin=<Min K1> and -fmax=<Max K1>",
79 " Enter the minimum and maximum perfusion value; defaults are",
80 " 0.005 and 4.0 mL/(mL*min), respectively.",
81 " -pmin=<Min p> and -pmax=<pmax>",
82 " Enter the minimum and maximum value for apparent partition coefficient",
83 " for water; defaults are 0.3 and 1.0 mL/mL, respectively.",
84 " -nr=<value>",
85 " Set number of basis functions; default is 500, minimum 100.",
86 " -bf=<filename>",
87 " Basis function curves are written in specified file.",
88 " -err=<filename>",
89 " Pixels with their k2 in its min or max value (calculated from min and",
90 " max K1 and p values) in the specified imagefile with values 1 and 2,",
91 " respectively, others with value 0.",
92 " -stdoptions", // List standard options like --help, -v, etc
93 " ",
94 "Example 1. Calculation of perfusion and arterial blood volume image,",
95 " stopping fit at 180 s:",
96 " @P -Va=s2345va.img s2345abfit.kbq s2345dy1.v 180 s2345flow.v",
97 " ",
98 "By default, the units of pixel values in the blood flow (K1) image is",
99 "(mL blood)/((mL tissue) * min), in Vd image (mL blood)/(mL tissue),",
100 "in k2 image 1/min, and in Va image (mL blood/mL tissue),",
101 "but the blood flow and Va units can be changed with above listed options.",
102 " ",
103 "References:",
104 "1. Lammertsma AA, Jones T. J Cereb Blood Flow Metab. 1983;3:416-424.",
105 "2. Koeppe RA et al. J Cereb Blood Flow metab. 1985;5:224-234.",
106 "3. Boellaard R et al. Mol Imaging Biol. 2005;7:273-285.",
107 "4. Ohta S, et al. J Cereb Blood Flow Metab. 1992;12:179-192.",
108 " ",
109 "See also: bfmh2o, imgflow, arlkup, fit_h2o, imgunit, fitdelay, imgcbv",
110 " ",
111 "Keywords: image, modelling, perfusion, radiowater, basis function method",
112 0};
113/*****************************************************************************/
114
115/*****************************************************************************/
116/* Turn on the globbing of the command line, since it is disabled by default in
117 mingw-w64 (_dowildcard=0); in MinGW32 define _CRT_glob instead, if necessary;
118 In Unix&Linux wildcard command line processing is enabled by default. */
119/*
120#undef _CRT_glob
121#define _CRT_glob -1
122*/
123int _dowildcard = -1;
124/*****************************************************************************/
125
126/*****************************************************************************/
130int main(int argc, char **argv)
131{
132 int ai, help=0, version=0, verbose=1;
133 int fitdimt;
134 char inpfile[FILENAME_MAX], petfile[FILENAME_MAX], flowfile[FILENAME_MAX];
135 char vdfile[FILENAME_MAX], k2file[FILENAME_MAX], wssfile[FILENAME_MAX];
136 char errfile[FILENAME_MAX], bfsfile[FILENAME_MAX];
137 char vafile[FILENAME_MAX], tmp[FILENAME_MAX+1], *cptr;
138 int bfNr=500, *bf_opt_nr;
139 float threshold, calcThreshold=0.0;
140 double fittime=0.0, f;
141 double flowMin=0.005, flowMax=4.0; // mL/(min*mL)
142 double pWaterMin=0.3, pWaterMax=1.0; // mL/mL
143 double k2min=-1.0, k2max=-1.0;
144 double wss, p1, p2, p3, rnorm_min;
145 int ret, fi, pi, xi, yi, bi, bi_min;
146 int nosolution_nr=0, thresholded_nr=0;
147 clock_t fitStart, fitFinish;
148 int fitVa=1; // 1=fitted, 0=Vb assumed to be zero
149 int per_dl=0; // 0 or 1
150 int per_gram=0; // 0 or 1
151 double density=1.04;
152
153 DFT blood, tac, bf;
154 IMG img, flowimg, k2img, vdimg, vaimg, wssimg, errimg;
155 /* qr */
156 int m, n, M, N;
157 double **mem, **A, *B, X[MAX_N], *tau, *residual, RNORM, *chain;
158 double *qrweight, **wws, *ws, *wwschain, *ct, *cti;
159
160
161 /*
162 * Get arguments
163 */
164 if(argc==1) {tpcPrintUsage(argv[0], info, stderr); return(1);}
165 inpfile[0]=petfile[0]=flowfile[0]=vdfile[0]=k2file[0]=wssfile[0]=(char)0;
166 errfile[0]=bfsfile[0]=vafile[0]=(char)0;
167 /* Get options */
168 for(ai=1; ai<argc; ai++) if(*argv[ai]=='-') {
169 if(tpcProcessStdOptions(argv[ai], &help, &version, &verbose)==0) continue;
170 cptr=argv[ai]+1; if(*cptr=='-') cptr++; if(cptr==NULL) continue;
171 if(strcasecmp(cptr, "DL")==0) {
172 per_dl=1; continue;
173 } else if(strcasecmp(cptr, "ML")==0) {
174 per_dl=0; continue;
175 } else if(strcasecmp(cptr, "DENSITY")==0) {
176 per_gram=1; continue; /* if plain -density, then use default density */
177 } else if(strncasecmp(cptr, "DENSITY=", 8)==0) {
178 per_gram=1; density=atof_dpi(cptr+8); if(density>0.0) continue;
179 } else if(strncasecmp(cptr, "VA=", 3)==0 || strncasecmp(cptr, "VB=", 3)==0) {
180 strlcpy(vafile, cptr+3, FILENAME_MAX); if(strlen(vafile)) continue;
181 } else if(strncasecmp(cptr, "k2min=", 6)==0) {
182 if(atof_with_check(cptr+6, &k2min)==0 && k2min>=0.0) continue;
183 } else if(strncasecmp(cptr, "k2max=", 6)==0) {
184 if(atof_with_check(cptr+6, &k2max)==0 && k2max>=0.0) continue;
185 } else if(strncasecmp(cptr, "fmin=", 5)==0 && strlen(cptr)>5) {
186 if(atof_with_check(cptr+5, &flowMin)==0 && flowMin>=0.0) continue;
187 } else if(strncasecmp(cptr, "fmax=", 5)==0 && strlen(cptr)>5) {
188 if(atof_with_check(cptr+5, &flowMax)==0 && flowMax>=0.0) continue;
189 } else if(strncasecmp(cptr, "pmin=", 5)==0 && strlen(cptr)>5) {
190 if(atof_with_check(cptr+5, &pWaterMin)==0 && pWaterMin>0.0) continue;
191 } else if(strncasecmp(cptr, "pmax=", 5)==0 && strlen(cptr)>5) {
192 if(atof_with_check(cptr+5, &pWaterMax)==0 && pWaterMax<1.25) continue;
193 } else if(strncasecmp(cptr, "VD=", 3)==0) {
194 strlcpy(vdfile, cptr+3, FILENAME_MAX); if(strlen(vdfile)>1) continue;
195 } else if(strncasecmp(cptr, "k2=", 3)==0) {
196 strlcpy(k2file, cptr+3, FILENAME_MAX); if(strlen(k2file)>1) continue;
197 } else if(strncasecmp(cptr, "NR=", 3)==0) {
198 bfNr=atoi(cptr+3); if(bfNr>5E+04) bfNr=5E+04;
199 if(bfNr>=100) continue;
200 } else if(strncasecmp(cptr, "BF=", 3)==0) {
201 strlcpy(bfsfile, cptr+3, FILENAME_MAX); if(strlen(bfsfile)>0) continue;
202 } else if(strncasecmp(cptr, "WSS=", 4)==0) {
203 strlcpy(wssfile, cptr+4, FILENAME_MAX); if(strlen(wssfile)>0) continue;
204 } else if(strncasecmp(cptr, "ERR=", 4)==0) {
205 strlcpy(errfile, cptr+4, FILENAME_MAX); if(strlen(errfile)>0) continue;
206 } else if(strncasecmp(cptr, "THR=", 4)==0 && strlen(cptr)>4) {
207 cptr+=4; if(isdigit(*cptr) || *cptr=='+' || *cptr=='-') {
208 calcThreshold=0.01*atof_dpi(cptr);
209 if(calcThreshold>=0.0 && calcThreshold<=2.0) continue;
210 }
211 }
212 fprintf(stderr, "Error: invalid option '%s'.\n", argv[ai]);
213 return(1);
214 } else break;
215
216 /* Print help or version? */
217 if(help==2) {tpcHtmlUsage(argv[0], info, ""); return(0);}
218 if(help) {tpcPrintUsage(argv[0], info, stdout); return(0);}
219 if(version) {tpcPrintBuild(argv[0], stdout); return(0);}
220
221 /* Process other arguments, starting from the first non-option */
222 if(ai<argc) strlcpy(inpfile, argv[ai++], FILENAME_MAX);
223 if(ai<argc) strlcpy(petfile, argv[ai++], FILENAME_MAX);
224 if(ai<argc) {
225 if(!atof_with_check(argv[ai], &fittime)) fittime/=60.0;
226 else {
227 fprintf(stderr, "Error: invalid fit time '%s'.\n", argv[ai]);
228 return(1);
229 }
230 ai++;
231 }
232 if(ai<argc) strlcpy(flowfile, argv[ai++], FILENAME_MAX);
233 if(ai<argc) {
234 fprintf(stderr, "Error: invalid argument '%s'.\n", argv[ai]);
235 return(1);
236 }
237 /* Did we get all the information that we need? */
238 if(!flowfile[0]) {
239 fprintf(stderr, "Error: missing command-line argument; use option --help\n");
240 return(1);
241 }
242 if(fittime<=0.0) fittime=1.0E+020;
243 /* Check if Va=0 is assumed */
244 if(strcasecmp(vafile, "NONE")==0 || strcasecmp(vafile, "ZERO")==0 ||
245 strcasecmp(vafile, "0")==0) {vafile[0]=(char)0; fitVa=0;}
246 /* In verbose mode print arguments and options */
247 if(verbose>1) {
248 printf("petfile := %s\n", petfile);
249 printf("inpfile := %s\n", inpfile);
250 printf("flowfile := %s\n", flowfile);
251 if(vdfile[0]) printf("vdfile := %s\n", vdfile);
252 if(k2file[0]) printf("k2file := %s\n", k2file);
253 if(vafile[0]) printf("vafile := %s\n", vafile);
254 if(wssfile[0]) printf("wssfile := %s\n", wssfile);
255 if(errfile[0]) printf("errfile := %s\n", errfile);
256 if(bfsfile[0]) printf("bfsfile := %s\n", bfsfile);
257 printf("fitVa := %d\n", fitVa);
258 printf("bfNr := %d\n", bfNr);
259 printf("per_dl := %d\n", per_dl);
260 printf("per_gram := %d\n", per_gram);
261 if(per_gram!=0) printf("density := %g\n", density);
262 printf("calcThreshold := %g\n", calcThreshold);
263 printf("requested_fittime := %g [min]\n", fittime);
264 if(k2min>0.0) printf("k2min := %g\n", k2min);
265 if(k2max>0.0) printf("k2max := %g\n", k2max);
266 printf("flowMax := %g\n", flowMax);
267 printf("flowMin := %g\n", flowMin);
268 printf("flowMax := %g\n", flowMax);
269 printf("pWaterMin := %g\n", pWaterMin);
270 printf("pWaterMax := %g\n", pWaterMax);
271 }
272 if(verbose>8) {IMG_TEST=verbose-8; SIF_TEST=verbose-8;} else IMG_TEST=SIF_TEST=0;
273
274 /* Check user-defined parameter ranges and calculate range of k2 */
275 if(flowMin>=flowMax) {
276 fprintf(stderr, "Error: invalid range for perfusion (%g - %g).\n",
277 flowMin, flowMax);
278 return(1);
279 }
280 if(pWaterMin>=pWaterMax) {
281 fprintf(stderr, "Error: invalid range for p (%g - %g).\n",
282 pWaterMin, pWaterMax);
283 return(1);
284 }
285 if(k2min<=0.) {
286 k2min=flowMin/pWaterMax;
287 if(verbose>1) printf("k2min := %g [1/min]\n", k2min);
288 }
289 if(k2max<=0.) {
290 k2max=flowMax/pWaterMin;
291 if(verbose>1) printf("k2max := %g [1/min]\n", k2max);
292 }
293 if(k2max<=k2min || k2min<=0.) {
294 fprintf(stderr, "Error: invalid range for k2 (%g - %g).\n", k2min, k2max);
295 return(1);
296 }
297
298
299 /*
300 * Read PET image and input TAC
301 */
302 if(verbose>1) printf("reading data files\n");
303 dftInit(&blood); dftInit(&tac); imgInit(&img);
305 petfile, NULL, inpfile, NULL, NULL, &fittime, &fitdimt, &img,
306 &blood, &tac, 1, stdout, verbose-2, tmp);
307 if(ret!=0) {
308 fprintf(stderr, "Error: %s.\n", tmp);
309 if(verbose>1) printf(" ret := %d\n", ret);
310 return(2);
311 }
312 if(imgNaNs(&img, 1)>0)
313 if(verbose>0) fprintf(stderr, "Warning: missing pixel values.\n");
314 //printf("last x2 = %g\n", tac.x2[tac.frameNr-1]);
315 /* Set time unit to min */
316 dftTimeunitConversion(&blood, TUNIT_MIN);
317 dftTimeunitConversion(&tac, TUNIT_MIN);
318 /* Check that the image is dynamic */
319 if(fitdimt<3) {
320 fprintf(stderr, "Error: too few time frames for fitting.\n");
321 if(verbose>0) imgInfo(&img);
322 imgEmpty(&img); dftEmpty(&blood); dftEmpty(&tac); return(2);
323 }
324 /* Allocate memory for tissue data and integrals */
325 ret=dftAddmem(&tac, 1);
326 if(ret!=0) {
327 fprintf(stderr, "Error (%d) in allocating memory.\n", ret); fflush(stderr);
328 imgEmpty(&img); dftEmpty(&blood); dftEmpty(&tac); return(3);
329 }
330 strcpy(tac.voi[0].voiname, "input");
331 strcpy(tac.voi[1].voiname, "tissue");
332 if(verbose>1) {
333 printf("fittimeFinal := %g min\n", fittime);
334 printf("fitdimt := %d\n", fitdimt);
335 fflush(stdout);
336 }
337
338 /* Determine the threshold */
339 if(verbose>50) dftPrint(&tac);
340 if(verbose>2) printf("input_AUC[%g] := %g\n",
341 tac.x2[tac.frameNr-1], tac.voi[0].y2[tac.frameNr-1]);
342 threshold=calcThreshold*tac.voi[0].y2[tac.frameNr-1]/60.0;
343 if(verbose>2) {printf("threshold_AUC = %g\n", threshold); fflush(stdout);}
344
345
346 /*
347 * Determine the weights for the fit
348 */
349 if(verbose>1) printf("setting weights\n");
350 ret=imgSetWeights(&img, 0, verbose-5);
351 if(ret) {
352 fprintf(stderr, "Warning: cannot calculate weights.\n"); fflush(stderr);
353 /* set weights to 1 */
354 for(fi=0; fi<img.dimt; fi++) img.weight[fi]=1.0;
355 img.isWeight=1;
356 }
357
358
359 /*
360 * Calculate the basis functions
361 */
362 if(verbose>1) {fprintf(stdout, "calculating basis functions\n"); fflush(stdout);}
363 dftInit(&bf);
364 ai=tac.frameNr; tac.frameNr=fitdimt;
365 ret=bfRadiowater(&blood, &tac, &bf, bfNr, k2min, k2max, tmp, verbose-2);
366 //ret=bfRadiowater(&blood, &tac, &bf, bfNr, k2min, k2max, tmp, 10);
367 tac.frameNr=ai;
368 if(ret) {
369 fprintf(stderr, "Error: cannot calculate basis functions (%d).\n", ret); fflush(stderr);
370 imgEmpty(&img); dftEmpty(&blood); dftEmpty(&tac); return(6);
371 }
372 /* Original sampling blood data not needed any more */
373 dftEmpty(&blood);
374 /* Note that basis functions are to be saved later (in bfsfile),
375 after it is known in how many image pixels each basis function was found
376 to give the best fit */
377
378
379 /*
380 * Allocate result images and fill the header info
381 */
382 if(verbose>1) {printf("allocating memory for parametric images\n"); fflush(stdout);}
383 imgInit(&flowimg); imgInit(&k2img); imgInit(&vdimg); imgInit(&vaimg);
384 imgInit(&wssimg); imgInit(&errimg);
385 ret=imgAllocateWithHeader(&flowimg, img.dimz, img.dimy, img.dimx, 1, &img);
386 if(ret!=0) {
387 fprintf(stderr, "Error: out of memory.\n");
388 imgEmpty(&img); dftEmpty(&tac); dftEmpty(&bf);
389 return(8);
390 }
391 flowimg.start[0]=0.0; flowimg.end[0]=fittime*60.0;
393 flowimg.isWeight=0;
394 flowimg.unit=IMGUNIT_ML_PER_ML_PER_MIN;
395 if(k2file[0]) {
396 ret=imgAllocateWithHeader(&k2img, img.dimz,img.dimy,img.dimx, 1, &flowimg);
397 k2img.unit=IMGUNIT_PER_MIN;
398 }
399 if(ret==0 && vdfile[0]) {
400 ret=imgAllocateWithHeader(&vdimg, img.dimz,img.dimy,img.dimx, 1, &flowimg);
401 vdimg.unit=IMGUNIT_UNITLESS;
402 }
403 if(ret==0 && vafile[0]) {
404 ret=imgAllocateWithHeader(&vaimg, img.dimz,img.dimy,img.dimx, 1, &flowimg);
405 vaimg.unit=IMGUNIT_UNITLESS;
406 }
407 if(ret==0 && wssfile[0]) {
408 ret=imgAllocateWithHeader(&wssimg, img.dimz,img.dimy,img.dimx, 1, &flowimg);
409 wssimg.unit=IMGUNIT_UNITLESS;
410 }
411 if(ret==0 && errfile[0]) {
412 ret=imgAllocateWithHeader(&errimg, img.dimz,img.dimy,img.dimx, 1, &flowimg);
413 errimg.unit=IMGUNIT_UNITLESS;
414 }
415 if(ret) {
416 fprintf(stderr, "Error: out of memory.\n"); fflush(stderr);
417 imgEmpty(&img); dftEmpty(&tac); dftEmpty(&bf);
418 imgEmpty(&flowimg); imgEmpty(&k2img); imgEmpty(&vaimg); imgEmpty(&vdimg);
419 imgEmpty(&wssimg); imgEmpty(&errimg);
420 return(8);
421 }
422
423
424 /*
425 * Allocate memory for QR
426 */
427 if(verbose>1) {fprintf(stdout, "allocating memory for QR\n"); fflush(stdout);}
428 M=tac.frameNr; N=2; if(fitVa==0) N--;
429 chain=(double*)malloc((M+1)*N*bf.voiNr * sizeof(double));
430 mem=(double**)malloc(bf.voiNr * sizeof(double*));
431 A=(double**)malloc(M * sizeof(double*));
432 B=(double*)malloc(M*sizeof(double));
433 residual=(double*)malloc(M*sizeof(double));
434 qrweight=(double*)malloc(M*sizeof(double));
435 wwschain=(double*)malloc((M*N+2*M)*sizeof(double));
436 wws=(double**)malloc(M * sizeof(double*));
437 if(chain==NULL || B==NULL || A==NULL || residual==NULL || qrweight==NULL || wwschain==NULL || wws==NULL)
438 {
439 fprintf(stderr, "Error: out of memory.\n"); fflush(stderr);
440 imgEmpty(&img); dftEmpty(&tac); dftEmpty(&bf);
441 imgEmpty(&flowimg); imgEmpty(&k2img); imgEmpty(&vaimg); imgEmpty(&vdimg);
442 imgEmpty(&wssimg); imgEmpty(&errimg);
443 return(8);
444 }
445 for(bi=0; bi<bf.voiNr; bi++) mem[bi]=chain+bi*(M+1)*N;
446 for(m=0; m<M; m++) wws[m]=wwschain+m*N;
447 ws=wwschain+M*N;
448
449 /* Pre-compute QR weights for faster execution */
450 for(m=0; m<M; m++) {
451 if(img.weight[m]<=1.0e-20) qrweight[m]=0.0;
452 else qrweight[m]=sqrt(img.weight[m]);
453 }
454
455 /* Make A matrix, and QR decomposition for it, for all pixels
456 beforehand for faster execution */
457 if(verbose>1) {fprintf(stdout, "calculating QR decomposition\n"); fflush(stdout);}
458 for(bi=0; bi<bf.voiNr; bi++) {
459
460 /* Define memory site for coefficient matrix and vector tau */
461 for(m=0; m<M; m++) A[m]=mem[bi]+m*N;
462 tau=mem[bi]+M*N;
463
464 /* Initiate matrix (A = mem[bi]) */
465 for(m=0; m<M; m++) {
466 A[m][0]=bf.voi[bi].y[m];
467 if(N>1) A[m][1]=tac.voi[0].y[m]; // blood TAC for Va estimation
468 }
469
470 /* Apply data weights */
471 for(m=0; m<M; m++) for(n=0; n<N; n++) A[m][n]*=qrweight[m];
472
473 /* Compute QR decomposition of the coefficient matrix */
474 ret=qr_decomp(A, M, N, tau, wws, ws);
475
476 if(ret>0) { /* Decomposition failed */
477 free(chain); free(B); free(residual);
478 free(A); free(wwschain); free(wws); free(qrweight); free(mem);
479 imgEmpty(&img); dftEmpty(&tac); dftEmpty(&bf);
480 imgEmpty(&flowimg); imgEmpty(&k2img); imgEmpty(&vaimg); imgEmpty(&vdimg);
481 imgEmpty(&wssimg); imgEmpty(&errimg);
482 return (9);
483 }
484 } /* next BF */
485
486
487 /*
488 * Compute pixel-by-pixel
489 */
490 if(verbose>0) {fprintf(stdout, "computing QR pixel-by-pixel\n"); fflush(stdout);}
491 thresholded_nr=0; nosolution_nr=0;
492 ct=tac.voi[1].y; cti=tac.voi[1].y2;
493 /* Allocate memory for BF counters on how often each BF is
494 found to provide the optimal fit */
495 bf_opt_nr=(int*)malloc(bfNr*sizeof(int));
496 for(bi=0; bi<bf.voiNr; bi++) bf_opt_nr[bi]=0.0;
497 fitStart=clock();
498 for(pi=0; pi<img.dimz; pi++) {
499 if(verbose>6) printf("computing plane %d\n", img.planeNumber[pi]);
500 else if(img.dimz>1 && verbose>0) {fprintf(stdout, "."); fflush(stdout);}
501 for(yi=0; yi<img.dimy; yi++) {
502 if(verbose>6 && yi==4*img.dimy/10) printf(" computing row %d\n", yi+1);
503 for(xi=0; xi<img.dimx; xi++) {
504 if(verbose>6 && yi==4*img.dimy/10 && xi==4*img.dimx/10)
505 printf(" computing column %d\n", xi+1);
506
507 /* if end AUC is less than threshold value, then set values to 0 */
508 /* Calculate pixel integral */
509 for(m=0; m<M; m++) {ct[m]=img.m[pi][yi][xi][m];}
510 ret=petintegral(tac.x1, tac.x2, ct, tac.frameNr, cti, NULL);
511//if(verbose>6 && yi==4*img.dimy/10 && xi==4*img.dimx/10)
512//printf("last x2 = %g\n", tac.x2[tac.frameNr-1]);
513 if(verbose>6 && yi==4*img.dimy/10 && xi==4*img.dimx/10) {
514 printf(" Pixel (%d,%d,%d), int= %f, threshold= %g\n",
515 pi, yi, xi, cti[M-1], threshold);
516 if(verbose>7) {
517 for(m=0; m<M; m++)
518 printf(" %02d:\t%g\t%g\n", m, img.m[pi][yi][xi][m], tac.voi[0].y[m]);
519 }
520 }
521 if(cti[M-1]<threshold) {
522 flowimg.m[pi][yi][xi][0]=0.0;
523 if(k2file[0]) k2img.m[pi][yi][xi][0]=0.0;
524 if(vafile[0]) vaimg.m[pi][yi][xi][0]=0.0;
525 if(vdfile[0]) vdimg.m[pi][yi][xi][0]=0.0;
526 if(wssfile[0]) wssimg.m[pi][yi][xi][0]=0.0;
527 if(errfile[0]) errimg.m[pi][yi][xi][0]=0.0;
528 thresholded_nr++;
529 continue;
530 }
531
532 /* Go through all basis functions */
533 bi_min=-1; rnorm_min=1.0E80; p1=p2=p3=0.0;
534 for(bi=0; bi<bf.voiNr; bi++) {
535
536 /* Define memory site for present coefficient matrix and vector tau */
537 for(m=0; m<M; m++) {A[m]=mem[bi]+ m*N;}
538 tau=mem[bi]+M*N;
539
540 /* Get data vector */
541 for(m=0; m<M; m++) {
542 B[m]=img.m[pi][yi][xi][m];
543 /* Apply data weights */
544 B[m]*=qrweight[m];
545 }
546
547 /* Compute solution */
548 ret=qr_solve(A, M, N, tau, B, X, residual, &RNORM, wws, ws);
549 if(ret>0) { /* no solution is possible */
550 for(n=0; n<N; n++) X[n]=0.0;
551 RNORM=1.0E80;
552 }
553 /* Check if this was best fit for now */
554 if(RNORM<rnorm_min) {
555 rnorm_min=RNORM; bi_min=bi;
556 /* K1 */ p1=X[0];
557 /* Va */ if(N>1) p2=X[1]; else p2=0.0;
558 /* k2 */ p3=bf.voi[bi_min].size;
559 }
560 } /* next basis function */
561 if(rnorm_min>=1.0E60) nosolution_nr++;
562 else bf_opt_nr[bi_min]+=1;
563
564 if(verbose>6 && yi==4*img.dimy/10 && xi==4*img.dimx/10) {
565 printf(" Pixel (%d,%d,%d), K1=%g Va=%g k2=%g\n",
566 pi, yi, xi, p1, p2, p3);
567 if(verbose>10) dftPrint(&tac);
568 }
569
570 /* Calculate WSS */
571 for(m=0, wss=0.0; m<M; m++) {
572 f=p1*bf.voi[bi_min].y[m];
573 if(N>1) f+=p2*tac.voi[0].y[m];
574 f-=img.m[pi][yi][xi][m];
575 wss+=img.weight[m]*f*f;
576 }
577 /* Put results to output images */
578 if(p1>flowMax) flowimg.m[pi][yi][xi][0]=flowMax;
579 else if(p1<0.0) flowimg.m[pi][yi][xi][0]=0.0;
580 else flowimg.m[pi][yi][xi][0]=p1;
581 if(vafile[0]) vaimg.m[pi][yi][xi][0]=p2;
582 if(k2file[0]) {
583 if(p3>k2max) k2img.m[pi][yi][xi][0]=k2max;
584 else k2img.m[pi][yi][xi][0]=p3;
585 }
586 if(vdfile[0]) {
587 f=p1/p3; if(f>pWaterMax) f=pWaterMax;
588 vdimg.m[pi][yi][xi][0]=f;
589 }
590 if(wssfile[0]) wssimg.m[pi][yi][xi][0]=wss;
591 if(errfile[0]) {
592 if(bi_min==0) ret=1; else if(bi_min==bf.voiNr-1) ret=2; else ret=0;
593 errimg.m[pi][yi][xi][0]=(float)ret;
594 }
595 } /* next column */
596 } /* next row */
597 } /* next plane */
598 fitFinish=clock();
599 if(verbose>0) {fprintf(stdout, "done.\n"); fflush(stdout);}
600 if(verbose>1 || thresholded_nr>0) {
601 double f;
602 f=(double)thresholded_nr/((double)(flowimg.dimx*flowimg.dimy*flowimg.dimz));
603 f*=100.; if(f<3.0) printf("%g%%", f); else printf("%.0f%%", f);
604 printf(" of pixels were not fitted due to threshold.\n");
605 if(verbose>2) printf("thresholded %d pixels\n", thresholded_nr);
606 fflush(stdout);
607 }
608 if(verbose>1 || nosolution_nr>0) {
609 fprintf(stdout, "no QR solution for %d pixels.\n", nosolution_nr); fflush(stdout);}
610
611 /* No need for dynamic image or input tac anymore */
612 imgEmpty(&img); dftEmpty(&tac);
613 /* Free memory of QR */
614 free(chain); free(B); free(residual); free(A); free(wwschain);
615 free(wws); free(qrweight); free(mem);
616
617 /*
618 * Save basis functions if required;
619 * this is done not before, so that also the number of optimal fits
620 * achieved with each BF can be saved as the "size".
621 */
622 if(bfsfile[0]) {
623 for(bi=0; bi<bf.voiNr; bi++) sprintf(bf.voi[bi].place, "%d", bf_opt_nr[bi]);
624 if(dftWrite(&bf, bfsfile)) {
625 fprintf(stderr, "Error in writing %s: %s\n", bfsfile, dfterrmsg); fflush(stderr);
626 imgEmpty(&flowimg); imgEmpty(&k2img); imgEmpty(&vaimg); imgEmpty(&vdimg);
627 imgEmpty(&wssimg); imgEmpty(&errimg);
628 dftEmpty(&bf); free(bf_opt_nr); return(11);
629 }
630 if(verbose>0) {fprintf(stdout, "basis functions were written in %s\n", bfsfile); fflush(stdout);}
631 }
632
633 /* No need for basis functions anymore */
634 dftEmpty(&bf); free(bf_opt_nr);
635
636
637 /*
638 * Convert units to per dL and/or gram if necessary
639 */
640 if(per_dl || per_gram) {
641 /* K1 */
642 float f;
643 if(per_dl && per_gram) {
644 f=100.0/density; flowimg.unit=IMGUNIT_ML_PER_DL_PER_MIN;
645 } else if(per_dl) {
646 f=100.0; flowimg.unit=IMGUNIT_ML_PER_DL_PER_MIN;
647 } else if(per_gram) {
648 f=1.0/density; flowimg.unit=IMGUNIT_ML_PER_ML_PER_MIN;
649 } else
650 f=1.0;
651 ret=imgArithmConst(&flowimg, f, '*', 1.0E+06, verbose-6);
652 }
653 if(per_dl && vafile[0]) {
654 /* Va */
655 vaimg.unit=IMGUNIT_ML_PER_DL;
656 ret=imgArithmConst(&vaimg, 100.0, '*', 1.0E+06, verbose-6);
657 }
658
659 /*
660 * Save parametric image(s)
661 */
662 if(verbose>1) {printf("Saving parametric images\n"); fflush(stdout);}
663 {
664 ret=imgWrite(flowfile, &flowimg);
665 if(ret) {fprintf(stderr, "Error: %s\n", flowimg.statmsg); fflush(stderr);}
666 else if(verbose>0) {fprintf(stdout, "Flow image %s saved.\n", flowfile); fflush(stdout);}
667 }
668 if(ret==0 && vafile[0]) {
669 ret=imgWrite(vafile, &vaimg);
670 if(ret) {fprintf(stderr, "Error: %s\n", vaimg.statmsg); fflush(stderr);}
671 else if(verbose>0) {fprintf(stdout, "Va image %s saved.\n", vafile); fflush(stdout);}
672 }
673 if(ret==0 && vdfile[0]) {
674 ret=imgWrite(vdfile, &vdimg);
675 if(ret) {fprintf(stderr, "Error: %s\n", vdimg.statmsg); fflush(stderr);}
676 else if(verbose>0) {fprintf(stdout, "Vd image %s saved.\n", vdfile); fflush(stdout);}
677 }
678 if(ret==0 && k2file[0]) {
679 ret=imgWrite(k2file, &k2img);
680 if(ret) {fprintf(stderr, "Error: %s\n", k2img.statmsg); fflush(stderr);}
681 else if(verbose>0) {fprintf(stdout, "k2 image %s saved.\n", k2file); fflush(stdout);}
682 }
683 if(ret==0 && wssfile[0]) {
684 ret=imgWrite(wssfile, &wssimg);
685 if(ret) {fprintf(stderr, "Error: %s\n", wssimg.statmsg); fflush(stderr);}
686 else if(verbose>0) {fprintf(stdout, "WSS image %s saved.\n", wssfile); fflush(stdout);}
687 }
688 if(ret==0 && errfile[0]) {
689 ret=imgWrite(errfile, &errimg);
690 if(ret) {fprintf(stderr, "Error: %s\n", errimg.statmsg); fflush(stderr);}
691 else if(verbose>0) {fprintf(stdout, "Error image %s saved.\n", errfile); fflush(stdout);}
692 }
693
694 imgEmpty(&flowimg); imgEmpty(&k2img); imgEmpty(&vaimg); imgEmpty(&vdimg);
695 imgEmpty(&wssimg); imgEmpty(&errimg);
696
697 /* How long did the fitting take */
698 if(verbose>0) fprintf(stdout, "parameter estimation time := %.1f [s]\n",
699 (double)(fitFinish-fitStart) / CLOCKS_PER_SEC );
700
701 return(ret);
702}
703/*****************************************************************************/
704
705/*****************************************************************************/
int bfRadiowater(DFT *input, DFT *tissue, DFT *bf, int bfNr, double k2min, double k2max, char *status, int verbose)
Definition bf_model.c:77
int atof_with_check(char *double_as_string, double *result_value)
Definition decpoint.c:107
double atof_dpi(char *str)
Definition decpoint.c:59
void dftInit(DFT *data)
Definition dft.c:38
char dfterrmsg[64]
Definition dft.c:6
int dftAddmem(DFT *dft, int voiNr)
Definition dft.c:107
void dftEmpty(DFT *data)
Definition dft.c:20
void dftPrint(DFT *data)
Definition dftio.c:538
int dftWrite(DFT *data, char *filename)
Definition dftio.c:594
int dftTimeunitConversion(DFT *dft, int tunit)
Definition dftunit.c:119
int IMG_TEST
Definition img.c:6
void imgInfo(IMG *image)
Definition img.c:359
unsigned long long imgNaNs(IMG *img, int fix)
Definition img.c:658
int imgAllocateWithHeader(IMG *image, int planes, int rows, int columns, int frames, IMG *image_from)
Definition img.c:279
void imgEmpty(IMG *image)
Definition img.c:121
void imgInit(IMG *image)
Definition img.c:60
int imgArithmConst(IMG *img, float operand, char operation, float ulimit, int verbose)
Definition imgarithm.c:100
int imgWrite(const char *fname, IMG *img)
Definition imgfile.c:136
int imgReadModelingData(char *petfile, char *siffile, char *inputfile1, char *inputfile2, char *inputfile3, double *fitdur, int *fitframeNr, IMG *img, DFT *inp, DFT *iinp, int verifypeak, FILE *loginfo, int verbose, char *status)
Definition imginput.c:24
int petintegral(double *x1, double *x2, double *y, int nr, double *ie, double *iie)
Integrate PET TAC data to frame mid times.
Definition integr.c:771
Header file for libtpccurveio.
Header file for libtpcimgio.
int SIF_TEST
Definition sif.c:6
#define IMG_DC_NONCORRECTED
Header file for libtpcimgp.
Header file for libtpcmisc.
int tpcProcessStdOptions(const char *s, int *print_usage, int *print_version, int *verbose_level)
Definition proginfo.c:40
size_t strlcpy(char *dst, const char *src, size_t dstsize)
Definition strext.c:245
int tpcHtmlUsage(const char *program, char *text[], const char *path)
Definition proginfo.c:213
void tpcPrintBuild(const char *program, FILE *fp)
Definition proginfo.c:383
void tpcPrintUsage(const char *program, char *text[], FILE *fp)
Definition proginfo.c:158
Header file for libtpcmodel.
int qr_solve(double **QR, int M, int N, double *tau, double *b, double *x, double *residual, double *resNorm, double **cchain, double *chain)
Definition qr.c:492
int qr_decomp(double **a, int M, int N, double *tau, double **cchain, double *chain)
Definition qr.c:427
Header file for libtpcmodext.
int imgSetWeights(IMG *img, int wmet, int verbose)
Voi * voi
double * x1
int voiNr
double * x2
int frameNr
unsigned short int dimx
float **** m
char decayCorrection
char unit
unsigned short int dimt
int * planeNumber
float * weight
float * start
unsigned short int dimz
unsigned short int dimy
float * end
const char * statmsg
char isWeight
double size
double * y2
char voiname[MAX_REGIONSUBNAME_LEN+1]
double * y
char place[MAX_REGIONSUBNAME_LEN+1]