File:  [Local Repository] / imach / src / imach.c
Revision 1.15: download - view: text, annotated - select for diffs
Wed Feb 20 17:08:52 2002 UTC (22 years, 4 months ago) by lievre
Branches: MAIN
CVS tags: HEAD
Moving average  on cross-sectional prevalences

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

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