![]() |
MAGiC
V5.0
Mailleurs Automatiques de Géometries intégrés à la Cao
|
Fonctions | |
| REAL | exactinit () |
| int | grow_expansion (int elen, REAL *e, REAL b, REAL *h) |
| int | grow_expansion_zeroelim (int elen, REAL *e, REAL b, REAL *h) |
| int | expansion_sum (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | expansion_sum_zeroelim1 (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | expansion_sum_zeroelim2 (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | fast_expansion_sum (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | fast_expansion_sum_zeroelim (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | linear_expansion_sum (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | linear_expansion_sum_zeroelim (int elen, REAL *e, int flen, REAL *f, REAL *h) |
| int | scale_expansion (int elen, REAL *e, REAL b, REAL *h) |
| int | scale_expansion_zeroelim (int elen, REAL *e, REAL b, REAL *h) |
| int | compress (int elen, REAL *e, REAL *h) |
| REAL | estimate (int elen, REAL *e) |
| REAL | orient2dfast (REAL *pa, REAL *pb, REAL *pc) |
| REAL | orient2dexact (REAL *pa, REAL *pb, REAL *pc) |
| REAL | orient2dslow (REAL *pa, REAL *pb, REAL *pc) |
| REAL | orient2dadapt (REAL *pa, REAL *pb, REAL *pc, REAL detsum) |
| REAL | orient2d (REAL *pa, REAL *pb, REAL *pc) |
| REAL | orient3dfast (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | orient3dexact (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | orient3dslow (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | orient3dadapt (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL permanent) |
| REAL | orient3d (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | incirclefast (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | incircleexact (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | incircleslow (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | incircleadapt (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL permanent) |
| REAL | incircle (REAL *pa, REAL *pb, REAL *pc, REAL *pd) |
| REAL | inspherefast (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe) |
| REAL | insphereexact (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe) |
| REAL | insphereslow (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe) |
| REAL | insphereadapt (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe, REAL permanent) |
| REAL | insphere (REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe) |
| double | incircle (double *pa, double *pb, double *pc, double *pd) |
| double | insphere (double *pa, double *pb, double *pc, double *pd, double *pe) |
| double | orient2d (double *pa, double *pb, double *pc) |
| double | orient3d (double *pa, double *pb, double *pc, double *pd) |
Variables | |
| static REAL | splitter |
| static REAL | epsilon |
| static REAL | resulterrbound |
| static REAL | ccwerrboundA |
| static REAL | ccwerrboundB |
| static REAL | ccwerrboundC |
| static REAL | o3derrboundA |
| static REAL | o3derrboundB |
| static REAL | o3derrboundC |
| static REAL | iccerrboundA |
| static REAL | iccerrboundB |
| static REAL | iccerrboundC |
| static REAL | isperrboundA |
| static REAL | isperrboundB |
| static REAL | isperrboundC |
Définition à la ligne 1349 du fichier robustpredicates.cc.
Références Fast_Two_Sum, INEXACT, et REAL.
Définition à la ligne 1392 du fichier robustpredicates.cc.
Références REAL.
Référencé par incircleadapt(), insphereadapt(), orient2dadapt(), et orient3dadapt().

| double robustPredicates::exactinit | ( | ) |
Définition à la ligne 664 du fichier robustpredicates.cc.
Références ccwerrboundA, ccwerrboundB, ccwerrboundC, epsilon, iccerrboundA, iccerrboundB, iccerrboundC, isperrboundA, isperrboundB, isperrboundC, o3derrboundA, o3derrboundB, o3derrboundC, REAL, resulterrbound, et splitter.
Référencé par MAILLEUR_DELAUNAY::MAILLEUR_DELAUNAY().

Définition à la ligne 966 du fichier robustpredicates.cc.
Références f(), Fast_Two_Sum, INEXACT, REAL, et Two_Sum.

| int robustPredicates::fast_expansion_sum_zeroelim | ( | int | elen, |
| REAL * | e, | ||
| int | flen, | ||
| REAL * | f, | ||
| REAL * | h | ||
| ) |
Définition à la ligne 1039 du fichier robustpredicates.cc.
Références f(), Fast_Two_Sum, INEXACT, REAL, et Two_Sum.
Référencé par incircleadapt(), incircleexact(), incircleslow(), insphereadapt(), insphereexact(), insphereslow(), orient2dadapt(), orient2dexact(), orient2dslow(), orient3dadapt(), orient3dexact(), et orient3dslow().


Définition à la ligne 742 du fichier robustpredicates.cc.
Définition à la ligne 776 du fichier robustpredicates.cc.
| double robustPredicates::incircle | ( | double * | pa, |
| double * | pb, | ||
| double * | pc, | ||
| double * | pd | ||
| ) |
Définition à la ligne 3259 du fichier robustpredicates.cc.
Références Absolute, iccerrboundA, incircleadapt(), et REAL.

Définition à la ligne 2642 du fichier robustpredicates.cc.
Références Absolute, estimate(), fast_expansion_sum_zeroelim(), iccerrboundB, iccerrboundC, INEXACT, REAL, resulterrbound, scale_expansion_zeroelim(), Square, Two_Diff_Tail, Two_Product, Two_Two_Diff, et Two_Two_Sum.
Référencé par incircle().


Définition à la ligne 2388 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, scale_expansion_zeroelim(), Two_Product, et Two_Two_Diff.

Définition à la ligne 2365 du fichier robustpredicates.cc.
Références REAL.
Définition à la ligne 2486 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, scale_expansion_zeroelim(), Two_Diff, et Two_Two_Product.

| double robustPredicates::insphere | ( | double * | pa, |
| double * | pb, | ||
| double * | pc, | ||
| double * | pd, | ||
| double * | pe | ||
| ) |
Définition à la ligne 4167 du fichier robustpredicates.cc.
Références Absolute, insphereadapt(), isperrboundA, et REAL.
Référencé par MAILLEUR2D_INS_NOEUD::inshper_point(), MAILLEUR2D_STL_REFINE_INS_NOEUD::inshper_point(), et DLY_TETRA::point_dans_la_sphere().


| REAL robustPredicates::insphereadapt | ( | REAL * | pa, |
| REAL * | pb, | ||
| REAL * | pc, | ||
| REAL * | pd, | ||
| REAL * | pe, | ||
| REAL | permanent | ||
| ) |
Définition à la ligne 3952 du fichier robustpredicates.cc.
Références Absolute, estimate(), fast_expansion_sum_zeroelim(), INEXACT, insphereexact(), isperrboundB, isperrboundC, REAL, resulterrbound, scale_expansion_zeroelim(), Two_Diff_Tail, Two_Product, et Two_Two_Diff.
Référencé par insphere().


Définition à la ligne 3371 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, scale_expansion_zeroelim(), Two_Product, et Two_Two_Diff.
Référencé par insphereadapt().


Définition à la ligne 3328 du fichier robustpredicates.cc.
Références REAL.
Définition à la ligne 3623 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, scale_expansion_zeroelim(), Two_Diff, et Two_Two_Product.

Définition à la ligne 1120 du fichier robustpredicates.cc.
Références f(), Fast_Two_Sum, INEXACT, REAL, et Two_Sum.

| int robustPredicates::linear_expansion_sum_zeroelim | ( | int | elen, |
| REAL * | e, | ||
| int | flen, | ||
| REAL * | f, | ||
| REAL * | h | ||
| ) |
Définition à la ligne 1180 du fichier robustpredicates.cc.
Références f(), Fast_Two_Sum, INEXACT, REAL, et Two_Sum.

| double robustPredicates::orient2d | ( | double * | pa, |
| double * | pb, | ||
| double * | pc | ||
| ) |
Définition à la ligne 1604 du fichier robustpredicates.cc.
Références ccwerrboundA, orient2dadapt(), et REAL.

Définition à la ligne 1524 du fichier robustpredicates.cc.
Références Absolute, ccwerrboundB, ccwerrboundC, estimate(), fast_expansion_sum_zeroelim(), INEXACT, REAL, resulterrbound, Two_Diff_Tail, Two_Product, et Two_Two_Diff.
Référencé par orient2d().


Définition à la ligne 1441 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, Two_Product, et Two_Two_Diff.

Définition à la ligne 1430 du fichier robustpredicates.cc.
Références REAL.
Définition à la ligne 1483 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, Two_Diff, et Two_Two_Product.

| double robustPredicates::orient3d | ( | double * | pa, |
| double * | pb, | ||
| double * | pc, | ||
| double * | pd | ||
| ) |
Définition à la ligne 2298 du fichier robustpredicates.cc.
Références Absolute, o3derrboundA, orient3dadapt(), et REAL.
Référencé par DLY_TETRA::point_dans_la_sphere().


Définition à la ligne 1856 du fichier robustpredicates.cc.
Références Absolute, estimate(), fast_expansion_sum_zeroelim(), INEXACT, o3derrboundB, o3derrboundC, REAL, resulterrbound, scale_expansion_zeroelim(), Two_Diff_Tail, Two_One_Product, Two_Product, et Two_Two_Diff.
Référencé par orient3d().


Définition à la ligne 1687 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, scale_expansion_zeroelim(), Two_Product, et Two_Two_Diff.

Définition à la ligne 1666 du fichier robustpredicates.cc.
Références REAL.
Définition à la ligne 1764 du fichier robustpredicates.cc.
Références fast_expansion_sum_zeroelim(), INEXACT, REAL, scale_expansion_zeroelim(), Two_Diff, et Two_Two_Product.

Définition à la ligne 1251 du fichier robustpredicates.cc.
Références INEXACT, REAL, Split, Two_Product_Presplit, et Two_Sum.
Définition à la ligne 1297 du fichier robustpredicates.cc.
Références Fast_Two_Sum, INEXACT, REAL, Split, Two_Product_Presplit, et Two_Sum.
Référencé par incircleadapt(), incircleexact(), incircleslow(), insphereadapt(), insphereexact(), insphereslow(), orient3dadapt(), orient3dexact(), et orient3dslow().

|
static |
Définition à la ligne 374 du fichier robustpredicates.cc.
Référencé par exactinit(), et orient2d().
|
static |
Définition à la ligne 374 du fichier robustpredicates.cc.
Référencé par exactinit(), et orient2dadapt().
|
static |
Définition à la ligne 374 du fichier robustpredicates.cc.
Référencé par exactinit(), et orient2dadapt().
|
static |
Définition à la ligne 371 du fichier robustpredicates.cc.
Référencé par MG_MAILLAGE_OUTILS::_FilterGeometricTolerance(), MGOPT_MVT_NORMAL::algorithme_gradient(), MG_CALCUL_FATIGUE::cycle_jump(), SLD_COURBE::deriver_seconde(), MG_MAILLAGE_OUTILS::est_delaunay_generalise(), exactinit(), MAILLEUR2D::examine_solution(), MSTRUCT_ANALYSE_ENERGIE_HILL::executer(), MSTRUCT_VES_DECOUP::generer_geometrie_virtuel(), OCC_IMPORT::importer(), SLD_IMPORT_TESSELLATION::importer_tessellation(), OCC_IMPORT::importer_triangulation_V2017(), MGOPT_POSTTRAITEMENT::lissage_chen2005(), MGOPT_POSTTRAITEMENT::lissage_chen2008(), MGOPT_POSTTRAITEMENT::lissage_McKenzie2016(), MGOPT_MVT_NORMAL::optimisation(), MGOPT_POSTTRAITEMENT::posttraite(), et MAILLEUR2D_STL_REFINE_INS_NOEUD::stlins_point_withbc().
|
static |
Définition à la ligne 376 du fichier robustpredicates.cc.
Référencé par exactinit(), et incircle().
|
static |
Définition à la ligne 376 du fichier robustpredicates.cc.
Référencé par exactinit(), et incircleadapt().
|
static |
Définition à la ligne 376 du fichier robustpredicates.cc.
Référencé par exactinit(), et incircleadapt().
|
static |
Définition à la ligne 377 du fichier robustpredicates.cc.
Référencé par exactinit(), et insphere().
|
static |
Définition à la ligne 377 du fichier robustpredicates.cc.
Référencé par exactinit(), et insphereadapt().
|
static |
Définition à la ligne 377 du fichier robustpredicates.cc.
Référencé par exactinit(), et insphereadapt().
|
static |
Définition à la ligne 375 du fichier robustpredicates.cc.
Référencé par exactinit(), et orient3d().
|
static |
Définition à la ligne 375 du fichier robustpredicates.cc.
Référencé par exactinit(), et orient3dadapt().
|
static |
Définition à la ligne 375 du fichier robustpredicates.cc.
Référencé par exactinit(), et orient3dadapt().
|
static |
Définition à la ligne 373 du fichier robustpredicates.cc.
Référencé par exactinit(), incircleadapt(), insphereadapt(), orient2dadapt(), et orient3dadapt().
|
static |
Définition à la ligne 370 du fichier robustpredicates.cc.
Référencé par exactinit().