/******************************************************** * ██████╗ ██████╗████████╗██╗ * ██╔════╝ ██╔════╝╚══██╔══╝██║ * ██║ ███╗██║ ██║ ██║ * ██║ ██║██║ ██║ ██║ * ╚██████╔╝╚██████╗ ██║ ███████╗ * ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝ * 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. ******************************************************/ #include "mkernel_tesseroid.h" #ifdef GCTL_POTENTIAL_MAGTESS #include "magtess/constants.h" #include "magtess/parsers.h" #include "magtess/linalg.h" #include "magtess/geometry.h" #include "magtess/grav_tess.h" void gctl::set_geomag_angles(array &out_para, const array &inclina, const array &declina) { if (out_para.size() != inclina.size() || out_para.size() != declina.size()) { throw std::runtime_error("[gctl::set_geomag_angles] Invalid sizes."); } for (size_t i = 0; i < out_para.size(); i++) { out_para[i].geo_iarc = inclina[i]*GCTL_Pi/180.0; out_para[i].geo_darc = declina[i]*GCTL_Pi/180.0; } return; } void gctl::callink_magnetic_para(array &in_tess, array &out_para, const point3dc &mag_B) { if (in_tess.size() != out_para.size()) { throw std::runtime_error("[gctl::callink_magnetic_para] Invalid sizes."); } for (size_t i = 0; i < in_tess.size(); i++) { // 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 out_para[i].bx = mag_B.x; out_para[i].by = mag_B.y; out_para[i].bz = mag_B.z; out_para[i].h = sqrt(power2(mag_B.x) + power2(mag_B.y)); out_para[i].f = sqrt(power2(mag_B.x) + power2(mag_B.y) + power2(mag_B.z)); in_tess[i].att = out_para.get(i); } return; } void gctl::callink_magnetic_para(array &in_tess, array &out_para, const array &mag_B) { if (in_tess.size() != out_para.size() || in_tess.size() != mag_B.size()) { throw std::runtime_error("[gctl::callink_magnetic_para] Invalid sizes."); } for (size_t i = 0; i < in_tess.size(); i++) { // 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 out_para[i].bx = mag_B[i].x; out_para[i].by = mag_B[i].y; out_para[i].bz = mag_B[i].z; out_para[i].h = sqrt(power2(mag_B[i].x) + power2(mag_B[i].y)); out_para[i].f = sqrt(power2(mag_B[i].x) + power2(mag_B[i].y) + power2(mag_B[i].z)); in_tess[i].att = out_para.get(i); } return; } void gctl::callink_magnetic_para_earth_sph(array &in_tess, array &out_para, double inclina_deg, double declina_deg, 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(...)"); } double I = inclina_deg*M_PI/180.0; double A = declina_deg*M_PI/180.0; point3dc mag_z; for (size_t i = 0; i < in_tess.size(); i++) { // 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); // magnetic susceptibility is taken as one here mag_z = field_tense * mag_z.normal(); out_para[i].bx = mag_z.x; out_para[i].by = mag_z.y; out_para[i].bz = mag_z.z; out_para[i].h = sqrt(power2(mag_z.x) + power2(mag_z.y)); out_para[i].f = sqrt(power2(mag_z.x) + power2(mag_z.y) + power2(mag_z.z)); in_tess[i].att = out_para.get(i); } return; } typedef void (*magkernel_tess_ptr)(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose); void magkernel_mag_tesseroid_hax(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose); void magkernel_mag_tesseroid_hay(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose); void magkernel_mag_tesseroid_za(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose); void magkernel_mag_tesseroid_deltaT(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose); void gctl::magkernel(matrix &out_kernel, const array &ele, const array &ops, magnetic_field_type_e comp_id, int gauss_order, verbose_type_e verbose) { magkernel_tess_ptr mag_tesseroid_kernel; switch (comp_id) { case Hax: mag_tesseroid_kernel = magkernel_mag_tesseroid_hax; break; case Hay: mag_tesseroid_kernel = magkernel_mag_tesseroid_hay; break; case Za: mag_tesseroid_kernel = magkernel_mag_tesseroid_za; break; case DeltaT: mag_tesseroid_kernel = magkernel_mag_tesseroid_deltaT; break; default: mag_tesseroid_kernel = magkernel_mag_tesseroid_za; break; } return mag_tesseroid_kernel(out_kernel, ele, ops, gauss_order, verbose); } typedef void (*magobser_tess_ptr)(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose); void magobser_mag_tesseroid_hax(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose); void magobser_mag_tesseroid_hay(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose); void magobser_mag_tesseroid_za(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose); void magobser_mag_tesseroid_deltaT(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose); void gctl::magobser(array &out_obs, const array &ele, const array &ops, const array &sus, magnetic_field_type_e comp_id, int gauss_order, verbose_type_e verbose) { magobser_tess_ptr mag_tesseroid_obser; switch (comp_id) { case Hax: mag_tesseroid_obser = magobser_mag_tesseroid_hax; break; case Hay: mag_tesseroid_obser = magobser_mag_tesseroid_hay; break; case Za: mag_tesseroid_obser = magobser_mag_tesseroid_za; break; case DeltaT: mag_tesseroid_obser = magobser_mag_tesseroid_deltaT; break; default: mag_tesseroid_obser = magobser_mag_tesseroid_za; break; } return mag_tesseroid_obser(out_obs, ele, ops, sus, gauss_order, verbose); } // 以下是具体的实现 void magkernel_mag_tesseroid_hax(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_kernel.resize(o_size, e_size); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[3]; double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magkernel_mag_tesseroid_hax(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magkernel_hax"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = 1.0; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; // We implemented the conversion in the function callink_magnetic_para_earth_sph() M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; #pragma omp parallel for private(i, gtt_v, M_vect_p, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxx, TESSEROID_GXX_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxy, TESSEROID_GXY_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxz, TESSEROID_GXZ_SIZE_RATIO); // This is in nT, and we reversed the positive direction here to the easting direction out_kernel[i][j] = -1.0*M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magkernel_mag_tesseroid_hay(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_kernel.resize(o_size, e_size); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[3]; double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magkernel_mag_tesseroid_hay(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magkernel_hay"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = 1.0; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; // We implemented the conversion in the function callink_magnetic_para_earth_sph() M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; #pragma omp parallel for private(i, gtt_v, M_vect_p, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxy, TESSEROID_GXY_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyy, TESSEROID_GYY_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyz, TESSEROID_GYZ_SIZE_RATIO); // This is in nT out_kernel[i][j] = M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magkernel_mag_tesseroid_za(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_kernel.resize(o_size, e_size); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[3]; double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magkernel_mag_tesseroid_za(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magkernel_za"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = 1.0; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; // We implemented the conversion in the function callink_magnetic_para_earth_sph() M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; #pragma omp parallel for private(i, gtt_v, M_vect_p, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxz, TESSEROID_GXZ_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyz, TESSEROID_GYZ_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gzz, TESSEROID_GZZ_SIZE_RATIO); // This is in nT, and we reversed the positive direction here to the downward direction out_kernel[i][j] = -1.0*M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magkernel_mag_tesseroid_deltaT(gctl::matrix &out_kernel, const gctl::array &ele, const gctl::array &ops, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_kernel.resize(o_size, e_size); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[6]; //double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; double Hax, Hay, Za; double hax_f, hay_f, za_f; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magkernel_mag_tesseroid_za(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magkernel_deltaT"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = 1.0; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); //B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; // We implemented the conversion in the function callink_magnetic_para_earth_sph() M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; hax_f = cos(mp->geo_iarc)*sin(GCTL_Pi*mp->geo_darc); hay_f = cos(mp->geo_iarc)*cos(GCTL_Pi*mp->geo_darc); za_f = sin(mp->geo_iarc); #pragma omp parallel for private(i, gtt_v, M_vect_p, Hax, Hay, Za, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxx, TESSEROID_GXX_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxy, TESSEROID_GXY_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxz, TESSEROID_GXZ_SIZE_RATIO); gtt_v[3] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyy, TESSEROID_GYY_SIZE_RATIO); gtt_v[4] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyz, TESSEROID_GYZ_SIZE_RATIO); gtt_v[5] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gzz, TESSEROID_GZZ_SIZE_RATIO); // This is in nT, and we reversed the positive direction here to the downward direction Hax = M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); Hay = M_0*(gtt_v[1]*M_vect_p[0] + gtt_v[3]*M_vect_p[1] + gtt_v[4]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); Za = -1.0*M_0*(gtt_v[2]*M_vect_p[0] + gtt_v[4]*M_vect_p[1] + gtt_v[5]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); out_kernel[i][j] = hax_f*Hax + hay_f*Hay + za_f*Za; glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magobser_mag_tesseroid_hax(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_obs.resize(o_size, 0.0); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[3]; //double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magobser_mag_tesseroid_hax(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magobser_hax"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = sus[j]; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); //B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; // We implemented the conversion in the function callink_magnetic_para_earth_sph() M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; #pragma omp parallel for private(i, gtt_v, M_vect_p, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxx, TESSEROID_GXX_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxy, TESSEROID_GXY_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxz, TESSEROID_GXZ_SIZE_RATIO); // This is in nT, and we reversed the positive direction here to the easting direction out_obs[i] -= sus[j]*M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magobser_mag_tesseroid_hay(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_obs.resize(o_size, 0.0); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[3]; double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magobser_mag_tesseroid_hay(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magobser_hay"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = sus[j]; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; // We implemented the conversion in the function callink_magnetic_para_earth_sph() M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; #pragma omp parallel for private(i, gtt_v, M_vect_p, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxy, TESSEROID_GXY_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyy, TESSEROID_GYY_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyz, TESSEROID_GYZ_SIZE_RATIO); // This is in nT out_obs[i] += sus[j]*M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magobser_mag_tesseroid_za(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_obs.resize(o_size, 0.0); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[3]; //double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magobser_mag_tesseroid_za(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magobser_za"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = sus[j]; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); // We implemented the conversion in the function callink_magnetic_para_earth_sph() //B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; #pragma omp parallel for private(i, gtt_v, M_vect_p, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxz, TESSEROID_GXZ_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyz, TESSEROID_GYZ_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gzz, TESSEROID_GZZ_SIZE_RATIO); // This is in nT, and we reversed the positive direction here to the downward direction out_obs[i] -= sus[j]*M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } void magobser_mag_tesseroid_deltaT(gctl::array &out_obs, const gctl::array &ele, const gctl::array &ops, const gctl::array &sus, int gauss_order, gctl::verbose_type_e verbose) { int i, j; int o_size = ops.size(); int e_size = ele.size(); out_obs.resize(o_size, 0.0); MAG_TESSEROID tess; GLQ *glqlon, *glqlat, *glqr; double gtt_v[6]; //double B_to_H; double M_vect[3] = {0, 0, 0}; double M_vect_p[3] = {0, 0, 0}; double Hax, Hay, Za; double hax_f, hay_f, za_f; for (j = 0; j < e_size; j++) { if (ele[j].att == nullptr) throw gctl::runtime_error("Magnetization parameter not set. from gctl::magobser_mag_tesseroid_deltaT(...)"); } gctl::magtess_para *mp; gctl::progress_bar bar(e_size, "magobser_deltaT"); for (j = 0; j < e_size; j++) { if (verbose == gctl::FullMsg) bar.progressed(j); else if (verbose == gctl::ShortMsg) bar.progressed_simple(j); mp = ele[j].att; tess.suscept = sus[j]; // unit susception tess.density = 1.0; // unit density tess.w = ele[j].dl->lon; tess.e = ele[j].ur->lon; tess.s = ele[j].dl->lat; tess.n = ele[j].ur->lat; tess.r1 = ele[j].dl->rad; tess.r2 = ele[j].ur->rad; tess.Bx = mp->bx/(400.0*GCTL_Pi); tess.By = mp->by/(400.0*GCTL_Pi); tess.Bz = mp->bz/(400.0*GCTL_Pi); tess.cos_a1 = cos(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.sin_a1 = sin(GCTL_Pi/2.0-GCTL_Pi*(tess.w+tess.e)*0.5/180.0); tess.cos_b1 = cos(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); tess.sin_b1 = sin(GCTL_Pi*(tess.s+tess.n)*0.5/180.0); // We implemented the conversion in the function callink_magnetic_para_earth_sph() //B_to_H = tess.suscept/(M_0); M_vect[0] = tess.By;// * B_to_H; M_vect[1] = tess.Bx;// * B_to_H; M_vect[2] = tess.Bz;// * B_to_H; hax_f = cos(mp->geo_iarc)*sin(GCTL_Pi*mp->geo_darc); hay_f = cos(mp->geo_iarc)*cos(GCTL_Pi*mp->geo_darc); za_f = sin(mp->geo_iarc); #pragma omp parallel for private(i, gtt_v, M_vect_p, Hax, Hay, Za, glqlon, glqlat, glqr) schedule(guided) for (i = 0; i < o_size; i++) { conv_vect_cblas(M_vect, (tess.w + tess.e)*0.5, (tess.s + tess.n)*0.5, ops[i].lon, ops[i].lat, M_vect_p); glqlon = glq_new(gauss_order, -1, 1); glqlat = glq_new(gauss_order, -1, 1); glqr = glq_new(gauss_order, -1, 1); gtt_v[0] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxx, TESSEROID_GXX_SIZE_RATIO); gtt_v[1] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxy, TESSEROID_GXY_SIZE_RATIO); gtt_v[2] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gxz, TESSEROID_GXZ_SIZE_RATIO); gtt_v[3] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyy, TESSEROID_GYY_SIZE_RATIO); gtt_v[4] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gyz, TESSEROID_GYZ_SIZE_RATIO); gtt_v[5] = calc_tess_model_adapt(&tess, 1, ops[i].lon, ops[i].lat, ops[i].rad, glqlon, glqlat, glqr, tess_gzz, TESSEROID_GZZ_SIZE_RATIO); // This is in nT, and we reversed the positive direction here to the downward direction Hax = M_0*(gtt_v[0]*M_vect_p[0] + gtt_v[1]*M_vect_p[1] + gtt_v[2]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); Hay = M_0*(gtt_v[1]*M_vect_p[0] + gtt_v[3]*M_vect_p[1] + gtt_v[4]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); Za = -1.0*M_0*(gtt_v[2]*M_vect_p[0] + gtt_v[4]*M_vect_p[1] + gtt_v[5]*M_vect_p[2])/(GCTL_G0*tess.density*4*GCTL_Pi); out_obs[i] += sus[j]*(hax_f*Hax + hay_f*Hay + za_f*Za); glq_free(glqlon); glq_free(glqlat); glq_free(glqr); } } return; } #endif // GCTL_POTENTIAL_MAGTESS