gctl_potential/example/mobser_block_gradient_ex.cpp
2024-09-10 19:56:41 +08:00

77 lines
3.1 KiB
C++

/********************************************************
* ██████╗ ██████╗████████╗██╗
* ██╔════╝ ██╔════╝╚══██╔══╝██║
* ██║ ███╗██║ ██║ ██║
* ██║ ██║██║ ██║ ██║
* ╚██████╔╝╚██████╗ ██║ ███████╗
* ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝
* 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 "../lib/potential.h"
using namespace gctl;
int main(int argc, char const *argv[])
{
// 设置块体模型
// 设置块体的磁化率与磁化参数
int m_num = 50, n_num = 50;
array<mag_block> bk(m_num*n_num);
array<double> sus(m_num*n_num, 0.1);
array<magblock_para> mag_para(m_num*n_num);
for (int i = 0; i < m_num; i++)
{
for (int j = 0; j < n_num; j++)
{
if (j >= 20 && j <= 29 && i >= 20 && i <= 29) bk[j + i*n_num].set(j*10.0, (j+1)*10.0, i*10.0, (i+1)*10.0, -10.0, 0.0);
else bk[j + i*n_num].set(j*10.0, (j+1)*10.0, i*10.0, (i+1)*10.0, -1e-6, 0.0);
mag_para[j + i*n_num].inclina_deg = 90.0;
mag_para[j + i*n_num].declina_deg = 0.0;
bk[j + i*n_num].att = mag_para.get(j + i*n_num);
}
}
// 设置观测点位
array<point3dc> obsp;
gridspace(point3dc(5, 0, 5), point3dc(495, 0, 5), point3dc(0, 5, 5), point3dc(0, 495, 5), 50, 50, obsp);
// 正演计算
array<double> za, za_grad;
magobser(za, bk, obsp, sus, Za, FullMsg);
for (int i = 0; i < bk.size(); i++)
{
bk[i].vert[0]->z = -1e-6;
bk[i].vert[1]->z = -1e-6;
bk[i].vert[2]->z = -1e-6;
bk[i].vert[3]->z = -1e-6;
}
loss_func lf(za, L2);
lf.set_uncertainty(1.0);
magobser_wrt_thickness(za_grad, lf, obsp, sus, bk, Za, ShortMsg);
// 保存网格
gctl::save_netcdf_grid("out_magobser_block_gradient", za, 50, 50, 5.0, 10.0, 5.0, 10.0, "x", "y", "Za");
gctl::append_netcdf_grid("out_magobser_block_gradient", za_grad, "x", "y", "Za_Grad");
return 0;
}