83 lines
3.2 KiB
C++
83 lines
3.2 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 "gkernel_sphere.h"
|
|
|
|
#ifdef GCTL_POTENTIAL_AUTODIFF
|
|
|
|
autodiff::var gobser_mass_depth_vz_sig(const gctl::point3dc &obsp, double loc_x,
|
|
double loc_y, const autodiff::ArrayXvar& depth_mass);
|
|
|
|
void gctl::gobser_wrt_mass_depth(double &depth_der, double &mass_der, const point3dc &mass_loc,
|
|
double mass, loss_func &lf, const array<point3dc> &obsp)
|
|
{
|
|
double tmp_val, tmp_grad;
|
|
autodiff::var tmp_obs;
|
|
autodiff::ArrayXvar depth_mass;
|
|
Eigen::VectorXd ders;
|
|
|
|
depth_der = 0.0;
|
|
mass_der = 0.0;
|
|
|
|
int i;
|
|
int o_size = obsp.size();
|
|
|
|
#pragma omp parallel for private (i, tmp_val, tmp_grad, tmp_obs, depth_mass, ders) schedule(guided)
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
depth_mass.resize(2);
|
|
depth_mass[0] = mass_loc.z;
|
|
depth_mass[1] = mass;
|
|
|
|
ders.setZero();
|
|
tmp_obs = gobser_mass_depth_vz_sig(obsp[i], mass_loc.x, mass_loc.y, depth_mass);
|
|
ders = autodiff::gradient(tmp_obs, depth_mass);
|
|
tmp_val = autodiff::val(tmp_obs);
|
|
tmp_grad = lf.gradient(tmp_val, i);
|
|
|
|
#pragma omp critical
|
|
{
|
|
lf.evaluate(tmp_val, i);
|
|
depth_der += tmp_grad*ders[0];
|
|
mass_der += tmp_grad*ders[1];
|
|
}
|
|
}
|
|
|
|
return;
|
|
}
|
|
|
|
|
|
autodiff::var gobser_mass_depth_vz_sig(const gctl::point3dc &obsp, double loc_x,
|
|
double loc_y, const autodiff::ArrayXvar& depth_mass)
|
|
{
|
|
autodiff::var tmp_obs = 1e+8*GCTL_G0*depth_mass[1]*(obsp.z - depth_mass[0])
|
|
/pow(pow(obsp.x - loc_x, 2) + pow(obsp.y - loc_y, 2) + pow(obsp.z - depth_mass[0], 2), 1.5);
|
|
|
|
return tmp_obs;
|
|
}
|
|
|
|
#endif // GCTL_POTENTIAL_AUTODIFF
|