gstlearn  1.0.0
CCC
dbtools.cpp File Reference
#include "geoslib_f.h"
#include "geoslib_old_f.h"
#include "geoslib_f_private.h"
#include "Enum/EJustify.hpp"
#include "Mesh/MeshETurbo.hpp"
#include "LinearOp/ShiftOpCs.hpp"
#include "LinearOp/PrecisionOp.hpp"
#include "LinearOp/ProjMatrix.hpp"
#include "LinearOp/OptimCostColored.hpp"
#include "Stats/Classical.hpp"
#include "Morpho/Morpho.hpp"
#include "Covariances/CovAniso.hpp"
#include "Model/ANoStat.hpp"
#include "Model/NoStatArray.hpp"
#include "Model/Model.hpp"
#include "Model/CovInternal.hpp"
#include "Db/Db.hpp"
#include "Basic/Law.hpp"
#include "Basic/NamingConvention.hpp"
#include "Basic/Utilities.hpp"
#include "Basic/String.hpp"
#include "Basic/File.hpp"
#include "Basic/VectorHelper.hpp"
#include "Basic/PolyLine2D.hpp"
#include "Basic/OptDbg.hpp"
#include "Polygon/Polygons.hpp"
#include "Skin/ISkinFunctions.hpp"
#include "Skin/Skin.hpp"
#include "Geometry/GeometryHelper.hpp"
#include "Tree/Ball.hpp"
#include "Tree/KNN.hpp"
#include <math.h>
#include <string.h>

Classes

class  LocalSkin
 

Functions

int migrate_grid_to_coor (const DbGrid *db_grid, int iatt, const VectorVectorDouble &coords, VectorDouble &tab)
 
int expand_point_to_coor (const Db *db1, int iatt, const VectorVectorDouble &coords, VectorDouble &tab)
 
int db_tool_duplicate (Db *db1, Db *db2, bool flag_same, bool verbose, int opt_code, double tolcode, const VectorDouble &dist, VectorDouble &sel)
 
int db_duplicate (Db *db, bool verbose, const VectorDouble &dist, int opt_code, double tolcode, const NamingConvention &namconv)
 
int surface (Db *db_point, DbGrid *db_grid, int, double dlim, double *dtab, double *gtab)
 
int db_edit (Db *db, int *flag_valid)
 
int db_normalize (Db *db, const char *oper, int ncol, int *cols, double center, double stdv)
 
int db_grid_fill (DbGrid *dbgrid, int mode, int seed, int radius, bool verbose, const NamingConvention &namconv)
 
int _db_category (Db *db, int iatt, const VectorDouble &mini, const VectorDouble &maxi, const VectorBool &incmini, const VectorBool &incmaxi, const NamingConvention &namconv)
 
int _db_indicator (Db *db, int iatt, int flag_indic, const VectorDouble &mini, const VectorDouble &maxi, const VectorBool &incmini, const VectorBool &incmaxi, bool flagBelow, bool flagAbove, const NamingConvention &namconv)
 
VectorDouble _db_limits_statistics (Db *db, int iatt, const VectorDouble &mini, const VectorDouble &maxi, const VectorBool &incmini, const VectorBool &incmaxi, int optionStat, bool flagBelow, bool flagAbove)
 
int interpolate_variable_to_point (DbGrid *db_grid, int iatt, int np, double *xp, double *yp, double *zp, double *tab)
 
void ut_trace_discretize (int nseg, double *trace, double disc, int *np_arg, double **xp_arg, double **yp_arg, double **dd_arg, double **del_arg, double *dist_arg)
 
void ut_trace_sample (Db *db, const ELoc &ptype, int np, double *xp, double *yp, double *dd, double radius, int *ns_arg, double **xs_arg, double **ys_arg, int **rks_arg, int **lys_arg, int **typ_arg)
 
int manage_external_info (int mode, const ELoc &locatorType, Db *dbin, Db *dbout, int *istart)
 
int manage_nostat_info (int mode, Model *model, Db *dbin, Db *dbout)
 
int db_center_point_to_grid (Db *db_point, DbGrid *db_grid, double eps_random)
 
DbGriddb_grid_sample (DbGrid *dbin, const VectorInt &nmult)
 
int st_next_sample (int ip0_init, int np, const VectorInt &order, const VectorDouble &xtab, double xtarget)
 
