1.4.0
CCC
 
spde.cpp File Reference
#include "geoslib_f_private.h"
#include "geoslib_old_f.h"
#include "Enum/ELoadBy.hpp"
#include "Matrix/MatrixFactory.hpp"
#include "Matrix/MatrixSquareGeneral.hpp"
#include "Matrix/MatrixSparse.hpp"
#include "Matrix/LinkMatrixSparse.hpp"
#include "Matrix/NF_Triplet.hpp"
#include "LinearOp/ProjMatrix.hpp"
#include "Mesh/MeshEStandard.hpp"
#include "Covariances/CovAniso.hpp"
#include "Covariances/ACovAnisoList.hpp"
#include "Covariances/CovFactory.hpp"
#include "Basic/AException.hpp"
#include "Basic/Utilities.hpp"
#include "Basic/Law.hpp"
#include "Basic/MathFunc.hpp"
#include "Basic/String.hpp"
#include "Basic/VectorHelper.hpp"
#include "Db/Db.hpp"
#include "Model/Model.hpp"
#include "Mesh/AMesh.hpp"
#include "Mesh/MeshETurbo.hpp"
#include "Calculators/CalcMigrate.hpp"
#include "Space/SpaceSN.hpp"
#include "Geometry/GeometryHelper.hpp"
#include <math.h>
#include <string.h>

Macros

#define NBLIN_TERMS   10
 
#define SPDE_MAX_NGRF   2
 

Functions

static int st_get_rank (int ivar, int jvar)
 
static void st_matelem_print (int icov)
 
void spde_option_update (SPDE_Option &s_option, const String &triswitch)
 
SPDE_Option spde_option_alloc (void)
 
SPDE_Matelemspde_get_current_matelem (int icov)
 
static void st_title (int flag_igrf, int flag_icov, int rank, const char *title)
 
static Modelst_get_model (void)
 
static int st_get_number_grf (void)
 
static void st_environ_init (void)
 
static cs_MGSst_mgs_manage (int mode, cs_MGS *mgs)
 
void simu_define_func_transf (void(*st_simu_transf)(Db *, int, int, int))
 
void simu_define_func_update (void(*st_simu_update)(Db *, int, int, int))
 
void simu_define_func_scale (void(*st_simu_scale)(Db *, int, int))
 
static bool st_is_model_nugget (void)
 
static CovAnisost_get_nugget (void)
 
static CovAnisost_get_cova (void)
 
static void st_set_model (Model *model)
 
static int st_get_nvs2 (void)
 
static int st_get_filnug (void)
 
static void st_set_filnug (int flag_filnug)
 
static double st_get_isill (int icov, int ivar, int jvar)
 
static void st_clean_Bhetero (void)
 
static void st_clean_Bnugget (void)
 
static void st_print_status (int auth)
 
static void st_qchol_filter (const char *title, int auth)
 
static void st_qchol_print (const char *title, QChol *QC)
 
static MatrixSparsest_extract_Q_from_Q (MatrixSparse *Q_in, int row_auth, int col_auth)
 
static QCholst_extract_QC_from_Q (const char *title, QChol *QC_in, int row_auth, int col_auth)
 
static QSimust_qsimu_manage (int mode, QSimu *qsimu)
 
QCholqchol_manage (int mode, QChol *QC)
 
static double st_get_nugget_sill (int ivar, int jvar)
 
static double st_get_cova_sill (int ivar, int jvar)
 
static double st_get_cova_param (void)
 
static int st_get_ncova_max (void)
 
static int st_get_ncova (void)
 
static int st_get_nvertex (int icov)
 
static int st_get_nvertex_max (void)
 
static int st_get_dimension (void)
 
static double st_get_model_mean (int ivar)
 
static double st_get_cova_range (void)
 
static double st_get_sill_total (int ivar, int jvar)
 
static void st_keypair_array (const char *name, int iter, double *tab)
 
static void st_print_all (const char *title)
 
static void st_compute_correc (void)
 
double spde_compute_correc (int ndim, double param)
 
static void st_compute_blin (void)
 
static void st_compute_hh ()
 
static void st_calcul_init (int ndim)
 
static void st_calcul_update (void)
 
static void st_convert_exponential2matern (CovAniso *cova)
 
int spde_attach_model (Model *model)
 
static int st_check_model (const Db *dbin, const Db *dbout, Model *model)
 
static int st_identify_nostat_param (const EConsElem &type0, int icov0=-1, int ivar0=-1, int jvar0=-1)
 
static double st_get_data_constraints (Db *db, int igrf, int iech)
 
static double st_simu_constraints (Db *db, int igrf, int iech, int iter, int ngibbs_burn, double yk, double sk)
 
static void st_gibbs (int igrf, int nout, int ngibbs_int, int iter0, int ngibbs_burn, Db *dbin, Db *dbout, double *zcur)
 
static void st_save_result (double *z, Db *dbout, const ELoc &locatorType, int iatt_simu)
 
static void st_merge (int ncur, double *z, double *zperm)
 
static void st_copy (int ncur, double *z, double *zperm)
 
static void st_simu_add_vertices (int number, const double *zsnc, double *zcur)
 
static void st_load_array (int number, const double *aux, double *perm)
 
static void st_init_array (int ncova, int nvar, int ncur, int flag_var, double *zperm)
 
static int st_simu_subtract_data (int ncur, const double *zsnc, const double *data, double *zerr)
 
static void st_kriging_one_rhs (QChol *QCtd, double *data, int ntarget, double *rhs)
 
static int st_kriging_cholesky (QChol *QC, double *rhs, VectorDouble &work, double *z)
 
static int st_filter (double *work, double *y)
 
double * _spde_get_mesh_dimension (AMesh *amesh)
 
static void st_calcul_update_nostat (AMesh *amesh, int imesh0)
 
static int st_fill_Isill (void)
 
static int st_fill_Csill (void)
 
static int st_fill_Bnugget (Db *dbin)
 
static int * st_get_vertex_ranks (AMesh *amesh, Db *dbin, Db *dbout)
 
static int st_fill_Bhetero (Db *dbin, Db *dbout)
 
static void st_triangle_center (AMesh *amesh, int ncorner, int imesh, double center[3], double xyz[3][3])
 
static void st_project_plane (double center[3], double axes[2][3], double xyz[3], double coeff[2])
 
static void st_tangent_calculate (double center[3], const double srot[2], double axes[2][3])
 
MatrixSparse_spde_fill_S (AMesh *amesh, Model *model, const double *units)
 
VectorDouble _spde_fill_TildeC (AMesh *amesh, const double *units)
 
VectorDouble _spde_fill_Lambda (Model *model, AMesh *amesh, const VectorDouble &TildeC)
 
static MatrixSparsest_extract_Q1_nugget (int row_var, int col_var, int *nrows, int *ncols)
 
static MatrixSparsest_extract_Q1_hetero (int row_var, int col_var, int row_oper, int col_oper, int *nrows, int *ncols)
 
static int st_build_QCov (SPDE_Matelem &Matelem)
 
MatrixSparse_spde_build_Q (MatrixSparse *S, const VectorDouble &Lambda, int nblin, double *blin)
 
static int st_build_Q (SPDE_Matelem &Matelem)
 
int spde_build_matrices (Model *model, int verbose)
 
static void st_load_data (AMesh *amesh, Db *dbin, Db *dbout, SPDE_Option &, int ivar0, double *data, double *zcur)
 
static double st_chebychev_function (double x, double power, const VectorDouble &blin)
 
static int st_chebychev_calculate_coeffs (Cheb_Elem *cheb_elem, int verbose, const VectorDouble &blin)
 
static void st_simulate_nugget (int ncur, VectorDouble &zsnc)
 
static int st_simulate_cholesky (QChol *QC, VectorDouble &work, VectorDouble &zsnc)
 
int spde_chebychev_operate (MatrixSparse *S, Cheb_Elem *cheb_elem, const VectorDouble &lambda, const VectorDouble &x, VectorDouble &y)
 
static int st_simulate_chebychev (VectorDouble &zsnc)
 
static int st_kriging_multigrid (QChol *QC, double *rhs, VectorDouble &work, double *z)
 
static int st_solve (QChol *QC, double *rhs, VectorDouble &work, double *z)
 
static int st_kriging_one (double *data, double *rhs, VectorDouble &work, double *z)
 
static int st_kriging_several_rhs (double *data, double *rhs, VectorDouble &work, double *ss_arg)
 
static int st_kriging_several_loop (int flag_crit, VectorDouble &work, double *rhscur, double *rhsloc, double *xcur, const double *rhs, double *crit)
 
static int st_kriging_several_results (const double *xcur, double *z)
 
static int st_kriging_several (double *data, double *rhs, VectorDouble &work, double *z)
 
static int st_kriging (AMesh *amesh, double *data, double *zkrig)
 
Cheb_Elemspde_cheb_manage (int mode, int verbose, double power, const VectorDouble &blin, MatrixSparse *S, Cheb_Elem *cheb_old)
 
Cheb_Elem_spde_cheb_duplicate (Cheb_Elem *cheb_in)
 
static void st_matelem_manage (int mode)
 
static int st_simulate (QChol *QC, double *zsnc)
 
int spde_process (Db *dbin, Db *dbout, SPDE_Option &s_option, int nbsimu, int ngibbs_nburn, int ngibbs_niter, int ngibbs_int)
 
static AMeshst_create_meshes (Db *dbin, Db *dbout, const VectorDouble &gext, SPDE_Option &s_option)
 
static int st_is_external_AQ_defined (int icov0)
 
int spde_external_copy (SPDE_Matelem &matelem, int icov0)
 
AMeshspde_mesh_load (Db *dbin, Db *dbout, const VectorDouble &gext, SPDE_Option &s_option, bool verbose)
 
