/********************************************************
* ██████╗ ███████╗████████╗██╗
* ██╔════╝ ██╔════╝╚══██╔══╝██║
* ██║ ███╗███████╗ ██║ ██║
* ██║ ██║╚════██║ ██║ ██║
* ╚██████╔╝███████║ ██║ ███████╗
* ╚═════╝ ╚══════╝ ╚═╝ ╚══════╝
* Generic Scientific Template Library
*
* Copyright (c) 2022 Yi Zhang (yizhang-geo@zju.edu.cn)
*
* The GSTL is distributed under a dual licensing scheme. You can redistribute
* it and/or modify it under the terms of the GNU Lesser General Public
* License as published by the Free Software Foundation, either version 2
* of the License, or (at your option) any later version. You should have
* received a copy of the GNU Lesser General Public License (LGPL) along with
* this program. If not, see .
*
* If the terms and conditions of the LGPL v.2. would prevent you from using
* the GSTL, please consider the option to obtain a commercial license for a
* fee. These licenses are offered by the GSTL's original author. As a rule,
* licenses are provided "as-is", unlimited in time for a one time fee. Please
* send corresponding requests to: yizhang-geo@zju.edu.cn. Please do not forget
* to include some description of your company and the realm of its activities.
* Also add information on how to contact you by electronic and paper mail.
******************************************************/
#include "../lib/optimization.h"
#define M 1000
#define N 800
double max_diff(const gctl::array &a, const gctl::array &b)
{
double max = -1.0;
for (size_t i = 0; i < a.size(); i++)
{
max = std::max(sqrt((a[i] - b[i])*(a[i] - b[i])), max);
}
return max;
}
class ex1 : public gctl::lcg_solver
{
public:
ex1();
virtual ~ex1();
// 计算共轭梯度的B项
void cal_partb(const gctl::array &x, gctl::array &B);
//定义共轭梯度中Ax的算法
virtual void LCG_Ax(const gctl::array &x, gctl::array &ax);
virtual void LCG_Mx(const gctl::array &x, gctl::array &mx);
private:
gctl::matrix kernel; // 普通二维数组做核矩阵
gctl::array tmp_arr; // 中间结果数组
gctl::array p; // 预优矩阵
};
ex1::ex1()
{
kernel.resize(M, N);
kernel.random(-1.0, 1.0, gctl::RdUniform);
tmp_arr.resize(M);
p.resize(N);
for (size_t i = 0; i < N; i++)
{
p[i] = 1.0;
}
double diag;
for (size_t i = 0; i < N; i++)
{
diag = 0.0;
for (size_t j = 0; j < M; j++)
{
diag += kernel[j][i]*kernel[j][i];
}
p[i] = 1.0/diag;
}
}
ex1::~ex1(){}
void ex1::cal_partb(const gctl::array &x, gctl::array &B)
{
LCG_Ax(x, B);
return;
}
void ex1::LCG_Ax(const gctl::array &x, gctl::array &ax)
{
matvec(tmp_arr, kernel, x);
matvec(ax, kernel, tmp_arr, gctl::Trans);
return;
}
void ex1::LCG_Mx(const gctl::array &x, gctl::array &mx)
{
vecmul(mx, p, x);
return;
}
int main(int argc, char const *argv[])
{
// 生成一组正演解
gctl::array fm(N);
gctl::random(fm, 1.0, 2.0, gctl::RdUniform);
ex1 test;
// 计算共轭梯度B项
gctl::array B(N);
test.cal_partb(fm, B);
// 声明一组解
gctl::array m(N, 0.0);
test.set_lcg_message(gctl::LCG_SOLUTION);
std::ofstream ofile("log.txt");
test.LCG_Minimize(m, B, gctl::LCG_CG, ofile);
ofile << "maximal difference: " << max_diff(fm, m) << std::endl;
m.assign_all(0.0);
test.LCG_Minimize(m, B, gctl::LCG_PCG, ofile);
ofile << "maximal difference: " << max_diff(fm, m) << std::endl;
m.assign_all(0.0);
test.LCG_Minimize(m, B, gctl::LCG_CGS, ofile);
ofile << "maximal difference: " << max_diff(fm, m) << std::endl;
ofile.close();
test.set_lcg_message(gctl::LCG_SOLUTION);
m.assign_all(0.0);
test.LCG_Minimize(m, B, gctl::LCG_BICGSTAB);
std::clog << "maximal difference: " << max_diff(fm, m) << std::endl;
m.assign_all(0.0);
test.LCG_Minimize(m, B, gctl::LCG_BICGSTAB2);
std::clog << "maximal difference: " << max_diff(fm, m) << std::endl;
gctl::array low(N, 1.0);
gctl::array hig(N, 2.0);
m.assign_all(0.0);
test.LCG_MinimizeConstrained(m, B, low, hig, gctl::LCG_PG);
std::clog << "maximal difference: " << max_diff(fm, m) << std::endl;
m.assign_all(0.0);
test.LCG_MinimizeConstrained(m, B, low, hig, gctl::LCG_SPG);
std::clog << "maximal difference: " << max_diff(fm, m) << std::endl;
return 0;
}