/**
\mainpage Comprehensive Documentation of C++ Functions and Classes (API)

TMB (Template Model Builder) is an R package for fitting statistical latent variable models to data. 
It is inspired by [ADMB](http://www.admb-project.org). Unlike most other R packages the model is formulated in C++. This
provides great flexibility, but requires some familiarity with the C/C++ programming language.

These pages describes the C++ part of TMB. The R part is documented via the ordinary R help system.


<ul>
  <li> \b Modules: Introduction and summary of TMB (start here if you are new).</li>
  <li> \b Namespaces: TMB's package system</li>
  <li> \b Classes: The C++ classes used in TMB</li>
  <li> \b Macros: The C++ Macros used in TMB</li>
  <li> \b Files: The C++ files that constitute TMB</li>
  <li> \b Examples: List of pre built models that comes with TMB</li>
  <li> \b Search \b field: Start typing a keyword, and all matching keywords will be listed. Very useful!</li>
</ul>

*/

/** \defgroup Tutorial Getting started
\brief A general introduction to TMB


A TMB project consists of an R file (*.R) and a C++ file (*.cpp). 
The R file does pre- and post processing of data in addition to maximizing
the log-likelihood contained in *.cpp. See \ref Examples for more details.
All R functions are documented within the standard help system in R.
This tutorial describes how to write the C++ file, and assumes familiarity with C++ and to some extent with R.

The purpose of the C++ program is to evaluate the objective function, i.e. the negative log-likelihood
of the model. The program is compiled and called from R, where
it can be fed to a function minimizer like \c nlminb().

The objective function should be of the following C++ type:
   \code
#include <TMB.hpp>

template<class Type>
Type objective_function<Type>::operator() ()
{
.... Here goes your C++ code ..... 
}
   \endcode
The first line includes the source code for the whole TMB package (and all its dependencies).
The objective function is a templated class where <tt><Type></tt> is
the data type of both the input values and the return value of the objective function.
This allows us to evaluate both the objective function and its derivatives using the same chunk of C++ code
(via the AD package [CppAD](http://www.coin-or.org/CppAD/Doc/doxydoc/html/index.html)). The technical aspects of this are hidden from the user. There is however one aspect 
that surprises the new TMB user. When a constant like "1.2" is used
in a calculation that affects the return value it must be "cast" to Type:
   \code
   Type nll;		// Define variable that holds the return value (neg. log. lik)
   nll = Type(1.2);	// Assign value 1.2; a cast is needed.
   \endcode

<h2>Obtaining data and parameter values from R</h2>
Obviously, we will need to pass both data and parameter values to the objective function.
This is done through a set of macros that TMB defines for us. To see
which macros are available start typing <tt>DATA_</tt> or <tt>PARAMETER_</tt> 
in the Doxygen search field of your browser (you may need to refresh the browser window between
each time you make a new search). A simple example if you want to read a vector
of numbers (doubles) is the following
   \code
  DATA_VECTOR(x);               // Vector x(0),x(1),...,x(n-1), where n is the length of x
   \endcode
Note that all vectors and matrices in TMB uses a \b zero-based indexing scheme.
It is not necessary to explicitly pass the dimension of \c x, as it can
be retrieved inside the C++ program:
   \code
  int n = x.size(); 
   \endcode

<h2>An extended C++ language</h2>
TMB extends C++ with functionality that is important for formulating
likelihood functions. You have different toolboxes available:
- Standard C++ used for infrastructure like loops etc.
- Vector, matrix and array library (see \ref matrix_arrays)
- Probability distributions (see \ref Densities and \ref R_style_distribution)

In addition to the variables defined through the <tt>DATA_</tt> or <tt>PARAMETER_</tt>
macros there can be "local" variables, for which ordinary C++ scoping rules apply.
There must also be a variable that holds the return value (neg. log. likelihood).
   \code
  DATA_VECTOR(x);               // Vector x(0),x(1),...,x(n-1)
  Type tmp = x(1);
  Type nll=tmp*tmp; 
   \endcode
As in ordinary C++ local variable tmp must be assigned a value before it can enter into
a calculation.

<h2>Statistical modelling</h2>
TMB can handle complex statistical problems with hierarchical structure (latent
random variables) and multiple data sources.  Latent random variables must be continuous 
(discrete distributions are not handled).  The <tt>PARAMETER_</tt> macros are used to pass 
two types of parameters.
- \b Parameters: to be estimated by maximum likelihood. These include fixed effects and variance
  components in the mixed model literature. They will also correspond to hyper parameters
  with non-informative priors in the Bayesian literature.
- \b Latent \b random \b variables: to be integrated out of the likelihood using a Laplace approximation. 

Which of these are chosen is controlled from R, via the <tt>random</tt> argument to the function <tt>MakeADFun</tt>. However,
on the C++ side it is usually necessary to assign a probability distribution to the parameter.

The purpose of the C++ program is to calculate the (negative) joint density of data and latent
random variables. Each datum and individual latent random effect gives a contribution
to log likelihood, which may be though of as a "distribution assignment" by users
familiar with software in the BUGS family.
   \code
  PARAMETER_VECTOR(u);          // Latent random variable 
  Type nll = Type(0);		// Return value
  nll -= dnorm(u(0),0,1,true)	// Distributional assignment: u(0) ~ N(0,1) 
   \endcode
The following rules apply:
- Distribution assignments do not need to take place before the latent variable
  is used in a calculation.
- More complicated distributional assignments are allowed, say u(0)-u(1) ~ N(0,1),
  but this requires the user to have a deeper understanding of the probabilistic aspects of the model.
- For latent variables only normal distributions should be used (otherwise
  the Laplace approximation will perform poorly). For response variables
  all probability distributions (discrete or continuous) are allowed.
  If a non-gaussian latent is needed the "transformation trick" can be used.
- The namespaces \ref R_style_distribution and \ref Densities contain many probability 
  distributions, including
  multivariate normal distributions. For probability distributions not available
  from these libraries, the user can use raw C++ code:
   \code
  DATA_VECTOR(y);          	// Data vector
  Type nll = Type(0);		// Return value
  nll -= sum(-log(Type(1.0)+y*y));	// y are i.i.d. Cauchy distributed
   \endcode

See \ref Toolbox for more about statistical modelling.


*/

