File:  [Local Repository] / imach / src / imach.c
Revision 1.21: download - view: text, annotated - select for diffs
Thu Feb 21 18:42:24 2002 UTC (22 years, 3 months ago) by lievre
Branches: MAIN
CVS tags: HEAD
Error 108 if stepm is not equal to 1 (month) because forecasting
is only possible with a former optimization with a step of one month.

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

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