File:  [Local Repository] / imach / src / imach.c
Revision 1.5: download - view: text, annotated - select for diffs
Wed May 2 17:42:45 2001 UTC (23 years ago) by lievre
Branches: MAIN
CVS tags: HEAD
suppressing cygwinpath. split directory instead.
version 15

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

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