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