TPCCLIB
Loading...
Searching...
No Matches
fitk3.c
Go to the documentation of this file.
1
7/*****************************************************************************/
8#include "tpcclibConfig.h"
9/*****************************************************************************/
10#include <stdio.h>
11#include <stdlib.h>
12#include <math.h>
13#include <unistd.h>
14#include <string.h>
15#include <time.h>
16/*****************************************************************************/
17#include "libtpcmisc.h"
18#include "libtpcmodel.h"
19#include "libtpccurveio.h"
20#include "libtpcsvg.h"
21#include "libtpcmodext.h"
22/*****************************************************************************/
23
24/*****************************************************************************/
25const int parNr=4;
26DFT input, dft;
27double *petmeas, *petsim;
28double fVb=-1.0;
29double pmin[MAX_PARAMETERS], pmax[MAX_PARAMETERS];
30double fk1k2;
31int fitframeNr;
32double wss_wo_penalty=0.0;
33/*****************************************************************************/
34/* Local functions */
35double cm3Func(int parNr, double *p, void*);
36/*****************************************************************************/
37
38/*****************************************************************************/
39static char *info[] = {
40 "Non-linear fitting of two-tissue compartment model to plasma input, blood,",
41 "and tissue time-activity curves (PTAC, BTAC, and TTAC) to estimate",
42 "parameters K1, k2, k3, and Vb. Sample times must be in minutes.",
43 " ",
44 " ______ ___________________ ",
45 " | | K1 | k3 | ",
46 " | | --> | -------> | ",
47 " | Ca | <-- | C1 <------- C2 | ",
48 " | | k2 | k4=0 | ",
49 " |______| |___________________| ",
50 " ",
51 "Usage: @P [Options] ptacfile btacfile ttacfile endtime resultfile",
52 " ",
53 "Options:",
54 " -lim[=<filename>]",
55 " Specify the constraints for model parameters;",
56 " This file with default values can be created by giving this",
57 " option as the only command-line argument to this program.",
58 " Without filename the default values are printed on screen.",
59 " -SD[=<y|N>]",
60 " Standard deviations are calculated and saved in results (y),",
61 " or not calculated (N, default).",
62 " Program runs a lot faster if SD and CL are not calculated.",
63 " -CL[=<y|N>]",
64 " 95% Confidence limits are calculated and saved in results (y), or",
65 " not calculated (N, default).",
66 " -Vb=<Vb(%)>",
67 " Enter a fixed Vb; fitted by default.",
68 " -r=<Reference region id or filename>",
69 " Optional reference region is used to constrain K1/k2 in other regions;",
70 " Also k3 is fitted to reference region data, thus any large region",
71 " (for example cortex) can be used here.",
72 " -fit=<Filename>",
73 " Fitted regional TACs are written in DFT format.",
74 " -svg=<Filename>",
75 " Fitted and measured TACs are plotted in specified SVG file.",
76 " -stdoptions", // List standard options like --help, -v, etc
77 " ",
78 "Example 1: estimate K1, K1/k2, k3 and Vb",
79 " @P ua919ap.bld ua919ab.bld ua919.tac 60 ua919k3.res",
80 " ",
81 "Example 2: estimate K1 and k3; Vb is constrained to 4.5% and K1/k2 is",
82 "constrained to K1/k2 estimated from region 'occip'",
83 " @P -Vb=4.5 -r=occip ua919ap.bld ua919ab.bld ua919.tac 60 ua919.res",
84 " ",
85 "See also: patlak, lhsol, p2t_v3c, tacweigh, taccbv, rescoll",
86 " ",
87 "Keywords: TAC, modelling, irreversible uptake, k3, Ki, 2TCM",
88 0};
89/*****************************************************************************/
90
91/*****************************************************************************/
92/* Turn on the globbing of the command line, since it is disabled by default in
93 mingw-w64 (_dowildcard=0); in MinGW32 define _CRT_glob instead, if necessary;
94 In Unix&Linux wildcard command line processing is enabled by default. */
95/*
96#undef _CRT_glob
97#define _CRT_glob -1
98*/
99int _dowildcard = -1;
100/*****************************************************************************/
101
102/*****************************************************************************/
106int main(int argc, char **argv)
107{
108 int ai, help=0, version=0, verbose=1;
109 int ri, fi, pi, m, n, ret;
110 int lambda_k3=0;
111 int ref=-1, refAdded=0, inputtype;
112 char dfile[FILENAME_MAX], pfile[FILENAME_MAX], bfile[FILENAME_MAX],
113 rfile[FILENAME_MAX], ffile[FILENAME_MAX], limfile[FILENAME_MAX];
114 char svgfile[FILENAME_MAX];
115 char *cptr, refname[FILENAME_MAX], tmp[FILENAME_MAX];
116 double fitdur, wss, aic, K1, k2, k3;
117 RES res;
118 IFT ift;
119 int doBootstrap=0, doSD=0, doCL=0;
120 double *sd, *cl1, *cl2;
121 double def_pmin[MAX_PARAMETERS], def_pmax[MAX_PARAMETERS];
122 int tgoNr=0, neighNr=0, iterNr=0, fittedparNr=0;
123
124
125 /* Set parameter initial values and constraints */
126 /* K1 */ def_pmin[0]=0.0; def_pmax[0]=5.0;
127 /* K1/k2 */ def_pmin[1]=0.00001; def_pmax[1]=10.0;
128 /* k3 */ def_pmin[2]=0.0; def_pmax[2]=2.0;
129 /* Vb */ def_pmin[3]=0.0; def_pmax[3]=0.08;
130
131
132 /*
133 * Get arguments
134 */
135 if(argc==1) {tpcPrintUsage(argv[0], info, stderr); return(1);}
136 dfile[0]=pfile[0]=bfile[0]=rfile[0]=ffile[0]=refname[0]=limfile[0]=(char)0;
137 svgfile[0]=(char)0;
138 fitdur=fVb=-1.0;
139 iftInit(&ift); resInit(&res); dftInit(&dft); dftInit(&input);
140 /* Get options first, because it affects what arguments are read */
141 for(ai=1; ai<argc; ai++) if(*argv[ai]=='-') {
142 cptr=argv[ai]+1; if(*cptr=='-') cptr++; if(cptr==NULL) continue;
143 if(tpcProcessStdOptions(argv[ai], &help, &version, &verbose)==0) continue;
144 if(strncasecmp(cptr, "CL", 2)==0) {
145 if(strlen(cptr)==2) {doCL=1; continue;}
146 cptr+=2; if(*cptr=='=') {
147 cptr++;
148 if(*cptr=='Y' || *cptr=='y') {doCL=1; continue;}
149 if(*cptr=='N' || *cptr=='n') {doCL=0; continue;}
150 }
151 } else if(strncasecmp(cptr, "SD", 2)==0) {
152 if(strlen(cptr)==2) {doSD=1; continue;}
153 cptr+=2; if(*cptr=='=') {
154 cptr++;
155 if(*cptr=='Y' || *cptr=='y') {doSD=1; continue;}
156 if(*cptr=='N' || *cptr=='n') {doSD=0; continue;}
157 }
158 } else if(strncasecmp(cptr, "LIM=", 4)==0 && strlen(cptr)>4) {
159 strlcpy(limfile, cptr+4, FILENAME_MAX); continue;
160 } else if(strcasecmp(cptr, "LIM")==0) {
161 strcpy(limfile, "stdout"); continue;
162 } else if(strcasecmp(cptr, "LK3")==0) {
163 lambda_k3=1; continue;
164 } else if(strncasecmp(cptr, "R=", 2)==0 && strlen(cptr)>2) {
165 strlcpy(refname, cptr+2, FILENAME_MAX); continue;
166 } else if(strncasecmp(cptr, "Vb=", 3)==0 && strlen(cptr)>3) {
167 fVb=0.01*atof_dpi(cptr+3);
168 if(fVb>=0.0 && fVb<1.0) {
169 if(fVb<0.01) fprintf(stderr, "Warning: Vb was set to %g%%\n", 100.*fVb);
170 def_pmin[3]=def_pmax[3]=fVb;
171 continue;
172 }
173 fVb=-1.0;
174 } else if(strncasecmp(cptr, "FIT=", 4)==0) {
175 strlcpy(ffile, cptr+4, FILENAME_MAX); if(strlen(ffile)>0) continue;
176 } else if(strncasecmp(cptr, "SVG=", 4)==0) {
177 strlcpy(svgfile, cptr+4, FILENAME_MAX); if(strlen(svgfile)>0) continue;
178 }
179 fprintf(stderr, "Error: invalid option '%s'.\n", argv[ai]);
180 return(1);
181 } else break;
182
183 /* Print help or version? */
184 if(help==2) {tpcHtmlUsage(argv[0], info, ""); return(0);}
185 if(help) {tpcPrintUsage(argv[0], info, stdout); return(0);}
186 if(version) {tpcPrintBuild(argv[0], stdout); return(0);}
187
188 /* Process other arguments, starting from the first non-option */
189 for(; ai<argc; ai++) {
190 if(!pfile[0]) {
191 strlcpy(pfile, argv[ai], FILENAME_MAX); continue;
192 } else if(!bfile[0]) {
193 strlcpy(bfile, argv[ai], FILENAME_MAX); continue;
194 } else if(!dfile[0]) {
195 strlcpy(dfile, argv[ai], FILENAME_MAX); continue;
196 } else if(fitdur<0) {
197 if(!atof_with_check(argv[ai], &fitdur) && fitdur>=0.0) continue;
198 fprintf(stderr, "Error: invalid fit time '%s'.\n", argv[ai]);
199 return(1);
200 } else if(!rfile[0]) {
201 strlcpy(rfile, argv[ai], FILENAME_MAX); continue;
202 }
203 /* we should never get this far */
204 fprintf(stderr, "Error: too many arguments: '%s'.\n", argv[ai]);
205 return(1);
206 }
207 if(doSD || doCL) doBootstrap=1; else doBootstrap=0;
208
209 /* In verbose mode print arguments and options */
210 if(verbose>1) {
211 printf("pfile := %s\n", pfile);
212 printf("dfile := %s\n", dfile);
213 printf("rfile := %s\n", rfile);
214 printf("ffile := %s\n", ffile);
215 printf("svgfile := %s\n", svgfile);
216 printf("limfile := %s\n", limfile);
217 printf("refname := %s\n", refname);
218 printf("lambda_k3 := %d\n", lambda_k3);
219 printf("fitdur := %g\n", fitdur);
220 printf("doBootstrap := %d\n", doBootstrap);
221 printf("doSD := %d\n", doSD);
222 printf("doCL := %d\n", doCL);
223 }
224
225 /* If only filename for parameter constraints was given, then write one
226 with default contents, and exit */
227 if(limfile[0] && !pfile[0]) {
228 /* Check that initial value file does not exist */
229 if(strcasecmp(limfile, "stdout")!=0 && access(limfile, 0) != -1) {
230 fprintf(stderr, "Error: parameter constraint file %s exists.\n", limfile);
231 return(9);
232 }
233 if(verbose>1 && strcasecmp(limfile, "stdout")!=0)
234 printf("writing parameter constraints file\n");
235 /* Create parameter file */
236 iftPutDouble(&ift, "K1_lower", def_pmin[0], NULL, 0);
237 iftPutDouble(&ift, "K1_upper", def_pmax[0], NULL, 0);
238 iftPutDouble(&ift, "K1k2_lower", def_pmin[1], NULL, 0);
239 iftPutDouble(&ift, "K1k2_upper", def_pmax[1], NULL, 0);
240 iftPutDouble(&ift, "k3_lower", def_pmin[2], NULL, 0);
241 iftPutDouble(&ift, "k3_upper", def_pmax[2], NULL, 0);
242 iftPutDouble(&ift, "Vb_lower", def_pmin[3], NULL, 0);
243 iftPutDouble(&ift, "Vb_upper", def_pmax[3], NULL, 0);
244 ret=iftWrite(&ift, limfile, 0);
245 if(ret) {
246 fprintf(stderr, "Error in writing '%s': %s\n", limfile, ift.status);
247 iftEmpty(&ift); return(9);
248 }
249 if(strcasecmp(limfile, "stdout")!=0)
250 fprintf(stdout, "Parameter file %s with initial values written.\n", limfile);
251 iftEmpty(&ift); return(0);
252 }
253
254 /* Did we get all the information from user that we need? */
255 if(fitdur==0) fitdur=1.0E+100;
256 else if(fitdur<0) {tpcPrintUsage(argv[0], info, stderr); return(1);}
257 if(!rfile[0]) {
258 fprintf(stderr, "Error: missing command-line argument; use option --help\n");
259 return(1);
260 }
261
262 /*
263 * Read model parameter upper and lower limits
264 * if file for that was given
265 */
266 if(limfile[0]) {
267 double v;
268 if(verbose>1) printf("reading %s\n", limfile);
269 ret=iftRead(&ift, limfile, 1, 0);
270 if(ret) {
271 fprintf(stderr, "Error in reading '%s': %s\n", limfile, ift.status);
272 return(9);
273 }
274 if(verbose>10) iftWrite(&ift, "stdout", 0);
275 int n=0;
276 /* K1 */
277 if(iftGetDoubleValue(&ift, 0, "K1_lower", &v, 0)>=0) {def_pmin[0]=v; n++;}
278 if(iftGetDoubleValue(&ift, 0, "K1_upper", &v, 0)>=0) {def_pmax[0]=v; n++;}
279 /* K1/k2 */
280 if(iftGetDoubleValue(&ift, 0, "K1k2_lower", &v, 0)>=0) {def_pmin[1]=v; n++;}
281 if(iftGetDoubleValue(&ift, 0, "K1k2_upper", &v, 0)>=0) {def_pmax[1]=v; n++;}
282 /* k3 */
283 if(iftGetDoubleValue(&ift, 0, "k3_lower", &v, 0)>=0) {def_pmin[2]=v; n++;}
284 if(iftGetDoubleValue(&ift, 0, "k3_upper", &v, 0)>=0) {def_pmax[2]=v; n++;}
285 /* Vb */
286 if(iftGetDoubleValue(&ift, 0, "Vb_lower", &v, 0)>=0) {def_pmin[3]=v; n++;}
287 if(iftGetDoubleValue(&ift, 0, "Vb_upper", &v, 0)>=0) {def_pmax[3]=v; n++;}
288 iftEmpty(&ift);
289 if(n==0) {fprintf(stderr, "Error: invalid parameter file.\n"); return(9);}
290 }
291 /* Check that these are ok */
292 for(pi=n=0, ret=0; pi<parNr; pi++) {
293 if(def_pmin[pi]<0.0) ret++;
294 if(def_pmax[pi]<def_pmin[pi]) ret++;
295 if(def_pmax[pi]>def_pmin[pi]) n++;
296 }
297 if(ret) {
298 fprintf(stderr, "Error: invalid parameter constraints.\n");
299 return(9);
300 }
301 if(n==0) {
302 fprintf(stderr, "Error: no model parameters left free for fitting.\n");
303 return(9);
304 }
305
306 /* Fixed/fitted Vb */
307 if(fVb>=0.0) def_pmin[3]=def_pmax[3]=fVb;
308 if(def_pmin[3]==def_pmax[3]) fVb=def_pmin[3];
309 if(fVb==0.0) strcpy(bfile, "");
310 if(verbose>1) {
311 printf("bfile := %s\n", bfile);
312 if(fVb>=0.0) printf("fVb := %g\n", fVb);
313 }
314
315
316 /*
317 * Read tissue and input data
318 */
319 if(verbose>1) printf("reading tissue and input data\n");
320 ret=dftReadModelingData(dfile, pfile, bfile, NULL, &fitdur,
321 &fitframeNr, &dft, &input, stdout, verbose-2, tmp);
322 if(ret!=0) {
323 fprintf(stderr, "Error: %s\n", tmp);
324 return(2);
325 }
326 if(fitframeNr<4 || input.frameNr<4) {
327 fprintf(stderr, "Error: too few samples in specified fit duration.\n");
328 dftEmpty(&input); dftEmpty(&dft); return(2);
329 }
330 /* If there is no blood TAC, then create a zero blood TAC */
331 if(input.voiNr<2) {
332 if(verbose>2) printf("setting blood tac to zero\n");
333 ret=dftAddmem(&input, 1);
334 if(ret) {
335 fprintf(stderr, "Error: cannot allocate more memory.\n");
336 dftEmpty(&dft); dftEmpty(&input); return(3);
337 }
338 strcpy(input.voi[1].voiname, "blood");
339 strcpy(input.voi[1].name, input.voi[1].voiname);
340 for(fi=0; fi<input.frameNr; fi++) input.voi[1].y[fi]=0.0;
341 input.voiNr=2;
342 }
343 if(verbose>10) dftPrint(&dft);
344 if(verbose>10) dftPrint(&input);
345 /* Print the weights */
346 if(verbose>2) {
347 fprintf(stdout, "common_data_weights := %g", dft.w[0]);
348 for(fi=1; fi<dft.frameNr; fi++) fprintf(stdout, ", %g", dft.w[fi]);
349 fprintf(stdout, "\n");
350 }
351
352
353 /*
354 * Read reference TAC
355 */
356 /* Check if user even wants any reference region */
357 if(!refname[0]) {
358 if(verbose>1) printf("no reference region data\n");
359 ref=-1;
360 } else {
361 if(verbose>1) printf("reading reference region data\n");
362 if((n=dftReadReference(&dft, refname, &inputtype, &ref, tmp, verbose-3))<1) {
363 fprintf(stderr, "Error in reading '%s': %s\n", refname, tmp);
364 if(verbose>2) printf("dftReadReference()=%d\n", n);
365 dftEmpty(&dft); dftEmpty(&input); return(6);
366 }
367 if(verbose>30) dftPrint(&dft);
368 if(n>1) {
369 fprintf(stderr, "Warning: %s selected of %d reference regions.\n",
370 dft.voi[ref].name, n);
371 if(verbose>2)
372 fprintf(stdout, "selected reference region := %s\n", dft.voi[ref].name);
373 }
374 if(inputtype==5) { // Reference region name was given
375 refAdded=0; strcpy(refname, "");
376 } else { // reference file was given; ref TACs may have to be deleted later
377 refAdded=1;
378 }
379 if(verbose>15) dftPrint(&dft);
380 if(verbose>1) printf("Reference region: %s\n", dft.voi[ref].name );
381 }
382
383 /* Allocate an extra TAC for the bootstrap */
384 if(doBootstrap) {
385 ret=dftAddmem(&dft, 1); if(ret) {
386 fprintf(stderr, "Error: cannot allocate more memory.\n");
387 dftEmpty(&dft); dftEmpty(&input); return(9);
388 }
389 strcpy(dft.voi[dft.voiNr].voiname, "BS");
390 strcpy(dft.voi[dft.voiNr].name, "BS");
391 }
392 if(verbose>10) dftPrint(&dft);
393
394
395 /*
396 * Prepare the room for the results
397 */
398 if(verbose>1) printf("initializing result data\n");
399 ret=res_allocate_with_dft(&res, &dft); if(ret!=0) {
400 fprintf(stderr, "Error: cannot setup memory for results.\n");
401 dftEmpty(&input); dftEmpty(&dft); return(7);
402 }
403 /* Copy titles & filenames */
404 tpcProgramName(argv[0], 1, 1, res.program, 256);
405 strcpy(res.datafile, dfile);
406 strcpy(res.plasmafile, pfile);
407 strcpy(res.bloodfile, bfile);
408 if(ref>=0) sprintf(res.refroi, "%s", dft.voi[ref].name);
409 if(refname[0]) strcpy(res.reffile, refname);
410 strcpy(res.fitmethod, "TGO");
411 /* Constants */
412 res.isweight=dft.isweight;
413 if(fVb>=0.0) res.Vb=100.0*fVb;
414 /* Set data range */
415 sprintf(res.datarange, "%g - %g %s", 0.0, fitdur, dftTimeunit(dft.timeunit));
416 res.datanr=fitframeNr;
417 /* Set current time to results */
418 res.time=time(NULL);
419 /* Set parameter number, including also the extra "parameters"
420 and the parameter names and units */
421 res.parNr=9;
422 pi=0; strcpy(res.parname[pi], "K1"); strcpy(res.parunit[pi], "ml/(min*ml)");
423 pi++; strcpy(res.parname[pi], "K1/k2"); strcpy(res.parunit[pi], "");
424 pi++; strcpy(res.parname[pi], "k3"); strcpy(res.parunit[pi], "1/min");
425 pi++; strcpy(res.parname[pi], "Vb"); strcpy(res.parunit[pi], "%");
426 pi++; strcpy(res.parname[pi], "Ki"); strcpy(res.parunit[pi], "ml/(min*ml)");
427 pi++; strcpy(res.parname[pi], "k3*K1/k2"); strcpy(res.parunit[pi], "1/min");
428 pi++; strcpy(res.parname[pi], "k3/(k2+k3)"); strcpy(res.parunit[pi], "");
429 pi++; strcpy(res.parname[pi], "WSS"); strcpy(res.parunit[pi], "");
430 pi++; strcpy(res.parname[pi], "AIC"); strcpy(res.parunit[pi], "");
431
432
433 /*
434 * Fit at first the reference ROI(s), if required
435 */
436 if(ref>=0) for(ri=0; ri<dft.voiNr; ri++) if(dft.voi[ri].sw>0) {
437 if(verbose>0)
438 fprintf(stdout, "fitting %s as reference region\n", dft.voi[ri].name);
439 /* Initiate values */
440 petmeas=dft.voi[ri].y; petsim=dft.voi[ri].y2;
441 /* Set constraints */
442 pmin[0]=def_pmin[0]; pmax[0]=def_pmax[0]; /* K1 */
443 pmin[1]=def_pmin[1]; pmax[1]=def_pmax[1]; /* K1/k2 */
444 pmin[2]=def_pmin[2]; pmax[2]=def_pmax[2]; /* k3 */
445 pmin[3]=def_pmin[3]; pmax[3]=def_pmax[3]; /* Vb */
446 for(pi=fittedparNr=0; pi<parNr; pi++) if(pmax[pi]>pmin[pi]) fittedparNr++;
447 if(verbose>3) {
448 printf(" ref_constraints :=");
449 for(pi=0; pi<parNr; pi++) printf(" [%g,%g]", pmin[pi], pmax[pi]);
450 printf("\n");
451 printf("fittedparNr := %d\n", fittedparNr);
452 }
453 /* Fit */
456 // neighNr=5; tgoNr=300; iterNr=0;
457 tgoNr=50+25*fittedparNr;
458 neighNr=6*fittedparNr;
459 ret=tgo(
460 pmin, pmax, cm3Func, NULL, parNr, neighNr,
461 &wss, res.voi[ri].parameter, tgoNr, iterNr, verbose-8);
462 if(ret>0) {
463 fprintf(stderr, "Error in optimization (%d).\n", ret);
464 dftEmpty(&input); dftEmpty(&dft); resEmpty(&res); return(8);
465 }
466 /* Correct fitted parameters to match constraints like inside the function */
467 (void)modelCheckParameters(parNr, pmin, pmax, res.voi[ri].parameter,
468 res.voi[ri].parameter, NULL);
469 wss=wss_wo_penalty;
470 /* Fix K1/k2 to non-reference regions, but only after all reference
471 regions are fitted */
472 if(ri==ref) {
473 fk1k2=res.voi[ri].parameter[1];
474 if(verbose>2) {fprintf(stdout, " K1/k2 := %g\n", fk1k2);}
475 }
476 /* Bootstrap */
477 if(doBootstrap) {
478 if(verbose>2) printf(" bootstrapping\n");
479 /* bootstrap changes measured and simulated data, therefore use copies */
480 petmeas=dft.voi[dft.voiNr].y2; petsim=dft.voi[dft.voiNr].y3;
481 if(doSD) sd=res.voi[ri].sd; else sd=NULL;
482 if(doCL) {cl1=res.voi[ri].cl1; cl2=res.voi[ri].cl2;} else cl1=cl2=NULL;
483 ret=bootstrap(
484 0, cl1, cl2, sd,
485 res.voi[ri].parameter, pmin, pmax, fitframeNr,
486 // measured original TAC, not modified
487 dft.voi[ri].y,
488 // fitted TAC, not modified
489 dft.voi[ri].y2,
490 // tissue TAC noisy data is written to be used by objf
491 petmeas,
492 parNr, dft.w, cm3Func, tmp, verbose-4
493 );
494 if(ret) {
495 fprintf(stderr, "Error in bootstrap: %s\n", tmp);
496 for(pi=0; pi<parNr; pi++) {
497 if(doSD) sd[pi]=nan("");
498 if(doCL) cl1[pi]=cl2[pi]=nan("");
499 }
500 }
501 }
502 /* Calculate AIC based on nr of parameters that actually are fitted */
503 for(pi=n=0; pi<parNr; pi++) if(pmax[pi]>pmin[pi]) n++;
504 if(verbose>2) printf("nr_of_fitted_parameters := %d\n", n);
505 for(fi=m=0; fi<fitframeNr; fi++) if(dft.w[fi]>0.0) m++;
506 if(verbose>2) printf("nr_of_fitted_samples := %d\n", m);
507 aic=aicSS(wss, m, n);
508
509 /* Set results wss and aic */
510 res.voi[ri].parameter[res.parNr-2]=wss;
511 res.voi[ri].parameter[res.parNr-1]=aic;
512
513 /* Convert Vb fraction to percentage */
514 res.voi[ri].parameter[3]*=100.;
515 if(doSD) res.voi[ri].sd[3]*=100.;
516 if(doCL) {res.voi[ri].cl1[3]*=100.; res.voi[ri].cl2[3]*=100.;}
517 /* Calculate Ki, lambda*k3, and k3/(k2+k3) */
518 K1=res.voi[ri].parameter[0];
519 k2=K1/res.voi[ri].parameter[1];
520 k3=res.voi[ri].parameter[2];
521 /* Ki */
522 res.voi[ri].parameter[4]=K1*k3/(k2+k3);
523 /* lambda*k3 */
524 res.voi[ri].parameter[5]=res.voi[ri].parameter[1]*res.voi[ri].parameter[2];
525 /* k3/(k2+k3) */
526 res.voi[ri].parameter[6]=k3/(k2+k3);
527
528 } // next reference region
529
530
531 /*
532 * Fit other than reference regions
533 */
534 if(verbose>0) {fprintf(stdout, "fitting regional TACs: "); fflush(stdout);}
535 if(verbose>1) printf("\n");
536 for(ri=0; ri<dft.voiNr; ri++) if(dft.voi[ri].sw==0) {
537 if(verbose>2) printf("\n %d %s:\n", ri, dft.voi[ri].name);
538
539 /* Initiate values */
540 petmeas=dft.voi[ri].y; petsim=dft.voi[ri].y2;
541
542 /* Set constraints */
543 pmin[0]=def_pmin[0]; pmax[0]=def_pmax[0]; /* K1 */
544 pmin[1]=def_pmin[1]; pmax[1]=def_pmax[1]; /* K1/k2 */
545 if(ref>=0) {
546 pmin[1]=pmax[1]=fk1k2; /* K1/k2 */
547 }
548 pmin[2]=def_pmin[2]; pmax[2]=def_pmax[2]; /* k3 */
549 pmin[3]=def_pmin[3]; pmax[3]=def_pmax[3]; /* Vb */
550 for(pi=fittedparNr=0; pi<parNr; pi++) if(pmax[pi]>pmin[pi]) fittedparNr++;
551 if(verbose>3) {
552 printf(" constraints :=");
553 for(pi=0; pi<parNr; pi++) printf(" [%g,%g]", pmin[pi], pmax[pi]);
554 printf("\n");
555 printf("fittedparNr := %d\n", fittedparNr);
556 }
557
558 /* Fit */
561 // tgoNr=300; iterNr=0; neighNr=5;
562 tgoNr=50+25*fittedparNr;
563 neighNr=6*fittedparNr;
564 iterNr=0;
565 ret=tgo(
566 pmin, pmax, cm3Func, NULL, parNr, 5,
567 &wss, res.voi[ri].parameter, 300, 0, verbose-8);
568 if(ret>0) {
569 fprintf(stderr, "\nError in optimization (%d).\n", ret);
570 dftEmpty(&input); dftEmpty(&dft); resEmpty(&res); return(8);
571 }
572 /* Correct fitted parameters to match constraints like inside the function */
573 (void)modelCheckParameters(parNr, pmin, pmax, res.voi[ri].parameter,
574 res.voi[ri].parameter, NULL);
575 wss=wss_wo_penalty;
576
577 /* Bootstrap */
578 if(doBootstrap) {
579 if(verbose>2) printf(" bootstrapping\n");
580 /* bootstrap changes measured and simulated data, therefore use copies */
581 petmeas=dft.voi[dft.voiNr].y2; petsim=dft.voi[dft.voiNr].y3;
582 if(doSD) sd=res.voi[ri].sd; else sd=NULL;
583 if(doCL) {cl1=res.voi[ri].cl1; cl2=res.voi[ri].cl2;} else cl1=cl2=NULL;
584 ret=bootstrap(
585 0, cl1, cl2, sd,
586 res.voi[ri].parameter, pmin, pmax, fitframeNr,
587 // measured original TAC, not modified
588 dft.voi[ri].y,
589 // fitted TAC, not modified
590 dft.voi[ri].y2,
591 // tissue TAC noisy data is written to be used by objf
592 petmeas,
593 parNr, dft.w, cm3Func, tmp, verbose-4
594 );
595 if(ret) {
596 fprintf(stderr, "Error in bootstrap: %s\n", tmp);
597 for(pi=0; pi<parNr; pi++) {
598 if(doSD) sd[pi]=nan("");
599 if(doCL) cl1[pi]=cl2[pi]=nan("");
600 }
601 }
602 }
603
604 /* Calculate AIC, based on nr of parameters that actually are fitted */
605 for(pi=n=0; pi<parNr; pi++) if(pmax[pi]>pmin[pi]) n++;
606 if(verbose>2) printf("nr_of_fitted_parameters := %d\n", n);
607 for(fi=m=0; fi<fitframeNr; fi++) if(dft.w[fi]>0.0) m++;
608 if(verbose>2) printf("nr_of_fitted_samples := %d\n", m);
609 aic=aicSS(wss, m, n);
610
611 /* Set results wss and aic */
612 res.voi[ri].parameter[res.parNr-2]=wss;
613 res.voi[ri].parameter[res.parNr-1]=aic;
614
615 /* Convert Vb fraction to percentage */
616 res.voi[ri].parameter[3]*=100.;
617 if(doSD) res.voi[ri].sd[3]*=100.;
618 if(doCL) {res.voi[ri].cl1[3]*=100.; res.voi[ri].cl2[3]*=100.;}
619 /* Calculate Ki, lambda*k3, and k3/(k2+k3) */
620 K1=res.voi[ri].parameter[0];
621 k2=K1/res.voi[ri].parameter[1];
622 k3=res.voi[ri].parameter[2];
623 /* Ki */
624 res.voi[ri].parameter[4]=K1*k3/(k2+k3);
625 /* lambda*k3 */
626 res.voi[ri].parameter[5]=res.voi[ri].parameter[1]*res.voi[ri].parameter[2];
627 /* k3/(k2+k3) */
628 res.voi[ri].parameter[6]=k3/(k2+k3);
629
630 /* done with this region */
631 if(dft.voiNr>2 && verbose==1) {fprintf(stdout, "."); fflush(stdout);}
632
633 } /* next region */
634 if(verbose>0) {fprintf(stdout, "\n"); fflush(stdout);}
635
636
637 /*
638 * Print results on screen
639 */
640 if(verbose>0) {resPrint(&res); fprintf(stdout, "\n");}
641
642
643 /*
644 * Save results
645 */
646 if(verbose>1) printf("saving results\n");
647 ret=resWrite(&res, rfile, verbose-3);
648 if(ret) {
649 fprintf(stderr, "Error in writing '%s': %s\n", rfile, reserrmsg);
650 dftEmpty(&dft); dftEmpty(&input); resEmpty(&res);
651 return(11);
652 }
653 if(verbose>0) fprintf(stdout, "Model parameters written in %s\n", rfile);
654
655 /*
656 * Saving and/or plotting of fitted TACs
657 */
658 if(svgfile[0] || ffile[0]) {
659
660 /* Create a DFT containing fitted TACs */
661 char tmp[64];
662 DFT dft2;
663 dftInit(&dft2); ret=dftdup(&dft, &dft2);
664 if(ret) {
665 fprintf(stderr, "Error: cannot save fitted curves.\n");
666 dftEmpty(&dft); dftEmpty(&input); resEmpty(&res);
667 return(21);
668 }
669 for(ri=0; ri<dft.voiNr; ri++) for(fi=0; fi<fitframeNr; fi++)
670 dft2.voi[ri].y[fi]=dft2.voi[ri].y2[fi];
671 dft2.frameNr=fitframeNr;
672
673 /* Save SVG plot of fitted and original data */
674 if(svgfile[0]) {
675 if(verbose>1) printf("saving SVG plot\n");
676 sprintf(tmp, "K1-k3 fit: ");
677 if(strlen(dft.studynr)>0) strcat(tmp, dft.studynr);
678 ret=plot_fitrange_svg(&dft, &dft2, tmp, 0.0, 1.02*dft.x[fitframeNr-1],
679 0.0, nan(""), svgfile, verbose-8);
680 if(ret) {
681 fprintf(stderr, "Error (%d) in writing '%s'.\n", ret, svgfile);
682 dftEmpty(&dft2); dftEmpty(&dft); dftEmpty(&input); resEmpty(&res);
683 return(30+ret);
684 }
685 if(verbose>0) printf("Plots written in %s\n", svgfile);
686 }
687
688 /* Delete reference region(s) from the data */
689 if(refAdded!=0) {
690 for(ri=dft2.voiNr-1; ri>=0; ri--) if(dft2.voi[ri].sw!=0)
691 dftDelete(&dft2, ri);
692 }
693
694 /* Save fitted TACs */
695 if(ffile[0]) {
696 if(verbose>1) printf("saving fitted curves\n");
697 tpcProgramName(argv[0], 1, 0, tmp, 128);
698 sprintf(dft2.comments, "# program := %s\n", tmp);
699 if(dftWrite(&dft2, ffile)) {
700 fprintf(stderr, "Error in writing '%s': %s\n", ffile, dfterrmsg);
701 dftEmpty(&dft2); dftEmpty(&dft); dftEmpty(&input); resEmpty(&res);
702 return(22);
703 }
704 if(verbose>0) printf("Fitted TACs written in %s\n", ffile);
705 }
706
707 dftEmpty(&dft2);
708 }
709
710 dftEmpty(&dft); dftEmpty(&input); resEmpty(&res);
711 return(0);
712}
713/*****************************************************************************/
714
715/*****************************************************************************
716 *
717 * Functions to be minimized
718 *
719 *****************************************************************************/
720double cm3Func(int parNr, double *p, void *fdata)
721{
722 int fi, ret;
723 double Vb, k2, k3, d, wss=0.0;
724 double pa[MAX_PARAMETERS], penalty=1.0;
725
726
727 /* Check parameters against the constraints */
728 ret=modelCheckParameters(parNr, pmin, pmax, p, pa, &penalty);
729 if(fdata) {}
730 /* Calculate k2 and k3 */
731 k2=pa[0]/pa[1]; k3=pa[2]; if(fVb>=0.0) Vb=fVb; else Vb=pa[3];
732
733 /* Simulate the tissue PET TAC */
734 ret=simC3vs(
735 input.x, input.voi[0].y, input.voi[1].y, input.frameNr,
736 pa[0], k2, k3, 0.0, 0.0, 0.0, 0.0, Vb, 1.0,
737 input.voi[0].y2, NULL, NULL, NULL, NULL, NULL);
738 if(ret) {
739 printf("error %d in simulation\n", ret);
740 return(nan(""));
741 }
742
743 /* Interpolate & integrate to measured PET frames */
745 ret=interpolate4pet(
746 input.x, input.voi[0].y2, input.frameNr,
747 dft.x1, dft.x2, petsim, NULL, NULL, fitframeNr);
748 else
749 ret=interpolate(
750 input.x, input.voi[0].y2, input.frameNr,
751 dft.x, petsim, NULL, NULL, fitframeNr);
752 if(ret) {
753 printf("error %d in interpolation\n", ret);
754 return(nan(""));
755 }
756
757 /* Calculate error */
758 for(fi=0, wss=0.0; fi<fitframeNr; fi++) if(dft.w[fi]>0.0) {
759 d=petmeas[fi]-petsim[fi];
760 wss+=dft.w[fi]*d*d;
761 }
762 wss_wo_penalty=wss;
763 wss*=penalty;
764 if(0) printf("K1=%g k2=%g k3=%g Vb=%g => %g\n",
765 pa[0], k2, pa[2], Vb, wss);
766
767 return(wss);
768}
769/*****************************************************************************/
770
771/*****************************************************************************/
double aicSS(double ss, const int n, const int k)
Definition aic.c:20
int bootstrap(int iterNr, double *cLim1, double *cLim2, double *SD, double *parameter, double *lowlim, double *uplim, int frameNr, double *origTac, double *fitTac, double *bsTac, int parNr, double *weight, double(*objf)(int, double *, void *), char *status, int verbose)
Definition bootstrap.c:55
int modelCheckParameters(int par_nr, double *lower_p, double *upper_p, double *test_p, double *accept_p, double *penalty)
Definition constraints.c:15
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
int dftdup(DFT *dft1, DFT *dft2)
Definition dft.c:655
char dfterrmsg[64]
Definition dft.c:6
int dftAddmem(DFT *dft, int voiNr)
Definition dft.c:107
int dftDelete(DFT *dft, int voi)
Definition dft.c:538
void dftEmpty(DFT *data)
Definition dft.c:20
int dftReadReference(DFT *tissue, char *filename, int *filetype, int *ref_index, char *status, int verbose)
Definition dftinput.c:204
int dftReadModelingData(char *tissuefile, char *inputfile1, char *inputfile2, char *inputfile3, double *fitdur, int *fitframeNr, DFT *tis, DFT *inp, FILE *loginfo, int verbose, char *status)
Definition dftinput.c:342
void dftPrint(DFT *data)
Definition dftio.c:538
int dftWrite(DFT *data, char *filename)
Definition dftio.c:594
int res_allocate_with_dft(RES *res, DFT *dft)
Definition dftres.c:14
int iftPutDouble(IFT *ift, char *key, double value, char *cmt_type, int verbose)
Definition ift.c:145
void iftEmpty(IFT *ift)
Definition ift.c:60
void iftInit(IFT *ift)
Definition ift.c:45
int iftRead(IFT *ift, char *filename, int is_key_required, int verbose)
Definition iftfile.c:24
int iftWrite(IFT *ift, char *filename, int verbose)
Definition iftfile.c:282
int iftGetDoubleValue(IFT *ift, int si, const char *key, double *value, int verbose)
Definition iftsrch.c:268
int interpolate(double *x, double *y, int nr, double *newx, double *newy, double *newyi, double *newyii, int newnr)
Linear interpolation and integration.
Definition integr.c:28
int interpolate4pet(double *x, double *y, int nr, double *newx1, double *newx2, double *newy, double *newyi, double *newyii, int newnr)
Interpolate and integrate TAC to PET frames.
Definition integr.c:510
Header file for libtpccurveio.
char reserrmsg[64]
Definition result.c:6
void resInit(RES *res)
Definition result.c:52
void resPrint(RES *res)
Definition result.c:186
int resWrite(RES *res, char *filename, int verbose)
Definition result.c:565
#define DFT_TIME_STARTEND
void resEmpty(RES *res)
Definition result.c:22
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
void tpcProgramName(const char *program, int version, int copyright, char *prname, int n)
Definition proginfo.c:101
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 simC3vs(double *t, double *ca, double *cb, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double f, double vb, double fa, double *cpet, double *cta, double *ctb, double *ctc, double *ctab, double *ctvb)
Definition simulate.c:243
int tgo(double *lowlim, double *uplim, double(*objf)(int, double *, void *), void *objfData, int dim, int neighNr, double *fmin, double *gmin, int samNr, int tgoNr, int verbose)
Definition tgo.c:39
int TGO_LOCAL_INSIDE
Definition tgo.c:29
#define MAX_PARAMETERS
Definition libtpcmodel.h:31
int TGO_SQUARED_TRANSF
Definition tgo.c:27
Header file for libtpcmodext.
int plot_fitrange_svg(DFT *dft1, DFT *dft2, char *main_title, double x1, double x2, double y1, double y2, char *fname, int verbose)
Definition plotfit.c:16
Header file for libtpcsvg.
int timetype
Voi * voi
int timeunit
char studynr[MAX_STUDYNR_LEN+1]
double * w
double * x1
char comments[_DFT_COMMENT_LEN+1]
int voiNr
double * x2
int frameNr
int isweight
double * x
const char * status
Definition libtpcmisc.h:277
int parNr
char parname[MAX_RESPARAMS][MAX_RESPARNAME_LEN+1]
int datanr
ResVOI * voi
double Vb
char program[1024]
char plasmafile[FILENAME_MAX]
char datarange[128]
int isweight
char datafile[FILENAME_MAX]
char reffile[FILENAME_MAX]
char fitmethod[128]
char refroi[64]
char parunit[MAX_RESPARAMS][MAX_RESPARNAME_LEN+1]
char bloodfile[FILENAME_MAX]
time_t time
double parameter[MAX_RESPARAMS]
double cl2[MAX_RESPARAMS]
double cl1[MAX_RESPARAMS]
double sd[MAX_RESPARAMS]
double * y2
char voiname[MAX_REGIONSUBNAME_LEN+1]
char sw
double * y
char name[MAX_REGIONNAME_LEN+1]
double * y3