void spde_external_mesh_define (int icov0, AMesh *mesh)
 
void spde_external_mesh_undefine (int icov0)
 
MatrixSparsespde_external_A_define (int icov0, MatrixSparse *A)
 
MatrixSparsespde_external_A_undefine (int icov0)
 
MatrixSparsespde_external_Q_define (int icov0, MatrixSparse *Q)
 
MatrixSparsespde_external_Q_undefine (int icov0)
 
void spde_mesh_assign (AMesh *amesh, int ndim, int ncorner, int nvertex, int nmesh, const VectorInt &arg_meshes, const VectorDouble &arg_points, int verbose)
 
int spde_prepar (Db *dbin, Db *dbout, const VectorDouble &gext, SPDE_Option &s_option)
 
int spde_posterior ()
 
static void st_environ_print (const Db *dbout, const VectorDouble &gext)
 
static bool is_in_mesh_neigh (AMesh *amesh, double *coor, double *caux, int ndim, int imesh, double radius)
 
static VectorDouble st_get_coords_3D (AMesh *amesh, const double *zcur)
 
void spde_free_all (void)
 
int spde_check (const Db *dbin, const Db *dbout, Model *model1, Model *model2, bool verbose, const VectorDouble &gext, bool mesh_dbin, bool mesh_dbout, bool flag_advanced, bool flag_est, bool flag_std, bool flag_gibbs, bool flag_modif)
 
int kriging2D_spde (Db *dbin, Model *model, SPDE_Option &s_option, int verbose, int *nmesh_arg, int *nvertex_arg, VectorInt &meshes_arg, VectorDouble &points_arg)
 
static void st_product_Q (const VectorDouble &blin, MatrixSparse *S, const VectorDouble &Lambda, const VectorDouble &TildeC, const VectorDouble &x, VectorDouble &y)
 
int spde_eval (const VectorDouble &blin, MatrixSparse *S, const VectorDouble &Lambda, const VectorDouble &TildeC, double power, VectorDouble &x, VectorDouble &y)
 
static int st_m2d_check_pinchout (Db *dbgrid, int icol_pinch)
 
static double st_m2d_draw_elevation (M2D_Environ *m2denv, int, int, double lower, double upper)
 
static void st_print_concatenate_interval (const char *title, double lower, double upper, int tail)
 
static void st_print_constraints_per_point (int ilayer, int iech, double value, double drift, double vgaus, double lower, double upper)
 
static int st_check_validity_MS (Db *db, int ilayer, int iech, int flag_positive, int flag_verbose, double M, double S)
 
static double st_m2d_get_M (M2D_Environ *m2denv, Db *db, int type, int ilayer, int iech)
 
static double st_m2d_get_S (M2D_Environ *m2denv, Db *, int, int, int)
 
static double st_m2d_external_drift_increment (M2D_Environ *m2denv, Db *db, int ilayer0, int iech0)
 
static double st_m2d_get_drift (M2D_Environ *m2denv, Db *db, int ilayer0, int iech0)
 
static void st_m2d_set_M (M2D_Environ *m2denv, int nlayer, int icol_pinch, Db *db, int iatt)
 
static int st_m2d_migrate_pinch_to_point (Db *dbout, Db *dbc, int icol_pinch)
 
static int st_m2d_drift_inc_manage (M2D_Environ *m2denv, int mode, int nlayer, int icol_pinch, Db *dbc, Db *dbout)
 
static void st_m2d_stats_init (M2D_Environ *m2denv, Db *dbin, int nlayer, int verbose)
 
static void st_m2d_stats_updt (M2D_Environ *m2denv, Db *dbc, int nlayer, int verbose)
 
static int st_m2d_initial_elevations (M2D_Environ *m2denv, Db *dbc, int nlayer, VectorDouble &work)
 
static int st_m2d_drift_manage (M2D_Environ *m2denv, Db *dbin, Db *dbout, int nlayer, int verbose, int *iatt_f)
 
static void st_print_details (Db *dbc, int nech, int ilayer)
 
static int st_m2d_drift_fitting (M2D_Environ *m2denv, Db *dbc, int nlayer, int number_hard, int verbose)
 
static void st_m2d_drift_save (M2D_Environ *m2denv, Db *dbout, int nlayer, double *gwork)
 
static int st_active_sample (Db *db, int ndim, int nlayer, int iech, int bypass)
 
static int st_record_sample (M2D_Environ *m2denv, Db *db, int iech, int ndim, int natt, int nlayer, int bypass, int *number_arg, double *tab)
 
static void st_define_locators (M2D_Environ *m2denv, Db *db, int ndim, int nvar, int nlayer)
 
static void st_m2d_print_environ (const char *title, M2D_Environ *m2denv)
 
static Dbst_m2d_create_constraints (M2D_Environ *m2denv, Db *dbin, Db *dbout, int ndim, int nlayer, int *number_hard)
 
static QCholst_derive_Qc (double s2, QChol *Qc, SPDE_Matelem &Matelem)
 
MatrixSparsedb_mesh_neigh (const Db *db, AMesh *amesh, double radius, int flag_exact, bool, int *nactive_arg, int **ranks_arg)
 
static double st_m2d_draw_gaussian (M2D_Environ *m2denv, Db *dbc, int verbose, int iter, int ilayer, int iech, double Zval, double Zcum, double Zmin, double Zmax, double Ymean, double Ysigma)
 
static void st_convert_Z2Y (M2D_Environ *m2denv, Db *dbc, int nlayer, int type, int iech, VectorDouble &tab)
 
static void st_convert_Y2Z (M2D_Environ *m2denv, Db *db, int nlayer, int type, int iech, VectorDouble &tab)
 
static void st_print_sample (const char *title, M2D_Environ *, Db *dbc, int nlayer, int iech, VectorDouble &work)
 
static int st_global_gibbs (M2D_Environ *m2denv, Db *dbc, int verbose, int iter, int nlayer, double sigma, VectorDouble &ymean, VectorDouble &ydat, VectorDouble &work)
 
static int st_check_gibbs_data (const char *title, M2D_Environ *m2denv, Db *dbc, int nlayer, int verbose, VectorDouble &ydat, VectorDouble &work)
 
static M2D_Environ * m2denv_manage (int mode, int flag_ed, double ystdv, M2D_Environ *m2denv_old)
 
static void st_m2d_vector_extract (M2D_Environ *m2denv, Db *dbc, int nlayer, VectorDouble &ydat, VectorDouble &work)
 
static void st_print_db_constraints (const char *title, Db *db, const VectorDouble &ydat, int nlayer, int verbose)
 
static void st_m2d_stats_gaus (const char *title, int nlayer, int nech, double *ydat)
 
int m2d_gibbs_spde (Db *dbin, Db *dbout, Model *model, int flag_ed, int nlayer, int niter, int seed, int nbsimu, int icol_pinch, int flag_drift, int flag_ce, int flag_cstd, int verbose)
 

Variables

static int DEBUG = 0
 
static int VERBOSE = 0
 
static int FLAG_KEYPAIR = 0
 
static double FACDIM [] = { 0., 1., 2., 6. }
 
static int SPDE_CURRENT_IGRF = 0
 
static int SPDE_CURRENT_ICOV = 0
 
static DbMEM_DBIN
 
static DbMEM_DBOUT
 
static AMeshS_EXTERNAL_MESH [3] = { NULL, NULL, NULL }
 
static MatrixSparseS_EXTERNAL_Q [3] = { NULL, NULL, NULL }
 
static MatrixSparseS_EXTERNAL_A [3] = { NULL, NULL, NULL }
 
static SPDE_Environ S_ENV
 
static SPDE_Decision S_DECIDE
 
static char NAME [100]
 
static char string_encode [100]
 
static SPDE_Calcul Calcul
 

Macro Definition Documentation

◆ NBLIN_TERMS

#define NBLIN_TERMS   10

◆ SPDE_MAX_NGRF

#define SPDE_MAX_NGRF   2

Function Documentation

◆ _spde_build_Q()

MatrixSparse* _spde_build_Q ( MatrixSparse S,
const VectorDouble Lambda,
int  nblin,
double *  blin 
)

Construct the final sparse matrix Q from the Model

Returns
Error return code
Parameters
[in]SShift operator
[in]LambdaLambda vector
[in]nblinNumber of blin coefficients
[in]blinArray of coefficients for Linear combinaison

◆ _spde_cheb_duplicate()

Cheb_Elem* _spde_cheb_duplicate ( Cheb_Elem cheb_in)

Duplicate a Cheb_Elem structure

Returns
Pointer to the newly created Cheb_Eleme structure
Parameters
[in]cheb_inInput Cheb_Eleme structure

◆ _spde_fill_Lambda()

VectorDouble _spde_fill_Lambda ( Model model,
AMesh amesh,
const VectorDouble TildeC 
)

Fill the vector for sill correction factors Works for both stationary and non-stationary cases

Parameters
[in]modelModel structure
[in]ameshMeshEStandard structure
[in]TildeCVector TildeC

◆ _spde_fill_S()

MatrixSparse* _spde_fill_S ( AMesh amesh,
Model model,
const double *  units 
)

Fill the sparse matrix S linked to mesh vertices

Returns
G sparse matrix
Parameters
[in]ameshMeshEStandard_Mesh structure
[in]modelModel structure
[in]unitsArray containing the mesh dimensions

◆ _spde_fill_TildeC()

VectorDouble _spde_fill_TildeC ( AMesh amesh,
const double *  units 
)

Fill the vector TildeC (Dimension: nvertex)

Returns
Error return code
Parameters
[in]ameshMeshEStandard structure
[in]unitsArray containing the element units

◆ _spde_get_mesh_dimension()

double* _spde_get_mesh_dimension ( AMesh amesh)

Calculate the array of dimensions of the meshes

