initial upload
This commit is contained in:
24
lki/CMakeLists.txt
Normal file
24
lki/CMakeLists.txt
Normal file
@@ -0,0 +1,24 @@
|
||||
set(TOOL_NAME lki)
|
||||
set(BIN_DIR bin)
|
||||
set(INSTALL_DIR sbin)
|
||||
|
||||
find_package(GCTL_OPTIMIZATION REQUIRED)
|
||||
include_directories(${GCTL_OPTIMIZATION_INC_DIR})
|
||||
|
||||
set(CMAKE_C_FLAGS "${CMAKE_C_FLAGS} -O3")
|
||||
if(WIN32)
|
||||
set(CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} -O3")
|
||||
else()
|
||||
set(CMAKE_CXX_FLAGS "${CMAKE_CXX_FLAGS} -O3")
|
||||
endif()
|
||||
|
||||
set(EXECUTABLE_OUTPUT_PATH ${PROJECT_BINARY_DIR}/${BIN_DIR})
|
||||
|
||||
aux_source_directory(. TOOL_SRC)
|
||||
add_executable(${TOOL_NAME} ${TOOL_SRC})
|
||||
set_target_properties(${TOOL_NAME} PROPERTIES INSTALL_RPATH /usr/local/lib)
|
||||
set_target_properties(${TOOL_NAME} PROPERTIES CXX_STANDARD 17 CXX_STANDARD_REQUIRED ON)
|
||||
|
||||
target_link_libraries(${TOOL_NAME} PUBLIC ${GCTL_OPTIMIZATION_LIB})
|
||||
|
||||
install(TARGETS ${TOOL_NAME} RUNTIME DESTINATION ${INSTALL_DIR})
|
96
lki/lki.cpp
Normal file
96
lki/lki.cpp
Normal file
@@ -0,0 +1,96 @@
|
||||
/********************************************************
|
||||
* ██████╗ ██████╗████████╗██╗
|
||||
* ██╔════╝ ██╔════╝╚══██╔══╝██║
|
||||
* ██║ ███╗██║ ██║ ██║
|
||||
* ██║ ██║██║ ██║ ██║
|
||||
* ╚██████╔╝╚██████╗ ██║ ███████╗
|
||||
* ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝
|
||||
* 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 "lki.h"
|
||||
|
||||
int main(int argc, char *argv[]) try
|
||||
{
|
||||
gctl::flags_parser fp;
|
||||
fp.set_proname("lki");
|
||||
fp.set_proinfo("Limited memory Kriging Interpolation of 2D arbitrary located data. \
|
||||
This program is a toolkit of the GCTL package. The GCTL comes with ABSOLUTE NO WARRANTY. \
|
||||
Please see instructions or contact the author for more information.");
|
||||
fp.add_opt('i', "in-table", required_argument, NULL, "Input table name of the constraint locations and data.", "<file-name>", true);
|
||||
fp.add_opt('t', "target-table", required_argument, NULL, "Input table name of the target locations or parameters of output grid locations.", "<file-name>|<xmin>/<dx>/<xmax>/<ymin>/<dx>/<ymax>", true);
|
||||
fp.add_opt('o', "out-table", required_argument, NULL, "Output table name of the interpolated data.", "<file-name>", true);
|
||||
fp.add_opt('k', "kernel-size", required_argument, NULL, "Kernel size for solving the linear system of the Variogram functions (default is 200).", "<size>", false);
|
||||
fp.add_opt('e', "epsilon", required_argument, NULL, "Threshold for solving the linear system of the Variogram functions (default is 1e-10).", "<value>", false);
|
||||
fp.add_opt('b', "box-dimension", required_argument, NULL, "2D boxes' dimensions used for pre-processing the constraint data (default is 10).", "<size>", false);
|
||||
fp.add_opt('v', "variogram", required_argument, NULL, "Variogram function used for interpolation. Available functions are spherical (default), gaussian, \
|
||||
wave, rational, circular and linear.", "<function>", false);
|
||||
fp.add_opt('p', "parameter", required_argument, NULL, "Parameters of the variogram function. The default values are 0.1, 10 and 100, respectively.", "<nugget>/<sill>/<range>", false);
|
||||
fp.add_opt('c', "column", required_argument, NULL, "Column index of x, y and data. The defaults are 0,1,2", "<col-x>,<col-y>,<col-data>", false);
|
||||
fp.add_opt('h', "head-record", required_argument, NULL, "Input file has [0] Header record(s).", "<size>", false);
|
||||
fp.add_opt('a', "annotate", required_argument, NULL, "The starting symbol of an annotate line.", "<sym>", false);
|
||||
fp.add_opt('d', "delimeter", required_argument, NULL, "Delimeter between different columns.", "<sym>", false);
|
||||
fp.configure(argc, argv);
|
||||
|
||||
if (argc == 1)
|
||||
{
|
||||
fp.show_help_page();
|
||||
return 0;
|
||||
}
|
||||
|
||||
std::string in_name, out_name, tar_name, kernel_str, eps_str, box_str, variogram_str, para_str, col_str, head_str, ant_str, del_str;
|
||||
fp.get_argv({'i', 't', 'o', 'k', 'e', 'b', 'v', 'p', 'c'},
|
||||
{&in_name, &tar_name, &out_name, &kernel_str, &eps_str, &box_str, &variogram_str, ¶_str, &col_str, &head_str, &ant_str, &del_str});
|
||||
// 查看是否通过强制参数检查
|
||||
if (!fp.pass_mandatory()) return 0;
|
||||
|
||||
int kernel_size, box_size;
|
||||
double epsilon;
|
||||
std::string variogram_type, variogram_para;
|
||||
if (kernel_str == "NULL") kernel_size = 200;
|
||||
else gctl::str2type(kernel_str, kernel_size);
|
||||
|
||||
if (eps_str == "NULL") epsilon = 1e-10;
|
||||
else gctl::str2type(eps_str, epsilon);
|
||||
|
||||
if (box_str == "NULL") box_size = 10;
|
||||
else gctl::str2type(box_str, box_size);
|
||||
|
||||
if (variogram_str == "NULL") variogram_type = "spherical";
|
||||
else variogram_type = variogram_str;
|
||||
|
||||
if (para_str == "NULL") variogram_para = "0.1/10/100";
|
||||
else variogram_para = para_str;
|
||||
|
||||
if (col_str == "NULL") col_str = "0,1,2";
|
||||
|
||||
gctl::text_descriptor desc;
|
||||
if (head_str != "NULL") gctl::str2type(head_str, desc.head_num_);
|
||||
if (ant_str != "NULL") gctl::str2type(ant_str, desc.att_sym_);
|
||||
if (del_str != "NULL") gctl::str2type(del_str, desc.delimiter_);
|
||||
|
||||
LKI ki;
|
||||
ki.Routine(in_name, tar_name, out_name, col_str, desc, kernel_size, box_size, epsilon, variogram_type, variogram_para);
|
||||
return 0;
|
||||
}
|
||||
catch (std::exception &e)
|
||||
{
|
||||
GCTL_ShowWhatError(e.what(), GCTL_ERROR_ERROR, 0, 0, 0);
|
||||
}
|
357
lki/lki.h
Normal file
357
lki/lki.h
Normal file
@@ -0,0 +1,357 @@
|
||||
/********************************************************
|
||||
* ██████╗ ██████╗████████╗██╗
|
||||
* ██╔════╝ ██╔════╝╚══██╔══╝██║
|
||||
* ██║ ███╗██║ ██║ ██║
|
||||
* ██║ ██║██║ ██║ ██║
|
||||
* ╚██████╔╝╚██████╗ ██║ ███████╗
|
||||
* ╚═════╝ ╚═════╝ ╚═╝ ╚══════╝
|
||||
* 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.
|
||||
******************************************************/
|
||||
|
||||
// linear conjugate gradient solver library
|
||||
//#include "lcg/solver.h"
|
||||
// GCTL library
|
||||
#include "gctl/core.h"
|
||||
#include "gctl/algorithm.h"
|
||||
#include "gctl/geometry.h"
|
||||
#include "gctl/io.h"
|
||||
#include "gctl/utility.h"
|
||||
#include "gctl/optimization.h"
|
||||
|
||||
#if defined _WINDOWS || __WIN32__
|
||||
#include "io.h"
|
||||
// Test for file existence
|
||||
#define F_OK 0
|
||||
#endif
|
||||
|
||||
class LKI : public gctl::lcg_solver
|
||||
{
|
||||
public:
|
||||
LKI(){}
|
||||
virtual ~LKI(){}
|
||||
void LCG_Ax(const gctl::array<double> &a, gctl::array<double> &b);
|
||||
void LCG_Mx(const gctl::array<double> &a, gctl::array<double> &b){}
|
||||
|
||||
void Routine(std::string inname, std::string tarname, std::string outname, std::string col_str,
|
||||
gctl::text_descriptor &desc, int kernel_size, int box_size, double epsilon,
|
||||
std::string variogram_type, std::string variogram_para);
|
||||
void ReadConstrainNodes(std::string filename, gctl::text_descriptor &desc);
|
||||
void WriteTargetNodes(std::string filename, const gctl::text_descriptor &desc);
|
||||
void InitTargetNodes(std::string para);
|
||||
void CalKernel();
|
||||
void CalKernel(const gctl::point3dc &tar_node);
|
||||
void set_kernel_size(unsigned int k){MatSize = k+1;}
|
||||
|
||||
public:
|
||||
gctl::array<gctl::point3dc> ConsNodes;
|
||||
gctl::array<gctl::point3dc> TargNodes;
|
||||
std::vector<gctl::point3dc*> LocalNodes;
|
||||
gctl::boxes2d<gctl::point3dc> PntBoxes;
|
||||
|
||||
int ConsSize, MatSize;
|
||||
gctl::matrix<double> Kernel;
|
||||
gctl::array<double> Wgts;
|
||||
gctl::array<double> MidPdt;
|
||||
gctl::array<double> B;
|
||||
|
||||
//gctl::array<double> GK;
|
||||
//gctl::array<double> DK;
|
||||
//gctl::array<double> ADK;
|
||||
|
||||
gctl::variogram_model VarModel;
|
||||
gctl::variogram_type_e VarType;
|
||||
|
||||
gctl::getoption gopt;
|
||||
|
||||
gctl::_1i_vector col_index;
|
||||
};
|
||||
|
||||
void LKI::Routine(std::string inname, std::string tarname, std::string outname, std::string col_str, gctl::text_descriptor &desc,
|
||||
int kernel_size, int box_size, double epsilon, std::string variogram_type, std::string variogram_para)
|
||||
{
|
||||
if (variogram_type == "spherical") VarType = gctl::SPHERICAL;
|
||||
else if (variogram_type == "exponential") VarType = gctl::EXPONENTIAL;
|
||||
else if (variogram_type == "gaussian") VarType = gctl::GAUSSIAN;
|
||||
else if (variogram_type == "wave") VarType = gctl::WAVE;
|
||||
else if (variogram_type == "rational") VarType = gctl::RATIONAL_Q;
|
||||
else if (variogram_type == "circular") VarType = gctl::CIRCULAR;
|
||||
else if (variogram_type == "linear") VarType = gctl::LINEAR;
|
||||
else throw gctl::runtime_error("Invalid variogram function type. From LKI::Routine(...)");
|
||||
|
||||
gctl::parse_string_to_value(variogram_para, '/', true, VarModel.nugget, VarModel.sill, VarModel.range);
|
||||
if (!VarModel.valid()) throw gctl::runtime_error("Invalid variogram parameters. From LKI::Routine(...)");
|
||||
|
||||
gctl::parse_string_to_vector(col_str, ',', col_index);
|
||||
if (col_index.size() < 3)
|
||||
{
|
||||
throw gctl::runtime_error("Invalid column index. From LKI::Routine(...)");
|
||||
}
|
||||
|
||||
ReadConstrainNodes(inname, desc);
|
||||
InitTargetNodes(tarname);
|
||||
|
||||
unsigned int k_size = kernel_size;
|
||||
if (k_size <= 1)
|
||||
{
|
||||
throw gctl::runtime_error("Invalid local size. From LBSI::Routine(...)");
|
||||
}
|
||||
|
||||
if (k_size >= ConsNodes.size())
|
||||
{
|
||||
GCTL_ShowWhatError("The local size is equal to or bigger than the input node's size. Reduced to the global algorithm.", GCTL_WARNING_ERROR, 0, 0, 0);
|
||||
|
||||
ConsSize = ConsNodes.size();
|
||||
MatSize = ConsSize + 1;
|
||||
Kernel.resize(MatSize, MatSize);
|
||||
Wgts.resize(MatSize);
|
||||
B.resize(MatSize);
|
||||
|
||||
CalKernel();
|
||||
|
||||
gctl::lu lu_k(Kernel);
|
||||
lu_k.decompose();
|
||||
|
||||
double dist, sum;
|
||||
for (int i = 0; i < TargNodes.size(); ++i)
|
||||
{
|
||||
for (int j = 0; j < ConsSize; ++j)
|
||||
{
|
||||
dist = sqrt((ConsNodes[j].x - TargNodes[i].x)*(ConsNodes[j].x - TargNodes[i].x)
|
||||
+ (ConsNodes[j].y - TargNodes[i].y)*(ConsNodes[j].y - TargNodes[i].y));
|
||||
|
||||
B[j] = variogram(dist, VarModel, VarType);
|
||||
}
|
||||
B[ConsSize] = 1.0;
|
||||
|
||||
lu_k.solve(B, Wgts);
|
||||
|
||||
sum = 0.0;
|
||||
for (int j = 0; j < ConsSize; ++j)
|
||||
{
|
||||
sum += ConsNodes[j].z*Wgts[j];
|
||||
}
|
||||
TargNodes[i].z = sum;
|
||||
}
|
||||
|
||||
WriteTargetNodes(outname, desc);
|
||||
return;
|
||||
}
|
||||
|
||||
set_kernel_size(k_size);
|
||||
LocalNodes.resize(MatSize-1);
|
||||
Kernel.resize(MatSize, MatSize);
|
||||
Wgts.resize(MatSize);
|
||||
MidPdt.resize(MatSize);
|
||||
B.resize(MatSize);
|
||||
//GK.resize(MatSize);
|
||||
//DK.resize(MatSize);
|
||||
//ADK.resize(MatSize);
|
||||
|
||||
gctl::lcg_para my_para = default_lcg_para();
|
||||
//my_para.max_iterations = 1000;
|
||||
my_para.epsilon = epsilon;
|
||||
set_lcg_para(my_para);
|
||||
|
||||
gctl::array<double> xs(ConsNodes.size());
|
||||
gctl::array<double> ys(ConsNodes.size());
|
||||
for (int i = 0; i < ConsNodes.size(); ++i)
|
||||
{
|
||||
xs[i] = ConsNodes[i].x;
|
||||
ys[i] = ConsNodes[i].y;
|
||||
}
|
||||
|
||||
PntBoxes.init(xs, ys, ConsNodes, box_size, box_size);
|
||||
|
||||
double dist, sum;
|
||||
gctl::progress_bar bar(TargNodes.size());
|
||||
for (int i = 0; i < TargNodes.size(); ++i)
|
||||
{
|
||||
bar.progressed(i);
|
||||
|
||||
Kernel.assign_all(0.0);
|
||||
Wgts.assign_all(0.0);
|
||||
|
||||
CalKernel(TargNodes[i]);
|
||||
|
||||
// run the inversion process in factory mode
|
||||
//lcg(_AxProduct, nullptr, Wgts.get(), B.get(), MatSize, &my_para, this, GK.get(), DK.get(), ADK.get());
|
||||
lcg(Wgts, B);
|
||||
|
||||
sum = 0.0;
|
||||
for (int j = 0; j < MatSize-1; ++j)
|
||||
{
|
||||
sum += LocalNodes[j]->z*Wgts[j];
|
||||
}
|
||||
TargNodes[i].z = sum;
|
||||
}
|
||||
|
||||
WriteTargetNodes(outname, desc);
|
||||
return;
|
||||
}
|
||||
|
||||
void LKI::ReadConstrainNodes(std::string filename, gctl::text_descriptor &desc)
|
||||
{
|
||||
//read_text2array(filename, ConsNodes);
|
||||
|
||||
gctl::_2d_vector table_data;
|
||||
gctl::read_text2vector2d(filename, table_data, desc);
|
||||
|
||||
if (table_data.size() <= 1)
|
||||
{
|
||||
throw gctl::runtime_error("Not enough constraint points. From LBSI::ReadConstrainNodes(...)");
|
||||
}
|
||||
|
||||
ConsNodes.resize(table_data.size());
|
||||
for (int i = 0; i < table_data.size(); ++i)
|
||||
{
|
||||
ConsNodes[i].x = table_data[i][col_index[0]];
|
||||
ConsNodes[i].y = table_data[i][col_index[1]];
|
||||
ConsNodes[i].z = table_data[i][col_index[2]];
|
||||
}
|
||||
|
||||
gctl::destroy_vector(table_data);
|
||||
return;
|
||||
}
|
||||
|
||||
void LKI::WriteTargetNodes(std::string filename, const gctl::text_descriptor &desc)
|
||||
{
|
||||
save_array2text(filename, TargNodes, desc);
|
||||
return;
|
||||
}
|
||||
|
||||
void LKI::InitTargetNodes(std::string para)
|
||||
{
|
||||
// try to use the para as a file name
|
||||
if (access(para.c_str(), F_OK) != -1)
|
||||
{
|
||||
std::vector<gctl::point2dc> tmp_vec;
|
||||
gctl::read_text2vector(para, tmp_vec);
|
||||
|
||||
TargNodes.resize(tmp_vec.size());
|
||||
for (int i = 0; i < tmp_vec.size(); ++i)
|
||||
{
|
||||
TargNodes[i].x = tmp_vec[i].x;
|
||||
TargNodes[i].y = tmp_vec[i].y;
|
||||
TargNodes[i].z = 0.0;
|
||||
}
|
||||
|
||||
gctl::destroy_vector(tmp_vec);
|
||||
return;
|
||||
}
|
||||
|
||||
double dx, dy, xmin, xmax, ymin, ymax, ele = 0.0;
|
||||
gctl::parse_string_to_value(para, '/', true, xmin, dx, xmax, ymin, dy, ymax);
|
||||
gctl::get_grid_point3c(TargNodes, xmin, xmax, ymin, ymax, dx, dy, ele);
|
||||
return;
|
||||
}
|
||||
|
||||
void LKI::CalKernel()
|
||||
{
|
||||
// 计算出所有成对的协方差函数值 成对的
|
||||
double dist;
|
||||
for (int i = 0; i < ConsSize; ++i)
|
||||
{
|
||||
Kernel[i][i] = variogram(0.0, VarModel, VarType);
|
||||
for (int j = i+1; j < ConsSize; ++j)
|
||||
{
|
||||
dist = sqrt((ConsNodes[i].x - ConsNodes[j].x)*(ConsNodes[i].x - ConsNodes[j].x)
|
||||
+ (ConsNodes[i].y - ConsNodes[j].y)*(ConsNodes[i].y - ConsNodes[j].y));
|
||||
|
||||
Kernel[i][j] = Kernel[j][i] = variogram(dist, VarModel, VarType);
|
||||
}
|
||||
}
|
||||
|
||||
// 最后一行系数全为1
|
||||
for (int i = 0; i < ConsSize; ++i)
|
||||
{
|
||||
Kernel[ConsSize][i] = Kernel[i][ConsSize] = 1.0;
|
||||
}
|
||||
|
||||
Kernel[ConsSize][ConsSize] = 0.0;
|
||||
return;
|
||||
}
|
||||
|
||||
void LKI::CalKernel(const gctl::point3dc &tar_node)
|
||||
{
|
||||
// 找出距离tar_node最近的一组控制点
|
||||
PntBoxes.get_by_number(tar_node.x, tar_node.y, MatSize-1, LocalNodes);
|
||||
|
||||
// 计算出所有成对的格林函数值
|
||||
double dist;
|
||||
for (int j = 0; j < MatSize-1; ++j)
|
||||
{
|
||||
Kernel[j][j] = variogram(0.0, VarModel, VarType);
|
||||
for (int k = j+1; k < MatSize-1; ++k)
|
||||
{
|
||||
dist = sqrt((LocalNodes[j]->x - LocalNodes[k]->x)*(LocalNodes[j]->x - LocalNodes[k]->x)
|
||||
+ (LocalNodes[j]->y - LocalNodes[k]->y)*(LocalNodes[j]->y - LocalNodes[k]->y));
|
||||
|
||||
Kernel[j][k] = Kernel[k][j] = variogram(dist, VarModel, VarType);
|
||||
}
|
||||
}
|
||||
|
||||
// 最后一行系数全为1
|
||||
for (int i = 0; i < MatSize-1; ++i)
|
||||
{
|
||||
Kernel[MatSize-1][i] = Kernel[i][MatSize-1] = 1.0;
|
||||
}
|
||||
Kernel[MatSize-1][MatSize-1] = 0.0;
|
||||
|
||||
for (int j = 0; j < MatSize-1; ++j)
|
||||
{
|
||||
dist = sqrt((LocalNodes[j]->x - tar_node.x)*(LocalNodes[j]->x - tar_node.x)
|
||||
+ (LocalNodes[j]->y - tar_node.y)*(LocalNodes[j]->y - tar_node.y));
|
||||
|
||||
MidPdt[j] = variogram(dist, VarModel, VarType);
|
||||
}
|
||||
MidPdt[MatSize-1] = 1.0;
|
||||
|
||||
for (int i = 0; i < MatSize; ++i)
|
||||
{
|
||||
B[i] = 0;
|
||||
for (int j = 0; j < MatSize; ++j)
|
||||
{
|
||||
B[i] += Kernel[j][i] * MidPdt[j];
|
||||
}
|
||||
}
|
||||
return;
|
||||
}
|
||||
|
||||
void LKI::LCG_Ax(const gctl::array<double> &a, gctl::array<double> &b)
|
||||
{
|
||||
for (int i = 0; i < MatSize; ++i)
|
||||
{
|
||||
MidPdt[i] = 0;
|
||||
for (int j = 0; j < MatSize; ++j)
|
||||
{
|
||||
MidPdt[i] += a[j] * Kernel[i][j];
|
||||
}
|
||||
}
|
||||
|
||||
for (int i = 0; i < MatSize; ++i)
|
||||
{
|
||||
b[i] = 0;
|
||||
for (int j = 0; j < MatSize; ++j)
|
||||
{
|
||||
b[i] += MidPdt[j] * Kernel[j][i];
|
||||
}
|
||||
}
|
||||
return;
|
||||
}
|
Reference in New Issue
Block a user