int expand_point_to_grid (Db *db_point, DbGrid *db_grid, int iatt, int iatt_time, int iatt_angle, int iatt_scaleu, int iatt_scalev, int iatt_scalew, int flag_index, int distType, const VectorDouble &dmax, VectorDouble &tab)
 
int db_compositional_transform (Db *db, int verbose, int mode, int type, int number, int *iatt_in, int *iatt_out, int *numout)
 
int points_to_block (Db *dbpoint, DbGrid *dbgrid, int option, int flag_size, int iatt_time, int iatt_size, int iatt_angle, int iatt_scaleu, int iatt_scalev, int iatt_scalew)
 
int db_resind (Db *db, int ivar, int ncut, double *zcut)
 
int db_gradient_components (DbGrid *dbgrid)
 
int db_streamline (DbGrid *dbgrid, Db *dbpoint, int niter, double step, int flag_norm, int use_grad, int save_grad, int *nbline_loc, int *npline_loc, double **line_loc)
 
int db_model_nostat (Db *db, Model *model, int icov, const NamingConvention &namconv)
 
int db_smooth_vpc (DbGrid *db, int width, double range)
 
Dbdb_extract (Db *db, int *ranks)
 
Dbdb_regularize (Db *db, DbGrid *dbgrid, int flag_center)
 
double * db_grid_sampling (DbGrid *dbgrid, double *x1, double *x2, int ndisc, int ncut, double *cuts, int *nval_ret)
 
int db_grid2point_sampling (DbGrid *dbgrid, int nvar, int *vars, int *npacks, int npcell, int nmini, int *nech_ret, double **coor_ret, double **data_ret)
 
Dbdb_point_init (int nech, const VectorDouble &coormin, const VectorDouble &coormax, DbGrid *dbgrid, bool flag_exact, bool flag_repulsion, double range, double beta, double extend, int seed, int flag_add_rank)
 
int db_grid1D_fill (DbGrid *dbgrid, int mode, int seed, const NamingConvention &namconv)
 
int _migrate (Db *db1, Db *db2, int iatt1, int iatt2, int distType, const VectorDouble &dmax, bool flag_fill, bool flag_inter, bool flag_ball)
 
int db_proportion_estimate (Db *dbin, DbGrid *dbout, Model *model, int niter, bool verbose, const NamingConvention &namconv)
 

Function Documentation

◆ _db_category()

int _db_category ( Db db,
int  iatt,
const VectorDouble mini,
const VectorDouble maxi,
const VectorBool incmini,
const VectorBool incmaxi,
const NamingConvention namconv 
)

Convert the contents of the ivar-th continuous variable into a categorical array

