File:  [Local Repository] / imach / src / imach.c
Revision 1.13: download - view: text, annotated - select for diffs
Wed Feb 20 17:02:08 2002 UTC (22 years, 3 months ago) by lievre
Branches: MAIN
CVS tags: HEAD
Variance of one-step transitions added and forecasting of prevalences

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

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