gctl_potential/lib/potential/gkernel_sphere.cpp
2024-09-10 19:56:41 +08:00

247 lines
9.0 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"
#include "cmath"
typedef void (*gobser_mass_ptr)(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype);
void gobser_mass_vz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype);
void gobser_mass_vzx(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype);
void gobser_mass_vzy(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype);
void gobser_mass_vzz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype);
void gctl::gobser_mass(array<double> &out_obs, const array<point3dc> &obsp, const point3dc &loc, double mass,
gravitational_field_type_e comp_id, value_operator_e valtype)
{
std::string err_str;
if (obsp.empty())
{
err_str = "Empty observation points. From gctl::gobser(...)";
throw runtime_error(err_str);
}
if((valtype == ReplaceVal || valtype == EraseVal) && out_obs.size() != obsp.size())
{
out_obs.resize(obsp.size());
}
else if (out_obs.size() != obsp.size())
{
err_str = "Arrays' size do not match. From gctl::gobser(...)";
throw runtime_error(err_str);
}
gobser_mass_ptr mass_obser;
switch (comp_id)
{
case Vz:
mass_obser = gobser_mass_vz;
break;
case Tzx:
mass_obser = gobser_mass_vzx;
break;
case Tzy:
mass_obser = gobser_mass_vzy;
break;
case Tzz:
mass_obser = gobser_mass_vzz;
break;
default:
mass_obser = gobser_mass_vz;
break;
}
mass_obser(out_obs, obsp, loc, mass, valtype);
return;
}
typedef void (*gobser_sphere_ptr)(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype);
void gobser_sphere_vz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype);
void gobser_sphere_vzx(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype);
void gobser_sphere_vzy(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype);
void gobser_sphere_vzz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype);
void gctl::gobser_sphere(array<double> &out_obs, const array<point3dc> &obsp, const gctl::sphere &gs,
gravitational_field_type_e comp_id, value_operator_e valtype)
{
std::string err_str;
if (obsp.empty())
{
err_str = "Empty observation points. From gctl::gobser(...)";
throw runtime_error(err_str);
}
if((valtype == ReplaceVal || valtype == EraseVal) && out_obs.size() != obsp.size())
{
out_obs.resize(obsp.size());
}
else if (out_obs.size() != obsp.size())
{
err_str = "Arrays' size do not match. From gctl::gobser(...)";
throw runtime_error(err_str);
}
gobser_sphere_ptr sphere_obser;
switch (comp_id)
{
case Vz:
sphere_obser = gobser_sphere_vz;
break;
case Tzx:
sphere_obser = gobser_sphere_vzx;
break;
case Tzy:
sphere_obser = gobser_sphere_vzy;
break;
case Tzz:
sphere_obser = gobser_sphere_vzz;
break;
default:
sphere_obser = gobser_sphere_vz;
break;
}
sphere_obser(out_obs, obsp, gs, valtype);
return;
}
void gobser_sphere_vz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype)
{
double *rho = (double*)gs.att;
double se_mass = 4.0*(*rho)*gs.rad*gs.rad*gs.rad*GCTL_Pi/3.0;
gctl::point3dc cen(gs.vert[0]->x, gs.vert[0]->y, gs.vert[0]->z);
gobser_mass_vz(out_obs, obsp, cen, se_mass, valtype);
return;
}
void gobser_sphere_vzx(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype)
{
double *rho = (double*)gs.att;
double se_mass = 4.0*(*rho)*gs.rad*gs.rad*gs.rad*GCTL_Pi/3.0;
gctl::point3dc cen(gs.vert[0]->x, gs.vert[0]->y, gs.vert[0]->z);
gobser_mass_vzx(out_obs, obsp, cen, se_mass, valtype);
return;
}
void gobser_sphere_vzy(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype)
{
double *rho = (double*)gs.att;
double se_mass = 4.0*(*rho)*gs.rad*gs.rad*gs.rad*GCTL_Pi/3.0;
gctl::point3dc cen(gs.vert[0]->x, gs.vert[0]->y, gs.vert[0]->z);
gobser_mass_vzy(out_obs, obsp, cen, se_mass, valtype);
return;
}
void gobser_sphere_vzz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::sphere &gs, gctl::value_operator_e valtype)
{
double *rho = (double*)gs.att;
double se_mass = 4.0*(*rho)*gs.rad*gs.rad*gs.rad*GCTL_Pi/3.0;
gctl::point3dc cen(gs.vert[0]->x, gs.vert[0]->y, gs.vert[0]->z);
gobser_mass_vzz(out_obs, obsp, cen, se_mass, valtype);
return;
}
void gobser_mass_vz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype)
{
double res;
for (int i = 0; i < obsp.size(); i++)
{
res = 1e+8*GCTL_G0*mass*(obsp[i].z-loc.z)
/pow(pow(obsp[i].x-loc.x,2)+pow(obsp[i].y-loc.y,2)+pow(obsp[i].z-loc.z,2),1.5);
if (valtype == gctl::ReplaceVal) out_obs[i] = res;
else if (valtype == gctl::AppendVal) out_obs[i] += res;
else if (valtype == gctl::RemoveVal) out_obs[i] -= res;
else out_obs[i] = 0.0;
}
return;
}
void gobser_mass_vzx(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype)
{
double res;
for (int i = 0; i < obsp.size(); i++)
{
res = -3e+8*GCTL_G0*mass*(obsp[i].z-loc.z)*(obsp[i].x-loc.x)
/pow(pow(obsp[i].x-loc.x,2)+pow(obsp[i].y-loc.y,2)+pow(obsp[i].z-loc.z,2),2.5);
if (valtype == gctl::ReplaceVal) out_obs[i] = res;
else if (valtype == gctl::AppendVal) out_obs[i] += res;
else if (valtype == gctl::RemoveVal) out_obs[i] -= res;
else out_obs[i] = 0.0;
}
return;
}
void gobser_mass_vzy(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype)
{
double res;
for (int i = 0; i < obsp.size(); i++)
{
res = -3e+8*GCTL_G0*mass*(obsp[i].z-loc.z)*(obsp[i].y-loc.y)
/pow(pow(obsp[i].x-loc.x,2)+pow(obsp[i].y-loc.y,2)+pow(obsp[i].z-loc.z,2),2.5);
if (valtype == gctl::ReplaceVal) out_obs[i] = res;
else if (valtype == gctl::AppendVal) out_obs[i] += res;
else if (valtype == gctl::RemoveVal) out_obs[i] -= res;
else out_obs[i] = 0.0;
}
return;
}
void gobser_mass_vzz(gctl::array<double> &out_obs, const gctl::array<gctl::point3dc> &obsp,
const gctl::point3dc &loc, double mass, gctl::value_operator_e valtype)
{
double res;
for (int i = 0; i < obsp.size(); i++)
{
res = 1e+8*GCTL_G0*mass*(2.0*pow(obsp[i].z-loc.z,2)-pow(obsp[i].x-loc.x,2)-pow(obsp[i].y-loc.y,2))
/pow(pow(obsp[i].x-loc.x,2)+pow(obsp[i].y-loc.y,2)+pow(obsp[i].z-loc.z,2),2.5);
if (valtype == gctl::ReplaceVal) out_obs[i] = res;
else if (valtype == gctl::AppendVal) out_obs[i] += res;
else if (valtype == gctl::RemoveVal) out_obs[i] -= res;
else out_obs[i] = 0.0;
}
return;
}