/** \defgroup Structure_TMB The structure of TMB
\brief An overview of the components that constitute TMB

This documentation only covers the TMB specific code, not
[CppAD](http://www.coin-or.org/CppAD/Doc/doxydoc/html/index.html)
or 
[Eigen](http://eigen.tuxfamily.org/dox/group__TutorialMatrixClass.html)
These packages have their own documentation, which may be relevant.
In particular, some of the standard functions like \c sin() and \c cos()
are part of CppAD, and are hence not documented through TMB.


![TMB components](TMB_components.png)
*/

/** \defgroup matrix_arrays Matrices and arrays
\brief Matrices, vectors and arrays in TMB versus R


\section Relation_R Relationship to R
In R you can apply both matrix multiplication ("%*%")
and elementwise multiplication ("*") to objects of
type "matrix", i.e. it is the operator that determines the operation.
In TMB we instead have two different types of objects, while
the multiplication operator "*" is the same: 
- \c matrix: linear algebra
- \c array: elementwise operations; () and [] style indexing.
- \c vector: can be used in linear algebra with \c matrix, but at the same
time admits R style element-wise operations.

See \ref matrix_arrays.cpp "matrix_arrays.cpp" for examples of use.


\section Relation_Eigen Relationship to Eigen
The TMB types \c matrix and \c array (in dimension 2) inherits directly from the the Eigen
types Matrix and Array. The advanced user of TMB will benefit from familiarity with
the <a href="http://eigen.tuxfamily.org/dox/group__TutorialMatrixClass.html">Eigen documentation</a>.
Note: Arrays of dimension 3 or higher are specially implemented in TMB, i.e. are 
not directly inherited from Eigen.

*/

