This file is indexed.

/usr/lib/ocaml/gsl/gsl_linalg.mli is in libocamlgsl-ocaml-dev 1.19.1-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
(* gsl-ocaml - OCaml interface to GSL                       *)
(* Copyright (©) 2002-2012 - Olivier Andrieu                *)
(* Distributed under the terms of the GPL version 3         *)

(** Simple linear algebra operations *)

open Gsl_vectmat
open Gsl_complex

(** {3 Simple matrix multiplication} *)

(** [matmult a ~transpa b ~transpb c] stores in matrix [c] the product
   of matrices [a] and [b]. [transpa] or [transpb] allow transposition
   of either matrix, so it can compute a.b or Trans(a).b or a.Trans(b)
   or Trans(a).Trans(b) . 

   See also {!Gsl.Blas.gemm}. *)
external matmult : 
  a:mat -> ?transpa:bool -> 
  b:mat -> ?transpb:bool -> mat -> unit
    = "ml_gsl_linalg_matmult_mod"


(** {3 LU decomposition} *)

(** {4 Low-level functions } *)

external _LU_decomp : mat -> Gsl_permut.permut -> int
    = "ml_gsl_linalg_LU_decomp"

external _LU_solve : mat -> Gsl_permut.permut -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_LU_solve"

external _LU_svx : mat -> Gsl_permut.permut -> vec -> unit
    = "ml_gsl_linalg_LU_svx"

external _LU_refine : a:mat -> lu:mat -> Gsl_permut.permut -> 
  b:vec -> x:vec -> res:vec -> unit
    = "ml_gsl_linalg_LU_refine_bc" "ml_gsl_linalg_LU_refine"

external _LU_invert : mat -> Gsl_permut.permut -> mat -> unit
    = "ml_gsl_linalg_LU_invert"

external _LU_det : mat -> int -> float
    = "ml_gsl_linalg_LU_det"

external _LU_lndet : mat -> float
    = "ml_gsl_linalg_LU_lndet"

external _LU_sgndet : mat -> int -> int
    = "ml_gsl_linalg_LU_sgndet"

(** {4 Higher-level functions} *)

(** With these, the arguments are protected (copied) and necessary
intermediate datastructures are allocated; *)

val decomp_LU : ?protect:bool ->
  [< `M of Gsl_matrix.matrix
   | `MF of Gsl_matrix_flat.matrix
   | `A of float array * int * int
   | `AA of float array array] -> mat * Gsl_permut.permut * int