Returns
Error return code
Parameters
[in]dbDb structure
[in]iattRank of the attribute
[in]miniArray containing the minima per class
[in]maxiArray containing the maxima per class
[in]incminiArray containing the inclusive flag for minima per class
[in]incmaxiArray containing the inclusive flag for maxima per class
[in]namconvNaming convention
Remarks
When the array mini or maxi are not provided, then
mini[iclass] = int(iclass-1)
maxi[iclass[ = int(iclass)

◆ _db_indicator()

int _db_indicator ( Db db,
int  iatt,
int  flag_indic,
const VectorDouble mini,
const VectorDouble maxi,
const VectorBool incmini,
const VectorBool incmaxi,
bool  flagBelow,
bool  flagAbove,
const NamingConvention namconv 
)

Create indicator variables

Returns
Error returned code
Parameters
[in]dbDb structure
[in]iattRank of the target variable
[in]flag_indicType of variable(s) to be stored: 1 the indicator variable 0 the mean variable per class
[in]miniArray containing the minima per class
[in]maxiArray containing the maxima per class
[in]incminiArray containing the inclusive flag for minima per class
[in]incmaxiArray containing the inclusive flag for maxima per class
[in]flagBelowIf True, consider the values below the lower bound
[in]flagAboveIf True, consider the values above the upper bound
[in]namconvNaming convention
Remarks
When both arrays mini and maxi are not provided, then:
mini[iclass] = iclass-1
maxi[iclass[ = iclass

◆ _db_limits_statistics()

VectorDouble _db_limits_statistics ( Db db,
int  iatt,
const VectorDouble mini,
const VectorDouble maxi,
const VectorBool incmini,
const VectorBool incmaxi,
int  optionStat,
bool  flagBelow,
bool  flagAbove 
)

Calculate statistics per class

Returns
Array of statistics
Parameters
[in]dbDb structure
[in]iattRank of the target variable
[in]miniArray containing the minima per class
[in]maxiArray containing the maxima per class
[in]incminiArray containing the inclusive flag for minima per class
[in]incmaxiArray containing the inclusive flag for maxima per class
[in]optionStat1 for Proportions; 2 for Mean
[in]flagBelowIf True, consider the values below the lower bound
[in]flagAboveIf True, consider the values above the upper bound

◆ _migrate()

int _migrate ( Db db1,
Db db2,
int  iatt1,
int  iatt2,
int  distType,
const VectorDouble dmax,
bool  flag_fill,
bool  flag_inter,
bool  flag_ball 
)

Migrates a variable from one Db to another one

Returns
Error return code
Parameters
[in]db1descriptor of the input Db
[in]db2descriptor of the output Db
[in]iatt1Attribute in Db1 to be migrated
[in]iatt2Attribute in Db2 where the result must be stored
[in]distTypeType of distance for calculating maximum distance 1 for L1 and 2 for L2 distance
[in]dmaxArray of maximum distances (optional)
[in]flag_fillFilling option
[in]flag_interInterpolation
[in]flag_ballUse the BallTree sort for speeding up calculations

◆ db_center_point_to_grid()

int db_center_point_to_grid ( Db db_point,
DbGrid db_grid,
double  eps_random 
)

Centers the samples of a Db to the center of blocks of a grid Db

Returns
Error return code
Parameters
[in]db_pointdescriptor of the point parameters
[in]db_griddescriptor of the grid parameters
[in]eps_randomRandomisation Epsilon
Remarks
The argument 'eps_random' allows perturbating the centered
coordinate so that it does not lie exactly on the node.
This possibility makes sense in order to identify centered data
from data actually located on the grid center (before migration)
The perturbation is calculated as DX(i) * eps

◆ db_compositional_transform()

int db_compositional_transform ( Db db,
int  verbose,
int  mode,
int  type,
int  number,
int *  iatt_in,
int *  iatt_out,
int *  numout 
)

Translate a set of compositional variables into auxiliary variables

Returns
Error return code
Parameters
[in]dbDb characteristics
[in]verbose1 for a Verbose option
[in]mode1 for Forward and -1 for Backward transformation
[in]typeType of conversion
  • 0 : Simple transformation in proportion
  • 1 : Additive logratio
  • 2 : Centered logratio
  • 3 : Isometric logratio
[in]numberNumber of input attributes
[in]iatt_inArray of the input attribute
[in]iatt_outArray of the output attribute
[out]numoutNumber of variables in output
Remarks
The additive and the isometric logratio transformations
transform N+1 compositional variables into N elements (Forward)
and from N transformed elements into N+1 compositional variables
(Backwards).
The zero-values are replaced by a conventional small value
This is defined by a variable that can be corrected using
the keypair facility with the keyword 'CompositionalEps'
The arguments 'iatt_in' and 'iatt_out' can coincide.
Outlier Detection for Compositional Data Using Robust Methods
Math Geosciences (2008) 40: 233-248

◆ db_duplicate()

int db_duplicate ( Db db,
bool  verbose,
const VectorDouble dist,
int  opt_code,
double  tolcode,
const NamingConvention namconv 
)

Look for duplicates within a Db

Returns
Error return code
Parameters
[in]dbDb Structure
[in]verboseTrue for verbose output
[in]distArray of the minimum distance
[in]opt_codecode selection option (if code is defined)
  • 0 : no use of the code selection
  • 1 : codes must be close enough
  • 2 : codes must be different
[in]tolcodeCode tolerance
[in]namconvNaming convention

◆ db_edit()

int db_edit ( Db db,
int *  flag_valid 
)

Edit the Data Base Db

Returns
Error return code
Parameters
[in]dbDb descriptor
[out]flag_valid1 for 'stop' and 0 for 'quit'

◆ db_extract()

Db* db_extract ( Db db,
int *  ranks 
)

Extract a new Db from an old Db using a sample seleciton vector

Returns
Pointer to the newly created Db
Parameters
[in]dbInitial Db
[in]ranksVector of selected ranks
Remarks
The 'ranks' array is dimensionned to the number of samples of 'db'
and ranks[i]=1 if the sample 'i' must be kept and 0 otherwise.
All the variables contained in the input Db are copied into the
output Db
This operations is limited to non-grid Db
Possible selection in the input Db is taken into account

◆ db_gradient_components()

int db_gradient_components ( DbGrid dbgrid)

Calculate the gradient over a grid

Returns
Rank of the newly created variables (or -1)
Parameters
[in]dbgridDb structure (grid organized)

◆ db_grid1D_fill()

int db_grid1D_fill ( DbGrid dbgrid,
int  mode,
int  seed,
const NamingConvention namconv 
)

Fill an incomplete 1-D grid

Returns
Error returned code
Parameters
[in]dbgridDb grid structure
[in]modeType of interpolation
  • 0 : Linear interpolation
  • 1 : Cubic Spline
[in]seedSeed used for the random number generation
[in]namconvNaming convention

◆ db_grid2point_sampling()

int db_grid2point_sampling ( DbGrid dbgrid,
int  nvar,
int *  vars,
int *  npacks,
int  npcell,
int  nmini,
int *  nech_ret,
double **  coor_ret,
double **  data_ret 
)

Sampling a fine grid in a coarser set of cells

Returns
Error returned code
Parameters
[in]dbgridreference Grid
[in]nvarNumber of variables
[in]varsArray of variable ranks
[in]npacksVector of packing factors
[in]npcellNumber of samples per cell
[in]nminiMinimum number of nodes before drawing
[out]nech_retNumber of selected samples
[out]coor_retArray of coordinates
[out]data_retArray of variables
Remarks
The returned arrays 'coor' and 'data' must be freed by
the calling function

◆ db_grid_fill()

int db_grid_fill ( DbGrid dbgrid,
int  mode,
int  seed,
int  radius,
bool  verbose,
const NamingConvention namconv 
)

Fill an incomplete grid

Returns
Error returned code
Parameters
[in]dbgridDb grid structure
[in]modeType of interpolation
  • 0 : Moving average
  • 1 : Inverse squared distance
  • 2 : Interpolation by a linear plane
  • 3 : Distance to the initial grains
[in]seedSeed used for the random number generation
[in]radiusRadius of the neighborhood
[in]verboseVerbose flag
[in]namconvNaming convention

◆ db_grid_sample()

DbGrid* db_grid_sample ( DbGrid dbin,
const VectorInt nmult 
)

Sample a grid into a finer subgrid (all variables)

Returns
Error return code
Parameters
[in]dbinDescriptor of the grid parameters
[in]nmultArray of multiplicity coefficients

◆ db_grid_sampling()

double* db_grid_sampling ( DbGrid dbgrid,
double *  x1,
double *  x2,
int  ndisc,
int  ncut,
double *  cuts,
int *  nval_ret 
)

Sampling vertices within a Grid between two points

Returns
Array of sampled vertices
Parameters
[in]dbgridreference Grid
[in]x1Array giving the coordinates of the first point
[in]x2Array giving the coordinates of the second point
[in]ndiscNumber of discretized points in the segment
[in]ncutNumber of cutoffs
[in]cutsArray of cutoffs
[out]nval_retNumber of samples in the output array
Remarks
This function considers the segment [x1,x2] and subdivises it
into 'ndisc' intervals. The endpoints of each interval correspond
to two points in the space
At each endpoint, the target variable is interpolated from the grid
If the target variable values cross a cutoff, the coordinates of
the intersection are calculated.
The program returns the list of all these intersection coordinates

◆ db_model_nostat()

int db_model_nostat ( Db db,
Model model,
int  icov,
const NamingConvention namconv 
)

Calculate and store new variables in the Db which contain the non-stationary Model component

Returns
Distance value
Parameters
[in]dbDb structure
[in]modelModel structure
[in]icovRank of the Basic structure
[in]namconvNaming convention
Remarks
This procedure automatically creates several fields:
ndim fields for storing the ranges
ndim fields for storing the angles
1 field for storing the sill

◆ db_normalize()

int db_normalize ( Db db,
const char *  oper,
int  ncol,
int *  cols,
double  center,
double  stdv 
)

Normalize a set of variables

Returns
Error returned code
Parameters
[in]dbDb structure
[in]operName of the operator
  • "mean" : Normalize the mean to the input value
  • "stdv" : Normalize the st. dev. to the input value
  • "scal" : Normalize the mean and st. dev.
  • "prop" : Normalize the variables to proportions
[in]ncolNumber of variables
[in]colsRanks of the variables
[in]centerTheoretical Mean value
[in]stdvTheoretical Standard Deviation value

◆ db_point_init()

Db* db_point_init ( int  nech,
const VectorDouble coormin,
const VectorDouble coormax,
DbGrid dbgrid,
bool  flag_exact,
bool  flag_repulsion,
double  range,
double  beta,
double  extend,
int  seed,
int  flag_add_rank 
)

Create a new Data Base with points generated at random

Returns
Pointer for the new Db structure
Parameters
[in]nechExpected number of samples
[in]coorminVector of lower coordinates
[in]coormaxVector of upper coordinates
[in]dbgridDescriptor of the Db grid parameters
[in]flag_exactTrue if the number of samples is dran from Poisson
[in]flag_repulsionUse repulsion (need: 'range' and 'beta')
[in]rangeRepulsion range
[in]betaBending coefficient
[in]extendExtension of the bounding box (when positive)
[in]seedSeed for the random number generator
[in]flag_add_rank1 if the Rank must be generated in the output Db
Remarks
Arguments 'extend' is only valid when 'dbgrid' is not defined

◆ db_proportion_estimate()

int db_proportion_estimate ( Db dbin,
DbGrid dbout,
Model model,
int  niter,
bool  verbose,
const NamingConvention namconv 
)

Standard Kriging

Returns
Error return code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel structure
[in]niterNumber of iterations
[in]verboseVerbose flag
[in]namconvNaming convention
Remarks
The procedure uses the FIRST covariance of the Model
to describe the spatial structure

◆ db_regularize()

Db* db_regularize ( Db db,
DbGrid dbgrid,
int  flag_center 
)

Regularize variables along vertical wells

Returns
Pointer to the newly created Db
Parameters
[in]dbInitial Db
[in]dbgridReference Grid
[in]flag_centerWhen TRUE, the sample is centered in the layer to which it belongs
Remarks
This function requires the input well ('db') and the grid to be
defined in space >= 3D
It requires a CODE variable to be defined in the input 'db'
This function regularizes all the variables marked with a Z-locator
This function takes a sample into account only if isotopic

◆ db_resind()

int db_resind ( Db db,
int  ivar,
int  ncut,
double *  zcut 
)

Create indicator residual variables

Returns
Error returned code
Parameters
[in]dbDb structure
[in]ivarIndex of the target variable
[in]ncutNumber of cutoffs
[in]zcutArray containing the cutoffs
Remarks
The array 'zcut' must be provided in increasing order

◆ db_smooth_vpc()

int db_smooth_vpc ( DbGrid db,
int  width,
double  range 
)

Smooth out the VPC

Returns
Distance value
Parameters
[in]db3-D Db structure containing the VPCs
[in]widthWidth of the Filter
[in]rangeRange of the Gaussian Weighting Function
Remarks
Work is performed IN PLACE

◆ db_streamline()

int db_streamline ( DbGrid dbgrid,
Db dbpoint,
int  niter,
double  step,
int  flag_norm,
int  use_grad,
int  save_grad,
int *  nbline_loc,
int *  npline_loc,
double **  line_loc 
)

Calculate the streamlines

Returns
Error return code
Parameters
[in]dbgridDb structure (grid organized)
[in]dbpointDb structure for control points
[in]niterMaximum number of iterations
[in]stepProgress step
[in]flag_norm1 if the gradients must be normalized
[in]use_grad1 if the existing gradients must be used 0 the gradients must be calculated here
[in]save_grad1 save the gradients generated in this function 0 delete gradients when calculated here
[out]nbline_locNumber of streamline steps
[out]npline_locNumber of values per line vertex
[out]line_locArray of streamline steps (Dimension: 5 * nbline)
Remarks
The returned array 'line_loc' must be freed by the calling function
Use get_keypone("Streamline_Skip",1) to define the skipping ratio

◆ db_tool_duplicate()

int db_tool_duplicate ( Db db1,
Db db2,
bool  flag_same,
bool  verbose,
int  opt_code,
double  tolcode,
const VectorDouble dist,
VectorDouble sel 
)

Look for duplicates

Returns
Error return code
Parameters
[in]db1First Db
[in]db2Second Db
[in]flag_sameTrue if the two Db files are the same
[in]verboseTrue for verbose output
[in]opt_codecode selection option (if code is defined)
  • 0 : no use of the code selection
  • 1 : codes must be close enough
  • 2 : codes must be different
[in]tolcodeCode tolerance
[in]distArray of the minimum distance (or NULL)
[out]selArray containing the selection

◆ expand_point_to_coor()

int expand_point_to_coor ( const Db db1,
int  iatt,
const VectorVectorDouble coords,
VectorDouble tab 
)

Expands a variable from one point Db into a variable at points defined by coordinate vectors (maximum 3D)

Returns
Error return code
Parameters
[in]db1descriptor of the input parameters
[in]iattrank of the input attribute
[in]coordsArray of coordinates
[out]tabOutput array (Dimension: number of discretized points)

◆ expand_point_to_grid()

int expand_point_to_grid ( Db db_point,
DbGrid db_grid,
int  iatt,
int  iatt_time,
int  iatt_angle,
int  iatt_scaleu,
int  iatt_scalev,
int  iatt_scalew,
int  flag_index,
int  distType,
const VectorDouble dmax,
VectorDouble tab 
)

Expands a variable from the point structure into a variable in the grid structure

Returns
Error return code
Parameters
[in]db_pointDescriptor of the point parameters
[in]db_gridDescriptor of the grid parameters
[in]iattRank of the point attribute
[in]iatt_timeOptional variable for Time shift
[in]iatt_angleOptional variable for anisotropy angle (around Z)
[in]iatt_scaleuOptional variable for anisotropy scale factor (U)
[in]iatt_scalevOptional variable for anisotropy scale factor (V)
[in]iatt_scalewOptional variable for anisotropy scale factor (W)
[in]flag_index1 if the Index must be assigned to grid node 0 the 'iatt' attribute is assigned instead
[in]distTypeType of distance for calculating maximum distance 1 for L1 and 2 for L2 distance
[in]dmaxArray of maximum distances (optional)
[out]tabOutput array
Remarks
When a Time Shift is present, this corresponds to Johnson-Mehl
The Time Shift is an optional variable which increases the
distance (the time-to-distance conversion is assumed to be 1)
Only positive Time Shifts are considered

◆ interpolate_variable_to_point()

int interpolate_variable_to_point ( DbGrid db_grid,
int  iatt,
int  np,
double *  xp,
double *  yp,
double *  zp,
double *  tab 
)

Interpolate a variable from a grid Db on discretization points

Parameters
[in]db_gridDescriptor of the grid parameters
[in]iattRank of the attribute in db_grid
[in]npNumber of discretized points
[in]xpArray of first coordinates
[in]ypArray of second coordinates
[in]zpArray of third coordinates
[out]tabOutput array
Remarks
The arguments 'xp', 'yp' and 'zp' must be defined in accordance
with the space dimension in argument 'db_grid'
A point which does not lie between two valuated grid nodes
(in all space dimensions) is always set to FFFF

◆ manage_external_info()

int manage_external_info ( int  mode,
const ELoc &  locatorType,
Db dbin,
Db dbout,
int *  istart 
)

Derive the external information(s) from the Output db (if Grid) to the Input Db

Returns
Error return code
Parameters
[in]mode1 for allocation; -1 for deallocation
[in]locatorTypeType of the pointer (ELoc)
[in]dbinDescriptor of the input Db
[in]dboutDescriptor of the output Db
[in,out]istartAddress of the first allocated external information
Remarks
This function only functions when the Output Db is a grid
However, in case of a Point output Db, this function should not
be used: the external drift functions should already be present
in the output Db.
If this is not the case, an error is issued.

◆ manage_nostat_info()

int manage_nostat_info ( int  mode,
Model model,
Db dbin,
Db dbout 
)

Derive the non-stationary information(s) from the Output db (if Grid) to the Input Db

Returns
Error return code
Parameters
[in]mode1 for allocation; -1 for deallocation
[in]modelDescriptor of the Model
[in]dbinDescriptor of the input Db
[in]dboutDescriptor of the output Db

◆ migrate_grid_to_coor()

int migrate_grid_to_coor ( const DbGrid db_grid,
int  iatt,
const VectorVectorDouble coords,
VectorDouble tab 
)

Migrates a variable from the grid structure into a variable at points defined by coordinate vectors

Returns
Error return code
Parameters
[in]db_griddescriptor of the grid parameters
[in]iattrank of the grid attribute
[in]coordsArray of coordinates (dimension: ndim, np)
[out]tabOutput array (Dimension: number of discretized points)

◆ points_to_block()

int points_to_block ( Db dbpoint,
DbGrid dbgrid,
int  option,
int  flag_size,
int  iatt_time,
int  iatt_size,
int  iatt_angle,
int  iatt_scaleu,
int  iatt_scalev,
int  iatt_scalew 
)

Plunge a set of isolated points within a discretization grid in order to compute the voronoi of the points and derive:

  • the statistics on the volume and perimeter of the cells
  • the edge between cells
Returns
Error return code
Parameters
[in]dbpointDescriptor of the point parameters
[in]dbgridDescriptor of the grid parameters
[in]optionConnectivity option (0 for cross and 1 for block)
[in]flag_sizeWhen 1, the border pixels report the border thickness When 0, the border pixels are painted in 1
[in]iatt_timeAttribute of 'dbpoint'for Time shift (optional)
[in]iatt_sizeAttribute of 'dbpoint' giving size (optional)
[in]iatt_angleOptional variable for anisotropy angle (around Z)
[in]iatt_scaleuOptional variable for anisotropy scale factor (U)
[in]iatt_scalevOptional variable for anisotropy scale factor (V)
[in]iatt_scalewOptional variable for anisotropy scale factor (W)
Remarks
The value of 'flag_index' can be turned on for assigning
the sample index to the grid cell (instead of the 'iatt' value)
using: set_keypair("PTB_flag_index")

◆ st_next_sample()

int st_next_sample ( int  ip0_init,
int  np,
const VectorInt order,
const VectorDouble xtab,
double  xtarget 
)

◆ surface()

int surface ( Db db_point,
DbGrid db_grid,
int  ,
double  dlim,
double *  dtab,
double *  gtab 
)

Calculate the (discretized) surface of influence

Returns
Error returned code
Parameters
[in]db_pointDb containing the data points
[in]db_gridDb containing the discretization grid
[in]dlimMaximum distance (TEST if not defined)
[out]dtabArray containing the surface of influence (Dimension = Number of samples in db_point)
[out]gtabArray containing the surface of influence of the polygon to which it belongs (or TEST) (Dimension = Number of samples in db_grid)

◆ ut_trace_discretize()

void ut_trace_discretize ( int  nseg,
double *  trace,
double  disc,
int *  np_arg,
double **  xp_arg,
double **  yp_arg,
double **  dd_arg,
double **  del_arg,
double *  dist_arg 
)

Generates the discretized points along the trace

Parameters
[in]nsegNumber of vertices along the trace
[in]traceArray defining the trace (Dimension: 2 * nseg)
[in]discDiscretization distance
[out]np_argNumber of discretized points
[out]xp_argArray of first coordinates
[out]yp_argArray of second coordinates
[out]dd_argArray of distances between discretized points
[out]del_argArray of distances between vertices
[out]dist_argTotal distance of the trace

◆ ut_trace_sample()

void ut_trace_sample ( Db db,
const ELoc &  ptype,
int  np,
double *  xp,
double *  yp,
double *  dd,
double  radius,
int *  ns_arg,
double **  xs_arg,
double **  ys_arg,
int **  rks_arg,
int **  lys_arg,
int **  typ_arg 
)

Sample the point Db close to discretized points of the trace

Parameters
[in]dbDb to be sampled
[in]ptypeType of locator
[in]npNumber of discretized points
[in]xpArray of first coordinates
[in]ypArray of second coordinates
[in]ddArray of distances
[in]radiusNeighborhood radius
[out]ns_argNumber of sampled points
[out]xs_argArray of first coordinates of sampled points
[out]ys_argArray of second coordinates of sampled points
[out]rks_argArray of sample indices (starting from 1)
[out]lys_argArray of layer indices of sampled points
[out]typ_argArray of data type 1 for hard data in Z or TIME 2 for lower bound 3 for upper bound