/** \defgroup R_style_distribution R style probability distributions

\brief Probability distributions (discrete and continuous) for use in the likelihood function. 

Attempts have been made to make the interface (function name and arguments) as close as possible 
to that of R. 

- The densities (\c d...) are provided both in the discrete and continuous case, 
	cumulative distributions (\c p...) and inverse cumulative distributions (\c q...) 
	are provided only for continuous distributions.
- Scalar and \c vector arguments (in combination) are supported, 
but not \c array or \c matrix arguments.
- The last argument (of type \c int) corresponds to the \c log argument 
  in R: 1=logaritm, 0=ordinary scale. \c true (logaritm) and \c false (ordinary scale) can 
  also be used.
- Vector arguments must all be of the same length (no recycling of elements). If vectors of different lengths are
  used an "out of range" error will occur, which can be picked up by the debugger.
- \c DATA_IVECTOR() and \c DATA_INTEGER() cannot be used with probability distributions, except possibly for the last (log) argument.
- Some example:
  \code 
  DATA_SCALAR(y);
  DATA_VECTOR(x);
  vector<Type> rate(10);
  matrix<Type> rate_matrix(10,10);
  dexp(y,rate,true);                      // OK, returns vector of length 10 of log-densities
  dexp(x,rate,true);                      // OK if x is length 10 
  dexp(x,rate_matrix,true);               // Error, matrix arguments not allowed
  \endcode 
- To sum over elements in the vector returned use \code sum(dexp(x,rate));\endcode 
*/

/** \defgroup special_functions Special mathematical functions
    \brief Special mathematical functions for which derivatives are implemented
*/

/** \defgroup matrix_functions Special matrix functions
    \brief Special matrix functions for which derivatives are implemented

    These matrix operations are designed to work with automatic
    differentiation for large dense matrices.
*/

