2024-09-10 20:04:47 +08:00
|
|
|
/********************************************************
|
|
|
|
* ██████╗ ██████╗████████╗██╗
|
|
|
|
* ██╔════╝ ██╔════╝╚══██╔══╝██║
|
|
|
|
* ██║ ███╗██║ ██║ ██║
|
|
|
|
* ██║ ██║██║ ██║ ██║
|
|
|
|
* ╚██████╔╝╚██████╗ ██║ ███████╗
|
|
|
|
* ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝
|
|
|
|
* Geophysical Computational Tools & Library (GCTL)
|
|
|
|
*
|
|
|
|
* Copyright (c) 2022 Yi Zhang (yizhang-geo@zju.edu.cn)
|
|
|
|
*
|
|
|
|
* GCTL 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 along with this
|
|
|
|
* program. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*
|
|
|
|
* If the terms and conditions of the LGPL v.2. would prevent you from using
|
|
|
|
* the GCTL, please consider the option to obtain a commercial license for a
|
|
|
|
* fee. These licenses are offered by the GCTL'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 "gctl/core.h"
|
|
|
|
#include "gctl/algorithm.h"
|
|
|
|
#include "../lib/optimization.h"
|
|
|
|
|
|
|
|
typedef gctl::array<std::complex<double>> cd_array;
|
|
|
|
|
|
|
|
#define N 1000
|
|
|
|
|
|
|
|
double max_diff(const cd_array &a, const cd_array &b)
|
|
|
|
{
|
|
|
|
double max = -1;
|
|
|
|
std::complex<double> t;
|
|
|
|
for (size_t i = 0; i < a.size(); i++)
|
|
|
|
{
|
|
|
|
t = a[i] - b[i];
|
|
|
|
max = std::max(std::norm(t), max);
|
|
|
|
}
|
|
|
|
return max;
|
|
|
|
}
|
|
|
|
|
|
|
|
class ex7 : public gctl::clcg_solver
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
ex7();
|
|
|
|
virtual ~ex7();
|
|
|
|
virtual void CLCG_Ax(const cd_array &x, cd_array &ax, gctl::matrix_layout_e layout, gctl::conjugate_type_e conj);
|
|
|
|
|
|
|
|
// 计算共轭梯度的B项
|
|
|
|
void cal_partb(const cd_array &x, cd_array &B);
|
|
|
|
|
|
|
|
private:
|
|
|
|
gctl::matrix<std::complex<double>> kernel; // 普通二维数组做核矩阵
|
|
|
|
};
|
|
|
|
|
|
|
|
ex7::ex7()
|
|
|
|
{
|
|
|
|
gctl::array<double> tmp(round(0.5*(N+1)*N));
|
2024-09-19 22:14:47 +08:00
|
|
|
tmp.random_float(1.0, 2.0, gctl::RdUniform);
|
2024-09-10 20:04:47 +08:00
|
|
|
|
|
|
|
size_t c = 0;
|
|
|
|
kernel.resize(N, N);
|
|
|
|
for (int i = 0; i < N; i++)
|
|
|
|
{
|
|
|
|
for (int j = i; j < N; j++)
|
|
|
|
{
|
|
|
|
kernel[i][j] = tmp[c];
|
|
|
|
kernel[j][i] = kernel[i][j];
|
|
|
|
c++;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
ex7::~ex7(){}
|
|
|
|
|
|
|
|
void ex7::cal_partb(const cd_array &x, cd_array &B)
|
|
|
|
{
|
|
|
|
gctl::matvec(B, kernel, x);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
void ex7::CLCG_Ax(const cd_array &x, cd_array &ax, gctl::matrix_layout_e layout, gctl::conjugate_type_e conj)
|
|
|
|
{
|
|
|
|
gctl::matvec(ax, kernel, x, layout, conj);
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
int main(int argc, char const *argv[])
|
|
|
|
{
|
|
|
|
// 生成一组正演解
|
|
|
|
gctl::array<double> tmp(2*N);
|
2024-09-19 22:14:47 +08:00
|
|
|
tmp.random_float(1.0, 2.0, gctl::RdUniform);
|
2024-09-10 20:04:47 +08:00
|
|
|
|
|
|
|
cd_array fm(N);
|
|
|
|
for (size_t i = 0; i < N; i++)
|
|
|
|
{
|
|
|
|
fm[i].real(tmp[2*i]);
|
|
|
|
fm[i].imag(tmp[2*i + 1]);
|
|
|
|
}
|
|
|
|
|
|
|
|
ex7 test;
|
|
|
|
|
|
|
|
// 计算共轭梯度B项
|
|
|
|
cd_array B(N);
|
|
|
|
test.cal_partb(fm, B);
|
|
|
|
|
|
|
|
// 声明一组解
|
|
|
|
cd_array m(N, std::complex<double>(0.0, 0.0));
|
|
|
|
|
|
|
|
gctl::clcg_para my_para = test.default_clcg_para();
|
|
|
|
my_para.abs_diff = 1;
|
|
|
|
|
|
|
|
test.set_clcg_para(my_para);
|
|
|
|
|
|
|
|
test.CLCG_Minimize(m, B, gctl::CLCG_BICG_SYM);
|
|
|
|
std::clog << "maximal difference: " << max_diff(fm, m) << std::endl;
|
|
|
|
return 0;
|
|
|
|
}
|