File:  [Local Repository] / imach / src / imach.c
Revision 1.165: download - view: text, annotated - select for diffs
Tue Dec 16 11:20:36 2014 UTC (9 years, 5 months ago) by brouard
Branches: MAIN
CVS tags: HEAD
Summary: After compiling on Visual C

* imach.c (Module): Merging 1.61 to 1.162

    1: /* $Id: imach.c,v 1.165 2014/12/16 11:20:36 brouard Exp $
    2:   $State: Exp $
    3:   $Log: imach.c,v $
    4:   Revision 1.165  2014/12/16 11:20:36  brouard
    5:   Summary: After compiling on Visual C
    6: 
    7:   * imach.c (Module): Merging 1.61 to 1.162
    8: 
    9:   Revision 1.164  2014/12/16 10:52:11  brouard
   10:   Summary: Merging with Visual C after suppressing some warnings for unused variables. Also fixing Saito's bug 0.98Xn
   11: 
   12:   * imach.c (Module): Merging 1.61 to 1.162
   13: 
   14:   Revision 1.163  2014/12/16 10:30:11  brouard
   15:   * imach.c (Module): Merging 1.61 to 1.162
   16: 
   17:   Revision 1.162  2014/09/25 11:43:39  brouard
   18:   Summary: temporary backup 0.99!
   19: 
   20:   Revision 1.1  2014/09/16 11:06:58  brouard
   21:   Summary: With some code (wrong) for nlopt
   22: 
   23:   Author:
   24: 
   25:   Revision 1.161  2014/09/15 20:41:41  brouard
   26:   Summary: Problem with macro SQR on Intel compiler
   27: 
   28:   Revision 1.160  2014/09/02 09:24:05  brouard
   29:   *** empty log message ***
   30: 
   31:   Revision 1.159  2014/09/01 10:34:10  brouard
   32:   Summary: WIN32
   33:   Author: Brouard
   34: 
   35:   Revision 1.158  2014/08/27 17:11:51  brouard
   36:   *** empty log message ***
   37: 
   38:   Revision 1.157  2014/08/27 16:26:55  brouard
   39:   Summary: Preparing windows Visual studio version
   40:   Author: Brouard
   41: 
   42:   In order to compile on Visual studio, time.h is now correct and time_t
   43:   and tm struct should be used. difftime should be used but sometimes I
   44:   just make the differences in raw time format (time(&now).
   45:   Trying to suppress #ifdef LINUX
   46:   Add xdg-open for __linux in order to open default browser.
   47: 
   48:   Revision 1.156  2014/08/25 20:10:10  brouard
   49:   *** empty log message ***
   50: 
   51:   Revision 1.155  2014/08/25 18:32:34  brouard
   52:   Summary: New compile, minor changes
   53:   Author: Brouard
   54: 
   55:   Revision 1.154  2014/06/20 17:32:08  brouard
   56:   Summary: Outputs now all graphs of convergence to period prevalence
   57: 
   58:   Revision 1.153  2014/06/20 16:45:46  brouard
   59:   Summary: If 3 live state, convergence to period prevalence on same graph
   60:   Author: Brouard
   61: 
   62:   Revision 1.152  2014/06/18 17:54:09  brouard
   63:   Summary: open browser, use gnuplot on same dir than imach if not found in the path
   64: 
   65:   Revision 1.151  2014/06/18 16:43:30  brouard
   66:   *** empty log message ***
   67: 
   68:   Revision 1.150  2014/06/18 16:42:35  brouard
   69:   Summary: If gnuplot is not in the path try on same directory than imach binary (OSX)
   70:   Author: brouard
   71: 
   72:   Revision 1.149  2014/06/18 15:51:14  brouard
   73:   Summary: Some fixes in parameter files errors
   74:   Author: Nicolas Brouard
   75: 
   76:   Revision 1.148  2014/06/17 17:38:48  brouard
   77:   Summary: Nothing new
   78:   Author: Brouard
   79: 
   80:   Just a new packaging for OS/X version 0.98nS
   81: 
   82:   Revision 1.147  2014/06/16 10:33:11  brouard
   83:   *** empty log message ***
   84: 
   85:   Revision 1.146  2014/06/16 10:20:28  brouard
   86:   Summary: Merge
   87:   Author: Brouard
   88: 
   89:   Merge, before building revised version.
   90: 
   91:   Revision 1.145  2014/06/10 21:23:15  brouard
   92:   Summary: Debugging with valgrind
   93:   Author: Nicolas Brouard
   94: 
   95:   Lot of changes in order to output the results with some covariates
   96:   After the Edimburgh REVES conference 2014, it seems mandatory to
   97:   improve the code.
   98:   No more memory valgrind error but a lot has to be done in order to
   99:   continue the work of splitting the code into subroutines.
  100:   Also, decodemodel has been improved. Tricode is still not
  101:   optimal. nbcode should be improved. Documentation has been added in
  102:   the source code.
  103: 
  104:   Revision 1.143  2014/01/26 09:45:38  brouard
  105:   Summary: Version 0.98nR (to be improved, but gives same optimization results as 0.98k. Nice, promising
  106: 
  107:   * imach.c (Module): Trying to merge old staffs together while being at Tokyo. Not tested...
  108:   (Module): Version 0.98nR Running ok, but output format still only works for three covariates.
  109: 
  110:   Revision 1.142  2014/01/26 03:57:36  brouard
  111:   Summary: gnuplot changed plot w l 1 has to be changed to plot w l lt 2
  112: 
  113:   * imach.c (Module): Trying to merge old staffs together while being at Tokyo. Not tested...
  114: 
  115:   Revision 1.141  2014/01/26 02:42:01  brouard
  116:   * imach.c (Module): Trying to merge old staffs together while being at Tokyo. Not tested...
  117: 
  118:   Revision 1.140  2011/09/02 10:37:54  brouard
  119:   Summary: times.h is ok with mingw32 now.
  120: 
  121:   Revision 1.139  2010/06/14 07:50:17  brouard
  122:   After the theft of my laptop, I probably lost some lines of codes which were not uploaded to the CVS tree.
  123:   I remember having already fixed agemin agemax which are pointers now but not cvs saved.
  124: 
  125:   Revision 1.138  2010/04/30 18:19:40  brouard
  126:   *** empty log message ***
  127: 
  128:   Revision 1.137  2010/04/29 18:11:38  brouard
  129:   (Module): Checking covariates for more complex models
  130:   than V1+V2. A lot of change to be done. Unstable.
  131: 
  132:   Revision 1.136  2010/04/26 20:30:53  brouard
  133:   (Module): merging some libgsl code. Fixing computation
  134:   of likelione (using inter/intrapolation if mle = 0) in order to
  135:   get same likelihood as if mle=1.
  136:   Some cleaning of code and comments added.
  137: 
  138:   Revision 1.135  2009/10/29 15:33:14  brouard
  139:   (Module): Now imach stops if date of birth, at least year of birth, is not given. Some cleaning of the code.
  140: 
  141:   Revision 1.134  2009/10/29 13:18:53  brouard
  142:   (Module): Now imach stops if date of birth, at least year of birth, is not given. Some cleaning of the code.
  143: 
  144:   Revision 1.133  2009/07/06 10:21:25  brouard
  145:   just nforces
  146: 
  147:   Revision 1.132  2009/07/06 08:22:05  brouard
  148:   Many tings
  149: 
  150:   Revision 1.131  2009/06/20 16:22:47  brouard
  151:   Some dimensions resccaled
  152: 
  153:   Revision 1.130  2009/05/26 06:44:34  brouard
  154:   (Module): Max Covariate is now set to 20 instead of 8. A
  155:   lot of cleaning with variables initialized to 0. Trying to make
  156:   V2+V3*age+V1+V4 strb=V3*age+V1+V4 working better.
  157: 
  158:   Revision 1.129  2007/08/31 13:49:27  lievre
  159:   Modification of the way of exiting when the covariate is not binary in order to see on the window the error message before exiting
  160: 
  161:   Revision 1.128  2006/06/30 13:02:05  brouard
  162:   (Module): Clarifications on computing e.j
  163: 
  164:   Revision 1.127  2006/04/28 18:11:50  brouard
  165:   (Module): Yes the sum of survivors was wrong since
  166:   imach-114 because nhstepm was no more computed in the age
  167:   loop. Now we define nhstepma in the age loop.
  168:   (Module): In order to speed up (in case of numerous covariates) we
  169:   compute health expectancies (without variances) in a first step
  170:   and then all the health expectancies with variances or standard
  171:   deviation (needs data from the Hessian matrices) which slows the
  172:   computation.
  173:   In the future we should be able to stop the program is only health
  174:   expectancies and graph are needed without standard deviations.
  175: 
  176:   Revision 1.126  2006/04/28 17:23:28  brouard
  177:   (Module): Yes the sum of survivors was wrong since
  178:   imach-114 because nhstepm was no more computed in the age
  179:   loop. Now we define nhstepma in the age loop.
  180:   Version 0.98h
  181: 
  182:   Revision 1.125  2006/04/04 15:20:31  lievre
  183:   Errors in calculation of health expectancies. Age was not initialized.
  184:   Forecasting file added.
  185: 
  186:   Revision 1.124  2006/03/22 17:13:53  lievre
  187:   Parameters are printed with %lf instead of %f (more numbers after the comma).
  188:   The log-likelihood is printed in the log file
  189: 
  190:   Revision 1.123  2006/03/20 10:52:43  brouard
  191:   * imach.c (Module): <title> changed, corresponds to .htm file
  192:   name. <head> headers where missing.
  193: 
  194:   * imach.c (Module): Weights can have a decimal point as for
  195:   English (a comma might work with a correct LC_NUMERIC environment,
  196:   otherwise the weight is truncated).
  197:   Modification of warning when the covariates values are not 0 or
  198:   1.
  199:   Version 0.98g
  200: 
  201:   Revision 1.122  2006/03/20 09:45:41  brouard
  202:   (Module): Weights can have a decimal point as for
  203:   English (a comma might work with a correct LC_NUMERIC environment,
  204:   otherwise the weight is truncated).
  205:   Modification of warning when the covariates values are not 0 or
  206:   1.
  207:   Version 0.98g
  208: 
  209:   Revision 1.121  2006/03/16 17:45:01  lievre
  210:   * imach.c (Module): Comments concerning covariates added
  211: 
  212:   * imach.c (Module): refinements in the computation of lli if
  213:   status=-2 in order to have more reliable computation if stepm is
  214:   not 1 month. Version 0.98f
  215: 
  216:   Revision 1.120  2006/03/16 15:10:38  lievre
  217:   (Module): refinements in the computation of lli if
  218:   status=-2 in order to have more reliable computation if stepm is
  219:   not 1 month. Version 0.98f
  220: 
  221:   Revision 1.119  2006/03/15 17:42:26  brouard
  222:   (Module): Bug if status = -2, the loglikelihood was
  223:   computed as likelihood omitting the logarithm. Version O.98e
  224: 
  225:   Revision 1.118  2006/03/14 18:20:07  brouard
  226:   (Module): varevsij Comments added explaining the second
  227:   table of variances if popbased=1 .
  228:   (Module): Covariances of eij, ekl added, graphs fixed, new html link.
  229:   (Module): Function pstamp added
  230:   (Module): Version 0.98d
  231: 
  232:   Revision 1.117  2006/03/14 17:16:22  brouard
  233:   (Module): varevsij Comments added explaining the second
  234:   table of variances if popbased=1 .
  235:   (Module): Covariances of eij, ekl added, graphs fixed, new html link.
  236:   (Module): Function pstamp added
  237:   (Module): Version 0.98d
  238: 
  239:   Revision 1.116  2006/03/06 10:29:27  brouard
  240:   (Module): Variance-covariance wrong links and
  241:   varian-covariance of ej. is needed (Saito).
  242: 
  243:   Revision 1.115  2006/02/27 12:17:45  brouard
  244:   (Module): One freematrix added in mlikeli! 0.98c
  245: 
  246:   Revision 1.114  2006/02/26 12:57:58  brouard
  247:   (Module): Some improvements in processing parameter
  248:   filename with strsep.
  249: 
  250:   Revision 1.113  2006/02/24 14:20:24  brouard
  251:   (Module): Memory leaks checks with valgrind and:
  252:   datafile was not closed, some imatrix were not freed and on matrix
  253:   allocation too.
  254: 
  255:   Revision 1.112  2006/01/30 09:55:26  brouard
  256:   (Module): Back to gnuplot.exe instead of wgnuplot.exe
  257: 
  258:   Revision 1.111  2006/01/25 20:38:18  brouard
  259:   (Module): Lots of cleaning and bugs added (Gompertz)
  260:   (Module): Comments can be added in data file. Missing date values
  261:   can be a simple dot '.'.
  262: 
  263:   Revision 1.110  2006/01/25 00:51:50  brouard
  264:   (Module): Lots of cleaning and bugs added (Gompertz)
  265: 
  266:   Revision 1.109  2006/01/24 19:37:15  brouard
  267:   (Module): Comments (lines starting with a #) are allowed in data.
  268: 
  269:   Revision 1.108  2006/01/19 18:05:42  lievre
  270:   Gnuplot problem appeared...
  271:   To be fixed
  272: 
  273:   Revision 1.107  2006/01/19 16:20:37  brouard
  274:   Test existence of gnuplot in imach path
  275: 
  276:   Revision 1.106  2006/01/19 13:24:36  brouard
  277:   Some cleaning and links added in html output
  278: 
  279:   Revision 1.105  2006/01/05 20:23:19  lievre
  280:   *** empty log message ***
  281: 
  282:   Revision 1.104  2005/09/30 16:11:43  lievre
  283:   (Module): sump fixed, loop imx fixed, and simplifications.
  284:   (Module): If the status is missing at the last wave but we know
  285:   that the person is alive, then we can code his/her status as -2
  286:   (instead of missing=-1 in earlier versions) and his/her
  287:   contributions to the likelihood is 1 - Prob of dying from last
  288:   health status (= 1-p13= p11+p12 in the easiest case of somebody in
  289:   the healthy state at last known wave). Version is 0.98
  290: 
  291:   Revision 1.103  2005/09/30 15:54:49  lievre
  292:   (Module): sump fixed, loop imx fixed, and simplifications.
  293: 
  294:   Revision 1.102  2004/09/15 17:31:30  brouard
  295:   Add the possibility to read data file including tab characters.
  296: 
  297:   Revision 1.101  2004/09/15 10:38:38  brouard
  298:   Fix on curr_time
  299: 
  300:   Revision 1.100  2004/07/12 18:29:06  brouard
  301:   Add version for Mac OS X. Just define UNIX in Makefile
  302: 
  303:   Revision 1.99  2004/06/05 08:57:40  brouard
  304:   *** empty log message ***
  305: 
  306:   Revision 1.98  2004/05/16 15:05:56  brouard
  307:   New version 0.97 . First attempt to estimate force of mortality
  308:   directly from the data i.e. without the need of knowing the health
  309:   state at each age, but using a Gompertz model: log u =a + b*age .
  310:   This is the basic analysis of mortality and should be done before any
  311:   other analysis, in order to test if the mortality estimated from the
  312:   cross-longitudinal survey is different from the mortality estimated
  313:   from other sources like vital statistic data.
  314: 
  315:   The same imach parameter file can be used but the option for mle should be -3.
  316: 
  317:   Agnès, who wrote this part of the code, tried to keep most of the
  318:   former routines in order to include the new code within the former code.
  319: 
  320:   The output is very simple: only an estimate of the intercept and of
  321:   the slope with 95% confident intervals.
  322: 
  323:   Current limitations:
  324:   A) Even if you enter covariates, i.e. with the
  325:   model= V1+V2 equation for example, the programm does only estimate a unique global model without covariates.
  326:   B) There is no computation of Life Expectancy nor Life Table.
  327: 
  328:   Revision 1.97  2004/02/20 13:25:42  lievre
  329:   Version 0.96d. Population forecasting command line is (temporarily)
  330:   suppressed.
  331: 
  332:   Revision 1.96  2003/07/15 15:38:55  brouard
  333:   * imach.c (Repository): Errors in subdirf, 2, 3 while printing tmpout is
  334:   rewritten within the same printf. Workaround: many printfs.
  335: 
  336:   Revision 1.95  2003/07/08 07:54:34  brouard
  337:   * imach.c (Repository):
  338:   (Repository): Using imachwizard code to output a more meaningful covariance
  339:   matrix (cov(a12,c31) instead of numbers.
  340: 
  341:   Revision 1.94  2003/06/27 13:00:02  brouard
  342:   Just cleaning
  343: 
  344:   Revision 1.93  2003/06/25 16:33:55  brouard
  345:   (Module): On windows (cygwin) function asctime_r doesn't
  346:   exist so I changed back to asctime which exists.
  347:   (Module): Version 0.96b
  348: 
  349:   Revision 1.92  2003/06/25 16:30:45  brouard
  350:   (Module): On windows (cygwin) function asctime_r doesn't
  351:   exist so I changed back to asctime which exists.
  352: 
  353:   Revision 1.91  2003/06/25 15:30:29  brouard
  354:   * imach.c (Repository): Duplicated warning errors corrected.
  355:   (Repository): Elapsed time after each iteration is now output. It
  356:   helps to forecast when convergence will be reached. Elapsed time
  357:   is stamped in powell.  We created a new html file for the graphs
  358:   concerning matrix of covariance. It has extension -cov.htm.
  359: 
  360:   Revision 1.90  2003/06/24 12:34:15  brouard
  361:   (Module): Some bugs corrected for windows. Also, when
  362:   mle=-1 a template is output in file "or"mypar.txt with the design
  363:   of the covariance matrix to be input.
  364: 
  365:   Revision 1.89  2003/06/24 12:30:52  brouard
  366:   (Module): Some bugs corrected for windows. Also, when
  367:   mle=-1 a template is output in file "or"mypar.txt with the design
  368:   of the covariance matrix to be input.
  369: 
  370:   Revision 1.88  2003/06/23 17:54:56  brouard
  371:   * imach.c (Repository): Create a sub-directory where all the secondary files are. Only imach, htm, gp and r(imach) are on the main directory. Correct time and other things.
  372: 
  373:   Revision 1.87  2003/06/18 12:26:01  brouard
  374:   Version 0.96
  375: 
  376:   Revision 1.86  2003/06/17 20:04:08  brouard
  377:   (Module): Change position of html and gnuplot routines and added
  378:   routine fileappend.
  379: 
  380:   Revision 1.85  2003/06/17 13:12:43  brouard
  381:   * imach.c (Repository): Check when date of death was earlier that
  382:   current date of interview. It may happen when the death was just
  383:   prior to the death. In this case, dh was negative and likelihood
  384:   was wrong (infinity). We still send an "Error" but patch by
  385:   assuming that the date of death was just one stepm after the
  386:   interview.
  387:   (Repository): Because some people have very long ID (first column)
  388:   we changed int to long in num[] and we added a new lvector for
  389:   memory allocation. But we also truncated to 8 characters (left
  390:   truncation)
  391:   (Repository): No more line truncation errors.
  392: 
  393:   Revision 1.84  2003/06/13 21:44:43  brouard
  394:   * imach.c (Repository): Replace "freqsummary" at a correct
  395:   place. It differs from routine "prevalence" which may be called
  396:   many times. Probs is memory consuming and must be used with
  397:   parcimony.
  398:   Version 0.95a3 (should output exactly the same maximization than 0.8a2)
  399: 
  400:   Revision 1.83  2003/06/10 13:39:11  lievre
  401:   *** empty log message ***
  402: 
  403:   Revision 1.82  2003/06/05 15:57:20  brouard
  404:   Add log in  imach.c and  fullversion number is now printed.
  405: 
  406: */
  407: /*
  408:    Interpolated Markov Chain
  409: 
  410:   Short summary of the programme:
  411:   
  412:   This program computes Healthy Life Expectancies from
  413:   cross-longitudinal data. Cross-longitudinal data consist in: -1- a
  414:   first survey ("cross") where individuals from different ages are
  415:   interviewed on their health status or degree of disability (in the
  416:   case of a health survey which is our main interest) -2- at least a
  417:   second wave of interviews ("longitudinal") which measure each change
  418:   (if any) in individual health status.  Health expectancies are
  419:   computed from the time spent in each health state according to a
  420:   model. More health states you consider, more time is necessary to reach the
  421:   Maximum Likelihood of the parameters involved in the model.  The
  422:   simplest model is the multinomial logistic model where pij is the
  423:   probability to be observed in state j at the second wave
  424:   conditional to be observed in state i at the first wave. Therefore
  425:   the model is: log(pij/pii)= aij + bij*age+ cij*sex + etc , where
  426:   'age' is age and 'sex' is a covariate. If you want to have a more
  427:   complex model than "constant and age", you should modify the program
  428:   where the markup *Covariates have to be included here again* invites
  429:   you to do it.  More covariates you add, slower the
  430:   convergence.
  431: 
  432:   The advantage of this computer programme, compared to a simple
  433:   multinomial logistic model, is clear when the delay between waves is not
  434:   identical for each individual. Also, if a individual missed an
  435:   intermediate interview, the information is lost, but taken into
  436:   account using an interpolation or extrapolation.  
  437: 
  438:   hPijx is the probability to be observed in state i at age x+h
  439:   conditional to the observed state i at age x. The delay 'h' can be
  440:   split into an exact number (nh*stepm) of unobserved intermediate
  441:   states. This elementary transition (by month, quarter,
  442:   semester or year) is modelled as a multinomial logistic.  The hPx
  443:   matrix is simply the matrix product of nh*stepm elementary matrices
  444:   and the contribution of each individual to the likelihood is simply
  445:   hPijx.
  446: 
  447:   Also this programme outputs the covariance matrix of the parameters but also
  448:   of the life expectancies. It also computes the period (stable) prevalence. 
  449:   
  450:   Authors: Nicolas Brouard (brouard@ined.fr) and Agnès Lièvre (lievre@ined.fr).
  451:            Institut national d'études démographiques, Paris.
  452:   This software have been partly granted by Euro-REVES, a concerted action
  453:   from the European Union.
  454:   It is copyrighted identically to a GNU software product, ie programme and
  455:   software can be distributed freely for non commercial use. Latest version
  456:   can be accessed at http://euroreves.ined.fr/imach .
  457: 
  458:   Help to debug: LD_PRELOAD=/usr/local/lib/libnjamd.so ./imach foo.imach
  459:   or better on gdb : set env LD_PRELOAD=/usr/local/lib/libnjamd.so
  460:   
  461:   **********************************************************************/
  462: /*
  463:   main
  464:   read parameterfile
  465:   read datafile
  466:   concatwav
  467:   freqsummary
  468:   if (mle >= 1)
  469:     mlikeli
  470:   print results files
  471:   if mle==1 
  472:      computes hessian
  473:   read end of parameter file: agemin, agemax, bage, fage, estepm
  474:       begin-prev-date,...
  475:   open gnuplot file
  476:   open html file
  477:   period (stable) prevalence      | pl_nom    1-1 2-2 etc by covariate
  478:    for age prevalim()             | #****** V1=0  V2=1  V3=1  V4=0 ******
  479:                                   | 65 1 0 2 1 3 1 4 0  0.96326 0.03674
  480:     freexexit2 possible for memory heap.
  481: 
  482:   h Pij x                         | pij_nom  ficrestpij
  483:    # Cov Agex agex+h hpijx with i,j= 1-1 1-2     1-3     2-1     2-2     2-3
  484:        1  85   85    1.00000             0.00000 0.00000 0.00000 1.00000 0.00000
  485:        1  85   86    0.68299             0.22291 0.09410 0.71093 0.00000 0.28907
  486: 
  487:        1  65   99    0.00364             0.00322 0.99314 0.00350 0.00310 0.99340
  488:        1  65  100    0.00214             0.00204 0.99581 0.00206 0.00196 0.99597
  489:   variance of p one-step probabilities varprob  | prob_nom   ficresprob #One-step probabilities and stand. devi in ()
  490:    Standard deviation of one-step probabilities | probcor_nom   ficresprobcor #One-step probabilities and correlation matrix
  491:    Matrix of variance covariance of one-step probabilities |  probcov_nom ficresprobcov #One-step probabilities and covariance matrix
  492: 
  493:   forecasting if prevfcast==1 prevforecast call prevalence()
  494:   health expectancies
  495:   Variance-covariance of DFLE
  496:   prevalence()
  497:    movingaverage()
  498:   varevsij() 
  499:   if popbased==1 varevsij(,popbased)
  500:   total life expectancies
  501:   Variance of period (stable) prevalence
  502:  end
  503: */
  504: 
  505: #define POWELL /* Instead of NLOPT */
  506: 
  507: #include <math.h>
  508: #include <stdio.h>
  509: #include <stdlib.h>
  510: #include <string.h>
  511: 
  512: #ifdef _WIN32
  513: #include <io.h>
  514: #else
  515: #include <unistd.h>
  516: #endif
  517: 
  518: #include <limits.h>
  519: #include <sys/types.h>
  520: #include <sys/stat.h>
  521: #include <errno.h>
  522: /* extern int errno; */
  523: 
  524: /* #ifdef LINUX */
  525: /* #include <time.h> */
  526: /* #include "timeval.h" */
  527: /* #else */
  528: /* #include <sys/time.h> */
  529: /* #endif */
  530: 
  531: #include <time.h>
  532: 
  533: #ifdef GSL
  534: #include <gsl/gsl_errno.h>
  535: #include <gsl/gsl_multimin.h>
  536: #endif
  537: 
  538: #ifdef NLOPT
  539: #include <nlopt.h>
  540: typedef struct {
  541:   double (* function)(double [] );
  542: } myfunc_data ;
  543: #endif
  544: 
  545: /* #include <libintl.h> */
  546: /* #define _(String) gettext (String) */
  547: 
  548: #define MAXLINE 1024 /* Was 256. Overflow with 312 with 2 states and 4 covariates. Should be ok */
  549: 
  550: #define GNUPLOTPROGRAM "gnuplot"
  551: /*#define GNUPLOTPROGRAM "..\\gp37mgw\\wgnuplot"*/
  552: #define FILENAMELENGTH 132
  553: 
  554: #define	GLOCK_ERROR_NOPATH		-1	/* empty path */
  555: #define	GLOCK_ERROR_GETCWD		-2	/* cannot get cwd */
  556: 
  557: #define MAXPARM 128 /**< Maximum number of parameters for the optimization */
  558: #define NPARMAX 64 /**< (nlstate+ndeath-1)*nlstate*ncovmodel */
  559: 
  560: #define NINTERVMAX 8
  561: #define NLSTATEMAX 8 /**< Maximum number of live states (for func) */
  562: #define NDEATHMAX 8 /**< Maximum number of dead states (for func) */
  563: #define NCOVMAX 20 /**< Maximum number of covariates, including generated covariates V1*V2 */
  564: #define codtabm(h,k)  1 & (h-1) >> (k-1) ;
  565: #define MAXN 20000
  566: #define YEARM 12. /**< Number of months per year */
  567: #define AGESUP 130
  568: #define AGEBASE 40
  569: #define AGEGOMP 10 /**< Minimal age for Gompertz adjustment */
  570: #ifdef _WIN32
  571: #define DIRSEPARATOR '\\'
  572: #define CHARSEPARATOR "\\"
  573: #define ODIRSEPARATOR '/'
  574: #else
  575: #define DIRSEPARATOR '/'
  576: #define CHARSEPARATOR "/"
  577: #define ODIRSEPARATOR '\\'
  578: #endif
  579: 
  580: /* $Id: imach.c,v 1.165 2014/12/16 11:20:36 brouard Exp $ */
  581: /* $State: Exp $ */
  582: 
  583: char version[]="Imach version 0.99, September 2014,INED-EUROREVES-Institut de longevite-Japan Society for the Promotion of Science (Grant-in-Aid for Scientific Research 25293121)";
  584: char fullversion[]="$Revision: 1.165 $ $Date: 2014/12/16 11:20:36 $"; 
  585: char strstart[80];
  586: char optionfilext[10], optionfilefiname[FILENAMELENGTH];
  587: int erreur=0, nberr=0, nbwarn=0; /* Error number, number of errors number of warnings  */
  588: int nvar=0, nforce=0; /* Number of variables, number of forces */
  589: /* Number of covariates model=V2+V1+ V3*age+V2*V4 */
  590: int cptcovn=0; /**< cptcovn number of covariates added in the model (excepting constant and age and age*product) */
  591: int cptcovt=0; /**< cptcovt number of covariates added in the model (excepting constant and age) */
  592: int cptcovs=0; /**< cptcovs number of simple covariates V2+V1 =2 */
  593: int cptcovage=0; /**< Number of covariates with age: V3*age only =1 */
  594: int cptcovprodnoage=0; /**< Number of covariate products without age */   
  595: int cptcoveff=0; /* Total number of covariates to vary for printing results */
  596: int cptcov=0; /* Working variable */
  597: int npar=NPARMAX;
  598: int nlstate=2; /* Number of live states */
  599: int ndeath=1; /* Number of dead states */
  600: int ncovmodel=0, ncovcol=0;     /* Total number of covariables including constant a12*1 +b12*x ncovmodel=2 */
  601: int popbased=0;
  602: 
  603: int *wav; /* Number of waves for this individuual 0 is possible */
  604: int maxwav=0; /* Maxim number of waves */
  605: int jmin=0, jmax=0; /* min, max spacing between 2 waves */
  606: int ijmin=0, ijmax=0; /* Individuals having jmin and jmax */ 
  607: int gipmx=0, gsw=0; /* Global variables on the number of contributions 
  608: 		   to the likelihood and the sum of weights (done by funcone)*/
  609: int mle=1, weightopt=0;
  610: int **mw; /* mw[mi][i] is number of the mi wave for this individual */
  611: int **dh; /* dh[mi][i] is number of steps between mi,mi+1 for this individual */
  612: int **bh; /* bh[mi][i] is the bias (+ or -) for this individual if the delay between
  613: 	   * wave mi and wave mi+1 is not an exact multiple of stepm. */
  614: int countcallfunc=0;  /* Count the number of calls to func */
  615: double jmean=1; /* Mean space between 2 waves */
  616: double **matprod2(); /* test */
  617: double **oldm, **newm, **savm; /* Working pointers to matrices */
  618: double **oldms, **newms, **savms; /* Fixed working pointers to matrices */
  619: /*FILE *fic ; */ /* Used in readdata only */
  620: FILE *ficpar, *ficparo,*ficres, *ficresp, *ficrespl, *ficrespij, *ficrest,*ficresf,*ficrespop;
  621: FILE *ficlog, *ficrespow;
  622: int globpr=0; /* Global variable for printing or not */
  623: double fretone; /* Only one call to likelihood */
  624: long ipmx=0; /* Number of contributions */
  625: double sw; /* Sum of weights */
  626: char filerespow[FILENAMELENGTH];
  627: char fileresilk[FILENAMELENGTH]; /* File of individual contributions to the likelihood */
  628: FILE *ficresilk;
  629: FILE *ficgp,*ficresprob,*ficpop, *ficresprobcov, *ficresprobcor;
  630: FILE *ficresprobmorprev;
  631: FILE *fichtm, *fichtmcov; /* Html File */
  632: FILE *ficreseij;
  633: char filerese[FILENAMELENGTH];
  634: FILE *ficresstdeij;
  635: char fileresstde[FILENAMELENGTH];
  636: FILE *ficrescveij;
  637: char filerescve[FILENAMELENGTH];
  638: FILE  *ficresvij;
  639: char fileresv[FILENAMELENGTH];
  640: FILE  *ficresvpl;
  641: char fileresvpl[FILENAMELENGTH];
  642: char title[MAXLINE];
  643: char optionfile[FILENAMELENGTH], datafile[FILENAMELENGTH],  filerespl[FILENAMELENGTH];
  644: char plotcmd[FILENAMELENGTH], pplotcmd[FILENAMELENGTH];
  645: char tmpout[FILENAMELENGTH],  tmpout2[FILENAMELENGTH]; 
  646: char command[FILENAMELENGTH];
  647: int  outcmd=0;
  648: 
  649: char fileres[FILENAMELENGTH], filerespij[FILENAMELENGTH], filereso[FILENAMELENGTH], rfileres[FILENAMELENGTH];
  650: 
  651: char filelog[FILENAMELENGTH]; /* Log file */
  652: char filerest[FILENAMELENGTH];
  653: char fileregp[FILENAMELENGTH];
  654: char popfile[FILENAMELENGTH];
  655: 
  656: char optionfilegnuplot[FILENAMELENGTH], optionfilehtm[FILENAMELENGTH], optionfilehtmcov[FILENAMELENGTH] ;
  657: 
  658: /* struct timeval start_time, end_time, curr_time, last_time, forecast_time; */
  659: /* struct timezone tzp; */
  660: /* extern int gettimeofday(); */
  661: struct tm tml, *gmtime(), *localtime();
  662: 
  663: extern time_t time();
  664: 
  665: struct tm start_time, end_time, curr_time, last_time, forecast_time;
  666: time_t  rstart_time, rend_time, rcurr_time, rlast_time, rforecast_time; /* raw time */
  667: struct tm tm;
  668: 
  669: char strcurr[80], strfor[80];
  670: 
  671: char *endptr;
  672: long lval;
  673: double dval;
  674: 
  675: #define NR_END 1
  676: #define FREE_ARG char*
  677: #define FTOL 1.0e-10
  678: 
  679: #define NRANSI 
  680: #define ITMAX 200 
  681: 
  682: #define TOL 2.0e-4 
  683: 
  684: #define CGOLD 0.3819660 
  685: #define ZEPS 1.0e-10 
  686: #define SHFT(a,b,c,d) (a)=(b);(b)=(c);(c)=(d); 
  687: 
  688: #define GOLD 1.618034 
  689: #define GLIMIT 100.0 
  690: #define TINY 1.0e-20 
  691: 
  692: static double maxarg1,maxarg2;
  693: #define FMAX(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)>(maxarg2)? (maxarg1):(maxarg2))
  694: #define FMIN(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)<(maxarg2)? (maxarg1):(maxarg2))
  695:   
  696: #define SIGN(a,b) ((b)>0.0 ? fabs(a) : -fabs(a))
  697: #define rint(a) floor(a+0.5)
  698: 
  699: static double sqrarg;
  700: #define SQR(a) ((sqrarg=(a)) == 0.0 ? 0.0 :sqrarg*sqrarg)
  701: #define SWAP(a,b) {temp=(a);(a)=(b);(b)=temp;} 
  702: int agegomp= AGEGOMP;
  703: 
  704: int imx; 
  705: int stepm=1;
  706: /* Stepm, step in month: minimum step interpolation*/
  707: 
  708: int estepm;
  709: /* Estepm, step in month to interpolate survival function in order to approximate Life Expectancy*/
  710: 
  711: int m,nb;
  712: long *num;
  713: int firstpass=0, lastpass=4,*cod, *ncodemax, *Tage,*cens;
  714: double **agev,*moisnais, *annais, *moisdc, *andc,**mint, **anint;
  715: double **pmmij, ***probs;
  716: double *ageexmed,*agecens;
  717: double dateintmean=0;
  718: 
  719: double *weight;
  720: int **s; /* Status */
  721: double *agedc;
  722: double  **covar; /**< covar[j,i], value of jth covariate for individual i,
  723: 		  * covar=matrix(0,NCOVMAX,1,n); 
  724: 		  * cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2]; */
  725: double  idx; 
  726: int **nbcode, *Tvar; /**< model=V2 => Tvar[1]= 2 */
  727: int *Ndum; /** Freq of modality (tricode */
  728: int **codtab; /**< codtab=imatrix(1,100,1,10); */
  729: int **Tvard, *Tprod, cptcovprod, *Tvaraff;
  730: double *lsurv, *lpop, *tpop;
  731: 
  732: double ftol=FTOL; /**< Tolerance for computing Max Likelihood */
  733: double ftolhess; /**< Tolerance for computing hessian */
  734: 
  735: /**************** split *************************/
  736: static	int split( char *path, char *dirc, char *name, char *ext, char *finame )
  737: {
  738:   /* From a file name with (full) path (either Unix or Windows) we extract the directory (dirc)
  739:      the name of the file (name), its extension only (ext) and its first part of the name (finame)
  740:   */ 
  741:   char	*ss;				/* pointer */
  742:   int	l1, l2;				/* length counters */
  743: 
  744:   l1 = strlen(path );			/* length of path */
  745:   if ( l1 == 0 ) return( GLOCK_ERROR_NOPATH );
  746:   ss= strrchr( path, DIRSEPARATOR );		/* find last / */
  747:   if ( ss == NULL ) {			/* no directory, so determine current directory */
  748:     strcpy( name, path );		/* we got the fullname name because no directory */
  749:     /*if(strrchr(path, ODIRSEPARATOR )==NULL)
  750:       printf("Warning you should use %s as a separator\n",DIRSEPARATOR);*/
  751:     /* get current working directory */
  752:     /*    extern  char* getcwd ( char *buf , int len);*/
  753:     if ( getcwd( dirc, FILENAME_MAX ) == NULL ) {
  754:       return( GLOCK_ERROR_GETCWD );
  755:     }
  756:     /* got dirc from getcwd*/
  757:     printf(" DIRC = %s \n",dirc);
  758:   } else {				/* strip direcotry from path */
  759:     ss++;				/* after this, the filename */
  760:     l2 = strlen( ss );			/* length of filename */
  761:     if ( l2 == 0 ) return( GLOCK_ERROR_NOPATH );
  762:     strcpy( name, ss );		/* save file name */
  763:     strncpy( dirc, path, l1 - l2 );	/* now the directory */
  764:     dirc[l1-l2] = 0;			/* add zero */
  765:     printf(" DIRC2 = %s \n",dirc);
  766:   }
  767:   /* We add a separator at the end of dirc if not exists */
  768:   l1 = strlen( dirc );			/* length of directory */
  769:   if( dirc[l1-1] != DIRSEPARATOR ){
  770:     dirc[l1] =  DIRSEPARATOR;
  771:     dirc[l1+1] = 0; 
  772:     printf(" DIRC3 = %s \n",dirc);
  773:   }
  774:   ss = strrchr( name, '.' );		/* find last / */
  775:   if (ss >0){
  776:     ss++;
  777:     strcpy(ext,ss);			/* save extension */
  778:     l1= strlen( name);
  779:     l2= strlen(ss)+1;
  780:     strncpy( finame, name, l1-l2);
  781:     finame[l1-l2]= 0;
  782:   }
  783: 
  784:   return( 0 );				/* we're done */
  785: }
  786: 
  787: 
  788: /******************************************/
  789: 
  790: void replace_back_to_slash(char *s, char*t)
  791: {
  792:   int i;
  793:   int lg=0;
  794:   i=0;
  795:   lg=strlen(t);
  796:   for(i=0; i<= lg; i++) {
  797:     (s[i] = t[i]);
  798:     if (t[i]== '\\') s[i]='/';
  799:   }
  800: }
  801: 
  802: char *trimbb(char *out, char *in)
  803: { /* Trim multiple blanks in line but keeps first blanks if line starts with blanks */
  804:   char *s;
  805:   s=out;
  806:   while (*in != '\0'){
  807:     while( *in == ' ' && *(in+1) == ' '){ /* && *(in+1) != '\0'){*/
  808:       in++;
  809:     }
  810:     *out++ = *in++;
  811:   }
  812:   *out='\0';
  813:   return s;
  814: }
  815: 
  816: char *cutl(char *blocc, char *alocc, char *in, char occ)
  817: {
  818:   /* cuts string in into blocc and alocc where blocc ends before first occurence of char 'occ' 
  819:      and alocc starts after first occurence of char 'occ' : ex cutv(blocc,alocc,"abcdef2ghi2j",'2')
  820:      gives blocc="abcdef2ghi" and alocc="j".
  821:      If occ is not found blocc is null and alocc is equal to in. Returns blocc
  822:   */
  823:   char *s, *t;
  824:   t=in;s=in;
  825:   while ((*in != occ) && (*in != '\0')){
  826:     *alocc++ = *in++;
  827:   }
  828:   if( *in == occ){
  829:     *(alocc)='\0';
  830:     s=++in;
  831:   }
  832:  
  833:   if (s == t) {/* occ not found */
  834:     *(alocc-(in-s))='\0';
  835:     in=s;
  836:   }
  837:   while ( *in != '\0'){
  838:     *blocc++ = *in++;
  839:   }
  840: 
  841:   *blocc='\0';
  842:   return t;
  843: }
  844: char *cutv(char *blocc, char *alocc, char *in, char occ)
  845: {
  846:   /* cuts string in into blocc and alocc where blocc ends before last occurence of char 'occ' 
  847:      and alocc starts after last occurence of char 'occ' : ex cutv(blocc,alocc,"abcdef2ghi2j",'2')
  848:      gives blocc="abcdef2ghi" and alocc="j".
  849:      If occ is not found blocc is null and alocc is equal to in. Returns alocc
  850:   */
  851:   char *s, *t;
  852:   t=in;s=in;
  853:   while (*in != '\0'){
  854:     while( *in == occ){
  855:       *blocc++ = *in++;
  856:       s=in;
  857:     }
  858:     *blocc++ = *in++;
  859:   }
  860:   if (s == t) /* occ not found */
  861:     *(blocc-(in-s))='\0';
  862:   else
  863:     *(blocc-(in-s)-1)='\0';
  864:   in=s;
  865:   while ( *in != '\0'){
  866:     *alocc++ = *in++;
  867:   }
  868: 
  869:   *alocc='\0';
  870:   return s;
  871: }
  872: 
  873: int nbocc(char *s, char occ)
  874: {
  875:   int i,j=0;
  876:   int lg=20;
  877:   i=0;
  878:   lg=strlen(s);
  879:   for(i=0; i<= lg; i++) {
  880:   if  (s[i] == occ ) j++;
  881:   }
  882:   return j;
  883: }
  884: 
  885: /* void cutv(char *u,char *v, char*t, char occ) */
  886: /* { */
  887: /*   /\* cuts string t into u and v where u ends before last occurence of char 'occ'  */
  888: /*      and v starts after last occurence of char 'occ' : ex cutv(u,v,"abcdef2ghi2j",'2') */
  889: /*      gives u="abcdef2ghi" and v="j" *\/ */
  890: /*   int i,lg,j,p=0; */
  891: /*   i=0; */
  892: /*   lg=strlen(t); */
  893: /*   for(j=0; j<=lg-1; j++) { */
  894: /*     if((t[j]!= occ) && (t[j+1]== occ)) p=j+1; */
  895: /*   } */
  896: 
  897: /*   for(j=0; j<p; j++) { */
  898: /*     (u[j] = t[j]); */
  899: /*   } */
  900: /*      u[p]='\0'; */
  901: 
  902: /*    for(j=0; j<= lg; j++) { */
  903: /*     if (j>=(p+1))(v[j-p-1] = t[j]); */
  904: /*   } */
  905: /* } */
  906: 
  907: #ifdef _WIN32
  908: char * strsep(char **pp, const char *delim)
  909: {
  910:   char *p, *q;
  911:          
  912:   if ((p = *pp) == NULL)
  913:     return 0;
  914:   if ((q = strpbrk (p, delim)) != NULL)
  915:   {
  916:     *pp = q + 1;
  917:     *q = '\0';
  918:   }
  919:   else
  920:     *pp = 0;
  921:   return p;
  922: }
  923: #endif
  924: 
  925: /********************** nrerror ********************/
  926: 
  927: void nrerror(char error_text[])
  928: {
  929:   fprintf(stderr,"ERREUR ...\n");
  930:   fprintf(stderr,"%s\n",error_text);
  931:   exit(EXIT_FAILURE);
  932: }
  933: /*********************** vector *******************/
  934: double *vector(int nl, int nh)
  935: {
  936:   double *v;
  937:   v=(double *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(double)));
  938:   if (!v) nrerror("allocation failure in vector");
  939:   return v-nl+NR_END;
  940: }
  941: 
  942: /************************ free vector ******************/
  943: void free_vector(double*v, int nl, int nh)
  944: {
  945:   free((FREE_ARG)(v+nl-NR_END));
  946: }
  947: 
  948: /************************ivector *******************************/
  949: int *ivector(long nl,long nh)
  950: {
  951:   int *v;
  952:   v=(int *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(int)));
  953:   if (!v) nrerror("allocation failure in ivector");
  954:   return v-nl+NR_END;
  955: }
  956: 
  957: /******************free ivector **************************/
  958: void free_ivector(int *v, long nl, long nh)
  959: {
  960:   free((FREE_ARG)(v+nl-NR_END));
  961: }
  962: 
  963: /************************lvector *******************************/
  964: long *lvector(long nl,long nh)
  965: {
  966:   long *v;
  967:   v=(long *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(long)));
  968:   if (!v) nrerror("allocation failure in ivector");
  969:   return v-nl+NR_END;
  970: }
  971: 
  972: /******************free lvector **************************/
  973: void free_lvector(long *v, long nl, long nh)
  974: {
  975:   free((FREE_ARG)(v+nl-NR_END));
  976: }
  977: 
  978: /******************* imatrix *******************************/
  979: int **imatrix(long nrl, long nrh, long ncl, long nch) 
  980:      /* allocate a int matrix with subscript range m[nrl..nrh][ncl..nch] */ 
  981: { 
  982:   long i, nrow=nrh-nrl+1,ncol=nch-ncl+1; 
  983:   int **m; 
  984:   
  985:   /* allocate pointers to rows */ 
  986:   m=(int **) malloc((size_t)((nrow+NR_END)*sizeof(int*))); 
  987:   if (!m) nrerror("allocation failure 1 in matrix()"); 
  988:   m += NR_END; 
  989:   m -= nrl; 
  990:   
  991:   
  992:   /* allocate rows and set pointers to them */ 
  993:   m[nrl]=(int *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(int))); 
  994:   if (!m[nrl]) nrerror("allocation failure 2 in matrix()"); 
  995:   m[nrl] += NR_END; 
  996:   m[nrl] -= ncl; 
  997:   
  998:   for(i=nrl+1;i<=nrh;i++) m[i]=m[i-1]+ncol; 
  999:   
 1000:   /* return pointer to array of pointers to rows */ 
 1001:   return m; 
 1002: } 
 1003: 
 1004: /****************** free_imatrix *************************/
 1005: void free_imatrix(m,nrl,nrh,ncl,nch)
 1006:       int **m;
 1007:       long nch,ncl,nrh,nrl; 
 1008:      /* free an int matrix allocated by imatrix() */ 
 1009: { 
 1010:   free((FREE_ARG) (m[nrl]+ncl-NR_END)); 
 1011:   free((FREE_ARG) (m+nrl-NR_END)); 
 1012: } 
 1013: 
 1014: /******************* matrix *******************************/
 1015: double **matrix(long nrl, long nrh, long ncl, long nch)
 1016: {
 1017:   long i, nrow=nrh-nrl+1, ncol=nch-ncl+1;
 1018:   double **m;
 1019: 
 1020:   m=(double **) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
 1021:   if (!m) nrerror("allocation failure 1 in matrix()");
 1022:   m += NR_END;
 1023:   m -= nrl;
 1024: 
 1025:   m[nrl]=(double *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
 1026:   if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
 1027:   m[nrl] += NR_END;
 1028:   m[nrl] -= ncl;
 1029: 
 1030:   for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
 1031:   return m;
 1032:   /* print *(*(m+1)+70) or print m[1][70]; print m+1 or print &(m[1]) or &(m[1][0])
 1033: m[i] = address of ith row of the table. &(m[i]) is its value which is another adress
 1034: that of m[i][0]. In order to get the value p m[i][0] but it is unitialized.
 1035:    */
 1036: }
 1037: 
 1038: /*************************free matrix ************************/
 1039: void free_matrix(double **m, long nrl, long nrh, long ncl, long nch)
 1040: {
 1041:   free((FREE_ARG)(m[nrl]+ncl-NR_END));
 1042:   free((FREE_ARG)(m+nrl-NR_END));
 1043: }
 1044: 
 1045: /******************* ma3x *******************************/
 1046: double ***ma3x(long nrl, long nrh, long ncl, long nch, long nll, long nlh)
 1047: {
 1048:   long i, j, nrow=nrh-nrl+1, ncol=nch-ncl+1, nlay=nlh-nll+1;
 1049:   double ***m;
 1050: 
 1051:   m=(double ***) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
 1052:   if (!m) nrerror("allocation failure 1 in matrix()");
 1053:   m += NR_END;
 1054:   m -= nrl;
 1055: 
 1056:   m[nrl]=(double **) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
 1057:   if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
 1058:   m[nrl] += NR_END;
 1059:   m[nrl] -= ncl;
 1060: 
 1061:   for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
 1062: 
 1063:   m[nrl][ncl]=(double *) malloc((size_t)((nrow*ncol*nlay+NR_END)*sizeof(double)));
 1064:   if (!m[nrl][ncl]) nrerror("allocation failure 3 in matrix()");
 1065:   m[nrl][ncl] += NR_END;
 1066:   m[nrl][ncl] -= nll;
 1067:   for (j=ncl+1; j<=nch; j++) 
 1068:     m[nrl][j]=m[nrl][j-1]+nlay;
 1069:   
 1070:   for (i=nrl+1; i<=nrh; i++) {
 1071:     m[i][ncl]=m[i-1l][ncl]+ncol*nlay;
 1072:     for (j=ncl+1; j<=nch; j++) 
 1073:       m[i][j]=m[i][j-1]+nlay;
 1074:   }
 1075:   return m; 
 1076:   /*  gdb: p *(m+1) <=> p m[1] and p (m+1) <=> p (m+1) <=> p &(m[1])
 1077:            &(m[i][j][k]) <=> *((*(m+i) + j)+k)
 1078:   */
 1079: }
 1080: 
 1081: /*************************free ma3x ************************/
 1082: void free_ma3x(double ***m, long nrl, long nrh, long ncl, long nch,long nll, long nlh)
 1083: {
 1084:   free((FREE_ARG)(m[nrl][ncl]+ nll-NR_END));
 1085:   free((FREE_ARG)(m[nrl]+ncl-NR_END));
 1086:   free((FREE_ARG)(m+nrl-NR_END));
 1087: }
 1088: 
 1089: /*************** function subdirf ***********/
 1090: char *subdirf(char fileres[])
 1091: {
 1092:   /* Caution optionfilefiname is hidden */
 1093:   strcpy(tmpout,optionfilefiname);
 1094:   strcat(tmpout,"/"); /* Add to the right */
 1095:   strcat(tmpout,fileres);
 1096:   return tmpout;
 1097: }
 1098: 
 1099: /*************** function subdirf2 ***********/
 1100: char *subdirf2(char fileres[], char *preop)
 1101: {
 1102:   
 1103:   /* Caution optionfilefiname is hidden */
 1104:   strcpy(tmpout,optionfilefiname);
 1105:   strcat(tmpout,"/");
 1106:   strcat(tmpout,preop);
 1107:   strcat(tmpout,fileres);
 1108:   return tmpout;
 1109: }
 1110: 
 1111: /*************** function subdirf3 ***********/
 1112: char *subdirf3(char fileres[], char *preop, char *preop2)
 1113: {
 1114:   
 1115:   /* Caution optionfilefiname is hidden */
 1116:   strcpy(tmpout,optionfilefiname);
 1117:   strcat(tmpout,"/");
 1118:   strcat(tmpout,preop);
 1119:   strcat(tmpout,preop2);
 1120:   strcat(tmpout,fileres);
 1121:   return tmpout;
 1122: }
 1123: 
 1124: char *asc_diff_time(long time_sec, char ascdiff[])
 1125: {
 1126:   long sec_left, days, hours, minutes;
 1127:   days = (time_sec) / (60*60*24);
 1128:   sec_left = (time_sec) % (60*60*24);
 1129:   hours = (sec_left) / (60*60) ;
 1130:   sec_left = (sec_left) %(60*60);
 1131:   minutes = (sec_left) /60;
 1132:   sec_left = (sec_left) % (60);
 1133:   sprintf(ascdiff,"%ld day(s) %ld hour(s) %ld minute(s) %ld second(s)",days, hours, minutes, sec_left);  
 1134:   return ascdiff;
 1135: }
 1136: 
 1137: /***************** f1dim *************************/
 1138: extern int ncom; 
 1139: extern double *pcom,*xicom;
 1140: extern double (*nrfunc)(double []); 
 1141:  
 1142: double f1dim(double x) 
 1143: { 
 1144:   int j; 
 1145:   double f;
 1146:   double *xt; 
 1147:  
 1148:   xt=vector(1,ncom); 
 1149:   for (j=1;j<=ncom;j++) xt[j]=pcom[j]+x*xicom[j]; 
 1150:   f=(*nrfunc)(xt); 
 1151:   free_vector(xt,1,ncom); 
 1152:   return f; 
 1153: } 
 1154: 
 1155: /*****************brent *************************/
 1156: double brent(double ax, double bx, double cx, double (*f)(double), double tol, 	double *xmin) 
 1157: { 
 1158:   int iter; 
 1159:   double a,b,d,etemp;
 1160:   double fu=0,fv,fw,fx;
 1161:   double ftemp=0.;
 1162:   double p,q,r,tol1,tol2,u,v,w,x,xm; 
 1163:   double e=0.0; 
 1164:  
 1165:   a=(ax < cx ? ax : cx); 
 1166:   b=(ax > cx ? ax : cx); 
 1167:   x=w=v=bx; 
 1168:   fw=fv=fx=(*f)(x); 
 1169:   for (iter=1;iter<=ITMAX;iter++) { 
 1170:     xm=0.5*(a+b); 
 1171:     tol2=2.0*(tol1=tol*fabs(x)+ZEPS); 
 1172:     /*		if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret)))*/
 1173:     printf(".");fflush(stdout);
 1174:     fprintf(ficlog,".");fflush(ficlog);
 1175: #ifdef DEBUGBRENT
 1176:     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);
 1177:     fprintf(ficlog,"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);
 1178:     /*		if ((fabs(x-xm) <= (tol2-0.5*(b-a)))||(2.0*fabs(fu-ftemp) <= ftol*1.e-2*(fabs(fu)+fabs(ftemp)))) { */
 1179: #endif
 1180:     if (fabs(x-xm) <= (tol2-0.5*(b-a))){ 
 1181:       *xmin=x; 
 1182:       return fx; 
 1183:     } 
 1184:     ftemp=fu;
 1185:     if (fabs(e) > tol1) { 
 1186:       r=(x-w)*(fx-fv); 
 1187:       q=(x-v)*(fx-fw); 
 1188:       p=(x-v)*q-(x-w)*r; 
 1189:       q=2.0*(q-r); 
 1190:       if (q > 0.0) p = -p; 
 1191:       q=fabs(q); 
 1192:       etemp=e; 
 1193:       e=d; 
 1194:       if (fabs(p) >= fabs(0.5*q*etemp) || p <= q*(a-x) || p >= q*(b-x)) 
 1195: 	d=CGOLD*(e=(x >= xm ? a-x : b-x)); 
 1196:       else { 
 1197: 	d=p/q; 
 1198: 	u=x+d; 
 1199: 	if (u-a < tol2 || b-u < tol2) 
 1200: 	  d=SIGN(tol1,xm-x); 
 1201:       } 
 1202:     } else { 
 1203:       d=CGOLD*(e=(x >= xm ? a-x : b-x)); 
 1204:     } 
 1205:     u=(fabs(d) >= tol1 ? x+d : x+SIGN(tol1,d)); 
 1206:     fu=(*f)(u); 
 1207:     if (fu <= fx) { 
 1208:       if (u >= x) a=x; else b=x; 
 1209:       SHFT(v,w,x,u) 
 1210: 	SHFT(fv,fw,fx,fu) 
 1211: 	} else { 
 1212: 	  if (u < x) a=u; else b=u; 
 1213: 	  if (fu <= fw || w == x) { 
 1214: 	    v=w; 
 1215: 	    w=u; 
 1216: 	    fv=fw; 
 1217: 	    fw=fu; 
 1218: 	  } else if (fu <= fv || v == x || v == w) { 
 1219: 	    v=u; 
 1220: 	    fv=fu; 
 1221: 	  } 
 1222: 	} 
 1223:   } 
 1224:   nrerror("Too many iterations in brent"); 
 1225:   *xmin=x; 
 1226:   return fx; 
 1227: } 
 1228: 
 1229: /****************** mnbrak ***********************/
 1230: 
 1231: void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb, double *fc, 
 1232: 	    double (*func)(double)) 
 1233: { 
 1234:   double ulim,u,r,q, dum;
 1235:   double fu; 
 1236:  
 1237:   *fa=(*func)(*ax); 
 1238:   *fb=(*func)(*bx); 
 1239:   if (*fb > *fa) { 
 1240:     SHFT(dum,*ax,*bx,dum) 
 1241:       SHFT(dum,*fb,*fa,dum) 
 1242:       } 
 1243:   *cx=(*bx)+GOLD*(*bx-*ax); 
 1244:   *fc=(*func)(*cx); 
 1245:   while (*fb > *fc) { /* Declining fa, fb, fc */
 1246:     r=(*bx-*ax)*(*fb-*fc); 
 1247:     q=(*bx-*cx)*(*fb-*fa); 
 1248:     u=(*bx)-((*bx-*cx)*q-(*bx-*ax)*r)/ 
 1249:       (2.0*SIGN(FMAX(fabs(q-r),TINY),q-r)); /* Minimum abscisse of a parabolic estimated from (a,fa), (b,fb) and (c,fc). */
 1250:     ulim=(*bx)+GLIMIT*(*cx-*bx); /* Maximum abscisse where function can be evaluated */
 1251:     if ((*bx-u)*(u-*cx) > 0.0) { /* if u between b and c */
 1252:       fu=(*func)(u); 
 1253: #ifdef DEBUG
 1254:       /* f(x)=A(x-u)**2+f(u) */
 1255:       double A, fparabu; 
 1256:       A= (*fb - *fa)/(*bx-*ax)/(*bx+*ax-2*u);
 1257:       fparabu= *fa - A*(*ax-u)*(*ax-u);
 1258:       printf("mnbrak (*ax=%.12f, *fa=%.12lf), (*bx=%.12f, *fb=%.12lf), (*cx=%.12f, *fc=%.12lf),  (*u=%.12f, fu=%.12lf, fparabu=%.12f)\n",*ax,*fa,*bx,*fb,*cx,*fc,u,fu, fparabu);
 1259:       fprintf(ficlog, "mnbrak (*ax=%.12f, *fa=%.12lf), (*bx=%.12f, *fb=%.12lf), (*cx=%.12f, *fc=%.12lf),  (*u=%.12f, fu=%.12lf, fparabu=%.12f)\n",*ax,*fa,*bx,*fb,*cx,*fc,u,fu, fparabu);
 1260: #endif 
 1261:     } else if ((*cx-u)*(u-ulim) > 0.0) { /* u is after c but before ulim */
 1262:       fu=(*func)(u); 
 1263:       if (fu < *fc) { 
 1264: 	SHFT(*bx,*cx,u,*cx+GOLD*(*cx-*bx)) 
 1265: 	  SHFT(*fb,*fc,fu,(*func)(u)) 
 1266: 	  } 
 1267:     } else if ((u-ulim)*(ulim-*cx) >= 0.0) { /* u outside ulim (verifying that ulim is beyond c) */
 1268:       u=ulim; 
 1269:       fu=(*func)(u); 
 1270:     } else { 
 1271:       u=(*cx)+GOLD*(*cx-*bx); 
 1272:       fu=(*func)(u); 
 1273:     } 
 1274:     SHFT(*ax,*bx,*cx,u) 
 1275:       SHFT(*fa,*fb,*fc,fu) 
 1276:       } 
 1277: } 
 1278: 
 1279: /*************** linmin ************************/
 1280: /* Given an n -dimensional point p[1..n] and an n -dimensional direction xi[1..n] , moves and
 1281: resets p to where the function func(p) takes on a minimum along the direction xi from p ,
 1282: and replaces xi by the actual vector displacement that p was moved. Also returns as fret
 1283: the value of func at the returned location p . This is actually all accomplished by calling the
 1284: routines mnbrak and brent .*/
 1285: int ncom; 
 1286: double *pcom,*xicom;
 1287: double (*nrfunc)(double []); 
 1288:  
 1289: void linmin(double p[], double xi[], int n, double *fret,double (*func)(double [])) 
 1290: { 
 1291:   double brent(double ax, double bx, double cx, 
 1292: 	       double (*f)(double), double tol, double *xmin); 
 1293:   double f1dim(double x); 
 1294:   void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb, 
 1295: 	      double *fc, double (*func)(double)); 
 1296:   int j; 
 1297:   double xx,xmin,bx,ax; 
 1298:   double fx,fb,fa;
 1299:  
 1300:   ncom=n; 
 1301:   pcom=vector(1,n); 
 1302:   xicom=vector(1,n); 
 1303:   nrfunc=func; 
 1304:   for (j=1;j<=n;j++) { 
 1305:     pcom[j]=p[j]; 
 1306:     xicom[j]=xi[j]; 
 1307:   } 
 1308:   ax=0.0; 
 1309:   xx=1.0; 
 1310:   mnbrak(&ax,&xx,&bx,&fa,&fx,&fb,f1dim); /* Find a bracket a,x,b in direction n=xi ie xicom */
 1311:   *fret=brent(ax,xx,bx,f1dim,TOL,&xmin); /* Find a minimum P+lambda n in that direction (lambdamin), with TOL between abscisses */
 1312: #ifdef DEBUG
 1313:   printf("retour brent fret=%.12e xmin=%.12e\n",*fret,xmin);
 1314:   fprintf(ficlog,"retour brent fret=%.12e xmin=%.12e\n",*fret,xmin);
 1315: #endif
 1316:   for (j=1;j<=n;j++) { 
 1317:     xi[j] *= xmin; 
 1318:     p[j] += xi[j]; 
 1319:   } 
 1320:   free_vector(xicom,1,n); 
 1321:   free_vector(pcom,1,n); 
 1322: } 
 1323: 
 1324: 
 1325: /*************** powell ************************/
 1326: /*
 1327: Minimization of a function func of n variables. Input consists of an initial starting point
 1328: p[1..n] ; an initial matrix xi[1..n][1..n] , whose columns contain the initial set of di-
 1329: rections (usually the n unit vectors); and ftol , the fractional tolerance in the function value
 1330: such that failure to decrease by more than this amount on one iteration signals doneness. On
 1331: output, p is set to the best point found, xi is the then-current direction set, fret is the returned
 1332: function value at p , and iter is the number of iterations taken. The routine linmin is used.
 1333:  */
 1334: void powell(double p[], double **xi, int n, double ftol, int *iter, double *fret, 
 1335: 	    double (*func)(double [])) 
 1336: { 
 1337:   void linmin(double p[], double xi[], int n, double *fret, 
 1338: 	      double (*func)(double [])); 
 1339:   int i,ibig,j; 
 1340:   double del,t,*pt,*ptt,*xit;
 1341:   double fp,fptt;
 1342:   double *xits;
 1343:   int niterf, itmp;
 1344: 
 1345:   pt=vector(1,n); 
 1346:   ptt=vector(1,n); 
 1347:   xit=vector(1,n); 
 1348:   xits=vector(1,n); 
 1349:   *fret=(*func)(p); 
 1350:   for (j=1;j<=n;j++) pt[j]=p[j]; 
 1351:     rcurr_time = time(NULL);  
 1352:   for (*iter=1;;++(*iter)) { 
 1353:     fp=(*fret); 
 1354:     ibig=0; 
 1355:     del=0.0; 
 1356:     rlast_time=rcurr_time;
 1357:     /* (void) gettimeofday(&curr_time,&tzp); */
 1358:     rcurr_time = time(NULL);  
 1359:     curr_time = *localtime(&rcurr_time);
 1360:     printf("\nPowell iter=%d -2*LL=%.12f %ld sec. %ld sec.",*iter,*fret, rcurr_time-rlast_time, rcurr_time-rstart_time);fflush(stdout);
 1361:     fprintf(ficlog,"\nPowell iter=%d -2*LL=%.12f %ld sec. %ld sec.",*iter,*fret,rcurr_time-rlast_time, rcurr_time-rstart_time); fflush(ficlog);
 1362: /*     fprintf(ficrespow,"%d %.12f %ld",*iter,*fret,curr_time.tm_sec-start_time.tm_sec); */
 1363:    for (i=1;i<=n;i++) {
 1364:       printf(" %d %.12f",i, p[i]);
 1365:       fprintf(ficlog," %d %.12lf",i, p[i]);
 1366:       fprintf(ficrespow," %.12lf", p[i]);
 1367:     }
 1368:     printf("\n");
 1369:     fprintf(ficlog,"\n");
 1370:     fprintf(ficrespow,"\n");fflush(ficrespow);
 1371:     if(*iter <=3){
 1372:       tml = *localtime(&rcurr_time);
 1373:       strcpy(strcurr,asctime(&tml));
 1374:       rforecast_time=rcurr_time; 
 1375:       itmp = strlen(strcurr);
 1376:       if(strcurr[itmp-1]=='\n')  /* Windows outputs with a new line */
 1377: 	strcurr[itmp-1]='\0';
 1378:       printf("\nConsidering the time needed for the last iteration #%d: %ld seconds,\n",*iter,rcurr_time-rlast_time);
 1379:       fprintf(ficlog,"\nConsidering the time needed for this last iteration #%d: %ld seconds,\n",*iter,rcurr_time-rlast_time);
 1380:       for(niterf=10;niterf<=30;niterf+=10){
 1381: 	rforecast_time=rcurr_time+(niterf-*iter)*(rcurr_time-rlast_time);
 1382: 	forecast_time = *localtime(&rforecast_time);
 1383: 	strcpy(strfor,asctime(&forecast_time));
 1384: 	itmp = strlen(strfor);
 1385: 	if(strfor[itmp-1]=='\n')
 1386: 	strfor[itmp-1]='\0';
 1387: 	printf("   - if your program needs %d iterations to converge, convergence will be \n   reached in %s i.e.\n   on %s (current time is %s);\n",niterf, asc_diff_time(rforecast_time-rcurr_time,tmpout),strfor,strcurr);
 1388: 	fprintf(ficlog,"   - if your program needs %d iterations to converge, convergence will be \n   reached in %s i.e.\n   on %s (current time is %s);\n",niterf, asc_diff_time(rforecast_time-rcurr_time,tmpout),strfor,strcurr);
 1389:       }
 1390:     }
 1391:     for (i=1;i<=n;i++) { 
 1392:       for (j=1;j<=n;j++) xit[j]=xi[j][i]; 
 1393:       fptt=(*fret); 
 1394: #ifdef DEBUG
 1395: 	  printf("fret=%lf, %lf, %lf \n", *fret, *fret, *fret);
 1396: 	  fprintf(ficlog, "fret=%lf, %lf, %lf \n", *fret, *fret, *fret);
 1397: #endif
 1398:       printf("%d",i);fflush(stdout);
 1399:       fprintf(ficlog,"%d",i);fflush(ficlog);
 1400:       linmin(p,xit,n,fret,func); 
 1401:       if (fabs(fptt-(*fret)) > del) { 
 1402: 	del=fabs(fptt-(*fret)); 
 1403: 	ibig=i; 
 1404:       } 
 1405: #ifdef DEBUG
 1406:       printf("%d %.12e",i,(*fret));
 1407:       fprintf(ficlog,"%d %.12e",i,(*fret));
 1408:       for (j=1;j<=n;j++) {
 1409: 	xits[j]=FMAX(fabs(p[j]-pt[j]),1.e-5);
 1410: 	printf(" x(%d)=%.12e",j,xit[j]);
 1411: 	fprintf(ficlog," x(%d)=%.12e",j,xit[j]);
 1412:       }
 1413:       for(j=1;j<=n;j++) {
 1414: 	printf(" p(%d)=%.12e",j,p[j]);
 1415: 	fprintf(ficlog," p(%d)=%.12e",j,p[j]);
 1416:       }
 1417:       printf("\n");
 1418:       fprintf(ficlog,"\n");
 1419: #endif
 1420:     } /* end i */
 1421:     if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret))) {
 1422: #ifdef DEBUG
 1423:       int k[2],l;
 1424:       k[0]=1;
 1425:       k[1]=-1;
 1426:       printf("Max: %.12e",(*func)(p));
 1427:       fprintf(ficlog,"Max: %.12e",(*func)(p));
 1428:       for (j=1;j<=n;j++) {
 1429: 	printf(" %.12e",p[j]);
 1430: 	fprintf(ficlog," %.12e",p[j]);
 1431:       }
 1432:       printf("\n");
 1433:       fprintf(ficlog,"\n");
 1434:       for(l=0;l<=1;l++) {
 1435: 	for (j=1;j<=n;j++) {
 1436: 	  ptt[j]=p[j]+(p[j]-pt[j])*k[l];
 1437: 	  printf("l=%d j=%d ptt=%.12e, xits=%.12e, p=%.12e, xit=%.12e", l,j,ptt[j],xits[j],p[j],xit[j]);
 1438: 	  fprintf(ficlog,"l=%d j=%d ptt=%.12e, xits=%.12e, p=%.12e, xit=%.12e", l,j,ptt[j],xits[j],p[j],xit[j]);
 1439: 	}
 1440: 	printf("func(ptt)=%.12e, deriv=%.12e\n",(*func)(ptt),(ptt[j]-p[j])/((*func)(ptt)-(*func)(p)));
 1441: 	fprintf(ficlog,"func(ptt)=%.12e, deriv=%.12e\n",(*func)(ptt),(ptt[j]-p[j])/((*func)(ptt)-(*func)(p)));
 1442:       }
 1443: #endif
 1444: 
 1445: 
 1446:       free_vector(xit,1,n); 
 1447:       free_vector(xits,1,n); 
 1448:       free_vector(ptt,1,n); 
 1449:       free_vector(pt,1,n); 
 1450:       return; 
 1451:     } 
 1452:     if (*iter == ITMAX) nrerror("powell exceeding maximum iterations."); 
 1453:     for (j=1;j<=n;j++) { /* Computes an extrapolated point */
 1454:       ptt[j]=2.0*p[j]-pt[j]; 
 1455:       xit[j]=p[j]-pt[j]; 
 1456:       pt[j]=p[j]; 
 1457:     } 
 1458:     fptt=(*func)(ptt); 
 1459:     if (fptt < fp) { /* If extrapolated point is better, decide if we keep that new direction or not */
 1460:       /* (x1 f1=fp), (x2 f2=*fret), (x3 f3=fptt), (xm fm) */
 1461:       /* From x1 (P0) distance of x2 is at h and x3 is 2h */
 1462:       /* Let f"(x2) be the 2nd derivative equal everywhere.  */
 1463:       /* Then the parabolic through (x1,f1), (x2,f2) and (x3,f3) */
 1464:       /* will reach at f3 = fm + h^2/2 f"m  ; f" = (f1 -2f2 +f3 ) / h**2 */
 1465:       /* f1-f3 = delta(2h) = 2 h**2 f'' = 2(f1- 2f2 +f3) */
 1466:       /* Thus we compare delta(2h) with observed f1-f3 */
 1467:       /* or best gain on one ancient line 'del' with total  */
 1468:       /* gain f1-f2 = f1 - f2 - 'del' with del  */
 1469:       /* t=2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del)-del*SQR(fp-fptt); */
 1470: 
 1471:       t=2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del);
 1472:       t= t- del*SQR(fp-fptt);
 1473:       printf("t1= %.12lf, t2= %.12lf, t=%.12lf\n", 2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del),del*SQR(fp-fptt),t);
 1474:       fprintf(ficlog,"t1= %.12lf, t2= %.12lf, t=%.12lf\n", 2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del),del*SQR(fp-fptt),t);
 1475: #ifdef DEBUG
 1476:       printf("t3= %.12lf, t4= %.12lf, t3*= %.12lf, t4*= %.12lf\n",SQR(fp-(*fret)-del),SQR(fp-fptt),
 1477: 	     (fp-(*fret)-del)*(fp-(*fret)-del),(fp-fptt)*(fp-fptt));
 1478:       fprintf(ficlog,"t3= %.12lf, t4= %.12lf, t3*= %.12lf, t4*= %.12lf\n",SQR(fp-(*fret)-del),SQR(fp-fptt),
 1479: 	     (fp-(*fret)-del)*(fp-(*fret)-del),(fp-fptt)*(fp-fptt));
 1480:       printf("tt= %.12lf, t=%.12lf\n",2.0*(fp-2.0*(*fret)+fptt)*(fp-(*fret)-del)*(fp-(*fret)-del)-del*(fp-fptt)*(fp-fptt),t);
 1481:       fprintf(ficlog, "tt= %.12lf, t=%.12lf\n",2.0*(fp-2.0*(*fret)+fptt)*(fp-(*fret)-del)*(fp-(*fret)-del)-del*(fp-fptt)*(fp-fptt),t);
 1482: #endif
 1483:       if (t < 0.0) { /* Then we use it for last direction */
 1484: 	linmin(p,xit,n,fret,func); /* computes mean on the extrapolated direction.*/
 1485: 	for (j=1;j<=n;j++) { 
 1486: 	  xi[j][ibig]=xi[j][n]; /* Replace the direction with biggest decrease by n */
 1487: 	  xi[j][n]=xit[j];      /* and nth direction by the extrapolated */
 1488: 	}
 1489: 	printf("Gaining to use average direction of P0 P%d instead of biggest increase direction %d :\n",n,ibig);
 1490: 	fprintf(ficlog,"Gaining to use average direction of P0 P%d instead of biggest increase direction :\n",n,ibig);
 1491: 
 1492: #ifdef DEBUG
 1493: 	printf("Direction changed  last moved %d in place of ibig=%d, new last is the average:\n",n,ibig);
 1494: 	fprintf(ficlog,"Direction changed  last moved %d in place of ibig=%d, new last is the average:\n",n,ibig);
 1495: 	for(j=1;j<=n;j++){
 1496: 	  printf(" %.12e",xit[j]);
 1497: 	  fprintf(ficlog," %.12e",xit[j]);
 1498: 	}
 1499: 	printf("\n");
 1500: 	fprintf(ficlog,"\n");
 1501: #endif
 1502:       } /* end of t negative */
 1503:     } /* end if (fptt < fp)  */
 1504:   } 
 1505: } 
 1506: 
 1507: /**** Prevalence limit (stable or period prevalence)  ****************/
 1508: 
 1509: double **prevalim(double **prlim, int nlstate, double x[], double age, double **oldm, double **savm, double ftolpl, int ij)
 1510: {
 1511:   /* Computes the prevalence limit in each live state at age x by left multiplying the unit
 1512:      matrix by transitions matrix until convergence is reached */
 1513: 
 1514:   int i, ii,j,k;
 1515:   double min, max, maxmin, maxmax,sumnew=0.;
 1516:   /* double **matprod2(); */ /* test */
 1517:   double **out, cov[NCOVMAX+1], **pmij();
 1518:   double **newm;
 1519:   double agefin, delaymax=50 ; /* Max number of years to converge */
 1520: 
 1521:   for (ii=1;ii<=nlstate+ndeath;ii++)
 1522:     for (j=1;j<=nlstate+ndeath;j++){
 1523:       oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1524:     }
 1525: 
 1526:    cov[1]=1.;
 1527:  
 1528:  /* Even if hstepm = 1, at least one multiplication by the unit matrix */
 1529:   for(agefin=age-stepm/YEARM; agefin>=age-delaymax; agefin=agefin-stepm/YEARM){
 1530:     newm=savm;
 1531:     /* Covariates have to be included here again */
 1532:     cov[2]=agefin;
 1533:     
 1534:     for (k=1; k<=cptcovn;k++) {
 1535:       cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]];
 1536:       /*printf("prevalim ij=%d k=%d Tvar[%d]=%d nbcode=%d cov=%lf codtab[%d][Tvar[%d]]=%d \n",ij,k, k, Tvar[k],nbcode[Tvar[k]][codtab[ij][Tvar[k]]],cov[2+k], ij, k, codtab[ij][Tvar[k]]);*/
 1537:     }
 1538:     /* for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2]; */
 1539:     /* for (k=1; k<=cptcovprod;k++) /\* Useless *\/ */
 1540:     /*   cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]] * nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]]; */
 1541:     
 1542:     /*printf("ij=%d cptcovprod=%d tvar=%d ", ij, cptcovprod, Tvar[1]);*/
 1543:     /*printf("ij=%d cov[3]=%lf cov[4]=%lf \n",ij, cov[3],cov[4]);*/
 1544:     /*printf("ij=%d cov[3]=%lf \n",ij, cov[3]);*/
 1545:     /* savm=pmij(pmmij,cov,ncovmodel,x,nlstate); */
 1546:     /* out=matprod2(newm, pmij(pmmij,cov,ncovmodel,x,nlstate),1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, oldm); /\* Bug Valgrind *\/ */
 1547:     out=matprod2(newm, pmij(pmmij,cov,ncovmodel,x,nlstate),1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, oldm); /* Bug Valgrind */
 1548:     
 1549:     savm=oldm;
 1550:     oldm=newm;
 1551:     maxmax=0.;
 1552:     for(j=1;j<=nlstate;j++){
 1553:       min=1.;
 1554:       max=0.;
 1555:       for(i=1; i<=nlstate; i++) {
 1556: 	sumnew=0;
 1557: 	for(k=1; k<=ndeath; k++) sumnew+=newm[i][nlstate+k];
 1558: 	prlim[i][j]= newm[i][j]/(1-sumnew);
 1559:         /*printf(" prevalim i=%d, j=%d, prmlim[%d][%d]=%f, agefin=%d \n", i, j, i, j, prlim[i][j],(int)agefin);*/
 1560: 	max=FMAX(max,prlim[i][j]);
 1561: 	min=FMIN(min,prlim[i][j]);
 1562:       }
 1563:       maxmin=max-min;
 1564:       maxmax=FMAX(maxmax,maxmin);
 1565:     }
 1566:     if(maxmax < ftolpl){
 1567:       return prlim;
 1568:     }
 1569:   }
 1570: }
 1571: 
 1572: /*************** transition probabilities ***************/ 
 1573: 
 1574: double **pmij(double **ps, double *cov, int ncovmodel, double *x, int nlstate )
 1575: {
 1576:   /* According to parameters values stored in x and the covariate's values stored in cov,
 1577:      computes the probability to be observed in state j being in state i by appying the
 1578:      model to the ncovmodel covariates (including constant and age).
 1579:      lnpijopii=ln(pij/pii)= aij+bij*age+cij*v1+dij*v2+... = sum_nc=1^ncovmodel xij(nc)*cov[nc]
 1580:      and, according on how parameters are entered, the position of the coefficient xij(nc) of the
 1581:      ncth covariate in the global vector x is given by the formula:
 1582:      j<i nc+((i-1)*(nlstate+ndeath-1)+j-1)*ncovmodel
 1583:      j>=i nc + ((i-1)*(nlstate+ndeath-1)+(j-2))*ncovmodel
 1584:      Computes ln(pij/pii) (lnpijopii), deduces pij/pii by exponentiation,
 1585:      sums on j different of i to get 1-pii/pii, deduces pii, and then all pij.
 1586:      Outputs ps[i][j] the probability to be observed in j being in j according to
 1587:      the values of the covariates cov[nc] and corresponding parameter values x[nc+shiftij]
 1588:   */
 1589:   double s1, lnpijopii;
 1590:   /*double t34;*/
 1591:   int i,j, nc, ii, jj;
 1592: 
 1593:     for(i=1; i<= nlstate; i++){
 1594:       for(j=1; j<i;j++){
 1595: 	for (nc=1, lnpijopii=0.;nc <=ncovmodel; nc++){
 1596: 	  /*lnpijopii += param[i][j][nc]*cov[nc];*/
 1597: 	  lnpijopii += x[nc+((i-1)*(nlstate+ndeath-1)+j-1)*ncovmodel]*cov[nc];
 1598: /* 	 printf("Int j<i s1=%.17e, lnpijopii=%.17e\n",s1,lnpijopii); */
 1599: 	}
 1600: 	ps[i][j]=lnpijopii; /* In fact ln(pij/pii) */
 1601: /* 	printf("s1=%.17e, lnpijopii=%.17e\n",s1,lnpijopii); */
 1602:       }
 1603:       for(j=i+1; j<=nlstate+ndeath;j++){
 1604: 	for (nc=1, lnpijopii=0.;nc <=ncovmodel; nc++){
 1605: 	  /*lnpijopii += x[(i-1)*nlstate*ncovmodel+(j-2)*ncovmodel+nc+(i-1)*(ndeath-1)*ncovmodel]*cov[nc];*/
 1606: 	  lnpijopii += x[nc + ((i-1)*(nlstate+ndeath-1)+(j-2))*ncovmodel]*cov[nc];
 1607: /* 	  printf("Int j>i s1=%.17e, lnpijopii=%.17e %lx %lx\n",s1,lnpijopii,s1,lnpijopii); */
 1608: 	}
 1609: 	ps[i][j]=lnpijopii; /* In fact ln(pij/pii) */
 1610:       }
 1611:     }
 1612:     
 1613:     for(i=1; i<= nlstate; i++){
 1614:       s1=0;
 1615:       for(j=1; j<i; j++){
 1616: 	s1+=exp(ps[i][j]); /* In fact sums pij/pii */
 1617: 	/*printf("debug1 %d %d ps=%lf exp(ps)=%lf s1+=%lf\n",i,j,ps[i][j],exp(ps[i][j]),s1); */
 1618:       }
 1619:       for(j=i+1; j<=nlstate+ndeath; j++){
 1620: 	s1+=exp(ps[i][j]); /* In fact sums pij/pii */
 1621: 	/*printf("debug2 %d %d ps=%lf exp(ps)=%lf s1+=%lf\n",i,j,ps[i][j],exp(ps[i][j]),s1); */
 1622:       }
 1623:       /* s1= sum_{j<>i} pij/pii=(1-pii)/pii and thus pii is known from s1 */
 1624:       ps[i][i]=1./(s1+1.);
 1625:       /* Computing other pijs */
 1626:       for(j=1; j<i; j++)
 1627: 	ps[i][j]= exp(ps[i][j])*ps[i][i];
 1628:       for(j=i+1; j<=nlstate+ndeath; j++)
 1629: 	ps[i][j]= exp(ps[i][j])*ps[i][i];
 1630:       /* ps[i][nlstate+1]=1.-s1- ps[i][i];*/ /* Sum should be 1 */
 1631:     } /* end i */
 1632:     
 1633:     for(ii=nlstate+1; ii<= nlstate+ndeath; ii++){
 1634:       for(jj=1; jj<= nlstate+ndeath; jj++){
 1635: 	ps[ii][jj]=0;
 1636: 	ps[ii][ii]=1;
 1637:       }
 1638:     }
 1639:     
 1640:     
 1641:     /* for(ii=1; ii<= nlstate+ndeath; ii++){ */
 1642:     /*   for(jj=1; jj<= nlstate+ndeath; jj++){ */
 1643:     /* 	printf(" pmij  ps[%d][%d]=%lf ",ii,jj,ps[ii][jj]); */
 1644:     /*   } */
 1645:     /*   printf("\n "); */
 1646:     /* } */
 1647:     /* printf("\n ");printf("%lf ",cov[2]);*/
 1648:     /*
 1649:       for(i=1; i<= npar; i++) printf("%f ",x[i]);
 1650:       goto end;*/
 1651:     return ps;
 1652: }
 1653: 
 1654: /**************** Product of 2 matrices ******************/
 1655: 
 1656: double **matprod2(double **out, double **in,int nrl, int nrh, int ncl, int nch, int ncolol, int ncoloh, double **b)
 1657: {
 1658:   /* Computes the matrix product of in(1,nrh-nrl+1)(1,nch-ncl+1) times
 1659:      b(1,nch-ncl+1)(1,ncoloh-ncolol+1) into out(...) */
 1660:   /* in, b, out are matrice of pointers which should have been initialized 
 1661:      before: only the contents of out is modified. The function returns
 1662:      a pointer to pointers identical to out */
 1663:   int i, j, k;
 1664:   for(i=nrl; i<= nrh; i++)
 1665:     for(k=ncolol; k<=ncoloh; k++){
 1666:       out[i][k]=0.;
 1667:       for(j=ncl; j<=nch; j++)
 1668:   	out[i][k] +=in[i][j]*b[j][k];
 1669:     }
 1670:   return out;
 1671: }
 1672: 
 1673: 
 1674: /************* Higher Matrix Product ***************/
 1675: 
 1676: double ***hpxij(double ***po, int nhstepm, double age, int hstepm, double *x, int nlstate, int stepm, double **oldm, double **savm, int ij )
 1677: {
 1678:   /* Computes the transition matrix starting at age 'age' over 
 1679:      'nhstepm*hstepm*stepm' months (i.e. until
 1680:      age (in years)  age+nhstepm*hstepm*stepm/12) by multiplying 
 1681:      nhstepm*hstepm matrices. 
 1682:      Output is stored in matrix po[i][j][h] for h every 'hstepm' step 
 1683:      (typically every 2 years instead of every month which is too big 
 1684:      for the memory).
 1685:      Model is determined by parameters x and covariates have to be 
 1686:      included manually here. 
 1687: 
 1688:      */
 1689: 
 1690:   int i, j, d, h, k;
 1691:   double **out, cov[NCOVMAX+1];
 1692:   double **newm;
 1693: 
 1694:   /* Hstepm could be zero and should return the unit matrix */
 1695:   for (i=1;i<=nlstate+ndeath;i++)
 1696:     for (j=1;j<=nlstate+ndeath;j++){
 1697:       oldm[i][j]=(i==j ? 1.0 : 0.0);
 1698:       po[i][j][0]=(i==j ? 1.0 : 0.0);
 1699:     }
 1700:   /* Even if hstepm = 1, at least one multiplication by the unit matrix */
 1701:   for(h=1; h <=nhstepm; h++){
 1702:     for(d=1; d <=hstepm; d++){
 1703:       newm=savm;
 1704:       /* Covariates have to be included here again */
 1705:       cov[1]=1.;
 1706:       cov[2]=age+((h-1)*hstepm + (d-1))*stepm/YEARM;
 1707:       for (k=1; k<=cptcovn;k++) 
 1708: 	cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]];
 1709:       for (k=1; k<=cptcovage;k++)
 1710: 	cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
 1711:       for (k=1; k<=cptcovprod;k++) /* Useless because included in cptcovn */
 1712: 	cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
 1713: 
 1714: 
 1715:       /*printf("hxi cptcov=%d cptcode=%d\n",cptcov,cptcode);*/
 1716:       /*printf("h=%d d=%d age=%f cov=%f\n",h,d,age,cov[2]);*/
 1717:       out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, 
 1718: 		   pmij(pmmij,cov,ncovmodel,x,nlstate));
 1719:       savm=oldm;
 1720:       oldm=newm;
 1721:     }
 1722:     for(i=1; i<=nlstate+ndeath; i++)
 1723:       for(j=1;j<=nlstate+ndeath;j++) {
 1724: 	po[i][j][h]=newm[i][j];
 1725: 	/*if(h==nhstepm) printf("po[%d][%d][%d]=%f ",i,j,h,po[i][j][h]);*/
 1726:       }
 1727:     /*printf("h=%d ",h);*/
 1728:   } /* end h */
 1729: /*     printf("\n H=%d \n",h); */
 1730:   return po;
 1731: }
 1732: 
 1733: #ifdef NLOPT
 1734:   double  myfunc(unsigned n, const double *p1, double *grad, void *pd){
 1735:   double fret;
 1736:   double *xt;
 1737:   int j;
 1738:   myfunc_data *d2 = (myfunc_data *) pd;
 1739: /* xt = (p1-1); */
 1740:   xt=vector(1,n); 
 1741:   for (j=1;j<=n;j++)   xt[j]=p1[j-1]; /* xt[1]=p1[0] */
 1742: 
 1743:   fret=(d2->function)(xt); /*  p xt[1]@8 is fine */
 1744:   /* fret=(*func)(xt); /\*  p xt[1]@8 is fine *\/ */
 1745:   printf("Function = %.12lf ",fret);
 1746:   for (j=1;j<=n;j++) printf(" %d %.8lf", j, xt[j]); 
 1747:   printf("\n");
 1748:  free_vector(xt,1,n);
 1749:   return fret;
 1750: }
 1751: #endif
 1752: 
 1753: /*************** log-likelihood *************/
 1754: double func( double *x)
 1755: {
 1756:   int i, ii, j, k, mi, d, kk;
 1757:   double l, ll[NLSTATEMAX+1], cov[NCOVMAX+1];
 1758:   double **out;
 1759:   double sw; /* Sum of weights */
 1760:   double lli; /* Individual log likelihood */
 1761:   int s1, s2;
 1762:   double bbh, survp;
 1763:   long ipmx;
 1764:   /*extern weight */
 1765:   /* We are differentiating ll according to initial status */
 1766:   /*  for (i=1;i<=npar;i++) printf("%f ", x[i]);*/
 1767:   /*for(i=1;i<imx;i++) 
 1768:     printf(" %d\n",s[4][i]);
 1769:   */
 1770: 
 1771:   ++countcallfunc;
 1772: 
 1773:   cov[1]=1.;
 1774: 
 1775:   for(k=1; k<=nlstate; k++) ll[k]=0.;
 1776: 
 1777:   if(mle==1){
 1778:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1779:       /* Computes the values of the ncovmodel covariates of the model
 1780: 	 depending if the covariates are fixed or variying (age dependent) and stores them in cov[]
 1781: 	 Then computes with function pmij which return a matrix p[i][j] giving the elementary probability
 1782: 	 to be observed in j being in i according to the model.
 1783:        */
 1784:       for (k=1; k<=cptcovn;k++){ /* Simple and product covariates without age* products */
 1785: 	cov[2+k]=covar[Tvar[k]][i];
 1786:       }
 1787:       /* In model V2+V1*V4+age*V3+V3*V2 Tvar[1] is V2, Tvar[2=V1*V4] 
 1788: 	 is 6, Tvar[3=age*V3] should not be computed because of age Tvar[4=V3*V2] 
 1789: 	 has been calculated etc */
 1790:       for(mi=1; mi<= wav[i]-1; mi++){
 1791: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1792: 	  for (j=1;j<=nlstate+ndeath;j++){
 1793: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1794: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1795: 	  }
 1796: 	for(d=0; d<dh[mi][i]; d++){
 1797: 	  newm=savm;
 1798: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1799: 	  for (kk=1; kk<=cptcovage;kk++) {
 1800: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2]; /* Tage[kk] gives the data-covariate associated with age */
 1801: 	  }
 1802: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1803: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1804: 	  savm=oldm;
 1805: 	  oldm=newm;
 1806: 	} /* end mult */
 1807:       
 1808: 	/*lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]);*/ /* Original formula */
 1809: 	/* But now since version 0.9 we anticipate for bias at large stepm.
 1810: 	 * If stepm is larger than one month (smallest stepm) and if the exact delay 
 1811: 	 * (in months) between two waves is not a multiple of stepm, we rounded to 
 1812: 	 * the nearest (and in case of equal distance, to the lowest) interval but now
 1813: 	 * we keep into memory the bias bh[mi][i] and also the previous matrix product
 1814: 	 * (i.e to dh[mi][i]-1) saved in 'savm'. Then we inter(extra)polate the
 1815: 	 * probability in order to take into account the bias as a fraction of the way
 1816: 	 * from savm to out if bh is negative or even beyond if bh is positive. bh varies
 1817: 	 * -stepm/2 to stepm/2 .
 1818: 	 * For stepm=1 the results are the same as for previous versions of Imach.
 1819: 	 * For stepm > 1 the results are less biased than in previous versions. 
 1820: 	 */
 1821: 	s1=s[mw[mi][i]][i];
 1822: 	s2=s[mw[mi+1][i]][i];
 1823: 	bbh=(double)bh[mi][i]/(double)stepm; 
 1824: 	/* bias bh is positive if real duration
 1825: 	 * is higher than the multiple of stepm and negative otherwise.
 1826: 	 */
 1827: 	/* lli= (savm[s1][s2]>1.e-8 ?(1.+bbh)*log(out[s1][s2])- bbh*log(savm[s1][s2]):log((1.+bbh)*out[s1][s2]));*/
 1828: 	if( s2 > nlstate){ 
 1829: 	  /* i.e. if s2 is a death state and if the date of death is known 
 1830: 	     then the contribution to the likelihood is the probability to 
 1831: 	     die between last step unit time and current  step unit time, 
 1832: 	     which is also equal to probability to die before dh 
 1833: 	     minus probability to die before dh-stepm . 
 1834: 	     In version up to 0.92 likelihood was computed
 1835: 	as if date of death was unknown. Death was treated as any other
 1836: 	health state: the date of the interview describes the actual state
 1837: 	and not the date of a change in health state. The former idea was
 1838: 	to consider that at each interview the state was recorded
 1839: 	(healthy, disable or death) and IMaCh was corrected; but when we
 1840: 	introduced the exact date of death then we should have modified
 1841: 	the contribution of an exact death to the likelihood. This new
 1842: 	contribution is smaller and very dependent of the step unit
 1843: 	stepm. It is no more the probability to die between last interview
 1844: 	and month of death but the probability to survive from last
 1845: 	interview up to one month before death multiplied by the
 1846: 	probability to die within a month. Thanks to Chris
 1847: 	Jackson for correcting this bug.  Former versions increased
 1848: 	mortality artificially. The bad side is that we add another loop
 1849: 	which slows down the processing. The difference can be up to 10%
 1850: 	lower mortality.
 1851: 	  */
 1852: 	  lli=log(out[s1][s2] - savm[s1][s2]);
 1853: 
 1854: 
 1855: 	} else if  (s2==-2) {
 1856: 	  for (j=1,survp=0. ; j<=nlstate; j++) 
 1857: 	    survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 1858: 	  /*survp += out[s1][j]; */
 1859: 	  lli= log(survp);
 1860: 	}
 1861: 	
 1862:  	else if  (s2==-4) { 
 1863: 	  for (j=3,survp=0. ; j<=nlstate; j++)  
 1864: 	    survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 1865:  	  lli= log(survp); 
 1866:  	} 
 1867: 
 1868:  	else if  (s2==-5) { 
 1869:  	  for (j=1,survp=0. ; j<=2; j++)  
 1870: 	    survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 1871:  	  lli= log(survp); 
 1872:  	} 
 1873: 	
 1874: 	else{
 1875: 	  lli= log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]); /* linear interpolation */
 1876: 	  /*  lli= (savm[s1][s2]>(double)1.e-8 ?log((1.+bbh)*out[s1][s2]- bbh*(savm[s1][s2])):log((1.+bbh)*out[s1][s2]));*/ /* linear interpolation */
 1877: 	} 
 1878: 	/*lli=(1.+bbh)*log(out[s1][s2])- bbh*log(savm[s1][s2]);*/
 1879: 	/*if(lli ==000.0)*/
 1880: 	/*printf("bbh= %f lli=%f savm=%f out=%f %d\n",bbh,lli,savm[s1][s2], out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]],i); */
 1881:   	ipmx +=1;
 1882: 	sw += weight[i];
 1883: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1884:       } /* end of wave */
 1885:     } /* end of individual */
 1886:   }  else if(mle==2){
 1887:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1888:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1889:       for(mi=1; mi<= wav[i]-1; mi++){
 1890: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1891: 	  for (j=1;j<=nlstate+ndeath;j++){
 1892: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1893: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1894: 	  }
 1895: 	for(d=0; d<=dh[mi][i]; d++){
 1896: 	  newm=savm;
 1897: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1898: 	  for (kk=1; kk<=cptcovage;kk++) {
 1899: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1900: 	  }
 1901: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1902: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1903: 	  savm=oldm;
 1904: 	  oldm=newm;
 1905: 	} /* end mult */
 1906:       
 1907: 	s1=s[mw[mi][i]][i];
 1908: 	s2=s[mw[mi+1][i]][i];
 1909: 	bbh=(double)bh[mi][i]/(double)stepm; 
 1910: 	lli= (savm[s1][s2]>(double)1.e-8 ?log((1.+bbh)*out[s1][s2]- bbh*(savm[s1][s2])):log((1.+bbh)*out[s1][s2])); /* linear interpolation */
 1911: 	ipmx +=1;
 1912: 	sw += weight[i];
 1913: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1914:       } /* end of wave */
 1915:     } /* end of individual */
 1916:   }  else if(mle==3){  /* exponential inter-extrapolation */
 1917:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1918:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1919:       for(mi=1; mi<= wav[i]-1; mi++){
 1920: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1921: 	  for (j=1;j<=nlstate+ndeath;j++){
 1922: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1923: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1924: 	  }
 1925: 	for(d=0; d<dh[mi][i]; d++){
 1926: 	  newm=savm;
 1927: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1928: 	  for (kk=1; kk<=cptcovage;kk++) {
 1929: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1930: 	  }
 1931: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1932: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1933: 	  savm=oldm;
 1934: 	  oldm=newm;
 1935: 	} /* end mult */
 1936:       
 1937: 	s1=s[mw[mi][i]][i];
 1938: 	s2=s[mw[mi+1][i]][i];
 1939: 	bbh=(double)bh[mi][i]/(double)stepm; 
 1940: 	lli= (savm[s1][s2]>1.e-8 ?(1.+bbh)*log(out[s1][s2])- bbh*log(savm[s1][s2]):log((1.+bbh)*out[s1][s2])); /* exponential inter-extrapolation */
 1941: 	ipmx +=1;
 1942: 	sw += weight[i];
 1943: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1944:       } /* end of wave */
 1945:     } /* end of individual */
 1946:   }else if (mle==4){  /* ml=4 no inter-extrapolation */
 1947:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1948:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1949:       for(mi=1; mi<= wav[i]-1; mi++){
 1950: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1951: 	  for (j=1;j<=nlstate+ndeath;j++){
 1952: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1953: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1954: 	  }
 1955: 	for(d=0; d<dh[mi][i]; d++){
 1956: 	  newm=savm;
 1957: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1958: 	  for (kk=1; kk<=cptcovage;kk++) {
 1959: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1960: 	  }
 1961: 	
 1962: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1963: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1964: 	  savm=oldm;
 1965: 	  oldm=newm;
 1966: 	} /* end mult */
 1967:       
 1968: 	s1=s[mw[mi][i]][i];
 1969: 	s2=s[mw[mi+1][i]][i];
 1970: 	if( s2 > nlstate){ 
 1971: 	  lli=log(out[s1][s2] - savm[s1][s2]);
 1972: 	}else{
 1973: 	  lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]); /* Original formula */
 1974: 	}
 1975: 	ipmx +=1;
 1976: 	sw += weight[i];
 1977: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1978: /* 	printf("i=%6d s1=%1d s2=%1d mi=%1d mw=%1d dh=%3d prob=%10.6f w=%6.4f out=%10.6f sav=%10.6f\n",i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],out[s1][s2],savm[s1][s2]); */
 1979:       } /* end of wave */
 1980:     } /* end of individual */
 1981:   }else{  /* ml=5 no inter-extrapolation no jackson =0.8a */
 1982:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1983:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1984:       for(mi=1; mi<= wav[i]-1; mi++){
 1985: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1986: 	  for (j=1;j<=nlstate+ndeath;j++){
 1987: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1988: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1989: 	  }
 1990: 	for(d=0; d<dh[mi][i]; d++){
 1991: 	  newm=savm;
 1992: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1993: 	  for (kk=1; kk<=cptcovage;kk++) {
 1994: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1995: 	  }
 1996: 	
 1997: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1998: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1999: 	  savm=oldm;
 2000: 	  oldm=newm;
 2001: 	} /* end mult */
 2002:       
 2003: 	s1=s[mw[mi][i]][i];
 2004: 	s2=s[mw[mi+1][i]][i];
 2005: 	lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]); /* Original formula */
 2006: 	ipmx +=1;
 2007: 	sw += weight[i];
 2008: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 2009: 	/*printf("i=%6d s1=%1d s2=%1d mi=%1d mw=%1d dh=%3d prob=%10.6f w=%6.4f out=%10.6f sav=%10.6f\n",i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],out[s1][s2],savm[s1][s2]);*/
 2010:       } /* end of wave */
 2011:     } /* end of individual */
 2012:   } /* End of if */
 2013:   for(k=1,l=0.; k<=nlstate; k++) l += ll[k];
 2014:   /* printf("l1=%f l2=%f ",ll[1],ll[2]); */
 2015:   l= l*ipmx/sw; /* To get the same order of magnitude as if weight=1 for every body */
 2016:   return -l;
 2017: }
 2018: 
 2019: /*************** log-likelihood *************/
 2020: double funcone( double *x)
 2021: {
 2022:   /* Same as likeli but slower because of a lot of printf and if */
 2023:   int i, ii, j, k, mi, d, kk;
 2024:   double l, ll[NLSTATEMAX+1], cov[NCOVMAX+1];
 2025:   double **out;
 2026:   double lli; /* Individual log likelihood */
 2027:   double llt;
 2028:   int s1, s2;
 2029:   double bbh, survp;
 2030:   /*extern weight */
 2031:   /* We are differentiating ll according to initial status */
 2032:   /*  for (i=1;i<=npar;i++) printf("%f ", x[i]);*/
 2033:   /*for(i=1;i<imx;i++) 
 2034:     printf(" %d\n",s[4][i]);
 2035:   */
 2036:   cov[1]=1.;
 2037: 
 2038:   for(k=1; k<=nlstate; k++) ll[k]=0.;
 2039: 
 2040:   for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 2041:     for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 2042:     for(mi=1; mi<= wav[i]-1; mi++){
 2043:       for (ii=1;ii<=nlstate+ndeath;ii++)
 2044: 	for (j=1;j<=nlstate+ndeath;j++){
 2045: 	  oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 2046: 	  savm[ii][j]=(ii==j ? 1.0 : 0.0);
 2047: 	}
 2048:       for(d=0; d<dh[mi][i]; d++){
 2049: 	newm=savm;
 2050: 	cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 2051: 	for (kk=1; kk<=cptcovage;kk++) {
 2052: 	  cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 2053: 	}
 2054: 	/* savm=pmij(pmmij,cov,ncovmodel,x,nlstate); */
 2055: 	out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 2056: 		     1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 2057: 	/* out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath, */
 2058: 	/* 	     1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate)); */
 2059: 	savm=oldm;
 2060: 	oldm=newm;
 2061:       } /* end mult */
 2062:       
 2063:       s1=s[mw[mi][i]][i];
 2064:       s2=s[mw[mi+1][i]][i];
 2065:       bbh=(double)bh[mi][i]/(double)stepm; 
 2066:       /* bias is positive if real duration
 2067:        * is higher than the multiple of stepm and negative otherwise.
 2068:        */
 2069:       if( s2 > nlstate && (mle <5) ){  /* Jackson */
 2070: 	lli=log(out[s1][s2] - savm[s1][s2]);
 2071:       } else if  (s2==-2) {
 2072: 	for (j=1,survp=0. ; j<=nlstate; j++) 
 2073: 	  survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 2074: 	lli= log(survp);
 2075:       }else if (mle==1){
 2076: 	lli= log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]); /* linear interpolation */
 2077:       } else if(mle==2){
 2078: 	lli= (savm[s1][s2]>(double)1.e-8 ?log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]):log((1.+bbh)*out[s1][s2])); /* linear interpolation */
 2079:       } else if(mle==3){  /* exponential inter-extrapolation */
 2080: 	lli= (savm[s1][s2]>(double)1.e-8 ?(1.+bbh)*log(out[s1][s2])- bbh*log(savm[s1][s2]):log((1.+bbh)*out[s1][s2])); /* exponential inter-extrapolation */
 2081:       } else if (mle==4){  /* mle=4 no inter-extrapolation */
 2082: 	lli=log(out[s1][s2]); /* Original formula */
 2083:       } else{  /* mle=0 back to 1 */
 2084: 	lli= log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]); /* linear interpolation */
 2085: 	/*lli=log(out[s1][s2]); */ /* Original formula */
 2086:       } /* End of if */
 2087:       ipmx +=1;
 2088:       sw += weight[i];
 2089:       ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 2090:       /*printf("i=%6d s1=%1d s2=%1d mi=%1d mw=%1d dh=%3d prob=%10.6f w=%6.4f out=%10.6f sav=%10.6f\n",i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],out[s1][s2],savm[s1][s2]); */
 2091:       if(globpr){
 2092: 	fprintf(ficresilk,"%9ld %6d %2d %2d %1d %1d %3d %11.6f %8.4f\
 2093:  %11.6f %11.6f %11.6f ", \
 2094: 		num[i],i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],
 2095: 		2*weight[i]*lli,out[s1][s2],savm[s1][s2]);
 2096: 	for(k=1,llt=0.,l=0.; k<=nlstate; k++){
 2097: 	  llt +=ll[k]*gipmx/gsw;
 2098: 	  fprintf(ficresilk," %10.6f",-ll[k]*gipmx/gsw);
 2099: 	}
 2100: 	fprintf(ficresilk," %10.6f\n", -llt);
 2101:       }
 2102:     } /* end of wave */
 2103:   } /* end of individual */
 2104:   for(k=1,l=0.; k<=nlstate; k++) l += ll[k];
 2105:   /* printf("l1=%f l2=%f ",ll[1],ll[2]); */
 2106:   l= l*ipmx/sw; /* To get the same order of magnitude as if weight=1 for every body */
 2107:   if(globpr==0){ /* First time we count the contributions and weights */
 2108:     gipmx=ipmx;
 2109:     gsw=sw;
 2110:   }
 2111:   return -l;
 2112: }
 2113: 
 2114: 
 2115: /*************** function likelione ***********/
 2116: void likelione(FILE *ficres,double p[], int npar, int nlstate, int *globpri, long *ipmx, double *sw, double *fretone, double (*funcone)(double []))
 2117: {
 2118:   /* This routine should help understanding what is done with 
 2119:      the selection of individuals/waves and
 2120:      to check the exact contribution to the likelihood.
 2121:      Plotting could be done.
 2122:    */
 2123:   int k;
 2124: 
 2125:   if(*globpri !=0){ /* Just counts and sums, no printings */
 2126:     strcpy(fileresilk,"ilk"); 
 2127:     strcat(fileresilk,fileres);
 2128:     if((ficresilk=fopen(fileresilk,"w"))==NULL) {
 2129:       printf("Problem with resultfile: %s\n", fileresilk);
 2130:       fprintf(ficlog,"Problem with resultfile: %s\n", fileresilk);
 2131:     }
 2132:     fprintf(ficresilk, "#individual(line's_record) s1 s2 wave# effective_wave# number_of_matrices_product pij weight -2ln(pij)*weight 0pij_x 0pij_(x-stepm) cumulating_loglikeli_by_health_state(reweighted=-2ll*weightXnumber_of_contribs/sum_of_weights) and_total\n");
 2133:     fprintf(ficresilk, "#num_i i s1 s2 mi mw dh likeli weight 2wlli out sav ");
 2134:     /* 	i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],2*weight[i]*lli,out[s1][s2],savm[s1][s2]); */
 2135:     for(k=1; k<=nlstate; k++) 
 2136:       fprintf(ficresilk," -2*gipw/gsw*weight*ll[%d]++",k);
 2137:     fprintf(ficresilk," -2*gipw/gsw*weight*ll(total)\n");
 2138:   }
 2139: 
 2140:   *fretone=(*funcone)(p);
 2141:   if(*globpri !=0){
 2142:     fclose(ficresilk);
 2143:     fprintf(fichtm,"\n<br>File of contributions to the likelihood: <a href=\"%s\">%s</a><br>\n",subdirf(fileresilk),subdirf(fileresilk));
 2144:     fflush(fichtm); 
 2145:   } 
 2146:   return;
 2147: }
 2148: 
 2149: 
 2150: /*********** Maximum Likelihood Estimation ***************/
 2151: 
 2152: void mlikeli(FILE *ficres,double p[], int npar, int ncovmodel, int nlstate, double ftol, double (*func)(double []))
 2153: {
 2154:   int i,j, iter=0;
 2155:   double **xi;
 2156:   double fret;
 2157:   double fretone; /* Only one call to likelihood */
 2158:   /*  char filerespow[FILENAMELENGTH];*/
 2159: 
 2160: #ifdef NLOPT
 2161:   int creturn;
 2162:   nlopt_opt opt;
 2163:   /* double lb[9] = { -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL }; /\* lower bounds *\/ */
 2164:   double *lb;
 2165:   double minf; /* the minimum objective value, upon return */
 2166:   double * p1; /* Shifted parameters from 0 instead of 1 */
 2167:   myfunc_data dinst, *d = &dinst;
 2168: #endif
 2169: 
 2170: 
 2171:   xi=matrix(1,npar,1,npar);
 2172:   for (i=1;i<=npar;i++)
 2173:     for (j=1;j<=npar;j++)
 2174:       xi[i][j]=(i==j ? 1.0 : 0.0);
 2175:   printf("Powell\n");  fprintf(ficlog,"Powell\n");
 2176:   strcpy(filerespow,"pow"); 
 2177:   strcat(filerespow,fileres);
 2178:   if((ficrespow=fopen(filerespow,"w"))==NULL) {
 2179:     printf("Problem with resultfile: %s\n", filerespow);
 2180:     fprintf(ficlog,"Problem with resultfile: %s\n", filerespow);
 2181:   }
 2182:   fprintf(ficrespow,"# Powell\n# iter -2*LL");
 2183:   for (i=1;i<=nlstate;i++)
 2184:     for(j=1;j<=nlstate+ndeath;j++)
 2185:       if(j!=i)fprintf(ficrespow," p%1d%1d",i,j);
 2186:   fprintf(ficrespow,"\n");
 2187: #ifdef POWELL
 2188:   powell(p,xi,npar,ftol,&iter,&fret,func);
 2189: #endif
 2190: 
 2191: #ifdef NLOPT
 2192: #ifdef NEWUOA
 2193:   opt = nlopt_create(NLOPT_LN_NEWUOA,npar);
 2194: #else
 2195:   opt = nlopt_create(NLOPT_LN_BOBYQA,npar);
 2196: #endif
 2197:   lb=vector(0,npar-1);
 2198:   for (i=0;i<npar;i++) lb[i]= -HUGE_VAL;
 2199:   nlopt_set_lower_bounds(opt, lb);
 2200:   nlopt_set_initial_step1(opt, 0.1);
 2201:   
 2202:   p1= (p+1); /*  p *(p+1)@8 and p *(p1)@8 are equal p1[0]=p[1] */
 2203:   d->function = func;
 2204:   printf(" Func %.12lf \n",myfunc(npar,p1,NULL,d));
 2205:   nlopt_set_min_objective(opt, myfunc, d);
 2206:   nlopt_set_xtol_rel(opt, ftol);
 2207:   if ((creturn=nlopt_optimize(opt, p1, &minf)) < 0) {
 2208:     printf("nlopt failed! %d\n",creturn); 
 2209:   }
 2210:   else {
 2211:     printf("found minimum after %d evaluations (NLOPT=%d)\n", countcallfunc ,NLOPT);
 2212:     printf("found minimum at f(%g,%g) = %0.10g\n", p[0], p[1], minf);
 2213:     iter=1; /* not equal */
 2214:   }
 2215:   nlopt_destroy(opt);
 2216: #endif
 2217:   free_matrix(xi,1,npar,1,npar);
 2218:   fclose(ficrespow);
 2219:   printf("\n#Number of iterations & function calls = %d & %d, -2 Log likelihood = %.12f\n",iter, countcallfunc,func(p));
 2220:   fprintf(ficlog,"\n#Number of iterations & function calls = %d & %d, -2 Log likelihood = %.12f\n",iter, countcallfunc,func(p));
 2221:   fprintf(ficres,"\n#Number of iterations & function calls = %d & %d, -2 Log likelihood = %.12f\n",iter, countcallfunc,func(p));
 2222: 
 2223: }
 2224: 
 2225: /**** Computes Hessian and covariance matrix ***/
 2226: void hesscov(double **matcov, double p[], int npar, double delti[], double ftolhess, double (*func)(double []))
 2227: {
 2228:   double  **a,**y,*x,pd;
 2229:   double **hess;
 2230:   int i, j;
 2231:   int *indx;
 2232: 
 2233:   double hessii(double p[], double delta, int theta, double delti[],double (*func)(double []),int npar);
 2234:   double hessij(double p[], double delti[], int i, int j,double (*func)(double []),int npar);
 2235:   void lubksb(double **a, int npar, int *indx, double b[]) ;
 2236:   void ludcmp(double **a, int npar, int *indx, double *d) ;
 2237:   double gompertz(double p[]);
 2238:   hess=matrix(1,npar,1,npar);
 2239: 
 2240:   printf("\nCalculation of the hessian matrix. Wait...\n");
 2241:   fprintf(ficlog,"\nCalculation of the hessian matrix. Wait...\n");
 2242:   for (i=1;i<=npar;i++){
 2243:     printf("%d",i);fflush(stdout);
 2244:     fprintf(ficlog,"%d",i);fflush(ficlog);
 2245:    
 2246:      hess[i][i]=hessii(p,ftolhess,i,delti,func,npar);
 2247:     
 2248:     /*  printf(" %f ",p[i]);
 2249: 	printf(" %lf %lf %lf",hess[i][i],ftolhess,delti[i]);*/
 2250:   }
 2251:   
 2252:   for (i=1;i<=npar;i++) {
 2253:     for (j=1;j<=npar;j++)  {
 2254:       if (j>i) { 
 2255: 	printf(".%d%d",i,j);fflush(stdout);
 2256: 	fprintf(ficlog,".%d%d",i,j);fflush(ficlog);
 2257: 	hess[i][j]=hessij(p,delti,i,j,func,npar);
 2258: 	
 2259: 	hess[j][i]=hess[i][j];    
 2260: 	/*printf(" %lf ",hess[i][j]);*/
 2261:       }
 2262:     }
 2263:   }
 2264:   printf("\n");
 2265:   fprintf(ficlog,"\n");
 2266: 
 2267:   printf("\nInverting the hessian to get the covariance matrix. Wait...\n");
 2268:   fprintf(ficlog,"\nInverting the hessian to get the covariance matrix. Wait...\n");
 2269:   
 2270:   a=matrix(1,npar,1,npar);
 2271:   y=matrix(1,npar,1,npar);
 2272:   x=vector(1,npar);
 2273:   indx=ivector(1,npar);
 2274:   for (i=1;i<=npar;i++)
 2275:     for (j=1;j<=npar;j++) a[i][j]=hess[i][j];
 2276:   ludcmp(a,npar,indx,&pd);
 2277: 
 2278:   for (j=1;j<=npar;j++) {
 2279:     for (i=1;i<=npar;i++) x[i]=0;
 2280:     x[j]=1;
 2281:     lubksb(a,npar,indx,x);
 2282:     for (i=1;i<=npar;i++){ 
 2283:       matcov[i][j]=x[i];
 2284:     }
 2285:   }
 2286: 
 2287:   printf("\n#Hessian matrix#\n");
 2288:   fprintf(ficlog,"\n#Hessian matrix#\n");
 2289:   for (i=1;i<=npar;i++) { 
 2290:     for (j=1;j<=npar;j++) { 
 2291:       printf("%.3e ",hess[i][j]);
 2292:       fprintf(ficlog,"%.3e ",hess[i][j]);
 2293:     }
 2294:     printf("\n");
 2295:     fprintf(ficlog,"\n");
 2296:   }
 2297: 
 2298:   /* Recompute Inverse */
 2299:   for (i=1;i<=npar;i++)
 2300:     for (j=1;j<=npar;j++) a[i][j]=matcov[i][j];
 2301:   ludcmp(a,npar,indx,&pd);
 2302: 
 2303:   /*  printf("\n#Hessian matrix recomputed#\n");
 2304: 
 2305:   for (j=1;j<=npar;j++) {
 2306:     for (i=1;i<=npar;i++) x[i]=0;
 2307:     x[j]=1;
 2308:     lubksb(a,npar,indx,x);
 2309:     for (i=1;i<=npar;i++){ 
 2310:       y[i][j]=x[i];
 2311:       printf("%.3e ",y[i][j]);
 2312:       fprintf(ficlog,"%.3e ",y[i][j]);
 2313:     }
 2314:     printf("\n");
 2315:     fprintf(ficlog,"\n");
 2316:   }
 2317:   */
 2318: 
 2319:   free_matrix(a,1,npar,1,npar);
 2320:   free_matrix(y,1,npar,1,npar);
 2321:   free_vector(x,1,npar);
 2322:   free_ivector(indx,1,npar);
 2323:   free_matrix(hess,1,npar,1,npar);
 2324: 
 2325: 
 2326: }
 2327: 
 2328: /*************** hessian matrix ****************/
 2329: double hessii(double x[], double delta, int theta, double delti[], double (*func)(double []), int npar)
 2330: {
 2331:   int i;
 2332:   int l=1, lmax=20;
 2333:   double k1,k2;
 2334:   double p2[MAXPARM+1]; /* identical to x */
 2335:   double res;
 2336:   double delt=0.0001, delts, nkhi=10.,nkhif=1., khi=1.e-4;
 2337:   double fx;
 2338:   int k=0,kmax=10;
 2339:   double l1;
 2340: 
 2341:   fx=func(x);
 2342:   for (i=1;i<=npar;i++) p2[i]=x[i];
 2343:   for(l=0 ; l <=lmax; l++){  /* Enlarging the zone around the Maximum */
 2344:     l1=pow(10,l);
 2345:     delts=delt;
 2346:     for(k=1 ; k <kmax; k=k+1){
 2347:       delt = delta*(l1*k);
 2348:       p2[theta]=x[theta] +delt;
 2349:       k1=func(p2)-fx;   /* Might be negative if too close to the theoretical maximum */
 2350:       p2[theta]=x[theta]-delt;
 2351:       k2=func(p2)-fx;
 2352:       /*res= (k1-2.0*fx+k2)/delt/delt; */
 2353:       res= (k1+k2)/delt/delt/2.; /* Divided by because L and not 2*L */
 2354:       
 2355: #ifdef DEBUGHESS
 2356:       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);
 2357:       fprintf(ficlog,"%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);
 2358: #endif
 2359:       /*if(fabs(k1-2.0*fx+k2) <1.e-13){ */
 2360:       if((k1 <khi/nkhi/2.) || (k2 <khi/nkhi/2.)){
 2361: 	k=kmax;
 2362:       }
 2363:       else if((k1 >khi/nkhif) || (k2 >khi/nkhif)){ /* Keeps lastvalue before 3.84/2 KHI2 5% 1d.f. */
 2364: 	k=kmax; l=lmax*10;
 2365:       }
 2366:       else if((k1 >khi/nkhi) || (k2 >khi/nkhi)){ 
 2367: 	delts=delt;
 2368:       }
 2369:     }
 2370:   }
 2371:   delti[theta]=delts;
 2372:   return res; 
 2373:   
 2374: }
 2375: 
 2376: double hessij( double x[], double delti[], int thetai,int thetaj,double (*func)(double []),int npar)
 2377: {
 2378:   int i;
 2379:   int l=1, lmax=20;
 2380:   double k1,k2,k3,k4,res,fx;
 2381:   double p2[MAXPARM+1];
 2382:   int k;
 2383: 
 2384:   fx=func(x);
 2385:   for (k=1; k<=2; k++) {
 2386:     for (i=1;i<=npar;i++) p2[i]=x[i];
 2387:     p2[thetai]=x[thetai]+delti[thetai]/k;
 2388:     p2[thetaj]=x[thetaj]+delti[thetaj]/k;
 2389:     k1=func(p2)-fx;
 2390:   
 2391:     p2[thetai]=x[thetai]+delti[thetai]/k;
 2392:     p2[thetaj]=x[thetaj]-delti[thetaj]/k;
 2393:     k2=func(p2)-fx;
 2394:   
 2395:     p2[thetai]=x[thetai]-delti[thetai]/k;
 2396:     p2[thetaj]=x[thetaj]+delti[thetaj]/k;
 2397:     k3=func(p2)-fx;
 2398:   
 2399:     p2[thetai]=x[thetai]-delti[thetai]/k;
 2400:     p2[thetaj]=x[thetaj]-delti[thetaj]/k;
 2401:     k4=func(p2)-fx;
 2402:     res=(k1-k2-k3+k4)/4.0/delti[thetai]*k/delti[thetaj]*k/2.; /* Because of L not 2*L */
 2403: #ifdef DEBUG
 2404:     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);
 2405:     fprintf(ficlog,"%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);
 2406: #endif
 2407:   }
 2408:   return res;
 2409: }
 2410: 
 2411: /************** Inverse of matrix **************/
 2412: void ludcmp(double **a, int n, int *indx, double *d) 
 2413: { 
 2414:   int i,imax,j,k; 
 2415:   double big,dum,sum,temp; 
 2416:   double *vv; 
 2417:  
 2418:   vv=vector(1,n); 
 2419:   *d=1.0; 
 2420:   for (i=1;i<=n;i++) { 
 2421:     big=0.0; 
 2422:     for (j=1;j<=n;j++) 
 2423:       if ((temp=fabs(a[i][j])) > big) big=temp; 
 2424:     if (big == 0.0) nrerror("Singular matrix in routine ludcmp"); 
 2425:     vv[i]=1.0/big; 
 2426:   } 
 2427:   for (j=1;j<=n;j++) { 
 2428:     for (i=1;i<j;i++) { 
 2429:       sum=a[i][j]; 
 2430:       for (k=1;k<i;k++) sum -= a[i][k]*a[k][j]; 
 2431:       a[i][j]=sum; 
 2432:     } 
 2433:     big=0.0; 
 2434:     for (i=j;i<=n;i++) { 
 2435:       sum=a[i][j]; 
 2436:       for (k=1;k<j;k++) 
 2437: 	sum -= a[i][k]*a[k][j]; 
 2438:       a[i][j]=sum; 
 2439:       if ( (dum=vv[i]*fabs(sum)) >= big) { 
 2440: 	big=dum; 
 2441: 	imax=i; 
 2442:       } 
 2443:     } 
 2444:     if (j != imax) { 
 2445:       for (k=1;k<=n;k++) { 
 2446: 	dum=a[imax][k]; 
 2447: 	a[imax][k]=a[j][k]; 
 2448: 	a[j][k]=dum; 
 2449:       } 
 2450:       *d = -(*d); 
 2451:       vv[imax]=vv[j]; 
 2452:     } 
 2453:     indx[j]=imax; 
 2454:     if (a[j][j] == 0.0) a[j][j]=TINY; 
 2455:     if (j != n) { 
 2456:       dum=1.0/(a[j][j]); 
 2457:       for (i=j+1;i<=n;i++) a[i][j] *= dum; 
 2458:     } 
 2459:   } 
 2460:   free_vector(vv,1,n);  /* Doesn't work */
 2461: ;
 2462: } 
 2463: 
 2464: void lubksb(double **a, int n, int *indx, double b[]) 
 2465: { 
 2466:   int i,ii=0,ip,j; 
 2467:   double sum; 
 2468:  
 2469:   for (i=1;i<=n;i++) { 
 2470:     ip=indx[i]; 
 2471:     sum=b[ip]; 
 2472:     b[ip]=b[i]; 
 2473:     if (ii) 
 2474:       for (j=ii;j<=i-1;j++) sum -= a[i][j]*b[j]; 
 2475:     else if (sum) ii=i; 
 2476:     b[i]=sum; 
 2477:   } 
 2478:   for (i=n;i>=1;i--) { 
 2479:     sum=b[i]; 
 2480:     for (j=i+1;j<=n;j++) sum -= a[i][j]*b[j]; 
 2481:     b[i]=sum/a[i][i]; 
 2482:   } 
 2483: } 
 2484: 
 2485: void pstamp(FILE *fichier)
 2486: {
 2487:   fprintf(fichier,"# %s.%s\n#%s\n#%s\n# %s", optionfilefiname,optionfilext,version,fullversion,strstart);
 2488: }
 2489: 
 2490: /************ Frequencies ********************/
 2491: void  freqsummary(char fileres[], int iagemin, int iagemax, int **s, double **agev, int nlstate, int imx, int *Tvaraff, int **nbcode, int *ncodemax,double **mint,double **anint, char strstart[])
 2492: {  /* Some frequencies */
 2493:   
 2494:   int i, m, jk, j1, bool, z1,j;
 2495:   int first;
 2496:   double ***freq; /* Frequencies */
 2497:   double *pp, **prop;
 2498:   double pos,posprop, k2, dateintsum=0,k2cpt=0;
 2499:   char fileresp[FILENAMELENGTH];
 2500:   
 2501:   pp=vector(1,nlstate);
 2502:   prop=matrix(1,nlstate,iagemin,iagemax+3);
 2503:   strcpy(fileresp,"p");
 2504:   strcat(fileresp,fileres);
 2505:   if((ficresp=fopen(fileresp,"w"))==NULL) {
 2506:     printf("Problem with prevalence resultfile: %s\n", fileresp);
 2507:     fprintf(ficlog,"Problem with prevalence resultfile: %s\n", fileresp);
 2508:     exit(0);
 2509:   }
 2510:   freq= ma3x(-5,nlstate+ndeath,-5,nlstate+ndeath,iagemin,iagemax+3);
 2511:   j1=0;
 2512:   
 2513:   j=cptcoveff;
 2514:   if (cptcovn<1) {j=1;ncodemax[1]=1;}
 2515: 
 2516:   first=1;
 2517: 
 2518:   /* for(k1=1; k1<=j ; k1++){   /* Loop on covariates */
 2519:   /*  for(i1=1; i1<=ncodemax[k1];i1++){ /* Now it is 2 */
 2520:   /*    j1++;
 2521: */
 2522:   for (j1 = 1; j1 <= (int) pow(2,cptcoveff); j1++){
 2523:       /*printf("cptcoveff=%d Tvaraff=%d", cptcoveff,Tvaraff[1]);
 2524: 	scanf("%d", i);*/
 2525:       for (i=-5; i<=nlstate+ndeath; i++)  
 2526: 	for (jk=-5; jk<=nlstate+ndeath; jk++)  
 2527: 	  for(m=iagemin; m <= iagemax+3; m++)
 2528: 	    freq[i][jk][m]=0;
 2529:       
 2530:       for (i=1; i<=nlstate; i++)  
 2531: 	for(m=iagemin; m <= iagemax+3; m++)
 2532: 	  prop[i][m]=0;
 2533:       
 2534:       dateintsum=0;
 2535:       k2cpt=0;
 2536:       for (i=1; i<=imx; i++) {
 2537: 	bool=1;
 2538: 	if  (cptcovn>0) { /* Filter is here: Must be looked at for model=V1+V2+V3+V4 */
 2539: 	  for (z1=1; z1<=cptcoveff; z1++)       
 2540:             if (covar[Tvaraff[z1]][i]!= nbcode[Tvaraff[z1]][codtab[j1][z1]]){
 2541:                 /* Tests if the value of each of the covariates of i is equal to filter j1 */
 2542:               bool=0;
 2543:               /* printf("bool=%d i=%d, z1=%d, Tvaraff[%d]=%d, covar[Tvarff][%d]=%2f, codtab[%d][%d]=%d, nbcode[Tvaraff][codtab[%d][%d]=%d, j1=%d\n", 
 2544:                 bool,i,z1, z1, Tvaraff[z1],i,covar[Tvaraff[z1]][i],j1,z1,codtab[j1][z1],
 2545:                 j1,z1,nbcode[Tvaraff[z1]][codtab[j1][z1]],j1);*/
 2546:               /* For j1=7 in V1+V2+V3+V4 = 0 1 1 0 and codtab[7][3]=1 and nbcde[3][?]=1*/
 2547:             } 
 2548: 	}
 2549:  
 2550: 	if (bool==1){
 2551: 	  for(m=firstpass; m<=lastpass; m++){
 2552: 	    k2=anint[m][i]+(mint[m][i]/12.);
 2553: 	    /*if ((k2>=dateprev1) && (k2<=dateprev2)) {*/
 2554: 	      if(agev[m][i]==0) agev[m][i]=iagemax+1;
 2555: 	      if(agev[m][i]==1) agev[m][i]=iagemax+2;
 2556: 	      if (s[m][i]>0 && s[m][i]<=nlstate) prop[s[m][i]][(int)agev[m][i]] += weight[i];
 2557: 	      if (m<lastpass) {
 2558: 		freq[s[m][i]][s[m+1][i]][(int)agev[m][i]] += weight[i];
 2559: 		freq[s[m][i]][s[m+1][i]][iagemax+3] += weight[i];
 2560: 	      }
 2561: 	      
 2562: 	      if ((agev[m][i]>1) && (agev[m][i]< (iagemax+3))) {
 2563: 		dateintsum=dateintsum+k2;
 2564: 		k2cpt++;
 2565: 	      }
 2566: 	      /*}*/
 2567: 	  }
 2568: 	}
 2569:       } /* end i */
 2570:        
 2571:       /*      fprintf(ficresp, "#Count between %.lf/%.lf/%.lf and %.lf/%.lf/%.lf\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2);*/
 2572:       pstamp(ficresp);
 2573:       if  (cptcovn>0) {
 2574: 	fprintf(ficresp, "\n#********** Variable "); 
 2575: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresp, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 2576: 	fprintf(ficresp, "**********\n#");
 2577: 	fprintf(ficlog, "\n#********** Variable "); 
 2578: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficlog, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 2579: 	fprintf(ficlog, "**********\n#");
 2580:       }
 2581:       for(i=1; i<=nlstate;i++) 
 2582: 	fprintf(ficresp, " Age Prev(%d) N(%d) N",i,i);
 2583:       fprintf(ficresp, "\n");
 2584:       
 2585:       for(i=iagemin; i <= iagemax+3; i++){
 2586: 	if(i==iagemax+3){
 2587: 	  fprintf(ficlog,"Total");
 2588: 	}else{
 2589: 	  if(first==1){
 2590: 	    first=0;
 2591: 	    printf("See log file for details...\n");
 2592: 	  }
 2593: 	  fprintf(ficlog,"Age %d", i);
 2594: 	}
 2595: 	for(jk=1; jk <=nlstate ; jk++){
 2596: 	  for(m=-1, pp[jk]=0; m <=nlstate+ndeath ; m++)
 2597: 	    pp[jk] += freq[jk][m][i]; 
 2598: 	}
 2599: 	for(jk=1; jk <=nlstate ; jk++){
 2600: 	  for(m=-1, pos=0; m <=0 ; m++)
 2601: 	    pos += freq[jk][m][i];
 2602: 	  if(pp[jk]>=1.e-10){
 2603: 	    if(first==1){
 2604: 	      printf(" %d.=%.0f loss[%d]=%.1f%%",jk,pp[jk],jk,100*pos/pp[jk]);
 2605: 	    }
 2606: 	    fprintf(ficlog," %d.=%.0f loss[%d]=%.1f%%",jk,pp[jk],jk,100*pos/pp[jk]);
 2607: 	  }else{
 2608: 	    if(first==1)
 2609: 	      printf(" %d.=%.0f loss[%d]=NaNQ%%",jk,pp[jk],jk);
 2610: 	    fprintf(ficlog," %d.=%.0f loss[%d]=NaNQ%%",jk,pp[jk],jk);
 2611: 	  }
 2612: 	}
 2613: 
 2614: 	for(jk=1; jk <=nlstate ; jk++){
 2615: 	  for(m=0, pp[jk]=0; m <=nlstate+ndeath; m++)
 2616: 	    pp[jk] += freq[jk][m][i];
 2617: 	}	
 2618: 	for(jk=1,pos=0,posprop=0; jk <=nlstate ; jk++){
 2619: 	  pos += pp[jk];
 2620: 	  posprop += prop[jk][i];
 2621: 	}
 2622: 	for(jk=1; jk <=nlstate ; jk++){
 2623: 	  if(pos>=1.e-5){
 2624: 	    if(first==1)
 2625: 	      printf(" %d.=%.0f prev[%d]=%.1f%%",jk,pp[jk],jk,100*pp[jk]/pos);
 2626: 	    fprintf(ficlog," %d.=%.0f prev[%d]=%.1f%%",jk,pp[jk],jk,100*pp[jk]/pos);
 2627: 	  }else{
 2628: 	    if(first==1)
 2629: 	      printf(" %d.=%.0f prev[%d]=NaNQ%%",jk,pp[jk],jk);
 2630: 	    fprintf(ficlog," %d.=%.0f prev[%d]=NaNQ%%",jk,pp[jk],jk);
 2631: 	  }
 2632: 	  if( i <= iagemax){
 2633: 	    if(pos>=1.e-5){
 2634: 	      fprintf(ficresp," %d %.5f %.0f %.0f",i,prop[jk][i]/posprop, prop[jk][i],posprop);
 2635: 	      /*probs[i][jk][j1]= pp[jk]/pos;*/
 2636: 	      /*printf("\ni=%d jk=%d j1=%d %.5f %.0f %.0f %f",i,jk,j1,pp[jk]/pos, pp[jk],pos,probs[i][jk][j1]);*/
 2637: 	    }
 2638: 	    else
 2639: 	      fprintf(ficresp," %d NaNq %.0f %.0f",i,prop[jk][i],posprop);
 2640: 	  }
 2641: 	}
 2642: 	
 2643: 	for(jk=-1; jk <=nlstate+ndeath; jk++)
 2644: 	  for(m=-1; m <=nlstate+ndeath; m++)
 2645: 	    if(freq[jk][m][i] !=0 ) {
 2646: 	    if(first==1)
 2647: 	      printf(" %d%d=%.0f",jk,m,freq[jk][m][i]);
 2648: 	      fprintf(ficlog," %d%d=%.0f",jk,m,freq[jk][m][i]);
 2649: 	    }
 2650: 	if(i <= iagemax)
 2651: 	  fprintf(ficresp,"\n");
 2652: 	if(first==1)
 2653: 	  printf("Others in log...\n");
 2654: 	fprintf(ficlog,"\n");
 2655:       }
 2656:       /*}*/
 2657:   }
 2658:   dateintmean=dateintsum/k2cpt; 
 2659:  
 2660:   fclose(ficresp);
 2661:   free_ma3x(freq,-5,nlstate+ndeath,-5,nlstate+ndeath, iagemin, iagemax+3);
 2662:   free_vector(pp,1,nlstate);
 2663:   free_matrix(prop,1,nlstate,iagemin, iagemax+3);
 2664:   /* End of Freq */
 2665: }
 2666: 
 2667: /************ Prevalence ********************/
 2668: void prevalence(double ***probs, double agemin, double agemax, int **s, double **agev, int nlstate, int imx, int *Tvar, int **nbcode, int *ncodemax,double **mint,double **anint, double dateprev1,double dateprev2, int firstpass, int lastpass)
 2669: {  
 2670:   /* Compute observed prevalence between dateprev1 and dateprev2 by counting the number of people
 2671:      in each health status at the date of interview (if between dateprev1 and dateprev2).
 2672:      We still use firstpass and lastpass as another selection.
 2673:   */
 2674:  
 2675:   int i, m, jk, j1, bool, z1,j;
 2676: 
 2677:   double **prop;
 2678:   double posprop; 
 2679:   double  y2; /* in fractional years */
 2680:   int iagemin, iagemax;
 2681:   int first; /** to stop verbosity which is redirected to log file */
 2682: 
 2683:   iagemin= (int) agemin;
 2684:   iagemax= (int) agemax;
 2685:   /*pp=vector(1,nlstate);*/
 2686:   prop=matrix(1,nlstate,iagemin,iagemax+3); 
 2687:   /*  freq=ma3x(-1,nlstate+ndeath,-1,nlstate+ndeath,iagemin,iagemax+3);*/
 2688:   j1=0;
 2689:   
 2690:   /*j=cptcoveff;*/
 2691:   if (cptcovn<1) {j=1;ncodemax[1]=1;}
 2692:   
 2693:   first=1;
 2694:   for(j1=1; j1<= (int) pow(2,cptcoveff);j1++){
 2695:     /*for(i1=1; i1<=ncodemax[k1];i1++){
 2696:       j1++;*/
 2697:       
 2698:       for (i=1; i<=nlstate; i++)  
 2699: 	for(m=iagemin; m <= iagemax+3; m++)
 2700: 	  prop[i][m]=0.0;
 2701:      
 2702:       for (i=1; i<=imx; i++) { /* Each individual */
 2703: 	bool=1;
 2704: 	if  (cptcovn>0) {
 2705: 	  for (z1=1; z1<=cptcoveff; z1++) 
 2706: 	    if (covar[Tvaraff[z1]][i]!= nbcode[Tvaraff[z1]][codtab[j1][z1]]) 
 2707: 	      bool=0;
 2708: 	} 
 2709: 	if (bool==1) { 
 2710: 	  for(m=firstpass; m<=lastpass; m++){/* Other selection (we can limit to certain interviews*/
 2711: 	    y2=anint[m][i]+(mint[m][i]/12.); /* Fractional date in year */
 2712: 	    if ((y2>=dateprev1) && (y2<=dateprev2)) { /* Here is the main selection (fractional years) */
 2713: 	      if(agev[m][i]==0) agev[m][i]=iagemax+1;
 2714: 	      if(agev[m][i]==1) agev[m][i]=iagemax+2;
 2715: 	      if((int)agev[m][i] <iagemin || (int)agev[m][i] >iagemax+3) printf("Error on individual =%d agev[m][i]=%f m=%d\n",i, agev[m][i],m); 
 2716:  	      if (s[m][i]>0 && s[m][i]<=nlstate) { 
 2717: 		/*if(i>4620) printf(" i=%d m=%d s[m][i]=%d (int)agev[m][i]=%d weight[i]=%f prop=%f\n",i,m,s[m][i],(int)agev[m][m],weight[i],prop[s[m][i]][(int)agev[m][i]]);*/
 2718:  		prop[s[m][i]][(int)agev[m][i]] += weight[i];
 2719:  		prop[s[m][i]][iagemax+3] += weight[i]; 
 2720:  	      } 
 2721: 	    }
 2722: 	  } /* end selection of waves */
 2723: 	}
 2724:       }
 2725:       for(i=iagemin; i <= iagemax+3; i++){  
 2726:  	for(jk=1,posprop=0; jk <=nlstate ; jk++) { 
 2727:  	  posprop += prop[jk][i]; 
 2728:  	} 
 2729: 	
 2730:  	for(jk=1; jk <=nlstate ; jk++){	    
 2731:  	  if( i <=  iagemax){ 
 2732:  	    if(posprop>=1.e-5){ 
 2733:  	      probs[i][jk][j1]= prop[jk][i]/posprop;
 2734:  	    } else{
 2735: 	      if(first==1){
 2736: 		first=0;
 2737: 		printf("Warning Observed prevalence probs[%d][%d][%d]=%lf because of lack of cases\nSee others on log file...\n",jk,i,j1,probs[i][jk][j1]);
 2738: 	      }
 2739: 	    }
 2740:  	  } 
 2741:  	}/* end jk */ 
 2742:       }/* end i */ 
 2743:     /*} *//* end i1 */
 2744:   } /* end j1 */
 2745:   
 2746:   /*  free_ma3x(freq,-1,nlstate+ndeath,-1,nlstate+ndeath, iagemin, iagemax+3);*/
 2747:   /*free_vector(pp,1,nlstate);*/
 2748:   free_matrix(prop,1,nlstate, iagemin,iagemax+3);
 2749: }  /* End of prevalence */
 2750: 
 2751: /************* Waves Concatenation ***************/
 2752: 
 2753: void  concatwav(int wav[], int **dh, int **bh,  int **mw, int **s, double *agedc, double **agev, int  firstpass, int lastpass, int imx, int nlstate, int stepm)
 2754: {
 2755:   /* Concatenates waves: wav[i] is the number of effective (useful waves) of individual i.
 2756:      Death is a valid wave (if date is known).
 2757:      mw[mi][i] is the mi (mi=1 to wav[i])  effective wave of individual i
 2758:      dh[m][i] or dh[mw[mi][i]][i] is the delay between two effective waves m=mw[mi][i]
 2759:      and mw[mi+1][i]. dh depends on stepm.
 2760:      */
 2761: 
 2762:   int i, mi, m;
 2763:   /* int j, k=0,jk, ju, jl,jmin=1e+5, jmax=-1;
 2764:      double sum=0., jmean=0.;*/
 2765:   int first;
 2766:   int j, k=0,jk, ju, jl;
 2767:   double sum=0.;
 2768:   first=0;
 2769:   jmin=100000;
 2770:   jmax=-1;
 2771:   jmean=0.;
 2772:   for(i=1; i<=imx; i++){
 2773:     mi=0;
 2774:     m=firstpass;
 2775:     while(s[m][i] <= nlstate){
 2776:       if(s[m][i]>=1 || s[m][i]==-2 || s[m][i]==-4 || s[m][i]==-5)
 2777: 	mw[++mi][i]=m;
 2778:       if(m >=lastpass)
 2779: 	break;
 2780:       else
 2781: 	m++;
 2782:     }/* end while */
 2783:     if (s[m][i] > nlstate){
 2784:       mi++;	/* Death is another wave */
 2785:       /* if(mi==0)  never been interviewed correctly before death */
 2786: 	 /* Only death is a correct wave */
 2787:       mw[mi][i]=m;
 2788:     }
 2789: 
 2790:     wav[i]=mi;
 2791:     if(mi==0){
 2792:       nbwarn++;
 2793:       if(first==0){
 2794: 	printf("Warning! No valid information for individual %ld line=%d (skipped) and may be others, see log file\n",num[i],i);
 2795: 	first=1;
 2796:       }
 2797:       if(first==1){
 2798: 	fprintf(ficlog,"Warning! No valid information for individual %ld line=%d (skipped)\n",num[i],i);
 2799:       }
 2800:     } /* end mi==0 */
 2801:   } /* End individuals */
 2802: 
 2803:   for(i=1; i<=imx; i++){
 2804:     for(mi=1; mi<wav[i];mi++){
 2805:       if (stepm <=0)
 2806: 	dh[mi][i]=1;
 2807:       else{
 2808: 	if (s[mw[mi+1][i]][i] > nlstate) { /* A death */
 2809: 	  if (agedc[i] < 2*AGESUP) {
 2810: 	    j= rint(agedc[i]*12-agev[mw[mi][i]][i]*12); 
 2811: 	    if(j==0) j=1;  /* Survives at least one month after exam */
 2812: 	    else if(j<0){
 2813: 	      nberr++;
 2814: 	      printf("Error! Negative delay (%d to death) between waves %d and %d of individual %ld at line %d who is aged %.1f with statuses from %d to %d\n ",j,mw[mi][i],mw[mi+1][i],num[i], i,agev[mw[mi][i]][i],s[mw[mi][i]][i] ,s[mw[mi+1][i]][i]);
 2815: 	      j=1; /* Temporary Dangerous patch */
 2816: 	      printf("   We assumed that the date of interview was correct (and not the date of death) and postponed the death %d month(s) (one stepm) after the interview. You MUST fix the contradiction between dates.\n",stepm);
 2817: 	      fprintf(ficlog,"Error! Negative delay (%d to death) between waves %d and %d of individual %ld at line %d who is aged %.1f with statuses from %d to %d\n ",j,mw[mi][i],mw[mi+1][i],num[i], i,agev[mw[mi][i]][i],s[mw[mi][i]][i] ,s[mw[mi+1][i]][i]);
 2818: 	      fprintf(ficlog,"   We assumed that the date of interview was correct (and not the date of death) and postponed the death %d month(s) (one stepm) after the interview. You MUST fix the contradiction between dates.\n",stepm);
 2819: 	    }
 2820: 	    k=k+1;
 2821: 	    if (j >= jmax){
 2822: 	      jmax=j;
 2823: 	      ijmax=i;
 2824: 	    }
 2825: 	    if (j <= jmin){
 2826: 	      jmin=j;
 2827: 	      ijmin=i;
 2828: 	    }
 2829: 	    sum=sum+j;
 2830: 	    /*if (j<0) printf("j=%d num=%d \n",j,i);*/
 2831: 	    /*	  printf("%d %d %d %d\n", s[mw[mi][i]][i] ,s[mw[mi+1][i]][i],j,i);*/
 2832: 	  }
 2833: 	}
 2834: 	else{
 2835: 	  j= rint( (agev[mw[mi+1][i]][i]*12 - agev[mw[mi][i]][i]*12));
 2836: /* 	  if (j<0) printf("%d %lf %lf %d %d %d\n", i,agev[mw[mi+1][i]][i], agev[mw[mi][i]][i],j,s[mw[mi][i]][i] ,s[mw[mi+1][i]][i]); */
 2837: 
 2838: 	  k=k+1;
 2839: 	  if (j >= jmax) {
 2840: 	    jmax=j;
 2841: 	    ijmax=i;
 2842: 	  }
 2843: 	  else if (j <= jmin){
 2844: 	    jmin=j;
 2845: 	    ijmin=i;
 2846: 	  }
 2847: 	  /*	    if (j<10) printf("j=%d jmin=%d num=%d ",j,jmin,i); */
 2848: 	  /*printf("%d %lf %d %d %d\n", i,agev[mw[mi][i]][i],j,s[mw[mi][i]][i] ,s[mw[mi+1][i]][i]);*/
 2849: 	  if(j<0){
 2850: 	    nberr++;
 2851: 	    printf("Error! Negative delay (%d) between waves %d and %d of individual %ld at line %d who is aged %.1f with statuses from %d to %d\n ",j,mw[mi][i],mw[mi+1][i],num[i], i,agev[mw[mi][i]][i],s[mw[mi][i]][i] ,s[mw[mi+1][i]][i]);
 2852: 	    fprintf(ficlog,"Error! Negative delay (%d) between waves %d and %d of individual %ld at line %d who is aged %.1f with statuses from %d to %d\n ",j,mw[mi][i],mw[mi+1][i],num[i], i,agev[mw[mi][i]][i],s[mw[mi][i]][i] ,s[mw[mi+1][i]][i]);
 2853: 	  }
 2854: 	  sum=sum+j;
 2855: 	}
 2856: 	jk= j/stepm;
 2857: 	jl= j -jk*stepm;
 2858: 	ju= j -(jk+1)*stepm;
 2859: 	if(mle <=1){ /* only if we use a the linear-interpoloation pseudo-likelihood */
 2860: 	  if(jl==0){
 2861: 	    dh[mi][i]=jk;
 2862: 	    bh[mi][i]=0;
 2863: 	  }else{ /* We want a negative bias in order to only have interpolation ie
 2864: 		  * to avoid the price of an extra matrix product in likelihood */
 2865: 	    dh[mi][i]=jk+1;
 2866: 	    bh[mi][i]=ju;
 2867: 	  }
 2868: 	}else{
 2869: 	  if(jl <= -ju){
 2870: 	    dh[mi][i]=jk;
 2871: 	    bh[mi][i]=jl;	/* bias is positive if real duration
 2872: 				 * is higher than the multiple of stepm and negative otherwise.
 2873: 				 */
 2874: 	  }
 2875: 	  else{
 2876: 	    dh[mi][i]=jk+1;
 2877: 	    bh[mi][i]=ju;
 2878: 	  }
 2879: 	  if(dh[mi][i]==0){
 2880: 	    dh[mi][i]=1; /* At least one step */
 2881: 	    bh[mi][i]=ju; /* At least one step */
 2882: 	    /*  printf(" bh=%d ju=%d jl=%d dh=%d jk=%d stepm=%d %d\n",bh[mi][i],ju,jl,dh[mi][i],jk,stepm,i);*/
 2883: 	  }
 2884: 	} /* end if mle */
 2885:       }
 2886:     } /* end wave */
 2887:   }
 2888:   jmean=sum/k;
 2889:   printf("Delay (in months) between two waves Min=%d (for indiviudal %ld) Max=%d (%ld) Mean=%f\n\n ",jmin, num[ijmin], jmax, num[ijmax], jmean);
 2890:   fprintf(ficlog,"Delay (in months) between two waves Min=%d (for indiviudal %d) Max=%d (%d) Mean=%f\n\n ",jmin, ijmin, jmax, ijmax, jmean);
 2891:  }
 2892: 
 2893: /*********** Tricode ****************************/
 2894: void tricode(int *Tvar, int **nbcode, int imx, int *Ndum)
 2895: {
 2896:   /**< Uses cptcovn+2*cptcovprod as the number of covariates */
 2897:   /*	  Tvar[i]=atoi(stre);  find 'n' in Vn and stores in Tvar. If model=V2+V1 Tvar[1]=2 and Tvar[2]=1 
 2898:   /* Boring subroutine which should only output nbcode[Tvar[j]][k]
 2899:    * Tvar[5] in V2+V1+V3*age+V2*V4 is 2 (V2)
 2900:   /* nbcode[Tvar[j]][1]= 
 2901:   */
 2902: 
 2903:   int ij=1, k=0, j=0, i=0, maxncov=NCOVMAX;
 2904:   int modmaxcovj=0; /* Modality max of covariates j */
 2905:   int cptcode=0; /* Modality max of covariates j */
 2906:   int modmincovj=0; /* Modality min of covariates j */
 2907: 
 2908: 
 2909:   cptcoveff=0; 
 2910:  
 2911:   for (k=-1; k < maxncov; k++) Ndum[k]=0;
 2912:   for (k=1; k <= maxncov; k++) ncodemax[k]=0; /* Horrible constant again replaced by NCOVMAX */
 2913: 
 2914:   /* Loop on covariates without age and products */
 2915:   for (j=1; j<=(cptcovs); j++) { /* model V1 + V2*age+ V3 + V3*V4 : V1 + V3 = 2 only */
 2916:     for (i=1; i<=imx; i++) { /* Lopp on individuals: reads the data file to get the maximum value of the 
 2917: 			       modality of this covariate Vj*/ 
 2918:       ij=(int)(covar[Tvar[j]][i]); /* ij=0 or 1 or -1. Value of the covariate Tvar[j] for individual i
 2919: 				    * If product of Vn*Vm, still boolean *:
 2920: 				    * If it was coded 1, 2, 3, 4 should be splitted into 3 boolean variables
 2921: 				    * 1 => 0 0 0, 2 => 0 0 1, 3 => 0 1 1, 4=1 0 0   */
 2922:       /* Finds for covariate j, n=Tvar[j] of Vn . ij is the
 2923: 				      modality of the nth covariate of individual i. */
 2924:       if (ij > modmaxcovj)
 2925:         modmaxcovj=ij; 
 2926:       else if (ij < modmincovj) 
 2927: 	modmincovj=ij; 
 2928:       if ((ij < -1) && (ij > NCOVMAX)){
 2929: 	printf( "Error: minimal is less than -1 or maximal is bigger than %d. Exiting. \n", NCOVMAX );
 2930: 	exit(1);
 2931:       }else
 2932:       Ndum[ij]++; /*counts and stores the occurence of this modality 0, 1, -1*/
 2933:       /*  If coded 1, 2, 3 , counts the number of 1 Ndum[1], number of 2, Ndum[2], etc */
 2934:       /*printf("i=%d ij=%d Ndum[ij]=%d imx=%d",i,ij,Ndum[ij],imx);*/
 2935:       /* getting the maximum value of the modality of the covariate
 2936: 	 (should be 0 or 1 now) Tvar[j]. If V=sex and male is coded 0 and
 2937: 	 female is 1, then modmaxcovj=1.*/
 2938:     }
 2939:     printf(" Minimal and maximal values of %d th covariate V%d: min=%d max=%d \n", j, Tvar[j], modmincovj, modmaxcovj);
 2940:     cptcode=modmaxcovj;
 2941:     /* Ndum[0] = frequency of 0 for model-covariate j, Ndum[1] frequency of 1 etc. */
 2942:    /*for (i=0; i<=cptcode; i++) {*/
 2943:     for (i=modmincovj;  i<=modmaxcovj; i++) { /* i=-1 ? 0 and 1*//* For each value of the modality of model-cov j */
 2944:       printf("Frequencies of covariates %d V%d %d\n", j, Tvar[j], Ndum[i]);
 2945:       if( Ndum[i] != 0 ){ /* Counts if nobody answered, empty modality */
 2946: 	ncodemax[j]++;  /* ncodemax[j]= Number of non-null modalities of the j th covariate. */
 2947:       }
 2948:       /* In fact  ncodemax[j]=2 (dichotom. variables only) but it could be more for
 2949: 	 historical reasons: 3 if coded 1, 2, 3 and 4 and Ndum[2]=0 */
 2950:     } /* Ndum[-1] number of undefined modalities */
 2951: 
 2952:     /* j is a covariate, n=Tvar[j] of Vn; Fills nbcode */
 2953:     /* For covariate j, modalities could be 1, 2, 3, 4. If Ndum[2]=0 ncodemax[j] is not 4 but 3 */
 2954:     /* If Ndum[3}= 635; Ndum[4]=0; Ndum[5]=0; Ndum[6]=27; Ndum[7]=125;
 2955:        modmincovj=3; modmaxcovj = 7;
 2956:        There are only 3 modalities non empty (or 2 if 27 is too few) : ncodemax[j]=3;
 2957:        which will be coded 0, 1, 2 which in binary on 3-1 digits are 0=00 1=01, 2=10; defining two dummy 
 2958:        variables V1_1 and V1_2.
 2959:        nbcode[Tvar[j]][ij]=k;
 2960:        nbcode[Tvar[j]][1]=0;
 2961:        nbcode[Tvar[j]][2]=1;
 2962:        nbcode[Tvar[j]][3]=2;
 2963:     */
 2964:     ij=1; /* ij is similar to i but can jumps over null modalities */
 2965:     for (i=modmincovj; i<=modmaxcovj; i++) { /* i= 1 to 2 for dichotomous, or from 1 to 3 */
 2966:       for (k=0; k<= cptcode; k++) { /* k=-1 ? k=0 to 1 *//* Could be 1 to 4 */
 2967: 	/*recode from 0 */
 2968: 	if (Ndum[k] != 0) { /* If at least one individual responded to this modality k */
 2969: 	  nbcode[Tvar[j]][ij]=k;  /* stores the modality in an array nbcode. 
 2970: 				     k is a modality. If we have model=V1+V1*sex 
 2971: 				     then: nbcode[1][1]=0 ; nbcode[1][2]=1; nbcode[2][1]=0 ; nbcode[2][2]=1; */
 2972: 	  ij++;
 2973: 	}
 2974: 	if (ij > ncodemax[j]) break; 
 2975:       }  /* end of loop on */
 2976:     } /* end of loop on modality */ 
 2977:   } /* end of loop on model-covariate j. nbcode[Tvarj][1]=0 and nbcode[Tvarj][2]=1 sets the value of covariate j*/  
 2978:   
 2979:  for (k=-1; k< maxncov; k++) Ndum[k]=0; 
 2980:   
 2981:   for (i=1; i<=ncovmodel-2; i++) { /* -2, cste and age */ 
 2982:    /* Listing of all covariables in statement model to see if some covariates appear twice. For example, V1 appears twice in V1+V1*V2.*/ 
 2983:    ij=Tvar[i]; /* Tvar might be -1 if status was unknown */ 
 2984:    Ndum[ij]++; 
 2985:  } 
 2986: 
 2987:  ij=1;
 2988:  for (i=0; i<=  maxncov-1; i++) { /* modmaxcovj is unknown here. Only Ndum[2(V2),3(age*V3), 5(V3*V2) 6(V1*V4) */
 2989:    /*printf("Ndum[%d]=%d\n",i, Ndum[i]);*/
 2990:    if((Ndum[i]!=0) && (i<=ncovcol)){
 2991:      /*printf("diff Ndum[%d]=%d\n",i, Ndum[i]);*/
 2992:      Tvaraff[ij]=i; /*For printing (unclear) */
 2993:      ij++;
 2994:    }else
 2995:        Tvaraff[ij]=0;
 2996:  }
 2997:  ij--;
 2998:  cptcoveff=ij; /*Number of total covariates*/
 2999: 
 3000: }
 3001: 
 3002: 
 3003: /*********** Health Expectancies ****************/
 3004: 
 3005: void evsij(double ***eij, double x[], int nlstate, int stepm, int bage, int fage, double **oldm, double **savm, int cij, int estepm,char strstart[] )
 3006: 
 3007: {
 3008:   /* Health expectancies, no variances */
 3009:   int i, j, nhstepm, hstepm, h, nstepm;
 3010:   int nhstepma, nstepma; /* Decreasing with age */
 3011:   double age, agelim, hf;
 3012:   double ***p3mat;
 3013:   double eip;
 3014: 
 3015:   pstamp(ficreseij);
 3016:   fprintf(ficreseij,"# (a) Life expectancies by health status at initial age and (b) health expectancies by health status at initial age\n");
 3017:   fprintf(ficreseij,"# Age");
 3018:   for(i=1; i<=nlstate;i++){
 3019:     for(j=1; j<=nlstate;j++){
 3020:       fprintf(ficreseij," e%1d%1d ",i,j);
 3021:     }
 3022:     fprintf(ficreseij," e%1d. ",i);
 3023:   }
 3024:   fprintf(ficreseij,"\n");
 3025: 
 3026:   
 3027:   if(estepm < stepm){
 3028:     printf ("Problem %d lower than %d\n",estepm, stepm);
 3029:   }
 3030:   else  hstepm=estepm;   
 3031:   /* We compute the life expectancy from trapezoids spaced every estepm months
 3032:    * This is mainly to measure the difference between two models: for example
 3033:    * if stepm=24 months pijx are given only every 2 years and by summing them
 3034:    * we are calculating an estimate of the Life Expectancy assuming a linear 
 3035:    * progression in between and thus overestimating or underestimating according
 3036:    * to the curvature of the survival function. If, for the same date, we 
 3037:    * estimate the model with stepm=1 month, we can keep estepm to 24 months
 3038:    * to compare the new estimate of Life expectancy with the same linear 
 3039:    * hypothesis. A more precise result, taking into account a more precise
 3040:    * curvature will be obtained if estepm is as small as stepm. */
 3041: 
 3042:   /* For example we decided to compute the life expectancy with the smallest unit */
 3043:   /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm. 
 3044:      nhstepm is the number of hstepm from age to agelim 
 3045:      nstepm is the number of stepm from age to agelin. 
 3046:      Look at hpijx to understand the reason of that which relies in memory size
 3047:      and note for a fixed period like estepm months */
 3048:   /* We decided (b) to get a life expectancy respecting the most precise curvature of the
 3049:      survival function given by stepm (the optimization length). Unfortunately it
 3050:      means that if the survival funtion is printed only each two years of age and if
 3051:      you sum them up and add 1 year (area under the trapezoids) you won't get the same 
 3052:      results. So we changed our mind and took the option of the best precision.
 3053:   */
 3054:   hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */ 
 3055: 
 3056:   agelim=AGESUP;
 3057:   /* If stepm=6 months */
 3058:     /* Computed by stepm unit matrices, product of hstepm matrices, stored
 3059:        in an array of nhstepm length: nhstepm=10, hstepm=4, stepm=6 months */
 3060:     
 3061: /* nhstepm age range expressed in number of stepm */
 3062:   nstepm=(int) rint((agelim-bage)*YEARM/stepm); /* Biggest nstepm */
 3063:   /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3064:   /* if (stepm >= YEARM) hstepm=1;*/
 3065:   nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
 3066:   p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3067: 
 3068:   for (age=bage; age<=fage; age ++){ 
 3069:     nstepma=(int) rint((agelim-bage)*YEARM/stepm); /* Biggest nstepm */
 3070:     /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3071:     /* if (stepm >= YEARM) hstepm=1;*/
 3072:     nhstepma = nstepma/hstepm;/* Expressed in hstepm, typically nhstepma=40/4=10 */
 3073: 
 3074:     /* If stepm=6 months */
 3075:     /* Computed by stepm unit matrices, product of hstepma matrices, stored
 3076:        in an array of nhstepma length: nhstepma=10, hstepm=4, stepm=6 months */
 3077:     
 3078:     hpxij(p3mat,nhstepma,age,hstepm,x,nlstate,stepm,oldm, savm, cij);  
 3079:     
 3080:     hf=hstepm*stepm/YEARM;  /* Duration of hstepm expressed in year unit. */
 3081:     
 3082:     printf("%d|",(int)age);fflush(stdout);
 3083:     fprintf(ficlog,"%d|",(int)age);fflush(ficlog);
 3084:     
 3085:     /* Computing expectancies */
 3086:     for(i=1; i<=nlstate;i++)
 3087:       for(j=1; j<=nlstate;j++)
 3088: 	for (h=0, eij[i][j][(int)age]=0; h<=nhstepm-1; h++){
 3089: 	  eij[i][j][(int)age] += (p3mat[i][j][h]+p3mat[i][j][h+1])/2.0*hf;
 3090: 	  
 3091: 	  /* if((int)age==70)printf("i=%2d,j=%2d,h=%2d,age=%3d,%9.4f,%9.4f,%9.4f\n",i,j,h,(int)age,p3mat[i][j][h],hf,eij[i][j][(int)age]);*/
 3092: 
 3093: 	}
 3094: 
 3095:     fprintf(ficreseij,"%3.0f",age );
 3096:     for(i=1; i<=nlstate;i++){
 3097:       eip=0;
 3098:       for(j=1; j<=nlstate;j++){
 3099: 	eip +=eij[i][j][(int)age];
 3100: 	fprintf(ficreseij,"%9.4f", eij[i][j][(int)age] );
 3101:       }
 3102:       fprintf(ficreseij,"%9.4f", eip );
 3103:     }
 3104:     fprintf(ficreseij,"\n");
 3105:     
 3106:   }
 3107:   free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3108:   printf("\n");
 3109:   fprintf(ficlog,"\n");
 3110:   
 3111: }
 3112: 
 3113: void cvevsij(double ***eij, double x[], int nlstate, int stepm, int bage, int fage, double **oldm, double **savm, int cij, int estepm,double delti[],double **matcov,char strstart[] )
 3114: 
 3115: {
 3116:   /* Covariances of health expectancies eij and of total life expectancies according
 3117:    to initial status i, ei. .
 3118:   */
 3119:   int i, j, nhstepm, hstepm, h, nstepm, k, cptj, cptj2, i2, j2, ij, ji;
 3120:   int nhstepma, nstepma; /* Decreasing with age */
 3121:   double age, agelim, hf;
 3122:   double ***p3matp, ***p3matm, ***varhe;
 3123:   double **dnewm,**doldm;
 3124:   double *xp, *xm;
 3125:   double **gp, **gm;
 3126:   double ***gradg, ***trgradg;
 3127:   int theta;
 3128: 
 3129:   double eip, vip;
 3130: 
 3131:   varhe=ma3x(1,nlstate*nlstate,1,nlstate*nlstate,(int) bage, (int) fage);
 3132:   xp=vector(1,npar);
 3133:   xm=vector(1,npar);
 3134:   dnewm=matrix(1,nlstate*nlstate,1,npar);
 3135:   doldm=matrix(1,nlstate*nlstate,1,nlstate*nlstate);
 3136:   
 3137:   pstamp(ficresstdeij);
 3138:   fprintf(ficresstdeij,"# Health expectancies with standard errors\n");
 3139:   fprintf(ficresstdeij,"# Age");
 3140:   for(i=1; i<=nlstate;i++){
 3141:     for(j=1; j<=nlstate;j++)
 3142:       fprintf(ficresstdeij," e%1d%1d (SE)",i,j);
 3143:     fprintf(ficresstdeij," e%1d. ",i);
 3144:   }
 3145:   fprintf(ficresstdeij,"\n");
 3146: 
 3147:   pstamp(ficrescveij);
 3148:   fprintf(ficrescveij,"# Subdiagonal matrix of covariances of health expectancies by age: cov(eij,ekl)\n");
 3149:   fprintf(ficrescveij,"# Age");
 3150:   for(i=1; i<=nlstate;i++)
 3151:     for(j=1; j<=nlstate;j++){
 3152:       cptj= (j-1)*nlstate+i;
 3153:       for(i2=1; i2<=nlstate;i2++)
 3154: 	for(j2=1; j2<=nlstate;j2++){
 3155: 	  cptj2= (j2-1)*nlstate+i2;
 3156: 	  if(cptj2 <= cptj)
 3157: 	    fprintf(ficrescveij,"  %1d%1d,%1d%1d",i,j,i2,j2);
 3158: 	}
 3159:     }
 3160:   fprintf(ficrescveij,"\n");
 3161:   
 3162:   if(estepm < stepm){
 3163:     printf ("Problem %d lower than %d\n",estepm, stepm);
 3164:   }
 3165:   else  hstepm=estepm;   
 3166:   /* We compute the life expectancy from trapezoids spaced every estepm months
 3167:    * This is mainly to measure the difference between two models: for example
 3168:    * if stepm=24 months pijx are given only every 2 years and by summing them
 3169:    * we are calculating an estimate of the Life Expectancy assuming a linear 
 3170:    * progression in between and thus overestimating or underestimating according
 3171:    * to the curvature of the survival function. If, for the same date, we 
 3172:    * estimate the model with stepm=1 month, we can keep estepm to 24 months
 3173:    * to compare the new estimate of Life expectancy with the same linear 
 3174:    * hypothesis. A more precise result, taking into account a more precise
 3175:    * curvature will be obtained if estepm is as small as stepm. */
 3176: 
 3177:   /* For example we decided to compute the life expectancy with the smallest unit */
 3178:   /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm. 
 3179:      nhstepm is the number of hstepm from age to agelim 
 3180:      nstepm is the number of stepm from age to agelin. 
 3181:      Look at hpijx to understand the reason of that which relies in memory size
 3182:      and note for a fixed period like estepm months */
 3183:   /* We decided (b) to get a life expectancy respecting the most precise curvature of the
 3184:      survival function given by stepm (the optimization length). Unfortunately it
 3185:      means that if the survival funtion is printed only each two years of age and if
 3186:      you sum them up and add 1 year (area under the trapezoids) you won't get the same 
 3187:      results. So we changed our mind and took the option of the best precision.
 3188:   */
 3189:   hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */ 
 3190: 
 3191:   /* If stepm=6 months */
 3192:   /* nhstepm age range expressed in number of stepm */
 3193:   agelim=AGESUP;
 3194:   nstepm=(int) rint((agelim-bage)*YEARM/stepm); 
 3195:   /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3196:   /* if (stepm >= YEARM) hstepm=1;*/
 3197:   nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
 3198:   
 3199:   p3matp=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3200:   p3matm=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3201:   gradg=ma3x(0,nhstepm,1,npar,1,nlstate*nlstate);
 3202:   trgradg =ma3x(0,nhstepm,1,nlstate*nlstate,1,npar);
 3203:   gp=matrix(0,nhstepm,1,nlstate*nlstate);
 3204:   gm=matrix(0,nhstepm,1,nlstate*nlstate);
 3205: 
 3206:   for (age=bage; age<=fage; age ++){ 
 3207:     nstepma=(int) rint((agelim-bage)*YEARM/stepm); /* Biggest nstepm */
 3208:     /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3209:     /* if (stepm >= YEARM) hstepm=1;*/
 3210:     nhstepma = nstepma/hstepm;/* Expressed in hstepm, typically nhstepma=40/4=10 */
 3211: 
 3212:     /* If stepm=6 months */
 3213:     /* Computed by stepm unit matrices, product of hstepma matrices, stored
 3214:        in an array of nhstepma length: nhstepma=10, hstepm=4, stepm=6 months */
 3215:     
 3216:     hf=hstepm*stepm/YEARM;  /* Duration of hstepm expressed in year unit. */
 3217: 
 3218:     /* Computing  Variances of health expectancies */
 3219:     /* Gradient is computed with plus gp and minus gm. Code is duplicated in order to
 3220:        decrease memory allocation */
 3221:     for(theta=1; theta <=npar; theta++){
 3222:       for(i=1; i<=npar; i++){ 
 3223: 	xp[i] = x[i] + (i==theta ?delti[theta]:0);
 3224: 	xm[i] = x[i] - (i==theta ?delti[theta]:0);
 3225:       }
 3226:       hpxij(p3matp,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, cij);  
 3227:       hpxij(p3matm,nhstepm,age,hstepm,xm,nlstate,stepm,oldm,savm, cij);  
 3228:   
 3229:       for(j=1; j<= nlstate; j++){
 3230: 	for(i=1; i<=nlstate; i++){
 3231: 	  for(h=0; h<=nhstepm-1; h++){
 3232: 	    gp[h][(j-1)*nlstate + i] = (p3matp[i][j][h]+p3matp[i][j][h+1])/2.;
 3233: 	    gm[h][(j-1)*nlstate + i] = (p3matm[i][j][h]+p3matm[i][j][h+1])/2.;
 3234: 	  }
 3235: 	}
 3236:       }
 3237:      
 3238:       for(ij=1; ij<= nlstate*nlstate; ij++)
 3239: 	for(h=0; h<=nhstepm-1; h++){
 3240: 	  gradg[h][theta][ij]= (gp[h][ij]-gm[h][ij])/2./delti[theta];
 3241: 	}
 3242:     }/* End theta */
 3243:     
 3244:     
 3245:     for(h=0; h<=nhstepm-1; h++)
 3246:       for(j=1; j<=nlstate*nlstate;j++)
 3247: 	for(theta=1; theta <=npar; theta++)
 3248: 	  trgradg[h][j][theta]=gradg[h][theta][j];
 3249:     
 3250: 
 3251:      for(ij=1;ij<=nlstate*nlstate;ij++)
 3252:       for(ji=1;ji<=nlstate*nlstate;ji++)
 3253: 	varhe[ij][ji][(int)age] =0.;
 3254: 
 3255:      printf("%d|",(int)age);fflush(stdout);
 3256:      fprintf(ficlog,"%d|",(int)age);fflush(ficlog);
 3257:      for(h=0;h<=nhstepm-1;h++){
 3258:       for(k=0;k<=nhstepm-1;k++){
 3259: 	matprod2(dnewm,trgradg[h],1,nlstate*nlstate,1,npar,1,npar,matcov);
 3260: 	matprod2(doldm,dnewm,1,nlstate*nlstate,1,npar,1,nlstate*nlstate,gradg[k]);
 3261: 	for(ij=1;ij<=nlstate*nlstate;ij++)
 3262: 	  for(ji=1;ji<=nlstate*nlstate;ji++)
 3263: 	    varhe[ij][ji][(int)age] += doldm[ij][ji]*hf*hf;
 3264:       }
 3265:     }
 3266: 
 3267:     /* Computing expectancies */
 3268:     hpxij(p3matm,nhstepm,age,hstepm,x,nlstate,stepm,oldm, savm, cij);  
 3269:     for(i=1; i<=nlstate;i++)
 3270:       for(j=1; j<=nlstate;j++)
 3271: 	for (h=0, eij[i][j][(int)age]=0; h<=nhstepm-1; h++){
 3272: 	  eij[i][j][(int)age] += (p3matm[i][j][h]+p3matm[i][j][h+1])/2.0*hf;
 3273: 	  
 3274: 	  /* if((int)age==70)printf("i=%2d,j=%2d,h=%2d,age=%3d,%9.4f,%9.4f,%9.4f\n",i,j,h,(int)age,p3mat[i][j][h],hf,eij[i][j][(int)age]);*/
 3275: 
 3276: 	}
 3277: 
 3278:     fprintf(ficresstdeij,"%3.0f",age );
 3279:     for(i=1; i<=nlstate;i++){
 3280:       eip=0.;
 3281:       vip=0.;
 3282:       for(j=1; j<=nlstate;j++){
 3283: 	eip += eij[i][j][(int)age];
 3284: 	for(k=1; k<=nlstate;k++) /* Sum on j and k of cov(eij,eik) */
 3285: 	  vip += varhe[(j-1)*nlstate+i][(k-1)*nlstate+i][(int)age];
 3286: 	fprintf(ficresstdeij," %9.4f (%.4f)", eij[i][j][(int)age], sqrt(varhe[(j-1)*nlstate+i][(j-1)*nlstate+i][(int)age]) );
 3287:       }
 3288:       fprintf(ficresstdeij," %9.4f (%.4f)", eip, sqrt(vip));
 3289:     }
 3290:     fprintf(ficresstdeij,"\n");
 3291: 
 3292:     fprintf(ficrescveij,"%3.0f",age );
 3293:     for(i=1; i<=nlstate;i++)
 3294:       for(j=1; j<=nlstate;j++){
 3295: 	cptj= (j-1)*nlstate+i;
 3296: 	for(i2=1; i2<=nlstate;i2++)
 3297: 	  for(j2=1; j2<=nlstate;j2++){
 3298: 	    cptj2= (j2-1)*nlstate+i2;
 3299: 	    if(cptj2 <= cptj)
 3300: 	      fprintf(ficrescveij," %.4f", varhe[cptj][cptj2][(int)age]);
 3301: 	  }
 3302:       }
 3303:     fprintf(ficrescveij,"\n");
 3304:    
 3305:   }
 3306:   free_matrix(gm,0,nhstepm,1,nlstate*nlstate);
 3307:   free_matrix(gp,0,nhstepm,1,nlstate*nlstate);
 3308:   free_ma3x(gradg,0,nhstepm,1,npar,1,nlstate*nlstate);
 3309:   free_ma3x(trgradg,0,nhstepm,1,nlstate*nlstate,1,npar);
 3310:   free_ma3x(p3matm,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3311:   free_ma3x(p3matp,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3312:   printf("\n");
 3313:   fprintf(ficlog,"\n");
 3314: 
 3315:   free_vector(xm,1,npar);
 3316:   free_vector(xp,1,npar);
 3317:   free_matrix(dnewm,1,nlstate*nlstate,1,npar);
 3318:   free_matrix(doldm,1,nlstate*nlstate,1,nlstate*nlstate);
 3319:   free_ma3x(varhe,1,nlstate*nlstate,1,nlstate*nlstate,(int) bage, (int)fage);
 3320: }
 3321: 
 3322: /************ Variance ******************/
 3323: void varevsij(char optionfilefiname[], 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, int estepm, int cptcov, int cptcod, int popbased, int mobilav, char strstart[])
 3324: {
 3325:   /* Variance of health expectancies */
 3326:   /*  double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double ** savm,double ftolpl);*/
 3327:   /* double **newm;*/
 3328:   double **dnewm,**doldm;
 3329:   double **dnewmp,**doldmp;
 3330:   int i, j, nhstepm, hstepm, h, nstepm ;
 3331:   int k;
 3332:   double *xp;
 3333:   double **gp, **gm;  /* for var eij */
 3334:   double ***gradg, ***trgradg; /*for var eij */
 3335:   double **gradgp, **trgradgp; /* for var p point j */
 3336:   double *gpp, *gmp; /* for var p point j */
 3337:   double **varppt; /* for var p point j nlstate to nlstate+ndeath */
 3338:   double ***p3mat;
 3339:   double age,agelim, hf;
 3340:   double ***mobaverage;
 3341:   int theta;
 3342:   char digit[4];
 3343:   char digitp[25];
 3344: 
 3345:   char fileresprobmorprev[FILENAMELENGTH];
 3346: 
 3347:   if(popbased==1){
 3348:     if(mobilav!=0)
 3349:       strcpy(digitp,"-populbased-mobilav-");
 3350:     else strcpy(digitp,"-populbased-nomobil-");
 3351:   }
 3352:   else 
 3353:     strcpy(digitp,"-stablbased-");
 3354: 
 3355:   if (mobilav!=0) {
 3356:     mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 3357:     if (movingaverage(probs, bage, fage, mobaverage,mobilav)!=0){
 3358:       fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 3359:       printf(" Error in movingaverage mobilav=%d\n",mobilav);
 3360:     }
 3361:   }
 3362: 
 3363:   strcpy(fileresprobmorprev,"prmorprev"); 
 3364:   sprintf(digit,"%-d",ij);
 3365:   /*printf("DIGIT=%s, ij=%d ijr=%-d|\n",digit, ij,ij);*/
 3366:   strcat(fileresprobmorprev,digit); /* Tvar to be done */
 3367:   strcat(fileresprobmorprev,digitp); /* Popbased or not, mobilav or not */
 3368:   strcat(fileresprobmorprev,fileres);
 3369:   if((ficresprobmorprev=fopen(fileresprobmorprev,"w"))==NULL) {
 3370:     printf("Problem with resultfile: %s\n", fileresprobmorprev);
 3371:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprobmorprev);
 3372:   }
 3373:   printf("Computing total mortality p.j=w1*p1j+w2*p2j+..: result on file '%s' \n",fileresprobmorprev);
 3374:  
 3375:   fprintf(ficlog,"Computing total mortality p.j=w1*p1j+w2*p2j+..: result on file '%s' \n",fileresprobmorprev);
 3376:   pstamp(ficresprobmorprev);
 3377:   fprintf(ficresprobmorprev,"# probabilities of dying before estepm=%d months for people of exact age and weighted probabilities w1*p1j+w2*p2j+... stand dev in()\n",estepm);
 3378:   fprintf(ficresprobmorprev,"# Age cov=%-d",ij);
 3379:   for(j=nlstate+1; j<=(nlstate+ndeath);j++){
 3380:     fprintf(ficresprobmorprev," p.%-d SE",j);
 3381:     for(i=1; i<=nlstate;i++)
 3382:       fprintf(ficresprobmorprev," w%1d p%-d%-d",i,i,j);
 3383:   }  
 3384:   fprintf(ficresprobmorprev,"\n");
 3385:   fprintf(ficgp,"\n# Routine varevsij");
 3386:   /* fprintf(fichtm, "#Local time at start: %s", strstart);*/
 3387:   fprintf(fichtm,"\n<li><h4> Computing probabilities of dying over estepm months as a weighted average (i.e global mortality independent of initial healh state)</h4></li>\n");
 3388:   fprintf(fichtm,"\n<br>%s  <br>\n",digitp);
 3389: /*   } */
 3390:   varppt = matrix(nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3391:   pstamp(ficresvij);
 3392:   fprintf(ficresvij,"# Variance and covariance of health expectancies e.j \n#  (weighted average of eij where weights are ");
 3393:   if(popbased==1)
 3394:     fprintf(ficresvij,"the age specific prevalence observed (cross-sectionally) in the population i.e cross-sectionally\n in each health state (popbased=1) (mobilav=%d\n",mobilav);
 3395:   else
 3396:     fprintf(ficresvij,"the age specific period (stable) prevalences in each health state \n");
 3397:   fprintf(ficresvij,"# Age");
 3398:   for(i=1; i<=nlstate;i++)
 3399:     for(j=1; j<=nlstate;j++)
 3400:       fprintf(ficresvij," Cov(e.%1d, e.%1d)",i,j);
 3401:   fprintf(ficresvij,"\n");
 3402: 
 3403:   xp=vector(1,npar);
 3404:   dnewm=matrix(1,nlstate,1,npar);
 3405:   doldm=matrix(1,nlstate,1,nlstate);
 3406:   dnewmp= matrix(nlstate+1,nlstate+ndeath,1,npar);
 3407:   doldmp= matrix(nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3408: 
 3409:   gradgp=matrix(1,npar,nlstate+1,nlstate+ndeath);
 3410:   gpp=vector(nlstate+1,nlstate+ndeath);
 3411:   gmp=vector(nlstate+1,nlstate+ndeath);
 3412:   trgradgp =matrix(nlstate+1,nlstate+ndeath,1,npar); /* mu or p point j*/
 3413:   
 3414:   if(estepm < stepm){
 3415:     printf ("Problem %d lower than %d\n",estepm, stepm);
 3416:   }
 3417:   else  hstepm=estepm;   
 3418:   /* For example we decided to compute the life expectancy with the smallest unit */
 3419:   /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm. 
 3420:      nhstepm is the number of hstepm from age to agelim 
 3421:      nstepm is the number of stepm from age to agelin. 
 3422:      Look at function hpijx to understand why (it is linked to memory size questions) */
 3423:   /* We decided (b) to get a life expectancy respecting the most precise curvature of the
 3424:      survival function given by stepm (the optimization length). Unfortunately it
 3425:      means that if the survival funtion is printed every two years of age and if
 3426:      you sum them up and add 1 year (area under the trapezoids) you won't get the same 
 3427:      results. So we changed our mind and took the option of the best precision.
 3428:   */
 3429:   hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */ 
 3430:   agelim = AGESUP;
 3431:   for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
 3432:     nstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */ 
 3433:     nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
 3434:     p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3435:     gradg=ma3x(0,nhstepm,1,npar,1,nlstate);
 3436:     gp=matrix(0,nhstepm,1,nlstate);
 3437:     gm=matrix(0,nhstepm,1,nlstate);
 3438: 
 3439: 
 3440:     for(theta=1; theta <=npar; theta++){
 3441:       for(i=1; i<=npar; i++){ /* Computes gradient x + delta*/
 3442: 	xp[i] = x[i] + (i==theta ?delti[theta]:0);
 3443:       }
 3444:       hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);  
 3445:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3446: 
 3447:       if (popbased==1) {
 3448: 	if(mobilav ==0){
 3449: 	  for(i=1; i<=nlstate;i++)
 3450: 	    prlim[i][i]=probs[(int)age][i][ij];
 3451: 	}else{ /* mobilav */ 
 3452: 	  for(i=1; i<=nlstate;i++)
 3453: 	    prlim[i][i]=mobaverage[(int)age][i][ij];
 3454: 	}
 3455:       }
 3456:   
 3457:       for(j=1; j<= nlstate; j++){
 3458: 	for(h=0; h<=nhstepm; h++){
 3459: 	  for(i=1, gp[h][j]=0.;i<=nlstate;i++)
 3460: 	    gp[h][j] += prlim[i][i]*p3mat[i][j][h];
 3461: 	}
 3462:       }
 3463:       /* This for computing probability of death (h=1 means
 3464:          computed over hstepm matrices product = hstepm*stepm months) 
 3465:          as a weighted average of prlim.
 3466:       */
 3467:       for(j=nlstate+1;j<=nlstate+ndeath;j++){
 3468: 	for(i=1,gpp[j]=0.; i<= nlstate; i++)
 3469: 	  gpp[j] += prlim[i][i]*p3mat[i][j][1];
 3470:       }    
 3471:       /* end probability of death */
 3472: 
 3473:       for(i=1; i<=npar; i++) /* Computes gradient x - delta */
 3474: 	xp[i] = x[i] - (i==theta ?delti[theta]:0);
 3475:       hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);  
 3476:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3477:  
 3478:       if (popbased==1) {
 3479: 	if(mobilav ==0){
 3480: 	  for(i=1; i<=nlstate;i++)
 3481: 	    prlim[i][i]=probs[(int)age][i][ij];
 3482: 	}else{ /* mobilav */ 
 3483: 	  for(i=1; i<=nlstate;i++)
 3484: 	    prlim[i][i]=mobaverage[(int)age][i][ij];
 3485: 	}
 3486:       }
 3487: 
 3488:       for(j=1; j<= nlstate; j++){  /* Sum of wi * eij = e.j */
 3489: 	for(h=0; h<=nhstepm; h++){
 3490: 	  for(i=1, gm[h][j]=0.;i<=nlstate;i++)
 3491: 	    gm[h][j] += prlim[i][i]*p3mat[i][j][h];
 3492: 	}
 3493:       }
 3494:       /* This for computing probability of death (h=1 means
 3495:          computed over hstepm matrices product = hstepm*stepm months) 
 3496:          as a weighted average of prlim.
 3497:       */
 3498:       for(j=nlstate+1;j<=nlstate+ndeath;j++){
 3499: 	for(i=1,gmp[j]=0.; i<= nlstate; i++)
 3500:          gmp[j] += prlim[i][i]*p3mat[i][j][1];
 3501:       }    
 3502:       /* end probability of death */
 3503: 
 3504:       for(j=1; j<= nlstate; j++) /* vareij */
 3505: 	for(h=0; h<=nhstepm; h++){
 3506: 	  gradg[h][theta][j]= (gp[h][j]-gm[h][j])/2./delti[theta];
 3507: 	}
 3508: 
 3509:       for(j=nlstate+1; j<= nlstate+ndeath; j++){ /* var mu */
 3510: 	gradgp[theta][j]= (gpp[j]-gmp[j])/2./delti[theta];
 3511:       }
 3512: 
 3513:     } /* End theta */
 3514: 
 3515:     trgradg =ma3x(0,nhstepm,1,nlstate,1,npar); /* veij */
 3516: 
 3517:     for(h=0; h<=nhstepm; h++) /* veij */
 3518:       for(j=1; j<=nlstate;j++)
 3519: 	for(theta=1; theta <=npar; theta++)
 3520: 	  trgradg[h][j][theta]=gradg[h][theta][j];
 3521: 
 3522:     for(j=nlstate+1; j<=nlstate+ndeath;j++) /* mu */
 3523:       for(theta=1; theta <=npar; theta++)
 3524: 	trgradgp[j][theta]=gradgp[theta][j];
 3525:   
 3526: 
 3527:     hf=hstepm*stepm/YEARM;  /* Duration of hstepm expressed in year unit. */
 3528:     for(i=1;i<=nlstate;i++)
 3529:       for(j=1;j<=nlstate;j++)
 3530: 	vareij[i][j][(int)age] =0.;
 3531: 
 3532:     for(h=0;h<=nhstepm;h++){
 3533:       for(k=0;k<=nhstepm;k++){
 3534: 	matprod2(dnewm,trgradg[h],1,nlstate,1,npar,1,npar,matcov);
 3535: 	matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg[k]);
 3536: 	for(i=1;i<=nlstate;i++)
 3537: 	  for(j=1;j<=nlstate;j++)
 3538: 	    vareij[i][j][(int)age] += doldm[i][j]*hf*hf;
 3539:       }
 3540:     }
 3541:   
 3542:     /* pptj */
 3543:     matprod2(dnewmp,trgradgp,nlstate+1,nlstate+ndeath,1,npar,1,npar,matcov);
 3544:     matprod2(doldmp,dnewmp,nlstate+1,nlstate+ndeath,1,npar,nlstate+1,nlstate+ndeath,gradgp);
 3545:     for(j=nlstate+1;j<=nlstate+ndeath;j++)
 3546:       for(i=nlstate+1;i<=nlstate+ndeath;i++)
 3547: 	varppt[j][i]=doldmp[j][i];
 3548:     /* end ppptj */
 3549:     /*  x centered again */
 3550:     hpxij(p3mat,nhstepm,age,hstepm,x,nlstate,stepm,oldm,savm, ij);  
 3551:     prevalim(prlim,nlstate,x,age,oldm,savm,ftolpl,ij);
 3552:  
 3553:     if (popbased==1) {
 3554:       if(mobilav ==0){
 3555: 	for(i=1; i<=nlstate;i++)
 3556: 	  prlim[i][i]=probs[(int)age][i][ij];
 3557:       }else{ /* mobilav */ 
 3558: 	for(i=1; i<=nlstate;i++)
 3559: 	  prlim[i][i]=mobaverage[(int)age][i][ij];
 3560:       }
 3561:     }
 3562:              
 3563:     /* This for computing probability of death (h=1 means
 3564:        computed over hstepm (estepm) matrices product = hstepm*stepm months) 
 3565:        as a weighted average of prlim.
 3566:     */
 3567:     for(j=nlstate+1;j<=nlstate+ndeath;j++){
 3568:       for(i=1,gmp[j]=0.;i<= nlstate; i++) 
 3569: 	gmp[j] += prlim[i][i]*p3mat[i][j][1]; 
 3570:     }    
 3571:     /* end probability of death */
 3572: 
 3573:     fprintf(ficresprobmorprev,"%3d %d ",(int) age, ij);
 3574:     for(j=nlstate+1; j<=(nlstate+ndeath);j++){
 3575:       fprintf(ficresprobmorprev," %11.3e %11.3e",gmp[j], sqrt(varppt[j][j]));
 3576:       for(i=1; i<=nlstate;i++){
 3577: 	fprintf(ficresprobmorprev," %11.3e %11.3e ",prlim[i][i],p3mat[i][j][1]);
 3578:       }
 3579:     } 
 3580:     fprintf(ficresprobmorprev,"\n");
 3581: 
 3582:     fprintf(ficresvij,"%.0f ",age );
 3583:     for(i=1; i<=nlstate;i++)
 3584:       for(j=1; j<=nlstate;j++){
 3585: 	fprintf(ficresvij," %.4f", vareij[i][j][(int)age]);
 3586:       }
 3587:     fprintf(ficresvij,"\n");
 3588:     free_matrix(gp,0,nhstepm,1,nlstate);
 3589:     free_matrix(gm,0,nhstepm,1,nlstate);
 3590:     free_ma3x(gradg,0,nhstepm,1,npar,1,nlstate);
 3591:     free_ma3x(trgradg,0,nhstepm,1,nlstate,1,npar);
 3592:     free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3593:   } /* End age */
 3594:   free_vector(gpp,nlstate+1,nlstate+ndeath);
 3595:   free_vector(gmp,nlstate+1,nlstate+ndeath);
 3596:   free_matrix(gradgp,1,npar,nlstate+1,nlstate+ndeath);
 3597:   free_matrix(trgradgp,nlstate+1,nlstate+ndeath,1,npar); /* mu or p point j*/
 3598:   fprintf(ficgp,"\nunset parametric;unset label; set ter png small size 320, 240");
 3599:   /* for(j=nlstate+1; j<= nlstate+ndeath; j++){ *//* Only the first actually */
 3600:   fprintf(ficgp,"\n set log y; unset log x;set xlabel \"Age\"; set ylabel \"Force of mortality (year-1)\";");
 3601: /*   fprintf(ficgp,"\n plot \"%s\"  u 1:($3*%6.3f) not w l 1 ",fileresprobmorprev,YEARM/estepm); */
 3602: /*   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3+1.96*$4)*%6.3f) t \"95\%% interval\" w l 2 ",fileresprobmorprev,YEARM/estepm); */
 3603: /*   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3-1.96*$4)*%6.3f) not w l 2 ",fileresprobmorprev,YEARM/estepm); */
 3604:   fprintf(ficgp,"\n plot \"%s\"  u 1:($3) not w l lt 1 ",subdirf(fileresprobmorprev));
 3605:   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3+1.96*$4)) t \"95\%% interval\" w l lt 2 ",subdirf(fileresprobmorprev));
 3606:   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3-1.96*$4)) not w l lt 2 ",subdirf(fileresprobmorprev));
 3607:   fprintf(fichtm,"\n<br> File (multiple files are possible if covariates are present): <A href=\"%s\">%s</a>\n",subdirf(fileresprobmorprev),subdirf(fileresprobmorprev));
 3608:   fprintf(fichtm,"\n<br> Probability is computed over estepm=%d months. <br> <img src=\"%s%s.png\"> <br>\n", estepm,subdirf3(optionfilefiname,"varmuptjgr",digitp),digit);
 3609:   /*  fprintf(fichtm,"\n<br> Probability is computed over estepm=%d months and then divided by estepm and multiplied by %.0f in order to have the probability to die over a year <br> <img src=\"varmuptjgr%s%s.png\"> <br>\n", stepm,YEARM,digitp,digit);
 3610: */
 3611: /*   fprintf(ficgp,"\nset out \"varmuptjgr%s%s%s.png\";replot;",digitp,optionfilefiname,digit); */
 3612:   fprintf(ficgp,"\nset out \"%s%s.png\";replot;\n",subdirf3(optionfilefiname,"varmuptjgr",digitp),digit);
 3613: 
 3614:   free_vector(xp,1,npar);
 3615:   free_matrix(doldm,1,nlstate,1,nlstate);
 3616:   free_matrix(dnewm,1,nlstate,1,npar);
 3617:   free_matrix(doldmp,nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3618:   free_matrix(dnewmp,nlstate+1,nlstate+ndeath,1,npar);
 3619:   free_matrix(varppt,nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3620:   if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 3621:   fclose(ficresprobmorprev);
 3622:   fflush(ficgp);
 3623:   fflush(fichtm); 
 3624: }  /* end varevsij */
 3625: 
 3626: /************ Variance of prevlim ******************/
 3627: 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, char strstart[])
 3628: {
 3629:   /* Variance of prevalence limit */
 3630:   /*  double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double **savm,double ftolpl);*/
 3631: 
 3632:   double **dnewm,**doldm;
 3633:   int i, j, nhstepm, hstepm;
 3634:   double *xp;
 3635:   double *gp, *gm;
 3636:   double **gradg, **trgradg;
 3637:   double age,agelim;
 3638:   int theta;
 3639:   
 3640:   pstamp(ficresvpl);
 3641:   fprintf(ficresvpl,"# Standard deviation of period (stable) prevalences \n");
 3642:   fprintf(ficresvpl,"# Age");
 3643:   for(i=1; i<=nlstate;i++)
 3644:       fprintf(ficresvpl," %1d-%1d",i,i);
 3645:   fprintf(ficresvpl,"\n");
 3646: 
 3647:   xp=vector(1,npar);
 3648:   dnewm=matrix(1,nlstate,1,npar);
 3649:   doldm=matrix(1,nlstate,1,nlstate);
 3650:   
 3651:   hstepm=1*YEARM; /* Every year of age */
 3652:   hstepm=hstepm/stepm; /* Typically in stepm units, if j= 2 years, = 2/6 months = 4 */ 
 3653:   agelim = AGESUP;
 3654:   for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
 3655:     nhstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */ 
 3656:     if (stepm >= YEARM) hstepm=1;
 3657:     nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
 3658:     gradg=matrix(1,npar,1,nlstate);
 3659:     gp=vector(1,nlstate);
 3660:     gm=vector(1,nlstate);
 3661: 
 3662:     for(theta=1; theta <=npar; theta++){
 3663:       for(i=1; i<=npar; i++){ /* Computes gradient */
 3664: 	xp[i] = x[i] + (i==theta ?delti[theta]:0);
 3665:       }
 3666:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3667:       for(i=1;i<=nlstate;i++)
 3668: 	gp[i] = prlim[i][i];
 3669:     
 3670:       for(i=1; i<=npar; i++) /* Computes gradient */
 3671: 	xp[i] = x[i] - (i==theta ?delti[theta]:0);
 3672:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3673:       for(i=1;i<=nlstate;i++)
 3674: 	gm[i] = prlim[i][i];
 3675: 
 3676:       for(i=1;i<=nlstate;i++)
 3677: 	gradg[theta][i]= (gp[i]-gm[i])/2./delti[theta];
 3678:     } /* End theta */
 3679: 
 3680:     trgradg =matrix(1,nlstate,1,npar);
 3681: 
 3682:     for(j=1; j<=nlstate;j++)
 3683:       for(theta=1; theta <=npar; theta++)
 3684: 	trgradg[j][theta]=gradg[theta][j];
 3685: 
 3686:     for(i=1;i<=nlstate;i++)
 3687:       varpl[i][(int)age] =0.;
 3688:     matprod2(dnewm,trgradg,1,nlstate,1,npar,1,npar,matcov);
 3689:     matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg);
 3690:     for(i=1;i<=nlstate;i++)
 3691:       varpl[i][(int)age] = doldm[i][i]; /* Covariances are useless */
 3692: 
 3693:     fprintf(ficresvpl,"%.0f ",age );
 3694:     for(i=1; i<=nlstate;i++)
 3695:       fprintf(ficresvpl," %.5f (%.5f)",prlim[i][i],sqrt(varpl[i][(int)age]));
 3696:     fprintf(ficresvpl,"\n");
 3697:     free_vector(gp,1,nlstate);
 3698:     free_vector(gm,1,nlstate);
 3699:     free_matrix(gradg,1,npar,1,nlstate);
 3700:     free_matrix(trgradg,1,nlstate,1,npar);
 3701:   } /* End age */
 3702: 
 3703:   free_vector(xp,1,npar);
 3704:   free_matrix(doldm,1,nlstate,1,npar);
 3705:   free_matrix(dnewm,1,nlstate,1,nlstate);
 3706: 
 3707: }
 3708: 
 3709: /************ Variance of one-step probabilities  ******************/
 3710: void varprob(char optionfilefiname[], double **matcov, double x[], double delti[], int nlstate, double bage, double fage, int ij, int *Tvar, int **nbcode, int *ncodemax, char strstart[])
 3711: {
 3712:   int i, j=0,  k1, l1, tj;
 3713:   int k2, l2, j1,  z1;
 3714:   int k=0, l;
 3715:   int first=1, first1, first2;
 3716:   double cv12, mu1, mu2, lc1, lc2, v12, v21, v11, v22,v1,v2, c12, tnalp;
 3717:   double **dnewm,**doldm;
 3718:   double *xp;
 3719:   double *gp, *gm;
 3720:   double **gradg, **trgradg;
 3721:   double **mu;
 3722:   double age, cov[NCOVMAX+1];
 3723:   double std=2.0; /* Number of standard deviation wide of confidence ellipsoids */
 3724:   int theta;
 3725:   char fileresprob[FILENAMELENGTH];
 3726:   char fileresprobcov[FILENAMELENGTH];
 3727:   char fileresprobcor[FILENAMELENGTH];
 3728:   double ***varpij;
 3729: 
 3730:   strcpy(fileresprob,"prob"); 
 3731:   strcat(fileresprob,fileres);
 3732:   if((ficresprob=fopen(fileresprob,"w"))==NULL) {
 3733:     printf("Problem with resultfile: %s\n", fileresprob);
 3734:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprob);
 3735:   }
 3736:   strcpy(fileresprobcov,"probcov"); 
 3737:   strcat(fileresprobcov,fileres);
 3738:   if((ficresprobcov=fopen(fileresprobcov,"w"))==NULL) {
 3739:     printf("Problem with resultfile: %s\n", fileresprobcov);
 3740:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprobcov);
 3741:   }
 3742:   strcpy(fileresprobcor,"probcor"); 
 3743:   strcat(fileresprobcor,fileres);
 3744:   if((ficresprobcor=fopen(fileresprobcor,"w"))==NULL) {
 3745:     printf("Problem with resultfile: %s\n", fileresprobcor);
 3746:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprobcor);
 3747:   }
 3748:   printf("Computing standard deviation of one-step probabilities: result on file '%s' \n",fileresprob);
 3749:   fprintf(ficlog,"Computing standard deviation of one-step probabilities: result on file '%s' \n",fileresprob);
 3750:   printf("Computing matrix of variance covariance of one-step probabilities: result on file '%s' \n",fileresprobcov);
 3751:   fprintf(ficlog,"Computing matrix of variance covariance of one-step probabilities: result on file '%s' \n",fileresprobcov);
 3752:   printf("and correlation matrix of one-step probabilities: result on file '%s' \n",fileresprobcor);
 3753:   fprintf(ficlog,"and correlation matrix of one-step probabilities: result on file '%s' \n",fileresprobcor);
 3754:   pstamp(ficresprob);
 3755:   fprintf(ficresprob,"#One-step probabilities and stand. devi in ()\n");
 3756:   fprintf(ficresprob,"# Age");
 3757:   pstamp(ficresprobcov);
 3758:   fprintf(ficresprobcov,"#One-step probabilities and covariance matrix\n");
 3759:   fprintf(ficresprobcov,"# Age");
 3760:   pstamp(ficresprobcor);
 3761:   fprintf(ficresprobcor,"#One-step probabilities and correlation matrix\n");
 3762:   fprintf(ficresprobcor,"# Age");
 3763: 
 3764: 
 3765:   for(i=1; i<=nlstate;i++)
 3766:     for(j=1; j<=(nlstate+ndeath);j++){
 3767:       fprintf(ficresprob," p%1d-%1d (SE)",i,j);
 3768:       fprintf(ficresprobcov," p%1d-%1d ",i,j);
 3769:       fprintf(ficresprobcor," p%1d-%1d ",i,j);
 3770:     }  
 3771:  /* fprintf(ficresprob,"\n");
 3772:   fprintf(ficresprobcov,"\n");
 3773:   fprintf(ficresprobcor,"\n");
 3774:  */
 3775:   xp=vector(1,npar);
 3776:   dnewm=matrix(1,(nlstate)*(nlstate+ndeath),1,npar);
 3777:   doldm=matrix(1,(nlstate)*(nlstate+ndeath),1,(nlstate)*(nlstate+ndeath));
 3778:   mu=matrix(1,(nlstate)*(nlstate+ndeath), (int) bage, (int)fage);
 3779:   varpij=ma3x(1,nlstate*(nlstate+ndeath),1,nlstate*(nlstate+ndeath),(int) bage, (int) fage);
 3780:   first=1;
 3781:   fprintf(ficgp,"\n# Routine varprob");
 3782:   fprintf(fichtm,"\n<li><h4> Computing and drawing one step probabilities with their confidence intervals</h4></li>\n");
 3783:   fprintf(fichtm,"\n");
 3784: 
 3785:   fprintf(fichtm,"\n<li><h4> <a href=\"%s\">Matrix of variance-covariance of pairs of step probabilities (drawings)</a></h4></li>\n",optionfilehtmcov);
 3786:   fprintf(fichtmcov,"\n<h4>Matrix of variance-covariance of pairs of step probabilities</h4>\n\
 3787:   file %s<br>\n",optionfilehtmcov);
 3788:   fprintf(fichtmcov,"\nEllipsoids of confidence centered on point (p<inf>ij</inf>, p<inf>kl</inf>) are estimated\
 3789: and drawn. It helps understanding how is the covariance between two incidences.\
 3790:  They are expressed in year<sup>-1</sup> in order to be less dependent of stepm.<br>\n");
 3791:   fprintf(fichtmcov,"\n<br> Contour plot corresponding to x'cov<sup>-1</sup>x = 4 (where x is the column vector (pij,pkl)) are drawn. \
 3792: It can be understood this way: if pij and pkl where uncorrelated the (2x2) matrix of covariance \
 3793: would have been (1/(var pij), 0 , 0, 1/(var pkl)), and the confidence interval would be 2 \
 3794: standard deviations wide on each axis. <br>\
 3795:  Now, if both incidences are correlated (usual case) we diagonalised the inverse of the covariance matrix\
 3796:  and made the appropriate rotation to look at the uncorrelated principal directions.<br>\
 3797: To be simple, these graphs help to understand the significativity of each parameter in relation to a second other one.<br> \n");
 3798: 
 3799:   cov[1]=1;
 3800:   /* tj=cptcoveff; */
 3801:   tj = (int) pow(2,cptcoveff);
 3802:   if (cptcovn<1) {tj=1;ncodemax[1]=1;}
 3803:   j1=0;
 3804:   for(j1=1; j1<=tj;j1++){
 3805:     /*for(i1=1; i1<=ncodemax[t];i1++){ */
 3806:     /*j1++;*/
 3807:       if  (cptcovn>0) {
 3808: 	fprintf(ficresprob, "\n#********** Variable "); 
 3809: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprob, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3810: 	fprintf(ficresprob, "**********\n#\n");
 3811: 	fprintf(ficresprobcov, "\n#********** Variable "); 
 3812: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprobcov, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3813: 	fprintf(ficresprobcov, "**********\n#\n");
 3814: 	
 3815: 	fprintf(ficgp, "\n#********** Variable "); 
 3816: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficgp, " V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3817: 	fprintf(ficgp, "**********\n#\n");
 3818: 	
 3819: 	
 3820: 	fprintf(fichtmcov, "\n<hr  size=\"2\" color=\"#EC5E5E\">********** Variable "); 
 3821: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(fichtm, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3822: 	fprintf(fichtmcov, "**********\n<hr size=\"2\" color=\"#EC5E5E\">");
 3823: 	
 3824: 	fprintf(ficresprobcor, "\n#********** Variable ");    
 3825: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprobcor, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3826: 	fprintf(ficresprobcor, "**********\n#");    
 3827:       }
 3828:       
 3829:       gradg=matrix(1,npar,1,(nlstate)*(nlstate+ndeath));
 3830:       trgradg=matrix(1,(nlstate)*(nlstate+ndeath),1,npar);
 3831:       gp=vector(1,(nlstate)*(nlstate+ndeath));
 3832:       gm=vector(1,(nlstate)*(nlstate+ndeath));
 3833:       for (age=bage; age<=fage; age ++){ 
 3834: 	cov[2]=age;
 3835: 	for (k=1; k<=cptcovn;k++) {
 3836: 	  cov[2+k]=nbcode[Tvar[k]][codtab[j1][Tvar[k]]];/* j1 1 2 3 4
 3837: 							 * 1  1 1 1 1
 3838: 							 * 2  2 1 1 1
 3839: 							 * 3  1 2 1 1
 3840: 							 */
 3841: 	  /* nbcode[1][1]=0 nbcode[1][2]=1;*/
 3842: 	}
 3843: 	for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
 3844: 	for (k=1; k<=cptcovprod;k++)
 3845: 	  cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
 3846: 	
 3847:     
 3848: 	for(theta=1; theta <=npar; theta++){
 3849: 	  for(i=1; i<=npar; i++)
 3850: 	    xp[i] = x[i] + (i==theta ?delti[theta]:(double)0);
 3851: 	  
 3852: 	  pmij(pmmij,cov,ncovmodel,xp,nlstate);
 3853: 	  
 3854: 	  k=0;
 3855: 	  for(i=1; i<= (nlstate); i++){
 3856: 	    for(j=1; j<=(nlstate+ndeath);j++){
 3857: 	      k=k+1;
 3858: 	      gp[k]=pmmij[i][j];
 3859: 	    }
 3860: 	  }
 3861: 	  
 3862: 	  for(i=1; i<=npar; i++)
 3863: 	    xp[i] = x[i] - (i==theta ?delti[theta]:(double)0);
 3864:     
 3865: 	  pmij(pmmij,cov,ncovmodel,xp,nlstate);
 3866: 	  k=0;
 3867: 	  for(i=1; i<=(nlstate); i++){
 3868: 	    for(j=1; j<=(nlstate+ndeath);j++){
 3869: 	      k=k+1;
 3870: 	      gm[k]=pmmij[i][j];
 3871: 	    }
 3872: 	  }
 3873:      
 3874: 	  for(i=1; i<= (nlstate)*(nlstate+ndeath); i++) 
 3875: 	    gradg[theta][i]=(gp[i]-gm[i])/(double)2./delti[theta];  
 3876: 	}
 3877: 
 3878: 	for(j=1; j<=(nlstate)*(nlstate+ndeath);j++)
 3879: 	  for(theta=1; theta <=npar; theta++)
 3880: 	    trgradg[j][theta]=gradg[theta][j];
 3881: 	
 3882: 	matprod2(dnewm,trgradg,1,(nlstate)*(nlstate+ndeath),1,npar,1,npar,matcov); 
 3883: 	matprod2(doldm,dnewm,1,(nlstate)*(nlstate+ndeath),1,npar,1,(nlstate)*(nlstate+ndeath),gradg);
 3884: 
 3885: 	pmij(pmmij,cov,ncovmodel,x,nlstate);
 3886: 	
 3887: 	k=0;
 3888: 	for(i=1; i<=(nlstate); i++){
 3889: 	  for(j=1; j<=(nlstate+ndeath);j++){
 3890: 	    k=k+1;
 3891: 	    mu[k][(int) age]=pmmij[i][j];
 3892: 	  }
 3893: 	}
 3894:      	for(i=1;i<=(nlstate)*(nlstate+ndeath);i++)
 3895: 	  for(j=1;j<=(nlstate)*(nlstate+ndeath);j++)
 3896: 	    varpij[i][j][(int)age] = doldm[i][j];
 3897: 
 3898: 	/*printf("\n%d ",(int)age);
 3899: 	  for (i=1; i<=(nlstate)*(nlstate+ndeath);i++){
 3900: 	  printf("%e [%e ;%e] ",gm[i],gm[i]-2*sqrt(doldm[i][i]),gm[i]+2*sqrt(doldm[i][i]));
 3901: 	  fprintf(ficlog,"%e [%e ;%e] ",gm[i],gm[i]-2*sqrt(doldm[i][i]),gm[i]+2*sqrt(doldm[i][i]));
 3902: 	  }*/
 3903: 
 3904: 	fprintf(ficresprob,"\n%d ",(int)age);
 3905: 	fprintf(ficresprobcov,"\n%d ",(int)age);
 3906: 	fprintf(ficresprobcor,"\n%d ",(int)age);
 3907: 
 3908: 	for (i=1; i<=(nlstate)*(nlstate+ndeath);i++)
 3909: 	  fprintf(ficresprob,"%11.3e (%11.3e) ",mu[i][(int) age],sqrt(varpij[i][i][(int)age]));
 3910: 	for (i=1; i<=(nlstate)*(nlstate+ndeath);i++){
 3911: 	  fprintf(ficresprobcov,"%11.3e ",mu[i][(int) age]);
 3912: 	  fprintf(ficresprobcor,"%11.3e ",mu[i][(int) age]);
 3913: 	}
 3914: 	i=0;
 3915: 	for (k=1; k<=(nlstate);k++){
 3916:  	  for (l=1; l<=(nlstate+ndeath);l++){ 
 3917:  	    i++;
 3918: 	    fprintf(ficresprobcov,"\n%d %d-%d",(int)age,k,l);
 3919: 	    fprintf(ficresprobcor,"\n%d %d-%d",(int)age,k,l);
 3920: 	    for (j=1; j<=i;j++){
 3921: 	      /* printf(" k=%d l=%d i=%d j=%d\n",k,l,i,j);fflush(stdout); */
 3922: 	      fprintf(ficresprobcov," %11.3e",varpij[i][j][(int)age]);
 3923: 	      fprintf(ficresprobcor," %11.3e",varpij[i][j][(int) age]/sqrt(varpij[i][i][(int) age])/sqrt(varpij[j][j][(int)age]));
 3924: 	    }
 3925: 	  }
 3926: 	}/* end of loop for state */
 3927:       } /* end of loop for age */
 3928:       free_vector(gp,1,(nlstate+ndeath)*(nlstate+ndeath));
 3929:       free_vector(gm,1,(nlstate+ndeath)*(nlstate+ndeath));
 3930:       free_matrix(trgradg,1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
 3931:       free_matrix(gradg,1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
 3932:       
 3933:       /* Confidence intervalle of pij  */
 3934:       /*
 3935: 	fprintf(ficgp,"\nunset parametric;unset label");
 3936: 	fprintf(ficgp,"\nset log y;unset log x; set xlabel \"Age\";set ylabel \"probability (year-1)\"");
 3937: 	fprintf(ficgp,"\nset ter png small\nset size 0.65,0.65");
 3938: 	fprintf(fichtm,"\n<br>Probability with  confidence intervals expressed in year<sup>-1</sup> :<a href=\"pijgr%s.png\">pijgr%s.png</A>, ",optionfilefiname,optionfilefiname);
 3939: 	fprintf(fichtm,"\n<br><img src=\"pijgr%s.png\"> ",optionfilefiname);
 3940: 	fprintf(ficgp,"\nset out \"pijgr%s.png\"",optionfilefiname);
 3941: 	fprintf(ficgp,"\nplot \"%s\" every :::%d::%d u 1:2 \"\%%lf",k1,k2,xfilevarprob);
 3942:       */
 3943: 
 3944:       /* Drawing ellipsoids of confidence of two variables p(k1-l1,k2-l2)*/
 3945:       first1=1;first2=2;
 3946:       for (k2=1; k2<=(nlstate);k2++){
 3947: 	for (l2=1; l2<=(nlstate+ndeath);l2++){ 
 3948: 	  if(l2==k2) continue;
 3949: 	  j=(k2-1)*(nlstate+ndeath)+l2;
 3950: 	  for (k1=1; k1<=(nlstate);k1++){
 3951: 	    for (l1=1; l1<=(nlstate+ndeath);l1++){ 
 3952: 	      if(l1==k1) continue;
 3953: 	      i=(k1-1)*(nlstate+ndeath)+l1;
 3954: 	      if(i<=j) continue;
 3955: 	      for (age=bage; age<=fage; age ++){ 
 3956: 		if ((int)age %5==0){
 3957: 		  v1=varpij[i][i][(int)age]/stepm*YEARM/stepm*YEARM;
 3958: 		  v2=varpij[j][j][(int)age]/stepm*YEARM/stepm*YEARM;
 3959: 		  cv12=varpij[i][j][(int)age]/stepm*YEARM/stepm*YEARM;
 3960: 		  mu1=mu[i][(int) age]/stepm*YEARM ;
 3961: 		  mu2=mu[j][(int) age]/stepm*YEARM;
 3962: 		  c12=cv12/sqrt(v1*v2);
 3963: 		  /* Computing eigen value of matrix of covariance */
 3964: 		  lc1=((v1+v2)+sqrt((v1+v2)*(v1+v2) - 4*(v1*v2-cv12*cv12)))/2.;
 3965: 		  lc2=((v1+v2)-sqrt((v1+v2)*(v1+v2) - 4*(v1*v2-cv12*cv12)))/2.;
 3966: 		  if ((lc2 <0) || (lc1 <0) ){
 3967: 		    if(first2==1){
 3968: 		      first1=0;
 3969: 		    printf("Strange: j1=%d One eigen value of 2x2 matrix of covariance is negative, lc1=%11.3e, lc2=%11.3e, v1=%11.3e, v2=%11.3e, cv12=%11.3e.\n It means that the matrix was not well estimated (varpij), for i=%2d, j=%2d, age=%4d .\n See files %s and %s. Probably WRONG RESULTS. See log file for details...\n", j1, lc1, lc2, v1, v2, cv12, i, j, (int)age,fileresprobcov, fileresprobcor);
 3970: 		    }
 3971: 		    fprintf(ficlog,"Strange: j1=%d One eigen value of 2x2 matrix of covariance is negative, lc1=%11.3e, lc2=%11.3e, v1=%11.3e, v2=%11.3e, cv12=%11.3e.\n It means that the matrix was not well estimated (varpij), for i=%2d, j=%2d, age=%4d .\n See files %s and %s. Probably WRONG RESULTS.\n", j1, lc1, lc2, v1, v2, cv12, i, j, (int)age,fileresprobcov, fileresprobcor);fflush(ficlog);
 3972: 		    /* lc1=fabs(lc1); */ /* If we want to have them positive */
 3973: 		    /* lc2=fabs(lc2); */
 3974: 		  }
 3975: 
 3976: 		  /* Eigen vectors */
 3977: 		  v11=(1./sqrt(1+(v1-lc1)*(v1-lc1)/cv12/cv12));
 3978: 		  /*v21=sqrt(1.-v11*v11); *//* error */
 3979: 		  v21=(lc1-v1)/cv12*v11;
 3980: 		  v12=-v21;
 3981: 		  v22=v11;
 3982: 		  tnalp=v21/v11;
 3983: 		  if(first1==1){
 3984: 		    first1=0;
 3985: 		    printf("%d %d%d-%d%d mu %.4e %.4e Var %.4e %.4e cor %.3f cov %.4e Eig %.3e %.3e 1stv %.3f %.3f tang %.3f\nOthers in log...\n",(int) age,k1,l1,k2,l2,mu1,mu2,v1,v2,c12,cv12,lc1,lc2,v11,v21,tnalp);
 3986: 		  }
 3987: 		  fprintf(ficlog,"%d %d%d-%d%d mu %.4e %.4e Var %.4e %.4e cor %.3f cov %.4e Eig %.3e %.3e 1stv %.3f %.3f tan %.3f\n",(int) age,k1,l1,k2,l2,mu1,mu2,v1,v2,c12,cv12,lc1,lc2,v11,v21,tnalp);
 3988: 		  /*printf(fignu*/
 3989: 		  /* mu1+ v11*lc1*cost + v12*lc2*sin(t) */
 3990: 		  /* mu2+ v21*lc1*cost + v22*lc2*sin(t) */
 3991: 		  if(first==1){
 3992: 		    first=0;
 3993:  		    fprintf(ficgp,"\nset parametric;unset label");
 3994: 		    fprintf(ficgp,"\nset log y;set log x; set xlabel \"p%1d%1d (year-1)\";set ylabel \"p%1d%1d (year-1)\"",k1,l1,k2,l2);
 3995: 		    fprintf(ficgp,"\nset ter png small size 320, 240");
 3996: 		    fprintf(fichtmcov,"\n<br>Ellipsoids of confidence cov(p%1d%1d,p%1d%1d) expressed in year<sup>-1</sup>\
 3997:  :<a href=\"%s%d%1d%1d-%1d%1d.png\">\
 3998: %s%d%1d%1d-%1d%1d.png</A>, ",k1,l1,k2,l2,\
 3999: 			    subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2,\
 4000: 			    subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4001: 		    fprintf(fichtmcov,"\n<br><img src=\"%s%d%1d%1d-%1d%1d.png\"> ",subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4002: 		    fprintf(fichtmcov,"\n<br> Correlation at age %d (%.3f),",(int) age, c12);
 4003: 		    fprintf(ficgp,"\nset out \"%s%d%1d%1d-%1d%1d.png\"",subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4004: 		    fprintf(ficgp,"\nset label \"%d\" at %11.3e,%11.3e center",(int) age, mu1,mu2);
 4005: 		    fprintf(ficgp,"\n# Age %d, p%1d%1d - p%1d%1d",(int) age, k1,l1,k2,l2);
 4006: 		    fprintf(ficgp,"\nplot [-pi:pi] %11.3e+ %.3f*(%11.3e*%11.3e*cos(t)+%11.3e*%11.3e*sin(t)), %11.3e +%.3f*(%11.3e*%11.3e*cos(t)+%11.3e*%11.3e*sin(t)) not",\
 4007: 			    mu1,std,v11,sqrt(lc1),v12,sqrt(lc2),\
 4008: 			    mu2,std,v21,sqrt(lc1),v22,sqrt(lc2));
 4009: 		  }else{
 4010: 		    first=0;
 4011: 		    fprintf(fichtmcov," %d (%.3f),",(int) age, c12);
 4012: 		    fprintf(ficgp,"\n# Age %d, p%1d%1d - p%1d%1d",(int) age, k1,l1,k2,l2);
 4013: 		    fprintf(ficgp,"\nset label \"%d\" at %11.3e,%11.3e center",(int) age, mu1,mu2);
 4014: 		    fprintf(ficgp,"\nreplot %11.3e+ %.3f*(%11.3e*%11.3e*cos(t)+%11.3e*%11.3e*sin(t)), %11.3e +%.3f*(%11.3e*%11.3e*cos(t)+%11.3e*%11.3e*sin(t)) not",\
 4015: 			    mu1,std,v11,sqrt(lc1),v12,sqrt(lc2),\
 4016: 			    mu2,std,v21,sqrt(lc1),v22,sqrt(lc2));
 4017: 		  }/* if first */
 4018: 		} /* age mod 5 */
 4019: 	      } /* end loop age */
 4020: 	      fprintf(ficgp,"\nset out \"%s%d%1d%1d-%1d%1d.png\";replot;",subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4021: 	      first=1;
 4022: 	    } /*l12 */
 4023: 	  } /* k12 */
 4024: 	} /*l1 */
 4025:       }/* k1 */
 4026:       /* } /* loop covariates */
 4027:   }
 4028:   free_ma3x(varpij,1,nlstate,1,nlstate+ndeath,(int) bage, (int)fage);
 4029:   free_matrix(mu,1,(nlstate+ndeath)*(nlstate+ndeath),(int) bage, (int)fage);
 4030:   free_matrix(doldm,1,(nlstate)*(nlstate+ndeath),1,(nlstate)*(nlstate+ndeath));
 4031:   free_matrix(dnewm,1,(nlstate)*(nlstate+ndeath),1,npar);
 4032:   free_vector(xp,1,npar);
 4033:   fclose(ficresprob);
 4034:   fclose(ficresprobcov);
 4035:   fclose(ficresprobcor);
 4036:   fflush(ficgp);
 4037:   fflush(fichtmcov);
 4038: }
 4039: 
 4040: 
 4041: /******************* Printing html file ***********/
 4042: void printinghtml(char fileres[], char title[], char datafile[], int firstpass, \
 4043: 		  int lastpass, int stepm, int weightopt, char model[],\
 4044: 		  int imx,int jmin, int jmax, double jmeanint,char rfileres[],\
 4045: 		  int popforecast, int estepm ,\
 4046: 		  double jprev1, double mprev1,double anprev1, \
 4047: 		  double jprev2, double mprev2,double anprev2){
 4048:   int jj1, k1, i1, cpt;
 4049: 
 4050:    fprintf(fichtm,"<ul><li><a href='#firstorder'>Result files (first order: no variance)</a>\n \
 4051:    <li><a href='#secondorder'>Result files (second order (variance)</a>\n \
 4052: </ul>");
 4053:    fprintf(fichtm,"<ul><li><h4><a name='firstorder'>Result files (first order: no variance)</a></h4>\n \
 4054:  - Observed prevalence in each state (during the period defined between %.lf/%.lf/%.lf and %.lf/%.lf/%.lf): <a href=\"%s\">%s</a> <br>\n ",
 4055: 	   jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,subdirf2(fileres,"p"),subdirf2(fileres,"p"));
 4056:    fprintf(fichtm,"\
 4057:  - Estimated transition probabilities over %d (stepm) months: <a href=\"%s\">%s</a><br>\n ",
 4058: 	   stepm,subdirf2(fileres,"pij"),subdirf2(fileres,"pij"));
 4059:    fprintf(fichtm,"\
 4060:  - Period (stable) prevalence in each health state: <a href=\"%s\">%s</a> <br>\n",
 4061: 	   subdirf2(fileres,"pl"),subdirf2(fileres,"pl"));
 4062:    fprintf(fichtm,"\
 4063:  - (a) Life expectancies by health status at initial age, ei. (b) health expectancies by health status at initial age, eij . If one or more covariates are included, specific tables for each value of the covariate are output in sequences within the same file (estepm=%2d months): \
 4064:    <a href=\"%s\">%s</a> <br>\n",
 4065: 	   estepm,subdirf2(fileres,"e"),subdirf2(fileres,"e"));
 4066:    fprintf(fichtm,"\
 4067:  - Population projections by age and states: \
 4068:    <a href=\"%s\">%s</a> <br>\n</li>", subdirf2(fileres,"f"),subdirf2(fileres,"f"));
 4069: 
 4070: fprintf(fichtm," \n<ul><li><b>Graphs</b></li><p>");
 4071: 
 4072:  m=pow(2,cptcoveff);
 4073:  if (cptcovn < 1) {m=1;ncodemax[1]=1;}
 4074: 
 4075:  jj1=0;
 4076:  for(k1=1; k1<=m;k1++){
 4077:    for(i1=1; i1<=ncodemax[k1];i1++){
 4078:      jj1++;
 4079:      if (cptcovn > 0) {
 4080:        fprintf(fichtm,"<hr  size=\"2\" color=\"#EC5E5E\">************ Results for covariates");
 4081:        for (cpt=1; cpt<=cptcoveff;cpt++) 
 4082: 	 fprintf(fichtm," V%d=%d ",Tvaraff[cpt],nbcode[Tvaraff[cpt]][codtab[jj1][cpt]]);
 4083:        fprintf(fichtm," ************\n<hr size=\"2\" color=\"#EC5E5E\">");
 4084:      }
 4085:      /* Pij */
 4086:      fprintf(fichtm,"<br>- Pij or Conditional probabilities to be observed in state j being in state i, %d (stepm) months before: <a href=\"%s%d_1.png\">%s%d_1.png</a><br> \
 4087: <img src=\"%s%d_1.png\">",stepm,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1);     
 4088:      /* Quasi-incidences */
 4089:      fprintf(fichtm,"<br>- Pij or Conditional probabilities to be observed in state j being in state i %d (stepm) months\
 4090:  before but expressed in per year i.e. quasi incidences if stepm is small and probabilities too: <a href=\"%s%d_2.png\">%s%d_2.png</a><br> \
 4091: <img src=\"%s%d_2.png\">",stepm,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1); 
 4092:        /* Period (stable) prevalence in each health state */
 4093:        for(cpt=1; cpt<=nlstate;cpt++){
 4094: 	 fprintf(fichtm,"<br>- Convergence from cross-sectional prevalence in each state (1 to %d) to period (stable) prevalence in specific state %d <a href=\"%s%d_%d.png\">%s%d_%d.png</a><br> \
 4095: <img src=\"%s%d_%d.png\">",nlstate, cpt, subdirf2(optionfilefiname,"p"),cpt,jj1,subdirf2(optionfilefiname,"p"),cpt,jj1,subdirf2(optionfilefiname,"p"),cpt,jj1);
 4096:        }
 4097:      for(cpt=1; cpt<=nlstate;cpt++) {
 4098:         fprintf(fichtm,"\n<br>- Life expectancy by health state (%d) at initial age and its decomposition into health expectancies in each alive state (1 to %d) : <a href=\"%s%d%d.png\">%s%d%d.png</a> <br> \
 4099: <img src=\"%s%d%d.png\">",cpt,nlstate,subdirf2(optionfilefiname,"exp"),cpt,jj1,subdirf2(optionfilefiname,"exp"),cpt,jj1,subdirf2(optionfilefiname,"exp"),cpt,jj1);
 4100:      }
 4101:    } /* end i1 */
 4102:  }/* End k1 */
 4103:  fprintf(fichtm,"</ul>");
 4104: 
 4105: 
 4106:  fprintf(fichtm,"\
 4107: \n<br><li><h4> <a name='secondorder'>Result files (second order: variances)</a></h4>\n\
 4108:  - Parameter file with estimated parameters and covariance matrix: <a href=\"%s\">%s</a> <br>\n", rfileres,rfileres);
 4109: 
 4110:  fprintf(fichtm," - Variance of one-step probabilities: <a href=\"%s\">%s</a> <br>\n",
 4111: 	 subdirf2(fileres,"prob"),subdirf2(fileres,"prob"));
 4112:  fprintf(fichtm,"\
 4113:  - Variance-covariance of one-step probabilities: <a href=\"%s\">%s</a> <br>\n",
 4114: 	 subdirf2(fileres,"probcov"),subdirf2(fileres,"probcov"));
 4115: 
 4116:  fprintf(fichtm,"\
 4117:  - Correlation matrix of one-step probabilities: <a href=\"%s\">%s</a> <br>\n",
 4118: 	 subdirf2(fileres,"probcor"),subdirf2(fileres,"probcor"));
 4119:  fprintf(fichtm,"\
 4120:  - Variances and covariances of health expectancies by age and <b>initial health status</b> (cov(e<sup>ij</sup>,e<sup>kl</sup>)(estepm=%2d months): \
 4121:    <a href=\"%s\">%s</a> <br>\n</li>",
 4122: 	   estepm,subdirf2(fileres,"cve"),subdirf2(fileres,"cve"));
 4123:  fprintf(fichtm,"\
 4124:  - (a) Health expectancies by health status at initial age (e<sup>ij</sup>) and standard errors (in parentheses) (b) life expectancies and standard errors (e<sup>i.</sup>=e<sup>i1</sup>+e<sup>i2</sup>+...)(estepm=%2d months): \
 4125:    <a href=\"%s\">%s</a> <br>\n</li>",
 4126: 	   estepm,subdirf2(fileres,"stde"),subdirf2(fileres,"stde"));
 4127:  fprintf(fichtm,"\
 4128:  - Variances and covariances of health expectancies by age. Status (i) based health expectancies (in state j), e<sup>ij</sup> are weighted by the period prevalences in each state i (if popbased=1, an additional computation is done using the cross-sectional prevalences, i.e population based) (estepm=%d months): <a href=\"%s\">%s</a><br>\n",
 4129: 	 estepm, subdirf2(fileres,"v"),subdirf2(fileres,"v"));
 4130:  fprintf(fichtm,"\
 4131:  - Total life expectancy and total health expectancies to be spent in each health state e<sup>.j</sup> with their standard errors (if popbased=1, an additional computation is done using the cross-sectional prevalences, i.e population based) (estepm=%d months): <a href=\"%s\">%s</a> <br>\n",
 4132: 	 estepm, subdirf2(fileres,"t"),subdirf2(fileres,"t"));
 4133:  fprintf(fichtm,"\
 4134:  - Standard deviation of period (stable) prevalences: <a href=\"%s\">%s</a> <br>\n",\
 4135: 	 subdirf2(fileres,"vpl"),subdirf2(fileres,"vpl"));
 4136: 
 4137: /*  if(popforecast==1) fprintf(fichtm,"\n */
 4138: /*  - Prevalences forecasting: <a href=\"f%s\">f%s</a> <br>\n */
 4139: /*  - Population forecasting (if popforecast=1): <a href=\"pop%s\">pop%s</a> <br>\n */
 4140: /* 	<br>",fileres,fileres,fileres,fileres); */
 4141: /*  else  */
 4142: /*    fprintf(fichtm,"\n No population forecast: popforecast = %d (instead of 1) or stepm = %d (instead of 1) or model=%s (instead of .)<br><br></li>\n",popforecast, stepm, model); */
 4143:  fflush(fichtm);
 4144:  fprintf(fichtm," <ul><li><b>Graphs</b></li><p>");
 4145: 
 4146:  m=pow(2,cptcoveff);
 4147:  if (cptcovn < 1) {m=1;ncodemax[1]=1;}
 4148: 
 4149:  jj1=0;
 4150:  for(k1=1; k1<=m;k1++){
 4151:    for(i1=1; i1<=ncodemax[k1];i1++){
 4152:      jj1++;
 4153:      if (cptcovn > 0) {
 4154:        fprintf(fichtm,"<hr  size=\"2\" color=\"#EC5E5E\">************ Results for covariates");
 4155:        for (cpt=1; cpt<=cptcoveff;cpt++) 
 4156: 	 fprintf(fichtm," V%d=%d ",Tvaraff[cpt],nbcode[Tvaraff[cpt]][codtab[jj1][cpt]]);
 4157:        fprintf(fichtm," ************\n<hr size=\"2\" color=\"#EC5E5E\">");
 4158:      }
 4159:      for(cpt=1; cpt<=nlstate;cpt++) {
 4160:        fprintf(fichtm,"<br>- Observed (cross-sectional) and period (incidence based) \
 4161: prevalence (with 95%% confidence interval) in state (%d): %s%d_%d.png <br>\
 4162: <img src=\"%s%d_%d.png\">",cpt,subdirf2(optionfilefiname,"v"),cpt,jj1,subdirf2(optionfilefiname,"v"),cpt,jj1);  
 4163:      }
 4164:      fprintf(fichtm,"\n<br>- Total life expectancy by age and \
 4165: health expectancies in states (1) and (2). If popbased=1 the smooth (due to the model) \
 4166: true period expectancies (those weighted with period prevalences are also\
 4167:  drawn in addition to the population based expectancies computed using\
 4168:  observed and cahotic prevalences: %s%d.png<br>\
 4169: <img src=\"%s%d.png\">",subdirf2(optionfilefiname,"e"),jj1,subdirf2(optionfilefiname,"e"),jj1);
 4170:    } /* end i1 */
 4171:  }/* End k1 */
 4172:  fprintf(fichtm,"</ul>");
 4173:  fflush(fichtm);
 4174: }
 4175: 
 4176: /******************* Gnuplot file **************/
 4177: void printinggnuplot(char fileres[], char optionfilefiname[], double ageminpar, double agemaxpar, double fage , char pathc[], double p[]){
 4178: 
 4179:   char dirfileres[132],optfileres[132];
 4180:   int cpt=0,k1=0,i=0,k=0,j=0,jk=0,k2=0,k3=0,ij=0,l=0;
 4181:   int ng=0;
 4182: /*   if((ficgp=fopen(optionfilegnuplot,"a"))==NULL) { */
 4183: /*     printf("Problem with file %s",optionfilegnuplot); */
 4184: /*     fprintf(ficlog,"Problem with file %s",optionfilegnuplot); */
 4185: /*   } */
 4186: 
 4187:   /*#ifdef windows */
 4188:   fprintf(ficgp,"cd \"%s\" \n",pathc);
 4189:     /*#endif */
 4190:   m=pow(2,cptcoveff);
 4191: 
 4192:   strcpy(dirfileres,optionfilefiname);
 4193:   strcpy(optfileres,"vpl");
 4194:  /* 1eme*/
 4195:   fprintf(ficgp,"\n# 1st: Period (stable) prevalence with CI: 'vpl' files\n");
 4196:   for (cpt=1; cpt<= nlstate ; cpt ++) {
 4197:     for (k1=1; k1<= m ; k1 ++) { /* plot [100000000000000000000:-100000000000000000000] "mysbiaspar/vplrmysbiaspar.txt to check */
 4198:      fprintf(ficgp,"\nset out \"%s%d_%d.png\" \n",subdirf2(optionfilefiname,"v"),cpt,k1);
 4199:      fprintf(ficgp,"\n#set out \"v%s%d_%d.png\" \n",optionfilefiname,cpt,k1);
 4200:      fprintf(ficgp,"set xlabel \"Age\" \n\
 4201: set ylabel \"Probability\" \n\
 4202: set ter png small size 320, 240\n\
 4203: plot [%.f:%.f] \"%s\" every :::%d::%d u 1:2 \"\%%lf",ageminpar,fage,subdirf2(fileres,"vpl"),k1-1,k1-1);
 4204: 
 4205:      for (i=1; i<= nlstate ; i ++) {
 4206:        if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
 4207:        else        fprintf(ficgp," \%%*lf (\%%*lf)");
 4208:      }
 4209:      fprintf(ficgp,"\" t\"Period (stable) prevalence\" w l lt 0,\"%s\" every :::%d::%d u 1:($2+1.96*$3) \"\%%lf",subdirf2(fileres,"vpl"),k1-1,k1-1);
 4210:      for (i=1; i<= nlstate ; i ++) {
 4211:        if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
 4212:        else fprintf(ficgp," \%%*lf (\%%*lf)");
 4213:      } 
 4214:      fprintf(ficgp,"\" t\"95\%% CI\" w l lt 1,\"%s\" every :::%d::%d u 1:($2-1.96*$3) \"\%%lf",subdirf2(fileres,"vpl"),k1-1,k1-1); 
 4215:      for (i=1; i<= nlstate ; i ++) {
 4216:        if (i==cpt) fprintf(ficgp," \%%lf (\%%lf)");
 4217:        else fprintf(ficgp," \%%*lf (\%%*lf)");
 4218:      }  
 4219:      fprintf(ficgp,"\" t\"\" w l lt 1,\"%s\" every :::%d::%d u 1:($%d) t\"Observed prevalence \" w l lt 2",subdirf2(fileres,"p"),k1-1,k1-1,2+4*(cpt-1));
 4220:    }
 4221:   }
 4222:   /*2 eme*/
 4223:   fprintf(ficgp,"\n# 2nd: Total life expectancy with CI: 't' files\n");
 4224:   for (k1=1; k1<= m ; k1 ++) { 
 4225:     fprintf(ficgp,"\nset out \"%s%d.png\" \n",subdirf2(optionfilefiname,"e"),k1);
 4226:     fprintf(ficgp,"set ylabel \"Years\" \nset ter png small size 320, 240\nplot [%.f:%.f] ",ageminpar,fage);
 4227:     
 4228:     for (i=1; i<= nlstate+1 ; i ++) {
 4229:       k=2*i;
 4230:       fprintf(ficgp,"\"%s\" every :::%d::%d u 1:2 \"\%%lf",subdirf2(fileres,"t"),k1-1,k1-1);
 4231:       for (j=1; j<= nlstate+1 ; j ++) {
 4232: 	if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
 4233: 	else fprintf(ficgp," \%%*lf (\%%*lf)");
 4234:       }   
 4235:       if (i== 1) fprintf(ficgp,"\" t\"TLE\" w l ,");
 4236:       else fprintf(ficgp,"\" t\"LE in state (%d)\" w l ,",i-1);
 4237:       fprintf(ficgp,"\"%s\" every :::%d::%d u 1:($2-$3*2) \"\%%lf",subdirf2(fileres,"t"),k1-1,k1-1);
 4238:       for (j=1; j<= nlstate+1 ; j ++) {
 4239: 	if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
 4240: 	else fprintf(ficgp," \%%*lf (\%%*lf)");
 4241:       }   
 4242:       fprintf(ficgp,"\" t\"\" w l lt 0,");
 4243:       fprintf(ficgp,"\"%s\" every :::%d::%d u 1:($2+$3*2) \"\%%lf",subdirf2(fileres,"t"),k1-1,k1-1);
 4244:       for (j=1; j<= nlstate+1 ; j ++) {
 4245: 	if (j==i) fprintf(ficgp," \%%lf (\%%lf)");
 4246: 	else fprintf(ficgp," \%%*lf (\%%*lf)");
 4247:       }   
 4248:       if (i== (nlstate+1)) fprintf(ficgp,"\" t\"\" w l lt 0");
 4249:       else fprintf(ficgp,"\" t\"\" w l lt 0,");
 4250:     }
 4251:   }
 4252:   
 4253:   /*3eme*/
 4254:   
 4255:   for (k1=1; k1<= m ; k1 ++) { 
 4256:     for (cpt=1; cpt<= nlstate ; cpt ++) {
 4257:       /*       k=2+nlstate*(2*cpt-2); */
 4258:       k=2+(nlstate+1)*(cpt-1);
 4259:       fprintf(ficgp,"\nset out \"%s%d%d.png\" \n",subdirf2(optionfilefiname,"exp"),cpt,k1);
 4260:       fprintf(ficgp,"set ter png small size 320, 240\n\
 4261: plot [%.f:%.f] \"%s\" every :::%d::%d u 1:%d t \"e%d1\" w l",ageminpar,fage,subdirf2(fileres,"e"),k1-1,k1-1,k,cpt);
 4262:       /*fprintf(ficgp,",\"e%s\" every :::%d::%d u 1:($%d-2*$%d) \"\%%lf ",fileres,k1-1,k1-1,k,k+1);
 4263: 	for (i=1; i<= nlstate*2 ; i ++) fprintf(ficgp,"\%%lf (\%%lf) ");
 4264: 	fprintf(ficgp,"\" t \"e%d1\" w l",cpt);
 4265: 	fprintf(ficgp,",\"e%s\" every :::%d::%d u 1:($%d+2*$%d) \"\%%lf ",fileres,k1-1,k1-1,k,k+1);
 4266: 	for (i=1; i<= nlstate*2 ; i ++) fprintf(ficgp,"\%%lf (\%%lf) ");
 4267: 	fprintf(ficgp,"\" t \"e%d1\" w l",cpt);
 4268: 	
 4269:       */
 4270:       for (i=1; i< nlstate ; i ++) {
 4271: 	fprintf(ficgp," ,\"%s\" every :::%d::%d u 1:%d t \"e%d%d\" w l",subdirf2(fileres,"e"),k1-1,k1-1,k+i,cpt,i+1);
 4272: 	/*	fprintf(ficgp," ,\"%s\" every :::%d::%d u 1:%d t \"e%d%d\" w l",subdirf2(fileres,"e"),k1-1,k1-1,k+2*i,cpt,i+1);*/
 4273: 	
 4274:       } 
 4275:       fprintf(ficgp," ,\"%s\" every :::%d::%d u 1:%d t \"e%d.\" w l",subdirf2(fileres,"e"),k1-1,k1-1,k+nlstate,cpt);
 4276:     }
 4277:   }
 4278:   
 4279:   /* CV preval stable (period) */
 4280:   for (k1=1; k1<= m ; k1 ++) { /* For each multivariate if any */
 4281:     for (cpt=1; cpt<=nlstate ; cpt ++) { /* For each life state */
 4282:       k=3;
 4283:       fprintf(ficgp,"\n#\n#\n#CV preval stable (period): 'pij' files, cov=%d state=%d",k1, cpt);
 4284:       fprintf(ficgp,"\nset out \"%s%d_%d.png\" \n",subdirf2(optionfilefiname,"p"),cpt,k1);
 4285:       fprintf(ficgp,"set xlabel \"Age\" \nset ylabel \"Probability\" \n\
 4286: set ter png small size 320, 240\n\
 4287: unset log y\n\
 4288: plot [%.f:%.f]  ", ageminpar, agemaxpar);
 4289:       for (i=1; i<= nlstate ; i ++){
 4290: 	if(i==1)
 4291: 	  fprintf(ficgp,"\"%s\"",subdirf2(fileres,"pij"));
 4292: 	else
 4293: 	  fprintf(ficgp,", '' ");
 4294: 	l=(nlstate+ndeath)*(i-1)+1;
 4295: 	fprintf(ficgp," u ($1==%d ? ($3):1/0):($%d/($%d",k1,k+l+(cpt-1),k+l);
 4296: 	for (j=1; j<= (nlstate-1) ; j ++)
 4297: 	  fprintf(ficgp,"+$%d",k+l+j);
 4298: 	fprintf(ficgp,")) t \"prev(%d,%d)\" w l",i,cpt);
 4299:       } /* nlstate */
 4300:       fprintf(ficgp,"\n");
 4301:     } /* end cpt state*/ 
 4302:   } /* end covariate */  
 4303:   
 4304:   /* proba elementaires */
 4305:   for(i=1,jk=1; i <=nlstate; i++){
 4306:     for(k=1; k <=(nlstate+ndeath); k++){
 4307:       if (k != i) {
 4308: 	for(j=1; j <=ncovmodel; j++){
 4309: 	  fprintf(ficgp,"p%d=%f ",jk,p[jk]);
 4310: 	  jk++; 
 4311: 	  fprintf(ficgp,"\n");
 4312: 	}
 4313:       }
 4314:     }
 4315:    }
 4316:   /*goto avoid;*/
 4317:    for(ng=1; ng<=2;ng++){ /* Number of graphics: first is probabilities second is incidence per year*/
 4318:      for(jk=1; jk <=m; jk++) {
 4319:        fprintf(ficgp,"\nset out \"%s%d_%d.png\" \n",subdirf2(optionfilefiname,"pe"),jk,ng); 
 4320:        if (ng==2)
 4321: 	 fprintf(ficgp,"\nset ylabel \"Quasi-incidence per year\"\n");
 4322:        else
 4323: 	 fprintf(ficgp,"\nset title \"Probability\"\n");
 4324:        fprintf(ficgp,"\nset ter png small size 320, 240\nset log y\nplot  [%.f:%.f] ",ageminpar,agemaxpar);
 4325:        i=1;
 4326:        for(k2=1; k2<=nlstate; k2++) {
 4327: 	 k3=i;
 4328: 	 for(k=1; k<=(nlstate+ndeath); k++) {
 4329: 	   if (k != k2){
 4330: 	     if(ng==2)
 4331: 	       fprintf(ficgp," %f*exp(p%d+p%d*x",YEARM/stepm,i,i+1);
 4332: 	     else
 4333: 	       fprintf(ficgp," exp(p%d+p%d*x",i,i+1);
 4334: 	     ij=1;/* To be checked else nbcode[0][0] wrong */
 4335: 	     for(j=3; j <=ncovmodel; j++) {
 4336: 	       /* if(((j-2)==Tage[ij]) &&(ij <=cptcovage)) { /\* Bug valgrind *\/ */
 4337: 	       /* 	 /\*fprintf(ficgp,"+p%d*%d*x",i+j-1,nbcode[Tvar[j-2]][codtab[jk][Tvar[j-2]]]);*\/ */
 4338: 	       /* 	 ij++; */
 4339: 	       /* } */
 4340: 	       /* else */
 4341: 		 fprintf(ficgp,"+p%d*%d",i+j-1,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
 4342: 	     }
 4343: 	     fprintf(ficgp,")/(1");
 4344: 	     
 4345: 	     for(k1=1; k1 <=nlstate; k1++){   
 4346: 	       fprintf(ficgp,"+exp(p%d+p%d*x",k3+(k1-1)*ncovmodel,k3+(k1-1)*ncovmodel+1);
 4347: 	       ij=1;
 4348: 	       for(j=3; j <=ncovmodel; j++){
 4349: 		 /* if(((j-2)==Tage[ij]) &&(ij <=cptcovage)) { */
 4350: 		 /*   fprintf(ficgp,"+p%d*%d*x",k3+(k1-1)*ncovmodel+1+j-2,nbcode[Tvar[j-2]][codtab[jk][Tvar[j-2]]]); */
 4351: 		 /*   ij++; */
 4352: 		 /* } */
 4353: 		 /* else */
 4354: 		   fprintf(ficgp,"+p%d*%d",k3+(k1-1)*ncovmodel+1+j-2,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
 4355: 	       }
 4356: 	       fprintf(ficgp,")");
 4357: 	     }
 4358: 	     fprintf(ficgp,") t \"p%d%d\" ", k2,k);
 4359: 	     if ((k+k2)!= (nlstate*2+ndeath)) fprintf(ficgp,",");
 4360: 	     i=i+ncovmodel;
 4361: 	   }
 4362: 	 } /* end k */
 4363:        } /* end k2 */
 4364:      } /* end jk */
 4365:    } /* end ng */
 4366:  /* avoid: */
 4367:    fflush(ficgp); 
 4368: }  /* end gnuplot */
 4369: 
 4370: 
 4371: /*************** Moving average **************/
 4372: int movingaverage(double ***probs, double bage,double fage, double ***mobaverage, int mobilav){
 4373: 
 4374:   int i, cpt, cptcod;
 4375:   int modcovmax =1;
 4376:   int mobilavrange, mob;
 4377:   double age;
 4378: 
 4379:   modcovmax=2*cptcoveff;/* Max number of modalities. We suppose 
 4380: 			   a covariate has 2 modalities */
 4381:   if (cptcovn<1) modcovmax=1; /* At least 1 pass */
 4382: 
 4383:   if(mobilav==1||mobilav ==3 ||mobilav==5 ||mobilav== 7){
 4384:     if(mobilav==1) mobilavrange=5; /* default */
 4385:     else mobilavrange=mobilav;
 4386:     for (age=bage; age<=fage; age++)
 4387:       for (i=1; i<=nlstate;i++)
 4388: 	for (cptcod=1;cptcod<=modcovmax;cptcod++)
 4389: 	  mobaverage[(int)age][i][cptcod]=probs[(int)age][i][cptcod];
 4390:     /* We keep the original values on the extreme ages bage, fage and for 
 4391:        fage+1 and bage-1 we use a 3 terms moving average; for fage+2 bage+2
 4392:        we use a 5 terms etc. until the borders are no more concerned. 
 4393:     */ 
 4394:     for (mob=3;mob <=mobilavrange;mob=mob+2){
 4395:       for (age=bage+(mob-1)/2; age<=fage-(mob-1)/2; age++){
 4396: 	for (i=1; i<=nlstate;i++){
 4397: 	  for (cptcod=1;cptcod<=modcovmax;cptcod++){
 4398: 	    mobaverage[(int)age][i][cptcod] =probs[(int)age][i][cptcod];
 4399: 	      for (cpt=1;cpt<=(mob-1)/2;cpt++){
 4400: 		mobaverage[(int)age][i][cptcod] +=probs[(int)age-cpt][i][cptcod];
 4401: 		mobaverage[(int)age][i][cptcod] +=probs[(int)age+cpt][i][cptcod];
 4402: 	      }
 4403: 	    mobaverage[(int)age][i][cptcod]=mobaverage[(int)age][i][cptcod]/mob;
 4404: 	  }
 4405: 	}
 4406:       }/* end age */
 4407:     }/* end mob */
 4408:   }else return -1;
 4409:   return 0;
 4410: }/* End movingaverage */
 4411: 
 4412: 
 4413: /************** Forecasting ******************/
 4414: prevforecast(char fileres[], double anproj1, double mproj1, double jproj1, double ageminpar, double agemax, double dateprev1, double dateprev2, int mobilav, double bage, double fage, int firstpass, int lastpass, double anproj2, double p[], int cptcoveff){
 4415:   /* proj1, year, month, day of starting projection 
 4416:      agemin, agemax range of age
 4417:      dateprev1 dateprev2 range of dates during which prevalence is computed
 4418:      anproj2 year of en of projection (same day and month as proj1).
 4419:   */
 4420:   int yearp, stepsize, hstepm, nhstepm, j, k, cptcod, i, h, i1;
 4421:   double agec; /* generic age */
 4422:   double agelim, ppij, yp,yp1,yp2,jprojmean,mprojmean,anprojmean;
 4423:   double *popeffectif,*popcount;
 4424:   double ***p3mat;
 4425:   double ***mobaverage;
 4426:   char fileresf[FILENAMELENGTH];
 4427: 
 4428:   agelim=AGESUP;
 4429:   prevalence(probs, ageminpar, agemax, s, agev, nlstate, imx, Tvar, nbcode, ncodemax, mint, anint, dateprev1, dateprev2, firstpass, lastpass);
 4430:  
 4431:   strcpy(fileresf,"f"); 
 4432:   strcat(fileresf,fileres);
 4433:   if((ficresf=fopen(fileresf,"w"))==NULL) {
 4434:     printf("Problem with forecast resultfile: %s\n", fileresf);
 4435:     fprintf(ficlog,"Problem with forecast resultfile: %s\n", fileresf);
 4436:   }
 4437:   printf("Computing forecasting: result on file '%s' \n", fileresf);
 4438:   fprintf(ficlog,"Computing forecasting: result on file '%s' \n", fileresf);
 4439: 
 4440:   if (cptcoveff==0) ncodemax[cptcoveff]=1;
 4441: 
 4442:   if (mobilav!=0) {
 4443:     mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4444:     if (movingaverage(probs, ageminpar, fage, mobaverage,mobilav)!=0){
 4445:       fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 4446:       printf(" Error in movingaverage mobilav=%d\n",mobilav);
 4447:     }
 4448:   }
 4449: 
 4450:   stepsize=(int) (stepm+YEARM-1)/YEARM;
 4451:   if (stepm<=12) stepsize=1;
 4452:   if(estepm < stepm){
 4453:     printf ("Problem %d lower than %d\n",estepm, stepm);
 4454:   }
 4455:   else  hstepm=estepm;   
 4456: 
 4457:   hstepm=hstepm/stepm; 
 4458:   yp1=modf(dateintmean,&yp);/* extracts integral of datemean in yp  and
 4459:                                fractional in yp1 */
 4460:   anprojmean=yp;
 4461:   yp2=modf((yp1*12),&yp);
 4462:   mprojmean=yp;
 4463:   yp1=modf((yp2*30.5),&yp);
 4464:   jprojmean=yp;
 4465:   if(jprojmean==0) jprojmean=1;
 4466:   if(mprojmean==0) jprojmean=1;
 4467: 
 4468:   i1=cptcoveff;
 4469:   if (cptcovn < 1){i1=1;}
 4470:   
 4471:   fprintf(ficresf,"# Mean day of interviews %.lf/%.lf/%.lf (%.2f) between %.2f and %.2f \n",jprojmean,mprojmean,anprojmean,dateintmean,dateprev1,dateprev2); 
 4472:   
 4473:   fprintf(ficresf,"#****** Routine prevforecast **\n");
 4474: 
 4475: /* 	      if (h==(int)(YEARM*yearp)){ */
 4476:   for(cptcov=1, k=0;cptcov<=i1;cptcov++){
 4477:     for(cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
 4478:       k=k+1;
 4479:       fprintf(ficresf,"\n#******");
 4480:       for(j=1;j<=cptcoveff;j++) {
 4481: 	fprintf(ficresf," V%d=%d, hpijx=probability over h years, hp.jx is weighted by observed prev ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 4482:       }
 4483:       fprintf(ficresf,"******\n");
 4484:       fprintf(ficresf,"# Covariate valuofcovar yearproj age");
 4485:       for(j=1; j<=nlstate+ndeath;j++){ 
 4486: 	for(i=1; i<=nlstate;i++) 	      
 4487:           fprintf(ficresf," p%d%d",i,j);
 4488: 	fprintf(ficresf," p.%d",j);
 4489:       }
 4490:       for (yearp=0; yearp<=(anproj2-anproj1);yearp +=stepsize) { 
 4491: 	fprintf(ficresf,"\n");
 4492: 	fprintf(ficresf,"\n# Forecasting at date %.lf/%.lf/%.lf ",jproj1,mproj1,anproj1+yearp);   
 4493: 
 4494:      	for (agec=fage; agec>=(ageminpar-1); agec--){ 
 4495: 	  nhstepm=(int) rint((agelim-agec)*YEARM/stepm); 
 4496: 	  nhstepm = nhstepm/hstepm; 
 4497: 	  p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4498: 	  oldm=oldms;savm=savms;
 4499: 	  hpxij(p3mat,nhstepm,agec,hstepm,p,nlstate,stepm,oldm,savm, k);  
 4500: 	
 4501: 	  for (h=0; h<=nhstepm; h++){
 4502: 	    if (h*hstepm/YEARM*stepm ==yearp) {
 4503:               fprintf(ficresf,"\n");
 4504:               for(j=1;j<=cptcoveff;j++) 
 4505:                 fprintf(ficresf,"%d %d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 4506: 	      fprintf(ficresf,"%.f %.f ",anproj1+yearp,agec+h*hstepm/YEARM*stepm);
 4507: 	    } 
 4508: 	    for(j=1; j<=nlstate+ndeath;j++) {
 4509: 	      ppij=0.;
 4510: 	      for(i=1; i<=nlstate;i++) {
 4511: 		if (mobilav==1) 
 4512: 		  ppij=ppij+p3mat[i][j][h]*mobaverage[(int)agec][i][cptcod];
 4513: 		else {
 4514: 		  ppij=ppij+p3mat[i][j][h]*probs[(int)(agec)][i][cptcod];
 4515: 		}
 4516: 		if (h*hstepm/YEARM*stepm== yearp) {
 4517: 		  fprintf(ficresf," %.3f", p3mat[i][j][h]);
 4518: 		}
 4519: 	      } /* end i */
 4520: 	      if (h*hstepm/YEARM*stepm==yearp) {
 4521: 		fprintf(ficresf," %.3f", ppij);
 4522: 	      }
 4523: 	    }/* end j */
 4524: 	  } /* end h */
 4525: 	  free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4526: 	} /* end agec */
 4527:       } /* end yearp */
 4528:     } /* end cptcod */
 4529:   } /* end  cptcov */
 4530:        
 4531:   if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4532: 
 4533:   fclose(ficresf);
 4534: }
 4535: 
 4536: /************** Forecasting *****not tested NB*************/
 4537: populforecast(char fileres[], double anpyram,double mpyram,double jpyram,double ageminpar, double agemax,double dateprev1, double dateprev2, int mobilav, double agedeb, double fage, int popforecast, char popfile[], double anpyram1,double p[], int i2){
 4538:   
 4539:   int cpt, stepsize, hstepm, nhstepm, j,k,c, cptcod, i,h;
 4540:   int *popage;
 4541:   double calagedatem, agelim, kk1, kk2;
 4542:   double *popeffectif,*popcount;
 4543:   double ***p3mat,***tabpop,***tabpopprev;
 4544:   double ***mobaverage;
 4545:   char filerespop[FILENAMELENGTH];
 4546: 
 4547:   tabpop= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4548:   tabpopprev= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4549:   agelim=AGESUP;
 4550:   calagedatem=(anpyram+mpyram/12.+jpyram/365.-dateintmean)*YEARM;
 4551:   
 4552:   prevalence(probs, ageminpar, agemax, s, agev, nlstate, imx, Tvar, nbcode, ncodemax, mint, anint, dateprev1, dateprev2, firstpass, lastpass);
 4553:   
 4554:   
 4555:   strcpy(filerespop,"pop"); 
 4556:   strcat(filerespop,fileres);
 4557:   if((ficrespop=fopen(filerespop,"w"))==NULL) {
 4558:     printf("Problem with forecast resultfile: %s\n", filerespop);
 4559:     fprintf(ficlog,"Problem with forecast resultfile: %s\n", filerespop);
 4560:   }
 4561:   printf("Computing forecasting: result on file '%s' \n", filerespop);
 4562:   fprintf(ficlog,"Computing forecasting: result on file '%s' \n", filerespop);
 4563: 
 4564:   if (cptcoveff==0) ncodemax[cptcoveff]=1;
 4565: 
 4566:   if (mobilav!=0) {
 4567:     mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4568:     if (movingaverage(probs, ageminpar, fage, mobaverage,mobilav)!=0){
 4569:       fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 4570:       printf(" Error in movingaverage mobilav=%d\n",mobilav);
 4571:     }
 4572:   }
 4573: 
 4574:   stepsize=(int) (stepm+YEARM-1)/YEARM;
 4575:   if (stepm<=12) stepsize=1;
 4576:   
 4577:   agelim=AGESUP;
 4578:   
 4579:   hstepm=1;
 4580:   hstepm=hstepm/stepm; 
 4581:   
 4582:   if (popforecast==1) {
 4583:     if((ficpop=fopen(popfile,"r"))==NULL) {
 4584:       printf("Problem with population file : %s\n",popfile);exit(0);
 4585:       fprintf(ficlog,"Problem with population file : %s\n",popfile);exit(0);
 4586:     } 
 4587:     popage=ivector(0,AGESUP);
 4588:     popeffectif=vector(0,AGESUP);
 4589:     popcount=vector(0,AGESUP);
 4590:     
 4591:     i=1;   
 4592:     while ((c=fscanf(ficpop,"%d %lf\n",&popage[i],&popcount[i])) != EOF) i=i+1;
 4593:    
 4594:     imx=i;
 4595:     for (i=1; i<imx;i++) popeffectif[popage[i]]=popcount[i];
 4596:   }
 4597: 
 4598:   for(cptcov=1,k=0;cptcov<=i2;cptcov++){
 4599:    for(cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
 4600:       k=k+1;
 4601:       fprintf(ficrespop,"\n#******");
 4602:       for(j=1;j<=cptcoveff;j++) {
 4603: 	fprintf(ficrespop," V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 4604:       }
 4605:       fprintf(ficrespop,"******\n");
 4606:       fprintf(ficrespop,"# Age");
 4607:       for(j=1; j<=nlstate+ndeath;j++) fprintf(ficrespop," P.%d",j);
 4608:       if (popforecast==1)  fprintf(ficrespop," [Population]");
 4609:       
 4610:       for (cpt=0; cpt<=0;cpt++) { 
 4611: 	fprintf(ficrespop,"\n\n# Forecasting at date %.lf/%.lf/%.lf ",jpyram,mpyram,anpyram+cpt);   
 4612: 	
 4613:      	for (agedeb=(fage-((int)calagedatem %12/12.)); agedeb>=(ageminpar-((int)calagedatem %12)/12.); agedeb--){ 
 4614: 	  nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm); 
 4615: 	  nhstepm = nhstepm/hstepm; 
 4616: 	  
 4617: 	  p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4618: 	  oldm=oldms;savm=savms;
 4619: 	  hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);  
 4620: 	
 4621: 	  for (h=0; h<=nhstepm; h++){
 4622: 	    if (h==(int) (calagedatem+YEARM*cpt)) {
 4623: 	      fprintf(ficrespop,"\n %3.f ",agedeb+h*hstepm/YEARM*stepm);
 4624: 	    } 
 4625: 	    for(j=1; j<=nlstate+ndeath;j++) {
 4626: 	      kk1=0.;kk2=0;
 4627: 	      for(i=1; i<=nlstate;i++) {	      
 4628: 		if (mobilav==1) 
 4629: 		  kk1=kk1+p3mat[i][j][h]*mobaverage[(int)agedeb+1][i][cptcod];
 4630: 		else {
 4631: 		  kk1=kk1+p3mat[i][j][h]*probs[(int)(agedeb+1)][i][cptcod];
 4632: 		}
 4633: 	      }
 4634: 	      if (h==(int)(calagedatem+12*cpt)){
 4635: 		tabpop[(int)(agedeb)][j][cptcod]=kk1;
 4636: 		  /*fprintf(ficrespop," %.3f", kk1);
 4637: 		    if (popforecast==1) fprintf(ficrespop," [%.f]", kk1*popeffectif[(int)agedeb+1]);*/
 4638: 	      }
 4639: 	    }
 4640: 	    for(i=1; i<=nlstate;i++){
 4641: 	      kk1=0.;
 4642: 		for(j=1; j<=nlstate;j++){
 4643: 		  kk1= kk1+tabpop[(int)(agedeb)][j][cptcod]; 
 4644: 		}
 4645: 		  tabpopprev[(int)(agedeb)][i][cptcod]=tabpop[(int)(agedeb)][i][cptcod]/kk1*popeffectif[(int)(agedeb+(calagedatem+12*cpt)*hstepm/YEARM*stepm-1)];
 4646: 	    }
 4647: 
 4648: 	    if (h==(int)(calagedatem+12*cpt)) for(j=1; j<=nlstate;j++) 
 4649: 	      fprintf(ficrespop," %15.2f",tabpopprev[(int)(agedeb+1)][j][cptcod]);
 4650: 	  }
 4651: 	  free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4652: 	}
 4653:       }
 4654:  
 4655:   /******/
 4656: 
 4657:       for (cpt=1; cpt<=(anpyram1-anpyram);cpt++) { 
 4658: 	fprintf(ficrespop,"\n\n# Forecasting at date %.lf/%.lf/%.lf ",jpyram,mpyram,anpyram+cpt);   
 4659: 	for (agedeb=(fage-((int)calagedatem %12/12.)); agedeb>=(ageminpar-((int)calagedatem %12)/12.); agedeb--){ 
 4660: 	  nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm); 
 4661: 	  nhstepm = nhstepm/hstepm; 
 4662: 	  
 4663: 	  p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4664: 	  oldm=oldms;savm=savms;
 4665: 	  hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);  
 4666: 	  for (h=0; h<=nhstepm; h++){
 4667: 	    if (h==(int) (calagedatem+YEARM*cpt)) {
 4668: 	      fprintf(ficresf,"\n %3.f ",agedeb+h*hstepm/YEARM*stepm);
 4669: 	    } 
 4670: 	    for(j=1; j<=nlstate+ndeath;j++) {
 4671: 	      kk1=0.;kk2=0;
 4672: 	      for(i=1; i<=nlstate;i++) {	      
 4673: 		kk1=kk1+p3mat[i][j][h]*tabpopprev[(int)agedeb+1][i][cptcod];	
 4674: 	      }
 4675: 	      if (h==(int)(calagedatem+12*cpt)) fprintf(ficresf," %15.2f", kk1);	
 4676: 	    }
 4677: 	  }
 4678: 	  free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4679: 	}
 4680:       }
 4681:    } 
 4682:   }
 4683:  
 4684:   if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4685: 
 4686:   if (popforecast==1) {
 4687:     free_ivector(popage,0,AGESUP);
 4688:     free_vector(popeffectif,0,AGESUP);
 4689:     free_vector(popcount,0,AGESUP);
 4690:   }
 4691:   free_ma3x(tabpop,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4692:   free_ma3x(tabpopprev,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4693:   fclose(ficrespop);
 4694: } /* End of popforecast */
 4695: 
 4696: int fileappend(FILE *fichier, char *optionfich)
 4697: {
 4698:   if((fichier=fopen(optionfich,"a"))==NULL) {
 4699:     printf("Problem with file: %s\n", optionfich);
 4700:     fprintf(ficlog,"Problem with file: %s\n", optionfich);
 4701:     return (0);
 4702:   }
 4703:   fflush(fichier);
 4704:   return (1);
 4705: }
 4706: 
 4707: 
 4708: /**************** function prwizard **********************/
 4709: void prwizard(int ncovmodel, int nlstate, int ndeath,  char model[], FILE *ficparo)
 4710: {
 4711: 
 4712:   /* Wizard to print covariance matrix template */
 4713: 
 4714:   char ca[32], cb[32];
 4715:   int i,j, k, li, lj, lk, ll, jj, npar, itimes;
 4716:   int numlinepar;
 4717: 
 4718:   printf("# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 4719:   fprintf(ficparo,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 4720:   for(i=1; i <=nlstate; i++){
 4721:     jj=0;
 4722:     for(j=1; j <=nlstate+ndeath; j++){
 4723:       if(j==i) continue;
 4724:       jj++;
 4725:       /*ca[0]= k+'a'-1;ca[1]='\0';*/
 4726:       printf("%1d%1d",i,j);
 4727:       fprintf(ficparo,"%1d%1d",i,j);
 4728:       for(k=1; k<=ncovmodel;k++){
 4729: 	/* 	  printf(" %lf",param[i][j][k]); */
 4730: 	/* 	  fprintf(ficparo," %lf",param[i][j][k]); */
 4731: 	printf(" 0.");
 4732: 	fprintf(ficparo," 0.");
 4733:       }
 4734:       printf("\n");
 4735:       fprintf(ficparo,"\n");
 4736:     }
 4737:   }
 4738:   printf("# Scales (for hessian or gradient estimation)\n");
 4739:   fprintf(ficparo,"# Scales (for hessian or gradient estimation)\n");
 4740:   npar= (nlstate+ndeath-1)*nlstate*ncovmodel; /* Number of parameters*/ 
 4741:   for(i=1; i <=nlstate; i++){
 4742:     jj=0;
 4743:     for(j=1; j <=nlstate+ndeath; j++){
 4744:       if(j==i) continue;
 4745:       jj++;
 4746:       fprintf(ficparo,"%1d%1d",i,j);
 4747:       printf("%1d%1d",i,j);
 4748:       fflush(stdout);
 4749:       for(k=1; k<=ncovmodel;k++){
 4750: 	/* 	printf(" %le",delti3[i][j][k]); */
 4751: 	/* 	fprintf(ficparo," %le",delti3[i][j][k]); */
 4752: 	printf(" 0.");
 4753: 	fprintf(ficparo," 0.");
 4754:       }
 4755:       numlinepar++;
 4756:       printf("\n");
 4757:       fprintf(ficparo,"\n");
 4758:     }
 4759:   }
 4760:   printf("# Covariance matrix\n");
 4761: /* # 121 Var(a12)\n\ */
 4762: /* # 122 Cov(b12,a12) Var(b12)\n\ */
 4763: /* # 131 Cov(a13,a12) Cov(a13,b12, Var(a13)\n\ */
 4764: /* # 132 Cov(b13,a12) Cov(b13,b12, Cov(b13,a13) Var(b13)\n\ */
 4765: /* # 212 Cov(a21,a12) Cov(a21,b12, Cov(a21,a13) Cov(a21,b13) Var(a21)\n\ */
 4766: /* # 212 Cov(b21,a12) Cov(b21,b12, Cov(b21,a13) Cov(b21,b13) Cov(b21,a21) Var(b21)\n\ */
 4767: /* # 232 Cov(a23,a12) Cov(a23,b12, Cov(a23,a13) Cov(a23,b13) Cov(a23,a21) Cov(a23,b21) Var(a23)\n\ */
 4768: /* # 232 Cov(b23,a12) Cov(b23,b12) ... Var (b23)\n" */
 4769:   fflush(stdout);
 4770:   fprintf(ficparo,"# Covariance matrix\n");
 4771:   /* # 121 Var(a12)\n\ */
 4772:   /* # 122 Cov(b12,a12) Var(b12)\n\ */
 4773:   /* #   ...\n\ */
 4774:   /* # 232 Cov(b23,a12)  Cov(b23,b12) ... Var (b23)\n" */
 4775:   
 4776:   for(itimes=1;itimes<=2;itimes++){
 4777:     jj=0;
 4778:     for(i=1; i <=nlstate; i++){
 4779:       for(j=1; j <=nlstate+ndeath; j++){
 4780: 	if(j==i) continue;
 4781: 	for(k=1; k<=ncovmodel;k++){
 4782: 	  jj++;
 4783: 	  ca[0]= k+'a'-1;ca[1]='\0';
 4784: 	  if(itimes==1){
 4785: 	    printf("#%1d%1d%d",i,j,k);
 4786: 	    fprintf(ficparo,"#%1d%1d%d",i,j,k);
 4787: 	  }else{
 4788: 	    printf("%1d%1d%d",i,j,k);
 4789: 	    fprintf(ficparo,"%1d%1d%d",i,j,k);
 4790: 	    /* 	printf(" %.5le",matcov[i][j]); */
 4791: 	  }
 4792: 	  ll=0;
 4793: 	  for(li=1;li <=nlstate; li++){
 4794: 	    for(lj=1;lj <=nlstate+ndeath; lj++){
 4795: 	      if(lj==li) continue;
 4796: 	      for(lk=1;lk<=ncovmodel;lk++){
 4797: 		ll++;
 4798: 		if(ll<=jj){
 4799: 		  cb[0]= lk +'a'-1;cb[1]='\0';
 4800: 		  if(ll<jj){
 4801: 		    if(itimes==1){
 4802: 		      printf(" Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 4803: 		      fprintf(ficparo," Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 4804: 		    }else{
 4805: 		      printf(" 0.");
 4806: 		      fprintf(ficparo," 0.");
 4807: 		    }
 4808: 		  }else{
 4809: 		    if(itimes==1){
 4810: 		      printf(" Var(%s%1d%1d)",ca,i,j);
 4811: 		      fprintf(ficparo," Var(%s%1d%1d)",ca,i,j);
 4812: 		    }else{
 4813: 		      printf(" 0.");
 4814: 		      fprintf(ficparo," 0.");
 4815: 		    }
 4816: 		  }
 4817: 		}
 4818: 	      } /* end lk */
 4819: 	    } /* end lj */
 4820: 	  } /* end li */
 4821: 	  printf("\n");
 4822: 	  fprintf(ficparo,"\n");
 4823: 	  numlinepar++;
 4824: 	} /* end k*/
 4825:       } /*end j */
 4826:     } /* end i */
 4827:   } /* end itimes */
 4828: 
 4829: } /* end of prwizard */
 4830: /******************* Gompertz Likelihood ******************************/
 4831: double gompertz(double x[])
 4832: { 
 4833:   double A,B,L=0.0,sump=0.,num=0.;
 4834:   int i,n=0; /* n is the size of the sample */
 4835: 
 4836:   for (i=0;i<=imx-1 ; i++) {
 4837:     sump=sump+weight[i];
 4838:     /*    sump=sump+1;*/
 4839:     num=num+1;
 4840:   }
 4841:  
 4842:  
 4843:   /* for (i=0; i<=imx; i++) 
 4844:      if (wav[i]>0) printf("i=%d ageex=%lf agecens=%lf agedc=%lf cens=%d %d\n" ,i,ageexmed[i],agecens[i],agedc[i],cens[i],wav[i]);*/
 4845: 
 4846:   for (i=1;i<=imx ; i++)
 4847:     {
 4848:       if (cens[i] == 1 && wav[i]>1)
 4849: 	A=-x[1]/(x[2])*(exp(x[2]*(agecens[i]-agegomp))-exp(x[2]*(ageexmed[i]-agegomp)));
 4850:       
 4851:       if (cens[i] == 0 && wav[i]>1)
 4852: 	A=-x[1]/(x[2])*(exp(x[2]*(agedc[i]-agegomp))-exp(x[2]*(ageexmed[i]-agegomp)))
 4853: 	     +log(x[1]/YEARM)+x[2]*(agedc[i]-agegomp)+log(YEARM);  
 4854:       
 4855:       /*if (wav[i] > 1 && agecens[i] > 15) {*/ /* ??? */
 4856:       if (wav[i] > 1 ) { /* ??? */
 4857: 	L=L+A*weight[i];
 4858: 	/* 	printf("\ni=%d A=%f L=%lf x[1]=%lf x[2]=%lf ageex=%lf agecens=%lf cens=%d agedc=%lf weight=%lf\n",i,A,L,x[1],x[2],ageexmed[i]*12,agecens[i]*12,cens[i],agedc[i]*12,weight[i]);*/
 4859:       }
 4860:     }
 4861: 
 4862:  /*printf("x1=%2.9f x2=%2.9f x3=%2.9f L=%f\n",x[1],x[2],x[3],L);*/
 4863:  
 4864:   return -2*L*num/sump;
 4865: }
 4866: 
 4867: #ifdef GSL
 4868: /******************* Gompertz_f Likelihood ******************************/
 4869: double gompertz_f(const gsl_vector *v, void *params)
 4870: { 
 4871:   double A,B,LL=0.0,sump=0.,num=0.;
 4872:   double *x= (double *) v->data;
 4873:   int i,n=0; /* n is the size of the sample */
 4874: 
 4875:   for (i=0;i<=imx-1 ; i++) {
 4876:     sump=sump+weight[i];
 4877:     /*    sump=sump+1;*/
 4878:     num=num+1;
 4879:   }
 4880:  
 4881:  
 4882:   /* for (i=0; i<=imx; i++) 
 4883:      if (wav[i]>0) printf("i=%d ageex=%lf agecens=%lf agedc=%lf cens=%d %d\n" ,i,ageexmed[i],agecens[i],agedc[i],cens[i],wav[i]);*/
 4884:   printf("x[0]=%lf x[1]=%lf\n",x[0],x[1]);
 4885:   for (i=1;i<=imx ; i++)
 4886:     {
 4887:       if (cens[i] == 1 && wav[i]>1)
 4888: 	A=-x[0]/(x[1])*(exp(x[1]*(agecens[i]-agegomp))-exp(x[1]*(ageexmed[i]-agegomp)));
 4889:       
 4890:       if (cens[i] == 0 && wav[i]>1)
 4891: 	A=-x[0]/(x[1])*(exp(x[1]*(agedc[i]-agegomp))-exp(x[1]*(ageexmed[i]-agegomp)))
 4892: 	     +log(x[0]/YEARM)+x[1]*(agedc[i]-agegomp)+log(YEARM);  
 4893:       
 4894:       /*if (wav[i] > 1 && agecens[i] > 15) {*/ /* ??? */
 4895:       if (wav[i] > 1 ) { /* ??? */
 4896: 	LL=LL+A*weight[i];
 4897: 	/* 	printf("\ni=%d A=%f L=%lf x[1]=%lf x[2]=%lf ageex=%lf agecens=%lf cens=%d agedc=%lf weight=%lf\n",i,A,L,x[1],x[2],ageexmed[i]*12,agecens[i]*12,cens[i],agedc[i]*12,weight[i]);*/
 4898:       }
 4899:     }
 4900: 
 4901:  /*printf("x1=%2.9f x2=%2.9f x3=%2.9f L=%f\n",x[1],x[2],x[3],L);*/
 4902:   printf("x[0]=%lf x[1]=%lf -2*LL*num/sump=%lf\n",x[0],x[1],-2*LL*num/sump);
 4903:  
 4904:   return -2*LL*num/sump;
 4905: }
 4906: #endif
 4907: 
 4908: /******************* Printing html file ***********/
 4909: void printinghtmlmort(char fileres[], char title[], char datafile[], int firstpass, \
 4910: 		  int lastpass, int stepm, int weightopt, char model[],\
 4911: 		  int imx,  double p[],double **matcov,double agemortsup){
 4912:   int i,k;
 4913: 
 4914:   fprintf(fichtm,"<ul><li><h4>Result files </h4>\n Force of mortality. Parameters of the Gompertz fit (with confidence interval in brackets):<br>");
 4915:   fprintf(fichtm,"  mu(age) =%lf*exp(%lf*(age-%d)) per year<br><br>",p[1],p[2],agegomp);
 4916:   for (i=1;i<=2;i++) 
 4917:     fprintf(fichtm," p[%d] = %lf [%f ; %f]<br>\n",i,p[i],p[i]-2*sqrt(matcov[i][i]),p[i]+2*sqrt(matcov[i][i]));
 4918:   fprintf(fichtm,"<br><br><img src=\"graphmort.png\">");
 4919:   fprintf(fichtm,"</ul>");
 4920: 
 4921: fprintf(fichtm,"<ul><li><h4>Life table</h4>\n <br>");
 4922: 
 4923:  fprintf(fichtm,"\nAge   l<inf>x</inf>     q<inf>x</inf> d(x,x+1)    L<inf>x</inf>     T<inf>x</inf>     e<infx</inf><br>");
 4924: 
 4925:  for (k=agegomp;k<(agemortsup-2);k++) 
 4926:    fprintf(fichtm,"%d %.0lf %lf %.0lf %.0lf %.0lf %lf<br>\n",k,lsurv[k],p[1]*exp(p[2]*(k-agegomp)),(p[1]*exp(p[2]*(k-agegomp)))*lsurv[k],lpop[k],tpop[k],tpop[k]/lsurv[k]);
 4927: 
 4928:  
 4929:   fflush(fichtm);
 4930: }
 4931: 
 4932: /******************* Gnuplot file **************/
 4933: void printinggnuplotmort(char fileres[], char optionfilefiname[], double ageminpar, double agemaxpar, double fage , char pathc[], double p[]){
 4934: 
 4935:   char dirfileres[132],optfileres[132];
 4936: 
 4937:   int ng;
 4938: 
 4939: 
 4940:   /*#ifdef windows */
 4941:   fprintf(ficgp,"cd \"%s\" \n",pathc);
 4942:     /*#endif */
 4943: 
 4944: 
 4945:   strcpy(dirfileres,optionfilefiname);
 4946:   strcpy(optfileres,"vpl");
 4947:   fprintf(ficgp,"set out \"graphmort.png\"\n "); 
 4948:   fprintf(ficgp,"set xlabel \"Age\"\n set ylabel \"Force of mortality (per year)\" \n "); 
 4949:   fprintf(ficgp, "set ter png small size 320, 240\n set log y\n"); 
 4950:   /* fprintf(ficgp, "set size 0.65,0.65\n"); */
 4951:   fprintf(ficgp,"plot [%d:100] %lf*exp(%lf*(x-%d))",agegomp,p[1],p[2],agegomp);
 4952: 
 4953: } 
 4954: 
 4955: int readdata(char datafile[], int firstobs, int lastobs, int *imax)
 4956: {
 4957: 
 4958:   /*-------- data file ----------*/
 4959:   FILE *fic;
 4960:   char dummy[]="                         ";
 4961:   int i=0, j=0, n=0;
 4962:   int linei, month, year,iout;
 4963:   char line[MAXLINE], linetmp[MAXLINE];
 4964:   char stra[MAXLINE], strb[MAXLINE];
 4965:   char *stratrunc;
 4966:   int lstra;
 4967: 
 4968: 
 4969:   if((fic=fopen(datafile,"r"))==NULL)    {
 4970:     printf("Problem while opening datafile: %s\n", datafile);return 1;
 4971:     fprintf(ficlog,"Problem while opening datafile: %s\n", datafile);return 1;
 4972:   }
 4973: 
 4974:   i=1;
 4975:   linei=0;
 4976:   while ((fgets(line, MAXLINE, fic) != NULL) &&((i >= firstobs) && (i <=lastobs))) {
 4977:     linei=linei+1;
 4978:     for(j=strlen(line); j>=0;j--){  /* Untabifies line */
 4979:       if(line[j] == '\t')
 4980: 	line[j] = ' ';
 4981:     }
 4982:     for(j=strlen(line)-1; (line[j]==' ')||(line[j]==10)||(line[j]==13);j--){
 4983:       ;
 4984:     };
 4985:     line[j+1]=0;  /* Trims blanks at end of line */
 4986:     if(line[0]=='#'){
 4987:       fprintf(ficlog,"Comment line\n%s\n",line);
 4988:       printf("Comment line\n%s\n",line);
 4989:       continue;
 4990:     }
 4991:     trimbb(linetmp,line); /* Trims multiple blanks in line */
 4992:     strcpy(line, linetmp);
 4993:   
 4994: 
 4995:     for (j=maxwav;j>=1;j--){
 4996:       cutv(stra, strb, line, ' '); 
 4997:       if(strb[0]=='.') { /* Missing status */
 4998: 	lval=-1;
 4999:       }else{
 5000: 	errno=0;
 5001: 	lval=strtol(strb,&endptr,10); 
 5002:       /*	if (errno == ERANGE && (lval == LONG_MAX || lval == LONG_MIN))*/
 5003: 	if( strb[0]=='\0' || (*endptr != '\0')){
 5004: 	  printf("Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a status of wave %d. Setting maxwav=%d might be wrong.  Exiting.\n", strb, linei,i,line,j,maxwav);
 5005: 	  fprintf(ficlog,"Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a status of wave %d. Setting maxwav=%d might be wrong.  Exiting.\n", strb, linei,i,line,j,maxwav);fflush(ficlog);
 5006: 	  return 1;
 5007: 	}
 5008:       }
 5009:       s[j][i]=lval;
 5010:       
 5011:       strcpy(line,stra);
 5012:       cutv(stra, strb,line,' ');
 5013:       if(iout=sscanf(strb,"%d/%d",&month, &year) != 0){
 5014:       }
 5015:       else  if(iout=sscanf(strb,"%s.",dummy) != 0){
 5016: 	month=99;
 5017: 	year=9999;
 5018:       }else{
 5019: 	printf("Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of interview (mm/yyyy or .) at wave %d.  Exiting.\n",strb, linei,i, line,j);
 5020: 	fprintf(ficlog,"Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of interview (mm/yyyy or .) at wave %d.  Exiting.\n",strb, linei,i, line,j);fflush(ficlog);
 5021: 	return 1;
 5022:       }
 5023:       anint[j][i]= (double) year; 
 5024:       mint[j][i]= (double)month; 
 5025:       strcpy(line,stra);
 5026:     } /* ENd Waves */
 5027:     
 5028:     cutv(stra, strb,line,' '); 
 5029:     if(iout=sscanf(strb,"%d/%d",&month, &year) != 0){
 5030:     }
 5031:     else  if(iout=sscanf(strb,"%s.",dummy) != 0){
 5032:       month=99;
 5033:       year=9999;
 5034:     }else{
 5035:       printf("Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of death (mm/yyyy or .).  Exiting.\n",strb, linei,i,line);
 5036: 	fprintf(ficlog,"Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of death (mm/yyyy or .).  Exiting.\n",strb, linei,i,line);fflush(ficlog);
 5037: 	return 1;
 5038:     }
 5039:     andc[i]=(double) year; 
 5040:     moisdc[i]=(double) month; 
 5041:     strcpy(line,stra);
 5042:     
 5043:     cutv(stra, strb,line,' '); 
 5044:     if(iout=sscanf(strb,"%d/%d",&month, &year) != 0){
 5045:     }
 5046:     else  if(iout=sscanf(strb,"%s.", dummy) != 0){
 5047:       month=99;
 5048:       year=9999;
 5049:     }else{
 5050:       printf("Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of birth (mm/yyyy or .).  Exiting.\n",strb, linei,i,line);
 5051:       fprintf(ficlog,"Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of birth (mm/yyyy or .).  Exiting.\n",strb, linei,i,line);fflush(ficlog);
 5052: 	return 1;
 5053:     }
 5054:     if (year==9999) {
 5055:       printf("Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of birth (mm/yyyy) but at least the year of birth should be given.  Exiting.\n",strb, linei,i,line);
 5056:       fprintf(ficlog,"Error reading data around '%s' at line number %d for individual %d, '%s'\nShould be a date of birth (mm/yyyy) but at least the year of birth should be given. Exiting.\n",strb, linei,i,line);fflush(ficlog);
 5057: 	return 1;
 5058: 
 5059:     }
 5060:     annais[i]=(double)(year);
 5061:     moisnais[i]=(double)(month); 
 5062:     strcpy(line,stra);
 5063:     
 5064:     cutv(stra, strb,line,' '); 
 5065:     errno=0;
 5066:     dval=strtod(strb,&endptr); 
 5067:     if( strb[0]=='\0' || (*endptr != '\0')){
 5068:       printf("Error reading data around '%f' at line number %d, \"%s\" for individual %d\nShould be a weight.  Exiting.\n",dval, i,line,linei);
 5069:       fprintf(ficlog,"Error reading data around '%f' at line number %d, \"%s\" for individual %d\nShould be a weight.  Exiting.\n",dval, i,line,linei);
 5070:       fflush(ficlog);
 5071:       return 1;
 5072:     }
 5073:     weight[i]=dval; 
 5074:     strcpy(line,stra);
 5075:     
 5076:     for (j=ncovcol;j>=1;j--){
 5077:       cutv(stra, strb,line,' '); 
 5078:       if(strb[0]=='.') { /* Missing status */
 5079: 	lval=-1;
 5080:       }else{
 5081: 	errno=0;
 5082: 	lval=strtol(strb,&endptr,10); 
 5083: 	if( strb[0]=='\0' || (*endptr != '\0')){
 5084: 	  printf("Error reading data around '%ld' at line number %d for individual %d, '%s'\nShould be a covariate value (=0 for the reference or 1 for alternative).  Exiting.\n",lval, linei,i, line);
 5085: 	  fprintf(ficlog,"Error reading data around '%ld' at line number %d for individual %d, '%s'\nShould be a covariate value (=0 for the reference or 1 for alternative).  Exiting.\n",lval, linei,i, line);fflush(ficlog);
 5086: 	  return 1;
 5087: 	}
 5088:       }
 5089:       if(lval <-1 || lval >1){
 5090: 	printf("Error reading data around '%ld' at line number %d for individual %d, '%s'\n \
 5091:  Should be a value of %d(nth) covariate (0 should be the value for the reference and 1\n \
 5092:  for the alternative. IMaCh does not build design variables automatically, do it yourself.\n \
 5093:  For example, for multinomial values like 1, 2 and 3,\n \
 5094:  build V1=0 V2=0 for the reference value (1),\n \
 5095:         V1=1 V2=0 for (2) \n \
 5096:  and V1=0 V2=1 for (3). V1=1 V2=1 should not exist and the corresponding\n \
 5097:  output of IMaCh is often meaningless.\n \
 5098:  Exiting.\n",lval,linei, i,line,j);
 5099: 	fprintf(ficlog,"Error reading data around '%ld' at line number %d for individual %d, '%s'\n \
 5100:  Should be a value of %d(nth) covariate (0 should be the value for the reference and 1\n \
 5101:  for the alternative. IMaCh does not build design variables automatically, do it yourself.\n \
 5102:  For example, for multinomial values like 1, 2 and 3,\n \
 5103:  build V1=0 V2=0 for the reference value (1),\n \
 5104:         V1=1 V2=0 for (2) \n \
 5105:  and V1=0 V2=1 for (3). V1=1 V2=1 should not exist and the corresponding\n \
 5106:  output of IMaCh is often meaningless.\n \
 5107:  Exiting.\n",lval,linei, i,line,j);fflush(ficlog);
 5108: 	return 1;
 5109:       }
 5110:       covar[j][i]=(double)(lval);
 5111:       strcpy(line,stra);
 5112:     }  
 5113:     lstra=strlen(stra);
 5114:      
 5115:     if(lstra > 9){ /* More than 2**32 or max of what printf can write with %ld */
 5116:       stratrunc = &(stra[lstra-9]);
 5117:       num[i]=atol(stratrunc);
 5118:     }
 5119:     else
 5120:       num[i]=atol(stra);
 5121:     /*if((s[2][i]==2) && (s[3][i]==-1)&&(s[4][i]==9)){
 5122:       printf("%ld %.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;}*/
 5123:     
 5124:     i=i+1;
 5125:   } /* End loop reading  data */
 5126: 
 5127:   *imax=i-1; /* Number of individuals */
 5128:   fclose(fic);
 5129:  
 5130:   return (0);
 5131:   /* endread: */
 5132:     printf("Exiting readdata: ");
 5133:     fclose(fic);
 5134:     return (1);
 5135: 
 5136: 
 5137: 
 5138: }
 5139: void removespace(char *str) {
 5140:   char *p1 = str, *p2 = str;
 5141:   do
 5142:     while (*p2 == ' ')
 5143:       p2++;
 5144:   while (*p1++ = *p2++);
 5145: }
 5146: 
 5147: int decodemodel ( char model[], int lastobs) /**< This routine decode the model and returns:
 5148:    * Model  V1+V2+V3+V8+V7*V8+V5*V6+V8*age+V3*age
 5149:    * - cptcovt total number of covariates of the model nbocc(+)+1 = 8
 5150:    * - cptcovn or number of covariates k of the models excluding age*products =6
 5151:    * - cptcovage number of covariates with age*products =2
 5152:    * - cptcovs number of simple covariates
 5153:    * - Tvar[k] is the id of the kth covariate Tvar[1]@12 {1, 2, 3, 8, 10, 11, 8, 3, 7, 8, 5, 6}, thus Tvar[5=V7*V8]=10
 5154:    *     which is a new column after the 9 (ncovcol) variables. 
 5155:    * - if k is a product Vn*Vm covar[k][i] is filled with correct values for each individual
 5156:    * - Tprod[l] gives the kth covariates of the product Vn*Vm l=1 to cptcovprod-cptcovage
 5157:    *    Tprod[1]@2 {5, 6}: position of first product V7*V8 is 5, and second V5*V6 is 6.
 5158:    * - Tvard[k]  p Tvard[1][1]@4 {7, 8, 5, 6} for V7*V8 and V5*V6 .
 5159:  */
 5160: {
 5161:   int i, j, k, ks;
 5162:   int  j1, k1, k2;
 5163:   char modelsav[80];
 5164:   char stra[80], strb[80], strc[80], strd[80],stre[80];
 5165: 
 5166:   /*removespace(model);*/
 5167:   if (strlen(model) >1){ /* If there is at least 1 covariate */
 5168:     j=0, j1=0, k1=0, k2=-1, ks=0, cptcovn=0;
 5169:     j=nbocc(model,'+'); /**< j=Number of '+' */
 5170:     j1=nbocc(model,'*'); /**< j1=Number of '*' */
 5171:     cptcovs=j+1-j1; /**<  Number of simple covariates V1+V2*age+V3 +V3*V4=> V1 + V3 =2  */
 5172:     cptcovt= j+1; /* Number of total covariates in the model V1 + V2*age+ V3 + V3*V4=> 4*/
 5173:                   /* including age products which are counted in cptcovage.
 5174: 		  /* but the covariates which are products must be treated separately: ncovn=4- 2=2 (V1+V3). */
 5175:     cptcovprod=j1; /**< Number of products  V1*V2 +v3*age = 2 */
 5176:     cptcovprodnoage=0; /**< Number of covariate products without age: V3*V4 =1  */
 5177:     strcpy(modelsav,model); 
 5178:     if (strstr(model,"AGE") !=0){
 5179:       printf("Error. AGE must be in lower case 'age' model=%s ",model);
 5180:       fprintf(ficlog,"Error. AGE must be in lower case model=%s ",model);fflush(ficlog);
 5181:       return 1;
 5182:     }
 5183:     if (strstr(model,"v") !=0){
 5184:       printf("Error. 'v' must be in upper case 'V' model=%s ",model);
 5185:       fprintf(ficlog,"Error. 'v' must be in upper case model=%s ",model);fflush(ficlog);
 5186:       return 1;
 5187:     }
 5188:     
 5189:     /*   Design
 5190:      *  V1   V2   V3   V4  V5  V6  V7  V8  V9 Weight
 5191:      *  <          ncovcol=8                >
 5192:      * Model V2 + V1 + V3*age + V3 + V5*V6 + V7*V8 + V8*age + V8
 5193:      *   k=  1    2      3       4     5       6      7        8
 5194:      *  cptcovn number of covariates (not including constant and age ) = # of + plus 1 = 7+1=8
 5195:      *  covar[k,i], value of kth covariate if not including age for individual i:
 5196:      *       covar[1][i]= (V2), covar[4][i]=(V3), covar[8][i]=(V8)
 5197:      *  Tvar[k] # of the kth covariate:  Tvar[1]=2  Tvar[4]=3 Tvar[8]=8
 5198:      *       if multiplied by age: V3*age Tvar[3=V3*age]=3 (V3) Tvar[7]=8 and 
 5199:      *  Tage[++cptcovage]=k
 5200:      *       if products, new covar are created after ncovcol with k1
 5201:      *  Tvar[k]=ncovcol+k1; # of the kth covariate product:  Tvar[5]=ncovcol+1=10  Tvar[6]=ncovcol+1=11
 5202:      *  Tprod[k1]=k; Tprod[1]=5 Tprod[2]= 6; gives the position of the k1th product
 5203:      *  Tvard[k1][1]=m Tvard[k1][2]=m; Tvard[1][1]=5 (V5) Tvard[1][2]=6 Tvard[2][1]=7 (V7) Tvard[2][2]=8
 5204:      *  Tvar[cptcovn+k2]=Tvard[k1][1];Tvar[cptcovn+k2+1]=Tvard[k1][2];
 5205:      *  Tvar[8+1]=5;Tvar[8+2]=6;Tvar[8+3]=7;Tvar[8+4]=8 inverted
 5206:      *  V1   V2   V3   V4  V5  V6  V7  V8  V9  V10  V11
 5207:      *  <          ncovcol=8                >
 5208:      *       Model V2 + V1 + V3*age + V3 + V5*V6 + V7*V8 + V8*age + V8    d1   d1   d2  d2
 5209:      *          k=  1    2      3       4     5       6      7        8    9   10   11  12
 5210:      *     Tvar[k]= 2    1      3       3    10      11      8        8    5    6    7   8
 5211:      * p Tvar[1]@12={2,   1,     3,      3,   11,     10,     8,       8,   7,   8,   5,  6}
 5212:      * p Tprod[1]@2={                         6, 5}
 5213:      *p Tvard[1][1]@4= {7, 8, 5, 6}
 5214:      * covar[k][i]= V2   V1      ?      V3    V5*V6?   V7*V8?  ?       V8   
 5215:      *  cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 5216:      *How to reorganize?
 5217:      * Model V1 + V2 + V3 + V8 + V5*V6 + V7*V8 + V3*age + V8*age
 5218:      * Tvars {2,   1,     3,      3,   11,     10,     8,       8,   7,   8,   5,  6}
 5219:      *       {2,   1,     4,      8,    5,      6,     3,       7}
 5220:      * Struct []
 5221:      */
 5222: 
 5223:     /* This loop fills the array Tvar from the string 'model'.*/
 5224:     /* j is the number of + signs in the model V1+V2+V3 j=2 i=3 to 1 */
 5225:     /*   modelsav=V2+V1+V4+age*V3 strb=age*V3 stra=V2+V1+V4  */
 5226:     /* 	k=4 (age*V3) Tvar[k=4]= 3 (from V3) Tage[cptcovage=1]=4 */
 5227:     /* 	k=3 V4 Tvar[k=3]= 4 (from V4) */
 5228:     /* 	k=2 V1 Tvar[k=2]= 1 (from V1) */
 5229:     /* 	k=1 Tvar[1]=2 (from V2) */
 5230:     /* 	k=5 Tvar[5] */
 5231:     /* for (k=1; k<=cptcovn;k++) { */
 5232:     /* 	cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]]; */
 5233:     /* 	} */
 5234:     /* for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2]; */
 5235:     /*
 5236:      * Treating invertedly V2+V1+V3*age+V2*V4 is as if written V2*V4 +V3*age + V1 + V2 */
 5237:     for(k=cptcovt; k>=1;k--) /**< Number of covariates */
 5238:         Tvar[k]=0;
 5239:     cptcovage=0;
 5240:     for(k=1; k<=cptcovt;k++){ /* Loop on total covariates of the model */
 5241:       cutl(stra,strb,modelsav,'+'); /* keeps in strb after the first '+' 
 5242: 				     modelsav==V2+V1+V4+V3*age strb=V3*age stra=V2+V1+V4 */ 
 5243:       if (nbocc(modelsav,'+')==0) strcpy(strb,modelsav); /* and analyzes it */
 5244:       /*      printf("i=%d a=%s b=%s sav=%s\n",i, stra,strb,modelsav);*/
 5245:       /*scanf("%d",i);*/
 5246:       if (strchr(strb,'*')) {  /**< Model includes a product V2+V1+V4+V3*age strb=V3*age */
 5247: 	cutl(strc,strd,strb,'*'); /**< strd*strc  Vm*Vn: strb=V3*age(input) strc=age strd=V3 ; V3*V2 strc=V2, strd=V3 */
 5248: 	if (strcmp(strc,"age")==0) { /**< Model includes age: Vn*age */
 5249: 	  /* covar is not filled and then is empty */
 5250: 	  cptcovprod--;
 5251: 	  cutl(stre,strb,strd,'V'); /* strd=V3(input): stre="3" */
 5252: 	  Tvar[k]=atoi(stre);  /* V2+V1+V4+V3*age Tvar[4]=3 ; V1+V2*age Tvar[2]=2 */
 5253: 	  cptcovage++; /* Sums the number of covariates which include age as a product */
 5254: 	  Tage[cptcovage]=k;  /* Tage[1] = 4 */
 5255: 	  /*printf("stre=%s ", stre);*/
 5256: 	} else if (strcmp(strd,"age")==0) { /* or age*Vn */
 5257: 	  cptcovprod--;
 5258: 	  cutl(stre,strb,strc,'V');
 5259: 	  Tvar[k]=atoi(stre);
 5260: 	  cptcovage++;
 5261: 	  Tage[cptcovage]=k;
 5262: 	} else {  /* Age is not in the model product V2+V1+V1*V4+V3*age+V3*V2  strb=V3*V2*/
 5263: 	  /* loops on k1=1 (V3*V2) and k1=2 V4*V3 */
 5264: 	  cptcovn++;
 5265: 	  cptcovprodnoage++;k1++;
 5266: 	  cutl(stre,strb,strc,'V'); /* strc= Vn, stre is n; strb=V3*V2 stre=3 strc=*/
 5267: 	  Tvar[k]=ncovcol+k1; /* For model-covariate k tells which data-covariate to use but
 5268: 				  because this model-covariate is a construction we invent a new column
 5269: 				  ncovcol + k1
 5270: 				  If already ncovcol=4 and model=V2+V1+V1*V4+age*V3+V3*V2
 5271: 				  Tvar[3=V1*V4]=4+1 Tvar[5=V3*V2]=4 + 2= 6, etc */
 5272: 	  cutl(strc,strb,strd,'V'); /* strd was Vm, strc is m */
 5273: 	  Tprod[k1]=k;  /* Tprod[1]=3(=V1*V4) for V2+V1+V1*V4+age*V3+V3*V2  */
 5274: 	  Tvard[k1][1] =atoi(strc); /* m 1 for V1*/
 5275: 	  Tvard[k1][2] =atoi(stre); /* n 4 for V4*/
 5276: 	  k2=k2+2;
 5277: 	  Tvar[cptcovt+k2]=Tvard[k1][1]; /* Tvar[(cptcovt=4+k2=1)=5]= 1 (V1) */
 5278: 	  Tvar[cptcovt+k2+1]=Tvard[k1][2];  /* Tvar[(cptcovt=4+(k2=1)+1)=6]= 4 (V4) */
 5279: 	  for (i=1; i<=lastobs;i++){
 5280: 	    /* Computes the new covariate which is a product of
 5281: 	       covar[n][i]* covar[m][i] and stores it at ncovol+k1 May not be defined */
 5282: 	    covar[ncovcol+k1][i]=covar[atoi(stre)][i]*covar[atoi(strc)][i];
 5283: 	  }
 5284: 	} /* End age is not in the model */
 5285:       } /* End if model includes a product */
 5286:       else { /* no more sum */
 5287: 	/*printf("d=%s c=%s b=%s\n", strd,strc,strb);*/
 5288:        /*  scanf("%d",i);*/
 5289: 	cutl(strd,strc,strb,'V');
 5290: 	ks++; /**< Number of simple covariates */
 5291: 	cptcovn++;
 5292: 	Tvar[k]=atoi(strd);
 5293:       }
 5294:       strcpy(modelsav,stra);  /* modelsav=V2+V1+V4 stra=V2+V1+V4 */ 
 5295:       /*printf("a=%s b=%s sav=%s\n", stra,strb,modelsav);
 5296: 	scanf("%d",i);*/
 5297:     } /* end of loop + */
 5298:   } /* end model */
 5299:   
 5300:   /*The number n of Vn is stored in Tvar. cptcovage =number of age covariate. Tage gives the position of age. cptcovprod= number of products.
 5301:     If model=V1+V1*age then Tvar[1]=1 Tvar[2]=1 cptcovage=1 Tage[1]=2 cptcovprod=0*/
 5302: 
 5303:   /* printf("tvar1=%d tvar2=%d tvar3=%d cptcovage=%d Tage=%d",Tvar[1],Tvar[2],Tvar[3],cptcovage,Tage[1]);
 5304:   printf("cptcovprod=%d ", cptcovprod);
 5305:   fprintf(ficlog,"cptcovprod=%d ", cptcovprod);
 5306: 
 5307:   scanf("%d ",i);*/
 5308: 
 5309: 
 5310:   return (0); /* with covar[new additional covariate if product] and Tage if age */ 
 5311:   /*endread:*/
 5312:     printf("Exiting decodemodel: ");
 5313:     return (1);
 5314: }
 5315: 
 5316: calandcheckages(int imx, int maxwav, double *agemin, double *agemax, int *nberr, int *nbwarn )
 5317: {
 5318:   int i, m;
 5319: 
 5320:   for (i=1; i<=imx; i++) {
 5321:     for(m=2; (m<= maxwav); m++) {
 5322:       if (((int)mint[m][i]== 99) && (s[m][i] <= nlstate)){
 5323: 	anint[m][i]=9999;
 5324: 	s[m][i]=-1;
 5325:       }
 5326:       if((int)moisdc[i]==99 && (int)andc[i]==9999 && s[m][i]>nlstate){
 5327: 	*nberr++;
 5328: 	printf("Error! Date of death (month %2d and year %4d) of individual %ld on line %d was unknown, you must set an arbitrary year of death or he/she is skipped and results are biased\n",(int)moisdc[i],(int)andc[i],num[i],i);
 5329: 	fprintf(ficlog,"Error! Date of death (month %2d and year %4d) of individual %ld on line %d was unknown, you must set an arbitrary year of death or he/she is skipped and results are biased\n",(int)moisdc[i],(int)andc[i],num[i],i);
 5330: 	s[m][i]=-1;
 5331:       }
 5332:       if((int)moisdc[i]==99 && (int)andc[i]!=9999 && s[m][i]>nlstate){
 5333: 	*nberr++;
 5334: 	printf("Error! Month of death of individual %ld on line %d was unknown %2d, you should set it otherwise the information on the death is skipped and results are biased.\n",num[i],i,(int)moisdc[i]); 
 5335: 	fprintf(ficlog,"Error! Month of death of individual %ld on line %d was unknown %f, you should set it otherwise the information on the death is skipped and results are biased.\n",num[i],i,moisdc[i]); 
 5336: 	s[m][i]=-1; /* We prefer to skip it (and to skip it in version 0.8a1 too */
 5337:       }
 5338:     }
 5339:   }
 5340: 
 5341:   for (i=1; i<=imx; i++)  {
 5342:     agedc[i]=(moisdc[i]/12.+andc[i])-(moisnais[i]/12.+annais[i]);
 5343:     for(m=firstpass; (m<= lastpass); m++){
 5344:       if(s[m][i] >0 || s[m][i]==-2 || s[m][i]==-4 || s[m][i]==-5){
 5345: 	if (s[m][i] >= nlstate+1) {
 5346: 	  if(agedc[i]>0)
 5347: 	    if((int)moisdc[i]!=99 && (int)andc[i]!=9999)
 5348: 	      agev[m][i]=agedc[i];
 5349: 	  /*if(moisdc[i]==99 && andc[i]==9999) s[m][i]=-1;*/
 5350: 	    else {
 5351: 	      if ((int)andc[i]!=9999){
 5352: 		nbwarn++;
 5353: 		printf("Warning negative age at death: %ld line:%d\n",num[i],i);
 5354: 		fprintf(ficlog,"Warning negative age at death: %ld line:%d\n",num[i],i);
 5355: 		agev[m][i]=-1;
 5356: 	      }
 5357: 	    }
 5358: 	}
 5359: 	else if(s[m][i] !=9){ /* Standard case, age in fractional
 5360: 				 years but with the precision of a month */
 5361: 	  agev[m][i]=(mint[m][i]/12.+1./24.+anint[m][i])-(moisnais[i]/12.+1./24.+annais[i]);
 5362: 	  if((int)mint[m][i]==99 || (int)anint[m][i]==9999)
 5363: 	    agev[m][i]=1;
 5364: 	  else if(agev[m][i] < *agemin){ 
 5365: 	    *agemin=agev[m][i];
 5366: 	    printf(" Min anint[%d][%d]=%.2f annais[%d]=%.2f, agemin=%.2f\n",m,i,anint[m][i], i,annais[i], *agemin);
 5367: 	  }
 5368: 	  else if(agev[m][i] >*agemax){
 5369: 	    *agemax=agev[m][i];
 5370: 	    /* printf(" Max anint[%d][%d]=%.0f annais[%d]=%.0f, agemax=%.2f\n",m,i,anint[m][i], i,annais[i], *agemax);*/
 5371: 	  }
 5372: 	  /*agev[m][i]=anint[m][i]-annais[i];*/
 5373: 	  /*	 agev[m][i] = age[i]+2*m;*/
 5374: 	}
 5375: 	else { /* =9 */
 5376: 	  agev[m][i]=1;
 5377: 	  s[m][i]=-1;
 5378: 	}
 5379:       }
 5380:       else /*= 0 Unknown */
 5381: 	agev[m][i]=1;
 5382:     }
 5383:     
 5384:   }
 5385:   for (i=1; i<=imx; i++)  {
 5386:     for(m=firstpass; (m<=lastpass); m++){
 5387:       if (s[m][i] > (nlstate+ndeath)) {
 5388: 	*nberr++;
 5389: 	printf("Error: on wave %d of individual %d status %d > (nlstate+ndeath)=(%d+%d)=%d\n",m,i,s[m][i],nlstate, ndeath, nlstate+ndeath);	
 5390: 	fprintf(ficlog,"Error: on wave %d of individual %d status %d > (nlstate+ndeath)=(%d+%d)=%d\n",m,i,s[m][i],nlstate, ndeath, nlstate+ndeath);	
 5391: 	return 1;
 5392:       }
 5393:     }
 5394:   }
 5395: 
 5396:   /*for (i=1; i<=imx; i++){
 5397:   for (m=firstpass; (m<lastpass); m++){
 5398:      printf("%ld %d %.lf %d %d\n", num[i],(covar[1][i]),agev[m][i],s[m][i],s[m+1][i]);
 5399: }
 5400: 
 5401: }*/
 5402: 
 5403: 
 5404:   printf("Total number of individuals= %d, Agemin = %.2f, Agemax= %.2f\n\n", imx, *agemin, *agemax);
 5405:   fprintf(ficlog,"Total number of individuals= %d, Agemin = %.2f, Agemax= %.2f\n\n", imx, *agemin, *agemax); 
 5406: 
 5407:   return (0);
 5408:  /* endread:*/
 5409:     printf("Exiting calandcheckages: ");
 5410:     return (1);
 5411: }
 5412: 
 5413: 
 5414: /***********************************************/
 5415: /**************** Main Program *****************/
 5416: /***********************************************/
 5417: 
 5418: int main(int argc, char *argv[])
 5419: {
 5420: #ifdef GSL
 5421:   const gsl_multimin_fminimizer_type *T;
 5422:   size_t iteri = 0, it;
 5423:   int rval = GSL_CONTINUE;
 5424:   int status = GSL_SUCCESS;
 5425:   double ssval;
 5426: #endif
 5427:   int movingaverage(double ***probs, double bage,double fage, double ***mobaverage, int mobilav);
 5428:   int i,j, k, n=MAXN,iter=0,m,size=100, cptcod;
 5429: 
 5430:   int jj, ll, li, lj, lk;
 5431:   int numlinepar=0; /* Current linenumber of parameter file */
 5432:   int itimes;
 5433:   int NDIM=2;
 5434:   int vpopbased=0;
 5435: 
 5436:   char ca[32], cb[32];
 5437:   /*  FILE *fichtm; *//* Html File */
 5438:   /* FILE *ficgp;*/ /*Gnuplot File */
 5439:   struct stat info;
 5440:   double agedeb;
 5441:   double ageminpar=1.e20,agemin=1.e20, agemaxpar=-1.e20, agemax=-1.e20;
 5442: 
 5443:   double fret;
 5444:   double dum; /* Dummy variable */
 5445:   double ***p3mat;
 5446:   double ***mobaverage;
 5447: 
 5448:   char line[MAXLINE];
 5449:   char path[MAXLINE],pathc[MAXLINE],pathcd[MAXLINE],pathtot[MAXLINE],model[MAXLINE];
 5450:   char pathr[MAXLINE], pathimach[MAXLINE]; 
 5451:   char *tok, *val; /* pathtot */
 5452:   int firstobs=1, lastobs=10;
 5453:   int c,  h , cpt;
 5454:   int jl;
 5455:   int i1, j1, jk, stepsize;
 5456:   int *tab; 
 5457:   int mobilavproj=0 , prevfcast=0 ; /* moving average of prev, If prevfcast=1 prevalence projection */
 5458:   int mobilav=0,popforecast=0;
 5459:   int hstepm, nhstepm;
 5460:   int agemortsup;
 5461:   float  sumlpop=0.;
 5462:   double jprev1=1, mprev1=1,anprev1=2000,jprev2=1, mprev2=1,anprev2=2000;
 5463:   double jpyram=1, mpyram=1,anpyram=2000,jpyram1=1, mpyram1=1,anpyram1=2000;
 5464: 
 5465:   double bage=0, fage=110, age, agelim, agebase;
 5466:   double ftolpl=FTOL;
 5467:   double **prlim;
 5468:   double ***param; /* Matrix of parameters */
 5469:   double  *p;
 5470:   double **matcov; /* Matrix of covariance */
 5471:   double ***delti3; /* Scale */
 5472:   double *delti; /* Scale */
 5473:   double ***eij, ***vareij;
 5474:   double **varpl; /* Variances of prevalence limits by age */
 5475:   double *epj, vepp;
 5476: 
 5477:   double dateprev1, dateprev2,jproj1=1,mproj1=1,anproj1=2000,jproj2=1,mproj2=1,anproj2=2000;
 5478:   double **ximort;
 5479:   char *alph[]={"a","a","b","c","d","e"}, str[4]="1234";
 5480:   int *dcwave;
 5481: 
 5482:   char z[1]="c";
 5483: 
 5484:   /*char  *strt;*/
 5485:   char strtend[80];
 5486: 
 5487: 
 5488: /*   setlocale (LC_ALL, ""); */
 5489: /*   bindtextdomain (PACKAGE, LOCALEDIR); */
 5490: /*   textdomain (PACKAGE); */
 5491: /*   setlocale (LC_CTYPE, ""); */
 5492: /*   setlocale (LC_MESSAGES, ""); */
 5493: 
 5494:   /*   gettimeofday(&start_time, (struct timezone*)0); */ /* at first time */
 5495:   rstart_time = time(NULL);  
 5496:   /*  (void) gettimeofday(&start_time,&tzp);*/
 5497:   start_time = *localtime(&rstart_time);
 5498:   curr_time=start_time;
 5499:   /*tml = *localtime(&start_time.tm_sec);*/
 5500:   /* strcpy(strstart,asctime(&tml)); */
 5501:   strcpy(strstart,asctime(&start_time));
 5502: 
 5503: /*  printf("Localtime (at start)=%s",strstart); */
 5504: /*  tp.tm_sec = tp.tm_sec +86400; */
 5505: /*  tm = *localtime(&start_time.tm_sec); */
 5506: /*   tmg.tm_year=tmg.tm_year +dsign*dyear; */
 5507: /*   tmg.tm_mon=tmg.tm_mon +dsign*dmonth; */
 5508: /*   tmg.tm_hour=tmg.tm_hour + 1; */
 5509: /*   tp.tm_sec = mktime(&tmg); */
 5510: /*   strt=asctime(&tmg); */
 5511: /*   printf("Time(after) =%s",strstart);  */
 5512: /*  (void) time (&time_value);
 5513: *  printf("time=%d,t-=%d\n",time_value,time_value-86400);
 5514: *  tm = *localtime(&time_value);
 5515: *  strstart=asctime(&tm);
 5516: *  printf("tim_value=%d,asctime=%s\n",time_value,strstart); 
 5517: */
 5518: 
 5519:   nberr=0; /* Number of errors and warnings */
 5520:   nbwarn=0;
 5521:   getcwd(pathcd, size);
 5522: 
 5523:   printf("\n%s\n%s",version,fullversion);
 5524:   if(argc <=1){
 5525:     printf("\nEnter the parameter file name: ");
 5526:     fgets(pathr,FILENAMELENGTH,stdin);
 5527:     i=strlen(pathr);
 5528:     if(pathr[i-1]=='\n')
 5529:       pathr[i-1]='\0';
 5530:     i=strlen(pathr);
 5531:     if(pathr[i-1]==' ') /* This may happen when dragging on oS/X! */
 5532:       pathr[i-1]='\0';
 5533:    for (tok = pathr; tok != NULL; ){
 5534:       printf("Pathr |%s|\n",pathr);
 5535:       while ((val = strsep(&tok, "\"" )) != NULL && *val == '\0');
 5536:       printf("val= |%s| pathr=%s\n",val,pathr);
 5537:       strcpy (pathtot, val);
 5538:       if(pathr[0] == '\0') break; /* Dirty */
 5539:     }
 5540:   }
 5541:   else{
 5542:     strcpy(pathtot,argv[1]);
 5543:   }
 5544:   /*if(getcwd(pathcd, MAXLINE)!= NULL)printf ("Error pathcd\n");*/
 5545:   /*cygwin_split_path(pathtot,path,optionfile);
 5546:     printf("pathtot=%s, path=%s, optionfile=%s\n",pathtot,path,optionfile);*/
 5547:   /* cutv(path,optionfile,pathtot,'\\');*/
 5548: 
 5549:   /* Split argv[0], imach program to get pathimach */
 5550:   printf("\nargv[0]=%s argv[1]=%s, \n",argv[0],argv[1]);
 5551:   split(argv[0],pathimach,optionfile,optionfilext,optionfilefiname);
 5552:   printf("\nargv[0]=%s pathimach=%s, \noptionfile=%s \noptionfilext=%s \noptionfilefiname=%s\n",argv[0],pathimach,optionfile,optionfilext,optionfilefiname);
 5553:  /*   strcpy(pathimach,argv[0]); */
 5554:   /* Split argv[1]=pathtot, parameter file name to get path, optionfile, extension and name */
 5555:   split(pathtot,path,optionfile,optionfilext,optionfilefiname);
 5556:   printf("\npathtot=%s,\npath=%s,\noptionfile=%s \noptionfilext=%s \noptionfilefiname=%s\n",pathtot,path,optionfile,optionfilext,optionfilefiname);
 5557:   chdir(path); /* Can be a relative path */
 5558:   if(getcwd(pathcd,MAXLINE) > 0) /* So pathcd is the full path */
 5559:     printf("Current directory %s!\n",pathcd);
 5560:   strcpy(command,"mkdir ");
 5561:   strcat(command,optionfilefiname);
 5562:   if((outcmd=system(command)) != 0){
 5563:     printf("Problem creating directory or it already exists %s%s, err=%d\n",path,optionfilefiname,outcmd);
 5564:     /* fprintf(ficlog,"Problem creating directory %s%s\n",path,optionfilefiname); */
 5565:     /* fclose(ficlog); */
 5566: /*     exit(1); */
 5567:   }
 5568: /*   if((imk=mkdir(optionfilefiname))<0){ */
 5569: /*     perror("mkdir"); */
 5570: /*   } */
 5571: 
 5572:   /*-------- arguments in the command line --------*/
 5573: 
 5574:   /* Log file */
 5575:   strcat(filelog, optionfilefiname);
 5576:   strcat(filelog,".log");    /* */
 5577:   if((ficlog=fopen(filelog,"w"))==NULL)    {
 5578:     printf("Problem with logfile %s\n",filelog);
 5579:     goto end;
 5580:   }
 5581:   fprintf(ficlog,"Log filename:%s\n",filelog);
 5582:   fprintf(ficlog,"\n%s\n%s",version,fullversion);
 5583:   fprintf(ficlog,"\nEnter the parameter file name: \n");
 5584:   fprintf(ficlog,"pathimach=%s\npathtot=%s\n\
 5585:  path=%s \n\
 5586:  optionfile=%s\n\
 5587:  optionfilext=%s\n\
 5588:  optionfilefiname='%s'\n",pathimach,pathtot,path,optionfile,optionfilext,optionfilefiname);
 5589: 
 5590:   printf("Local time (at start):%s",strstart);
 5591:   fprintf(ficlog,"Local time (at start): %s",strstart);
 5592:   fflush(ficlog);
 5593: /*   (void) gettimeofday(&curr_time,&tzp); */
 5594: /*   printf("Elapsed time %d\n", asc_diff_time(curr_time.tm_sec-start_time.tm_sec,tmpout)); */
 5595: 
 5596:   /* */
 5597:   strcpy(fileres,"r");
 5598:   strcat(fileres, optionfilefiname);
 5599:   strcat(fileres,".txt");    /* Other files have txt extension */
 5600: 
 5601:   /*---------arguments file --------*/
 5602: 
 5603:   if((ficpar=fopen(optionfile,"r"))==NULL)    {
 5604:     printf("Problem with optionfile '%s' with errno='%s'\n",optionfile,strerror(errno));
 5605:     fprintf(ficlog,"Problem with optionfile '%s' with errno='%s'\n",optionfile,strerror(errno));
 5606:     fflush(ficlog);
 5607:     /* goto end; */
 5608:     exit(70); 
 5609:   }
 5610: 
 5611: 
 5612: 
 5613:   strcpy(filereso,"o");
 5614:   strcat(filereso,fileres);
 5615:   if((ficparo=fopen(filereso,"w"))==NULL) { /* opened on subdirectory */
 5616:     printf("Problem with Output resultfile: %s\n", filereso);
 5617:     fprintf(ficlog,"Problem with Output resultfile: %s\n", filereso);
 5618:     fflush(ficlog);
 5619:     goto end;
 5620:   }
 5621: 
 5622:   /* Reads comments: lines beginning with '#' */
 5623:   numlinepar=0;
 5624:   while((c=getc(ficpar))=='#' && c!= EOF){
 5625:     ungetc(c,ficpar);
 5626:     fgets(line, MAXLINE, ficpar);
 5627:     numlinepar++;
 5628:     fputs(line,stdout);
 5629:     fputs(line,ficparo);
 5630:     fputs(line,ficlog);
 5631:   }
 5632:   ungetc(c,ficpar);
 5633: 
 5634:   fscanf(ficpar,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%lf stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d model=%s\n",title, datafile, &lastobs, &firstpass,&lastpass,&ftol, &stepm, &ncovcol, &nlstate,&ndeath, &maxwav, &mle, &weightopt,model);
 5635:   numlinepar++;
 5636:   printf("title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol, stepm, ncovcol, nlstate,ndeath, maxwav, mle, weightopt,model);
 5637:   fprintf(ficparo,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol,stepm,ncovcol,nlstate,ndeath,maxwav, mle, weightopt,model);
 5638:   fprintf(ficlog,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle=%d weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol,stepm,ncovcol,nlstate,ndeath,maxwav, mle, weightopt,model);
 5639:   fflush(ficlog);
 5640:   while((c=getc(ficpar))=='#' && c!= EOF){
 5641:     ungetc(c,ficpar);
 5642:     fgets(line, MAXLINE, ficpar);
 5643:     numlinepar++;
 5644:     fputs(line, stdout);
 5645:     //puts(line);
 5646:     fputs(line,ficparo);
 5647:     fputs(line,ficlog);
 5648:   }
 5649:   ungetc(c,ficpar);
 5650: 
 5651:    
 5652:   covar=matrix(0,NCOVMAX,1,n);  /**< used in readdata */
 5653:   cptcovn=0; /*Number of covariates, i.e. number of '+' in model statement plus one, indepently of n in Vn*/
 5654:   /* v1+v2+v3+v2*v4+v5*age makes cptcovn = 5
 5655:      v1+v2*age+v2*v3 makes cptcovn = 3
 5656:   */
 5657:   if (strlen(model)>1) 
 5658:     ncovmodel=2+nbocc(model,'+')+1; /*Number of variables including intercept and age = cptcovn + intercept + age : v1+v2+v3+v2*v4+v5*age makes 5+2=7*/
 5659:   else
 5660:     ncovmodel=2;
 5661:   nvar=ncovmodel-1; /* Suppressing age as a basic covariate */
 5662:   nforce= (nlstate+ndeath-1)*nlstate; /* Number of forces ij from state i to j */
 5663:   npar= nforce*ncovmodel; /* Number of parameters like aij*/
 5664:   if(npar >MAXPARM || nlstate >NLSTATEMAX || ndeath >NDEATHMAX || ncovmodel>NCOVMAX){
 5665:     printf("Too complex model for current IMaCh: npar=(nlstate+ndeath-1)*nlstate*ncovmodel=%d >= %d(MAXPARM) or nlstate=%d >= %d(NLSTATEMAX) or ndeath=%d >= %d(NDEATHMAX) or ncovmodel=(k+age+#of+signs)=%d(NCOVMAX) >= %d\n",npar, MAXPARM, nlstate, NLSTATEMAX, ndeath, NDEATHMAX, ncovmodel, NCOVMAX);
 5666:     fprintf(ficlog,"Too complex model for current IMaCh: %d >=%d(MAXPARM) or %d >=%d(NLSTATEMAX) or %d >=%d(NDEATHMAX) or %d(NCOVMAX) >=%d\n",npar, MAXPARM, nlstate, NLSTATEMAX, ndeath, NDEATHMAX, ncovmodel, NCOVMAX);
 5667:     fflush(stdout);
 5668:     fclose (ficlog);
 5669:     goto end;
 5670:   }
 5671:   delti3= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
 5672:   delti=delti3[1][1];
 5673:   /*delti=vector(1,npar); *//* Scale of each paramater (output from hesscov)*/
 5674:   if(mle==-1){ /* Print a wizard for help writing covariance matrix */
 5675:     prwizard(ncovmodel, nlstate, ndeath, model, ficparo);
 5676:     printf(" You choose mle=-1, look at file %s for a template of covariance matrix \n",filereso);
 5677:     fprintf(ficlog," You choose mle=-1, look at file %s for a template of covariance matrix \n",filereso);
 5678:     free_ma3x(delti3,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel); 
 5679:     fclose (ficparo);
 5680:     fclose (ficlog);
 5681:     goto end;
 5682:     exit(0);
 5683:   }
 5684:   else if(mle==-3) {
 5685:     prwizard(ncovmodel, nlstate, ndeath, model, ficparo);
 5686:     printf(" You choose mle=-3, look at file %s for a template of covariance matrix \n",filereso);
 5687:     fprintf(ficlog," You choose mle=-3, look at file %s for a template of covariance matrix \n",filereso);
 5688:     param= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
 5689:     matcov=matrix(1,npar,1,npar);
 5690:   }
 5691:   else{
 5692:     /* Read guessed parameters */
 5693:     /* Reads comments: lines beginning with '#' */
 5694:     while((c=getc(ficpar))=='#' && c!= EOF){
 5695:       ungetc(c,ficpar);
 5696:       fgets(line, MAXLINE, ficpar);
 5697:       numlinepar++;
 5698:       fputs(line,stdout);
 5699:       fputs(line,ficparo);
 5700:       fputs(line,ficlog);
 5701:     }
 5702:     ungetc(c,ficpar);
 5703:     
 5704:     param= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
 5705:     for(i=1; i <=nlstate; i++){
 5706:       j=0;
 5707:       for(jj=1; jj <=nlstate+ndeath; jj++){
 5708: 	if(jj==i) continue;
 5709: 	j++;
 5710: 	fscanf(ficpar,"%1d%1d",&i1,&j1);
 5711: 	if ((i1 != i) && (j1 != j)){
 5712: 	  printf("Error in line parameters number %d, %1d%1d instead of %1d%1d \n \
 5713: It might be a problem of design; if ncovcol and the model are correct\n \
 5714: run imach with mle=-1 to get a correct template of the parameter file.\n",numlinepar, i,j, i1, j1);
 5715: 	  exit(1);
 5716: 	}
 5717: 	fprintf(ficparo,"%1d%1d",i1,j1);
 5718: 	if(mle==1)
 5719: 	  printf("%1d%1d",i,j);
 5720: 	fprintf(ficlog,"%1d%1d",i,j);
 5721: 	for(k=1; k<=ncovmodel;k++){
 5722: 	  fscanf(ficpar," %lf",&param[i][j][k]);
 5723: 	  if(mle==1){
 5724: 	    printf(" %lf",param[i][j][k]);
 5725: 	    fprintf(ficlog," %lf",param[i][j][k]);
 5726: 	  }
 5727: 	  else
 5728: 	    fprintf(ficlog," %lf",param[i][j][k]);
 5729: 	  fprintf(ficparo," %lf",param[i][j][k]);
 5730: 	}
 5731: 	fscanf(ficpar,"\n");
 5732: 	numlinepar++;
 5733: 	if(mle==1)
 5734: 	  printf("\n");
 5735: 	fprintf(ficlog,"\n");
 5736: 	fprintf(ficparo,"\n");
 5737:       }
 5738:     }  
 5739:     fflush(ficlog);
 5740: 
 5741:     /* Reads scales values */
 5742:     p=param[1][1];
 5743:     
 5744:     /* Reads comments: lines beginning with '#' */
 5745:     while((c=getc(ficpar))=='#' && c!= EOF){
 5746:       ungetc(c,ficpar);
 5747:       fgets(line, MAXLINE, ficpar);
 5748:       numlinepar++;
 5749:       fputs(line,stdout);
 5750:       fputs(line,ficparo);
 5751:       fputs(line,ficlog);
 5752:     }
 5753:     ungetc(c,ficpar);
 5754: 
 5755:     for(i=1; i <=nlstate; i++){
 5756:       for(j=1; j <=nlstate+ndeath-1; j++){
 5757: 	fscanf(ficpar,"%1d%1d",&i1,&j1);
 5758: 	if ( (i1-i) * (j1-j) != 0){
 5759: 	  printf("Error in line parameters number %d, %1d%1d instead of %1d%1d \n",numlinepar, i,j, i1, j1);
 5760: 	  exit(1);
 5761: 	}
 5762: 	printf("%1d%1d",i,j);
 5763: 	fprintf(ficparo,"%1d%1d",i1,j1);
 5764: 	fprintf(ficlog,"%1d%1d",i1,j1);
 5765: 	for(k=1; k<=ncovmodel;k++){
 5766: 	  fscanf(ficpar,"%le",&delti3[i][j][k]);
 5767: 	  printf(" %le",delti3[i][j][k]);
 5768: 	  fprintf(ficparo," %le",delti3[i][j][k]);
 5769: 	  fprintf(ficlog," %le",delti3[i][j][k]);
 5770: 	}
 5771: 	fscanf(ficpar,"\n");
 5772: 	numlinepar++;
 5773: 	printf("\n");
 5774: 	fprintf(ficparo,"\n");
 5775: 	fprintf(ficlog,"\n");
 5776:       }
 5777:     }
 5778:     fflush(ficlog);
 5779: 
 5780:     /* Reads covariance matrix */
 5781:     delti=delti3[1][1];
 5782: 
 5783: 
 5784:     /* free_ma3x(delti3,1,nlstate,1,nlstate+ndeath-1,1,ncovmodel); */ /* Hasn't to to freed here otherwise delti is no more allocated */
 5785:   
 5786:     /* Reads comments: lines beginning with '#' */
 5787:     while((c=getc(ficpar))=='#' && c!= EOF){
 5788:       ungetc(c,ficpar);
 5789:       fgets(line, MAXLINE, ficpar);
 5790:       numlinepar++;
 5791:       fputs(line,stdout);
 5792:       fputs(line,ficparo);
 5793:       fputs(line,ficlog);
 5794:     }
 5795:     ungetc(c,ficpar);
 5796:   
 5797:     matcov=matrix(1,npar,1,npar);
 5798:     for(i=1; i <=npar; i++)
 5799:       for(j=1; j <=npar; j++) matcov[i][j]=0.;
 5800:       
 5801:     for(i=1; i <=npar; i++){
 5802:       fscanf(ficpar,"%s",str);
 5803:       if(mle==1)
 5804: 	printf("%s",str);
 5805:       fprintf(ficlog,"%s",str);
 5806:       fprintf(ficparo,"%s",str);
 5807:       for(j=1; j <=i; j++){
 5808: 	fscanf(ficpar," %le",&matcov[i][j]);
 5809: 	if(mle==1){
 5810: 	  printf(" %.5le",matcov[i][j]);
 5811: 	}
 5812: 	fprintf(ficlog," %.5le",matcov[i][j]);
 5813: 	fprintf(ficparo," %.5le",matcov[i][j]);
 5814:       }
 5815:       fscanf(ficpar,"\n");
 5816:       numlinepar++;
 5817:       if(mle==1)
 5818: 	printf("\n");
 5819:       fprintf(ficlog,"\n");
 5820:       fprintf(ficparo,"\n");
 5821:     }
 5822:     for(i=1; i <=npar; i++)
 5823:       for(j=i+1;j<=npar;j++)
 5824: 	matcov[i][j]=matcov[j][i];
 5825:     
 5826:     if(mle==1)
 5827:       printf("\n");
 5828:     fprintf(ficlog,"\n");
 5829:     
 5830:     fflush(ficlog);
 5831:     
 5832:     /*-------- Rewriting parameter file ----------*/
 5833:     strcpy(rfileres,"r");    /* "Rparameterfile */
 5834:     strcat(rfileres,optionfilefiname);    /* Parameter file first name*/
 5835:     strcat(rfileres,".");    /* */
 5836:     strcat(rfileres,optionfilext);    /* Other files have txt extension */
 5837:     if((ficres =fopen(rfileres,"w"))==NULL) {
 5838:       printf("Problem writing new parameter file: %s\n", fileres);goto end;
 5839:       fprintf(ficlog,"Problem writing new parameter file: %s\n", fileres);goto end;
 5840:     }
 5841:     fprintf(ficres,"#%s\n",version);
 5842:   }    /* End of mle != -3 */
 5843: 
 5844: 
 5845:   n= lastobs;
 5846:   num=lvector(1,n);
 5847:   moisnais=vector(1,n);
 5848:   annais=vector(1,n);
 5849:   moisdc=vector(1,n);
 5850:   andc=vector(1,n);
 5851:   agedc=vector(1,n);
 5852:   cod=ivector(1,n);
 5853:   weight=vector(1,n);
 5854:   for(i=1;i<=n;i++) weight[i]=1.0; /* Equal weights, 1 by default */
 5855:   mint=matrix(1,maxwav,1,n);
 5856:   anint=matrix(1,maxwav,1,n);
 5857:   s=imatrix(1,maxwav+1,1,n); /* s[i][j] health state for wave i and individual j */ 
 5858:   tab=ivector(1,NCOVMAX);
 5859:   ncodemax=ivector(1,NCOVMAX); /* Number of code per covariate; if O and 1 only, 2**ncov; V1+V2+V3+V4=>16 */
 5860: 
 5861:   /* Reads data from file datafile */
 5862:   if (readdata(datafile, firstobs, lastobs, &imx)==1)
 5863:     goto end;
 5864: 
 5865:   /* Calculation of the number of parameters from char model */
 5866:     /*    modelsav=V2+V1+V4+age*V3 strb=age*V3 stra=V2+V1+V4 
 5867: 	k=4 (age*V3) Tvar[k=4]= 3 (from V3) Tag[cptcovage=1]=4
 5868: 	k=3 V4 Tvar[k=3]= 4 (from V4)
 5869: 	k=2 V1 Tvar[k=2]= 1 (from V1)
 5870: 	k=1 Tvar[1]=2 (from V2)
 5871:     */
 5872:   Tvar=ivector(1,NCOVMAX); /* Was 15 changed to NCOVMAX. */
 5873:   /*  V2+V1+V4+age*V3 is a model with 4 covariates (3 plus signs). 
 5874:       For each model-covariate stores the data-covariate id. Tvar[1]=2, Tvar[2]=1, Tvar[3]=4, 
 5875:       Tvar[4=age*V3] is 3 and 'age' is recorded in Tage.
 5876:   */
 5877:   /* For model-covariate k tells which data-covariate to use but
 5878:     because this model-covariate is a construction we invent a new column
 5879:     ncovcol + k1
 5880:     If already ncovcol=4 and model=V2+V1+V1*V4+age*V3
 5881:     Tvar[3=V1*V4]=4+1 etc */
 5882:   Tprod=ivector(1,NCOVMAX); /* Gives the position of a product */
 5883:   /* Tprod[k1=1]=3(=V1*V4) for V2+V1+V1*V4+age*V3
 5884:      if  V2+V1+V1*V4+age*V3+V3*V2   TProd[k1=2]=5 (V3*V2)
 5885:   */
 5886:   Tvaraff=ivector(1,NCOVMAX); /* Unclear */
 5887:   Tvard=imatrix(1,NCOVMAX,1,2); /* n=Tvard[k1][1]  and m=Tvard[k1][2] gives the couple n,m of the k1 th product Vn*Vm
 5888: 			    * For V3*V2 (in V2+V1+V1*V4+age*V3+V3*V2), V3*V2 position is 2nd. 
 5889: 			    * Tvard[k1=2][1]=3 (V3) Tvard[k1=2][2]=2(V2) */
 5890:   Tage=ivector(1,NCOVMAX); /* Gives the covariate id of covariates associated with age: V2 + V1 + age*V4 + V3*age
 5891: 			 4 covariates (3 plus signs)
 5892: 			 Tage[1=V3*age]= 4; Tage[2=age*V4] = 3
 5893: 		      */  
 5894: 
 5895:   if(decodemodel(model, lastobs) == 1)
 5896:     goto end;
 5897: 
 5898:   if((double)(lastobs-imx)/(double)imx > 1.10){
 5899:     nbwarn++;
 5900:     printf("Warning: The value of parameter lastobs=%d is big compared to the \n  effective number of cases imx=%d, please adjust, \n  otherwise you are allocating more memory than necessary.\n",lastobs, imx); 
 5901:     fprintf(ficlog,"Warning: The value of parameter lastobs=%d is big compared to the \n  effective number of cases imx=%d, please adjust, \n  otherwise you are allocating more memory than necessary.\n",lastobs, imx); 
 5902:   }
 5903:     /*  if(mle==1){*/
 5904:   if (weightopt != 1) { /* Maximisation without weights. We can have weights different from 1 but want no weight*/
 5905:     for(i=1;i<=imx;i++) weight[i]=1.0; /* changed to imx */
 5906:   }
 5907: 
 5908:     /*-calculation of age at interview from date of interview and age at death -*/
 5909:   agev=matrix(1,maxwav,1,imx);
 5910: 
 5911:   if(calandcheckages(imx, maxwav, &agemin, &agemax, &nberr, &nbwarn) == 1)
 5912:     goto end;
 5913: 
 5914: 
 5915:   agegomp=(int)agemin;
 5916:   free_vector(moisnais,1,n);
 5917:   free_vector(annais,1,n);
 5918:   /* free_matrix(mint,1,maxwav,1,n);
 5919:      free_matrix(anint,1,maxwav,1,n);*/
 5920:   free_vector(moisdc,1,n);
 5921:   free_vector(andc,1,n);
 5922:   /* */
 5923:   
 5924:   wav=ivector(1,imx);
 5925:   dh=imatrix(1,lastpass-firstpass+1,1,imx);
 5926:   bh=imatrix(1,lastpass-firstpass+1,1,imx);
 5927:   mw=imatrix(1,lastpass-firstpass+1,1,imx);
 5928:    
 5929:   /* Concatenates waves */
 5930:   concatwav(wav, dh, bh, mw, s, agedc, agev,  firstpass, lastpass, imx, nlstate, stepm);
 5931:   /* */
 5932:  
 5933:   /* Routine tricode is to calculate cptcoveff (real number of unique covariates) and to associate covariable number and modality */
 5934: 
 5935:   nbcode=imatrix(0,NCOVMAX,0,NCOVMAX); 
 5936:   ncodemax[1]=1;
 5937:   Ndum =ivector(-1,NCOVMAX);  
 5938:   if (ncovmodel > 2)
 5939:     tricode(Tvar,nbcode,imx, Ndum); /**< Fills nbcode[Tvar[j]][l]; */
 5940: 
 5941:   codtab=imatrix(1,100,1,10); /* codtab[h,k]=( (h-1) - mod(k-1,2**(k-1) )/2**(k-1) */
 5942:   /*printf(" codtab[1,1],codtab[100,10]=%d,%d\n", codtab[1][1],codtab[100][10]);*/
 5943:   h=0;
 5944: 
 5945: 
 5946:   /*if (cptcovn > 0) */
 5947:       
 5948:  
 5949:   m=pow(2,cptcoveff);
 5950:  
 5951:   for(k=1;k<=cptcoveff; k++){ /* scans any effective covariate */
 5952:     for(i=1; i <=pow(2,cptcoveff-k);i++){ /* i=1 to 8/1=8; i=1 to 8/2=4; i=1 to 8/8=1 */ 
 5953:       for(j=1; j <= ncodemax[k]; j++){ /* For each modality of this covariate ncodemax=2*/
 5954: 	for(cpt=1; cpt <=pow(2,k-1); cpt++){  /* cpt=1 to 8/2**(3+1-1 or 3+1-3) =1 or 4 */ 
 5955: 	  h++;
 5956: 	  if (h>m) 
 5957: 	    h=1;
 5958: 	  /**< codtab(h,k)  k   = codtab[h,k]=( (h-1) - mod(k-1,2**(k-1) )/2**(k-1) + 1
 5959: 	   *     h     1     2     3     4
 5960: 	   *______________________________  
 5961: 	   *     1 i=1 1 i=1 1 i=1 1 i=1 1
 5962: 	   *     2     2     1     1     1
 5963: 	   *     3 i=2 1     2     1     1
 5964: 	   *     4     2     2     1     1
 5965: 	   *     5 i=3 1 i=2 1     2     1
 5966: 	   *     6     2     1     2     1
 5967: 	   *     7 i=4 1     2     2     1
 5968: 	   *     8     2     2     2     1
 5969: 	   *     9 i=5 1 i=3 1 i=2 1     1
 5970: 	   *    10     2     1     1     1
 5971: 	   *    11 i=6 1     2     1     1
 5972: 	   *    12     2     2     1     1
 5973: 	   *    13 i=7 1 i=4 1     2     1    
 5974: 	   *    14     2     1     2     1
 5975: 	   *    15 i=8 1     2     2     1
 5976: 	   *    16     2     2     2     1
 5977: 	   */
 5978: 	  codtab[h][k]=j;
 5979: 	  /*codtab[h][Tvar[k]]=j;*/
 5980: 	  printf("h=%d k=%d j=%d codtab[h][k]=%d Tvar[k]=%d codtab[h][Tvar[k]]=%d \n",h, k,j,codtab[h][k],Tvar[k],codtab[h][Tvar[k]]);
 5981: 	} 
 5982:       }
 5983:     }
 5984:   } 
 5985:   /* printf("codtab[1][2]=%d codtab[2][2]=%d",codtab[1][2],codtab[2][2]); 
 5986:      codtab[1][2]=1;codtab[2][2]=2; */
 5987:   /* for(i=1; i <=m ;i++){ 
 5988:      for(k=1; k <=cptcovn; k++){
 5989:        printf("i=%d k=%d %d %d ",i,k,codtab[i][k], cptcoveff);
 5990:      }
 5991:      printf("\n");
 5992:      }
 5993:      scanf("%d",i);*/
 5994: 
 5995:  free_ivector(Ndum,-1,NCOVMAX);
 5996: 
 5997: 
 5998:     
 5999:   /*------------ gnuplot -------------*/
 6000:   strcpy(optionfilegnuplot,optionfilefiname);
 6001:   if(mle==-3)
 6002:     strcat(optionfilegnuplot,"-mort");
 6003:   strcat(optionfilegnuplot,".gp");
 6004: 
 6005:   if((ficgp=fopen(optionfilegnuplot,"w"))==NULL) {
 6006:     printf("Problem with file %s",optionfilegnuplot);
 6007:   }
 6008:   else{
 6009:     fprintf(ficgp,"\n# %s\n", version); 
 6010:     fprintf(ficgp,"# %s\n", optionfilegnuplot); 
 6011:     //fprintf(ficgp,"set missing 'NaNq'\n");
 6012:     fprintf(ficgp,"set datafile missing 'NaNq'\n");
 6013:   }
 6014:   /*  fclose(ficgp);*/
 6015:   /*--------- index.htm --------*/
 6016: 
 6017:   strcpy(optionfilehtm,optionfilefiname); /* Main html file */
 6018:   if(mle==-3)
 6019:     strcat(optionfilehtm,"-mort");
 6020:   strcat(optionfilehtm,".htm");
 6021:   if((fichtm=fopen(optionfilehtm,"w"))==NULL)    {
 6022:     printf("Problem with %s \n",optionfilehtm);
 6023:     exit(0);
 6024:   }
 6025: 
 6026:   strcpy(optionfilehtmcov,optionfilefiname); /* Only for matrix of covariance */
 6027:   strcat(optionfilehtmcov,"-cov.htm");
 6028:   if((fichtmcov=fopen(optionfilehtmcov,"w"))==NULL)    {
 6029:     printf("Problem with %s \n",optionfilehtmcov), exit(0);
 6030:   }
 6031:   else{
 6032:   fprintf(fichtmcov,"<html><head>\n<title>IMaCh Cov %s</title></head>\n <body><font size=\"2\">%s <br> %s</font> \
 6033: <hr size=\"2\" color=\"#EC5E5E\"> \n\
 6034: Title=%s <br>Datafile=%s Firstpass=%d Lastpass=%d Stepm=%d Weight=%d Model=%s<br>\n",\
 6035: 	  optionfilehtmcov,version,fullversion,title,datafile,firstpass,lastpass,stepm, weightopt, model);
 6036:   }
 6037: 
 6038:   fprintf(fichtm,"<html><head>\n<title>IMaCh %s</title></head>\n <body><font size=\"2\">%s <br> %s</font> \
 6039: <hr size=\"2\" color=\"#EC5E5E\"> \n\
 6040: Title=%s <br>Datafile=%s Firstpass=%d Lastpass=%d Stepm=%d Weight=%d Model=%s<br>\n\
 6041: \n\
 6042: <hr  size=\"2\" color=\"#EC5E5E\">\
 6043:  <ul><li><h4>Parameter files</h4>\n\
 6044:  - Parameter file: <a href=\"%s.%s\">%s.%s</a><br>\n\
 6045:  - Copy of the parameter file: <a href=\"o%s\">o%s</a><br>\n\
 6046:  - Log file of the run: <a href=\"%s\">%s</a><br>\n\
 6047:  - Gnuplot file name: <a href=\"%s\">%s</a><br>\n\
 6048:  - Date and time at start: %s</ul>\n",\
 6049: 	  optionfilehtm,version,fullversion,title,datafile,firstpass,lastpass,stepm, weightopt, model,\
 6050: 	  optionfilefiname,optionfilext,optionfilefiname,optionfilext,\
 6051: 	  fileres,fileres,\
 6052: 	  filelog,filelog,optionfilegnuplot,optionfilegnuplot,strstart);
 6053:   fflush(fichtm);
 6054: 
 6055:   strcpy(pathr,path);
 6056:   strcat(pathr,optionfilefiname);
 6057:   chdir(optionfilefiname); /* Move to directory named optionfile */
 6058:   
 6059:   /* Calculates basic frequencies. Computes observed prevalence at single age
 6060:      and prints on file fileres'p'. */
 6061:   freqsummary(fileres, agemin, agemax, s, agev, nlstate, imx,Tvaraff,nbcode, ncodemax,mint,anint,strstart);
 6062: 
 6063:   fprintf(fichtm,"\n");
 6064:   fprintf(fichtm,"<br>Total number of observations=%d <br>\n\
 6065: Youngest age at first (selected) pass %.2f, oldest age %.2f<br>\n\
 6066: Interval (in months) between two waves: Min=%d Max=%d Mean=%.2lf<br>\n",\
 6067: 	  imx,agemin,agemax,jmin,jmax,jmean);
 6068:   pmmij= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6069:     oldms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6070:     newms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6071:     savms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6072:     oldm=oldms; newm=newms; savm=savms; /* Keeps fixed addresses to free */
 6073:     
 6074:    
 6075:   /* For Powell, parameters are in a vector p[] starting at p[1]
 6076:      so we point p on param[1][1] so that p[1] maps on param[1][1][1] */
 6077:   p=param[1][1]; /* *(*(*(param +1)+1)+0) */
 6078: 
 6079:   globpr=0; /* To get the number ipmx of contributions and the sum of weights*/
 6080: 
 6081:   if (mle==-3){
 6082:     ximort=matrix(1,NDIM,1,NDIM); 
 6083: /*     ximort=gsl_matrix_alloc(1,NDIM,1,NDIM); */
 6084:     cens=ivector(1,n);
 6085:     ageexmed=vector(1,n);
 6086:     agecens=vector(1,n);
 6087:     dcwave=ivector(1,n);
 6088:  
 6089:     for (i=1; i<=imx; i++){
 6090:       dcwave[i]=-1;
 6091:       for (m=firstpass; m<=lastpass; m++)
 6092: 	if (s[m][i]>nlstate) {
 6093: 	  dcwave[i]=m;
 6094: 	  /*	printf("i=%d j=%d s=%d dcwave=%d\n",i,j, s[j][i],dcwave[i]);*/
 6095: 	  break;
 6096: 	}
 6097:     }
 6098: 
 6099:     for (i=1; i<=imx; i++) {
 6100:       if (wav[i]>0){
 6101: 	ageexmed[i]=agev[mw[1][i]][i];
 6102: 	j=wav[i];
 6103: 	agecens[i]=1.; 
 6104: 
 6105: 	if (ageexmed[i]> 1 && wav[i] > 0){
 6106: 	  agecens[i]=agev[mw[j][i]][i];
 6107: 	  cens[i]= 1;
 6108: 	}else if (ageexmed[i]< 1) 
 6109: 	  cens[i]= -1;
 6110: 	if (agedc[i]< AGESUP && agedc[i]>1 && dcwave[i]>firstpass && dcwave[i]<=lastpass)
 6111: 	  cens[i]=0 ;
 6112:       }
 6113:       else cens[i]=-1;
 6114:     }
 6115:     
 6116:     for (i=1;i<=NDIM;i++) {
 6117:       for (j=1;j<=NDIM;j++)
 6118: 	ximort[i][j]=(i == j ? 1.0 : 0.0);
 6119:     }
 6120:     
 6121:     /*p[1]=0.0268; p[NDIM]=0.083;*/
 6122:     /*printf("%lf %lf", p[1], p[2]);*/
 6123:     
 6124:     
 6125: #ifdef GSL
 6126:     printf("GSL optimization\n");  fprintf(ficlog,"Powell\n");
 6127: #else
 6128:     printf("Powell\n");  fprintf(ficlog,"Powell\n");
 6129: #endif
 6130:     strcpy(filerespow,"pow-mort"); 
 6131:     strcat(filerespow,fileres);
 6132:     if((ficrespow=fopen(filerespow,"w"))==NULL) {
 6133:       printf("Problem with resultfile: %s\n", filerespow);
 6134:       fprintf(ficlog,"Problem with resultfile: %s\n", filerespow);
 6135:     }
 6136: #ifdef GSL
 6137:     fprintf(ficrespow,"# GSL optimization\n# iter -2*LL");
 6138: #else
 6139:     fprintf(ficrespow,"# Powell\n# iter -2*LL");
 6140: #endif
 6141:     /*  for (i=1;i<=nlstate;i++)
 6142: 	for(j=1;j<=nlstate+ndeath;j++)
 6143: 	if(j!=i)fprintf(ficrespow," p%1d%1d",i,j);
 6144:     */
 6145:     fprintf(ficrespow,"\n");
 6146: #ifdef GSL
 6147:     /* gsl starts here */ 
 6148:     T = gsl_multimin_fminimizer_nmsimplex;
 6149:     gsl_multimin_fminimizer *sfm = NULL;
 6150:     gsl_vector *ss, *x;
 6151:     gsl_multimin_function minex_func;
 6152: 
 6153:     /* Initial vertex size vector */
 6154:     ss = gsl_vector_alloc (NDIM);
 6155:     
 6156:     if (ss == NULL){
 6157:       GSL_ERROR_VAL ("failed to allocate space for ss", GSL_ENOMEM, 0);
 6158:     }
 6159:     /* Set all step sizes to 1 */
 6160:     gsl_vector_set_all (ss, 0.001);
 6161: 
 6162:     /* Starting point */
 6163:     
 6164:     x = gsl_vector_alloc (NDIM);
 6165:     
 6166:     if (x == NULL){
 6167:       gsl_vector_free(ss);
 6168:       GSL_ERROR_VAL ("failed to allocate space for x", GSL_ENOMEM, 0);
 6169:     }
 6170:   
 6171:     /* Initialize method and iterate */
 6172:     /*     p[1]=0.0268; p[NDIM]=0.083; */
 6173: /*     gsl_vector_set(x, 0, 0.0268); */
 6174: /*     gsl_vector_set(x, 1, 0.083); */
 6175:     gsl_vector_set(x, 0, p[1]);
 6176:     gsl_vector_set(x, 1, p[2]);
 6177: 
 6178:     minex_func.f = &gompertz_f;
 6179:     minex_func.n = NDIM;
 6180:     minex_func.params = (void *)&p; /* ??? */
 6181:     
 6182:     sfm = gsl_multimin_fminimizer_alloc (T, NDIM);
 6183:     gsl_multimin_fminimizer_set (sfm, &minex_func, x, ss);
 6184:     
 6185:     printf("Iterations beginning .....\n\n");
 6186:     printf("Iter. #    Intercept       Slope     -Log Likelihood     Simplex size\n");
 6187: 
 6188:     iteri=0;
 6189:     while (rval == GSL_CONTINUE){
 6190:       iteri++;
 6191:       status = gsl_multimin_fminimizer_iterate(sfm);
 6192:       
 6193:       if (status) printf("error: %s\n", gsl_strerror (status));
 6194:       fflush(0);
 6195:       
 6196:       if (status) 
 6197:         break;
 6198:       
 6199:       rval = gsl_multimin_test_size (gsl_multimin_fminimizer_size (sfm), 1e-6);
 6200:       ssval = gsl_multimin_fminimizer_size (sfm);
 6201:       
 6202:       if (rval == GSL_SUCCESS)
 6203:         printf ("converged to a local maximum at\n");
 6204:       
 6205:       printf("%5d ", iteri);
 6206:       for (it = 0; it < NDIM; it++){
 6207: 	printf ("%10.5f ", gsl_vector_get (sfm->x, it));
 6208:       }
 6209:       printf("f() = %-10.5f ssize = %.7f\n", sfm->fval, ssval);
 6210:     }
 6211:     
 6212:     printf("\n\n Please note: Program should be run many times with varying starting points to detemine global maximum\n\n");
 6213:     
 6214:     gsl_vector_free(x); /* initial values */
 6215:     gsl_vector_free(ss); /* inital step size */
 6216:     for (it=0; it<NDIM; it++){
 6217:       p[it+1]=gsl_vector_get(sfm->x,it);
 6218:       fprintf(ficrespow," %.12lf", p[it]);
 6219:     }
 6220:     gsl_multimin_fminimizer_free (sfm); /* p *(sfm.x.data) et p *(sfm.x.data+1)  */
 6221: #endif
 6222: #ifdef POWELL
 6223:      powell(p,ximort,NDIM,ftol,&iter,&fret,gompertz);
 6224: #endif  
 6225:     fclose(ficrespow);
 6226:     
 6227:     hesscov(matcov, p, NDIM, delti, 1e-4, gompertz); 
 6228: 
 6229:     for(i=1; i <=NDIM; i++)
 6230:       for(j=i+1;j<=NDIM;j++)
 6231: 	matcov[i][j]=matcov[j][i];
 6232:     
 6233:     printf("\nCovariance matrix\n ");
 6234:     for(i=1; i <=NDIM; i++) {
 6235:       for(j=1;j<=NDIM;j++){ 
 6236: 	printf("%f ",matcov[i][j]);
 6237:       }
 6238:       printf("\n ");
 6239:     }
 6240:     
 6241:     printf("iter=%d MLE=%f Eq=%lf*exp(%lf*(age-%d))\n",iter,-gompertz(p),p[1],p[2],agegomp);
 6242:     for (i=1;i<=NDIM;i++) 
 6243:       printf("%f [%f ; %f]\n",p[i],p[i]-2*sqrt(matcov[i][i]),p[i]+2*sqrt(matcov[i][i]));
 6244: 
 6245:     lsurv=vector(1,AGESUP);
 6246:     lpop=vector(1,AGESUP);
 6247:     tpop=vector(1,AGESUP);
 6248:     lsurv[agegomp]=100000;
 6249:     
 6250:     for (k=agegomp;k<=AGESUP;k++) {
 6251:       agemortsup=k;
 6252:       if (p[1]*exp(p[2]*(k-agegomp))>1) break;
 6253:     }
 6254:     
 6255:     for (k=agegomp;k<agemortsup;k++)
 6256:       lsurv[k+1]=lsurv[k]-lsurv[k]*(p[1]*exp(p[2]*(k-agegomp)));
 6257:     
 6258:     for (k=agegomp;k<agemortsup;k++){
 6259:       lpop[k]=(lsurv[k]+lsurv[k+1])/2.;
 6260:       sumlpop=sumlpop+lpop[k];
 6261:     }
 6262:     
 6263:     tpop[agegomp]=sumlpop;
 6264:     for (k=agegomp;k<(agemortsup-3);k++){
 6265:       /*  tpop[k+1]=2;*/
 6266:       tpop[k+1]=tpop[k]-lpop[k];
 6267:     }
 6268:     
 6269:     
 6270:     printf("\nAge   lx     qx    dx    Lx     Tx     e(x)\n");
 6271:     for (k=agegomp;k<(agemortsup-2);k++) 
 6272:       printf("%d %.0lf %lf %.0lf %.0lf %.0lf %lf\n",k,lsurv[k],p[1]*exp(p[2]*(k-agegomp)),(p[1]*exp(p[2]*(k-agegomp)))*lsurv[k],lpop[k],tpop[k],tpop[k]/lsurv[k]);
 6273:     
 6274:     
 6275:     replace_back_to_slash(pathc,pathcd); /* Even gnuplot wants a / */
 6276:     printinggnuplotmort(fileres, optionfilefiname,ageminpar,agemaxpar,fage, pathc,p);
 6277:     
 6278:     printinghtmlmort(fileres,title,datafile, firstpass, lastpass, \
 6279: 		     stepm, weightopt,\
 6280: 		     model,imx,p,matcov,agemortsup);
 6281:     
 6282:     free_vector(lsurv,1,AGESUP);
 6283:     free_vector(lpop,1,AGESUP);
 6284:     free_vector(tpop,1,AGESUP);
 6285: #ifdef GSL
 6286:     free_ivector(cens,1,n);
 6287:     free_vector(agecens,1,n);
 6288:     free_ivector(dcwave,1,n);
 6289:     free_matrix(ximort,1,NDIM,1,NDIM);
 6290: #endif
 6291:   } /* Endof if mle==-3 */
 6292:   
 6293:   else{ /* For mle >=1 */
 6294:     globpr=0;/* debug */
 6295:     likelione(ficres, p, npar, nlstate, &globpr, &ipmx, &sw, &fretone, funcone); /* Prints the contributions to the likelihood */
 6296:     printf("First Likeli=%12.6f ipmx=%ld sw=%12.6f",fretone,ipmx,sw);
 6297:     for (k=1; k<=npar;k++)
 6298:       printf(" %d %8.5f",k,p[k]);
 6299:     printf("\n");
 6300:     globpr=1; /* to print the contributions */
 6301:     likelione(ficres, p, npar, nlstate, &globpr, &ipmx, &sw, &fretone, funcone); /* Prints the contributions to the likelihood */
 6302:     printf("Second Likeli=%12.6f ipmx=%ld sw=%12.6f",fretone,ipmx,sw);
 6303:     for (k=1; k<=npar;k++)
 6304:       printf(" %d %8.5f",k,p[k]);
 6305:     printf("\n");
 6306:     if(mle>=1){ /* Could be 1 or 2 */
 6307:       mlikeli(ficres,p, npar, ncovmodel, nlstate, ftol, func);
 6308:     }
 6309:     
 6310:     /*--------- results files --------------*/
 6311:     fprintf(ficres,"title=%s datafile=%s lastobs=%d firstpass=%d lastpass=%d\nftol=%e stepm=%d ncovcol=%d nlstate=%d ndeath=%d maxwav=%d mle= 0 weight=%d\nmodel=%s\n", title, datafile, lastobs, firstpass,lastpass,ftol, stepm, ncovcol, nlstate, ndeath, maxwav, weightopt,model);
 6312:     
 6313:     
 6314:     fprintf(ficres,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 6315:     printf("# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 6316:     fprintf(ficlog,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 6317:     for(i=1,jk=1; i <=nlstate; i++){
 6318:       for(k=1; k <=(nlstate+ndeath); k++){
 6319: 	if (k != i) {
 6320: 	  printf("%d%d ",i,k);
 6321: 	  fprintf(ficlog,"%d%d ",i,k);
 6322: 	  fprintf(ficres,"%1d%1d ",i,k);
 6323: 	  for(j=1; j <=ncovmodel; j++){
 6324: 	    printf("%lf ",p[jk]);
 6325: 	    fprintf(ficlog,"%lf ",p[jk]);
 6326: 	    fprintf(ficres,"%lf ",p[jk]);
 6327: 	    jk++; 
 6328: 	  }
 6329: 	  printf("\n");
 6330: 	  fprintf(ficlog,"\n");
 6331: 	  fprintf(ficres,"\n");
 6332: 	}
 6333:       }
 6334:     }
 6335:     if(mle!=0){
 6336:       /* Computing hessian and covariance matrix */
 6337:       ftolhess=ftol; /* Usually correct */
 6338:       hesscov(matcov, p, npar, delti, ftolhess, func);
 6339:     }
 6340:     fprintf(ficres,"# Scales (for hessian or gradient estimation)\n");
 6341:     printf("# Scales (for hessian or gradient estimation)\n");
 6342:     fprintf(ficlog,"# Scales (for hessian or gradient estimation)\n");
 6343:     for(i=1,jk=1; i <=nlstate; i++){
 6344:       for(j=1; j <=nlstate+ndeath; j++){
 6345: 	if (j!=i) {
 6346: 	  fprintf(ficres,"%1d%1d",i,j);
 6347: 	  printf("%1d%1d",i,j);
 6348: 	  fprintf(ficlog,"%1d%1d",i,j);
 6349: 	  for(k=1; k<=ncovmodel;k++){
 6350: 	    printf(" %.5e",delti[jk]);
 6351: 	    fprintf(ficlog," %.5e",delti[jk]);
 6352: 	    fprintf(ficres," %.5e",delti[jk]);
 6353: 	    jk++;
 6354: 	  }
 6355: 	  printf("\n");
 6356: 	  fprintf(ficlog,"\n");
 6357: 	  fprintf(ficres,"\n");
 6358: 	}
 6359:       }
 6360:     }
 6361:     
 6362:     fprintf(ficres,"# Covariance matrix \n# 121 Var(a12)\n# 122 Cov(b12,a12) Var(b12)\n#   ...\n# 232 Cov(b23,a12)  Cov(b23,b12) ... Var (b23)\n");
 6363:     if(mle>=1)
 6364:       printf("# Covariance matrix \n# 121 Var(a12)\n# 122 Cov(b12,a12) Var(b12)\n#   ...\n# 232 Cov(b23,a12)  Cov(b23,b12) ... Var (b23)\n");
 6365:     fprintf(ficlog,"# Covariance matrix \n# 121 Var(a12)\n# 122 Cov(b12,a12) Var(b12)\n#   ...\n# 232 Cov(b23,a12)  Cov(b23,b12) ... Var (b23)\n");
 6366:     /* # 121 Var(a12)\n\ */
 6367:     /* # 122 Cov(b12,a12) Var(b12)\n\ */
 6368:     /* # 131 Cov(a13,a12) Cov(a13,b12, Var(a13)\n\ */
 6369:     /* # 132 Cov(b13,a12) Cov(b13,b12, Cov(b13,a13) Var(b13)\n\ */
 6370:     /* # 212 Cov(a21,a12) Cov(a21,b12, Cov(a21,a13) Cov(a21,b13) Var(a21)\n\ */
 6371:     /* # 212 Cov(b21,a12) Cov(b21,b12, Cov(b21,a13) Cov(b21,b13) Cov(b21,a21) Var(b21)\n\ */
 6372:     /* # 232 Cov(a23,a12) Cov(a23,b12, Cov(a23,a13) Cov(a23,b13) Cov(a23,a21) Cov(a23,b21) Var(a23)\n\ */
 6373:     /* # 232 Cov(b23,a12) Cov(b23,b12) ... Var (b23)\n" */
 6374:     
 6375:     
 6376:     /* Just to have a covariance matrix which will be more understandable
 6377:        even is we still don't want to manage dictionary of variables
 6378:     */
 6379:     for(itimes=1;itimes<=2;itimes++){
 6380:       jj=0;
 6381:       for(i=1; i <=nlstate; i++){
 6382: 	for(j=1; j <=nlstate+ndeath; j++){
 6383: 	  if(j==i) continue;
 6384: 	  for(k=1; k<=ncovmodel;k++){
 6385: 	    jj++;
 6386: 	    ca[0]= k+'a'-1;ca[1]='\0';
 6387: 	    if(itimes==1){
 6388: 	      if(mle>=1)
 6389: 		printf("#%1d%1d%d",i,j,k);
 6390: 	      fprintf(ficlog,"#%1d%1d%d",i,j,k);
 6391: 	      fprintf(ficres,"#%1d%1d%d",i,j,k);
 6392: 	    }else{
 6393: 	      if(mle>=1)
 6394: 		printf("%1d%1d%d",i,j,k);
 6395: 	      fprintf(ficlog,"%1d%1d%d",i,j,k);
 6396: 	      fprintf(ficres,"%1d%1d%d",i,j,k);
 6397: 	    }
 6398: 	    ll=0;
 6399: 	    for(li=1;li <=nlstate; li++){
 6400: 	      for(lj=1;lj <=nlstate+ndeath; lj++){
 6401: 		if(lj==li) continue;
 6402: 		for(lk=1;lk<=ncovmodel;lk++){
 6403: 		  ll++;
 6404: 		  if(ll<=jj){
 6405: 		    cb[0]= lk +'a'-1;cb[1]='\0';
 6406: 		    if(ll<jj){
 6407: 		      if(itimes==1){
 6408: 			if(mle>=1)
 6409: 			  printf(" Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 6410: 			fprintf(ficlog," Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 6411: 			fprintf(ficres," Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 6412: 		      }else{
 6413: 			if(mle>=1)
 6414: 			  printf(" %.5e",matcov[jj][ll]); 
 6415: 			fprintf(ficlog," %.5e",matcov[jj][ll]); 
 6416: 			fprintf(ficres," %.5e",matcov[jj][ll]); 
 6417: 		      }
 6418: 		    }else{
 6419: 		      if(itimes==1){
 6420: 			if(mle>=1)
 6421: 			  printf(" Var(%s%1d%1d)",ca,i,j);
 6422: 			fprintf(ficlog," Var(%s%1d%1d)",ca,i,j);
 6423: 			fprintf(ficres," Var(%s%1d%1d)",ca,i,j);
 6424: 		      }else{
 6425: 			if(mle>=1)
 6426: 			  printf(" %.5e",matcov[jj][ll]); 
 6427: 			fprintf(ficlog," %.5e",matcov[jj][ll]); 
 6428: 			fprintf(ficres," %.5e",matcov[jj][ll]); 
 6429: 		      }
 6430: 		    }
 6431: 		  }
 6432: 		} /* end lk */
 6433: 	      } /* end lj */
 6434: 	    } /* end li */
 6435: 	    if(mle>=1)
 6436: 	      printf("\n");
 6437: 	    fprintf(ficlog,"\n");
 6438: 	    fprintf(ficres,"\n");
 6439: 	    numlinepar++;
 6440: 	  } /* end k*/
 6441: 	} /*end j */
 6442:       } /* end i */
 6443:     } /* end itimes */
 6444:     
 6445:     fflush(ficlog);
 6446:     fflush(ficres);
 6447:     
 6448:     while((c=getc(ficpar))=='#' && c!= EOF){
 6449:       ungetc(c,ficpar);
 6450:       fgets(line, MAXLINE, ficpar);
 6451:       fputs(line,stdout);
 6452:       fputs(line,ficparo);
 6453:     }
 6454:     ungetc(c,ficpar);
 6455:     
 6456:     estepm=0;
 6457:     fscanf(ficpar,"agemin=%lf agemax=%lf bage=%lf fage=%lf estepm=%d\n",&ageminpar,&agemaxpar, &bage, &fage, &estepm);
 6458:     if (estepm==0 || estepm < stepm) estepm=stepm;
 6459:     if (fage <= 2) {
 6460:       bage = ageminpar;
 6461:       fage = agemaxpar;
 6462:     }
 6463:     
 6464:     fprintf(ficres,"# agemin agemax for life expectancy, bage fage (if mle==0 ie no data nor Max likelihood).\n");
 6465:     fprintf(ficres,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f estepm=%d\n",ageminpar,agemaxpar,bage,fage, estepm);
 6466:     fprintf(ficparo,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f estepm=%d\n",ageminpar,agemaxpar,bage,fage, estepm);
 6467:     
 6468:     while((c=getc(ficpar))=='#' && c!= EOF){
 6469:       ungetc(c,ficpar);
 6470:       fgets(line, MAXLINE, ficpar);
 6471:       fputs(line,stdout);
 6472:       fputs(line,ficparo);
 6473:     }
 6474:     ungetc(c,ficpar);
 6475:     
 6476:     fscanf(ficpar,"begin-prev-date=%lf/%lf/%lf end-prev-date=%lf/%lf/%lf mov_average=%d\n",&jprev1, &mprev1,&anprev1,&jprev2, &mprev2,&anprev2,&mobilav);
 6477:     fprintf(ficparo,"begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf mov_average=%d\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,mobilav);
 6478:     fprintf(ficres,"begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf mov_average=%d\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,mobilav);
 6479:     printf("begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf mov_average=%d\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,mobilav);
 6480:     fprintf(ficlog,"begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf mov_average=%d\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,mobilav);
 6481:     
 6482:     while((c=getc(ficpar))=='#' && c!= EOF){
 6483:       ungetc(c,ficpar);
 6484:       fgets(line, MAXLINE, ficpar);
 6485:       fputs(line,stdout);
 6486:       fputs(line,ficparo);
 6487:     }
 6488:     ungetc(c,ficpar);
 6489:     
 6490:     
 6491:     dateprev1=anprev1+(mprev1-1)/12.+(jprev1-1)/365.;
 6492:     dateprev2=anprev2+(mprev2-1)/12.+(jprev2-1)/365.;
 6493:     
 6494:     fscanf(ficpar,"pop_based=%d\n",&popbased);
 6495:     fprintf(ficparo,"pop_based=%d\n",popbased);   
 6496:     fprintf(ficres,"pop_based=%d\n",popbased);   
 6497:     
 6498:     while((c=getc(ficpar))=='#' && c!= EOF){
 6499:       ungetc(c,ficpar);
 6500:       fgets(line, MAXLINE, ficpar);
 6501:       fputs(line,stdout);
 6502:       fputs(line,ficparo);
 6503:     }
 6504:     ungetc(c,ficpar);
 6505:     
 6506:     fscanf(ficpar,"prevforecast=%d starting-proj-date=%lf/%lf/%lf final-proj-date=%lf/%lf/%lf mobil_average=%d\n",&prevfcast,&jproj1,&mproj1,&anproj1,&jproj2,&mproj2,&anproj2,&mobilavproj);
 6507:     fprintf(ficparo,"prevforecast=%d starting-proj-date=%.lf/%.lf/%.lf final-proj-date=%.lf/%.lf/%.lf mobil_average=%d\n",prevfcast,jproj1,mproj1,anproj1,jproj2,mproj2,anproj2,mobilavproj);
 6508:     printf("prevforecast=%d starting-proj-date=%.lf/%.lf/%.lf final-proj-date=%.lf/%.lf/%.lf mobil_average=%d\n",prevfcast,jproj1,mproj1,anproj1,jproj2,mproj2,anproj2,mobilavproj);
 6509:     fprintf(ficlog,"prevforecast=%d starting-proj-date=%.lf/%.lf/%.lf final-proj-date=%.lf/%.lf/%.lf mobil_average=%d\n",prevfcast,jproj1,mproj1,anproj1,jproj2,mproj2,anproj2,mobilavproj);
 6510:     fprintf(ficres,"prevforecast=%d starting-proj-date=%.lf/%.lf/%.lf final-proj-date=%.lf/%.lf/%.lf mobil_average=%d\n",prevfcast,jproj1,mproj1,anproj1,jproj2,mproj2,anproj2,mobilavproj);
 6511:     /* day and month of proj2 are not used but only year anproj2.*/
 6512:     
 6513:     
 6514:     
 6515:      /* freqsummary(fileres, agemin, agemax, s, agev, nlstate, imx,Tvaraff,nbcode, ncodemax,mint,anint); */
 6516:     /* ,dateprev1,dateprev2,jprev1, mprev1,anprev1,jprev2, mprev2,anprev2); */
 6517:     
 6518:     replace_back_to_slash(pathc,pathcd); /* Even gnuplot wants a / */
 6519:     printinggnuplot(fileres, optionfilefiname,ageminpar,agemaxpar,fage, pathc,p);
 6520:     
 6521:     printinghtml(fileres,title,datafile, firstpass, lastpass, stepm, weightopt,\
 6522: 		 model,imx,jmin,jmax,jmean,rfileres,popforecast,estepm,\
 6523: 		 jprev1,mprev1,anprev1,jprev2,mprev2,anprev2);
 6524:       
 6525:    /*------------ free_vector  -------------*/
 6526:    /*  chdir(path); */
 6527:  
 6528:     free_ivector(wav,1,imx);
 6529:     free_imatrix(dh,1,lastpass-firstpass+1,1,imx);
 6530:     free_imatrix(bh,1,lastpass-firstpass+1,1,imx);
 6531:     free_imatrix(mw,1,lastpass-firstpass+1,1,imx);   
 6532:     free_lvector(num,1,n);
 6533:     free_vector(agedc,1,n);
 6534:     /*free_matrix(covar,0,NCOVMAX,1,n);*/
 6535:     /*free_matrix(covar,1,NCOVMAX,1,n);*/
 6536:     fclose(ficparo);
 6537:     fclose(ficres);
 6538: 
 6539: 
 6540:     /*--------------- Prevalence limit  (period or stable prevalence) --------------*/
 6541: #include "prevlim.h"  /* Use ficrespl, ficlog */
 6542:     fclose(ficrespl);
 6543: 
 6544: #ifdef FREEEXIT2
 6545: #include "freeexit2.h"
 6546: #endif
 6547: 
 6548:     /*------------- h Pij x at various ages ------------*/
 6549: #include "hpijx.h"
 6550:     fclose(ficrespij);
 6551: 
 6552:   /*-------------- Variance of one-step probabilities---*/
 6553:     k=1;
 6554:     varprob(optionfilefiname, matcov, p, delti, nlstate, bage, fage,k,Tvar,nbcode, ncodemax,strstart);
 6555: 
 6556: 
 6557:     probs= ma3x(1,AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6558:     for(i=1;i<=AGESUP;i++)
 6559:       for(j=1;j<=NCOVMAX;j++)
 6560: 	for(k=1;k<=NCOVMAX;k++)
 6561: 	  probs[i][j][k]=0.;
 6562: 
 6563:     /*---------- Forecasting ------------------*/
 6564:     /*if((stepm == 1) && (strcmp(model,".")==0)){*/
 6565:     if(prevfcast==1){
 6566:       /*    if(stepm ==1){*/
 6567:       prevforecast(fileres, anproj1, mproj1, jproj1, agemin, agemax, dateprev1, dateprev2, mobilavproj, bage, fage, firstpass, lastpass, anproj2, p, cptcoveff);
 6568:       /* (popforecast==1) populforecast(fileres, anpyram,mpyram,jpyram, agemin,agemax, dateprev1, dateprev2,mobilav, agedeb, fage, popforecast, popfile, anpyram1,p, i1);*/
 6569:       /*      }  */
 6570:       /*      else{ */
 6571:       /*        erreur=108; */
 6572:       /*        printf("Warning %d!! You can only forecast the prevalences if the optimization\n  has been performed with stepm = 1 (month) instead of %d or model=. instead of '%s'\n", erreur, stepm, model); */
 6573:       /*        fprintf(ficlog,"Warning %d!! You can only forecast the prevalences if the optimization\n  has been performed with stepm = 1 (month) instead of %d or model=. instead of '%s'\n", erreur, stepm, model); */
 6574:       /*      } */
 6575:     }
 6576:   
 6577: 
 6578:     /* Computes prevalence between agemin (i.e minimal age computed) and no more ageminpar */
 6579: 
 6580:     prevalence(probs, agemin, agemax, s, agev, nlstate, imx, Tvar, nbcode, ncodemax, mint, anint, dateprev1, dateprev2, firstpass, lastpass);
 6581:     /*  printf("ageminpar=%f, agemax=%f, s[lastpass][imx]=%d, agev[lastpass][imx]=%f, nlstate=%d, imx=%d,  mint[lastpass][imx]=%f, anint[lastpass][imx]=%f,dateprev1=%f, dateprev2=%f, firstpass=%d, lastpass=%d\n",\
 6582: 	ageminpar, agemax, s[lastpass][imx], agev[lastpass][imx], nlstate, imx, mint[lastpass][imx],anint[lastpass][imx], dateprev1, dateprev2, firstpass, lastpass);
 6583:     */
 6584: 
 6585:     if (mobilav!=0) {
 6586:       mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6587:       if (movingaverage(probs, bage, fage, mobaverage,mobilav)!=0){
 6588: 	fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 6589: 	printf(" Error in movingaverage mobilav=%d\n",mobilav);
 6590:       }
 6591:     }
 6592: 
 6593: 
 6594:     /*---------- Health expectancies, no variances ------------*/
 6595: 
 6596:     strcpy(filerese,"e");
 6597:     strcat(filerese,fileres);
 6598:     if((ficreseij=fopen(filerese,"w"))==NULL) {
 6599:       printf("Problem with Health Exp. resultfile: %s\n", filerese); exit(0);
 6600:       fprintf(ficlog,"Problem with Health Exp. resultfile: %s\n", filerese); exit(0);
 6601:     }
 6602:     printf("Computing Health Expectancies: result on file '%s' \n", filerese);
 6603:     fprintf(ficlog,"Computing Health Expectancies: result on file '%s' \n", filerese);
 6604:     /*for(cptcov=1,k=0;cptcov<=i1;cptcov++){
 6605:       for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){*/
 6606:           
 6607:     for (k=1; k <= (int) pow(2,cptcoveff); k++){
 6608: 	fprintf(ficreseij,"\n#****** ");
 6609: 	for(j=1;j<=cptcoveff;j++) {
 6610: 	  fprintf(ficreseij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6611: 	}
 6612: 	fprintf(ficreseij,"******\n");
 6613: 
 6614: 	eij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
 6615: 	oldm=oldms;savm=savms;
 6616: 	evsij(eij, p, nlstate, stepm, (int) bage, (int)fage, oldm, savm, k, estepm, strstart);  
 6617:       
 6618: 	free_ma3x(eij,1,nlstate,1,nlstate,(int) bage, (int)fage);
 6619:       /*}*/
 6620:     }
 6621:     fclose(ficreseij);
 6622: 
 6623: 
 6624:     /*---------- Health expectancies and variances ------------*/
 6625: 
 6626: 
 6627:     strcpy(filerest,"t");
 6628:     strcat(filerest,fileres);
 6629:     if((ficrest=fopen(filerest,"w"))==NULL) {
 6630:       printf("Problem with total LE resultfile: %s\n", filerest);goto end;
 6631:       fprintf(ficlog,"Problem with total LE resultfile: %s\n", filerest);goto end;
 6632:     }
 6633:     printf("Computing Total Life expectancies with their standard errors: file '%s' \n", filerest); 
 6634:     fprintf(ficlog,"Computing Total Life expectancies with their standard errors: file '%s' \n", filerest); 
 6635: 
 6636: 
 6637:     strcpy(fileresstde,"stde");
 6638:     strcat(fileresstde,fileres);
 6639:     if((ficresstdeij=fopen(fileresstde,"w"))==NULL) {
 6640:       printf("Problem with Health Exp. and std errors resultfile: %s\n", fileresstde); exit(0);
 6641:       fprintf(ficlog,"Problem with Health Exp. and std errors resultfile: %s\n", fileresstde); exit(0);
 6642:     }
 6643:     printf("Computing Health Expectancies and standard errors: result on file '%s' \n", fileresstde);
 6644:     fprintf(ficlog,"Computing Health Expectancies and standard errors: result on file '%s' \n", fileresstde);
 6645: 
 6646:     strcpy(filerescve,"cve");
 6647:     strcat(filerescve,fileres);
 6648:     if((ficrescveij=fopen(filerescve,"w"))==NULL) {
 6649:       printf("Problem with Covar. Health Exp. resultfile: %s\n", filerescve); exit(0);
 6650:       fprintf(ficlog,"Problem with Covar. Health Exp. resultfile: %s\n", filerescve); exit(0);
 6651:     }
 6652:     printf("Computing Covar. of Health Expectancies: result on file '%s' \n", filerescve);
 6653:     fprintf(ficlog,"Computing Covar. of Health Expectancies: result on file '%s' \n", filerescve);
 6654: 
 6655:     strcpy(fileresv,"v");
 6656:     strcat(fileresv,fileres);
 6657:     if((ficresvij=fopen(fileresv,"w"))==NULL) {
 6658:       printf("Problem with variance resultfile: %s\n", fileresv);exit(0);
 6659:       fprintf(ficlog,"Problem with variance resultfile: %s\n", fileresv);exit(0);
 6660:     }
 6661:     printf("Computing Variance-covariance of DFLEs: file '%s' \n", fileresv);
 6662:     fprintf(ficlog,"Computing Variance-covariance of DFLEs: file '%s' \n", fileresv);
 6663: 
 6664:     /*for(cptcov=1,k=0;cptcov<=i1;cptcov++){
 6665:       for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){*/
 6666:           
 6667:     for (k=1; k <= (int) pow(2,cptcoveff); k++){
 6668:     	fprintf(ficrest,"\n#****** ");
 6669: 	for(j=1;j<=cptcoveff;j++) 
 6670: 	  fprintf(ficrest,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6671: 	fprintf(ficrest,"******\n");
 6672: 
 6673: 	fprintf(ficresstdeij,"\n#****** ");
 6674: 	fprintf(ficrescveij,"\n#****** ");
 6675: 	for(j=1;j<=cptcoveff;j++) {
 6676: 	  fprintf(ficresstdeij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6677: 	  fprintf(ficrescveij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6678: 	}
 6679: 	fprintf(ficresstdeij,"******\n");
 6680: 	fprintf(ficrescveij,"******\n");
 6681: 
 6682: 	fprintf(ficresvij,"\n#****** ");
 6683: 	for(j=1;j<=cptcoveff;j++) 
 6684: 	  fprintf(ficresvij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6685: 	fprintf(ficresvij,"******\n");
 6686: 
 6687: 	eij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
 6688: 	oldm=oldms;savm=savms;
 6689: 	cvevsij(eij, p, nlstate, stepm, (int) bage, (int)fage, oldm, savm, k, estepm, delti, matcov, strstart);  
 6690: 	/*
 6691: 	 */
 6692: 	/* goto endfree; */
 6693:  
 6694: 	vareij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
 6695: 	pstamp(ficrest);
 6696: 
 6697: 
 6698: 	for(vpopbased=0; vpopbased <= popbased; vpopbased++){ /* Done for vpopbased=0 and vpopbased=1 if popbased==1*/
 6699: 	  oldm=oldms;savm=savms; /* Segmentation fault */
 6700: 	  cptcod= 0; /* To be deleted */
 6701: 	  varevsij(optionfilefiname, vareij, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k, estepm, cptcov,cptcod,vpopbased,mobilav, strstart); /* cptcod not initialized Intel */
 6702: 	  fprintf(ficrest,"# Total life expectancy with std error and decomposition into time to be expected in each health state\n#  (weighted average of eij where weights are ");
 6703: 	  if(vpopbased==1)
 6704: 	    fprintf(ficrest,"the age specific prevalence observed (cross-sectionally) in the population i.e cross-sectionally\n in each health state (popbased=1) (mobilav=%d)\n",mobilav);
 6705: 	  else
 6706: 	    fprintf(ficrest,"the age specific period (stable) prevalences in each health state \n");
 6707: 	  fprintf(ficrest,"# Age e.. (std) ");
 6708: 	  for (i=1;i<=nlstate;i++) fprintf(ficrest,"e.%d (std) ",i);
 6709: 	  fprintf(ficrest,"\n");
 6710: 
 6711: 	  epj=vector(1,nlstate+1);
 6712: 	  for(age=bage; age <=fage ;age++){
 6713: 	    prevalim(prlim, nlstate, p, age, oldm, savm,ftolpl,k);
 6714: 	    if (vpopbased==1) {
 6715: 	      if(mobilav ==0){
 6716: 		for(i=1; i<=nlstate;i++)
 6717: 		  prlim[i][i]=probs[(int)age][i][k];
 6718: 	      }else{ /* mobilav */ 
 6719: 		for(i=1; i<=nlstate;i++)
 6720: 		  prlim[i][i]=mobaverage[(int)age][i][k];
 6721: 	      }
 6722: 	    }
 6723: 	
 6724: 	    fprintf(ficrest," %4.0f",age);
 6725: 	    for(j=1, epj[nlstate+1]=0.;j <=nlstate;j++){
 6726: 	      for(i=1, epj[j]=0.;i <=nlstate;i++) {
 6727: 		epj[j] += prlim[i][i]*eij[i][j][(int)age];
 6728: 		/*  printf("%lf %lf ", prlim[i][i] ,eij[i][j][(int)age]);*/
 6729: 	      }
 6730: 	      epj[nlstate+1] +=epj[j];
 6731: 	    }
 6732: 
 6733: 	    for(i=1, vepp=0.;i <=nlstate;i++)
 6734: 	      for(j=1;j <=nlstate;j++)
 6735: 		vepp += vareij[i][j][(int)age];
 6736: 	    fprintf(ficrest," %7.3f (%7.3f)", epj[nlstate+1],sqrt(vepp));
 6737: 	    for(j=1;j <=nlstate;j++){
 6738: 	      fprintf(ficrest," %7.3f (%7.3f)", epj[j],sqrt(vareij[j][j][(int)age]));
 6739: 	    }
 6740: 	    fprintf(ficrest,"\n");
 6741: 	  }
 6742: 	}
 6743: 	free_ma3x(eij,1,nlstate,1,nlstate,(int) bage, (int)fage);
 6744: 	free_ma3x(vareij,1,nlstate,1,nlstate,(int) bage, (int)fage);
 6745: 	free_vector(epj,1,nlstate+1);
 6746:       /*}*/
 6747:     }
 6748:     free_vector(weight,1,n);
 6749:     free_imatrix(Tvard,1,NCOVMAX,1,2);
 6750:     free_imatrix(s,1,maxwav+1,1,n);
 6751:     free_matrix(anint,1,maxwav,1,n); 
 6752:     free_matrix(mint,1,maxwav,1,n);
 6753:     free_ivector(cod,1,n);
 6754:     free_ivector(tab,1,NCOVMAX);
 6755:     fclose(ficresstdeij);
 6756:     fclose(ficrescveij);
 6757:     fclose(ficresvij);
 6758:     fclose(ficrest);
 6759:     fclose(ficpar);
 6760:   
 6761:     /*------- Variance of period (stable) prevalence------*/   
 6762: 
 6763:     strcpy(fileresvpl,"vpl");
 6764:     strcat(fileresvpl,fileres);
 6765:     if((ficresvpl=fopen(fileresvpl,"w"))==NULL) {
 6766:       printf("Problem with variance of period (stable) prevalence  resultfile: %s\n", fileresvpl);
 6767:       exit(0);
 6768:     }
 6769:     printf("Computing Variance-covariance of period (stable) prevalence: file '%s' \n", fileresvpl);
 6770: 
 6771:     /*for(cptcov=1,k=0;cptcov<=i1;cptcov++){
 6772:       for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){*/
 6773:           
 6774:     for (k=1; k <= (int) pow(2,cptcoveff); k++){
 6775:     	fprintf(ficresvpl,"\n#****** ");
 6776: 	for(j=1;j<=cptcoveff;j++) 
 6777: 	  fprintf(ficresvpl,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6778: 	fprintf(ficresvpl,"******\n");
 6779:       
 6780: 	varpl=matrix(1,nlstate,(int) bage, (int) fage);
 6781: 	oldm=oldms;savm=savms;
 6782: 	varprevlim(fileres, varpl, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k,strstart);
 6783: 	free_matrix(varpl,1,nlstate,(int) bage, (int)fage);
 6784:       /*}*/
 6785:     }
 6786: 
 6787:     fclose(ficresvpl);
 6788: 
 6789:     /*---------- End : free ----------------*/
 6790:     if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6791:     free_ma3x(probs,1,AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6792:   }  /* mle==-3 arrives here for freeing */
 6793:  /* endfree:*/
 6794:     free_matrix(prlim,1,nlstate,1,nlstate); /*here or after loop ? */
 6795:     free_matrix(pmmij,1,nlstate+ndeath,1,nlstate+ndeath);
 6796:     free_matrix(oldms, 1,nlstate+ndeath,1,nlstate+ndeath);
 6797:     free_matrix(newms, 1,nlstate+ndeath,1,nlstate+ndeath);
 6798:     free_matrix(savms, 1,nlstate+ndeath,1,nlstate+ndeath);
 6799:     free_matrix(covar,0,NCOVMAX,1,n);
 6800:     free_matrix(matcov,1,npar,1,npar);
 6801:     /*free_vector(delti,1,npar);*/
 6802:     free_ma3x(delti3,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel); 
 6803:     free_matrix(agev,1,maxwav,1,imx);
 6804:     free_ma3x(param,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel);
 6805: 
 6806:     free_ivector(ncodemax,1,NCOVMAX);
 6807:     free_ivector(Tvar,1,NCOVMAX);
 6808:     free_ivector(Tprod,1,NCOVMAX);
 6809:     free_ivector(Tvaraff,1,NCOVMAX);
 6810:     free_ivector(Tage,1,NCOVMAX);
 6811: 
 6812:     free_imatrix(nbcode,0,NCOVMAX,0,NCOVMAX);
 6813:     free_imatrix(codtab,1,100,1,10);
 6814:   fflush(fichtm);
 6815:   fflush(ficgp);
 6816:   
 6817: 
 6818:   if((nberr >0) || (nbwarn>0)){
 6819:     printf("End of Imach with %d errors and/or %d warnings\n",nberr,nbwarn);
 6820:     fprintf(ficlog,"End of Imach with %d errors and/or warnings %d\n",nberr,nbwarn);
 6821:   }else{
 6822:     printf("End of Imach\n");
 6823:     fprintf(ficlog,"End of Imach\n");
 6824:   }
 6825:   printf("See log file on %s\n",filelog);
 6826:   /*  gettimeofday(&end_time, (struct timezone*)0);*/  /* after time */
 6827:   /*(void) gettimeofday(&end_time,&tzp);*/
 6828:   rend_time = time(NULL);  
 6829:   end_time = *localtime(&rend_time);
 6830:   /* tml = *localtime(&end_time.tm_sec); */
 6831:   strcpy(strtend,asctime(&end_time));
 6832:   printf("Local time at start %s\nLocal time at end   %s",strstart, strtend); 
 6833:   fprintf(ficlog,"Local time at start %s\nLocal time at end   %s\n",strstart, strtend); 
 6834:   printf("Total time used %s\n", asc_diff_time(rend_time -rstart_time,tmpout));
 6835: 
 6836:   printf("Total time was %.0lf Sec.\n", difftime(rend_time,rstart_time));
 6837:   fprintf(ficlog,"Total time used %s\n", asc_diff_time(rend_time -rstart_time,tmpout));
 6838:   fprintf(ficlog,"Total time was %.0lf Sec.\n", difftime(rend_time,rstart_time));
 6839:   /*  printf("Total time was %d uSec.\n", total_usecs);*/
 6840: /*   if(fileappend(fichtm,optionfilehtm)){ */
 6841:   fprintf(fichtm,"<br>Local time at start %s<br>Local time at end   %s<br>\n</body></html>",strstart, strtend);
 6842:   fclose(fichtm);
 6843:   fprintf(fichtmcov,"<br>Local time at start %s<br>Local time at end   %s<br>\n</body></html>",strstart, strtend);
 6844:   fclose(fichtmcov);
 6845:   fclose(ficgp);
 6846:   fclose(ficlog);
 6847:   /*------ End -----------*/
 6848: 
 6849: 
 6850:    printf("Before Current directory %s!\n",pathcd);
 6851:    if(chdir(pathcd) != 0)
 6852:     printf("Can't move to directory %s!\n",path);
 6853:   if(getcwd(pathcd,MAXLINE) > 0)
 6854:     printf("Current directory %s!\n",pathcd);
 6855:   /*strcat(plotcmd,CHARSEPARATOR);*/
 6856:   sprintf(plotcmd,"gnuplot");
 6857: #ifdef _WIN32
 6858:   sprintf(plotcmd,"\"%sgnuplot.exe\"",pathimach);
 6859: #endif
 6860:   if(!stat(plotcmd,&info)){
 6861:     printf("Error or gnuplot program not found: '%s'\n",plotcmd);fflush(stdout);
 6862:     if(!stat(getenv("GNUPLOTBIN"),&info)){
 6863:       printf("Error or gnuplot program not found: '%s' Environment GNUPLOTBIN not set.\n",plotcmd);fflush(stdout);
 6864:     }else
 6865:       strcpy(pplotcmd,plotcmd);
 6866: #ifdef __unix
 6867:     strcpy(plotcmd,GNUPLOTPROGRAM);
 6868:     if(!stat(plotcmd,&info)){
 6869:       printf("Error gnuplot program not found: '%s'\n",plotcmd);fflush(stdout);
 6870:     }else
 6871:       strcpy(pplotcmd,plotcmd);
 6872: #endif
 6873:   }else
 6874:     strcpy(pplotcmd,plotcmd);
 6875:   
 6876:   sprintf(plotcmd,"%s %s",pplotcmd, optionfilegnuplot);
 6877:   printf("Starting graphs with: '%s'\n",plotcmd);fflush(stdout);
 6878: 
 6879:   if((outcmd=system(plotcmd)) != 0){
 6880:     printf("gnuplot command might not be in your path: '%s', err=%d\n", plotcmd, outcmd);
 6881:     printf("\n Trying if gnuplot resides on the same directory that IMaCh\n");
 6882:     sprintf(plotcmd,"%sgnuplot %s", pathimach, optionfilegnuplot);
 6883:     if((outcmd=system(plotcmd)) != 0)
 6884:       printf("\n Still a problem with gnuplot command %s, err=%d\n", plotcmd, outcmd);
 6885:   }
 6886:   printf(" Successful, please wait...");
 6887:   while (z[0] != 'q') {
 6888:     /* chdir(path); */
 6889:     printf("\nType e to edit results with your browser, g to graph again and q for exit: ");
 6890:     scanf("%s",z);
 6891: /*     if (z[0] == 'c') system("./imach"); */
 6892:     if (z[0] == 'e') {
 6893: #ifdef __APPLE__
 6894:       sprintf(pplotcmd, "open %s", optionfilehtm);
 6895: #elif __linux
 6896:       sprintf(pplotcmd, "xdg-open %s", optionfilehtm);
 6897: #else
 6898:       sprintf(pplotcmd, "%s", optionfilehtm);
 6899: #endif
 6900:       printf("Starting browser with: %s",pplotcmd);fflush(stdout);
 6901:       system(pplotcmd);
 6902:     }
 6903:     else if (z[0] == 'g') system(plotcmd);
 6904:     else if (z[0] == 'q') exit(0);
 6905:   }
 6906:   end:
 6907:   while (z[0] != 'q') {
 6908:     printf("\nType  q for exiting: ");
 6909:     scanf("%s",z);
 6910:   }
 6911: }

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