File:  [Local Repository] / imach / src / imach.c
Revision 1.22: download - view: text, annotated - select for diffs
Fri Feb 22 17:54:20 2002 UTC (22 years, 2 months ago) by brouard
Branches: MAIN
CVS tags: HEAD
First merge between Agnès and Nicolas!!

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

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