val solve_LU : ?protect:bool -> 
  [< `M of Gsl_matrix.matrix
   | `MF of Gsl_matrix_flat.matrix
   | `A of float array * int * int
   | `AA of float array array] ->
  [< `A of float array 
   | `VF of Gsl_vector_flat.vector 
   | `V of Gsl_vector.vector] -> float array

val det_LU : ?protect:bool ->
  [< `M of Gsl_matrix.matrix
   | `MF of Gsl_matrix_flat.matrix
   | `A of float array * int * int
   | `AA of float array array] -> float

val invert_LU : ?protect:bool ->
  ?result:mat ->
  [< `M of Gsl_matrix.matrix
   | `MF of Gsl_matrix_flat.matrix
   | `A of float array * int * int
   | `AA of float array array] -> mat


(** {3 Complex LU decomposition} *)

external complex_LU_decomp : cmat -> Gsl_permut.permut -> int
    = "ml_gsl_linalg_complex_LU_decomp"

external complex_LU_solve : cmat -> Gsl_permut.permut -> b:cvec -> x:cvec -> unit
    = "ml_gsl_linalg_complex_LU_solve"

external complex_LU_svx : cmat -> Gsl_permut.permut -> cvec -> unit
    = "ml_gsl_linalg_complex_LU_svx"

external complex_LU_refine : a:cmat -> lu:cmat -> Gsl_permut.permut -> 
  b:cvec -> x:cvec -> res:cvec -> unit
    = "ml_gsl_linalg_complex_LU_refine_bc" "ml_gsl_linalg_complex_LU_refine"

external complex_LU_invert : cmat -> Gsl_permut.permut -> cmat -> unit
    = "ml_gsl_linalg_complex_LU_invert"

external complex_LU_det : cmat -> int -> complex
    = "ml_gsl_linalg_complex_LU_det"

external complex_LU_lndet : cmat -> float
    = "ml_gsl_linalg_complex_LU_lndet"

external complex_LU_sgndet : cmat -> int -> complex
    = "ml_gsl_linalg_complex_LU_sgndet"


(** {3 QR decomposition} *)

external _QR_decomp : mat -> vec -> unit
    = "ml_gsl_linalg_QR_decomp"

external _QR_solve : mat -> vec -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_QR_solve"

external _QR_svx : mat -> vec -> x:vec -> unit
    = "ml_gsl_linalg_QR_svx"

external _QR_lssolve : mat -> vec -> b:vec -> x:vec -> res:vec -> unit
    = "ml_gsl_linalg_QR_lssolve"

external _QR_QTvec : mat -> vec -> v:vec -> unit
    = "ml_gsl_linalg_QR_QTvec"

external _QR_Qvec : mat -> vec -> v:vec -> unit
    = "ml_gsl_linalg_QR_Qvec"

external _QR_Rsolve : mat -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_QR_Rsolve"

external _QR_Rsvx : mat -> x:vec -> unit
    = "ml_gsl_linalg_QR_Rsvx"

external _QR_unpack : mat -> tau:vec -> q:mat -> r:mat -> unit
    = "ml_gsl_linalg_QR_unpack"

external _QR_QRsolve : mat -> r:mat -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_QR_QRsolve"

external _QR_update : mat -> r:mat -> w:vec -> v:vec -> unit
    = "ml_gsl_linalg_QR_update"

external _R_solve : r:mat -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_R_solve"

(* external _R_svx : r:mat -> x:vec -> unit*)
(*     = "ml_gsl_linalg_R_svx"*)



(** {3 QR Decomposition with Column Pivoting} *)

external _QRPT_decomp : a:mat -> tau:vec -> p:Gsl_permut.permut -> norm:vec -> int
    = "ml_gsl_linalg_QRPT_decomp"

external _QRPT_decomp2 : a:mat -> q:mat -> r:mat -> tau:vec -> p:Gsl_permut.permut -> norm:vec -> int
    = "ml_gsl_linalg_QRPT_decomp2_bc" "ml_gsl_linalg_QRPT_decomp2"

external _QRPT_solve : qr:mat -> tau:vec -> p:Gsl_permut.permut -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_QRPT_solve"

external _QRPT_svx : qr:mat -> tau:vec -> p:Gsl_permut.permut -> x:vec -> unit
    = "ml_gsl_linalg_QRPT_svx"

external _QRPT_QRsolve : q:mat -> r:mat -> p:Gsl_permut.permut -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_QRPT_QRsolve"

external _QRPT_update : q:mat -> r:mat -> p:Gsl_permut.permut -> u:vec -> v:vec -> unit
    = "ml_gsl_linalg_QRPT_update"

external _QRPT_Rsolve : qr:mat -> p:Gsl_permut.permut -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_QRPT_Rsolve"

external _QRPT_Rsvx : qr:mat -> p:Gsl_permut.permut -> x:vec -> unit
    = "ml_gsl_linalg_QRPT_Rsolve"



(** {3 Singular Value Decomposition} *)

external _SV_decomp : a:mat -> v:mat -> s:vec -> work:vec -> unit
    = "ml_gsl_linalg_SV_decomp"

external _SV_decomp_mod : a:mat -> x:mat -> v:mat -> 
  s:vec -> work:vec -> unit
    = "ml_gsl_linalg_SV_decomp_mod"

external _SV_decomp_jacobi : a:mat -> v:mat -> s:vec -> unit
    = "ml_gsl_linalg_SV_decomp_jacobi"

external _SV_solve : u:mat -> v:mat -> s:vec -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_SV_solve"



(** {3 LQ decomposition} *)

external _LQ_decomp : a:mat -> tau:vec -> unit = "ml_gsl_linalg_LQ_decomp"
external _LQ_solve_T : lq:mat -> tau:vec -> b:vec -> x:vec -> unit = "ml_gsl_linalg_LQ_solve_T"
external _LQ_svx_T : lq:mat -> tau:vec -> x:vec -> unit = "ml_gsl_linalg_LQ_svx_T"
external _LQ_lssolve_T : lq:mat -> tau:vec -> b:vec -> x:vec -> res:vec -> unit = "ml_gsl_linalg_LQ_lssolve_T"
external _LQ_Lsolve_T : lq:mat -> b:vec -> x:vec -> unit = "ml_gsl_linalg_LQ_Lsolve_T"
external _LQ_Lsvx_T : lq:mat -> x:vec -> unit = "ml_gsl_linalg_LQ_Lsvx_T"
external _L_solve_T : l:mat -> b:vec -> x:vec -> unit = "ml_gsl_linalg_L_solve_T"
external _LQ_vecQ : lq:mat -> tau:vec -> v:vec -> unit = "ml_gsl_linalg_LQ_vecQ"
external _LQ_vecQT : lq:mat -> tau:vec -> v:vec -> unit = "ml_gsl_linalg_LQ_vecQT"
external _LQ_unpack : lq:mat -> tau:vec -> q:mat -> l:mat -> unit = "ml_gsl_linalg_LQ_unpack"
external _LQ_update : q:mat -> r:mat -> v:vec -> w:vec -> unit  = "ml_gsl_linalg_LQ_update"
external _LQ_LQsolve : q:mat -> l:mat -> b:vec -> x:vec -> unit = "ml_gsl_linalg_LQ_LQsolve"



(** {3 P^T L Q decomposition} *)

external _PTLQ_decomp : a:mat -> tau:vec -> Gsl_permut.permut -> norm:vec -> int = "ml_gsl_linalg_PTLQ_decomp"
external _PTLQ_decomp2 : a:mat -> q:mat -> r:mat -> tau:vec -> Gsl_permut.permut -> norm:vec -> int = "ml_gsl_linalg_PTLQ_decomp2_bc" "ml_gsl_linalg_PTLQ_decomp2"
external _PTLQ_solve_T : qr:mat -> tau:vec -> Gsl_permut.permut -> b:vec -> x:vec -> unit = "ml_gsl_linalg_PTLQ_solve_T"
external _PTLQ_svx_T : lq:mat -> tau:vec -> Gsl_permut.permut -> x:vec -> unit = "ml_gsl_linalg_PTLQ_svx_T"
external _PTLQ_LQsolve_T : q:mat -> l:mat -> Gsl_permut.permut -> b:vec -> x:vec -> unit = "ml_gsl_linalg_PTLQ_LQsolve_T"
external _PTLQ_Lsolve_T : lq:mat -> Gsl_permut.permut -> b:vec -> x:vec -> unit = "ml_gsl_linalg_PTLQ_Lsolve_T"
external _PTLQ_Lsvx_T : lq:mat -> Gsl_permut.permut -> x:vec -> unit = "ml_gsl_linalg_PTLQ_Lsvx_T"
external _PTLQ_update : q:mat -> l:mat -> Gsl_permut.permut -> v:vec -> w:vec -> unit  = "ml_gsl_linalg_PTLQ_update"



(** {3 Cholesky decomposition} *)

external cho_decomp : mat -> unit
    = "ml_gsl_linalg_cholesky_decomp"

external cho_solve : mat -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_cholesky_solve"

external cho_svx : mat -> vec -> unit
    = "ml_gsl_linalg_cholesky_svx"

external cho_decomp_unit : mat -> vec -> unit
    = "ml_gsl_linalg_cholesky_decomp_unit"


(** {3 Tridiagonal Decomposition of Real Symmetric Matrices} *)

external symmtd_decomp : a:mat -> tau:vec -> unit
    = "ml_gsl_linalg_symmtd_decomp"

external symmtd_unpack : a:mat -> tau:vec -> 
  q:mat -> diag:vec -> subdiag:vec -> unit
    = "ml_gsl_linalg_symmtd_unpack"

external symmtd_unpack_T : a:mat -> diag:vec -> subdiag:vec -> unit
    = "ml_gsl_linalg_symmtd_unpack_T"



(** {3 Tridiagonal Decomposition of Hermitian Matrices} *)

external hermtd_decomp : a:cmat -> tau:cvec -> unit
    = "ml_gsl_linalg_hermtd_decomp"

external hermtd_unpack : a:cmat -> tau:cvec -> 
  q:cmat -> diag:vec -> subdiag:vec -> unit
    = "ml_gsl_linalg_hermtd_unpack"

external hermtd_unpack_T : a:cmat -> diag:vec -> subdiag:vec -> unit
    = "ml_gsl_linalg_hermtd_unpack_T"


(** {3 Bidiagonalization} *)

external bidiag_decomp : a:mat -> tau_u:vec -> tau_v:vec -> unit
    = "ml_gsl_linalg_bidiag_decomp"

external bidiag_unpack : a:mat -> tau_u:vec -> u:mat -> tau_v:vec -> v:mat -> diag:vec -> superdiag:vec -> unit
    = "ml_gsl_linalg_bidiag_unpack_bc" "ml_gsl_linalg_bidiag_unpack"

external bidiag_unpack2 : a:mat -> tau_u:vec -> tau_v:vec -> v:mat -> unit
    = "ml_gsl_linalg_bidiag_unpack2"

external bidiag_unpack_B : a:mat -> diag:vec -> superdiag:vec -> unit
    = "ml_gsl_linalg_bidiag_unpack_B"



(** {3 Householder solver} *)

external _HH_solve : mat -> b:vec -> x:vec -> unit
      = "ml_gsl_linalg_HH_solve"

external _HH_svx : mat -> vec -> unit
    = "ml_gsl_linalg_HH_svx"

val solve_HH : ?protect:bool -> 
  [< `M of Gsl_matrix.matrix
   | `MF of Gsl_matrix_flat.matrix
   | `A of float array * int * int
   | `AA of float array array] ->
  [< `A of float array 
   | `VF of Gsl_vector_flat.vector 
   | `V of Gsl_vector.vector] -> float array



(** {3 Tridiagonal Systems} *)

external solve_symm_tridiag : diag:vec -> offdiag:vec -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_solve_symm_tridiag"

external solve_tridiag : diag:vec -> abovediag:vec -> belowdiag:vec -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_solve_tridiag"

external solve_symm_cyc_tridiag : diag:vec -> offdiag:vec -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_solve_symm_cyc_tridiag"

external solve_cyc_tridiag : diag:vec -> abovediag:vec -> belowdiag:vec -> b:vec -> x:vec -> unit
    = "ml_gsl_linalg_solve_cyc_tridiag"



(** {3 Exponential} *)

external _exponential : mat -> mat -> Gsl_fun.mode -> unit
    = "ml_gsl_linalg_exponential_ss"

val exponential : ?mode:Gsl_fun.mode -> 
    [< `M of Gsl_matrix.matrix
     | `MF of Gsl_matrix_flat.matrix
     | `A of float array * int * int] -> [ `M of Gsl_matrix.matrix]