Returns
Pointer to the newly allocated array containing mesh dimensions
Parameters
[in]ameshMeshEStandard structure
Remarks
The array returned by this function must be deallocated

◆ db_mesh_neigh()

MatrixSparse* db_mesh_neigh ( const Db db,
AMesh amesh,
double  radius,
int  flag_exact,
bool  ,
int *  nactive_arg,
int **  ranks_arg 
)

Returns the projection matrix of a set of points (contained in a Db) onto a meshing

Returns
Pointer to the newly created sparse matrix (or NULL)
Parameters
[in]dbDb structure
[in]ameshAMesh structure
[in]flag_exactType of test for intersection (See remarks)
[in]radiusNeighborhood radius
[out]nactive_argNumber of active samples from the Db
[out]ranks_argRanks of the active samples
Remarks
The calling function must free the argument 'ranks'
When flag_exact is TRUE, for each active sample of Db, a vertex
of the mesh is active as soon as it lies within the vicinity
of the sample.
If flag_exact is FALSE, all vertices of a mesh are considered as
active as soon as the mesh intersects the ball around a sample.
The vicinity is defined as any point located at a distance
from the sample smaller than 'radius'. The distance is calculated
as the Euclidean distance over the space whose dimension is
if the smallest value between the Db et Mesh space dimensions.

◆ is_in_mesh_neigh()

static bool is_in_mesh_neigh ( AMesh amesh,
double *  coor,
double *  caux,
int  ndim,
int  imesh,
double  radius 
)
static

True if a (dilated) point intersects a mesh

Returns
1 if the intersection is not empty; 0 otherwise
Parameters
[in]ameshAMesh structure
[in]coorPoint coordinates (Dimension: ndim)
[in]cauxWorking array (Dimension: ndim)
[in]ndimComon space dimension between coor and s_mesh
[in]imeshMesh Index
[in]radiusDilation radius
Remarks
This function must not be called if on a sphere

◆ kriging2D_spde()

int kriging2D_spde ( Db dbin,
Model model,
SPDE_Option s_option,
int  verbose,
int *  nmesh_arg,
int *  nvertex_arg,
VectorInt meshes_arg,
VectorDouble points_arg 
)

Perform Kriging using SPDE on the set of constructed 2-D triangles

Returns
Error return code
Parameters
[in]dbinDb input structure
[in]modelModel structure
[in]s_optionSPDE_Option structure
[in]verboseVerbose option
[out]nmesh_argNumber of triangles generated
[out]nvertex_argNumber of vertices
[out]meshes_argArray of triangulate vertex indices
[out]points_argArray containing the 2-D vertices
Remarks
It is the responsability of the calling function to free
the returned arrays 'points_arg' and 'meshes_arg'.
This function assumes that there is only one Model defined
Hence the test on the number of models

◆ m2d_gibbs_spde()

int m2d_gibbs_spde ( Db dbin,
Db dbout,
Model model,
int  flag_ed,
int  nlayer,
int  niter,
int  seed,
int  nbsimu,
int  icol_pinch,
int  flag_drift,
int  flag_ce,
int  flag_cstd,
int  verbose 
)

Perform Gibbs on a multilayer setup

