8.9 KiB
C++ Library of the Linear Conjugate Gradient Methods (LibLCG) 说明文档
张壹(yizhang-geo@zju.edu.cn)
浙江大学地球科学学院·地球物理研究所
此说明仅覆盖算法库的简单介绍及使用,更详细的内容请查看代码注释。如果还有问题,请发邮件联系我。同时也欢迎有兴趣的同学加入开发团队!
简介
liblcg 是一个高效的、可扩展的 C++ 线性共轭梯度算法库,在原生数据结构接口的基础上,同时提供基于Eigen3和CUDA的算法接口,可以方便的实现基于CPU或GPU并行的加速计算,其中基于Eigen3的算法包含了稠密与稀疏矩阵的实现,而基于CUDA的算法主要为稀疏矩阵的实现。liblcg 包含多种实数与复数域的共轭梯度算法与其他一些迭代求解方法。目前已有得方法包括共轭梯度法、预优的共轭梯度算法、共轭梯度平方算法、双稳共轭梯度算法、BB步共轭梯度投影法与SPG共轭梯度投影法;复数域的双共轭梯度法、共轭梯度平方法、预优的共轭梯度法与TFQMR法。共轭梯度法广泛应用于无约束与不等式约束的线性最优化问题,拥有优良的收敛与计算效率。
共轭梯度算法可用于求解如下形式的线性方程组:
Ax = B
其中,A 是一个 N 阶的方阵、x 为 N*1 大小的待求解的模型向量,B 为 N*1 大小的需拟合的目标向量。需要注意的是,不同种类的共轭梯度算法对A可能有不同的要求,比如必须是正定的,或者对称的。不同算法的具体要求可以查阅其他参考文献或者查看代码中的注释。
安装
算法库使用 CMake 工具进行汇编,可在不同操作平台生成相应的Makefile或工程文件。
编译选项
算法库目前可用的编译选项有:
- LibLCG_OPENMP:是否使用OpenMP进行加速,需要安装OpeMP。默认为ON。
- LibLCG_EIGEN:是否编译基于Eigen的算法与借口,需要安装Eigen。默认为ON。
- LibLCG_STD_COMPLEX:是否使用std::complex<double>作为复数的默认类型。默认为ON。
- LibLCG_CUDA:是否编译基于CUDA的算法与借口,需要安装CUDA。默认为ON。
用户可以使用cmake命令中的-D选项对编译选项进行设置,比如关闭LibLCG_Eigen:
cmake -DLibLCG_EIGEN=OFF
Linux 与 MacOS
liblcg的默认安装路径为 /usr/local。头文件与动态库分别安装于 include 与 lib 文件夹。具体的编译与安装步骤如下:
- 下载安装CMake软件;
- 下载安装GCC编译器(常见系统已内置);
- 在源文件路径内使用如下命令进行编译与安装:
mkdir build && cd build && cmake .. && make install
Windows
MinGW 和 GCC
Windows系统不包含GNU编译环境,用户需自行下载并配置。方法如下:
- 下载MinGW安装文件,并选择gcc、pthreads与make相关软件包安装;
- 下载安装CMake软件;
- 添加CMake与MinGW可执行文件路径至Windows环境变量;
- 在源文件路径内使用如下命令进行编译与安装:
mkdir build && cd build && cmake .. -G "MinGW Makefiles" && make install
默认的安装路径为C:/Program\ Files。头文件与动态库分别安装于 include 与 lib 文件夹。
注意:用户需要手动添加头文件与动态库地址到计算机的环境变量中。
Visual Studio
用户可使用CMake工具构建VS工程文件并编译使用动态库。方法如下:
- 下载安装 Visual Studio 软件;
- 下载安装CMake软件;
- 在源文件路径内使用如下命令生成VS工程文件:
mkdir build && cd build && cmake .. -G "Visual Studio 16 2019"
注:如需生成其他版本的VS工程文件,请使用-G命令查看相应的识别码。
- 使用 Visual Studio 打开.sln工程文件并编译动态库。
使用与编译
用户使用库函数时需在源文件中引入相应的头文件,如:
#include "lcg/lcg.h"
编译可执行文件时需链接lcg动态库。以g++为例:
g++ example.cpp -llcg -o example_out
快速开始
要使用liblcg求解线性方程组Ax=B,用户需要定义Ax乘积的计算函数(回调函数),该函数的功能为计算不同的x所对应的乘积Ax。以实数类型的共轭梯度算法为例,其回调函数的接口定义为:
typedef void (*lcg_axfunc_ptr)(void* instance, const lcg_float* x, lcg_float* prod_Ax, const int n_size);
其中,x
为输入的向量,prod_Ax
为返回的乘积向量,n
为这两个向量的长度。注意此处参数列表中并不包含矩阵A,这意味这A必须为全局或者类变量。这样设计的主要原因是在某些复杂最优化问题的编程中,计算并存储A并不实际或者划算,此时一般采用的策略是存储相关变量且仅计算Ax的乘积,所以矩阵A并不总是存在。
用户在定义Ax计算函数后即可调用求解函数 lcg_solver() 对线性方程组进行求解。以无约束的求解函数为例,其声明如下:
int lcg_solver(lcg_axfunc_ptr Afp, lcg_progress_ptr Pfp, lcg_float* m, const lcg_float* B, const int n_size,
const lcg_para* param, void* instance, lcg_solver_enum solver_id = LCG_CGS);
其中:
lcg_axfunc_ptr Afp
为正演计算的回调函数;lcg_progress_ptr Pfp
监控迭代过程的回调函数(非必须,无需监控时使用 nullptr 参数即可);lcg_float* m
初始解向量,迭代取得的解也保存与此数组;const lcg_float* B
Ax = B 中的 B 项;const int n_size
解向量的大小;const lcg_para* param
迭代使用的参数,此参数为 nullptr 即使用默认参数;void* instance
传入的实例对象, 此函数在类中使用即为类的 this 指针, 在普通函数中使用时即为 nullptr;int solver_id
求解函数使用的求解方法,具体的方法代号可查看对应的头文件;
一个简单的例子
#include "cmath"
#include "iostream"
#include "lcg/lcg.h"
#define M 100
#define N 80
// 返回两个数组元素之间的最大差值
lcg_float max_diff(const lcg_float *a, const lcg_float *b, int size)
{
lcg_float max = -1;
for (int i = 0; i < size; i++)
{
max = lcg_max(sqrt((a[i] - b[i])*(a[i] - b[i])), max);
}
return max;
}
// 普通二维数组做核矩阵
lcg_float **kernel;
// 中间结果数组
lcg_float *tmp_arr;
// 计算核矩阵乘向量的乘积 lcg_solver的回调函数
void CalAx(void* instance, const lcg_float* x, lcg_float* prod_Ax, const int n_s)
{
// 注意核矩阵实际为 kernel^T * kernel,大小为N*N
lcg_matvec(kernel, x, tmp_arr, M, n_s, MatNormal); // tmp_tar = kernel * x
lcg_matvec(kernel, tmp_arr, prod_Ax, M, n_s, MatTranspose); // prod_Ax = kernel^T * tmp_tar
return;
}
// 定义监控函数 lcg_solver的回调函数
// 这个函数显示当前的迭代次数与收敛值
int Prog(void* instance, const lcg_float* m, const lcg_float converge, const lcg_para* param, const int n_s, const int k)
{
std::clog << "\rIteration-times: " << k << "\tconvergence: " << converge;
return 0;
}
int main(int argc, char const *argv[])
{
// 开辟数组空间
kernel = lcg_malloc(M, N);
tmp_arr = lcg_malloc(M);
// 为核矩阵赋初值
lcg_vecrnd(kernel, -1.0, 1.0, M, N);
// 生成一组理论解
lcg_float *fm = lcg_malloc(N);
lcg_vecrnd(fm, 1.0, 2.0, N);
// 计算共轭梯度B项
lcg_float *B = lcg_malloc(N);
lcg_matvec(kernel, fm, tmp_arr, M, N, MatNormal);
lcg_matvec(kernel, tmp_arr, B, M, N, MatTranspose);
// 设置共轭梯度参数
lcg_para self_para = lcg_default_parameters();
self_para.epsilon = 1e-5;
self_para.abs_diff = 0;
// 声明一组解
lcg_float *m = lcg_malloc(N);
lcg_vecset(m, 0.0, N);
// 使用标准共轭梯度方法(LCG_CG)求解线性方程组
// 将回调函数传递给solver
// 由于回调函数为全局函数,因此instance变量的值为NULL
int ret = lcg_solver(CalAx, Prog, m, B, N, &self_para, NULL, LCG_CG);
std::clog << std::endl; lcg_error_str(ret);
std::clog << "maximal difference: " << max_diff(fm, m, N) << std::endl << std::endl;
// 销毁数组
lcg_free(kernel, M);
lcg_free(tmp_arr);
lcg_free(fm);
lcg_free(B);
lcg_free(m);
return 0;
}
完整的例子储存在sample文件夹内。
类模版
liblcg为不同类型的共轭梯度算法定义了通用的求解类模版,包含了类中函数的指针代理及通用的监控函数实现,用户可直接继承并使用。需要注意的是这些类模版中定义了纯虚的函数接口,用户需要全部实现。其中没用到的定义成空函数就行了。以实数的求解类模版为例,需要实现的接口函数包括:
void AxProduct(const lcg_float* a, lcg_float* b, const int num) = 0
void MxProduct(const lcg_float* a, lcg_float* b, const int num) = 0
其中AxProduct
是Ax的计算函数,MxProduct
是预优过程的计算函数,即M^-1x。