55 double b, c, d, w, z, dt2;
56 double cai, ca_last, t_last;
57 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
58 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
63 if(ct==NULL)
return 2;
67 if(k3<=0.0) {k3=0.0;
if(k2<=0.0) {k2=0.0;}}
68 else if(k5<=0.0) {k5=0.0;
if(k4<=0.0) {k4=0.0;}}
69 else if(k6<=0.0) {k6=0.0;}
72 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
74 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=0.0;
75 ct1=ct2=ct3=ct1i=ct2i=ct3i=0.0;
76 for(
int i=0; i<nr; i++) {
78 dt2=0.5*(t[i]-t_last);
84 cai+=(ca[i]+ca_last)*dt2;
86 b=ct1i_last+dt2*ct1_last;
87 c=ct2i_last+dt2*ct2_last;
88 d=ct3i_last+dt2*ct3_last;
89 w=k4 + k5 - (k5*k6*dt2)/(1.0+k6*dt2);
93 + k1*z*cai + (k3*k4*dt2 - (k2+k3)*z)*b
94 + k4*c + k4*k6*dt2*d/(1.0+k6*dt2)
95 ) / ( z*(1.0 + dt2*(k2+k3)) - k3*k4*dt2*dt2 );
96 ct1i = ct1i_last + dt2*(ct1_last+ct1);
98 ct2 = (k3*ct1i - w*c + k6*d/(1.0+k6*dt2)) / z;
99 ct2i = ct2i_last + dt2*(ct2_last+ct2);
101 ct3 = (k5*ct2i - k6*d) / (1.0 + k6*dt2);
102 ct3i = ct3i_last + dt2*(ct3_last+ct3);
105 ct[i]=ct1+ct2+ct3;
if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
106 if(cta!=NULL) {cta[i]=ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
107 if(ctb!=NULL) {ctb[i]=ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
108 if(ctc!=NULL) {ctc[i]=ct3;
if(fabs(ctc[i])<1.0e-12) ctc[i]=0.0;}
110 t_last=t[i]; ca_last=ca[i];
111 ct1_last=ct1; ct1i_last=ct1i;
112 ct2_last=ct2; ct2i_last=ct2i;
113 ct3_last=ct3; ct3i_last=ct3i;
165 double dt2, r, s, u, v, w;
166 double cai, ca_last, t_last;
167 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
168 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
173 if(ct==NULL)
return 2;
179 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
181 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=0.0;
182 ct1=ct2=ct3=ct1i=ct2i=ct3i=0.0;
183 for(i=0; i<nr; i++) {
185 dt2=0.5*(t[i]-t_last);
191 cai+=(ca[i]+ca_last)*dt2;
195 u=ct1i_last+dt2*ct1_last;
196 v=ct2i_last+dt2*ct2_last;
197 w=ct3i_last+dt2*ct3_last;
199 ct1 = ( k1*cai - (k2 + (k3/r) + (k5/s))*u + (k4/r)*v + (k6/s)*w )
200 / ( 1.0 + dt2*(k2 + (k3/r) + (k5/s)) );
201 ct1i = ct1i_last + dt2*(ct1_last+ct1);
203 ct2 = (k3*ct1i - k4*v) / r;
204 ct2i = ct2i_last + dt2*(ct2_last+ct2);
206 ct3 = (k5*ct1i - k6*w) / s;
207 ct3i = ct3i_last + dt2*(ct3_last+ct3);
210 ct[i]=ct1+ct2+ct3;
if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
211 if(cta!=NULL) {cta[i]=ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
212 if(ctb!=NULL) {ctb[i]=ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
213 if(ctc!=NULL) {ctc[i]=ct3;
if(fabs(ctc[i])<1.0e-12) ctc[i]=0.0;}
215 t_last=t[i]; ca_last=ca[i];
216 ct1_last=ct1; ct1i_last=ct1i;
217 ct2_last=ct2; ct2i_last=ct2i;
218 ct3_last=ct3; ct3i_last=ct3i;
284 double b, c, d, w, z, dt2, va, vv;
285 double cai, ca_last, t_last, dct, cvb;
286 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
287 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
292 if(cpet==NULL)
return 2;
296 if(vb<0.0 || vb>=1.0)
return 4;
297 if(fa<=0.0 || fa>1.0)
return 5;
298 va=fa*vb; vv=(1.0-fa)*vb;
301 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
303 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=0.0;
304 ct1=ct2=ct3=ct1i=ct2i=ct3i=0.0;
305 for(i=0; i<nr; i++) {
307 dt2=0.5*(t[i]-t_last);
313 cai+=(ca[i]+ca_last)*dt2;
315 b=ct1i_last+dt2*ct1_last;
316 c=ct2i_last+dt2*ct2_last;
317 d=ct3i_last+dt2*ct3_last;
318 w=k4 + k5 - (k5*k6*dt2)/(1.0+k6*dt2);
322 + k1*z*cai + (k3*k4*dt2 - (k2+k3)*z)*b
323 + k4*c + k4*k6*dt2*d/(1.0+k6*dt2)
324 ) / ( z*(1.0 + dt2*(k2+k3)) - k3*k4*dt2*dt2 );
325 ct1i = ct1i_last + dt2*(ct1_last+ct1);
327 ct2 = (k3*ct1i - w*c + k6*d/(1.0+k6*dt2)) / z;
328 ct2i = ct2i_last + dt2*(ct2_last+ct2);
330 ct3 = (k5*ct2i - k6*d) / (1.0 + k6*dt2);
331 ct3i = ct3i_last + dt2*(ct3_last+ct3);
334 if(f>0.) {dct = k1*ca[i] - k2*ct1; cvb = cb[i] - dct/f;}
else cvb=cb[i];
336 cpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2+ct3);
337 if(fabs(cpet[i])<1.0e-12) cpet[i]=0.0;
338 if(cta!=NULL) {cta[i]=(1.0-vb)*ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
339 if(ctb!=NULL) {ctb[i]=(1.0-vb)*ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
340 if(ctc!=NULL) {ctc[i]=(1.0-vb)*ct3;
if(fabs(ctc[i])<1.0e-12) ctc[i]=0.0;}
341 if(ctab!=NULL) {ctab[i]=va*cb[i];
if(fabs(ctab[i])<1.0e-12) ctab[i]=0.0;}
342 if(ctvb!=NULL) {ctvb[i]=vv*cvb;
if(fabs(ctvb[i])<1.0e-12) ctvb[i]=0.0;}
344 t_last=t[i]; ca_last=ca[i];
345 ct1_last=ct1; ct1i_last=ct1i;
346 ct2_last=ct2; ct2i_last=ct2i;
347 ct3_last=ct3; ct3i_last=ct3i;
414 double dt2, r, s, u, v, w, va, vv;
415 double cai, ca_last, t_last, dct, cvb;
416 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
417 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
422 if(cpet==NULL)
return 2;
426 if(vb<0.0 || vb>=1.0)
return 4;
427 if(fa<=0.0 || fa>1.0)
return 5;
428 va=fa*vb; vv=(1.0-fa)*vb;
431 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
433 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=0.0;
434 ct1=ct2=ct3=ct1i=ct2i=ct3i=0.0;
435 for(i=0; i<nr; i++) {
437 dt2=0.5*(t[i]-t_last);
443 cai+=(ca[i]+ca_last)*dt2;
447 u=ct1i_last+dt2*ct1_last;
448 v=ct2i_last+dt2*ct2_last;
449 w=ct3i_last+dt2*ct3_last;
451 ct1 = ( k1*cai - (k2 + (k3/r) + (k5/s))*u + (k4/r)*v + (k6/s)*w )
452 / ( 1.0 + dt2*(k2 + (k3/r) + (k5/s)) );
453 ct1i = ct1i_last + dt2*(ct1_last+ct1);
455 ct2 = (k3*ct1i - k4*v) / r;
456 ct2i = ct2i_last + dt2*(ct2_last+ct2);
458 ct3 = (k5*ct1i - k6*w) / s;
459 ct3i = ct3i_last + dt2*(ct3_last+ct3);
462 if(f>0.) {dct = k1*ca[i] - k2*ct1; cvb = cb[i] - dct/f;}
else cvb=cb[i];
464 cpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2+ct3);
465 if(fabs(cpet[i])<1.0e-12) cpet[i]=0.0;
466 if(cta!=NULL) {cta[i]=(1.0-vb)*ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
467 if(ctb!=NULL) {ctb[i]=(1.0-vb)*ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
468 if(ctc!=NULL) {ctc[i]=(1.0-vb)*ct3;
if(fabs(ctc[i])<1.0e-12) ctc[i]=0.0;}
469 if(ctab!=NULL) {ctab[i]=va*cb[i];
if(fabs(ctab[i])<1.0e-12) ctab[i]=0.0;}
470 if(ctvb!=NULL) {ctvb[i]=vv*cvb;
if(fabs(ctvb[i])<1.0e-12) ctvb[i]=0.0;}
472 t_last=t[i]; ca_last=ca[i];
473 ct1_last=ct1; ct1i_last=ct1i;
474 ct2_last=ct2; ct2i_last=ct2i;
475 ct3_last=ct3; ct3i_last=ct3i;
524 double cai, ca_last, t_last;
525 double ct1, ct1_last, ct2, ct2_last;
526 double ct1i, ct1i_last, ct2i, ct2i_last;
531 if(ct==NULL)
return 2;
537 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
539 ct1_last=ct2_last=ct1i_last=ct2i_last=ct1=ct2=ct1i=ct2i=0.0;
540 for(i=0; i<nr; i++) {
542 dt2=0.5*(t[i]-t_last);
548 cai+=(ca[i]+ca_last)*dt2;
550 b=ct1i_last+dt2*ct1_last;
551 c=ct2i_last+dt2*ct2_last;
553 ct1 = (k1*cai - (k2+k3)*b) / (1.0 + (k2+k3)*dt2 );
554 ct1i = ct1i_last + dt2*(ct1_last+ct1);
556 ct2 = (k3*ct1i - kLoss*c) / (1.0 + kLoss*dt2);
557 ct2i = ct2i_last + dt2*(ct2_last+ct2);
560 ct[i]=ct1+ct2;
if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
561 if(cta!=NULL) {cta[i]=ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
562 if(ctb!=NULL) {ctb[i]=ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
564 t_last=t[i]; ca_last=ca[i];
565 ct1_last=ct1; ct1i_last=ct1i;
566 ct2_last=ct2; ct2i_last=ct2i;
627 double dt2, b, c, va, vv;
628 double cai, ca_last, t_last, dct, cvb;
629 double ct1, ct1_last, ct2, ct2_last;
630 double ct1i, ct1i_last, ct2i, ct2i_last;
635 if(cpet==NULL)
return 2;
639 if(vb<0.0 || vb>=1.0)
return 4;
640 if(fa<=0.0 || fa>1.0)
return 5;
641 va=fa*vb; vv=(1.0-fa)*vb;
644 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
646 ct1_last=ct2_last=ct1i_last=ct2i_last=ct1=ct2=ct1i=ct2i=0.0;
647 for(i=0; i<nr; i++) {
649 dt2=0.5*(t[i]-t_last);
655 cai+=(ca[i]+ca_last)*dt2;
657 b=ct1i_last+dt2*ct1_last;
658 c=ct2i_last+dt2*ct2_last;
660 ct1 = (k1*cai - (k2+k3)*b) / (1.0 + (k2+k3)*dt2 );
661 ct1i = ct1i_last + dt2*(ct1_last+ct1);
663 ct2 = (k3*ct1i - kL*c) / (1.0 + kL*dt2);
664 ct2i = ct2i_last + dt2*(ct2_last+ct2);
667 if(f>0.) {dct = k1*ca[i] - k2*ct1 - kL*ct2; cvb = cb[i] - dct/f;}
670 cpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2);
671 if(fabs(cpet[i])<1.0e-12) cpet[i]=0.0;
672 if(cta!=NULL) {cta[i]=(1.0-vb)*ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
673 if(ctb!=NULL) {ctb[i]=(1.0-vb)*ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
674 if(ctab!=NULL) {ctab[i]=va*cb[i];
if(fabs(ctab[i])<1.0e-12) ctab[i]=0.0;}
675 if(ctvb!=NULL) {ctvb[i]=vv*cvb;
if(fabs(ctvb[i])<1.0e-12) ctvb[i]=0.0;}
677 t_last=t[i]; ca_last=ca[i];
678 ct1_last=ct1; ct1i_last=ct1i;
679 ct2_last=ct2; ct2i_last=ct2i;
750 double dt2, b, c, d, w, z, u, va, vv;
751 double cai, ca_last, t_last, dct, cvb;
752 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
753 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
758 if(cpet==NULL)
return 2;
762 if(vb<0.0 || vb>=1.0)
return 4;
763 if(fa<=0.0 || fa>1.0)
return 5;
764 va=fa*vb; vv=(1.0-fa)*vb;
767 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
769 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=0.0;
770 ct1=ct2=ct3=ct1i=ct2i=ct3i=0.0;
771 for(i=0; i<nr; i++) {
773 dt2=0.5*(t[i]-t_last);
779 cai+=(ca[i]+ca_last)*dt2;
782 z = 1.0 + dt2*(k6 + kLoss);
783 u = k2 + k3 + k5 - k3*k4*dt2/w - k5*k6*dt2/z;
784 b = ct1i_last+dt2*ct1_last;
785 c = ct2i_last+dt2*ct2_last;
786 d = ct3i_last+dt2*ct3_last;
788 ct1 = ( k1*cai - u*b + k4*c/w + k6*d/z ) / ( 1.0 + dt2*u);
789 ct1i = ct1i_last + dt2*(ct1_last+ct1);
791 ct2 = (k3*ct1i - k4*c) / w;
792 ct2i = ct2i_last + dt2*(ct2_last+ct2);
794 ct3 = (k5*ct1i - (k6 + kLoss)*d) / z;
795 ct3i = ct3i_last + dt2*(ct3_last+ct3);
798 if(f>0.) {dct = k1*ca[i] - k2*ct1 - kLoss*ct3; cvb = cb[i] - dct/f;}
801 cpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2+ct3);
802 if(fabs(cpet[i])<1.0e-12) cpet[i]=0.0;
803 if(cta!=NULL) {cta[i]=(1.0-vb)*ct1;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
804 if(ctb!=NULL) {ctb[i]=(1.0-vb)*ct2;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
805 if(ctc!=NULL) {ctc[i]=(1.0-vb)*ct3;
if(fabs(ctc[i])<1.0e-12) ctc[i]=0.0;}
806 if(ctab!=NULL) {ctab[i]=va*cb[i];
if(fabs(ctab[i])<1.0e-12) ctab[i]=0.0;}
807 if(ctvb!=NULL) {ctvb[i]=vv*cvb;
if(fabs(ctvb[i])<1.0e-12) ctvb[i]=0.0;}
809 t_last=t[i]; ca_last=ca[i];
810 ct1_last=ct1; ct1i_last=ct1i;
811 ct2_last=ct2; ct2i_last=ct2i;
812 ct3_last=ct3; ct3i_last=ct3i;
859 double cri, cr_last, t_last;
860 double cf, cf_last, cb, cb_last;
861 double cfi, cfi_last, cbi, cbi_last;
866 if(ct==NULL)
return 2;
869 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
870 cri=cr_last=0.0; cf_last=cb_last=cfi_last=cbi_last=cf=cb=cfi=cbi=0.0;
871 for(i=0; i<nr; i++) {
873 dt2=0.5*(t[i]-t_last);
879 cri+=(cr[i]+cr_last)*dt2;
881 f=cfi_last+dt2*cf_last;
882 b=cbi_last+dt2*cb_last;
883 w=k2 + k3 + k2*k4*dt2;
885 cf = ( (1.0 + k4*dt2)*(R1*cr[i] + k2*cri) + k4*b - w*f ) /
886 ( 1.0 + dt2*(w+k4) );
887 cfi = cfi_last + dt2*(cf_last+cf);
889 cb = (k3*cfi - k4*b) / (1.0 + k4*dt2);
890 cbi = cbi_last + dt2*(cb_last+cb);
893 ct[i]=cf+cb;
if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
894 if(cta!=NULL) {cta[i]=cf;
if(fabs(cta[i])<1.0e-12) cta[i]=0.0;}
895 if(ctb!=NULL) {ctb[i]=cb;
if(fabs(ctb[i])<1.0e-12) ctb[i]=0.0;}
897 t_last=t[i]; cr_last=cr[i];
898 cf_last=cf; cfi_last=cfi;
899 cb_last=cb; cbi_last=cbi;
937 double cri, cr_last, t_last;
938 double ct_last, cti, cti_last;
943 if(ct==NULL)
return 2;
946 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
947 cri=cr_last=0.0; cti=ct_last=cti_last=0.0;
948 for(i=0; i<nr; i++) {
950 dt2=0.5*(t[i]-t_last);
956 cri+=(cr[i]+cr_last)*dt2;
958 ct[i] = ( R1*cr[i] + k2*cri - (k2/(1.0+BP))*(cti_last+dt2*ct_last) ) /
959 ( 1.0 + dt2*(k2/(1.0+BP)) );
960 cti = cti_last + dt2*(ct_last+ct[i]);
963 if(!(fabs(ct[i])>=1.0e-12)) ct[i]=0.0;
965 t_last=t[i]; cr_last=cr[i];
966 ct_last=ct[i]; cti_last=cti;
1004 double cri, cr_last, t_last;
1005 double ct_last, cti, cti_last;
1010 if(ct==NULL)
return 2;
1013 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1014 cri=cr_last=0.0; cti=ct_last=cti_last=0.0;
1015 for(i=0; i<nr; i++) {
1017 dt2=0.5*(t[i]-t_last);
1021 }
else if(dt2>0.0) {
1023 cri+=(cr[i]+cr_last)*dt2;
1025 ct[i] = ( R1*cr[i] + R1*k3*cri - (k2+k3)*(cti_last+dt2*ct_last) ) /
1026 ( 1.0 + dt2*(k2+k3) );
1027 cti = cti_last + dt2*(ct_last+ct[i]);
1030 if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
1032 t_last=t[i]; cr_last=cr[i];
1033 ct_last=ct[i]; cti_last=cti;
1081 double dt2, t_last, ictot, ctot_last;
1082 double c0_, c1_, c2_, c3_, c4_;
1083 double ic0_, ic1_, ic2_, ic3_, ic4_;
1084 double c0_last, c1_last, c2_last, c3_last, c4_last;
1085 double ic0_last, ic1_last, ic2_last, ic3_last, ic4_last;
1086 double a, b, am1, am2, am3, am4;
1090 if(t==NULL || ctot==NULL || nr<2)
return(1);
1091 if(k01<0 || k12<0 || k21<0 || k03<0 || k34<0 || k43<0)
return(2);
1092 if(t[0]<0.0)
return(3);
1095 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1096 ictot=ctot_last=0.0;
1097 c0_=c0_last=ic0_=ic0_last=0.0;
1098 c1_=c1_last=ic1_=ic1_last=0.0;
1099 c2_=c2_last=ic2_=ic2_last=0.0;
1100 c3_=c3_last=ic3_=ic3_last=0.0;
1101 c4_=c4_last=ic4_=ic4_last=0.0;
1102 for(i=0; i<nr; i++) {
1104 dt2=0.5*(t[i]-t_last);
if(dt2<0.0)
return(5);
1107 a=k01+k12-(k12*k21*dt2/(1.0+dt2*k21));
1108 b=k03+k34-(k34*k43*dt2/(1.0+dt2*k43));
1109 am1=ic1_last+dt2*c1_last;
1110 am2=ic2_last+dt2*c2_last;
1111 am3=ic3_last+dt2*c3_last;
1112 am4=ic4_last+dt2*c4_last;
1115 ictot+=(ctot[i]+ctot_last)*dt2;
1117 c1_= ( k01*(1.0-k03*dt2/(1.0+dt2*b))*ictot
1118 -(a-k01*k03*dt2/(1.0+dt2*b))*am1
1119 +(k21/(1.0+dt2*k21))*am2
1120 -(k01/(1.0+dt2*b))*am3
1121 -(k01*k43*dt2/((1.0+dt2*b)*(1.0+dt2*k43)))*am4
1122 ) / ( 1.0+dt2*(a-k01*k03*dt2/(1.0+dt2*b)) );
1123 ic1_= ic1_last + dt2*(c1_+c1_last);
1125 c2_= (k12*ic1_-k21*am2)/(1.0+dt2*k21);
1126 ic2_= ic2_last + dt2*(c2_+c2_last);
1129 c3_= ( k03*(1.0-k01*dt2/(1.0+dt2*a))*ictot
1130 -(b-k01*k03*dt2/(1.0+dt2*a))*am3
1131 +(k43/(1.0+dt2*k43))*am4
1132 -(k03/(1.0+dt2*a))*am1
1133 -(k03*k21*dt2/((1.0+dt2*a)*(1.0+dt2*k21)))*am2
1134 ) / ( 1.0+dt2*(b-k01*k03*dt2/(1.0+dt2*a)) );
1135 ic3_= ic3_last + dt2*(c3_+c3_last);
1137 c4_= (k34*ic3_-k43*am4)/(1.0+dt2*k43);
1138 ic4_= ic4_last + dt2*(c4_+c4_last);
1141 c0_=ctot[i]-c1_-c3_;
1148 c0_last=c0_; c1_last=c1_; c2_last=c2_; c3_last=c3_; c4_last=c4_;
1149 ic0_last=ic0_; ic1_last=ic1_; ic2_last=ic2_; ic3_last=ic3_; ic4_last=ic4_;
1150 ctot_last=ctot[i]; t_last=t[i];
1189 double ctoti, ctot_last;
1191 double ct1m, ct1mi, ct1m_last, ct1mi_last;
1192 double ct2m, ct2mi, ct2m_last, ct2mi_last;
1193 double cpm, cpmi, cpm_last, cpmi_last;
1198 if(t==NULL || ctot==NULL || ca==NULL)
return 2;
1201 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1202 ctoti=ctot_last=0.0;
1203 ct1m_last=ct1mi_last=ct2m_last=ct2mi_last=cpm_last=cpmi_last=0.0;
1205 for(i=0; i<nr; i++) {
1207 dt2=0.5*(t[i]-t_last);
1213 ctoti+=(ctot[i]+ctot_last)*dt2;
1215 a=k4m/(1.0+k4m*dt2);
1216 b=(k1m-km)/(1.0+dt2*(k1m+k3m-k3m*a*dt2));
1218 ct1m = (km*ctoti - k2m*(1.0-dt2*b)*(ct1mi_last+dt2*ct1m_last)
1219 + b*(cpmi_last+dt2*cpm_last) + a*b*dt2*(ct2mi_last+dt2*ct2m_last)
1220 ) / (1.0 + dt2*k2m*(1.0-dt2*b));
1221 ct1mi = ct1mi_last + dt2*(ct1m_last+ct1m);
1223 cpm = (k2m*ct1mi + a*(ct2mi_last+dt2*ct2m_last)
1224 - (k1m+k3m-k3m*dt2*a)*(cpmi_last+dt2*cpm_last)
1225 ) / (1.0 + dt2*(k1m+k3m-k3m*dt2*a));
1226 cpmi = cpmi_last + dt2*(cpm_last+cpm);
1228 ct2m = (k3m*cpmi - k4m*(ct2mi_last+dt2*ct2m_last)) / (1.0 + dt2*k4m);
1229 ct2mi = ct2mi_last + dt2*(ct2m_last+ct2m);
1232 if(ca!=NULL) {ca[i]=ctot[i]-cpm;
if(fabs(ca[i])<1.0e-12) ca[i]=0.0;}
1233 if(cm!=NULL) {cm[i]=cpm;
if(fabs(cm[i])<1.0e-12) cm[i]=0.0;}
1235 t_last=t[i]; ctot_last=ctot[i];
1236 ct1m_last=ct1m; ct1mi_last=ct1mi;
1237 ct2m_last=ct2m; ct2mi_last=ct2mi;
1238 cpm_last=cpm; cpmi_last=cpmi;
1271 double cii, ci_last, t_last;
1272 double ct_last, cti, cti_last;
1277 if(ct==NULL)
return(2);
1280 t_last=0.0; cii=ci_last=0.0;
1281 cti=ct_last=cti_last=0.0;
1282 for(i=0; i<nr; i++) {
1284 dt2=0.5*(t[i]-t_last);
1288 }
else if(dt2>0.0) {
1290 cii+=(ci[i]+ci_last)*dt2;
1292 ct[i] = (Vfit*ci[i] + k1*cii - k2*(cti_last+dt2*ct_last)) / (1.0 + dt2*k2);
1293 cti = cti_last + dt2*(ct_last+ct[i]);
1296 if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
1298 t_last=t[i]; ci_last=ci[i];
1299 ct_last=ct[i]; cti_last=cti;
1333 double cai, ca_last, t_last;
1334 double ct1, ct1_last;
1335 double ct1i, ct1i_last;
1340 if(ct==NULL)
return 2;
1343 if(k1<0.0)
return 3;
1346 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1348 ct1_last=ct1i_last=0.0;
1350 for(i=0; i<nr; i++) {
1352 dt2=0.5*(t[i]-t_last);
1356 }
else if(dt2>0.0) {
1358 cai+=(ca[i]+ca_last)*dt2;
1360 ct1 = (k1*cai - k2*(ct1i_last+dt2*ct1_last)) / (1.0 + dt2*k2);
1361 ct1i = ct1i_last + dt2*(ct1_last+ct1);
1364 ct[i]=ct1;
if(fabs(ct[i])<1.0e-12) ct[i]=0.0;
1366 t_last=t[i]; ca_last=ca[i];
1367 ct1_last=ct1; ct1i_last=ct1i;
1436 double b, c, d, e, w, z, dt2, va, vv;
1437 double ca1i, ca1_last, ca2i, ca2_last, t_last, dct, cvb;
1438 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
1439 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
1440 double ct1b, ct1b_last, ct1bi, ct1bi_last;
1445 if(scpet==NULL)
return 2;
1448 if(k1<0.0)
return 3;
1449 if(vb<0.0 || vb>=1.0)
return 4;
1450 if(fa<=0.0 || fa>1.0)
return 5;
1451 va=fa*vb; vv=(1.0-fa)*vb;
1454 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1455 ca1i=ca1_last=ca2i=ca2_last=0.0;
1456 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=0.0;
1457 ct1b_last=ct1bi_last=0.0;
1458 ct1=ct2=ct3=ct1b=ct1i=ct2i=ct3i=ct1bi=0.0;
1459 for(i=0; i<nr; i++) {
1461 dt2=0.5*(t[i]-t_last);
1465 }
else if(dt2>0.0) {
1467 ca1i+=(ca1[i]+ca1_last)*dt2;
1468 ca2i+=(ca2[i]+ca2_last)*dt2;
1470 b=ct1i_last+dt2*ct1_last;
1471 c=ct2i_last+dt2*ct2_last;
1472 d=ct3i_last+dt2*ct3_last;
1473 e=ct1bi_last+dt2*ct1b_last;
1474 w=k4 + k5 - (k5*k6*dt2)/(1.0+k6*dt2);
1478 + k1*z*ca1i + (k3*k4*dt2 - (k2+k3)*z)*b
1479 + k4*c + k4*k6*dt2*d/(1.0+k6*dt2)
1480 ) / ( z*(1.0 + dt2*(k2+k3)) - k3*k4*dt2*dt2 );
1481 ct1i = ct1i_last + dt2*(ct1_last+ct1);
1482 ct1b = (k1b*ca2i - k2b*e) / (1.0 + dt2*k2b);
1483 ct1bi = ct1bi_last + dt2*(ct1b_last+ct1b);
1485 ct2 = (k3*ct1i - w*c + k6*d/(1.0+k6*dt2)) / z;
1486 ct2i = ct2i_last + dt2*(ct2_last+ct2);
1488 ct3 = (k5*ct2i - k6*d) / (1.0 + k6*dt2);
1489 ct3i = ct3i_last + dt2*(ct3_last+ct3);
1493 dct = k1*ca1[i] - k2*ct1 + k1b*ca2[i] - k2b*ct1b;
1494 cvb = cb[i] - dct/f;
1497 scpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2+ct3+ct1b);
1498 if(fabs(scpet[i])<1.0e-12) scpet[i]=0.0;
1499 if(sct1!=NULL) {sct1[i]=(1.0-vb)*ct1;
if(fabs(sct1[i])<1.0e-12) sct1[i]=0.0;}
1500 if(sct2!=NULL) {sct2[i]=(1.0-vb)*ct2;
if(fabs(sct2[i])<1.0e-12) sct2[i]=0.0;}
1501 if(sct3!=NULL) {sct3[i]=(1.0-vb)*ct3;
if(fabs(sct3[i])<1.0e-12) sct3[i]=0.0;}
1503 sct1b[i]=(1.0-vb)*ct1b;
if(fabs(sct1b[i])<1.0e-12) sct1b[i]=0.0;}
1504 if(sctab!=NULL) {sctab[i]=va*cb[i];
if(fabs(sctab[i])<1.0e-12) sctab[i]=0.0;}
1505 if(sctvb!=NULL) {sctvb[i]=vv*cvb;
if(fabs(sctvb[i])<1.0e-12) sctvb[i]=0.0;}
1507 t_last=t[i]; ca1_last=ca1[i]; ca2_last=ca2[i];
1508 ct1_last=ct1; ct1i_last=ct1i;
1509 ct2_last=ct2; ct2i_last=ct2i;
1510 ct3_last=ct3; ct3i_last=ct3i;
1511 ct1b_last=ct1b; ct1bi_last=ct1bi;
1588 double b, c, d, e, pt, qt, dt2, va, vv;
1589 double ca1i, ca1_last, ca2i, ca2_last, t_last, dct, cvb;
1590 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
1591 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
1592 double ct1b, ct1b_last, ct1bi, ct1bi_last;
1596 printf(
"simC4DIvp()\n");
1598 printf(
" k1 := %g\n", k1);
1599 printf(
" k2 := %g\n", k2);
1600 printf(
" k3 := %g\n", k3);
1601 printf(
" k4 := %g\n", k4);
1602 printf(
" k5 := %g\n", k5);
1603 printf(
" k6 := %g\n", k6);
1604 printf(
" k7 := %g\n", k7);
1605 printf(
" km := %g\n", km);
1606 printf(
" k1b := %g\n", k1b);
1607 printf(
" k2b := %g\n", k2b);
1608 printf(
" vb := %g\n", vb);
1609 printf(
" fa := %g\n", fa);
1610 printf(
" f := %g\n", f);
1616 if(scpet==NULL)
return 2;
1619 if(k1<0.0 || k1b<0.0)
return 3;
1620 if(vb<0.0 || vb>=1.0)
return 4;
1621 if(fa<0.0 || fa>1.0)
return 5;
1622 va=fa*vb; vv=(1.0-fa)*vb;
1625 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1626 ca1i=ca1_last=ca2i=ca2_last=0.0;
1627 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=ct1b_last=
1629 ct1=ct2=ct3=ct1b=ct1i=ct2i=ct3i=ct1bi=0.0;
1630 for(i=0; i<nr; i++) {
1632 dt2=0.5*(t[i]-t_last);
1636 }
else if(dt2>0.0) {
1638 ca1i+=(ca1[i]+ca1_last)*dt2;
1639 ca2i+=(ca2[i]+ca2_last)*dt2;
1641 b=ct1i_last+dt2*ct1_last;
1642 c=ct2i_last+dt2*ct2_last;
1643 d=ct3i_last+dt2*ct3_last;
1644 e=ct1bi_last+dt2*ct1b_last;
1646 qt=k2+k3+k5+km-(k3*k4*dt2)/(1.0+k4*dt2)-(k5*k6*dt2)/(1.0+pt*dt2);
1648 ct1 = (k1/(1.0+qt*dt2))*ca1i
1649 - (qt/(1.0+qt*dt2))*b
1650 + (k4/((1.0+qt*dt2)*(1.0+k4*dt2)))*c
1651 + (k6/((1.0+qt*dt2)*(1.0+pt*dt2)))*d;
1652 ct1i = ct1i_last + dt2*(ct1_last+ct1);
1654 ct2 = (k3/(1.0+k4*dt2))*ct1i
1655 - (k4/(1.0+k4*dt2))*c;
1656 ct2i = ct2i_last + dt2*(ct2_last+ct2);
1658 ct3 = (k5/(1.0+pt*dt2))*ct1i
1659 - (pt/(1.0+pt*dt2))*d;
1660 ct3i = ct3i_last + dt2*(ct3_last+ct3);
1662 ct1b = (k1b/(1.0+k2b*dt2))*ca2i
1663 - (k2b/(1.0+k2b*dt2))*e
1664 + (km/(1.0+k2b*dt2))*ct1i;
1665 ct1bi = ct1bi_last + dt2*(ct1b_last+ct1b);
1669 dct = k1*ca1[i] - k2*ct1 - k7*ct3 + k1b*ca2[i] - k2b*ct1b;
1670 cvb = cb[i] - dct/f;
1673 scpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2+ct3+ct1b);
1674 if(fabs(scpet[i])<1.0e-12) scpet[i]=0.0;
1676 sct1[i]=(1.0-vb)*ct1;
if(fabs(sct1[i])<1.0e-12) sct1[i]=0.0;}
1678 sct2[i]=(1.0-vb)*ct2;
if(fabs(sct2[i])<1.0e-12) sct2[i]=0.0;}
1680 sct3[i]=(1.0-vb)*ct3;
if(fabs(sct3[i])<1.0e-12) sct3[i]=0.0;}
1682 sct1b[i]=(1.0-vb)*ct1b;
if(fabs(sct1b[i])<1.0e-12) sct1b[i]=0.0;}
1684 sctab[i]=va*cb[i];
if(fabs(sctab[i])<1.0e-12) sctab[i]=0.0;}
1686 sctvb[i]=vv*cvb;
if(fabs(sctvb[i])<1.0e-12) sctvb[i]=0.0;}
1688 t_last=t[i]; ca1_last=ca1[i]; ca2_last=ca2[i];
1689 ct1_last=ct1; ct1i_last=ct1i;
1690 ct2_last=ct2; ct2i_last=ct2i;
1691 ct3_last=ct3; ct3i_last=ct3i;
1692 ct1b_last=ct1b; ct1bi_last=ct1bi;
1696 printf(
"AUC 0-%g:\n", t_last);
1697 printf(
" ca1i := %g\n", ca1i);
1698 printf(
" ca2i := %g\n", ca2i);
1699 printf(
" ct1i := %g\n", ct1i_last);
1700 printf(
" ct2i := %g\n", ct2i_last);
1701 printf(
" ct3i := %g\n", ct3i_last);
1702 printf(
" ct1bi := %g\n", ct1bi_last);
1779 double b, c, d, e, pt, qt, rt, dt2, va, vv;
1780 double ca1i, ca1_last, ca2i, ca2_last, t_last, dct, cvb;
1781 double ct1, ct1_last, ct2, ct2_last, ct3, ct3_last;
1782 double ct1i, ct1i_last, ct2i, ct2i_last, ct3i, ct3i_last;
1783 double ct1b, ct1b_last, ct1bi, ct1bi_last;
1787 printf(
"simC4DIvs()\n");
1789 printf(
" k1 := %g\n", k1);
1790 printf(
" k2 := %g\n", k2);
1791 printf(
" k3 := %g\n", k3);
1792 printf(
" k4 := %g\n", k4);
1793 printf(
" k5 := %g\n", k5);
1794 printf(
" k6 := %g\n", k6);
1795 printf(
" k7 := %g\n", k7);
1796 printf(
" km := %g\n", km);
1797 printf(
" k1b := %g\n", k1b);
1798 printf(
" k2b := %g\n", k2b);
1799 printf(
" vb := %g\n", vb);
1800 printf(
" fa := %g\n", fa);
1801 printf(
" f := %g\n", f);
1807 if(scpet==NULL)
return 2;
1810 if(k1<0.0 || k1b<0.0)
return 3;
1811 if(vb<0.0 || vb>=1.0)
return 4;
1812 if(fa<0.0 || fa>1.0)
return 5;
1813 va=fa*vb; vv=(1.0-fa)*vb;
1816 t_last=0.0;
if(t[0]<t_last) t_last=t[0];
1817 ca1i=ca1_last=ca2i=ca2_last=0.0;
1818 ct1_last=ct2_last=ct3_last=ct1i_last=ct2i_last=ct3i_last=ct1b_last=
1820 ct1=ct2=ct3=ct1b=ct1i=ct2i=ct3i=ct1bi=0.0;
1821 for(i=0; i<nr; i++) {
1823 dt2=0.5*(t[i]-t_last);
1827 }
else if(dt2>0.0) {
1829 ca1i+=(ca1[i]+ca1_last)*dt2;
1830 ca2i+=(ca2[i]+ca2_last)*dt2;
1834 b=ct1i_last+dt2*ct1_last;
1835 c=ct2i_last+dt2*ct2_last;
1836 d=ct3i_last+dt2*ct3_last;
1837 e=ct1bi_last+dt2*ct1b_last;
1839 qt=k4+k5-(k5*k6*dt2)/(1.0+pt*dt2);
1840 rt=k2+k3+km-(k3*k4*dt2)/(1.0+qt*dt2);
1842 ct1 = (k1/(1.0+rt*dt2))*ca1i
1843 - (rt/(1.0+rt*dt2))*b
1844 + (k4/((1.0+qt*dt2)*(1.0+rt*dt2)))*c
1845 + ((k4*k6*dt2)/((1.0+pt*dt2)*(1.0+qt*dt2)*(1.0+rt*dt2)))*d;
1846 ct1i = ct1i_last + dt2*(ct1_last+ct1);
1848 ct2 = (k3/(1.0+qt*dt2))*ct1i
1849 - (qt/(1.0+qt*dt2))*c
1850 + (k6/((1.0+pt*dt2)*(1.0+qt*dt2)))*d;
1851 ct2i = ct2i_last + dt2*(ct2_last+ct2);
1853 ct3 = (k5/(1.0+pt*dt2))*ct2i
1854 - (pt/(1.0+pt*dt2))*d;
1855 ct3i = ct3i_last + dt2*(ct3_last+ct3);
1857 ct1b = (k1b/(1.0+k2b*dt2))*ca2i
1858 - (k2b/(1.0+k2b*dt2))*e
1859 + (km/(1.0+k2b*dt2))*ct1i;
1860 ct1bi = ct1bi_last + dt2*(ct1b_last+ct1b);
1864 dct = k1*ca1[i] - k2*ct1 - k7*ct3 + k1b*ca2[i] - k2b*ct1b;
1865 cvb = cb[i] - dct/f;
1868 scpet[i]= va*cb[i] + vv*cvb + (1.0-vb)*(ct1+ct2+ct3+ct1b);
1869 if(fabs(scpet[i])<1.0e-12) scpet[i]=0.0;
1871 sct1[i]=(1.0-vb)*ct1;
if(fabs(sct1[i])<1.0e-12) sct1[i]=0.0;}
1873 sct2[i]=(1.0-vb)*ct2;
if(fabs(sct2[i])<1.0e-12) sct2[i]=0.0;}
1875 sct3[i]=(1.0-vb)*ct3;
if(fabs(sct3[i])<1.0e-12) sct3[i]=0.0;}
1877 sct1b[i]=(1.0-vb)*ct1b;
if(fabs(sct1b[i])<1.0e-12) sct1b[i]=0.0;}
1879 sctab[i]=va*cb[i];
if(fabs(sctab[i])<1.0e-12) sctab[i]=0.0;}
1881 sctvb[i]=vv*cvb;
if(fabs(sctvb[i])<1.0e-12) sctvb[i]=0.0;}
1883 t_last=t[i]; ca1_last=ca1[i]; ca2_last=ca2[i];
1884 ct1_last=ct1; ct1i_last=ct1i;
1885 ct2_last=ct2; ct2i_last=ct2i;
1886 ct3_last=ct3; ct3i_last=ct3i;
1887 ct1b_last=ct1b; ct1bi_last=ct1bi;
1891 printf(
"AUC 0-%g:\n", t_last);
1892 printf(
" ca1i := %g\n", ca1i);
1893 printf(
" ca2i := %g\n", ca2i);
1894 printf(
" ct1i := %g\n", ct1i_last);
1895 printf(
" ct2i := %g\n", ct2i_last);
1896 printf(
" ct3i := %g\n", ct3i_last);
1897 printf(
" ct1bi := %g\n", ct1bi_last);
1927 if(x==NULL || y==NULL || n<2)
return 1;
1928 if(tau1<0.0 || tau2<0.0)
return 2;
1932 if(tmp!=NULL) buf=tmp;
else buf=(
double*)malloc(n*
sizeof(
double));
1933 if(buf==NULL)
return 3;
1938 int ret=
simC1(x, y, n, k, k, buf);
1940 if(tmp==NULL) free(buf);
1943 for(
int i=0; i<n; i++) y[i]=buf[i];
1949 int ret=
simC1(x, y, n, k, k, buf);
1951 if(tmp==NULL) free(buf);
1954 for(
int i=0; i<n; i++) y[i]=buf[i];
1957 if(tmp==NULL) free(buf);
2031 printf(
"simOxygen()\n");
2033 printf(
" k1a := %g\n", k1a);
2034 printf(
" k2a := %g\n", k2a);
2035 printf(
" km := %g\n", km);
2036 printf(
" k1b := %g\n", k1b);
2037 printf(
" k2b := %g\n", k2b);
2038 printf(
" vb := %g\n", vb);
2039 printf(
" fa := %g\n", fa);
2040 printf(
" n := %d\n", n);
2046 if(t==NULL)
return 2;
2047 if(ca1==NULL || ca2==NULL)
return 3;
2050 if(k1a<0.0 || k1b<0.0 || k2a<0.0 || k2b<0.0)
return(4);
2051 if(vb<0.0 || vb>=1.0)
return(5);
2052 if(fa<0.0 || fa>1.0)
return(6);
2054 double vv=(1.0-fa)*vb;
2058 if(t[0]<t_last) t_last=t[0];
2060 double cba1i=0.0, cba1_last=0.0;
2061 double cba2i=0.0, cba2_last=0.0;
2062 double ct1=0.0, ct1_last=0.0, ct1i=0.0, ct1i_last=0.0;
2063 double ct2=0.0, ct2_last=0.0, ct2i=0.0, ct2i_last=0.0;
2064 double cvb1=0.0, cvb2=0.0;
2069 for(
int i=0; i<n; i++) {
2070 dt2=0.5*(t[i]-t_last);
2073 if(ca1i!=NULL) cba1i=ca1i[i];
else cba1i+=(ca1[i]+cba1_last)*dt2;
2074 if(ca2i!=NULL) cba2i=ca2i[i];
else cba2i+=(ca2[i]+cba2_last)*dt2;
2076 p=ct1i_last+dt2*ct1_last;
2077 q=ct2i_last+dt2*ct2_last;
2079 ct1 = (k1a*cba1i - (k2a+km)*p) / (1.0 + dt2*(k2a+km));
2080 ct1i = ct1i_last + dt2*(ct1_last+ct1);
2082 ct2 = (km*ct1i + k1b*cba2i - k2b*q) / (1.0 + dt2*k2b);
2083 ct2i = ct2i_last + dt2*(ct2_last+ct2);
2086 if(k1a>0.0 && k2a>0.0) cvb1=ct1/(k1a/k2a);
2087 else if(k2a>0.0) cvb1=0.0;
else cvb1=ca1[i];
2088 if(k1b>0.0 && k2b>0.0) cvb2=ct2/(k1b/k2b);
2089 else if(k2b>0.0) cvb2=0.0;
else cvb2=ca2[i];
2092 scpet[i]= va*(ca1[i]+ca2[i]) + vv*(cvb1+cvb2) + (1.0-vb)*(ct1+ct2);
2093 if(fabs(scpet[i])<1.0e-12) scpet[i]=0.0;
2096 sct1[i]=(1.0-vb)*ct1;
2097 if(fabs(sct1[i])<1.0e-12) sct1[i]=0.0;
2100 sct2[i]=(1.0-vb)*ct2;
2101 if(fabs(sct2[i])<1.0e-12) sct2[i]=0.0;
2104 sctab[i]=va*(ca1[i]+ca2[i]);
2105 if(fabs(sctab[i])<1.0e-12) sctab[i]=0.0;
2109 if(fabs(sctvb1[i])<1.0e-12) sctvb1[i]=0.0;
2113 if(fabs(sctvb2[i])<1.0e-12) sctvb2[i]=0.0;
2115 if(scvb1!=NULL) scvb1[i]=cvb1;
2116 if(scvb2!=NULL) scvb2[i]=cvb2;
2119 cba1_last=ca1[i]; cba2_last=ca2[i];
2120 ct1_last=ct1; ct1i_last=ct1i;
2121 ct2_last=ct2; ct2i_last=ct2i;
2125 printf(
"AUC 0-%g:\n", t_last);
2126 printf(
" cba1i := %g\n", cba1i);
2127 printf(
" cba2i := %g\n", cba2i);
2128 printf(
" ct1i := %g\n", ct1i_last);
2129 printf(
" ct2i := %g\n", ct2i_last);
Header file for libtpcmodel.
int simTRTM(double *t, double *cr, int nr, double R1, double k2, double k3, double *ct)
int simC3vp(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)
int simC4DIvp(double *t, double *ca1, double *ca2, double *cb, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double k7, double km, double k1b, double k2b, double f, double vb, double fa, double *scpet, double *sct1, double *sct2, double *sct3, double *sct1b, double *sctab, double *sctvb, int verbose)
int simTPCMOD0009c(double *t, double *ctot, int nr, double km, double k1m, double k2m, double k3m, double k4m, double *ca, double *cm)
int simOxygen(double *t, double *ca1, double *ca2, double *ca1i, double *ca2i, const int n, const double k1a, const double k2a, const double km, const double k1b, const double k2b, const double vb, const double fa, double *scpet, double *sct1, double *sct2, double *sctab, double *sctvb1, double *sctvb2, double *scvb1, double *scvb2, const int verbose)
int simHuangmet(double *t, double *ctot, int nr, double k01, double k12, double k21, double k03, double k34, double k43, double *c0, double *c1, double *c3)
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)
int simDispersion(double *x, double *y, int n, double tau1, double tau2, double *tmp)
int simC4DIvs(double *t, double *ca1, double *ca2, double *cb, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double k7, double km, double k1b, double k2b, double f, double vb, double fa, double *scpet, double *sct1, double *sct2, double *sct3, double *sct1b, double *sctab, double *sctvb, int verbose)
int simSRTM(double *t, double *cr, int nr, double R1, double k2, double BP, double *ct)
int simC1(double *t, double *ca, int nr, double k1, double k2, double *ct)
int simRTCM(double *t, double *cr, int nr, double R1, double k2, double k3, double k4, double *ct, double *cta, double *ctb)
int simMBF(double *t, double *ci, int nr, double k1, double k2, double Vfit, double *ct)
int simC3vpKLoss(double *t, double *ca, double *cb, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double kLoss, double f, double vb, double fa, double *cpet, double *cta, double *ctb, double *ctc, double *ctab, double *ctvb)
int simC2vl(double *t, double *ca, double *cb, int nr, double k1, double k2, double k3, double kL, double f, double vb, double fa, double *cpet, double *cta, double *ctb, double *ctab, double *ctvb)
int simC3p(double *t, double *ca, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double *ct, double *cta, double *ctb, double *ctc)
int simC2l(double *t, double *ca, int nr, double k1, double k2, double k3, double kLoss, double *ct, double *cta, double *ctb)
int simC3DIvs(double *t, double *ca1, double *ca2, double *cb, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double k1b, double k2b, double f, double vb, double fa, double *scpet, double *sct1, double *sct2, double *sct3, double *sct1b, double *sctab, double *sctvb)
int simC3s(double *t, double *ca, int nr, double k1, double k2, double k3, double k4, double k5, double k6, double *ct, double *cta, double *ctb, double *ctc)