Returns
Error return code
Parameters
[in]dbinDb input structure
[in]dboutDb output structure
[in]modelModel structure
[in]flag_ed1 if External Drit is used
[in]nlayerNumber of layers
[in]niterNumber of iterations
[in]seedSeed for random number generator
[in]nbsimuNumber of simulaations
[in]icol_pinchAddress of the variable containing the pinchout
[in]flag_drift1 to return the drift only 0 the simulations
[in]flag_ce1 if the conditional expectation should be returned instead of simulations
[in]flag_cstd1 if the conditional standard deviation should be returned instead of simulations
[in]verboseVerbose option
Remarks
In 'dbin':
- the lower and upper bounds must be defined for each datum
(set to the locator ELoc::L and ELoc::U
In 'dbout':
- the trend (if flag_ed is 1) must be defined and set to
the locator ELoc::F
When defined, the pinchout should be defined as a grid variable
with values ranging between 0 and 1 (FFFF are admitted).
It will serve as a multiplier to the Mean thickness maps.

◆ m2denv_manage()

static M2D_Environ* m2denv_manage ( int  mode,
int  flag_ed,
double  ystdv,
M2D_Environ *  m2denv_old 
)
static

Manage the M2D_Environ structure

Returns
Pointer to the M2D_Environ structure
Parameters
[in]mode1 for allocation; -1 for deallocation
[in]flag_ed1 if external drift is used; 0 otherwise
[in]ystdvStamdard deviation of the Gaussian Transformed
[in]m2denv_oldPointer to the already existing M2D_Environ (only used when mode==-1)

◆ qchol_manage()

QChol* qchol_manage ( int  mode,
QChol QC 
)

Manage the Sparse matrix and its cholesky decomposition

Returns
Pointer to the QChol structure
Parameters
[in]modeManagement mode 1 : Allocation -1 : Total deallocation (NULL is returned)
[in]QCPointer to the QChol structure (when mode < 0)

◆ simu_define_func_scale()

void simu_define_func_scale ( void(*)(Db *, int, int)  st_simu_scale)

Define the function to scale the Modification arrays

Parameters
[in]st_simu_scalePointer to the scaling function

◆ simu_define_func_transf()

void simu_define_func_transf ( void(*)(Db *, int, int, int)  st_simu_transf)

Define the function to transform a simulation

Parameters
[in]st_simu_transfPointer to the transformation function

◆ simu_define_func_update()

void simu_define_func_update ( void(*)(Db *, int, int, int)  st_simu_update)

Define the function to account for the current simulation outcome in the calculation of the Modification arrays

Parameters
[in]st_simu_updatePointer to the update function

◆ spde_attach_model()

int spde_attach_model ( Model model)

Attach the model

Returns
Error returned code
Parameters
[in]modelModel structure

◆ spde_build_matrices()

int spde_build_matrices ( Model model,
int  verbose 
)

Build all matrices needed for establishing the Q sparse matrix

Returns
Error return code
Parameters
[in]modelModel structure
[in]verboseVerbose option
Remarks
Contents of SP_MAT (sparse matrices or vectors) is allocated here
It must be freed by the calling functions

◆ spde_cheb_manage()

Cheb_Elem* spde_cheb_manage ( int  mode,
int  verbose,
double  power,
const VectorDouble blin,
MatrixSparse S,
Cheb_Elem cheb_old 
)

Manage Cheb_Elem structure

Returns
Error return code
Parameters
[in]mode1 for allocation; -1 for deallocation
[in]verboseVerbose flag
[in]powerParameter passed to Chebychev function
[in]blinArray of coefficients for Linear combinaison
[in]SShift operator
[in]cheb_oldCheb_Elem to be freed (only for mode=-1)
Remarks
Arguments 'power', 'nblin', 'blin' and 'B' are used if mode=1
Argument 'cheb_old' is used if mode=-1

◆ spde_chebychev_operate()

int spde_chebychev_operate ( MatrixSparse S,
Cheb_Elem cheb_elem,
const VectorDouble lambda,
const VectorDouble x,
VectorDouble y 
)

Perform the Chebychev polynomial procedure on an input vector

Returns
Error return code
Parameters
[in]SShift operator
[in]cheb_elemCheb_Elem structure
[in]lambdaScaling vector
[in]xInput array (Dimension: nvertex)
[out]yOutput array (Dimension: nvertex)

◆ spde_check()

int spde_check ( const Db dbin,
const Db dbout,
Model model1,
Model model2,
bool  verbose,
const VectorDouble gext,
bool  mesh_dbin,
bool  mesh_dbout,
bool  flag_advanced,
bool  flag_est,
bool  flag_std,
bool  flag_gibbs,
bool  flag_modif 
)

Define the main options

Returns
Error return code
Parameters
[in]dbinPointer to the input Db
[in]dboutPointer to the output Db
[in]model1Model structure (first)
[in]model2Model structure (second)
[in]verboseVerbose flag
[in]gextArray of domain dilation
[in]mesh_dbinTrue if Input points must participate to meshing
[in]mesh_dboutTrue if Output points must participate to meshing
[in]flag_advancedTrue for advanced calculus (estimation or simulation) False if only matrices are required
[in]flag_estTrue for estimation
[in]flag_stdTrue for standard deviation
[in]flag_gibbsTrue for Gibbs sampler
[in]flag_modifTrue for post-processing simulations
Remarks
This function initiates the Matelem structures

◆ spde_compute_correc()

double spde_compute_correc ( int  ndim,
double  param 
)

◆ spde_eval()

int spde_eval ( const VectorDouble blin,
MatrixSparse S,
const VectorDouble Lambda,
const VectorDouble TildeC,
double  power,
VectorDouble x,
VectorDouble y 
)

Perform the product of a vector by the inverse of the power of a sparse matrix using the Chebychev Polynomial procedure

Returns
Error return code
Parameters
[in]blinArray of coefficients for Linear combinaison
[in]SShift operator
[in]LambdaVector Lambda
[in]TildeCVector TildeC
[in]powerParameter used in the Chebychev approximation
[in]xInput array
[out]yOutput array

◆ spde_external_A_define()

MatrixSparse* spde_external_A_define ( int  icov0,
MatrixSparse A 
)

◆ spde_external_A_undefine()

MatrixSparse* spde_external_A_undefine ( int  icov0)

◆ spde_external_copy()

int spde_external_copy ( SPDE_Matelem matelem,
int  icov0 
)

Copy the contents of the internal S_EXTERNAL_AQ into an output Matelem

Error return code

Parameters
[in]matelemOutput SPDE_Matelem structure
[in]icov0Rank of the current Covariance

◆ spde_external_mesh_define()

void spde_external_mesh_define ( int  icov0,
AMesh mesh 
)

Manage the contents of the External AMesh structure used for storing external meshing information

Parameters
[in]icov0Rank of the current covariance (from 0 to 2)
[in]meshAMesh pointer

◆ spde_external_mesh_undefine()

void spde_external_mesh_undefine ( int  icov0)

◆ spde_external_Q_define()

MatrixSparse* spde_external_Q_define ( int  icov0,
MatrixSparse Q 
)

◆ spde_external_Q_undefine()

MatrixSparse* spde_external_Q_undefine ( int  icov0)

◆ spde_free_all()

void spde_free_all ( void  )

Free all memory used in SPDE

◆ spde_get_current_matelem()

SPDE_Matelem& spde_get_current_matelem ( int  icov)

Get the pointer to the current SPDE_Matelem structure

Parameters
[in]icovRank of the target Covariance (or -1)typedef struct

◆ spde_mesh_assign()

void spde_mesh_assign ( AMesh amesh,
int  ndim,
int  ncorner,
int  nvertex,
int  nmesh,
const VectorInt arg_meshes,
const VectorDouble arg_points,
int  verbose 
)

Assign fields of the AMesh structure which would have been calculated elsewhere

Parameters
[in]ndimSpace dimension
[in]ncornerNumber of vertices per element
[in]nvertexNumber of points
[in]nmeshNumber of meshes
[in]arg_meshesArray containing the meshes
[in]arg_pointsArray containing the vertex coordinates
[in]verboseVerbose flag
[in,out]ameshPointer to AMesh to be assigned

◆ spde_mesh_load()

AMesh* spde_mesh_load ( Db dbin,
Db dbout,
const VectorDouble gext,
SPDE_Option s_option,
bool  verbose 
)

Load the AMesh structure

Returns
Pointer on the newly allocated AMesh
Parameters
[in]dbinDb structure for the conditioning data
[in]dboutDb structure of the grid
[in]gextArray of domain dilation
[in]s_optionSPDE_Option structure
[in]verboseVerbose option

◆ spde_option_alloc()

SPDE_Option spde_option_alloc ( void  )

Manage the SPDE_Option structure

◆ spde_option_update()

void spde_option_update ( SPDE_Option s_option,
const String triswitch 
)

Add a new option item for an additional covariance structure

Parameters
[in]s_optionPointer to SPDE_Option to be freed (if mode<0)
[in]triswitchString defining the meshing characteristics

◆ spde_posterior()

int spde_posterior ( )

Cleaning operation after SPDE

Returns
Error return code

◆ spde_prepar()

int spde_prepar ( Db dbin,
Db dbout,
const VectorDouble gext,
SPDE_Option s_option 
)

Preparation using SPDE (for all GRF and COV)

Returns
Error return code
Parameters
[in]dbinDb structure for the conditioning data
[in]dboutDb structure of the grid
[in]gextArray of domain dilation
[in]s_optionSPDE_Option structure

◆ spde_process()

int spde_process ( Db dbin,
Db dbout,
SPDE_Option s_option,
int  nbsimu,
int  ngibbs_nburn,
int  ngibbs_niter,
int  ngibbs_int 
)

Perform the main Simulations steps after Q has been constructed

Returns
Error return code
Parameters
[in]dbinInput Db grid structure (optional)
[in]dboutOutput Db grid structure
[in]s_optionSPDE_Option structure
[in]nbsimuNumber of simulations
[in]ngibbs_nburnNumber of iterations (Burning step)
[in]ngibbs_niterNumber of iterations
[in]ngibbs_intNumber of iterations internal to Gibbs (SPDE)
Remarks
If the number of simulations 'nbsimu' is set to 0,
the simulation algorithm is turned into a kriging one

◆ st_active_sample()

static int st_active_sample ( Db db,
int  ndim,
int  nlayer,
int  iech,
int  bypass 
)
static

Check if a sample must be considered as an active constraint

Returns
1 if the sample is active; 0 otherwise
Parameters
[in]dbDb input structure
[in]ndimSpace dimension
[in]nlayerNumber of layers
[in]iechRank of the sample
[in]bypass1 to bypass check that at least one bound is defined
Remarks
A sample is an active constraint if at least one constraint
is defined

◆ st_build_Q()

static int st_build_Q ( SPDE_Matelem Matelem)
static

Construct the final sparse matrix Q from the Model

Returns
Error return code
Parameters
[in]MatelemSPDE_Matelem structure

◆ st_build_QCov()

static int st_build_QCov ( SPDE_Matelem Matelem)
static

Construct the sparse matrix QCov (used in multistructure - multivariable)

Returns
Error return code
Parameters
[in]MatelemSPDE_Matelem structure
Remarks
This function requires the Q matrices to be established already,
as well as the Aproj matrices.
In case of presence of nugget effect, we also need 'Bnugget'
Otherwise, we need 'BheteroD' and 'BheteroT'

◆ st_calcul_init()

static void st_calcul_init ( int  ndim)
static

Initialize the contents of SPDE_Calcul structure

◆ st_calcul_update()

static void st_calcul_update ( void  )
static

Update the contents of SPDE_Calcul structure

◆ st_calcul_update_nostat()

static void st_calcul_update_nostat ( AMesh amesh,
int  imesh0 
)
static

Update parameters in S_ENV structure in non-stationary case

Parameters
[in]ameshMeshEStandard structure
[in]imesh0Rank of the current mesh

◆ st_chebychev_calculate_coeffs()

static int st_chebychev_calculate_coeffs ( Cheb_Elem cheb_elem,
int  verbose,
const VectorDouble blin 
)
static

Evaluate the number of coefficients necessary to evaluate a function (at a sample location) at a given approximation

Returns
Error return code
Parameters
[in]cheb_elemCheb_Elem structure to be filled
[in]verboseVerbose flag
[in]blinArray of coefficients for Linear combination

◆ st_chebychev_function()

static double st_chebychev_function ( double  x,
double  power,
const VectorDouble blin 
)
static

Internal function used for the Chebychev approximation

Returns
Returned value
Parameters
[in]xInput value
[in]powerParameter used in the Chebychev approximation
[in]blinArray of coefficients for Linear combination

◆ st_check_gibbs_data()

static int st_check_gibbs_data ( const char *  title,
M2D_Environ *  m2denv,
Db dbc,
int  nlayer,
int  verbose,
VectorDouble ydat,
VectorDouble work 
)
static

Check that the Gibbs constraints are fullfilled at datum locations

Returns
Error return code
Parameters
[in]titleTitle for the printout (if error)
[in]m2denvM2D_Environ structure
[in]dbcDb structure containing the constraints
[in]nlayerNumber of layers
[in]verboseVerbose flag
[in]ydatArray of simulations on the data
[out]workArray of tentative values (Dimension: nlayer)

◆ st_check_model()

static int st_check_model ( const Db dbin,
const Db dbout,
Model model 
)
static

Check that the Model is authorized for SPDE

Returns
Error returned code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure
[in]modelModel structure

◆ st_check_validity_MS()

static int st_check_validity_MS ( Db db,
int  ilayer,
int  iech,
int  flag_positive,
int  flag_verbose,
double  M,
double  S 
)
static

Check the validity of the Mean and Variance values

Returns
Error return code
Parameters
[in]dbDb structure containing the constraints
[in]ilayerRank of the layer of interest
[in]iechRank of the sample of interest
[in]flag_positivePositivity check
[in]flag_verboseVerbose output
[in]MValue for the Mean
[in]SValue for the Variance

◆ st_clean_Bhetero()

static void st_clean_Bhetero ( void  )
static

Clean the Bhetero sparse matrices

◆ st_clean_Bnugget()

static void st_clean_Bnugget ( void  )
static

Clean the Bnugget sparse matrices

◆ st_compute_blin()

static void st_compute_blin ( void  )
static

Compute the coefficients of the linear combination

Remarks
This function stores the coefficients 'blin' in SPDE_Calcul

◆ st_compute_correc()

static void st_compute_correc ( void  )
static

Compute the variance correction term Store in the SPDE_Calcul structure

◆ st_compute_hh()

static void st_compute_hh ( )
static

Compute H matrix for anisotropic case and the square root of determinant Requires the knowledge of the actual parameters of the current Covariance Fills the SPDE_Calcul structure

◆ st_convert_exponential2matern()

static void st_convert_exponential2matern ( CovAniso cova)
static

Modify the Exponential into a Matern

Parameters
[in]covaCovariance sructure

◆ st_convert_Y2Z()

static void st_convert_Y2Z ( M2D_Environ *  m2denv,
Db db,
int  nlayer,
int  type,
int  iech,
VectorDouble tab 
)
static

Convert a layer-pile at a datum from the true to the working domain

Parameters
[in]m2denvM2D_Environ structure
[in]dbDb structure
[in]nlayerNumber of layers
[in]type1 for the constraining Db 2 for the grid output Db
[in]iechRank of the sample
[in,out]tabInput/Ouput array of Y-values (Dimension: nlayer)

◆ st_convert_Z2Y()

static void st_convert_Z2Y ( M2D_Environ *  m2denv,
Db dbc,
int  nlayer,
int  type,
int  iech,
VectorDouble tab 
)
static

Convert a layer-pile at a datum from the working to the true domain

Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb structure
[in]nlayerNumber of layers
[in]type1 for the constraining Db 2 for the grid output Db
[in]iechRank of the sample
[in,out]tabInput/Output array of Z-values (Dimension: nlayer)

◆ st_copy()

static void st_copy ( int  ncur,
double *  z,
double *  zperm 
)
static

Copy an auxiliary array into an permanent array

Parameters
[in]ncurNumber of elements per variable
[in]zArray of values to be merged
[out]zpermPermanent array
Remarks
This operation takes place for all samples located within the Part

◆ st_create_meshes()

static AMesh* st_create_meshes ( Db dbin,
Db dbout,
const VectorDouble gext,
SPDE_Option s_option 
)
static

Load the meshes

Returns
Pointer to the newly created AMesh structure
Parameters
[in]dbinDb structure for the conditioning data
[in]dboutDb structure of the grid
[in]gextArray of domain dilation
[in]s_optionSPDE_Option structure
Remarks
The option 'flag_force' forces to use the regular meshing rather
than the Turbo one

◆ st_define_locators()

static void st_define_locators ( M2D_Environ *  m2denv,
Db db,
int  ndim,
int  nvar,
int  nlayer 
)
static

Define the locators on the newly created Db

Parameters
[in]m2denvM2D_Environ structure
[in]dbDb constraints structure
[in]ndimNumber of coodinates
[in]nlayerNumber of layers
[in]nvarNumber of variables

◆ st_derive_Qc()

static QChol* st_derive_Qc ( double  s2,
QChol Qc,
SPDE_Matelem Matelem 
)
static

Calculate the inverse of (s2 * Q + B %*% Bt) and store it into a new QChol object

Returns
The calculated sparse matrix
Parameters
[in]s2Nugget effect value
[in]QcQc structure (already existing)
[in]MatelemMatelem structure

◆ st_environ_init()

static void st_environ_init ( void  )
static

Initialize the S_ENV Environment structure

◆ st_environ_print()

static void st_environ_print ( const Db dbout,
const VectorDouble gext 
)
static

Print the environment

Parameters
[in]dboutDb output structure
[in]gextArray of domain dilation

◆ st_extract_Q1_hetero()

static MatrixSparse* st_extract_Q1_hetero ( int  row_var,
int  col_var,
int  row_oper,
int  col_oper,
int *  nrows,
int *  ncols 
)
static

Extract the sparse matrix from the Q matrix (coninuous structure)

Returns
The extracted sparse matrix or NULL
Parameters
[in]row_varRank of the variable for the row
[in]col_varRank of the variable for the column
[in]row_operOperator type for row (1:Data or 2:Target)
[in]col_operOperator type for column (1:Data or 2:Target)
[out]nrowsNumber of rows
[out]ncolsNumber of columns
Remarks
Extracts a part of Q matrix (for the first structure) for:
- a given pair of variables
- a given pair of operators (Data or Target)
The returned matrix is multipled by the inverse of the Sill

◆ st_extract_Q1_nugget()

static MatrixSparse* st_extract_Q1_nugget ( int  row_var,
int  col_var,
int *  nrows,
int *  ncols 
)
static

Extract the sparse matrix from the Q matrix (case of nugget effect)

Returns
The extracted sparse matrix or NULL
Parameters
[in]row_varRank of the variable for the row
[in]col_varRank of the variable for the column
[out]nrowsNumber of rows
[out]ncolsNumber of columns
Remarks
Extracts a part of Bnugget matrix for:
- a given pair of variables
- for Data-Data operators

◆ st_extract_Q_from_Q()

static MatrixSparse* st_extract_Q_from_Q ( MatrixSparse Q_in,
int  row_auth,
int  col_auth 
)
static

Construct the sparse matrix Q from another sparse matrix

Returns
The Q structure or NULL
Parameters
[in]Q_inInput sparse matrix
[in]row_authSpecification for rows extraction
[in]col_authSpecification for columns extraction
Remarks
The Cholesky decomposition is performed (if possible)

◆ st_extract_QC_from_Q()

static QChol* st_extract_QC_from_Q ( const char *  title,
QChol QC_in,
int  row_auth,
int  col_auth 
)
static

Construct the QChol sub-structure from another QChol structure

Returns
The QChol structure or NULL
Parameters
[in]titleName of the QChol item
[in]QC_inInput QChol structure
[in]row_authSpecification for rows extraction
[in]col_authSpecification for columns extraction

◆ st_fill_Bhetero()

static int st_fill_Bhetero ( Db dbin,
Db dbout 
)
static

Fill some matrices for Kriging in the case of a model without nugget effect Constitute the Bhetero sparse matrices

Returns
Error returned code
Parameters
[in]dbinInput Db structure
[in]dboutOutput Db structure

◆ st_fill_Bnugget()

static int st_fill_Bnugget ( Db dbin)
static

Fill the Bnugget sparse matrix linked to nugget effect

Returns
Error returned code
Parameters
[in]dbinDb structure
Remarks
This function allocates 'nvs2' sparse matrices of dimension 'ndata'.
where nvs2 is the product nvar * (nvar+1) / 2

◆ st_fill_Csill()

static int st_fill_Csill ( void  )
static

Fill the Csill matrix linked to the continuous parts of the Model

Returns
Error returned code
Remarks
The matrix 'Csill' is dimensioned to ncova * nvar * (nvar+1)/2 where
- ncova designates the number of continuous structures of the Model

◆ st_fill_Isill()

static int st_fill_Isill ( void  )
static

Fill the Isill matrix linked to the covariance of the Model

Returns
Error returned code
Remarks
The matrix 'Isill' is dimensioned to nvar * nvar where

◆ st_filter()

static int st_filter ( double *  work,
double *  y 
)
static

Perform the Filtering of Nugget Effect

Returns
Error return code
Parameters
[out]workWorking array (Dimension: nvertex)
[in]yOutput Array (Dimension: nvertex)

◆ st_get_coords_3D()

static VectorDouble st_get_coords_3D ( AMesh amesh,
const double *  zcur 
)
static

Construct the array of coordinates of triangles in 3-D space by gluing the original 2-D coordinates (turned back in the original space) with the information provided by the third dimension array

Returns
Array of 3-D coordinates
Parameters
[in]ameshAMesh structure
[in]zcurArray of results

◆ st_get_cova()

static CovAniso* st_get_cova ( void  )
static

Returns the pointer to the structure

◆ st_get_cova_param()

static double st_get_cova_param ( void  )
static

Return the param of the model

◆ st_get_cova_range()

static double st_get_cova_range ( void  )
static

Get the normalized range

◆ st_get_cova_sill()

static double st_get_cova_sill ( int  ivar,
int  jvar 
)
static

Return the sill of the model (or TEST)

Parameters
[in]ivarRank of the first variable
[in]jvarRank of the second variable
Remarks
To save time, no check is performed with respect to the rank
of the structure or of the variables

◆ st_get_data_constraints()

static double st_get_data_constraints ( Db db,
int  igrf,
int  iech 
)
static

Return the initial value given to a data under constraints

Returns
Initial value
Parameters
[in]dbDb structure
[in]igrfRank of the GRF
[in]iechSample rank
Remarks
The test that Interval exists has already been performed

◆ st_get_dimension()

static int st_get_dimension ( void  )
static

Return the dimension for array allocation This dimension is the aximum of the maximum number of vertices and the number of data

◆ st_get_filnug()

static int st_get_filnug ( void  )
static

Return if a nugget effect component must be filtered

◆ st_get_isill()

static double st_get_isill ( int  icov,
int  ivar,
int  jvar 
)
static

Get the value of the Inverse of the sill for a given covariance and a pair of variables

◆ st_get_model()

static Model* st_get_model ( void  )
static

Return the global non-stationary characteristics

◆ st_get_model_mean()

static double st_get_model_mean ( int  ivar)
static

Return the mean of the model

Parameters
[in]ivarRank of the target variable (for its mean)

◆ st_get_ncova()

static int st_get_ncova ( void  )
static

Returns the number of structures in the Model (nugget excluded)

◆ st_get_ncova_max()

static int st_get_ncova_max ( void  )
static

Return the maximum number of covariances in the different Models for all GRFS (nugget included)

◆ st_get_nugget()

static CovAniso* st_get_nugget ( void  )
static

Returns the pointer to structure containing the Nugget Effect (or NULL)

◆ st_get_nugget_sill()

static double st_get_nugget_sill ( int  ivar,
int  jvar 
)
static

Return the sill of the Nugget Effect (or TEST)

Parameters
[in]ivarRank of the first variable
[in]jvarRank of the second variable
Remarks
To save time, no check is performed with respect to the rank
of the variables

◆ st_get_number_grf()

static int st_get_number_grf ( void  )
static

Returns the number of GRFs of the environment

◆ st_get_nvertex()

static int st_get_nvertex ( int  icov)
static

Return the number of vertices for the current Matelem

Parameters
[in]icovRank of the target Covariance (or -1)

◆ st_get_nvertex_max()

static int st_get_nvertex_max ( void  )
static

Return the maximum number of vertices for all meshes for all GRFs

◆ st_get_nvs2()

static int st_get_nvs2 ( void  )
static

Return the number of variables of the environment

◆ st_get_rank()

static int st_get_rank ( int  ivar,
int  jvar 
)
static

Returns the index of a pair of variable ranks within the triangle

Returns
Absolute index
Parameters
[in]ivarRank of the first variable
[in]jvarRank of the second variable
Remarks
The calling function does not have to bother of the relative
order between 'ivar' and 'jvar'

◆ st_get_sill_total()

static double st_get_sill_total ( int  ivar,
int  jvar 
)
static

Return the total sill of the model

Parameters
[in]ivarRank of the first variable
[in]jvarRank of the second variable
Remarks
To save time, no check is performed with respect to the rank
of the variables

◆ st_get_vertex_ranks()

static int* st_get_vertex_ranks ( AMesh amesh,
Db dbin,
Db dbout 
)
static

Return the list of (target+data) indices for a given mesh

Returns
An array of vertex identification (Dimension: nvertex) or NULL
Parameters
[in]ameshAMesh structure
[in]dbinDb structure for input (optional)
[in]dboutDb structure for the output
Remarks
The array ranks is filled as follows:
- Its contents follows the mesh numbering
- If positive, its value provides the rank of the data
- If negative, its absolute value provides the rank of the target
- If zero, theses are Steiner points
Warning: Ranks are counted from 1
The returned array must be freed by the calling function
Dimension: nvertex

◆ st_gibbs()

static void st_gibbs ( int  igrf,
int  nout,
int  ngibbs_int,
int  iter0,
int  ngibbs_burn,
Db dbin,
Db dbout,
double *  zcur 
)
static

Perform one iteration of the Gibbs sampler

Parameters
[in]igrfRank of the GRF
[in]noutNumber of samples on which Gibbs should be run
[in]ngibbs_intNumber of internal Gibbs iterations
[in]iter0Gibbs iteration rank (starting from 0)
[in]ngibbs_burnNumber of iterations (Burning step)
[in]dbinInput Db
[in]dboutOutput Db
[out]zcurVector of simulation (VT_FREE|VT_GIBBS|VT_HARD)
Remarks
The argument 'iter0' refers to the iteration rank, including
'ngibbs_burn' and 'ngibbs_iter'

◆ st_global_gibbs()

static int st_global_gibbs ( M2D_Environ *  m2denv,
Db dbc,
int  verbose,
int  iter,
int  nlayer,
double  sigma,
VectorDouble ymean,
VectorDouble ydat,
VectorDouble work 
)
static

Perform the Gibbs iterations

Returns
Error returned code
Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb structure containing the constraints
[in]verboseVerbose flag
[in]iterRank of the iteration
[in]nlayerNumber of layers
[in]sigmaStandard deviation of the nugget value
[in]ymeanArray of mean values at constraints
[in,out]ydatArray of values at constraints samples
[out]workArray of tentative values (Dimension: nlayer)
Remarks
The [in,out] argument 'ydat' is expressed in working domain
It needs to be locally transformed in real domain in order to
be compared to the bounds.

◆ st_identify_nostat_param()

static int st_identify_nostat_param ( const EConsElem &  type0,
int  icov0 = -1,
int  ivar0 = -1,
int  jvar0 = -1 
)
static

Identify a parameter among the non-stationary ones

Returns
The rank of the parameter of -1 (not found)
Parameters
[in]type0Type of parameter (EConsElem)
[in]icov0Rank of the target covariance
[in]ivar0Rank of the target variable (only when type=EConsElem::SILL)
[in]jvar0Rank of the target variable (only when type=EConsElem::SILL)
Remarks
The covariance are ranked from 0 for non-nugget ones
The nugget effect corresponds to rank (-1)

◆ st_init_array()

static void st_init_array ( int  ncova,
int  nvar,
int  ncur,
int  flag_var,
double *  zperm 
)
static

Initialize an array with a constant mean value

Parameters
[in]ncovaNumber of structures
[in]nvarNumber of variables
[in]ncurNumber of terms in arrays 'perm' and 'aux'
[in]flag_var1 : to initialize a variable; 0 : to initialize variance
[out]zpermOutput array

◆ st_is_external_AQ_defined()

static int st_is_external_AQ_defined ( int  icov0)
static

Check if External Q has been defined

Returns
1 if External A-Q has been defined (for 'icov'); 0 otherwise
Parameters
[in]icov0Rank of the (non-nugget) covariance

◆ st_is_model_nugget()

static bool st_is_model_nugget ( void  )
static

Checks if there is a nugget component in the Model

Returns
true if a Nugget component is present; false otherwise

◆ st_keypair_array()

static void st_keypair_array ( const char *  name,
int  iter,
double *  tab 
)
static

Use the keypair mechanism to output some arrays

Parameters
[in]nameName of the output
[in]iterConcatenated in name if >= 0; ignored otherwise
[in]tabTab to be output

◆ st_kriging()

static int st_kriging ( AMesh amesh,
double *  data,
double *  zkrig 
)
static

Perform the Kriging operation

Returns
Error return code
Parameters
[in]ameshAMesh structure
[in]dataVector of data
[out]zkrigOutput Array

◆ st_kriging_cholesky()

static int st_kriging_cholesky ( QChol QC,
double *  rhs,
VectorDouble work,
double *  z 
)
static

Perform the Calculation of the Kriging estimate

Returns
Error return code
Parameters
[in]QCPointer to QChol structure
[in]rhsR.H.S. array (Dimension: ntarget)
[out]workWorking array (Dimension: ntarget)
[out]zOutput array (Dimension: ntarget)

◆ st_kriging_multigrid()

static int st_kriging_multigrid ( QChol QC,
double *  rhs,
VectorDouble work,
double *  z 
)
static

Perform the Calculation of the Kriging estimate (Multigrid case)

Returns
Error return code
Parameters
[in]QCPointer to QChol structure (target-target)(finalized)
[in]rhsR.H.S. array (Dimension: ntarget)
[out]workWorking array (Dimension: ntarget)
[out]zOutput array (Dimension: ntarget)

◆ st_kriging_one()

static int st_kriging_one ( double *  data,
double *  rhs,
VectorDouble work,
double *  z 
)
static

Perform the Kriging operation (monovariate - monostructure)

Returns
Error return code
Parameters
[in]dataVector of data
[out]rhsR.H.S.
[out]workWorking array
[out]zOutput Array

◆ st_kriging_one_rhs()

static void st_kriging_one_rhs ( QChol QCtd,
double *  data,
int  ntarget,
double *  rhs 
)
static

Prepare the RHS for Kriging

Parameters
[in]QCtdPointer to QChol structure (target-data)
[in]dataInput array (Dimension: ndata)
[in]ntargetNumber of target values
[out]rhsOutput array (Dimension: ntarget)

◆ st_kriging_several()

static int st_kriging_several ( double *  data,
double *  rhs,
VectorDouble work,
double *  z 
)
static

Perform the Kriging operation (multivariate and/or multistructure)

Returns
Error return code
Parameters
[in]dataVector of data
[out]rhsR.H.S.
[out]workWorking array
[out]zOutput Array
Remarks
This code is only provided for the Cholesky solver

◆ st_kriging_several_loop()

static int st_kriging_several_loop ( int  flag_crit,
VectorDouble work,
double *  rhscur,
double *  rhsloc,
double *  xcur,
const double *  rhs,
double *  crit 
)
static

Perform the Loop for the Kriging operation (several)

Returns
Error return code
Parameters
[in]flag_crit1 if used for criterion evaluation; 0 when used for calculation loop
[out]workWorking array
[out]rhscurWorking array
[out]rhslocWorking array
[in,out]xcurCurrent solution vector
[in,out]rhsR.H.S.
[in,out]critCriterion

◆ st_kriging_several_results()

static int st_kriging_several_results ( const double *  xcur,
double *  z 
)
static

Perform the final reconstitution for multivariate or multistructure Kriging

Returns
Error return code
Parameters
[in]xcurSolution vector (Dimension: ncur * ncova * nvar)
[out]zOutput Array

◆ st_kriging_several_rhs()

static int st_kriging_several_rhs ( double *  data,
double *  rhs,
VectorDouble work,
double *  ss_arg 
)
static

Prepare the R.H.S. for the Kriging operation (several)

Returns
Error return code
Parameters
[in]dataVector of data
[out]rhsR.H.S.
[out]workWorking array
[out]ss_argGlobal returned criterion

◆ st_load_array()

static void st_load_array ( int  number,
const double *  aux,
double *  perm 
)
static

Load an auxiliary array into a permanent array: perm <- aux

Parameters
[in]number: Number of terms in arrays 'perm' and 'aux'
[in]aux: Auxiliary array
[in,out]perm: Input/Output array

◆ st_load_data()

static void st_load_data ( AMesh amesh,
Db dbin,
Db dbout,
SPDE_Option ,
int  ivar0,
double *  data,
double *  zcur 
)
static

Load the data information within the complete output array

Parameters
[in]ameshAMesh structure
[in]dbinDb input structure
[in]dboutDb output structure
[in]ivar0Rank of the variable or -1
[out]dataVector of active data
[out]zcurOutput array (Dimension: number of meshes)
Remarks
When 'ivar0' is defined, this corresponds to the flag_mult_data
flag where the target is the simulation
Otherwise the Z-variables must all be loaded

◆ st_m2d_check_pinchout()

static int st_m2d_check_pinchout ( Db dbgrid,
int  icol_pinch 
)
static

Check the pinchout variable

Returns
Error returned code
Parameters
[in]dbgridGrid structure
[in]icol_pinchPointer to the pinchout variabe

◆ st_m2d_create_constraints()

static Db* st_m2d_create_constraints ( M2D_Environ *  m2denv,
Db dbin,
Db dbout,
int  ndim,
int  nlayer,
int *  number_hard 
)
static

Create a Db containing all the constraining information

Returns
Pointer to the newly created Db or NULL
Parameters
[in]m2denvM2D_Environ structure
[in]dbinDb input structure
[in]dboutDb output structure
[in]ndimSpace dimension
[in]nlayerNumber of layers
[out]number_hardNumber of hard data which will serve for seting the optimal drift
Remarks
Note that the file constructed here contains as many samples
as the number of ACTIVE samples of the input Db

◆ st_m2d_draw_elevation()

static double st_m2d_draw_elevation ( M2D_Environ *  m2denv,
int  ,
int  ,
double  lower,
double  upper 
)
static

Get the elevation within bounds

Returns
The value assigned to this inequality
Parameters
[in]m2denvM2D_Environ structure
[in]lowerLower bound
[in]upperUpper bound

◆ st_m2d_draw_gaussian()

static double st_m2d_draw_gaussian ( M2D_Environ *  m2denv,
Db dbc,
int  verbose,
int  iter,
int  ilayer,
int  iech,
double  Zval,
double  Zcum,
double  Zmin,
double  Zmax,
double  Ymean,
double  Ysigma 
)
static

Draw a Z-value within bounds

Returns
Z-Value in the working domain
Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb structure
[in]verboseVerbose flag
[in]iterRank of the iteration
[in]ilayerRank of the layer
[in]iechRank of the sample
[in]ZvalInput value
[in]ZcumCumulated Z value of layers above
[in]ZminLower bound in Z
[in]ZmaxUpper bound in Z
[in]YmeanMean of the Y Law
[in]YsigmaStandard deviation of the Y Law

◆ st_m2d_drift_fitting()

static int st_m2d_drift_fitting ( M2D_Environ *  m2denv,
Db dbc,
int  nlayer,
int  number_hard,
int  verbose 
)
static

Fit the coefficients of the trend terms for each layer

Returns
Error returned code
Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb structure for constraints
[in]nlayerNumber of layers
[in]number_hardNumber of hard data used to fit the drift
[in]verboseVerbose flag
Remarks
The drift is only established on data where lower and upper bounds
are both defined. The drift coefficients are assumed to be the same
for all layers
The impact of areal constraints being to important, it has been
chosen to base the drift fitting only on the first 'number_hard'
samples (which correspond to constraints coming from 'dbin'.

◆ st_m2d_drift_inc_manage()

static int st_m2d_drift_inc_manage ( M2D_Environ *  m2denv,
int  mode,
int  nlayer,
int  icol_pinch,
Db dbc,
Db dbout 
)
static

Calculate and store drift value per point in constraints and output Db Check the validity of the drift at points

Returns
Error returned code
Parameters
[in]m2denvM2D_Environ structure
[in]mode1 adding; -1 deleting
[in]nlayerNumber of layers
[in]icol_pinchPointer to the pinchout variabe
[in]dbcDb constraints structure
[in]dboutDb output structure

◆ st_m2d_drift_manage()

static int st_m2d_drift_manage ( M2D_Environ *  m2denv,
Db dbin,
Db dbout,
int  nlayer,
int  verbose,
int *  iatt_f 
)
static

Fit the coefficients of the trend terms for each layer

Returns
Error returned code
Parameters
[in]m2denvM2D_Environ structure
[in]dbinDb input structure
[in]dboutDb output structure
[in]nlayerNumber of layers
[in]verboseVerbose flag
[out]iatt_fPointer in dbin to the added variables ELoc::F
Remarks
This function also add the attributes to 'dbin' per layer:
- the external drift values (ELoc::F)

◆ st_m2d_drift_save()

static void st_m2d_drift_save ( M2D_Environ *  m2denv,
Db dbout,
int  nlayer,
double *  gwork 
)
static

Save the drift at the grid nodes

Parameters
[in]m2denvM2D_Environ structure
[in]dboutDb otput structure
[in]nlayerNumber of layers
[out]gworkWorking array
Remarks
The drift returned as a surface uses directly the coefficients
of the linear combinaison.

◆ st_m2d_external_drift_increment()

static double st_m2d_external_drift_increment ( M2D_Environ *  m2denv,
Db db,
int  ilayer0,
int  iech0 
)
static

At a point, returns the external drift increment from previous layer

Returns
The external drift increment
Parameters
[in]m2denvM2D_Environ structure
[in]dbDb structure containing the constraints
[in]ilayer0Rank of the layer of interest
[in]iech0Rank of the sample of interest

◆ st_m2d_get_drift()

static double st_m2d_get_drift ( M2D_Environ *  m2denv,
Db db,
int  ilayer0,
int  iech0 
)
static

Returns the value of the drift contribution at a sample This value is a weighted combinaison of constant and external drift term

Returns
The drift interval value
Parameters
[in]m2denvM2D_Environ structure
[in]dbDb structure containing the constraints
[in]ilayer0Rank of the layer of interest
[in]iech0Rank of the sample of interest

◆ st_m2d_get_M()

static double st_m2d_get_M ( M2D_Environ *  m2denv,
Db db,
int  type,
int  ilayer,
int  iech 
)
static

Returns the value of the drift increment at a sample (mean)

Returns
The mean value or TEST value
Parameters
[in]m2denvM2D_Environ structure
[in]dbDb structure containing the constraints
[in]type1 for the constraining Db 2 for the grid output Db
[in]ilayerRank of the layer of interest
[in]iechRank of the sample of interest

◆ st_m2d_get_S()

static double st_m2d_get_S ( M2D_Environ *  m2denv,
Db ,
int  ,
int  ,
int   
)
static

Returns the value of the gaussian standard deviation

Returns
The mean value or TEST value
Parameters
[in]m2denvM2D_Environ structure

◆ st_m2d_initial_elevations()

static int st_m2d_initial_elevations ( M2D_Environ *  m2denv,
Db dbc,
int  nlayer,
VectorDouble work 
)
static

Set the initial elevations at the constraining information

Returns
Error returned code
Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb structure for constraints
[in]nlayerNumber of layers
[out]workArray of tentative values (Dimension: nlayer)
Remarks
This function also add the attributes to 'dbin' per layer:
- the initial value (ELoc::Z)

◆ st_m2d_migrate_pinch_to_point()

static int st_m2d_migrate_pinch_to_point ( Db dbout,
Db dbc,
int  icol_pinch 
)
static

Locally migrate the pinchout distance from grid to point

Returns
Address of the newly added vector in 'dbc'
Parameters
[in]dboutDb output structure
[in]dbcDb constraints structure
[in]icol_pinchPointer to the pinchout variabe

◆ st_m2d_print_environ()

static void st_m2d_print_environ ( const char *  title,
M2D_Environ *  m2denv 
)
static

Print the Environnement

◆ st_m2d_set_M()

static void st_m2d_set_M ( M2D_Environ *  m2denv,
int  nlayer,
int  icol_pinch,
Db db,
int  iatt 
)
static

Calculate the drift increment in a Db

Parameters
[in]m2denvM2D_Environ structure
[in]nlayerNumber of layers
[in]icol_pinchPointer to the pinchout variabe
[in]dbDb structure
[in]iattPointer to the drift vector

◆ st_m2d_stats_gaus()

static void st_m2d_stats_gaus ( const char *  title,
int  nlayer,
int  nech,
double *  ydat 
)
static

Print the statistics on the current array

Parameters
[in]titleTitle attache to the printou
[in]nlayerNumber of layers
[in]nechNumber of samples per layer in the target array
[in]ydatTarget (generic) array
Remarks
This function can be used for any array

◆ st_m2d_stats_init()

static void st_m2d_stats_init ( M2D_Environ *  m2denv,
Db dbin,
int  nlayer,
int  verbose 
)
static

Calculate global statistics on elevations

Parameters
[in]m2denvM2D_Environ structure
[in]dbinDb input structure
[in]nlayerNumber of layers
[in]verboseVerbose flag

◆ st_m2d_stats_updt()

static void st_m2d_stats_updt ( M2D_Environ *  m2denv,
Db dbc,
int  nlayer,
int  verbose 
)
static

Update global statistics on the raw information

Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb structure for constraints
[in]nlayerNumber of layers
[in]verboseVerbose flag

◆ st_m2d_vector_extract()

static void st_m2d_vector_extract ( M2D_Environ *  m2denv,
Db dbc,
int  nlayer,
VectorDouble ydat,
VectorDouble work 
)
static

Extract a vector containing the constraints

Parameters
[in]m2denvM2D_Environ structure
[in]dbcDb constraints structure
[in]nlayerNumber of layers
[out]ydatArray of values at constraints samples (Dimension: nech * nlayer)
[out]workArray of tentative values (Dimension: nlayer)

◆ st_matelem_manage()

static void st_matelem_manage ( int  mode)
static

Initialize one SP_Mat structure

Parameters
[in]modeType of the action 1 for allocation; 0 for partial deallocation (of current Matelem) -1 for deallocation
Remarks
This function is called when the current IGRF has been chosen

◆ st_matelem_print()

static void st_matelem_print ( int  icov)
static

Print the contents of one SP_Mat structure

Parameters
[in]icovRank of the covariance

◆ st_merge()

static void st_merge ( int  ncur,
double *  z,
double *  zperm 
)
static

Merge an auxiliary array into an permanent array

Parameters
[in]ncurNumber of elements per variable
[in]zArray of values to be merged
[out]zpermFill permanent array

◆ st_mgs_manage()

static cs_MGS* st_mgs_manage ( int  mode,
cs_MGS mgs 
)
static

Manage the Multigrid solving operations

Returns
Error return code
Parameters
[in]mode1 for allocation; -1 for deallocation
[in]mgscs_MGS to be freed (only for mode=-1)

◆ st_print_all()

static void st_print_all ( const char *  title)
static

Print the Matelem characteristics (for given GRF and COV)

Parameters
[in]titleTitle to be printed

◆ st_print_concatenate_interval()

static void st_print_concatenate_interval ( const char *  title,
double  lower,
double  upper,
int  tail 
)
static

Print (concatenate) the printout of an interval

Parameters
[in]titleOptional title
[in]lowerLower bound or FFFF
[in]upperUpper bound or FFFF
[in]tail0: blank character; 1: "\n" character
Remarks
The printed string starts with a blank character
It ends with either a blank or a <SR/LF> character (see 'tail')

◆ st_print_constraints_per_point()

static void st_print_constraints_per_point ( int  ilayer,
int  iech,
double  value,
double  drift,
double  vgaus,
double  lower,
double  upper 
)
static

Print the constraints information for a single point

Parameters
[in]ilayerRank of the layer
[in]iechRank of the sample
[in]valueCurrent value
[in]driftDrift value (or TEST)
[in]vgausCurrent Gaussian value (or TEST)
[in]lowerLower bound or FFFF
[in]upperUpper bound or FFFF

◆ st_print_db_constraints()

static void st_print_db_constraints ( const char *  title,
Db db,
const VectorDouble ydat,
int  nlayer,
int  verbose 
)
static

Print the set of constraints

Parameters
[in]titleTitle
[in]dbDb constraints structure
[in]ydatArray of gaussian values at constraints (optional)
[in]nlayerNumber of layers
[in]verboseVerbose flag
Remarks
This function tends to produce verbose outputs.
This is the reason why it has been conditioned to print only
the values of the first samples. This is controled by the
internal parameter 'nprint' which can be ruled by keypair:
set.keypair("Print_Data",0)
The default number of samples i s 0 (no printout)

◆ st_print_details()

static void st_print_details ( Db dbc,
int  nech,
int  ilayer 
)
static

Print the details of the constraints

Parameters
[in]dbcDb structure for constraints
[in]nechNumber of hard data
[in]ilayerRank of the target layer

◆ st_print_sample()

static void st_print_sample ( const char *  title,
M2D_Environ *  ,
Db dbc,
int  nlayer,
int  iech,
VectorDouble work 
)
static

Print the values at a sample location

Parameters
[in]titleTitle
[in]dbcDb structure containing the constraints
[in]nlayerNumber of layers
[in]iechSample rank
[in]workArray of values (defined in Z)

◆ st_print_status()

static void st_print_status ( int  auth)
static

Encode the status of the variable

Parameters
[in]authStatus option

◆ st_product_Q()

static void st_product_Q ( const VectorDouble blin,
MatrixSparse S,
const VectorDouble Lambda,
const VectorDouble TildeC,
const VectorDouble x,
VectorDouble y 
)
static

Perform the product of x by Q using the blin decomposition

Parameters
[in]blinArray of coefficients for Linear combinaison
[in]SShift operator
[in]LambdaVector Lambda
[in]TildeCVector TildeC
[in]xInput array
[out]yOutput array

◆ st_project_plane()

static void st_project_plane ( double  center[3],
double  axes[2][3],
double  xyz[3],
double  coeff[2] 
)
static

Project a point on the tangent plane

Parameters
[in]centerCoordinates of the reference point (Dimension: 3)
[in]axesCoordinates of the endpoints (Dimension: 2 * 3)
[in]xyzCoordinates of the target point (Dimension: 3)
[out]coeffCoordinate of point in the local system (Dimension: 2)

◆ st_qchol_filter()

static void st_qchol_filter ( const char *  title,
int  auth 
)
static

Print information on the Filter

Parameters
[in]titleOptional title
[in]authFilter option

◆ st_qchol_print()

static void st_qchol_print ( const char *  title,
QChol QC 
)
static

Print information on the Sparse matrix and its cholesky decomposition

Parameters
[in]titleOptional title
[in]QCPointer to the QChol structure

◆ st_qsimu_manage()

static QSimu* st_qsimu_manage ( int  mode,
QSimu qsimu 
)
static

Manage the QSimu structure

Returns
Pointer to the QSimu structure
Parameters
[in]modeType of operation 1 : Allocation -1 : Deallocation
[in]qsimuQSimu structure

◆ st_record_sample()

static int st_record_sample ( M2D_Environ *  m2denv,
Db db,
int  iech,
int  ndim,
int  natt,
int  nlayer,
int  bypass,
int *  number_arg,
double *  tab 
)
static

Record a new active point

Returns
Error returned code
Parameters
[in]m2denvM2D_Environ structure
[in]dbDb structure
[in]iechSample rank in 'db'
[in]ndimSpace dimension
[in]nattNumber of attributes
[in]nlayerNumber of layers
[in]bypass1 to bypass check that at least one bound is defined
[in,out]number_argNumber of samples
[in,out]tabArray of samples

◆ st_save_result()

static void st_save_result ( double *  z,
Db dbout,
const ELoc &  locatorType,
int  iatt_simu 
)
static

Copy an array into the Db

Parameters
[in]zArray of values
[in]dboutOutput Db structure
[in]locatorTypeRank of the pointer
[in]iatt_simuPointer to the current simulation

◆ st_set_filnug()

static void st_set_filnug ( int  flag_filnug)
static

Defines if a nugget effect component must be filtered

Parameters
[in]flag_filnugFlag to define if a nugget effect must be filtered

◆ st_set_model()

static void st_set_model ( Model model)
static

Set the pointer to the model of the environment

Parameters
[in]modelPointer to the Model structure
Remarks
The pointer is copied. The contents is NOT duplicated

◆ st_simu_add_vertices()

static void st_simu_add_vertices ( int  number,
const double *  zsnc,
double *  zcur 
)
static

Restores the conditional simulation array at mesh vertices

Parameters
[in]number: Number of terms in arrays 'perm' and 'aux'
[in]zsnc: Array containing the non-conditional simulation
[in,out]zcur: Array containing the simulation error in input and the conditional simulation in output

◆ st_simu_constraints()

static double st_simu_constraints ( Db db,
int  igrf,
int  iech,
int  iter,
int  ngibbs_burn,
double  yk,
double  sk 
)
static

Return the simulated value under constraints

Returns
Simulated value
Parameters
[in]dbDb structure
[in]igrfRank of the GRF
[in]iechRank of the sample
[in]iterGibbs iteration rank (starting from 0)
[in]ngibbs_burnNumber of iterations (Burning step)
[in]ykKriged value
[in]skStandard deviation of the kriged error
Remarks
In the gradual case, the constraints is calculated as a function
of the iteration rank
Otherwise the constant bounds are used

◆ st_simu_subtract_data()

static int st_simu_subtract_data ( int  ncur,
const double *  zsnc,
const double *  data,
double *  zerr 
)
static

Constitutes the array containing the simulation error

Returns
Error return code
Parameters
[in]ncurNumber of elements per variable
[in]zsncPermanent array
[in]dataArray of datas
[in]zerrArray of extracted values

◆ st_simulate()

static int st_simulate ( QChol QC,
double *  zsnc 
)
static

Perform the Simulations

Returns
Error return code
Parameters
[in]QCPointer to the QChol structure (finalized)
[out]zsncOutput simulation array (Dimension: nvertex * nvar)

◆ st_simulate_chebychev()

static int st_simulate_chebychev ( VectorDouble zsnc)
static

Perform the basic non-conditional Simulation using the Chebychev Polynomial procedure

Returns
Error return code
Parameters
[out]zsncOutput array (Dimension: nvertex)

◆ st_simulate_cholesky()

static int st_simulate_cholesky ( QChol QC,
VectorDouble work,
VectorDouble zsnc 
)
static

Perform the basic non-conditional Simulation using the Cholesky decomposition method

Returns
Error return code
Parameters
[in]QCPointer to the QChol structure (finalized)
[out]workWorking array (Dimension: nvertex)
[out]zsncOutput array (Dimension: nvertex)

◆ st_simulate_nugget()

static void st_simulate_nugget ( int  ncur,
VectorDouble zsnc 
)
static

Simulate the nugget effect component

Parameters
[in]ncurNumber of target to be simulated
[out]zsncOutput array (Dimension: nvertex)

◆ st_solve()

static int st_solve ( QChol QC,
double *  rhs,
VectorDouble work,
double *  z 
)
static

Solve the linear system (either using Cholesky or Multigrid)

Returns
Error return code
Parameters
[in]QCQchol structure
[in]rhsR.H.S.
[out]workWorking array (Dimension: ntarget)
[out]zOutput Array (Dimension: ntarget)

◆ st_tangent_calculate()

static void st_tangent_calculate ( double  center[3],
const double  srot[2],
double  axes[2][3] 
)
static

Get the coordinates of the axis endpoints in the tangent plane

Parameters
[in]centerCoordinates of the reference point (Dimension: 3)
[in]srotRotation angles on sphere (Dimension: 2)
[out]axesCoordinates of the endpoints (Dimension: 2 * 3)

◆ st_title()

static void st_title ( int  flag_igrf,
int  flag_icov,
int  rank,
const char *  title 
)
static

Update a string to include the rank of the current GRF and Covariance

Parameters
[in]flag_igrfTo add current GRF
[in]flag_icovTo add current COV
[in]rankRank of the highlight (see mestitle or -1)
[in]titleInput title

◆ st_triangle_center()

static void st_triangle_center ( AMesh amesh,
int  ncorner,
int  imesh,
double  center[3],
double  xyz[3][3] 
)
static

Get the 3-D coordinates of the center of a triangle on the sphere

Parameters
[in]ameshMeshEStandard structure
[in]ncornerNumber of vertices per element
[in]imeshRank of the current mesh
[out]centerCoordinates of the center point (Dimension: 3)
[out]xyzCoordinate of the point (Dimension: 3x3)

Variable Documentation

◆ Calcul

SPDE_Calcul Calcul
static

◆ DEBUG

int DEBUG = 0
static

◆ FACDIM

double FACDIM[] = { 0., 1., 2., 6. }
static

◆ FLAG_KEYPAIR

int FLAG_KEYPAIR = 0
static

◆ MEM_DBIN

Db* MEM_DBIN
static

◆ MEM_DBOUT

Db * MEM_DBOUT
static

◆ NAME

char NAME[100]
static

◆ S_DECIDE

SPDE_Decision S_DECIDE
static

◆ S_ENV

SPDE_Environ S_ENV
static

◆ S_EXTERNAL_A

MatrixSparse* S_EXTERNAL_A[3] = { NULL, NULL, NULL }
static

◆ S_EXTERNAL_MESH

AMesh* S_EXTERNAL_MESH[3] = { NULL, NULL, NULL }
static

◆ S_EXTERNAL_Q

MatrixSparse* S_EXTERNAL_Q[3] = { NULL, NULL, NULL }
static

◆ SPDE_CURRENT_ICOV

int SPDE_CURRENT_ICOV = 0
static

◆ SPDE_CURRENT_IGRF

int SPDE_CURRENT_IGRF = 0
static

◆ string_encode

char string_encode[100]
static

◆ VERBOSE

int VERBOSE = 0
static