/** \defgroup Densities Multivariate distributions

\brief Multivariate normal distributions, AR-processes, GMRF's, separable covariance functions, etc.

The namespace
  \code 
  using namespace density;
  \endcode 
gives access to a variety of multivariate normal distributions:
 - Multivariate normal distributions specified via a covariance matrix (structured or unstructured).
 - Autoregressive (AR) processes.
 - Gaussian Markov random fields (GMRF) defined on regular grids or defined via
  	    a (sparse) precision matrix.
 - Separable covariance functions, i.e. time-space separability.
  
These seemingly unrelated concepts are all implemented via the notion of a \c distribution, 
which explains why they are placed in the same namespace. You can combine two
\c distributions, and this lets you build up complex multivariate distributions using extremely
compact notation. Due to the flexibility of the approach it is more abstract than
other parts of TMB, but here it will be explained from scratch. Before looking 
at the different categories of multivariate distributions we note the following
which is of practical importance:
- All members in the \c density namespace return the \b negative log density,
  opposed to the univariate densities in \ref R_style_distribution.
  
<h2>Multivariate normal distributions</h2>
Consider a zero-mean multivariate normal distribution
with covariance matrix *Sigma* (symmetric positive definite),
that we want to evaluate at _x_:
  \code 
  int n = 10;
  vector<Type> x(n);                     // Evaluation point           
  x.fill(0.0);                           // Point of evaluation: x = (0,0,...,0)
  \endcode


The negative log-normal density is evaluated as follows:
\code 
  using namespace density;
  matrix<Type> Sigma(n,n);     // Covariance matrix
  // ..... User must assign value to Sigma here
  res = MVNORM(Sigma)(x);      // Evaluate negative log likelihod
\endcode
In the last line \c MVNORM(Sigma) should be interpreted as
a multivariate density, which via the last parenthesis \c (x) is evaluated at \c x.
A less compact way of expressing this is
\code 
  MVNORM_t<Type> your_dmnorm(Sigma);
  res = your_dmnorm(x);                
\endcode
in which \c your_dmnorm is a variable that holds the "density". 
    
_Sigma_ can be parameterized in different ways. Due to the symmetry of *Sigma* there
are at most *n(n+1)/2* free parameters (*n* variances and *n(n-1)/2*
correlation parameters). If you want to estimate all of
these freely (modulo the positive definite constraint) 
you can use \ref UNSTRUCTURED_CORR to specify the correlation matrix, 
and \ref VECSCALE to specify variances.
\ref UNSTRUCTURED_CORR takes as input a vector a dummy parameters that internally 
is used to build the correlation matrix via its cholesky factor. 
  \code 
  using namespace density;
  int n = 10;
  vector<Type> unconstrained_params(n*(n-1)/2);  // Dummy parameterization of correlation matrix         
  vector<Type> sds(n);                   // Standard deviations
  res = VECSCALE(UNSTRUCTURED_CORR(unconstrained_params),sds)(x);
  \endcode
If all elements of \c dummy_params are estimated we are in effect
estimating a full correlation matrix without any constraints
on its elements (except for the mandatory positive definiteness).
The actual value of the correlation matrix, but not the full covariance matrix,
can easily be assessed
using the \c ".cov()" operator
  \code 
  matrix<Type> Sigma(n,n);      
  Sigma = UNSTRUCTURED_CORR(unconstrained_params).cov();
  REPORT(Sigma);                                         // Report back to R session
  \endcode


<h2>Autoregressive processes</h2>
Consider a stationary univariate Gaussian AR1 process *x(t), t=0,...,n-1*. The stationary
distribution is choosen so that:
- *x(t)* has mean 0 and variance 1 (for all *t*).

The multivariate density of the vector *x* can be evaluated as follows
  \code 
  int n = 10;
  using namespace density;
 
  vector<Type> x(n);           // Evaluation point           
  x.fill(0.0);                    // Point of evaluation: x = (0,0,...,0)
  Type rho = 0.2;              // Correlation parameter
  res = AR1(rho)(x);           // Evaluate negative log-density of AR1 process at point x 
  \endcode	

Due to the assumed stationarity the correlation parameter must satisfy: 
- \b Stationarity constraint: -1 < *rho* < 1

Note that *cor[x(t),x(t-1)] = rho*. 

The \c SCALE function can be used to set the standard deviation.
  \code 
  Type sigma = 2.1;            // standard deviation of x
  res = SCALE(AR1(rho),sigma)(x);  
  \endcode
Now, *var[x(t)] = sigma^2*. Because all elements of \c x are scaled by the
same constant we use SCALE rather than VECSCALE. 

<h4>Multivariate AR1 processes</h4>
This is the first real illustrations of how distributions can be used as building blocks
to obtain more complex distributions. Consider the *p* dimensional AR1 process
  \code
  int n = 10;                   // Number of time steps
  int p=3;                      // dim(x)
  array<Type> x(p,n);           // Evaluation point       
  \endcode
The columns in \c x refer to the different time points.   
We then evaluate the (negative log) joint density of the time series.
  \code
  MVNORM_t<Type> your_dmnorm(Sigma);  // Density of x(t) 
  Type phi;                           // Correlation parameter
  res = AR1(phi,your_dmnorm)(x);
  \endcode
Note the following:
- We have introduced an intermediate variable \c your_dmnorm, which
  holds the p-dim density marginal density of *x(t)*. 
  This is a zero-mean normal density with covariance matrix \c Sigma.
- All *p* univarite time series have the same serial correlation \phi. 
- The multivariate process *x(t)* is stationary in the same sense as the 
   univariate AR1 process described above.

<h4>Higher order AR processes</h4>
  
There also exists \ref ARk_t of arbitrary autoregressive order. 

<h2>Gaussian Markov random fields (GMRF)</h2>
 GMRF may be defined in two ways:
 1. Via a (sparse) precision matrix Q 
 2. Via a d-dimensional lattice
 
 For further details please see \ref GMRF_t. Under 1) a sparse Q 
 corresponding to a Matern covariance function can be
 obtained via the \ref R_inla namespace.

<h2>Separable construction of covariance (precision) matrices</h2>
A typical use of separability is to create space-time models with
a sparse precision matrix. Details are given in \ref SEPARABLE_t. Here
we give a simple example.

Assume that we study a quantity \c x that changes both in space
and time. For simplicity we consider only a one-dimensional space.
We discretize space and time using equidistant grids, and assume
that the distance between grid points is 1 in both dimensions. We then
define an AR1(rho_s) process in space and one in time AR1(rho_t).
The separable assumption is that two points \c x1 and \c x2, separated
in space by a distance \c ds and in time by a distance \c dt, have correlation
given by

\c rho_s^ds\*rho_t^dt

This is implemented as
\code
  using namespace density;
  int n_s = 10;                   // Number of grid points in space
  int n_t = 10;                   // Number of grid points in time
  Type rho_s = 0.2;               // Correlation in space
  Type rho_t = 0.4;               // Correlation in time
 
  array<Type> x(n_s,n_t);         
  x.setZero();                    // x = 0

  res = SEPARABLE(AR1(rho_t),AR1(rho_s))(x);
\endcode

Note that the arguments to \c SEPARABLE() are given in the
opposite order to the dimensions of \c x.

*/

