/******************************************************** * ██████╗ ██████╗████████╗██╗ * ██╔════╝ ██╔════╝╚══██╔══╝██║ * ██║ ███╗██║ ██║ ██║ * ██║ ██║██║ ██║ ██║ * ╚██████╔╝╚██████╗ ██║ ███████╗ * ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝ * 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 "gm_regular_grid.h" #ifdef GCTL_POTENTIAL_FFTW3 // include fftw head files #include "fftw3.h" #endif //GCTL_POTENTIAL_FFTW3 gctl::gm_regular_grid::gm_regular_grid(){} gctl::gm_regular_grid::~gm_regular_grid(){} gctl::gm_regular_grid::gm_regular_grid(std::string in_name, std::string in_info, int xnum, int ynum, double xmin, double ymin, double dx, double dy) : regular_grid::regular_grid(in_name, in_info, xnum, ynum, xmin, ymin, dx, dy){} #ifdef GCTL_POTENTIAL_FFTW3 void gctl::gm_regular_grid::gradient(std::string datname, std::string gradname, gravitational_field_type_e c_type, int order) { meshdata &data = get_data(datname); mesh_data_value_e value_type = data.valtype_; mesh_data_type_e data_type = data.loctype_; if(value_type != Scalar) { throw std::runtime_error("[gctl::gm_regular_grid::gradient] Invalid value type."); } if (c_type != Tzz && c_type != Tzx && c_type != Tzy) { throw std::runtime_error("[gctl::gm_regular_grid::gradient] The gradient type must be Tzz, Tzx and Tzy."); } // 检查是否存在与梯度数据同名的数据 若有则检查是否需要重新建立数据 meshdata &new_data = add_data(data_type, value_type, gradname, 0.0); gctl::array ingrid_ex; int M = -1, N = -1; int ori_row, ori_col; if (data_type == NodeData) gctl::cosine_extrapolate_2d(data.datval_, ingrid_ex, rg_ynum, rg_xnum, M, N, ori_row, ori_col); else gctl::cosine_extrapolate_2d(data.datval_, ingrid_ex, rg_ynum-1, rg_xnum-1, M, N, ori_row, ori_col); fftw_complex *in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*M*N); fftw_complex *out= (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*M*N); for (int i = 0; i < M*N; i++) { in[i][0] = ingrid_ex[i]; in[i][1] = 0.0; } fftw_plan p = fftw_plan_dft_2d(M, N, in, out, FFTW_FORWARD, FFTW_ESTIMATE); fftw_execute(p); int half_M = ceil(M/2); int half_N = ceil(N/2); double du = 2.0*GCTL_Pi/((M-1)*rg_dy); double dv = 2.0*GCTL_Pi/((N-1)*rg_dx); int idx1, idx2, idx3, idx4; double u, v; if (c_type == Tzx) { for (int i = 0; i < M; i++) { for (int j = 0; j < half_N; j++) { v = dv*(j+1); idx1 = i*N+j; idx2 = i*N+N-1-j; in[idx1][0] = -1.0*out[idx1][1]*pow(v, order); in[idx1][1] = out[idx1][0]*pow(v, order); in[idx2][0] = -1.0*out[idx2][1]*pow(-1.0*v, order); in[idx2][1] = out[idx2][0]*pow(-1.0*v, order); } } } else if (c_type == Tzy) { for (int i = 0; i < half_M; i++) { u = du*(i+1); for (int j = 0; j < N; j++) { idx1 = i*N+j; idx3 = (M-1-i)*N+j; in[idx1][0] = -1.0*out[idx1][1]*pow(u, order); in[idx1][1] = out[idx1][0]*pow(u, order); in[idx3][0] = -1.0*out[idx3][1]*pow(-1.0*u, order); in[idx3][1] = out[idx3][0]*pow(-1.0*u, order); } } } else { for (int i = 0; i < half_M; i++) { u = du*(i+1); for (int j = 0; j < half_N; j++) { v = dv*(j+1); idx1 = i*N+j; idx2 = i*N+N-1-j; idx3 = (M-1-i)*N+j; idx4 = (M-1-i)*N+N-1-j; in[idx1][0] = out[idx1][0]*pow(sqrt(u*u+v*v), order); in[idx1][1] = out[idx1][1]*pow(sqrt(u*u+v*v), order); in[idx2][0] = out[idx2][0]*pow(sqrt(u*u+v*v), order); in[idx2][1] = out[idx2][1]*pow(sqrt(u*u+v*v), order); in[idx3][0] = out[idx3][0]*pow(sqrt(u*u+v*v), order); in[idx3][1] = out[idx3][1]*pow(sqrt(u*u+v*v), order); in[idx4][0] = out[idx4][0]*pow(sqrt(u*u+v*v), order); in[idx4][1] = out[idx4][1]*pow(sqrt(u*u+v*v), order); } } } p = fftw_plan_dft_2d(M, N, in, out, FFTW_BACKWARD, FFTW_ESTIMATE); fftw_execute(p); for (int i = 0; i < M; i++) { for (int j = 0; j < N; j++) { ingrid_ex[i*N+j] = out[i*N+j][0]/(M*N); } } if (data_type == NodeData) { for (int i = 0; i < rg_ynum; i++) { for (int j = 0; j < rg_xnum; j++) { new_data.datval_[i*rg_xnum+j] = ingrid_ex[(i+ori_row)*N+j+ori_col]; } } } else { for (int i = 0; i < rg_ynum-1; i++) { for (int j = 0; j < rg_xnum-1; j++) { new_data.datval_[i*(rg_xnum-1)+j] = ingrid_ex[(i+ori_row)*N+j+ori_col]; } } } fftw_destroy_plan(p); fftw_free(in); fftw_free(out); return; } void gctl::gm_regular_grid::rtp(std::string datname, std::string rtpname, double inc, double dec) { meshdata &data = get_data(datname); mesh_data_value_e value_type = data.valtype_; mesh_data_type_e data_type = data.loctype_; if(value_type != Scalar) { std::string err_str = "Invalid value type. From gctl::gm_regular_grid::rtp(...)"; throw runtime_error(err_str); } // 检查是否存在与梯度数据同名的数据 若有则检查是否需要重新建立数据 meshdata &new_data = add_data(data_type, value_type, rtpname, 0.0); gctl::array ingrid_ex; int M = -1, N = -1; int ori_row, ori_col; if (data_type == NodeData) gctl::cosine_extrapolate_2d(data.datval_, ingrid_ex, rg_ynum, rg_xnum, M, N, ori_row, ori_col); else gctl::cosine_extrapolate_2d(data.datval_, ingrid_ex, rg_ynum-1, rg_xnum-1, M, N, ori_row, ori_col); fftw_complex *in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*M*N); fftw_complex *out= (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*M*N); for (int i = 0; i < M*N; i++) { in[i][0] = ingrid_ex[i]; in[i][1] = 0.0; } fftw_plan p = fftw_plan_dft_2d(M, N, in, out, FFTW_FORWARD, FFTW_ESTIMATE); fftw_execute(p); int half_M = ceil(M/2); int half_N = ceil(N/2); double du = 2.0*GCTL_Pi/((M-1)*rg_dy); double dv = 2.0*GCTL_Pi/((N-1)*rg_dx); double cc = cosd(dec)*cosd(inc); double sc = sind(dec)*cosd(inc); double s = sind(inc); int idx; double U, V, Rpart, Ipart; double temp, temp1, temp2, temp3, temp4, base; for (int i = 0; i < half_M; i++) { for (int j = 0; j < half_N; j++) { U = (i+1)*du; V = (j+1)*dv; temp = U*U+V*V; idx = i*N+j; temp1 = cc*U+sc*V; temp2 = s*s*temp-temp1*temp1; base = temp2*temp2+4.0*s*s*temp1*temp1*temp; Rpart = temp*temp2/base; Ipart = -2.0*s*temp1*pow(temp,1.5)/base; temp3 = out[idx][0]; temp4 = out[idx][1]; in[idx][0] = Rpart*temp3-temp4*Ipart; in[idx][1] = Rpart*temp4+Ipart*temp3; idx = (M-1-i)*N+j; temp1 = -cc*U+sc*V; temp2 = s*s*temp-temp1*temp1; base = temp2*temp2+4.0*s*s*temp1*temp1*temp; Rpart = temp*temp2/base; Ipart = -2.0*s*temp1*pow(temp,1.5)/base; temp3 = out[idx][0]; temp4 = out[idx][1]; in[idx][0] = Rpart*temp3-temp4*Ipart; in[idx][1] = Rpart*temp4+Ipart*temp3; idx = i*N+N-1-j; temp1 = cc*U-sc*V; temp2 = s*s*temp-temp1*temp1; base = temp2*temp2+4.0*s*s*temp1*temp1*temp; Rpart = temp*temp2/base; Ipart = -2.0*s*temp1*pow(temp,1.5)/base; temp3 = out[idx][0]; temp4 = out[idx][1]; in[idx][0] = Rpart*temp3-temp4*Ipart; in[idx][1] = Rpart*temp4+Ipart*temp3; idx = (M-1-i)*N+N-1-j; temp1 = -cc*U-sc*V; temp2 = s*s*temp-temp1*temp1; base = temp2*temp2+4.0*s*s*temp1*temp1*temp; Rpart = temp*temp2/base; Ipart = -2.0*s*temp1*pow(temp,1.5)/base; temp3 = out[idx][0]; temp4 = out[idx][1]; in[idx][0] = Rpart*temp3-temp4*Ipart; in[idx][1] = Rpart*temp4+Ipart*temp3; } } p = fftw_plan_dft_2d(M, N, in, out, FFTW_BACKWARD, FFTW_ESTIMATE); fftw_execute(p); for (int i = 0; i < M; i++) { for (int j = 0; j < N; j++) { ingrid_ex[i*N+j] = out[i*N+j][0]/(M*N); } } if (data_type == NodeData) { for (int i = 0; i < rg_ynum; i++) { for (int j = 0; j < rg_xnum; j++) { new_data.datval_[i*rg_xnum+j] = ingrid_ex[(i+ori_row)*N+j+ori_col]; } } } else { for (int i = 0; i < rg_ynum-1; i++) { for (int j = 0; j < rg_xnum-1; j++) { new_data.datval_[i*(rg_xnum-1)+j] = ingrid_ex[(i+ori_row)*N+j+ori_col]; } } } fftw_destroy_plan(p); fftw_free(in); fftw_free(out); return; } void gctl::gm_regular_grid::drtp(std::string datname, std::string incname, std::string decname, std::string drtpname, int order) { // 首先获取数据并确定所有数据类型和格式都是一致的 meshdata *md_ptr = get_data_ptr(datname); mesh_data_value_e dat_valtype = md_ptr->valtype_; mesh_data_type_e dat_datype = md_ptr->loctype_; md_ptr = get_data_ptr(incname); mesh_data_value_e inc_valtype = md_ptr->valtype_; mesh_data_type_e inc_datype = md_ptr->loctype_; md_ptr = get_data_ptr(decname); mesh_data_value_e dec_valtype = md_ptr->valtype_; mesh_data_type_e dec_datype = md_ptr->loctype_; std::string err_str; if (dat_valtype != Scalar) { err_str = "Invalid value type. From gctl::gm_regular_grid::drtp(...)"; throw runtime_error(err_str); } if (dat_datype != inc_datype || dat_datype != dec_datype) { err_str = "Incompatible data type. From gctl::gm_regular_grid::drtp(...)"; throw runtime_error(err_str); } // 获取数据指针 md_ptr = get_data_ptr(datname); array *dat_ptr = &md_ptr->datval_; md_ptr = get_data_ptr(incname); array *inc_ptr = &md_ptr->datval_; md_ptr = get_data_ptr(decname); array *dec_ptr = &md_ptr->datval_; // 检查是否存在同名的数据 若有则检查是否需要重新建立数据 add_data(dat_datype, dat_valtype, drtpname, 0.0); md_ptr = get_data_ptr(drtpname); // 获取数据指针 array *drtp_ptr = &md_ptr->datval_; // 新建临时数据并获取数据指针 std::string tmp_name = drtpname+"tmp"; add_data(dat_datype, dat_valtype, tmp_name, 0.0); md_ptr = get_data_ptr(tmp_name); array *tmp_ptr = &md_ptr->datval_; std::string diff1_name = drtpname+"diff1"; add_data(dat_datype, dat_valtype, diff1_name, 0.0); md_ptr = get_data_ptr(diff1_name); array *diff1_ptr = &md_ptr->datval_; std::string diff2_name = drtpname+"diff2"; add_data(dat_datype, dat_valtype, diff2_name, 0.0); array *diff2_ptr = &md_ptr->datval_; // 计算平均磁化参数 double meaninc = 0.0, meandec = 0.0; for (int i = 0; i < inc_ptr->size(); i++) { meaninc += inc_ptr->at(i); meandec += dec_ptr->at(i); } meaninc /= inc_ptr->size(); meandec /= dec_ptr->size(); // 计算平均化极结果 rtp(datname, tmp_name, meaninc, meandec); for (int i = 0; i < drtp_ptr->size(); i++) { drtp_ptr->at(i) = tmp_ptr->at(i); } // 使用微小扰动计算RTP相对与倾角与偏角的导数 double dinc = 0.01, ddec = 0.01; // 计算一阶变角度化极结果 rtp(datname, diff1_name, meaninc+dinc, meandec); rtp(datname, diff2_name, meaninc, meandec+ddec); for (int i = 0; i < drtp_ptr->size(); i++) { diff1_ptr->at(i) = (inc_ptr->at(i) - meaninc) *(diff1_ptr->at(i) - tmp_ptr->at(i))/dinc; drtp_ptr->at(i) += diff1_ptr->at(i); diff2_ptr->at(i) = (dec_ptr->at(i) - meandec) *(diff2_ptr->at(i) - tmp_ptr->at(i))/ddec; drtp_ptr->at(i) += diff2_ptr->at(i); } double factor; for (int o = 2; o <= order; o++) { factor = 1.0; for (int n = o; n > 1; n--) { factor *= n; } factor = 1.0/factor; rtp(diff1_name, tmp_name, meaninc, meandec); rtp(diff1_name, diff1_name, meaninc+dinc, meandec); for (int i = 0; i < drtp_ptr->size(); i++) { diff1_ptr->at(i) = factor*pow(inc_ptr->at(i) - meaninc, o) *(diff1_ptr->at(i) - tmp_ptr->at(i))/dinc; drtp_ptr->at(i) += diff1_ptr->at(i); } rtp(diff2_name, tmp_name, meaninc, meandec); rtp(diff2_name, diff2_name, meaninc, meandec+ddec); for (int i = 0; i < drtp_ptr->size(); i++) { diff2_ptr->at(i) = factor*pow(dec_ptr->at(i) - meandec, o) *(diff2_ptr->at(i) - tmp_ptr->at(i))/ddec; drtp_ptr->at(i) += diff2_ptr->at(i); } } // 删除临时数据 remove_data(tmp_name); remove_data(diff1_name); remove_data(diff2_name); return; } void gctl::gm_regular_grid::conti(std::string datname, std::string retname, double height) { meshdata &data = get_data(datname); mesh_data_value_e value_type = data.valtype_; mesh_data_type_e data_type = data.loctype_; if(value_type != Scalar) { throw std::runtime_error("[gctl::gm_regular_grid::conti] Invalid value type."); } // 检查是否存在与梯度数据同名的数据 若有则检查是否需要重新建立数据 meshdata &new_data = add_data(data_type, value_type, retname, 0.0); gctl::array ingrid_ex; int M = -1, N = -1; int ori_row, ori_col; if (data_type == NodeData) gctl::cosine_extrapolate_2d(data.datval_, ingrid_ex, rg_ynum, rg_xnum, M, N, ori_row, ori_col); else gctl::cosine_extrapolate_2d(data.datval_, ingrid_ex, rg_ynum-1, rg_xnum-1, M, N, ori_row, ori_col); fftw_complex *in = (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*M*N); fftw_complex *out= (fftw_complex*) fftw_malloc(sizeof(fftw_complex)*M*N); for (int i = 0; i < M*N; i++) { in[i][0] = ingrid_ex[i]; in[i][1] = 0.0; } fftw_plan p = fftw_plan_dft_2d(M, N, in, out, FFTW_FORWARD, FFTW_ESTIMATE); fftw_execute(p); int half_M = ceil(M/2); int half_N = ceil(N/2); double du = 2.0*GCTL_Pi/((M-1)*rg_dy); double dv = 2.0*GCTL_Pi/((N-1)*rg_dx); int idx1, idx2, idx3, idx4; double u, v; for (int i = 0; i < half_M; i++) { u = du*(i+1); for (int j = 0; j < half_N; j++) { v = dv*(j+1); idx1 = i*N+j; idx2 = i*N+N-1-j; idx3 = (M-1-i)*N+j; idx4 = (M-1-i)*N+N-1-j; in[idx1][0] = out[idx1][0]*exp(sqrt(u*u+v*v)*height); in[idx1][1] = out[idx1][1]*exp(sqrt(u*u+v*v)*height); in[idx2][0] = out[idx2][0]*exp(sqrt(u*u+v*v)*height); in[idx2][1] = out[idx2][1]*exp(sqrt(u*u+v*v)*height); in[idx3][0] = out[idx3][0]*exp(sqrt(u*u+v*v)*height); in[idx3][1] = out[idx3][1]*exp(sqrt(u*u+v*v)*height); in[idx4][0] = out[idx4][0]*exp(sqrt(u*u+v*v)*height); in[idx4][1] = out[idx4][1]*exp(sqrt(u*u+v*v)*height); } } p = fftw_plan_dft_2d(M, N, in, out, FFTW_BACKWARD, FFTW_ESTIMATE); fftw_execute(p); for (int i = 0; i < M; i++) { for (int j = 0; j < N; j++) { ingrid_ex[i*N+j] = out[i*N+j][0]/(M*N); } } if (data_type == NodeData) { for (int i = 0; i < rg_ynum; i++) { for (int j = 0; j < rg_xnum; j++) { new_data.datval_[i*rg_xnum+j] = ingrid_ex[(i+ori_row)*N+j+ori_col]; } } } else { for (int i = 0; i < rg_ynum-1; i++) { for (int j = 0; j < rg_xnum-1; j++) { new_data.datval_[i*(rg_xnum-1)+j] = ingrid_ex[(i+ori_row)*N+j+ori_col]; } } } fftw_destroy_plan(p); fftw_free(in); fftw_free(out); return; } #endif //GCTL_POTENTIAL_FFTW3 void gctl::gm_regular_grid::trend(std::string datname, std::string regname, std::string resname, int wx_size, int wy_size, int x_order, int y_order) { meshdata &data = get_data(datname); mesh_data_value_e value_type = data.valtype_; mesh_data_type_e data_type = data.loctype_; if(value_type != Scalar) { throw std::runtime_error("Invalid value type. From gctl::gm_regular_grid::trend(...)"); } // 检查是否存在与趋势场局部场数据同名的数据 若有则检查是否需要重新建立数据 meshdata ®_data = add_data(data_type, value_type, regname, 0.0); meshdata &res_data = add_data(data_type, value_type, resname, 0.0); matrix data_mat; if (data_type == NodeData) { data_mat.resize(rg_ynum, rg_xnum); for (size_t i = 0; i < rg_ynum; i++) { for (size_t j = 0; j < rg_xnum; j++) { data_mat[i][j] = data.datval_[i*rg_xnum+j]; } } } else { data_mat.resize(rg_ynum-1, rg_xnum-1); for (size_t i = 0; i < rg_ynum-1; i++) { for (size_t j = 0; j < rg_xnum-1; j++) { data_mat[i][j] = data.datval_[i*(rg_xnum-1)+j]; } } } trend_2d(data_mat, wy_size, wx_size, y_order, x_order); if (data_type == NodeData) { for (size_t i = 0; i < rg_ynum; i++) { for (size_t j = 0; j < rg_xnum; j++) { reg_data.datval_[i*rg_xnum+j] = data_mat[i][j]; res_data.datval_[i*rg_xnum+j] = data.datval_[i*rg_xnum+j] - data_mat[i][j]; } } } else { for (size_t i = 0; i < rg_ynum-1; i++) { for (size_t j = 0; j < rg_xnum-1; j++) { reg_data.datval_[i*(rg_xnum-1)+j] = data_mat[i][j]; res_data.datval_[i*(rg_xnum-1)+j] = data.datval_[i*(rg_xnum-1)+j] - data_mat[i][j]; } } } return; }