Annotation of imach/src/imach.c, revision 1.40
1.40 ! lievre 1: /* $Id: imach.c,v 1.39 2002/04/05 15:45:00 lievre Exp $
1.24 lievre 2: Interpolated Markov Chain
1.22 brouard 3:
4: Short summary of the programme:
5:
6: This program computes Healthy Life Expectancies from
7: cross-longitudinal data. Cross-longitudinal data consist in: -1- a
8: first survey ("cross") where individuals from different ages are
9: interviewed on their health status or degree of disability (in the
10: case of a health survey which is our main interest) -2- at least a
11: second wave of interviews ("longitudinal") which measure each change
12: (if any) in individual health status. Health expectancies are
13: computed from the time spent in each health state according to a
14: model. More health states you consider, more time is necessary to reach the
15: Maximum Likelihood of the parameters involved in the model. The
16: simplest model is the multinomial logistic model where pij is the
1.39 lievre 17: probability to be observed in state j at the second wave
1.22 brouard 18: conditional to be observed in state i at the first wave. Therefore
19: the model is: log(pij/pii)= aij + bij*age+ cij*sex + etc , where
20: 'age' is age and 'sex' is a covariate. If you want to have a more
21: complex model than "constant and age", you should modify the program
22: where the markup *Covariates have to be included here again* invites
23: you to do it. More covariates you add, slower the
24: convergence.
25:
26: The advantage of this computer programme, compared to a simple
27: multinomial logistic model, is clear when the delay between waves is not
28: identical for each individual. Also, if a individual missed an
29: intermediate interview, the information is lost, but taken into
30: account using an interpolation or extrapolation.
31:
32: hPijx is the probability to be observed in state i at age x+h
33: conditional to the observed state i at age x. The delay 'h' can be
34: split into an exact number (nh*stepm) of unobserved intermediate
35: states. This elementary transition (by month or quarter trimester,
36: semester or year) is model as a multinomial logistic. The hPx
37: matrix is simply the matrix product of nh*stepm elementary matrices
38: and the contribution of each individual to the likelihood is simply
39: hPijx.
1.2 lievre 40:
41: Also this programme outputs the covariance matrix of the parameters but also
42: of the life expectancies. It also computes the prevalence limits.
43:
44: Authors: Nicolas Brouard (brouard@ined.fr) and Agnès Lièvre (lievre@ined.fr).
45: Institut national d'études démographiques, Paris.
46: This software have been partly granted by Euro-REVES, a concerted action
47: from the European Union.
48: It is copyrighted identically to a GNU software product, ie programme and
49: software can be distributed freely for non commercial use. Latest version
50: can be accessed at http://euroreves.ined.fr/imach .
51: **********************************************************************/
52:
53: #include <math.h>
54: #include <stdio.h>
55: #include <stdlib.h>
56: #include <unistd.h>
57:
58: #define MAXLINE 256
1.35 lievre 59: #define GNUPLOTPROGRAM "wgnuplot"
60: /*#define GNUPLOTPROGRAM "..\\gp37mgw\\wgnuplot"*/
1.2 lievre 61: #define FILENAMELENGTH 80
62: /*#define DEBUG*/
63: #define windows
1.5 lievre 64: #define GLOCK_ERROR_NOPATH -1 /* empty path */
65: #define GLOCK_ERROR_GETCWD -2 /* cannot get cwd */
66:
1.2 lievre 67: #define MAXPARM 30 /* Maximum number of parameters for the optimization */
68: #define NPARMAX 64 /* (nlstate+ndeath-1)*nlstate*ncovmodel */
69:
70: #define NINTERVMAX 8
71: #define NLSTATEMAX 8 /* Maximum number of live states (for func) */
72: #define NDEATHMAX 8 /* Maximum number of dead states (for func) */
73: #define NCOVMAX 8 /* Maximum number of covariates */
1.3 lievre 74: #define MAXN 20000
1.2 lievre 75: #define YEARM 12. /* Number of months per year */
76: #define AGESUP 130
77: #define AGEBASE 40
78:
79:
1.21 lievre 80: int erreur; /* Error number */
1.2 lievre 81: int nvar;
1.8 lievre 82: int cptcovn, cptcovage=0, cptcoveff=0,cptcov;
1.2 lievre 83: int npar=NPARMAX;
84: int nlstate=2; /* Number of live states */
85: int ndeath=1; /* Number of dead states */
1.34 brouard 86: int ncovmodel, ncovcol; /* Total number of covariables including constant a12*1 +b12*x ncovmodel=2 */
1.15 lievre 87: int popbased=0;
1.2 lievre 88:
89: int *wav; /* Number of waves for this individuual 0 is possible */
90: int maxwav; /* Maxim number of waves */
1.8 lievre 91: int jmin, jmax; /* min, max spacing between 2 waves */
1.2 lievre 92: int mle, weightopt;
93: int **mw; /* mw[mi][i] is number of the mi wave for this individual */
94: int **dh; /* dh[mi][i] is number of steps between mi,mi+1 for this individual */
1.8 lievre 95: double jmean; /* Mean space between 2 waves */
1.2 lievre 96: double **oldm, **newm, **savm; /* Working pointers to matrices */
97: double **oldms, **newms, **savms; /* Fixed working pointers to matrices */
1.27 lievre 98: FILE *fic,*ficpar, *ficparo,*ficres, *ficrespl, *ficrespij, *ficrest,*ficresf,*ficrespop;
1.25 lievre 99: FILE *ficgp,*ficresprob,*ficpop;
1.2 lievre 100: FILE *ficreseij;
101: char filerese[FILENAMELENGTH];
102: FILE *ficresvij;
103: char fileresv[FILENAMELENGTH];
104: FILE *ficresvpl;
105: char fileresvpl[FILENAMELENGTH];
106:
107: #define NR_END 1
108: #define FREE_ARG char*
109: #define FTOL 1.0e-10
110:
111: #define NRANSI
112: #define ITMAX 200
113:
114: #define TOL 2.0e-4
115:
116: #define CGOLD 0.3819660
117: #define ZEPS 1.0e-10
118: #define SHFT(a,b,c,d) (a)=(b);(b)=(c);(c)=(d);
119:
120: #define GOLD 1.618034
121: #define GLIMIT 100.0
122: #define TINY 1.0e-20
123:
124: static double maxarg1,maxarg2;
125: #define FMAX(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)>(maxarg2)? (maxarg1):(maxarg2))
126: #define FMIN(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)<(maxarg2)? (maxarg1):(maxarg2))
1.25 lievre 127:
1.2 lievre 128: #define SIGN(a,b) ((b)>0.0 ? fabs(a) : -fabs(a))
129: #define rint(a) floor(a+0.5)
130:
131: static double sqrarg;
132: #define SQR(a) ((sqrarg=(a)) == 0.0 ? 0.0 :sqrarg*sqrarg)
133: #define SWAP(a,b) {temp=(a);(a)=(b);(b)=temp;}
134:
135: int imx;
136: int stepm;
137: /* Stepm, step in month: minimum step interpolation*/
138:
1.36 brouard 139: int estepm;
140: /* Estepm, step in month to interpolate survival function in order to approximate Life Expectancy*/
141:
1.2 lievre 142: int m,nb;
1.6 lievre 143: int *num, firstpass=0, lastpass=4,*cod, *ncodemax, *Tage;
1.2 lievre 144: double **agev,*moisnais, *annais, *moisdc, *andc,**mint, **anint;
1.13 lievre 145: double **pmmij, ***probs, ***mobaverage;
1.19 lievre 146: double dateintmean=0;
1.2 lievre 147:
148: double *weight;
149: int **s; /* Status */
150: double *agedc, **covar, idx;
1.7 lievre 151: int **nbcode, *Tcode, *Tvar, **codtab, **Tvard, *Tprod, cptcovprod, *Tvaraff;
1.2 lievre 152:
153: double ftol=FTOL; /* Tolerance for computing Max Likelihood */
154: double ftolhess; /* Tolerance for computing hessian */
155:
1.7 lievre 156: /**************** split *************************/
1.22 brouard 157: static int split( char *path, char *dirc, char *name, char *ext, char *finame )
1.5 lievre 158: {
159: char *s; /* pointer */
160: int l1, l2; /* length counters */
161:
162: l1 = strlen( path ); /* length of path */
163: if ( l1 == 0 ) return( GLOCK_ERROR_NOPATH );
1.22 brouard 164: #ifdef windows
1.5 lievre 165: s = strrchr( path, '\\' ); /* find last / */
1.22 brouard 166: #else
167: s = strrchr( path, '/' ); /* find last / */
168: #endif
1.5 lievre 169: if ( s == NULL ) { /* no directory, so use current */
170: #if defined(__bsd__) /* get current working directory */
171: extern char *getwd( );
172:
173: if ( getwd( dirc ) == NULL ) {
174: #else
175: extern char *getcwd( );
176:
177: if ( getcwd( dirc, FILENAME_MAX ) == NULL ) {
178: #endif
179: return( GLOCK_ERROR_GETCWD );
180: }
181: strcpy( name, path ); /* we've got it */
182: } else { /* strip direcotry from path */
183: s++; /* after this, the filename */
184: l2 = strlen( s ); /* length of filename */
185: if ( l2 == 0 ) return( GLOCK_ERROR_NOPATH );
186: strcpy( name, s ); /* save file name */
187: strncpy( dirc, path, l1 - l2 ); /* now the directory */
188: dirc[l1-l2] = 0; /* add zero */
189: }
190: l1 = strlen( dirc ); /* length of directory */
1.22 brouard 191: #ifdef windows
1.5 lievre 192: if ( dirc[l1-1] != '\\' ) { dirc[l1] = '\\'; dirc[l1+1] = 0; }
1.22 brouard 193: #else
194: if ( dirc[l1-1] != '/' ) { dirc[l1] = '/'; dirc[l1+1] = 0; }
195: #endif
196: s = strrchr( name, '.' ); /* find last / */
197: s++;
198: strcpy(ext,s); /* save extension */
199: l1= strlen( name);
200: l2= strlen( s)+1;
201: strncpy( finame, name, l1-l2);
202: finame[l1-l2]= 0;
1.5 lievre 203: return( 0 ); /* we're done */
204: }
205:
206:
1.2 lievre 207: /******************************************/
208:
209: void replace(char *s, char*t)
210: {
211: int i;
212: int lg=20;
213: i=0;
214: lg=strlen(t);
215: for(i=0; i<= lg; i++) {
216: (s[i] = t[i]);
217: if (t[i]== '\\') s[i]='/';
218: }
219: }
220:
221: int nbocc(char *s, char occ)
222: {
223: int i,j=0;
224: int lg=20;
225: i=0;
226: lg=strlen(s);
227: for(i=0; i<= lg; i++) {
228: if (s[i] == occ ) j++;
229: }
230: return j;
231: }
232:
233: void cutv(char *u,char *v, char*t, char occ)
234: {
1.6 lievre 235: int i,lg,j,p=0;
1.2 lievre 236: i=0;
237: for(j=0; j<=strlen(t)-1; j++) {
1.3 lievre 238: if((t[j]!= occ) && (t[j+1]== occ)) p=j+1;
1.2 lievre 239: }
240:
241: lg=strlen(t);
242: for(j=0; j<p; j++) {
243: (u[j] = t[j]);
244: }
1.6 lievre 245: u[p]='\0';
1.2 lievre 246:
247: for(j=0; j<= lg; j++) {
248: if (j>=(p+1))(v[j-p-1] = t[j]);
249: }
250: }
251:
252: /********************** nrerror ********************/
253:
254: void nrerror(char error_text[])
255: {
256: fprintf(stderr,"ERREUR ...\n");
257: fprintf(stderr,"%s\n",error_text);
258: exit(1);
259: }
260: /*********************** vector *******************/
261: double *vector(int nl, int nh)
262: {
263: double *v;
264: v=(double *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(double)));
265: if (!v) nrerror("allocation failure in vector");
266: return v-nl+NR_END;
267: }
268:
269: /************************ free vector ******************/
270: void free_vector(double*v, int nl, int nh)
271: {
272: free((FREE_ARG)(v+nl-NR_END));
273: }
274:
275: /************************ivector *******************************/
276: int *ivector(long nl,long nh)
277: {
278: int *v;
279: v=(int *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(int)));
280: if (!v) nrerror("allocation failure in ivector");
281: return v-nl+NR_END;
282: }
283:
284: /******************free ivector **************************/
285: void free_ivector(int *v, long nl, long nh)
286: {
287: free((FREE_ARG)(v+nl-NR_END));
288: }
289:
290: /******************* imatrix *******************************/
291: int **imatrix(long nrl, long nrh, long ncl, long nch)
292: /* allocate a int matrix with subscript range m[nrl..nrh][ncl..nch] */
293: {
294: long i, nrow=nrh-nrl+1,ncol=nch-ncl+1;
295: int **m;
296:
297: /* allocate pointers to rows */
298: m=(int **) malloc((size_t)((nrow+NR_END)*sizeof(int*)));
299: if (!m) nrerror("allocation failure 1 in matrix()");
300: m += NR_END;
301: m -= nrl;
302:
303:
304: /* allocate rows and set pointers to them */
305: m[nrl]=(int *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(int)));
306: if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
307: m[nrl] += NR_END;
308: m[nrl] -= ncl;
309:
310: for(i=nrl+1;i<=nrh;i++) m[i]=m[i-1]+ncol;
311:
312: /* return pointer to array of pointers to rows */
313: return m;
314: }
315:
316: /****************** free_imatrix *************************/
317: void free_imatrix(m,nrl,nrh,ncl,nch)
318: int **m;
319: long nch,ncl,nrh,nrl;
320: /* free an int matrix allocated by imatrix() */
321: {
322: free((FREE_ARG) (m[nrl]+ncl-NR_END));
323: free((FREE_ARG) (m+nrl-NR_END));
324: }
325:
326: /******************* matrix *******************************/
327: double **matrix(long nrl, long nrh, long ncl, long nch)
328: {
329: long i, nrow=nrh-nrl+1, ncol=nch-ncl+1;
330: double **m;
331:
332: m=(double **) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
333: if (!m) nrerror("allocation failure 1 in matrix()");
334: m += NR_END;
335: m -= nrl;
336:
337: m[nrl]=(double *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
338: if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
339: m[nrl] += NR_END;
340: m[nrl] -= ncl;
341:
342: for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
343: return m;
344: }
345:
346: /*************************free matrix ************************/
347: void free_matrix(double **m, long nrl, long nrh, long ncl, long nch)
348: {
349: free((FREE_ARG)(m[nrl]+ncl-NR_END));
350: free((FREE_ARG)(m+nrl-NR_END));
351: }
352:
353: /******************* ma3x *******************************/
354: double ***ma3x(long nrl, long nrh, long ncl, long nch, long nll, long nlh)
355: {
356: long i, j, nrow=nrh-nrl+1, ncol=nch-ncl+1, nlay=nlh-nll+1;
357: double ***m;
358:
359: m=(double ***) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
360: if (!m) nrerror("allocation failure 1 in matrix()");
361: m += NR_END;
362: m -= nrl;
363:
364: m[nrl]=(double **) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
365: if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
366: m[nrl] += NR_END;
367: m[nrl] -= ncl;
368:
369: for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
370:
371: m[nrl][ncl]=(double *) malloc((size_t)((nrow*ncol*nlay+NR_END)*sizeof(double)));
372: if (!m[nrl][ncl]) nrerror("allocation failure 3 in matrix()");
373: m[nrl][ncl] += NR_END;
374: m[nrl][ncl] -= nll;
375: for (j=ncl+1; j<=nch; j++)
376: m[nrl][j]=m[nrl][j-1]+nlay;
377:
378: for (i=nrl+1; i<=nrh; i++) {
379: m[i][ncl]=m[i-1l][ncl]+ncol*nlay;
380: for (j=ncl+1; j<=nch; j++)
381: m[i][j]=m[i][j-1]+nlay;
382: }
383: return m;
384: }
385:
386: /*************************free ma3x ************************/
387: void free_ma3x(double ***m, long nrl, long nrh, long ncl, long nch,long nll, long nlh)
388: {
389: free((FREE_ARG)(m[nrl][ncl]+ nll-NR_END));
390: free((FREE_ARG)(m[nrl]+ncl-NR_END));
391: free((FREE_ARG)(m+nrl-NR_END));
392: }
393:
394: /***************** f1dim *************************/
395: extern int ncom;
396: extern double *pcom,*xicom;
397: extern double (*nrfunc)(double []);
398:
399: double f1dim(double x)
400: {
401: int j;
402: double f;
403: double *xt;
404:
405: xt=vector(1,ncom);
406: for (j=1;j<=ncom;j++) xt[j]=pcom[j]+x*xicom[j];
407: f=(*nrfunc)(xt);
408: free_vector(xt,1,ncom);
409: return f;
410: }
411:
412: /*****************brent *************************/
413: double brent(double ax, double bx, double cx, double (*f)(double), double tol, double *xmin)
414: {
415: int iter;
416: double a,b,d,etemp;
417: double fu,fv,fw,fx;
418: double ftemp;
419: double p,q,r,tol1,tol2,u,v,w,x,xm;
420: double e=0.0;
421:
422: a=(ax < cx ? ax : cx);
423: b=(ax > cx ? ax : cx);
424: x=w=v=bx;
425: fw=fv=fx=(*f)(x);
426: for (iter=1;iter<=ITMAX;iter++) {
427: xm=0.5*(a+b);
428: tol2=2.0*(tol1=tol*fabs(x)+ZEPS);
429: /* if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret)))*/
430: printf(".");fflush(stdout);
431: #ifdef DEBUG
432: 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);
433: /* if ((fabs(x-xm) <= (tol2-0.5*(b-a)))||(2.0*fabs(fu-ftemp) <= ftol*1.e-2*(fabs(fu)+fabs(ftemp)))) { */
434: #endif
435: if (fabs(x-xm) <= (tol2-0.5*(b-a))){
436: *xmin=x;
437: return fx;
438: }
439: ftemp=fu;
440: if (fabs(e) > tol1) {
441: r=(x-w)*(fx-fv);
442: q=(x-v)*(fx-fw);
443: p=(x-v)*q-(x-w)*r;
444: q=2.0*(q-r);
445: if (q > 0.0) p = -p;
446: q=fabs(q);
447: etemp=e;
448: e=d;
449: if (fabs(p) >= fabs(0.5*q*etemp) || p <= q*(a-x) || p >= q*(b-x))
450: d=CGOLD*(e=(x >= xm ? a-x : b-x));
451: else {
452: d=p/q;
453: u=x+d;
454: if (u-a < tol2 || b-u < tol2)
455: d=SIGN(tol1,xm-x);
456: }
457: } else {
458: d=CGOLD*(e=(x >= xm ? a-x : b-x));
459: }
460: u=(fabs(d) >= tol1 ? x+d : x+SIGN(tol1,d));
461: fu=(*f)(u);
462: if (fu <= fx) {
463: if (u >= x) a=x; else b=x;
464: SHFT(v,w,x,u)
465: SHFT(fv,fw,fx,fu)
466: } else {
467: if (u < x) a=u; else b=u;
468: if (fu <= fw || w == x) {
469: v=w;
470: w=u;
471: fv=fw;
472: fw=fu;
473: } else if (fu <= fv || v == x || v == w) {
474: v=u;
475: fv=fu;
476: }
477: }
478: }
479: nrerror("Too many iterations in brent");
480: *xmin=x;
481: return fx;
482: }
483:
484: /****************** mnbrak ***********************/
485:
486: void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb, double *fc,
487: double (*func)(double))
488: {
489: double ulim,u,r,q, dum;
490: double fu;
491:
492: *fa=(*func)(*ax);
493: *fb=(*func)(*bx);
494: if (*fb > *fa) {
495: SHFT(dum,*ax,*bx,dum)
496: SHFT(dum,*fb,*fa,dum)
497: }
498: *cx=(*bx)+GOLD*(*bx-*ax);
499: *fc=(*func)(*cx);
500: while (*fb > *fc) {
501: r=(*bx-*ax)*(*fb-*fc);
502: q=(*bx-*cx)*(*fb-*fa);
503: u=(*bx)-((*bx-*cx)*q-(*bx-*ax)*r)/
504: (2.0*SIGN(FMAX(fabs(q-r),TINY),q-r));
505: ulim=(*bx)+GLIMIT*(*cx-*bx);
506: if ((*bx-u)*(u-*cx) > 0.0) {
507: fu=(*func)(u);
508: } else if ((*cx-u)*(u-ulim) > 0.0) {
509: fu=(*func)(u);
510: if (fu < *fc) {
511: SHFT(*bx,*cx,u,*cx+GOLD*(*cx-*bx))
512: SHFT(*fb,*fc,fu,(*func)(u))
513: }
514: } else if ((u-ulim)*(ulim-*cx) >= 0.0) {
515: u=ulim;
516: fu=(*func)(u);
517: } else {
518: u=(*cx)+GOLD*(*cx-*bx);
519: fu=(*func)(u);
520: }
521: SHFT(*ax,*bx,*cx,u)
522: SHFT(*fa,*fb,*fc,fu)
523: }
524: }
525:
526: /*************** linmin ************************/
527:
528: int ncom;
529: double *pcom,*xicom;
530: double (*nrfunc)(double []);
531:
532: void linmin(double p[], double xi[], int n, double *fret,double (*func)(double []))
533: {
534: double brent(double ax, double bx, double cx,
535: double (*f)(double), double tol, double *xmin);
536: double f1dim(double x);
537: void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb,
538: double *fc, double (*func)(double));
539: int j;
540: double xx,xmin,bx,ax;
541: double fx,fb,fa;
542:
543: ncom=n;
544: pcom=vector(1,n);
545: xicom=vector(1,n);
546: nrfunc=func;
547: for (j=1;j<=n;j++) {
548: pcom[j]=p[j];
549: xicom[j]=xi[j];
550: }
551: ax=0.0;
552: xx=1.0;
553: mnbrak(&ax,&xx,&bx,&fa,&fx,&fb,f1dim);
554: *fret=brent(ax,xx,bx,f1dim,TOL,&xmin);
555: #ifdef DEBUG
556: printf("retour brent fret=%.12e xmin=%.12e\n",*fret,xmin);
557: #endif
558: for (j=1;j<=n;j++) {
559: xi[j] *= xmin;
560: p[j] += xi[j];
561: }
562: free_vector(xicom,1,n);
563: free_vector(pcom,1,n);
564: }
565:
566: /*************** powell ************************/
567: void powell(double p[], double **xi, int n, double ftol, int *iter, double *fret,
568: double (*func)(double []))
569: {
570: void linmin(double p[], double xi[], int n, double *fret,
571: double (*func)(double []));
572: int i,ibig,j;
573: double del,t,*pt,*ptt,*xit;
574: double fp,fptt;
575: double *xits;
576: pt=vector(1,n);
577: ptt=vector(1,n);
578: xit=vector(1,n);
579: xits=vector(1,n);
580: *fret=(*func)(p);
581: for (j=1;j<=n;j++) pt[j]=p[j];
582: for (*iter=1;;++(*iter)) {
583: fp=(*fret);
584: ibig=0;
585: del=0.0;
586: printf("\nPowell iter=%d -2*LL=%.12f",*iter,*fret);
587: for (i=1;i<=n;i++)
588: printf(" %d %.12f",i, p[i]);
589: printf("\n");
590: for (i=1;i<=n;i++) {
591: for (j=1;j<=n;j++) xit[j]=xi[j][i];
592: fptt=(*fret);
593: #ifdef DEBUG
594: printf("fret=%lf \n",*fret);
595: #endif
596: printf("%d",i);fflush(stdout);
597: linmin(p,xit,n,fret,func);
598: if (fabs(fptt-(*fret)) > del) {
599: del=fabs(fptt-(*fret));
600: ibig=i;
601: }
602: #ifdef DEBUG
603: printf("%d %.12e",i,(*fret));
604: for (j=1;j<=n;j++) {
605: xits[j]=FMAX(fabs(p[j]-pt[j]),1.e-5);
606: printf(" x(%d)=%.12e",j,xit[j]);
607: }
608: for(j=1;j<=n;j++)
609: printf(" p=%.12e",p[j]);
610: printf("\n");
611: #endif
612: }
613: if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret))) {
614: #ifdef DEBUG
615: int k[2],l;
616: k[0]=1;
617: k[1]=-1;
618: printf("Max: %.12e",(*func)(p));
619: for (j=1;j<=n;j++)
620: printf(" %.12e",p[j]);
621: printf("\n");
622: for(l=0;l<=1;l++) {
623: for (j=1;j<=n;j++) {
624: ptt[j]=p[j]+(p[j]-pt[j])*k[l];
625: printf("l=%d j=%d ptt=%.12e, xits=%.12e, p=%.12e, xit=%.12e", l,j,ptt[j],xits[j],p[j],xit[j]);
626: }
627: printf("func(ptt)=%.12e, deriv=%.12e\n",(*func)(ptt),(ptt[j]-p[j])/((*func)(ptt)-(*func)(p)));
628: }
629: #endif
630:
631:
632: free_vector(xit,1,n);
633: free_vector(xits,1,n);
634: free_vector(ptt,1,n);
635: free_vector(pt,1,n);
636: return;
637: }
638: if (*iter == ITMAX) nrerror("powell exceeding maximum iterations.");
639: for (j=1;j<=n;j++) {
640: ptt[j]=2.0*p[j]-pt[j];
641: xit[j]=p[j]-pt[j];
642: pt[j]=p[j];
643: }
644: fptt=(*func)(ptt);
645: if (fptt < fp) {
646: t=2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del)-del*SQR(fp-fptt);
647: if (t < 0.0) {
648: linmin(p,xit,n,fret,func);
649: for (j=1;j<=n;j++) {
650: xi[j][ibig]=xi[j][n];
651: xi[j][n]=xit[j];
652: }
653: #ifdef DEBUG
654: printf("Direction changed last moved %d in place of ibig=%d, new last is the average:\n",n,ibig);
655: for(j=1;j<=n;j++)
656: printf(" %.12e",xit[j]);
657: printf("\n");
658: #endif
659: }
660: }
661: }
662: }
663:
664: /**** Prevalence limit ****************/
665:
666: double **prevalim(double **prlim, int nlstate, double x[], double age, double **oldm, double **savm, double ftolpl, int ij)
667: {
668: /* Computes the prevalence limit in each live state at age x by left multiplying the unit
669: matrix by transitions matrix until convergence is reached */
670:
671: int i, ii,j,k;
672: double min, max, maxmin, maxmax,sumnew=0.;
673: double **matprod2();
674: double **out, cov[NCOVMAX], **pmij();
675: double **newm;
676: double agefin, delaymax=50 ; /* Max number of years to converge */
677:
678: for (ii=1;ii<=nlstate+ndeath;ii++)
679: for (j=1;j<=nlstate+ndeath;j++){
680: oldm[ii][j]=(ii==j ? 1.0 : 0.0);
681: }
1.6 lievre 682:
683: cov[1]=1.;
684:
685: /* Even if hstepm = 1, at least one multiplication by the unit matrix */
1.2 lievre 686: for(agefin=age-stepm/YEARM; agefin>=age-delaymax; agefin=agefin-stepm/YEARM){
687: newm=savm;
688: /* Covariates have to be included here again */
1.6 lievre 689: cov[2]=agefin;
690:
691: for (k=1; k<=cptcovn;k++) {
1.7 lievre 692: cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]];
1.35 lievre 693: /* printf("ij=%d k=%d Tvar[k]=%d nbcode=%d cov=%lf codtab[ij][Tvar[k]]=%d \n",ij,k, Tvar[k],nbcode[Tvar[k]][codtab[ij][Tvar[k]]],cov[2+k], codtab[ij][Tvar[k]]);*/
1.6 lievre 694: }
1.35 lievre 695: for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
1.7 lievre 696: for (k=1; k<=cptcovprod;k++)
697: cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
698:
699: /*printf("ij=%d cptcovprod=%d tvar=%d ", ij, cptcovprod, Tvar[1]);*/
700: /*printf("ij=%d cov[3]=%lf cov[4]=%lf \n",ij, cov[3],cov[4]);*/
1.35 lievre 701: /*printf("ij=%d cov[3]=%lf \n",ij, cov[3]);*/
1.2 lievre 702: out=matprod2(newm, pmij(pmmij,cov,ncovmodel,x,nlstate),1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, oldm);
703:
704: savm=oldm;
705: oldm=newm;
706: maxmax=0.;
707: for(j=1;j<=nlstate;j++){
708: min=1.;
709: max=0.;
710: for(i=1; i<=nlstate; i++) {
711: sumnew=0;
712: for(k=1; k<=ndeath; k++) sumnew+=newm[i][nlstate+k];
713: prlim[i][j]= newm[i][j]/(1-sumnew);
714: max=FMAX(max,prlim[i][j]);
715: min=FMIN(min,prlim[i][j]);
716: }
717: maxmin=max-min;
718: maxmax=FMAX(maxmax,maxmin);
719: }
720: if(maxmax < ftolpl){
721: return prlim;
722: }
723: }
724: }
725:
1.12 lievre 726: /*************** transition probabilities ***************/
1.2 lievre 727:
728: double **pmij(double **ps, double *cov, int ncovmodel, double *x, int nlstate )
729: {
730: double s1, s2;
731: /*double t34;*/
732: int i,j,j1, nc, ii, jj;
733:
734: for(i=1; i<= nlstate; i++){
735: for(j=1; j<i;j++){
736: for (nc=1, s2=0.;nc <=ncovmodel; nc++){
737: /*s2 += param[i][j][nc]*cov[nc];*/
738: s2 += x[(i-1)*nlstate*ncovmodel+(j-1)*ncovmodel+nc+(i-1)*(ndeath-1)*ncovmodel]*cov[nc];
739: /*printf("Int j<i s1=%.17e, s2=%.17e\n",s1,s2);*/
740: }
741: ps[i][j]=s2;
742: /*printf("s1=%.17e, s2=%.17e\n",s1,s2);*/
743: }
744: for(j=i+1; j<=nlstate+ndeath;j++){
745: for (nc=1, s2=0.;nc <=ncovmodel; nc++){
746: s2 += x[(i-1)*nlstate*ncovmodel+(j-2)*ncovmodel+nc+(i-1)*(ndeath-1)*ncovmodel]*cov[nc];
747: /*printf("Int j>i s1=%.17e, s2=%.17e %lx %lx\n",s1,s2,s1,s2);*/
748: }
1.22 brouard 749: ps[i][j]=s2;
1.2 lievre 750: }
751: }
1.12 lievre 752: /*ps[3][2]=1;*/
753:
1.2 lievre 754: for(i=1; i<= nlstate; i++){
755: s1=0;
756: for(j=1; j<i; j++)
757: s1+=exp(ps[i][j]);
758: for(j=i+1; j<=nlstate+ndeath; j++)
759: s1+=exp(ps[i][j]);
760: ps[i][i]=1./(s1+1.);
761: for(j=1; j<i; j++)
762: ps[i][j]= exp(ps[i][j])*ps[i][i];
763: for(j=i+1; j<=nlstate+ndeath; j++)
764: ps[i][j]= exp(ps[i][j])*ps[i][i];
765: /* ps[i][nlstate+1]=1.-s1- ps[i][i];*/ /* Sum should be 1 */
766: } /* end i */
767:
768: for(ii=nlstate+1; ii<= nlstate+ndeath; ii++){
769: for(jj=1; jj<= nlstate+ndeath; jj++){
770: ps[ii][jj]=0;
771: ps[ii][ii]=1;
772: }
773: }
774:
1.12 lievre 775:
1.2 lievre 776: /* for(ii=1; ii<= nlstate+ndeath; ii++){
777: for(jj=1; jj<= nlstate+ndeath; jj++){
778: printf("%lf ",ps[ii][jj]);
779: }
780: printf("\n ");
781: }
782: printf("\n ");printf("%lf ",cov[2]);*/
783: /*
784: for(i=1; i<= npar; i++) printf("%f ",x[i]);
785: goto end;*/
786: return ps;
787: }
788:
789: /**************** Product of 2 matrices ******************/
790:
791: double **matprod2(double **out, double **in,long nrl, long nrh, long ncl, long nch, long ncolol, long ncoloh, double **b)
792: {
1.13 lievre 793: /* Computes the matrix product of in(1,nrh-nrl+1)(1,nch-ncl+1) times
1.2 lievre 794: b(1,nch-ncl+1)(1,ncoloh-ncolol+1) into out(...) */
795: /* in, b, out are matrice of pointers which should have been initialized
796: before: only the contents of out is modified. The function returns
797: a pointer to pointers identical to out */
798: long i, j, k;
799: for(i=nrl; i<= nrh; i++)
800: for(k=ncolol; k<=ncoloh; k++)
801: for(j=ncl,out[i][k]=0.; j<=nch; j++)
802: out[i][k] +=in[i][j]*b[j][k];
803:
804: return out;
805: }
806:
807:
808: /************* Higher Matrix Product ***************/
809:
810: double ***hpxij(double ***po, int nhstepm, double age, int hstepm, double *x, int nlstate, int stepm, double **oldm, double **savm, int ij )
811: {
812: /* Computes the transition matrix starting at age 'age' over 'nhstepm*hstepm*stepm' month
813: duration (i.e. until
814: age (in years) age+nhstepm*stepm/12) by multiplying nhstepm*hstepm matrices.
815: Output is stored in matrix po[i][j][h] for h every 'hstepm' step
816: (typically every 2 years instead of every month which is too big).
817: Model is determined by parameters x and covariates have to be
818: included manually here.
819:
820: */
821:
822: int i, j, d, h, k;
823: double **out, cov[NCOVMAX];
824: double **newm;
825:
826: /* Hstepm could be zero and should return the unit matrix */
827: for (i=1;i<=nlstate+ndeath;i++)
828: for (j=1;j<=nlstate+ndeath;j++){
829: oldm[i][j]=(i==j ? 1.0 : 0.0);
830: po[i][j][0]=(i==j ? 1.0 : 0.0);
831: }
832: /* Even if hstepm = 1, at least one multiplication by the unit matrix */
833: for(h=1; h <=nhstepm; h++){
834: for(d=1; d <=hstepm; d++){
835: newm=savm;
836: /* Covariates have to be included here again */
837: cov[1]=1.;
838: cov[2]=age+((h-1)*hstepm + (d-1))*stepm/YEARM;
1.7 lievre 839: for (k=1; k<=cptcovn;k++) cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]];
1.12 lievre 840: for (k=1; k<=cptcovage;k++)
1.7 lievre 841: cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
1.12 lievre 842: for (k=1; k<=cptcovprod;k++)
1.7 lievre 843: cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
844:
845:
1.2 lievre 846: /*printf("hxi cptcov=%d cptcode=%d\n",cptcov,cptcode);*/
847: /*printf("h=%d d=%d age=%f cov=%f\n",h,d,age,cov[2]);*/
848: out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath,
849: pmij(pmmij,cov,ncovmodel,x,nlstate));
850: savm=oldm;
851: oldm=newm;
852: }
853: for(i=1; i<=nlstate+ndeath; i++)
854: for(j=1;j<=nlstate+ndeath;j++) {
855: po[i][j][h]=newm[i][j];
856: /*printf("i=%d j=%d h=%d po[i][j][h]=%f ",i,j,h,po[i][j][h]);
857: */
858: }
859: } /* end h */
860: return po;
861: }
862:
863:
864: /*************** log-likelihood *************/
865: double func( double *x)
866: {
1.6 lievre 867: int i, ii, j, k, mi, d, kk;
1.2 lievre 868: double l, ll[NLSTATEMAX], cov[NCOVMAX];
869: double **out;
870: double sw; /* Sum of weights */
871: double lli; /* Individual log likelihood */
872: long ipmx;
873: /*extern weight */
874: /* We are differentiating ll according to initial status */
875: /* for (i=1;i<=npar;i++) printf("%f ", x[i]);*/
876: /*for(i=1;i<imx;i++)
1.8 lievre 877: printf(" %d\n",s[4][i]);
1.2 lievre 878: */
1.6 lievre 879: cov[1]=1.;
1.2 lievre 880:
881: for(k=1; k<=nlstate; k++) ll[k]=0.;
882: for (i=1,ipmx=0, sw=0.; i<=imx; i++){
1.6 lievre 883: for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
1.8 lievre 884: for(mi=1; mi<= wav[i]-1; mi++){
1.2 lievre 885: for (ii=1;ii<=nlstate+ndeath;ii++)
886: for (j=1;j<=nlstate+ndeath;j++) oldm[ii][j]=(ii==j ? 1.0 : 0.0);
1.8 lievre 887: for(d=0; d<dh[mi][i]; d++){
888: newm=savm;
889: cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
890: for (kk=1; kk<=cptcovage;kk++) {
891: cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
892: }
893:
894: out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
895: 1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
896: savm=oldm;
897: oldm=newm;
898:
899:
1.2 lievre 900: } /* end mult */
1.8 lievre 901:
1.2 lievre 902: lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]);
903: /* printf(" %f ",out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]);*/
904: ipmx +=1;
905: sw += weight[i];
906: ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
907: } /* end of wave */
908: } /* end of individual */
909:
910: for(k=1,l=0.; k<=nlstate; k++) l += ll[k];
911: /* printf("l1=%f l2=%f ",ll[1],ll[2]); */
912: l= l*ipmx/sw; /* To get the same order of magnitude as if weight=1 for every body */
913: return -l;
914: }
915:
916:
917: /*********** Maximum Likelihood Estimation ***************/
918:
919: void mlikeli(FILE *ficres,double p[], int npar, int ncovmodel, int nlstate, double ftol, double (*func)(double []))
920: {
921: int i,j, iter;
922: double **xi,*delti;
923: double fret;
924: xi=matrix(1,npar,1,npar);
925: for (i=1;i<=npar;i++)
926: for (j=1;j<=npar;j++)
927: xi[i][j]=(i==j ? 1.0 : 0.0);
928: printf("Powell\n");
929: powell(p,xi,npar,ftol,&iter,&fret,func);
930:
931: printf("\n#Number of iterations = %d, -2 Log likelihood = %.12f\n",iter,func(p));
1.21 lievre 932: fprintf(ficres,"#Number of iterations = %d, -2 Log likelihood = %.12f \n",iter,func(p));
1.2 lievre 933:
934: }
935:
936: /**** Computes Hessian and covariance matrix ***/
937: void hesscov(double **matcov, double p[], int npar, double delti[], double ftolhess, double (*func)(double []))
938: {
939: double **a,**y,*x,pd;
940: double **hess;
941: int i, j,jk;
942: int *indx;
943:
944: double hessii(double p[], double delta, int theta, double delti[]);
945: double hessij(double p[], double delti[], int i, int j);
946: void lubksb(double **a, int npar, int *indx, double b[]) ;
947: void ludcmp(double **a, int npar, int *indx, double *d) ;
948:
949: hess=matrix(1,npar,1,npar);
950:
951: printf("\nCalculation of the hessian matrix. Wait...\n");
952: for (i=1;i<=npar;i++){
953: printf("%d",i);fflush(stdout);
954: hess[i][i]=hessii(p,ftolhess,i,delti);
955: /*printf(" %f ",p[i]);*/
1.12 lievre 956: /*printf(" %lf ",hess[i][i]);*/
1.2 lievre 957: }
1.12 lievre 958:
1.2 lievre 959: for (i=1;i<=npar;i++) {
960: for (j=1;j<=npar;j++) {
961: if (j>i) {
962: printf(".%d%d",i,j);fflush(stdout);
963: hess[i][j]=hessij(p,delti,i,j);
1.12 lievre 964: hess[j][i]=hess[i][j];
965: /*printf(" %lf ",hess[i][j]);*/
1.2 lievre 966: }
967: }
968: }
969: printf("\n");
970:
971: printf("\nInverting the hessian to get the covariance matrix. Wait...\n");
972:
973: a=matrix(1,npar,1,npar);
974: y=matrix(1,npar,1,npar);
975: x=vector(1,npar);
976: indx=ivector(1,npar);
977: for (i=1;i<=npar;i++)
978: for (j=1;j<=npar;j++) a[i][j]=hess[i][j];
979: ludcmp(a,npar,indx,&pd);
980:
981: for (j=1;j<=npar;j++) {
982: for (i=1;i<=npar;i++) x[i]=0;
983: x[j]=1;
984: lubksb(a,npar,indx,x);
985: for (i=1;i<=npar;i++){
986: matcov[i][j]=x[i];
987: }
988: }
989:
990: printf("\n#Hessian matrix#\n");
991: for (i=1;i<=npar;i++) {
992: for (j=1;j<=npar;j++) {
993: printf("%.3e ",hess[i][j]);
994: }
995: printf("\n");
996: }
997:
998: /* Recompute Inverse */
999: for (i=1;i<=npar;i++)
1000: for (j=1;j<=npar;j++) a[i][j]=matcov[i][j];
1001: ludcmp(a,npar,indx,&pd);
1002:
1003: /* printf("\n#Hessian matrix recomputed#\n");
1004:
1005: for (j=1;j<=npar;j++) {
1006: for (i=1;i<=npar;i++) x[i]=0;
1007: x[j]=1;
1008: lubksb(a,npar,indx,x);
1009: for (i=1;i<=npar;i++){
1010: y[i][j]=x[i];
1011: printf("%.3e ",y[i][j]);
1012: }
1013: printf("\n");
1014: }
1015: */
1016:
1017: free_matrix(a,1,npar,1,npar);
1018: free_matrix(y,1,npar,1,npar);
1019: free_vector(x,1,npar);
1020: free_ivector(indx,1,npar);
1021: free_matrix(hess,1,npar,1,npar);
1022:
1023:
1024: }
1025:
1026: /*************** hessian matrix ****************/
1027: double hessii( double x[], double delta, int theta, double delti[])
1028: {
1029: int i;
1030: int l=1, lmax=20;
1031: double k1,k2;
1032: double p2[NPARMAX+1];
1033: double res;
1034: double delt, delts, nkhi=10.,nkhif=1., khi=1.e-4;
1035: double fx;
1036: int k=0,kmax=10;
1037: double l1;
1038:
1039: fx=func(x);
1040: for (i=1;i<=npar;i++) p2[i]=x[i];
1041: for(l=0 ; l <=lmax; l++){
1042: l1=pow(10,l);
1043: delts=delt;
1044: for(k=1 ; k <kmax; k=k+1){
1045: delt = delta*(l1*k);
1046: p2[theta]=x[theta] +delt;
1047: k1=func(p2)-fx;
1048: p2[theta]=x[theta]-delt;
1049: k2=func(p2)-fx;
1050: /*res= (k1-2.0*fx+k2)/delt/delt; */
1051: res= (k1+k2)/delt/delt/2.; /* Divided by because L and not 2*L */
1052:
1053: #ifdef DEBUG
1054: 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);
1055: #endif
1056: /*if(fabs(k1-2.0*fx+k2) <1.e-13){ */
1057: if((k1 <khi/nkhi/2.) || (k2 <khi/nkhi/2.)){
1058: k=kmax;
1059: }
1060: else if((k1 >khi/nkhif) || (k2 >khi/nkhif)){ /* Keeps lastvalue before 3.84/2 KHI2 5% 1d.f. */
1061: k=kmax; l=lmax*10.;
1062: }
1063: else if((k1 >khi/nkhi) || (k2 >khi/nkhi)){
1064: delts=delt;
1065: }
1066: }
1067: }
1068: delti[theta]=delts;
1.12 lievre 1069: return res;
1.3 lievre 1070:
1.2 lievre 1071: }
1072:
1073: double hessij( double x[], double delti[], int thetai,int thetaj)
1074: {
1075: int i;
1076: int l=1, l1, lmax=20;
1077: double k1,k2,k3,k4,res,fx;
1078: double p2[NPARMAX+1];
1079: int k;
1080:
1081: fx=func(x);
1082: for (k=1; k<=2; k++) {
1083: for (i=1;i<=npar;i++) p2[i]=x[i];
1084: p2[thetai]=x[thetai]+delti[thetai]/k;
1085: p2[thetaj]=x[thetaj]+delti[thetaj]/k;
1086: k1=func(p2)-fx;
1087:
1088: p2[thetai]=x[thetai]+delti[thetai]/k;
1089: p2[thetaj]=x[thetaj]-delti[thetaj]/k;
1090: k2=func(p2)-fx;
1091:
1092: p2[thetai]=x[thetai]-delti[thetai]/k;
1093: p2[thetaj]=x[thetaj]+delti[thetaj]/k;
1094: k3=func(p2)-fx;
1095:
1096: p2[thetai]=x[thetai]-delti[thetai]/k;
1097: p2[thetaj]=x[thetaj]-delti[thetaj]/k;
1098: k4=func(p2)-fx;
1099: res=(k1-k2-k3+k4)/4.0/delti[thetai]*k/delti[thetaj]*k/2.; /* Because of L not 2*L */
1100: #ifdef DEBUG
1101: 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);
1102: #endif
1103: }
1104: return res;
1105: }
1106:
1107: /************** Inverse of matrix **************/
1108: void ludcmp(double **a, int n, int *indx, double *d)
1109: {
1110: int i,imax,j,k;
1111: double big,dum,sum,temp;
1112: double *vv;
1113:
1114: vv=vector(1,n);
1115: *d=1.0;
1116: for (i=1;i<=n;i++) {
1117: big=0.0;
1118: for (j=1;j<=n;j++)
1119: if ((temp=fabs(a[i][j])) > big) big=temp;
1120: if (big == 0.0) nrerror("Singular matrix in routine ludcmp");
1121: vv[i]=1.0/big;
1122: }
1123: for (j=1;j<=n;j++) {
1124: for (i=1;i<j;i++) {
1125: sum=a[i][j];
1126: for (k=1;k<i;k++) sum -= a[i][k]*a[k][j];
1127: a[i][j]=sum;
1128: }
1129: big=0.0;
1130: for (i=j;i<=n;i++) {
1131: sum=a[i][j];
1132: for (k=1;k<j;k++)
1133: sum -= a[i][k]*a[k][j];
1134: a[i][j]=sum;
1135: if ( (dum=vv[i]*fabs(sum)) >= big) {
1136: big=dum;
1137: imax=i;
1138: }
1139: }
1140: if (j != imax) {
1141: for (k=1;k<=n;k++) {
1142: dum=a[imax][k];
1143: a[imax][k]=a[j][k];
1144: a[j][k]=dum;
1145: }
1146: *d = -(*d);
1147: vv[imax]=vv[j];
1148: }
1149: indx[j]=imax;
1150: if (a[j][j] == 0.0) a[j][j]=TINY;
1151: if (j != n) {
1152: dum=1.0/(a[j][j]);
1153: for (i=j+1;i<=n;i++) a[i][j] *= dum;
1154: }
1155: }
1156: free_vector(vv,1,n); /* Doesn't work */
1157: ;
1158: }
1159:
1160: void lubksb(double **a, int n, int *indx, double b[])
1161: {
1162: int i,ii=0,ip,j;
1163: double sum;
1164:
1165: for (i=1;i<=n;i++) {
1166: ip=indx[i];
1167: sum=b[ip];
1168: b[ip]=b[i];
1169: if (ii)
1170: for (j=ii;j<=i-1;j++) sum -= a[i][j]*b[j];
1171: else if (sum) ii=i;
1172: b[i]=sum;
1173: }
1174: for (i=n;i>=1;i--) {
1175: sum=b[i];
1176: for (j=i+1;j<=n;j++) sum -= a[i][j]*b[j];
1177: b[i]=sum/a[i][i];
1178: }
1179: }
1180:
1181: /************ Frequencies ********************/
1.26 lievre 1182: void freqsummary(char fileres[], int agemin, int agemax, int **s, double **agev, int nlstate, int imx, int *Tvar, int **nbcode, int *ncodemax,double **mint,double **anint, double dateprev1,double dateprev2,double jprev1, double mprev1,double anprev1,double jprev2, double mprev2,double anprev2)
1.2 lievre 1183: { /* Some frequencies */
1.35 lievre 1184:
1.18 lievre 1185: int i, m, jk, k1,i1, j1, bool, z1,z2,j;
1.2 lievre 1186: double ***freq; /* Frequencies */
1187: double *pp;
1.19 lievre 1188: double pos, k2, dateintsum=0,k2cpt=0;
1.2 lievre 1189: FILE *ficresp;
1190: char fileresp[FILENAMELENGTH];
1.35 lievre 1191:
1.2 lievre 1192: pp=vector(1,nlstate);
1.19 lievre 1193: probs= ma3x(1,AGESUP,1,NCOVMAX, 1,NCOVMAX);
1.2 lievre 1194: strcpy(fileresp,"p");
1195: strcat(fileresp,fileres);
1196: if((ficresp=fopen(fileresp,"w"))==NULL) {
1197: printf("Problem with prevalence resultfile: %s\n", fileresp);
1198: exit(0);
1199: }
1200: freq= ma3x(-1,nlstate+ndeath,-1,nlstate+ndeath,agemin,agemax+3);
1201: j1=0;
1.35 lievre 1202:
1.7 lievre 1203: j=cptcoveff;
1.2 lievre 1204: if (cptcovn<1) {j=1;ncodemax[1]=1;}
1.35 lievre 1205:
1.2 lievre 1206: for(k1=1; k1<=j;k1++){
1.35 lievre 1207: for(i1=1; i1<=ncodemax[k1];i1++){
1208: j1++;
1209: /*printf("cptcoveff=%d Tvaraff=%d", cptcoveff,Tvaraff[1]);
1210: scanf("%d", i);*/
1211: for (i=-1; i<=nlstate+ndeath; i++)
1212: for (jk=-1; jk<=nlstate+ndeath; jk++)
1213: for(m=agemin; m <= agemax+3; m++)
1214: freq[i][jk][m]=0;
1215:
1216: dateintsum=0;
1217: k2cpt=0;
1218: for (i=1; i<=imx; i++) {
1219: bool=1;
1220: if (cptcovn>0) {
1221: for (z1=1; z1<=cptcoveff; z1++)
1222: if (covar[Tvaraff[z1]][i]!= nbcode[Tvaraff[z1]][codtab[j1][z1]])
1223: bool=0;
1224: }
1225: if (bool==1) {
1226: for(m=firstpass; m<=lastpass; m++){
1227: k2=anint[m][i]+(mint[m][i]/12.);
1228: if ((k2>=dateprev1) && (k2<=dateprev2)) {
1229: if(agev[m][i]==0) agev[m][i]=agemax+1;
1230: if(agev[m][i]==1) agev[m][i]=agemax+2;
1231: if (m<lastpass) {
1232: freq[s[m][i]][s[m+1][i]][(int)agev[m][i]] += weight[i];
1233: freq[s[m][i]][s[m+1][i]][(int) agemax+3] += weight[i];
1234: }
1235:
1236: if ((agev[m][i]>1) && (agev[m][i]< (agemax+3))) {
1237: dateintsum=dateintsum+k2;
1238: k2cpt++;
1239: }
1240: }
1241: }
1242: }
1243: }
1.26 lievre 1244:
1.35 lievre 1245: fprintf(ficresp, "#Count between %.lf/%.lf/%.lf and %.lf/%.lf/%.lf\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2);
1.26 lievre 1246:
1.35 lievre 1247: if (cptcovn>0) {
1248: fprintf(ficresp, "\n#********** Variable ");
1249: for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresp, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
1250: fprintf(ficresp, "**********\n#");
1251: }
1252: for(i=1; i<=nlstate;i++)
1253: fprintf(ficresp, " Age Prev(%d) N(%d) N",i,i);
1254: fprintf(ficresp, "\n");
1255:
1256: for(i=(int)agemin; i <= (int)agemax+3; i++){
1257: if(i==(int)agemax+3)
1258: printf("Total");
1259: else
1260: printf("Age %d", i);
1261: for(jk=1; jk <=nlstate ; jk++){
1262: for(m=-1, pp[jk]=0; m <=nlstate+ndeath ; m++)
1263: pp[jk] += freq[jk][m][i];
1264: }
1265: for(jk=1; jk <=nlstate ; jk++){
1266: for(m=-1, pos=0; m <=0 ; m++)
1267: pos += freq[jk][m][i];
1268: if(pp[jk]>=1.e-10)
1269: printf(" %d.=%.0f loss[%d]=%.1f%%",jk,pp[jk],jk,100*pos/pp[jk]);
1270: else
1271: printf(" %d.=%.0f loss[%d]=NaNQ%%",jk,pp[jk],jk);
1272: }
1.14 lievre 1273:
1.35 lievre 1274: for(jk=1; jk <=nlstate ; jk++){
1275: for(m=0, pp[jk]=0; m <=nlstate+ndeath; m++)
1276: pp[jk] += freq[jk][m][i];
1277: }
1.14 lievre 1278:
1.35 lievre 1279: for(jk=1,pos=0; jk <=nlstate ; jk++)
1280: pos += pp[jk];
1281: for(jk=1; jk <=nlstate ; jk++){
1282: if(pos>=1.e-5)
1283: printf(" %d.=%.0f prev[%d]=%.1f%%",jk,pp[jk],jk,100*pp[jk]/pos);
1284: else
1285: printf(" %d.=%.0f prev[%d]=NaNQ%%",jk,pp[jk],jk);
1286: if( i <= (int) agemax){
1287: if(pos>=1.e-5){
1288: fprintf(ficresp," %d %.5f %.0f %.0f",i,pp[jk]/pos, pp[jk],pos);
1289: probs[i][jk][j1]= pp[jk]/pos;
1290: /*printf("\ni=%d jk=%d j1=%d %.5f %.0f %.0f %f",i,jk,j1,pp[jk]/pos, pp[jk],pos,probs[i][jk][j1]);*/
1291: }
1292: else
1293: fprintf(ficresp," %d NaNq %.0f %.0f",i,pp[jk],pos);
1294: }
1.13 lievre 1295: }
1.35 lievre 1296:
1297: for(jk=-1; jk <=nlstate+ndeath; jk++)
1298: for(m=-1; m <=nlstate+ndeath; m++)
1299: if(freq[jk][m][i] !=0 ) printf(" %d%d=%.0f",jk,m,freq[jk][m][i]);
1300: if(i <= (int) agemax)
1301: fprintf(ficresp,"\n");
1302: printf("\n");
1.2 lievre 1303: }
1304: }
1.35 lievre 1305: }
1.19 lievre 1306: dateintmean=dateintsum/k2cpt;
1.2 lievre 1307:
1308: fclose(ficresp);
1309: free_ma3x(freq,-1,nlstate+ndeath,-1,nlstate+ndeath,(int) agemin,(int) agemax+3);
1310: free_vector(pp,1,nlstate);
1.35 lievre 1311:
1.19 lievre 1312: /* End of Freq */
1313: }
1.2 lievre 1314:
1.15 lievre 1315: /************ Prevalence ********************/
1.28 lievre 1316: void prevalence(int agemin, float agemax, int **s, double **agev, int nlstate, int imx, int *Tvar, int **nbcode, int *ncodemax,double **mint,double **anint, double dateprev1,double dateprev2, double calagedate)
1.15 lievre 1317: { /* Some frequencies */
1318:
1319: int i, m, jk, k1, i1, j1, bool, z1,z2,j;
1320: double ***freq; /* Frequencies */
1321: double *pp;
1.18 lievre 1322: double pos, k2;
1.15 lievre 1323:
1324: pp=vector(1,nlstate);
1.19 lievre 1325: probs= ma3x(1,AGESUP,1,NCOVMAX, 1,NCOVMAX);
1.15 lievre 1326:
1327: freq=ma3x(-1,nlstate+ndeath,-1,nlstate+ndeath,agemin,agemax+3);
1328: j1=0;
1329:
1330: j=cptcoveff;
1331: if (cptcovn<1) {j=1;ncodemax[1]=1;}
1332:
1333: for(k1=1; k1<=j;k1++){
1334: for(i1=1; i1<=ncodemax[k1];i1++){
1335: j1++;
1336:
1337: for (i=-1; i<=nlstate+ndeath; i++)
1338: for (jk=-1; jk<=nlstate+ndeath; jk++)
1339: for(m=agemin; m <= agemax+3; m++)
1.19 lievre 1340: freq[i][jk][m]=0;
1.28 lievre 1341:
1.15 lievre 1342: for (i=1; i<=imx; i++) {
1343: bool=1;
1344: if (cptcovn>0) {
1345: for (z1=1; z1<=cptcoveff; z1++)
1346: if (covar[Tvaraff[z1]][i]!= nbcode[Tvaraff[z1]][codtab[j1][z1]])
1347: bool=0;
1.28 lievre 1348: }
1.19 lievre 1349: if (bool==1) {
1350: for(m=firstpass; m<=lastpass; m++){
1351: k2=anint[m][i]+(mint[m][i]/12.);
1352: if ((k2>=dateprev1) && (k2<=dateprev2)) {
1.18 lievre 1353: if(agev[m][i]==0) agev[m][i]=agemax+1;
1354: if(agev[m][i]==1) agev[m][i]=agemax+2;
1.35 lievre 1355: if (m<lastpass) freq[s[m][i]][s[m+1][i]][(int)(agev[m][i]+1-((int)calagedate %12)/12.)] += weight[i];
1.28 lievre 1356: /* freq[s[m][i]][s[m+1][i]][(int)(agemax+3+1)] += weight[i]; */
1.18 lievre 1357: }
1.15 lievre 1358: }
1359: }
1360: }
1.18 lievre 1361: for(i=(int)agemin; i <= (int)agemax+3; i++){
1362: for(jk=1; jk <=nlstate ; jk++){
1363: for(m=-1, pp[jk]=0; m <=nlstate+ndeath ; m++)
1364: pp[jk] += freq[jk][m][i];
1365: }
1366: for(jk=1; jk <=nlstate ; jk++){
1367: for(m=-1, pos=0; m <=0 ; m++)
1.15 lievre 1368: pos += freq[jk][m][i];
1369: }
1370:
1371: for(jk=1; jk <=nlstate ; jk++){
1372: for(m=0, pp[jk]=0; m <=nlstate+ndeath; m++)
1373: pp[jk] += freq[jk][m][i];
1374: }
1375:
1376: for(jk=1,pos=0; jk <=nlstate ; jk++) pos += pp[jk];
1377:
1378: for(jk=1; jk <=nlstate ; jk++){
1379: if( i <= (int) agemax){
1380: if(pos>=1.e-5){
1381: probs[i][jk][j1]= pp[jk]/pos;
1382: }
1383: }
1384: }
1385:
1.18 lievre 1386: }
1.15 lievre 1387: }
1388: }
1389:
1390:
1391: free_ma3x(freq,-1,nlstate+ndeath,-1,nlstate+ndeath,(int) agemin,(int) agemax+3);
1392: free_vector(pp,1,nlstate);
1393:
1394: } /* End of Freq */
1.19 lievre 1395:
1.2 lievre 1396: /************* Waves Concatenation ***************/
1397:
1398: void concatwav(int wav[], int **dh, int **mw, int **s, double *agedc, double **agev, int firstpass, int lastpass, int imx, int nlstate, int stepm)
1399: {
1400: /* Concatenates waves: wav[i] is the number of effective (useful waves) of individual i.
1401: Death is a valid wave (if date is known).
1402: mw[mi][i] is the mi (mi=1 to wav[i]) effective wave of individual i
1403: dh[m][i] of dh[mw[mi][i][i] is the delay between two effective waves m=mw[mi][i]
1404: and mw[mi+1][i]. dh depends on stepm.
1405: */
1406:
1407: int i, mi, m;
1.8 lievre 1408: /* int j, k=0,jk, ju, jl,jmin=1e+5, jmax=-1;
1409: double sum=0., jmean=0.;*/
1.2 lievre 1410:
1.11 lievre 1411: int j, k=0,jk, ju, jl;
1412: double sum=0.;
1413: jmin=1e+5;
1414: jmax=-1;
1415: jmean=0.;
1.2 lievre 1416: for(i=1; i<=imx; i++){
1417: mi=0;
1418: m=firstpass;
1419: while(s[m][i] <= nlstate){
1420: if(s[m][i]>=1)
1421: mw[++mi][i]=m;
1422: if(m >=lastpass)
1423: break;
1424: else
1425: m++;
1426: }/* end while */
1427: if (s[m][i] > nlstate){
1428: mi++; /* Death is another wave */
1429: /* if(mi==0) never been interviewed correctly before death */
1430: /* Only death is a correct wave */
1431: mw[mi][i]=m;
1432: }
1433:
1434: wav[i]=mi;
1435: if(mi==0)
1436: printf("Warning, no any valid information for:%d line=%d\n",num[i],i);
1437: }
1438:
1439: for(i=1; i<=imx; i++){
1440: for(mi=1; mi<wav[i];mi++){
1441: if (stepm <=0)
1442: dh[mi][i]=1;
1443: else{
1444: if (s[mw[mi+1][i]][i] > nlstate) {
1.10 lievre 1445: if (agedc[i] < 2*AGESUP) {
1.2 lievre 1446: j= rint(agedc[i]*12-agev[mw[mi][i]][i]*12);
1.8 lievre 1447: if(j==0) j=1; /* Survives at least one month after exam */
1448: k=k+1;
1449: if (j >= jmax) jmax=j;
1.11 lievre 1450: if (j <= jmin) jmin=j;
1.8 lievre 1451: sum=sum+j;
1.30 lievre 1452: /*if (j<0) printf("j=%d num=%d \n",j,i); */
1.10 lievre 1453: }
1.2 lievre 1454: }
1455: else{
1456: j= rint( (agev[mw[mi+1][i]][i]*12 - agev[mw[mi][i]][i]*12));
1457: k=k+1;
1458: if (j >= jmax) jmax=j;
1459: else if (j <= jmin)jmin=j;
1.30 lievre 1460: /* if (j<10) printf("j=%d jmin=%d num=%d ",j,jmin,i); */
1.2 lievre 1461: sum=sum+j;
1462: }
1463: jk= j/stepm;
1464: jl= j -jk*stepm;
1465: ju= j -(jk+1)*stepm;
1466: if(jl <= -ju)
1467: dh[mi][i]=jk;
1468: else
1469: dh[mi][i]=jk+1;
1470: if(dh[mi][i]==0)
1471: dh[mi][i]=1; /* At least one step */
1472: }
1473: }
1474: }
1.8 lievre 1475: jmean=sum/k;
1476: printf("Delay (in months) between two waves Min=%d Max=%d Mean=%f\n\n ",jmin, jmax,jmean);
1.12 lievre 1477: }
1.2 lievre 1478: /*********** Tricode ****************************/
1479: void tricode(int *Tvar, int **nbcode, int imx)
1480: {
1.7 lievre 1481: int Ndum[20],ij=1, k, j, i;
1.2 lievre 1482: int cptcode=0;
1.7 lievre 1483: cptcoveff=0;
1484:
1485: for (k=0; k<19; k++) Ndum[k]=0;
1.2 lievre 1486: for (k=1; k<=7; k++) ncodemax[k]=0;
1.6 lievre 1487:
1.7 lievre 1488: for (j=1; j<=(cptcovn+2*cptcovprod); j++) {
1.2 lievre 1489: for (i=1; i<=imx; i++) {
1490: ij=(int)(covar[Tvar[j]][i]);
1491: Ndum[ij]++;
1.8 lievre 1492: /*printf("i=%d ij=%d Ndum[ij]=%d imx=%d",i,ij,Ndum[ij],imx);*/
1.2 lievre 1493: if (ij > cptcode) cptcode=ij;
1494: }
1.7 lievre 1495:
1.2 lievre 1496: for (i=0; i<=cptcode; i++) {
1497: if(Ndum[i]!=0) ncodemax[j]++;
1498: }
1499: ij=1;
1.7 lievre 1500:
1.8 lievre 1501:
1.2 lievre 1502: for (i=1; i<=ncodemax[j]; i++) {
1.7 lievre 1503: for (k=0; k<=19; k++) {
1.2 lievre 1504: if (Ndum[k] != 0) {
1505: nbcode[Tvar[j]][ij]=k;
1.39 lievre 1506:
1.2 lievre 1507: ij++;
1508: }
1509: if (ij > ncodemax[j]) break;
1510: }
1511: }
1.7 lievre 1512: }
1.8 lievre 1513:
1514: for (k=0; k<19; k++) Ndum[k]=0;
1515:
1.12 lievre 1516: for (i=1; i<=ncovmodel-2; i++) {
1.7 lievre 1517: ij=Tvar[i];
1518: Ndum[ij]++;
1519: }
1.8 lievre 1520:
1.7 lievre 1521: ij=1;
1.8 lievre 1522: for (i=1; i<=10; i++) {
1.34 brouard 1523: if((Ndum[i]!=0) && (i<=ncovcol)){
1.8 lievre 1524: Tvaraff[ij]=i;
1525: ij++;
1.7 lievre 1526: }
1527: }
1528:
1.8 lievre 1529: cptcoveff=ij-1;
1.6 lievre 1530: }
1.2 lievre 1531:
1532: /*********** Health Expectancies ****************/
1533:
1.36 brouard 1534: void evsij(char fileres[], double ***eij, double x[], int nlstate, int stepm, int bage, int fage, double **oldm, double **savm, int ij, int estepm)
1.2 lievre 1535: {
1536: /* Health expectancies */
1.36 brouard 1537: int i, j, nhstepm, hstepm, h, nstepm;
1.35 lievre 1538: double age, agelim, hf;
1.2 lievre 1539: double ***p3mat;
1540:
1541: fprintf(ficreseij,"# Health expectancies\n");
1542: fprintf(ficreseij,"# Age");
1543: for(i=1; i<=nlstate;i++)
1544: for(j=1; j<=nlstate;j++)
1545: fprintf(ficreseij," %1d-%1d",i,j);
1546: fprintf(ficreseij,"\n");
1547:
1.36 brouard 1548: if(estepm < stepm){
1549: printf ("Problem %d lower than %d\n",estepm, stepm);
1550: }
1551: else hstepm=estepm;
1552: /* We compute the life expectancy from trapezoids spaced every estepm months
1553: * This is mainly to measure the difference between two models: for example
1554: * if stepm=24 months pijx are given only every 2 years and by summing them
1555: * we are calculating an estimate of the Life Expectancy assuming a linear
1556: * progression inbetween and thus overestimating or underestimating according
1557: * to the curvature of the survival function. If, for the same date, we
1558: * estimate the model with stepm=1 month, we can keep estepm to 24 months
1559: * to compare the new estimate of Life expectancy with the same linear
1560: * hypothesis. A more precise result, taking into account a more precise
1561: * curvature will be obtained if estepm is as small as stepm. */
1562:
1563: /* For example we decided to compute the life expectancy with the smallest unit */
1.31 brouard 1564: /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm.
1565: nhstepm is the number of hstepm from age to agelim
1566: nstepm is the number of stepm from age to agelin.
1567: Look at hpijx to understand the reason of that which relies in memory size
1.36 brouard 1568: and note for a fixed period like estepm months */
1.31 brouard 1569: /* We decided (b) to get a life expectancy respecting the most precise curvature of the
1.32 brouard 1570: survival function given by stepm (the optimization length). Unfortunately it
1.31 brouard 1571: means that if the survival funtion is printed only each two years of age and if
1572: you sum them up and add 1 year (area under the trapezoids) you won't get the same
1573: results. So we changed our mind and took the option of the best precision.
1574: */
1.36 brouard 1575: hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */
1.2 lievre 1576:
1577: agelim=AGESUP;
1578: for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
1579: /* nhstepm age range expressed in number of stepm */
1.31 brouard 1580: nstepm=(int) rint((agelim-age)*YEARM/stepm);
1581: /* Typically if 20 years nstepm = 20*12/6=40 stepm */
1.33 brouard 1582: /* if (stepm >= YEARM) hstepm=1;*/
1.31 brouard 1583: nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
1.2 lievre 1584: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1585: /* Computed by stepm unit matrices, product of hstepm matrices, stored
1586: in an array of nhstepm length: nhstepm=10, hstepm=4, stepm=6 months */
1587: hpxij(p3mat,nhstepm,age,hstepm,x,nlstate,stepm,oldm, savm, ij);
1.39 lievre 1588:
1589: /*for (h=0, eij[i][j][(int)age]=0; h<=nhstepm-1; h++) printf("%f %.5f\n", age*12+h, p3mat[1][1][h]);*/
1590:
1.32 brouard 1591: hf=hstepm*stepm/YEARM; /* Duration of hstepm expressed in year unit. */
1.2 lievre 1592: for(i=1; i<=nlstate;i++)
1593: for(j=1; j<=nlstate;j++)
1.29 lievre 1594: for (h=0, eij[i][j][(int)age]=0; h<=nhstepm-1; h++){
1.31 brouard 1595: eij[i][j][(int)age] += (p3mat[i][j][h]+p3mat[i][j][h+1])/2.0*hf;
1596: /* if((int)age==70)printf("i=%2d,j=%2d,h=%2d,age=%3d,%9.4f,%9.4f,%9.4f\n",i,j,h,(int)age,p3mat[i][j][h],hf,eij[i][j][(int)age]);*/
1.2 lievre 1597: }
1.29 lievre 1598: fprintf(ficreseij,"%3.0f",age );
1599: for(i=1; i<=nlstate;i++)
1600: for(j=1; j<=nlstate;j++){
1.30 lievre 1601: fprintf(ficreseij," %9.4f", eij[i][j][(int)age]);
1.2 lievre 1602: }
1603: fprintf(ficreseij,"\n");
1604: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1605: }
1606: }
1607:
1608: /************ Variance ******************/
1.36 brouard 1609: 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, int estepm)
1.2 lievre 1610: {
1611: /* Variance of health expectancies */
1612: /* double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double ** savm,double ftolpl);*/
1613: double **newm;
1614: double **dnewm,**doldm;
1.36 brouard 1615: int i, j, nhstepm, hstepm, h, nstepm ;
1.2 lievre 1616: int k, cptcode;
1.12 lievre 1617: double *xp;
1.2 lievre 1618: double **gp, **gm;
1619: double ***gradg, ***trgradg;
1620: double ***p3mat;
1.35 lievre 1621: double age,agelim, hf;
1.2 lievre 1622: int theta;
1623:
1624: fprintf(ficresvij,"# Covariances of life expectancies\n");
1625: fprintf(ficresvij,"# Age");
1626: for(i=1; i<=nlstate;i++)
1627: for(j=1; j<=nlstate;j++)
1628: fprintf(ficresvij," Cov(e%1d, e%1d)",i,j);
1629: fprintf(ficresvij,"\n");
1630:
1631: xp=vector(1,npar);
1632: dnewm=matrix(1,nlstate,1,npar);
1633: doldm=matrix(1,nlstate,1,nlstate);
1634:
1.36 brouard 1635: if(estepm < stepm){
1636: printf ("Problem %d lower than %d\n",estepm, stepm);
1637: }
1638: else hstepm=estepm;
1639: /* For example we decided to compute the life expectancy with the smallest unit */
1.35 lievre 1640: /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm.
1641: nhstepm is the number of hstepm from age to agelim
1642: nstepm is the number of stepm from age to agelin.
1643: Look at hpijx to understand the reason of that which relies in memory size
1644: and note for a fixed period like k years */
1645: /* We decided (b) to get a life expectancy respecting the most precise curvature of the
1646: survival function given by stepm (the optimization length). Unfortunately it
1647: means that if the survival funtion is printed only each two years of age and if
1648: you sum them up and add 1 year (area under the trapezoids) you won't get the same
1649: results. So we changed our mind and took the option of the best precision.
1650: */
1.36 brouard 1651: hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */
1.2 lievre 1652: agelim = AGESUP;
1653: for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
1.35 lievre 1654: nstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */
1655: nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
1.2 lievre 1656: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1657: gradg=ma3x(0,nhstepm,1,npar,1,nlstate);
1658: gp=matrix(0,nhstepm,1,nlstate);
1659: gm=matrix(0,nhstepm,1,nlstate);
1660:
1661: for(theta=1; theta <=npar; theta++){
1662: for(i=1; i<=npar; i++){ /* Computes gradient */
1663: xp[i] = x[i] + (i==theta ?delti[theta]:0);
1664: }
1665: hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);
1666: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1.14 lievre 1667:
1668: if (popbased==1) {
1669: for(i=1; i<=nlstate;i++)
1670: prlim[i][i]=probs[(int)age][i][ij];
1671: }
1.26 lievre 1672:
1.2 lievre 1673: for(j=1; j<= nlstate; j++){
1674: for(h=0; h<=nhstepm; h++){
1675: for(i=1, gp[h][j]=0.;i<=nlstate;i++)
1676: gp[h][j] += prlim[i][i]*p3mat[i][j][h];
1677: }
1678: }
1679:
1680: for(i=1; i<=npar; i++) /* Computes gradient */
1681: xp[i] = x[i] - (i==theta ?delti[theta]:0);
1682: hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);
1683: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1.26 lievre 1684:
1.14 lievre 1685: if (popbased==1) {
1686: for(i=1; i<=nlstate;i++)
1687: prlim[i][i]=probs[(int)age][i][ij];
1688: }
1689:
1.2 lievre 1690: for(j=1; j<= nlstate; j++){
1691: for(h=0; h<=nhstepm; h++){
1692: for(i=1, gm[h][j]=0.;i<=nlstate;i++)
1693: gm[h][j] += prlim[i][i]*p3mat[i][j][h];
1694: }
1695: }
1.14 lievre 1696:
1.2 lievre 1697: for(j=1; j<= nlstate; j++)
1698: for(h=0; h<=nhstepm; h++){
1699: gradg[h][theta][j]= (gp[h][j]-gm[h][j])/2./delti[theta];
1700: }
1701: } /* End theta */
1702:
1703: trgradg =ma3x(0,nhstepm,1,nlstate,1,npar);
1704:
1705: for(h=0; h<=nhstepm; h++)
1706: for(j=1; j<=nlstate;j++)
1707: for(theta=1; theta <=npar; theta++)
1708: trgradg[h][j][theta]=gradg[h][theta][j];
1709:
1.35 lievre 1710: hf=hstepm*stepm/YEARM; /* Duration of hstepm expressed in year unit. */
1.2 lievre 1711: for(i=1;i<=nlstate;i++)
1712: for(j=1;j<=nlstate;j++)
1713: vareij[i][j][(int)age] =0.;
1.35 lievre 1714:
1.2 lievre 1715: for(h=0;h<=nhstepm;h++){
1716: for(k=0;k<=nhstepm;k++){
1717: matprod2(dnewm,trgradg[h],1,nlstate,1,npar,1,npar,matcov);
1718: matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg[k]);
1719: for(i=1;i<=nlstate;i++)
1720: for(j=1;j<=nlstate;j++)
1.35 lievre 1721: vareij[i][j][(int)age] += doldm[i][j]*hf*hf;
1.2 lievre 1722: }
1723: }
1.35 lievre 1724:
1.2 lievre 1725: fprintf(ficresvij,"%.0f ",age );
1726: for(i=1; i<=nlstate;i++)
1727: for(j=1; j<=nlstate;j++){
1.35 lievre 1728: fprintf(ficresvij," %.4f", vareij[i][j][(int)age]);
1.2 lievre 1729: }
1730: fprintf(ficresvij,"\n");
1731: free_matrix(gp,0,nhstepm,1,nlstate);
1732: free_matrix(gm,0,nhstepm,1,nlstate);
1733: free_ma3x(gradg,0,nhstepm,1,npar,1,nlstate);
1734: free_ma3x(trgradg,0,nhstepm,1,nlstate,1,npar);
1735: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
1736: } /* End age */
1.26 lievre 1737:
1.2 lievre 1738: free_vector(xp,1,npar);
1739: free_matrix(doldm,1,nlstate,1,npar);
1740: free_matrix(dnewm,1,nlstate,1,nlstate);
1741:
1742: }
1743:
1744: /************ Variance of prevlim ******************/
1745: 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)
1746: {
1747: /* Variance of prevalence limit */
1748: /* double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double ** savm,double ftolpl);*/
1749: double **newm;
1750: double **dnewm,**doldm;
1751: int i, j, nhstepm, hstepm;
1752: int k, cptcode;
1753: double *xp;
1754: double *gp, *gm;
1755: double **gradg, **trgradg;
1756: double age,agelim;
1757: int theta;
1758:
1759: fprintf(ficresvpl,"# Standard deviation of prevalences limit\n");
1760: fprintf(ficresvpl,"# Age");
1761: for(i=1; i<=nlstate;i++)
1762: fprintf(ficresvpl," %1d-%1d",i,i);
1763: fprintf(ficresvpl,"\n");
1764:
1765: xp=vector(1,npar);
1766: dnewm=matrix(1,nlstate,1,npar);
1767: doldm=matrix(1,nlstate,1,nlstate);
1768:
1769: hstepm=1*YEARM; /* Every year of age */
1770: hstepm=hstepm/stepm; /* Typically in stepm units, if j= 2 years, = 2/6 months = 4 */
1771: agelim = AGESUP;
1772: for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
1773: nhstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */
1774: if (stepm >= YEARM) hstepm=1;
1775: nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
1776: gradg=matrix(1,npar,1,nlstate);
1777: gp=vector(1,nlstate);
1778: gm=vector(1,nlstate);
1779:
1780: for(theta=1; theta <=npar; theta++){
1781: for(i=1; i<=npar; i++){ /* Computes gradient */
1782: xp[i] = x[i] + (i==theta ?delti[theta]:0);
1783: }
1784: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1785: for(i=1;i<=nlstate;i++)
1786: gp[i] = prlim[i][i];
1787:
1788: for(i=1; i<=npar; i++) /* Computes gradient */
1789: xp[i] = x[i] - (i==theta ?delti[theta]:0);
1790: prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
1791: for(i=1;i<=nlstate;i++)
1792: gm[i] = prlim[i][i];
1793:
1794: for(i=1;i<=nlstate;i++)
1795: gradg[theta][i]= (gp[i]-gm[i])/2./delti[theta];
1796: } /* End theta */
1797:
1798: trgradg =matrix(1,nlstate,1,npar);
1799:
1800: for(j=1; j<=nlstate;j++)
1801: for(theta=1; theta <=npar; theta++)
1802: trgradg[j][theta]=gradg[theta][j];
1803:
1804: for(i=1;i<=nlstate;i++)
1805: varpl[i][(int)age] =0.;
1806: matprod2(dnewm,trgradg,1,nlstate,1,npar,1,npar,matcov);
1807: matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg);
1808: for(i=1;i<=nlstate;i++)
1809: varpl[i][(int)age] = doldm[i][i]; /* Covariances are useless */
1810:
1811: fprintf(ficresvpl,"%.0f ",age );
1812: for(i=1; i<=nlstate;i++)
1813: fprintf(ficresvpl," %.5f (%.5f)",prlim[i][i],sqrt(varpl[i][(int)age]));
1814: fprintf(ficresvpl,"\n");
1815: free_vector(gp,1,nlstate);
1816: free_vector(gm,1,nlstate);
1817: free_matrix(gradg,1,npar,1,nlstate);
1818: free_matrix(trgradg,1,nlstate,1,npar);
1819: } /* End age */
1820:
1821: free_vector(xp,1,npar);
1822: free_matrix(doldm,1,nlstate,1,npar);
1823: free_matrix(dnewm,1,nlstate,1,nlstate);
1824:
1825: }
1826:
1.13 lievre 1827: /************ Variance of one-step probabilities ******************/
1.39 lievre 1828: void varprob(char fileres[], double **matcov, double x[], double delti[], int nlstate, double bage, double fage, int ij, int *Tvar, int **nbcode, int *ncodemax)
1.13 lievre 1829: {
1.39 lievre 1830: int i, j, i1, k1, j1, z1;
1.13 lievre 1831: int k=0, cptcode;
1832: double **dnewm,**doldm;
1833: double *xp;
1834: double *gp, *gm;
1835: double **gradg, **trgradg;
1836: double age,agelim, cov[NCOVMAX];
1837: int theta;
1838: char fileresprob[FILENAMELENGTH];
1839:
1840: strcpy(fileresprob,"prob");
1841: strcat(fileresprob,fileres);
1842: if((ficresprob=fopen(fileresprob,"w"))==NULL) {
1843: printf("Problem with resultfile: %s\n", fileresprob);
1844: }
1845: printf("Computing variance of one-step probabilities: result on file '%s' \n",fileresprob);
1846:
1.2 lievre 1847:
1.13 lievre 1848: xp=vector(1,npar);
1849: dnewm=matrix(1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
1850: doldm=matrix(1,(nlstate+ndeath)*(nlstate+ndeath),1,(nlstate+ndeath)*(nlstate+ndeath));
1851:
1852: cov[1]=1;
1.39 lievre 1853: j=cptcoveff;
1854: if (cptcovn<1) {j=1;ncodemax[1]=1;}
1855: j1=0;
1856: for(k1=1; k1<=1;k1++){
1857: for(i1=1; i1<=ncodemax[k1];i1++){
1858: j1++;
1859:
1860: if (cptcovn>0) {
1861: fprintf(ficresprob, "\n#********** Variable ");
1862: for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprob, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
1863: fprintf(ficresprob, "**********\n#");
1864: }
1.13 lievre 1865:
1.39 lievre 1866: for (age=bage; age<=fage; age ++){
1867: cov[2]=age;
1868: for (k=1; k<=cptcovn;k++) {
1869: cov[2+k]=nbcode[Tvar[k]][codtab[j1][Tvar[k]]];
1870:
1.13 lievre 1871: }
1.39 lievre 1872: for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
1873: for (k=1; k<=cptcovprod;k++)
1874: cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
1875:
1876: gradg=matrix(1,npar,1,9);
1877: trgradg=matrix(1,9,1,npar);
1878: gp=vector(1,(nlstate+ndeath)*(nlstate+ndeath));
1879: gm=vector(1,(nlstate+ndeath)*(nlstate+ndeath));
1880:
1881: for(theta=1; theta <=npar; theta++){
1882: for(i=1; i<=npar; i++)
1883: xp[i] = x[i] + (i==theta ?delti[theta]:0);
1884:
1885: pmij(pmmij,cov,ncovmodel,xp,nlstate);
1886:
1887: k=0;
1888: for(i=1; i<= (nlstate+ndeath); i++){
1889: for(j=1; j<=(nlstate+ndeath);j++){
1890: k=k+1;
1891: gp[k]=pmmij[i][j];
1892: }
1893: }
1894:
1895: for(i=1; i<=npar; i++)
1896: xp[i] = x[i] - (i==theta ?delti[theta]:0);
1.13 lievre 1897:
1.39 lievre 1898: pmij(pmmij,cov,ncovmodel,xp,nlstate);
1899: k=0;
1900: for(i=1; i<=(nlstate+ndeath); i++){
1901: for(j=1; j<=(nlstate+ndeath);j++){
1902: k=k+1;
1903: gm[k]=pmmij[i][j];
1904: }
1905: }
1906:
1907: for(i=1; i<= (nlstate+ndeath)*(nlstate+ndeath); i++)
1908: gradg[theta][i]=(gp[i]-gm[i])/2./delti[theta];
1.13 lievre 1909: }
1910:
1.39 lievre 1911: for(j=1; j<=(nlstate+ndeath)*(nlstate+ndeath);j++)
1912: for(theta=1; theta <=npar; theta++)
1913: trgradg[j][theta]=gradg[theta][j];
1914:
1915: matprod2(dnewm,trgradg,1,9,1,npar,1,npar,matcov);
1916: matprod2(doldm,dnewm,1,9,1,npar,1,9,gradg);
1917:
1918: pmij(pmmij,cov,ncovmodel,x,nlstate);
1919:
1920: k=0;
1921: for(i=1; i<=(nlstate+ndeath); i++){
1922: for(j=1; j<=(nlstate+ndeath);j++){
1923: k=k+1;
1924: gm[k]=pmmij[i][j];
1925: }
1.13 lievre 1926: }
1927:
1928: /*printf("\n%d ",(int)age);
1929: for (i=1; i<=(nlstate+ndeath)*(nlstate+ndeath-1);i++){
1930: printf("%e [%e ;%e] ",gm[i],gm[i]-2*sqrt(doldm[i][i]),gm[i]+2*sqrt(doldm[i][i]));
1931: }*/
1932:
1.39 lievre 1933: fprintf(ficresprob,"\n%d ",(int)age);
1.13 lievre 1934:
1.39 lievre 1935: for (i=1; i<=(nlstate+ndeath)*(nlstate+ndeath-1);i++)
1936: fprintf(ficresprob,"%.3e (%.3e) ",gm[i],doldm[i][i]);
1937:
1938: }
1939: }
1.13 lievre 1940: free_vector(gp,1,(nlstate+ndeath)*(nlstate+ndeath));
1941: free_vector(gm,1,(nlstate+ndeath)*(nlstate+ndeath));
1942: free_matrix(trgradg,1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
1943: free_matrix(gradg,1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
1.39 lievre 1944: }
1945: free_vector(xp,1,npar);
1946: fclose(ficresprob);
1947:
1.13 lievre 1948: }
1.2 lievre 1949:
1.25 lievre 1950: /******************* Printing html file ***********/
1.35 lievre 1951: void printinghtml(char fileres[], char title[], char datafile[], int firstpass, \
1952: int lastpass, int stepm, int weightopt, char model[],\
1953: int imx,int jmin, int jmax, double jmeanint,char optionfile[], \
1954: char optionfilehtm[],char rfileres[], char optionfilegnuplot[],\
1.36 brouard 1955: char version[], int popforecast, int estepm ){
1.25 lievre 1956: int jj1, k1, i1, cpt;
1957: FILE *fichtm;
1958: /*char optionfilehtm[FILENAMELENGTH];*/
1959:
1960: strcpy(optionfilehtm,optionfile);
1961: strcat(optionfilehtm,".htm");
1962: if((fichtm=fopen(optionfilehtm,"w"))==NULL) {
1963: printf("Problem with %s \n",optionfilehtm), exit(0);
1964: }
1965:
1.37 brouard 1966: fprintf(fichtm,"<body> <font size=\"2\">%s </font> <hr size=\"2\" color=\"#EC5E5E\"> \n
1.35 lievre 1967: Title=%s <br>Datafile=%s Firstpass=%d Lastpass=%d Stepm=%d Weight=%d Model=%s<br>\n
1968: \n
1969: Total number of observations=%d <br>\n
1970: Interval (in months) between two waves: Min=%d Max=%d Mean=%.2lf<br>\n
1.25 lievre 1971: <hr size=\"2\" color=\"#EC5E5E\">
1.35 lievre 1972: <ul><li>Outputs files<br>\n
1973: - Copy of the parameter file: <a href=\"o%s\">o%s</a><br>\n
1974: - Gnuplot file name: <a href=\"%s\">%s</a><br>\n
1975: - Observed prevalence in each state: <a href=\"p%s\">p%s</a> <br>\n
1976: - Stationary prevalence in each state: <a href=\"pl%s\">pl%s</a> <br>\n
1977: - Transition probabilities: <a href=\"pij%s\">pij%s</a><br>\n
1.38 lievre 1978: - Life expectancies by age and initial health status (estepm=%2d months): <a href=\"e%s\">e%s</a> <br>\n",version,title,datafile,firstpass,lastpass,stepm, weightopt,model,imx,jmin,jmax,jmean,fileres,fileres,optionfilegnuplot,optionfilegnuplot,fileres,fileres,fileres,fileres,fileres,fileres,estepm,fileres,fileres);
1.35 lievre 1979:
1980: fprintf(fichtm,"\n
1981: - Parameter file with estimated parameters and the covariance matrix: <a href=\"%s\">%s</a> <br>\n
1.36 brouard 1982: - Variances of life expectancies by age and initial health status (estepm=%d months): <a href=\"v%s\">v%s</a><br>\n
1.35 lievre 1983: - Health expectancies with their variances: <a href=\"t%s\">t%s</a> <br>\n
1.36 brouard 1984: - Standard deviation of stationary prevalences: <a href=\"vpl%s\">vpl%s</a> <br>\n",rfileres,rfileres, estepm, fileres,fileres,fileres,fileres,fileres,fileres);
1.35 lievre 1985:
1986: if(popforecast==1) fprintf(fichtm,"\n
1987: - Prevalences forecasting: <a href=\"f%s\">f%s</a> <br>\n
1988: - Population forecasting (if popforecast=1): <a href=\"pop%s\">pop%s</a> <br>\n
1989: <br>",fileres,fileres,fileres,fileres);
1990: else
1991: fprintf(fichtm,"\n No population forecast: popforecast = %d (instead of 1) or stepm = %d (instead of 1) or model=%s (instead of .)<br><br></li>\n",popforecast, stepm, model);
1.25 lievre 1992: fprintf(fichtm," <li>Graphs</li><p>");
1993:
1994: m=cptcoveff;
1995: if (cptcovn < 1) {m=1;ncodemax[1]=1;}
1996:
1997: jj1=0;
1998: for(k1=1; k1<=m;k1++){
1999: for(i1=1; i1<=ncodemax[k1];i1++){
2000: jj1++;
2001: if (cptcovn > 0) {
2002: fprintf(fichtm,"<hr size=\"2\" color=\"#EC5E5E\">************ Results for covariates");
2003: for (cpt=1; cpt<=cptcoveff;cpt++)
2004: fprintf(fichtm," V%d=%d ",Tvaraff[cpt],nbcode[Tvaraff[cpt]][codtab[jj1][cpt]]);
2005: fprintf(fichtm," ************\n<hr size=\"2\" color=\"#EC5E5E\">");
2006: }
2007: fprintf(fichtm,"<br>- Probabilities: pe%s%d.gif<br>
2008: <img src=\"pe%s%d.gif\">",strtok(optionfile, "."),jj1,strtok(optionfile, "."),jj1);
2009: for(cpt=1; cpt<nlstate;cpt++){
2010: fprintf(fichtm,"<br>- Prevalence of disability : p%s%d%d.gif<br>
2011: <img src=\"p%s%d%d.gif\">",strtok(optionfile, "."),cpt,jj1,strtok(optionfile, "."),cpt,jj1);
2012: }
2013: for(cpt=1; cpt<=nlstate;cpt++) {
2014: fprintf(fichtm,"<br>- Observed and stationary prevalence (with confident
2015: interval) in state (%d): v%s%d%d.gif <br>
2016: <img src=\"v%s%d%d.gif\">",cpt,strtok(optionfile, "."),cpt,jj1,strtok(optionfile, "."),cpt,jj1);
2017: }
2018: for(cpt=1; cpt<=nlstate;cpt++) {
2019: fprintf(fichtm,"\n<br>- Health life expectancies by age and initial health state (%d): exp%s%d%d.gif <br>
2020: <img src=\"exp%s%d%d.gif\">",cpt,strtok(optionfile, "."),cpt,jj1,strtok(optionfile, "."),cpt,jj1);
2021: }
2022: fprintf(fichtm,"\n<br>- Total life expectancy by age and
2023: health expectancies in states (1) and (2): e%s%d.gif<br>
2024: <img src=\"e%s%d.gif\">",strtok(optionfile, "."),jj1,strtok(optionfile, "."),jj1);
2025: fprintf(fichtm,"\n</body>");
2026: }
2027: }
2028: fclose(fichtm);
2029: }
2030:
2031: /******************* Gnuplot file **************/
1.35 lievre 2032: void printinggnuplot(char fileres[],char optionfilefiname[],char optionfile[],char optionfilegnuplot[], double ageminpar, double agemaxpar, double fage , char pathc[], double p[]){
1.25 lievre 2033:
2034: int m,cpt,k1,i,k,j,jk,k2,k3,ij,l;
2035:
2036: strcpy(optionfilegnuplot,optionfilefiname);
1.35 lievre 2037: strcat(optionfilegnuplot,".gp.txt");
1.25 lievre 2038: if((ficgp=fopen(optionfilegnuplot,"w"))==NULL) {
2039: printf("Problem with file %s",optionfilegnuplot);
2040: }
2041:
2042: #ifdef windows
2043: fprintf(ficgp,"cd \"%s\" \n",pathc);
2044: #endif
2045: m=pow(2,cptcoveff);
2046:
2047: /* 1eme*/
2048: for (cpt=1; cpt<= nlstate ; cpt ++) {
2049: for (k1=1; k1<= m ; k1 ++) {
2050:
2051: #ifdef windows
1.35 lievre 2052: 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",ageminpar,fage,fileres,k1-1,k1-1);
1.25 lievre 2053: #endif
2054: #ifdef unix
1.35 lievre 2055: fprintf(ficgp,"set xlabel \"Age\" \nset ylabel \"Probability\" \nplot [%.f:%.f] \"vpl%s\" u 1:2 \"\%%lf",ageminpar,fage,fileres);
1.25 lievre 2056: #endif
2057:
2058: for (i=1; i<= nlstate ; i ++) {
2059: if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
2060: else fprintf(ficgp," \%%*lf (\%%*lf)");
2061: }
2062: fprintf(ficgp,"\" t\"Stationary prevalence\" w l 0,\"vpl%s\" every :::%d::%d u 1:($2+2*$3) \"\%%lf",fileres,k1-1,k1-1);
2063: for (i=1; i<= nlstate ; i ++) {
2064: if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
2065: else fprintf(ficgp," \%%*lf (\%%*lf)");
2066: }
2067: fprintf(ficgp,"\" t\"95\%% CI\" w l 1,\"vpl%s\" every :::%d::%d u 1:($2-2*$3) \"\%%lf",fileres,k1-1,k1-1);
2068: for (i=1; i<= nlstate ; i ++) {
2069: if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
2070: else fprintf(ficgp," \%%*lf (\%%*lf)");
2071: }
2072: 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));
2073: #ifdef unix
2074: fprintf(ficgp,"\nset ter gif small size 400,300");
2075: #endif
2076: fprintf(ficgp,"\nset out \"v%s%d%d.gif\" \nreplot\n\n",strtok(optionfile, "."),cpt,k1);
2077: }
2078: }
2079: /*2 eme*/
2080:
2081: for (k1=1; k1<= m ; k1 ++) {
1.35 lievre 2082: fprintf(ficgp,"set ylabel \"Years\" \nset ter gif small size 400,300\nplot [%.f:%.f] ",ageminpar,fage);
1.25 lievre 2083:
2084: for (i=1; i<= nlstate+1 ; i ++) {
2085: k=2*i;
2086: fprintf(ficgp,"\"t%s\" every :::%d::%d u 1:2 \"\%%lf",fileres,k1-1,k1-1);
2087: for (j=1; j<= nlstate+1 ; j ++) {
2088: if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
2089: else fprintf(ficgp," \%%*lf (\%%*lf)");
2090: }
2091: if (i== 1) fprintf(ficgp,"\" t\"TLE\" w l ,");
2092: else fprintf(ficgp,"\" t\"LE in state (%d)\" w l ,",i-1);
2093: fprintf(ficgp,"\"t%s\" every :::%d::%d u 1:($2-$3*2) \"\%%lf",fileres,k1-1,k1-1);
2094: for (j=1; j<= nlstate+1 ; j ++) {
2095: if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
2096: else fprintf(ficgp," \%%*lf (\%%*lf)");
2097: }
2098: fprintf(ficgp,"\" t\"\" w l 0,");
2099: fprintf(ficgp,"\"t%s\" every :::%d::%d u 1:($2+$3*2) \"\%%lf",fileres,k1-1,k1-1);
2100: for (j=1; j<= nlstate+1 ; j ++) {
2101: if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
2102: else fprintf(ficgp," \%%*lf (\%%*lf)");
2103: }
2104: if (i== (nlstate+1)) fprintf(ficgp,"\" t\"\" w l 0");
2105: else fprintf(ficgp,"\" t\"\" w l 0,");
2106: }
2107: fprintf(ficgp,"\nset out \"e%s%d.gif\" \nreplot\n\n",strtok(optionfile, "."),k1);
2108: }
2109:
2110: /*3eme*/
2111:
2112: for (k1=1; k1<= m ; k1 ++) {
2113: for (cpt=1; cpt<= nlstate ; cpt ++) {
2114: k=2+nlstate*(cpt-1);
1.35 lievre 2115: 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",ageminpar,fage,fileres,k1-1,k1-1,k,cpt);
1.25 lievre 2116: for (i=1; i< nlstate ; i ++) {
2117: 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);
2118: }
2119: fprintf(ficgp,"\nset out \"exp%s%d%d.gif\" \nreplot\n\n",strtok(optionfile, "."),cpt,k1);
2120: }
2121: }
2122:
2123: /* CV preval stat */
2124: for (k1=1; k1<= m ; k1 ++) {
2125: for (cpt=1; cpt<nlstate ; cpt ++) {
2126: k=3;
1.35 lievre 2127: 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",ageminpar,agemaxpar,fileres,k1,k+cpt+1,k+1);
1.25 lievre 2128:
2129: for (i=1; i< nlstate ; i ++)
2130: fprintf(ficgp,"+$%d",k+i+1);
2131: fprintf(ficgp,")) t\"prev(%d,%d)\" w l",cpt,cpt+1);
2132:
2133: l=3+(nlstate+ndeath)*cpt;
2134: fprintf(ficgp,",\"pij%s\" u ($1==%d ? ($3):1/0):($%d/($%d",fileres,k1,l+cpt+1,l+1);
2135: for (i=1; i< nlstate ; i ++) {
2136: l=3+(nlstate+ndeath)*cpt;
2137: fprintf(ficgp,"+$%d",l+i+1);
2138: }
2139: fprintf(ficgp,")) t\"prev(%d,%d)\" w l\n",cpt+1,cpt+1);
2140: fprintf(ficgp,"set out \"p%s%d%d.gif\" \nreplot\n\n",strtok(optionfile, "."),cpt,k1);
2141: }
2142: }
2143:
2144: /* proba elementaires */
2145: for(i=1,jk=1; i <=nlstate; i++){
2146: for(k=1; k <=(nlstate+ndeath); k++){
2147: if (k != i) {
2148: for(j=1; j <=ncovmodel; j++){
2149:
2150: fprintf(ficgp,"p%d=%f ",jk,p[jk]);
2151: jk++;
2152: fprintf(ficgp,"\n");
2153: }
2154: }
2155: }
2156: }
2157:
2158: for(jk=1; jk <=m; jk++) {
1.35 lievre 2159: fprintf(ficgp,"\nset ter gif small size 400,300\nset log y\nplot [%.f:%.f] ",ageminpar,agemaxpar);
1.25 lievre 2160: i=1;
2161: for(k2=1; k2<=nlstate; k2++) {
2162: k3=i;
2163: for(k=1; k<=(nlstate+ndeath); k++) {
2164: if (k != k2){
2165: fprintf(ficgp," exp(p%d+p%d*x",i,i+1);
2166: ij=1;
2167: for(j=3; j <=ncovmodel; j++) {
2168: if(((j-2)==Tage[ij]) &&(ij <=cptcovage)) {
2169: fprintf(ficgp,"+p%d*%d*x",i+j-1,nbcode[Tvar[j-2]][codtab[jk][Tvar[j-2]]]);
2170: ij++;
2171: }
2172: else
2173: fprintf(ficgp,"+p%d*%d",i+j-1,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
2174: }
2175: fprintf(ficgp,")/(1");
2176:
2177: for(k1=1; k1 <=nlstate; k1++){
2178: fprintf(ficgp,"+exp(p%d+p%d*x",k3+(k1-1)*ncovmodel,k3+(k1-1)*ncovmodel+1);
2179: ij=1;
2180: for(j=3; j <=ncovmodel; j++){
2181: if(((j-2)==Tage[ij]) &&(ij <=cptcovage)) {
2182: fprintf(ficgp,"+p%d*%d*x",k3+(k1-1)*ncovmodel+1+j-2,nbcode[Tvar[j-2]][codtab[jk][Tvar[j-2]]]);
2183: ij++;
2184: }
2185: else
2186: fprintf(ficgp,"+p%d*%d",k3+(k1-1)*ncovmodel+1+j-2,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
2187: }
2188: fprintf(ficgp,")");
2189: }
2190: fprintf(ficgp,") t \"p%d%d\" ", k2,k);
2191: if ((k+k2)!= (nlstate*2+ndeath)) fprintf(ficgp,",");
2192: i=i+ncovmodel;
2193: }
2194: }
2195: }
2196: fprintf(ficgp,"\nset out \"pe%s%d.gif\" \nreplot\n\n",strtok(optionfile, "."),jk);
2197: }
2198:
2199: fclose(ficgp);
2200: } /* end gnuplot */
2201:
2202:
2203: /*************** Moving average **************/
1.35 lievre 2204: void movingaverage(double agedeb, double fage,double ageminpar, double ***mobaverage){
1.25 lievre 2205:
2206: int i, cpt, cptcod;
1.35 lievre 2207: for (agedeb=ageminpar; agedeb<=fage; agedeb++)
1.25 lievre 2208: for (i=1; i<=nlstate;i++)
2209: for (cptcod=1;cptcod<=ncodemax[cptcov];cptcod++)
2210: mobaverage[(int)agedeb][i][cptcod]=0.;
2211:
1.35 lievre 2212: for (agedeb=ageminpar+4; agedeb<=fage; agedeb++){
1.25 lievre 2213: for (i=1; i<=nlstate;i++){
2214: for (cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
2215: for (cpt=0;cpt<=4;cpt++){
2216: mobaverage[(int)agedeb-2][i][cptcod]=mobaverage[(int)agedeb-2][i][cptcod]+probs[(int)agedeb-cpt][i][cptcod];
2217: }
2218: mobaverage[(int)agedeb-2][i][cptcod]=mobaverage[(int)agedeb-2][i][cptcod]/5;
2219: }
2220: }
2221: }
2222:
2223: }
2224:
1.27 lievre 2225:
2226: /************** Forecasting ******************/
1.35 lievre 2227: prevforecast(char fileres[], double anproj1,double mproj1,double jproj1,double ageminpar, double agemax,double dateprev1, double dateprev2, int mobilav, double agedeb, double fage, int popforecast, char popfile[], double anproj2,double p[], int i2){
1.27 lievre 2228:
2229: int cpt, stepsize, hstepm, nhstepm, j,k,c, cptcod, i,h;
2230: int *popage;
2231: double calagedate, agelim, kk1, kk2, yp,yp1,yp2,jprojmean,mprojmean,anprojmean;
2232: double *popeffectif,*popcount;
2233: double ***p3mat;
2234: char fileresf[FILENAMELENGTH];
2235:
2236: agelim=AGESUP;
2237: calagedate=(anproj1+mproj1/12.+jproj1/365.-dateintmean)*YEARM;
1.28 lievre 2238:
1.35 lievre 2239: prevalence(ageminpar, agemax, s, agev, nlstate, imx,Tvar,nbcode, ncodemax,mint,anint,dateprev1,dateprev2, calagedate);
1.28 lievre 2240:
1.27 lievre 2241:
2242: strcpy(fileresf,"f");
2243: strcat(fileresf,fileres);
2244: if((ficresf=fopen(fileresf,"w"))==NULL) {
2245: printf("Problem with forecast resultfile: %s\n", fileresf);
2246: }
2247: printf("Computing forecasting: result on file '%s' \n", fileresf);
2248:
2249: if (cptcoveff==0) ncodemax[cptcoveff]=1;
2250:
2251: if (mobilav==1) {
2252: mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
1.35 lievre 2253: movingaverage(agedeb, fage, ageminpar, mobaverage);
1.27 lievre 2254: }
2255:
2256: stepsize=(int) (stepm+YEARM-1)/YEARM;
2257: if (stepm<=12) stepsize=1;
2258:
2259: agelim=AGESUP;
2260:
2261: hstepm=1;
2262: hstepm=hstepm/stepm;
2263: yp1=modf(dateintmean,&yp);
2264: anprojmean=yp;
2265: yp2=modf((yp1*12),&yp);
2266: mprojmean=yp;
2267: yp1=modf((yp2*30.5),&yp);
2268: jprojmean=yp;
2269: if(jprojmean==0) jprojmean=1;
2270: if(mprojmean==0) jprojmean=1;
2271:
2272: fprintf(ficresf,"# Estimated date of observed prevalence: %.lf/%.lf/%.lf ",jprojmean,mprojmean,anprojmean);
2273:
2274: for(cptcov=1;cptcov<=i2;cptcov++){
2275: for(cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
2276: k=k+1;
2277: fprintf(ficresf,"\n#******");
2278: for(j=1;j<=cptcoveff;j++) {
2279: fprintf(ficresf," V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
2280: }
2281: fprintf(ficresf,"******\n");
2282: fprintf(ficresf,"# StartingAge FinalAge");
2283: for(j=1; j<=nlstate+ndeath;j++) fprintf(ficresf," P.%d",j);
2284:
2285:
2286: for (cpt=0; cpt<=(anproj2-anproj1);cpt++) {
2287: fprintf(ficresf,"\n");
2288: fprintf(ficresf,"\n# Forecasting at date %.lf/%.lf/%.lf ",jproj1,mproj1,anproj1+cpt);
1.28 lievre 2289:
1.35 lievre 2290: for (agedeb=(fage-((int)calagedate %12/12.)); agedeb>=(ageminpar-((int)calagedate %12)/12.); agedeb--){
1.27 lievre 2291: nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm);
2292: nhstepm = nhstepm/hstepm;
2293:
2294: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2295: oldm=oldms;savm=savms;
2296: hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);
2297:
2298: for (h=0; h<=nhstepm; h++){
2299: if (h==(int) (calagedate+YEARM*cpt)) {
1.35 lievre 2300: fprintf(ficresf,"\n %.f %.f ",anproj1+cpt,agedeb+h*hstepm/YEARM*stepm);
1.27 lievre 2301: }
2302: for(j=1; j<=nlstate+ndeath;j++) {
2303: kk1=0.;kk2=0;
2304: for(i=1; i<=nlstate;i++) {
2305: if (mobilav==1)
2306: kk1=kk1+p3mat[i][j][h]*mobaverage[(int)agedeb+1][i][cptcod];
2307: else {
2308: kk1=kk1+p3mat[i][j][h]*probs[(int)(agedeb+1)][i][cptcod];
2309: }
2310:
2311: }
2312: if (h==(int)(calagedate+12*cpt)){
2313: fprintf(ficresf," %.3f", kk1);
2314:
2315: }
2316: }
2317: }
2318: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2319: }
2320: }
2321: }
2322: }
2323:
2324: if (mobilav==1) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
2325:
2326: fclose(ficresf);
2327: }
2328: /************** Forecasting ******************/
1.35 lievre 2329: populforecast(char fileres[], double anpyram,double mpyram,double jpyram,double ageminpar, double agemax,double dateprev1, double dateprev2, int mobilav, double agedeb, double fage, int popforecast, char popfile[], double anpyram1,double p[], int i2){
1.27 lievre 2330:
2331: int cpt, stepsize, hstepm, nhstepm, j,k,c, cptcod, i,h;
2332: int *popage;
2333: double calagedate, agelim, kk1, kk2, yp,yp1,yp2,jprojmean,mprojmean,anprojmean;
2334: double *popeffectif,*popcount;
2335: double ***p3mat,***tabpop,***tabpopprev;
2336: char filerespop[FILENAMELENGTH];
2337:
1.28 lievre 2338: tabpop= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
2339: tabpopprev= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
2340: agelim=AGESUP;
2341: calagedate=(anpyram+mpyram/12.+jpyram/365.-dateintmean)*YEARM;
1.27 lievre 2342:
1.35 lievre 2343: prevalence(ageminpar, agemax, s, agev, nlstate, imx,Tvar,nbcode, ncodemax,mint,anint,dateprev1,dateprev2, calagedate);
1.28 lievre 2344:
2345:
1.27 lievre 2346: strcpy(filerespop,"pop");
2347: strcat(filerespop,fileres);
2348: if((ficrespop=fopen(filerespop,"w"))==NULL) {
2349: printf("Problem with forecast resultfile: %s\n", filerespop);
2350: }
2351: printf("Computing forecasting: result on file '%s' \n", filerespop);
2352:
2353: if (cptcoveff==0) ncodemax[cptcoveff]=1;
2354:
2355: if (mobilav==1) {
2356: mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
1.35 lievre 2357: movingaverage(agedeb, fage, ageminpar, mobaverage);
1.27 lievre 2358: }
2359:
2360: stepsize=(int) (stepm+YEARM-1)/YEARM;
2361: if (stepm<=12) stepsize=1;
2362:
2363: agelim=AGESUP;
2364:
2365: hstepm=1;
2366: hstepm=hstepm/stepm;
2367:
2368: if (popforecast==1) {
2369: if((ficpop=fopen(popfile,"r"))==NULL) {
2370: printf("Problem with population file : %s\n",popfile);exit(0);
2371: }
2372: popage=ivector(0,AGESUP);
2373: popeffectif=vector(0,AGESUP);
2374: popcount=vector(0,AGESUP);
2375:
2376: i=1;
2377: while ((c=fscanf(ficpop,"%d %lf\n",&popage[i],&popcount[i])) != EOF) i=i+1;
2378:
2379: imx=i;
2380: for (i=1; i<imx;i++) popeffectif[popage[i]]=popcount[i];
2381: }
2382:
2383: for(cptcov=1;cptcov<=i2;cptcov++){
2384: for(cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
2385: k=k+1;
2386: fprintf(ficrespop,"\n#******");
2387: for(j=1;j<=cptcoveff;j++) {
2388: fprintf(ficrespop," V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
2389: }
2390: fprintf(ficrespop,"******\n");
1.28 lievre 2391: fprintf(ficrespop,"# Age");
1.27 lievre 2392: for(j=1; j<=nlstate+ndeath;j++) fprintf(ficrespop," P.%d",j);
2393: if (popforecast==1) fprintf(ficrespop," [Population]");
2394:
2395: for (cpt=0; cpt<=0;cpt++) {
1.28 lievre 2396: fprintf(ficrespop,"\n\n# Forecasting at date %.lf/%.lf/%.lf ",jpyram,mpyram,anpyram+cpt);
1.27 lievre 2397:
1.35 lievre 2398: for (agedeb=(fage-((int)calagedate %12/12.)); agedeb>=(ageminpar-((int)calagedate %12)/12.); agedeb--){
1.27 lievre 2399: nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm);
2400: nhstepm = nhstepm/hstepm;
2401:
2402: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2403: oldm=oldms;savm=savms;
2404: hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);
2405:
2406: for (h=0; h<=nhstepm; h++){
2407: if (h==(int) (calagedate+YEARM*cpt)) {
1.28 lievre 2408: fprintf(ficrespop,"\n %3.f ",agedeb+h*hstepm/YEARM*stepm);
1.27 lievre 2409: }
2410: for(j=1; j<=nlstate+ndeath;j++) {
2411: kk1=0.;kk2=0;
2412: for(i=1; i<=nlstate;i++) {
2413: if (mobilav==1)
2414: kk1=kk1+p3mat[i][j][h]*mobaverage[(int)agedeb+1][i][cptcod];
2415: else {
2416: kk1=kk1+p3mat[i][j][h]*probs[(int)(agedeb+1)][i][cptcod];
2417: }
2418: }
2419: if (h==(int)(calagedate+12*cpt)){
2420: tabpop[(int)(agedeb)][j][cptcod]=kk1;
2421: /*fprintf(ficrespop," %.3f", kk1);
2422: if (popforecast==1) fprintf(ficrespop," [%.f]", kk1*popeffectif[(int)agedeb+1]);*/
2423: }
2424: }
2425: for(i=1; i<=nlstate;i++){
2426: kk1=0.;
2427: for(j=1; j<=nlstate;j++){
1.28 lievre 2428: kk1= kk1+tabpop[(int)(agedeb)][j][cptcod];
1.27 lievre 2429: }
2430: tabpopprev[(int)(agedeb)][i][cptcod]=tabpop[(int)(agedeb)][i][cptcod]/kk1*popeffectif[(int)(agedeb+(calagedate+12*cpt)*hstepm/YEARM*stepm-1)];
2431: }
2432:
1.28 lievre 2433: if (h==(int)(calagedate+12*cpt)) for(j=1; j<=nlstate;j++)
2434: fprintf(ficrespop," %15.2f",tabpopprev[(int)(agedeb+1)][j][cptcod]);
1.27 lievre 2435: }
2436: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2437: }
2438: }
2439:
2440: /******/
2441:
1.28 lievre 2442: for (cpt=1; cpt<=(anpyram1-anpyram);cpt++) {
2443: fprintf(ficrespop,"\n\n# Forecasting at date %.lf/%.lf/%.lf ",jpyram,mpyram,anpyram+cpt);
1.35 lievre 2444: for (agedeb=(fage-((int)calagedate %12/12.)); agedeb>=(ageminpar-((int)calagedate %12)/12.); agedeb--){
1.27 lievre 2445: nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm);
2446: nhstepm = nhstepm/hstepm;
2447:
2448: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2449: oldm=oldms;savm=savms;
2450: hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);
2451: for (h=0; h<=nhstepm; h++){
2452: if (h==(int) (calagedate+YEARM*cpt)) {
1.28 lievre 2453: fprintf(ficresf,"\n %3.f ",agedeb+h*hstepm/YEARM*stepm);
1.27 lievre 2454: }
2455: for(j=1; j<=nlstate+ndeath;j++) {
2456: kk1=0.;kk2=0;
2457: for(i=1; i<=nlstate;i++) {
2458: kk1=kk1+p3mat[i][j][h]*tabpopprev[(int)agedeb+1][i][cptcod];
2459: }
1.28 lievre 2460: if (h==(int)(calagedate+12*cpt)) fprintf(ficresf," %15.2f", kk1);
1.27 lievre 2461: }
2462: }
2463: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
2464: }
2465: }
2466: }
2467: }
2468:
2469: if (mobilav==1) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
2470:
2471: if (popforecast==1) {
2472: free_ivector(popage,0,AGESUP);
2473: free_vector(popeffectif,0,AGESUP);
2474: free_vector(popcount,0,AGESUP);
2475: }
2476: free_ma3x(tabpop,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
2477: free_ma3x(tabpopprev,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
2478: fclose(ficrespop);
2479: }
2480:
1.2 lievre 2481: /***********************************************/
2482: /**************** Main Program *****************/
2483: /***********************************************/
2484:
1.22 brouard 2485: int main(int argc, char *argv[])
1.2 lievre 2486: {
2487:
1.8 lievre 2488: int i,j, k, n=MAXN,iter,m,size,cptcode, cptcod;
1.2 lievre 2489: double agedeb, agefin,hf;
1.35 lievre 2490: double ageminpar=1.e20,agemin=1.e20, agemaxpar=-1.e20, agemax=-1.e20;
1.2 lievre 2491:
2492: double fret;
2493: double **xi,tmp,delta;
2494:
2495: double dum; /* Dummy variable */
2496: double ***p3mat;
2497: int *indx;
2498: char line[MAXLINE], linepar[MAXLINE];
2499: char title[MAXLINE];
1.25 lievre 2500: char optionfile[FILENAMELENGTH], datafile[FILENAMELENGTH], filerespl[FILENAMELENGTH];
2501: char optionfilext[10], optionfilefiname[FILENAMELENGTH], optionfilehtm[FILENAMELENGTH], optionfilegnuplot[FILENAMELENGTH], plotcmd[FILENAMELENGTH];
1.22 brouard 2502:
1.29 lievre 2503: char fileres[FILENAMELENGTH], filerespij[FILENAMELENGTH], filereso[FILENAMELENGTH], rfileres[FILENAMELENGTH];
1.22 brouard 2504:
1.2 lievre 2505: char filerest[FILENAMELENGTH];
2506: char fileregp[FILENAMELENGTH];
1.16 lievre 2507: char popfile[FILENAMELENGTH];
1.2 lievre 2508: char path[80],pathc[80],pathcd[80],pathtot[80],model[20];
2509: int firstobs=1, lastobs=10;
2510: int sdeb, sfin; /* Status at beginning and end */
2511: int c, h , cpt,l;
2512: int ju,jl, mi;
1.7 lievre 2513: int i1,j1, k1,k2,k3,jk,aa,bb, stepsize, ij;
1.14 lievre 2514: int jnais,jdc,jint4,jint1,jint2,jint3,**outcome,**adl,*tab;
1.19 lievre 2515: int mobilav=0,popforecast=0;
1.2 lievre 2516: int hstepm, nhstepm;
1.28 lievre 2517: double jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,jpyram, mpyram,anpyram,jpyram1, mpyram1,anpyram1;
1.14 lievre 2518:
1.2 lievre 2519: double bage, fage, age, agelim, agebase;
2520: double ftolpl=FTOL;
2521: double **prlim;
2522: double *severity;
2523: double ***param; /* Matrix of parameters */
2524: double *p;
2525: double **matcov; /* Matrix of covariance */
2526: double ***delti3; /* Scale */
2527: double *delti; /* Scale */
2528: double ***eij, ***vareij;
2529: double **varpl; /* Variances of prevalence limits by age */
2530: double *epj, vepp;
1.16 lievre 2531: double kk1, kk2;
1.27 lievre 2532: double dateprev1, dateprev2,jproj1,mproj1,anproj1,jproj2,mproj2,anproj2;
2533:
1.13 lievre 2534:
1.37 brouard 2535: char version[80]="Imach version 0.8b, March 2002, INED-EUROREVES ";
1.2 lievre 2536: char *alph[]={"a","a","b","c","d","e"}, str[4];
1.5 lievre 2537:
1.13 lievre 2538:
1.2 lievre 2539: char z[1]="c", occ;
2540: #include <sys/time.h>
2541: #include <time.h>
2542: char stra[80], strb[80], strc[80], strd[80],stre[80],modelsav[80];
1.19 lievre 2543:
1.2 lievre 2544: /* long total_usecs;
2545: struct timeval start_time, end_time;
2546:
2547: gettimeofday(&start_time, (struct timezone*)0); */ /* at first time */
1.35 lievre 2548: getcwd(pathcd, size);
1.2 lievre 2549:
1.22 brouard 2550: printf("\n%s",version);
2551: if(argc <=1){
2552: printf("\nEnter the parameter file name: ");
2553: scanf("%s",pathtot);
2554: }
2555: else{
2556: strcpy(pathtot,argv[1]);
2557: }
2558: /*if(getcwd(pathcd, 80)!= NULL)printf ("Error pathcd\n");*/
1.5 lievre 2559: /*cygwin_split_path(pathtot,path,optionfile);
2560: printf("pathtot=%s, path=%s, optionfile=%s\n",pathtot,path,optionfile);*/
2561: /* cutv(path,optionfile,pathtot,'\\');*/
2562:
1.22 brouard 2563: split(pathtot,path,optionfile,optionfilext,optionfilefiname);
2564: printf("pathtot=%s, path=%s, optionfile=%s optionfilext=%s optionfilefiname=%s\n",pathtot,path,optionfile,optionfilext,optionfilefiname);
1.2 lievre 2565: chdir(path);
2566: replace(pathc,path);
2567:
2568: /*-------- arguments in the command line --------*/
2569:
2570: strcpy(fileres,"r");
1.22 brouard 2571: strcat(fileres, optionfilefiname);
2572: strcat(fileres,".txt"); /* Other files have txt extension */
1.2 lievre 2573:
2574: /*---------arguments file --------*/
2575:
2576: if((ficpar=fopen(optionfile,"r"))==NULL) {
2577: printf("Problem with optionfile %s\n",optionfile);
2578: goto end;
2579: }
2580:
2581: strcpy(filereso,"o");
2582: strcat(filereso,fileres);
2583: if((ficparo=fopen(filereso,"w"))==NULL) {
2584: printf("Problem with Output resultfile: %s\n", filereso);goto end;
2585: }
2586:
2587: /* Reads comments: lines beginning with '#' */
2588: while((c=getc(ficpar))=='#' && c!= EOF){
2589: ungetc(c,ficpar);
2590: fgets(line, MAXLINE, ficpar);
2591: puts(line);
2592: fputs(line,ficparo);
2593: }
2594: ungetc(c,ficpar);
2595:
1.34 brouard 2596: fscanf(ficpar,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%lf stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n",title, datafile, &lastobs, &firstpass,&lastpass,&ftol, &stepm, &ncovcol, &nlstate,&ndeath, &maxwav, &mle, &weightopt,model);
2597: printf("title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol, stepm, ncovcol, nlstate,ndeath, maxwav, mle, weightopt,model);
2598: fprintf(ficparo,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol,stepm,ncovcol,nlstate,ndeath,maxwav, mle, weightopt,model);
1.14 lievre 2599: while((c=getc(ficpar))=='#' && c!= EOF){
2600: ungetc(c,ficpar);
2601: fgets(line, MAXLINE, ficpar);
2602: puts(line);
2603: fputs(line,ficparo);
2604: }
2605: ungetc(c,ficpar);
2606:
1.19 lievre 2607:
1.8 lievre 2608: covar=matrix(0,NCOVMAX,1,n);
2609: cptcovn=0;
2610: if (strlen(model)>1) cptcovn=nbocc(model,'+')+1;
1.2 lievre 2611:
2612: ncovmodel=2+cptcovn;
2613: nvar=ncovmodel-1; /* Suppressing age as a basic covariate */
2614:
2615: /* Read guess parameters */
2616: /* Reads comments: lines beginning with '#' */
2617: while((c=getc(ficpar))=='#' && c!= EOF){
2618: ungetc(c,ficpar);
2619: fgets(line, MAXLINE, ficpar);
2620: puts(line);
2621: fputs(line,ficparo);
2622: }
2623: ungetc(c,ficpar);
2624:
2625: param= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
2626: for(i=1; i <=nlstate; i++)
2627: for(j=1; j <=nlstate+ndeath-1; j++){
2628: fscanf(ficpar,"%1d%1d",&i1,&j1);
2629: fprintf(ficparo,"%1d%1d",i1,j1);
2630: printf("%1d%1d",i,j);
2631: for(k=1; k<=ncovmodel;k++){
2632: fscanf(ficpar," %lf",¶m[i][j][k]);
2633: printf(" %lf",param[i][j][k]);
2634: fprintf(ficparo," %lf",param[i][j][k]);
2635: }
2636: fscanf(ficpar,"\n");
2637: printf("\n");
2638: fprintf(ficparo,"\n");
2639: }
2640:
1.12 lievre 2641: npar= (nlstate+ndeath-1)*nlstate*ncovmodel;
2642:
1.2 lievre 2643: p=param[1][1];
2644:
2645: /* Reads comments: lines beginning with '#' */
2646: while((c=getc(ficpar))=='#' && c!= EOF){
2647: ungetc(c,ficpar);
2648: fgets(line, MAXLINE, ficpar);
2649: puts(line);
2650: fputs(line,ficparo);
2651: }
2652: ungetc(c,ficpar);
2653:
2654: delti3= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
2655: delti=vector(1,npar); /* Scale of each paramater (output from hesscov) */
2656: for(i=1; i <=nlstate; i++){
2657: for(j=1; j <=nlstate+ndeath-1; j++){
2658: fscanf(ficpar,"%1d%1d",&i1,&j1);
2659: printf("%1d%1d",i,j);
2660: fprintf(ficparo,"%1d%1d",i1,j1);
2661: for(k=1; k<=ncovmodel;k++){
2662: fscanf(ficpar,"%le",&delti3[i][j][k]);
2663: printf(" %le",delti3[i][j][k]);
2664: fprintf(ficparo," %le",delti3[i][j][k]);
2665: }
2666: fscanf(ficpar,"\n");
2667: printf("\n");
2668: fprintf(ficparo,"\n");
2669: }
2670: }
2671: delti=delti3[1][1];
2672:
2673: /* Reads comments: lines beginning with '#' */
2674: while((c=getc(ficpar))=='#' && c!= EOF){
2675: ungetc(c,ficpar);
2676: fgets(line, MAXLINE, ficpar);
2677: puts(line);
2678: fputs(line,ficparo);
2679: }
2680: ungetc(c,ficpar);
2681:
2682: matcov=matrix(1,npar,1,npar);
2683: for(i=1; i <=npar; i++){
2684: fscanf(ficpar,"%s",&str);
2685: printf("%s",str);
2686: fprintf(ficparo,"%s",str);
2687: for(j=1; j <=i; j++){
2688: fscanf(ficpar," %le",&matcov[i][j]);
2689: printf(" %.5le",matcov[i][j]);
2690: fprintf(ficparo," %.5le",matcov[i][j]);
2691: }
2692: fscanf(ficpar,"\n");
2693: printf("\n");
2694: fprintf(ficparo,"\n");
2695: }
2696: for(i=1; i <=npar; i++)
2697: for(j=i+1;j<=npar;j++)
2698: matcov[i][j]=matcov[j][i];
2699:
2700: printf("\n");
2701:
2702:
1.29 lievre 2703: /*-------- Rewriting paramater file ----------*/
2704: strcpy(rfileres,"r"); /* "Rparameterfile */
2705: strcat(rfileres,optionfilefiname); /* Parameter file first name*/
2706: strcat(rfileres,"."); /* */
2707: strcat(rfileres,optionfilext); /* Other files have txt extension */
2708: if((ficres =fopen(rfileres,"w"))==NULL) {
2709: printf("Problem writing new parameter file: %s\n", fileres);goto end;
1.2 lievre 2710: }
2711: fprintf(ficres,"#%s\n",version);
2712:
1.29 lievre 2713: /*-------- data file ----------*/
1.2 lievre 2714: if((fic=fopen(datafile,"r"))==NULL) {
2715: printf("Problem with datafile: %s\n", datafile);goto end;
2716: }
2717:
2718: n= lastobs;
2719: severity = vector(1,maxwav);
2720: outcome=imatrix(1,maxwav+1,1,n);
2721: num=ivector(1,n);
2722: moisnais=vector(1,n);
2723: annais=vector(1,n);
2724: moisdc=vector(1,n);
2725: andc=vector(1,n);
2726: agedc=vector(1,n);
2727: cod=ivector(1,n);
2728: weight=vector(1,n);
2729: for(i=1;i<=n;i++) weight[i]=1.0; /* Equal weights, 1 by default */
2730: mint=matrix(1,maxwav,1,n);
2731: anint=matrix(1,maxwav,1,n);
2732: s=imatrix(1,maxwav+1,1,n);
2733: adl=imatrix(1,maxwav+1,1,n);
2734: tab=ivector(1,NCOVMAX);
1.3 lievre 2735: ncodemax=ivector(1,8);
1.2 lievre 2736:
1.12 lievre 2737: i=1;
1.2 lievre 2738: while (fgets(line, MAXLINE, fic) != NULL) {
2739: if ((i >= firstobs) && (i <=lastobs)) {
2740:
2741: for (j=maxwav;j>=1;j--){
2742: cutv(stra, strb,line,' '); s[j][i]=atoi(strb);
2743: strcpy(line,stra);
2744: cutv(stra, strb,line,'/'); anint[j][i]=(double)(atoi(strb)); strcpy(line,stra);
2745: cutv(stra, strb,line,' '); mint[j][i]=(double)(atoi(strb)); strcpy(line,stra);
2746: }
2747:
2748: cutv(stra, strb,line,'/'); andc[i]=(double)(atoi(strb)); strcpy(line,stra);
2749: cutv(stra, strb,line,' '); moisdc[i]=(double)(atoi(strb)); strcpy(line,stra);
2750:
2751: cutv(stra, strb,line,'/'); annais[i]=(double)(atoi(strb)); strcpy(line,stra);
2752: cutv(stra, strb,line,' '); moisnais[i]=(double)(atoi(strb)); strcpy(line,stra);
2753:
2754: cutv(stra, strb,line,' '); weight[i]=(double)(atoi(strb)); strcpy(line,stra);
1.34 brouard 2755: for (j=ncovcol;j>=1;j--){
1.2 lievre 2756: cutv(stra, strb,line,' '); covar[j][i]=(double)(atoi(strb)); strcpy(line,stra);
2757: }
2758: num[i]=atol(stra);
1.12 lievre 2759:
2760: /*if((s[2][i]==2) && (s[3][i]==-1)&&(s[4][i]==9)){
2761: 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])); ij=ij+1;}*/
1.2 lievre 2762:
2763: i=i+1;
2764: }
2765: }
1.12 lievre 2766: /* printf("ii=%d", ij);
2767: scanf("%d",i);*/
2768: imx=i-1; /* Number of individuals */
1.3 lievre 2769:
1.12 lievre 2770: /* for (i=1; i<=imx; i++){
2771: if ((s[1][i]==3) && (s[2][i]==2)) s[2][i]=3;
2772: if ((s[2][i]==3) && (s[3][i]==2)) s[3][i]=3;
2773: if ((s[3][i]==3) && (s[4][i]==2)) s[4][i]=3;
1.35 lievre 2774: }*/
1.39 lievre 2775: /* for (i=1; i<=imx; i++){
1.35 lievre 2776: if (s[4][i]==9) s[4][i]=-1;
1.39 lievre 2777: 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]));}*/
2778:
1.35 lievre 2779:
1.2 lievre 2780: /* Calculation of the number of parameter from char model*/
1.7 lievre 2781: Tvar=ivector(1,15);
2782: Tprod=ivector(1,15);
2783: Tvaraff=ivector(1,15);
2784: Tvard=imatrix(1,15,1,2);
1.6 lievre 2785: Tage=ivector(1,15);
1.2 lievre 2786:
2787: if (strlen(model) >1){
1.7 lievre 2788: j=0, j1=0, k1=1, k2=1;
1.2 lievre 2789: j=nbocc(model,'+');
1.6 lievre 2790: j1=nbocc(model,'*');
1.2 lievre 2791: cptcovn=j+1;
1.7 lievre 2792: cptcovprod=j1;
1.3 lievre 2793:
1.2 lievre 2794: strcpy(modelsav,model);
1.8 lievre 2795: if ((strcmp(model,"age")==0) || (strcmp(model,"age*age")==0)){
2796: printf("Error. Non available option model=%s ",model);
2797: goto end;
2798: }
2799:
2800: for(i=(j+1); i>=1;i--){
2801: cutv(stra,strb,modelsav,'+');
2802: if (nbocc(modelsav,'+')==0) strcpy(strb,modelsav);
2803: /* printf("i=%d a=%s b=%s sav=%s\n",i, stra,strb,modelsav);*/
2804: /*scanf("%d",i);*/
2805: if (strchr(strb,'*')) {
2806: cutv(strd,strc,strb,'*');
2807: if (strcmp(strc,"age")==0) {
1.7 lievre 2808: cptcovprod--;
1.8 lievre 2809: cutv(strb,stre,strd,'V');
2810: Tvar[i]=atoi(stre);
2811: cptcovage++;
2812: Tage[cptcovage]=i;
2813: /*printf("stre=%s ", stre);*/
1.7 lievre 2814: }
1.8 lievre 2815: else if (strcmp(strd,"age")==0) {
1.7 lievre 2816: cptcovprod--;
1.8 lievre 2817: cutv(strb,stre,strc,'V');
2818: Tvar[i]=atoi(stre);
2819: cptcovage++;
2820: Tage[cptcovage]=i;
1.7 lievre 2821: }
2822: else {
1.8 lievre 2823: cutv(strb,stre,strc,'V');
1.34 brouard 2824: Tvar[i]=ncovcol+k1;
1.8 lievre 2825: cutv(strb,strc,strd,'V');
2826: Tprod[k1]=i;
2827: Tvard[k1][1]=atoi(strc);
2828: Tvard[k1][2]=atoi(stre);
2829: Tvar[cptcovn+k2]=Tvard[k1][1];
2830: Tvar[cptcovn+k2+1]=Tvard[k1][2];
1.7 lievre 2831: for (k=1; k<=lastobs;k++)
1.34 brouard 2832: covar[ncovcol+k1][k]=covar[atoi(stre)][k]*covar[atoi(strc)][k];
1.8 lievre 2833: k1++;
2834: k2=k2+2;
1.7 lievre 2835: }
1.2 lievre 2836: }
1.8 lievre 2837: else {
2838: /*printf("d=%s c=%s b=%s\n", strd,strc,strb);*/
2839: /* scanf("%d",i);*/
2840: cutv(strd,strc,strb,'V');
2841: Tvar[i]=atoi(strc);
2842: }
2843: strcpy(modelsav,stra);
2844: /*printf("a=%s b=%s sav=%s\n", stra,strb,modelsav);
2845: scanf("%d",i);*/
1.2 lievre 2846: }
1.8 lievre 2847: }
2848:
1.35 lievre 2849: /* printf("tvar1=%d tvar2=%d tvar3=%d cptcovage=%d Tage=%d",Tvar[1],Tvar[2],Tvar[3],cptcovage,Tage[1]);
1.8 lievre 2850: printf("cptcovprod=%d ", cptcovprod);
2851: scanf("%d ",i);*/
1.2 lievre 2852: fclose(fic);
2853:
1.7 lievre 2854: /* if(mle==1){*/
1.2 lievre 2855: if (weightopt != 1) { /* Maximisation without weights*/
2856: for(i=1;i<=n;i++) weight[i]=1.0;
2857: }
2858: /*-calculation of age at interview from date of interview and age at death -*/
2859: agev=matrix(1,maxwav,1,imx);
1.12 lievre 2860:
1.35 lievre 2861: for (i=1; i<=imx; i++) {
2862: for(m=2; (m<= maxwav); m++) {
1.12 lievre 2863: if ((mint[m][i]== 99) && (s[m][i] <= nlstate)){
2864: anint[m][i]=9999;
2865: s[m][i]=-1;
2866: }
1.35 lievre 2867: if(moisdc[i]==99 && andc[i]==9999 & s[m][i]>nlstate) s[m][i]=-1;
2868: }
2869: }
2870:
1.2 lievre 2871: for (i=1; i<=imx; i++) {
2872: agedc[i]=(moisdc[i]/12.+andc[i])-(moisnais[i]/12.+annais[i]);
2873: for(m=1; (m<= maxwav); m++){
2874: if(s[m][i] >0){
1.35 lievre 2875: if (s[m][i] >= nlstate+1) {
1.2 lievre 2876: if(agedc[i]>0)
2877: if(moisdc[i]!=99 && andc[i]!=9999)
1.35 lievre 2878: agev[m][i]=agedc[i];
2879: /*if(moisdc[i]==99 && andc[i]==9999) s[m][i]=-1;*/
2880: else {
1.8 lievre 2881: if (andc[i]!=9999){
1.2 lievre 2882: printf("Warning negative age at death: %d line:%d\n",num[i],i);
2883: agev[m][i]=-1;
1.8 lievre 2884: }
1.2 lievre 2885: }
2886: }
2887: else if(s[m][i] !=9){ /* Should no more exist */
2888: agev[m][i]=(mint[m][i]/12.+1./24.+anint[m][i])-(moisnais[i]/12.+1./24.+annais[i]);
1.3 lievre 2889: if(mint[m][i]==99 || anint[m][i]==9999)
1.2 lievre 2890: agev[m][i]=1;
2891: else if(agev[m][i] <agemin){
2892: agemin=agev[m][i];
2893: /*printf(" Min anint[%d][%d]=%.2f annais[%d]=%.2f, agemin=%.2f\n",m,i,anint[m][i], i,annais[i], agemin);*/
2894: }
2895: else if(agev[m][i] >agemax){
2896: agemax=agev[m][i];
2897: /* printf(" anint[%d][%d]=%.0f annais[%d]=%.0f, agemax=%.0f\n",m,i,anint[m][i], i,annais[i], agemax);*/
2898: }
2899: /*agev[m][i]=anint[m][i]-annais[i];*/
2900: /* agev[m][i] = age[i]+2*m;*/
2901: }
2902: else { /* =9 */
2903: agev[m][i]=1;
2904: s[m][i]=-1;
2905: }
2906: }
2907: else /*= 0 Unknown */
2908: agev[m][i]=1;
2909: }
2910:
2911: }
2912: for (i=1; i<=imx; i++) {
2913: for(m=1; (m<= maxwav); m++){
2914: if (s[m][i] > (nlstate+ndeath)) {
2915: printf("Error: Wrong value in nlstate or ndeath\n");
2916: goto end;
2917: }
2918: }
2919: }
2920:
2921: printf("Total number of individuals= %d, Agemin = %.2f, Agemax= %.2f\n\n", imx, agemin, agemax);
2922:
2923: free_vector(severity,1,maxwav);
2924: free_imatrix(outcome,1,maxwav+1,1,n);
2925: free_vector(moisnais,1,n);
2926: free_vector(annais,1,n);
1.17 lievre 2927: /* free_matrix(mint,1,maxwav,1,n);
2928: free_matrix(anint,1,maxwav,1,n);*/
1.2 lievre 2929: free_vector(moisdc,1,n);
2930: free_vector(andc,1,n);
2931:
2932:
2933: wav=ivector(1,imx);
2934: dh=imatrix(1,lastpass-firstpass+1,1,imx);
2935: mw=imatrix(1,lastpass-firstpass+1,1,imx);
2936:
2937: /* Concatenates waves */
2938: concatwav(wav, dh, mw, s, agedc, agev, firstpass, lastpass, imx, nlstate, stepm);
2939:
2940:
1.6 lievre 2941: Tcode=ivector(1,100);
1.8 lievre 2942: nbcode=imatrix(0,NCOVMAX,0,NCOVMAX);
1.7 lievre 2943: ncodemax[1]=1;
2944: if (cptcovn > 0) tricode(Tvar,nbcode,imx);
2945:
1.2 lievre 2946: codtab=imatrix(1,100,1,10);
2947: h=0;
1.7 lievre 2948: m=pow(2,cptcoveff);
1.2 lievre 2949:
1.7 lievre 2950: for(k=1;k<=cptcoveff; k++){
1.2 lievre 2951: for(i=1; i <=(m/pow(2,k));i++){
2952: for(j=1; j <= ncodemax[k]; j++){
1.7 lievre 2953: for(cpt=1; cpt <=(m/pow(2,cptcoveff+1-k)); cpt++){
1.2 lievre 2954: h++;
1.35 lievre 2955: if (h>m) h=1;codtab[h][k]=j;codtab[h][Tvar[k]]=j;
2956: /* printf("h=%d k=%d j=%d codtab[h][k]=%d tvar[k]=%d \n",h, k,j,codtab[h][k],Tvar[k]);*/
1.2 lievre 2957: }
2958: }
2959: }
2960: }
1.35 lievre 2961: /* printf("codtab[1][2]=%d codtab[2][2]=%d",codtab[1][2],codtab[2][2]);
2962: codtab[1][2]=1;codtab[2][2]=2; */
2963: /* for(i=1; i <=m ;i++){
2964: for(k=1; k <=cptcovn; k++){
2965: printf("i=%d k=%d %d %d ",i,k,codtab[i][k], cptcoveff);
2966: }
2967: printf("\n");
2968: }
2969: scanf("%d",i);*/
1.2 lievre 2970:
2971: /* Calculates basic frequencies. Computes observed prevalence at single age
2972: and prints on file fileres'p'. */
1.18 lievre 2973:
1.19 lievre 2974:
1.18 lievre 2975:
1.19 lievre 2976: pmmij= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
1.2 lievre 2977: oldms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2978: newms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2979: savms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
2980: oldm=oldms; newm=newms; savm=savms; /* Keeps fixed addresses to free */
1.12 lievre 2981:
1.2 lievre 2982: /* For Powell, parameters are in a vector p[] starting at p[1]
2983: so we point p on param[1][1] so that p[1] maps on param[1][1][1] */
2984: p=param[1][1]; /* *(*(*(param +1)+1)+0) */
1.7 lievre 2985:
2986: if(mle==1){
1.2 lievre 2987: mlikeli(ficres,p, npar, ncovmodel, nlstate, ftol, func);
1.7 lievre 2988: }
1.2 lievre 2989:
2990: /*--------- results files --------------*/
1.34 brouard 2991: fprintf(ficres,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle= 0 weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol, stepm, ncovcol, nlstate, ndeath, maxwav, weightopt,model);
1.19 lievre 2992:
1.16 lievre 2993:
1.2 lievre 2994: jk=1;
1.34 brouard 2995: fprintf(ficres,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
2996: printf("# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
1.2 lievre 2997: for(i=1,jk=1; i <=nlstate; i++){
2998: for(k=1; k <=(nlstate+ndeath); k++){
2999: if (k != i)
3000: {
3001: printf("%d%d ",i,k);
3002: fprintf(ficres,"%1d%1d ",i,k);
3003: for(j=1; j <=ncovmodel; j++){
3004: printf("%f ",p[jk]);
3005: fprintf(ficres,"%f ",p[jk]);
3006: jk++;
3007: }
3008: printf("\n");
3009: fprintf(ficres,"\n");
3010: }
3011: }
3012: }
1.7 lievre 3013: if(mle==1){
1.2 lievre 3014: /* Computing hessian and covariance matrix */
3015: ftolhess=ftol; /* Usually correct */
3016: hesscov(matcov, p, npar, delti, ftolhess, func);
1.7 lievre 3017: }
1.34 brouard 3018: fprintf(ficres,"# Scales (for hessian or gradient estimation)\n");
3019: printf("# Scales (for hessian or gradient estimation)\n");
1.2 lievre 3020: for(i=1,jk=1; i <=nlstate; i++){
3021: for(j=1; j <=nlstate+ndeath; j++){
3022: if (j!=i) {
3023: fprintf(ficres,"%1d%1d",i,j);
3024: printf("%1d%1d",i,j);
3025: for(k=1; k<=ncovmodel;k++){
3026: printf(" %.5e",delti[jk]);
3027: fprintf(ficres," %.5e",delti[jk]);
3028: jk++;
3029: }
3030: printf("\n");
3031: fprintf(ficres,"\n");
3032: }
3033: }
1.18 lievre 3034: }
1.2 lievre 3035:
3036: k=1;
1.34 brouard 3037: fprintf(ficres,"# Covariance matrix \n# 121 Var(a12)\n# 122 Cov(b12,a12) Var(b12)\n# ...\n# 232 Cov(b23,a12) Cov(b23,b12) ... Var (b23)\n");
3038: printf("# Covariance matrix \n# 121 Var(a12)\n# 122 Cov(b12,a12) Var(b12)\n# ...\n# 232 Cov(b23,a12) Cov(b23,b12) ... Var (b23)\n");
1.2 lievre 3039: for(i=1;i<=npar;i++){
3040: /* if (k>nlstate) k=1;
3041: i1=(i-1)/(ncovmodel*nlstate)+1;
3042: fprintf(ficres,"%s%d%d",alph[k],i1,tab[i]);
3043: printf("%s%d%d",alph[k],i1,tab[i]);*/
3044: fprintf(ficres,"%3d",i);
3045: printf("%3d",i);
3046: for(j=1; j<=i;j++){
3047: fprintf(ficres," %.5e",matcov[i][j]);
3048: printf(" %.5e",matcov[i][j]);
3049: }
3050: fprintf(ficres,"\n");
3051: printf("\n");
3052: k++;
3053: }
3054:
3055: while((c=getc(ficpar))=='#' && c!= EOF){
3056: ungetc(c,ficpar);
3057: fgets(line, MAXLINE, ficpar);
3058: puts(line);
3059: fputs(line,ficparo);
3060: }
3061: ungetc(c,ficpar);
1.36 brouard 3062: estepm=0;
3063: fscanf(ficpar,"agemin=%lf agemax=%lf bage=%lf fage=%lf estepm=%d\n",&ageminpar,&agemaxpar, &bage, &fage, &estepm);
3064: if (estepm==0 || estepm < stepm) estepm=stepm;
1.2 lievre 3065: if (fage <= 2) {
1.35 lievre 3066: bage = ageminpar;
1.28 lievre 3067: fage = agemaxpar;
1.2 lievre 3068: }
1.22 brouard 3069:
3070: fprintf(ficres,"# agemin agemax for life expectancy, bage fage (if mle==0 ie no data nor Max likelihood).\n");
1.36 brouard 3071: fprintf(ficres,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f estepm=%d\n",ageminpar,agemaxpar,bage,fage, estepm);
3072: fprintf(ficparo,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f estepm=%d\n",ageminpar,agemaxpar,bage,fage, estepm);
1.19 lievre 3073:
3074: while((c=getc(ficpar))=='#' && c!= EOF){
3075: ungetc(c,ficpar);
3076: fgets(line, MAXLINE, ficpar);
3077: puts(line);
3078: fputs(line,ficparo);
3079: }
3080: ungetc(c,ficpar);
3081:
1.25 lievre 3082: fscanf(ficpar,"begin-prev-date=%lf/%lf/%lf end-prev-date=%lf/%lf/%lf\n",&jprev1, &mprev1,&anprev1,&jprev2, &mprev2,&anprev2);
3083: fprintf(ficparo,"begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2);
3084: fprintf(ficres,"begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2);
1.19 lievre 3085:
3086: while((c=getc(ficpar))=='#' && c!= EOF){
3087: ungetc(c,ficpar);
3088: fgets(line, MAXLINE, ficpar);
3089: puts(line);
3090: fputs(line,ficparo);
3091: }
3092: ungetc(c,ficpar);
3093:
1.7 lievre 3094:
1.19 lievre 3095: dateprev1=anprev1+mprev1/12.+jprev1/365.;
3096: dateprev2=anprev2+mprev2/12.+jprev2/365.;
3097:
3098: fscanf(ficpar,"pop_based=%d\n",&popbased);
1.28 lievre 3099: fprintf(ficparo,"pop_based=%d\n",popbased);
3100: fprintf(ficres,"pop_based=%d\n",popbased);
3101:
3102: while((c=getc(ficpar))=='#' && c!= EOF){
3103: ungetc(c,ficpar);
3104: fgets(line, MAXLINE, ficpar);
3105: puts(line);
3106: fputs(line,ficparo);
3107: }
3108: ungetc(c,ficpar);
1.19 lievre 3109:
1.28 lievre 3110: fscanf(ficpar,"starting-proj-date=%lf/%lf/%lf final-proj-date=%lf/%lf/%lf mov_average=%d\n",&jproj1,&mproj1,&anproj1,&jproj2,&mproj2,&anproj2,&mobilav);
3111: fprintf(ficparo,"starting-proj-date=%.lf/%.lf/%.lf final-proj-date=%.lf/%.lf/%.lf mov_average=%d\n",jproj1,mproj1,anproj1,jproj2,mproj2,anproj2,mobilav);
3112: fprintf(ficres,"starting-proj-date=%.lf/%.lf/%.lf final-proj-date=%.lf/%.lf/%.lf mov_average=%d\n",jproj1,mproj1,anproj1,jproj2,mproj2,anproj2,mobilav);
3113:
3114:
3115: while((c=getc(ficpar))=='#' && c!= EOF){
1.19 lievre 3116: ungetc(c,ficpar);
3117: fgets(line, MAXLINE, ficpar);
3118: puts(line);
3119: fputs(line,ficparo);
3120: }
3121: ungetc(c,ficpar);
1.28 lievre 3122:
3123: fscanf(ficpar,"popforecast=%d popfile=%s popfiledate=%lf/%lf/%lf last-popfiledate=%lf/%lf/%lf\n",&popforecast,popfile,&jpyram,&mpyram,&anpyram,&jpyram1,&mpyram1,&anpyram1);
3124: fprintf(ficparo,"popforecast=%d popfile=%s popfiledate=%.lf/%.lf/%.lf last-popfiledate=%.lf/%.lf/%.lf\n",popforecast,popfile,jpyram,mpyram,anpyram,jpyram1,mpyram1,anpyram1);
3125: fprintf(ficres,"popforecast=%d popfile=%s popfiledate=%.lf/%.lf/%.lf last-popfiledate=%.lf/%.lf/%.lf\n",popforecast,popfile,jpyram,mpyram,anpyram,jpyram1,mpyram1,anpyram1);
1.19 lievre 3126:
1.26 lievre 3127: freqsummary(fileres, agemin, agemax, s, agev, nlstate, imx,Tvar,nbcode, ncodemax,mint,anint,dateprev1,dateprev2,jprev1, mprev1,anprev1,jprev2, mprev2,anprev2);
1.19 lievre 3128:
1.25 lievre 3129: /*------------ gnuplot -------------*/
1.35 lievre 3130: printinggnuplot(fileres,optionfilefiname,optionfile,optionfilegnuplot, ageminpar,agemaxpar,fage, pathc,p);
1.25 lievre 3131:
3132: /*------------ free_vector -------------*/
3133: chdir(path);
1.2 lievre 3134:
1.25 lievre 3135: free_ivector(wav,1,imx);
3136: free_imatrix(dh,1,lastpass-firstpass+1,1,imx);
3137: free_imatrix(mw,1,lastpass-firstpass+1,1,imx);
3138: free_ivector(num,1,n);
3139: free_vector(agedc,1,n);
3140: /*free_matrix(covar,1,NCOVMAX,1,n);*/
3141: fclose(ficparo);
3142: fclose(ficres);
1.28 lievre 3143:
1.2 lievre 3144: /*--------- index.htm --------*/
3145:
1.36 brouard 3146: printinghtml(fileres,title,datafile, firstpass, lastpass, stepm, weightopt,model,imx,jmin,jmax,jmean,optionfile,optionfilehtm,rfileres,optionfilegnuplot,version,popforecast,estepm);
1.2 lievre 3147:
1.25 lievre 3148:
1.2 lievre 3149: /*--------------- Prevalence limit --------------*/
3150:
3151: strcpy(filerespl,"pl");
3152: strcat(filerespl,fileres);
3153: if((ficrespl=fopen(filerespl,"w"))==NULL) {
3154: printf("Problem with Prev limit resultfile: %s\n", filerespl);goto end;
3155: }
3156: printf("Computing prevalence limit: result on file '%s' \n", filerespl);
3157: fprintf(ficrespl,"#Prevalence limit\n");
3158: fprintf(ficrespl,"#Age ");
3159: for(i=1; i<=nlstate;i++) fprintf(ficrespl,"%d-%d ",i,i);
3160: fprintf(ficrespl,"\n");
3161:
3162: prlim=matrix(1,nlstate,1,nlstate);
3163: pmmij= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
3164: oldms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
3165: newms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
3166: savms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
3167: oldm=oldms; newm=newms; savm=savms; /* Keeps fixed addresses to free */
3168: k=0;
1.35 lievre 3169: agebase=ageminpar;
1.28 lievre 3170: agelim=agemaxpar;
1.2 lievre 3171: ftolpl=1.e-10;
1.7 lievre 3172: i1=cptcoveff;
1.2 lievre 3173: if (cptcovn < 1){i1=1;}
3174:
3175: for(cptcov=1;cptcov<=i1;cptcov++){
3176: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
3177: k=k+1;
3178: /*printf("cptcov=%d cptcod=%d codtab=%d nbcode=%d\n",cptcov, cptcod,Tcode[cptcode],codtab[cptcod][cptcov]);*/
1.6 lievre 3179: fprintf(ficrespl,"\n#******");
1.7 lievre 3180: for(j=1;j<=cptcoveff;j++)
3181: fprintf(ficrespl," V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
1.2 lievre 3182: fprintf(ficrespl,"******\n");
3183:
3184: for (age=agebase; age<=agelim; age++){
3185: prevalim(prlim, nlstate, p, age, oldm, savm,ftolpl,k);
3186: fprintf(ficrespl,"%.0f",age );
3187: for(i=1; i<=nlstate;i++)
3188: fprintf(ficrespl," %.5f", prlim[i][i]);
3189: fprintf(ficrespl,"\n");
3190: }
3191: }
3192: }
3193: fclose(ficrespl);
1.13 lievre 3194:
1.2 lievre 3195: /*------------- h Pij x at various ages ------------*/
3196:
3197: strcpy(filerespij,"pij"); strcat(filerespij,fileres);
3198: if((ficrespij=fopen(filerespij,"w"))==NULL) {
3199: printf("Problem with Pij resultfile: %s\n", filerespij);goto end;
3200: }
3201: printf("Computing pij: result on file '%s' \n", filerespij);
3202:
3203: stepsize=(int) (stepm+YEARM-1)/YEARM;
1.13 lievre 3204: /*if (stepm<=24) stepsize=2;*/
1.2 lievre 3205:
3206: agelim=AGESUP;
3207: hstepm=stepsize*YEARM; /* Every year of age */
3208: hstepm=hstepm/stepm; /* Typically 2 years, = 2/6 months = 4 */
3209:
3210: k=0;
3211: for(cptcov=1;cptcov<=i1;cptcov++){
3212: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
3213: k=k+1;
3214: fprintf(ficrespij,"\n#****** ");
1.7 lievre 3215: for(j=1;j<=cptcoveff;j++)
3216: fprintf(ficrespij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
1.2 lievre 3217: fprintf(ficrespij,"******\n");
3218:
3219: for (agedeb=fage; agedeb>=bage; agedeb--){ /* If stepm=6 months */
3220: nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */
3221: nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
3222: p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
3223: oldm=oldms;savm=savms;
3224: hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);
3225: fprintf(ficrespij,"# Age");
3226: for(i=1; i<=nlstate;i++)
3227: for(j=1; j<=nlstate+ndeath;j++)
3228: fprintf(ficrespij," %1d-%1d",i,j);
3229: fprintf(ficrespij,"\n");
1.40 ! lievre 3230: for (h=0; h<=nhstepm; h++){
1.2 lievre 3231: fprintf(ficrespij,"%d %.0f %.0f",k,agedeb, agedeb+ h*hstepm/YEARM*stepm );
3232: for(i=1; i<=nlstate;i++)
3233: for(j=1; j<=nlstate+ndeath;j++)
3234: fprintf(ficrespij," %.5f", p3mat[i][j][h]);
3235: fprintf(ficrespij,"\n");
1.40 ! lievre 3236: }
1.2 lievre 3237: free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
3238: fprintf(ficrespij,"\n");
3239: }
3240: }
3241: }
3242:
1.39 lievre 3243: varprob(fileres, matcov, p, delti, nlstate, (int) bage, (int) fage,k,Tvar,nbcode, ncodemax);
1.13 lievre 3244:
1.2 lievre 3245: fclose(ficrespij);
3246:
1.27 lievre 3247:
3248: /*---------- Forecasting ------------------*/
1.32 brouard 3249: if((stepm == 1) && (strcmp(model,".")==0)){
1.27 lievre 3250: prevforecast(fileres, anproj1,mproj1,jproj1, agemin,agemax, dateprev1, dateprev2,mobilav, agedeb, fage, popforecast, popfile, anproj2,p, i1);
1.32 brouard 3251: if (popforecast==1) populforecast(fileres, anpyram,mpyram,jpyram, agemin,agemax, dateprev1, dateprev2,mobilav, agedeb, fage, popforecast, popfile, anpyram1,p, i1);
1.27 lievre 3252: free_matrix(mint,1,maxwav,1,n);
3253: free_matrix(anint,1,maxwav,1,n); free_imatrix(s,1,maxwav+1,1,n);
3254: free_vector(weight,1,n);}
1.21 lievre 3255: else{
3256: erreur=108;
1.32 brouard 3257: printf("Warning %d!! You can only forecast the prevalences if the optimization\n has been performed with stepm = 1 (month) instead of %d or model=. instead of '%s'\n", erreur, stepm, model);
1.21 lievre 3258: }
1.26 lievre 3259:
1.27 lievre 3260:
1.2 lievre 3261: /*---------- Health expectancies and variances ------------*/
3262:
3263: strcpy(filerest,"t");
3264: strcat(filerest,fileres);
3265: if((ficrest=fopen(filerest,"w"))==NULL) {
3266: printf("Problem with total LE resultfile: %s\n", filerest);goto end;
3267: }
3268: printf("Computing Total LEs with variances: file '%s' \n", filerest);
3269:
3270:
3271: strcpy(filerese,"e");
3272: strcat(filerese,fileres);
3273: if((ficreseij=fopen(filerese,"w"))==NULL) {
3274: printf("Problem with Health Exp. resultfile: %s\n", filerese); exit(0);
3275: }
3276: printf("Computing Health Expectancies: result on file '%s' \n", filerese);
3277:
3278: strcpy(fileresv,"v");
3279: strcat(fileresv,fileres);
3280: if((ficresvij=fopen(fileresv,"w"))==NULL) {
3281: printf("Problem with variance resultfile: %s\n", fileresv);exit(0);
3282: }
3283: printf("Computing Variance-covariance of DFLEs: file '%s' \n", fileresv);
3284:
3285: k=0;
3286: for(cptcov=1;cptcov<=i1;cptcov++){
3287: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
3288: k=k+1;
3289: fprintf(ficrest,"\n#****** ");
1.7 lievre 3290: for(j=1;j<=cptcoveff;j++)
3291: fprintf(ficrest,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
1.2 lievre 3292: fprintf(ficrest,"******\n");
3293:
3294: fprintf(ficreseij,"\n#****** ");
1.7 lievre 3295: for(j=1;j<=cptcoveff;j++)
1.35 lievre 3296: fprintf(ficreseij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
1.2 lievre 3297: fprintf(ficreseij,"******\n");
3298:
3299: fprintf(ficresvij,"\n#****** ");
1.7 lievre 3300: for(j=1;j<=cptcoveff;j++)
1.35 lievre 3301: fprintf(ficresvij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
1.2 lievre 3302: fprintf(ficresvij,"******\n");
3303:
3304: eij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
3305: oldm=oldms;savm=savms;
1.36 brouard 3306: evsij(fileres, eij, p, nlstate, stepm, (int) bage, (int)fage, oldm, savm, k, estepm);
1.2 lievre 3307: vareij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
3308: oldm=oldms;savm=savms;
1.36 brouard 3309: varevsij(fileres, vareij, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k, estepm);
1.26 lievre 3310:
3311:
3312:
1.2 lievre 3313: fprintf(ficrest,"#Total LEs with variances: e.. (std) ");
3314: for (i=1;i<=nlstate;i++) fprintf(ficrest,"e.%d (std) ",i);
3315: fprintf(ficrest,"\n");
1.26 lievre 3316:
1.2 lievre 3317: epj=vector(1,nlstate+1);
3318: for(age=bage; age <=fage ;age++){
3319: prevalim(prlim, nlstate, p, age, oldm, savm,ftolpl,k);
1.14 lievre 3320: if (popbased==1) {
3321: for(i=1; i<=nlstate;i++)
3322: prlim[i][i]=probs[(int)age][i][k];
3323: }
3324:
1.33 brouard 3325: fprintf(ficrest," %4.0f",age);
1.2 lievre 3326: for(j=1, epj[nlstate+1]=0.;j <=nlstate;j++){
3327: for(i=1, epj[j]=0.;i <=nlstate;i++) {
1.33 brouard 3328: epj[j] += prlim[i][i]*eij[i][j][(int)age];
1.2 lievre 3329: }
3330: epj[nlstate+1] +=epj[j];
3331: }
3332: for(i=1, vepp=0.;i <=nlstate;i++)
3333: for(j=1;j <=nlstate;j++)
3334: vepp += vareij[i][j][(int)age];
1.38 lievre 3335: fprintf(ficrest," %7.3f (%7.3f)", epj[nlstate+1],sqrt(vepp));
1.2 lievre 3336: for(j=1;j <=nlstate;j++){
1.38 lievre 3337: fprintf(ficrest," %7.3f (%7.3f)", epj[j],sqrt(vareij[j][j][(int)age]));
1.2 lievre 3338: }
3339: fprintf(ficrest,"\n");
3340: }
3341: }
3342: }
1.13 lievre 3343:
1.27 lievre 3344: fclose(ficreseij);
3345: fclose(ficresvij);
1.2 lievre 3346: fclose(ficrest);
3347: fclose(ficpar);
3348: free_vector(epj,1,nlstate+1);
1.26 lievre 3349:
1.2 lievre 3350: /*------- Variance limit prevalence------*/
3351:
1.27 lievre 3352: strcpy(fileresvpl,"vpl");
1.2 lievre 3353: strcat(fileresvpl,fileres);
3354: if((ficresvpl=fopen(fileresvpl,"w"))==NULL) {
3355: printf("Problem with variance prev lim resultfile: %s\n", fileresvpl);
3356: exit(0);
3357: }
3358: printf("Computing Variance-covariance of Prevalence limit: file '%s' \n", fileresvpl);
3359:
1.27 lievre 3360: k=0;
3361: for(cptcov=1;cptcov<=i1;cptcov++){
3362: for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){
3363: k=k+1;
3364: fprintf(ficresvpl,"\n#****** ");
3365: for(j=1;j<=cptcoveff;j++)
3366: fprintf(ficresvpl,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
3367: fprintf(ficresvpl,"******\n");
3368:
3369: varpl=matrix(1,nlstate,(int) bage, (int) fage);
3370: oldm=oldms;savm=savms;
1.2 lievre 3371: varprevlim(fileres, varpl, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k);
1.27 lievre 3372: }
1.2 lievre 3373: }
3374:
3375: fclose(ficresvpl);
3376:
3377: /*---------- End : free ----------------*/
3378: free_matrix(varpl,1,nlstate,(int) bage, (int)fage);
3379:
3380: free_ma3x(vareij,1,nlstate,1,nlstate,(int) bage, (int)fage);
3381: free_ma3x(eij,1,nlstate,1,nlstate,(int) bage, (int)fage);
3382:
3383:
3384: free_matrix(pmmij,1,nlstate+ndeath,1,nlstate+ndeath);
3385: free_matrix(oldms, 1,nlstate+ndeath,1,nlstate+ndeath);
3386: free_matrix(newms, 1,nlstate+ndeath,1,nlstate+ndeath);
3387: free_matrix(savms, 1,nlstate+ndeath,1,nlstate+ndeath);
1.13 lievre 3388:
1.2 lievre 3389: free_matrix(matcov,1,npar,1,npar);
3390: free_vector(delti,1,npar);
1.26 lievre 3391: free_matrix(agev,1,maxwav,1,imx);
1.2 lievre 3392: free_ma3x(param,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel);
3393:
1.21 lievre 3394: if(erreur >0)
1.34 brouard 3395: printf("End of Imach with error or warning %d\n",erreur);
1.21 lievre 3396: else printf("End of Imach\n");
1.2 lievre 3397: /* gettimeofday(&end_time, (struct timezone*)0);*/ /* after time */
3398:
3399: /* 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);*/
3400: /*printf("Total time was %d uSec.\n", total_usecs);*/
3401: /*------ End -----------*/
1.12 lievre 3402:
1.2 lievre 3403:
3404: end:
3405: #ifdef windows
1.22 brouard 3406: /* chdir(pathcd);*/
1.2 lievre 3407: #endif
1.22 brouard 3408: /*system("wgnuplot graph.plt");*/
3409: /*system("../gp37mgw/wgnuplot graph.plt");*/
3410: /*system("cd ../gp37mgw");*/
3411: /* system("..\\gp37mgw\\wgnuplot graph.plt");*/
3412: strcpy(plotcmd,GNUPLOTPROGRAM);
3413: strcat(plotcmd," ");
3414: strcat(plotcmd,optionfilegnuplot);
3415: system(plotcmd);
1.2 lievre 3416:
3417: #ifdef windows
3418: while (z[0] != 'q') {
1.35 lievre 3419: /* chdir(path); */
3420: printf("\nType e to edit output files, g to graph again, c to start again, and q for exiting: ");
1.2 lievre 3421: scanf("%s",z);
3422: if (z[0] == 'c') system("./imach");
1.35 lievre 3423: else if (z[0] == 'e') system(optionfilehtm);
3424: else if (z[0] == 'g') system(plotcmd);
1.2 lievre 3425: else if (z[0] == 'q') exit(0);
3426: }
3427: #endif
3428: }
3429:
3430:
FreeBSD-CVSweb <freebsd-cvsweb@FreeBSD.org>