File:  [Local Repository] / imach / src / imach.c
Revision 1.3: download - view: text, annotated - select for diffs
Wed May 2 17:21:42 2001 UTC (23 years, 1 month ago) by lievre
Branches: MAIN
CVS tags: HEAD
version 13

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

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