File:  [Local Repository] / imach / src / imach.c
Revision 1.171: download - view: text, annotated - select for diffs
Tue Dec 23 13:26:59 2014 UTC (9 years, 5 months ago) by brouard
Branches: MAIN
CVS tags: HEAD
Summary: Back from Visual C

Still problem with utsname.h on Windows

    1: /* $Id: imach.c,v 1.171 2014/12/23 13:26:59 brouard Exp $
    2:   $State: Exp $
    3:   $Log: imach.c,v $
    4:   Revision 1.171  2014/12/23 13:26:59  brouard
    5:   Summary: Back from Visual C
    6: 
    7:   Still problem with utsname.h on Windows
    8: 
    9:   Revision 1.170  2014/12/23 11:17:12  brouard
   10:   Summary: Cleaning some \%% back to %%
   11: 
   12:   The escape was mandatory for a specific compiler (which one?), but too many warnings.
   13: 
   14:   Revision 1.169  2014/12/22 23:08:31  brouard
   15:   Summary: 0.98p
   16: 
   17:   Outputs some informations on compiler used, OS etc. Testing on different platforms.
   18: 
   19:   Revision 1.168  2014/12/22 15:17:42  brouard
   20:   Summary: update
   21: 
   22:   Revision 1.167  2014/12/22 13:50:56  brouard
   23:   Summary: Testing uname and compiler version and if compiled 32 or 64
   24: 
   25:   Testing on Linux 64
   26: 
   27:   Revision 1.166  2014/12/22 11:40:47  brouard
   28:   *** empty log message ***
   29: 
   30:   Revision 1.165  2014/12/16 11:20:36  brouard
   31:   Summary: After compiling on Visual C
   32: 
   33:   * imach.c (Module): Merging 1.61 to 1.162
   34: 
   35:   Revision 1.164  2014/12/16 10:52:11  brouard
   36:   Summary: Merging with Visual C after suppressing some warnings for unused variables. Also fixing Saito's bug 0.98Xn
   37: 
   38:   * imach.c (Module): Merging 1.61 to 1.162
   39: 
   40:   Revision 1.163  2014/12/16 10:30:11  brouard
   41:   * imach.c (Module): Merging 1.61 to 1.162
   42: 
   43:   Revision 1.162  2014/09/25 11:43:39  brouard
   44:   Summary: temporary backup 0.99!
   45: 
   46:   Revision 1.1  2014/09/16 11:06:58  brouard
   47:   Summary: With some code (wrong) for nlopt
   48: 
   49:   Author:
   50: 
   51:   Revision 1.161  2014/09/15 20:41:41  brouard
   52:   Summary: Problem with macro SQR on Intel compiler
   53: 
   54:   Revision 1.160  2014/09/02 09:24:05  brouard
   55:   *** empty log message ***
   56: 
   57:   Revision 1.159  2014/09/01 10:34:10  brouard
   58:   Summary: WIN32
   59:   Author: Brouard
   60: 
   61:   Revision 1.158  2014/08/27 17:11:51  brouard
   62:   *** empty log message ***
   63: 
   64:   Revision 1.157  2014/08/27 16:26:55  brouard
   65:   Summary: Preparing windows Visual studio version
   66:   Author: Brouard
   67: 
   68:   In order to compile on Visual studio, time.h is now correct and time_t
   69:   and tm struct should be used. difftime should be used but sometimes I
   70:   just make the differences in raw time format (time(&now).
   71:   Trying to suppress #ifdef LINUX
   72:   Add xdg-open for __linux in order to open default browser.
   73: 
   74:   Revision 1.156  2014/08/25 20:10:10  brouard
   75:   *** empty log message ***
   76: 
   77:   Revision 1.155  2014/08/25 18:32:34  brouard
   78:   Summary: New compile, minor changes
   79:   Author: Brouard
   80: 
   81:   Revision 1.154  2014/06/20 17:32:08  brouard
   82:   Summary: Outputs now all graphs of convergence to period prevalence
   83: 
   84:   Revision 1.153  2014/06/20 16:45:46  brouard
   85:   Summary: If 3 live state, convergence to period prevalence on same graph
   86:   Author: Brouard
   87: 
   88:   Revision 1.152  2014/06/18 17:54:09  brouard
   89:   Summary: open browser, use gnuplot on same dir than imach if not found in the path
   90: 
   91:   Revision 1.151  2014/06/18 16:43:30  brouard
   92:   *** empty log message ***
   93: 
   94:   Revision 1.150  2014/06/18 16:42:35  brouard
   95:   Summary: If gnuplot is not in the path try on same directory than imach binary (OSX)
   96:   Author: brouard
   97: 
   98:   Revision 1.149  2014/06/18 15:51:14  brouard
   99:   Summary: Some fixes in parameter files errors
  100:   Author: Nicolas Brouard
  101: 
  102:   Revision 1.148  2014/06/17 17:38:48  brouard
  103:   Summary: Nothing new
  104:   Author: Brouard
  105: 
  106:   Just a new packaging for OS/X version 0.98nS
  107: 
  108:   Revision 1.147  2014/06/16 10:33:11  brouard
  109:   *** empty log message ***
  110: 
  111:   Revision 1.146  2014/06/16 10:20:28  brouard
  112:   Summary: Merge
  113:   Author: Brouard
  114: 
  115:   Merge, before building revised version.
  116: 
  117:   Revision 1.145  2014/06/10 21:23:15  brouard
  118:   Summary: Debugging with valgrind
  119:   Author: Nicolas Brouard
  120: 
  121:   Lot of changes in order to output the results with some covariates
  122:   After the Edimburgh REVES conference 2014, it seems mandatory to
  123:   improve the code.
  124:   No more memory valgrind error but a lot has to be done in order to
  125:   continue the work of splitting the code into subroutines.
  126:   Also, decodemodel has been improved. Tricode is still not
  127:   optimal. nbcode should be improved. Documentation has been added in
  128:   the source code.
  129: 
  130:   Revision 1.143  2014/01/26 09:45:38  brouard
  131:   Summary: Version 0.98nR (to be improved, but gives same optimization results as 0.98k. Nice, promising
  132: 
  133:   * imach.c (Module): Trying to merge old staffs together while being at Tokyo. Not tested...
  134:   (Module): Version 0.98nR Running ok, but output format still only works for three covariates.
  135: 
  136:   Revision 1.142  2014/01/26 03:57:36  brouard
  137:   Summary: gnuplot changed plot w l 1 has to be changed to plot w l lt 2
  138: 
  139:   * imach.c (Module): Trying to merge old staffs together while being at Tokyo. Not tested...
  140: 
  141:   Revision 1.141  2014/01/26 02:42:01  brouard
  142:   * imach.c (Module): Trying to merge old staffs together while being at Tokyo. Not tested...
  143: 
  144:   Revision 1.140  2011/09/02 10:37:54  brouard
  145:   Summary: times.h is ok with mingw32 now.
  146: 
  147:   Revision 1.139  2010/06/14 07:50:17  brouard
  148:   After the theft of my laptop, I probably lost some lines of codes which were not uploaded to the CVS tree.
  149:   I remember having already fixed agemin agemax which are pointers now but not cvs saved.
  150: 
  151:   Revision 1.138  2010/04/30 18:19:40  brouard
  152:   *** empty log message ***
  153: 
  154:   Revision 1.137  2010/04/29 18:11:38  brouard
  155:   (Module): Checking covariates for more complex models
  156:   than V1+V2. A lot of change to be done. Unstable.
  157: 
  158:   Revision 1.136  2010/04/26 20:30:53  brouard
  159:   (Module): merging some libgsl code. Fixing computation
  160:   of likelione (using inter/intrapolation if mle = 0) in order to
  161:   get same likelihood as if mle=1.
  162:   Some cleaning of code and comments added.
  163: 
  164:   Revision 1.135  2009/10/29 15:33:14  brouard
  165:   (Module): Now imach stops if date of birth, at least year of birth, is not given. Some cleaning of the code.
  166: 
  167:   Revision 1.134  2009/10/29 13:18:53  brouard
  168:   (Module): Now imach stops if date of birth, at least year of birth, is not given. Some cleaning of the code.
  169: 
  170:   Revision 1.133  2009/07/06 10:21:25  brouard
  171:   just nforces
  172: 
  173:   Revision 1.132  2009/07/06 08:22:05  brouard
  174:   Many tings
  175: 
  176:   Revision 1.131  2009/06/20 16:22:47  brouard
  177:   Some dimensions resccaled
  178: 
  179:   Revision 1.130  2009/05/26 06:44:34  brouard
  180:   (Module): Max Covariate is now set to 20 instead of 8. A
  181:   lot of cleaning with variables initialized to 0. Trying to make
  182:   V2+V3*age+V1+V4 strb=V3*age+V1+V4 working better.
  183: 
  184:   Revision 1.129  2007/08/31 13:49:27  lievre
  185:   Modification of the way of exiting when the covariate is not binary in order to see on the window the error message before exiting
  186: 
  187:   Revision 1.128  2006/06/30 13:02:05  brouard
  188:   (Module): Clarifications on computing e.j
  189: 
  190:   Revision 1.127  2006/04/28 18:11:50  brouard
  191:   (Module): Yes the sum of survivors was wrong since
  192:   imach-114 because nhstepm was no more computed in the age
  193:   loop. Now we define nhstepma in the age loop.
  194:   (Module): In order to speed up (in case of numerous covariates) we
  195:   compute health expectancies (without variances) in a first step
  196:   and then all the health expectancies with variances or standard
  197:   deviation (needs data from the Hessian matrices) which slows the
  198:   computation.
  199:   In the future we should be able to stop the program is only health
  200:   expectancies and graph are needed without standard deviations.
  201: 
  202:   Revision 1.126  2006/04/28 17:23:28  brouard
  203:   (Module): Yes the sum of survivors was wrong since
  204:   imach-114 because nhstepm was no more computed in the age
  205:   loop. Now we define nhstepma in the age loop.
  206:   Version 0.98h
  207: 
  208:   Revision 1.125  2006/04/04 15:20:31  lievre
  209:   Errors in calculation of health expectancies. Age was not initialized.
  210:   Forecasting file added.
  211: 
  212:   Revision 1.124  2006/03/22 17:13:53  lievre
  213:   Parameters are printed with %lf instead of %f (more numbers after the comma).
  214:   The log-likelihood is printed in the log file
  215: 
  216:   Revision 1.123  2006/03/20 10:52:43  brouard
  217:   * imach.c (Module): <title> changed, corresponds to .htm file
  218:   name. <head> headers where missing.
  219: 
  220:   * imach.c (Module): Weights can have a decimal point as for
  221:   English (a comma might work with a correct LC_NUMERIC environment,
  222:   otherwise the weight is truncated).
  223:   Modification of warning when the covariates values are not 0 or
  224:   1.
  225:   Version 0.98g
  226: 
  227:   Revision 1.122  2006/03/20 09:45:41  brouard
  228:   (Module): Weights can have a decimal point as for
  229:   English (a comma might work with a correct LC_NUMERIC environment,
  230:   otherwise the weight is truncated).
  231:   Modification of warning when the covariates values are not 0 or
  232:   1.
  233:   Version 0.98g
  234: 
  235:   Revision 1.121  2006/03/16 17:45:01  lievre
  236:   * imach.c (Module): Comments concerning covariates added
  237: 
  238:   * imach.c (Module): refinements in the computation of lli if
  239:   status=-2 in order to have more reliable computation if stepm is
  240:   not 1 month. Version 0.98f
  241: 
  242:   Revision 1.120  2006/03/16 15:10:38  lievre
  243:   (Module): refinements in the computation of lli if
  244:   status=-2 in order to have more reliable computation if stepm is
  245:   not 1 month. Version 0.98f
  246: 
  247:   Revision 1.119  2006/03/15 17:42:26  brouard
  248:   (Module): Bug if status = -2, the loglikelihood was
  249:   computed as likelihood omitting the logarithm. Version O.98e
  250: 
  251:   Revision 1.118  2006/03/14 18:20:07  brouard
  252:   (Module): varevsij Comments added explaining the second
  253:   table of variances if popbased=1 .
  254:   (Module): Covariances of eij, ekl added, graphs fixed, new html link.
  255:   (Module): Function pstamp added
  256:   (Module): Version 0.98d
  257: 
  258:   Revision 1.117  2006/03/14 17:16:22  brouard
  259:   (Module): varevsij Comments added explaining the second
  260:   table of variances if popbased=1 .
  261:   (Module): Covariances of eij, ekl added, graphs fixed, new html link.
  262:   (Module): Function pstamp added
  263:   (Module): Version 0.98d
  264: 
  265:   Revision 1.116  2006/03/06 10:29:27  brouard
  266:   (Module): Variance-covariance wrong links and
  267:   varian-covariance of ej. is needed (Saito).
  268: 
  269:   Revision 1.115  2006/02/27 12:17:45  brouard
  270:   (Module): One freematrix added in mlikeli! 0.98c
  271: 
  272:   Revision 1.114  2006/02/26 12:57:58  brouard
  273:   (Module): Some improvements in processing parameter
  274:   filename with strsep.
  275: 
  276:   Revision 1.113  2006/02/24 14:20:24  brouard
  277:   (Module): Memory leaks checks with valgrind and:
  278:   datafile was not closed, some imatrix were not freed and on matrix
  279:   allocation too.
  280: 
  281:   Revision 1.112  2006/01/30 09:55:26  brouard
  282:   (Module): Back to gnuplot.exe instead of wgnuplot.exe
  283: 
  284:   Revision 1.111  2006/01/25 20:38:18  brouard
  285:   (Module): Lots of cleaning and bugs added (Gompertz)
  286:   (Module): Comments can be added in data file. Missing date values
  287:   can be a simple dot '.'.
  288: 
  289:   Revision 1.110  2006/01/25 00:51:50  brouard
  290:   (Module): Lots of cleaning and bugs added (Gompertz)
  291: 
  292:   Revision 1.109  2006/01/24 19:37:15  brouard
  293:   (Module): Comments (lines starting with a #) are allowed in data.
  294: 
  295:   Revision 1.108  2006/01/19 18:05:42  lievre
  296:   Gnuplot problem appeared...
  297:   To be fixed
  298: 
  299:   Revision 1.107  2006/01/19 16:20:37  brouard
  300:   Test existence of gnuplot in imach path
  301: 
  302:   Revision 1.106  2006/01/19 13:24:36  brouard
  303:   Some cleaning and links added in html output
  304: 
  305:   Revision 1.105  2006/01/05 20:23:19  lievre
  306:   *** empty log message ***
  307: 
  308:   Revision 1.104  2005/09/30 16:11:43  lievre
  309:   (Module): sump fixed, loop imx fixed, and simplifications.
  310:   (Module): If the status is missing at the last wave but we know
  311:   that the person is alive, then we can code his/her status as -2
  312:   (instead of missing=-1 in earlier versions) and his/her
  313:   contributions to the likelihood is 1 - Prob of dying from last
  314:   health status (= 1-p13= p11+p12 in the easiest case of somebody in
  315:   the healthy state at last known wave). Version is 0.98
  316: 
  317:   Revision 1.103  2005/09/30 15:54:49  lievre
  318:   (Module): sump fixed, loop imx fixed, and simplifications.
  319: 
  320:   Revision 1.102  2004/09/15 17:31:30  brouard
  321:   Add the possibility to read data file including tab characters.
  322: 
  323:   Revision 1.101  2004/09/15 10:38:38  brouard
  324:   Fix on curr_time
  325: 
  326:   Revision 1.100  2004/07/12 18:29:06  brouard
  327:   Add version for Mac OS X. Just define UNIX in Makefile
  328: 
  329:   Revision 1.99  2004/06/05 08:57:40  brouard
  330:   *** empty log message ***
  331: 
  332:   Revision 1.98  2004/05/16 15:05:56  brouard
  333:   New version 0.97 . First attempt to estimate force of mortality
  334:   directly from the data i.e. without the need of knowing the health
  335:   state at each age, but using a Gompertz model: log u =a + b*age .
  336:   This is the basic analysis of mortality and should be done before any
  337:   other analysis, in order to test if the mortality estimated from the
  338:   cross-longitudinal survey is different from the mortality estimated
  339:   from other sources like vital statistic data.
  340: 
  341:   The same imach parameter file can be used but the option for mle should be -3.
  342: 
  343:   Agnès, who wrote this part of the code, tried to keep most of the
  344:   former routines in order to include the new code within the former code.
  345: 
  346:   The output is very simple: only an estimate of the intercept and of
  347:   the slope with 95% confident intervals.
  348: 
  349:   Current limitations:
  350:   A) Even if you enter covariates, i.e. with the
  351:   model= V1+V2 equation for example, the programm does only estimate a unique global model without covariates.
  352:   B) There is no computation of Life Expectancy nor Life Table.
  353: 
  354:   Revision 1.97  2004/02/20 13:25:42  lievre
  355:   Version 0.96d. Population forecasting command line is (temporarily)
  356:   suppressed.
  357: 
  358:   Revision 1.96  2003/07/15 15:38:55  brouard
  359:   * imach.c (Repository): Errors in subdirf, 2, 3 while printing tmpout is
  360:   rewritten within the same printf. Workaround: many printfs.
  361: 
  362:   Revision 1.95  2003/07/08 07:54:34  brouard
  363:   * imach.c (Repository):
  364:   (Repository): Using imachwizard code to output a more meaningful covariance
  365:   matrix (cov(a12,c31) instead of numbers.
  366: 
  367:   Revision 1.94  2003/06/27 13:00:02  brouard
  368:   Just cleaning
  369: 
  370:   Revision 1.93  2003/06/25 16:33:55  brouard
  371:   (Module): On windows (cygwin) function asctime_r doesn't
  372:   exist so I changed back to asctime which exists.
  373:   (Module): Version 0.96b
  374: 
  375:   Revision 1.92  2003/06/25 16:30:45  brouard
  376:   (Module): On windows (cygwin) function asctime_r doesn't
  377:   exist so I changed back to asctime which exists.
  378: 
  379:   Revision 1.91  2003/06/25 15:30:29  brouard
  380:   * imach.c (Repository): Duplicated warning errors corrected.
  381:   (Repository): Elapsed time after each iteration is now output. It
  382:   helps to forecast when convergence will be reached. Elapsed time
  383:   is stamped in powell.  We created a new html file for the graphs
  384:   concerning matrix of covariance. It has extension -cov.htm.
  385: 
  386:   Revision 1.90  2003/06/24 12:34:15  brouard
  387:   (Module): Some bugs corrected for windows. Also, when
  388:   mle=-1 a template is output in file "or"mypar.txt with the design
  389:   of the covariance matrix to be input.
  390: 
  391:   Revision 1.89  2003/06/24 12:30:52  brouard
  392:   (Module): Some bugs corrected for windows. Also, when
  393:   mle=-1 a template is output in file "or"mypar.txt with the design
  394:   of the covariance matrix to be input.
  395: 
  396:   Revision 1.88  2003/06/23 17:54:56  brouard
  397:   * 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.
  398: 
  399:   Revision 1.87  2003/06/18 12:26:01  brouard
  400:   Version 0.96
  401: 
  402:   Revision 1.86  2003/06/17 20:04:08  brouard
  403:   (Module): Change position of html and gnuplot routines and added
  404:   routine fileappend.
  405: 
  406:   Revision 1.85  2003/06/17 13:12:43  brouard
  407:   * imach.c (Repository): Check when date of death was earlier that
  408:   current date of interview. It may happen when the death was just
  409:   prior to the death. In this case, dh was negative and likelihood
  410:   was wrong (infinity). We still send an "Error" but patch by
  411:   assuming that the date of death was just one stepm after the
  412:   interview.
  413:   (Repository): Because some people have very long ID (first column)
  414:   we changed int to long in num[] and we added a new lvector for
  415:   memory allocation. But we also truncated to 8 characters (left
  416:   truncation)
  417:   (Repository): No more line truncation errors.
  418: 
  419:   Revision 1.84  2003/06/13 21:44:43  brouard
  420:   * imach.c (Repository): Replace "freqsummary" at a correct
  421:   place. It differs from routine "prevalence" which may be called
  422:   many times. Probs is memory consuming and must be used with
  423:   parcimony.
  424:   Version 0.95a3 (should output exactly the same maximization than 0.8a2)
  425: 
  426:   Revision 1.83  2003/06/10 13:39:11  lievre
  427:   *** empty log message ***
  428: 
  429:   Revision 1.82  2003/06/05 15:57:20  brouard
  430:   Add log in  imach.c and  fullversion number is now printed.
  431: 
  432: */
  433: /*
  434:    Interpolated Markov Chain
  435: 
  436:   Short summary of the programme:
  437:   
  438:   This program computes Healthy Life Expectancies from
  439:   cross-longitudinal data. Cross-longitudinal data consist in: -1- a
  440:   first survey ("cross") where individuals from different ages are
  441:   interviewed on their health status or degree of disability (in the
  442:   case of a health survey which is our main interest) -2- at least a
  443:   second wave of interviews ("longitudinal") which measure each change
  444:   (if any) in individual health status.  Health expectancies are
  445:   computed from the time spent in each health state according to a
  446:   model. More health states you consider, more time is necessary to reach the
  447:   Maximum Likelihood of the parameters involved in the model.  The
  448:   simplest model is the multinomial logistic model where pij is the
  449:   probability to be observed in state j at the second wave
  450:   conditional to be observed in state i at the first wave. Therefore
  451:   the model is: log(pij/pii)= aij + bij*age+ cij*sex + etc , where
  452:   'age' is age and 'sex' is a covariate. If you want to have a more
  453:   complex model than "constant and age", you should modify the program
  454:   where the markup *Covariates have to be included here again* invites
  455:   you to do it.  More covariates you add, slower the
  456:   convergence.
  457: 
  458:   The advantage of this computer programme, compared to a simple
  459:   multinomial logistic model, is clear when the delay between waves is not
  460:   identical for each individual. Also, if a individual missed an
  461:   intermediate interview, the information is lost, but taken into
  462:   account using an interpolation or extrapolation.  
  463: 
  464:   hPijx is the probability to be observed in state i at age x+h
  465:   conditional to the observed state i at age x. The delay 'h' can be
  466:   split into an exact number (nh*stepm) of unobserved intermediate
  467:   states. This elementary transition (by month, quarter,
  468:   semester or year) is modelled as a multinomial logistic.  The hPx
  469:   matrix is simply the matrix product of nh*stepm elementary matrices
  470:   and the contribution of each individual to the likelihood is simply
  471:   hPijx.
  472: 
  473:   Also this programme outputs the covariance matrix of the parameters but also
  474:   of the life expectancies. It also computes the period (stable) prevalence. 
  475:   
  476:   Authors: Nicolas Brouard (brouard@ined.fr) and Agnès Lièvre (lievre@ined.fr).
  477:            Institut national d'études démographiques, Paris.
  478:   This software have been partly granted by Euro-REVES, a concerted action
  479:   from the European Union.
  480:   It is copyrighted identically to a GNU software product, ie programme and
  481:   software can be distributed freely for non commercial use. Latest version
  482:   can be accessed at http://euroreves.ined.fr/imach .
  483: 
  484:   Help to debug: LD_PRELOAD=/usr/local/lib/libnjamd.so ./imach foo.imach
  485:   or better on gdb : set env LD_PRELOAD=/usr/local/lib/libnjamd.so
  486:   
  487:   **********************************************************************/
  488: /*
  489:   main
  490:   read parameterfile
  491:   read datafile
  492:   concatwav
  493:   freqsummary
  494:   if (mle >= 1)
  495:     mlikeli
  496:   print results files
  497:   if mle==1 
  498:      computes hessian
  499:   read end of parameter file: agemin, agemax, bage, fage, estepm
  500:       begin-prev-date,...
  501:   open gnuplot file
  502:   open html file
  503:   period (stable) prevalence      | pl_nom    1-1 2-2 etc by covariate
  504:    for age prevalim()             | #****** V1=0  V2=1  V3=1  V4=0 ******
  505:                                   | 65 1 0 2 1 3 1 4 0  0.96326 0.03674
  506:     freexexit2 possible for memory heap.
  507: 
  508:   h Pij x                         | pij_nom  ficrestpij
  509:    # Cov Agex agex+h hpijx with i,j= 1-1 1-2     1-3     2-1     2-2     2-3
  510:        1  85   85    1.00000             0.00000 0.00000 0.00000 1.00000 0.00000
  511:        1  85   86    0.68299             0.22291 0.09410 0.71093 0.00000 0.28907
  512: 
  513:        1  65   99    0.00364             0.00322 0.99314 0.00350 0.00310 0.99340
  514:        1  65  100    0.00214             0.00204 0.99581 0.00206 0.00196 0.99597
  515:   variance of p one-step probabilities varprob  | prob_nom   ficresprob #One-step probabilities and stand. devi in ()
  516:    Standard deviation of one-step probabilities | probcor_nom   ficresprobcor #One-step probabilities and correlation matrix
  517:    Matrix of variance covariance of one-step probabilities |  probcov_nom ficresprobcov #One-step probabilities and covariance matrix
  518: 
  519:   forecasting if prevfcast==1 prevforecast call prevalence()
  520:   health expectancies
  521:   Variance-covariance of DFLE
  522:   prevalence()
  523:    movingaverage()
  524:   varevsij() 
  525:   if popbased==1 varevsij(,popbased)
  526:   total life expectancies
  527:   Variance of period (stable) prevalence
  528:  end
  529: */
  530: 
  531: #define POWELL /* Instead of NLOPT */
  532: 
  533: #include <math.h>
  534: #include <stdio.h>
  535: #include <stdlib.h>
  536: #include <string.h>
  537: 
  538: #ifdef _WIN32
  539: #include <io.h>
  540: #else
  541: #include <unistd.h>
  542: #endif
  543: 
  544: #include <limits.h>
  545: #include <sys/types.h>
  546: 
  547: #if defined(__GNUC__)
  548: #include <sys/utsname.h> /* Doesn't work on Windows */
  549: #endif
  550: 
  551: #include <sys/stat.h>
  552: #include <errno.h>
  553: /* extern int errno; */
  554: 
  555: /* #ifdef LINUX */
  556: /* #include <time.h> */
  557: /* #include "timeval.h" */
  558: /* #else */
  559: /* #include <sys/time.h> */
  560: /* #endif */
  561: 
  562: #include <time.h>
  563: 
  564: #ifdef GSL
  565: #include <gsl/gsl_errno.h>
  566: #include <gsl/gsl_multimin.h>
  567: #endif
  568: 
  569: 
  570: #ifdef NLOPT
  571: #include <nlopt.h>
  572: typedef struct {
  573:   double (* function)(double [] );
  574: } myfunc_data ;
  575: #endif
  576: 
  577: /* #include <libintl.h> */
  578: /* #define _(String) gettext (String) */
  579: 
  580: #define MAXLINE 1024 /* Was 256. Overflow with 312 with 2 states and 4 covariates. Should be ok */
  581: 
  582: #define GNUPLOTPROGRAM "gnuplot"
  583: /*#define GNUPLOTPROGRAM "..\\gp37mgw\\wgnuplot"*/
  584: #define FILENAMELENGTH 132
  585: 
  586: #define	GLOCK_ERROR_NOPATH		-1	/* empty path */
  587: #define	GLOCK_ERROR_GETCWD		-2	/* cannot get cwd */
  588: 
  589: #define MAXPARM 128 /**< Maximum number of parameters for the optimization */
  590: #define NPARMAX 64 /**< (nlstate+ndeath-1)*nlstate*ncovmodel */
  591: 
  592: #define NINTERVMAX 8
  593: #define NLSTATEMAX 8 /**< Maximum number of live states (for func) */
  594: #define NDEATHMAX 8 /**< Maximum number of dead states (for func) */
  595: #define NCOVMAX 20 /**< Maximum number of covariates, including generated covariates V1*V2 */
  596: #define codtabm(h,k)  1 & (h-1) >> (k-1) ;
  597: #define MAXN 20000
  598: #define YEARM 12. /**< Number of months per year */
  599: #define AGESUP 130
  600: #define AGEBASE 40
  601: #define AGEGOMP 10 /**< Minimal age for Gompertz adjustment */
  602: #ifdef _WIN32
  603: #define DIRSEPARATOR '\\'
  604: #define CHARSEPARATOR "\\"
  605: #define ODIRSEPARATOR '/'
  606: #else
  607: #define DIRSEPARATOR '/'
  608: #define CHARSEPARATOR "/"
  609: #define ODIRSEPARATOR '\\'
  610: #endif
  611: 
  612: /* $Id: imach.c,v 1.171 2014/12/23 13:26:59 brouard Exp $ */
  613: /* $State: Exp $ */
  614: 
  615: char version[]="Imach version 0.98p, December 2014,INED-EUROREVES-Institut de longevite-Japan Society for the Promotion of Science (Grant-in-Aid for Scientific Research 25293121), Intel Software 2015";
  616: char fullversion[]="$Revision: 1.171 $ $Date: 2014/12/23 13:26:59 $"; 
  617: char strstart[80];
  618: char optionfilext[10], optionfilefiname[FILENAMELENGTH];
  619: int erreur=0, nberr=0, nbwarn=0; /* Error number, number of errors number of warnings  */
  620: int nvar=0, nforce=0; /* Number of variables, number of forces */
  621: /* Number of covariates model=V2+V1+ V3*age+V2*V4 */
  622: int cptcovn=0; /**< cptcovn number of covariates added in the model (excepting constant and age and age*product) */
  623: int cptcovt=0; /**< cptcovt number of covariates added in the model (excepting constant and age) */
  624: int cptcovs=0; /**< cptcovs number of simple covariates V2+V1 =2 */
  625: int cptcovage=0; /**< Number of covariates with age: V3*age only =1 */
  626: int cptcovprodnoage=0; /**< Number of covariate products without age */   
  627: int cptcoveff=0; /* Total number of covariates to vary for printing results */
  628: int cptcov=0; /* Working variable */
  629: int npar=NPARMAX;
  630: int nlstate=2; /* Number of live states */
  631: int ndeath=1; /* Number of dead states */
  632: int ncovmodel=0, ncovcol=0;     /* Total number of covariables including constant a12*1 +b12*x ncovmodel=2 */
  633: int popbased=0;
  634: 
  635: int *wav; /* Number of waves for this individuual 0 is possible */
  636: int maxwav=0; /* Maxim number of waves */
  637: int jmin=0, jmax=0; /* min, max spacing between 2 waves */
  638: int ijmin=0, ijmax=0; /* Individuals having jmin and jmax */ 
  639: int gipmx=0, gsw=0; /* Global variables on the number of contributions 
  640: 		   to the likelihood and the sum of weights (done by funcone)*/
  641: int mle=1, weightopt=0;
  642: int **mw; /* mw[mi][i] is number of the mi wave for this individual */
  643: int **dh; /* dh[mi][i] is number of steps between mi,mi+1 for this individual */
  644: int **bh; /* bh[mi][i] is the bias (+ or -) for this individual if the delay between
  645: 	   * wave mi and wave mi+1 is not an exact multiple of stepm. */
  646: int countcallfunc=0;  /* Count the number of calls to func */
  647: double jmean=1; /* Mean space between 2 waves */
  648: double **matprod2(); /* test */
  649: double **oldm, **newm, **savm; /* Working pointers to matrices */
  650: double **oldms, **newms, **savms; /* Fixed working pointers to matrices */
  651: /*FILE *fic ; */ /* Used in readdata only */
  652: FILE *ficpar, *ficparo,*ficres, *ficresp, *ficrespl, *ficrespij, *ficrest,*ficresf,*ficrespop;
  653: FILE *ficlog, *ficrespow;
  654: int globpr=0; /* Global variable for printing or not */
  655: double fretone; /* Only one call to likelihood */
  656: long ipmx=0; /* Number of contributions */
  657: double sw; /* Sum of weights */
  658: char filerespow[FILENAMELENGTH];
  659: char fileresilk[FILENAMELENGTH]; /* File of individual contributions to the likelihood */
  660: FILE *ficresilk;
  661: FILE *ficgp,*ficresprob,*ficpop, *ficresprobcov, *ficresprobcor;
  662: FILE *ficresprobmorprev;
  663: FILE *fichtm, *fichtmcov; /* Html File */
  664: FILE *ficreseij;
  665: char filerese[FILENAMELENGTH];
  666: FILE *ficresstdeij;
  667: char fileresstde[FILENAMELENGTH];
  668: FILE *ficrescveij;
  669: char filerescve[FILENAMELENGTH];
  670: FILE  *ficresvij;
  671: char fileresv[FILENAMELENGTH];
  672: FILE  *ficresvpl;
  673: char fileresvpl[FILENAMELENGTH];
  674: char title[MAXLINE];
  675: char optionfile[FILENAMELENGTH], datafile[FILENAMELENGTH],  filerespl[FILENAMELENGTH];
  676: char plotcmd[FILENAMELENGTH], pplotcmd[FILENAMELENGTH];
  677: char tmpout[FILENAMELENGTH],  tmpout2[FILENAMELENGTH]; 
  678: char command[FILENAMELENGTH];
  679: int  outcmd=0;
  680: 
  681: char fileres[FILENAMELENGTH], filerespij[FILENAMELENGTH], filereso[FILENAMELENGTH], rfileres[FILENAMELENGTH];
  682: 
  683: char filelog[FILENAMELENGTH]; /* Log file */
  684: char filerest[FILENAMELENGTH];
  685: char fileregp[FILENAMELENGTH];
  686: char popfile[FILENAMELENGTH];
  687: 
  688: char optionfilegnuplot[FILENAMELENGTH], optionfilehtm[FILENAMELENGTH], optionfilehtmcov[FILENAMELENGTH] ;
  689: 
  690: /* struct timeval start_time, end_time, curr_time, last_time, forecast_time; */
  691: /* struct timezone tzp; */
  692: /* extern int gettimeofday(); */
  693: struct tm tml, *gmtime(), *localtime();
  694: 
  695: extern time_t time();
  696: 
  697: struct tm start_time, end_time, curr_time, last_time, forecast_time;
  698: time_t  rstart_time, rend_time, rcurr_time, rlast_time, rforecast_time; /* raw time */
  699: struct tm tm;
  700: 
  701: char strcurr[80], strfor[80];
  702: 
  703: char *endptr;
  704: long lval;
  705: double dval;
  706: 
  707: #define NR_END 1
  708: #define FREE_ARG char*
  709: #define FTOL 1.0e-10
  710: 
  711: #define NRANSI 
  712: #define ITMAX 200 
  713: 
  714: #define TOL 2.0e-4 
  715: 
  716: #define CGOLD 0.3819660 
  717: #define ZEPS 1.0e-10 
  718: #define SHFT(a,b,c,d) (a)=(b);(b)=(c);(c)=(d); 
  719: 
  720: #define GOLD 1.618034 
  721: #define GLIMIT 100.0 
  722: #define TINY 1.0e-20 
  723: 
  724: static double maxarg1,maxarg2;
  725: #define FMAX(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)>(maxarg2)? (maxarg1):(maxarg2))
  726: #define FMIN(a,b) (maxarg1=(a),maxarg2=(b),(maxarg1)<(maxarg2)? (maxarg1):(maxarg2))
  727:   
  728: #define SIGN(a,b) ((b)>0.0 ? fabs(a) : -fabs(a))
  729: #define rint(a) floor(a+0.5)
  730: /* http://www.thphys.uni-heidelberg.de/~robbers/cmbeasy/doc/html/myutils_8h-source.html */
  731: /* #define mytinydouble 1.0e-16 */
  732: /* #define DEQUAL(a,b) (fabs((a)-(b))<mytinydouble) */
  733: /* http://www.thphys.uni-heidelberg.de/~robbers/cmbeasy/doc/html/mynrutils_8h-source.html */
  734: /* static double dsqrarg; */
  735: /* #define DSQR(a) (DEQUAL((dsqrarg=(a)),0.0) ? 0.0 : dsqrarg*dsqrarg) */
  736: static double sqrarg;
  737: #define SQR(a) ((sqrarg=(a)) == 0.0 ? 0.0 :sqrarg*sqrarg)
  738: #define SWAP(a,b) {temp=(a);(a)=(b);(b)=temp;} 
  739: int agegomp= AGEGOMP;
  740: 
  741: int imx; 
  742: int stepm=1;
  743: /* Stepm, step in month: minimum step interpolation*/
  744: 
  745: int estepm;
  746: /* Estepm, step in month to interpolate survival function in order to approximate Life Expectancy*/
  747: 
  748: int m,nb;
  749: long *num;
  750: int firstpass=0, lastpass=4,*cod, *ncodemax, *Tage,*cens;
  751: double **agev,*moisnais, *annais, *moisdc, *andc,**mint, **anint;
  752: double **pmmij, ***probs;
  753: double *ageexmed,*agecens;
  754: double dateintmean=0;
  755: 
  756: double *weight;
  757: int **s; /* Status */
  758: double *agedc;
  759: double  **covar; /**< covar[j,i], value of jth covariate for individual i,
  760: 		  * covar=matrix(0,NCOVMAX,1,n); 
  761: 		  * cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2]; */
  762: double  idx; 
  763: int **nbcode, *Tvar; /**< model=V2 => Tvar[1]= 2 */
  764: int *Ndum; /** Freq of modality (tricode */
  765: int **codtab; /**< codtab=imatrix(1,100,1,10); */
  766: int **Tvard, *Tprod, cptcovprod, *Tvaraff;
  767: double *lsurv, *lpop, *tpop;
  768: 
  769: double ftol=FTOL; /**< Tolerance for computing Max Likelihood */
  770: double ftolhess; /**< Tolerance for computing hessian */
  771: 
  772: /**************** split *************************/
  773: static	int split( char *path, char *dirc, char *name, char *ext, char *finame )
  774: {
  775:   /* From a file name with (full) path (either Unix or Windows) we extract the directory (dirc)
  776:      the name of the file (name), its extension only (ext) and its first part of the name (finame)
  777:   */ 
  778:   char	*ss;				/* pointer */
  779:   int	l1, l2;				/* length counters */
  780: 
  781:   l1 = strlen(path );			/* length of path */
  782:   if ( l1 == 0 ) return( GLOCK_ERROR_NOPATH );
  783:   ss= strrchr( path, DIRSEPARATOR );		/* find last / */
  784:   if ( ss == NULL ) {			/* no directory, so determine current directory */
  785:     strcpy( name, path );		/* we got the fullname name because no directory */
  786:     /*if(strrchr(path, ODIRSEPARATOR )==NULL)
  787:       printf("Warning you should use %s as a separator\n",DIRSEPARATOR);*/
  788:     /* get current working directory */
  789:     /*    extern  char* getcwd ( char *buf , int len);*/
  790:     if ( getcwd( dirc, FILENAME_MAX ) == NULL ) {
  791:       return( GLOCK_ERROR_GETCWD );
  792:     }
  793:     /* got dirc from getcwd*/
  794:     printf(" DIRC = %s \n",dirc);
  795:   } else {				/* strip direcotry from path */
  796:     ss++;				/* after this, the filename */
  797:     l2 = strlen( ss );			/* length of filename */
  798:     if ( l2 == 0 ) return( GLOCK_ERROR_NOPATH );
  799:     strcpy( name, ss );		/* save file name */
  800:     strncpy( dirc, path, l1 - l2 );	/* now the directory */
  801:     dirc[l1-l2] = 0;			/* add zero */
  802:     printf(" DIRC2 = %s \n",dirc);
  803:   }
  804:   /* We add a separator at the end of dirc if not exists */
  805:   l1 = strlen( dirc );			/* length of directory */
  806:   if( dirc[l1-1] != DIRSEPARATOR ){
  807:     dirc[l1] =  DIRSEPARATOR;
  808:     dirc[l1+1] = 0; 
  809:     printf(" DIRC3 = %s \n",dirc);
  810:   }
  811:   ss = strrchr( name, '.' );		/* find last / */
  812:   if (ss >0){
  813:     ss++;
  814:     strcpy(ext,ss);			/* save extension */
  815:     l1= strlen( name);
  816:     l2= strlen(ss)+1;
  817:     strncpy( finame, name, l1-l2);
  818:     finame[l1-l2]= 0;
  819:   }
  820: 
  821:   return( 0 );				/* we're done */
  822: }
  823: 
  824: 
  825: /******************************************/
  826: 
  827: void replace_back_to_slash(char *s, char*t)
  828: {
  829:   int i;
  830:   int lg=0;
  831:   i=0;
  832:   lg=strlen(t);
  833:   for(i=0; i<= lg; i++) {
  834:     (s[i] = t[i]);
  835:     if (t[i]== '\\') s[i]='/';
  836:   }
  837: }
  838: 
  839: char *trimbb(char *out, char *in)
  840: { /* Trim multiple blanks in line but keeps first blanks if line starts with blanks */
  841:   char *s;
  842:   s=out;
  843:   while (*in != '\0'){
  844:     while( *in == ' ' && *(in+1) == ' '){ /* && *(in+1) != '\0'){*/
  845:       in++;
  846:     }
  847:     *out++ = *in++;
  848:   }
  849:   *out='\0';
  850:   return s;
  851: }
  852: 
  853: char *cutl(char *blocc, char *alocc, char *in, char occ)
  854: {
  855:   /* cuts string in into blocc and alocc where blocc ends before first occurence of char 'occ' 
  856:      and alocc starts after first occurence of char 'occ' : ex cutv(blocc,alocc,"abcdef2ghi2j",'2')
  857:      gives blocc="abcdef2ghi" and alocc="j".
  858:      If occ is not found blocc is null and alocc is equal to in. Returns blocc
  859:   */
  860:   char *s, *t;
  861:   t=in;s=in;
  862:   while ((*in != occ) && (*in != '\0')){
  863:     *alocc++ = *in++;
  864:   }
  865:   if( *in == occ){
  866:     *(alocc)='\0';
  867:     s=++in;
  868:   }
  869:  
  870:   if (s == t) {/* occ not found */
  871:     *(alocc-(in-s))='\0';
  872:     in=s;
  873:   }
  874:   while ( *in != '\0'){
  875:     *blocc++ = *in++;
  876:   }
  877: 
  878:   *blocc='\0';
  879:   return t;
  880: }
  881: char *cutv(char *blocc, char *alocc, char *in, char occ)
  882: {
  883:   /* cuts string in into blocc and alocc where blocc ends before last occurence of char 'occ' 
  884:      and alocc starts after last occurence of char 'occ' : ex cutv(blocc,alocc,"abcdef2ghi2j",'2')
  885:      gives blocc="abcdef2ghi" and alocc="j".
  886:      If occ is not found blocc is null and alocc is equal to in. Returns alocc
  887:   */
  888:   char *s, *t;
  889:   t=in;s=in;
  890:   while (*in != '\0'){
  891:     while( *in == occ){
  892:       *blocc++ = *in++;
  893:       s=in;
  894:     }
  895:     *blocc++ = *in++;
  896:   }
  897:   if (s == t) /* occ not found */
  898:     *(blocc-(in-s))='\0';
  899:   else
  900:     *(blocc-(in-s)-1)='\0';
  901:   in=s;
  902:   while ( *in != '\0'){
  903:     *alocc++ = *in++;
  904:   }
  905: 
  906:   *alocc='\0';
  907:   return s;
  908: }
  909: 
  910: int nbocc(char *s, char occ)
  911: {
  912:   int i,j=0;
  913:   int lg=20;
  914:   i=0;
  915:   lg=strlen(s);
  916:   for(i=0; i<= lg; i++) {
  917:   if  (s[i] == occ ) j++;
  918:   }
  919:   return j;
  920: }
  921: 
  922: /* void cutv(char *u,char *v, char*t, char occ) */
  923: /* { */
  924: /*   /\* cuts string t into u and v where u ends before last occurence of char 'occ'  */
  925: /*      and v starts after last occurence of char 'occ' : ex cutv(u,v,"abcdef2ghi2j",'2') */
  926: /*      gives u="abcdef2ghi" and v="j" *\/ */
  927: /*   int i,lg,j,p=0; */
  928: /*   i=0; */
  929: /*   lg=strlen(t); */
  930: /*   for(j=0; j<=lg-1; j++) { */
  931: /*     if((t[j]!= occ) && (t[j+1]== occ)) p=j+1; */
  932: /*   } */
  933: 
  934: /*   for(j=0; j<p; j++) { */
  935: /*     (u[j] = t[j]); */
  936: /*   } */
  937: /*      u[p]='\0'; */
  938: 
  939: /*    for(j=0; j<= lg; j++) { */
  940: /*     if (j>=(p+1))(v[j-p-1] = t[j]); */
  941: /*   } */
  942: /* } */
  943: 
  944: #ifdef _WIN32
  945: char * strsep(char **pp, const char *delim)
  946: {
  947:   char *p, *q;
  948:          
  949:   if ((p = *pp) == NULL)
  950:     return 0;
  951:   if ((q = strpbrk (p, delim)) != NULL)
  952:   {
  953:     *pp = q + 1;
  954:     *q = '\0';
  955:   }
  956:   else
  957:     *pp = 0;
  958:   return p;
  959: }
  960: #endif
  961: 
  962: /********************** nrerror ********************/
  963: 
  964: void nrerror(char error_text[])
  965: {
  966:   fprintf(stderr,"ERREUR ...\n");
  967:   fprintf(stderr,"%s\n",error_text);
  968:   exit(EXIT_FAILURE);
  969: }
  970: /*********************** vector *******************/
  971: double *vector(int nl, int nh)
  972: {
  973:   double *v;
  974:   v=(double *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(double)));
  975:   if (!v) nrerror("allocation failure in vector");
  976:   return v-nl+NR_END;
  977: }
  978: 
  979: /************************ free vector ******************/
  980: void free_vector(double*v, int nl, int nh)
  981: {
  982:   free((FREE_ARG)(v+nl-NR_END));
  983: }
  984: 
  985: /************************ivector *******************************/
  986: int *ivector(long nl,long nh)
  987: {
  988:   int *v;
  989:   v=(int *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(int)));
  990:   if (!v) nrerror("allocation failure in ivector");
  991:   return v-nl+NR_END;
  992: }
  993: 
  994: /******************free ivector **************************/
  995: void free_ivector(int *v, long nl, long nh)
  996: {
  997:   free((FREE_ARG)(v+nl-NR_END));
  998: }
  999: 
 1000: /************************lvector *******************************/
 1001: long *lvector(long nl,long nh)
 1002: {
 1003:   long *v;
 1004:   v=(long *) malloc((size_t)((nh-nl+1+NR_END)*sizeof(long)));
 1005:   if (!v) nrerror("allocation failure in ivector");
 1006:   return v-nl+NR_END;
 1007: }
 1008: 
 1009: /******************free lvector **************************/
 1010: void free_lvector(long *v, long nl, long nh)
 1011: {
 1012:   free((FREE_ARG)(v+nl-NR_END));
 1013: }
 1014: 
 1015: /******************* imatrix *******************************/
 1016: int **imatrix(long nrl, long nrh, long ncl, long nch) 
 1017:      /* allocate a int matrix with subscript range m[nrl..nrh][ncl..nch] */ 
 1018: { 
 1019:   long i, nrow=nrh-nrl+1,ncol=nch-ncl+1; 
 1020:   int **m; 
 1021:   
 1022:   /* allocate pointers to rows */ 
 1023:   m=(int **) malloc((size_t)((nrow+NR_END)*sizeof(int*))); 
 1024:   if (!m) nrerror("allocation failure 1 in matrix()"); 
 1025:   m += NR_END; 
 1026:   m -= nrl; 
 1027:   
 1028:   
 1029:   /* allocate rows and set pointers to them */ 
 1030:   m[nrl]=(int *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(int))); 
 1031:   if (!m[nrl]) nrerror("allocation failure 2 in matrix()"); 
 1032:   m[nrl] += NR_END; 
 1033:   m[nrl] -= ncl; 
 1034:   
 1035:   for(i=nrl+1;i<=nrh;i++) m[i]=m[i-1]+ncol; 
 1036:   
 1037:   /* return pointer to array of pointers to rows */ 
 1038:   return m; 
 1039: } 
 1040: 
 1041: /****************** free_imatrix *************************/
 1042: void free_imatrix(m,nrl,nrh,ncl,nch)
 1043:       int **m;
 1044:       long nch,ncl,nrh,nrl; 
 1045:      /* free an int matrix allocated by imatrix() */ 
 1046: { 
 1047:   free((FREE_ARG) (m[nrl]+ncl-NR_END)); 
 1048:   free((FREE_ARG) (m+nrl-NR_END)); 
 1049: } 
 1050: 
 1051: /******************* matrix *******************************/
 1052: double **matrix(long nrl, long nrh, long ncl, long nch)
 1053: {
 1054:   long i, nrow=nrh-nrl+1, ncol=nch-ncl+1;
 1055:   double **m;
 1056: 
 1057:   m=(double **) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
 1058:   if (!m) nrerror("allocation failure 1 in matrix()");
 1059:   m += NR_END;
 1060:   m -= nrl;
 1061: 
 1062:   m[nrl]=(double *) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
 1063:   if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
 1064:   m[nrl] += NR_END;
 1065:   m[nrl] -= ncl;
 1066: 
 1067:   for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
 1068:   return m;
 1069:   /* print *(*(m+1)+70) or print m[1][70]; print m+1 or print &(m[1]) or &(m[1][0])
 1070: m[i] = address of ith row of the table. &(m[i]) is its value which is another adress
 1071: that of m[i][0]. In order to get the value p m[i][0] but it is unitialized.
 1072:    */
 1073: }
 1074: 
 1075: /*************************free matrix ************************/
 1076: void free_matrix(double **m, long nrl, long nrh, long ncl, long nch)
 1077: {
 1078:   free((FREE_ARG)(m[nrl]+ncl-NR_END));
 1079:   free((FREE_ARG)(m+nrl-NR_END));
 1080: }
 1081: 
 1082: /******************* ma3x *******************************/
 1083: double ***ma3x(long nrl, long nrh, long ncl, long nch, long nll, long nlh)
 1084: {
 1085:   long i, j, nrow=nrh-nrl+1, ncol=nch-ncl+1, nlay=nlh-nll+1;
 1086:   double ***m;
 1087: 
 1088:   m=(double ***) malloc((size_t)((nrow+NR_END)*sizeof(double*)));
 1089:   if (!m) nrerror("allocation failure 1 in matrix()");
 1090:   m += NR_END;
 1091:   m -= nrl;
 1092: 
 1093:   m[nrl]=(double **) malloc((size_t)((nrow*ncol+NR_END)*sizeof(double)));
 1094:   if (!m[nrl]) nrerror("allocation failure 2 in matrix()");
 1095:   m[nrl] += NR_END;
 1096:   m[nrl] -= ncl;
 1097: 
 1098:   for (i=nrl+1; i<=nrh; i++) m[i]=m[i-1]+ncol;
 1099: 
 1100:   m[nrl][ncl]=(double *) malloc((size_t)((nrow*ncol*nlay+NR_END)*sizeof(double)));
 1101:   if (!m[nrl][ncl]) nrerror("allocation failure 3 in matrix()");
 1102:   m[nrl][ncl] += NR_END;
 1103:   m[nrl][ncl] -= nll;
 1104:   for (j=ncl+1; j<=nch; j++) 
 1105:     m[nrl][j]=m[nrl][j-1]+nlay;
 1106:   
 1107:   for (i=nrl+1; i<=nrh; i++) {
 1108:     m[i][ncl]=m[i-1l][ncl]+ncol*nlay;
 1109:     for (j=ncl+1; j<=nch; j++) 
 1110:       m[i][j]=m[i][j-1]+nlay;
 1111:   }
 1112:   return m; 
 1113:   /*  gdb: p *(m+1) <=> p m[1] and p (m+1) <=> p (m+1) <=> p &(m[1])
 1114:            &(m[i][j][k]) <=> *((*(m+i) + j)+k)
 1115:   */
 1116: }
 1117: 
 1118: /*************************free ma3x ************************/
 1119: void free_ma3x(double ***m, long nrl, long nrh, long ncl, long nch,long nll, long nlh)
 1120: {
 1121:   free((FREE_ARG)(m[nrl][ncl]+ nll-NR_END));
 1122:   free((FREE_ARG)(m[nrl]+ncl-NR_END));
 1123:   free((FREE_ARG)(m+nrl-NR_END));
 1124: }
 1125: 
 1126: /*************** function subdirf ***********/
 1127: char *subdirf(char fileres[])
 1128: {
 1129:   /* Caution optionfilefiname is hidden */
 1130:   strcpy(tmpout,optionfilefiname);
 1131:   strcat(tmpout,"/"); /* Add to the right */
 1132:   strcat(tmpout,fileres);
 1133:   return tmpout;
 1134: }
 1135: 
 1136: /*************** function subdirf2 ***********/
 1137: char *subdirf2(char fileres[], char *preop)
 1138: {
 1139:   
 1140:   /* Caution optionfilefiname is hidden */
 1141:   strcpy(tmpout,optionfilefiname);
 1142:   strcat(tmpout,"/");
 1143:   strcat(tmpout,preop);
 1144:   strcat(tmpout,fileres);
 1145:   return tmpout;
 1146: }
 1147: 
 1148: /*************** function subdirf3 ***********/
 1149: char *subdirf3(char fileres[], char *preop, char *preop2)
 1150: {
 1151:   
 1152:   /* Caution optionfilefiname is hidden */
 1153:   strcpy(tmpout,optionfilefiname);
 1154:   strcat(tmpout,"/");
 1155:   strcat(tmpout,preop);
 1156:   strcat(tmpout,preop2);
 1157:   strcat(tmpout,fileres);
 1158:   return tmpout;
 1159: }
 1160: 
 1161: char *asc_diff_time(long time_sec, char ascdiff[])
 1162: {
 1163:   long sec_left, days, hours, minutes;
 1164:   days = (time_sec) / (60*60*24);
 1165:   sec_left = (time_sec) % (60*60*24);
 1166:   hours = (sec_left) / (60*60) ;
 1167:   sec_left = (sec_left) %(60*60);
 1168:   minutes = (sec_left) /60;
 1169:   sec_left = (sec_left) % (60);
 1170:   sprintf(ascdiff,"%ld day(s) %ld hour(s) %ld minute(s) %ld second(s)",days, hours, minutes, sec_left);  
 1171:   return ascdiff;
 1172: }
 1173: 
 1174: /***************** f1dim *************************/
 1175: extern int ncom; 
 1176: extern double *pcom,*xicom;
 1177: extern double (*nrfunc)(double []); 
 1178:  
 1179: double f1dim(double x) 
 1180: { 
 1181:   int j; 
 1182:   double f;
 1183:   double *xt; 
 1184:  
 1185:   xt=vector(1,ncom); 
 1186:   for (j=1;j<=ncom;j++) xt[j]=pcom[j]+x*xicom[j]; 
 1187:   f=(*nrfunc)(xt); 
 1188:   free_vector(xt,1,ncom); 
 1189:   return f; 
 1190: } 
 1191: 
 1192: /*****************brent *************************/
 1193: double brent(double ax, double bx, double cx, double (*f)(double), double tol, 	double *xmin) 
 1194: { 
 1195:   int iter; 
 1196:   double a,b,d,etemp;
 1197:   double fu=0,fv,fw,fx;
 1198:   double ftemp=0.;
 1199:   double p,q,r,tol1,tol2,u,v,w,x,xm; 
 1200:   double e=0.0; 
 1201:  
 1202:   a=(ax < cx ? ax : cx); 
 1203:   b=(ax > cx ? ax : cx); 
 1204:   x=w=v=bx; 
 1205:   fw=fv=fx=(*f)(x); 
 1206:   for (iter=1;iter<=ITMAX;iter++) { 
 1207:     xm=0.5*(a+b); 
 1208:     tol2=2.0*(tol1=tol*fabs(x)+ZEPS); 
 1209:     /*		if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret)))*/
 1210:     printf(".");fflush(stdout);
 1211:     fprintf(ficlog,".");fflush(ficlog);
 1212: #ifdef DEBUGBRENT
 1213:     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);
 1214:     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);
 1215:     /*		if ((fabs(x-xm) <= (tol2-0.5*(b-a)))||(2.0*fabs(fu-ftemp) <= ftol*1.e-2*(fabs(fu)+fabs(ftemp)))) { */
 1216: #endif
 1217:     if (fabs(x-xm) <= (tol2-0.5*(b-a))){ 
 1218:       *xmin=x; 
 1219:       return fx; 
 1220:     } 
 1221:     ftemp=fu;
 1222:     if (fabs(e) > tol1) { 
 1223:       r=(x-w)*(fx-fv); 
 1224:       q=(x-v)*(fx-fw); 
 1225:       p=(x-v)*q-(x-w)*r; 
 1226:       q=2.0*(q-r); 
 1227:       if (q > 0.0) p = -p; 
 1228:       q=fabs(q); 
 1229:       etemp=e; 
 1230:       e=d; 
 1231:       if (fabs(p) >= fabs(0.5*q*etemp) || p <= q*(a-x) || p >= q*(b-x)) 
 1232: 	d=CGOLD*(e=(x >= xm ? a-x : b-x)); 
 1233:       else { 
 1234: 	d=p/q; 
 1235: 	u=x+d; 
 1236: 	if (u-a < tol2 || b-u < tol2) 
 1237: 	  d=SIGN(tol1,xm-x); 
 1238:       } 
 1239:     } else { 
 1240:       d=CGOLD*(e=(x >= xm ? a-x : b-x)); 
 1241:     } 
 1242:     u=(fabs(d) >= tol1 ? x+d : x+SIGN(tol1,d)); 
 1243:     fu=(*f)(u); 
 1244:     if (fu <= fx) { 
 1245:       if (u >= x) a=x; else b=x; 
 1246:       SHFT(v,w,x,u) 
 1247: 	SHFT(fv,fw,fx,fu) 
 1248: 	} else { 
 1249: 	  if (u < x) a=u; else b=u; 
 1250: 	  if (fu <= fw || w == x) { 
 1251: 	    v=w; 
 1252: 	    w=u; 
 1253: 	    fv=fw; 
 1254: 	    fw=fu; 
 1255: 	  } else if (fu <= fv || v == x || v == w) { 
 1256: 	    v=u; 
 1257: 	    fv=fu; 
 1258: 	  } 
 1259: 	} 
 1260:   } 
 1261:   nrerror("Too many iterations in brent"); 
 1262:   *xmin=x; 
 1263:   return fx; 
 1264: } 
 1265: 
 1266: /****************** mnbrak ***********************/
 1267: 
 1268: void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb, double *fc, 
 1269: 	    double (*func)(double)) 
 1270: { 
 1271:   double ulim,u,r,q, dum;
 1272:   double fu; 
 1273:  
 1274:   *fa=(*func)(*ax); 
 1275:   *fb=(*func)(*bx); 
 1276:   if (*fb > *fa) { 
 1277:     SHFT(dum,*ax,*bx,dum) 
 1278:       SHFT(dum,*fb,*fa,dum) 
 1279:       } 
 1280:   *cx=(*bx)+GOLD*(*bx-*ax); 
 1281:   *fc=(*func)(*cx); 
 1282:   while (*fb > *fc) { /* Declining fa, fb, fc */
 1283:     r=(*bx-*ax)*(*fb-*fc); 
 1284:     q=(*bx-*cx)*(*fb-*fa); 
 1285:     u=(*bx)-((*bx-*cx)*q-(*bx-*ax)*r)/ 
 1286:       (2.0*SIGN(FMAX(fabs(q-r),TINY),q-r)); /* Minimum abscisse of a parabolic estimated from (a,fa), (b,fb) and (c,fc). */
 1287:     ulim=(*bx)+GLIMIT*(*cx-*bx); /* Maximum abscisse where function can be evaluated */
 1288:     if ((*bx-u)*(u-*cx) > 0.0) { /* if u between b and c */
 1289:       fu=(*func)(u); 
 1290: #ifdef DEBUG
 1291:       /* f(x)=A(x-u)**2+f(u) */
 1292:       double A, fparabu; 
 1293:       A= (*fb - *fa)/(*bx-*ax)/(*bx+*ax-2*u);
 1294:       fparabu= *fa - A*(*ax-u)*(*ax-u);
 1295:       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);
 1296:       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);
 1297: #endif 
 1298:     } else if ((*cx-u)*(u-ulim) > 0.0) { /* u is after c but before ulim */
 1299:       fu=(*func)(u); 
 1300:       if (fu < *fc) { 
 1301: 	SHFT(*bx,*cx,u,*cx+GOLD*(*cx-*bx)) 
 1302: 	  SHFT(*fb,*fc,fu,(*func)(u)) 
 1303: 	  } 
 1304:     } else if ((u-ulim)*(ulim-*cx) >= 0.0) { /* u outside ulim (verifying that ulim is beyond c) */
 1305:       u=ulim; 
 1306:       fu=(*func)(u); 
 1307:     } else { 
 1308:       u=(*cx)+GOLD*(*cx-*bx); 
 1309:       fu=(*func)(u); 
 1310:     } 
 1311:     SHFT(*ax,*bx,*cx,u) 
 1312:       SHFT(*fa,*fb,*fc,fu) 
 1313:       } 
 1314: } 
 1315: 
 1316: /*************** linmin ************************/
 1317: /* Given an n -dimensional point p[1..n] and an n -dimensional direction xi[1..n] , moves and
 1318: resets p to where the function func(p) takes on a minimum along the direction xi from p ,
 1319: and replaces xi by the actual vector displacement that p was moved. Also returns as fret
 1320: the value of func at the returned location p . This is actually all accomplished by calling the
 1321: routines mnbrak and brent .*/
 1322: int ncom; 
 1323: double *pcom,*xicom;
 1324: double (*nrfunc)(double []); 
 1325:  
 1326: void linmin(double p[], double xi[], int n, double *fret,double (*func)(double [])) 
 1327: { 
 1328:   double brent(double ax, double bx, double cx, 
 1329: 	       double (*f)(double), double tol, double *xmin); 
 1330:   double f1dim(double x); 
 1331:   void mnbrak(double *ax, double *bx, double *cx, double *fa, double *fb, 
 1332: 	      double *fc, double (*func)(double)); 
 1333:   int j; 
 1334:   double xx,xmin,bx,ax; 
 1335:   double fx,fb,fa;
 1336:  
 1337:   ncom=n; 
 1338:   pcom=vector(1,n); 
 1339:   xicom=vector(1,n); 
 1340:   nrfunc=func; 
 1341:   for (j=1;j<=n;j++) { 
 1342:     pcom[j]=p[j]; 
 1343:     xicom[j]=xi[j]; 
 1344:   } 
 1345:   ax=0.0; 
 1346:   xx=1.0; 
 1347:   mnbrak(&ax,&xx,&bx,&fa,&fx,&fb,f1dim); /* Find a bracket a,x,b in direction n=xi ie xicom */
 1348:   *fret=brent(ax,xx,bx,f1dim,TOL,&xmin); /* Find a minimum P+lambda n in that direction (lambdamin), with TOL between abscisses */
 1349: #ifdef DEBUG
 1350:   printf("retour brent fret=%.12e xmin=%.12e\n",*fret,xmin);
 1351:   fprintf(ficlog,"retour brent fret=%.12e xmin=%.12e\n",*fret,xmin);
 1352: #endif
 1353:   for (j=1;j<=n;j++) { 
 1354:     xi[j] *= xmin; 
 1355:     p[j] += xi[j]; 
 1356:   } 
 1357:   free_vector(xicom,1,n); 
 1358:   free_vector(pcom,1,n); 
 1359: } 
 1360: 
 1361: 
 1362: /*************** powell ************************/
 1363: /*
 1364: Minimization of a function func of n variables. Input consists of an initial starting point
 1365: p[1..n] ; an initial matrix xi[1..n][1..n] , whose columns contain the initial set of di-
 1366: rections (usually the n unit vectors); and ftol , the fractional tolerance in the function value
 1367: such that failure to decrease by more than this amount on one iteration signals doneness. On
 1368: output, p is set to the best point found, xi is the then-current direction set, fret is the returned
 1369: function value at p , and iter is the number of iterations taken. The routine linmin is used.
 1370:  */
 1371: void powell(double p[], double **xi, int n, double ftol, int *iter, double *fret, 
 1372: 	    double (*func)(double [])) 
 1373: { 
 1374:   void linmin(double p[], double xi[], int n, double *fret, 
 1375: 	      double (*func)(double [])); 
 1376:   int i,ibig,j; 
 1377:   double del,t,*pt,*ptt,*xit;
 1378:   double fp,fptt;
 1379:   double *xits;
 1380:   int niterf, itmp;
 1381: 
 1382:   pt=vector(1,n); 
 1383:   ptt=vector(1,n); 
 1384:   xit=vector(1,n); 
 1385:   xits=vector(1,n); 
 1386:   *fret=(*func)(p); 
 1387:   for (j=1;j<=n;j++) pt[j]=p[j]; 
 1388:     rcurr_time = time(NULL);  
 1389:   for (*iter=1;;++(*iter)) { 
 1390:     fp=(*fret); 
 1391:     ibig=0; 
 1392:     del=0.0; 
 1393:     rlast_time=rcurr_time;
 1394:     /* (void) gettimeofday(&curr_time,&tzp); */
 1395:     rcurr_time = time(NULL);  
 1396:     curr_time = *localtime(&rcurr_time);
 1397:     printf("\nPowell iter=%d -2*LL=%.12f %ld sec. %ld sec.",*iter,*fret, rcurr_time-rlast_time, rcurr_time-rstart_time);fflush(stdout);
 1398:     fprintf(ficlog,"\nPowell iter=%d -2*LL=%.12f %ld sec. %ld sec.",*iter,*fret,rcurr_time-rlast_time, rcurr_time-rstart_time); fflush(ficlog);
 1399: /*     fprintf(ficrespow,"%d %.12f %ld",*iter,*fret,curr_time.tm_sec-start_time.tm_sec); */
 1400:    for (i=1;i<=n;i++) {
 1401:       printf(" %d %.12f",i, p[i]);
 1402:       fprintf(ficlog," %d %.12lf",i, p[i]);
 1403:       fprintf(ficrespow," %.12lf", p[i]);
 1404:     }
 1405:     printf("\n");
 1406:     fprintf(ficlog,"\n");
 1407:     fprintf(ficrespow,"\n");fflush(ficrespow);
 1408:     if(*iter <=3){
 1409:       tml = *localtime(&rcurr_time);
 1410:       strcpy(strcurr,asctime(&tml));
 1411:       rforecast_time=rcurr_time; 
 1412:       itmp = strlen(strcurr);
 1413:       if(strcurr[itmp-1]=='\n')  /* Windows outputs with a new line */
 1414: 	strcurr[itmp-1]='\0';
 1415:       printf("\nConsidering the time needed for the last iteration #%d: %ld seconds,\n",*iter,rcurr_time-rlast_time);
 1416:       fprintf(ficlog,"\nConsidering the time needed for this last iteration #%d: %ld seconds,\n",*iter,rcurr_time-rlast_time);
 1417:       for(niterf=10;niterf<=30;niterf+=10){
 1418: 	rforecast_time=rcurr_time+(niterf-*iter)*(rcurr_time-rlast_time);
 1419: 	forecast_time = *localtime(&rforecast_time);
 1420: 	strcpy(strfor,asctime(&forecast_time));
 1421: 	itmp = strlen(strfor);
 1422: 	if(strfor[itmp-1]=='\n')
 1423: 	strfor[itmp-1]='\0';
 1424: 	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);
 1425: 	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);
 1426:       }
 1427:     }
 1428:     for (i=1;i<=n;i++) { 
 1429:       for (j=1;j<=n;j++) xit[j]=xi[j][i]; 
 1430:       fptt=(*fret); 
 1431: #ifdef DEBUG
 1432: 	  printf("fret=%lf, %lf, %lf \n", *fret, *fret, *fret);
 1433: 	  fprintf(ficlog, "fret=%lf, %lf, %lf \n", *fret, *fret, *fret);
 1434: #endif
 1435:       printf("%d",i);fflush(stdout);
 1436:       fprintf(ficlog,"%d",i);fflush(ficlog);
 1437:       linmin(p,xit,n,fret,func); 
 1438:       if (fabs(fptt-(*fret)) > del) { 
 1439: 	del=fabs(fptt-(*fret)); 
 1440: 	ibig=i; 
 1441:       } 
 1442: #ifdef DEBUG
 1443:       printf("%d %.12e",i,(*fret));
 1444:       fprintf(ficlog,"%d %.12e",i,(*fret));
 1445:       for (j=1;j<=n;j++) {
 1446: 	xits[j]=FMAX(fabs(p[j]-pt[j]),1.e-5);
 1447: 	printf(" x(%d)=%.12e",j,xit[j]);
 1448: 	fprintf(ficlog," x(%d)=%.12e",j,xit[j]);
 1449:       }
 1450:       for(j=1;j<=n;j++) {
 1451: 	printf(" p(%d)=%.12e",j,p[j]);
 1452: 	fprintf(ficlog," p(%d)=%.12e",j,p[j]);
 1453:       }
 1454:       printf("\n");
 1455:       fprintf(ficlog,"\n");
 1456: #endif
 1457:     } /* end i */
 1458:     if (2.0*fabs(fp-(*fret)) <= ftol*(fabs(fp)+fabs(*fret))) {
 1459: #ifdef DEBUG
 1460:       int k[2],l;
 1461:       k[0]=1;
 1462:       k[1]=-1;
 1463:       printf("Max: %.12e",(*func)(p));
 1464:       fprintf(ficlog,"Max: %.12e",(*func)(p));
 1465:       for (j=1;j<=n;j++) {
 1466: 	printf(" %.12e",p[j]);
 1467: 	fprintf(ficlog," %.12e",p[j]);
 1468:       }
 1469:       printf("\n");
 1470:       fprintf(ficlog,"\n");
 1471:       for(l=0;l<=1;l++) {
 1472: 	for (j=1;j<=n;j++) {
 1473: 	  ptt[j]=p[j]+(p[j]-pt[j])*k[l];
 1474: 	  printf("l=%d j=%d ptt=%.12e, xits=%.12e, p=%.12e, xit=%.12e", l,j,ptt[j],xits[j],p[j],xit[j]);
 1475: 	  fprintf(ficlog,"l=%d j=%d ptt=%.12e, xits=%.12e, p=%.12e, xit=%.12e", l,j,ptt[j],xits[j],p[j],xit[j]);
 1476: 	}
 1477: 	printf("func(ptt)=%.12e, deriv=%.12e\n",(*func)(ptt),(ptt[j]-p[j])/((*func)(ptt)-(*func)(p)));
 1478: 	fprintf(ficlog,"func(ptt)=%.12e, deriv=%.12e\n",(*func)(ptt),(ptt[j]-p[j])/((*func)(ptt)-(*func)(p)));
 1479:       }
 1480: #endif
 1481: 
 1482: 
 1483:       free_vector(xit,1,n); 
 1484:       free_vector(xits,1,n); 
 1485:       free_vector(ptt,1,n); 
 1486:       free_vector(pt,1,n); 
 1487:       return; 
 1488:     } 
 1489:     if (*iter == ITMAX) nrerror("powell exceeding maximum iterations."); 
 1490:     for (j=1;j<=n;j++) { /* Computes an extrapolated point */
 1491:       ptt[j]=2.0*p[j]-pt[j]; 
 1492:       xit[j]=p[j]-pt[j]; 
 1493:       pt[j]=p[j]; 
 1494:     } 
 1495:     fptt=(*func)(ptt); 
 1496:     if (fptt < fp) { /* If extrapolated point is better, decide if we keep that new direction or not */
 1497:       /* (x1 f1=fp), (x2 f2=*fret), (x3 f3=fptt), (xm fm) */
 1498:       /* From x1 (P0) distance of x2 is at h and x3 is 2h */
 1499:       /* Let f"(x2) be the 2nd derivative equal everywhere.  */
 1500:       /* Then the parabolic through (x1,f1), (x2,f2) and (x3,f3) */
 1501:       /* will reach at f3 = fm + h^2/2 f"m  ; f" = (f1 -2f2 +f3 ) / h**2 */
 1502:       /* f1-f3 = delta(2h) = 2 h**2 f'' = 2(f1- 2f2 +f3) */
 1503:       /* Thus we compare delta(2h) with observed f1-f3 */
 1504:       /* or best gain on one ancient line 'del' with total  */
 1505:       /* gain f1-f2 = f1 - f2 - 'del' with del  */
 1506:       /* t=2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del)-del*SQR(fp-fptt); */
 1507: 
 1508:       t=2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del);
 1509:       t= t- del*SQR(fp-fptt);
 1510:       printf("t1= %.12lf, t2= %.12lf, t=%.12lf\n", 2.0*(fp-2.0*(*fret)+fptt)*SQR(fp-(*fret)-del),del*SQR(fp-fptt),t);
 1511:       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);
 1512: #ifdef DEBUG
 1513:       printf("t3= %.12lf, t4= %.12lf, t3*= %.12lf, t4*= %.12lf\n",SQR(fp-(*fret)-del),SQR(fp-fptt),
 1514: 	     (fp-(*fret)-del)*(fp-(*fret)-del),(fp-fptt)*(fp-fptt));
 1515:       fprintf(ficlog,"t3= %.12lf, t4= %.12lf, t3*= %.12lf, t4*= %.12lf\n",SQR(fp-(*fret)-del),SQR(fp-fptt),
 1516: 	     (fp-(*fret)-del)*(fp-(*fret)-del),(fp-fptt)*(fp-fptt));
 1517:       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);
 1518:       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);
 1519: #endif
 1520:       if (t < 0.0) { /* Then we use it for last direction */
 1521: 	linmin(p,xit,n,fret,func); /* computes mean on the extrapolated direction.*/
 1522: 	for (j=1;j<=n;j++) { 
 1523: 	  xi[j][ibig]=xi[j][n]; /* Replace the direction with biggest decrease by n */
 1524: 	  xi[j][n]=xit[j];      /* and nth direction by the extrapolated */
 1525: 	}
 1526: 	printf("Gaining to use average direction of P0 P%d instead of biggest increase direction %d :\n",n,ibig);
 1527: 	fprintf(ficlog,"Gaining to use average direction of P0 P%d instead of biggest increase direction %d :\n",n,ibig);
 1528: 
 1529: #ifdef DEBUG
 1530: 	printf("Direction changed  last moved %d in place of ibig=%d, new last is the average:\n",n,ibig);
 1531: 	fprintf(ficlog,"Direction changed  last moved %d in place of ibig=%d, new last is the average:\n",n,ibig);
 1532: 	for(j=1;j<=n;j++){
 1533: 	  printf(" %.12e",xit[j]);
 1534: 	  fprintf(ficlog," %.12e",xit[j]);
 1535: 	}
 1536: 	printf("\n");
 1537: 	fprintf(ficlog,"\n");
 1538: #endif
 1539:       } /* end of t negative */
 1540:     } /* end if (fptt < fp)  */
 1541:   } 
 1542: } 
 1543: 
 1544: /**** Prevalence limit (stable or period prevalence)  ****************/
 1545: 
 1546: double **prevalim(double **prlim, int nlstate, double x[], double age, double **oldm, double **savm, double ftolpl, int ij)
 1547: {
 1548:   /* Computes the prevalence limit in each live state at age x by left multiplying the unit
 1549:      matrix by transitions matrix until convergence is reached */
 1550:   
 1551:   int i, ii,j,k;
 1552:   double min, max, maxmin, maxmax,sumnew=0.;
 1553:   /* double **matprod2(); */ /* test */
 1554:   double **out, cov[NCOVMAX+1], **pmij();
 1555:   double **newm;
 1556:   double agefin, delaymax=50 ; /* Max number of years to converge */
 1557:   
 1558:   for (ii=1;ii<=nlstate+ndeath;ii++)
 1559:     for (j=1;j<=nlstate+ndeath;j++){
 1560:       oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1561:     }
 1562:   
 1563:   cov[1]=1.;
 1564:   
 1565:   /* Even if hstepm = 1, at least one multiplication by the unit matrix */
 1566:   for(agefin=age-stepm/YEARM; agefin>=age-delaymax; agefin=agefin-stepm/YEARM){
 1567:     newm=savm;
 1568:     /* Covariates have to be included here again */
 1569:     cov[2]=agefin;
 1570:     
 1571:     for (k=1; k<=cptcovn;k++) {
 1572:       cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]];
 1573:       /*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]]);*/
 1574:     }
 1575:     /* for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2]; */
 1576:     /* for (k=1; k<=cptcovprod;k++) /\* Useless *\/ */
 1577:     /*   cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]] * nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]]; */
 1578:     
 1579:     /*printf("ij=%d cptcovprod=%d tvar=%d ", ij, cptcovprod, Tvar[1]);*/
 1580:     /*printf("ij=%d cov[3]=%lf cov[4]=%lf \n",ij, cov[3],cov[4]);*/
 1581:     /*printf("ij=%d cov[3]=%lf \n",ij, cov[3]);*/
 1582:     /* savm=pmij(pmmij,cov,ncovmodel,x,nlstate); */
 1583:     /* out=matprod2(newm, pmij(pmmij,cov,ncovmodel,x,nlstate),1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, oldm); /\* Bug Valgrind *\/ */
 1584:     out=matprod2(newm, pmij(pmmij,cov,ncovmodel,x,nlstate),1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, oldm); /* Bug Valgrind */
 1585:     
 1586:     savm=oldm;
 1587:     oldm=newm;
 1588:     maxmax=0.;
 1589:     for(j=1;j<=nlstate;j++){
 1590:       min=1.;
 1591:       max=0.;
 1592:       for(i=1; i<=nlstate; i++) {
 1593: 	sumnew=0;
 1594: 	for(k=1; k<=ndeath; k++) sumnew+=newm[i][nlstate+k];
 1595: 	prlim[i][j]= newm[i][j]/(1-sumnew);
 1596:         /*printf(" prevalim i=%d, j=%d, prmlim[%d][%d]=%f, agefin=%d \n", i, j, i, j, prlim[i][j],(int)agefin);*/
 1597: 	max=FMAX(max,prlim[i][j]);
 1598: 	min=FMIN(min,prlim[i][j]);
 1599:       }
 1600:       maxmin=max-min;
 1601:       maxmax=FMAX(maxmax,maxmin);
 1602:     } /* j loop */
 1603:     if(maxmax < ftolpl){
 1604:       return prlim;
 1605:     }
 1606:   } /* age loop */
 1607:   return prlim; /* should not reach here */
 1608: }
 1609: 
 1610: /*************** transition probabilities ***************/ 
 1611: 
 1612: double **pmij(double **ps, double *cov, int ncovmodel, double *x, int nlstate )
 1613: {
 1614:   /* According to parameters values stored in x and the covariate's values stored in cov,
 1615:      computes the probability to be observed in state j being in state i by appying the
 1616:      model to the ncovmodel covariates (including constant and age).
 1617:      lnpijopii=ln(pij/pii)= aij+bij*age+cij*v1+dij*v2+... = sum_nc=1^ncovmodel xij(nc)*cov[nc]
 1618:      and, according on how parameters are entered, the position of the coefficient xij(nc) of the
 1619:      ncth covariate in the global vector x is given by the formula:
 1620:      j<i nc+((i-1)*(nlstate+ndeath-1)+j-1)*ncovmodel
 1621:      j>=i nc + ((i-1)*(nlstate+ndeath-1)+(j-2))*ncovmodel
 1622:      Computes ln(pij/pii) (lnpijopii), deduces pij/pii by exponentiation,
 1623:      sums on j different of i to get 1-pii/pii, deduces pii, and then all pij.
 1624:      Outputs ps[i][j] the probability to be observed in j being in j according to
 1625:      the values of the covariates cov[nc] and corresponding parameter values x[nc+shiftij]
 1626:   */
 1627:   double s1, lnpijopii;
 1628:   /*double t34;*/
 1629:   int i,j, nc, ii, jj;
 1630: 
 1631:     for(i=1; i<= nlstate; i++){
 1632:       for(j=1; j<i;j++){
 1633: 	for (nc=1, lnpijopii=0.;nc <=ncovmodel; nc++){
 1634: 	  /*lnpijopii += param[i][j][nc]*cov[nc];*/
 1635: 	  lnpijopii += x[nc+((i-1)*(nlstate+ndeath-1)+j-1)*ncovmodel]*cov[nc];
 1636: /* 	 printf("Int j<i s1=%.17e, lnpijopii=%.17e\n",s1,lnpijopii); */
 1637: 	}
 1638: 	ps[i][j]=lnpijopii; /* In fact ln(pij/pii) */
 1639: /* 	printf("s1=%.17e, lnpijopii=%.17e\n",s1,lnpijopii); */
 1640:       }
 1641:       for(j=i+1; j<=nlstate+ndeath;j++){
 1642: 	for (nc=1, lnpijopii=0.;nc <=ncovmodel; nc++){
 1643: 	  /*lnpijopii += x[(i-1)*nlstate*ncovmodel+(j-2)*ncovmodel+nc+(i-1)*(ndeath-1)*ncovmodel]*cov[nc];*/
 1644: 	  lnpijopii += x[nc + ((i-1)*(nlstate+ndeath-1)+(j-2))*ncovmodel]*cov[nc];
 1645: /* 	  printf("Int j>i s1=%.17e, lnpijopii=%.17e %lx %lx\n",s1,lnpijopii,s1,lnpijopii); */
 1646: 	}
 1647: 	ps[i][j]=lnpijopii; /* In fact ln(pij/pii) */
 1648:       }
 1649:     }
 1650:     
 1651:     for(i=1; i<= nlstate; i++){
 1652:       s1=0;
 1653:       for(j=1; j<i; j++){
 1654: 	s1+=exp(ps[i][j]); /* In fact sums pij/pii */
 1655: 	/*printf("debug1 %d %d ps=%lf exp(ps)=%lf s1+=%lf\n",i,j,ps[i][j],exp(ps[i][j]),s1); */
 1656:       }
 1657:       for(j=i+1; j<=nlstate+ndeath; j++){
 1658: 	s1+=exp(ps[i][j]); /* In fact sums pij/pii */
 1659: 	/*printf("debug2 %d %d ps=%lf exp(ps)=%lf s1+=%lf\n",i,j,ps[i][j],exp(ps[i][j]),s1); */
 1660:       }
 1661:       /* s1= sum_{j<>i} pij/pii=(1-pii)/pii and thus pii is known from s1 */
 1662:       ps[i][i]=1./(s1+1.);
 1663:       /* Computing other pijs */
 1664:       for(j=1; j<i; j++)
 1665: 	ps[i][j]= exp(ps[i][j])*ps[i][i];
 1666:       for(j=i+1; j<=nlstate+ndeath; j++)
 1667: 	ps[i][j]= exp(ps[i][j])*ps[i][i];
 1668:       /* ps[i][nlstate+1]=1.-s1- ps[i][i];*/ /* Sum should be 1 */
 1669:     } /* end i */
 1670:     
 1671:     for(ii=nlstate+1; ii<= nlstate+ndeath; ii++){
 1672:       for(jj=1; jj<= nlstate+ndeath; jj++){
 1673: 	ps[ii][jj]=0;
 1674: 	ps[ii][ii]=1;
 1675:       }
 1676:     }
 1677:     
 1678:     
 1679:     /* for(ii=1; ii<= nlstate+ndeath; ii++){ */
 1680:     /*   for(jj=1; jj<= nlstate+ndeath; jj++){ */
 1681:     /* 	printf(" pmij  ps[%d][%d]=%lf ",ii,jj,ps[ii][jj]); */
 1682:     /*   } */
 1683:     /*   printf("\n "); */
 1684:     /* } */
 1685:     /* printf("\n ");printf("%lf ",cov[2]);*/
 1686:     /*
 1687:       for(i=1; i<= npar; i++) printf("%f ",x[i]);
 1688:       goto end;*/
 1689:     return ps;
 1690: }
 1691: 
 1692: /**************** Product of 2 matrices ******************/
 1693: 
 1694: double **matprod2(double **out, double **in,int nrl, int nrh, int ncl, int nch, int ncolol, int ncoloh, double **b)
 1695: {
 1696:   /* Computes the matrix product of in(1,nrh-nrl+1)(1,nch-ncl+1) times
 1697:      b(1,nch-ncl+1)(1,ncoloh-ncolol+1) into out(...) */
 1698:   /* in, b, out are matrice of pointers which should have been initialized 
 1699:      before: only the contents of out is modified. The function returns
 1700:      a pointer to pointers identical to out */
 1701:   int i, j, k;
 1702:   for(i=nrl; i<= nrh; i++)
 1703:     for(k=ncolol; k<=ncoloh; k++){
 1704:       out[i][k]=0.;
 1705:       for(j=ncl; j<=nch; j++)
 1706:   	out[i][k] +=in[i][j]*b[j][k];
 1707:     }
 1708:   return out;
 1709: }
 1710: 
 1711: 
 1712: /************* Higher Matrix Product ***************/
 1713: 
 1714: double ***hpxij(double ***po, int nhstepm, double age, int hstepm, double *x, int nlstate, int stepm, double **oldm, double **savm, int ij )
 1715: {
 1716:   /* Computes the transition matrix starting at age 'age' over 
 1717:      'nhstepm*hstepm*stepm' months (i.e. until
 1718:      age (in years)  age+nhstepm*hstepm*stepm/12) by multiplying 
 1719:      nhstepm*hstepm matrices. 
 1720:      Output is stored in matrix po[i][j][h] for h every 'hstepm' step 
 1721:      (typically every 2 years instead of every month which is too big 
 1722:      for the memory).
 1723:      Model is determined by parameters x and covariates have to be 
 1724:      included manually here. 
 1725: 
 1726:      */
 1727: 
 1728:   int i, j, d, h, k;
 1729:   double **out, cov[NCOVMAX+1];
 1730:   double **newm;
 1731: 
 1732:   /* Hstepm could be zero and should return the unit matrix */
 1733:   for (i=1;i<=nlstate+ndeath;i++)
 1734:     for (j=1;j<=nlstate+ndeath;j++){
 1735:       oldm[i][j]=(i==j ? 1.0 : 0.0);
 1736:       po[i][j][0]=(i==j ? 1.0 : 0.0);
 1737:     }
 1738:   /* Even if hstepm = 1, at least one multiplication by the unit matrix */
 1739:   for(h=1; h <=nhstepm; h++){
 1740:     for(d=1; d <=hstepm; d++){
 1741:       newm=savm;
 1742:       /* Covariates have to be included here again */
 1743:       cov[1]=1.;
 1744:       cov[2]=age+((h-1)*hstepm + (d-1))*stepm/YEARM;
 1745:       for (k=1; k<=cptcovn;k++) 
 1746: 	cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]];
 1747:       for (k=1; k<=cptcovage;k++)
 1748: 	cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
 1749:       for (k=1; k<=cptcovprod;k++) /* Useless because included in cptcovn */
 1750: 	cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
 1751: 
 1752: 
 1753:       /*printf("hxi cptcov=%d cptcode=%d\n",cptcov,cptcode);*/
 1754:       /*printf("h=%d d=%d age=%f cov=%f\n",h,d,age,cov[2]);*/
 1755:       out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,1,nlstate+ndeath, 
 1756: 		   pmij(pmmij,cov,ncovmodel,x,nlstate));
 1757:       savm=oldm;
 1758:       oldm=newm;
 1759:     }
 1760:     for(i=1; i<=nlstate+ndeath; i++)
 1761:       for(j=1;j<=nlstate+ndeath;j++) {
 1762: 	po[i][j][h]=newm[i][j];
 1763: 	/*if(h==nhstepm) printf("po[%d][%d][%d]=%f ",i,j,h,po[i][j][h]);*/
 1764:       }
 1765:     /*printf("h=%d ",h);*/
 1766:   } /* end h */
 1767: /*     printf("\n H=%d \n",h); */
 1768:   return po;
 1769: }
 1770: 
 1771: #ifdef NLOPT
 1772:   double  myfunc(unsigned n, const double *p1, double *grad, void *pd){
 1773:   double fret;
 1774:   double *xt;
 1775:   int j;
 1776:   myfunc_data *d2 = (myfunc_data *) pd;
 1777: /* xt = (p1-1); */
 1778:   xt=vector(1,n); 
 1779:   for (j=1;j<=n;j++)   xt[j]=p1[j-1]; /* xt[1]=p1[0] */
 1780: 
 1781:   fret=(d2->function)(xt); /*  p xt[1]@8 is fine */
 1782:   /* fret=(*func)(xt); /\*  p xt[1]@8 is fine *\/ */
 1783:   printf("Function = %.12lf ",fret);
 1784:   for (j=1;j<=n;j++) printf(" %d %.8lf", j, xt[j]); 
 1785:   printf("\n");
 1786:  free_vector(xt,1,n);
 1787:   return fret;
 1788: }
 1789: #endif
 1790: 
 1791: /*************** log-likelihood *************/
 1792: double func( double *x)
 1793: {
 1794:   int i, ii, j, k, mi, d, kk;
 1795:   double l, ll[NLSTATEMAX+1], cov[NCOVMAX+1];
 1796:   double **out;
 1797:   double sw; /* Sum of weights */
 1798:   double lli; /* Individual log likelihood */
 1799:   int s1, s2;
 1800:   double bbh, survp;
 1801:   long ipmx;
 1802:   /*extern weight */
 1803:   /* We are differentiating ll according to initial status */
 1804:   /*  for (i=1;i<=npar;i++) printf("%f ", x[i]);*/
 1805:   /*for(i=1;i<imx;i++) 
 1806:     printf(" %d\n",s[4][i]);
 1807:   */
 1808: 
 1809:   ++countcallfunc;
 1810: 
 1811:   cov[1]=1.;
 1812: 
 1813:   for(k=1; k<=nlstate; k++) ll[k]=0.;
 1814: 
 1815:   if(mle==1){
 1816:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1817:       /* Computes the values of the ncovmodel covariates of the model
 1818: 	 depending if the covariates are fixed or variying (age dependent) and stores them in cov[]
 1819: 	 Then computes with function pmij which return a matrix p[i][j] giving the elementary probability
 1820: 	 to be observed in j being in i according to the model.
 1821:        */
 1822:       for (k=1; k<=cptcovn;k++){ /* Simple and product covariates without age* products */
 1823: 	cov[2+k]=covar[Tvar[k]][i];
 1824:       }
 1825:       /* In model V2+V1*V4+age*V3+V3*V2 Tvar[1] is V2, Tvar[2=V1*V4] 
 1826: 	 is 6, Tvar[3=age*V3] should not be computed because of age Tvar[4=V3*V2] 
 1827: 	 has been calculated etc */
 1828:       for(mi=1; mi<= wav[i]-1; mi++){
 1829: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1830: 	  for (j=1;j<=nlstate+ndeath;j++){
 1831: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1832: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1833: 	  }
 1834: 	for(d=0; d<dh[mi][i]; d++){
 1835: 	  newm=savm;
 1836: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1837: 	  for (kk=1; kk<=cptcovage;kk++) {
 1838: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2]; /* Tage[kk] gives the data-covariate associated with age */
 1839: 	  }
 1840: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1841: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1842: 	  savm=oldm;
 1843: 	  oldm=newm;
 1844: 	} /* end mult */
 1845:       
 1846: 	/*lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]);*/ /* Original formula */
 1847: 	/* But now since version 0.9 we anticipate for bias at large stepm.
 1848: 	 * If stepm is larger than one month (smallest stepm) and if the exact delay 
 1849: 	 * (in months) between two waves is not a multiple of stepm, we rounded to 
 1850: 	 * the nearest (and in case of equal distance, to the lowest) interval but now
 1851: 	 * we keep into memory the bias bh[mi][i] and also the previous matrix product
 1852: 	 * (i.e to dh[mi][i]-1) saved in 'savm'. Then we inter(extra)polate the
 1853: 	 * probability in order to take into account the bias as a fraction of the way
 1854: 	 * from savm to out if bh is negative or even beyond if bh is positive. bh varies
 1855: 	 * -stepm/2 to stepm/2 .
 1856: 	 * For stepm=1 the results are the same as for previous versions of Imach.
 1857: 	 * For stepm > 1 the results are less biased than in previous versions. 
 1858: 	 */
 1859: 	s1=s[mw[mi][i]][i];
 1860: 	s2=s[mw[mi+1][i]][i];
 1861: 	bbh=(double)bh[mi][i]/(double)stepm; 
 1862: 	/* bias bh is positive if real duration
 1863: 	 * is higher than the multiple of stepm and negative otherwise.
 1864: 	 */
 1865: 	/* lli= (savm[s1][s2]>1.e-8 ?(1.+bbh)*log(out[s1][s2])- bbh*log(savm[s1][s2]):log((1.+bbh)*out[s1][s2]));*/
 1866: 	if( s2 > nlstate){ 
 1867: 	  /* i.e. if s2 is a death state and if the date of death is known 
 1868: 	     then the contribution to the likelihood is the probability to 
 1869: 	     die between last step unit time and current  step unit time, 
 1870: 	     which is also equal to probability to die before dh 
 1871: 	     minus probability to die before dh-stepm . 
 1872: 	     In version up to 0.92 likelihood was computed
 1873: 	as if date of death was unknown. Death was treated as any other
 1874: 	health state: the date of the interview describes the actual state
 1875: 	and not the date of a change in health state. The former idea was
 1876: 	to consider that at each interview the state was recorded
 1877: 	(healthy, disable or death) and IMaCh was corrected; but when we
 1878: 	introduced the exact date of death then we should have modified
 1879: 	the contribution of an exact death to the likelihood. This new
 1880: 	contribution is smaller and very dependent of the step unit
 1881: 	stepm. It is no more the probability to die between last interview
 1882: 	and month of death but the probability to survive from last
 1883: 	interview up to one month before death multiplied by the
 1884: 	probability to die within a month. Thanks to Chris
 1885: 	Jackson for correcting this bug.  Former versions increased
 1886: 	mortality artificially. The bad side is that we add another loop
 1887: 	which slows down the processing. The difference can be up to 10%
 1888: 	lower mortality.
 1889: 	  */
 1890: 	  lli=log(out[s1][s2] - savm[s1][s2]);
 1891: 
 1892: 
 1893: 	} else if  (s2==-2) {
 1894: 	  for (j=1,survp=0. ; j<=nlstate; j++) 
 1895: 	    survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 1896: 	  /*survp += out[s1][j]; */
 1897: 	  lli= log(survp);
 1898: 	}
 1899: 	
 1900:  	else if  (s2==-4) { 
 1901: 	  for (j=3,survp=0. ; j<=nlstate; j++)  
 1902: 	    survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 1903:  	  lli= log(survp); 
 1904:  	} 
 1905: 
 1906:  	else if  (s2==-5) { 
 1907:  	  for (j=1,survp=0. ; j<=2; j++)  
 1908: 	    survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 1909:  	  lli= log(survp); 
 1910:  	} 
 1911: 	
 1912: 	else{
 1913: 	  lli= log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]); /* linear interpolation */
 1914: 	  /*  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 */
 1915: 	} 
 1916: 	/*lli=(1.+bbh)*log(out[s1][s2])- bbh*log(savm[s1][s2]);*/
 1917: 	/*if(lli ==000.0)*/
 1918: 	/*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); */
 1919:   	ipmx +=1;
 1920: 	sw += weight[i];
 1921: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1922:       } /* end of wave */
 1923:     } /* end of individual */
 1924:   }  else if(mle==2){
 1925:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1926:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1927:       for(mi=1; mi<= wav[i]-1; mi++){
 1928: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1929: 	  for (j=1;j<=nlstate+ndeath;j++){
 1930: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1931: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1932: 	  }
 1933: 	for(d=0; d<=dh[mi][i]; d++){
 1934: 	  newm=savm;
 1935: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1936: 	  for (kk=1; kk<=cptcovage;kk++) {
 1937: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1938: 	  }
 1939: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1940: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1941: 	  savm=oldm;
 1942: 	  oldm=newm;
 1943: 	} /* end mult */
 1944:       
 1945: 	s1=s[mw[mi][i]][i];
 1946: 	s2=s[mw[mi+1][i]][i];
 1947: 	bbh=(double)bh[mi][i]/(double)stepm; 
 1948: 	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 */
 1949: 	ipmx +=1;
 1950: 	sw += weight[i];
 1951: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1952:       } /* end of wave */
 1953:     } /* end of individual */
 1954:   }  else if(mle==3){  /* exponential inter-extrapolation */
 1955:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1956:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1957:       for(mi=1; mi<= wav[i]-1; mi++){
 1958: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1959: 	  for (j=1;j<=nlstate+ndeath;j++){
 1960: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1961: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1962: 	  }
 1963: 	for(d=0; d<dh[mi][i]; d++){
 1964: 	  newm=savm;
 1965: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1966: 	  for (kk=1; kk<=cptcovage;kk++) {
 1967: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1968: 	  }
 1969: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 1970: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 1971: 	  savm=oldm;
 1972: 	  oldm=newm;
 1973: 	} /* end mult */
 1974:       
 1975: 	s1=s[mw[mi][i]][i];
 1976: 	s2=s[mw[mi+1][i]][i];
 1977: 	bbh=(double)bh[mi][i]/(double)stepm; 
 1978: 	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 */
 1979: 	ipmx +=1;
 1980: 	sw += weight[i];
 1981: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 1982:       } /* end of wave */
 1983:     } /* end of individual */
 1984:   }else if (mle==4){  /* ml=4 no inter-extrapolation */
 1985:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 1986:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 1987:       for(mi=1; mi<= wav[i]-1; mi++){
 1988: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 1989: 	  for (j=1;j<=nlstate+ndeath;j++){
 1990: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 1991: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 1992: 	  }
 1993: 	for(d=0; d<dh[mi][i]; d++){
 1994: 	  newm=savm;
 1995: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 1996: 	  for (kk=1; kk<=cptcovage;kk++) {
 1997: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 1998: 	  }
 1999: 	
 2000: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 2001: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 2002: 	  savm=oldm;
 2003: 	  oldm=newm;
 2004: 	} /* end mult */
 2005:       
 2006: 	s1=s[mw[mi][i]][i];
 2007: 	s2=s[mw[mi+1][i]][i];
 2008: 	if( s2 > nlstate){ 
 2009: 	  lli=log(out[s1][s2] - savm[s1][s2]);
 2010: 	}else{
 2011: 	  lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]); /* Original formula */
 2012: 	}
 2013: 	ipmx +=1;
 2014: 	sw += weight[i];
 2015: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 2016: /* 	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]); */
 2017:       } /* end of wave */
 2018:     } /* end of individual */
 2019:   }else{  /* ml=5 no inter-extrapolation no jackson =0.8a */
 2020:     for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 2021:       for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 2022:       for(mi=1; mi<= wav[i]-1; mi++){
 2023: 	for (ii=1;ii<=nlstate+ndeath;ii++)
 2024: 	  for (j=1;j<=nlstate+ndeath;j++){
 2025: 	    oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 2026: 	    savm[ii][j]=(ii==j ? 1.0 : 0.0);
 2027: 	  }
 2028: 	for(d=0; d<dh[mi][i]; d++){
 2029: 	  newm=savm;
 2030: 	  cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 2031: 	  for (kk=1; kk<=cptcovage;kk++) {
 2032: 	    cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 2033: 	  }
 2034: 	
 2035: 	  out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 2036: 		       1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 2037: 	  savm=oldm;
 2038: 	  oldm=newm;
 2039: 	} /* end mult */
 2040:       
 2041: 	s1=s[mw[mi][i]][i];
 2042: 	s2=s[mw[mi+1][i]][i];
 2043: 	lli=log(out[s[mw[mi][i]][i]][s[mw[mi+1][i]][i]]); /* Original formula */
 2044: 	ipmx +=1;
 2045: 	sw += weight[i];
 2046: 	ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 2047: 	/*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]);*/
 2048:       } /* end of wave */
 2049:     } /* end of individual */
 2050:   } /* End of if */
 2051:   for(k=1,l=0.; k<=nlstate; k++) l += ll[k];
 2052:   /* printf("l1=%f l2=%f ",ll[1],ll[2]); */
 2053:   l= l*ipmx/sw; /* To get the same order of magnitude as if weight=1 for every body */
 2054:   return -l;
 2055: }
 2056: 
 2057: /*************** log-likelihood *************/
 2058: double funcone( double *x)
 2059: {
 2060:   /* Same as likeli but slower because of a lot of printf and if */
 2061:   int i, ii, j, k, mi, d, kk;
 2062:   double l, ll[NLSTATEMAX+1], cov[NCOVMAX+1];
 2063:   double **out;
 2064:   double lli; /* Individual log likelihood */
 2065:   double llt;
 2066:   int s1, s2;
 2067:   double bbh, survp;
 2068:   /*extern weight */
 2069:   /* We are differentiating ll according to initial status */
 2070:   /*  for (i=1;i<=npar;i++) printf("%f ", x[i]);*/
 2071:   /*for(i=1;i<imx;i++) 
 2072:     printf(" %d\n",s[4][i]);
 2073:   */
 2074:   cov[1]=1.;
 2075: 
 2076:   for(k=1; k<=nlstate; k++) ll[k]=0.;
 2077: 
 2078:   for (i=1,ipmx=0, sw=0.; i<=imx; i++){
 2079:     for (k=1; k<=cptcovn;k++) cov[2+k]=covar[Tvar[k]][i];
 2080:     for(mi=1; mi<= wav[i]-1; mi++){
 2081:       for (ii=1;ii<=nlstate+ndeath;ii++)
 2082: 	for (j=1;j<=nlstate+ndeath;j++){
 2083: 	  oldm[ii][j]=(ii==j ? 1.0 : 0.0);
 2084: 	  savm[ii][j]=(ii==j ? 1.0 : 0.0);
 2085: 	}
 2086:       for(d=0; d<dh[mi][i]; d++){
 2087: 	newm=savm;
 2088: 	cov[2]=agev[mw[mi][i]][i]+d*stepm/YEARM;
 2089: 	for (kk=1; kk<=cptcovage;kk++) {
 2090: 	  cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 2091: 	}
 2092: 	/* savm=pmij(pmmij,cov,ncovmodel,x,nlstate); */
 2093: 	out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath,
 2094: 		     1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate));
 2095: 	/* out=matprod2(newm,oldm,1,nlstate+ndeath,1,nlstate+ndeath, */
 2096: 	/* 	     1,nlstate+ndeath,pmij(pmmij,cov,ncovmodel,x,nlstate)); */
 2097: 	savm=oldm;
 2098: 	oldm=newm;
 2099:       } /* end mult */
 2100:       
 2101:       s1=s[mw[mi][i]][i];
 2102:       s2=s[mw[mi+1][i]][i];
 2103:       bbh=(double)bh[mi][i]/(double)stepm; 
 2104:       /* bias is positive if real duration
 2105:        * is higher than the multiple of stepm and negative otherwise.
 2106:        */
 2107:       if( s2 > nlstate && (mle <5) ){  /* Jackson */
 2108: 	lli=log(out[s1][s2] - savm[s1][s2]);
 2109:       } else if  (s2==-2) {
 2110: 	for (j=1,survp=0. ; j<=nlstate; j++) 
 2111: 	  survp += (1.+bbh)*out[s1][j]- bbh*savm[s1][j];
 2112: 	lli= log(survp);
 2113:       }else if (mle==1){
 2114: 	lli= log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]); /* linear interpolation */
 2115:       } else if(mle==2){
 2116: 	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 */
 2117:       } else if(mle==3){  /* exponential inter-extrapolation */
 2118: 	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 */
 2119:       } else if (mle==4){  /* mle=4 no inter-extrapolation */
 2120: 	lli=log(out[s1][s2]); /* Original formula */
 2121:       } else{  /* mle=0 back to 1 */
 2122: 	lli= log((1.+bbh)*out[s1][s2]- bbh*savm[s1][s2]); /* linear interpolation */
 2123: 	/*lli=log(out[s1][s2]); */ /* Original formula */
 2124:       } /* End of if */
 2125:       ipmx +=1;
 2126:       sw += weight[i];
 2127:       ll[s[mw[mi][i]][i]] += 2*weight[i]*lli;
 2128:       /*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]); */
 2129:       if(globpr){
 2130: 	fprintf(ficresilk,"%9ld %6d %2d %2d %1d %1d %3d %11.6f %8.4f\
 2131:  %11.6f %11.6f %11.6f ", \
 2132: 		num[i],i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],
 2133: 		2*weight[i]*lli,out[s1][s2],savm[s1][s2]);
 2134: 	for(k=1,llt=0.,l=0.; k<=nlstate; k++){
 2135: 	  llt +=ll[k]*gipmx/gsw;
 2136: 	  fprintf(ficresilk," %10.6f",-ll[k]*gipmx/gsw);
 2137: 	}
 2138: 	fprintf(ficresilk," %10.6f\n", -llt);
 2139:       }
 2140:     } /* end of wave */
 2141:   } /* end of individual */
 2142:   for(k=1,l=0.; k<=nlstate; k++) l += ll[k];
 2143:   /* printf("l1=%f l2=%f ",ll[1],ll[2]); */
 2144:   l= l*ipmx/sw; /* To get the same order of magnitude as if weight=1 for every body */
 2145:   if(globpr==0){ /* First time we count the contributions and weights */
 2146:     gipmx=ipmx;
 2147:     gsw=sw;
 2148:   }
 2149:   return -l;
 2150: }
 2151: 
 2152: 
 2153: /*************** function likelione ***********/
 2154: void likelione(FILE *ficres,double p[], int npar, int nlstate, int *globpri, long *ipmx, double *sw, double *fretone, double (*funcone)(double []))
 2155: {
 2156:   /* This routine should help understanding what is done with 
 2157:      the selection of individuals/waves and
 2158:      to check the exact contribution to the likelihood.
 2159:      Plotting could be done.
 2160:    */
 2161:   int k;
 2162: 
 2163:   if(*globpri !=0){ /* Just counts and sums, no printings */
 2164:     strcpy(fileresilk,"ilk"); 
 2165:     strcat(fileresilk,fileres);
 2166:     if((ficresilk=fopen(fileresilk,"w"))==NULL) {
 2167:       printf("Problem with resultfile: %s\n", fileresilk);
 2168:       fprintf(ficlog,"Problem with resultfile: %s\n", fileresilk);
 2169:     }
 2170:     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");
 2171:     fprintf(ficresilk, "#num_i i s1 s2 mi mw dh likeli weight 2wlli out sav ");
 2172:     /* 	i,s1,s2,mi,mw[mi][i],dh[mi][i],exp(lli),weight[i],2*weight[i]*lli,out[s1][s2],savm[s1][s2]); */
 2173:     for(k=1; k<=nlstate; k++) 
 2174:       fprintf(ficresilk," -2*gipw/gsw*weight*ll[%d]++",k);
 2175:     fprintf(ficresilk," -2*gipw/gsw*weight*ll(total)\n");
 2176:   }
 2177: 
 2178:   *fretone=(*funcone)(p);
 2179:   if(*globpri !=0){
 2180:     fclose(ficresilk);
 2181:     fprintf(fichtm,"\n<br>File of contributions to the likelihood: <a href=\"%s\">%s</a><br>\n",subdirf(fileresilk),subdirf(fileresilk));
 2182:     fflush(fichtm); 
 2183:   } 
 2184:   return;
 2185: }
 2186: 
 2187: 
 2188: /*********** Maximum Likelihood Estimation ***************/
 2189: 
 2190: void mlikeli(FILE *ficres,double p[], int npar, int ncovmodel, int nlstate, double ftol, double (*func)(double []))
 2191: {
 2192:   int i,j, iter=0;
 2193:   double **xi;
 2194:   double fret;
 2195:   double fretone; /* Only one call to likelihood */
 2196:   /*  char filerespow[FILENAMELENGTH];*/
 2197: 
 2198: #ifdef NLOPT
 2199:   int creturn;
 2200:   nlopt_opt opt;
 2201:   /* double lb[9] = { -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL, -HUGE_VAL }; /\* lower bounds *\/ */
 2202:   double *lb;
 2203:   double minf; /* the minimum objective value, upon return */
 2204:   double * p1; /* Shifted parameters from 0 instead of 1 */
 2205:   myfunc_data dinst, *d = &dinst;
 2206: #endif
 2207: 
 2208: 
 2209:   xi=matrix(1,npar,1,npar);
 2210:   for (i=1;i<=npar;i++)
 2211:     for (j=1;j<=npar;j++)
 2212:       xi[i][j]=(i==j ? 1.0 : 0.0);
 2213:   printf("Powell\n");  fprintf(ficlog,"Powell\n");
 2214:   strcpy(filerespow,"pow"); 
 2215:   strcat(filerespow,fileres);
 2216:   if((ficrespow=fopen(filerespow,"w"))==NULL) {
 2217:     printf("Problem with resultfile: %s\n", filerespow);
 2218:     fprintf(ficlog,"Problem with resultfile: %s\n", filerespow);
 2219:   }
 2220:   fprintf(ficrespow,"# Powell\n# iter -2*LL");
 2221:   for (i=1;i<=nlstate;i++)
 2222:     for(j=1;j<=nlstate+ndeath;j++)
 2223:       if(j!=i)fprintf(ficrespow," p%1d%1d",i,j);
 2224:   fprintf(ficrespow,"\n");
 2225: #ifdef POWELL
 2226:   powell(p,xi,npar,ftol,&iter,&fret,func);
 2227: #endif
 2228: 
 2229: #ifdef NLOPT
 2230: #ifdef NEWUOA
 2231:   opt = nlopt_create(NLOPT_LN_NEWUOA,npar);
 2232: #else
 2233:   opt = nlopt_create(NLOPT_LN_BOBYQA,npar);
 2234: #endif
 2235:   lb=vector(0,npar-1);
 2236:   for (i=0;i<npar;i++) lb[i]= -HUGE_VAL;
 2237:   nlopt_set_lower_bounds(opt, lb);
 2238:   nlopt_set_initial_step1(opt, 0.1);
 2239:   
 2240:   p1= (p+1); /*  p *(p+1)@8 and p *(p1)@8 are equal p1[0]=p[1] */
 2241:   d->function = func;
 2242:   printf(" Func %.12lf \n",myfunc(npar,p1,NULL,d));
 2243:   nlopt_set_min_objective(opt, myfunc, d);
 2244:   nlopt_set_xtol_rel(opt, ftol);
 2245:   if ((creturn=nlopt_optimize(opt, p1, &minf)) < 0) {
 2246:     printf("nlopt failed! %d\n",creturn); 
 2247:   }
 2248:   else {
 2249:     printf("found minimum after %d evaluations (NLOPT=%d)\n", countcallfunc ,NLOPT);
 2250:     printf("found minimum at f(%g,%g) = %0.10g\n", p[0], p[1], minf);
 2251:     iter=1; /* not equal */
 2252:   }
 2253:   nlopt_destroy(opt);
 2254: #endif
 2255:   free_matrix(xi,1,npar,1,npar);
 2256:   fclose(ficrespow);
 2257:   printf("\n#Number of iterations & function calls = %d & %d, -2 Log likelihood = %.12f\n",iter, countcallfunc,func(p));
 2258:   fprintf(ficlog,"\n#Number of iterations & function calls = %d & %d, -2 Log likelihood = %.12f\n",iter, countcallfunc,func(p));
 2259:   fprintf(ficres,"\n#Number of iterations & function calls = %d & %d, -2 Log likelihood = %.12f\n",iter, countcallfunc,func(p));
 2260: 
 2261: }
 2262: 
 2263: /**** Computes Hessian and covariance matrix ***/
 2264: void hesscov(double **matcov, double p[], int npar, double delti[], double ftolhess, double (*func)(double []))
 2265: {
 2266:   double  **a,**y,*x,pd;
 2267:   double **hess;
 2268:   int i, j;
 2269:   int *indx;
 2270: 
 2271:   double hessii(double p[], double delta, int theta, double delti[],double (*func)(double []),int npar);
 2272:   double hessij(double p[], double delti[], int i, int j,double (*func)(double []),int npar);
 2273:   void lubksb(double **a, int npar, int *indx, double b[]) ;
 2274:   void ludcmp(double **a, int npar, int *indx, double *d) ;
 2275:   double gompertz(double p[]);
 2276:   hess=matrix(1,npar,1,npar);
 2277: 
 2278:   printf("\nCalculation of the hessian matrix. Wait...\n");
 2279:   fprintf(ficlog,"\nCalculation of the hessian matrix. Wait...\n");
 2280:   for (i=1;i<=npar;i++){
 2281:     printf("%d",i);fflush(stdout);
 2282:     fprintf(ficlog,"%d",i);fflush(ficlog);
 2283:    
 2284:      hess[i][i]=hessii(p,ftolhess,i,delti,func,npar);
 2285:     
 2286:     /*  printf(" %f ",p[i]);
 2287: 	printf(" %lf %lf %lf",hess[i][i],ftolhess,delti[i]);*/
 2288:   }
 2289:   
 2290:   for (i=1;i<=npar;i++) {
 2291:     for (j=1;j<=npar;j++)  {
 2292:       if (j>i) { 
 2293: 	printf(".%d%d",i,j);fflush(stdout);
 2294: 	fprintf(ficlog,".%d%d",i,j);fflush(ficlog);
 2295: 	hess[i][j]=hessij(p,delti,i,j,func,npar);
 2296: 	
 2297: 	hess[j][i]=hess[i][j];    
 2298: 	/*printf(" %lf ",hess[i][j]);*/
 2299:       }
 2300:     }
 2301:   }
 2302:   printf("\n");
 2303:   fprintf(ficlog,"\n");
 2304: 
 2305:   printf("\nInverting the hessian to get the covariance matrix. Wait...\n");
 2306:   fprintf(ficlog,"\nInverting the hessian to get the covariance matrix. Wait...\n");
 2307:   
 2308:   a=matrix(1,npar,1,npar);
 2309:   y=matrix(1,npar,1,npar);
 2310:   x=vector(1,npar);
 2311:   indx=ivector(1,npar);
 2312:   for (i=1;i<=npar;i++)
 2313:     for (j=1;j<=npar;j++) a[i][j]=hess[i][j];
 2314:   ludcmp(a,npar,indx,&pd);
 2315: 
 2316:   for (j=1;j<=npar;j++) {
 2317:     for (i=1;i<=npar;i++) x[i]=0;
 2318:     x[j]=1;
 2319:     lubksb(a,npar,indx,x);
 2320:     for (i=1;i<=npar;i++){ 
 2321:       matcov[i][j]=x[i];
 2322:     }
 2323:   }
 2324: 
 2325:   printf("\n#Hessian matrix#\n");
 2326:   fprintf(ficlog,"\n#Hessian matrix#\n");
 2327:   for (i=1;i<=npar;i++) { 
 2328:     for (j=1;j<=npar;j++) { 
 2329:       printf("%.3e ",hess[i][j]);
 2330:       fprintf(ficlog,"%.3e ",hess[i][j]);
 2331:     }
 2332:     printf("\n");
 2333:     fprintf(ficlog,"\n");
 2334:   }
 2335: 
 2336:   /* Recompute Inverse */
 2337:   for (i=1;i<=npar;i++)
 2338:     for (j=1;j<=npar;j++) a[i][j]=matcov[i][j];
 2339:   ludcmp(a,npar,indx,&pd);
 2340: 
 2341:   /*  printf("\n#Hessian matrix recomputed#\n");
 2342: 
 2343:   for (j=1;j<=npar;j++) {
 2344:     for (i=1;i<=npar;i++) x[i]=0;
 2345:     x[j]=1;
 2346:     lubksb(a,npar,indx,x);
 2347:     for (i=1;i<=npar;i++){ 
 2348:       y[i][j]=x[i];
 2349:       printf("%.3e ",y[i][j]);
 2350:       fprintf(ficlog,"%.3e ",y[i][j]);
 2351:     }
 2352:     printf("\n");
 2353:     fprintf(ficlog,"\n");
 2354:   }
 2355:   */
 2356: 
 2357:   free_matrix(a,1,npar,1,npar);
 2358:   free_matrix(y,1,npar,1,npar);
 2359:   free_vector(x,1,npar);
 2360:   free_ivector(indx,1,npar);
 2361:   free_matrix(hess,1,npar,1,npar);
 2362: 
 2363: 
 2364: }
 2365: 
 2366: /*************** hessian matrix ****************/
 2367: double hessii(double x[], double delta, int theta, double delti[], double (*func)(double []), int npar)
 2368: {
 2369:   int i;
 2370:   int l=1, lmax=20;
 2371:   double k1,k2;
 2372:   double p2[MAXPARM+1]; /* identical to x */
 2373:   double res;
 2374:   double delt=0.0001, delts, nkhi=10.,nkhif=1., khi=1.e-4;
 2375:   double fx;
 2376:   int k=0,kmax=10;
 2377:   double l1;
 2378: 
 2379:   fx=func(x);
 2380:   for (i=1;i<=npar;i++) p2[i]=x[i];
 2381:   for(l=0 ; l <=lmax; l++){  /* Enlarging the zone around the Maximum */
 2382:     l1=pow(10,l);
 2383:     delts=delt;
 2384:     for(k=1 ; k <kmax; k=k+1){
 2385:       delt = delta*(l1*k);
 2386:       p2[theta]=x[theta] +delt;
 2387:       k1=func(p2)-fx;   /* Might be negative if too close to the theoretical maximum */
 2388:       p2[theta]=x[theta]-delt;
 2389:       k2=func(p2)-fx;
 2390:       /*res= (k1-2.0*fx+k2)/delt/delt; */
 2391:       res= (k1+k2)/delt/delt/2.; /* Divided by because L and not 2*L */
 2392:       
 2393: #ifdef DEBUGHESS
 2394:       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);
 2395:       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);
 2396: #endif
 2397:       /*if(fabs(k1-2.0*fx+k2) <1.e-13){ */
 2398:       if((k1 <khi/nkhi/2.) || (k2 <khi/nkhi/2.)){
 2399: 	k=kmax;
 2400:       }
 2401:       else if((k1 >khi/nkhif) || (k2 >khi/nkhif)){ /* Keeps lastvalue before 3.84/2 KHI2 5% 1d.f. */
 2402: 	k=kmax; l=lmax*10;
 2403:       }
 2404:       else if((k1 >khi/nkhi) || (k2 >khi/nkhi)){ 
 2405: 	delts=delt;
 2406:       }
 2407:     }
 2408:   }
 2409:   delti[theta]=delts;
 2410:   return res; 
 2411:   
 2412: }
 2413: 
 2414: double hessij( double x[], double delti[], int thetai,int thetaj,double (*func)(double []),int npar)
 2415: {
 2416:   int i;
 2417:   int l=1, lmax=20;
 2418:   double k1,k2,k3,k4,res,fx;
 2419:   double p2[MAXPARM+1];
 2420:   int k;
 2421: 
 2422:   fx=func(x);
 2423:   for (k=1; k<=2; k++) {
 2424:     for (i=1;i<=npar;i++) p2[i]=x[i];
 2425:     p2[thetai]=x[thetai]+delti[thetai]/k;
 2426:     p2[thetaj]=x[thetaj]+delti[thetaj]/k;
 2427:     k1=func(p2)-fx;
 2428:   
 2429:     p2[thetai]=x[thetai]+delti[thetai]/k;
 2430:     p2[thetaj]=x[thetaj]-delti[thetaj]/k;
 2431:     k2=func(p2)-fx;
 2432:   
 2433:     p2[thetai]=x[thetai]-delti[thetai]/k;
 2434:     p2[thetaj]=x[thetaj]+delti[thetaj]/k;
 2435:     k3=func(p2)-fx;
 2436:   
 2437:     p2[thetai]=x[thetai]-delti[thetai]/k;
 2438:     p2[thetaj]=x[thetaj]-delti[thetaj]/k;
 2439:     k4=func(p2)-fx;
 2440:     res=(k1-k2-k3+k4)/4.0/delti[thetai]*k/delti[thetaj]*k/2.; /* Because of L not 2*L */
 2441: #ifdef DEBUG
 2442:     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);
 2443:     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);
 2444: #endif
 2445:   }
 2446:   return res;
 2447: }
 2448: 
 2449: /************** Inverse of matrix **************/
 2450: void ludcmp(double **a, int n, int *indx, double *d) 
 2451: { 
 2452:   int i,imax,j,k; 
 2453:   double big,dum,sum,temp; 
 2454:   double *vv; 
 2455:  
 2456:   vv=vector(1,n); 
 2457:   *d=1.0; 
 2458:   for (i=1;i<=n;i++) { 
 2459:     big=0.0; 
 2460:     for (j=1;j<=n;j++) 
 2461:       if ((temp=fabs(a[i][j])) > big) big=temp; 
 2462:     if (big == 0.0) nrerror("Singular matrix in routine ludcmp"); 
 2463:     vv[i]=1.0/big; 
 2464:   } 
 2465:   for (j=1;j<=n;j++) { 
 2466:     for (i=1;i<j;i++) { 
 2467:       sum=a[i][j]; 
 2468:       for (k=1;k<i;k++) sum -= a[i][k]*a[k][j]; 
 2469:       a[i][j]=sum; 
 2470:     } 
 2471:     big=0.0; 
 2472:     for (i=j;i<=n;i++) { 
 2473:       sum=a[i][j]; 
 2474:       for (k=1;k<j;k++) 
 2475: 	sum -= a[i][k]*a[k][j]; 
 2476:       a[i][j]=sum; 
 2477:       if ( (dum=vv[i]*fabs(sum)) >= big) { 
 2478: 	big=dum; 
 2479: 	imax=i; 
 2480:       } 
 2481:     } 
 2482:     if (j != imax) { 
 2483:       for (k=1;k<=n;k++) { 
 2484: 	dum=a[imax][k]; 
 2485: 	a[imax][k]=a[j][k]; 
 2486: 	a[j][k]=dum; 
 2487:       } 
 2488:       *d = -(*d); 
 2489:       vv[imax]=vv[j]; 
 2490:     } 
 2491:     indx[j]=imax; 
 2492:     if (a[j][j] == 0.0) a[j][j]=TINY; 
 2493:     if (j != n) { 
 2494:       dum=1.0/(a[j][j]); 
 2495:       for (i=j+1;i<=n;i++) a[i][j] *= dum; 
 2496:     } 
 2497:   } 
 2498:   free_vector(vv,1,n);  /* Doesn't work */
 2499: ;
 2500: } 
 2501: 
 2502: void lubksb(double **a, int n, int *indx, double b[]) 
 2503: { 
 2504:   int i,ii=0,ip,j; 
 2505:   double sum; 
 2506:  
 2507:   for (i=1;i<=n;i++) { 
 2508:     ip=indx[i]; 
 2509:     sum=b[ip]; 
 2510:     b[ip]=b[i]; 
 2511:     if (ii) 
 2512:       for (j=ii;j<=i-1;j++) sum -= a[i][j]*b[j]; 
 2513:     else if (sum) ii=i; 
 2514:     b[i]=sum; 
 2515:   } 
 2516:   for (i=n;i>=1;i--) { 
 2517:     sum=b[i]; 
 2518:     for (j=i+1;j<=n;j++) sum -= a[i][j]*b[j]; 
 2519:     b[i]=sum/a[i][i]; 
 2520:   } 
 2521: } 
 2522: 
 2523: void pstamp(FILE *fichier)
 2524: {
 2525:   fprintf(fichier,"# %s.%s\n#%s\n#%s\n# %s", optionfilefiname,optionfilext,version,fullversion,strstart);
 2526: }
 2527: 
 2528: /************ Frequencies ********************/
 2529: 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[])
 2530: {  /* Some frequencies */
 2531:   
 2532:   int i, m, jk, j1, bool, z1,j;
 2533:   int first;
 2534:   double ***freq; /* Frequencies */
 2535:   double *pp, **prop;
 2536:   double pos,posprop, k2, dateintsum=0,k2cpt=0;
 2537:   char fileresp[FILENAMELENGTH];
 2538:   
 2539:   pp=vector(1,nlstate);
 2540:   prop=matrix(1,nlstate,iagemin,iagemax+3);
 2541:   strcpy(fileresp,"p");
 2542:   strcat(fileresp,fileres);
 2543:   if((ficresp=fopen(fileresp,"w"))==NULL) {
 2544:     printf("Problem with prevalence resultfile: %s\n", fileresp);
 2545:     fprintf(ficlog,"Problem with prevalence resultfile: %s\n", fileresp);
 2546:     exit(0);
 2547:   }
 2548:   freq= ma3x(-5,nlstate+ndeath,-5,nlstate+ndeath,iagemin,iagemax+3);
 2549:   j1=0;
 2550:   
 2551:   j=cptcoveff;
 2552:   if (cptcovn<1) {j=1;ncodemax[1]=1;}
 2553: 
 2554:   first=1;
 2555: 
 2556:   /* for(k1=1; k1<=j ; k1++){ */  /* Loop on covariates */
 2557:   /*  for(i1=1; i1<=ncodemax[k1];i1++){ */ /* Now it is 2 */
 2558:   /*    j1++; */
 2559:   for (j1 = 1; j1 <= (int) pow(2,cptcoveff); j1++){
 2560:       /*printf("cptcoveff=%d Tvaraff=%d", cptcoveff,Tvaraff[1]);
 2561: 	scanf("%d", i);*/
 2562:       for (i=-5; i<=nlstate+ndeath; i++)  
 2563: 	for (jk=-5; jk<=nlstate+ndeath; jk++)  
 2564: 	  for(m=iagemin; m <= iagemax+3; m++)
 2565: 	    freq[i][jk][m]=0;
 2566:       
 2567:       for (i=1; i<=nlstate; i++)  
 2568: 	for(m=iagemin; m <= iagemax+3; m++)
 2569: 	  prop[i][m]=0;
 2570:       
 2571:       dateintsum=0;
 2572:       k2cpt=0;
 2573:       for (i=1; i<=imx; i++) {
 2574: 	bool=1;
 2575: 	if  (cptcovn>0) { /* Filter is here: Must be looked at for model=V1+V2+V3+V4 */
 2576: 	  for (z1=1; z1<=cptcoveff; z1++)       
 2577:             if (covar[Tvaraff[z1]][i]!= nbcode[Tvaraff[z1]][codtab[j1][z1]]){
 2578:                 /* Tests if the value of each of the covariates of i is equal to filter j1 */
 2579:               bool=0;
 2580:               /* 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", 
 2581:                 bool,i,z1, z1, Tvaraff[z1],i,covar[Tvaraff[z1]][i],j1,z1,codtab[j1][z1],
 2582:                 j1,z1,nbcode[Tvaraff[z1]][codtab[j1][z1]],j1);*/
 2583:               /* For j1=7 in V1+V2+V3+V4 = 0 1 1 0 and codtab[7][3]=1 and nbcde[3][?]=1*/
 2584:             } 
 2585: 	}
 2586:  
 2587: 	if (bool==1){
 2588: 	  for(m=firstpass; m<=lastpass; m++){
 2589: 	    k2=anint[m][i]+(mint[m][i]/12.);
 2590: 	    /*if ((k2>=dateprev1) && (k2<=dateprev2)) {*/
 2591: 	      if(agev[m][i]==0) agev[m][i]=iagemax+1;
 2592: 	      if(agev[m][i]==1) agev[m][i]=iagemax+2;
 2593: 	      if (s[m][i]>0 && s[m][i]<=nlstate) prop[s[m][i]][(int)agev[m][i]] += weight[i];
 2594: 	      if (m<lastpass) {
 2595: 		freq[s[m][i]][s[m+1][i]][(int)agev[m][i]] += weight[i];
 2596: 		freq[s[m][i]][s[m+1][i]][iagemax+3] += weight[i];
 2597: 	      }
 2598: 	      
 2599: 	      if ((agev[m][i]>1) && (agev[m][i]< (iagemax+3))) {
 2600: 		dateintsum=dateintsum+k2;
 2601: 		k2cpt++;
 2602: 	      }
 2603: 	      /*}*/
 2604: 	  }
 2605: 	}
 2606:       } /* end i */
 2607:        
 2608:       /*      fprintf(ficresp, "#Count between %.lf/%.lf/%.lf and %.lf/%.lf/%.lf\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2);*/
 2609:       pstamp(ficresp);
 2610:       if  (cptcovn>0) {
 2611: 	fprintf(ficresp, "\n#********** Variable "); 
 2612: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresp, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 2613: 	fprintf(ficresp, "**********\n#");
 2614: 	fprintf(ficlog, "\n#********** Variable "); 
 2615: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficlog, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 2616: 	fprintf(ficlog, "**********\n#");
 2617:       }
 2618:       for(i=1; i<=nlstate;i++) 
 2619: 	fprintf(ficresp, " Age Prev(%d) N(%d) N",i,i);
 2620:       fprintf(ficresp, "\n");
 2621:       
 2622:       for(i=iagemin; i <= iagemax+3; i++){
 2623: 	if(i==iagemax+3){
 2624: 	  fprintf(ficlog,"Total");
 2625: 	}else{
 2626: 	  if(first==1){
 2627: 	    first=0;
 2628: 	    printf("See log file for details...\n");
 2629: 	  }
 2630: 	  fprintf(ficlog,"Age %d", i);
 2631: 	}
 2632: 	for(jk=1; jk <=nlstate ; jk++){
 2633: 	  for(m=-1, pp[jk]=0; m <=nlstate+ndeath ; m++)
 2634: 	    pp[jk] += freq[jk][m][i]; 
 2635: 	}
 2636: 	for(jk=1; jk <=nlstate ; jk++){
 2637: 	  for(m=-1, pos=0; m <=0 ; m++)
 2638: 	    pos += freq[jk][m][i];
 2639: 	  if(pp[jk]>=1.e-10){
 2640: 	    if(first==1){
 2641: 	      printf(" %d.=%.0f loss[%d]=%.1f%%",jk,pp[jk],jk,100*pos/pp[jk]);
 2642: 	    }
 2643: 	    fprintf(ficlog," %d.=%.0f loss[%d]=%.1f%%",jk,pp[jk],jk,100*pos/pp[jk]);
 2644: 	  }else{
 2645: 	    if(first==1)
 2646: 	      printf(" %d.=%.0f loss[%d]=NaNQ%%",jk,pp[jk],jk);
 2647: 	    fprintf(ficlog," %d.=%.0f loss[%d]=NaNQ%%",jk,pp[jk],jk);
 2648: 	  }
 2649: 	}
 2650: 
 2651: 	for(jk=1; jk <=nlstate ; jk++){
 2652: 	  for(m=0, pp[jk]=0; m <=nlstate+ndeath; m++)
 2653: 	    pp[jk] += freq[jk][m][i];
 2654: 	}	
 2655: 	for(jk=1,pos=0,posprop=0; jk <=nlstate ; jk++){
 2656: 	  pos += pp[jk];
 2657: 	  posprop += prop[jk][i];
 2658: 	}
 2659: 	for(jk=1; jk <=nlstate ; jk++){
 2660: 	  if(pos>=1.e-5){
 2661: 	    if(first==1)
 2662: 	      printf(" %d.=%.0f prev[%d]=%.1f%%",jk,pp[jk],jk,100*pp[jk]/pos);
 2663: 	    fprintf(ficlog," %d.=%.0f prev[%d]=%.1f%%",jk,pp[jk],jk,100*pp[jk]/pos);
 2664: 	  }else{
 2665: 	    if(first==1)
 2666: 	      printf(" %d.=%.0f prev[%d]=NaNQ%%",jk,pp[jk],jk);
 2667: 	    fprintf(ficlog," %d.=%.0f prev[%d]=NaNQ%%",jk,pp[jk],jk);
 2668: 	  }
 2669: 	  if( i <= iagemax){
 2670: 	    if(pos>=1.e-5){
 2671: 	      fprintf(ficresp," %d %.5f %.0f %.0f",i,prop[jk][i]/posprop, prop[jk][i],posprop);
 2672: 	      /*probs[i][jk][j1]= pp[jk]/pos;*/
 2673: 	      /*printf("\ni=%d jk=%d j1=%d %.5f %.0f %.0f %f",i,jk,j1,pp[jk]/pos, pp[jk],pos,probs[i][jk][j1]);*/
 2674: 	    }
 2675: 	    else
 2676: 	      fprintf(ficresp," %d NaNq %.0f %.0f",i,prop[jk][i],posprop);
 2677: 	  }
 2678: 	}
 2679: 	
 2680: 	for(jk=-1; jk <=nlstate+ndeath; jk++)
 2681: 	  for(m=-1; m <=nlstate+ndeath; m++)
 2682: 	    if(freq[jk][m][i] !=0 ) {
 2683: 	    if(first==1)
 2684: 	      printf(" %d%d=%.0f",jk,m,freq[jk][m][i]);
 2685: 	      fprintf(ficlog," %d%d=%.0f",jk,m,freq[jk][m][i]);
 2686: 	    }
 2687: 	if(i <= iagemax)
 2688: 	  fprintf(ficresp,"\n");
 2689: 	if(first==1)
 2690: 	  printf("Others in log...\n");
 2691: 	fprintf(ficlog,"\n");
 2692:       }
 2693:       /*}*/
 2694:   }
 2695:   dateintmean=dateintsum/k2cpt; 
 2696:  
 2697:   fclose(ficresp);
 2698:   free_ma3x(freq,-5,nlstate+ndeath,-5,nlstate+ndeath, iagemin, iagemax+3);
 2699:   free_vector(pp,1,nlstate);
 2700:   free_matrix(prop,1,nlstate,iagemin, iagemax+3);
 2701:   /* End of Freq */
 2702: }
 2703: 
 2704: /************ Prevalence ********************/
 2705: 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)
 2706: {  
 2707:   /* Compute observed prevalence between dateprev1 and dateprev2 by counting the number of people
 2708:      in each health status at the date of interview (if between dateprev1 and dateprev2).
 2709:      We still use firstpass and lastpass as another selection.
 2710:   */
 2711:  
 2712:   int i, m, jk, j1, bool, z1,j;
 2713: 
 2714:   double **prop;
 2715:   double posprop; 
 2716:   double  y2; /* in fractional years */
 2717:   int iagemin, iagemax;
 2718:   int first; /** to stop verbosity which is redirected to log file */
 2719: 
 2720:   iagemin= (int) agemin;
 2721:   iagemax= (int) agemax;
 2722:   /*pp=vector(1,nlstate);*/
 2723:   prop=matrix(1,nlstate,iagemin,iagemax+3); 
 2724:   /*  freq=ma3x(-1,nlstate+ndeath,-1,nlstate+ndeath,iagemin,iagemax+3);*/
 2725:   j1=0;
 2726:   
 2727:   /*j=cptcoveff;*/
 2728:   if (cptcovn<1) {j=1;ncodemax[1]=1;}
 2729:   
 2730:   first=1;
 2731:   for(j1=1; j1<= (int) pow(2,cptcoveff);j1++){
 2732:     /*for(i1=1; i1<=ncodemax[k1];i1++){
 2733:       j1++;*/
 2734:       
 2735:       for (i=1; i<=nlstate; i++)  
 2736: 	for(m=iagemin; m <= iagemax+3; m++)
 2737: 	  prop[i][m]=0.0;
 2738:      
 2739:       for (i=1; i<=imx; i++) { /* Each individual */
 2740: 	bool=1;
 2741: 	if  (cptcovn>0) {
 2742: 	  for (z1=1; z1<=cptcoveff; z1++) 
 2743: 	    if (covar[Tvaraff[z1]][i]!= nbcode[Tvaraff[z1]][codtab[j1][z1]]) 
 2744: 	      bool=0;
 2745: 	} 
 2746: 	if (bool==1) { 
 2747: 	  for(m=firstpass; m<=lastpass; m++){/* Other selection (we can limit to certain interviews*/
 2748: 	    y2=anint[m][i]+(mint[m][i]/12.); /* Fractional date in year */
 2749: 	    if ((y2>=dateprev1) && (y2<=dateprev2)) { /* Here is the main selection (fractional years) */
 2750: 	      if(agev[m][i]==0) agev[m][i]=iagemax+1;
 2751: 	      if(agev[m][i]==1) agev[m][i]=iagemax+2;
 2752: 	      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); 
 2753:  	      if (s[m][i]>0 && s[m][i]<=nlstate) { 
 2754: 		/*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]]);*/
 2755:  		prop[s[m][i]][(int)agev[m][i]] += weight[i];
 2756:  		prop[s[m][i]][iagemax+3] += weight[i]; 
 2757:  	      } 
 2758: 	    }
 2759: 	  } /* end selection of waves */
 2760: 	}
 2761:       }
 2762:       for(i=iagemin; i <= iagemax+3; i++){  
 2763:  	for(jk=1,posprop=0; jk <=nlstate ; jk++) { 
 2764:  	  posprop += prop[jk][i]; 
 2765:  	} 
 2766: 	
 2767:  	for(jk=1; jk <=nlstate ; jk++){	    
 2768:  	  if( i <=  iagemax){ 
 2769:  	    if(posprop>=1.e-5){ 
 2770:  	      probs[i][jk][j1]= prop[jk][i]/posprop;
 2771:  	    } else{
 2772: 	      if(first==1){
 2773: 		first=0;
 2774: 		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]);
 2775: 	      }
 2776: 	    }
 2777:  	  } 
 2778:  	}/* end jk */ 
 2779:       }/* end i */ 
 2780:     /*} *//* end i1 */
 2781:   } /* end j1 */
 2782:   
 2783:   /*  free_ma3x(freq,-1,nlstate+ndeath,-1,nlstate+ndeath, iagemin, iagemax+3);*/
 2784:   /*free_vector(pp,1,nlstate);*/
 2785:   free_matrix(prop,1,nlstate, iagemin,iagemax+3);
 2786: }  /* End of prevalence */
 2787: 
 2788: /************* Waves Concatenation ***************/
 2789: 
 2790: 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)
 2791: {
 2792:   /* Concatenates waves: wav[i] is the number of effective (useful waves) of individual i.
 2793:      Death is a valid wave (if date is known).
 2794:      mw[mi][i] is the mi (mi=1 to wav[i])  effective wave of individual i
 2795:      dh[m][i] or dh[mw[mi][i]][i] is the delay between two effective waves m=mw[mi][i]
 2796:      and mw[mi+1][i]. dh depends on stepm.
 2797:      */
 2798: 
 2799:   int i, mi, m;
 2800:   /* int j, k=0,jk, ju, jl,jmin=1e+5, jmax=-1;
 2801:      double sum=0., jmean=0.;*/
 2802:   int first;
 2803:   int j, k=0,jk, ju, jl;
 2804:   double sum=0.;
 2805:   first=0;
 2806:   jmin=100000;
 2807:   jmax=-1;
 2808:   jmean=0.;
 2809:   for(i=1; i<=imx; i++){
 2810:     mi=0;
 2811:     m=firstpass;
 2812:     while(s[m][i] <= nlstate){
 2813:       if(s[m][i]>=1 || s[m][i]==-2 || s[m][i]==-4 || s[m][i]==-5)
 2814: 	mw[++mi][i]=m;
 2815:       if(m >=lastpass)
 2816: 	break;
 2817:       else
 2818: 	m++;
 2819:     }/* end while */
 2820:     if (s[m][i] > nlstate){
 2821:       mi++;	/* Death is another wave */
 2822:       /* if(mi==0)  never been interviewed correctly before death */
 2823: 	 /* Only death is a correct wave */
 2824:       mw[mi][i]=m;
 2825:     }
 2826: 
 2827:     wav[i]=mi;
 2828:     if(mi==0){
 2829:       nbwarn++;
 2830:       if(first==0){
 2831: 	printf("Warning! No valid information for individual %ld line=%d (skipped) and may be others, see log file\n",num[i],i);
 2832: 	first=1;
 2833:       }
 2834:       if(first==1){
 2835: 	fprintf(ficlog,"Warning! No valid information for individual %ld line=%d (skipped)\n",num[i],i);
 2836:       }
 2837:     } /* end mi==0 */
 2838:   } /* End individuals */
 2839: 
 2840:   for(i=1; i<=imx; i++){
 2841:     for(mi=1; mi<wav[i];mi++){
 2842:       if (stepm <=0)
 2843: 	dh[mi][i]=1;
 2844:       else{
 2845: 	if (s[mw[mi+1][i]][i] > nlstate) { /* A death */
 2846: 	  if (agedc[i] < 2*AGESUP) {
 2847: 	    j= rint(agedc[i]*12-agev[mw[mi][i]][i]*12); 
 2848: 	    if(j==0) j=1;  /* Survives at least one month after exam */
 2849: 	    else if(j<0){
 2850: 	      nberr++;
 2851: 	      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]);
 2852: 	      j=1; /* Temporary Dangerous patch */
 2853: 	      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);
 2854: 	      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]);
 2855: 	      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);
 2856: 	    }
 2857: 	    k=k+1;
 2858: 	    if (j >= jmax){
 2859: 	      jmax=j;
 2860: 	      ijmax=i;
 2861: 	    }
 2862: 	    if (j <= jmin){
 2863: 	      jmin=j;
 2864: 	      ijmin=i;
 2865: 	    }
 2866: 	    sum=sum+j;
 2867: 	    /*if (j<0) printf("j=%d num=%d \n",j,i);*/
 2868: 	    /*	  printf("%d %d %d %d\n", s[mw[mi][i]][i] ,s[mw[mi+1][i]][i],j,i);*/
 2869: 	  }
 2870: 	}
 2871: 	else{
 2872: 	  j= rint( (agev[mw[mi+1][i]][i]*12 - agev[mw[mi][i]][i]*12));
 2873: /* 	  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]); */
 2874: 
 2875: 	  k=k+1;
 2876: 	  if (j >= jmax) {
 2877: 	    jmax=j;
 2878: 	    ijmax=i;
 2879: 	  }
 2880: 	  else if (j <= jmin){
 2881: 	    jmin=j;
 2882: 	    ijmin=i;
 2883: 	  }
 2884: 	  /*	    if (j<10) printf("j=%d jmin=%d num=%d ",j,jmin,i); */
 2885: 	  /*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]);*/
 2886: 	  if(j<0){
 2887: 	    nberr++;
 2888: 	    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]);
 2889: 	    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]);
 2890: 	  }
 2891: 	  sum=sum+j;
 2892: 	}
 2893: 	jk= j/stepm;
 2894: 	jl= j -jk*stepm;
 2895: 	ju= j -(jk+1)*stepm;
 2896: 	if(mle <=1){ /* only if we use a the linear-interpoloation pseudo-likelihood */
 2897: 	  if(jl==0){
 2898: 	    dh[mi][i]=jk;
 2899: 	    bh[mi][i]=0;
 2900: 	  }else{ /* We want a negative bias in order to only have interpolation ie
 2901: 		  * to avoid the price of an extra matrix product in likelihood */
 2902: 	    dh[mi][i]=jk+1;
 2903: 	    bh[mi][i]=ju;
 2904: 	  }
 2905: 	}else{
 2906: 	  if(jl <= -ju){
 2907: 	    dh[mi][i]=jk;
 2908: 	    bh[mi][i]=jl;	/* bias is positive if real duration
 2909: 				 * is higher than the multiple of stepm and negative otherwise.
 2910: 				 */
 2911: 	  }
 2912: 	  else{
 2913: 	    dh[mi][i]=jk+1;
 2914: 	    bh[mi][i]=ju;
 2915: 	  }
 2916: 	  if(dh[mi][i]==0){
 2917: 	    dh[mi][i]=1; /* At least one step */
 2918: 	    bh[mi][i]=ju; /* At least one step */
 2919: 	    /*  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);*/
 2920: 	  }
 2921: 	} /* end if mle */
 2922:       }
 2923:     } /* end wave */
 2924:   }
 2925:   jmean=sum/k;
 2926:   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);
 2927:   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);
 2928:  }
 2929: 
 2930: /*********** Tricode ****************************/
 2931: void tricode(int *Tvar, int **nbcode, int imx, int *Ndum)
 2932: {
 2933:   /**< Uses cptcovn+2*cptcovprod as the number of covariates */
 2934:   /*	  Tvar[i]=atoi(stre);  find 'n' in Vn and stores in Tvar. If model=V2+V1 Tvar[1]=2 and Tvar[2]=1 
 2935:    * Boring subroutine which should only output nbcode[Tvar[j]][k]
 2936:    * Tvar[5] in V2+V1+V3*age+V2*V4 is 2 (V2)
 2937:    * nbcode[Tvar[j]][1]= 
 2938:   */
 2939: 
 2940:   int ij=1, k=0, j=0, i=0, maxncov=NCOVMAX;
 2941:   int modmaxcovj=0; /* Modality max of covariates j */
 2942:   int cptcode=0; /* Modality max of covariates j */
 2943:   int modmincovj=0; /* Modality min of covariates j */
 2944: 
 2945: 
 2946:   cptcoveff=0; 
 2947:  
 2948:   for (k=-1; k < maxncov; k++) Ndum[k]=0;
 2949:   for (k=1; k <= maxncov; k++) ncodemax[k]=0; /* Horrible constant again replaced by NCOVMAX */
 2950: 
 2951:   /* Loop on covariates without age and products */
 2952:   for (j=1; j<=(cptcovs); j++) { /* model V1 + V2*age+ V3 + V3*V4 : V1 + V3 = 2 only */
 2953:     for (i=1; i<=imx; i++) { /* Lopp on individuals: reads the data file to get the maximum value of the 
 2954: 			       modality of this covariate Vj*/ 
 2955:       ij=(int)(covar[Tvar[j]][i]); /* ij=0 or 1 or -1. Value of the covariate Tvar[j] for individual i
 2956: 				    * If product of Vn*Vm, still boolean *:
 2957: 				    * If it was coded 1, 2, 3, 4 should be splitted into 3 boolean variables
 2958: 				    * 1 => 0 0 0, 2 => 0 0 1, 3 => 0 1 1, 4=1 0 0   */
 2959:       /* Finds for covariate j, n=Tvar[j] of Vn . ij is the
 2960: 				      modality of the nth covariate of individual i. */
 2961:       if (ij > modmaxcovj)
 2962:         modmaxcovj=ij; 
 2963:       else if (ij < modmincovj) 
 2964: 	modmincovj=ij; 
 2965:       if ((ij < -1) && (ij > NCOVMAX)){
 2966: 	printf( "Error: minimal is less than -1 or maximal is bigger than %d. Exiting. \n", NCOVMAX );
 2967: 	exit(1);
 2968:       }else
 2969:       Ndum[ij]++; /*counts and stores the occurence of this modality 0, 1, -1*/
 2970:       /*  If coded 1, 2, 3 , counts the number of 1 Ndum[1], number of 2, Ndum[2], etc */
 2971:       /*printf("i=%d ij=%d Ndum[ij]=%d imx=%d",i,ij,Ndum[ij],imx);*/
 2972:       /* getting the maximum value of the modality of the covariate
 2973: 	 (should be 0 or 1 now) Tvar[j]. If V=sex and male is coded 0 and
 2974: 	 female is 1, then modmaxcovj=1.*/
 2975:     }
 2976:     printf(" Minimal and maximal values of %d th covariate V%d: min=%d max=%d \n", j, Tvar[j], modmincovj, modmaxcovj);
 2977:     cptcode=modmaxcovj;
 2978:     /* Ndum[0] = frequency of 0 for model-covariate j, Ndum[1] frequency of 1 etc. */
 2979:    /*for (i=0; i<=cptcode; i++) {*/
 2980:     for (i=modmincovj;  i<=modmaxcovj; i++) { /* i=-1 ? 0 and 1*//* For each value of the modality of model-cov j */
 2981:       printf("Frequencies of covariates %d V%d %d\n", j, Tvar[j], Ndum[i]);
 2982:       if( Ndum[i] != 0 ){ /* Counts if nobody answered, empty modality */
 2983: 	ncodemax[j]++;  /* ncodemax[j]= Number of non-null modalities of the j th covariate. */
 2984:       }
 2985:       /* In fact  ncodemax[j]=2 (dichotom. variables only) but it could be more for
 2986: 	 historical reasons: 3 if coded 1, 2, 3 and 4 and Ndum[2]=0 */
 2987:     } /* Ndum[-1] number of undefined modalities */
 2988: 
 2989:     /* j is a covariate, n=Tvar[j] of Vn; Fills nbcode */
 2990:     /* For covariate j, modalities could be 1, 2, 3, 4. If Ndum[2]=0 ncodemax[j] is not 4 but 3 */
 2991:     /* If Ndum[3}= 635; Ndum[4]=0; Ndum[5]=0; Ndum[6]=27; Ndum[7]=125;
 2992:        modmincovj=3; modmaxcovj = 7;
 2993:        There are only 3 modalities non empty (or 2 if 27 is too few) : ncodemax[j]=3;
 2994:        which will be coded 0, 1, 2 which in binary on 3-1 digits are 0=00 1=01, 2=10; defining two dummy 
 2995:        variables V1_1 and V1_2.
 2996:        nbcode[Tvar[j]][ij]=k;
 2997:        nbcode[Tvar[j]][1]=0;
 2998:        nbcode[Tvar[j]][2]=1;
 2999:        nbcode[Tvar[j]][3]=2;
 3000:     */
 3001:     ij=1; /* ij is similar to i but can jumps over null modalities */
 3002:     for (i=modmincovj; i<=modmaxcovj; i++) { /* i= 1 to 2 for dichotomous, or from 1 to 3 */
 3003:       for (k=0; k<= cptcode; k++) { /* k=-1 ? k=0 to 1 *//* Could be 1 to 4 */
 3004: 	/*recode from 0 */
 3005: 	if (Ndum[k] != 0) { /* If at least one individual responded to this modality k */
 3006: 	  nbcode[Tvar[j]][ij]=k;  /* stores the modality in an array nbcode. 
 3007: 				     k is a modality. If we have model=V1+V1*sex 
 3008: 				     then: nbcode[1][1]=0 ; nbcode[1][2]=1; nbcode[2][1]=0 ; nbcode[2][2]=1; */
 3009: 	  ij++;
 3010: 	}
 3011: 	if (ij > ncodemax[j]) break; 
 3012:       }  /* end of loop on */
 3013:     } /* end of loop on modality */ 
 3014:   } /* end of loop on model-covariate j. nbcode[Tvarj][1]=0 and nbcode[Tvarj][2]=1 sets the value of covariate j*/  
 3015:   
 3016:  for (k=-1; k< maxncov; k++) Ndum[k]=0; 
 3017:   
 3018:   for (i=1; i<=ncovmodel-2; i++) { /* -2, cste and age */ 
 3019:    /* Listing of all covariables in statement model to see if some covariates appear twice. For example, V1 appears twice in V1+V1*V2.*/ 
 3020:    ij=Tvar[i]; /* Tvar might be -1 if status was unknown */ 
 3021:    Ndum[ij]++; 
 3022:  } 
 3023: 
 3024:  ij=1;
 3025:  for (i=0; i<=  maxncov-1; i++) { /* modmaxcovj is unknown here. Only Ndum[2(V2),3(age*V3), 5(V3*V2) 6(V1*V4) */
 3026:    /*printf("Ndum[%d]=%d\n",i, Ndum[i]);*/
 3027:    if((Ndum[i]!=0) && (i<=ncovcol)){
 3028:      /*printf("diff Ndum[%d]=%d\n",i, Ndum[i]);*/
 3029:      Tvaraff[ij]=i; /*For printing (unclear) */
 3030:      ij++;
 3031:    }else
 3032:        Tvaraff[ij]=0;
 3033:  }
 3034:  ij--;
 3035:  cptcoveff=ij; /*Number of total covariates*/
 3036: 
 3037: }
 3038: 
 3039: 
 3040: /*********** Health Expectancies ****************/
 3041: 
 3042: void evsij(double ***eij, double x[], int nlstate, int stepm, int bage, int fage, double **oldm, double **savm, int cij, int estepm,char strstart[] )
 3043: 
 3044: {
 3045:   /* Health expectancies, no variances */
 3046:   int i, j, nhstepm, hstepm, h, nstepm;
 3047:   int nhstepma, nstepma; /* Decreasing with age */
 3048:   double age, agelim, hf;
 3049:   double ***p3mat;
 3050:   double eip;
 3051: 
 3052:   pstamp(ficreseij);
 3053:   fprintf(ficreseij,"# (a) Life expectancies by health status at initial age and (b) health expectancies by health status at initial age\n");
 3054:   fprintf(ficreseij,"# Age");
 3055:   for(i=1; i<=nlstate;i++){
 3056:     for(j=1; j<=nlstate;j++){
 3057:       fprintf(ficreseij," e%1d%1d ",i,j);
 3058:     }
 3059:     fprintf(ficreseij," e%1d. ",i);
 3060:   }
 3061:   fprintf(ficreseij,"\n");
 3062: 
 3063:   
 3064:   if(estepm < stepm){
 3065:     printf ("Problem %d lower than %d\n",estepm, stepm);
 3066:   }
 3067:   else  hstepm=estepm;   
 3068:   /* We compute the life expectancy from trapezoids spaced every estepm months
 3069:    * This is mainly to measure the difference between two models: for example
 3070:    * if stepm=24 months pijx are given only every 2 years and by summing them
 3071:    * we are calculating an estimate of the Life Expectancy assuming a linear 
 3072:    * progression in between and thus overestimating or underestimating according
 3073:    * to the curvature of the survival function. If, for the same date, we 
 3074:    * estimate the model with stepm=1 month, we can keep estepm to 24 months
 3075:    * to compare the new estimate of Life expectancy with the same linear 
 3076:    * hypothesis. A more precise result, taking into account a more precise
 3077:    * curvature will be obtained if estepm is as small as stepm. */
 3078: 
 3079:   /* For example we decided to compute the life expectancy with the smallest unit */
 3080:   /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm. 
 3081:      nhstepm is the number of hstepm from age to agelim 
 3082:      nstepm is the number of stepm from age to agelin. 
 3083:      Look at hpijx to understand the reason of that which relies in memory size
 3084:      and note for a fixed period like estepm months */
 3085:   /* We decided (b) to get a life expectancy respecting the most precise curvature of the
 3086:      survival function given by stepm (the optimization length). Unfortunately it
 3087:      means that if the survival funtion is printed only each two years of age and if
 3088:      you sum them up and add 1 year (area under the trapezoids) you won't get the same 
 3089:      results. So we changed our mind and took the option of the best precision.
 3090:   */
 3091:   hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */ 
 3092: 
 3093:   agelim=AGESUP;
 3094:   /* If stepm=6 months */
 3095:     /* Computed by stepm unit matrices, product of hstepm matrices, stored
 3096:        in an array of nhstepm length: nhstepm=10, hstepm=4, stepm=6 months */
 3097:     
 3098: /* nhstepm age range expressed in number of stepm */
 3099:   nstepm=(int) rint((agelim-bage)*YEARM/stepm); /* Biggest nstepm */
 3100:   /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3101:   /* if (stepm >= YEARM) hstepm=1;*/
 3102:   nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
 3103:   p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3104: 
 3105:   for (age=bage; age<=fage; age ++){ 
 3106:     nstepma=(int) rint((agelim-bage)*YEARM/stepm); /* Biggest nstepm */
 3107:     /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3108:     /* if (stepm >= YEARM) hstepm=1;*/
 3109:     nhstepma = nstepma/hstepm;/* Expressed in hstepm, typically nhstepma=40/4=10 */
 3110: 
 3111:     /* If stepm=6 months */
 3112:     /* Computed by stepm unit matrices, product of hstepma matrices, stored
 3113:        in an array of nhstepma length: nhstepma=10, hstepm=4, stepm=6 months */
 3114:     
 3115:     hpxij(p3mat,nhstepma,age,hstepm,x,nlstate,stepm,oldm, savm, cij);  
 3116:     
 3117:     hf=hstepm*stepm/YEARM;  /* Duration of hstepm expressed in year unit. */
 3118:     
 3119:     printf("%d|",(int)age);fflush(stdout);
 3120:     fprintf(ficlog,"%d|",(int)age);fflush(ficlog);
 3121:     
 3122:     /* Computing expectancies */
 3123:     for(i=1; i<=nlstate;i++)
 3124:       for(j=1; j<=nlstate;j++)
 3125: 	for (h=0, eij[i][j][(int)age]=0; h<=nhstepm-1; h++){
 3126: 	  eij[i][j][(int)age] += (p3mat[i][j][h]+p3mat[i][j][h+1])/2.0*hf;
 3127: 	  
 3128: 	  /* 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]);*/
 3129: 
 3130: 	}
 3131: 
 3132:     fprintf(ficreseij,"%3.0f",age );
 3133:     for(i=1; i<=nlstate;i++){
 3134:       eip=0;
 3135:       for(j=1; j<=nlstate;j++){
 3136: 	eip +=eij[i][j][(int)age];
 3137: 	fprintf(ficreseij,"%9.4f", eij[i][j][(int)age] );
 3138:       }
 3139:       fprintf(ficreseij,"%9.4f", eip );
 3140:     }
 3141:     fprintf(ficreseij,"\n");
 3142:     
 3143:   }
 3144:   free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3145:   printf("\n");
 3146:   fprintf(ficlog,"\n");
 3147:   
 3148: }
 3149: 
 3150: 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[] )
 3151: 
 3152: {
 3153:   /* Covariances of health expectancies eij and of total life expectancies according
 3154:    to initial status i, ei. .
 3155:   */
 3156:   int i, j, nhstepm, hstepm, h, nstepm, k, cptj, cptj2, i2, j2, ij, ji;
 3157:   int nhstepma, nstepma; /* Decreasing with age */
 3158:   double age, agelim, hf;
 3159:   double ***p3matp, ***p3matm, ***varhe;
 3160:   double **dnewm,**doldm;
 3161:   double *xp, *xm;
 3162:   double **gp, **gm;
 3163:   double ***gradg, ***trgradg;
 3164:   int theta;
 3165: 
 3166:   double eip, vip;
 3167: 
 3168:   varhe=ma3x(1,nlstate*nlstate,1,nlstate*nlstate,(int) bage, (int) fage);
 3169:   xp=vector(1,npar);
 3170:   xm=vector(1,npar);
 3171:   dnewm=matrix(1,nlstate*nlstate,1,npar);
 3172:   doldm=matrix(1,nlstate*nlstate,1,nlstate*nlstate);
 3173:   
 3174:   pstamp(ficresstdeij);
 3175:   fprintf(ficresstdeij,"# Health expectancies with standard errors\n");
 3176:   fprintf(ficresstdeij,"# Age");
 3177:   for(i=1; i<=nlstate;i++){
 3178:     for(j=1; j<=nlstate;j++)
 3179:       fprintf(ficresstdeij," e%1d%1d (SE)",i,j);
 3180:     fprintf(ficresstdeij," e%1d. ",i);
 3181:   }
 3182:   fprintf(ficresstdeij,"\n");
 3183: 
 3184:   pstamp(ficrescveij);
 3185:   fprintf(ficrescveij,"# Subdiagonal matrix of covariances of health expectancies by age: cov(eij,ekl)\n");
 3186:   fprintf(ficrescveij,"# Age");
 3187:   for(i=1; i<=nlstate;i++)
 3188:     for(j=1; j<=nlstate;j++){
 3189:       cptj= (j-1)*nlstate+i;
 3190:       for(i2=1; i2<=nlstate;i2++)
 3191: 	for(j2=1; j2<=nlstate;j2++){
 3192: 	  cptj2= (j2-1)*nlstate+i2;
 3193: 	  if(cptj2 <= cptj)
 3194: 	    fprintf(ficrescveij,"  %1d%1d,%1d%1d",i,j,i2,j2);
 3195: 	}
 3196:     }
 3197:   fprintf(ficrescveij,"\n");
 3198:   
 3199:   if(estepm < stepm){
 3200:     printf ("Problem %d lower than %d\n",estepm, stepm);
 3201:   }
 3202:   else  hstepm=estepm;   
 3203:   /* We compute the life expectancy from trapezoids spaced every estepm months
 3204:    * This is mainly to measure the difference between two models: for example
 3205:    * if stepm=24 months pijx are given only every 2 years and by summing them
 3206:    * we are calculating an estimate of the Life Expectancy assuming a linear 
 3207:    * progression in between and thus overestimating or underestimating according
 3208:    * to the curvature of the survival function. If, for the same date, we 
 3209:    * estimate the model with stepm=1 month, we can keep estepm to 24 months
 3210:    * to compare the new estimate of Life expectancy with the same linear 
 3211:    * hypothesis. A more precise result, taking into account a more precise
 3212:    * curvature will be obtained if estepm is as small as stepm. */
 3213: 
 3214:   /* For example we decided to compute the life expectancy with the smallest unit */
 3215:   /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm. 
 3216:      nhstepm is the number of hstepm from age to agelim 
 3217:      nstepm is the number of stepm from age to agelin. 
 3218:      Look at hpijx to understand the reason of that which relies in memory size
 3219:      and note for a fixed period like estepm months */
 3220:   /* We decided (b) to get a life expectancy respecting the most precise curvature of the
 3221:      survival function given by stepm (the optimization length). Unfortunately it
 3222:      means that if the survival funtion is printed only each two years of age and if
 3223:      you sum them up and add 1 year (area under the trapezoids) you won't get the same 
 3224:      results. So we changed our mind and took the option of the best precision.
 3225:   */
 3226:   hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */ 
 3227: 
 3228:   /* If stepm=6 months */
 3229:   /* nhstepm age range expressed in number of stepm */
 3230:   agelim=AGESUP;
 3231:   nstepm=(int) rint((agelim-bage)*YEARM/stepm); 
 3232:   /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3233:   /* if (stepm >= YEARM) hstepm=1;*/
 3234:   nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
 3235:   
 3236:   p3matp=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3237:   p3matm=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3238:   gradg=ma3x(0,nhstepm,1,npar,1,nlstate*nlstate);
 3239:   trgradg =ma3x(0,nhstepm,1,nlstate*nlstate,1,npar);
 3240:   gp=matrix(0,nhstepm,1,nlstate*nlstate);
 3241:   gm=matrix(0,nhstepm,1,nlstate*nlstate);
 3242: 
 3243:   for (age=bage; age<=fage; age ++){ 
 3244:     nstepma=(int) rint((agelim-bage)*YEARM/stepm); /* Biggest nstepm */
 3245:     /* Typically if 20 years nstepm = 20*12/6=40 stepm */ 
 3246:     /* if (stepm >= YEARM) hstepm=1;*/
 3247:     nhstepma = nstepma/hstepm;/* Expressed in hstepm, typically nhstepma=40/4=10 */
 3248: 
 3249:     /* If stepm=6 months */
 3250:     /* Computed by stepm unit matrices, product of hstepma matrices, stored
 3251:        in an array of nhstepma length: nhstepma=10, hstepm=4, stepm=6 months */
 3252:     
 3253:     hf=hstepm*stepm/YEARM;  /* Duration of hstepm expressed in year unit. */
 3254: 
 3255:     /* Computing  Variances of health expectancies */
 3256:     /* Gradient is computed with plus gp and minus gm. Code is duplicated in order to
 3257:        decrease memory allocation */
 3258:     for(theta=1; theta <=npar; theta++){
 3259:       for(i=1; i<=npar; i++){ 
 3260: 	xp[i] = x[i] + (i==theta ?delti[theta]:0);
 3261: 	xm[i] = x[i] - (i==theta ?delti[theta]:0);
 3262:       }
 3263:       hpxij(p3matp,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, cij);  
 3264:       hpxij(p3matm,nhstepm,age,hstepm,xm,nlstate,stepm,oldm,savm, cij);  
 3265:   
 3266:       for(j=1; j<= nlstate; j++){
 3267: 	for(i=1; i<=nlstate; i++){
 3268: 	  for(h=0; h<=nhstepm-1; h++){
 3269: 	    gp[h][(j-1)*nlstate + i] = (p3matp[i][j][h]+p3matp[i][j][h+1])/2.;
 3270: 	    gm[h][(j-1)*nlstate + i] = (p3matm[i][j][h]+p3matm[i][j][h+1])/2.;
 3271: 	  }
 3272: 	}
 3273:       }
 3274:      
 3275:       for(ij=1; ij<= nlstate*nlstate; ij++)
 3276: 	for(h=0; h<=nhstepm-1; h++){
 3277: 	  gradg[h][theta][ij]= (gp[h][ij]-gm[h][ij])/2./delti[theta];
 3278: 	}
 3279:     }/* End theta */
 3280:     
 3281:     
 3282:     for(h=0; h<=nhstepm-1; h++)
 3283:       for(j=1; j<=nlstate*nlstate;j++)
 3284: 	for(theta=1; theta <=npar; theta++)
 3285: 	  trgradg[h][j][theta]=gradg[h][theta][j];
 3286:     
 3287: 
 3288:      for(ij=1;ij<=nlstate*nlstate;ij++)
 3289:       for(ji=1;ji<=nlstate*nlstate;ji++)
 3290: 	varhe[ij][ji][(int)age] =0.;
 3291: 
 3292:      printf("%d|",(int)age);fflush(stdout);
 3293:      fprintf(ficlog,"%d|",(int)age);fflush(ficlog);
 3294:      for(h=0;h<=nhstepm-1;h++){
 3295:       for(k=0;k<=nhstepm-1;k++){
 3296: 	matprod2(dnewm,trgradg[h],1,nlstate*nlstate,1,npar,1,npar,matcov);
 3297: 	matprod2(doldm,dnewm,1,nlstate*nlstate,1,npar,1,nlstate*nlstate,gradg[k]);
 3298: 	for(ij=1;ij<=nlstate*nlstate;ij++)
 3299: 	  for(ji=1;ji<=nlstate*nlstate;ji++)
 3300: 	    varhe[ij][ji][(int)age] += doldm[ij][ji]*hf*hf;
 3301:       }
 3302:     }
 3303: 
 3304:     /* Computing expectancies */
 3305:     hpxij(p3matm,nhstepm,age,hstepm,x,nlstate,stepm,oldm, savm, cij);  
 3306:     for(i=1; i<=nlstate;i++)
 3307:       for(j=1; j<=nlstate;j++)
 3308: 	for (h=0, eij[i][j][(int)age]=0; h<=nhstepm-1; h++){
 3309: 	  eij[i][j][(int)age] += (p3matm[i][j][h]+p3matm[i][j][h+1])/2.0*hf;
 3310: 	  
 3311: 	  /* 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]);*/
 3312: 
 3313: 	}
 3314: 
 3315:     fprintf(ficresstdeij,"%3.0f",age );
 3316:     for(i=1; i<=nlstate;i++){
 3317:       eip=0.;
 3318:       vip=0.;
 3319:       for(j=1; j<=nlstate;j++){
 3320: 	eip += eij[i][j][(int)age];
 3321: 	for(k=1; k<=nlstate;k++) /* Sum on j and k of cov(eij,eik) */
 3322: 	  vip += varhe[(j-1)*nlstate+i][(k-1)*nlstate+i][(int)age];
 3323: 	fprintf(ficresstdeij," %9.4f (%.4f)", eij[i][j][(int)age], sqrt(varhe[(j-1)*nlstate+i][(j-1)*nlstate+i][(int)age]) );
 3324:       }
 3325:       fprintf(ficresstdeij," %9.4f (%.4f)", eip, sqrt(vip));
 3326:     }
 3327:     fprintf(ficresstdeij,"\n");
 3328: 
 3329:     fprintf(ficrescveij,"%3.0f",age );
 3330:     for(i=1; i<=nlstate;i++)
 3331:       for(j=1; j<=nlstate;j++){
 3332: 	cptj= (j-1)*nlstate+i;
 3333: 	for(i2=1; i2<=nlstate;i2++)
 3334: 	  for(j2=1; j2<=nlstate;j2++){
 3335: 	    cptj2= (j2-1)*nlstate+i2;
 3336: 	    if(cptj2 <= cptj)
 3337: 	      fprintf(ficrescveij," %.4f", varhe[cptj][cptj2][(int)age]);
 3338: 	  }
 3339:       }
 3340:     fprintf(ficrescveij,"\n");
 3341:    
 3342:   }
 3343:   free_matrix(gm,0,nhstepm,1,nlstate*nlstate);
 3344:   free_matrix(gp,0,nhstepm,1,nlstate*nlstate);
 3345:   free_ma3x(gradg,0,nhstepm,1,npar,1,nlstate*nlstate);
 3346:   free_ma3x(trgradg,0,nhstepm,1,nlstate*nlstate,1,npar);
 3347:   free_ma3x(p3matm,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3348:   free_ma3x(p3matp,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3349:   printf("\n");
 3350:   fprintf(ficlog,"\n");
 3351: 
 3352:   free_vector(xm,1,npar);
 3353:   free_vector(xp,1,npar);
 3354:   free_matrix(dnewm,1,nlstate*nlstate,1,npar);
 3355:   free_matrix(doldm,1,nlstate*nlstate,1,nlstate*nlstate);
 3356:   free_ma3x(varhe,1,nlstate*nlstate,1,nlstate*nlstate,(int) bage, (int)fage);
 3357: }
 3358: 
 3359: /************ Variance ******************/
 3360: 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[])
 3361: {
 3362:   /* Variance of health expectancies */
 3363:   /*  double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double ** savm,double ftolpl);*/
 3364:   /* double **newm;*/
 3365:   /* int movingaverage(double ***probs, double bage,double fage, double ***mobaverage, int mobilav)*/
 3366:   
 3367:   int movingaverage();
 3368:   double **dnewm,**doldm;
 3369:   double **dnewmp,**doldmp;
 3370:   int i, j, nhstepm, hstepm, h, nstepm ;
 3371:   int k;
 3372:   double *xp;
 3373:   double **gp, **gm;  /* for var eij */
 3374:   double ***gradg, ***trgradg; /*for var eij */
 3375:   double **gradgp, **trgradgp; /* for var p point j */
 3376:   double *gpp, *gmp; /* for var p point j */
 3377:   double **varppt; /* for var p point j nlstate to nlstate+ndeath */
 3378:   double ***p3mat;
 3379:   double age,agelim, hf;
 3380:   double ***mobaverage;
 3381:   int theta;
 3382:   char digit[4];
 3383:   char digitp[25];
 3384: 
 3385:   char fileresprobmorprev[FILENAMELENGTH];
 3386: 
 3387:   if(popbased==1){
 3388:     if(mobilav!=0)
 3389:       strcpy(digitp,"-populbased-mobilav-");
 3390:     else strcpy(digitp,"-populbased-nomobil-");
 3391:   }
 3392:   else 
 3393:     strcpy(digitp,"-stablbased-");
 3394: 
 3395:   if (mobilav!=0) {
 3396:     mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 3397:     if (movingaverage(probs, bage, fage, mobaverage,mobilav)!=0){
 3398:       fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 3399:       printf(" Error in movingaverage mobilav=%d\n",mobilav);
 3400:     }
 3401:   }
 3402: 
 3403:   strcpy(fileresprobmorprev,"prmorprev"); 
 3404:   sprintf(digit,"%-d",ij);
 3405:   /*printf("DIGIT=%s, ij=%d ijr=%-d|\n",digit, ij,ij);*/
 3406:   strcat(fileresprobmorprev,digit); /* Tvar to be done */
 3407:   strcat(fileresprobmorprev,digitp); /* Popbased or not, mobilav or not */
 3408:   strcat(fileresprobmorprev,fileres);
 3409:   if((ficresprobmorprev=fopen(fileresprobmorprev,"w"))==NULL) {
 3410:     printf("Problem with resultfile: %s\n", fileresprobmorprev);
 3411:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprobmorprev);
 3412:   }
 3413:   printf("Computing total mortality p.j=w1*p1j+w2*p2j+..: result on file '%s' \n",fileresprobmorprev);
 3414:  
 3415:   fprintf(ficlog,"Computing total mortality p.j=w1*p1j+w2*p2j+..: result on file '%s' \n",fileresprobmorprev);
 3416:   pstamp(ficresprobmorprev);
 3417:   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);
 3418:   fprintf(ficresprobmorprev,"# Age cov=%-d",ij);
 3419:   for(j=nlstate+1; j<=(nlstate+ndeath);j++){
 3420:     fprintf(ficresprobmorprev," p.%-d SE",j);
 3421:     for(i=1; i<=nlstate;i++)
 3422:       fprintf(ficresprobmorprev," w%1d p%-d%-d",i,i,j);
 3423:   }  
 3424:   fprintf(ficresprobmorprev,"\n");
 3425:   fprintf(ficgp,"\n# Routine varevsij");
 3426:   /* fprintf(fichtm, "#Local time at start: %s", strstart);*/
 3427:   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");
 3428:   fprintf(fichtm,"\n<br>%s  <br>\n",digitp);
 3429: /*   } */
 3430:   varppt = matrix(nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3431:   pstamp(ficresvij);
 3432:   fprintf(ficresvij,"# Variance and covariance of health expectancies e.j \n#  (weighted average of eij where weights are ");
 3433:   if(popbased==1)
 3434:     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);
 3435:   else
 3436:     fprintf(ficresvij,"the age specific period (stable) prevalences in each health state \n");
 3437:   fprintf(ficresvij,"# Age");
 3438:   for(i=1; i<=nlstate;i++)
 3439:     for(j=1; j<=nlstate;j++)
 3440:       fprintf(ficresvij," Cov(e.%1d, e.%1d)",i,j);
 3441:   fprintf(ficresvij,"\n");
 3442: 
 3443:   xp=vector(1,npar);
 3444:   dnewm=matrix(1,nlstate,1,npar);
 3445:   doldm=matrix(1,nlstate,1,nlstate);
 3446:   dnewmp= matrix(nlstate+1,nlstate+ndeath,1,npar);
 3447:   doldmp= matrix(nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3448: 
 3449:   gradgp=matrix(1,npar,nlstate+1,nlstate+ndeath);
 3450:   gpp=vector(nlstate+1,nlstate+ndeath);
 3451:   gmp=vector(nlstate+1,nlstate+ndeath);
 3452:   trgradgp =matrix(nlstate+1,nlstate+ndeath,1,npar); /* mu or p point j*/
 3453:   
 3454:   if(estepm < stepm){
 3455:     printf ("Problem %d lower than %d\n",estepm, stepm);
 3456:   }
 3457:   else  hstepm=estepm;   
 3458:   /* For example we decided to compute the life expectancy with the smallest unit */
 3459:   /* hstepm beeing the number of stepms, if hstepm=1 the length of hstepm is stepm. 
 3460:      nhstepm is the number of hstepm from age to agelim 
 3461:      nstepm is the number of stepm from age to agelin. 
 3462:      Look at function hpijx to understand why (it is linked to memory size questions) */
 3463:   /* We decided (b) to get a life expectancy respecting the most precise curvature of the
 3464:      survival function given by stepm (the optimization length). Unfortunately it
 3465:      means that if the survival funtion is printed every two years of age and if
 3466:      you sum them up and add 1 year (area under the trapezoids) you won't get the same 
 3467:      results. So we changed our mind and took the option of the best precision.
 3468:   */
 3469:   hstepm=hstepm/stepm; /* Typically in stepm units, if stepm=6 & estepm=24 , = 24/6 months = 4 */ 
 3470:   agelim = AGESUP;
 3471:   for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
 3472:     nstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */ 
 3473:     nhstepm = nstepm/hstepm;/* Expressed in hstepm, typically nhstepm=40/4=10 */
 3474:     p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3475:     gradg=ma3x(0,nhstepm,1,npar,1,nlstate);
 3476:     gp=matrix(0,nhstepm,1,nlstate);
 3477:     gm=matrix(0,nhstepm,1,nlstate);
 3478: 
 3479: 
 3480:     for(theta=1; theta <=npar; theta++){
 3481:       for(i=1; i<=npar; i++){ /* Computes gradient x + delta*/
 3482: 	xp[i] = x[i] + (i==theta ?delti[theta]:0);
 3483:       }
 3484:       hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);  
 3485:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3486: 
 3487:       if (popbased==1) {
 3488: 	if(mobilav ==0){
 3489: 	  for(i=1; i<=nlstate;i++)
 3490: 	    prlim[i][i]=probs[(int)age][i][ij];
 3491: 	}else{ /* mobilav */ 
 3492: 	  for(i=1; i<=nlstate;i++)
 3493: 	    prlim[i][i]=mobaverage[(int)age][i][ij];
 3494: 	}
 3495:       }
 3496:   
 3497:       for(j=1; j<= nlstate; j++){
 3498: 	for(h=0; h<=nhstepm; h++){
 3499: 	  for(i=1, gp[h][j]=0.;i<=nlstate;i++)
 3500: 	    gp[h][j] += prlim[i][i]*p3mat[i][j][h];
 3501: 	}
 3502:       }
 3503:       /* This for computing probability of death (h=1 means
 3504:          computed over hstepm matrices product = hstepm*stepm months) 
 3505:          as a weighted average of prlim.
 3506:       */
 3507:       for(j=nlstate+1;j<=nlstate+ndeath;j++){
 3508: 	for(i=1,gpp[j]=0.; i<= nlstate; i++)
 3509: 	  gpp[j] += prlim[i][i]*p3mat[i][j][1];
 3510:       }    
 3511:       /* end probability of death */
 3512: 
 3513:       for(i=1; i<=npar; i++) /* Computes gradient x - delta */
 3514: 	xp[i] = x[i] - (i==theta ?delti[theta]:0);
 3515:       hpxij(p3mat,nhstepm,age,hstepm,xp,nlstate,stepm,oldm,savm, ij);  
 3516:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3517:  
 3518:       if (popbased==1) {
 3519: 	if(mobilav ==0){
 3520: 	  for(i=1; i<=nlstate;i++)
 3521: 	    prlim[i][i]=probs[(int)age][i][ij];
 3522: 	}else{ /* mobilav */ 
 3523: 	  for(i=1; i<=nlstate;i++)
 3524: 	    prlim[i][i]=mobaverage[(int)age][i][ij];
 3525: 	}
 3526:       }
 3527: 
 3528:       for(j=1; j<= nlstate; j++){  /* Sum of wi * eij = e.j */
 3529: 	for(h=0; h<=nhstepm; h++){
 3530: 	  for(i=1, gm[h][j]=0.;i<=nlstate;i++)
 3531: 	    gm[h][j] += prlim[i][i]*p3mat[i][j][h];
 3532: 	}
 3533:       }
 3534:       /* This for computing probability of death (h=1 means
 3535:          computed over hstepm matrices product = hstepm*stepm months) 
 3536:          as a weighted average of prlim.
 3537:       */
 3538:       for(j=nlstate+1;j<=nlstate+ndeath;j++){
 3539: 	for(i=1,gmp[j]=0.; i<= nlstate; i++)
 3540:          gmp[j] += prlim[i][i]*p3mat[i][j][1];
 3541:       }    
 3542:       /* end probability of death */
 3543: 
 3544:       for(j=1; j<= nlstate; j++) /* vareij */
 3545: 	for(h=0; h<=nhstepm; h++){
 3546: 	  gradg[h][theta][j]= (gp[h][j]-gm[h][j])/2./delti[theta];
 3547: 	}
 3548: 
 3549:       for(j=nlstate+1; j<= nlstate+ndeath; j++){ /* var mu */
 3550: 	gradgp[theta][j]= (gpp[j]-gmp[j])/2./delti[theta];
 3551:       }
 3552: 
 3553:     } /* End theta */
 3554: 
 3555:     trgradg =ma3x(0,nhstepm,1,nlstate,1,npar); /* veij */
 3556: 
 3557:     for(h=0; h<=nhstepm; h++) /* veij */
 3558:       for(j=1; j<=nlstate;j++)
 3559: 	for(theta=1; theta <=npar; theta++)
 3560: 	  trgradg[h][j][theta]=gradg[h][theta][j];
 3561: 
 3562:     for(j=nlstate+1; j<=nlstate+ndeath;j++) /* mu */
 3563:       for(theta=1; theta <=npar; theta++)
 3564: 	trgradgp[j][theta]=gradgp[theta][j];
 3565:   
 3566: 
 3567:     hf=hstepm*stepm/YEARM;  /* Duration of hstepm expressed in year unit. */
 3568:     for(i=1;i<=nlstate;i++)
 3569:       for(j=1;j<=nlstate;j++)
 3570: 	vareij[i][j][(int)age] =0.;
 3571: 
 3572:     for(h=0;h<=nhstepm;h++){
 3573:       for(k=0;k<=nhstepm;k++){
 3574: 	matprod2(dnewm,trgradg[h],1,nlstate,1,npar,1,npar,matcov);
 3575: 	matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg[k]);
 3576: 	for(i=1;i<=nlstate;i++)
 3577: 	  for(j=1;j<=nlstate;j++)
 3578: 	    vareij[i][j][(int)age] += doldm[i][j]*hf*hf;
 3579:       }
 3580:     }
 3581:   
 3582:     /* pptj */
 3583:     matprod2(dnewmp,trgradgp,nlstate+1,nlstate+ndeath,1,npar,1,npar,matcov);
 3584:     matprod2(doldmp,dnewmp,nlstate+1,nlstate+ndeath,1,npar,nlstate+1,nlstate+ndeath,gradgp);
 3585:     for(j=nlstate+1;j<=nlstate+ndeath;j++)
 3586:       for(i=nlstate+1;i<=nlstate+ndeath;i++)
 3587: 	varppt[j][i]=doldmp[j][i];
 3588:     /* end ppptj */
 3589:     /*  x centered again */
 3590:     hpxij(p3mat,nhstepm,age,hstepm,x,nlstate,stepm,oldm,savm, ij);  
 3591:     prevalim(prlim,nlstate,x,age,oldm,savm,ftolpl,ij);
 3592:  
 3593:     if (popbased==1) {
 3594:       if(mobilav ==0){
 3595: 	for(i=1; i<=nlstate;i++)
 3596: 	  prlim[i][i]=probs[(int)age][i][ij];
 3597:       }else{ /* mobilav */ 
 3598: 	for(i=1; i<=nlstate;i++)
 3599: 	  prlim[i][i]=mobaverage[(int)age][i][ij];
 3600:       }
 3601:     }
 3602:              
 3603:     /* This for computing probability of death (h=1 means
 3604:        computed over hstepm (estepm) matrices product = hstepm*stepm months) 
 3605:        as a weighted average of prlim.
 3606:     */
 3607:     for(j=nlstate+1;j<=nlstate+ndeath;j++){
 3608:       for(i=1,gmp[j]=0.;i<= nlstate; i++) 
 3609: 	gmp[j] += prlim[i][i]*p3mat[i][j][1]; 
 3610:     }    
 3611:     /* end probability of death */
 3612: 
 3613:     fprintf(ficresprobmorprev,"%3d %d ",(int) age, ij);
 3614:     for(j=nlstate+1; j<=(nlstate+ndeath);j++){
 3615:       fprintf(ficresprobmorprev," %11.3e %11.3e",gmp[j], sqrt(varppt[j][j]));
 3616:       for(i=1; i<=nlstate;i++){
 3617: 	fprintf(ficresprobmorprev," %11.3e %11.3e ",prlim[i][i],p3mat[i][j][1]);
 3618:       }
 3619:     } 
 3620:     fprintf(ficresprobmorprev,"\n");
 3621: 
 3622:     fprintf(ficresvij,"%.0f ",age );
 3623:     for(i=1; i<=nlstate;i++)
 3624:       for(j=1; j<=nlstate;j++){
 3625: 	fprintf(ficresvij," %.4f", vareij[i][j][(int)age]);
 3626:       }
 3627:     fprintf(ficresvij,"\n");
 3628:     free_matrix(gp,0,nhstepm,1,nlstate);
 3629:     free_matrix(gm,0,nhstepm,1,nlstate);
 3630:     free_ma3x(gradg,0,nhstepm,1,npar,1,nlstate);
 3631:     free_ma3x(trgradg,0,nhstepm,1,nlstate,1,npar);
 3632:     free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 3633:   } /* End age */
 3634:   free_vector(gpp,nlstate+1,nlstate+ndeath);
 3635:   free_vector(gmp,nlstate+1,nlstate+ndeath);
 3636:   free_matrix(gradgp,1,npar,nlstate+1,nlstate+ndeath);
 3637:   free_matrix(trgradgp,nlstate+1,nlstate+ndeath,1,npar); /* mu or p point j*/
 3638:   fprintf(ficgp,"\nunset parametric;unset label; set ter png small size 320, 240");
 3639:   /* for(j=nlstate+1; j<= nlstate+ndeath; j++){ *//* Only the first actually */
 3640:   fprintf(ficgp,"\n set log y; unset log x;set xlabel \"Age\"; set ylabel \"Force of mortality (year-1)\";");
 3641: /*   fprintf(ficgp,"\n plot \"%s\"  u 1:($3*%6.3f) not w l 1 ",fileresprobmorprev,YEARM/estepm); */
 3642: /*   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3+1.96*$4)*%6.3f) t \"95\%% interval\" w l 2 ",fileresprobmorprev,YEARM/estepm); */
 3643: /*   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3-1.96*$4)*%6.3f) not w l 2 ",fileresprobmorprev,YEARM/estepm); */
 3644:   fprintf(ficgp,"\n plot \"%s\"  u 1:($3) not w l lt 1 ",subdirf(fileresprobmorprev));
 3645:   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3+1.96*$4)) t \"95%% interval\" w l lt 2 ",subdirf(fileresprobmorprev));
 3646:   fprintf(ficgp,"\n replot \"%s\"  u 1:(($3-1.96*$4)) not w l lt 2 ",subdirf(fileresprobmorprev));
 3647:   fprintf(fichtm,"\n<br> File (multiple files are possible if covariates are present): <A href=\"%s\">%s</a>\n",subdirf(fileresprobmorprev),subdirf(fileresprobmorprev));
 3648:   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);
 3649:   /*  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);
 3650: */
 3651: /*   fprintf(ficgp,"\nset out \"varmuptjgr%s%s%s.png\";replot;",digitp,optionfilefiname,digit); */
 3652:   fprintf(ficgp,"\nset out \"%s%s.png\";replot;\n",subdirf3(optionfilefiname,"varmuptjgr",digitp),digit);
 3653: 
 3654:   free_vector(xp,1,npar);
 3655:   free_matrix(doldm,1,nlstate,1,nlstate);
 3656:   free_matrix(dnewm,1,nlstate,1,npar);
 3657:   free_matrix(doldmp,nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3658:   free_matrix(dnewmp,nlstate+1,nlstate+ndeath,1,npar);
 3659:   free_matrix(varppt,nlstate+1,nlstate+ndeath,nlstate+1,nlstate+ndeath);
 3660:   if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 3661:   fclose(ficresprobmorprev);
 3662:   fflush(ficgp);
 3663:   fflush(fichtm); 
 3664: }  /* end varevsij */
 3665: 
 3666: /************ Variance of prevlim ******************/
 3667: 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[])
 3668: {
 3669:   /* Variance of prevalence limit */
 3670:   /*  double **prevalim(double **prlim, int nlstate, double *xp, double age, double **oldm, double **savm,double ftolpl);*/
 3671: 
 3672:   double **dnewm,**doldm;
 3673:   int i, j, nhstepm, hstepm;
 3674:   double *xp;
 3675:   double *gp, *gm;
 3676:   double **gradg, **trgradg;
 3677:   double age,agelim;
 3678:   int theta;
 3679:   
 3680:   pstamp(ficresvpl);
 3681:   fprintf(ficresvpl,"# Standard deviation of period (stable) prevalences \n");
 3682:   fprintf(ficresvpl,"# Age");
 3683:   for(i=1; i<=nlstate;i++)
 3684:       fprintf(ficresvpl," %1d-%1d",i,i);
 3685:   fprintf(ficresvpl,"\n");
 3686: 
 3687:   xp=vector(1,npar);
 3688:   dnewm=matrix(1,nlstate,1,npar);
 3689:   doldm=matrix(1,nlstate,1,nlstate);
 3690:   
 3691:   hstepm=1*YEARM; /* Every year of age */
 3692:   hstepm=hstepm/stepm; /* Typically in stepm units, if j= 2 years, = 2/6 months = 4 */ 
 3693:   agelim = AGESUP;
 3694:   for (age=bage; age<=fage; age ++){ /* If stepm=6 months */
 3695:     nhstepm=(int) rint((agelim-age)*YEARM/stepm); /* Typically 20 years = 20*12/6=40 */ 
 3696:     if (stepm >= YEARM) hstepm=1;
 3697:     nhstepm = nhstepm/hstepm; /* Typically 40/4=10 */
 3698:     gradg=matrix(1,npar,1,nlstate);
 3699:     gp=vector(1,nlstate);
 3700:     gm=vector(1,nlstate);
 3701: 
 3702:     for(theta=1; theta <=npar; theta++){
 3703:       for(i=1; i<=npar; i++){ /* Computes gradient */
 3704: 	xp[i] = x[i] + (i==theta ?delti[theta]:0);
 3705:       }
 3706:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3707:       for(i=1;i<=nlstate;i++)
 3708: 	gp[i] = prlim[i][i];
 3709:     
 3710:       for(i=1; i<=npar; i++) /* Computes gradient */
 3711: 	xp[i] = x[i] - (i==theta ?delti[theta]:0);
 3712:       prevalim(prlim,nlstate,xp,age,oldm,savm,ftolpl,ij);
 3713:       for(i=1;i<=nlstate;i++)
 3714: 	gm[i] = prlim[i][i];
 3715: 
 3716:       for(i=1;i<=nlstate;i++)
 3717: 	gradg[theta][i]= (gp[i]-gm[i])/2./delti[theta];
 3718:     } /* End theta */
 3719: 
 3720:     trgradg =matrix(1,nlstate,1,npar);
 3721: 
 3722:     for(j=1; j<=nlstate;j++)
 3723:       for(theta=1; theta <=npar; theta++)
 3724: 	trgradg[j][theta]=gradg[theta][j];
 3725: 
 3726:     for(i=1;i<=nlstate;i++)
 3727:       varpl[i][(int)age] =0.;
 3728:     matprod2(dnewm,trgradg,1,nlstate,1,npar,1,npar,matcov);
 3729:     matprod2(doldm,dnewm,1,nlstate,1,npar,1,nlstate,gradg);
 3730:     for(i=1;i<=nlstate;i++)
 3731:       varpl[i][(int)age] = doldm[i][i]; /* Covariances are useless */
 3732: 
 3733:     fprintf(ficresvpl,"%.0f ",age );
 3734:     for(i=1; i<=nlstate;i++)
 3735:       fprintf(ficresvpl," %.5f (%.5f)",prlim[i][i],sqrt(varpl[i][(int)age]));
 3736:     fprintf(ficresvpl,"\n");
 3737:     free_vector(gp,1,nlstate);
 3738:     free_vector(gm,1,nlstate);
 3739:     free_matrix(gradg,1,npar,1,nlstate);
 3740:     free_matrix(trgradg,1,nlstate,1,npar);
 3741:   } /* End age */
 3742: 
 3743:   free_vector(xp,1,npar);
 3744:   free_matrix(doldm,1,nlstate,1,npar);
 3745:   free_matrix(dnewm,1,nlstate,1,nlstate);
 3746: 
 3747: }
 3748: 
 3749: /************ Variance of one-step probabilities  ******************/
 3750: 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[])
 3751: {
 3752:   int i, j=0,  k1, l1, tj;
 3753:   int k2, l2, j1,  z1;
 3754:   int k=0, l;
 3755:   int first=1, first1, first2;
 3756:   double cv12, mu1, mu2, lc1, lc2, v12, v21, v11, v22,v1,v2, c12, tnalp;
 3757:   double **dnewm,**doldm;
 3758:   double *xp;
 3759:   double *gp, *gm;
 3760:   double **gradg, **trgradg;
 3761:   double **mu;
 3762:   double age, cov[NCOVMAX+1];
 3763:   double std=2.0; /* Number of standard deviation wide of confidence ellipsoids */
 3764:   int theta;
 3765:   char fileresprob[FILENAMELENGTH];
 3766:   char fileresprobcov[FILENAMELENGTH];
 3767:   char fileresprobcor[FILENAMELENGTH];
 3768:   double ***varpij;
 3769: 
 3770:   strcpy(fileresprob,"prob"); 
 3771:   strcat(fileresprob,fileres);
 3772:   if((ficresprob=fopen(fileresprob,"w"))==NULL) {
 3773:     printf("Problem with resultfile: %s\n", fileresprob);
 3774:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprob);
 3775:   }
 3776:   strcpy(fileresprobcov,"probcov"); 
 3777:   strcat(fileresprobcov,fileres);
 3778:   if((ficresprobcov=fopen(fileresprobcov,"w"))==NULL) {
 3779:     printf("Problem with resultfile: %s\n", fileresprobcov);
 3780:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprobcov);
 3781:   }
 3782:   strcpy(fileresprobcor,"probcor"); 
 3783:   strcat(fileresprobcor,fileres);
 3784:   if((ficresprobcor=fopen(fileresprobcor,"w"))==NULL) {
 3785:     printf("Problem with resultfile: %s\n", fileresprobcor);
 3786:     fprintf(ficlog,"Problem with resultfile: %s\n", fileresprobcor);
 3787:   }
 3788:   printf("Computing standard deviation of one-step probabilities: result on file '%s' \n",fileresprob);
 3789:   fprintf(ficlog,"Computing standard deviation of one-step probabilities: result on file '%s' \n",fileresprob);
 3790:   printf("Computing matrix of variance covariance of one-step probabilities: result on file '%s' \n",fileresprobcov);
 3791:   fprintf(ficlog,"Computing matrix of variance covariance of one-step probabilities: result on file '%s' \n",fileresprobcov);
 3792:   printf("and correlation matrix of one-step probabilities: result on file '%s' \n",fileresprobcor);
 3793:   fprintf(ficlog,"and correlation matrix of one-step probabilities: result on file '%s' \n",fileresprobcor);
 3794:   pstamp(ficresprob);
 3795:   fprintf(ficresprob,"#One-step probabilities and stand. devi in ()\n");
 3796:   fprintf(ficresprob,"# Age");
 3797:   pstamp(ficresprobcov);
 3798:   fprintf(ficresprobcov,"#One-step probabilities and covariance matrix\n");
 3799:   fprintf(ficresprobcov,"# Age");
 3800:   pstamp(ficresprobcor);
 3801:   fprintf(ficresprobcor,"#One-step probabilities and correlation matrix\n");
 3802:   fprintf(ficresprobcor,"# Age");
 3803: 
 3804: 
 3805:   for(i=1; i<=nlstate;i++)
 3806:     for(j=1; j<=(nlstate+ndeath);j++){
 3807:       fprintf(ficresprob," p%1d-%1d (SE)",i,j);
 3808:       fprintf(ficresprobcov," p%1d-%1d ",i,j);
 3809:       fprintf(ficresprobcor," p%1d-%1d ",i,j);
 3810:     }  
 3811:  /* fprintf(ficresprob,"\n");
 3812:   fprintf(ficresprobcov,"\n");
 3813:   fprintf(ficresprobcor,"\n");
 3814:  */
 3815:   xp=vector(1,npar);
 3816:   dnewm=matrix(1,(nlstate)*(nlstate+ndeath),1,npar);
 3817:   doldm=matrix(1,(nlstate)*(nlstate+ndeath),1,(nlstate)*(nlstate+ndeath));
 3818:   mu=matrix(1,(nlstate)*(nlstate+ndeath), (int) bage, (int)fage);
 3819:   varpij=ma3x(1,nlstate*(nlstate+ndeath),1,nlstate*(nlstate+ndeath),(int) bage, (int) fage);
 3820:   first=1;
 3821:   fprintf(ficgp,"\n# Routine varprob");
 3822:   fprintf(fichtm,"\n<li><h4> Computing and drawing one step probabilities with their confidence intervals</h4></li>\n");
 3823:   fprintf(fichtm,"\n");
 3824: 
 3825:   fprintf(fichtm,"\n<li><h4> <a href=\"%s\">Matrix of variance-covariance of pairs of step probabilities (drawings)</a></h4></li>\n",optionfilehtmcov);
 3826:   fprintf(fichtmcov,"\n<h4>Matrix of variance-covariance of pairs of step probabilities</h4>\n\
 3827:   file %s<br>\n",optionfilehtmcov);
 3828:   fprintf(fichtmcov,"\nEllipsoids of confidence centered on point (p<inf>ij</inf>, p<inf>kl</inf>) are estimated\
 3829: and drawn. It helps understanding how is the covariance between two incidences.\
 3830:  They are expressed in year<sup>-1</sup> in order to be less dependent of stepm.<br>\n");
 3831:   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. \
 3832: It can be understood this way: if pij and pkl where uncorrelated the (2x2) matrix of covariance \
 3833: would have been (1/(var pij), 0 , 0, 1/(var pkl)), and the confidence interval would be 2 \
 3834: standard deviations wide on each axis. <br>\
 3835:  Now, if both incidences are correlated (usual case) we diagonalised the inverse of the covariance matrix\
 3836:  and made the appropriate rotation to look at the uncorrelated principal directions.<br>\
 3837: To be simple, these graphs help to understand the significativity of each parameter in relation to a second other one.<br> \n");
 3838: 
 3839:   cov[1]=1;
 3840:   /* tj=cptcoveff; */
 3841:   tj = (int) pow(2,cptcoveff);
 3842:   if (cptcovn<1) {tj=1;ncodemax[1]=1;}
 3843:   j1=0;
 3844:   for(j1=1; j1<=tj;j1++){
 3845:     /*for(i1=1; i1<=ncodemax[t];i1++){ */
 3846:     /*j1++;*/
 3847:       if  (cptcovn>0) {
 3848: 	fprintf(ficresprob, "\n#********** Variable "); 
 3849: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprob, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3850: 	fprintf(ficresprob, "**********\n#\n");
 3851: 	fprintf(ficresprobcov, "\n#********** Variable "); 
 3852: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprobcov, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3853: 	fprintf(ficresprobcov, "**********\n#\n");
 3854: 	
 3855: 	fprintf(ficgp, "\n#********** Variable "); 
 3856: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficgp, " V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3857: 	fprintf(ficgp, "**********\n#\n");
 3858: 	
 3859: 	
 3860: 	fprintf(fichtmcov, "\n<hr  size=\"2\" color=\"#EC5E5E\">********** Variable "); 
 3861: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(fichtm, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3862: 	fprintf(fichtmcov, "**********\n<hr size=\"2\" color=\"#EC5E5E\">");
 3863: 	
 3864: 	fprintf(ficresprobcor, "\n#********** Variable ");    
 3865: 	for (z1=1; z1<=cptcoveff; z1++) fprintf(ficresprobcor, "V%d=%d ",Tvaraff[z1],nbcode[Tvaraff[z1]][codtab[j1][z1]]);
 3866: 	fprintf(ficresprobcor, "**********\n#");    
 3867:       }
 3868:       
 3869:       gradg=matrix(1,npar,1,(nlstate)*(nlstate+ndeath));
 3870:       trgradg=matrix(1,(nlstate)*(nlstate+ndeath),1,npar);
 3871:       gp=vector(1,(nlstate)*(nlstate+ndeath));
 3872:       gm=vector(1,(nlstate)*(nlstate+ndeath));
 3873:       for (age=bage; age<=fage; age ++){ 
 3874: 	cov[2]=age;
 3875: 	for (k=1; k<=cptcovn;k++) {
 3876: 	  cov[2+k]=nbcode[Tvar[k]][codtab[j1][Tvar[k]]];/* j1 1 2 3 4
 3877: 							 * 1  1 1 1 1
 3878: 							 * 2  2 1 1 1
 3879: 							 * 3  1 2 1 1
 3880: 							 */
 3881: 	  /* nbcode[1][1]=0 nbcode[1][2]=1;*/
 3882: 	}
 3883: 	for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2];
 3884: 	for (k=1; k<=cptcovprod;k++)
 3885: 	  cov[2+Tprod[k]]=nbcode[Tvard[k][1]][codtab[ij][Tvard[k][1]]]*nbcode[Tvard[k][2]][codtab[ij][Tvard[k][2]]];
 3886: 	
 3887:     
 3888: 	for(theta=1; theta <=npar; theta++){
 3889: 	  for(i=1; i<=npar; i++)
 3890: 	    xp[i] = x[i] + (i==theta ?delti[theta]:(double)0);
 3891: 	  
 3892: 	  pmij(pmmij,cov,ncovmodel,xp,nlstate);
 3893: 	  
 3894: 	  k=0;
 3895: 	  for(i=1; i<= (nlstate); i++){
 3896: 	    for(j=1; j<=(nlstate+ndeath);j++){
 3897: 	      k=k+1;
 3898: 	      gp[k]=pmmij[i][j];
 3899: 	    }
 3900: 	  }
 3901: 	  
 3902: 	  for(i=1; i<=npar; i++)
 3903: 	    xp[i] = x[i] - (i==theta ?delti[theta]:(double)0);
 3904:     
 3905: 	  pmij(pmmij,cov,ncovmodel,xp,nlstate);
 3906: 	  k=0;
 3907: 	  for(i=1; i<=(nlstate); i++){
 3908: 	    for(j=1; j<=(nlstate+ndeath);j++){
 3909: 	      k=k+1;
 3910: 	      gm[k]=pmmij[i][j];
 3911: 	    }
 3912: 	  }
 3913:      
 3914: 	  for(i=1; i<= (nlstate)*(nlstate+ndeath); i++) 
 3915: 	    gradg[theta][i]=(gp[i]-gm[i])/(double)2./delti[theta];  
 3916: 	}
 3917: 
 3918: 	for(j=1; j<=(nlstate)*(nlstate+ndeath);j++)
 3919: 	  for(theta=1; theta <=npar; theta++)
 3920: 	    trgradg[j][theta]=gradg[theta][j];
 3921: 	
 3922: 	matprod2(dnewm,trgradg,1,(nlstate)*(nlstate+ndeath),1,npar,1,npar,matcov); 
 3923: 	matprod2(doldm,dnewm,1,(nlstate)*(nlstate+ndeath),1,npar,1,(nlstate)*(nlstate+ndeath),gradg);
 3924: 
 3925: 	pmij(pmmij,cov,ncovmodel,x,nlstate);
 3926: 	
 3927: 	k=0;
 3928: 	for(i=1; i<=(nlstate); i++){
 3929: 	  for(j=1; j<=(nlstate+ndeath);j++){
 3930: 	    k=k+1;
 3931: 	    mu[k][(int) age]=pmmij[i][j];
 3932: 	  }
 3933: 	}
 3934:      	for(i=1;i<=(nlstate)*(nlstate+ndeath);i++)
 3935: 	  for(j=1;j<=(nlstate)*(nlstate+ndeath);j++)
 3936: 	    varpij[i][j][(int)age] = doldm[i][j];
 3937: 
 3938: 	/*printf("\n%d ",(int)age);
 3939: 	  for (i=1; i<=(nlstate)*(nlstate+ndeath);i++){
 3940: 	  printf("%e [%e ;%e] ",gm[i],gm[i]-2*sqrt(doldm[i][i]),gm[i]+2*sqrt(doldm[i][i]));
 3941: 	  fprintf(ficlog,"%e [%e ;%e] ",gm[i],gm[i]-2*sqrt(doldm[i][i]),gm[i]+2*sqrt(doldm[i][i]));
 3942: 	  }*/
 3943: 
 3944: 	fprintf(ficresprob,"\n%d ",(int)age);
 3945: 	fprintf(ficresprobcov,"\n%d ",(int)age);
 3946: 	fprintf(ficresprobcor,"\n%d ",(int)age);
 3947: 
 3948: 	for (i=1; i<=(nlstate)*(nlstate+ndeath);i++)
 3949: 	  fprintf(ficresprob,"%11.3e (%11.3e) ",mu[i][(int) age],sqrt(varpij[i][i][(int)age]));
 3950: 	for (i=1; i<=(nlstate)*(nlstate+ndeath);i++){
 3951: 	  fprintf(ficresprobcov,"%11.3e ",mu[i][(int) age]);
 3952: 	  fprintf(ficresprobcor,"%11.3e ",mu[i][(int) age]);
 3953: 	}
 3954: 	i=0;
 3955: 	for (k=1; k<=(nlstate);k++){
 3956:  	  for (l=1; l<=(nlstate+ndeath);l++){ 
 3957:  	    i++;
 3958: 	    fprintf(ficresprobcov,"\n%d %d-%d",(int)age,k,l);
 3959: 	    fprintf(ficresprobcor,"\n%d %d-%d",(int)age,k,l);
 3960: 	    for (j=1; j<=i;j++){
 3961: 	      /* printf(" k=%d l=%d i=%d j=%d\n",k,l,i,j);fflush(stdout); */
 3962: 	      fprintf(ficresprobcov," %11.3e",varpij[i][j][(int)age]);
 3963: 	      fprintf(ficresprobcor," %11.3e",varpij[i][j][(int) age]/sqrt(varpij[i][i][(int) age])/sqrt(varpij[j][j][(int)age]));
 3964: 	    }
 3965: 	  }
 3966: 	}/* end of loop for state */
 3967:       } /* end of loop for age */
 3968:       free_vector(gp,1,(nlstate+ndeath)*(nlstate+ndeath));
 3969:       free_vector(gm,1,(nlstate+ndeath)*(nlstate+ndeath));
 3970:       free_matrix(trgradg,1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
 3971:       free_matrix(gradg,1,(nlstate+ndeath)*(nlstate+ndeath),1,npar);
 3972:       
 3973:       /* Confidence intervalle of pij  */
 3974:       /*
 3975: 	fprintf(ficgp,"\nunset parametric;unset label");
 3976: 	fprintf(ficgp,"\nset log y;unset log x; set xlabel \"Age\";set ylabel \"probability (year-1)\"");
 3977: 	fprintf(ficgp,"\nset ter png small\nset size 0.65,0.65");
 3978: 	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);
 3979: 	fprintf(fichtm,"\n<br><img src=\"pijgr%s.png\"> ",optionfilefiname);
 3980: 	fprintf(ficgp,"\nset out \"pijgr%s.png\"",optionfilefiname);
 3981: 	fprintf(ficgp,"\nplot \"%s\" every :::%d::%d u 1:2 \"\%%lf",k1,k2,xfilevarprob);
 3982:       */
 3983: 
 3984:       /* Drawing ellipsoids of confidence of two variables p(k1-l1,k2-l2)*/
 3985:       first1=1;first2=2;
 3986:       for (k2=1; k2<=(nlstate);k2++){
 3987: 	for (l2=1; l2<=(nlstate+ndeath);l2++){ 
 3988: 	  if(l2==k2) continue;
 3989: 	  j=(k2-1)*(nlstate+ndeath)+l2;
 3990: 	  for (k1=1; k1<=(nlstate);k1++){
 3991: 	    for (l1=1; l1<=(nlstate+ndeath);l1++){ 
 3992: 	      if(l1==k1) continue;
 3993: 	      i=(k1-1)*(nlstate+ndeath)+l1;
 3994: 	      if(i<=j) continue;
 3995: 	      for (age=bage; age<=fage; age ++){ 
 3996: 		if ((int)age %5==0){
 3997: 		  v1=varpij[i][i][(int)age]/stepm*YEARM/stepm*YEARM;
 3998: 		  v2=varpij[j][j][(int)age]/stepm*YEARM/stepm*YEARM;
 3999: 		  cv12=varpij[i][j][(int)age]/stepm*YEARM/stepm*YEARM;
 4000: 		  mu1=mu[i][(int) age]/stepm*YEARM ;
 4001: 		  mu2=mu[j][(int) age]/stepm*YEARM;
 4002: 		  c12=cv12/sqrt(v1*v2);
 4003: 		  /* Computing eigen value of matrix of covariance */
 4004: 		  lc1=((v1+v2)+sqrt((v1+v2)*(v1+v2) - 4*(v1*v2-cv12*cv12)))/2.;
 4005: 		  lc2=((v1+v2)-sqrt((v1+v2)*(v1+v2) - 4*(v1*v2-cv12*cv12)))/2.;
 4006: 		  if ((lc2 <0) || (lc1 <0) ){
 4007: 		    if(first2==1){
 4008: 		      first1=0;
 4009: 		    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);
 4010: 		    }
 4011: 		    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);
 4012: 		    /* lc1=fabs(lc1); */ /* If we want to have them positive */
 4013: 		    /* lc2=fabs(lc2); */
 4014: 		  }
 4015: 
 4016: 		  /* Eigen vectors */
 4017: 		  v11=(1./sqrt(1+(v1-lc1)*(v1-lc1)/cv12/cv12));
 4018: 		  /*v21=sqrt(1.-v11*v11); *//* error */
 4019: 		  v21=(lc1-v1)/cv12*v11;
 4020: 		  v12=-v21;
 4021: 		  v22=v11;
 4022: 		  tnalp=v21/v11;
 4023: 		  if(first1==1){
 4024: 		    first1=0;
 4025: 		    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);
 4026: 		  }
 4027: 		  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);
 4028: 		  /*printf(fignu*/
 4029: 		  /* mu1+ v11*lc1*cost + v12*lc2*sin(t) */
 4030: 		  /* mu2+ v21*lc1*cost + v22*lc2*sin(t) */
 4031: 		  if(first==1){
 4032: 		    first=0;
 4033:  		    fprintf(ficgp,"\nset parametric;unset label");
 4034: 		    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);
 4035: 		    fprintf(ficgp,"\nset ter png small size 320, 240");
 4036: 		    fprintf(fichtmcov,"\n<br>Ellipsoids of confidence cov(p%1d%1d,p%1d%1d) expressed in year<sup>-1</sup>\
 4037:  :<a href=\"%s%d%1d%1d-%1d%1d.png\">\
 4038: %s%d%1d%1d-%1d%1d.png</A>, ",k1,l1,k2,l2,\
 4039: 			    subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2,\
 4040: 			    subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4041: 		    fprintf(fichtmcov,"\n<br><img src=\"%s%d%1d%1d-%1d%1d.png\"> ",subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4042: 		    fprintf(fichtmcov,"\n<br> Correlation at age %d (%.3f),",(int) age, c12);
 4043: 		    fprintf(ficgp,"\nset out \"%s%d%1d%1d-%1d%1d.png\"",subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4044: 		    fprintf(ficgp,"\nset label \"%d\" at %11.3e,%11.3e center",(int) age, mu1,mu2);
 4045: 		    fprintf(ficgp,"\n# Age %d, p%1d%1d - p%1d%1d",(int) age, k1,l1,k2,l2);
 4046: 		    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",\
 4047: 			    mu1,std,v11,sqrt(lc1),v12,sqrt(lc2),\
 4048: 			    mu2,std,v21,sqrt(lc1),v22,sqrt(lc2));
 4049: 		  }else{
 4050: 		    first=0;
 4051: 		    fprintf(fichtmcov," %d (%.3f),",(int) age, c12);
 4052: 		    fprintf(ficgp,"\n# Age %d, p%1d%1d - p%1d%1d",(int) age, k1,l1,k2,l2);
 4053: 		    fprintf(ficgp,"\nset label \"%d\" at %11.3e,%11.3e center",(int) age, mu1,mu2);
 4054: 		    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",\
 4055: 			    mu1,std,v11,sqrt(lc1),v12,sqrt(lc2),\
 4056: 			    mu2,std,v21,sqrt(lc1),v22,sqrt(lc2));
 4057: 		  }/* if first */
 4058: 		} /* age mod 5 */
 4059: 	      } /* end loop age */
 4060: 	      fprintf(ficgp,"\nset out \"%s%d%1d%1d-%1d%1d.png\";replot;",subdirf2(optionfilefiname,"varpijgr"), j1,k1,l1,k2,l2);
 4061: 	      first=1;
 4062: 	    } /*l12 */
 4063: 	  } /* k12 */
 4064: 	} /*l1 */
 4065:       }/* k1 */
 4066:       /* } */ /* loop covariates */
 4067:   }
 4068:   free_ma3x(varpij,1,nlstate,1,nlstate+ndeath,(int) bage, (int)fage);
 4069:   free_matrix(mu,1,(nlstate+ndeath)*(nlstate+ndeath),(int) bage, (int)fage);
 4070:   free_matrix(doldm,1,(nlstate)*(nlstate+ndeath),1,(nlstate)*(nlstate+ndeath));
 4071:   free_matrix(dnewm,1,(nlstate)*(nlstate+ndeath),1,npar);
 4072:   free_vector(xp,1,npar);
 4073:   fclose(ficresprob);
 4074:   fclose(ficresprobcov);
 4075:   fclose(ficresprobcor);
 4076:   fflush(ficgp);
 4077:   fflush(fichtmcov);
 4078: }
 4079: 
 4080: 
 4081: /******************* Printing html file ***********/
 4082: void printinghtml(char fileres[], char title[], char datafile[], int firstpass, \
 4083: 		  int lastpass, int stepm, int weightopt, char model[],\
 4084: 		  int imx,int jmin, int jmax, double jmeanint,char rfileres[],\
 4085: 		  int popforecast, int estepm ,\
 4086: 		  double jprev1, double mprev1,double anprev1, \
 4087: 		  double jprev2, double mprev2,double anprev2){
 4088:   int jj1, k1, i1, cpt;
 4089: 
 4090:    fprintf(fichtm,"<ul><li><a href='#firstorder'>Result files (first order: no variance)</a>\n \
 4091:    <li><a href='#secondorder'>Result files (second order (variance)</a>\n \
 4092: </ul>");
 4093:    fprintf(fichtm,"<ul><li><h4><a name='firstorder'>Result files (first order: no variance)</a></h4>\n \
 4094:  - Observed prevalence in each state (during the period defined between %.lf/%.lf/%.lf and %.lf/%.lf/%.lf): <a href=\"%s\">%s</a> <br>\n ",
 4095: 	   jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,subdirf2(fileres,"p"),subdirf2(fileres,"p"));
 4096:    fprintf(fichtm,"\
 4097:  - Estimated transition probabilities over %d (stepm) months: <a href=\"%s\">%s</a><br>\n ",
 4098: 	   stepm,subdirf2(fileres,"pij"),subdirf2(fileres,"pij"));
 4099:    fprintf(fichtm,"\
 4100:  - Period (stable) prevalence in each health state: <a href=\"%s\">%s</a> <br>\n",
 4101: 	   subdirf2(fileres,"pl"),subdirf2(fileres,"pl"));
 4102:    fprintf(fichtm,"\
 4103:  - (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): \
 4104:    <a href=\"%s\">%s</a> <br>\n",
 4105: 	   estepm,subdirf2(fileres,"e"),subdirf2(fileres,"e"));
 4106:    fprintf(fichtm,"\
 4107:  - Population projections by age and states: \
 4108:    <a href=\"%s\">%s</a> <br>\n</li>", subdirf2(fileres,"f"),subdirf2(fileres,"f"));
 4109: 
 4110: fprintf(fichtm," \n<ul><li><b>Graphs</b></li><p>");
 4111: 
 4112:  m=pow(2,cptcoveff);
 4113:  if (cptcovn < 1) {m=1;ncodemax[1]=1;}
 4114: 
 4115:  jj1=0;
 4116:  for(k1=1; k1<=m;k1++){
 4117:    for(i1=1; i1<=ncodemax[k1];i1++){
 4118:      jj1++;
 4119:      if (cptcovn > 0) {
 4120:        fprintf(fichtm,"<hr  size=\"2\" color=\"#EC5E5E\">************ Results for covariates");
 4121:        for (cpt=1; cpt<=cptcoveff;cpt++) 
 4122: 	 fprintf(fichtm," V%d=%d ",Tvaraff[cpt],nbcode[Tvaraff[cpt]][codtab[jj1][cpt]]);
 4123:        fprintf(fichtm," ************\n<hr size=\"2\" color=\"#EC5E5E\">");
 4124:      }
 4125:      /* Pij */
 4126:      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> \
 4127: <img src=\"%s%d_1.png\">",stepm,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1);     
 4128:      /* Quasi-incidences */
 4129:      fprintf(fichtm,"<br>- Pij or Conditional probabilities to be observed in state j being in state i %d (stepm) months\
 4130:  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> \
 4131: <img src=\"%s%d_2.png\">",stepm,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1,subdirf2(optionfilefiname,"pe"),jj1); 
 4132:        /* Period (stable) prevalence in each health state */
 4133:        for(cpt=1; cpt<=nlstate;cpt++){
 4134: 	 fprintf(fichtm,"<br>- Convergence to period (stable) prevalence in state %d. Or probability to be in state %d being in state (1 to %d) at different ages. <a href=\"%s%d_%d.png\">%s%d_%d.png</a><br> \
 4135: <img src=\"%s%d_%d.png\">", cpt, cpt, nlstate, subdirf2(optionfilefiname,"p"),cpt,jj1,subdirf2(optionfilefiname,"p"),cpt,jj1,subdirf2(optionfilefiname,"p"),cpt,jj1);
 4136:        }
 4137:      for(cpt=1; cpt<=nlstate;cpt++) {
 4138:         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> \
 4139: <img src=\"%s%d%d.png\">",cpt,nlstate,subdirf2(optionfilefiname,"exp"),cpt,jj1,subdirf2(optionfilefiname,"exp"),cpt,jj1,subdirf2(optionfilefiname,"exp"),cpt,jj1);
 4140:      }
 4141:    } /* end i1 */
 4142:  }/* End k1 */
 4143:  fprintf(fichtm,"</ul>");
 4144: 
 4145: 
 4146:  fprintf(fichtm,"\
 4147: \n<br><li><h4> <a name='secondorder'>Result files (second order: variances)</a></h4>\n\
 4148:  - Parameter file with estimated parameters and covariance matrix: <a href=\"%s\">%s</a> <br>\n", rfileres,rfileres);
 4149: 
 4150:  fprintf(fichtm," - Variance of one-step probabilities: <a href=\"%s\">%s</a> <br>\n",
 4151: 	 subdirf2(fileres,"prob"),subdirf2(fileres,"prob"));
 4152:  fprintf(fichtm,"\
 4153:  - Variance-covariance of one-step probabilities: <a href=\"%s\">%s</a> <br>\n",
 4154: 	 subdirf2(fileres,"probcov"),subdirf2(fileres,"probcov"));
 4155: 
 4156:  fprintf(fichtm,"\
 4157:  - Correlation matrix of one-step probabilities: <a href=\"%s\">%s</a> <br>\n",
 4158: 	 subdirf2(fileres,"probcor"),subdirf2(fileres,"probcor"));
 4159:  fprintf(fichtm,"\
 4160:  - 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): \
 4161:    <a href=\"%s\">%s</a> <br>\n</li>",
 4162: 	   estepm,subdirf2(fileres,"cve"),subdirf2(fileres,"cve"));
 4163:  fprintf(fichtm,"\
 4164:  - (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): \
 4165:    <a href=\"%s\">%s</a> <br>\n</li>",
 4166: 	   estepm,subdirf2(fileres,"stde"),subdirf2(fileres,"stde"));
 4167:  fprintf(fichtm,"\
 4168:  - 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",
 4169: 	 estepm, subdirf2(fileres,"v"),subdirf2(fileres,"v"));
 4170:  fprintf(fichtm,"\
 4171:  - 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",
 4172: 	 estepm, subdirf2(fileres,"t"),subdirf2(fileres,"t"));
 4173:  fprintf(fichtm,"\
 4174:  - Standard deviation of period (stable) prevalences: <a href=\"%s\">%s</a> <br>\n",\
 4175: 	 subdirf2(fileres,"vpl"),subdirf2(fileres,"vpl"));
 4176: 
 4177: /*  if(popforecast==1) fprintf(fichtm,"\n */
 4178: /*  - Prevalences forecasting: <a href=\"f%s\">f%s</a> <br>\n */
 4179: /*  - Population forecasting (if popforecast=1): <a href=\"pop%s\">pop%s</a> <br>\n */
 4180: /* 	<br>",fileres,fileres,fileres,fileres); */
 4181: /*  else  */
 4182: /*    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); */
 4183:  fflush(fichtm);
 4184:  fprintf(fichtm," <ul><li><b>Graphs</b></li><p>");
 4185: 
 4186:  m=pow(2,cptcoveff);
 4187:  if (cptcovn < 1) {m=1;ncodemax[1]=1;}
 4188: 
 4189:  jj1=0;
 4190:  for(k1=1; k1<=m;k1++){
 4191:    for(i1=1; i1<=ncodemax[k1];i1++){
 4192:      jj1++;
 4193:      if (cptcovn > 0) {
 4194:        fprintf(fichtm,"<hr  size=\"2\" color=\"#EC5E5E\">************ Results for covariates");
 4195:        for (cpt=1; cpt<=cptcoveff;cpt++) 
 4196: 	 fprintf(fichtm," V%d=%d ",Tvaraff[cpt],nbcode[Tvaraff[cpt]][codtab[jj1][cpt]]);
 4197:        fprintf(fichtm," ************\n<hr size=\"2\" color=\"#EC5E5E\">");
 4198:      }
 4199:      for(cpt=1; cpt<=nlstate;cpt++) {
 4200:        fprintf(fichtm,"<br>- Observed (cross-sectional) and period (incidence based) \
 4201: prevalence (with 95%% confidence interval) in state (%d): %s%d_%d.png <br>\
 4202: <img src=\"%s%d_%d.png\">",cpt,subdirf2(optionfilefiname,"v"),cpt,jj1,subdirf2(optionfilefiname,"v"),cpt,jj1);  
 4203:      }
 4204:      fprintf(fichtm,"\n<br>- Total life expectancy by age and \
 4205: health expectancies in states (1) and (2). If popbased=1 the smooth (due to the model) \
 4206: true period expectancies (those weighted with period prevalences are also\
 4207:  drawn in addition to the population based expectancies computed using\
 4208:  observed and cahotic prevalences: %s%d.png<br>\
 4209: <img src=\"%s%d.png\">",subdirf2(optionfilefiname,"e"),jj1,subdirf2(optionfilefiname,"e"),jj1);
 4210:    } /* end i1 */
 4211:  }/* End k1 */
 4212:  fprintf(fichtm,"</ul>");
 4213:  fflush(fichtm);
 4214: }
 4215: 
 4216: /******************* Gnuplot file **************/
 4217: void printinggnuplot(char fileres[], char optionfilefiname[], double ageminpar, double agemaxpar, double fage , char pathc[], double p[]){
 4218: 
 4219:   char dirfileres[132],optfileres[132];
 4220:   int cpt=0,k1=0,i=0,k=0,j=0,jk=0,k2=0,k3=0,ij=0,l=0;
 4221:   int ng=0;
 4222: /*   if((ficgp=fopen(optionfilegnuplot,"a"))==NULL) { */
 4223: /*     printf("Problem with file %s",optionfilegnuplot); */
 4224: /*     fprintf(ficlog,"Problem with file %s",optionfilegnuplot); */
 4225: /*   } */
 4226: 
 4227:   /*#ifdef windows */
 4228:   fprintf(ficgp,"cd \"%s\" \n",pathc);
 4229:     /*#endif */
 4230:   m=pow(2,cptcoveff);
 4231: 
 4232:   strcpy(dirfileres,optionfilefiname);
 4233:   strcpy(optfileres,"vpl");
 4234:  /* 1eme*/
 4235:   fprintf(ficgp,"\n# 1st: Period (stable) prevalence with CI: 'vpl' files\n");
 4236:   for (cpt=1; cpt<= nlstate ; cpt ++) {
 4237:     for (k1=1; k1<= m ; k1 ++) { /* plot [100000000000000000000:-100000000000000000000] "mysbiaspar/vplrmysbiaspar.txt to check */
 4238:      fprintf(ficgp,"\nset out \"%s%d_%d.png\" \n",subdirf2(optionfilefiname,"v"),cpt,k1);
 4239:      fprintf(ficgp,"\n#set out \"v%s%d_%d.png\" \n",optionfilefiname,cpt,k1);
 4240:      fprintf(ficgp,"set xlabel \"Age\" \n\
 4241: set ylabel \"Probability\" \n\
 4242: set ter png small size 320, 240\n\
 4243: plot [%.f:%.f] \"%s\" every :::%d::%d u 1:2 \"%%lf",ageminpar,fage,subdirf2(fileres,"vpl"),k1-1,k1-1);
 4244: 
 4245:      for (i=1; i<= nlstate ; i ++) {
 4246:        if (i==cpt) fprintf(ficgp," %%lf (%%lf)");
 4247:        else        fprintf(ficgp," %%*lf (%%*lf)");
 4248:      }
 4249:      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);
 4250:      for (i=1; i<= nlstate ; i ++) {
 4251:        if (i==cpt) fprintf(ficgp," %%lf (%%lf)");
 4252:        else fprintf(ficgp," %%*lf (%%*lf)");
 4253:      } 
 4254:      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); 
 4255:      for (i=1; i<= nlstate ; i ++) {
 4256:        if (i==cpt) fprintf(ficgp," %%lf (%%lf)");
 4257:        else fprintf(ficgp," %%*lf (%%*lf)");
 4258:      }  
 4259:      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));
 4260:    }
 4261:   }
 4262:   /*2 eme*/
 4263:   fprintf(ficgp,"\n# 2nd: Total life expectancy with CI: 't' files\n");
 4264:   for (k1=1; k1<= m ; k1 ++) { 
 4265:     fprintf(ficgp,"\nset out \"%s%d.png\" \n",subdirf2(optionfilefiname,"e"),k1);
 4266:     fprintf(ficgp,"set ylabel \"Years\" \nset ter png small size 320, 240\nplot [%.f:%.f] ",ageminpar,fage);
 4267:     
 4268:     for (i=1; i<= nlstate+1 ; i ++) {
 4269:       k=2*i;
 4270:       fprintf(ficgp,"\"%s\" every :::%d::%d u 1:2 \"%%lf",subdirf2(fileres,"t"),k1-1,k1-1);
 4271:       for (j=1; j<= nlstate+1 ; j ++) {
 4272: 	if (j==i) fprintf(ficgp," %%lf (%%lf)");
 4273: 	else fprintf(ficgp," %%*lf (%%*lf)");
 4274:       }   
 4275:       if (i== 1) fprintf(ficgp,"\" t\"TLE\" w l ,");
 4276:       else fprintf(ficgp,"\" t\"LE in state (%d)\" w l ,",i-1);
 4277:       fprintf(ficgp,"\"%s\" every :::%d::%d u 1:($2-$3*2) \"%%lf",subdirf2(fileres,"t"),k1-1,k1-1);
 4278:       for (j=1; j<= nlstate+1 ; j ++) {
 4279: 	if (j==i) fprintf(ficgp," %%lf (%%lf)");
 4280: 	else fprintf(ficgp," %%*lf (%%*lf)");
 4281:       }   
 4282:       fprintf(ficgp,"\" t\"\" w l lt 0,");
 4283:       fprintf(ficgp,"\"%s\" every :::%d::%d u 1:($2+$3*2) \"%%lf",subdirf2(fileres,"t"),k1-1,k1-1);
 4284:       for (j=1; j<= nlstate+1 ; j ++) {
 4285: 	if (j==i) fprintf(ficgp," %%lf (%%lf)");
 4286: 	else fprintf(ficgp," %%*lf (%%*lf)");
 4287:       }   
 4288:       if (i== (nlstate+1)) fprintf(ficgp,"\" t\"\" w l lt 0");
 4289:       else fprintf(ficgp,"\" t\"\" w l lt 0,");
 4290:     }
 4291:   }
 4292:   
 4293:   /*3eme*/
 4294:   
 4295:   for (k1=1; k1<= m ; k1 ++) { 
 4296:     for (cpt=1; cpt<= nlstate ; cpt ++) {
 4297:       /*       k=2+nlstate*(2*cpt-2); */
 4298:       k=2+(nlstate+1)*(cpt-1);
 4299:       fprintf(ficgp,"\nset out \"%s%d%d.png\" \n",subdirf2(optionfilefiname,"exp"),cpt,k1);
 4300:       fprintf(ficgp,"set ter png small size 320, 240\n\
 4301: 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);
 4302:       /*fprintf(ficgp,",\"e%s\" every :::%d::%d u 1:($%d-2*$%d) \"\%%lf ",fileres,k1-1,k1-1,k,k+1);
 4303: 	for (i=1; i<= nlstate*2 ; i ++) fprintf(ficgp,"\%%lf (\%%lf) ");
 4304: 	fprintf(ficgp,"\" t \"e%d1\" w l",cpt);
 4305: 	fprintf(ficgp,",\"e%s\" every :::%d::%d u 1:($%d+2*$%d) \"\%%lf ",fileres,k1-1,k1-1,k,k+1);
 4306: 	for (i=1; i<= nlstate*2 ; i ++) fprintf(ficgp,"\%%lf (\%%lf) ");
 4307: 	fprintf(ficgp,"\" t \"e%d1\" w l",cpt);
 4308: 	
 4309:       */
 4310:       for (i=1; i< nlstate ; i ++) {
 4311: 	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);
 4312: 	/*	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);*/
 4313: 	
 4314:       } 
 4315:       fprintf(ficgp," ,\"%s\" every :::%d::%d u 1:%d t \"e%d.\" w l",subdirf2(fileres,"e"),k1-1,k1-1,k+nlstate,cpt);
 4316:     }
 4317:   }
 4318:   
 4319:   /* CV preval stable (period) */
 4320:   for (k1=1; k1<= m ; k1 ++) { /* For each multivariate if any */
 4321:     for (cpt=1; cpt<=nlstate ; cpt ++) { /* For each life state */
 4322:       k=3;
 4323:       fprintf(ficgp,"\n#\n#\n#CV preval stable (period): 'pij' files, cov=%d state=%d",k1, cpt);
 4324:       fprintf(ficgp,"\nset out \"%s%d_%d.png\" \n",subdirf2(optionfilefiname,"p"),cpt,k1);
 4325:       fprintf(ficgp,"set xlabel \"Age\" \nset ylabel \"Probability\" \n\
 4326: set ter png small size 320, 240\n\
 4327: unset log y\n\
 4328: plot [%.f:%.f]  ", ageminpar, agemaxpar);
 4329:       for (i=1; i<= nlstate ; i ++){
 4330: 	if(i==1)
 4331: 	  fprintf(ficgp,"\"%s\"",subdirf2(fileres,"pij"));
 4332: 	else
 4333: 	  fprintf(ficgp,", '' ");
 4334: 	l=(nlstate+ndeath)*(i-1)+1;
 4335: 	fprintf(ficgp," u ($1==%d ? ($3):1/0):($%d/($%d",k1,k+l+(cpt-1),k+l);
 4336: 	for (j=1; j<= (nlstate-1) ; j ++)
 4337: 	  fprintf(ficgp,"+$%d",k+l+j);
 4338: 	fprintf(ficgp,")) t \"prev(%d,%d)\" w l",i,cpt);
 4339:       } /* nlstate */
 4340:       fprintf(ficgp,"\n");
 4341:     } /* end cpt state*/ 
 4342:   } /* end covariate */  
 4343:   
 4344:   /* proba elementaires */
 4345:   for(i=1,jk=1; i <=nlstate; i++){
 4346:     for(k=1; k <=(nlstate+ndeath); k++){
 4347:       if (k != i) {
 4348: 	for(j=1; j <=ncovmodel; j++){
 4349: 	  fprintf(ficgp,"p%d=%f ",jk,p[jk]);
 4350: 	  jk++; 
 4351: 	  fprintf(ficgp,"\n");
 4352: 	}
 4353:       }
 4354:     }
 4355:    }
 4356:   /*goto avoid;*/
 4357:    for(ng=1; ng<=2;ng++){ /* Number of graphics: first is probabilities second is incidence per year*/
 4358:      for(jk=1; jk <=m; jk++) {
 4359:        fprintf(ficgp,"\nset out \"%s%d_%d.png\" \n",subdirf2(optionfilefiname,"pe"),jk,ng); 
 4360:        if (ng==2)
 4361: 	 fprintf(ficgp,"\nset ylabel \"Quasi-incidence per year\"\n");
 4362:        else
 4363: 	 fprintf(ficgp,"\nset title \"Probability\"\n");
 4364:        fprintf(ficgp,"\nset ter png small size 320, 240\nset log y\nplot  [%.f:%.f] ",ageminpar,agemaxpar);
 4365:        i=1;
 4366:        for(k2=1; k2<=nlstate; k2++) {
 4367: 	 k3=i;
 4368: 	 for(k=1; k<=(nlstate+ndeath); k++) {
 4369: 	   if (k != k2){
 4370: 	     if(ng==2)
 4371: 	       fprintf(ficgp," %f*exp(p%d+p%d*x",YEARM/stepm,i,i+1);
 4372: 	     else
 4373: 	       fprintf(ficgp," exp(p%d+p%d*x",i,i+1);
 4374: 	     ij=1;/* To be checked else nbcode[0][0] wrong */
 4375: 	     for(j=3; j <=ncovmodel; j++) {
 4376: 	       /* if(((j-2)==Tage[ij]) &&(ij <=cptcovage)) { /\* Bug valgrind *\/ */
 4377: 	       /* 	 /\*fprintf(ficgp,"+p%d*%d*x",i+j-1,nbcode[Tvar[j-2]][codtab[jk][Tvar[j-2]]]);*\/ */
 4378: 	       /* 	 ij++; */
 4379: 	       /* } */
 4380: 	       /* else */
 4381: 		 fprintf(ficgp,"+p%d*%d",i+j-1,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
 4382: 	     }
 4383: 	     fprintf(ficgp,")/(1");
 4384: 	     
 4385: 	     for(k1=1; k1 <=nlstate; k1++){   
 4386: 	       fprintf(ficgp,"+exp(p%d+p%d*x",k3+(k1-1)*ncovmodel,k3+(k1-1)*ncovmodel+1);
 4387: 	       ij=1;
 4388: 	       for(j=3; j <=ncovmodel; j++){
 4389: 		 /* if(((j-2)==Tage[ij]) &&(ij <=cptcovage)) { */
 4390: 		 /*   fprintf(ficgp,"+p%d*%d*x",k3+(k1-1)*ncovmodel+1+j-2,nbcode[Tvar[j-2]][codtab[jk][Tvar[j-2]]]); */
 4391: 		 /*   ij++; */
 4392: 		 /* } */
 4393: 		 /* else */
 4394: 		   fprintf(ficgp,"+p%d*%d",k3+(k1-1)*ncovmodel+1+j-2,nbcode[Tvar[j-2]][codtab[jk][j-2]]);
 4395: 	       }
 4396: 	       fprintf(ficgp,")");
 4397: 	     }
 4398: 	     fprintf(ficgp,") t \"p%d%d\" ", k2,k);
 4399: 	     if ((k+k2)!= (nlstate*2+ndeath)) fprintf(ficgp,",");
 4400: 	     i=i+ncovmodel;
 4401: 	   }
 4402: 	 } /* end k */
 4403:        } /* end k2 */
 4404:      } /* end jk */
 4405:    } /* end ng */
 4406:  /* avoid: */
 4407:    fflush(ficgp); 
 4408: }  /* end gnuplot */
 4409: 
 4410: 
 4411: /*************** Moving average **************/
 4412: int movingaverage(double ***probs, double bage,double fage, double ***mobaverage, int mobilav){
 4413: 
 4414:   int i, cpt, cptcod;
 4415:   int modcovmax =1;
 4416:   int mobilavrange, mob;
 4417:   double age;
 4418: 
 4419:   modcovmax=2*cptcoveff;/* Max number of modalities. We suppose 
 4420: 			   a covariate has 2 modalities */
 4421:   if (cptcovn<1) modcovmax=1; /* At least 1 pass */
 4422: 
 4423:   if(mobilav==1||mobilav ==3 ||mobilav==5 ||mobilav== 7){
 4424:     if(mobilav==1) mobilavrange=5; /* default */
 4425:     else mobilavrange=mobilav;
 4426:     for (age=bage; age<=fage; age++)
 4427:       for (i=1; i<=nlstate;i++)
 4428: 	for (cptcod=1;cptcod<=modcovmax;cptcod++)
 4429: 	  mobaverage[(int)age][i][cptcod]=probs[(int)age][i][cptcod];
 4430:     /* We keep the original values on the extreme ages bage, fage and for 
 4431:        fage+1 and bage-1 we use a 3 terms moving average; for fage+2 bage+2
 4432:        we use a 5 terms etc. until the borders are no more concerned. 
 4433:     */ 
 4434:     for (mob=3;mob <=mobilavrange;mob=mob+2){
 4435:       for (age=bage+(mob-1)/2; age<=fage-(mob-1)/2; age++){
 4436: 	for (i=1; i<=nlstate;i++){
 4437: 	  for (cptcod=1;cptcod<=modcovmax;cptcod++){
 4438: 	    mobaverage[(int)age][i][cptcod] =probs[(int)age][i][cptcod];
 4439: 	      for (cpt=1;cpt<=(mob-1)/2;cpt++){
 4440: 		mobaverage[(int)age][i][cptcod] +=probs[(int)age-cpt][i][cptcod];
 4441: 		mobaverage[(int)age][i][cptcod] +=probs[(int)age+cpt][i][cptcod];
 4442: 	      }
 4443: 	    mobaverage[(int)age][i][cptcod]=mobaverage[(int)age][i][cptcod]/mob;
 4444: 	  }
 4445: 	}
 4446:       }/* end age */
 4447:     }/* end mob */
 4448:   }else return -1;
 4449:   return 0;
 4450: }/* End movingaverage */
 4451: 
 4452: 
 4453: /************** Forecasting ******************/
 4454: void 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){
 4455:   /* proj1, year, month, day of starting projection 
 4456:      agemin, agemax range of age
 4457:      dateprev1 dateprev2 range of dates during which prevalence is computed
 4458:      anproj2 year of en of projection (same day and month as proj1).
 4459:   */
 4460:   int yearp, stepsize, hstepm, nhstepm, j, k, cptcod, i, h, i1;
 4461:   double agec; /* generic age */
 4462:   double agelim, ppij, yp,yp1,yp2,jprojmean,mprojmean,anprojmean;
 4463:   double *popeffectif,*popcount;
 4464:   double ***p3mat;
 4465:   double ***mobaverage;
 4466:   char fileresf[FILENAMELENGTH];
 4467: 
 4468:   agelim=AGESUP;
 4469:   prevalence(probs, ageminpar, agemax, s, agev, nlstate, imx, Tvar, nbcode, ncodemax, mint, anint, dateprev1, dateprev2, firstpass, lastpass);
 4470:  
 4471:   strcpy(fileresf,"f"); 
 4472:   strcat(fileresf,fileres);
 4473:   if((ficresf=fopen(fileresf,"w"))==NULL) {
 4474:     printf("Problem with forecast resultfile: %s\n", fileresf);
 4475:     fprintf(ficlog,"Problem with forecast resultfile: %s\n", fileresf);
 4476:   }
 4477:   printf("Computing forecasting: result on file '%s' \n", fileresf);
 4478:   fprintf(ficlog,"Computing forecasting: result on file '%s' \n", fileresf);
 4479: 
 4480:   if (cptcoveff==0) ncodemax[cptcoveff]=1;
 4481: 
 4482:   if (mobilav!=0) {
 4483:     mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4484:     if (movingaverage(probs, ageminpar, fage, mobaverage,mobilav)!=0){
 4485:       fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 4486:       printf(" Error in movingaverage mobilav=%d\n",mobilav);
 4487:     }
 4488:   }
 4489: 
 4490:   stepsize=(int) (stepm+YEARM-1)/YEARM;
 4491:   if (stepm<=12) stepsize=1;
 4492:   if(estepm < stepm){
 4493:     printf ("Problem %d lower than %d\n",estepm, stepm);
 4494:   }
 4495:   else  hstepm=estepm;   
 4496: 
 4497:   hstepm=hstepm/stepm; 
 4498:   yp1=modf(dateintmean,&yp);/* extracts integral of datemean in yp  and
 4499:                                fractional in yp1 */
 4500:   anprojmean=yp;
 4501:   yp2=modf((yp1*12),&yp);
 4502:   mprojmean=yp;
 4503:   yp1=modf((yp2*30.5),&yp);
 4504:   jprojmean=yp;
 4505:   if(jprojmean==0) jprojmean=1;
 4506:   if(mprojmean==0) jprojmean=1;
 4507: 
 4508:   i1=cptcoveff;
 4509:   if (cptcovn < 1){i1=1;}
 4510:   
 4511:   fprintf(ficresf,"# Mean day of interviews %.lf/%.lf/%.lf (%.2f) between %.2f and %.2f \n",jprojmean,mprojmean,anprojmean,dateintmean,dateprev1,dateprev2); 
 4512:   
 4513:   fprintf(ficresf,"#****** Routine prevforecast **\n");
 4514: 
 4515: /* 	      if (h==(int)(YEARM*yearp)){ */
 4516:   for(cptcov=1, k=0;cptcov<=i1;cptcov++){
 4517:     for(cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
 4518:       k=k+1;
 4519:       fprintf(ficresf,"\n#******");
 4520:       for(j=1;j<=cptcoveff;j++) {
 4521: 	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]]);
 4522:       }
 4523:       fprintf(ficresf,"******\n");
 4524:       fprintf(ficresf,"# Covariate valuofcovar yearproj age");
 4525:       for(j=1; j<=nlstate+ndeath;j++){ 
 4526: 	for(i=1; i<=nlstate;i++) 	      
 4527:           fprintf(ficresf," p%d%d",i,j);
 4528: 	fprintf(ficresf," p.%d",j);
 4529:       }
 4530:       for (yearp=0; yearp<=(anproj2-anproj1);yearp +=stepsize) { 
 4531: 	fprintf(ficresf,"\n");
 4532: 	fprintf(ficresf,"\n# Forecasting at date %.lf/%.lf/%.lf ",jproj1,mproj1,anproj1+yearp);   
 4533: 
 4534:      	for (agec=fage; agec>=(ageminpar-1); agec--){ 
 4535: 	  nhstepm=(int) rint((agelim-agec)*YEARM/stepm); 
 4536: 	  nhstepm = nhstepm/hstepm; 
 4537: 	  p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4538: 	  oldm=oldms;savm=savms;
 4539: 	  hpxij(p3mat,nhstepm,agec,hstepm,p,nlstate,stepm,oldm,savm, k);  
 4540: 	
 4541: 	  for (h=0; h<=nhstepm; h++){
 4542: 	    if (h*hstepm/YEARM*stepm ==yearp) {
 4543:               fprintf(ficresf,"\n");
 4544:               for(j=1;j<=cptcoveff;j++) 
 4545:                 fprintf(ficresf,"%d %d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 4546: 	      fprintf(ficresf,"%.f %.f ",anproj1+yearp,agec+h*hstepm/YEARM*stepm);
 4547: 	    } 
 4548: 	    for(j=1; j<=nlstate+ndeath;j++) {
 4549: 	      ppij=0.;
 4550: 	      for(i=1; i<=nlstate;i++) {
 4551: 		if (mobilav==1) 
 4552: 		  ppij=ppij+p3mat[i][j][h]*mobaverage[(int)agec][i][cptcod];
 4553: 		else {
 4554: 		  ppij=ppij+p3mat[i][j][h]*probs[(int)(agec)][i][cptcod];
 4555: 		}
 4556: 		if (h*hstepm/YEARM*stepm== yearp) {
 4557: 		  fprintf(ficresf," %.3f", p3mat[i][j][h]);
 4558: 		}
 4559: 	      } /* end i */
 4560: 	      if (h*hstepm/YEARM*stepm==yearp) {
 4561: 		fprintf(ficresf," %.3f", ppij);
 4562: 	      }
 4563: 	    }/* end j */
 4564: 	  } /* end h */
 4565: 	  free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4566: 	} /* end agec */
 4567:       } /* end yearp */
 4568:     } /* end cptcod */
 4569:   } /* end  cptcov */
 4570:        
 4571:   if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4572: 
 4573:   fclose(ficresf);
 4574: }
 4575: 
 4576: /************** Forecasting *****not tested NB*************/
 4577: void 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){
 4578:   
 4579:   int cpt, stepsize, hstepm, nhstepm, j,k,c, cptcod, i,h;
 4580:   int *popage;
 4581:   double calagedatem, agelim, kk1, kk2;
 4582:   double *popeffectif,*popcount;
 4583:   double ***p3mat,***tabpop,***tabpopprev;
 4584:   double ***mobaverage;
 4585:   char filerespop[FILENAMELENGTH];
 4586: 
 4587:   tabpop= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4588:   tabpopprev= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4589:   agelim=AGESUP;
 4590:   calagedatem=(anpyram+mpyram/12.+jpyram/365.-dateintmean)*YEARM;
 4591:   
 4592:   prevalence(probs, ageminpar, agemax, s, agev, nlstate, imx, Tvar, nbcode, ncodemax, mint, anint, dateprev1, dateprev2, firstpass, lastpass);
 4593:   
 4594:   
 4595:   strcpy(filerespop,"pop"); 
 4596:   strcat(filerespop,fileres);
 4597:   if((ficrespop=fopen(filerespop,"w"))==NULL) {
 4598:     printf("Problem with forecast resultfile: %s\n", filerespop);
 4599:     fprintf(ficlog,"Problem with forecast resultfile: %s\n", filerespop);
 4600:   }
 4601:   printf("Computing forecasting: result on file '%s' \n", filerespop);
 4602:   fprintf(ficlog,"Computing forecasting: result on file '%s' \n", filerespop);
 4603: 
 4604:   if (cptcoveff==0) ncodemax[cptcoveff]=1;
 4605: 
 4606:   if (mobilav!=0) {
 4607:     mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4608:     if (movingaverage(probs, ageminpar, fage, mobaverage,mobilav)!=0){
 4609:       fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 4610:       printf(" Error in movingaverage mobilav=%d\n",mobilav);
 4611:     }
 4612:   }
 4613: 
 4614:   stepsize=(int) (stepm+YEARM-1)/YEARM;
 4615:   if (stepm<=12) stepsize=1;
 4616:   
 4617:   agelim=AGESUP;
 4618:   
 4619:   hstepm=1;
 4620:   hstepm=hstepm/stepm; 
 4621:   
 4622:   if (popforecast==1) {
 4623:     if((ficpop=fopen(popfile,"r"))==NULL) {
 4624:       printf("Problem with population file : %s\n",popfile);exit(0);
 4625:       fprintf(ficlog,"Problem with population file : %s\n",popfile);exit(0);
 4626:     } 
 4627:     popage=ivector(0,AGESUP);
 4628:     popeffectif=vector(0,AGESUP);
 4629:     popcount=vector(0,AGESUP);
 4630:     
 4631:     i=1;   
 4632:     while ((c=fscanf(ficpop,"%d %lf\n",&popage[i],&popcount[i])) != EOF) i=i+1;
 4633:    
 4634:     imx=i;
 4635:     for (i=1; i<imx;i++) popeffectif[popage[i]]=popcount[i];
 4636:   }
 4637: 
 4638:   for(cptcov=1,k=0;cptcov<=i2;cptcov++){
 4639:    for(cptcod=1;cptcod<=ncodemax[cptcoveff];cptcod++){
 4640:       k=k+1;
 4641:       fprintf(ficrespop,"\n#******");
 4642:       for(j=1;j<=cptcoveff;j++) {
 4643: 	fprintf(ficrespop," V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 4644:       }
 4645:       fprintf(ficrespop,"******\n");
 4646:       fprintf(ficrespop,"# Age");
 4647:       for(j=1; j<=nlstate+ndeath;j++) fprintf(ficrespop," P.%d",j);
 4648:       if (popforecast==1)  fprintf(ficrespop," [Population]");
 4649:       
 4650:       for (cpt=0; cpt<=0;cpt++) { 
 4651: 	fprintf(ficrespop,"\n\n# Forecasting at date %.lf/%.lf/%.lf ",jpyram,mpyram,anpyram+cpt);   
 4652: 	
 4653:      	for (agedeb=(fage-((int)calagedatem %12/12.)); agedeb>=(ageminpar-((int)calagedatem %12)/12.); agedeb--){ 
 4654: 	  nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm); 
 4655: 	  nhstepm = nhstepm/hstepm; 
 4656: 	  
 4657: 	  p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4658: 	  oldm=oldms;savm=savms;
 4659: 	  hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);  
 4660: 	
 4661: 	  for (h=0; h<=nhstepm; h++){
 4662: 	    if (h==(int) (calagedatem+YEARM*cpt)) {
 4663: 	      fprintf(ficrespop,"\n %3.f ",agedeb+h*hstepm/YEARM*stepm);
 4664: 	    } 
 4665: 	    for(j=1; j<=nlstate+ndeath;j++) {
 4666: 	      kk1=0.;kk2=0;
 4667: 	      for(i=1; i<=nlstate;i++) {	      
 4668: 		if (mobilav==1) 
 4669: 		  kk1=kk1+p3mat[i][j][h]*mobaverage[(int)agedeb+1][i][cptcod];
 4670: 		else {
 4671: 		  kk1=kk1+p3mat[i][j][h]*probs[(int)(agedeb+1)][i][cptcod];
 4672: 		}
 4673: 	      }
 4674: 	      if (h==(int)(calagedatem+12*cpt)){
 4675: 		tabpop[(int)(agedeb)][j][cptcod]=kk1;
 4676: 		  /*fprintf(ficrespop," %.3f", kk1);
 4677: 		    if (popforecast==1) fprintf(ficrespop," [%.f]", kk1*popeffectif[(int)agedeb+1]);*/
 4678: 	      }
 4679: 	    }
 4680: 	    for(i=1; i<=nlstate;i++){
 4681: 	      kk1=0.;
 4682: 		for(j=1; j<=nlstate;j++){
 4683: 		  kk1= kk1+tabpop[(int)(agedeb)][j][cptcod]; 
 4684: 		}
 4685: 		  tabpopprev[(int)(agedeb)][i][cptcod]=tabpop[(int)(agedeb)][i][cptcod]/kk1*popeffectif[(int)(agedeb+(calagedatem+12*cpt)*hstepm/YEARM*stepm-1)];
 4686: 	    }
 4687: 
 4688: 	    if (h==(int)(calagedatem+12*cpt)) for(j=1; j<=nlstate;j++) 
 4689: 	      fprintf(ficrespop," %15.2f",tabpopprev[(int)(agedeb+1)][j][cptcod]);
 4690: 	  }
 4691: 	  free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4692: 	}
 4693:       }
 4694:  
 4695:   /******/
 4696: 
 4697:       for (cpt=1; cpt<=(anpyram1-anpyram);cpt++) { 
 4698: 	fprintf(ficrespop,"\n\n# Forecasting at date %.lf/%.lf/%.lf ",jpyram,mpyram,anpyram+cpt);   
 4699: 	for (agedeb=(fage-((int)calagedatem %12/12.)); agedeb>=(ageminpar-((int)calagedatem %12)/12.); agedeb--){ 
 4700: 	  nhstepm=(int) rint((agelim-agedeb)*YEARM/stepm); 
 4701: 	  nhstepm = nhstepm/hstepm; 
 4702: 	  
 4703: 	  p3mat=ma3x(1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4704: 	  oldm=oldms;savm=savms;
 4705: 	  hpxij(p3mat,nhstepm,agedeb,hstepm,p,nlstate,stepm,oldm,savm, k);  
 4706: 	  for (h=0; h<=nhstepm; h++){
 4707: 	    if (h==(int) (calagedatem+YEARM*cpt)) {
 4708: 	      fprintf(ficresf,"\n %3.f ",agedeb+h*hstepm/YEARM*stepm);
 4709: 	    } 
 4710: 	    for(j=1; j<=nlstate+ndeath;j++) {
 4711: 	      kk1=0.;kk2=0;
 4712: 	      for(i=1; i<=nlstate;i++) {	      
 4713: 		kk1=kk1+p3mat[i][j][h]*tabpopprev[(int)agedeb+1][i][cptcod];	
 4714: 	      }
 4715: 	      if (h==(int)(calagedatem+12*cpt)) fprintf(ficresf," %15.2f", kk1);	
 4716: 	    }
 4717: 	  }
 4718: 	  free_ma3x(p3mat,1,nlstate+ndeath,1, nlstate+ndeath, 0,nhstepm);
 4719: 	}
 4720:       }
 4721:    } 
 4722:   }
 4723:  
 4724:   if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4725: 
 4726:   if (popforecast==1) {
 4727:     free_ivector(popage,0,AGESUP);
 4728:     free_vector(popeffectif,0,AGESUP);
 4729:     free_vector(popcount,0,AGESUP);
 4730:   }
 4731:   free_ma3x(tabpop,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4732:   free_ma3x(tabpopprev,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 4733:   fclose(ficrespop);
 4734: } /* End of popforecast */
 4735: 
 4736: int fileappend(FILE *fichier, char *optionfich)
 4737: {
 4738:   if((fichier=fopen(optionfich,"a"))==NULL) {
 4739:     printf("Problem with file: %s\n", optionfich);
 4740:     fprintf(ficlog,"Problem with file: %s\n", optionfich);
 4741:     return (0);
 4742:   }
 4743:   fflush(fichier);
 4744:   return (1);
 4745: }
 4746: 
 4747: 
 4748: /**************** function prwizard **********************/
 4749: void prwizard(int ncovmodel, int nlstate, int ndeath,  char model[], FILE *ficparo)
 4750: {
 4751: 
 4752:   /* Wizard to print covariance matrix template */
 4753: 
 4754:   char ca[32], cb[32];
 4755:   int i,j, k, li, lj, lk, ll, jj, npar, itimes;
 4756:   int numlinepar;
 4757: 
 4758:   printf("# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 4759:   fprintf(ficparo,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 4760:   for(i=1; i <=nlstate; i++){
 4761:     jj=0;
 4762:     for(j=1; j <=nlstate+ndeath; j++){
 4763:       if(j==i) continue;
 4764:       jj++;
 4765:       /*ca[0]= k+'a'-1;ca[1]='\0';*/
 4766:       printf("%1d%1d",i,j);
 4767:       fprintf(ficparo,"%1d%1d",i,j);
 4768:       for(k=1; k<=ncovmodel;k++){
 4769: 	/* 	  printf(" %lf",param[i][j][k]); */
 4770: 	/* 	  fprintf(ficparo," %lf",param[i][j][k]); */
 4771: 	printf(" 0.");
 4772: 	fprintf(ficparo," 0.");
 4773:       }
 4774:       printf("\n");
 4775:       fprintf(ficparo,"\n");
 4776:     }
 4777:   }
 4778:   printf("# Scales (for hessian or gradient estimation)\n");
 4779:   fprintf(ficparo,"# Scales (for hessian or gradient estimation)\n");
 4780:   npar= (nlstate+ndeath-1)*nlstate*ncovmodel; /* Number of parameters*/ 
 4781:   for(i=1; i <=nlstate; i++){
 4782:     jj=0;
 4783:     for(j=1; j <=nlstate+ndeath; j++){
 4784:       if(j==i) continue;
 4785:       jj++;
 4786:       fprintf(ficparo,"%1d%1d",i,j);
 4787:       printf("%1d%1d",i,j);
 4788:       fflush(stdout);
 4789:       for(k=1; k<=ncovmodel;k++){
 4790: 	/* 	printf(" %le",delti3[i][j][k]); */
 4791: 	/* 	fprintf(ficparo," %le",delti3[i][j][k]); */
 4792: 	printf(" 0.");
 4793: 	fprintf(ficparo," 0.");
 4794:       }
 4795:       numlinepar++;
 4796:       printf("\n");
 4797:       fprintf(ficparo,"\n");
 4798:     }
 4799:   }
 4800:   printf("# Covariance matrix\n");
 4801: /* # 121 Var(a12)\n\ */
 4802: /* # 122 Cov(b12,a12) Var(b12)\n\ */
 4803: /* # 131 Cov(a13,a12) Cov(a13,b12, Var(a13)\n\ */
 4804: /* # 132 Cov(b13,a12) Cov(b13,b12, Cov(b13,a13) Var(b13)\n\ */
 4805: /* # 212 Cov(a21,a12) Cov(a21,b12, Cov(a21,a13) Cov(a21,b13) Var(a21)\n\ */
 4806: /* # 212 Cov(b21,a12) Cov(b21,b12, Cov(b21,a13) Cov(b21,b13) Cov(b21,a21) Var(b21)\n\ */
 4807: /* # 232 Cov(a23,a12) Cov(a23,b12, Cov(a23,a13) Cov(a23,b13) Cov(a23,a21) Cov(a23,b21) Var(a23)\n\ */
 4808: /* # 232 Cov(b23,a12) Cov(b23,b12) ... Var (b23)\n" */
 4809:   fflush(stdout);
 4810:   fprintf(ficparo,"# Covariance matrix\n");
 4811:   /* # 121 Var(a12)\n\ */
 4812:   /* # 122 Cov(b12,a12) Var(b12)\n\ */
 4813:   /* #   ...\n\ */
 4814:   /* # 232 Cov(b23,a12)  Cov(b23,b12) ... Var (b23)\n" */
 4815:   
 4816:   for(itimes=1;itimes<=2;itimes++){
 4817:     jj=0;
 4818:     for(i=1; i <=nlstate; i++){
 4819:       for(j=1; j <=nlstate+ndeath; j++){
 4820: 	if(j==i) continue;
 4821: 	for(k=1; k<=ncovmodel;k++){
 4822: 	  jj++;
 4823: 	  ca[0]= k+'a'-1;ca[1]='\0';
 4824: 	  if(itimes==1){
 4825: 	    printf("#%1d%1d%d",i,j,k);
 4826: 	    fprintf(ficparo,"#%1d%1d%d",i,j,k);
 4827: 	  }else{
 4828: 	    printf("%1d%1d%d",i,j,k);
 4829: 	    fprintf(ficparo,"%1d%1d%d",i,j,k);
 4830: 	    /* 	printf(" %.5le",matcov[i][j]); */
 4831: 	  }
 4832: 	  ll=0;
 4833: 	  for(li=1;li <=nlstate; li++){
 4834: 	    for(lj=1;lj <=nlstate+ndeath; lj++){
 4835: 	      if(lj==li) continue;
 4836: 	      for(lk=1;lk<=ncovmodel;lk++){
 4837: 		ll++;
 4838: 		if(ll<=jj){
 4839: 		  cb[0]= lk +'a'-1;cb[1]='\0';
 4840: 		  if(ll<jj){
 4841: 		    if(itimes==1){
 4842: 		      printf(" Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 4843: 		      fprintf(ficparo," Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 4844: 		    }else{
 4845: 		      printf(" 0.");
 4846: 		      fprintf(ficparo," 0.");
 4847: 		    }
 4848: 		  }else{
 4849: 		    if(itimes==1){
 4850: 		      printf(" Var(%s%1d%1d)",ca,i,j);
 4851: 		      fprintf(ficparo," Var(%s%1d%1d)",ca,i,j);
 4852: 		    }else{
 4853: 		      printf(" 0.");
 4854: 		      fprintf(ficparo," 0.");
 4855: 		    }
 4856: 		  }
 4857: 		}
 4858: 	      } /* end lk */
 4859: 	    } /* end lj */
 4860: 	  } /* end li */
 4861: 	  printf("\n");
 4862: 	  fprintf(ficparo,"\n");
 4863: 	  numlinepar++;
 4864: 	} /* end k*/
 4865:       } /*end j */
 4866:     } /* end i */
 4867:   } /* end itimes */
 4868: 
 4869: } /* end of prwizard */
 4870: /******************* Gompertz Likelihood ******************************/
 4871: double gompertz(double x[])
 4872: { 
 4873:   double A,B,L=0.0,sump=0.,num=0.;
 4874:   int i,n=0; /* n is the size of the sample */
 4875: 
 4876:   for (i=0;i<=imx-1 ; i++) {
 4877:     sump=sump+weight[i];
 4878:     /*    sump=sump+1;*/
 4879:     num=num+1;
 4880:   }
 4881:  
 4882:  
 4883:   /* for (i=0; i<=imx; i++) 
 4884:      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]);*/
 4885: 
 4886:   for (i=1;i<=imx ; i++)
 4887:     {
 4888:       if (cens[i] == 1 && wav[i]>1)
 4889: 	A=-x[1]/(x[2])*(exp(x[2]*(agecens[i]-agegomp))-exp(x[2]*(ageexmed[i]-agegomp)));
 4890:       
 4891:       if (cens[i] == 0 && wav[i]>1)
 4892: 	A=-x[1]/(x[2])*(exp(x[2]*(agedc[i]-agegomp))-exp(x[2]*(ageexmed[i]-agegomp)))
 4893: 	     +log(x[1]/YEARM)+x[2]*(agedc[i]-agegomp)+log(YEARM);  
 4894:       
 4895:       /*if (wav[i] > 1 && agecens[i] > 15) {*/ /* ??? */
 4896:       if (wav[i] > 1 ) { /* ??? */
 4897: 	L=L+A*weight[i];
 4898: 	/* 	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]);*/
 4899:       }
 4900:     }
 4901: 
 4902:  /*printf("x1=%2.9f x2=%2.9f x3=%2.9f L=%f\n",x[1],x[2],x[3],L);*/
 4903:  
 4904:   return -2*L*num/sump;
 4905: }
 4906: 
 4907: #ifdef GSL
 4908: /******************* Gompertz_f Likelihood ******************************/
 4909: double gompertz_f(const gsl_vector *v, void *params)
 4910: { 
 4911:   double A,B,LL=0.0,sump=0.,num=0.;
 4912:   double *x= (double *) v->data;
 4913:   int i,n=0; /* n is the size of the sample */
 4914: 
 4915:   for (i=0;i<=imx-1 ; i++) {
 4916:     sump=sump+weight[i];
 4917:     /*    sump=sump+1;*/
 4918:     num=num+1;
 4919:   }
 4920:  
 4921:  
 4922:   /* for (i=0; i<=imx; i++) 
 4923:      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]);*/
 4924:   printf("x[0]=%lf x[1]=%lf\n",x[0],x[1]);
 4925:   for (i=1;i<=imx ; i++)
 4926:     {
 4927:       if (cens[i] == 1 && wav[i]>1)
 4928: 	A=-x[0]/(x[1])*(exp(x[1]*(agecens[i]-agegomp))-exp(x[1]*(ageexmed[i]-agegomp)));
 4929:       
 4930:       if (cens[i] == 0 && wav[i]>1)
 4931: 	A=-x[0]/(x[1])*(exp(x[1]*(agedc[i]-agegomp))-exp(x[1]*(ageexmed[i]-agegomp)))
 4932: 	     +log(x[0]/YEARM)+x[1]*(agedc[i]-agegomp)+log(YEARM);  
 4933:       
 4934:       /*if (wav[i] > 1 && agecens[i] > 15) {*/ /* ??? */
 4935:       if (wav[i] > 1 ) { /* ??? */
 4936: 	LL=LL+A*weight[i];
 4937: 	/* 	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]);*/
 4938:       }
 4939:     }
 4940: 
 4941:  /*printf("x1=%2.9f x2=%2.9f x3=%2.9f L=%f\n",x[1],x[2],x[3],L);*/
 4942:   printf("x[0]=%lf x[1]=%lf -2*LL*num/sump=%lf\n",x[0],x[1],-2*LL*num/sump);
 4943:  
 4944:   return -2*LL*num/sump;
 4945: }
 4946: #endif
 4947: 
 4948: /******************* Printing html file ***********/
 4949: void printinghtmlmort(char fileres[], char title[], char datafile[], int firstpass, \
 4950: 		  int lastpass, int stepm, int weightopt, char model[],\
 4951: 		  int imx,  double p[],double **matcov,double agemortsup){
 4952:   int i,k;
 4953: 
 4954:   fprintf(fichtm,"<ul><li><h4>Result files </h4>\n Force of mortality. Parameters of the Gompertz fit (with confidence interval in brackets):<br>");
 4955:   fprintf(fichtm,"  mu(age) =%lf*exp(%lf*(age-%d)) per year<br><br>",p[1],p[2],agegomp);
 4956:   for (i=1;i<=2;i++) 
 4957:     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]));
 4958:   fprintf(fichtm,"<br><br><img src=\"graphmort.png\">");
 4959:   fprintf(fichtm,"</ul>");
 4960: 
 4961: fprintf(fichtm,"<ul><li><h4>Life table</h4>\n <br>");
 4962: 
 4963:  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>");
 4964: 
 4965:  for (k=agegomp;k<(agemortsup-2);k++) 
 4966:    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]);
 4967: 
 4968:  
 4969:   fflush(fichtm);
 4970: }
 4971: 
 4972: /******************* Gnuplot file **************/
 4973: void printinggnuplotmort(char fileres[], char optionfilefiname[], double ageminpar, double agemaxpar, double fage , char pathc[], double p[]){
 4974: 
 4975:   char dirfileres[132],optfileres[132];
 4976: 
 4977:   int ng;
 4978: 
 4979: 
 4980:   /*#ifdef windows */
 4981:   fprintf(ficgp,"cd \"%s\" \n",pathc);
 4982:     /*#endif */
 4983: 
 4984: 
 4985:   strcpy(dirfileres,optionfilefiname);
 4986:   strcpy(optfileres,"vpl");
 4987:   fprintf(ficgp,"set out \"graphmort.png\"\n "); 
 4988:   fprintf(ficgp,"set xlabel \"Age\"\n set ylabel \"Force of mortality (per year)\" \n "); 
 4989:   fprintf(ficgp, "set ter png small size 320, 240\n set log y\n"); 
 4990:   /* fprintf(ficgp, "set size 0.65,0.65\n"); */
 4991:   fprintf(ficgp,"plot [%d:100] %lf*exp(%lf*(x-%d))",agegomp,p[1],p[2],agegomp);
 4992: 
 4993: } 
 4994: 
 4995: int readdata(char datafile[], int firstobs, int lastobs, int *imax)
 4996: {
 4997: 
 4998:   /*-------- data file ----------*/
 4999:   FILE *fic;
 5000:   char dummy[]="                         ";
 5001:   int i=0, j=0, n=0;
 5002:   int linei, month, year,iout;
 5003:   char line[MAXLINE], linetmp[MAXLINE];
 5004:   char stra[MAXLINE], strb[MAXLINE];
 5005:   char *stratrunc;
 5006:   int lstra;
 5007: 
 5008: 
 5009:   if((fic=fopen(datafile,"r"))==NULL)    {
 5010:     printf("Problem while opening datafile: %s\n", datafile);return 1;
 5011:     fprintf(ficlog,"Problem while opening datafile: %s\n", datafile);return 1;
 5012:   }
 5013: 
 5014:   i=1;
 5015:   linei=0;
 5016:   while ((fgets(line, MAXLINE, fic) != NULL) &&((i >= firstobs) && (i <=lastobs))) {
 5017:     linei=linei+1;
 5018:     for(j=strlen(line); j>=0;j--){  /* Untabifies line */
 5019:       if(line[j] == '\t')
 5020: 	line[j] = ' ';
 5021:     }
 5022:     for(j=strlen(line)-1; (line[j]==' ')||(line[j]==10)||(line[j]==13);j--){
 5023:       ;
 5024:     };
 5025:     line[j+1]=0;  /* Trims blanks at end of line */
 5026:     if(line[0]=='#'){
 5027:       fprintf(ficlog,"Comment line\n%s\n",line);
 5028:       printf("Comment line\n%s\n",line);
 5029:       continue;
 5030:     }
 5031:     trimbb(linetmp,line); /* Trims multiple blanks in line */
 5032:     strcpy(line, linetmp);
 5033:   
 5034: 
 5035:     for (j=maxwav;j>=1;j--){
 5036:       cutv(stra, strb, line, ' '); 
 5037:       if(strb[0]=='.') { /* Missing status */
 5038: 	lval=-1;
 5039:       }else{
 5040: 	errno=0;
 5041: 	lval=strtol(strb,&endptr,10); 
 5042:       /*	if (errno == ERANGE && (lval == LONG_MAX || lval == LONG_MIN))*/
 5043: 	if( strb[0]=='\0' || (*endptr != '\0')){
 5044: 	  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);
 5045: 	  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);
 5046: 	  return 1;
 5047: 	}
 5048:       }
 5049:       s[j][i]=lval;
 5050:       
 5051:       strcpy(line,stra);
 5052:       cutv(stra, strb,line,' ');
 5053:       if( (iout=sscanf(strb,"%d/%d",&month, &year)) != 0){
 5054:       }
 5055:       else  if( (iout=sscanf(strb,"%s.",dummy)) != 0){
 5056: 	month=99;
 5057: 	year=9999;
 5058:       }else{
 5059: 	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);
 5060: 	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);
 5061: 	return 1;
 5062:       }
 5063:       anint[j][i]= (double) year; 
 5064:       mint[j][i]= (double)month; 
 5065:       strcpy(line,stra);
 5066:     } /* ENd Waves */
 5067:     
 5068:     cutv(stra, strb,line,' '); 
 5069:     if( (iout=sscanf(strb,"%d/%d",&month, &year)) != 0){
 5070:     }
 5071:     else  if( (iout=sscanf(strb,"%s.",dummy)) != 0){
 5072:       month=99;
 5073:       year=9999;
 5074:     }else{
 5075:       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);
 5076: 	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);
 5077: 	return 1;
 5078:     }
 5079:     andc[i]=(double) year; 
 5080:     moisdc[i]=(double) month; 
 5081:     strcpy(line,stra);
 5082:     
 5083:     cutv(stra, strb,line,' '); 
 5084:     if( (iout=sscanf(strb,"%d/%d",&month, &year)) != 0){
 5085:     }
 5086:     else  if( (iout=sscanf(strb,"%s.", dummy)) != 0){
 5087:       month=99;
 5088:       year=9999;
 5089:     }else{
 5090:       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);
 5091:       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);
 5092: 	return 1;
 5093:     }
 5094:     if (year==9999) {
 5095:       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);
 5096:       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);
 5097: 	return 1;
 5098: 
 5099:     }
 5100:     annais[i]=(double)(year);
 5101:     moisnais[i]=(double)(month); 
 5102:     strcpy(line,stra);
 5103:     
 5104:     cutv(stra, strb,line,' '); 
 5105:     errno=0;
 5106:     dval=strtod(strb,&endptr); 
 5107:     if( strb[0]=='\0' || (*endptr != '\0')){
 5108:       printf("Error reading data around '%f' at line number %d, \"%s\" for individual %d\nShould be a weight.  Exiting.\n",dval, i,line,linei);
 5109:       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);
 5110:       fflush(ficlog);
 5111:       return 1;
 5112:     }
 5113:     weight[i]=dval; 
 5114:     strcpy(line,stra);
 5115:     
 5116:     for (j=ncovcol;j>=1;j--){
 5117:       cutv(stra, strb,line,' '); 
 5118:       if(strb[0]=='.') { /* Missing status */
 5119: 	lval=-1;
 5120:       }else{
 5121: 	errno=0;
 5122: 	lval=strtol(strb,&endptr,10); 
 5123: 	if( strb[0]=='\0' || (*endptr != '\0')){
 5124: 	  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);
 5125: 	  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);
 5126: 	  return 1;
 5127: 	}
 5128:       }
 5129:       if(lval <-1 || lval >1){
 5130: 	printf("Error reading data around '%ld' at line number %d for individual %d, '%s'\n \
 5131:  Should be a value of %d(nth) covariate (0 should be the value for the reference and 1\n \
 5132:  for the alternative. IMaCh does not build design variables automatically, do it yourself.\n \
 5133:  For example, for multinomial values like 1, 2 and 3,\n \
 5134:  build V1=0 V2=0 for the reference value (1),\n \
 5135:         V1=1 V2=0 for (2) \n \
 5136:  and V1=0 V2=1 for (3). V1=1 V2=1 should not exist and the corresponding\n \
 5137:  output of IMaCh is often meaningless.\n \
 5138:  Exiting.\n",lval,linei, i,line,j);
 5139: 	fprintf(ficlog,"Error reading data around '%ld' at line number %d for individual %d, '%s'\n \
 5140:  Should be a value of %d(nth) covariate (0 should be the value for the reference and 1\n \
 5141:  for the alternative. IMaCh does not build design variables automatically, do it yourself.\n \
 5142:  For example, for multinomial values like 1, 2 and 3,\n \
 5143:  build V1=0 V2=0 for the reference value (1),\n \
 5144:         V1=1 V2=0 for (2) \n \
 5145:  and V1=0 V2=1 for (3). V1=1 V2=1 should not exist and the corresponding\n \
 5146:  output of IMaCh is often meaningless.\n \
 5147:  Exiting.\n",lval,linei, i,line,j);fflush(ficlog);
 5148: 	return 1;
 5149:       }
 5150:       covar[j][i]=(double)(lval);
 5151:       strcpy(line,stra);
 5152:     }  
 5153:     lstra=strlen(stra);
 5154:      
 5155:     if(lstra > 9){ /* More than 2**32 or max of what printf can write with %ld */
 5156:       stratrunc = &(stra[lstra-9]);
 5157:       num[i]=atol(stratrunc);
 5158:     }
 5159:     else
 5160:       num[i]=atol(stra);
 5161:     /*if((s[2][i]==2) && (s[3][i]==-1)&&(s[4][i]==9)){
 5162:       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;}*/
 5163:     
 5164:     i=i+1;
 5165:   } /* End loop reading  data */
 5166: 
 5167:   *imax=i-1; /* Number of individuals */
 5168:   fclose(fic);
 5169:  
 5170:   return (0);
 5171:   /* endread: */
 5172:     printf("Exiting readdata: ");
 5173:     fclose(fic);
 5174:     return (1);
 5175: 
 5176: 
 5177: 
 5178: }
 5179: void removespace(char *str) {
 5180:   char *p1 = str, *p2 = str;
 5181:   do
 5182:     while (*p2 == ' ')
 5183:       p2++;
 5184:   while (*p1++ == *p2++);
 5185: }
 5186: 
 5187: int decodemodel ( char model[], int lastobs) /**< This routine decode the model and returns:
 5188:    * Model  V1+V2+V3+V8+V7*V8+V5*V6+V8*age+V3*age
 5189:    * - cptcovt total number of covariates of the model nbocc(+)+1 = 8
 5190:    * - cptcovn or number of covariates k of the models excluding age*products =6
 5191:    * - cptcovage number of covariates with age*products =2
 5192:    * - cptcovs number of simple covariates
 5193:    * - 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
 5194:    *     which is a new column after the 9 (ncovcol) variables. 
 5195:    * - if k is a product Vn*Vm covar[k][i] is filled with correct values for each individual
 5196:    * - Tprod[l] gives the kth covariates of the product Vn*Vm l=1 to cptcovprod-cptcovage
 5197:    *    Tprod[1]@2 {5, 6}: position of first product V7*V8 is 5, and second V5*V6 is 6.
 5198:    * - Tvard[k]  p Tvard[1][1]@4 {7, 8, 5, 6} for V7*V8 and V5*V6 .
 5199:  */
 5200: {
 5201:   int i, j, k, ks;
 5202:   int  j1, k1, k2;
 5203:   char modelsav[80];
 5204:   char stra[80], strb[80], strc[80], strd[80],stre[80];
 5205: 
 5206:   /*removespace(model);*/
 5207:   if (strlen(model) >1){ /* If there is at least 1 covariate */
 5208:     j=0, j1=0, k1=0, k2=-1, ks=0, cptcovn=0;
 5209:     j=nbocc(model,'+'); /**< j=Number of '+' */
 5210:     j1=nbocc(model,'*'); /**< j1=Number of '*' */
 5211:     cptcovs=j+1-j1; /**<  Number of simple covariates V1+V2*age+V3 +V3*V4=> V1 + V3 =2  */
 5212:     cptcovt= j+1; /* Number of total covariates in the model V1 + V2*age+ V3 + V3*V4=> 4*/
 5213:                   /* including age products which are counted in cptcovage.
 5214: 		  * but the covariates which are products must be treated separately: ncovn=4- 2=2 (V1+V3). */
 5215:     cptcovprod=j1; /**< Number of products  V1*V2 +v3*age = 2 */
 5216:     cptcovprodnoage=0; /**< Number of covariate products without age: V3*V4 =1  */
 5217:     strcpy(modelsav,model); 
 5218:     if (strstr(model,"AGE") !=0){
 5219:       printf("Error. AGE must be in lower case 'age' model=%s ",model);
 5220:       fprintf(ficlog,"Error. AGE must be in lower case model=%s ",model);fflush(ficlog);
 5221:       return 1;
 5222:     }
 5223:     if (strstr(model,"v") !=0){
 5224:       printf("Error. 'v' must be in upper case 'V' model=%s ",model);
 5225:       fprintf(ficlog,"Error. 'v' must be in upper case model=%s ",model);fflush(ficlog);
 5226:       return 1;
 5227:     }
 5228:     
 5229:     /*   Design
 5230:      *  V1   V2   V3   V4  V5  V6  V7  V8  V9 Weight
 5231:      *  <          ncovcol=8                >
 5232:      * Model V2 + V1 + V3*age + V3 + V5*V6 + V7*V8 + V8*age + V8
 5233:      *   k=  1    2      3       4     5       6      7        8
 5234:      *  cptcovn number of covariates (not including constant and age ) = # of + plus 1 = 7+1=8
 5235:      *  covar[k,i], value of kth covariate if not including age for individual i:
 5236:      *       covar[1][i]= (V2), covar[4][i]=(V3), covar[8][i]=(V8)
 5237:      *  Tvar[k] # of the kth covariate:  Tvar[1]=2  Tvar[4]=3 Tvar[8]=8
 5238:      *       if multiplied by age: V3*age Tvar[3=V3*age]=3 (V3) Tvar[7]=8 and 
 5239:      *  Tage[++cptcovage]=k
 5240:      *       if products, new covar are created after ncovcol with k1
 5241:      *  Tvar[k]=ncovcol+k1; # of the kth covariate product:  Tvar[5]=ncovcol+1=10  Tvar[6]=ncovcol+1=11
 5242:      *  Tprod[k1]=k; Tprod[1]=5 Tprod[2]= 6; gives the position of the k1th product
 5243:      *  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
 5244:      *  Tvar[cptcovn+k2]=Tvard[k1][1];Tvar[cptcovn+k2+1]=Tvard[k1][2];
 5245:      *  Tvar[8+1]=5;Tvar[8+2]=6;Tvar[8+3]=7;Tvar[8+4]=8 inverted
 5246:      *  V1   V2   V3   V4  V5  V6  V7  V8  V9  V10  V11
 5247:      *  <          ncovcol=8                >
 5248:      *       Model V2 + V1 + V3*age + V3 + V5*V6 + V7*V8 + V8*age + V8    d1   d1   d2  d2
 5249:      *          k=  1    2      3       4     5       6      7        8    9   10   11  12
 5250:      *     Tvar[k]= 2    1      3       3    10      11      8        8    5    6    7   8
 5251:      * p Tvar[1]@12={2,   1,     3,      3,   11,     10,     8,       8,   7,   8,   5,  6}
 5252:      * p Tprod[1]@2={                         6, 5}
 5253:      *p Tvard[1][1]@4= {7, 8, 5, 6}
 5254:      * covar[k][i]= V2   V1      ?      V3    V5*V6?   V7*V8?  ?       V8   
 5255:      *  cov[Tage[kk]+2]=covar[Tvar[Tage[kk]]][i]*cov[2];
 5256:      *How to reorganize?
 5257:      * Model V1 + V2 + V3 + V8 + V5*V6 + V7*V8 + V3*age + V8*age
 5258:      * Tvars {2,   1,     3,      3,   11,     10,     8,       8,   7,   8,   5,  6}
 5259:      *       {2,   1,     4,      8,    5,      6,     3,       7}
 5260:      * Struct []
 5261:      */
 5262: 
 5263:     /* This loop fills the array Tvar from the string 'model'.*/
 5264:     /* j is the number of + signs in the model V1+V2+V3 j=2 i=3 to 1 */
 5265:     /*   modelsav=V2+V1+V4+age*V3 strb=age*V3 stra=V2+V1+V4  */
 5266:     /* 	k=4 (age*V3) Tvar[k=4]= 3 (from V3) Tage[cptcovage=1]=4 */
 5267:     /* 	k=3 V4 Tvar[k=3]= 4 (from V4) */
 5268:     /* 	k=2 V1 Tvar[k=2]= 1 (from V1) */
 5269:     /* 	k=1 Tvar[1]=2 (from V2) */
 5270:     /* 	k=5 Tvar[5] */
 5271:     /* for (k=1; k<=cptcovn;k++) { */
 5272:     /* 	cov[2+k]=nbcode[Tvar[k]][codtab[ij][Tvar[k]]]; */
 5273:     /* 	} */
 5274:     /* for (k=1; k<=cptcovage;k++) cov[2+Tage[k]]=cov[2+Tage[k]]*cov[2]; */
 5275:     /*
 5276:      * Treating invertedly V2+V1+V3*age+V2*V4 is as if written V2*V4 +V3*age + V1 + V2 */
 5277:     for(k=cptcovt; k>=1;k--) /**< Number of covariates */
 5278:         Tvar[k]=0;
 5279:     cptcovage=0;
 5280:     for(k=1; k<=cptcovt;k++){ /* Loop on total covariates of the model */
 5281:       cutl(stra,strb,modelsav,'+'); /* keeps in strb after the first '+' 
 5282: 				     modelsav==V2+V1+V4+V3*age strb=V3*age stra=V2+V1+V4 */ 
 5283:       if (nbocc(modelsav,'+')==0) strcpy(strb,modelsav); /* and analyzes it */
 5284:       /*      printf("i=%d a=%s b=%s sav=%s\n",i, stra,strb,modelsav);*/
 5285:       /*scanf("%d",i);*/
 5286:       if (strchr(strb,'*')) {  /**< Model includes a product V2+V1+V4+V3*age strb=V3*age */
 5287: 	cutl(strc,strd,strb,'*'); /**< strd*strc  Vm*Vn: strb=V3*age(input) strc=age strd=V3 ; V3*V2 strc=V2, strd=V3 */
 5288: 	if (strcmp(strc,"age")==0) { /**< Model includes age: Vn*age */
 5289: 	  /* covar is not filled and then is empty */
 5290: 	  cptcovprod--;
 5291: 	  cutl(stre,strb,strd,'V'); /* strd=V3(input): stre="3" */
 5292: 	  Tvar[k]=atoi(stre);  /* V2+V1+V4+V3*age Tvar[4]=3 ; V1+V2*age Tvar[2]=2 */
 5293: 	  cptcovage++; /* Sums the number of covariates which include age as a product */
 5294: 	  Tage[cptcovage]=k;  /* Tage[1] = 4 */
 5295: 	  /*printf("stre=%s ", stre);*/
 5296: 	} else if (strcmp(strd,"age")==0) { /* or age*Vn */
 5297: 	  cptcovprod--;
 5298: 	  cutl(stre,strb,strc,'V');
 5299: 	  Tvar[k]=atoi(stre);
 5300: 	  cptcovage++;
 5301: 	  Tage[cptcovage]=k;
 5302: 	} else {  /* Age is not in the model product V2+V1+V1*V4+V3*age+V3*V2  strb=V3*V2*/
 5303: 	  /* loops on k1=1 (V3*V2) and k1=2 V4*V3 */
 5304: 	  cptcovn++;
 5305: 	  cptcovprodnoage++;k1++;
 5306: 	  cutl(stre,strb,strc,'V'); /* strc= Vn, stre is n; strb=V3*V2 stre=3 strc=*/
 5307: 	  Tvar[k]=ncovcol+k1; /* For model-covariate k tells which data-covariate to use but
 5308: 				  because this model-covariate is a construction we invent a new column
 5309: 				  ncovcol + k1
 5310: 				  If already ncovcol=4 and model=V2+V1+V1*V4+age*V3+V3*V2
 5311: 				  Tvar[3=V1*V4]=4+1 Tvar[5=V3*V2]=4 + 2= 6, etc */
 5312: 	  cutl(strc,strb,strd,'V'); /* strd was Vm, strc is m */
 5313: 	  Tprod[k1]=k;  /* Tprod[1]=3(=V1*V4) for V2+V1+V1*V4+age*V3+V3*V2  */
 5314: 	  Tvard[k1][1] =atoi(strc); /* m 1 for V1*/
 5315: 	  Tvard[k1][2] =atoi(stre); /* n 4 for V4*/
 5316: 	  k2=k2+2;
 5317: 	  Tvar[cptcovt+k2]=Tvard[k1][1]; /* Tvar[(cptcovt=4+k2=1)=5]= 1 (V1) */
 5318: 	  Tvar[cptcovt+k2+1]=Tvard[k1][2];  /* Tvar[(cptcovt=4+(k2=1)+1)=6]= 4 (V4) */
 5319: 	  for (i=1; i<=lastobs;i++){
 5320: 	    /* Computes the new covariate which is a product of
 5321: 	       covar[n][i]* covar[m][i] and stores it at ncovol+k1 May not be defined */
 5322: 	    covar[ncovcol+k1][i]=covar[atoi(stre)][i]*covar[atoi(strc)][i];
 5323: 	  }
 5324: 	} /* End age is not in the model */
 5325:       } /* End if model includes a product */
 5326:       else { /* no more sum */
 5327: 	/*printf("d=%s c=%s b=%s\n", strd,strc,strb);*/
 5328:        /*  scanf("%d",i);*/
 5329: 	cutl(strd,strc,strb,'V');
 5330: 	ks++; /**< Number of simple covariates */
 5331: 	cptcovn++;
 5332: 	Tvar[k]=atoi(strd);
 5333:       }
 5334:       strcpy(modelsav,stra);  /* modelsav=V2+V1+V4 stra=V2+V1+V4 */ 
 5335:       /*printf("a=%s b=%s sav=%s\n", stra,strb,modelsav);
 5336: 	scanf("%d",i);*/
 5337:     } /* end of loop + */
 5338:   } /* end model */
 5339:   
 5340:   /*The number n of Vn is stored in Tvar. cptcovage =number of age covariate. Tage gives the position of age. cptcovprod= number of products.
 5341:     If model=V1+V1*age then Tvar[1]=1 Tvar[2]=1 cptcovage=1 Tage[1]=2 cptcovprod=0*/
 5342: 
 5343:   /* printf("tvar1=%d tvar2=%d tvar3=%d cptcovage=%d Tage=%d",Tvar[1],Tvar[2],Tvar[3],cptcovage,Tage[1]);
 5344:   printf("cptcovprod=%d ", cptcovprod);
 5345:   fprintf(ficlog,"cptcovprod=%d ", cptcovprod);
 5346: 
 5347:   scanf("%d ",i);*/
 5348: 
 5349: 
 5350:   return (0); /* with covar[new additional covariate if product] and Tage if age */ 
 5351:   /*endread:*/
 5352:     printf("Exiting decodemodel: ");
 5353:     return (1);
 5354: }
 5355: 
 5356: int calandcheckages(int imx, int maxwav, double *agemin, double *agemax, int *nberr, int *nbwarn )
 5357: {
 5358:   int i, m;
 5359: 
 5360:   for (i=1; i<=imx; i++) {
 5361:     for(m=2; (m<= maxwav); m++) {
 5362:       if (((int)mint[m][i]== 99) && (s[m][i] <= nlstate)){
 5363: 	anint[m][i]=9999;
 5364: 	s[m][i]=-1;
 5365:       }
 5366:       if((int)moisdc[i]==99 && (int)andc[i]==9999 && s[m][i]>nlstate){
 5367: 	*nberr = *nberr + 1;
 5368: 	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 (%d)\n",(int)moisdc[i],(int)andc[i],num[i],i, *nberr);
 5369: 	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 (%d)\n",(int)moisdc[i],(int)andc[i],num[i],i, *nberr);
 5370: 	s[m][i]=-1;
 5371:       }
 5372:       if((int)moisdc[i]==99 && (int)andc[i]!=9999 && s[m][i]>nlstate){
 5373: 	(*nberr)++;
 5374: 	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]); 
 5375: 	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]); 
 5376: 	s[m][i]=-1; /* We prefer to skip it (and to skip it in version 0.8a1 too */
 5377:       }
 5378:     }
 5379:   }
 5380: 
 5381:   for (i=1; i<=imx; i++)  {
 5382:     agedc[i]=(moisdc[i]/12.+andc[i])-(moisnais[i]/12.+annais[i]);
 5383:     for(m=firstpass; (m<= lastpass); m++){
 5384:       if(s[m][i] >0 || s[m][i]==-2 || s[m][i]==-4 || s[m][i]==-5){
 5385: 	if (s[m][i] >= nlstate+1) {
 5386: 	  if(agedc[i]>0){
 5387: 	    if((int)moisdc[i]!=99 && (int)andc[i]!=9999){
 5388: 	      agev[m][i]=agedc[i];
 5389: 	  /*if(moisdc[i]==99 && andc[i]==9999) s[m][i]=-1;*/
 5390: 	    }else {
 5391: 	      if ((int)andc[i]!=9999){
 5392: 		nbwarn++;
 5393: 		printf("Warning negative age at death: %ld line:%d\n",num[i],i);
 5394: 		fprintf(ficlog,"Warning negative age at death: %ld line:%d\n",num[i],i);
 5395: 		agev[m][i]=-1;
 5396: 	      }
 5397: 	    }
 5398: 	  } /* agedc > 0 */
 5399: 	}
 5400: 	else if(s[m][i] !=9){ /* Standard case, age in fractional
 5401: 				 years but with the precision of a month */
 5402: 	  agev[m][i]=(mint[m][i]/12.+1./24.+anint[m][i])-(moisnais[i]/12.+1./24.+annais[i]);
 5403: 	  if((int)mint[m][i]==99 || (int)anint[m][i]==9999)
 5404: 	    agev[m][i]=1;
 5405: 	  else if(agev[m][i] < *agemin){ 
 5406: 	    *agemin=agev[m][i];
 5407: 	    printf(" Min anint[%d][%d]=%.2f annais[%d]=%.2f, agemin=%.2f\n",m,i,anint[m][i], i,annais[i], *agemin);
 5408: 	  }
 5409: 	  else if(agev[m][i] >*agemax){
 5410: 	    *agemax=agev[m][i];
 5411: 	    /* printf(" Max anint[%d][%d]=%.0f annais[%d]=%.0f, agemax=%.2f\n",m,i,anint[m][i], i,annais[i], *agemax);*/
 5412: 	  }
 5413: 	  /*agev[m][i]=anint[m][i]-annais[i];*/
 5414: 	  /*	 agev[m][i] = age[i]+2*m;*/
 5415: 	}
 5416: 	else { /* =9 */
 5417: 	  agev[m][i]=1;
 5418: 	  s[m][i]=-1;
 5419: 	}
 5420:       }
 5421:       else /*= 0 Unknown */
 5422: 	agev[m][i]=1;
 5423:     }
 5424:     
 5425:   }
 5426:   for (i=1; i<=imx; i++)  {
 5427:     for(m=firstpass; (m<=lastpass); m++){
 5428:       if (s[m][i] > (nlstate+ndeath)) {
 5429: 	(*nberr)++;
 5430: 	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);	
 5431: 	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);	
 5432: 	return 1;
 5433:       }
 5434:     }
 5435:   }
 5436: 
 5437:   /*for (i=1; i<=imx; i++){
 5438:   for (m=firstpass; (m<lastpass); m++){
 5439:      printf("%ld %d %.lf %d %d\n", num[i],(covar[1][i]),agev[m][i],s[m][i],s[m+1][i]);
 5440: }
 5441: 
 5442: }*/
 5443: 
 5444: 
 5445:   printf("Total number of individuals= %d, Agemin = %.2f, Agemax= %.2f\n\n", imx, *agemin, *agemax);
 5446:   fprintf(ficlog,"Total number of individuals= %d, Agemin = %.2f, Agemax= %.2f\n\n", imx, *agemin, *agemax); 
 5447: 
 5448:   return (0);
 5449:  /* endread:*/
 5450:     printf("Exiting calandcheckages: ");
 5451:     return (1);
 5452: }
 5453: 
 5454: void syscompilerinfo()
 5455:  {
 5456:    /* #include "syscompilerinfo.h"*/
 5457:    /* #include <gnu/libc-version.h> */ /* Only on gnu */
 5458: #include <stdint.h>
 5459:    printf("Compiled with:");fprintf(ficlog,"Compiled with:");
 5460: #if defined(__clang__)
 5461:    printf(" Clang/LLVM");fprintf(ficlog," Clang/LLVM");	/* Clang/LLVM. ---------------------------------------------- */
 5462: #endif
 5463: #if defined(__ICC) || defined(__INTEL_COMPILER)
 5464:    printf(" Intel ICC/ICPC");fprintf(ficlog," Intel ICC/ICPC");/* Intel ICC/ICPC. ------------------------------------------ */
 5465: #endif
 5466: #if defined(__GNUC__) || defined(__GNUG__)
 5467:    printf(" GNU GCC/G++");fprintf(ficlog," GNU GCC/G++");/* GNU GCC/G++. --------------------------------------------- */
 5468: #endif
 5469: #if defined(__HP_cc) || defined(__HP_aCC)
 5470:    printf(" Hewlett-Packard C/aC++");fprintf(fcilog," Hewlett-Packard C/aC++"); /* Hewlett-Packard C/aC++. ---------------------------------- */
 5471: #endif
 5472: #if defined(__IBMC__) || defined(__IBMCPP__)
 5473:    printf(" IBM XL C/C++"); fprintf(ficlog," IBM XL C/C++");/* IBM XL C/C++. -------------------------------------------- */
 5474: #endif
 5475: #if defined(_MSC_VER)
 5476:    printf(" Microsoft Visual Studio");fprintf(ficlog," Microsoft Visual Studio");/* Microsoft Visual Studio. --------------------------------- */
 5477: #endif
 5478: #if defined(__PGI)
 5479:    printf(" Portland Group PGCC/PGCPP");fprintf(ficlog," Portland Group PGCC/PGCPP");/* Portland Group PGCC/PGCPP. ------------------------------- */
 5480: #endif
 5481: #if defined(__SUNPRO_C) || defined(__SUNPRO_CC)
 5482:    printf(" Oracle Solaris Studio");fprintf(ficlog," Oracle Solaris Studio\n");/* Oracle Solaris Studio. ----------------------------------- */
 5483: #endif
 5484:    printf(". ");fprintf(ficlog,". ");
 5485:    
 5486: // http://stackoverflow.com/questions/4605842/how-to-identify-platform-compiler-from-preprocessor-macros
 5487: #ifdef _WIN32 // note the underscore: without it, it's not msdn official!
 5488:     // Windows (x64 and x86)
 5489: #elif __unix__ // all unices, not all compilers
 5490:     // Unix
 5491: #elif __linux__
 5492:     // linux
 5493: #elif __APPLE__
 5494:     // Mac OS, not sure if this is covered by __posix__ and/or __unix__ though...
 5495: #endif
 5496: 
 5497: /*  __MINGW32__	  */
 5498: /*  __CYGWIN__	 */
 5499: /* __MINGW64__  */
 5500: // http://msdn.microsoft.com/en-us/library/b0084kay.aspx
 5501: /* _MSC_VER  //the Visual C++ compiler is 17.00.51106.1, the _MSC_VER macro evaluates to 1700. Type cl /?  */
 5502: /* _MSC_FULL_VER //the Visual C++ compiler is 15.00.20706.01, the _MSC_FULL_VER macro evaluates to 150020706 */
 5503: /* _WIN64  // Defined for applications for Win64. */
 5504: /* _M_X64 // Defined for compilations that target x64 processors. */
 5505: /* _DEBUG // Defined when you compile with /LDd, /MDd, and /MTd. */
 5506: 
 5507: #if UINTPTR_MAX == 0xffffffff
 5508:    printf(" 32-bit.\n"); fprintf(ficlog," 32-bit.\n");/* 32-bit */
 5509: #elif UINTPTR_MAX == 0xffffffffffffffff
 5510:    printf(" 64-bit.\n"); fprintf(ficlog," 64-bit.\n");/* 64-bit */
 5511: #else
 5512:    printf(" wtf-bit.\n"); fprintf(ficlog," wtf-bit.\n");/* wtf */
 5513: #endif
 5514: 
 5515: /* struct utsname sysInfo;
 5516: 
 5517:    if (uname(&sysInfo) != -1) {
 5518:      printf(" %s %s %s %s %s\n",sysInfo.sysname, sysInfo.nodename, sysInfo.release, sysInfo.version, sysInfo.machine);
 5519:      fprintf(ficlog," %s %s %s %s %s\n ",sysInfo.sysname, sysInfo.nodename, sysInfo.release, sysInfo.version, sysInfo.machine);
 5520:    }
 5521:    else
 5522:       perror("uname() error");
 5523: 	  */
 5524: #if defined(__GNUC__)
 5525: # if defined(__GNUC_PATCHLEVEL__)
 5526: #  define __GNUC_VERSION__ (__GNUC__ * 10000 \
 5527:                             + __GNUC_MINOR__ * 100 \
 5528:                             + __GNUC_PATCHLEVEL__)
 5529: # else
 5530: #  define __GNUC_VERSION__ (__GNUC__ * 10000 \
 5531:                             + __GNUC_MINOR__ * 100)
 5532: # endif
 5533:    printf("GNU C version %d.\n", __GNUC_VERSION__);
 5534:    fprintf(ficlog, "GNU C version %d.\n", __GNUC_VERSION__);
 5535: #endif
 5536: #if defined(_MSC_VER)
 5537:    /*printf("Visual C++ compiler: %s \n;", _MSC_FULL_VER);*/
 5538:    /*fprintf(ficlog, "Visual C++ compiler: %s \n;", _MSC_FULL_VER);*/
 5539: #endif
 5540:    
 5541:   /* printf("GNU libc version: %s\n", gnu_get_libc_version()); */
 5542: 
 5543:  }
 5544: 
 5545: /***********************************************/
 5546: /**************** Main Program *****************/
 5547: /***********************************************/
 5548: 
 5549: int main(int argc, char *argv[])
 5550: {
 5551: #ifdef GSL
 5552:   const gsl_multimin_fminimizer_type *T;
 5553:   size_t iteri = 0, it;
 5554:   int rval = GSL_CONTINUE;
 5555:   int status = GSL_SUCCESS;
 5556:   double ssval;
 5557: #endif
 5558:   int movingaverage(double ***probs, double bage,double fage, double ***mobaverage, int mobilav);
 5559:   int i,j, k, n=MAXN,iter=0,m,size=100, cptcod;
 5560: 
 5561:   int jj, ll, li, lj, lk;
 5562:   int numlinepar=0; /* Current linenumber of parameter file */
 5563:   int itimes;
 5564:   int NDIM=2;
 5565:   int vpopbased=0;
 5566: 
 5567:   char ca[32], cb[32];
 5568:   /*  FILE *fichtm; *//* Html File */
 5569:   /* FILE *ficgp;*/ /*Gnuplot File */
 5570:   struct stat info;
 5571:   double agedeb;
 5572:   double ageminpar=1.e20,agemin=1.e20, agemaxpar=-1.e20, agemax=-1.e20;
 5573: 
 5574:   double fret;
 5575:   double dum; /* Dummy variable */
 5576:   double ***p3mat;
 5577:   double ***mobaverage;
 5578: 
 5579:   char line[MAXLINE];
 5580:   char path[MAXLINE],pathc[MAXLINE],pathcd[MAXLINE],pathtot[MAXLINE],model[MAXLINE];
 5581:   char pathr[MAXLINE], pathimach[MAXLINE]; 
 5582:   char *tok, *val; /* pathtot */
 5583:   int firstobs=1, lastobs=10;
 5584:   int c,  h , cpt;
 5585:   int jl;
 5586:   int i1, j1, jk, stepsize;
 5587:   int *tab; 
 5588:   int mobilavproj=0 , prevfcast=0 ; /* moving average of prev, If prevfcast=1 prevalence projection */
 5589:   int mobilav=0,popforecast=0;
 5590:   int hstepm, nhstepm;
 5591:   int agemortsup;
 5592:   float  sumlpop=0.;
 5593:   double jprev1=1, mprev1=1,anprev1=2000,jprev2=1, mprev2=1,anprev2=2000;
 5594:   double jpyram=1, mpyram=1,anpyram=2000,jpyram1=1, mpyram1=1,anpyram1=2000;
 5595: 
 5596:   double bage=0, fage=110, age, agelim, agebase;
 5597:   double ftolpl=FTOL;
 5598:   double **prlim;
 5599:   double ***param; /* Matrix of parameters */
 5600:   double  *p;
 5601:   double **matcov; /* Matrix of covariance */
 5602:   double ***delti3; /* Scale */
 5603:   double *delti; /* Scale */
 5604:   double ***eij, ***vareij;
 5605:   double **varpl; /* Variances of prevalence limits by age */
 5606:   double *epj, vepp;
 5607: 
 5608:   double dateprev1, dateprev2,jproj1=1,mproj1=1,anproj1=2000,jproj2=1,mproj2=1,anproj2=2000;
 5609:   double **ximort;
 5610:   char *alph[]={"a","a","b","c","d","e"}, str[4]="1234";
 5611:   int *dcwave;
 5612: 
 5613:   char z[1]="c";
 5614: 
 5615:   /*char  *strt;*/
 5616:   char strtend[80];
 5617: 
 5618: 
 5619: /*   setlocale (LC_ALL, ""); */
 5620: /*   bindtextdomain (PACKAGE, LOCALEDIR); */
 5621: /*   textdomain (PACKAGE); */
 5622: /*   setlocale (LC_CTYPE, ""); */
 5623: /*   setlocale (LC_MESSAGES, ""); */
 5624: 
 5625:   /*   gettimeofday(&start_time, (struct timezone*)0); */ /* at first time */
 5626:   rstart_time = time(NULL);  
 5627:   /*  (void) gettimeofday(&start_time,&tzp);*/
 5628:   start_time = *localtime(&rstart_time);
 5629:   curr_time=start_time;
 5630:   /*tml = *localtime(&start_time.tm_sec);*/
 5631:   /* strcpy(strstart,asctime(&tml)); */
 5632:   strcpy(strstart,asctime(&start_time));
 5633: 
 5634: /*  printf("Localtime (at start)=%s",strstart); */
 5635: /*  tp.tm_sec = tp.tm_sec +86400; */
 5636: /*  tm = *localtime(&start_time.tm_sec); */
 5637: /*   tmg.tm_year=tmg.tm_year +dsign*dyear; */
 5638: /*   tmg.tm_mon=tmg.tm_mon +dsign*dmonth; */
 5639: /*   tmg.tm_hour=tmg.tm_hour + 1; */
 5640: /*   tp.tm_sec = mktime(&tmg); */
 5641: /*   strt=asctime(&tmg); */
 5642: /*   printf("Time(after) =%s",strstart);  */
 5643: /*  (void) time (&time_value);
 5644: *  printf("time=%d,t-=%d\n",time_value,time_value-86400);
 5645: *  tm = *localtime(&time_value);
 5646: *  strstart=asctime(&tm);
 5647: *  printf("tim_value=%d,asctime=%s\n",time_value,strstart); 
 5648: */
 5649: 
 5650:   nberr=0; /* Number of errors and warnings */
 5651:   nbwarn=0;
 5652:   getcwd(pathcd, size);
 5653: 
 5654:   printf("\n%s\n%s",version,fullversion);
 5655:   if(argc <=1){
 5656:     printf("\nEnter the parameter file name: ");
 5657:     fgets(pathr,FILENAMELENGTH,stdin);
 5658:     i=strlen(pathr);
 5659:     if(pathr[i-1]=='\n')
 5660:       pathr[i-1]='\0';
 5661:     i=strlen(pathr);
 5662:     if(pathr[i-1]==' ') /* This may happen when dragging on oS/X! */
 5663:       pathr[i-1]='\0';
 5664:    for (tok = pathr; tok != NULL; ){
 5665:       printf("Pathr |%s|\n",pathr);
 5666:       while ((val = strsep(&tok, "\"" )) != NULL && *val == '\0');
 5667:       printf("val= |%s| pathr=%s\n",val,pathr);
 5668:       strcpy (pathtot, val);
 5669:       if(pathr[0] == '\0') break; /* Dirty */
 5670:     }
 5671:   }
 5672:   else{
 5673:     strcpy(pathtot,argv[1]);
 5674:   }
 5675:   /*if(getcwd(pathcd, MAXLINE)!= NULL)printf ("Error pathcd\n");*/
 5676:   /*cygwin_split_path(pathtot,path,optionfile);
 5677:     printf("pathtot=%s, path=%s, optionfile=%s\n",pathtot,path,optionfile);*/
 5678:   /* cutv(path,optionfile,pathtot,'\\');*/
 5679: 
 5680:   /* Split argv[0], imach program to get pathimach */
 5681:   printf("\nargv[0]=%s argv[1]=%s, \n",argv[0],argv[1]);
 5682:   split(argv[0],pathimach,optionfile,optionfilext,optionfilefiname);
 5683:   printf("\nargv[0]=%s pathimach=%s, \noptionfile=%s \noptionfilext=%s \noptionfilefiname=%s\n",argv[0],pathimach,optionfile,optionfilext,optionfilefiname);
 5684:  /*   strcpy(pathimach,argv[0]); */
 5685:   /* Split argv[1]=pathtot, parameter file name to get path, optionfile, extension and name */
 5686:   split(pathtot,path,optionfile,optionfilext,optionfilefiname);
 5687:   printf("\npathtot=%s,\npath=%s,\noptionfile=%s \noptionfilext=%s \noptionfilefiname=%s\n",pathtot,path,optionfile,optionfilext,optionfilefiname);
 5688:   chdir(path); /* Can be a relative path */
 5689:   if(getcwd(pathcd,MAXLINE) > 0) /* So pathcd is the full path */
 5690:     printf("Current directory %s!\n",pathcd);
 5691:   strcpy(command,"mkdir ");
 5692:   strcat(command,optionfilefiname);
 5693:   if((outcmd=system(command)) != 0){
 5694:     printf("Directory already exists (or can't create it) %s%s, err=%d\n",path,optionfilefiname,outcmd);
 5695:     /* fprintf(ficlog,"Problem creating directory %s%s\n",path,optionfilefiname); */
 5696:     /* fclose(ficlog); */
 5697: /*     exit(1); */
 5698:   }
 5699: /*   if((imk=mkdir(optionfilefiname))<0){ */
 5700: /*     perror("mkdir"); */
 5701: /*   } */
 5702: 
 5703:   /*-------- arguments in the command line --------*/
 5704: 
 5705:   /* Log file */
 5706:   strcat(filelog, optionfilefiname);
 5707:   strcat(filelog,".log");    /* */
 5708:   if((ficlog=fopen(filelog,"w"))==NULL)    {
 5709:     printf("Problem with logfile %s\n",filelog);
 5710:     goto end;
 5711:   }
 5712:   fprintf(ficlog,"Log filename:%s\n",filelog);
 5713:   fprintf(ficlog,"\n%s\n%s",version,fullversion);
 5714:   fprintf(ficlog,"\nEnter the parameter file name: \n");
 5715:   fprintf(ficlog,"pathimach=%s\npathtot=%s\n\
 5716:  path=%s \n\
 5717:  optionfile=%s\n\
 5718:  optionfilext=%s\n\
 5719:  optionfilefiname='%s'\n",pathimach,pathtot,path,optionfile,optionfilext,optionfilefiname);
 5720: 
 5721:   syscompilerinfo();
 5722: 
 5723:   printf("Local time (at start):%s",strstart);
 5724:   fprintf(ficlog,"Local time (at start): %s",strstart);
 5725:   fflush(ficlog);
 5726: /*   (void) gettimeofday(&curr_time,&tzp); */
 5727: /*   printf("Elapsed time %d\n", asc_diff_time(curr_time.tm_sec-start_time.tm_sec,tmpout)); */
 5728: 
 5729:   /* */
 5730:   strcpy(fileres,"r");
 5731:   strcat(fileres, optionfilefiname);
 5732:   strcat(fileres,".txt");    /* Other files have txt extension */
 5733: 
 5734:   /*---------arguments file --------*/
 5735: 
 5736:   if((ficpar=fopen(optionfile,"r"))==NULL)    {
 5737:     printf("Problem with optionfile '%s' with errno='%s'\n",optionfile,strerror(errno));
 5738:     fprintf(ficlog,"Problem with optionfile '%s' with errno='%s'\n",optionfile,strerror(errno));
 5739:     fflush(ficlog);
 5740:     /* goto end; */
 5741:     exit(70); 
 5742:   }
 5743: 
 5744: 
 5745: 
 5746:   strcpy(filereso,"o");
 5747:   strcat(filereso,fileres);
 5748:   if((ficparo=fopen(filereso,"w"))==NULL) { /* opened on subdirectory */
 5749:     printf("Problem with Output resultfile: %s\n", filereso);
 5750:     fprintf(ficlog,"Problem with Output resultfile: %s\n", filereso);
 5751:     fflush(ficlog);
 5752:     goto end;
 5753:   }
 5754: 
 5755:   /* Reads comments: lines beginning with '#' */
 5756:   numlinepar=0;
 5757:   while((c=getc(ficpar))=='#' && c!= EOF){
 5758:     ungetc(c,ficpar);
 5759:     fgets(line, MAXLINE, ficpar);
 5760:     numlinepar++;
 5761:     fputs(line,stdout);
 5762:     fputs(line,ficparo);
 5763:     fputs(line,ficlog);
 5764:   }
 5765:   ungetc(c,ficpar);
 5766: 
 5767:   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);
 5768:   numlinepar++;
 5769:   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);
 5770:   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);
 5771:   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);
 5772:   fflush(ficlog);
 5773:   while((c=getc(ficpar))=='#' && c!= EOF){
 5774:     ungetc(c,ficpar);
 5775:     fgets(line, MAXLINE, ficpar);
 5776:     numlinepar++;
 5777:     fputs(line, stdout);
 5778:     //puts(line);
 5779:     fputs(line,ficparo);
 5780:     fputs(line,ficlog);
 5781:   }
 5782:   ungetc(c,ficpar);
 5783: 
 5784:    
 5785:   covar=matrix(0,NCOVMAX,1,n);  /**< used in readdata */
 5786:   cptcovn=0; /*Number of covariates, i.e. number of '+' in model statement plus one, indepently of n in Vn*/
 5787:   /* v1+v2+v3+v2*v4+v5*age makes cptcovn = 5
 5788:      v1+v2*age+v2*v3 makes cptcovn = 3
 5789:   */
 5790:   if (strlen(model)>1) 
 5791:     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*/
 5792:   else
 5793:     ncovmodel=2;
 5794:   nvar=ncovmodel-1; /* Suppressing age as a basic covariate */
 5795:   nforce= (nlstate+ndeath-1)*nlstate; /* Number of forces ij from state i to j */
 5796:   npar= nforce*ncovmodel; /* Number of parameters like aij*/
 5797:   if(npar >MAXPARM || nlstate >NLSTATEMAX || ndeath >NDEATHMAX || ncovmodel>NCOVMAX){
 5798:     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);
 5799:     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);
 5800:     fflush(stdout);
 5801:     fclose (ficlog);
 5802:     goto end;
 5803:   }
 5804:   delti3= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
 5805:   delti=delti3[1][1];
 5806:   /*delti=vector(1,npar); *//* Scale of each paramater (output from hesscov)*/
 5807:   if(mle==-1){ /* Print a wizard for help writing covariance matrix */
 5808:     prwizard(ncovmodel, nlstate, ndeath, model, ficparo);
 5809:     printf(" You choose mle=-1, look at file %s for a template of covariance matrix \n",filereso);
 5810:     fprintf(ficlog," You choose mle=-1, look at file %s for a template of covariance matrix \n",filereso);
 5811:     free_ma3x(delti3,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel); 
 5812:     fclose (ficparo);
 5813:     fclose (ficlog);
 5814:     goto end;
 5815:     exit(0);
 5816:   }
 5817:   else if(mle==-3) {
 5818:     prwizard(ncovmodel, nlstate, ndeath, model, ficparo);
 5819:     printf(" You choose mle=-3, look at file %s for a template of covariance matrix \n",filereso);
 5820:     fprintf(ficlog," You choose mle=-3, look at file %s for a template of covariance matrix \n",filereso);
 5821:     param= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
 5822:     matcov=matrix(1,npar,1,npar);
 5823:   }
 5824:   else{
 5825:     /* Read guessed parameters */
 5826:     /* Reads comments: lines beginning with '#' */
 5827:     while((c=getc(ficpar))=='#' && c!= EOF){
 5828:       ungetc(c,ficpar);
 5829:       fgets(line, MAXLINE, ficpar);
 5830:       numlinepar++;
 5831:       fputs(line,stdout);
 5832:       fputs(line,ficparo);
 5833:       fputs(line,ficlog);
 5834:     }
 5835:     ungetc(c,ficpar);
 5836:     
 5837:     param= ma3x(1,nlstate,1,nlstate+ndeath-1,1,ncovmodel);
 5838:     for(i=1; i <=nlstate; i++){
 5839:       j=0;
 5840:       for(jj=1; jj <=nlstate+ndeath; jj++){
 5841: 	if(jj==i) continue;
 5842: 	j++;
 5843: 	fscanf(ficpar,"%1d%1d",&i1,&j1);
 5844: 	if ((i1 != i) && (j1 != j)){
 5845: 	  printf("Error in line parameters number %d, %1d%1d instead of %1d%1d \n \
 5846: It might be a problem of design; if ncovcol and the model are correct\n \
 5847: run imach with mle=-1 to get a correct template of the parameter file.\n",numlinepar, i,j, i1, j1);
 5848: 	  exit(1);
 5849: 	}
 5850: 	fprintf(ficparo,"%1d%1d",i1,j1);
 5851: 	if(mle==1)
 5852: 	  printf("%1d%1d",i,j);
 5853: 	fprintf(ficlog,"%1d%1d",i,j);
 5854: 	for(k=1; k<=ncovmodel;k++){
 5855: 	  fscanf(ficpar," %lf",&param[i][j][k]);
 5856: 	  if(mle==1){
 5857: 	    printf(" %lf",param[i][j][k]);
 5858: 	    fprintf(ficlog," %lf",param[i][j][k]);
 5859: 	  }
 5860: 	  else
 5861: 	    fprintf(ficlog," %lf",param[i][j][k]);
 5862: 	  fprintf(ficparo," %lf",param[i][j][k]);
 5863: 	}
 5864: 	fscanf(ficpar,"\n");
 5865: 	numlinepar++;
 5866: 	if(mle==1)
 5867: 	  printf("\n");
 5868: 	fprintf(ficlog,"\n");
 5869: 	fprintf(ficparo,"\n");
 5870:       }
 5871:     }  
 5872:     fflush(ficlog);
 5873: 
 5874:     /* Reads scales values */
 5875:     p=param[1][1];
 5876:     
 5877:     /* Reads comments: lines beginning with '#' */
 5878:     while((c=getc(ficpar))=='#' && c!= EOF){
 5879:       ungetc(c,ficpar);
 5880:       fgets(line, MAXLINE, ficpar);
 5881:       numlinepar++;
 5882:       fputs(line,stdout);
 5883:       fputs(line,ficparo);
 5884:       fputs(line,ficlog);
 5885:     }
 5886:     ungetc(c,ficpar);
 5887: 
 5888:     for(i=1; i <=nlstate; i++){
 5889:       for(j=1; j <=nlstate+ndeath-1; j++){
 5890: 	fscanf(ficpar,"%1d%1d",&i1,&j1);
 5891: 	if ( (i1-i) * (j1-j) != 0){
 5892: 	  printf("Error in line parameters number %d, %1d%1d instead of %1d%1d \n",numlinepar, i,j, i1, j1);
 5893: 	  exit(1);
 5894: 	}
 5895: 	printf("%1d%1d",i,j);
 5896: 	fprintf(ficparo,"%1d%1d",i1,j1);
 5897: 	fprintf(ficlog,"%1d%1d",i1,j1);
 5898: 	for(k=1; k<=ncovmodel;k++){
 5899: 	  fscanf(ficpar,"%le",&delti3[i][j][k]);
 5900: 	  printf(" %le",delti3[i][j][k]);
 5901: 	  fprintf(ficparo," %le",delti3[i][j][k]);
 5902: 	  fprintf(ficlog," %le",delti3[i][j][k]);
 5903: 	}
 5904: 	fscanf(ficpar,"\n");
 5905: 	numlinepar++;
 5906: 	printf("\n");
 5907: 	fprintf(ficparo,"\n");
 5908: 	fprintf(ficlog,"\n");
 5909:       }
 5910:     }
 5911:     fflush(ficlog);
 5912: 
 5913:     /* Reads covariance matrix */
 5914:     delti=delti3[1][1];
 5915: 
 5916: 
 5917:     /* free_ma3x(delti3,1,nlstate,1,nlstate+ndeath-1,1,ncovmodel); */ /* Hasn't to to freed here otherwise delti is no more allocated */
 5918:   
 5919:     /* Reads comments: lines beginning with '#' */
 5920:     while((c=getc(ficpar))=='#' && c!= EOF){
 5921:       ungetc(c,ficpar);
 5922:       fgets(line, MAXLINE, ficpar);
 5923:       numlinepar++;
 5924:       fputs(line,stdout);
 5925:       fputs(line,ficparo);
 5926:       fputs(line,ficlog);
 5927:     }
 5928:     ungetc(c,ficpar);
 5929:   
 5930:     matcov=matrix(1,npar,1,npar);
 5931:     for(i=1; i <=npar; i++)
 5932:       for(j=1; j <=npar; j++) matcov[i][j]=0.;
 5933:       
 5934:     for(i=1; i <=npar; i++){
 5935:       fscanf(ficpar,"%s",str);
 5936:       if(mle==1)
 5937: 	printf("%s",str);
 5938:       fprintf(ficlog,"%s",str);
 5939:       fprintf(ficparo,"%s",str);
 5940:       for(j=1; j <=i; j++){
 5941: 	fscanf(ficpar," %le",&matcov[i][j]);
 5942: 	if(mle==1){
 5943: 	  printf(" %.5le",matcov[i][j]);
 5944: 	}
 5945: 	fprintf(ficlog," %.5le",matcov[i][j]);
 5946: 	fprintf(ficparo," %.5le",matcov[i][j]);
 5947:       }
 5948:       fscanf(ficpar,"\n");
 5949:       numlinepar++;
 5950:       if(mle==1)
 5951: 	printf("\n");
 5952:       fprintf(ficlog,"\n");
 5953:       fprintf(ficparo,"\n");
 5954:     }
 5955:     for(i=1; i <=npar; i++)
 5956:       for(j=i+1;j<=npar;j++)
 5957: 	matcov[i][j]=matcov[j][i];
 5958:     
 5959:     if(mle==1)
 5960:       printf("\n");
 5961:     fprintf(ficlog,"\n");
 5962:     
 5963:     fflush(ficlog);
 5964:     
 5965:     /*-------- Rewriting parameter file ----------*/
 5966:     strcpy(rfileres,"r");    /* "Rparameterfile */
 5967:     strcat(rfileres,optionfilefiname);    /* Parameter file first name*/
 5968:     strcat(rfileres,".");    /* */
 5969:     strcat(rfileres,optionfilext);    /* Other files have txt extension */
 5970:     if((ficres =fopen(rfileres,"w"))==NULL) {
 5971:       printf("Problem writing new parameter file: %s\n", fileres);goto end;
 5972:       fprintf(ficlog,"Problem writing new parameter file: %s\n", fileres);goto end;
 5973:     }
 5974:     fprintf(ficres,"#%s\n",version);
 5975:   }    /* End of mle != -3 */
 5976: 
 5977: 
 5978:   n= lastobs;
 5979:   num=lvector(1,n);
 5980:   moisnais=vector(1,n);
 5981:   annais=vector(1,n);
 5982:   moisdc=vector(1,n);
 5983:   andc=vector(1,n);
 5984:   agedc=vector(1,n);
 5985:   cod=ivector(1,n);
 5986:   weight=vector(1,n);
 5987:   for(i=1;i<=n;i++) weight[i]=1.0; /* Equal weights, 1 by default */
 5988:   mint=matrix(1,maxwav,1,n);
 5989:   anint=matrix(1,maxwav,1,n);
 5990:   s=imatrix(1,maxwav+1,1,n); /* s[i][j] health state for wave i and individual j */ 
 5991:   tab=ivector(1,NCOVMAX);
 5992:   ncodemax=ivector(1,NCOVMAX); /* Number of code per covariate; if O and 1 only, 2**ncov; V1+V2+V3+V4=>16 */
 5993: 
 5994:   /* Reads data from file datafile */
 5995:   if (readdata(datafile, firstobs, lastobs, &imx)==1)
 5996:     goto end;
 5997: 
 5998:   /* Calculation of the number of parameters from char model */
 5999:     /*    modelsav=V2+V1+V4+age*V3 strb=age*V3 stra=V2+V1+V4 
 6000: 	k=4 (age*V3) Tvar[k=4]= 3 (from V3) Tag[cptcovage=1]=4
 6001: 	k=3 V4 Tvar[k=3]= 4 (from V4)
 6002: 	k=2 V1 Tvar[k=2]= 1 (from V1)
 6003: 	k=1 Tvar[1]=2 (from V2)
 6004:     */
 6005:   Tvar=ivector(1,NCOVMAX); /* Was 15 changed to NCOVMAX. */
 6006:   /*  V2+V1+V4+age*V3 is a model with 4 covariates (3 plus signs). 
 6007:       For each model-covariate stores the data-covariate id. Tvar[1]=2, Tvar[2]=1, Tvar[3]=4, 
 6008:       Tvar[4=age*V3] is 3 and 'age' is recorded in Tage.
 6009:   */
 6010:   /* For model-covariate k tells which data-covariate to use but
 6011:     because this model-covariate is a construction we invent a new column
 6012:     ncovcol + k1
 6013:     If already ncovcol=4 and model=V2+V1+V1*V4+age*V3
 6014:     Tvar[3=V1*V4]=4+1 etc */
 6015:   Tprod=ivector(1,NCOVMAX); /* Gives the position of a product */
 6016:   /* Tprod[k1=1]=3(=V1*V4) for V2+V1+V1*V4+age*V3
 6017:      if  V2+V1+V1*V4+age*V3+V3*V2   TProd[k1=2]=5 (V3*V2)
 6018:   */
 6019:   Tvaraff=ivector(1,NCOVMAX); /* Unclear */
 6020:   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
 6021: 			    * For V3*V2 (in V2+V1+V1*V4+age*V3+V3*V2), V3*V2 position is 2nd. 
 6022: 			    * Tvard[k1=2][1]=3 (V3) Tvard[k1=2][2]=2(V2) */
 6023:   Tage=ivector(1,NCOVMAX); /* Gives the covariate id of covariates associated with age: V2 + V1 + age*V4 + V3*age
 6024: 			 4 covariates (3 plus signs)
 6025: 			 Tage[1=V3*age]= 4; Tage[2=age*V4] = 3
 6026: 		      */  
 6027: 
 6028:   if(decodemodel(model, lastobs) == 1)
 6029:     goto end;
 6030: 
 6031:   if((double)(lastobs-imx)/(double)imx > 1.10){
 6032:     nbwarn++;
 6033:     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); 
 6034:     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); 
 6035:   }
 6036:     /*  if(mle==1){*/
 6037:   if (weightopt != 1) { /* Maximisation without weights. We can have weights different from 1 but want no weight*/
 6038:     for(i=1;i<=imx;i++) weight[i]=1.0; /* changed to imx */
 6039:   }
 6040: 
 6041:     /*-calculation of age at interview from date of interview and age at death -*/
 6042:   agev=matrix(1,maxwav,1,imx);
 6043: 
 6044:   if(calandcheckages(imx, maxwav, &agemin, &agemax, &nberr, &nbwarn) == 1)
 6045:     goto end;
 6046: 
 6047: 
 6048:   agegomp=(int)agemin;
 6049:   free_vector(moisnais,1,n);
 6050:   free_vector(annais,1,n);
 6051:   /* free_matrix(mint,1,maxwav,1,n);
 6052:      free_matrix(anint,1,maxwav,1,n);*/
 6053:   free_vector(moisdc,1,n);
 6054:   free_vector(andc,1,n);
 6055:   /* */
 6056:   
 6057:   wav=ivector(1,imx);
 6058:   dh=imatrix(1,lastpass-firstpass+1,1,imx);
 6059:   bh=imatrix(1,lastpass-firstpass+1,1,imx);
 6060:   mw=imatrix(1,lastpass-firstpass+1,1,imx);
 6061:    
 6062:   /* Concatenates waves */
 6063:   concatwav(wav, dh, bh, mw, s, agedc, agev,  firstpass, lastpass, imx, nlstate, stepm);
 6064:   /* */
 6065:  
 6066:   /* Routine tricode is to calculate cptcoveff (real number of unique covariates) and to associate covariable number and modality */
 6067: 
 6068:   nbcode=imatrix(0,NCOVMAX,0,NCOVMAX); 
 6069:   ncodemax[1]=1;
 6070:   Ndum =ivector(-1,NCOVMAX);  
 6071:   if (ncovmodel > 2)
 6072:     tricode(Tvar,nbcode,imx, Ndum); /**< Fills nbcode[Tvar[j]][l]; */
 6073: 
 6074:   codtab=imatrix(1,100,1,10); /* codtab[h,k]=( (h-1) - mod(k-1,2**(k-1) )/2**(k-1) */
 6075:   /*printf(" codtab[1,1],codtab[100,10]=%d,%d\n", codtab[1][1],codtab[100][10]);*/
 6076:   h=0;
 6077: 
 6078: 
 6079:   /*if (cptcovn > 0) */
 6080:       
 6081:  
 6082:   m=pow(2,cptcoveff);
 6083:  
 6084:   for(k=1;k<=cptcoveff; k++){ /* scans any effective covariate */
 6085:     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 */ 
 6086:       for(j=1; j <= ncodemax[k]; j++){ /* For each modality of this covariate ncodemax=2*/
 6087: 	for(cpt=1; cpt <=pow(2,k-1); cpt++){  /* cpt=1 to 8/2**(3+1-1 or 3+1-3) =1 or 4 */ 
 6088: 	  h++;
 6089: 	  if (h>m) 
 6090: 	    h=1;
 6091: 	  /**< codtab(h,k)  k   = codtab[h,k]=( (h-1) - mod(k-1,2**(k-1) )/2**(k-1) + 1
 6092: 	   *     h     1     2     3     4
 6093: 	   *______________________________  
 6094: 	   *     1 i=1 1 i=1 1 i=1 1 i=1 1
 6095: 	   *     2     2     1     1     1
 6096: 	   *     3 i=2 1     2     1     1
 6097: 	   *     4     2     2     1     1
 6098: 	   *     5 i=3 1 i=2 1     2     1
 6099: 	   *     6     2     1     2     1
 6100: 	   *     7 i=4 1     2     2     1
 6101: 	   *     8     2     2     2     1
 6102: 	   *     9 i=5 1 i=3 1 i=2 1     1
 6103: 	   *    10     2     1     1     1
 6104: 	   *    11 i=6 1     2     1     1
 6105: 	   *    12     2     2     1     1
 6106: 	   *    13 i=7 1 i=4 1     2     1    
 6107: 	   *    14     2     1     2     1
 6108: 	   *    15 i=8 1     2     2     1
 6109: 	   *    16     2     2     2     1
 6110: 	   */
 6111: 	  codtab[h][k]=j;
 6112: 	  /*codtab[h][Tvar[k]]=j;*/
 6113: 	  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]]);
 6114: 	} 
 6115:       }
 6116:     }
 6117:   } 
 6118:   /* printf("codtab[1][2]=%d codtab[2][2]=%d",codtab[1][2],codtab[2][2]); 
 6119:      codtab[1][2]=1;codtab[2][2]=2; */
 6120:   /* for(i=1; i <=m ;i++){ 
 6121:      for(k=1; k <=cptcovn; k++){
 6122:        printf("i=%d k=%d %d %d ",i,k,codtab[i][k], cptcoveff);
 6123:      }
 6124:      printf("\n");
 6125:      }
 6126:      scanf("%d",i);*/
 6127: 
 6128:  free_ivector(Ndum,-1,NCOVMAX);
 6129: 
 6130: 
 6131:     
 6132:   /*------------ gnuplot -------------*/
 6133:   strcpy(optionfilegnuplot,optionfilefiname);
 6134:   if(mle==-3)
 6135:     strcat(optionfilegnuplot,"-mort");
 6136:   strcat(optionfilegnuplot,".gp");
 6137: 
 6138:   if((ficgp=fopen(optionfilegnuplot,"w"))==NULL) {
 6139:     printf("Problem with file %s",optionfilegnuplot);
 6140:   }
 6141:   else{
 6142:     fprintf(ficgp,"\n# %s\n", version); 
 6143:     fprintf(ficgp,"# %s\n", optionfilegnuplot); 
 6144:     //fprintf(ficgp,"set missing 'NaNq'\n");
 6145:     fprintf(ficgp,"set datafile missing 'NaNq'\n");
 6146:   }
 6147:   /*  fclose(ficgp);*/
 6148:   /*--------- index.htm --------*/
 6149: 
 6150:   strcpy(optionfilehtm,optionfilefiname); /* Main html file */
 6151:   if(mle==-3)
 6152:     strcat(optionfilehtm,"-mort");
 6153:   strcat(optionfilehtm,".htm");
 6154:   if((fichtm=fopen(optionfilehtm,"w"))==NULL)    {
 6155:     printf("Problem with %s \n",optionfilehtm);
 6156:     exit(0);
 6157:   }
 6158: 
 6159:   strcpy(optionfilehtmcov,optionfilefiname); /* Only for matrix of covariance */
 6160:   strcat(optionfilehtmcov,"-cov.htm");
 6161:   if((fichtmcov=fopen(optionfilehtmcov,"w"))==NULL)    {
 6162:     printf("Problem with %s \n",optionfilehtmcov), exit(0);
 6163:   }
 6164:   else{
 6165:   fprintf(fichtmcov,"<html><head>\n<title>IMaCh Cov %s</title></head>\n <body><font size=\"2\">%s <br> %s</font> \
 6166: <hr size=\"2\" color=\"#EC5E5E\"> \n\
 6167: Title=%s <br>Datafile=%s Firstpass=%d Lastpass=%d Stepm=%d Weight=%d Model=%s<br>\n",\
 6168: 	  optionfilehtmcov,version,fullversion,title,datafile,firstpass,lastpass,stepm, weightopt, model);
 6169:   }
 6170: 
 6171:   fprintf(fichtm,"<html><head>\n<title>IMaCh %s</title></head>\n <body><font size=\"2\">%s <br> %s</font> \
 6172: <hr size=\"2\" color=\"#EC5E5E\"> \n\
 6173: Title=%s <br>Datafile=%s Firstpass=%d Lastpass=%d Stepm=%d Weight=%d Model=%s<br>\n\
 6174: \n\
 6175: <hr  size=\"2\" color=\"#EC5E5E\">\
 6176:  <ul><li><h4>Parameter files</h4>\n\
 6177:  - Parameter file: <a href=\"%s.%s\">%s.%s</a><br>\n\
 6178:  - Copy of the parameter file: <a href=\"o%s\">o%s</a><br>\n\
 6179:  - Log file of the run: <a href=\"%s\">%s</a><br>\n\
 6180:  - Gnuplot file name: <a href=\"%s\">%s</a><br>\n\
 6181:  - Date and time at start: %s</ul>\n",\
 6182: 	  optionfilehtm,version,fullversion,title,datafile,firstpass,lastpass,stepm, weightopt, model,\
 6183: 	  optionfilefiname,optionfilext,optionfilefiname,optionfilext,\
 6184: 	  fileres,fileres,\
 6185: 	  filelog,filelog,optionfilegnuplot,optionfilegnuplot,strstart);
 6186:   fflush(fichtm);
 6187: 
 6188:   strcpy(pathr,path);
 6189:   strcat(pathr,optionfilefiname);
 6190:   chdir(optionfilefiname); /* Move to directory named optionfile */
 6191:   
 6192:   /* Calculates basic frequencies. Computes observed prevalence at single age
 6193:      and prints on file fileres'p'. */
 6194:   freqsummary(fileres, agemin, agemax, s, agev, nlstate, imx,Tvaraff,nbcode, ncodemax,mint,anint,strstart);
 6195: 
 6196:   fprintf(fichtm,"\n");
 6197:   fprintf(fichtm,"<br>Total number of observations=%d <br>\n\
 6198: Youngest age at first (selected) pass %.2f, oldest age %.2f<br>\n\
 6199: Interval (in months) between two waves: Min=%d Max=%d Mean=%.2lf<br>\n",\
 6200: 	  imx,agemin,agemax,jmin,jmax,jmean);
 6201:   pmmij= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6202:     oldms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6203:     newms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6204:     savms= matrix(1,nlstate+ndeath,1,nlstate+ndeath); /* creation */
 6205:     oldm=oldms; newm=newms; savm=savms; /* Keeps fixed addresses to free */
 6206:     
 6207:    
 6208:   /* For Powell, parameters are in a vector p[] starting at p[1]
 6209:      so we point p on param[1][1] so that p[1] maps on param[1][1][1] */
 6210:   p=param[1][1]; /* *(*(*(param +1)+1)+0) */
 6211: 
 6212:   globpr=0; /* To get the number ipmx of contributions and the sum of weights*/
 6213: 
 6214:   if (mle==-3){
 6215:     ximort=matrix(1,NDIM,1,NDIM); 
 6216: /*     ximort=gsl_matrix_alloc(1,NDIM,1,NDIM); */
 6217:     cens=ivector(1,n);
 6218:     ageexmed=vector(1,n);
 6219:     agecens=vector(1,n);
 6220:     dcwave=ivector(1,n);
 6221:  
 6222:     for (i=1; i<=imx; i++){
 6223:       dcwave[i]=-1;
 6224:       for (m=firstpass; m<=lastpass; m++)
 6225: 	if (s[m][i]>nlstate) {
 6226: 	  dcwave[i]=m;
 6227: 	  /*	printf("i=%d j=%d s=%d dcwave=%d\n",i,j, s[j][i],dcwave[i]);*/
 6228: 	  break;
 6229: 	}
 6230:     }
 6231: 
 6232:     for (i=1; i<=imx; i++) {
 6233:       if (wav[i]>0){
 6234: 	ageexmed[i]=agev[mw[1][i]][i];
 6235: 	j=wav[i];
 6236: 	agecens[i]=1.; 
 6237: 
 6238: 	if (ageexmed[i]> 1 && wav[i] > 0){
 6239: 	  agecens[i]=agev[mw[j][i]][i];
 6240: 	  cens[i]= 1;
 6241: 	}else if (ageexmed[i]< 1) 
 6242: 	  cens[i]= -1;
 6243: 	if (agedc[i]< AGESUP && agedc[i]>1 && dcwave[i]>firstpass && dcwave[i]<=lastpass)
 6244: 	  cens[i]=0 ;
 6245:       }
 6246:       else cens[i]=-1;
 6247:     }
 6248:     
 6249:     for (i=1;i<=NDIM;i++) {
 6250:       for (j=1;j<=NDIM;j++)
 6251: 	ximort[i][j]=(i == j ? 1.0 : 0.0);
 6252:     }
 6253:     
 6254:     /*p[1]=0.0268; p[NDIM]=0.083;*/
 6255:     /*printf("%lf %lf", p[1], p[2]);*/
 6256:     
 6257:     
 6258: #ifdef GSL
 6259:     printf("GSL optimization\n");  fprintf(ficlog,"Powell\n");
 6260: #else
 6261:     printf("Powell\n");  fprintf(ficlog,"Powell\n");
 6262: #endif
 6263:     strcpy(filerespow,"pow-mort"); 
 6264:     strcat(filerespow,fileres);
 6265:     if((ficrespow=fopen(filerespow,"w"))==NULL) {
 6266:       printf("Problem with resultfile: %s\n", filerespow);
 6267:       fprintf(ficlog,"Problem with resultfile: %s\n", filerespow);
 6268:     }
 6269: #ifdef GSL
 6270:     fprintf(ficrespow,"# GSL optimization\n# iter -2*LL");
 6271: #else
 6272:     fprintf(ficrespow,"# Powell\n# iter -2*LL");
 6273: #endif
 6274:     /*  for (i=1;i<=nlstate;i++)
 6275: 	for(j=1;j<=nlstate+ndeath;j++)
 6276: 	if(j!=i)fprintf(ficrespow," p%1d%1d",i,j);
 6277:     */
 6278:     fprintf(ficrespow,"\n");
 6279: #ifdef GSL
 6280:     /* gsl starts here */ 
 6281:     T = gsl_multimin_fminimizer_nmsimplex;
 6282:     gsl_multimin_fminimizer *sfm = NULL;
 6283:     gsl_vector *ss, *x;
 6284:     gsl_multimin_function minex_func;
 6285: 
 6286:     /* Initial vertex size vector */
 6287:     ss = gsl_vector_alloc (NDIM);
 6288:     
 6289:     if (ss == NULL){
 6290:       GSL_ERROR_VAL ("failed to allocate space for ss", GSL_ENOMEM, 0);
 6291:     }
 6292:     /* Set all step sizes to 1 */
 6293:     gsl_vector_set_all (ss, 0.001);
 6294: 
 6295:     /* Starting point */
 6296:     
 6297:     x = gsl_vector_alloc (NDIM);
 6298:     
 6299:     if (x == NULL){
 6300:       gsl_vector_free(ss);
 6301:       GSL_ERROR_VAL ("failed to allocate space for x", GSL_ENOMEM, 0);
 6302:     }
 6303:   
 6304:     /* Initialize method and iterate */
 6305:     /*     p[1]=0.0268; p[NDIM]=0.083; */
 6306: /*     gsl_vector_set(x, 0, 0.0268); */
 6307: /*     gsl_vector_set(x, 1, 0.083); */
 6308:     gsl_vector_set(x, 0, p[1]);
 6309:     gsl_vector_set(x, 1, p[2]);
 6310: 
 6311:     minex_func.f = &gompertz_f;
 6312:     minex_func.n = NDIM;
 6313:     minex_func.params = (void *)&p; /* ??? */
 6314:     
 6315:     sfm = gsl_multimin_fminimizer_alloc (T, NDIM);
 6316:     gsl_multimin_fminimizer_set (sfm, &minex_func, x, ss);
 6317:     
 6318:     printf("Iterations beginning .....\n\n");
 6319:     printf("Iter. #    Intercept       Slope     -Log Likelihood     Simplex size\n");
 6320: 
 6321:     iteri=0;
 6322:     while (rval == GSL_CONTINUE){
 6323:       iteri++;
 6324:       status = gsl_multimin_fminimizer_iterate(sfm);
 6325:       
 6326:       if (status) printf("error: %s\n", gsl_strerror (status));
 6327:       fflush(0);
 6328:       
 6329:       if (status) 
 6330:         break;
 6331:       
 6332:       rval = gsl_multimin_test_size (gsl_multimin_fminimizer_size (sfm), 1e-6);
 6333:       ssval = gsl_multimin_fminimizer_size (sfm);
 6334:       
 6335:       if (rval == GSL_SUCCESS)
 6336:         printf ("converged to a local maximum at\n");
 6337:       
 6338:       printf("%5d ", iteri);
 6339:       for (it = 0; it < NDIM; it++){
 6340: 	printf ("%10.5f ", gsl_vector_get (sfm->x, it));
 6341:       }
 6342:       printf("f() = %-10.5f ssize = %.7f\n", sfm->fval, ssval);
 6343:     }
 6344:     
 6345:     printf("\n\n Please note: Program should be run many times with varying starting points to detemine global maximum\n\n");
 6346:     
 6347:     gsl_vector_free(x); /* initial values */
 6348:     gsl_vector_free(ss); /* inital step size */
 6349:     for (it=0; it<NDIM; it++){
 6350:       p[it+1]=gsl_vector_get(sfm->x,it);
 6351:       fprintf(ficrespow," %.12lf", p[it]);
 6352:     }
 6353:     gsl_multimin_fminimizer_free (sfm); /* p *(sfm.x.data) et p *(sfm.x.data+1)  */
 6354: #endif
 6355: #ifdef POWELL
 6356:      powell(p,ximort,NDIM,ftol,&iter,&fret,gompertz);
 6357: #endif  
 6358:     fclose(ficrespow);
 6359:     
 6360:     hesscov(matcov, p, NDIM, delti, 1e-4, gompertz); 
 6361: 
 6362:     for(i=1; i <=NDIM; i++)
 6363:       for(j=i+1;j<=NDIM;j++)
 6364: 	matcov[i][j]=matcov[j][i];
 6365:     
 6366:     printf("\nCovariance matrix\n ");
 6367:     for(i=1; i <=NDIM; i++) {
 6368:       for(j=1;j<=NDIM;j++){ 
 6369: 	printf("%f ",matcov[i][j]);
 6370:       }
 6371:       printf("\n ");
 6372:     }
 6373:     
 6374:     printf("iter=%d MLE=%f Eq=%lf*exp(%lf*(age-%d))\n",iter,-gompertz(p),p[1],p[2],agegomp);
 6375:     for (i=1;i<=NDIM;i++) 
 6376:       printf("%f [%f ; %f]\n",p[i],p[i]-2*sqrt(matcov[i][i]),p[i]+2*sqrt(matcov[i][i]));
 6377: 
 6378:     lsurv=vector(1,AGESUP);
 6379:     lpop=vector(1,AGESUP);
 6380:     tpop=vector(1,AGESUP);
 6381:     lsurv[agegomp]=100000;
 6382:     
 6383:     for (k=agegomp;k<=AGESUP;k++) {
 6384:       agemortsup=k;
 6385:       if (p[1]*exp(p[2]*(k-agegomp))>1) break;
 6386:     }
 6387:     
 6388:     for (k=agegomp;k<agemortsup;k++)
 6389:       lsurv[k+1]=lsurv[k]-lsurv[k]*(p[1]*exp(p[2]*(k-agegomp)));
 6390:     
 6391:     for (k=agegomp;k<agemortsup;k++){
 6392:       lpop[k]=(lsurv[k]+lsurv[k+1])/2.;
 6393:       sumlpop=sumlpop+lpop[k];
 6394:     }
 6395:     
 6396:     tpop[agegomp]=sumlpop;
 6397:     for (k=agegomp;k<(agemortsup-3);k++){
 6398:       /*  tpop[k+1]=2;*/
 6399:       tpop[k+1]=tpop[k]-lpop[k];
 6400:     }
 6401:     
 6402:     
 6403:     printf("\nAge   lx     qx    dx    Lx     Tx     e(x)\n");
 6404:     for (k=agegomp;k<(agemortsup-2);k++) 
 6405:       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]);
 6406:     
 6407:     
 6408:     replace_back_to_slash(pathc,pathcd); /* Even gnuplot wants a / */
 6409:     printinggnuplotmort(fileres, optionfilefiname,ageminpar,agemaxpar,fage, pathc,p);
 6410:     
 6411:     printinghtmlmort(fileres,title,datafile, firstpass, lastpass, \
 6412: 		     stepm, weightopt,\
 6413: 		     model,imx,p,matcov,agemortsup);
 6414:     
 6415:     free_vector(lsurv,1,AGESUP);
 6416:     free_vector(lpop,1,AGESUP);
 6417:     free_vector(tpop,1,AGESUP);
 6418: #ifdef GSL
 6419:     free_ivector(cens,1,n);
 6420:     free_vector(agecens,1,n);
 6421:     free_ivector(dcwave,1,n);
 6422:     free_matrix(ximort,1,NDIM,1,NDIM);
 6423: #endif
 6424:   } /* Endof if mle==-3 */
 6425:   
 6426:   else{ /* For mle >=1 */
 6427:     globpr=0;/* debug */
 6428:     likelione(ficres, p, npar, nlstate, &globpr, &ipmx, &sw, &fretone, funcone); /* Prints the contributions to the likelihood */
 6429:     printf("First Likeli=%12.6f ipmx=%ld sw=%12.6f",fretone,ipmx,sw);
 6430:     for (k=1; k<=npar;k++)
 6431:       printf(" %d %8.5f",k,p[k]);
 6432:     printf("\n");
 6433:     globpr=1; /* to print the contributions */
 6434:     likelione(ficres, p, npar, nlstate, &globpr, &ipmx, &sw, &fretone, funcone); /* Prints the contributions to the likelihood */
 6435:     printf("Second Likeli=%12.6f ipmx=%ld sw=%12.6f",fretone,ipmx,sw);
 6436:     for (k=1; k<=npar;k++)
 6437:       printf(" %d %8.5f",k,p[k]);
 6438:     printf("\n");
 6439:     if(mle>=1){ /* Could be 1 or 2 */
 6440:       mlikeli(ficres,p, npar, ncovmodel, nlstate, ftol, func);
 6441:     }
 6442:     
 6443:     /*--------- results files --------------*/
 6444:     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);
 6445:     
 6446:     
 6447:     fprintf(ficres,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 6448:     printf("# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 6449:     fprintf(ficlog,"# Parameters nlstate*nlstate*ncov a12*1 + b12 * age + ...\n");
 6450:     for(i=1,jk=1; i <=nlstate; i++){
 6451:       for(k=1; k <=(nlstate+ndeath); k++){
 6452: 	if (k != i) {
 6453: 	  printf("%d%d ",i,k);
 6454: 	  fprintf(ficlog,"%d%d ",i,k);
 6455: 	  fprintf(ficres,"%1d%1d ",i,k);
 6456: 	  for(j=1; j <=ncovmodel; j++){
 6457: 	    printf("%lf ",p[jk]);
 6458: 	    fprintf(ficlog,"%lf ",p[jk]);
 6459: 	    fprintf(ficres,"%lf ",p[jk]);
 6460: 	    jk++; 
 6461: 	  }
 6462: 	  printf("\n");
 6463: 	  fprintf(ficlog,"\n");
 6464: 	  fprintf(ficres,"\n");
 6465: 	}
 6466:       }
 6467:     }
 6468:     if(mle!=0){
 6469:       /* Computing hessian and covariance matrix */
 6470:       ftolhess=ftol; /* Usually correct */
 6471:       hesscov(matcov, p, npar, delti, ftolhess, func);
 6472:     }
 6473:     fprintf(ficres,"# Scales (for hessian or gradient estimation)\n");
 6474:     printf("# Scales (for hessian or gradient estimation)\n");
 6475:     fprintf(ficlog,"# Scales (for hessian or gradient estimation)\n");
 6476:     for(i=1,jk=1; i <=nlstate; i++){
 6477:       for(j=1; j <=nlstate+ndeath; j++){
 6478: 	if (j!=i) {
 6479: 	  fprintf(ficres,"%1d%1d",i,j);
 6480: 	  printf("%1d%1d",i,j);
 6481: 	  fprintf(ficlog,"%1d%1d",i,j);
 6482: 	  for(k=1; k<=ncovmodel;k++){
 6483: 	    printf(" %.5e",delti[jk]);
 6484: 	    fprintf(ficlog," %.5e",delti[jk]);
 6485: 	    fprintf(ficres," %.5e",delti[jk]);
 6486: 	    jk++;
 6487: 	  }
 6488: 	  printf("\n");
 6489: 	  fprintf(ficlog,"\n");
 6490: 	  fprintf(ficres,"\n");
 6491: 	}
 6492:       }
 6493:     }
 6494:     
 6495:     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");
 6496:     if(mle>=1)
 6497:       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");
 6498:     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");
 6499:     /* # 121 Var(a12)\n\ */
 6500:     /* # 122 Cov(b12,a12) Var(b12)\n\ */
 6501:     /* # 131 Cov(a13,a12) Cov(a13,b12, Var(a13)\n\ */
 6502:     /* # 132 Cov(b13,a12) Cov(b13,b12, Cov(b13,a13) Var(b13)\n\ */
 6503:     /* # 212 Cov(a21,a12) Cov(a21,b12, Cov(a21,a13) Cov(a21,b13) Var(a21)\n\ */
 6504:     /* # 212 Cov(b21,a12) Cov(b21,b12, Cov(b21,a13) Cov(b21,b13) Cov(b21,a21) Var(b21)\n\ */
 6505:     /* # 232 Cov(a23,a12) Cov(a23,b12, Cov(a23,a13) Cov(a23,b13) Cov(a23,a21) Cov(a23,b21) Var(a23)\n\ */
 6506:     /* # 232 Cov(b23,a12) Cov(b23,b12) ... Var (b23)\n" */
 6507:     
 6508:     
 6509:     /* Just to have a covariance matrix which will be more understandable
 6510:        even is we still don't want to manage dictionary of variables
 6511:     */
 6512:     for(itimes=1;itimes<=2;itimes++){
 6513:       jj=0;
 6514:       for(i=1; i <=nlstate; i++){
 6515: 	for(j=1; j <=nlstate+ndeath; j++){
 6516: 	  if(j==i) continue;
 6517: 	  for(k=1; k<=ncovmodel;k++){
 6518: 	    jj++;
 6519: 	    ca[0]= k+'a'-1;ca[1]='\0';
 6520: 	    if(itimes==1){
 6521: 	      if(mle>=1)
 6522: 		printf("#%1d%1d%d",i,j,k);
 6523: 	      fprintf(ficlog,"#%1d%1d%d",i,j,k);
 6524: 	      fprintf(ficres,"#%1d%1d%d",i,j,k);
 6525: 	    }else{
 6526: 	      if(mle>=1)
 6527: 		printf("%1d%1d%d",i,j,k);
 6528: 	      fprintf(ficlog,"%1d%1d%d",i,j,k);
 6529: 	      fprintf(ficres,"%1d%1d%d",i,j,k);
 6530: 	    }
 6531: 	    ll=0;
 6532: 	    for(li=1;li <=nlstate; li++){
 6533: 	      for(lj=1;lj <=nlstate+ndeath; lj++){
 6534: 		if(lj==li) continue;
 6535: 		for(lk=1;lk<=ncovmodel;lk++){
 6536: 		  ll++;
 6537: 		  if(ll<=jj){
 6538: 		    cb[0]= lk +'a'-1;cb[1]='\0';
 6539: 		    if(ll<jj){
 6540: 		      if(itimes==1){
 6541: 			if(mle>=1)
 6542: 			  printf(" Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 6543: 			fprintf(ficlog," Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 6544: 			fprintf(ficres," Cov(%s%1d%1d,%s%1d%1d)",ca,i,j,cb, li,lj);
 6545: 		      }else{
 6546: 			if(mle>=1)
 6547: 			  printf(" %.5e",matcov[jj][ll]); 
 6548: 			fprintf(ficlog," %.5e",matcov[jj][ll]); 
 6549: 			fprintf(ficres," %.5e",matcov[jj][ll]); 
 6550: 		      }
 6551: 		    }else{
 6552: 		      if(itimes==1){
 6553: 			if(mle>=1)
 6554: 			  printf(" Var(%s%1d%1d)",ca,i,j);
 6555: 			fprintf(ficlog," Var(%s%1d%1d)",ca,i,j);
 6556: 			fprintf(ficres," Var(%s%1d%1d)",ca,i,j);
 6557: 		      }else{
 6558: 			if(mle>=1)
 6559: 			  printf(" %.5e",matcov[jj][ll]); 
 6560: 			fprintf(ficlog," %.5e",matcov[jj][ll]); 
 6561: 			fprintf(ficres," %.5e",matcov[jj][ll]); 
 6562: 		      }
 6563: 		    }
 6564: 		  }
 6565: 		} /* end lk */
 6566: 	      } /* end lj */
 6567: 	    } /* end li */
 6568: 	    if(mle>=1)
 6569: 	      printf("\n");
 6570: 	    fprintf(ficlog,"\n");
 6571: 	    fprintf(ficres,"\n");
 6572: 	    numlinepar++;
 6573: 	  } /* end k*/
 6574: 	} /*end j */
 6575:       } /* end i */
 6576:     } /* end itimes */
 6577:     
 6578:     fflush(ficlog);
 6579:     fflush(ficres);
 6580:     
 6581:     while((c=getc(ficpar))=='#' && c!= EOF){
 6582:       ungetc(c,ficpar);
 6583:       fgets(line, MAXLINE, ficpar);
 6584:       fputs(line,stdout);
 6585:       fputs(line,ficparo);
 6586:     }
 6587:     ungetc(c,ficpar);
 6588:     
 6589:     estepm=0;
 6590:     fscanf(ficpar,"agemin=%lf agemax=%lf bage=%lf fage=%lf estepm=%d\n",&ageminpar,&agemaxpar, &bage, &fage, &estepm);
 6591:     if (estepm==0 || estepm < stepm) estepm=stepm;
 6592:     if (fage <= 2) {
 6593:       bage = ageminpar;
 6594:       fage = agemaxpar;
 6595:     }
 6596:     
 6597:     fprintf(ficres,"# agemin agemax for life expectancy, bage fage (if mle==0 ie no data nor Max likelihood).\n");
 6598:     fprintf(ficres,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f estepm=%d\n",ageminpar,agemaxpar,bage,fage, estepm);
 6599:     fprintf(ficparo,"agemin=%.0f agemax=%.0f bage=%.0f fage=%.0f estepm=%d\n",ageminpar,agemaxpar,bage,fage, estepm);
 6600:     
 6601:     while((c=getc(ficpar))=='#' && c!= EOF){
 6602:       ungetc(c,ficpar);
 6603:       fgets(line, MAXLINE, ficpar);
 6604:       fputs(line,stdout);
 6605:       fputs(line,ficparo);
 6606:     }
 6607:     ungetc(c,ficpar);
 6608:     
 6609:     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);
 6610:     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);
 6611:     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);
 6612:     printf("begin-prev-date=%.lf/%.lf/%.lf end-prev-date=%.lf/%.lf/%.lf mov_average=%d\n",jprev1, mprev1,anprev1,jprev2, mprev2,anprev2,mobilav);
 6613:     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);
 6614:     
 6615:     while((c=getc(ficpar))=='#' && c!= EOF){
 6616:       ungetc(c,ficpar);
 6617:       fgets(line, MAXLINE, ficpar);
 6618:       fputs(line,stdout);
 6619:       fputs(line,ficparo);
 6620:     }
 6621:     ungetc(c,ficpar);
 6622:     
 6623:     
 6624:     dateprev1=anprev1+(mprev1-1)/12.+(jprev1-1)/365.;
 6625:     dateprev2=anprev2+(mprev2-1)/12.+(jprev2-1)/365.;
 6626:     
 6627:     fscanf(ficpar,"pop_based=%d\n",&popbased);
 6628:     fprintf(ficparo,"pop_based=%d\n",popbased);   
 6629:     fprintf(ficres,"pop_based=%d\n",popbased);   
 6630:     
 6631:     while((c=getc(ficpar))=='#' && c!= EOF){
 6632:       ungetc(c,ficpar);
 6633:       fgets(line, MAXLINE, ficpar);
 6634:       fputs(line,stdout);
 6635:       fputs(line,ficparo);
 6636:     }
 6637:     ungetc(c,ficpar);
 6638:     
 6639:     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);
 6640:     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);
 6641:     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);
 6642:     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);
 6643:     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);
 6644:     /* day and month of proj2 are not used but only year anproj2.*/
 6645:     
 6646:     
 6647:     
 6648:      /* freqsummary(fileres, agemin, agemax, s, agev, nlstate, imx,Tvaraff,nbcode, ncodemax,mint,anint); */
 6649:     /* ,dateprev1,dateprev2,jprev1, mprev1,anprev1,jprev2, mprev2,anprev2); */
 6650:     
 6651:     replace_back_to_slash(pathc,pathcd); /* Even gnuplot wants a / */
 6652:     printinggnuplot(fileres, optionfilefiname,ageminpar,agemaxpar,fage, pathc,p);
 6653:     
 6654:     printinghtml(fileres,title,datafile, firstpass, lastpass, stepm, weightopt,\
 6655: 		 model,imx,jmin,jmax,jmean,rfileres,popforecast,estepm,\
 6656: 		 jprev1,mprev1,anprev1,jprev2,mprev2,anprev2);
 6657:       
 6658:    /*------------ free_vector  -------------*/
 6659:    /*  chdir(path); */
 6660:  
 6661:     free_ivector(wav,1,imx);
 6662:     free_imatrix(dh,1,lastpass-firstpass+1,1,imx);
 6663:     free_imatrix(bh,1,lastpass-firstpass+1,1,imx);
 6664:     free_imatrix(mw,1,lastpass-firstpass+1,1,imx);   
 6665:     free_lvector(num,1,n);
 6666:     free_vector(agedc,1,n);
 6667:     /*free_matrix(covar,0,NCOVMAX,1,n);*/
 6668:     /*free_matrix(covar,1,NCOVMAX,1,n);*/
 6669:     fclose(ficparo);
 6670:     fclose(ficres);
 6671: 
 6672: 
 6673:     /*--------------- Prevalence limit  (period or stable prevalence) --------------*/
 6674: #include "prevlim.h"  /* Use ficrespl, ficlog */
 6675:     fclose(ficrespl);
 6676: 
 6677: #ifdef FREEEXIT2
 6678: #include "freeexit2.h"
 6679: #endif
 6680: 
 6681:     /*------------- h Pij x at various ages ------------*/
 6682: #include "hpijx.h"
 6683:     fclose(ficrespij);
 6684: 
 6685:   /*-------------- Variance of one-step probabilities---*/
 6686:     k=1;
 6687:     varprob(optionfilefiname, matcov, p, delti, nlstate, bage, fage,k,Tvar,nbcode, ncodemax,strstart);
 6688: 
 6689: 
 6690:     probs= ma3x(1,AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6691:     for(i=1;i<=AGESUP;i++)
 6692:       for(j=1;j<=NCOVMAX;j++)
 6693: 	for(k=1;k<=NCOVMAX;k++)
 6694: 	  probs[i][j][k]=0.;
 6695: 
 6696:     /*---------- Forecasting ------------------*/
 6697:     /*if((stepm == 1) && (strcmp(model,".")==0)){*/
 6698:     if(prevfcast==1){
 6699:       /*    if(stepm ==1){*/
 6700:       prevforecast(fileres, anproj1, mproj1, jproj1, agemin, agemax, dateprev1, dateprev2, mobilavproj, bage, fage, firstpass, lastpass, anproj2, p, cptcoveff);
 6701:       /* (popforecast==1) populforecast(fileres, anpyram,mpyram,jpyram, agemin,agemax, dateprev1, dateprev2,mobilav, agedeb, fage, popforecast, popfile, anpyram1,p, i1);*/
 6702:       /*      }  */
 6703:       /*      else{ */
 6704:       /*        erreur=108; */
 6705:       /*        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); */
 6706:       /*        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); */
 6707:       /*      } */
 6708:     }
 6709:   
 6710: 
 6711:     /* Computes prevalence between agemin (i.e minimal age computed) and no more ageminpar */
 6712: 
 6713:     prevalence(probs, agemin, agemax, s, agev, nlstate, imx, Tvar, nbcode, ncodemax, mint, anint, dateprev1, dateprev2, firstpass, lastpass);
 6714:     /*  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",\
 6715: 	ageminpar, agemax, s[lastpass][imx], agev[lastpass][imx], nlstate, imx, mint[lastpass][imx],anint[lastpass][imx], dateprev1, dateprev2, firstpass, lastpass);
 6716:     */
 6717: 
 6718:     if (mobilav!=0) {
 6719:       mobaverage= ma3x(1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6720:       if (movingaverage(probs, bage, fage, mobaverage,mobilav)!=0){
 6721: 	fprintf(ficlog," Error in movingaverage mobilav=%d\n",mobilav);
 6722: 	printf(" Error in movingaverage mobilav=%d\n",mobilav);
 6723:       }
 6724:     }
 6725: 
 6726: 
 6727:     /*---------- Health expectancies, no variances ------------*/
 6728: 
 6729:     strcpy(filerese,"e");
 6730:     strcat(filerese,fileres);
 6731:     if((ficreseij=fopen(filerese,"w"))==NULL) {
 6732:       printf("Problem with Health Exp. resultfile: %s\n", filerese); exit(0);
 6733:       fprintf(ficlog,"Problem with Health Exp. resultfile: %s\n", filerese); exit(0);
 6734:     }
 6735:     printf("Computing Health Expectancies: result on file '%s' \n", filerese);
 6736:     fprintf(ficlog,"Computing Health Expectancies: result on file '%s' \n", filerese);
 6737:     /*for(cptcov=1,k=0;cptcov<=i1;cptcov++){
 6738:       for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){*/
 6739:           
 6740:     for (k=1; k <= (int) pow(2,cptcoveff); k++){
 6741: 	fprintf(ficreseij,"\n#****** ");
 6742: 	for(j=1;j<=cptcoveff;j++) {
 6743: 	  fprintf(ficreseij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6744: 	}
 6745: 	fprintf(ficreseij,"******\n");
 6746: 
 6747: 	eij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
 6748: 	oldm=oldms;savm=savms;
 6749: 	evsij(eij, p, nlstate, stepm, (int) bage, (int)fage, oldm, savm, k, estepm, strstart);  
 6750:       
 6751: 	free_ma3x(eij,1,nlstate,1,nlstate,(int) bage, (int)fage);
 6752:       /*}*/
 6753:     }
 6754:     fclose(ficreseij);
 6755: 
 6756: 
 6757:     /*---------- Health expectancies and variances ------------*/
 6758: 
 6759: 
 6760:     strcpy(filerest,"t");
 6761:     strcat(filerest,fileres);
 6762:     if((ficrest=fopen(filerest,"w"))==NULL) {
 6763:       printf("Problem with total LE resultfile: %s\n", filerest);goto end;
 6764:       fprintf(ficlog,"Problem with total LE resultfile: %s\n", filerest);goto end;
 6765:     }
 6766:     printf("Computing Total Life expectancies with their standard errors: file '%s' \n", filerest); 
 6767:     fprintf(ficlog,"Computing Total Life expectancies with their standard errors: file '%s' \n", filerest); 
 6768: 
 6769: 
 6770:     strcpy(fileresstde,"stde");
 6771:     strcat(fileresstde,fileres);
 6772:     if((ficresstdeij=fopen(fileresstde,"w"))==NULL) {
 6773:       printf("Problem with Health Exp. and std errors resultfile: %s\n", fileresstde); exit(0);
 6774:       fprintf(ficlog,"Problem with Health Exp. and std errors resultfile: %s\n", fileresstde); exit(0);
 6775:     }
 6776:     printf("Computing Health Expectancies and standard errors: result on file '%s' \n", fileresstde);
 6777:     fprintf(ficlog,"Computing Health Expectancies and standard errors: result on file '%s' \n", fileresstde);
 6778: 
 6779:     strcpy(filerescve,"cve");
 6780:     strcat(filerescve,fileres);
 6781:     if((ficrescveij=fopen(filerescve,"w"))==NULL) {
 6782:       printf("Problem with Covar. Health Exp. resultfile: %s\n", filerescve); exit(0);
 6783:       fprintf(ficlog,"Problem with Covar. Health Exp. resultfile: %s\n", filerescve); exit(0);
 6784:     }
 6785:     printf("Computing Covar. of Health Expectancies: result on file '%s' \n", filerescve);
 6786:     fprintf(ficlog,"Computing Covar. of Health Expectancies: result on file '%s' \n", filerescve);
 6787: 
 6788:     strcpy(fileresv,"v");
 6789:     strcat(fileresv,fileres);
 6790:     if((ficresvij=fopen(fileresv,"w"))==NULL) {
 6791:       printf("Problem with variance resultfile: %s\n", fileresv);exit(0);
 6792:       fprintf(ficlog,"Problem with variance resultfile: %s\n", fileresv);exit(0);
 6793:     }
 6794:     printf("Computing Variance-covariance of DFLEs: file '%s' \n", fileresv);
 6795:     fprintf(ficlog,"Computing Variance-covariance of DFLEs: file '%s' \n", fileresv);
 6796: 
 6797:     /*for(cptcov=1,k=0;cptcov<=i1;cptcov++){
 6798:       for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){*/
 6799:           
 6800:     for (k=1; k <= (int) pow(2,cptcoveff); k++){
 6801:     	fprintf(ficrest,"\n#****** ");
 6802: 	for(j=1;j<=cptcoveff;j++) 
 6803: 	  fprintf(ficrest,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6804: 	fprintf(ficrest,"******\n");
 6805: 
 6806: 	fprintf(ficresstdeij,"\n#****** ");
 6807: 	fprintf(ficrescveij,"\n#****** ");
 6808: 	for(j=1;j<=cptcoveff;j++) {
 6809: 	  fprintf(ficresstdeij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6810: 	  fprintf(ficrescveij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6811: 	}
 6812: 	fprintf(ficresstdeij,"******\n");
 6813: 	fprintf(ficrescveij,"******\n");
 6814: 
 6815: 	fprintf(ficresvij,"\n#****** ");
 6816: 	for(j=1;j<=cptcoveff;j++) 
 6817: 	  fprintf(ficresvij,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6818: 	fprintf(ficresvij,"******\n");
 6819: 
 6820: 	eij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
 6821: 	oldm=oldms;savm=savms;
 6822: 	cvevsij(eij, p, nlstate, stepm, (int) bage, (int)fage, oldm, savm, k, estepm, delti, matcov, strstart);  
 6823: 	/*
 6824: 	 */
 6825: 	/* goto endfree; */
 6826:  
 6827: 	vareij=ma3x(1,nlstate,1,nlstate,(int) bage, (int) fage);
 6828: 	pstamp(ficrest);
 6829: 
 6830: 
 6831: 	for(vpopbased=0; vpopbased <= popbased; vpopbased++){ /* Done for vpopbased=0 and vpopbased=1 if popbased==1*/
 6832: 	  oldm=oldms;savm=savms; /* Segmentation fault */
 6833: 	  cptcod= 0; /* To be deleted */
 6834: 	  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 */
 6835: 	  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 ");
 6836: 	  if(vpopbased==1)
 6837: 	    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);
 6838: 	  else
 6839: 	    fprintf(ficrest,"the age specific period (stable) prevalences in each health state \n");
 6840: 	  fprintf(ficrest,"# Age e.. (std) ");
 6841: 	  for (i=1;i<=nlstate;i++) fprintf(ficrest,"e.%d (std) ",i);
 6842: 	  fprintf(ficrest,"\n");
 6843: 
 6844: 	  epj=vector(1,nlstate+1);
 6845: 	  for(age=bage; age <=fage ;age++){
 6846: 	    prevalim(prlim, nlstate, p, age, oldm, savm,ftolpl,k);
 6847: 	    if (vpopbased==1) {
 6848: 	      if(mobilav ==0){
 6849: 		for(i=1; i<=nlstate;i++)
 6850: 		  prlim[i][i]=probs[(int)age][i][k];
 6851: 	      }else{ /* mobilav */ 
 6852: 		for(i=1; i<=nlstate;i++)
 6853: 		  prlim[i][i]=mobaverage[(int)age][i][k];
 6854: 	      }
 6855: 	    }
 6856: 	
 6857: 	    fprintf(ficrest," %4.0f",age);
 6858: 	    for(j=1, epj[nlstate+1]=0.;j <=nlstate;j++){
 6859: 	      for(i=1, epj[j]=0.;i <=nlstate;i++) {
 6860: 		epj[j] += prlim[i][i]*eij[i][j][(int)age];
 6861: 		/*  printf("%lf %lf ", prlim[i][i] ,eij[i][j][(int)age]);*/
 6862: 	      }
 6863: 	      epj[nlstate+1] +=epj[j];
 6864: 	    }
 6865: 
 6866: 	    for(i=1, vepp=0.;i <=nlstate;i++)
 6867: 	      for(j=1;j <=nlstate;j++)
 6868: 		vepp += vareij[i][j][(int)age];
 6869: 	    fprintf(ficrest," %7.3f (%7.3f)", epj[nlstate+1],sqrt(vepp));
 6870: 	    for(j=1;j <=nlstate;j++){
 6871: 	      fprintf(ficrest," %7.3f (%7.3f)", epj[j],sqrt(vareij[j][j][(int)age]));
 6872: 	    }
 6873: 	    fprintf(ficrest,"\n");
 6874: 	  }
 6875: 	}
 6876: 	free_ma3x(eij,1,nlstate,1,nlstate,(int) bage, (int)fage);
 6877: 	free_ma3x(vareij,1,nlstate,1,nlstate,(int) bage, (int)fage);
 6878: 	free_vector(epj,1,nlstate+1);
 6879:       /*}*/
 6880:     }
 6881:     free_vector(weight,1,n);
 6882:     free_imatrix(Tvard,1,NCOVMAX,1,2);
 6883:     free_imatrix(s,1,maxwav+1,1,n);
 6884:     free_matrix(anint,1,maxwav,1,n); 
 6885:     free_matrix(mint,1,maxwav,1,n);
 6886:     free_ivector(cod,1,n);
 6887:     free_ivector(tab,1,NCOVMAX);
 6888:     fclose(ficresstdeij);
 6889:     fclose(ficrescveij);
 6890:     fclose(ficresvij);
 6891:     fclose(ficrest);
 6892:     fclose(ficpar);
 6893:   
 6894:     /*------- Variance of period (stable) prevalence------*/   
 6895: 
 6896:     strcpy(fileresvpl,"vpl");
 6897:     strcat(fileresvpl,fileres);
 6898:     if((ficresvpl=fopen(fileresvpl,"w"))==NULL) {
 6899:       printf("Problem with variance of period (stable) prevalence  resultfile: %s\n", fileresvpl);
 6900:       exit(0);
 6901:     }
 6902:     printf("Computing Variance-covariance of period (stable) prevalence: file '%s' \n", fileresvpl);
 6903: 
 6904:     /*for(cptcov=1,k=0;cptcov<=i1;cptcov++){
 6905:       for(cptcod=1;cptcod<=ncodemax[cptcov];cptcod++){*/
 6906:           
 6907:     for (k=1; k <= (int) pow(2,cptcoveff); k++){
 6908:     	fprintf(ficresvpl,"\n#****** ");
 6909: 	for(j=1;j<=cptcoveff;j++) 
 6910: 	  fprintf(ficresvpl,"V%d=%d ",Tvaraff[j],nbcode[Tvaraff[j]][codtab[k][j]]);
 6911: 	fprintf(ficresvpl,"******\n");
 6912:       
 6913: 	varpl=matrix(1,nlstate,(int) bage, (int) fage);
 6914: 	oldm=oldms;savm=savms;
 6915: 	varprevlim(fileres, varpl, matcov, p, delti, nlstate, stepm, (int) bage, (int) fage, oldm, savm, prlim, ftolpl,k,strstart);
 6916: 	free_matrix(varpl,1,nlstate,(int) bage, (int)fage);
 6917:       /*}*/
 6918:     }
 6919: 
 6920:     fclose(ficresvpl);
 6921: 
 6922:     /*---------- End : free ----------------*/
 6923:     if (mobilav!=0) free_ma3x(mobaverage,1, AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6924:     free_ma3x(probs,1,AGESUP,1,NCOVMAX, 1,NCOVMAX);
 6925:   }  /* mle==-3 arrives here for freeing */
 6926:  /* endfree:*/
 6927:     free_matrix(prlim,1,nlstate,1,nlstate); /*here or after loop ? */
 6928:     free_matrix(pmmij,1,nlstate+ndeath,1,nlstate+ndeath);
 6929:     free_matrix(oldms, 1,nlstate+ndeath,1,nlstate+ndeath);
 6930:     free_matrix(newms, 1,nlstate+ndeath,1,nlstate+ndeath);
 6931:     free_matrix(savms, 1,nlstate+ndeath,1,nlstate+ndeath);
 6932:     free_matrix(covar,0,NCOVMAX,1,n);
 6933:     free_matrix(matcov,1,npar,1,npar);
 6934:     /*free_vector(delti,1,npar);*/
 6935:     free_ma3x(delti3,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel); 
 6936:     free_matrix(agev,1,maxwav,1,imx);
 6937:     free_ma3x(param,1,nlstate,1, nlstate+ndeath-1,1,ncovmodel);
 6938: 
 6939:     free_ivector(ncodemax,1,NCOVMAX);
 6940:     free_ivector(Tvar,1,NCOVMAX);
 6941:     free_ivector(Tprod,1,NCOVMAX);
 6942:     free_ivector(Tvaraff,1,NCOVMAX);
 6943:     free_ivector(Tage,1,NCOVMAX);
 6944: 
 6945:     free_imatrix(nbcode,0,NCOVMAX,0,NCOVMAX);
 6946:     free_imatrix(codtab,1,100,1,10);
 6947:   fflush(fichtm);
 6948:   fflush(ficgp);
 6949:   
 6950: 
 6951:   if((nberr >0) || (nbwarn>0)){
 6952:     printf("End of Imach with %d errors and/or %d warnings\n",nberr,nbwarn);
 6953:     fprintf(ficlog,"End of Imach with %d errors and/or warnings %d\n",nberr,nbwarn);
 6954:   }else{
 6955:     printf("End of Imach\n");
 6956:     fprintf(ficlog,"End of Imach\n");
 6957:   }
 6958:   printf("See log file on %s\n",filelog);
 6959:   /*  gettimeofday(&end_time, (struct timezone*)0);*/  /* after time */
 6960:   /*(void) gettimeofday(&end_time,&tzp);*/
 6961:   rend_time = time(NULL);  
 6962:   end_time = *localtime(&rend_time);
 6963:   /* tml = *localtime(&end_time.tm_sec); */
 6964:   strcpy(strtend,asctime(&end_time));
 6965:   printf("Local time at start %s\nLocal time at end   %s",strstart, strtend); 
 6966:   fprintf(ficlog,"Local time at start %s\nLocal time at end   %s\n",strstart, strtend); 
 6967:   printf("Total time used %s\n", asc_diff_time(rend_time -rstart_time,tmpout));
 6968: 
 6969:   printf("Total time was %.0lf Sec.\n", difftime(rend_time,rstart_time));
 6970:   fprintf(ficlog,"Total time used %s\n", asc_diff_time(rend_time -rstart_time,tmpout));
 6971:   fprintf(ficlog,"Total time was %.0lf Sec.\n", difftime(rend_time,rstart_time));
 6972:   /*  printf("Total time was %d uSec.\n", total_usecs);*/
 6973: /*   if(fileappend(fichtm,optionfilehtm)){ */
 6974:   fprintf(fichtm,"<br>Local time at start %s<br>Local time at end   %s<br>\n</body></html>",strstart, strtend);
 6975:   fclose(fichtm);
 6976:   fprintf(fichtmcov,"<br>Local time at start %s<br>Local time at end   %s<br>\n</body></html>",strstart, strtend);
 6977:   fclose(fichtmcov);
 6978:   fclose(ficgp);
 6979:   fclose(ficlog);
 6980:   /*------ End -----------*/
 6981: 
 6982: 
 6983:    printf("Before Current directory %s!\n",pathcd);
 6984:    if(chdir(pathcd) != 0)
 6985:     printf("Can't move to directory %s!\n",path);
 6986:   if(getcwd(pathcd,MAXLINE) > 0)
 6987:     printf("Current directory %s!\n",pathcd);
 6988:   /*strcat(plotcmd,CHARSEPARATOR);*/
 6989:   sprintf(plotcmd,"gnuplot");
 6990: #ifdef _WIN32
 6991:   sprintf(plotcmd,"\"%sgnuplot.exe\"",pathimach);
 6992: #endif
 6993:   if(!stat(plotcmd,&info)){
 6994:     printf("Error or gnuplot program not found: '%s'\n",plotcmd);fflush(stdout);
 6995:     if(!stat(getenv("GNUPLOTBIN"),&info)){
 6996:       printf("Error or gnuplot program not found: '%s' Environment GNUPLOTBIN not set.\n",plotcmd);fflush(stdout);
 6997:     }else
 6998:       strcpy(pplotcmd,plotcmd);
 6999: #ifdef __unix
 7000:     strcpy(plotcmd,GNUPLOTPROGRAM);
 7001:     if(!stat(plotcmd,&info)){
 7002:       printf("Error gnuplot program not found: '%s'\n",plotcmd);fflush(stdout);
 7003:     }else
 7004:       strcpy(pplotcmd,plotcmd);
 7005: #endif
 7006:   }else
 7007:     strcpy(pplotcmd,plotcmd);
 7008:   
 7009:   sprintf(plotcmd,"%s %s",pplotcmd, optionfilegnuplot);
 7010:   printf("Starting graphs with: '%s'\n",plotcmd);fflush(stdout);
 7011: 
 7012:   if((outcmd=system(plotcmd)) != 0){
 7013:     printf("gnuplot command might not be in your path: '%s', err=%d\n", plotcmd, outcmd);
 7014:     printf("\n Trying if gnuplot resides on the same directory that IMaCh\n");
 7015:     sprintf(plotcmd,"%sgnuplot %s", pathimach, optionfilegnuplot);
 7016:     if((outcmd=system(plotcmd)) != 0)
 7017:       printf("\n Still a problem with gnuplot command %s, err=%d\n", plotcmd, outcmd);
 7018:   }
 7019:   printf(" Successful, please wait...");
 7020:   while (z[0] != 'q') {
 7021:     /* chdir(path); */
 7022:     printf("\nType e to edit results with your browser, g to graph again and q for exit: ");
 7023:     scanf("%s",z);
 7024: /*     if (z[0] == 'c') system("./imach"); */
 7025:     if (z[0] == 'e') {
 7026: #ifdef __APPLE__
 7027:       sprintf(pplotcmd, "open %s", optionfilehtm);
 7028: #elif __linux
 7029:       sprintf(pplotcmd, "xdg-open %s", optionfilehtm);
 7030: #else
 7031:       sprintf(pplotcmd, "%s", optionfilehtm);
 7032: #endif
 7033:       printf("Starting browser with: %s",pplotcmd);fflush(stdout);
 7034:       system(pplotcmd);
 7035:     }
 7036:     else if (z[0] == 'g') system(plotcmd);
 7037:     else if (z[0] == 'q') exit(0);
 7038:   }
 7039:   end:
 7040:   while (z[0] != 'q') {
 7041:     printf("\nType  q for exiting: ");
 7042:     scanf("%s",z);
 7043:   }
 7044: }

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