/********************************************************
* ██████╗ ██████╗████████╗██╗
* ██╔════╝ ██╔════╝╚══██╔══╝██║
* ██║ ███╗██║ ██║ ██║
* ██║ ██║██║ ██║ ██║
* ╚██████╔╝╚██████╗ ██║ ███████╗
* ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝
* 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 .
*
* 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.
******************************************************/
#ifndef _GCTL_GRAV_KERNEL_SPHERE_H
#define _GCTL_GRAV_KERNEL_SPHERE_H
#include "gctl/core/array.h"
#include "gctl/geometry/point3c.h"
#include "gctl/geometry/sphere.h"
#include "gctl/optimization/loss_func.h"
#include "gm_data.h"
#ifdef GCTL_POTENTIAL_AUTODIFF
#include "autodiff/reverse/var.hpp"
#include "autodiff/reverse/var/eigen.hpp"
#endif // GCTL_POTENTIAL_AUTODIFF
namespace gctl
{
/**
* @brief 计算质量点的重力异常值
*
* @param out_obs 输出的重力数据
* @param[in] obsp 直角坐标下的观测点数组
* @param[in] loc 质量点的位置
* @param[in] mass 质量点的重量
* @param[in] comp_id 待计算的重力场分量 默认为重力值
* @param[in] valtype 当前计算的重力观测值在输出数据中的赋值类型
*/
void gobser_mass(array &out_obs, const array &obsp, const point3dc &loc, double mass,
gravitational_field_type_e comp_id = Vz, value_operator_e valtype = ReplaceVal);
/**
* @brief 计算球体的重力异常值
*
* @param out_obs 输出的重力数据
* @param[in] obsp 直角坐标下的观测点数组
* @param[in] gs 正演的球体
* @param[in] comp_id 待计算的重力场分量 默认为重力值
* @param[in] valtype 当前计算的重力观测值在输出数据中的赋值类型
*/
void gobser_sphere(array &out_obs, const array &obsp, const sphere &gs,
gravitational_field_type_e comp_id = Vz, value_operator_e valtype = ReplaceVal);
#ifdef GCTL_POTENTIAL_AUTODIFF
/**
* @brief 计算球体的正演重力异常与观测数据的数据拟合差,以及其相对于球体质量与埋深的偏导数
*
* @param depth_der 相对于埋深的偏导数
* @param mass_der 相对于质量的偏导数
* @param mass_loc 球体的位置
* @param mass 球体的质量
* @param lf 数据拟合差对象(包含了观测数据与不确定度等信息)
* @param obsp 观测点位
*/
void gobser_wrt_mass_depth(double &depth_der, double &mass_der, const point3dc &mass_loc,
double mass, loss_func &lf, const array &obsp);
#endif // GCTL_POTENTIAL_AUTODIFF
}
#endif // _GCTL_GRAV_KERNEL_SPHERE_H