2025-01-08 19:27:51 +08:00
|
|
|
/********************************************************
|
|
|
|
* ██████╗ ██████╗████████╗██╗
|
|
|
|
* ██╔════╝ ██╔════╝╚══██╔══╝██║
|
|
|
|
* ██║ ███╗██║ ██║ ██║
|
|
|
|
* ██║ ██║██║ ██║ ██║
|
|
|
|
* ╚██████╔╝╚██████╗ ██║ ███████╗
|
|
|
|
* ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝
|
|
|
|
* 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_Ren2017.h"
|
|
|
|
#include "cmath"
|
|
|
|
|
2025-01-09 00:08:55 +08:00
|
|
|
double mag_tolerance = 1e-16;
|
|
|
|
|
|
|
|
void gctl::set_magcone_ren17_tolerance(double tol)
|
|
|
|
{
|
|
|
|
if (tol > 0) mag_tolerance = tol;
|
|
|
|
else throw invalid_argument("Invalid tolerance value. From gctl::set_magcone_ren17_tolerance(...)");
|
|
|
|
return;
|
|
|
|
}
|
2025-01-08 19:27:51 +08:00
|
|
|
|
|
|
|
void gctl::callink_magnetic_para(array<magcone_ren17> &in_cone, array<magcone_para_ren17> &out_para, const array<point3dc> &mag_B)
|
|
|
|
{
|
|
|
|
point3dc v1, v2, v3, ne, nf;
|
|
|
|
|
|
|
|
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(...)");
|
|
|
|
}
|
|
|
|
|
|
|
|
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].nf[f] = nf;
|
|
|
|
|
|
|
|
for (int e = 0; e < 3; ++e)
|
|
|
|
{
|
|
|
|
v3 = *in_cone[i].fget(f, (e+1)%3) - *in_cone[i].fget(f, e);
|
|
|
|
ne = cross(v3, nf).normal();
|
|
|
|
out_para[i].ne[e+f*3] = ne;
|
|
|
|
out_para[i].te[e+f*3] = cross(nf, ne).normal();
|
|
|
|
}
|
|
|
|
|
|
|
|
out_para[i].mag_amp[f] = dot(mag_B[i], nf);
|
|
|
|
}
|
|
|
|
|
|
|
|
in_cone[i].att = out_para.get(i);
|
|
|
|
}
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
void gctl::callink_magnetic_para_earth_sph(array<magcone_ren17> &in_cone, array<magcone_para_ren17> &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*M_PI/180.0;
|
|
|
|
double A = declina_deg*M_PI/180.0;
|
|
|
|
|
|
|
|
if (mag_vec != nullptr) mag_vec->resize(in_cone.size());
|
|
|
|
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_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_cone[i].vert[0] + *in_cone[i].vert[1] + *in_cone[i].vert[2]);
|
|
|
|
s = c.c2s();
|
|
|
|
// magnetic susceptibility is taken as one here
|
|
|
|
// rotate the local coordinate system to the regular status
|
|
|
|
mag_z = field_tense * mag_z.rotate((90.0 - s.lat)*M_PI/180.0, 0.0, (90.0 + s.lon)*M_PI/180.0).normal();
|
|
|
|
if (mag_vec != nullptr) mag_vec->at(i) = mag_z;
|
|
|
|
|
|
|
|
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].nf[f] = nf;
|
|
|
|
|
|
|
|
for (int e = 0; e < 3; ++e)
|
|
|
|
{
|
|
|
|
v3 = *in_cone[i].fget(f, (e+1)%3) - *in_cone[i].fget(f, e);
|
|
|
|
ne = cross(v3, nf).normal();
|
|
|
|
out_para[i].ne[e+f*3] = ne;
|
|
|
|
out_para[i].te[e+f*3] = cross(nf, ne).normal();
|
|
|
|
}
|
|
|
|
|
|
|
|
out_para[i].mag_amp[f] = dot(mag_z, nf);
|
|
|
|
}
|
|
|
|
|
|
|
|
in_cone[i].att = out_para.get(i);
|
|
|
|
}
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
gctl::point3dc gctl::magkernel_single(const magcone_ren17 &a_ele, const point3ds &a_op, tensor *R_ptr)
|
|
|
|
{
|
|
|
|
double Rij_minus, Rij_plus, Sij_plus, Sij_minus, Rij0, mij0, wi0;
|
|
|
|
double beta, Aij, sig, absw;
|
|
|
|
gctl::point3dc oi, k1, part1, part2;
|
|
|
|
|
|
|
|
// get attribute pointer
|
|
|
|
if (a_ele.att == nullptr) throw gctl::runtime_error("[gctl_potential::magcone_ren17] Magnetization parameter not set.");
|
|
|
|
gctl::magcone_para_ren17 *mpara = a_ele.att;
|
|
|
|
|
|
|
|
gctl::tensor R;
|
|
|
|
if (R_ptr != nullptr) R = *R_ptr;
|
|
|
|
else R = transform_matrix(a_op);
|
|
|
|
|
|
|
|
// get the observation site in the local coordinates
|
|
|
|
gctl::point3dc site = a_op.s2c();
|
|
|
|
gctl::point3dc out_grad(0.0, 0.0, 0.0);
|
|
|
|
for (int f = 0; f < 4; ++f)
|
|
|
|
{
|
|
|
|
k1.set(0.0, 0.0, 0.0);
|
|
|
|
for (int j = 0; j < 3; ++j)
|
|
|
|
{
|
|
|
|
Rij_minus = (site - *a_ele.fget(f, j)).module();
|
|
|
|
Rij_plus = (site - *a_ele.fget(f, (j+1)%3)).module();
|
|
|
|
|
|
|
|
if (j == 0)
|
|
|
|
{
|
|
|
|
wi0 = gctl::dot(site - *a_ele.fget(f, j), mpara->nf[f]);
|
|
|
|
sig = gctl::sign(wi0);
|
|
|
|
absw = std::abs(wi0);
|
|
|
|
}
|
|
|
|
|
|
|
|
oi = site - wi0*mpara->nf[f];
|
|
|
|
Sij_minus = gctl::dot(*a_ele.fget(f, j) - oi, mpara->te[3*f+j]);
|
|
|
|
Sij_plus = gctl::dot(*a_ele.fget(f, (j+1)%3) - oi, mpara->te[3*f+j]);
|
|
|
|
mij0 = gctl::dot(*a_ele.fget(f, j) - oi, mpara->ne[3*f+j]);
|
|
|
|
Rij0 = std::sqrt(wi0*wi0 + mij0*mij0);
|
|
|
|
|
|
|
|
part2.set(0.0, 0.0, 0.0);
|
2025-01-09 00:08:55 +08:00
|
|
|
if (absw > mag_tolerance)
|
2025-01-08 19:27:51 +08:00
|
|
|
{
|
|
|
|
beta = atan((mij0*Sij_plus)/(Rij0*Rij0 + absw*Rij_plus))
|
|
|
|
- atan((mij0*Sij_minus)/(Rij0*Rij0 + absw*Rij_minus));
|
|
|
|
|
|
|
|
part2 = sig*beta*mpara->nf[f];
|
|
|
|
}
|
|
|
|
|
2025-01-09 00:08:55 +08:00
|
|
|
if (std::abs(Rij0) > mag_tolerance)
|
2025-01-08 19:27:51 +08:00
|
|
|
{
|
|
|
|
Aij = std::log((long double)(Rij_plus+Sij_plus)) - std::log((long double)(Rij_minus+Sij_minus));
|
|
|
|
}
|
|
|
|
else if (Sij_plus > 0.0 && Sij_minus > 0.0)
|
|
|
|
{
|
|
|
|
Aij = std::log(Sij_plus) - std::log(Sij_minus);
|
|
|
|
}
|
|
|
|
else if (Sij_plus < 0.0 && Sij_minus < 0.0)
|
|
|
|
{
|
|
|
|
Aij = std::log(-1.0*Sij_minus) - std::log(-1.0*Sij_plus);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
throw gctl::runtime_error("[gctl_potential::magcone_ren17] Observation site on edge.");
|
|
|
|
}
|
|
|
|
|
|
|
|
part1 = Aij * mpara->ne[3*f+j];
|
|
|
|
k1 = k1 - (part1 + part2);
|
|
|
|
}
|
|
|
|
|
|
|
|
// nT -> T 1e-9 / 4*pi*e-7 = 1/(400*pi)*100 = 1/(4*pi) -> A/m
|
|
|
|
out_grad = out_grad - 1.0/(4.0*GCTL_Pi)*k1*mpara->mag_amp[f];
|
|
|
|
}
|
|
|
|
|
|
|
|
out_grad = R*out_grad;
|
|
|
|
out_grad.x *= -1.0;
|
|
|
|
out_grad.y *= -1.0;
|
|
|
|
return out_grad;
|
|
|
|
}
|
|
|
|
|
|
|
|
gctl::tensor gctl::magkernel_single_tensor(const magcone_ren17 &a_ele, const point3ds &a_op, tensor *R_ptr)
|
|
|
|
{
|
|
|
|
double Rij_minus, Rij_plus, Sij_plus, Sij_minus, Rij0, mij0, wi0;
|
|
|
|
double beta, sig, absw;
|
|
|
|
double factor_n_mij, factor_tij;
|
|
|
|
gctl::point3dc oi, grad_Aij;
|
|
|
|
gctl::point3dc grad_Rij_plus, grad_Rij_minus, grad_Sij_plus, grad_Sij_minus;
|
|
|
|
gctl::point3dc grad_mij0, grad_Rij0, grad_abs_wi0;
|
|
|
|
gctl::point3dc grad_a_plus, grad_b_plus, grad_a_minus, grad_b_minus;
|
|
|
|
gctl::point3dc grad_betaij_plus, grad_betaij_minus, grad_betaij;
|
|
|
|
double a_plus, b_plus, a_minus, b_minus;
|
|
|
|
double k3;
|
|
|
|
gctl::tensor tmp_k, k2;
|
|
|
|
|
|
|
|
// get attribute pointer
|
|
|
|
if (a_ele.att == nullptr) throw gctl::runtime_error("[gctl_potential::magcone_ren17] Magnetization parameter not set.");
|
|
|
|
gctl::magcone_para_ren17 *mpara = a_ele.att;
|
|
|
|
|
|
|
|
gctl::tensor R, R_T;
|
|
|
|
if (R_ptr != nullptr) R = *R_ptr;
|
|
|
|
else R = transform_matrix(a_op);
|
|
|
|
R_T = R.transpose();
|
|
|
|
|
|
|
|
// get the observation site in the local coordinates
|
|
|
|
gctl::point3dc site = a_op.s2c();
|
|
|
|
gctl::tensor out_tensor(0.0);
|
|
|
|
for (int f = 0; f < 4; ++f)
|
|
|
|
{
|
|
|
|
k2.set(0.0);
|
|
|
|
k3 = 0.0;
|
|
|
|
for (int j = 0; j < 3; ++j)
|
|
|
|
{
|
|
|
|
Rij_minus = (site - *a_ele.fget(f, j)).module();
|
|
|
|
Rij_plus = (site - *a_ele.fget(f, (j+1)%3)).module();
|
|
|
|
|
|
|
|
if (j == 0)
|
|
|
|
{
|
|
|
|
wi0 = gctl::dot(site - *a_ele.fget(f, j), mpara->nf[f]);
|
|
|
|
sig = gctl::sign(wi0);
|
|
|
|
absw = std::abs(wi0);
|
|
|
|
}
|
|
|
|
|
|
|
|
oi = site - wi0*mpara->nf[f];
|
|
|
|
Sij_minus = gctl::dot(*a_ele.fget(f, j) - oi, mpara->te[3*f+j]);
|
|
|
|
Sij_plus = gctl::dot(*a_ele.fget(f, (j+1)%3) - oi, mpara->te[3*f+j]);
|
|
|
|
mij0 = gctl::dot(*a_ele.fget(f, j) - oi, mpara->ne[3*f+j]);
|
|
|
|
Rij0 = std::sqrt(wi0*wi0 + mij0*mij0);
|
|
|
|
|
2025-01-09 00:08:55 +08:00
|
|
|
if (std::abs(Rij0) > mag_tolerance)
|
2025-01-08 19:27:51 +08:00
|
|
|
{
|
|
|
|
factor_n_mij = -1.0*(Sij_plus/(Rij0*Rij0*Rij_plus) - Sij_minus/(Rij0*Rij0*Rij_minus));
|
|
|
|
factor_tij = -1.0/Rij_plus + 1.0/Rij_minus;
|
|
|
|
grad_Aij = (wi0*factor_n_mij)*mpara->nf[f] + factor_tij*mpara->te[3*f+j]
|
|
|
|
- (mij0*factor_n_mij)*mpara->ne[3*f+j];
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
factor_tij = -1.0/Rij_plus + 1.0/Rij_minus;
|
|
|
|
grad_Aij = factor_tij*mpara->te[3*f+j];
|
|
|
|
}
|
|
|
|
|
|
|
|
//tmp_k = gctl::kron(grad_Aij, ne_[e][3*f+j]);
|
|
|
|
tmp_k = gctl::kron(mpara->ne[3*f+j], grad_Aij);
|
|
|
|
k2 = k2 - tmp_k;
|
|
|
|
|
2025-01-09 00:08:55 +08:00
|
|
|
if (absw > mag_tolerance)
|
2025-01-08 19:27:51 +08:00
|
|
|
{
|
|
|
|
grad_Rij_plus = (1.0/Rij_plus)*(site - *a_ele.fget(f, (j+1)%3));
|
|
|
|
grad_Rij_minus = (1.0/Rij_minus)*(site - *a_ele.fget(f, j));
|
|
|
|
grad_Sij_plus = -1.0*mpara->te[3*f+j];
|
|
|
|
grad_Sij_minus = -1.0*mpara->te[3*f+j];
|
|
|
|
grad_mij0 = -1.0*mpara->ne[3*f+j];
|
|
|
|
grad_Rij0 = (1.0/Rij0)*(wi0*mpara->nf[f] - mij0*mpara->ne[3*f+j]);
|
|
|
|
grad_abs_wi0 = sig*mpara->nf[f];
|
|
|
|
a_plus = Rij0*Rij0 + absw*Rij_plus;
|
|
|
|
b_plus = mij0*Sij_plus;
|
|
|
|
grad_a_plus = (2.0*Rij0)*grad_Rij0 + Rij_plus*grad_abs_wi0 + absw*grad_Rij_plus;
|
|
|
|
grad_b_plus = Sij_plus*grad_mij0 + mij0*grad_Sij_plus;
|
|
|
|
a_minus = Rij0*Rij0 + absw*Rij_minus;
|
|
|
|
b_minus = mij0*Sij_minus;
|
|
|
|
grad_a_minus = (2.0*Rij0)*grad_Rij0 + Rij_minus*grad_abs_wi0 + absw*grad_Rij_minus;
|
|
|
|
grad_b_minus = Sij_minus*grad_mij0 + mij0*grad_Sij_minus;
|
|
|
|
grad_betaij_plus = (1.0/(a_plus*a_plus + b_plus*b_plus))*(a_plus*grad_b_plus - b_plus*grad_a_plus);
|
|
|
|
grad_betaij_minus = (1.0/(a_minus*a_minus + b_minus*b_minus))*(a_minus*grad_b_minus - b_minus*grad_a_minus);
|
|
|
|
|
|
|
|
grad_betaij = grad_betaij_plus - grad_betaij_minus;
|
|
|
|
|
|
|
|
tmp_k = gctl::kron(grad_betaij, mpara->nf[f]);
|
|
|
|
k2 = k2 - sig*tmp_k;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
2025-01-09 00:08:55 +08:00
|
|
|
if (std::abs(mij0) > mag_tolerance)
|
2025-01-08 19:27:51 +08:00
|
|
|
{
|
|
|
|
k3 += (-1.0/mij0)*(Sij_plus/Rij_plus - Sij_minus/Rij_minus);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
if (k3 != 0.0)
|
|
|
|
{
|
|
|
|
tmp_k = gctl::kron(mpara->nf[f], mpara->nf[f]);
|
|
|
|
k2 = k2 - k3*tmp_k;
|
|
|
|
}
|
|
|
|
|
|
|
|
// nT -> T 1e-9 / 4*pi*e-7 = 1/(400*pi)*100 = 1/(4*pi) -> A/m
|
|
|
|
out_tensor = out_tensor - 1.0/(4.0*GCTL_Pi)*k2*mpara->mag_amp[f];
|
|
|
|
}
|
|
|
|
|
|
|
|
out_tensor = R*out_tensor*R_T;
|
|
|
|
out_tensor[0][0] *= -1.0;
|
|
|
|
out_tensor[0][1] *= -1.0;
|
|
|
|
out_tensor[0][2] *= -1.0;
|
|
|
|
out_tensor[1][0] *= -1.0;
|
|
|
|
out_tensor[2][0] *= -1.0;
|
|
|
|
return out_tensor;
|
|
|
|
}
|
|
|
|
|
|
|
|
void gctl::magkernel(matrix<double> &kernel, const array<magcone_ren17> &top_ele,
|
|
|
|
const array<magcone_ren17> &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_potential::magcone_ren17] Invalid magnetic componment type. Must be Bx, By or Bz.");
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (top_ele.size() != btm_ele.size())
|
|
|
|
{
|
|
|
|
throw std::runtime_error("[gctl_potential::magcone_ren17] 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<magcone_ren17> &top_ele,
|
|
|
|
const array<magcone_ren17> &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_potential::magcone_ren17] Invalid magnetic componment type. Must be Bx, By or Bz.");
|
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (top_ele.size() != btm_ele.size())
|
|
|
|
{
|
|
|
|
throw std::runtime_error("[gctl_potential::magcone_ren17] 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<magcone_ren17> &top_ele,
|
|
|
|
const array<magcone_ren17> &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_potential::magcone_ren17] 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;
|
|
|
|
}
|
|
|
|
|
|
|
|
void gctl::magobser(array<tensor> &out_obs, const array<magcone_ren17> &top_ele,
|
|
|
|
const array<magcone_ren17> &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_potential::magcone_ren17] 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, tensor(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_tensor(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_tensor(btm_ele[j], obsp[i], Rs.get(i));
|
|
|
|
}
|
|
|
|
}
|
|
|
|
return;
|
|
|
|
}
|