File:  [Local Repository] / imach / src / imach.c
Revision 1.8: download - view: text, annotated - select for diffs
Wed May 2 17:54:31 2001 UTC (23 years, 1 month ago) by lievre
Branches: MAIN
CVS tags: HEAD
Changes index.htm: title, model, ...

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

FreeBSD-CVSweb <freebsd-cvsweb@FreeBSD.org>