/** \defgroup Examples Example collection
\brief Overview of available examples.

- A list of all examples is found on the "Examples" tab on the top of the page. 
- Locations of example files: \c adcomp/tmb_examples and \c adcomp/TMB/inst/examples.
- For each example there is both a \c .cpp is a \c .R
- To run an example, e.g. \c ar1_4D, use the R command \code source("ar1_4D.R") \endcode.

\example nan_error_ex.cpp
\example matrix_arrays.cpp
\example multivariate_distributions.cpp

\example adaptive_integration.cpp
\example ar1_4D.cpp
\example hmm.cpp
\example hmm_filter.hpp
\example laplace.cpp
\example linreg.cpp
\example linreg_parallel.cpp
\example longlinreg.cpp
\example lr_test.cpp
\example matern.cpp
\example mvrw.cpp
\example mvrw_sparse.cpp
\example nmix.cpp
\example orange_big.cpp
\example sam.cpp
\example sde_linear.cpp
\example sdv_multi_compact.cpp
\example sdv_multi.cpp
\example socatt.cpp
\example spatial.cpp
\example spde_aniso.cpp
\example spde.cpp
\example thetalog.cpp
\example transform2.cpp
\example transform.cpp
\example transform_parallel.cpp
\example tweedie.cpp

\example MVRandomWalkValidation.cpp
\example randomwalkvalidation.cpp
\example rickervalidation.cpp

*/



/** \defgroup Errors Compilation and run time errors
\brief How to debug TMB programs

The R interface to the debugger (gdb) is documented as part of the R help system,
i.e. you can type \c "?gdbsource" in R to get info. The current document
only adresses isses that the relate to C++.

<h2>Compilation errors</h2>
It may be hard to understand the compilation errors for the following reasons
- The Eigen libraries use templated C++ which generate non-intuitive error messages.

<h2>Run time errors</h2>
Run time errors are broadly speaking of two types:
-  Out-of-bounds (you are "walking out of an array")
-  Floating point exceptions 

You can use the debugger to locate both types of errors, but the procedure is
a little bit different in the two cases. The following assumes
that you have the GNU debugger \c gdb installed.

<h3>Out-of-bounds error</h3>
An example is:
   \code
        vector<Type> y(4);
        y(5);                // 5 is not a valid index value here
   \endcode
This will cause TMB and R to crash with the following error message:
   \code
TMB has received an error from Eigen. The following condition was not met:
index >= 0 && index < size()
Please check your matrix-vector bounds etc., or run your program through a debugger.
Aborted (core dumped)
   \endcode
So, you must restart R and give the commands
\code
> library(TMB)
> gdbsource("my_project.R")
#5  objective_function<double>::operator() (this=<optimized out>) at nan_error_ex.cpp:11 
\endcode
and you can see that the debugger points to line number 11 in the .cpp file.
\c gdbsource() is an R function that is part of TMB.

<h3>Floating point execption</h3>
If you on the other hand perform an illegal mathematical operation, such as
   \code
        Type f = sqrt(-1.);
   \endcode
R will not crash, but the objective function will return a NaN value.
However, you will not know in which part of your C++ code the error occured.
By including the \c fenv.h library (part of many C++ compilers, 
but can otherwise be downloaded from http://www.scs.stanford.edu/histar/src/uinc/fenv.h)
 
\c nan_error_ex.cpp:
\include nan_error_ex.cpp

a floating point exception will be turned into an actual error that can be picked up by the debugger.
There are only two extra lines that need to be included ("//Extra line needed" in the above example).
 

When we try to run this program in the usual way, the program crashes:
\code
> source("nan_error_ex.R")
Floating point exception (core dumped)
tmp3>   
\endcode

At this stage you should run the debugger to find out that the floating point exception occurs at line number 14:
\code
> library(TMB)
Loading required package: Matrix
> gdbsource("nan_error_ex.R")
#1  0x00007ffff0e7eb09 in objective_function<double>::operator() (this=<optimized out>) at nan_error_ex.cpp:14
\endcode


*/

