1:
2: /*********************** Imach **************************************
3: This program computes Healthy Life Expectancies from cross-longitudinal
4: data. Cross-longitudinal consist in a first survey ("cross") where
5: individuals from different ages are interviewed on their health status
6: or degree of disability. At least a second wave of interviews
7: ("longitudinal") should measure each new individual health status.
8: Health expectancies are computed from the transistions observed between
9: waves and are computed for each degree of severity of disability (number
10: of life states). More degrees you consider, more time is necessary to
11: reach the Maximum Likelihood of the parameters involved in the model.
12: The simplest model is the multinomial logistic model where pij is
13: the probabibility to be observed in state j at the second wave conditional
14: to be observed in state i at the first wave. Therefore the model is:
15: log(pij/pii)= aij + bij*age+ cij*sex + etc , where 'age' is age and 'sex'
16: is a covariate. If you want to have a more complex model than "constant and
17: age", you should modify the program where the markup
18: *Covariates have to be included here again* invites you to do it.
19: More covariates you add, less is the speed of the convergence.
20:
21: The advantage that this computer programme claims, comes from that if the
22: delay between waves is not identical for each individual, or if some
23: individual missed an interview, the information is not rounded or lost, but
24: taken into account using an interpolation or extrapolation.
25: hPijx is the probability to be
26: observed in state i at age x+h conditional to the observed state i at age
27: x. The delay 'h' can be split into an exact number (nh*stepm) of
28: unobserved intermediate states. This elementary transition (by month or
29: quarter trimester, semester or year) is model as a multinomial logistic.
30: The hPx matrix is simply the matrix product of nh*stepm elementary matrices
31: and the contribution of each individual to the likelihood is simply hPijx.
32:
33: Also this programme outputs the covariance matrix of the parameters but also
34: of the life expectancies. It also computes the prevalence limits.
35:
36: Authors: Nicolas Brouard (brouard@ined.fr) and Agnès Lièvre (lievre@ined.fr).
37: Institut national d'études démographiques, Paris.
38: This software have been partly granted by Euro-REVES, a concerted action
39: from the European Union.
40: It is copyrighted identically to a GNU software product, ie programme and
41: software can be distributed freely for non commercial use. Latest version
42: can be accessed at http://euroreves.ined.fr/imach .
43: **********************************************************************/
44:
45: #include <math.h>
46: #include <stdio.h>
47: #include <stdlib.h>
48: #include <unistd.h>
49:
50: #define MAXLINE 256
51: #define FILENAMELENGTH 80
52: /*#define DEBUG*/
53: #define windows
54:
55: #define MAXPARM 30 /* Maximum number of parameters for the optimization */
56: #define NPARMAX 64 /* (nlstate+ndeath-1)*nlstate*ncovmodel */
57:
58: #define NINTERVMAX 8
59: #define NLSTATEMAX 8 /* Maximum number of live states (for func) */
60: #define NDEATHMAX 8 /* Maximum number of dead states (for func) */
61: #define NCOVMAX 8 /* Maximum number of covariates */
62: #define MAXN 80000
63: #define YEARM 12. /* Number of months per year */
64: #define AGESUP 130
65: #define AGEBASE 40
66:
67:
68: int nvar;
69: static int cptcov;
70: int cptcovn;
71: int npar=NPARMAX;
72: int nlstate=2; /* Number of live states */
73: int ndeath=1; /* Number of dead states */
74: int ncovmodel, ncov; /* Total number of covariables including constant a12*1 +b12*x ncovmodel=2 */
75:
76: int *wav; /* Number of waves for this individuual 0 is possible */
77: int maxwav; /* Maxim number of waves */
78: int mle, weightopt;
79: int **mw; /* mw[mi][i] is number of the mi wave for this individual */
80: int **dh; /* dh[mi][i] is number of steps between mi,mi+1 for this individual */
81: double **oldm, **newm, **savm; /* Working pointers to matrices */
82: double **oldms, **newms, **savms; /* Fixed working pointers to matrices */
83: FILE *fic,*ficpar, *ficparo,*ficres, *ficrespl, *ficrespij, *ficrest;
84: FILE *ficgp, *fichtm;
85: FILE *ficreseij;
86: char filerese[FILENAMELENGTH];
87: FILE *ficresvij;
88: char fileresv[FILENAMELENGTH];
89: FILE *ficresvpl;
90: char fileresvpl[FILENAMELENGTH];
91:
92:
93:
94:
95: #define NR_END 1
96: #define FREE_ARG char*
97: #define FTOL 1.0e-10
98:
99: #define NRANSI
100: #define ITMAX 200
101:
102: #define TOL 2.0e-4
103:
104: #define CGOLD 0.3819660
105: #define ZEPS 1.0e-10
106: #define SHFT(a,b,c,d) (a)=(b);(b)=(c);(c)=(d);
107:
108: #define GOLD 1.618034
109: #define GLIMIT 100.0
110: #define TINY 1.0e-20
111:
112: static double maxarg1,maxarg2;
113: #define FMAX(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)>(maxarg2)? (maxarg1):(maxarg2))
114: #define FMIN(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)<(maxarg2)? (maxarg1):(maxarg2))
115:
116: #define SIGN(a,b) ((b)>0.0 ? fabs(a) : -fabs(a))
117: #define rint(a) floor(a+0.5)
118:
119: static double sqrarg;
120: #define SQR(a) ((sqrarg=(a)) == 0.0 ? 0.0 :sqrarg*sqrarg)
121: #define SWAP(a,b) {temp=(a);(a)=(b);(b)=temp;}
122:
123: int imx;
124: int stepm;
125: /* Stepm, step in month: minimum step interpolation*/
126:
127: int m,nb;
128: int *num, firstpass=0, lastpass=4,*cod, *ncodemax;
129: double **agev,*moisnais, *annais, *moisdc, *andc,**mint, **anint;
130: double **pmmij;
131:
132: double *weight;
133: int **s; /* Status */
134: double *agedc, **covar, idx;
135: int **nbcode, *Tcode, *Tvar, **codtab;
136:
137: double ftol=FTOL; /* Tolerance for computing Max Likelihood */
138: double ftolhess; /* Tolerance for computing hessian */
139:
140:
141: /******************************************/
142:
143: void replace(char *s, char*t)
144: {
145: int i;
146: int lg=20;
147: i=0;
148: lg=strlen(t);
149: for(i=0; i<= lg; i++) {
150: (s[i] = t[i]);
151: if (t[i]== '\\') s[i]='/';
152: }
153: }
154:
155: int nbocc(char *s, char occ)
156: {
157: int i,j=0;
158: int lg=20;
159: i=0;
160: lg=strlen(s);
161: for(i=0; i<= lg; i++) {
162: if (s[i] == occ ) j++;
163: }
164: return j;
165: }
166:
167: void cutv(char *u,char *v, char*t, char occ)
168: {
169: int i,lg,j,p;
170: i=0;
171: if (t[0]== occ) p=0;
172: for(j=0; j<=strlen(t)-1; j++) {
173: if((t[j]!= occ) && (t[j+1]==occ)) p=j+1;
174: }
175:
176: lg=strlen(t);
177: for(j=0; j<p; j++) {
178: (u[j] = t[j]);
179: u[p]='\0';
180: }
181:
182: for(j=0; j<= lg; j++) {
183: if (j>=(p+1))(v[j-p-1] = t[j]);
184: }
185: }
186:
187:
188: /********************** nrerror ********************/
189:
190: void nrerror(char error_text[])
191: {
192: fprintf(stderr,"ERREUR ...\n");
193: fprintf(stderr,"%s\n",error_text);
194: exit(1);
195: }
196: /*********************** vector *******************/
197: double *vector(int nl, int nh)
198: {
199: double *v;
200: v=(double *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(double)));
201: if (!v) nrerror("allocation failure in vector");
202: return v-nl+NR_END;
203: }
204:
205: /************************ free vector ******************/
206: void free_vector(double*v, int nl, int nh)
207: {
208: free((FREE_ARG)(v+nl-NR_END));
209: }
210:
211: /************************ivector *******************************/
212: int *ivector(long nl,long nh)
213: {
214: int *v;
215: v=(int *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(int)));
216: if (!v) nrerror("allocation failure in ivector");
217: return v-nl+NR_END;
218: }
219:
220: /******************free ivector **************************/
221: void free_ivector(int *v, long nl, long nh)
222: {
223: free((FREE_ARG)(v+nl-NR_END));
224: }
225:
226: /******************* imatrix *******************************/
227: int **imatrix(long nrl, long nrh, long ncl, long nch)
228: /* allocate a int matrix with subscript range m[nrl..nrh][ncl..nch] */
229: {
230: long i, nrow=nrh-nrl+1,ncol=nch-ncl+1;
231: int **m;
232:
233: /* allocate pointers to rows */
234: m=(int **) malloc((size_t)((nrow+NR_END)*sizeof(int*)));
235: if (!m) nrerror("allocation failure 1 in matrix()");
236: m += NR_END;
237: m -= nrl;
238:
239:
240: /* allocate rows and set pointers to them */
241: m[nrl]=(int *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(int)));
242: if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
243: m[nrl] += NR_END;
244: m[nrl] -= ncl;
245:
246: for(i=nrl+1;i<=nrh;i++) m[i]=m[i-1]+ncol;
247:
248: /* return pointer to array of pointers to rows */
249: return m;
250: }
251:
252: /****************** free_imatrix *************************/
253: void free_imatrix(m,nrl,nrh,ncl,nch)
254: int **m;
255: long nch,ncl,nrh,nrl;
256: /* free an int matrix allocated by imatrix() */
257: {
258: free((FREE_ARG) (m[nrl]+ncl-NR_END));
259: free((FREE_ARG) (m+nrl-NR_END));
260: }
261:
262: /******************* matrix *******************************/
263: double **matrix(long nrl, long nrh, long ncl, long nch)
264: {
265: long i, nrow=nrh-nrl+1, ncol=nch-ncl+1;
266: double **m;
267:
268: m=(double **) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
269: if (!m) nrerror("allocation failure 1 in matrix()");
270: m += NR_END;
271: m -= nrl;
272:
273: m[nrl]=(double *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
274: if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
275: m[nrl] += NR_END;
276: m[nrl] -= ncl;
277:
278: for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
279: return m;
280: }
281:
282: /*************************free matrix ************************/
283: void free_matrix(double **m, long nrl, long nrh, long ncl, long nch)
284: {
285: free((FREE_ARG)(m[nrl]+ncl-NR_END));
286: free((FREE_ARG)(m+nrl-NR_END));
287: }
288:
289: /******************* ma3x *******************************/
290: double ***ma3x(long nrl, long nrh, long ncl, long nch, long nll, long nlh)
291: {
292: long i, j, nrow=nrh-nrl+1, ncol=nch-ncl+1, nlay=nlh-nll+1;
293: double ***m;
294:
295: m=(double ***) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
296: if (!m) nrerror("allocation failure 1 in matrix()");
297: m += NR_END;
298: m -= nrl;
299:
300: m[nrl]=(double **) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
301: if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
302: m[nrl] += NR_END;
303: m[nrl] -= ncl;
304:
305: for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
306:
307: m[nrl][ncl]=(double *) malloc((size_t)((nrow*ncol*nlay+NR_END)*sizeof(double)));
308: if (!m[nrl][ncl]) nrerror("allocation failure 3 in matrix()");
309: m[nrl][ncl] += NR_END;
310: m[nrl][ncl] -= nll;
311: for (j=ncl+1; j<=nch; j++)
312: m[nrl][j]=m[nrl][j-1]+nlay;
313:
314: for (i=nrl+1; i<=nrh; i++) {
315: m[i][ncl]=m[i-1l][ncl]+ncol*nlay;
316: for (j=ncl+1; j<=nch; j++)
317: m[i][j]=m[i][j-1]+nlay;
318: }
319: return m;
320: }
321:
322: /*************************free ma3x ************************/
323: void free_ma3x(double ***m, long nrl, long nrh, long ncl, long nch,long nll, long nlh)
324: {
325: free((FREE_ARG)(m[nrl][ncl]+ nll-NR_END));
326: free((FREE_ARG)(m[nrl]+ncl-NR_END));
327: free((FREE_ARG)(m+nrl-NR_END));
328: }
329:
330: /***************** f1dim *************************/
331: extern int ncom;
332: extern double *pcom,*xicom;
333: extern double (*nrfunc)(double []);
334:
335: double f1dim(double x)
336: {
337: int j;
338: double f;
339: double *xt;
340:
341: xt=vector(1,ncom);
342: for (j=1;j<=ncom;j++) xt[j]=pcom[j]+x*xicom[j];
343: f=(*nrfunc)(xt);
344: free_vector(xt,1,ncom);
345: return f;
346: }
347:
348: /*****************brent *************************/
349: double brent(double ax, double bx, double cx, double (*f)(double), double tol, double *xmin)
350: {
351: int iter;
352: double a,b,d,etemp;
353: double fu,fv,fw,fx;
354: double ftemp;
355: double p,q,r,tol1,tol2,u,v,w,x,xm;
356: double e=0.0;
357:
358: a=(ax < cx ? ax : cx);
359: b=(ax > cx ? ax : cx);
360: x=w=v=bx;
361: fw=fv=fx=(*f)(x);
362: for (iter=1;iter<=ITMAX;iter++) {
363: xm=0.5*(a+b);
364: tol2=2.0*(tol1=tol*fabs(x)+ZEPS);
365: /* if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret)))*/
366: printf(".");fflush(stdout);
367: #ifdef DEBUG
368: printf("br %d,x=%.10e xm=%.10e b=%.10e a=%.10e tol=%.10e tol1=%.10e tol2=%.10e x-xm=%.10e fx=%.12e fu=%.12e,fw=%.12e,ftemp=%.12e,ftol=%.12e\n",iter,x,xm,b,a,tol,tol1,tol2,(x-xm),fx,fu,fw,ftemp,ftol);
369: /* if ((fabs(x-xm) <= (tol2-0.5*(b-a)))||(2.0*fabs(fu-ftemp) <= ftol*1.e-2*(fabs(fu)+fabs(ftemp)))) { */
370: #endif
371: if (fabs(x-xm) <= (tol2-0.5*(b-a))){
372: *xmin=x;
373: return fx;
374: }
375: ftemp=fu;
376: if (fabs(e) > tol1) {
377: r=(x-w)*(fx-fv);
378: q=(x-v)*(fx-fw);
379: p=(x-v)*q-(x-w)*r;
380: q=2.0*(q-r);
381: if (q > 0.0) p = -p;
382: q=fabs(q);
383: etemp=e;
384: e=d;
385: if (fabs(p) >= fabs(0.5*q*etemp) || p <= q*(a-x) || p >= q*(b-x))
386: d=CGOLD*(e=(x >= xm ? a-x : b-x));
387: else {
388: d=p/q;
389: u=x+d;
390: if (u-a < tol2 || b-u < tol2)
391: d=SIGN(tol1,xm-x);
392: }
393: } else {
394: d=CGOLD*(e=(x >= xm ? a-x : b-x));
395: }
396: u=(fabs(d) >= tol1 ? x+d : x+SIGN(tol1,d));
397: fu=(*f)(u);
398: if (fu <= fx) {
399: if (u >= x) a=x; else b=x;
400: SHFT(v,w,x,u)
401: SHFT(fv,fw,fx,fu)
402: } else {
403: if (u < x) a=u; else b=u;
404: if (fu <= fw || w == x) {
405: v=w;
406: w=u;
407: fv=fw;
408: fw=fu;
409: } else if (fu <= fv || v == x || v == w) {
410: v=u;
411: fv=fu;
412: }
413: }
414: }
415: nrerror("Too many iterations in brent");
416: *xmin=x;
417: return fx;
418: }
419:
420: /****************** mnbrak ***********************/
421:
422: void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb, double *fc,
423: double (*func)(double))
424: {
425: double ulim,u,r,q, dum;
426: double fu;
427:
428: *fa=(*func)(*ax);
429: *fb=(*func)(*bx);
430: if (*fb > *fa) {
431: SHFT(dum,*ax,*bx,dum)
432: SHFT(dum,*fb,*fa,dum)
433: }
434: *cx=(*bx)+GOLD*(*bx-*ax);
435: *fc=(*func)(*cx);
436: while (*fb > *fc) {
437: r=(*bx-*ax)*(*fb-*fc);
438: q=(*bx-*cx)*(*fb-*fa);
439: u=(*bx)-((*bx-*cx)*q-(*bx-*ax)*r)/
440: (2.0*SIGN(FMAX(fabs(q-r),TINY),q-r));
441: ulim=(*bx)+GLIMIT*(*cx-*bx);
442: if ((*bx-u)*(u-*cx) > 0.0) {
443: fu=(*func)(u);
444: } else if ((*cx-u)*(u-ulim) > 0.0) {
445: fu=(*func)(u);
446: if (fu < *fc) {
447: SHFT(*bx,*cx,u,*cx+GOLD*(*cx-*bx))
448: SHFT(*fb,*fc,fu,(*func)(u))
449: }
450: } else if ((u-ulim)*(ulim-*cx) >= 0.0) {
451: u=ulim;
452: fu=(*func)(u);
453: } else {
454: u=(*cx)+GOLD*(*cx-*bx);
455: fu=(*func)(u);
456: }
457: SHFT(*ax,*bx,*cx,u)
458: SHFT(*fa,*fb,*fc,fu)
459: }
460: }
461:
462: /*************** linmin ************************/
463:
464: int ncom;
465: double *pcom,*xicom;
466: double (*nrfunc)(double []);
467:
468: void linmin(double p[], double xi[], int n, double *fret,double (*func)(double []))
469: {
470: double brent(double ax, double bx, double cx,
471: double (*f)(double), double tol, double *xmin);
472: double f1dim(double x);
473: void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb,
474: double *fc, double (*func)(double));
475: int j;
476: double xx,xmin,bx,ax;
477: double fx,fb,fa;
478:
479: ncom=n;
480: pcom=vector(1,n);
481: xicom=vector(1,n);
482: nrfunc=func;
483: for (j=1;j<=n;j++) {
484: pcom[j]=p[j];
485: xicom[j]=xi[j];
486: }
487: ax=0.0;
488: xx=1.0;
489: mnbrak(&ax,&xx,&bx,&fa,&fx,&fb,f1dim);
490: *fret=brent(ax,xx,bx,f1dim,TOL,&xmin);
491: #ifdef DEBUG
492: printf("retour brent fret=%.12e xmin=%.12e\n",*fret,xmin);
493: #endif
494: for (j=1;j<=n;j++) {
495: xi[j] *= xmin;
496: p[j] += xi[j];
497: }
498: free_vector(xicom,1,n);
499: free_vector(pcom,1,n);
500: }
501:
502: /*************** powell ************************/
503: void powell(double p[], double **xi, int n, double ftol, int *iter, double *fret,
504: double (*func)(double []))
505: {
506: void linmin(double p[], double xi[], int n, double *fret,
507: double (*func)(double []));
508: int i,ibig,j;
509: double del,t,*pt,*ptt,*xit;
510: double fp,fptt;
511: double *xits;
512: pt=vector(1,n);
513: ptt=vector(1,n);
514: xit=vector(1,n);
515: xits=vector(1,n);
516: *fret=(*func)(p);
517: for (j=1;j<=n;j++) pt[j]=p[j];
518: for (*iter=1;;++(*iter)) {
519: fp=(*fret);
520: ibig=0;
521: del=0.0;
522: printf("\nPowell iter=%d -2*LL=%.12f",*iter,*fret);
523: for (i=1;i<=n;i++)
524: printf(" %d %.12f",i, p[i]);
525: printf("\n");
526: for (i=1;i<=n;i++) {
527: for (j=1;j<=n;j++) xit[j]=xi[j][i];
528: fptt=(*fret);
529: #ifdef DEBUG
530: printf("fret=%lf \n",*fret);
531: #endif
532: printf("%d",i);fflush(stdout);
533: linmin(p,xit,n,fret,func);
534: if (fabs(fptt-(*fret)) > del) {
535: del=fabs(fptt-(*fret));
536: ibig=i;
537: }
538: #ifdef DEBUG
539: printf("%d %.12e",i,(*fret));
540: for (j=1;j<=n;j++) {
541: xits[j]=FMAX(fabs(p[j]-pt[j]),1.e-5);
542: printf(" x(%d)=%.12e",j,xit[j]);
543: }
544: for(j=1;j<=n;j++)
545: printf(" p=%.12e",p[j]);
546: printf("\n");
547: #endif
548: }
549: if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret))) {
550: #ifdef DEBUG
551: int k[2],l;
552: k[0]=1;
553: k[1]=-1;
554: printf("Max: %.12e",(*func)(p));
555: for (j=1;j<=n;j++)
556: printf(" %.12e",p[j]);
557: printf("\n");
558: for(l=0;l<=1;l++) {
559: for (j=1;j<=n;j++) {
560: ptt[j]=p[j]+(p[j]-pt[j])*k[l];
561: printf("l=%d j=%d ptt=%.12e, xits=%.12e, p=%.12e, xit=%.12e", l,j,ptt[j],xits[j],p[j],xit[j]);
562: }
563: printf("func(ptt)=%.12e, deriv=%.12e\n",(*func)(ptt),(ptt[j]-p[j])/((*func)(ptt)-(*func)(p)));
564: }
565: #endif
566:
567:
568: free_vector(xit,1,n);
569: free_vector(xits,1,n);
570: free_vector(ptt,1,n);
571: free_vector(pt,1,n);
572: return;
573: }
574: if (*iter == ITMAX) nrerror("powell exceeding maximum iterations.");
575: for (j=1;j<=n;j++) {
576: ptt[j]=2.0*p[j]-pt[j];
577: xit[j]=p[j]-pt[j];
578: pt[j]=p[j];
579: }
580: fptt=(*func)(ptt);
581: if (fptt < fp) {
582: t=2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del)-del*SQR(fp-fptt);
583: if (t < 0.0) {
584: linmin(p,xit,n,fret,func);
585: for (j=1;j<=n;j++) {
586: xi[j][ibig]=xi[j][n];
587: xi[j][n]=xit[j];
588: }
589: #ifdef DEBUG
590: printf("Direction changed last moved %d in place of ibig=%d, new last is the average:\n",n,ibig);
591: for(j=1;j<=n;j++)
592: printf(" %.12e",xit[j]);
593: printf("\n");
594: #endif
595: }
596: }
597: }
598: }
599:
600: /**** Prevalence limit ****************/
601:
602: double **prevalim(double **prlim, int nlstate, double x[], double age, double **oldm, double **savm, double ftolpl, int ij)
603: {
604: /* Computes the prevalence limit in each live state at age x by left multiplying the unit
605: matrix by transitions matrix until convergence is reached */
606:
607: int i, ii,j,k;
608: double min, max, maxmin, maxmax,sumnew=0.;
609: double **matprod2();
610: double **out, cov[NCOVMAX], **pmij();
611: double **newm;
612: double agefin, delaymax=50 ; /* Max number of years to converge */
613:
614: for (ii=1;ii<=nlstate+ndeath;ii++)
615: for (j=1;j<=nlstate+ndeath;j++){
616: oldm[ii][j]=(ii==j ? 1.0 : 0.0);
617: }
618: /* Even if hstepm = 1, at least one multiplication by the unit matrix */
619: for(agefin=age-stepm/YEARM; agefin>=age-delaymax; agefin=agefin-stepm/YEARM){
620: newm=savm;
621: /* Covariates have to be included here again */
622: cov[1]=1.;
623: cov[2]=agefin;
624: if (cptcovn>0){
625: for (k=1; k<=cptcovn;k++) {cov[2+k]=nbcode[Tvar[k]][codtab[ij][k]];/*printf("Tcode[ij]=%d nbcode=%d\n",Tcode[ij],nbcode[k][Tcode[ij]]);*/}
626: }
627: out=matprod2(newm, pmij(pmmij,cov,ncovmodel,x,nlstate),1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, oldm);
628:
629: savm=oldm;
630: oldm=newm;
631: maxmax=0.;
632: for(j=1;j<=nlstate;j++){
633: min=1.;
634: max=0.;
635: for(i=1; i<=nlstate; i++) {
636: sumnew=0;
637: for(k=1; k<=ndeath; k++) sumnew+=newm[i][nlstate+k];
638: prlim[i][j]= newm[i][j]/(1-sumnew);
639: max=FMAX(max,prlim[i][j]);
640: min=FMIN(min,prlim[i][j]);
641: }
642: maxmin=max-min;
643: maxmax=FMAX(maxmax,maxmin);
644: }
645: if(maxmax < ftolpl){
646: return prlim;
647: }
648: }
649: }
650:
651: /*************** transition probabilities **********/
652:
653: double **pmij(double **ps, double *cov, int ncovmodel, double *x, int nlstate )
654: {
655: double s1, s2;
656: /*double t34;*/
657: int i,j,j1, nc, ii, jj;
658:
659: for(i=1; i<= nlstate; i++){
660: for(j=1; j<i;j++){
661: for (nc=1, s2=0.;nc <=ncovmodel; nc++){
662: /*s2 += param[i][j][nc]*cov[nc];*/
663: s2 += x[(i-1)*nlstate*ncovmodel+(j-1)*ncovmodel+nc+(i-1)*(ndeath-1)*ncovmodel]*cov[nc];
664: /*printf("Int j<i s1=%.17e, s2=%.17e\n",s1,s2);*/
665: }
666: ps[i][j]=s2;
667: /*printf("s1=%.17e, s2=%.17e\n",s1,s2);*/
668: }
669: for(j=i+1; j<=nlstate+ndeath;j++){
670: for (nc=1, s2=0.;nc <=ncovmodel; nc++){
671: s2 += x[(i-1)*nlstate*ncovmodel+(j-2)*ncovmodel+nc+(i-1)*(ndeath-1)*ncovmodel]*cov[nc];
672: /*printf("Int j>i s1=%.17e, s2=%.17e %lx %lx\n",s1,s2,s1,s2);*/
673: }
674: ps[i][j]=s2;
675: }
676: }
677: for(i=1; i<= nlstate; i++){
678: s1=0;
679: for(j=1; j<i; j++)
680: s1+=exp(ps[i][j]);
681: for(j=i+1; j<=nlstate+ndeath; j++)
682: s1+=exp(ps[i][j]);
683: ps[i][i]=1./(s1+1.);
684: for(j=1; j<i; j++)
685: ps[i][j]= exp(ps[i][j])*ps[i][i];
686: for(j=i+1; j<=nlstate+ndeath; j++)
687: ps[i][j]= exp(ps[i][j])*ps[i][i];
688: /* ps[i][nlstate+1]=1.-s1- ps[i][i];*/ /* Sum should be 1 */
689: } /* end i */
690:
691: for(ii=nlstate+1; ii<= nlstate+ndeath; ii++){
692: for(jj=1; jj<= nlstate+ndeath; jj++){
693: ps[ii][jj]=0;
694: ps[ii][ii]=1;
695: }
696: }
697:
698: /* for(ii=1; ii<= nlstate+ndeath; ii++){
699: for(jj=1; jj<= nlstate+ndeath; jj++){
700: printf("%lf ",ps[ii][jj]);
701: }
702: printf("\n ");
703: }
704: printf("\n ");printf("%lf ",cov[2]);*/
705: /*
706: for(i=1; i<= npar; i++) printf("%f ",x[i]);
707: goto end;*/
708: return ps;
709: }
710:
711: /**************** Product of 2 matrices ******************/
712:
713: double **matprod2(double **out, double **in,long nrl, long nrh, long ncl, long nch, long ncolol, long ncoloh, double **b)
714: {
715: /* Computes the matric product of in(1,nrh-nrl+1)(1,nch-ncl+1) times
716: b(1,nch-ncl+1)(1,ncoloh-ncolol+1) into out(...) */
717: /* in, b, out are matrice of pointers which should have been initialized
718: before: only the contents of out is modified. The function returns
719: a pointer to pointers identical to out */
720: long i, j, k;
721: for(i=nrl; i<= nrh; i++)
722: for(k=ncolol; k<=ncoloh; k++)
723: for(j=ncl,out[i][k]=0.; j<=nch; j++)
724: out[i][k] +=in[i][j]*b[j][k];
725:
726: return out;
727: }
728:
729:
730: /************* Higher Matrix Product ***************/
731:
732: double ***hpxij(double ***po, int nhstepm, double age, int hstepm, double *x, int nlstate, int stepm, double **oldm, double **savm, int ij )
733: {
734: /* Computes the transition matrix starting at age 'age' over 'nhstepm*hstepm*stepm' month
735: duration (i.e. until
736: age (in years) age+nhstepm*stepm/12) by multiplying nhstepm*hstepm matrices.
737: Output is stored in matrix po[i][j][h] for h every 'hstepm' step
738: (typically every 2 years instead of every month which is too big).
739: Model is determined by parameters x and covariates have to be
740: included manually here.
741:
742: */
743:
744: int i, j, d, h, k;
745: double **out, cov[NCOVMAX];
746: double **newm;
747:
748: /* Hstepm could be zero and should return the unit matrix */
749: for (i=1;i<=nlstate+ndeath;i++)
750: for (j=1;j<=nlstate+ndeath;j++){
751: oldm[i][j]=(i==j ? 1.0 : 0.0);
752: po[i][j][0]=(i==j ? 1.0 : 0.0);
753: }
754: /* Even if hstepm = 1, at least one multiplication by the unit matrix */
755: for(h=1; h <=nhstepm; h++){
756: for(d=1; d <=hstepm; d++){
757: newm=savm;
758: /* Covariates have to be included here again */
759: cov[1]=1.;
760: cov[2]=age+((h-1)*hstepm + (d-1))*stepm/YEARM;
761: if (cptcovn>0){
762: for (k=1; k<=cptcovn;k++) cov[2+k]=nbcode[Tvar[k]][codtab[ij][k]];
763: }
764: /*printf("hxi cptcov=%d cptcode=%d\n",cptcov,cptcode);*/
765: /*printf("h=%d d=%d age=%f cov=%f\n",h,d,age,cov[2]);*/
766: out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath,
767: pmij(pmmij,cov,ncovmodel,x,nlstate));
768: savm=oldm;
769: oldm=newm;
770: }
771: for(i=1; i<=nlstate+ndeath; i++)
772: for(j=1;j<=nlstate+ndeath;j++) {
773: po[i][j][h]=newm[i][j];
774: /*printf("i=%d j=%d h=%d po[i][j][h]=%f ",i,j,h,po[i][j][h]);
775: */
776: }
777: } /* end h */
778: return po;
779: }
780:
781:
782: /*************** log-likelihood *************/
783: double func( double *x)
784: {
785: int i, ii, j, k, mi, d;
786: double l, ll[NLSTATEMAX], cov[NCOVMAX];
787: double **out;
788: double sw; /* Sum of weights */
789: double lli; /* Individual log likelihood */
790: long ipmx;
791: /*extern weight */
792: /* We are differentiating ll according to initial status */
793: /* for (i=1;i<=npar;i++) printf("%f ", x[i]);*/
794: /*for(i=1;i<imx;i++)
795: printf(" %d\n",s[4][i]);
796: */
797:
798: for(k=1; k<=nlstate; k++) ll[k]=0.;
799: for (i=1,ipmx=0, sw=0.; i<=imx; i++){
800: for(mi=1; mi<= wav[i]-1; mi++){
801: for (ii=1;ii<=nlstate+ndeath;ii++)
802: for (j=1;j<=nlstate+ndeath;j++) oldm[ii][j]=(ii==j ? 1.0 : 0.0);
803: for(d=0; d<dh[mi][i]; d++){
804: newm=savm;
805: cov[1]=1.;
806: cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
807: if (cptcovn>0){
808: for (k=1; k<=cptcovn;k++) {
809: cov[2+k]=covar[Tvar[k]][i];
810: /* printf("k=%d cptcovn=%d %lf\n",k,cptcovn,covar[Tvar[k]][i]);*/
811: }
812: }
813: out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
814: 1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
815: savm=oldm;
816: oldm=newm;
817: } /* end mult */
818:
819: lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]);
820: /* printf(" %f ",out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]);*/
821: ipmx +=1;
822: sw += weight[i];
823: ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
824: } /* end of wave */
825: } /* end of individual */
826:
827: for(k=1,l=0.; k<=nlstate; k++) l += ll[k];
828: /* printf("l1=%f l2=%f ",ll[1],ll[2]); */
829: l= l*ipmx/sw; /* To get the same order of magnitude as if weight=1 for every body */
830:
831: return -l;
832: }
833:
834:
835: /*********** Maximum Likelihood Estimation ***************/
836:
837: void mlikeli(FILE *ficres,double p[], int npar, int ncovmodel, int nlstate, double ftol, double (*func)(double []))
838: {
839: int i,j, iter;
840: double **xi,*delti;
841: double fret;
842: xi=matrix(1,npar,1,npar);
843: for (i=1;i<=npar;i++)
844: for (j=1;j<=npar;j++)
845: xi[i][j]=(i==j ? 1.0 : 0.0);
846: printf("Powell\n");
847: powell(p,xi,npar,ftol,&iter,&fret,func);
848:
849: printf("\n#Number of iterations = %d, -2 Log likelihood = %.12f\n",iter,func(p));
850: fprintf(ficres,"#Number of iterations = %d, -2 Log likelihood = %.12f ",iter,func(p));
851:
852: }
853:
854: /**** Computes Hessian and covariance matrix ***/
855: void hesscov(double **matcov, double p[], int npar, double delti[], double ftolhess, double (*func)(double []))
856: {
857: double **a,**y,*x,pd;
858: double **hess;
859: int i, j,jk;
860: int *indx;
861:
862: double hessii(double p[], double delta, int theta, double delti[]);
863: double hessij(double p[], double delti[], int i, int j);
864: void lubksb(double **a, int npar, int *indx, double b[]) ;
865: void ludcmp(double **a, int npar, int *indx, double *d) ;
866:
867:
868: hess=matrix(1,npar,1,npar);
869:
870: printf("\nCalculation of the hessian matrix. Wait...\n");
871: for (i=1;i<=npar;i++){
872: printf("%d",i);fflush(stdout);
873: hess[i][i]=hessii(p,ftolhess,i,delti);
874: /*printf(" %f ",p[i]);*/
875: }
876:
877: for (i=1;i<=npar;i++) {
878: for (j=1;j<=npar;j++) {
879: if (j>i) {
880: printf(".%d%d",i,j);fflush(stdout);
881: hess[i][j]=hessij(p,delti,i,j);
882: hess[j][i]=hess[i][j];
883: }
884: }
885: }
886: printf("\n");
887:
888: printf("\nInverting the hessian to get the covariance matrix. Wait...\n");
889:
890: a=matrix(1,npar,1,npar);
891: y=matrix(1,npar,1,npar);
892: x=vector(1,npar);
893: indx=ivector(1,npar);
894: for (i=1;i<=npar;i++)
895: for (j=1;j<=npar;j++) a[i][j]=hess[i][j];
896: ludcmp(a,npar,indx,&pd);
897:
898: for (j=1;j<=npar;j++) {
899: for (i=1;i<=npar;i++) x[i]=0;
900: x[j]=1;
901: lubksb(a,npar,indx,x);
902: for (i=1;i<=npar;i++){
903: matcov[i][j]=x[i];
904: }
905: }
906:
907: printf("\n#Hessian matrix#\n");
908: for (i=1;i<=npar;i++) {
909: for (j=1;j<=npar;j++) {
910: printf("%.3e ",hess[i][j]);
911: }
912: printf("\n");
913: }
914:
915: /* Recompute Inverse */
916: for (i=1;i<=npar;i++)
917: for (j=1;j<=npar;j++) a[i][j]=matcov[i][j];
918: ludcmp(a,npar,indx,&pd);
919:
920: /* printf("\n#Hessian matrix recomputed#\n");
921:
922: for (j=1;j<=npar;j++) {
923: for (i=1;i<=npar;i++) x[i]=0;
924: x[j]=1;
925: lubksb(a,npar,indx,x);
926: for (i=1;i<=npar;i++){
927: y[i][j]=x[i];
928: printf("%.3e ",y[i][j]);
929: }
930: printf("\n");
931: }
932: */
933:
934: free_matrix(a,1,npar,1,npar);
935: free_matrix(y,1,npar,1,npar);
936: free_vector(x,1,npar);
937: free_ivector(indx,1,npar);
938: free_matrix(hess,1,npar,1,npar);
939:
940:
941: }
942:
943: /*************** hessian matrix ****************/
944: double hessii( double x[], double delta, int theta, double delti[])
945: {
946: int i;
947: int l=1, lmax=20;
948: double k1,k2;
949: double p2[NPARMAX+1];
950: double res;
951: double delt, delts, nkhi=10.,nkhif=1., khi=1.e-4;
952: double fx;
953: int k=0,kmax=10;
954: double l1;
955:
956: fx=func(x);
957: for (i=1;i<=npar;i++) p2[i]=x[i];
958: for(l=0 ; l <=lmax; l++){
959: l1=pow(10,l);
960: delts=delt;
961: for(k=1 ; k <kmax; k=k+1){
962: delt = delta*(l1*k);
963: p2[theta]=x[theta] +delt;
964: k1=func(p2)-fx;
965: p2[theta]=x[theta]-delt;
966: k2=func(p2)-fx;
967: /*res= (k1-2.0*fx+k2)/delt/delt; */
968: res= (k1+k2)/delt/delt/2.; /* Divided by because L and not 2*L */
969:
970: #ifdef DEBUG
971: printf("%d %d k1=%.12e k2=%.12e xk1=%.12e xk2=%.12e delt=%.12e res=%.12e l=%d k=%d,fx=%.12e\n",theta,theta,k1,k2,x[theta]+delt,x[theta]-delt,delt,res, l, k,fx);
972: #endif
973: /*if(fabs(k1-2.0*fx+k2) <1.e-13){ */
974: if((k1 <khi/nkhi/2.) || (k2 <khi/nkhi/2.)){
975: k=kmax;
976: }
977: else if((k1 >khi/nkhif) || (k2 >khi/nkhif)){ /* Keeps lastvalue before 3.84/2 KHI2 5% 1d.f. */
978: k=kmax; l=lmax*10.;
979: }
980: else if((k1 >khi/nkhi) || (k2 >khi/nkhi)){
981: delts=delt;
982: }
983: }
984: }
985: delti[theta]=delts;
986: return res;
987: }
988:
989: double hessij( double x[], double delti[], int thetai,int thetaj)
990: {
991: int i;
992: int l=1, l1, lmax=20;
993: double k1,k2,k3,k4,res,fx;
994: double p2[NPARMAX+1];
995: int k;
996:
997: fx=func(x);
998: for (k=1; k<=2; k++) {
999: for (i=1;i<=npar;i++) p2[i]=x[i];
1000: p2[thetai]=x[thetai]+delti[thetai]/k;
1001: p2[thetaj]=x[thetaj]+delti[thetaj]/k;
1002: k1=func(p2)-fx;
1003:
1004: p2[thetai]=x[thetai]+delti[thetai]/k;
1005: p2[thetaj]=x[thetaj]-delti[thetaj]/k;
1006: k2=func(p2)-fx;
1007:
1008: p2[thetai]=x[thetai]-delti[thetai]/k;
1009: p2[thetaj]=x[thetaj]+delti[thetaj]/k;
1010: k3=func(p2)-fx;
1011:
1012: p2[thetai]=x[thetai]-delti[thetai]/k;
1013: p2[thetaj]=x[thetaj]-delti[thetaj]/k;
1014: k4=func(p2)-fx;
1015: res=(k1-k2-k3+k4)/4.0/delti[thetai]*k/delti[thetaj]*k/2.; /* Because of L not 2*L */
1016: #ifdef DEBUG
1017: printf("%d %d k=%d, k1=%.12e k2=%.12e k3=%.12e k4=%.12e delti/k=%.12e deltj/k=%.12e, xi-de/k=%.12e xj-de/k=%.12e res=%.12e k1234=%.12e,k1-2=%.12e,k3-4=%.12e\n",thetai,thetaj,k,k1,k2,k3,k4,delti[thetai]/k,delti[thetaj]/k,x[thetai]-delti[thetai]/k,x[thetaj]-delti[thetaj]/k, res,k1-k2-k3+k4,k1-k2,k3-k4);
1018: #endif
1019: }
1020: return res;
1021: }
1022:
1023: /************** Inverse of matrix **************/
1024: void ludcmp(double **a, int n, int *indx, double *d)
1025: {
1026: int i,imax,j,k;
1027: double big,dum,sum,temp;
1028: double *vv;
1029:
1030: vv=vector(1,n);
1031: *d=1.0;
1032: for (i=1;i<=n;i++) {
1033: big=0.0;
1034: for (j=1;j<=n;j++)
1035: if ((temp=fabs(a[i][j])) > big) big=temp;
1036: if (big == 0.0) nrerror("Singular matrix in routine ludcmp");
1037: vv[i]=1.0/big;
1038: }
1039: for (j=1;j<=n;j++) {
1040: for (i=1;i<j;i++) {
1041: sum=a[i][j];
1042: for (k=1;k<i;k++) sum -= a[i][k]*a[k][j];
1043: a[i][j]=sum;
1044: }
1045: big=0.0;
1046: for (i=j;i<=n;i++) {
1047: sum=a[i][j];
1048: for (k=1;k<j;k++)
1049: sum -= a[i][k]*a[k][j];
1050: a[i][j]=sum;
1051: if ( (dum=vv[i]*fabs(sum)) >= big) {
1052: big=dum;
1053: imax=i;
1054: }
1055: }
1056: if (j != imax) {
1057: for (k=1;k<=n;k++) {
1058: dum=a[imax][k];
1059: a[imax][k]=a[j][k];
1060: a[j][k]=dum;
1061: }
1062: *d = -(*d);
1063: vv[imax]=vv[j];
1064: }
1065: indx[j]=imax;
1066: if (a[j][j] == 0.0) a[j][j]=TINY;
1067: if (j != n) {
1068: dum=1.0/(a[j][j]);
1069: for (i=j+1;i<=n;i++) a[i][j] *= dum;
1070: }
1071: }
1072: free_vector(vv,1,n); /* Doesn't work */
1073: ;
1074: }
1075:
1076: void lubksb(double **a, int n, int *indx, double b[])
1077: {
1078: int i,ii=0,ip,j;
1079: double sum;
1080:
1081: for (i=1;i<=n;i++) {
1082: ip=indx[i];
1083: sum=b[ip];
1084: b[ip]=b[i];
1085: if (ii)
1086: for (j=ii;j<=i-1;j++) sum -= a[i][j]*b[j];
1087: else if (sum) ii=i;
1088: b[i]=sum;
1089: }
1090: for (i=n;i>=1;i--) {
1091: sum=b[i];
1092: for (j=i+1;j<=n;j++) sum -= a[i][j]*b[j];
1093: b[i]=sum/a[i][i];
1094: }
1095: }
1096:
1097: /************ Frequencies ********************/
1098: void freqsummary(char fileres[], int agemin, int agemax, int **s, double **agev, int nlstate, int imx, int *Tvar, int **nbcode, int *ncodemax)
1099: { /* Some frequencies */
1100:
1101: int i, m, jk, k1, i1, j1, bool, z1,z2,j;
1102: double ***freq; /* Frequencies */
1103: double *pp;
1104: double pos;
1105: FILE *ficresp;
1106: char fileresp[FILENAMELENGTH];
1107:
1108: pp=vector(1,nlstate);
1109:
1110: strcpy(fileresp,"p");
1111: strcat(fileresp,fileres);
1112: if((ficresp=fopen(fileresp,"w"))==NULL) {
1113: printf("Problem with prevalence resultfile: %s\n", fileresp);
1114: exit(0);
1115: }
1116: freq= ma3x(-1,nlstate+ndeath,-1,nlstate+ndeath,agemin,agemax+3);
1117: j1=0;
1118:
1119: j=cptcovn;
1120: if (cptcovn<1) {j=1;ncodemax[1]=1;}
1121:
1122: for(k1=1; k1<=j;k1++){
1123: for(i1=1; i1<=ncodemax[k1];i1++){
1124: j1++;
1125:
1126: for (i=-1; i<=nlstate+ndeath; i++)
1127: for (jk=-1; jk<=nlstate+ndeath; jk++)
1128: for(m=agemin; m <= agemax+3; m++)
1129: freq[i][jk][m]=0;
1130:
1131: for (i=1; i<=imx; i++) {
1132: bool=1;
1133: if (cptcovn>0) {
1134: for (z1=1; z1<=cptcovn; z1++)
1135: if (covar[Tvar[z1]][i]!= nbcode[Tvar[z1]][codtab[j1][z1]]) bool=0;
1136: }
1137: if (bool==1) {
1138: for(m=firstpass; m<=lastpass-1; m++){
1139: if(agev[m][i]==0) agev[m][i]=agemax+1;
1140: if(agev[m][i]==1) agev[m][i]=agemax+2;
1141: freq[s[m][i]][s[m+1][i]][(int)agev[m][i]] += weight[i];
1142: freq[s[m][i]][s[m+1][i]][(int) agemax+3] += weight[i];
1143: }
1144: }
1145: }
1146: if (cptcovn>0) {
1147: fprintf(ficresp, "\n#Variable");
1148: for (z1=1; z1<=cptcovn; z1++) fprintf(ficresp, " V%d=%d",Tvar[z1],nbcode[Tvar[z1]][codtab[j1][z1]]);
1149: }
1150: fprintf(ficresp, "\n#");
1151: for(i=1; i<=nlstate;i++)
1152: fprintf(ficresp, " Age Prev(%d) N(%d) N",i,i);
1153: fprintf(ficresp, "\n");
1154:
1155: for(i=(int)agemin; i <= (int)agemax+3; i++){
1156: if(i==(int)agemax+3)
1157: printf("Total");
1158: else
1159: printf("Age %d", i);
1160: for(jk=1; jk <=nlstate ; jk++){
1161: for(m=-1, pp[jk]=0; m <=nlstate+ndeath ; m++)
1162: pp[jk] += freq[jk][m][i];
1163: }
1164: for(jk=1; jk <=nlstate ; jk++){
1165: for(m=-1, pos=0; m <=0 ; m++)
1166: pos += freq[jk][m][i];
1167: if(pp[jk]>=1.e-10)
1168: printf(" %d.=%.0f loss[%d]=%.1f%%",jk,pp[jk],jk,100*pos/pp[jk]);
1169: else
1170: printf(" %d.=%.0f loss[%d]=NaNQ%%",jk,pp[jk],jk);
1171: }
1172: for(jk=1; jk <=nlstate ; jk++){
1173: for(m=1, pp[jk]=0; m <=nlstate+ndeath; m++)
1174: pp[jk] += freq[jk][m][i];
1175: }
1176: for(jk=1,pos=0; jk <=nlstate ; jk++)
1177: pos += pp[jk];
1178: for(jk=1; jk <=nlstate ; jk++){
1179: if(pos>=1.e-5)
1180: printf(" %d.=%.0f prev[%d]=%.1f%%",jk,pp[jk],jk,100*pp[jk]/pos);
1181: else
1182: printf(" %d.=%.0f prev[%d]=NaNQ%%",jk,pp[jk],jk);
1183: if( i <= (int) agemax){
1184: if(pos>=1.e-5)
1185: fprintf(ficresp," %d %.5f %.0f %.0f",i,pp[jk]/pos, pp[jk],pos);
1186: else
1187: fprintf(ficresp," %d NaNq %.0f %.0f",i,pp[jk],pos);
1188: }
1189: }
1190: for(jk=-1; jk <=nlstate+ndeath; jk++)
1191: for(m=-1; m <=nlstate+ndeath; m++)
1192: if(freq[jk][m][i] !=0 ) printf(" %d%d=%.0f",jk,m,freq[jk][m][i]);
1193: if(i <= (int) agemax)
1194: fprintf(ficresp,"\n");
1195: printf("\n");
1196: }
1197: }
1198: }
1199:
1200: fclose(ficresp);
1201: free_ma3x(freq,-1,nlstate+ndeath,-1,nlstate+ndeath,(int) agemin,(int) agemax+3);
1202: free_vector(pp,1,nlstate);
1203:
1204: } /* End of Freq */
1205:
1206: /************* Waves Concatenation ***************/
1207:
1208: void concatwav(int wav[], int **dh, int **mw, int **s, double *agedc, double **agev, int firstpass, int lastpass, int imx, int nlstate, int stepm)
1209: {
1210: /* Concatenates waves: wav[i] is the number of effective (useful waves) of individual i.
1211: Death is a valid wave (if date is known).
1212: mw[mi][i] is the mi (mi=1 to wav[i]) effective wave of individual i
1213: dh[m][i] of dh[mw[mi][i][i] is the delay between two effective waves m=mw[mi][i]
1214: and mw[mi+1][i]. dh depends on stepm.
1215: */
1216:
1217: int i, mi, m;
1218: int j, k=0,jk, ju, jl,jmin=1e+5, jmax=-1;
1219: float sum=0.;
1220:
1221: for(i=1; i<=imx; i++){
1222: mi=0;
1223: m=firstpass;
1224: while(s[m][i] <= nlstate){
1225: if(s[m][i]>=1)
1226: mw[++mi][i]=m;
1227: if(m >=lastpass)
1228: break;
1229: else
1230: m++;
1231: }/* end while */
1232: if (s[m][i] > nlstate){
1233: mi++; /* Death is another wave */
1234: /* if(mi==0) never been interviewed correctly before death */
1235: /* Only death is a correct wave */
1236: mw[mi][i]=m;
1237: }
1238:
1239: wav[i]=mi;
1240: if(mi==0)
1241: printf("Warning, no any valid information for:%d line=%d\n",num[i],i);
1242: }
1243:
1244: for(i=1; i<=imx; i++){
1245: for(mi=1; mi<wav[i];mi++){
1246: if (stepm <=0)
1247: dh[mi][i]=1;
1248: else{
1249: if (s[mw[mi+1][i]][i] > nlstate) {
1250: j= rint(agedc[i]*12-agev[mw[mi][i]][i]*12);
1251: if(j=0) j=1; /* Survives at least one month after exam */
1252: }
1253: else{
1254: j= rint( (agev[mw[mi+1][i]][i]*12 - agev[mw[mi][i]][i]*12));
1255: /*printf("i=%d agevi+1=%lf agevi=%lf j=%d\n", i,agev[mw[mi+1][i]][i],agev[mw[mi][i]][i],j);*/
1256:
1257: k=k+1;
1258: if (j >= jmax) jmax=j;
1259: else if (j <= jmin)jmin=j;
1260: sum=sum+j;
1261: }
1262: jk= j/stepm;
1263: jl= j -jk*stepm;
1264: ju= j -(jk+1)*stepm;
1265: if(jl <= -ju)
1266: dh[mi][i]=jk;
1267: else
1268: dh[mi][i]=jk+1;
1269: if(dh[mi][i]==0)
1270: dh[mi][i]=1; /* At least one step */
1271: }
1272: }
1273: }
1274: printf("Delay (in months) between two waves Min=%d Max=%d Mean=%f\n\n ",jmin, jmax,sum/k);
1275: }
1276: /*********** Tricode ****************************/
1277: void tricode(int *Tvar, int **nbcode, int imx)
1278: {
1279: int Ndum[80],ij, k, j, i;
1280: int cptcode=0;
1281: for (k=0; k<79; k++) Ndum[k]=0;
1282: for (k=1; k<=7; k++) ncodemax[k]=0;
1283:
1284: for (j=1; j<=cptcovn; j++) {
1285: for (i=1; i<=imx; i++) {
1286: ij=(int)(covar[Tvar[j]][i]);
1287: Ndum[ij]++;
1288: if (ij > cptcode) cptcode=ij;
1289: }
1290: /*printf("cptcode=%d cptcovn=%d ",cptcode,cptcovn);*/
1291: for (i=0; i<=cptcode; i++) {
1292: if(Ndum[i]!=0) ncodemax[j]++;
1293: }
1294:
1295: ij=1;
1296: for (i=1; i<=ncodemax[j]; i++) {
1297: for (k=0; k<=79; k++) {
1298: if (Ndum[k] != 0) {
1299: nbcode[Tvar[j]][ij]=k;
1300: ij++;
1301: }
1302: if (ij > ncodemax[j]) break;
1303: }
1304: }
1305: }
1306:
1307: }
1308:
1309: /*********** Health Expectancies ****************/
1310:
1311: void evsij(char fileres[], double ***eij, double x[], int nlstate, int stepm, int bage, int fage, double **oldm, double **savm, int ij)
1312: {
1313: /* Health expectancies */
1314: int i, j, nhstepm, hstepm, h;
1315: double age, agelim,hf;
1316: double ***p3mat;
1317:
1318: fprintf(ficreseij,"# Health expectancies\n");
1319: fprintf(ficreseij,"# Age");
1320: for(i=1; i<=nlstate;i++)
1321: for(j=1; j<=nlstate;j++)
1322: fprintf(ficreseij," %1d-%1d",i,j);
1323: fprintf(ficreseij,"\n");
1324:
1325: hstepm=1*YEARM; /* Every j years of age (in month) */
1326: hstepm=hstepm/stepm; /* Typically in stepm units, if j= 2 years, = 2/6 months = 4 */
1327:
1328: agelim=AGESUP;
1329: for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
1330: /* nhstepm age range expressed in number of stepm */
1331: nhstepm=(int) rint((agelim-age)*YEARM/stepm);
1332: /* Typically if 20 years = 20*12/6=40 stepm */
1333: if (stepm >= YEARM) hstepm=1;
1334: nhstepm = nhstepm/hstepm;/* Expressed in hstepm, typically 40/4=10 */
1335: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1336: /* Computed by stepm unit matrices, product of hstepm matrices, stored
1337: in an array of nhstepm length: nhstepm=10, hstepm=4, stepm=6 months */
1338: hpxij(p3mat,nhstepm,age,hstepm,x,nlstate,stepm,oldm, savm, ij);
1339:
1340:
1341: for(i=1; i<=nlstate;i++)
1342: for(j=1; j<=nlstate;j++)
1343: for (h=0, eij[i][j][(int)age]=0; h<=nhstepm; h++){
1344: eij[i][j][(int)age] +=p3mat[i][j][h];
1345: }
1346:
1347: hf=1;
1348: if (stepm >= YEARM) hf=stepm/YEARM;
1349: fprintf(ficreseij,"%.0f",age );
1350: for(i=1; i<=nlstate;i++)
1351: for(j=1; j<=nlstate;j++){
1352: fprintf(ficreseij," %.4f", hf*eij[i][j][(int)age]);
1353: }
1354: fprintf(ficreseij,"\n");
1355: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1356: }
1357: }
1358:
1359: /************ Variance ******************/
1360: void varevsij(char fileres[], double ***vareij, double **matcov, double x[], double delti[], int nlstate, int stepm, double bage, double fage, double **oldm, double **savm, double **prlim, double ftolpl, int ij)
1361: {
1362: /* Variance of health expectancies */
1363: /* double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double ** savm,double ftolpl);*/
1364: double **newm;
1365: double **dnewm,**doldm;
1366: int i, j, nhstepm, hstepm, h;
1367: int k, cptcode;
1368: double *xp;
1369: double **gp, **gm;
1370: double ***gradg, ***trgradg;
1371: double ***p3mat;
1372: double age,agelim;
1373: int theta;
1374:
1375: fprintf(ficresvij,"# Covariances of life expectancies\n");
1376: fprintf(ficresvij,"# Age");
1377: for(i=1; i<=nlstate;i++)
1378: for(j=1; j<=nlstate;j++)
1379: fprintf(ficresvij," Cov(e%1d, e%1d)",i,j);
1380: fprintf(ficresvij,"\n");
1381:
1382: xp=vector(1,npar);
1383: dnewm=matrix(1,nlstate,1,npar);
1384: doldm=matrix(1,nlstate,1,nlstate);
1385:
1386: hstepm=1*YEARM; /* Every year of age */
1387: hstepm=hstepm/stepm; /* Typically in stepm units, if j= 2 years, = 2/6 months = 4 */
1388: agelim = AGESUP;
1389: for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
1390: nhstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */
1391: if (stepm >= YEARM) hstepm=1;
1392: nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
1393: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1394: gradg=ma3x(0,nhstepm,1,npar,1,nlstate);
1395: gp=matrix(0,nhstepm,1,nlstate);
1396: gm=matrix(0,nhstepm,1,nlstate);
1397:
1398: for(theta=1; theta <=npar; theta++){
1399: for(i=1; i<=npar; i++){ /* Computes gradient */
1400: xp[i] = x[i] + (i==theta ?delti[theta]:0);
1401: }
1402: hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);
1403: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1404: for(j=1; j<= nlstate; j++){
1405: for(h=0; h<=nhstepm; h++){
1406: for(i=1, gp[h][j]=0.;i<=nlstate;i++)
1407: gp[h][j] += prlim[i][i]*p3mat[i][j][h];
1408: }
1409: }
1410:
1411: for(i=1; i<=npar; i++) /* Computes gradient */
1412: xp[i] = x[i] - (i==theta ?delti[theta]:0);
1413: hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);
1414: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1415: for(j=1; j<= nlstate; j++){
1416: for(h=0; h<=nhstepm; h++){
1417: for(i=1, gm[h][j]=0.;i<=nlstate;i++)
1418: gm[h][j] += prlim[i][i]*p3mat[i][j][h];
1419: }
1420: }
1421: for(j=1; j<= nlstate; j++)
1422: for(h=0; h<=nhstepm; h++){
1423: gradg[h][theta][j]= (gp[h][j]-gm[h][j])/2./delti[theta];
1424: }
1425: } /* End theta */
1426:
1427: trgradg =ma3x(0,nhstepm,1,nlstate,1,npar);
1428:
1429: for(h=0; h<=nhstepm; h++)
1430: for(j=1; j<=nlstate;j++)
1431: for(theta=1; theta <=npar; theta++)
1432: trgradg[h][j][theta]=gradg[h][theta][j];
1433:
1434: for(i=1;i<=nlstate;i++)
1435: for(j=1;j<=nlstate;j++)
1436: vareij[i][j][(int)age] =0.;
1437: for(h=0;h<=nhstepm;h++){
1438: for(k=0;k<=nhstepm;k++){
1439: matprod2(dnewm,trgradg[h],1,nlstate,1,npar,1,npar,matcov);
1440: matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg[k]);
1441: for(i=1;i<=nlstate;i++)
1442: for(j=1;j<=nlstate;j++)
1443: vareij[i][j][(int)age] += doldm[i][j];
1444: }
1445: }
1446: h=1;
1447: if (stepm >= YEARM) h=stepm/YEARM;
1448: fprintf(ficresvij,"%.0f ",age );
1449: for(i=1; i<=nlstate;i++)
1450: for(j=1; j<=nlstate;j++){
1451: fprintf(ficresvij," %.4f", h*vareij[i][j][(int)age]);
1452: }
1453: fprintf(ficresvij,"\n");
1454: free_matrix(gp,0,nhstepm,1,nlstate);
1455: free_matrix(gm,0,nhstepm,1,nlstate);
1456: free_ma3x(gradg,0,nhstepm,1,npar,1,nlstate);
1457: free_ma3x(trgradg,0,nhstepm,1,nlstate,1,npar);
1458: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1459: } /* End age */
1460:
1461: free_vector(xp,1,npar);
1462: free_matrix(doldm,1,nlstate,1,npar);
1463: free_matrix(dnewm,1,nlstate,1,nlstate);
1464:
1465: }
1466:
1467: /************ Variance of prevlim ******************/
1468: void varprevlim(char fileres[], double **varpl, double **matcov, double x[], double delti[], int nlstate, int stepm, double bage, double fage, double **oldm, double **savm, double **prlim, double ftolpl, int ij)
1469: {
1470: /* Variance of prevalence limit */
1471: /* double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double ** savm,double ftolpl);*/
1472: double **newm;
1473: double **dnewm,**doldm;
1474: int i, j, nhstepm, hstepm;
1475: int k, cptcode;
1476: double *xp;
1477: double *gp, *gm;
1478: double **gradg, **trgradg;
1479: double age,agelim;
1480: int theta;
1481:
1482: fprintf(ficresvpl,"# Standard deviation of prevalences limit\n");
1483: fprintf(ficresvpl,"# Age");
1484: for(i=1; i<=nlstate;i++)
1485: fprintf(ficresvpl," %1d-%1d",i,i);
1486: fprintf(ficresvpl,"\n");
1487:
1488: xp=vector(1,npar);
1489: dnewm=matrix(1,nlstate,1,npar);
1490: doldm=matrix(1,nlstate,1,nlstate);
1491:
1492: hstepm=1*YEARM; /* Every year of age */
1493: hstepm=hstepm/stepm; /* Typically in stepm units, if j= 2 years, = 2/6 months = 4 */
1494: agelim = AGESUP;
1495: for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
1496: nhstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */
1497: if (stepm >= YEARM) hstepm=1;
1498: nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
1499: gradg=matrix(1,npar,1,nlstate);
1500: gp=vector(1,nlstate);
1501: gm=vector(1,nlstate);
1502:
1503: for(theta=1; theta <=npar; theta++){
1504: for(i=1; i<=npar; i++){ /* Computes gradient */
1505: xp[i] = x[i] + (i==theta ?delti[theta]:0);
1506: }
1507: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1508: for(i=1;i<=nlstate;i++)
1509: gp[i] = prlim[i][i];
1510:
1511: for(i=1; i<=npar; i++) /* Computes gradient */
1512: xp[i] = x[i] - (i==theta ?delti[theta]:0);
1513: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1514: for(i=1;i<=nlstate;i++)
1515: gm[i] = prlim[i][i];
1516:
1517: for(i=1;i<=nlstate;i++)
1518: gradg[theta][i]= (gp[i]-gm[i])/2./delti[theta];
1519: } /* End theta */
1520:
1521: trgradg =matrix(1,nlstate,1,npar);
1522:
1523: for(j=1; j<=nlstate;j++)
1524: for(theta=1; theta <=npar; theta++)
1525: trgradg[j][theta]=gradg[theta][j];
1526:
1527: for(i=1;i<=nlstate;i++)
1528: varpl[i][(int)age] =0.;
1529: matprod2(dnewm,trgradg,1,nlstate,1,npar,1,npar,matcov);
1530: matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg);
1531: for(i=1;i<=nlstate;i++)
1532: varpl[i][(int)age] = doldm[i][i]; /* Covariances are useless */
1533:
1534: fprintf(ficresvpl,"%.0f ",age );
1535: for(i=1; i<=nlstate;i++)
1536: fprintf(ficresvpl," %.5f (%.5f)",prlim[i][i],sqrt(varpl[i][(int)age]));
1537: fprintf(ficresvpl,"\n");
1538: free_vector(gp,1,nlstate);
1539: free_vector(gm,1,nlstate);
1540: free_matrix(gradg,1,npar,1,nlstate);
1541: free_matrix(trgradg,1,nlstate,1,npar);
1542: } /* End age */
1543:
1544: free_vector(xp,1,npar);
1545: free_matrix(doldm,1,nlstate,1,npar);
1546: free_matrix(dnewm,1,nlstate,1,nlstate);
1547:
1548: }
1549:
1550:
1551:
1552: /***********************************************/
1553: /**************** Main Program *****************/
1554: /***********************************************/
1555:
1556: /*int main(int argc, char *argv[])*/
1557: int main()
1558: {
1559:
1560: int i,j, k, n=MAXN,iter,m,size,cptcode, aaa, cptcod;
1561: double agedeb, agefin,hf;
1562: double agemin=1.e20, agemax=-1.e20;
1563:
1564: double fret;
1565: double **xi,tmp,delta;
1566:
1567: double dum; /* Dummy variable */
1568: double ***p3mat;
1569: int *indx;
1570: char line[MAXLINE], linepar[MAXLINE];
1571: char title[MAXLINE];
1572: char optionfile[FILENAMELENGTH], datafile[FILENAMELENGTH], filerespl[FILENAMELENGTH];
1573: char fileres[FILENAMELENGTH], filerespij[FILENAMELENGTH], filereso[FILENAMELENGTH];
1574: char filerest[FILENAMELENGTH];
1575: char fileregp[FILENAMELENGTH];
1576: char path[80],pathc[80],pathcd[80],pathtot[80],model[20];
1577: int firstobs=1, lastobs=10;
1578: int sdeb, sfin; /* Status at beginning and end */
1579: int c, h , cpt,l;
1580: int ju,jl, mi;
1581: int i1,j1, k1,jk,aa,bb, stepsize;
1582: int jnais,jdc,jint4,jint1,jint2,jint3,**outcome,**adl,*tab;
1583:
1584: int hstepm, nhstepm;
1585: double bage, fage, age, agelim, agebase;
1586: double ftolpl=FTOL;
1587: double **prlim;
1588: double *severity;
1589: double ***param; /* Matrix of parameters */
1590: double *p;
1591: double **matcov; /* Matrix of covariance */
1592: double ***delti3; /* Scale */
1593: double *delti; /* Scale */
1594: double ***eij, ***vareij;
1595: double **varpl; /* Variances of prevalence limits by age */
1596: double *epj, vepp;
1597: char version[80]="Imach version 62c, May 1999, INED-EUROREVES ";
1598: char *alph[]={"a","a","b","c","d","e"}, str[4];
1599: char z[1]="c", occ;
1600: #include <sys/time.h>
1601: #include <time.h>
1602: char stra[80], strb[80], strc[80], strd[80],stre[80],modelsav[80];
1603: /* long total_usecs;
1604: struct timeval start_time, end_time;
1605:
1606: gettimeofday(&start_time, (struct timezone*)0); */ /* at first time */
1607:
1608:
1609: printf("\nIMACH, Version 0.63");
1610: printf("\nEnter the parameter file name: ");
1611:
1612: #ifdef windows
1613: scanf("%s",pathtot);
1614: getcwd(pathcd, size);
1615: cutv(path,optionfile,pathtot,'\\');
1616: chdir(path);
1617: replace(pathc,path);
1618: #endif
1619: #ifdef unix
1620: scanf("%s",optionfile);
1621: #endif
1622:
1623: /*-------- arguments in the command line --------*/
1624:
1625: strcpy(fileres,"r");
1626: strcat(fileres, optionfile);
1627:
1628: /*---------arguments file --------*/
1629:
1630: if((ficpar=fopen(optionfile,"r"))==NULL) {
1631: printf("Problem with optionfile %s\n",optionfile);
1632: goto end;
1633: }
1634:
1635: strcpy(filereso,"o");
1636: strcat(filereso,fileres);
1637: if((ficparo=fopen(filereso,"w"))==NULL) {
1638: printf("Problem with Output resultfile: %s\n", filereso);goto end;
1639: }
1640:
1641: /* Reads comments: lines beginning with '#' */
1642: while((c=getc(ficpar))=='#' && c!= EOF){
1643: ungetc(c,ficpar);
1644: fgets(line, MAXLINE, ficpar);
1645: puts(line);
1646: fputs(line,ficparo);
1647: }
1648: ungetc(c,ficpar);
1649:
1650: fscanf(ficpar,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%lf stepm=%d ncov=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n",title, datafile, &lastobs, &firstpass,&lastpass,&ftol, &stepm, &ncov, &nlstate,&ndeath, &maxwav, &mle, &weightopt,model);
1651: printf("title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncov=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol, stepm, ncov, nlstate,ndeath, maxwav, mle, weightopt,model);
1652: fprintf(ficparo,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncov=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol,stepm,ncov,nlstate,ndeath,maxwav, mle, weightopt,model);
1653:
1654: covar=matrix(1,NCOVMAX,1,n);
1655: if (strlen(model)<=1) cptcovn=0;
1656: else {
1657: j=0;
1658: j=nbocc(model,'+');
1659: cptcovn=j+1;
1660: }
1661:
1662: ncovmodel=2+cptcovn;
1663: nvar=ncovmodel-1; /* Suppressing age as a basic covariate */
1664:
1665: /* Read guess parameters */
1666: /* Reads comments: lines beginning with '#' */
1667: while((c=getc(ficpar))=='#' && c!= EOF){
1668: ungetc(c,ficpar);
1669: fgets(line, MAXLINE, ficpar);
1670: puts(line);
1671: fputs(line,ficparo);
1672: }
1673: ungetc(c,ficpar);
1674:
1675: param= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
1676: for(i=1; i <=nlstate; i++)
1677: for(j=1; j <=nlstate+ndeath-1; j++){
1678: fscanf(ficpar,"%1d%1d",&i1,&j1);
1679: fprintf(ficparo,"%1d%1d",i1,j1);
1680: printf("%1d%1d",i,j);
1681: for(k=1; k<=ncovmodel;k++){
1682: fscanf(ficpar," %lf",¶m[i][j][k]);
1683: printf(" %lf",param[i][j][k]);
1684: fprintf(ficparo," %lf",param[i][j][k]);
1685: }
1686: fscanf(ficpar,"\n");
1687: printf("\n");
1688: fprintf(ficparo,"\n");
1689: }
1690:
1691: npar= (nlstate+ndeath-1)*nlstate*ncovmodel;
1692: p=param[1][1];
1693:
1694: /* Reads comments: lines beginning with '#' */
1695: while((c=getc(ficpar))=='#' && c!= EOF){
1696: ungetc(c,ficpar);
1697: fgets(line, MAXLINE, ficpar);
1698: puts(line);
1699: fputs(line,ficparo);
1700: }
1701: ungetc(c,ficpar);
1702:
1703: delti3= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
1704: delti=vector(1,npar); /* Scale of each paramater (output from hesscov) */
1705: for(i=1; i <=nlstate; i++){
1706: for(j=1; j <=nlstate+ndeath-1; j++){
1707: fscanf(ficpar,"%1d%1d",&i1,&j1);
1708: printf("%1d%1d",i,j);
1709: fprintf(ficparo,"%1d%1d",i1,j1);
1710: for(k=1; k<=ncovmodel;k++){
1711: fscanf(ficpar,"%le",&delti3[i][j][k]);
1712: printf(" %le",delti3[i][j][k]);
1713: fprintf(ficparo," %le",delti3[i][j][k]);
1714: }
1715: fscanf(ficpar,"\n");
1716: printf("\n");
1717: fprintf(ficparo,"\n");
1718: }
1719: }
1720: delti=delti3[1][1];
1721:
1722: /* Reads comments: lines beginning with '#' */
1723: while((c=getc(ficpar))=='#' && c!= EOF){
1724: ungetc(c,ficpar);
1725: fgets(line, MAXLINE, ficpar);
1726: puts(line);
1727: fputs(line,ficparo);
1728: }
1729: ungetc(c,ficpar);
1730:
1731: matcov=matrix(1,npar,1,npar);
1732: for(i=1; i <=npar; i++){
1733: fscanf(ficpar,"%s",&str);
1734: printf("%s",str);
1735: fprintf(ficparo,"%s",str);
1736: for(j=1; j <=i; j++){
1737: fscanf(ficpar," %le",&matcov[i][j]);
1738: printf(" %.5le",matcov[i][j]);
1739: fprintf(ficparo," %.5le",matcov[i][j]);
1740: }
1741: fscanf(ficpar,"\n");
1742: printf("\n");
1743: fprintf(ficparo,"\n");
1744: }
1745: for(i=1; i <=npar; i++)
1746: for(j=i+1;j<=npar;j++)
1747: matcov[i][j]=matcov[j][i];
1748:
1749: printf("\n");
1750:
1751:
1752: if(mle==1){
1753: /*-------- data file ----------*/
1754: if((ficres =fopen(fileres,"w"))==NULL) {
1755: printf("Problem with resultfile: %s\n", fileres);goto end;
1756: }
1757: fprintf(ficres,"#%s\n",version);
1758:
1759: if((fic=fopen(datafile,"r"))==NULL) {
1760: printf("Problem with datafile: %s\n", datafile);goto end;
1761: }
1762:
1763: n= lastobs;
1764: severity = vector(1,maxwav);
1765: outcome=imatrix(1,maxwav+1,1,n);
1766: num=ivector(1,n);
1767: moisnais=vector(1,n);
1768: annais=vector(1,n);
1769: moisdc=vector(1,n);
1770: andc=vector(1,n);
1771: agedc=vector(1,n);
1772: cod=ivector(1,n);
1773: weight=vector(1,n);
1774: for(i=1;i<=n;i++) weight[i]=1.0; /* Equal weights, 1 by default */
1775: mint=matrix(1,maxwav,1,n);
1776: anint=matrix(1,maxwav,1,n);
1777: s=imatrix(1,maxwav+1,1,n);
1778: adl=imatrix(1,maxwav+1,1,n);
1779: tab=ivector(1,NCOVMAX);
1780: ncodemax=ivector(1,NCOVMAX);
1781:
1782: i=1;
1783: while (fgets(line, MAXLINE, fic) != NULL) {
1784: if ((i >= firstobs) && (i <=lastobs)) {
1785:
1786: for (j=maxwav;j>=1;j--){
1787: cutv(stra, strb,line,' '); s[j][i]=atoi(strb);
1788: strcpy(line,stra);
1789: cutv(stra, strb,line,'/'); anint[j][i]=(double)(atoi(strb)); strcpy(line,stra);
1790: cutv(stra, strb,line,' '); mint[j][i]=(double)(atoi(strb)); strcpy(line,stra);
1791: }
1792:
1793: cutv(stra, strb,line,'/'); andc[i]=(double)(atoi(strb)); strcpy(line,stra);
1794: cutv(stra, strb,line,' '); moisdc[i]=(double)(atoi(strb)); strcpy(line,stra);
1795:
1796: cutv(stra, strb,line,'/'); annais[i]=(double)(atoi(strb)); strcpy(line,stra);
1797: cutv(stra, strb,line,' '); moisnais[i]=(double)(atoi(strb)); strcpy(line,stra);
1798:
1799: cutv(stra, strb,line,' '); weight[i]=(double)(atoi(strb)); strcpy(line,stra);
1800: for (j=ncov;j>=1;j--){
1801: cutv(stra, strb,line,' '); covar[j][i]=(double)(atoi(strb)); strcpy(line,stra);
1802: }
1803: num[i]=atol(stra);
1804:
1805: /* printf("%d %.lf %.lf %.lf %.lf/%.lf %.lf/%.lf %.lf/%.lf %d %.lf/%.lf %d %.lf/%.lf %d %.lf/%.lf %d\n",num[i],(covar[1][i]), (covar[2][i]), (weight[i]), (moisnais[i]), (annais[i]), (moisdc[i]), (andc[i]), (mint[1][i]), (anint[1][i]), (s[1][i]), (mint[2][i]), (anint[2][i]), (s[2][i]), (mint[3][i]), (anint[3][i]), (s[3][i]), (mint[4][i]), (anint[4][i]), (s[4][i]));*/
1806:
1807: /*printf("%d %.lf %.lf %.lf %.lf %.lf %.lf/%.lf %.lf/%.lf %.lf/%.lf %d %.lf/%.lf %d\n",num[i],(covar[1][i]), (covar[2][i]),(covar[3][i]), (covar[4][i]), (weight[i]), (moisnais[i]), (annais[i]), (moisdc[i]), (andc[i]), (mint[1][i]), (anint[1][i]), (s[1][i]), (mint[2][i]), (anint[2][i]), (s[2][i]));*/
1808:
1809: i=i+1;
1810: }
1811: }
1812: /*scanf("%d",i);*/
1813:
1814: imx=i-1; /* Number of individuals */
1815:
1816: /* Calculation of the number of parameter from char model*/
1817: Tvar=ivector(1,8);
1818:
1819: if (strlen(model) >1){
1820: j=0;
1821: j=nbocc(model,'+');
1822: cptcovn=j+1;
1823:
1824: strcpy(modelsav,model);
1825: if (j==0) {
1826: cutv(stra,strb,modelsav,'V'); Tvar[1]=atoi(strb);
1827: }
1828: else {
1829: for(i=j; i>=1;i--){
1830: cutv(stra,strb,modelsav,'+');
1831: if (strchr(strb,'*')) {
1832: cutv(strd,strc,strb,'*');
1833: cutv(strb,stre,strc,'V');Tvar[i+1]=ncov+1;
1834: cutv(strb,strc,strd,'V');
1835: for (k=1; k<=lastobs;k++)
1836: covar[ncov+1][k]=covar[atoi(stre)][k]*covar[atoi(strc)][k];
1837: }
1838: else {
1839: cutv(strd,strc,strb,'V');
1840: Tvar[i+1]=atoi(strc);
1841: }
1842: strcpy(modelsav,stra);
1843: }
1844: /*cutv(strd,strc,stra,'V');*/
1845: Tvar[1]=atoi(strc);
1846: }
1847: }
1848: /*printf("tvar=%d ",Tvar[1]);*/
1849: /*scanf("%d ",i);*/
1850: fclose(fic);
1851:
1852: if (weightopt != 1) { /* Maximisation without weights*/
1853: for(i=1;i<=n;i++) weight[i]=1.0;
1854: }
1855: /*-calculation of age at interview from date of interview and age at death -*/
1856: agev=matrix(1,maxwav,1,imx);
1857:
1858: for (i=1; i<=imx; i++) {
1859: agedc[i]=(moisdc[i]/12.+andc[i])-(moisnais[i]/12.+annais[i]);
1860: for(m=1; (m<= maxwav); m++){
1861: if (mint[m][i]==99 || anint[m][i]==9999) s[m][i]=-1;
1862: if(s[m][i] >0){
1863: if (s[m][i] == nlstate+1) {
1864: if(agedc[i]>0)
1865: if(moisdc[i]!=99 && andc[i]!=9999)
1866: agev[m][i]=agedc[i];
1867: else{
1868: printf("Warning negative age at death: %d line:%d\n",num[i],i);
1869: agev[m][i]=-1;
1870: }
1871: }
1872: else if(s[m][i] !=9){ /* Should no more exist */
1873: agev[m][i]=(mint[m][i]/12.+1./24.+anint[m][i])-(moisnais[i]/12.+1./24.+annais[i]);
1874: if(mint[m][i]==99 || anint[m][i]==9999){
1875: agev[m][i]=1;
1876: /* printf("i=%d m=%d agev=%lf \n",i,m, agev[m][i]); */
1877: }
1878: else if(agev[m][i] <agemin){
1879: agemin=agev[m][i];
1880: /*printf(" Min anint[%d][%d]=%.2f annais[%d]=%.2f, agemin=%.2f\n",m,i,anint[m][i], i,annais[i], agemin);*/
1881: }
1882: else if(agev[m][i] >agemax){
1883: agemax=agev[m][i];
1884: /* printf(" anint[%d][%d]=%.0f annais[%d]=%.0f, agemax=%.0f\n",m,i,anint[m][i], i,annais[i], agemax);*/
1885: }
1886: /*agev[m][i]=anint[m][i]-annais[i];*/
1887: /* agev[m][i] = age[i]+2*m;*/
1888: }
1889: else { /* =9 */
1890: agev[m][i]=1;
1891: s[m][i]=-1;
1892: }
1893: }
1894: else /*= 0 Unknown */
1895: agev[m][i]=1;
1896: }
1897:
1898: }
1899: for (i=1; i<=imx; i++) {
1900: for(m=1; (m<= maxwav); m++){
1901: if (s[m][i] > (nlstate+ndeath)) {
1902: printf("Error: Wrong value in nlstate or ndeath\n");
1903: goto end;
1904: }
1905: }
1906: }
1907:
1908: printf("Total number of individuals= %d, Agemin = %.2f, Agemax= %.2f\n\n", imx, agemin, agemax);
1909:
1910: free_vector(severity,1,maxwav);
1911: free_imatrix(outcome,1,maxwav+1,1,n);
1912: free_vector(moisnais,1,n);
1913: free_vector(annais,1,n);
1914: free_matrix(mint,1,maxwav,1,n);
1915: free_matrix(anint,1,maxwav,1,n);
1916: free_vector(moisdc,1,n);
1917: free_vector(andc,1,n);
1918:
1919:
1920: wav=ivector(1,imx);
1921: dh=imatrix(1,lastpass-firstpass+1,1,imx);
1922: mw=imatrix(1,lastpass-firstpass+1,1,imx);
1923:
1924: /* Concatenates waves */
1925: concatwav(wav, dh, mw, s, agedc, agev, firstpass, lastpass, imx, nlstate, stepm);
1926:
1927:
1928: Tcode=ivector(1,100);
1929: nbcode=imatrix(1,nvar,1,8);
1930: ncodemax[1]=1;
1931: if (cptcovn > 0) tricode(Tvar,nbcode,imx);
1932:
1933: codtab=imatrix(1,100,1,10);
1934: h=0;
1935: m=pow(2,cptcovn);
1936:
1937: for(k=1;k<=cptcovn; k++){
1938: for(i=1; i <=(m/pow(2,k));i++){
1939: for(j=1; j <= ncodemax[k]; j++){
1940: for(cpt=1; cpt <=(m/pow(2,cptcovn+1-k)); cpt++){
1941: h++;
1942: if (h>m) h=1;codtab[h][k]=j;
1943: }
1944: }
1945: }
1946: }
1947:
1948: /*for(i=1; i <=m ;i++){
1949: for(k=1; k <=cptcovn; k++){
1950: printf("i=%d k=%d %d ",i,k,codtab[i][k]);
1951: }
1952: printf("\n");
1953: }
1954: scanf("%d",i);*/
1955:
1956: /* Calculates basic frequencies. Computes observed prevalence at single age
1957: and prints on file fileres'p'. */
1958: freqsummary(fileres, agemin, agemax, s, agev, nlstate, imx,Tvar,nbcode, ncodemax);
1959:
1960: pmmij= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
1961: oldms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
1962: newms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
1963: savms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
1964: oldm=oldms; newm=newms; savm=savms; /* Keeps fixed addresses to free */
1965:
1966: /* For Powell, parameters are in a vector p[] starting at p[1]
1967: so we point p on param[1][1] so that p[1] maps on param[1][1][1] */
1968: p=param[1][1]; /* *(*(*(param +1)+1)+0) */
1969: /*scanf("%d",i);*/
1970: mlikeli(ficres,p, npar, ncovmodel, nlstate, ftol, func);
1971:
1972:
1973: /*--------- results files --------------*/
1974: fprintf(ficres,"\ntitle=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncov=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol, stepm, ncov, nlstate, ndeath, maxwav, mle,weightopt,model);
1975:
1976: jk=1;
1977: fprintf(ficres,"# Parameters\n");
1978: printf("# Parameters\n");
1979: for(i=1,jk=1; i <=nlstate; i++){
1980: for(k=1; k <=(nlstate+ndeath); k++){
1981: if (k != i)
1982: {
1983: printf("%d%d ",i,k);
1984: fprintf(ficres,"%1d%1d ",i,k);
1985: for(j=1; j <=ncovmodel; j++){
1986: printf("%f ",p[jk]);
1987: fprintf(ficres,"%f ",p[jk]);
1988: jk++;
1989: }
1990: printf("\n");
1991: fprintf(ficres,"\n");
1992: }
1993: }
1994: }
1995:
1996: /* Computing hessian and covariance matrix */
1997: ftolhess=ftol; /* Usually correct */
1998: hesscov(matcov, p, npar, delti, ftolhess, func);
1999: fprintf(ficres,"# Scales\n");
2000: printf("# Scales\n");
2001: for(i=1,jk=1; i <=nlstate; i++){
2002: for(j=1; j <=nlstate+ndeath; j++){
2003: if (j!=i) {
2004: fprintf(ficres,"%1d%1d",i,j);
2005: printf("%1d%1d",i,j);
2006: for(k=1; k<=ncovmodel;k++){
2007: printf(" %.5e",delti[jk]);
2008: fprintf(ficres," %.5e",delti[jk]);
2009: jk++;
2010: }
2011: printf("\n");
2012: fprintf(ficres,"\n");
2013: }
2014: }
2015: }
2016:
2017: k=1;
2018: fprintf(ficres,"# Covariance\n");
2019: printf("# Covariance\n");
2020: for(i=1;i<=npar;i++){
2021: /* if (k>nlstate) k=1;
2022: i1=(i-1)/(ncovmodel*nlstate)+1;
2023: fprintf(ficres,"%s%d%d",alph[k],i1,tab[i]);
2024: printf("%s%d%d",alph[k],i1,tab[i]);*/
2025: fprintf(ficres,"%3d",i);
2026: printf("%3d",i);
2027: for(j=1; j<=i;j++){
2028: fprintf(ficres," %.5e",matcov[i][j]);
2029: printf(" %.5e",matcov[i][j]);
2030: }
2031: fprintf(ficres,"\n");
2032: printf("\n");
2033: k++;
2034: }
2035:
2036: while((c=getc(ficpar))=='#' && c!= EOF){
2037: ungetc(c,ficpar);
2038: fgets(line, MAXLINE, ficpar);
2039: puts(line);
2040: fputs(line,ficparo);
2041: }
2042: ungetc(c,ficpar);
2043:
2044: fscanf(ficpar,"agemin=%lf agemax=%lf bage=%lf fage=%lf\n",&agemin,&agemax, &bage, &fage);
2045:
2046: if (fage <= 2) {
2047: bage = agemin;
2048: fage = agemax;
2049: }
2050:
2051: fprintf(ficres,"# agemin agemax for life expectancy, bage fage (if mle==0 ie no data nor Max likelihood).\n");
2052: fprintf(ficres,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f\n",agemin,agemax,bage,fage);
2053: /*------------ gnuplot -------------*/
2054: chdir(pathcd);
2055: if((ficgp=fopen("graph.plt","w"))==NULL) {
2056: printf("Problem with file graph.gp");goto end;
2057: }
2058: #ifdef windows
2059: fprintf(ficgp,"cd \"%s\" \n",pathc);
2060: #endif
2061: m=pow(2,cptcovn);
2062:
2063: /* 1eme*/
2064: for (cpt=1; cpt<= nlstate ; cpt ++) {
2065: for (k1=1; k1<= m ; k1 ++) {
2066:
2067: #ifdef windows
2068: fprintf(ficgp,"set xlabel \"Age\" \nset ylabel \"Probability\" \nset ter gif small size 400,300\nplot [%.f:%.f] \"vpl%s\" every :::%d::%d u 1:2 \"\%%lf",agemin,fage,fileres,k1-1,k1-1);
2069: #endif
2070: #ifdef unix
2071: fprintf(ficgp,"set xlabel \"Age\" \nset ylabel \"Probability\" \nplot [%.f:%.f] \"vpl%s\" u 1:2 \"\%%lf",agemin,fage,fileres);
2072: #endif
2073:
2074: for (i=1; i<= nlstate ; i ++) {
2075: if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
2076: else fprintf(ficgp," \%%*lf (\%%*lf)");
2077: }
2078: fprintf(ficgp,"\" t\"Stationary prevalence\" w l 0,\"vpl%s\" every :::%d::%d u 1:($2+2*$3) \"\%%lf",fileres,k1-1,k1-1);
2079: for (i=1; i<= nlstate ; i ++) {
2080: if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
2081: else fprintf(ficgp," \%%*lf (\%%*lf)");
2082: }
2083: fprintf(ficgp,"\" t\"95\%% CI\" w l 1,\"vpl%s\" every :::%d::%d u 1:($2-2*$3) \"\%%lf",fileres,k1-1,k1-1);
2084: for (i=1; i<= nlstate ; i ++) {
2085: if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
2086: else fprintf(ficgp," \%%*lf (\%%*lf)");
2087: }
2088: fprintf(ficgp,"\" t\"\" w l 1,\"p%s\" every :::%d::%d u 1:($%d) t\"Observed prevalence \" w l 2",fileres,k1-1,k1-1,2+4*(cpt-1));
2089: #ifdef unix
2090: fprintf(ficgp,"\nset ter gif small size 400,300");
2091: #endif
2092: fprintf(ficgp,"\nset out \"v%s%d%d.gif\" \nreplot\n\n",strtok(optionfile, "."),cpt,k1);
2093: }
2094: }
2095: /*2 eme*/
2096:
2097: for (k1=1; k1<= m ; k1 ++) {
2098: fprintf(ficgp,"set ylabel \"Years\" \nset ter gif small size 400,300\nplot [%.f:%.f] ",agemin,fage);
2099:
2100: for (i=1; i<= nlstate+1 ; i ++) {
2101: k=2*i;
2102: fprintf(ficgp,"\"t%s\" every :::%d::%d u 1:2 \"\%%lf",fileres,k1-1,k1-1);
2103: for (j=1; j<= nlstate+1 ; j ++) {
2104: if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
2105: else fprintf(ficgp," \%%*lf (\%%*lf)");
2106: }
2107: if (i== 1) fprintf(ficgp,"\" t\"TLE\" w l ,");
2108: else fprintf(ficgp,"\" t\"LE in state (%d)\" w l ,",i-1);
2109: fprintf(ficgp,"\"t%s\" every :::%d::%d u 1:($2-$3*2) \"\%%lf",fileres,k1-1,k1-1);
2110: for (j=1; j<= nlstate+1 ; j ++) {
2111: if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
2112: else fprintf(ficgp," \%%*lf (\%%*lf)");
2113: }
2114: fprintf(ficgp,"\" t\"\" w l 0,");
2115: fprintf(ficgp,"\"t%s\" every :::%d::%d u 1:($2+$3*2) \"\%%lf",fileres,k1-1,k1-1);
2116: for (j=1; j<= nlstate+1 ; j ++) {
2117: if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
2118: else fprintf(ficgp," \%%*lf (\%%*lf)");
2119: }
2120: if (i== (nlstate+1)) fprintf(ficgp,"\" t\"\" w l 0");
2121: else fprintf(ficgp,"\" t\"\" w l 0,");
2122: }
2123: fprintf(ficgp,"\nset out \"e%s%d.gif\" \nreplot\n\n",strtok(optionfile, "."),k1);
2124: }
2125:
2126: /*3eme*/
2127:
2128: for (k1=1; k1<= m ; k1 ++) {
2129: for (cpt=1; cpt<= nlstate ; cpt ++) {
2130: k=2+nlstate*(cpt-1);
2131: fprintf(ficgp,"set ter gif small size 400,300\nplot [%.f:%.f] \"e%s\" every :::%d::%d u 1:%d t \"e%d1\" w l",agemin,fage,fileres,k1-1,k1-1,k,cpt);
2132: for (i=1; i< nlstate ; i ++) {
2133: fprintf(ficgp,",\"e%s\" every :::%d::%d u 1:%d t \"e%d%d\" w l",fileres,k1-1,k1-1,k+i,cpt,i+1);
2134: }
2135: fprintf(ficgp,"\nset out \"exp%s%d%d.gif\" \nreplot\n\n",strtok(optionfile, "."),cpt,k1);
2136: }
2137: }
2138:
2139: /* CV preval stat */
2140: for (k1=1; k1<= m ; k1 ++) {
2141: for (cpt=1; cpt<nlstate ; cpt ++) {
2142: k=3;
2143: fprintf(ficgp,"set xlabel \"Age\" \nset ylabel \"Probability\" \nset ter gif small size 400,300\nplot [%.f:%.f] \"pij%s\" u ($1==%d ? ($3):1/0):($%d/($%d",agemin,agemax,fileres,k1,k+cpt+1,k+1);
2144: for (i=1; i< nlstate ; i ++)
2145: fprintf(ficgp,"+$%d",k+i+1);
2146: fprintf(ficgp,")) t\"prev(%d,%d)\" w l",cpt,cpt+1);
2147:
2148: l=3+(nlstate+ndeath)*cpt;
2149: fprintf(ficgp,",\"pij%s\" u ($1==%d ? ($3):1/0):($%d/($%d",fileres,k1,l+cpt+1,l+1);
2150: for (i=1; i< nlstate ; i ++) {
2151: l=3+(nlstate+ndeath)*cpt;
2152: fprintf(ficgp,"+$%d",l+i+1);
2153: }
2154: fprintf(ficgp,")) t\"prev(%d,%d)\" w l\n",cpt+1,cpt+1);
2155: fprintf(ficgp,"set out \"p%s%d%d.gif\" \nreplot\n\n",strtok(optionfile, "."),cpt,k1);
2156: }
2157: }
2158:
2159: /* proba elementaires */
2160: for(i=1,jk=1; i <=nlstate; i++){
2161: for(k=1; k <=(nlstate+ndeath); k++){
2162: if (k != i) {
2163: /* fprintf(ficgp,"%1d%1d ",i,k);*/
2164: for(j=1; j <=ncovmodel; j++){
2165: fprintf(ficgp,"%s%1d%1d=%f ",alph[j],i,k,p[jk]);
2166: jk++;
2167: fprintf(ficgp,"\n");
2168: }
2169: }
2170: }
2171: }
2172: for(jk=1; jk <=m; jk++) {
2173: fprintf(ficgp,"\nset ter gif small size 400,300\nset log y\nplot [%.f:%.f] ",agemin,agemax);
2174: for(i=1; i <=nlstate; i++) {
2175: for(k=1; k <=(nlstate+ndeath); k++){
2176: if (k != i) {
2177: fprintf(ficgp," exp(a%d%d+b%d%d*x",i,k,i,k);
2178: for(j=3; j <=ncovmodel; j++)
2179: fprintf(ficgp,"+%s%d%d*%d",alph[j],i,k,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
2180: fprintf(ficgp,")/(1");
2181: for(k1=1; k1 <=(nlstate+ndeath); k1++)
2182: if (k1 != i) {
2183: fprintf(ficgp,"+exp(a%d%d+b%d%d*x",i,k1,i,k1);
2184: for(j=3; j <=ncovmodel; j++)
2185: fprintf(ficgp,"+%s%d%d*%d",alph[j],i,k,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
2186: fprintf(ficgp,")");
2187: }
2188: fprintf(ficgp,") t \"p%d%d\" ", i,k);
2189: if ((i+k)!= (nlstate*2+ndeath)) fprintf(ficgp,",");
2190: }
2191: }
2192: }
2193: fprintf(ficgp,"\nset out \"pe%s%d.gif\" \nreplot\n\n",strtok(optionfile, "."),jk);
2194: }
2195: fclose(ficgp);
2196:
2197: chdir(path);
2198: free_matrix(agev,1,maxwav,1,imx);
2199: free_ivector(wav,1,imx);
2200: free_imatrix(dh,1,lastpass-firstpass+1,1,imx);
2201: free_imatrix(mw,1,lastpass-firstpass+1,1,imx);
2202:
2203: free_imatrix(s,1,maxwav+1,1,n);
2204:
2205:
2206: free_ivector(num,1,n);
2207: free_vector(agedc,1,n);
2208: free_vector(weight,1,n);
2209: /*free_matrix(covar,1,NCOVMAX,1,n);*/
2210: fclose(ficparo);
2211: fclose(ficres);
2212: }
2213:
2214: /*________fin mle=1_________*/
2215:
2216:
2217:
2218: /* No more information from the sample is required now */
2219: /* Reads comments: lines beginning with '#' */
2220: while((c=getc(ficpar))=='#' && c!= EOF){
2221: ungetc(c,ficpar);
2222: fgets(line, MAXLINE, ficpar);
2223: puts(line);
2224: fputs(line,ficparo);
2225: }
2226: ungetc(c,ficpar);
2227:
2228: fscanf(ficpar,"agemin=%lf agemax=%lf bage=%lf fage=%lf\n",&agemin,&agemax, &bage, &fage);
2229: printf("agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f\n",agemin,agemax, bage, fage);
2230: fprintf(ficparo,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f\n",agemin,agemax,bage,fage);
2231: /*--------- index.htm --------*/
2232:
2233: if((fichtm=fopen("index.htm","w"))==NULL) {
2234: printf("Problem with index.htm \n");goto end;
2235: }
2236:
2237: fprintf(fichtm,"<body><ul> Imach, Version 0.63<hr> <li>Outputs files<br><br>\n
2238: - Observed prevalence in each state: <a href=\"p%s\">p%s</a> <br>\n
2239: - Estimated parameters and the covariance matrix: <a href=\"%s\">%s</a> <br>
2240: - Stationary prevalence in each state: <a href=\"pl%s\">pl%s</a> <br>
2241: - Transition probabilities: <a href=\"pij%s\">pij%s</a><br>
2242: - Copy of the parameter file: <a href=\"o%s\">o%s</a><br>
2243: - Life expectancies by age and initial health status: <a href=\"e%s\">e%s</a> <br>
2244: - Variances of life expectancies by age and initial health status: <a href=\"v%s\">v%s</a><br>
2245: - Health expectancies with their variances: <a href=\"t%s\">t%s</a> <br>
2246: - Standard deviation of stationary prevalences: <a href=\"vpl%s\">vpl%s</a> <br><br>",fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres,fileres);
2247:
2248: fprintf(fichtm," <li>Graphs</li>\n<p>");
2249:
2250: m=cptcovn;
2251: if (cptcovn < 1) {m=1;ncodemax[1]=1;}
2252:
2253: j1=0;
2254: for(k1=1; k1<=m;k1++){
2255: for(i1=1; i1<=ncodemax[k1];i1++){
2256: j1++;
2257: if (cptcovn > 0) {
2258: fprintf(fichtm,"<hr>************ Results for covariates");
2259: for (cpt=1; cpt<=cptcovn;cpt++)
2260: fprintf(fichtm," V%d=%d ",Tvar[cpt],nbcode[Tvar[cpt]][codtab[j1][cpt]]);
2261: fprintf(fichtm," ************\n<hr>");
2262: }
2263: fprintf(fichtm,"<br>- Probabilities: pe%s%d.gif<br>
2264: <img src=\"pe%s%d.gif\">",strtok(optionfile, "."),j1,strtok(optionfile, "."),j1);
2265: for(cpt=1; cpt<nlstate;cpt++){
2266: fprintf(fichtm,"<br>- Prevalence of disability : p%s%d%d.gif<br>
2267: <img src=\"p%s%d%d.gif\">",strtok(optionfile, "."),cpt,j1,strtok(optionfile, "."),cpt,j1);
2268: }
2269: for(cpt=1; cpt<=nlstate;cpt++) {
2270: fprintf(fichtm,"<br>- Observed and stationary prevalence (with confident
2271: interval) in state (%d): v%s%d%d.gif <br>
2272: <img src=\"v%s%d%d.gif\">",cpt,strtok(optionfile, "."),cpt,j1,strtok(optionfile, "."),cpt,j1);
2273: }
2274: for(cpt=1; cpt<=nlstate;cpt++) {
2275: fprintf(fichtm,"\n<br>- Health life expectancies by age and initial health state (%d): exp%s%d%d.gif <br>
2276: <img src=\"exp%s%d%d.gif\">",cpt,strtok(optionfile, "."),cpt,j1,strtok(optionfile, "."),cpt,j1);
2277: }
2278: fprintf(fichtm,"\n<br>- Total life expectancy by age and
2279: health expectancies in states (1) and (2): e%s%d.gif<br>
2280: <img src=\"e%s%d.gif\">",strtok(optionfile, "."),j1,strtok(optionfile, "."),j1);
2281: fprintf(fichtm,"\n</body>");
2282: }
2283: }
2284: fclose(fichtm);
2285:
2286: /*--------------- Prevalence limit --------------*/
2287:
2288: strcpy(filerespl,"pl");
2289: strcat(filerespl,fileres);
2290: if((ficrespl=fopen(filerespl,"w"))==NULL) {
2291: printf("Problem with Prev limit resultfile: %s\n", filerespl);goto end;
2292: }
2293: printf("Computing prevalence limit: result on file '%s' \n", filerespl);
2294: fprintf(ficrespl,"#Prevalence limit\n");
2295: fprintf(ficrespl,"#Age ");
2296: for(i=1; i<=nlstate;i++) fprintf(ficrespl,"%d-%d ",i,i);
2297: fprintf(ficrespl,"\n");
2298:
2299: prlim=matrix(1,nlstate,1,nlstate);
2300: pmmij= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2301: oldms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2302: newms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2303: savms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2304: oldm=oldms; newm=newms; savm=savms; /* Keeps fixed addresses to free */
2305: k=0;
2306: agebase=agemin;
2307: agelim=agemax;
2308: ftolpl=1.e-10;
2309: i1=cptcovn;
2310: if (cptcovn < 1){i1=1;}
2311:
2312: for(cptcov=1;cptcov<=i1;cptcov++){
2313: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
2314: k=k+1;
2315: /*printf("cptcov=%d cptcod=%d codtab=%d nbcode=%d\n",cptcov, cptcod,Tcode[cptcode],codtab[cptcod][cptcov]);*/
2316: fprintf(ficrespl,"\n#****** ");
2317: for(j=1;j<=cptcovn;j++)
2318: fprintf(ficrespl,"V%d=%d ",Tvar[j],nbcode[Tvar[j]][codtab[k][j]]);
2319: fprintf(ficrespl,"******\n");
2320:
2321: for (age=agebase; age<=agelim; age++){
2322: prevalim(prlim, nlstate, p, age, oldm, savm,ftolpl,k);
2323: fprintf(ficrespl,"%.0f",age );
2324: for(i=1; i<=nlstate;i++)
2325: fprintf(ficrespl," %.5f", prlim[i][i]);
2326: fprintf(ficrespl,"\n");
2327: }
2328: }
2329: }
2330: fclose(ficrespl);
2331: /*------------- h Pij x at various ages ------------*/
2332:
2333: strcpy(filerespij,"pij"); strcat(filerespij,fileres);
2334: if((ficrespij=fopen(filerespij,"w"))==NULL) {
2335: printf("Problem with Pij resultfile: %s\n", filerespij);goto end;
2336: }
2337: printf("Computing pij: result on file '%s' \n", filerespij);
2338:
2339: stepsize=(int) (stepm+YEARM-1)/YEARM;
2340: if (stepm<=24) stepsize=2;
2341:
2342: agelim=AGESUP;
2343: hstepm=stepsize*YEARM; /* Every year of age */
2344: hstepm=hstepm/stepm; /* Typically 2 years, = 2/6 months = 4 */
2345:
2346: k=0;
2347: for(cptcov=1;cptcov<=i1;cptcov++){
2348: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
2349: k=k+1;
2350: fprintf(ficrespij,"\n#****** ");
2351: for(j=1;j<=cptcovn;j++)
2352: fprintf(ficrespij,"V%d=%d ",Tvar[j],nbcode[Tvar[j]][codtab[k][j]]);
2353: fprintf(ficrespij,"******\n");
2354:
2355: for (agedeb=fage; agedeb>=bage; agedeb--){ /* If stepm=6 months */
2356: nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */
2357: nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
2358: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2359: oldm=oldms;savm=savms;
2360: hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);
2361: fprintf(ficrespij,"# Age");
2362: for(i=1; i<=nlstate;i++)
2363: for(j=1; j<=nlstate+ndeath;j++)
2364: fprintf(ficrespij," %1d-%1d",i,j);
2365: fprintf(ficrespij,"\n");
2366: for (h=0; h<=nhstepm; h++){
2367: fprintf(ficrespij,"%d %.0f %.0f",k,agedeb, agedeb+ h*hstepm/YEARM*stepm );
2368: for(i=1; i<=nlstate;i++)
2369: for(j=1; j<=nlstate+ndeath;j++)
2370: fprintf(ficrespij," %.5f", p3mat[i][j][h]);
2371: fprintf(ficrespij,"\n");
2372: }
2373: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2374: fprintf(ficrespij,"\n");
2375: }
2376: }
2377: }
2378:
2379: fclose(ficrespij);
2380:
2381: /*---------- Health expectancies and variances ------------*/
2382:
2383: strcpy(filerest,"t");
2384: strcat(filerest,fileres);
2385: if((ficrest=fopen(filerest,"w"))==NULL) {
2386: printf("Problem with total LE resultfile: %s\n", filerest);goto end;
2387: }
2388: printf("Computing Total LEs with variances: file '%s' \n", filerest);
2389:
2390:
2391: strcpy(filerese,"e");
2392: strcat(filerese,fileres);
2393: if((ficreseij=fopen(filerese,"w"))==NULL) {
2394: printf("Problem with Health Exp. resultfile: %s\n", filerese); exit(0);
2395: }
2396: printf("Computing Health Expectancies: result on file '%s' \n", filerese);
2397:
2398: strcpy(fileresv,"v");
2399: strcat(fileresv,fileres);
2400: if((ficresvij=fopen(fileresv,"w"))==NULL) {
2401: printf("Problem with variance resultfile: %s\n", fileresv);exit(0);
2402: }
2403: printf("Computing Variance-covariance of DFLEs: file '%s' \n", fileresv);
2404:
2405: k=0;
2406: for(cptcov=1;cptcov<=i1;cptcov++){
2407: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
2408: k=k+1;
2409: fprintf(ficrest,"\n#****** ");
2410: for(j=1;j<=cptcovn;j++)
2411: fprintf(ficrest,"V%d=%d ",Tvar[j],nbcode[Tvar[j]][codtab[k][j]]);
2412: fprintf(ficrest,"******\n");
2413:
2414: fprintf(ficreseij,"\n#****** ");
2415: for(j=1;j<=cptcovn;j++)
2416: fprintf(ficreseij,"V%d=%d ",j,nbcode[j][codtab[k][j]]);
2417: fprintf(ficreseij,"******\n");
2418:
2419: fprintf(ficresvij,"\n#****** ");
2420: for(j=1;j<=cptcovn;j++)
2421: fprintf(ficresvij,"V%d=%d ",j,nbcode[j][codtab[k][j]]);
2422: fprintf(ficresvij,"******\n");
2423:
2424: eij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
2425: oldm=oldms;savm=savms;
2426: evsij(fileres, eij, p, nlstate, stepm, (int) bage, (int)fage, oldm, savm, k);
2427: vareij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
2428: oldm=oldms;savm=savms;
2429: varevsij(fileres, vareij, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k);
2430:
2431: fprintf(ficrest,"#Total LEs with variances: e.. (std) ");
2432: for (i=1;i<=nlstate;i++) fprintf(ficrest,"e.%d (std) ",i);
2433: fprintf(ficrest,"\n");
2434:
2435: hf=1;
2436: if (stepm >= YEARM) hf=stepm/YEARM;
2437: epj=vector(1,nlstate+1);
2438: for(age=bage; age <=fage ;age++){
2439: prevalim(prlim, nlstate, p, age, oldm, savm,ftolpl,k);
2440: fprintf(ficrest," %.0f",age);
2441: for(j=1, epj[nlstate+1]=0.;j <=nlstate;j++){
2442: for(i=1, epj[j]=0.;i <=nlstate;i++) {
2443: epj[j] += prlim[i][i]*hf*eij[i][j][(int)age];
2444: }
2445: epj[nlstate+1] +=epj[j];
2446: }
2447: for(i=1, vepp=0.;i <=nlstate;i++)
2448: for(j=1;j <=nlstate;j++)
2449: vepp += vareij[i][j][(int)age];
2450: fprintf(ficrest," %.2f (%.2f)", epj[nlstate+1],hf*sqrt(vepp));
2451: for(j=1;j <=nlstate;j++){
2452: fprintf(ficrest," %.2f (%.2f)", epj[j],hf*sqrt(vareij[j][j][(int)age]));
2453: }
2454: fprintf(ficrest,"\n");
2455: }
2456: }
2457: }
2458:
2459: fclose(ficreseij);
2460: fclose(ficresvij);
2461: fclose(ficrest);
2462: fclose(ficpar);
2463: free_vector(epj,1,nlstate+1);
2464: /* scanf("%d ",i); */
2465:
2466: /*------- Variance limit prevalence------*/
2467:
2468: strcpy(fileresvpl,"vpl");
2469: strcat(fileresvpl,fileres);
2470: if((ficresvpl=fopen(fileresvpl,"w"))==NULL) {
2471: printf("Problem with variance prev lim resultfile: %s\n", fileresvpl);
2472: exit(0);
2473: }
2474: printf("Computing Variance-covariance of Prevalence limit: file '%s' \n", fileresvpl);
2475:
2476: k=0;
2477: for(cptcov=1;cptcov<=i1;cptcov++){
2478: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
2479: k=k+1;
2480: fprintf(ficresvpl,"\n#****** ");
2481: for(j=1;j<=cptcovn;j++)
2482: fprintf(ficresvpl,"V%d=%d ",Tvar[j],nbcode[Tvar[j]][codtab[k][j]]);
2483: fprintf(ficresvpl,"******\n");
2484:
2485: varpl=matrix(1,nlstate,(int) bage, (int) fage);
2486: oldm=oldms;savm=savms;
2487: varprevlim(fileres, varpl, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k);
2488: }
2489: }
2490:
2491: fclose(ficresvpl);
2492:
2493: /*---------- End : free ----------------*/
2494: free_matrix(varpl,1,nlstate,(int) bage, (int)fage);
2495:
2496: free_ma3x(vareij,1,nlstate,1,nlstate,(int) bage, (int)fage);
2497: free_ma3x(eij,1,nlstate,1,nlstate,(int) bage, (int)fage);
2498:
2499:
2500: free_matrix(pmmij,1,nlstate+ndeath,1,nlstate+ndeath);
2501: free_matrix(oldms, 1,nlstate+ndeath,1,nlstate+ndeath);
2502: free_matrix(newms, 1,nlstate+ndeath,1,nlstate+ndeath);
2503: free_matrix(savms, 1,nlstate+ndeath,1,nlstate+ndeath);
2504:
2505: free_matrix(matcov,1,npar,1,npar);
2506: free_vector(delti,1,npar);
2507:
2508: free_ma3x(param,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel);
2509:
2510: printf("End of Imach\n");
2511: /* gettimeofday(&end_time, (struct timezone*)0);*/ /* after time */
2512:
2513: /* printf("Total time was %d Sec. %d uSec.\n", end_time.tv_sec -start_time.tv_sec, end_time.tv_usec -start_time.tv_usec);*/
2514: /*printf("Total time was %d uSec.\n", total_usecs);*/
2515: /*------ End -----------*/
2516:
2517: end:
2518: #ifdef windows
2519: chdir(pathcd);
2520: #endif
2521: system("wgnuplot ../gp37mgw/graph.plt");
2522:
2523: #ifdef windows
2524: while (z[0] != 'q') {
2525: chdir(pathcd);
2526: printf("\nType e to edit output files, c to start again, and q for exiting: ");
2527: scanf("%s",z);
2528: if (z[0] == 'c') system("./imach");
2529: else if (z[0] == 'e') {
2530: chdir(path);
2531: system("index.htm");
2532: }
2533: else if (z[0] == 'q') exit(0);
2534: }
2535: #endif
2536: }
2537:
2538:
FreeBSD-CVSweb <freebsd-cvsweb@FreeBSD.org>