383 lines
12 KiB
C++
383 lines
12 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 "mkernel_tricone.h"
|
|
|
|
void gctl::callink_magnetic_para(array<mag_tricone> &in_cone, array<magcone_para> &out_para, const array<point3dc> &mag_B)
|
|
{
|
|
point3dc v1, v2, v3, ne, nf, mag_z;
|
|
|
|
out_para.resize(in_cone.size());
|
|
for (int i = 0; i < in_cone.size(); ++i)
|
|
{
|
|
if (in_cone[i].vert[0] == nullptr || in_cone[i].vert[1] == nullptr ||
|
|
in_cone[i].vert[2] == nullptr || in_cone[i].vert[3] == nullptr)
|
|
{
|
|
throw runtime_error("Invalid vertex pointer. From callink_magnetic_para(...)");
|
|
}
|
|
|
|
// magnetic susceptibility is taken as one here
|
|
out_para[i].B = mag_B[i];
|
|
|
|
for (int f = 0; f < 4; ++f)
|
|
{
|
|
v1 = *in_cone[i].fget(f, 1) - *in_cone[i].fget(f, 0);
|
|
v2 = *in_cone[i].fget(f, 2) - *in_cone[i].fget(f, 0);
|
|
nf = cross(v1, v2).normal();
|
|
out_para[i].F[f] = kron(nf, nf);
|
|
|
|
for (int e = 0; e < 3; ++e)
|
|
{
|
|
v3 = *in_cone[i].fget(f, (e+1)%3) - *in_cone[i].fget(f, e);
|
|
out_para[i].edglen[e+f*3] = v3.module();
|
|
out_para[i].E[e+f*3] = kron(nf, cross(v3, nf).normal());
|
|
}
|
|
}
|
|
|
|
in_cone[i].att = out_para.get(i);
|
|
}
|
|
return;
|
|
}
|
|
|
|
void gctl::callink_magnetic_para_earth_sph(array<mag_tricone> &in_tet, array<magcone_para> &out_para,
|
|
double inclina_deg, double declina_deg, array<point3dc> *mag_vec, double field_tense)
|
|
{
|
|
if (declina_deg < -180.0 || declina_deg > 180.0 ||
|
|
inclina_deg < -90 || inclina_deg > 90 || field_tense < 0.0)
|
|
{
|
|
throw invalid_argument("Invalid parameters. From gctl::callink_magnetic_para_earth_sph(...)");
|
|
}
|
|
|
|
point3dc v1, v2, v3, ne, nf, mag_z, c;
|
|
point3ds s;
|
|
|
|
double I = inclina_deg*GCTL_Pi/180.0;
|
|
double A = declina_deg*GCTL_Pi/180.0;
|
|
|
|
if (mag_vec != nullptr) mag_vec->resize(in_tet.size());
|
|
out_para.resize(in_tet.size());
|
|
for (int i = 0; i < in_tet.size(); ++i)
|
|
{
|
|
if (in_tet[i].vert[0] == nullptr || in_tet[i].vert[1] == nullptr ||
|
|
in_tet[i].vert[2] == nullptr || in_tet[i].vert[3] == nullptr)
|
|
{
|
|
throw runtime_error("Invalid vertex pointer. From callink_magnetic_para_earth_sph(...)");
|
|
}
|
|
|
|
// magnetization vector at the local coordinates
|
|
// note here the postive direction of the inclination angle is downward
|
|
// note here the postive direction of the declination angle is clockwise
|
|
mag_z.x = cos(I)*sin(A);
|
|
mag_z.y = cos(I)*cos(A);
|
|
mag_z.z = -1.0*sin(I);
|
|
|
|
c = 1.0/3.0*(*in_tet[i].vert[0] + *in_tet[i].vert[1] + *in_tet[i].vert[2]);
|
|
s = c.c2s();
|
|
// magnetic susceptibility is taken as one here
|
|
// rotate the local coordinate system to the regular status
|
|
out_para[i].B = field_tense * mag_z.rotate((90.0 - s.lat)*GCTL_Pi/180.0, 0.0, (90.0 + s.lon)*GCTL_Pi/180.0).normal();
|
|
if (mag_vec != nullptr) mag_vec->at(i) = out_para[i].B;
|
|
|
|
for (int f = 0; f < 4; ++f)
|
|
{
|
|
v1 = *in_tet[i].fget(f, 1) - *in_tet[i].fget(f, 0);
|
|
v2 = *in_tet[i].fget(f, 2) - *in_tet[i].fget(f, 0);
|
|
nf = cross(v1, v2).normal();
|
|
out_para[i].F[f] = kron(nf, nf);
|
|
|
|
for (int e = 0; e < 3; ++e)
|
|
{
|
|
v3 = *in_tet[i].fget(f, (e+1)%3) - *in_tet[i].fget(f, e);
|
|
out_para[i].edglen[e+f*3] = v3.module();
|
|
out_para[i].E[e+f*3] = kron(nf, cross(v3, nf).normal());
|
|
}
|
|
}
|
|
|
|
in_tet[i].att = out_para.get(i);
|
|
}
|
|
return;
|
|
}
|
|
|
|
gctl::point3dc gctl::magkernel_single(const mag_tricone &a_ele, const point3ds &a_op, tensor *R_ptr)
|
|
{
|
|
int f,e;
|
|
double Le,wf;
|
|
double dv1,dv2;
|
|
gctl::point3dc r_ijk[3];
|
|
gctl::tensor R;
|
|
gctl::tensor face_sum(0.0);
|
|
gctl::tensor edge_sum(0.0);
|
|
double L_ijk[3];
|
|
|
|
if (R_ptr != nullptr) R = *R_ptr;
|
|
else R = transform_matrix(a_op);
|
|
|
|
gctl::magcone_para* gp = a_ele.att;
|
|
gctl::point3dc pc = a_op.s2c();
|
|
|
|
for (f = 0; f < 4; f++)
|
|
{
|
|
r_ijk[0] = *a_ele.fget(f, 0) - pc;
|
|
r_ijk[1] = *a_ele.fget(f, 1) - pc;
|
|
r_ijk[2] = *a_ele.fget(f, 2) - pc;
|
|
|
|
L_ijk[0] = r_ijk[0].module();
|
|
L_ijk[1] = r_ijk[1].module();
|
|
L_ijk[2] = r_ijk[2].module();
|
|
|
|
wf = 2*atan2(dot(r_ijk[0], cross(r_ijk[1],r_ijk[2])),
|
|
L_ijk[0]*L_ijk[1]*L_ijk[2] + L_ijk[0]*dot(r_ijk[1],r_ijk[2]) +
|
|
L_ijk[1]*dot(r_ijk[2],r_ijk[0]) + L_ijk[2]*dot(r_ijk[0],r_ijk[1]));
|
|
|
|
face_sum = face_sum + wf * gp->F[f];
|
|
|
|
for (e = 0; e < 3; e++)
|
|
{
|
|
dv1 = distance(*a_ele.fget(f, e), pc);
|
|
dv2 = distance(*a_ele.fget(f, (e+1)%3), pc);
|
|
|
|
Le = log((dv1+dv2+gp->edglen[e+3*f])/(dv1+dv2-gp->edglen[e+3*f]));
|
|
|
|
edge_sum = edge_sum + Le * gp->E[e+3*f];
|
|
}
|
|
}
|
|
|
|
// nT -> T 1e-9 / 4*pi*e-7 = 1/(400*pi)*100 = 1/(4*pi) -> A/m
|
|
gctl::tensor gt = -1.0/(4.0*GCTL_Pi)*(face_sum - edge_sum);
|
|
gctl::point3dc out_c = R*(gt*gp->B);
|
|
out_c.x *= -1.0; // x对应半径向分量
|
|
out_c.y *= -1.0; // y对应纬向分量
|
|
// z对应经向分量
|
|
return out_c;
|
|
}
|
|
|
|
void gctl::magkernel(matrix<double> &kernel, const array<mag_tricone> &top_ele,
|
|
const array<mag_tricone> &btm_ele, const array<point3ds> &obsp,
|
|
magnetic_field_type_e comp_type, verbose_type_e verbose)
|
|
{
|
|
if (comp_type != Bx && comp_type != By && comp_type != Bz)
|
|
{
|
|
throw std::runtime_error("[gctl::magkernel] Invalid magnetic componment type. Must be Bx, By or Bz.");
|
|
return;
|
|
}
|
|
|
|
if (top_ele.size() != btm_ele.size())
|
|
{
|
|
throw std::runtime_error("[gctl::magkernel] Elements' size don't equal.");
|
|
return;
|
|
}
|
|
|
|
int i, j;
|
|
int o_size = obsp.size();
|
|
int e_size = top_ele.size();
|
|
kernel.resize(o_size, e_size);
|
|
|
|
array<tensor> Rs(o_size);
|
|
#pragma omp parallel for private (i) schedule(guided)
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
Rs[i] = transform_matrix(obsp[i]);
|
|
}
|
|
|
|
point3dc mag_b;
|
|
gctl::progress_bar bar(o_size*2, "magkernel_componments");
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar.progressed(i);
|
|
else if (verbose == gctl::ShortMsg) bar.progressed_simple(i);
|
|
|
|
#pragma omp parallel for private (j, mag_b) shared(i) schedule(guided)
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
mag_b = magkernel_single(top_ele[j], obsp[i], Rs.get(i));
|
|
if (comp_type == Bx) kernel[i][j] = mag_b.y;
|
|
if (comp_type == By) kernel[i][j] = mag_b.z;
|
|
if (comp_type == Bz) kernel[i][j] = mag_b.x;
|
|
}
|
|
}
|
|
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar.progressed(i + o_size);
|
|
else if (verbose == gctl::ShortMsg) bar.progressed_simple(i + o_size);
|
|
|
|
#pragma omp parallel for private (j, mag_b) shared(i) schedule(guided)
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
mag_b = magkernel_single(btm_ele[j], obsp[i], Rs.get(i));
|
|
if (comp_type == Bx) kernel[i][j] -= mag_b.y;
|
|
if (comp_type == By) kernel[i][j] -= mag_b.z;
|
|
if (comp_type == Bz) kernel[i][j] -= mag_b.x;
|
|
}
|
|
}
|
|
return;
|
|
}
|
|
|
|
void gctl::magkernel(spmat<double> &kernel, const array<mag_tricone> &top_ele, const array<mag_tricone> &btm_ele,
|
|
const array<point3ds> &obsp, double cut_angle, magnetic_field_type_e comp_type, verbose_type_e verbose)
|
|
{
|
|
if (comp_type != Bx && comp_type != By && comp_type != Bz)
|
|
{
|
|
throw std::runtime_error("[gctl::magkernel] Invalid magnetic componment type. Must be Bx, By or Bz.");
|
|
return;
|
|
}
|
|
|
|
if (top_ele.size() != btm_ele.size())
|
|
{
|
|
throw std::runtime_error("[gctl::magkernel] Elements' size don't equal.");
|
|
return;
|
|
}
|
|
|
|
int i, j;
|
|
int o_size = obsp.size();
|
|
int e_size = top_ele.size();
|
|
|
|
array<tensor> Rs(o_size);
|
|
#pragma omp parallel for private (i) schedule(guided)
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
Rs[i] = transform_matrix(obsp[i]);
|
|
}
|
|
|
|
point3dc cen, mag_b;
|
|
mat_node<double> tmp_plt;
|
|
std::map<int, int> tri_idx;
|
|
std::vector<mat_node<double>> triplts;
|
|
|
|
gctl::progress_bar bar2(o_size, "initializing triplts");
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar2.progressed(i);
|
|
else if (verbose == gctl::ShortMsg) bar2.progressed_simple(i);
|
|
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
cen = 1.0/3.0*(*top_ele[j].vert[0] + *top_ele[j].vert[1] + *top_ele[j].vert[2]);
|
|
if (geometry3d::angle(obsp[i].s2c(), cen) < cut_angle*GCTL_Pi/180.0)
|
|
{
|
|
tmp_plt.r_id = i;
|
|
tmp_plt.c_id = j;
|
|
tmp_plt.val = 1.0;
|
|
|
|
tri_idx[j + i*e_size] = triplts.size();
|
|
triplts.push_back(tmp_plt);
|
|
}
|
|
}
|
|
}
|
|
|
|
gctl::progress_bar bar(o_size*2, "magkernel_componments");
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar.progressed(i);
|
|
else if (verbose == gctl::ShortMsg) bar.progressed_simple(i);
|
|
|
|
#pragma omp parallel for private (j, mag_b, cen) shared(i, cut_angle) schedule(guided)
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
cen = 1.0/3.0*(*top_ele[j].vert[0] + *top_ele[j].vert[1] + *top_ele[j].vert[2]);
|
|
if (geometry3d::angle(obsp[i].s2c(), cen) < cut_angle*GCTL_Pi/180.0)
|
|
{
|
|
mag_b = magkernel_single(top_ele[j], obsp[i], Rs.get(i));
|
|
if (comp_type == Bx) triplts[tri_idx[j + i*e_size]].val = mag_b.y;
|
|
if (comp_type == By) triplts[tri_idx[j + i*e_size]].val = mag_b.z;
|
|
if (comp_type == Bz) triplts[tri_idx[j + i*e_size]].val = mag_b.x;
|
|
}
|
|
}
|
|
}
|
|
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar.progressed(i + o_size);
|
|
else if (verbose == gctl::ShortMsg) bar.progressed_simple(i + o_size);
|
|
|
|
#pragma omp parallel for private (j, mag_b, cen) shared(i, cut_angle) schedule(guided)
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
cen = 1.0/3.0*(*btm_ele[j].vert[0] + *btm_ele[j].vert[1] + *btm_ele[j].vert[2]);
|
|
if (geometry3d::angle(obsp[i].s2c(), cen) < cut_angle*GCTL_Pi/180.0)
|
|
{
|
|
mag_b = magkernel_single(btm_ele[j], obsp[i], Rs.get(i));
|
|
if (comp_type == Bx) triplts[tri_idx[j + i*e_size]].val -= mag_b.y;
|
|
if (comp_type == By) triplts[tri_idx[j + i*e_size]].val -= mag_b.z;
|
|
if (comp_type == Bz) triplts[tri_idx[j + i*e_size]].val -= mag_b.x;
|
|
}
|
|
}
|
|
}
|
|
|
|
kernel.malloc(o_size, e_size, 0.0);
|
|
kernel.set_triplts(triplts);
|
|
return;
|
|
}
|
|
|
|
void gctl::magobser(array<point3dc> &out_obs, const array<mag_tricone> &top_ele,
|
|
const array<mag_tricone> &btm_ele, const array<point3ds> &obsp, const array<double> &sus,
|
|
verbose_type_e verbose)
|
|
{
|
|
if (top_ele.size() != btm_ele.size())
|
|
{
|
|
throw std::runtime_error("[gctl::magobser] Elements' size don't equal.");
|
|
return;
|
|
}
|
|
|
|
int i, j;
|
|
int o_size = obsp.size();
|
|
int e_size = top_ele.size();
|
|
|
|
out_obs.resize(o_size, point3dc(0.0, 0.0, 0.0));
|
|
array<tensor> Rs(o_size);
|
|
#pragma omp parallel for private (i) schedule(guided)
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
Rs[i] = transform_matrix(obsp[i]);
|
|
}
|
|
|
|
gctl::progress_bar bar(e_size*2, "magobser_componments");
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar.progressed(j);
|
|
else if (verbose == gctl::ShortMsg) bar.progressed_simple(j);
|
|
|
|
#pragma omp parallel for private (i) shared (j) schedule(guided)
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
out_obs[i] = out_obs[i] + sus[j] * magkernel_single(top_ele[j], obsp[i], Rs.get(i));
|
|
}
|
|
}
|
|
|
|
for (j = 0; j < e_size; j++)
|
|
{
|
|
if (verbose == gctl::FullMsg) bar.progressed(j + e_size);
|
|
else if (verbose == gctl::ShortMsg) bar.progressed_simple(j + e_size);
|
|
|
|
#pragma omp parallel for private (i) shared (j) schedule(guided)
|
|
for (i = 0; i < o_size; i++)
|
|
{
|
|
out_obs[i] = out_obs[i] - sus[j] * magkernel_single(btm_ele[j], obsp[i], Rs.get(i));
|
|
}
|
|
}
|
|
return;
|
|
} |