/** \defgroup Toolbox Statistical toolbox
\brief Collection of useful statistical techniques


First read the *Statistical Modelling* section of \ref Tutorial.



<h3>Non-normal latent variables (random effects)</h3>

The underlying latent random variables in TMB must be Gaussian for the Laplace
approximation to be accurate. To obtain other distributions,
say a gamma distribution, the "transformation trick" can be used. We start
out with normally distributed variables \c u and transform these into
new variables \c w via the \c pnorm and \c qgamma functions as follows:
\code
  PARAMETER_VECTOR(u);                             // Underlying latent random variables 
  Type nll=Type(0.0);
  nll -= sum(dnorm(u,Type(0),Type(1),true));       // Assign N(0,1) distribution u 
  vector<Type> v = pnorm(u,Type(0),Type(1),true);  // Uniformly distributed variables (on [0,1])
  vector<Type> w = qgamma(v,shape,scale,true);
\endcode
\c w now has a gamma distribution.

<h3>Discrete latent variables</h3>

The Laplace approximation can not be applied to discrete latent variables that
occur in mixture models and HMMs (Hidden Markov models). However,
such likelihoods have analytic expressions, and may be coded up in TMB.
TMB would still calculate the exact gradient of the HMM likelihood.

<h3>Mixture models </h3>

Althought mixture mode are a special case of discrete latent variable models,
they do diserve special attention. Consider the case that we want a mixture
of two zero-mean normal distributions (with different standard deviations). 
This can be implemented as:
\code
  DATA_VECTOR(x);              			
  PARAMETER_VECTOR(sigma);      // sigma0 og sigma1
  PARAMETER(p);                 // Mixture proportion of model 0
  Type nll=Type(0.0);
  nll -= sum(log(p*dnorm(x,Type(0),sigma(0),false)
  					+(Type(1)-p)*dnorm(x,Type(0),sigma(1),false)));
\endcode

<h3>Time series</h3>
Autoregressive (AR) processes may be implemented using the compact notation
of section \ref  Densities. The resulting AR process may be applied
both in the observational part and in the distribution of a latent variable.

Nonlinear time must be implemented from scratch, as in the example \ref thetalog.cpp

<h3>Spatial models</h3>
TMB has strong support for spatial model and space-time models via the \ref GMRF
and \ref SEPARABLE functions, and the notion of a \c distribution. The reader
is referred to section \ref  Densities for details and examples.

*/

/** \defgroup parallel Parallel computations
\brief How to parallelize the likelihood function using TMB.
*/

/** \defgroup CppTutorial C++ tutorial
\brief Help for R users to write the C++ user template \b and for C++ programmers about
  TMB specifics

<h3>I know R but not C++</h3>

Summary of how syntax differs between R and C++:
\code
              R code             C++/TMB code                
              
Comments      #                  //                          // Comment symbol
Constants     3.4                Type(3.4);                  // Explicit casting recommended in TMB              
Scalar        x = 5.2            Type x = Type(5.2);         // Variables must have type
Arrays        x = numeric(10)    vector<Type> x(10);         // C++ code here does NOT initialize to 0     
Indexing      x[1]+x[10]         x(0)+x(9);                  // C++ indexing is zero-based
Loops         for(i in 1:10)     for(i=1;i<=10;i++)          // Integer i must be declared in C++    
Increments    x[1] = x[1] + 3    x(0) += 3.0;                // += -= *= /= incremental operators in C++

\endcode

<h3>I know C++</h3>
TMB specific C++ include:
- You should not use \c if(x) statements where \c x is 
   a \c PARAMETER, or is derived from a variable of type
   \c PARAMETER. (It is OK to use \c if on \c DATA types, however.)
   TMB will remove the \c if(x) statement, so the code will
   produce unexpected results.
 
*/

/** \defgroup macros C++ Macros
\brief C++ macros used to enable interchange of data structures 
between R and C++.
*/
