max_inter_ll.h.html mathcode2html   
 Source file:   max_inter_ll.h
 Converted:   Sat May 9 2015 at 14:46:13
 This documentation file will not reflect any later changes in the source file.

$$\phantom{******** If you see this on the webpage then the browser could not locate *********}$$
$$\phantom{******** jsMath/easy/load.js or the variable root is set wrong in this file *********}$$
$$\newcommand{\vector}[1]{\left[\begin{array}{c} #1 \end{array}\right]}$$ $$\newenvironment{matrix}{\left[\begin{array}{cccccccccc}} {\end{array}\right]}$$ $$\newcommand{\A}{{\cal A}}$$ $$\newcommand{\W}{{\cal W}}$$

/* Copyright 2008-2012 Research Foundation State University of New York   */

/* This file is part of QUB Express.                                      */

/* QUB Express is free software; you can redistribute it and/or modify    */
/* it under the terms of the GNU General Public License as published by   */
/* the Free Software Foundation, either version 3 of the License, or      */
/* (at your option) any later version.                                    */

/* QUB Express is distributed in the hope that it will be useful,         */
/* but WITHOUT ANY WARRANTY; without even the implied warranty of         */
/* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the          */
/* GNU General Public License for more details.                           */

/* You should have received a copy of the GNU General Public License,     */
/* named LICENSE.txt, in the QUB Express program directory.  If not, see  */
/* <http://www.gnu.org/licenses/>.                                        */

#ifndef MAX_INTER_LL_H
#define MAX_INTER_LL_H

#include "maxill.h"

/*
  

Maximum Interval Likelihood (Qin et al, 1996)

This function computes the forward log-likelihood of data given a model.

See also:

Up: Index

The Model

Our Markov model is a graph with colored vertices. A vertex is called a "state," and its color is a nonnegative integer called its "class." States in the same class are indistinguishable (same amp). To describe the vertices, provide the array

  int clazz[Ns] = [class of each state]
  int Ns        = number of states
  

Each edge is labeled with its transition rate (probability per second). These form the matrix

Q, a Ns x Ns matrix with
\(Q_{a,b}\) = rate from state a to state b
\(Q_{a,a} = - \sum_i Q_{a,i}\) where \(i \neq a\)

Each \(Q_{a,b} = K0_{a,b} * Ligand_{a,b} * e^{K1_{a,b} * Voltage_{a,b}}\). You provide the Ns x Ns matrices

  double **K0, **K1            of kinetic parameters
  int    **Ligand, **Voltage   index of the ligand or voltage signal influencing each rate, or 0
  

with the diagonals undefined.

A pair of states is either connected in both directions or neither. To indicate un-connectedness, set

\(K0_{a,b} = K0_{b,a} = 0\)

For constant state entry probabilities, as described in (Qin...1996), provide the array

  int P0[Ns]
  

For equilibrium entry probabilities, as described in (Milescu...2005), provide

  P0 = NULL
  

The Data

For data, you provide one or more idealized segments. A segment is described by

  int   DwellCount            number of events
  int   Classes[DwellCount]   class index of each event
  float Durations[DwellCount] length of each dwell in seconds, less tdead
  

and you give the segments together as

  int   Nseg                  number of segments
  int   **dwellCounts
  int   **classeses
  float **durationses
  

tdead is the longest duration of events that can't be reliably detected. MIL assumes you have deleted any such events, by merging them with their prior. Also, for computational reasons, tdead should be subtracted from each event. We provide a utility which merges and shortens events, in-place:

*/

extern "C" MAXILL_API void mil_prep_events( int *dwellCount, int *classes, float *durations, double tdead_ms );
extern "C" MAXILL_API void mil_prep_segments( int Nseg, int *dwellCounts, int **classeses, float **durationses, double tdead_ms );

/*
  

If any \(Ligand_{a,b}\ \neq 0\) or \(Voltage_{a,b}\ \neq 0\), you must provide an array

  double *Constants
  

where e.g. if \(Ligand_{a,b} = 2\) then Constants[2] holds the ligand concentration. We assume all segments were recorded at the same constant conditions. For global fitting, call max_inter_ll separately for each dataset, with different constants, and sum the LL.

double *LL: upon return, contains log(prob. that model generated this data)

Returns: 0 on success, error codes to be defined.

Limitations

  • no gradients
  • no maximization
  • no support for multiple identical models (channel count > 1)

*/

extern "C" MAXILL_API int inter_ll_arr(
						 int Ns, int *clazz, double *P0_in,  double *K0_in, double *K1_in, double *K2_in,
						 int *Ligand_in, int *Voltage_in, int *Pressure_in,
						 double *Constants, double tdead,
						 int Nseg, int *dwellCounts, int **classeses, float **durationses,
						 double *ll);

extern "C" MAXILL_API int inter_ll(
					     int Ns, int *clazz, double *P0,  double **K0, double **K1, double **K2,
					     int **Ligand, int **Voltage, int **Pressure,
					     double *Constants, double tdead,
					     int Nseg, int *dwellCounts, int **classeses, float **durationses,
					     double *LL);


#endif