/* linalg/gsl_linalg.h
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000, 2006, 2007, 2019 Gerard Jungman, Brian Gough, Patrick Alken
*
* This program 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.
*
* This program 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
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
#ifndef __GSL_LINALG_H__
#define __GSL_LINALG_H__
#include <stdlib.h>
#include <gsl/gsl_mode.h>
#include <gsl/gsl_permutation.h>
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_inline.h>
#include <gsl/gsl_blas.h>
#undef __BEGIN_DECLS
#undef __END_DECLS
#ifdef __cplusplus
#define __BEGIN_DECLS extern "C" {
#define __END_DECLS }
#else
#define __BEGIN_DECLS /* empty */
#define __END_DECLS /* empty */
#endif
__BEGIN_DECLS
typedef enum
{
GSL_LINALG_MOD_NONE = 0,
GSL_LINALG_MOD_TRANSPOSE = 1,
GSL_LINALG_MOD_CONJUGATE = 2
}
gsl_linalg_matrix_mod_t;
/* Note: You can now use the gsl_blas_dgemm function instead of matmult */
/* Simple implementation of matrix multiply.
* Calculates C = A.B
*
* exceptions: GSL_EBADLEN
*/
int gsl_linalg_matmult (const gsl_matrix * A,
const gsl_matrix * B,
gsl_matrix * C);
/* Simple implementation of matrix multiply.
* Allows transposition of either matrix, so it
* can compute A.B or Trans(A).B or A.Trans(B) or Trans(A).Trans(B)
*
* exceptions: GSL_EBADLEN
*/
int gsl_linalg_matmult_mod (const gsl_matrix * A,
gsl_linalg_matrix_mod_t modA,
const gsl_matrix * B,
gsl_linalg_matrix_mod_t modB,
gsl_matrix * C);
/* Calculate the matrix exponential by the scaling and
* squaring method described in Moler + Van Loan,
* SIAM Rev 20, 801 (1978). The mode argument allows
* choosing an optimal strategy, from the table
* given in the paper, for a given precision.
*
* exceptions: GSL_ENOTSQR, GSL_EBADLEN
*/
int gsl_linalg_exponential_ss(
const gsl_matrix * A,
gsl_matrix * eA,
gsl_mode_t mode
);
/* Householder Transformations */
double gsl_linalg_householder_transform (gsl_vector * v);
double gsl_linalg_householder_transform2 (double * alpha, gsl_vector * v);
gsl_complex gsl_linalg_complex_householder_transform (gsl_vector_complex * v);
int gsl_linalg_householder_hm (double tau,
const gsl_vector * v,
gsl_matrix * A);
int gsl_linalg_householder_mh (double tau,
const gsl_vector * v,
gsl_matrix * A);
int gsl_linalg_householder_hv (double tau,
const gsl_vector * v,
gsl_vector * w);
int gsl_linalg_householder_left(const double tau,
const gsl_vector * v,
gsl_matrix * A,
gsl_vector * work);
int gsl_linalg_householder_right(const double tau,
const gsl_vector * v,
gsl_matrix * A,
gsl_vector * work);
int gsl_linalg_householder_hm1 (double tau,
gsl_matrix * A);
int gsl_linalg_complex_householder_hm (gsl_complex tau,
const gsl_vector_complex * v,
gsl_matrix_complex * A);
int gsl_linalg_complex_householder_mh (gsl_complex tau,
const gsl_vector_complex * v,
gsl_matrix_complex * A);
int gsl_linalg_complex_householder_hv (gsl_complex tau,
const gsl_vector_complex * v,
gsl_vector_complex * w);
int gsl_linalg_complex_householder_left (const gsl_complex tau,
const gsl_vector_complex * v,
gsl_matrix_complex * A,
gsl_vector_complex * work);
/* Hessenberg reduction */
int gsl_linalg_hessenberg_decomp(gsl_matrix *A, gsl_vector *tau);
int gsl_linalg_hessenberg_unpack(gsl_matrix * H, gsl_vector * tau,
gsl_matrix * U);
int gsl_linalg_hessenberg_unpack_accum(gsl_matrix * H, gsl_vector * tau,
gsl_matrix * U);
int gsl_linalg_hessenberg_set_zero(gsl_matrix * H);
int gsl_linalg_hessenberg_submatrix(gsl_matrix *M, gsl_matrix *A,
size_t top, gsl_vector *tau);
/* Hessenberg-Triangular reduction */
int gsl_linalg_hesstri_decomp(gsl_matrix * A, gsl_matrix * B,
gsl_matrix * U, gsl_matrix * V,
gsl_vector * work);
/* Singular Value Decomposition
* exceptions:
*/
int
gsl_linalg_SV_decomp (gsl_matrix * A,
gsl_matrix * V,
gsl_vector * S,
gsl_vector * work);
int
gsl_linalg_SV_decomp_mod (gsl_matrix * A,
gsl_matrix * X,
gsl_matrix * V,
gsl_vector * S,
gsl_vector * work);
int gsl_linalg_SV_decomp_jacobi (gsl_matrix * A,
gsl_matrix * Q,
gsl_vector * S);
int
gsl_linalg_SV_solve (const gsl_matrix * U,
const gsl_matrix * Q,
const gsl_vector * S,
const gsl_vector * b,
gsl_vector * x);
int gsl_linalg_SV_leverage(const gsl_matrix *U, gsl_vector *h);
/* LU Decomposition, Gaussian elimination with partial pivoting
*/
int gsl_linalg_LU_decomp (gsl_matrix * A, gsl_permutation * p, int *signum);
int gsl_linalg_LU_solve (const gsl_matrix * LU,
const gsl_permutation * p,
const gsl_vector * b,
gsl_vector * x);
int gsl_linalg_LU_svx (const gsl_matrix * LU,
const gsl_permutation * p,
gsl_vector * x);
int gsl_linalg_LU_refine (const gsl_matrix * A,
const gsl_matrix * LU,
const gsl_permutation * p,
const gsl_vector * b,
gsl_vector * x,
gsl_vector * work);
int gsl_linalg_LU_invert (const gsl_matrix * LU,
const gsl_permutation * p,
gsl_matrix * inverse);
int gsl_linalg_LU_invx (gsl_matrix * LU, const gsl_permutation * p);
double gsl_linalg_LU_det (gsl_matrix * LU, int signum);
double gsl_linalg_LU_lndet (gsl_matrix * LU);
int gsl_linalg_LU_sgndet (gsl_matrix * lu, int signum);
/* Banded LU decomposition */
int gsl_linalg_LU_band_decomp (const size_t M, const size_t lb, const size_t ub, gsl_matrix * AB, gsl_vector_uint * piv);
int gsl_linalg_LU_band_solve (const size_t lb, const size_t ub, const gsl_matrix * LUB,
const gsl_vector_uint * piv, const gsl_vector *
没有合适的资源?快使用搜索试试~ 我知道了~
GSL库VS2015 x64编译结果,静态和动态两种都有。
共557个文件
h:530个
cmake:9个
lib:8个
需积分: 21 7 下载量 100 浏览量
2022-03-05
16:58:13
上传
评论 1
收藏 10.67MB ZIP 举报
温馨提示
GSL库VS2015 x64编译结果,静态和动态两种都有。
资源详情
资源评论
资源推荐
收起资源包目录
GSL库VS2015 x64编译结果,静态和动态两种都有。 (557个子文件)
gsl-targets.cmake 3KB
gsl-targets.cmake 3KB
gsl-config-version.cmake 3KB
gsl-config-version.cmake 3KB
gsl-targets-release.cmake 1KB
gsl-targets-debug.cmake 1KB
gsl-targets-release.cmake 1KB
gsl-config.cmake 372B
gsl-config.cmake 372B
gsl.dll 5.9MB
gsl.dll 2.53MB
gslcblas.dll 570KB
gslcblas.dll 440KB
gsl-config 1KB
gsl-config 1KB
gsl_linalg.h 36KB
gsl_linalg.h 36KB
gsl_cblas.h 33KB
gsl_cblas.h 33KB
gsl_blas.h 22KB
gsl_blas.h 22KB
gsl_matrix_complex_long_double.h 16KB
gsl_matrix_complex_long_double.h 16KB
gsl_multifit.h 15KB
gsl_multifit.h 15KB
gsl_matrix_complex_float.h 14KB
gsl_matrix_complex_float.h 14KB
gsl_sf_bessel.h 14KB
gsl_sf_bessel.h 14KB
gsl_matrix_long_double.h 14KB
gsl_matrix_long_double.h 14KB
gsl_integration.h 14KB
gsl_integration.h 14KB
gsl_odeiv2.h 14KB
gsl_odeiv2.h 14KB
gsl_eigen.h 13KB
gsl_eigen.h 13KB
gsl_matrix_ushort.h 13KB
gsl_matrix_ushort.h 13KB
gsl_matrix_complex_double.h 13KB
gsl_matrix_complex_double.h 13KB
gsl_multilarge_nlinear.h 13KB
gsl_multilarge_nlinear.h 13KB
gsl_matrix_ulong.h 13KB
gsl_matrix_uchar.h 13KB
gsl_matrix_ulong.h 13KB
gsl_matrix_uchar.h 13KB
gsl_matrix_float.h 13KB
gsl_matrix_short.h 13KB
gsl_matrix_float.h 13KB
gsl_matrix_short.h 13KB
gsl_matrix_uint.h 13KB
gsl_matrix_uint.h 13KB
gsl_matrix_long.h 13KB
gsl_matrix_char.h 13KB
gsl_matrix_long.h 13KB
gsl_matrix_char.h 13KB
gsl_matrix_int.h 12KB
gsl_matrix_int.h 12KB
gsl_sf_legendre.h 12KB
gsl_sf_legendre.h 12KB
gsl_multifit_nlinear.h 12KB
gsl_multifit_nlinear.h 12KB
gsl_matrix_double.h 12KB
gsl_matrix_double.h 12KB
gsl_multifit_nlin.h 11KB
gsl_multifit_nlin.h 11KB
gsl_randist.h 10KB
gsl_randist.h 10KB
gsl_vector_complex_long_double.h 10KB
gsl_vector_complex_long_double.h 10KB
gsl_interp2d.h 9KB
gsl_interp2d.h 9KB
gsl_vector_complex_float.h 9KB
gsl_vector_complex_float.h 9KB
gsl_statistics_long_double.h 9KB
gsl_statistics_long_double.h 9KB
gsl_vector_long_double.h 9KB
gsl_vector_long_double.h 9KB
gsl_vector_complex_double.h 8KB
gsl_vector_complex_double.h 8KB
gsl_statistics_float.h 8KB
gsl_statistics_float.h 8KB
gsl_vector_ushort.h 8KB
gsl_vector_ushort.h 8KB
gsl_vector_uchar.h 8KB
gsl_vector_uchar.h 8KB
gsl_vector_ulong.h 8KB
gsl_vector_ulong.h 8KB
gsl_statistics_double.h 8KB
gsl_statistics_double.h 8KB
gsl_sf_gamma.h 8KB
gsl_sf_gamma.h 8KB
gsl_odeiv.h 8KB
gsl_odeiv.h 8KB
gsl_vector_uint.h 8KB
gsl_vector_uint.h 8KB
gsl_vector_float.h 8KB
gsl_vector_short.h 8KB
gsl_vector_float.h 8KB
共 557 条
- 1
- 2
- 3
- 4
- 5
- 6
仟人斩
- 粉丝: 4805
- 资源: 37
上传资源 快速赚钱
- 我的内容管理 展开
- 我的资源 快来上传第一个资源
- 我的收益 登录查看自己的收益
- 我的积分 登录查看自己的积分
- 我的C币 登录后查看C币余额
- 我的收藏
- 我的下载
- 下载帮助
安全验证
文档复制为VIP权益,开通VIP直接复制
信息提交成功
评论0