mmg3d
mmg3d.h File Reference
#include "libmmg3d.h"
#include "libmmgcommon.h"
Include dependency graph for mmg3d.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  _MMG3D_octree_s
 
struct  _MMG3D_octree
 

Macros

#define _MMG5_RETURN_AND_FREE(mesh, met, disp, val)
 
#define _MMG5_POINT_REALLOC(mesh, sol, ip, wantedGap, law, o, tag, retval)
 
#define _MMG5_TETRA_REALLOC(mesh, jel, wantedGap, law, retval)
 
#define _MMG3D_ALPHAD   20.7846096908265
 
#define _MMG3D_LLONG   2.5
 
#define _MMG3D_LSHRT   0.3
 
#define _MMG3D_LOPTL   1.3
 
#define _MMG3D_LOPTS   0.6
 
#define _MMG3D_BADKAL   0.2
 
#define _MMG3D_NPMAX   1000000
 
#define _MMG3D_NAMAX   200000
 
#define _MMG3D_NTMAX   2000000
 
#define _MMG3D_NEMAX   6000000
 
#define _MMG3D_SHORTMAX   0x7fff
 

Typedefs

typedef struct _MMG3D_octree_s _MMG3D_octree_s
 
typedef _MMG3D_octree_MMG3D_pOctree
 

Functions

void _MMG3D_initOctree_s (_MMG3D_octree_s *q)
 
int _MMG3D_initOctree (MMG5_pMesh, _MMG3D_pOctree *q, int nv)
 
void _MMG3D_freeOctree_s (MMG5_pMesh, _MMG3D_octree_s *q, int nv)
 
void _MMG3D_freeOctree (MMG5_pMesh, _MMG3D_octree **q)
 
int _MMG3D_isCellIncluded (double *cellCenter, double l, double *zoneCenter, double l0)
 
void _MMG3D_placeInListDouble (double *, double, int, int)
 
void _MMG3D_placeInListOctree (_MMG3D_octree_s **, _MMG3D_octree_s *, int, int)
 
int _MMG3D_seekIndex (double *distList, double dist, int indexMin, int indexMax)
 
int _MMG3D_intersectRect (double *rectin, double *rectinout)
 
int _MMG3D_getListSquareRec (_MMG3D_octree_s *, double *, double *, _MMG3D_octree_s ***, double *, double *, double, int, int, int *)
 
int _MMG3D_getListSquare (MMG5_pMesh, double *, _MMG3D_octree *, double *, _MMG3D_octree_s ***)
 
int _MMG3D_addOctreeRec (MMG5_pMesh, _MMG3D_octree_s *, double *, const int, int)
 
int _MMG3D_addOctree (MMG5_pMesh mesh, _MMG3D_octree *q, const int no)
 
int _MMG3D_delOctreeVertex (MMG5_pMesh, _MMG3D_octree_s *q, int no)
 
int _MMG3D_moveOctree (MMG5_pMesh, _MMG3D_pOctree, int, double *, double *)
 
void _MMG3D_mergeBranchesRec (_MMG3D_octree_s *, _MMG3D_octree_s *, int, int, int *)
 
void _MMG3D_mergeBranches (MMG5_pMesh mesh, _MMG3D_octree_s *q, int dim, int nv)
 
int _MMG3D_delOctreeRec (MMG5_pMesh, _MMG3D_octree_s *, double *, const int, const int)
 
int _MMG3D_delOctree (MMG5_pMesh mesh, _MMG3D_pOctree q, const int no)
 
void _MMG3D_printArbreDepth (_MMG3D_octree_s *q, int depth, int nv, int dim)
 
void _MMG3D_printArbre (_MMG3D_octree *q)
 
void _MMG3D_sizeArbreRec (_MMG3D_octree_s *q, int nv, int dim, int *, int *)
 
int * _MMG3D_sizeArbre (_MMG3D_octree *q, int dim)
 
int _MMG3D_octreein_iso (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int, double)
 
int _MMG3D_octreein_ani (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int, double)
 
int _MMG3D_tetraQual (MMG5_pMesh mesh, MMG5_pSol met, char metRidTyp)
 
void _MMG3D_solTruncature (MMG5_pMesh mesh, MMG5_pSol met)
 
int _MMG5_directsurfball (MMG5_pMesh mesh, int ip, int *list, int ilist, double n[3])
 
int _MMG3D_Init_mesh_var (va_list argptr)
 
int _MMG3D_Free_all_var (va_list argptr)
 
int _MMG3D_Free_structures_var (va_list argptr)
 
int _MMG3D_Free_names_var (va_list argptr)
 
int _MMG3D_newPt (MMG5_pMesh mesh, double c[3], int16_t tag)
 
int _MMG3D_newElt (MMG5_pMesh mesh)
 
int _MMG3D_delElt (MMG5_pMesh mesh, int iel)
 
void _MMG3D_delPt (MMG5_pMesh mesh, int ip)
 
int _MMG3D_zaldy (MMG5_pMesh mesh)
 
void _MMG5_freeXTets (MMG5_pMesh mesh)
 
void _MMG5_freeXPrisms (MMG5_pMesh mesh)
 
void _MMG3D_Free_topoTables (MMG5_pMesh mesh)
 
char _MMG5_chkedg (MMG5_pMesh mesh, MMG5_pTria pt, char ori, double, double, int)
 
int _MMG5_chkBdryTria (MMG5_pMesh mesh)
 
void _MMG5_tet2tri (MMG5_pMesh mesh, int k, char ie, MMG5_Tria *ptt)
 
int _MMG5_mmg3dBezierCP (MMG5_pMesh mesh, MMG5_Tria *pt, _MMG5_pBezier pb, char ori)
 
int _MMG5_BezierTgt (double c1[3], double c2[3], double n1[3], double n2[3], double t1[3], double t2[3])
 
double _MMG5_BezierGeod (double c1[3], double c2[3], double t1[3], double t2[3])
 
int _MMG3D_bezierInt (_MMG5_pBezier pb, double uv[2], double o[3], double no[3], double to[3])
 
int _MMG5_BezierReg (MMG5_pMesh mesh, int ip0, int ip1, double s, double v[3], double *o, double *no)
 
int _MMG5_BezierRef (MMG5_pMesh mesh, int ip0, int ip1, double s, double *o, double *no, double *to)
 
int _MMG5_BezierEdge (MMG5_pMesh mesh, int ip0, int ip1, double b0[3], double b1[3], char isrid, double v[3])
 
int _MMG5_BezierRidge (MMG5_pMesh mesh, int ip0, int ip1, double s, double *o, double *no1, double *no2, double *to)
 
int _MMG5_BezierNom (MMG5_pMesh mesh, int ip0, int ip1, double s, double *o, double *no, double *to)
 
int _MMG5_norface (MMG5_pMesh mesh, int k, int iface, double v[3])
 
int _MMG5_boulernm (MMG5_pMesh mesh, int start, int ip, int *ng, int *nr)
 
int _MMG5_boulenm (MMG5_pMesh mesh, int start, int ip, int iface, double n[3], double t[3])
 
int _MMG5_boulevolp (MMG5_pMesh mesh, int start, int ip, int *list)
 
int _MMG5_boulesurfvolp (MMG5_pMesh mesh, int start, int ip, int iface, int *listv, int *ilistv, int *lists, int *ilists, int isnm)
 
int _MMG5_bouletrid (MMG5_pMesh, int, int, int, int *, int *, int *, int *, int *, int *)
 
int _MMG5_startedgsurfball (MMG5_pMesh mesh, int nump, int numq, int *list, int ilist)
 
int _MMG5_srcbdy (MMG5_pMesh mesh, int start, int ia)
 
int _MMG5_coquil (MMG5_pMesh mesh, int start, int ia, int *list)
 
int _MMG5_coquilface (MMG5_pMesh mesh, int start, char iface, int, int *, int *, int *, int)
 
int _MMG3D_coquilFaceFirstLoop (MMG5_pMesh mesh, int start, int na, int nb, char iface, char ia, int *list, int *ilist, int *it1, int *it2, int *piv, int *adj, char *hasadja, int *nbdy, int silent)
 
void _MMG3D_coquilFaceSecondLoopInit (MMG5_pMesh mesh, int piv, char *iface, int *i, int *list, int *ilist, int *it1, int *pradj, int *adj)
 
void _MMG5_coquilFaceErrorMessage (MMG5_pMesh mesh, int k1, int k2)
 
int16_t _MMG5_coquilTravel (MMG5_pMesh, int, int, int *, int *, char *, int *)
 
void _MMG5_openCoquilTravel (MMG5_pMesh, int, int, int *, int *, char *, int *)
 
int _MMG5_settag (MMG5_pMesh, int, int, int16_t, int)
 
int _MMG5_deltag (MMG5_pMesh, int, int, int16_t)
 
int _MMG5_setNmTag (MMG5_pMesh mesh, _MMG5_Hash *hash)
 
int _MMG5_chkcol_int (MMG5_pMesh, MMG5_pSol, int, char, char, int *, int, char)
 
int _MMG5_chkcol_bdy (MMG5_pMesh, MMG5_pSol, int, char, char, int *, int, int *, int, char)
 
int _MMG5_chkmanicoll (MMG5_pMesh, int, int, int, int, int, char, char)
 
int _MMG5_chkmani (MMG5_pMesh mesh)
 
int _MMG5_colver (MMG5_pMesh, MMG5_pSol, int *, int, char, char)
 
int _MMG3D_analys (MMG5_pMesh mesh)
 
int _MMG3D_hashTria (MMG5_pMesh mesh, _MMG5_Hash *)
 
int MMG3D_hashPrism (MMG5_pMesh mesh)
 
int _MMG5_hashPop (_MMG5_Hash *hash, int a, int b)
 
int _MMG5_hPop (MMG5_HGeom *hash, int a, int b, int *ref, int16_t *tag)
 
int _MMG5_hTag (MMG5_HGeom *hash, int a, int b, int ref, int16_t tag)
 
int _MMG5_hGet (MMG5_HGeom *hash, int a, int b, int *ref, int16_t *tag)
 
int _MMG5_hEdge (MMG5_pMesh mesh, MMG5_HGeom *hash, int a, int b, int ref, int16_t tag)
 
int _MMG5_hNew (MMG5_pMesh mesh, MMG5_HGeom *hash, int hsiz, int hmax)
 
int _MMG5_hGeom (MMG5_pMesh mesh)
 
int _MMG5_bdryIso (MMG5_pMesh)
 
int _MMG5_bdrySet (MMG5_pMesh)
 
int _MMG5_bdryUpdate (MMG5_pMesh)
 
int _MMG5_bdryPerm (MMG5_pMesh)
 
int _MMG5_chkfemtopo (MMG5_pMesh mesh)
 
int _MMG5_cntbdypt (MMG5_pMesh mesh, int nump)
 
long long _MMG5_memSize (void)
 
int _MMG3D_memOption (MMG5_pMesh mesh)
 
int _MMG5_mmg3d1_pattern (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_mmg3d1_delone (MMG5_pMesh, MMG5_pSol)
 
int _MMG3D_mmg3d2 (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_mmg3dChkmsh (MMG5_pMesh, int, int)
 
int _MMG3D_split1_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split1 (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char metRidTyp)
 
int _MMG5_split1b (MMG5_pMesh, MMG5_pSol, int *, int, int, int, char, char)
 
int _MMG5_splitedg (MMG5_pMesh mesh, MMG5_pSol met, int iel, int iar, double crit)
 
int _MMG3D_split2sf_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split2sf (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split2_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split2 (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split3_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split3 (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split3cone_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split3cone (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split3op_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split3op (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split4sf_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split4sf (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split4op_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split4op (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split5_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split5 (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG3D_split6_sim (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6])
 
int _MMG5_split6 (MMG5_pMesh mesh, MMG5_pSol met, int k, int vx[6], char)
 
int _MMG5_split4bar (MMG5_pMesh mesh, MMG5_pSol met, int k, char)
 
int _MMG3D_simbulgept (MMG5_pMesh mesh, MMG5_pSol met, int *list, int ilist, int)
 
void _MMG5_nsort (int, double *, char *)
 
int _MMG3D_optlap (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_movintpt_iso (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int)
 
int _MMG3D_movnormal_iso (MMG5_pMesh, MMG5_pSol, int, int)
 
int _MMG5_movintptLES_iso (MMG5_pMesh mesh, MMG5_pSol met, _MMG3D_pOctree, int *, int, int)
 
int _MMG5_movintpt_ani (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int)
 
int _MMG5_movbdyregpt_iso (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int, int)
 
int _MMG5_movbdyregpt_ani (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int, int)
 
int _MMG5_movbdyrefpt_iso (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int _MMG5_movbdyrefpt_ani (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int _MMG5_movbdynompt_iso (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int _MMG5_movbdynompt_ani (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int _MMG5_movbdyridpt_iso (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int _MMG5_movbdyridpt_ani (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int _MMG3D_movv_ani (MMG5_pMesh, MMG5_pSol, int, int)
 
int _MMG3D_movv_iso (MMG5_pMesh, MMG5_pSol, int, int)
 
int _MMG3D_normalAdjaTri (MMG5_pMesh, int, char, int, double n[3])
 
int _MMG5_chkswpbdy (MMG5_pMesh, MMG5_pSol, int *, int, int, int, char)
 
int _MMG5_swpbdy (MMG5_pMesh, MMG5_pSol, int *, int, int, _MMG3D_pOctree, char)
 
int _MMG5_swpgen (MMG5_pMesh, MMG5_pSol, int, int, int *, _MMG3D_pOctree, char)
 
int _MMG5_chkswpgen (MMG5_pMesh, MMG5_pSol, int, int, int *, int *, double, char)
 
int MMG3D_swap23 (MMG5_pMesh mesh, MMG5_pSol met, int k, char metRidTyp)
 
int _MMG5_srcface (MMG5_pMesh mesh, int n0, int n1, int n2)
 
int _MMG5_chkptonbdy (MMG5_pMesh, int)
 
double _MMG5_orcal_poi (double a[3], double b[3], double c[3], double d[3])
 
int _MMG5_countelt (MMG5_pMesh mesh, MMG5_pSol sol, double *weightelt, long *npcible)
 
int MMG3D_opttyp (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree)
 
int _MMG3D_swpItem (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int, int)
 
int _MMG3D_splitItem (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int, int, double)
 
int MMG3D_optbdry (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int)
 
int MMG3D_movetetrapoints (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int)
 
int _MMG5_trydisp (MMG5_pMesh, double *, short)
 
int _MMG5_dichodisp (MMG5_pMesh, double *)
 
int _MMG5_lapantilap (MMG5_pMesh, double *)
 
int _MMG5_ppgdisp (MMG5_pMesh, double *)
 
int _MMG5_denoisbdy (MMG5_pMesh)
 
int _MMG3D_inqua (MMG5_pMesh mesh, MMG5_pSol met)
 
int _MMG3D_outqua (MMG5_pMesh mesh, MMG5_pSol met)
 
int _MMG3D_prilen (MMG5_pMesh mesh, MMG5_pSol met, char)
 
void _MMG3D_solTruncatureForOptim (MMG5_pMesh mesh, MMG5_pSol met)
 
void _MMG5_defaultValues (MMG5_pMesh)
 
int _MMG5_intridmet (MMG5_pMesh, MMG5_pSol, int, int, double, double *, double *)
 
int _MMG5_intregmet (MMG5_pMesh, MMG5_pSol, int, char, double, double *)
 
int _MMG5_intvolmet (MMG5_pMesh, MMG5_pSol, int, char, double, double *)
 
int _MMG3D_localParamReg (MMG5_pMesh, int, int *, int, int *, int, double *, double *, double *)
 
int _MMG3D_localParamNm (MMG5_pMesh, int, int, int, double *, double *, double *)
 
int _MMG3D_packMesh (MMG5_pMesh, MMG5_pSol, MMG5_pSol)
 
int _MMG3D_bdryBuild (MMG5_pMesh)
 
int _MMG3D_indElt (MMG5_pMesh mesh, int kel)
 
int _MMG3D_indPt (MMG5_pMesh mesh, int kp)
 
void _MMG5_printTetra (MMG5_pMesh mesh, char *fileName)
 
int _MMG5_meancur (MMG5_pMesh mesh, int np, double c[3], int ilist, int *list, double h[3])
 
double _MMG5_surftri (MMG5_pMesh, int, int)
 
double _MMG5_timestepMCF (MMG5_pMesh, double)
 
int _MMG5_bdyMCF (MMG5_pMesh)
 
double _MMG5_volint (MMG5_pMesh)
 
double _MMG5_estavglen (MMG5_pMesh)
 
int _MMG5_stiffelt (MMG5_pMesh, int, double *, double *)
 
int _MMG5_mmg3d3 (MMG5_pMesh, MMG5_pSol, MMG5_pSol)
 
int _MMG5_velextLS (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_saveDisp (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_delone (MMG5_pMesh mesh, MMG5_pSol sol, int ip, int *list, int ilist)
 
int _MMG5_cavity_iso (MMG5_pMesh mesh, MMG5_pSol sol, int iel, int ip, int *list, int lon, double volmin)
 
int _MMG5_cavity_ani (MMG5_pMesh mesh, MMG5_pSol sol, int iel, int ip, int *list, int lon, double volmin)
 
int _MMG5_cenrad_iso (MMG5_pMesh mesh, double *ct, double *c, double *rad)
 
int _MMG5_cenrad_ani (MMG5_pMesh mesh, double *ct, double *m, double *c, double *rad)
 
int _MMG3D_dichoto (MMG5_pMesh mesh, MMG5_pSol met, int k, int *vx)
 
int _MMG3D_dichoto1b (MMG5_pMesh mesh, MMG5_pSol met, int *list, int ret, int)
 
char _MMG5_chkedg (MMG5_pMesh mesh, MMG5_Tria *pt, char ori, double, double, int)
 
int _MMG5_anatet (MMG5_pMesh mesh, MMG5_pSol met, char typchk, int patternMode)
 
int _MMG5_movtet (MMG5_pMesh mesh, MMG5_pSol met, _MMG3D_pOctree octree, double clickSurf, double clickVol, int moveVol, int improveSurf, int improveVolSurf, int improveVol, int maxit)
 
int _MMG5_swpmsh (MMG5_pMesh mesh, MMG5_pSol met, _MMG3D_pOctree octree, int)
 
int _MMG5_swptet (MMG5_pMesh mesh, MMG5_pSol met, double, double, _MMG3D_pOctree, int)
 
void _MMG5_Init_parameters (MMG5_pMesh mesh)
 
double _MMG5_caltet33_ani (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)
 
double _MMG5_lenedgCoor_iso (double *, double *, double *, double *)
 Compute edge length from edge's coordinates. More...
 
int _MMG5_intmet_iso (MMG5_pMesh, MMG5_pSol, int, char, int, double)
 
int _MMG5_intmet_ani (MMG5_pMesh, MMG5_pSol, int, char, int, double)
 
int _MMG3D_intmet33_ani (MMG5_pMesh, MMG5_pSol, int, char, int, double)
 
int _MMG5_interp4bar_ani (MMG5_pMesh, MMG5_pSol, int, int, double *)
 
int _MMG5_interp4bar33_ani (MMG5_pMesh, MMG5_pSol, int, int, double *)
 
int _MMG5_interp4bar_iso (MMG5_pMesh, MMG5_pSol, int, int, double *)
 
int _MMG3D_defsiz_iso (MMG5_pMesh, MMG5_pSol)
 
int _MMG3D_defsiz_ani (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_gradsiz_iso (MMG5_pMesh, MMG5_pSol)
 
int _MMG5_gradsiz_ani (MMG5_pMesh, MMG5_pSol)
 
double _MMG5_meansizreg_iso (MMG5_pMesh, MMG5_pSol, int, int *, int, double, double)
 
int _MMG3D_chk4ridVertices (MMG5_pMesh mesh, MMG5_pTetra pt)
 
int _MMG5_moymet (MMG5_pMesh, MMG5_pSol, MMG5_pTetra, double *)
 
static void _MMG5_warnOrientation (MMG5_pMesh mesh)
 
static void _MMG3D_Set_commonFunc ()
 

Variables

static const unsigned char _MMG5_inxt3 [7] = { 1,2,3,0,1,2,3 }
 next vertex of tetra: {1,2,3,0,1,2,3} More...
 
static const unsigned char _MMG5_iprv3 [7] = { 3,0,1,2,3,0,1 }
 previous vertex of tetra: {3,0,1,2,3,0,1} More...
 
static const unsigned char _MMG5_idir [4][3] = { {1,2,3}, {0,3,2}, {0,1,3}, {0,2,1} }
 idir[i]: vertices of face opposite to vertex i More...
 
static const char _MMG5_idirinv [4][4] = {{-1,0,1,2},{0,-1,2,1},{0,1,-1,2},{0,2,1,-1}}
 
static const unsigned char _MMG5_iarf [4][3] = { {5,4,3}, {5,1,2}, {4,2,0}, {3,0,1} }
 iarf[i]: edges of face opposite to vertex i More...
 
static const unsigned char _MMG5_iarfinv [4][6] = { {-1,-1,-1,2,1,0}, {-1,1,2,-1,-1,0},{2,-1,1,-1,0,-1},{1,2,-1,0,-1,-1}}
 num of the j^th edge in the i^th face More...
 
static const unsigned char _MMG5_iare [6][2] = { {0,1}, {0,2}, {0,3}, {1,2}, {1,3}, {2,3} }
 vertices of extremities of the edges of the tetra More...
 
static const unsigned char _MMG5_ifar [6][2] = { {2,3}, {1,3}, {1,2}, {0,3}, {0,2}, {0,1} }
 ifar[i][]: faces sharing the ith edge of the tetra More...
 
static const unsigned char _MMG5_isar [6][2] = { {2,3}, {3,1}, {1,2}, {0,3}, {2,0}, {0,1} }
 isar[i][]: vertices of extremities of the edge opposite to the ith edge More...
 
static const unsigned char _MMG5_arpt [4][3] = { {0,1,2}, {0,4,3}, {1,3,5}, {2,5,4} }
 arpt[i]: edges passing through vertex i More...
 
static const unsigned char _MMG5_idir_pr [5][4] = { {0,1,2,0},{3,5,4,3},{1,4,5,2},{0,2,5,3},{0,3,4,1} }
 idir[i]: vertices of face i for a prism More...
 
static const unsigned char _MMG5_iarf_pr [5][5] = { {0,1,3,0}, {6,8,7,6}, {3,5,8,4}, {5,1,2,7},{0,4,6,2} }
 iarf[i]: edges of face i for a prism More...
 
static const unsigned char MMG5_permedge [12][6]
 
double(* _MMG5_lenedg )(MMG5_pMesh, MMG5_pSol, int, MMG5_pTetra)
 
double(* _MMG5_lenedgspl )(MMG5_pMesh, MMG5_pSol, int, MMG5_pTetra)
 
double(* _MMG5_caltet )(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)
 
double(* _MMG5_caltri )(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
 
int(* _MMG5_defsiz )(MMG5_pMesh, MMG5_pSol)
 
int(* _MMG5_gradsiz )(MMG5_pMesh, MMG5_pSol)
 
int(* _MMG5_intmet )(MMG5_pMesh, MMG5_pSol, int, char, int, double)
 
int(* _MMG5_interp4bar )(MMG5_pMesh, MMG5_pSol, int, int, double *)
 
int(* _MMG5_movintpt )(MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int)
 
int(* _MMG5_movbdyregpt )(MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int, int)
 
int(* _MMG5_movbdyrefpt )(MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int(* _MMG5_movbdynompt )(MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int(* _MMG5_movbdyridpt )(MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int *, int, int *, int, int)
 
int(* _MMG5_cavity )(MMG5_pMesh, MMG5_pSol, int, int, int *, int, double)
 
int(* _MMG3D_octreein )(MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree, int, double)
 

Macro Definition Documentation

◆ _MMG3D_ALPHAD

#define _MMG3D_ALPHAD   20.7846096908265

◆ _MMG3D_BADKAL

#define _MMG3D_BADKAL   0.2

◆ _MMG3D_LLONG

#define _MMG3D_LLONG   2.5

◆ _MMG3D_LOPTL

#define _MMG3D_LOPTL   1.3

◆ _MMG3D_LOPTS

#define _MMG3D_LOPTS   0.6

◆ _MMG3D_LSHRT

#define _MMG3D_LSHRT   0.3

◆ _MMG3D_NAMAX

#define _MMG3D_NAMAX   200000

◆ _MMG3D_NEMAX

#define _MMG3D_NEMAX   6000000

◆ _MMG3D_NPMAX

#define _MMG3D_NPMAX   1000000

◆ _MMG3D_NTMAX

#define _MMG3D_NTMAX   2000000

◆ _MMG3D_SHORTMAX

#define _MMG3D_SHORTMAX   0x7fff

◆ _MMG5_POINT_REALLOC

#define _MMG5_POINT_REALLOC (   mesh,
  sol,
  ip,
  wantedGap,
  law,
  o,
  tag,
  retval 
)
Value:
do \
{ \
int klink; \
"larger point table",law,retval); \
\
mesh->npnil = mesh->np+1; \
for (klink=mesh->npnil; klink<mesh->npmax-1; klink++) \
mesh->point[klink].tmp = klink+1; \
\
/* solution */ \
if ( sol->m ) { \
_MMG5_ADD_MEM(mesh,(sol->size*(mesh->npmax-sol->npmax))*sizeof(double), \
"larger solution",law); \
_MMG5_SAFE_REALLOC(sol->m,sol->size*(mesh->npmax+1), \
double,"larger solution",retval); \
} \
sol->npmax = mesh->npmax; \
\
/* We try again to add the point */ \
ip = _MMG3D_newPt(mesh,o,tag); \
if ( !ip ) {law;} \
}while(0)
int _MMG3D_newPt(MMG5_pMesh mesh, double c[3], int16_t tag)
Definition: zaldy_3d.c:39
MMG5_pPoint point
Definition: libmmgtypes.h:505
#define _MMG5_TAB_RECALLOC(mesh, ptr, initSize, wantedGap, type, message, law, retval)
Definition: mmgcommon.h:254
int npmax
Definition: libmmgtypes.h:480
int npnil
Definition: libmmgtypes.h:490
! int npmax
Definition: libmmgtypesf.h:530
int np
Definition: libmmgtypes.h:480
MMG5_pMesh * mesh
Definition: API_functionsf_3d.c:66
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:205
! int16_t tag
Definition: libmmgtypesf.h:252
MMG5_pMesh char int int * retval
Definition: API_functionsf_3d.c:775

Reallocation of point table and sol table and creation of point ip with coordinates o and tag tag

◆ _MMG5_RETURN_AND_FREE

#define _MMG5_RETURN_AND_FREE (   mesh,
  met,
  disp,
  val 
)
Value:
do \
{ \
MMG5_ARG_ppDisp,&disp, \
MMG5_ARG_end) ) { \
return(MMG5_LOWFAILURE); \
} \
return(val); \
}while(0)
#define MMG5_ARG_start
Definition: libmmgtypes.h:73
double val
Definition: mmgcommon.h:428
#define MMG5_ARG_ppDisp
Definition: libmmgtypes.h:112
int MMG3D_Free_all(const int starter,...)
Definition: API_functions_3d.c:1976
MMG5_pMesh * mesh
Definition: API_functionsf_3d.c:66
#define MMG5_LOWFAILURE
Definition: libmmgtypes.h:48
#define MMG5_ARG_end
Definition: libmmgtypes.h:159
#define MMG5_ARG_ppMet
Definition: libmmgtypes.h:102
#define MMG5_ARG_ppMesh
Definition: libmmgtypes.h:82

Free allocated pointers of mesh and sol structure and return value val

◆ _MMG5_TETRA_REALLOC

#define _MMG5_TETRA_REALLOC (   mesh,
  jel,
  wantedGap,
  law,
  retval 
)
Value:
do \
{ \
int klink,oldSiz; \
\
oldSiz = mesh->nemax; \
_MMG5_TAB_RECALLOC(mesh,mesh->tetra,mesh->nemax,wantedGap,MMG5_Tetra, \
"larger tetra table",law,retval); \
\
mesh->nenil = mesh->ne+1; \
for (klink=mesh->nenil; klink<mesh->nemax-1; klink++) \
mesh->tetra[klink].v[3] = klink+1; \
if ( mesh->adja ) { \
/* adja table */ \
_MMG5_ADD_MEM(mesh,4*(mesh->nemax-oldSiz)*sizeof(int), \
"larger adja table",law); \
_MMG5_SAFE_RECALLOC(mesh->adja,4*mesh->ne+5,4*mesh->nemax+5,int \
,"larger adja table",retval); \
} \
\
/* We try again to add the point */ \
jel = _MMG3D_newElt(mesh); \
if ( !jel ) {law;} \
}while(0)
int _MMG3D_newElt(MMG5_pMesh mesh)
Definition: zaldy_3d.c:94
if(!ier) exit(EXIT_FAILURE)
Definition: libmmgtypes.h:330
MMG5_pTetra tetra
Definition: libmmgtypes.h:507
int nenil
Definition: libmmgtypes.h:491
MMG5_pMesh * mesh
Definition: API_functionsf_3d.c:66
int nemax
Definition: libmmgtypes.h:480
int * adja
Definition: libmmgtypes.h:493
int ne
Definition: libmmgtypes.h:480
MMG5_pMesh char int int * retval
Definition: API_functionsf_3d.c:775
! int nemax
Definition: libmmgtypesf.h:530

Reallocation of tetra table and creation of tetra jel

Typedef Documentation

◆ _MMG3D_octree_s

Octree cell.

◆ _MMG3D_pOctree

Function Documentation

◆ _MMG3D_addOctree()

int _MMG3D_addOctree ( MMG5_pMesh  mesh,
_MMG3D_octree q,
const int  no 
)
Here is the caller graph for this function:

◆ _MMG3D_addOctreeRec()

int _MMG3D_addOctreeRec ( MMG5_pMesh  mesh,
_MMG3D_octree_s q,
double *  ver,
const int  no,
int  nv 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward an octree cell.
ververtex coordinates scaled such that the quadrant is [0;1]x[0;1]x[0;1]
novertex index in the mesh.
nvmaximum number of points in an octree cell.
Returns
1 if ok 0 if memory saturated

Add vertex in the suitable quadrant of the octree. This function is recursively called until we reach the last one. At each step, the vertex coordinates are scaled such as the quadrant is the [0;1]x[0;1]x[0;1] box.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_analys()

int _MMG3D_analys ( MMG5_pMesh  mesh)

preprocessing stage: mesh analysis

— stage 1: data structures for surface
— stage 2: surface analysis

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_bdryBuild()

int _MMG3D_bdryBuild ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure (unused).
Returns
-1 if fail, the number of detected ridges otherwise

Create the boundary entities of the mesh (triangles and edges).

Warning
mesh must be packed and hashed
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_bezierInt()

int _MMG3D_bezierInt ( _MMG5_pBezier  pb,
double  uv[2],
double  o[3],
double  no[3],
double  to[3] 
)
Parameters
pbpointer toward the Bezier structure.
uvcoordinates of the point in the parametric space.
ocomputed coordinates of the point in the real space.
nocomputed normal.
tocomputed tangent.
Returns
1.

Compute o, no and to at $(u,v)$ in Bezier patch.

Here is the caller graph for this function:

◆ _MMG3D_chk4ridVertices()

int _MMG3D_chk4ridVertices ( MMG5_pMesh  mesh,
MMG5_pTetra  pt 
)
Here is the caller graph for this function:

◆ _MMG3D_coquilFaceFirstLoop()

int _MMG3D_coquilFaceFirstLoop ( MMG5_pMesh  mesh,
int  start,
int  na,
int  nb,
char  iface,
char  ia,
int *  list,
int *  ilist,
int *  it1,
int *  it2,
int *  piv,
int *  adj,
char *  hasadja,
int *  nbdy,
int  silent 
)
Parameters
meshpointer toward the mesh structure.
startindex of the starting tetrahedron.
naglobal index of the 1st extremity of the edge whose shell is computed
nbglobal index of the 2d extremity of the edge whose shell is computed
ifaceindex of the face from which we come.
iaindex of edge whose shell is computed (in tetra).
listpointer toward the list of tetra in the shell (to fill).
ilistpointer toward the number of tetra in the shell (to fill).
it1pointer toward the index of the 1st boundary face sharing ia
it2pointer toward the index of the 2d boundary face sharing ia (to fill).
adjpointer toward the adjacent to treat in the shell (to update)
hasadjapointer toward 0 if we don't have adja through iface, 0 otherwise (to fill)
nbdypointer toward the number of boundaries found minus 1 (to update)
silentif 1, print error message for more than 2 boundary triangles in the shell
Returns
-1 if fail, 1 otherwise

Travel in the shell of the edge until meeting the first tetra or reaching a tetra without adjacent. Fill it2 and list.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_coquilFaceSecondLoopInit()

void _MMG3D_coquilFaceSecondLoopInit ( MMG5_pMesh  mesh,
int  piv,
char *  iface,
int *  i,
int *  list,
int *  ilist,
int *  it1,
int *  pradj,
int *  adj 
)
Parameters
meshpointer toward the mesh structure.
pivglobal index of the pivot.
ifaceindex of the face from which we come.
iindex of edge whose shell is computed (in tetra).
listpointer toward the list of tetra in the shell (to fill).
ilistpointer toward the number of tetra in the shell (to fill).
it1pointer toward the index of the 1st boundary face sharing ia
pradjpointer toward the first tetra of the shell (to fill).
adjpointer toward the adjacent to treat in the shell (to update)

Initialize the travel in the shell of the edge in reverse direction than in the coquilFaceFirstLoop function.

Here is the caller graph for this function:

◆ _MMG3D_defsiz_ani()

int _MMG3D_defsiz_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric stucture.
Returns
0 if fail, 1 otherwise.

Define size at points by intersecting the surfacic metric and the physical metric.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_defsiz_iso()

int _MMG3D_defsiz_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
0 if fail, 1 otherwise.

Define isotropic size map at all boundary vertices of the mesh, associated with geometric approx, and prescribe hmax at the internal vertices Field h of Point is used, to store the prescribed size (not inverse, squared,...)

1) Size at internal points











First step: search for local parameters










Second step: set the metric









Set size at points that cannot be reached from the tetra








First step: search for local parameters







Second step: set the metric






Set size at points that cannot be reached from the tetra





2) size at regular surface points




First step: search for local parameters



Second step: set the metric


3) Travel all boundary faces to update size prescription for points on ridges/edges

First step: search for local parameters
Second step: set metric

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_delElt()

int _MMG3D_delElt ( MMG5_pMesh  mesh,
int  iel 
)
Parameters
meshpointer toward the mesh
ielindex of the element to delete
Returns
1 if success, 0 if fail

Delete the element iel

Here is the caller graph for this function:

◆ _MMG3D_delOctree()

int _MMG3D_delOctree ( MMG5_pMesh  mesh,
_MMG3D_pOctree  q,
const int  no 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward the global octree.
noreference of the vertex to be deleted.
Returns
1 if ok 0 if memory saturated

Delete the vertex no from the octree structure.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_delOctreeRec()

int _MMG3D_delOctreeRec ( MMG5_pMesh  mesh,
_MMG3D_octree_s q,
double *  ver,
const int  no,
const int  nv 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward an octree cell.
ververtex coordinates scaled such that the quadrant is [0;1]x[0;1]x[0;1]
novertex index in the mesh.
nvmaximum number of points in an octree cell.
Returns
1 if ok 0 if memory saturated

Delete vertex no from the octree. This function is recursively called until we reach the terminal octree cell containing the vertex no. At each step, the vertex coordinates are scaled such as the quadrant is the [0;1]x[0;1]x[0;1] box.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_delOctreeVertex()

int _MMG3D_delOctreeVertex ( MMG5_pMesh  mesh,
_MMG3D_octree_s q,
int  indNo 
)
Parameters
qpointer toward a terminal octree cell (containing vertex)
noindex of the point to delete from the octree
Returns
1 if ok 0 if memory saturated

Delete the vertex of index no from the terminal octree cell, merge the cells if necessary.

Here is the caller graph for this function:

◆ _MMG3D_delPt()

void _MMG3D_delPt ( MMG5_pMesh  mesh,
int  ip 
)
Here is the caller graph for this function:

◆ _MMG3D_dichoto()

int _MMG3D_dichoto ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int *  vx 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ktetrahedron index.
vxpointer toward table of edges to split.
Returns
1 if success, 0 if fail.

Find acceptable position for splitting.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_dichoto1b()

int _MMG3D_dichoto1b ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  list,
int  ret,
int  ip 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
listpointer toward the shell of edge.
retdouble of the number of tetrahedra in the shell.
ipnew point index.
Returns
1 if success, 0 if fail

Find acceptable position for _MMG5_split1b, passing the shell of considered edge, starting from o point.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_Free_all_var()

int _MMG3D_Free_all_var ( va_list  argptr)
Parameters
argptrlist of the mmg structures that must be deallocated. Each structure must follow one of the MMG5_ARG* preprocessor variable that allow to identify it.

argptr contains at least a pointer toward a MMG5_pMesh structure (that will contain the mesh and identified by the MMG5_ARG_ppMesh keyword)

To call the MMG3D_mmg3dlib function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the ouput metric (and the input one, if provided) and identified by the MMG5_ARG_ppMet keyword).

To call the MMG3D_mmg3dls function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the level-set function and identified by the MMG5_ARG_ppLs keyword).

To call the MMG3D_mmg3dmov library, you must also provide a pointer toward a MMG5_pSol structure storing the displacement (and identified by the MMG5_ARG_ppDisp keyword).

Returns
0 if fail, 1 if success

Internal function for deallocations before return (taking a va_list as argument).

Remarks
we pass the structures by reference in order to have argument compatibility between the library call from a Fortran code and a C code.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_Free_names_var()

int _MMG3D_Free_names_var ( va_list  argptr)
Parameters
argptrlist of the mmg structures for whose we want to deallocate the name. Each structure must follow one of the MMG5_ARG* preprocessor variable that allow to identify it.

argptr contains at least a pointer toward a MMG5_pMesh structure (that will contain the mesh and identified by the MMG5_ARG_ppMesh keyword)

To call the MMG3D_mmg3dlib function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the ouput metric (and the input one, if provided) and identified by the MMG5_ARG_ppMet keyword).

To call the MMG3D_mmg3dls function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the level-set function and identified by the MMG5_ARG_ppLs keyword).

To call the MMG3D_mmg3dmov library, you must also provide a pointer toward a MMG5_pSol structure storing the displacement (and identified by the MMG5_ARG_ppDisp keyword).

Returns
0 if fail, 1 if success

Internal function for name deallocations before return (taking a va_list as argument).

Remarks
we pass the structures by reference in order to have argument compatibility between the library call from a Fortran code and a C code.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_Free_structures_var()

int _MMG3D_Free_structures_var ( va_list  argptr)
Parameters
argptrlist of the mmg structures that must be deallocated. Each structure must follow one of the MMG5_ARG* preprocessor variable that allow to identify it.

argptr contains at least a pointer toward a MMG5_pMesh structure (that will contain the mesh and identified by the MMG5_ARG_ppMesh keyword)

To call the MMG3D_mmg3dlib function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the ouput metric (and the input one, if provided) and identified by the MMG5_ARG_ppMet keyword).

To call the MMG3D_mmg3dls function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the level-set function and identified by the MMG5_ARG_ppLs keyword).

To call the MMG3D_mmg3dmov library, you must also provide a pointer toward a MMG5_pSol structure storing the displacement (and identified by the MMG5_ARG_ppDisp keyword).

Returns
0 if fail, 1 if success

Internal function for structures deallocations before return (taking a va_list as argument).

Remarks
we pass the structures by reference in order to have argument compatibility between the library call from a Fortran code and a C code.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_Free_topoTables()

void _MMG3D_Free_topoTables ( MMG5_pMesh  mesh)
inline

Free adja, xtetra and xpoint tables

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_freeOctree()

void _MMG3D_freeOctree ( MMG5_pMesh  ,
_MMG3D_octree **  q 
)
Here is the caller graph for this function:

◆ _MMG3D_freeOctree_s()

void _MMG3D_freeOctree_s ( MMG5_pMesh  mesh,
_MMG3D_octree_s q,
int  nv 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward the octree cell
nvnumber of vertices in the cell subtree

Free the octree cell.

Here is the caller graph for this function:

◆ _MMG3D_getListSquare()

int _MMG3D_getListSquare ( MMG5_pMesh  ,
double *  ,
_MMG3D_octree ,
double *  ,
_MMG3D_octree_s ***   
)

◆ _MMG3D_getListSquareRec()

int _MMG3D_getListSquareRec ( _MMG3D_octree_s q,
double *  center,
double *  rect,
_MMG3D_octree_s ***  qlist,
double *  dist,
double *  ani,
double  l0,
int  nc,
int  dim,
int *  index 
)
Parameters
qpointer toward the octree cell.
centercoordinates of the centre of the current subtree.
rectrectangle that we want to intersect with the subtree. We define it given: the coordinates of one corner of the rectange and the length of the rectangle in each dimension.
qlistpointer toward the list of pointer over the sub octrees that intersect rect.
distpointer toward the list of distances between center of the octree cells in qlist and the last 3 elements are the coordinates of the center of the whole recangle.
animetric of the point.
l0radius of the search zone.
ncnumber max of cell in the list +3 (the three last.
dimdimension =3.
indexnumber of octree cells that intersect rect
Returns
0 if the rectangle doesn't intersect the octree (possible due to the surface reconstruction), 1 otherwise.

List the number of octree cells that intersect the rectangle rect. To avoid counting of the cells, a maximum is set.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_hashTria()

int _MMG3D_hashTria ( MMG5_pMesh  mesh,
_MMG5_Hash hash 
)
Parameters
meshpointer toward the mesh structure.
hashEdges hash table.
Returns
1 if success, 0 if failed.

Create surface adjacency table. Allocate the edge hash table hash but don't free it.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_indElt()

int _MMG3D_indElt ( MMG5_pMesh  mesh,
int  kel 
)

find the element number in packed numerotation

Here is the caller graph for this function:

◆ _MMG3D_indPt()

int _MMG3D_indPt ( MMG5_pMesh  mesh,
int  kp 
)

find the point number in packed numerotation

Here is the caller graph for this function:

◆ _MMG3D_Init_mesh_var()

int _MMG3D_Init_mesh_var ( va_list  argptr)
Parameters
argptrlist of the mmg structures that must be initialized. Each structure must follow one of the MMG5_ARG* preprocessor variable that allow to identify it.

argptr contains at least a pointer toward a MMG5_pMesh structure (that will contain the mesh and identified by the MMG5_ARG_ppMesh keyword)

To call the MMG3D_mmg3dlib function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the ouput metric (and the input one, if provided) and identified by the MMG5_ARG_ppMet keyword).

To call the MMG3D_mmg3dls function, you must also provide a pointer toward a MMG5_pSol structure (that will contain the level-set function and identified by the MMG5_ARG_ppLs keyword).

To call the MMG3D_mmg3dmov library, you must also provide a pointer toward a MMG5_pSol structure storing the displacement (and identified by the MMG5_ARG_ppDisp keyword).

Returns
1 if success, 0 if fail

Internal function for structure allocations (taking a va_list argument).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_initOctree()

int _MMG3D_initOctree ( MMG5_pMesh  mesh,
_MMG3D_pOctree q,
int  nv 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward the global octree
nvmaximum number of vertices in each cell before subdivision
Returns
1 if ok 0 if memory saturated

Initialisation of the octree cell.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_initOctree_s()

void _MMG3D_initOctree_s ( _MMG3D_octree_s q)
Parameters
qpointer toward the octree cell

Initialisation of the octree cell.

Here is the caller graph for this function:

◆ _MMG3D_inqua()

int _MMG3D_inqua ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
0 if the worst element has a nul quality, 1 otherwise.

Print histogram of mesh qualities for classic storage of metric at ridges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_intersectRect()

int _MMG3D_intersectRect ( double *  rectin,
double *  rectinout 
)
Parameters
rectinrectangle to intersect, is not modified.
rectinoutrectangle to intersect, is set to the intersection.
Returns
1 if rectinout intersect rectin, 0 otherwise (possible because the surface reconstruction may leads to point outside the [0;1]x[0;1]x[0;1] bounding box)

Set rectinout to the intersection of the two rectangles. Rectangles are defined by: the coordinates of the lower left corner of the rectange and the length of the rectangle in each dimension.

Here is the caller graph for this function:

◆ _MMG3D_intmet33_ani()

int _MMG3D_intmet33_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  i,
int  ip,
double  s 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kelement index.
ilocal index of edge in k.
ipglobal index of the new point in which we want to compute the metric.
sinterpolation parameter (between 0 and 1).
Returns
0 if fail, 1 otherwise.

Interpolation of anisotropic sizemap at parameter s along edge i of elt k for a classic storage of ridges metrics (before defsiz call).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_isCellIncluded()

int _MMG3D_isCellIncluded ( double *  cellCenter,
double  l,
double *  zoneCenter,
double  l0 
)
Parameters
cellCenter3 coordinates of the center of the octree cell to test.
lsize of the cell
zoneCenter3 coordinates of the center of the search zone
radiusof the search zone
Returns
wether the cell is included in the search zone.

◆ _MMG3D_localParamNm()

int _MMG3D_localParamNm ( MMG5_pMesh  mesh,
int  iel,
int  iface,
int  ia,
double *  hausd_ip,
double *  hmin_ip,
double *  hmax_ip 
)
Parameters
meshpointer toward the mesh structure.
ielindex of tetra in which we work
ifaceindex of face in iel
iaindex of edge in iel along which we want to compute the local parameters
hausd_ippointer toward the local hausdorff parameter to compute
hmin_ippointer toward the local minimal edge size to compute
hmax_ippointer toward the local maximal edge size to compute
Returns
1 if success, 0 if fail

Compute the local parameters at non manifold point ip.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_localParamReg()

int _MMG3D_localParamReg ( MMG5_pMesh  mesh,
int  ip,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
double *  hausd_ip,
double *  hmin_ip,
double *  hmax_ip 
)
Parameters
meshpointer toward the mesh structure.
ipglobal index of point in which we want to compute the local parameters
listvpointer toward the ball of ip
ilistvnumber of tetra in the ball of ip
listspointer toward the surface ball of ip
ilistsnumber of tetra in the surface ball of ip
hausd_ippointer toward the local hausdorff parameter to compute
hmin_ippointer toward the local minimal edge size to compute
hmax_ippointer toward the local maximal edge size to compute
Returns
1 if success, 0 if fail

Compute the local parameters at point ip (the volume and surface ball of point must be provided).

Here is the caller graph for this function:

◆ _MMG3D_memOption()

int _MMG3D_memOption ( MMG5_pMesh  mesh)

memory repartition for the -m option

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_mergeBranches()

void _MMG3D_mergeBranches ( MMG5_pMesh  mesh,
_MMG3D_octree_s q,
int  dim,
int  nv 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward an octree cell.
dimdimension of the space (=3)
nvmaximum number of points in an octree cell.

Merge branches that have a parent counting less than nv vertices.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_mergeBranchesRec()

void _MMG3D_mergeBranchesRec ( _MMG3D_octree_s q0,
_MMG3D_octree_s q,
int  dim,
int  nv,
int *  index 
)
Parameters
q0pointer toward an octree cell.
qpointer toward an octree cell.
dimdimension of the space (=3).
nvmaximum number of points in an octree cell.
indexnext index in the array to be filled.

Merge sub-branches q of q0, in their parent q0. q0 should contain no more than nv vertices.

Here is the caller graph for this function:

◆ _MMG3D_mmg3d2()

int _MMG3D_mmg3d2 ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the solution structure
Returns
0 if fail, 1 otherwise.

Create implicit surface in mesh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_moveOctree()

int _MMG3D_moveOctree ( MMG5_pMesh  mesh,
_MMG3D_pOctree  q,
int  no,
double *  newVer,
double *  oldVer 
)
Parameters
meshpointer toward the mesh structure.
qpointer toward the global octree.
noindex of the moved point.
newVernew coordinates for the moved point.
oldVerold coordinates for the moved point.
Returns
1 if ok 0 if memory saturated

Move one point in the octree structure. /!\ the vertex of index no can have either the new or the old coordinates in the mesh but all other vertice should have the same coordinates as when they were inserted into the octree. (ie: one move at a time in the mesh and the octree)

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_movnormal_iso()

int _MMG3D_movnormal_iso ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
int  k,
int  ib 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the metric structure.
ktetra index.
iblocal index of the point inside the tetra k.
Returns
0 if fail, 1 if success.

Move internal point according to the normal at the opposite face Try to increase the volume of the tetra.

Remarks
the metric is not interpolated at the new position.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_movv_ani()

int _MMG3D_movv_ani ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
int   
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_movv_iso()

int _MMG3D_movv_iso ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
int   
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_newElt()

int _MMG3D_newElt ( MMG5_pMesh  mesh)

get new elt address

Here is the caller graph for this function:

◆ _MMG3D_newPt()

int _MMG3D_newPt ( MMG5_pMesh  mesh,
double  c[3],
int16_t  tag 
)

get new point address

Here is the caller graph for this function:

◆ _MMG3D_normalAdjaTri()

int _MMG3D_normalAdjaTri ( MMG5_pMesh  mesh,
int  start,
char  iface,
int  ia,
double  n[3] 
)
Parameters
meshpointer toward the mesh structure
startindex of the working tetra
ifacelocal index of the boundary face of the tetra start
ialocal index on face iface of the edge through which we seek the adjacent triangle of the triangle iface of start.
nnormal of the new boundary face in the tetra idx.
Returns
1 if success, 0 if we want to refuse the collapse, -1 if fail.

Compute the normal of the adjacent triangle of the triangle iface of the tetra start through the edge ia (in local numbering of the face).

Store the adjacent boundary triangle (triangle adjacent to iface through the edge ia
Compute the normal of the second triangle

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_octreein_ani()

int _MMG3D_octreein_ani ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
_MMG3D_pOctree  octree,
int  ip,
double  lmax 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the solution structure.
octreepointer toward the octree structure.
ipindex of point to check.
Returns
1 if we can insert ip, 0 otherwise

Check if the vertex ip is not too close from another one (for an anisotropic metric).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_octreein_iso()

int _MMG3D_octreein_iso ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
_MMG3D_pOctree  octree,
int  ip,
double  lmax 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the solution structure.
octreepointer toward the octree structure.
ipindex of point to check.
Returns
1 if we can insert ip, 0 otherwise

Check if the vertex ip is not too close from another one (for an isotropic metric).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_optlap()

int _MMG3D_optlap ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the sol structure
Returns
0 if fail, 1 otherwise.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_outqua()

int _MMG3D_outqua ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
0 if the worst element has a nul quality, 1 otherwise.

Print histogram of mesh qualities for special storage of metric at ridges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_packMesh()

int _MMG3D_packMesh ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  disp 
)
Parameters
meshpointer toward the mesh structure (unused).
metpointer toward the solution (metric or level-set) structure.
disppointer toward the solution (displacement) structure.
Returns
1 if success, 0 if chkmsh fail or if we are unable to build triangles.

Pack the sparse mesh and create triangles and edges before getting out of library

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_placeInListDouble()

void _MMG3D_placeInListDouble ( double *  distList,
double  dist,
int  index,
int  size 
)
Parameters
distListlist of values.
distvalue to insert in the list.
indexposition of the element before the place where dist should be inserted.
sizesize of the list before insertion.

Insert the value dist in the list distList at position index+1. Moves other data so nothing is lost. No memory check performed, this function should be called with coherent parameters.

Here is the caller graph for this function:

◆ _MMG3D_placeInListOctree()

void _MMG3D_placeInListOctree ( _MMG3D_octree_s **  qlist,
_MMG3D_octree_s q,
int  index,
int  size 
)
Parameters
qListlist of pointer on octree.
qpointer on octree to be inserted in the list.
indexposition of the element before the place where q should be inserted.
sizesize of the list before insertion.

Insert the pointer q in the list qList at position index+1. Moves other data so nothing is lost. No memory check performed, this function should be called with coherent parameters.

Here is the caller graph for this function:

◆ _MMG3D_prilen()

int _MMG3D_prilen ( MMG5_pMesh  mesh,
MMG5_pSol  met,
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
metRidTypType of storage of ridges metrics: 0 for classic storage, 1 for special storage.
Returns
0 if fail, 1 otherwise.

Compute sizes of edges of the mesh, and displays histo.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_printArbre()

void _MMG3D_printArbre ( _MMG3D_octree q)

◆ _MMG3D_printArbreDepth()

void _MMG3D_printArbreDepth ( _MMG3D_octree_s q,
int  depth,
int  nv,
int  dim 
)
Parameters
qpointer toward an octree cell
depthdepth of the subtree
nvnumber of vertices in the subtree
dimdimension in which we work

Print the depth depth of the subtree of q.

Warning
debug function, not safe
Here is the caller graph for this function:

◆ _MMG3D_seekIndex()

int _MMG3D_seekIndex ( double *  distList,
double  dist,
int  indexMin,
int  indexMax 
)
Parameters
distListordered list of value from smallest to largest.
distvalue to be compared to elements in the list.
indexMinminimum index of the list.
indexMaxmaximum index of the list.
Returns
Index of the biggest value of disList that is strictly smaller than dist. Only search in the bounds of indexMin and indexMax.
Here is the caller graph for this function:

◆ _MMG3D_Set_commonFunc()

static void _MMG3D_Set_commonFunc ( )
inlinestatic

Set common pointer functions between mmgs and mmg3d to the matching mmg3d functions.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_simbulgept()

int _MMG3D_simbulgept ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  list,
int  ret,
int  ip 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric.
listpointer toward the edge shell.
retsize of the edge shell.
ipnew point index.
Returns
0 if final position is invalid, 1 if all checks are ok.

Simulate at the same time creation and bulging of one point, with new position o and tag tag, to be inserted at an edge, whose shell is passed.

Check the deviation for new triangles

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_sizeArbre()

int* _MMG3D_sizeArbre ( _MMG3D_octree q,
int  dim 
)

◆ _MMG3D_sizeArbreRec()

void _MMG3D_sizeArbreRec ( _MMG3D_octree_s q,
int  nv,
int  dim,
int *  s1,
int *  s2 
)
Parameters
qpointer toward an octree cell
nvmaximum number of vertices in an octree leaf
dimdimension in which we work
ssize of the octree

Print the memory size of the octree.

Warning
debug function, not safe
Here is the caller graph for this function:

◆ _MMG3D_solTruncature()

void _MMG3D_solTruncature ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)

◆ _MMG3D_solTruncatureForOptim()

void _MMG3D_solTruncatureForOptim ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the solution structure.

Truncate the metric computed by the DoSol function by hmax and hmin values (if setted by the user). Set hmin and hmax if they are not setted.

Here is the caller graph for this function:

◆ _MMG3D_split1_sim()

int _MMG3D_split1_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if split leads to invalid situation, else 1.

Simulate the splitting of 1 edge of element

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split2_sim()

int _MMG3D_split2_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of two opposite edges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split2sf_sim()

int _MMG3D_split2sf_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of two edges that belong to a common face

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split3_sim()

int _MMG3D_split3_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)

Simulate split of 1 face (3 edges)

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split3cone_sim()

int _MMG3D_split3cone_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of 3 edges in cone configuration.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split3op_sim()

int _MMG3D_split3op_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of 3 edges in opposite configuration.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split4op_sim()

int _MMG3D_split4op_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of 4 edges in opposite configuration.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split4sf_sim()

int _MMG3D_split4sf_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of 4 edges in a configuration when 3 lie on the same face.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split5_sim()

int _MMG3D_split5_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of 5 edges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_split6_sim()

int _MMG3D_split6_sim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
Returns
0 if the split fail, 1 otherwise

Simulate split of 6 edges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_splitItem()

int _MMG3D_splitItem ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int  k,
int  iar,
double  OCRIT 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
kelt index.
iarindex of edge to split.
OCRITquality threshold.
Returns
1 if success, 0 otherwise

Try to split edge number iar of tetra k

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_swpItem()

int _MMG3D_swpItem ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int  k,
int  iar 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
kelt index.
iarindex of edge to not try to swap.
Returns
-1 if fail, 0 if we don't swap anything, 1 otherwise.

Try to swap edge iar of tetra k.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_tetraQual()

int _MMG3D_tetraQual ( MMG5_pMesh  mesh,
MMG5_pSol  met,
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the meric structure.
metRidTypmetric storage (classic or special)
Returns
1 if success, 0 if fail.

Compute the quality of the tetras over the mesh.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG3D_zaldy()

int _MMG3D_zaldy ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh
Returns
1 if success, 0 if fail

allocate main structure

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_anatet()

int _MMG5_anatet ( MMG5_pMesh  mesh,
MMG5_pSol  met,
char  typchk,
int  patternMode 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
typchktype of checking for edges length.
patternModeflag to say if we perform vertex insertion by patterns or by delaunay kernel.
Returns
0 if fail, number of new points otherwise.

Analyze tetrahedra and split if needed.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_bdryIso()

int _MMG5_bdryIso ( MMG5_pMesh  )

◆ _MMG5_bdryPerm()

int _MMG5_bdryPerm ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
0 if failed, 1 otherwise.

Make orientation of triangles compatible with tetra faces.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_bdrySet()

int _MMG5_bdrySet ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
0 if failed, 1 if success.

Set the triangles references to the tetrahedra faces and edges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_bdryUpdate()

int _MMG5_bdryUpdate ( MMG5_pMesh  mesh)

Update tag and refs of tetra edges. If tetra is required, set the faces/edges to required

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_bdyMCF()

int _MMG5_bdyMCF ( MMG5_pMesh  )

◆ _MMG5_BezierEdge()

int _MMG5_BezierEdge ( MMG5_pMesh  mesh,
int  ip0,
int  ip1,
double  b0[3],
double  b1[3],
char  ised,
double  v[3] 
)
inline
Parameters
meshpointer toward the mesh structure.
ip0index of the first point of the curve.
ip1index of the second point of the curve.
b0the first computed extrapolated control point.
b1the second computed extrapolated control point.
isedflag for special edge.
vdirection for normal vectors.

Compute control points associated to the underlying curve to $[p0;p1]$. ised = 1 if $[p0;p1]$ must be considered as a special edge. Provide a direction v which will be considered as reference when dealing with choice of normal vectors.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_BezierGeod()

double _MMG5_BezierGeod ( double  c1[3],
double  c2[3],
double  t1[3],
double  t2[3] 
)
inline
Parameters
c1coordinates of the first point of the curve.
c2coordinates of the second point of the curve.
t1normal at the first point of the curve.
t2normal at the second point of the curve.
Returns
The parameter value.

Compute value of the parameter that makes the underlying Bezier curve with 'constant speed'

Here is the caller graph for this function:

◆ _MMG5_BezierNom()

int _MMG5_BezierNom ( MMG5_pMesh  mesh,
int  ip0,
int  ip1,
double  s,
double *  o,
double *  no,
double *  to 
)
inline

Compute point located at parameter value step from point ip0, as well as interpolate of normals, tangent for a NOM edge

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_BezierRef()

int _MMG5_BezierRef ( MMG5_pMesh  mesh,
int  ip0,
int  ip1,
double  s,
double *  o,
double *  no,
double *  to 
)
inline

Compute point located at parameter value step from point ip0, as well as interpolate of normals, tangent for a REF edge

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_BezierReg()

int _MMG5_BezierReg ( MMG5_pMesh  mesh,
int  ip0,
int  ip1,
double  s,
double  v[3],
double *  o,
double *  no 
)
inline

Compute point located at parameter value step from point ip0, as well as interpolate of normals, tangent for a regular edge ; v = ref vector (normal) for choice of normals if need be

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_BezierRidge()

int _MMG5_BezierRidge ( MMG5_pMesh  mesh,
int  ip0,
int  ip1,
double  s,
double *  o,
double *  no1,
double *  no2,
double *  to 
)
inline

Compute point located at parameter value step from point ip0, as well as interpolate of normals, tangent for a RIDGE edge

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_BezierTgt()

int _MMG5_BezierTgt ( double  c1[3],
double  c2[3],
double  n1[3],
double  n2[3],
double  t1[3],
double  t2[3] 
)
inline
Parameters
c1coordinates of the first point of the curve.
c2coordinates of the second point of the curve.
n1normal at the first point of the curve.
n2normal at the second point of the curve.
t1computed normal at the first point of the curve.
t2computed normal at the second point of the curve.
Returns
0 if failed, 1 otherwise.

Compute tangent to geometric support curve passing through c1,c2, with normals n1,n2

Here is the caller graph for this function:

◆ _MMG5_boulenm()

int _MMG5_boulenm ( MMG5_pMesh  mesh,
int  start,
int  ip,
int  iface,
double  n[3],
double  t[3] 
)
Parameters
meshpointer toward the mesh structure.
starttetra index.
ippoint index.
ifaceface index.
ncomputed normal vector.
tcomputed tangent vector.
Returns
0 if point is singular, 1 otherwise.

Define normal and tangent vectors at a non manifold point (ip in start, supported by face iface), enumerating its (outer)surfacic ball.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_boulernm()

int _MMG5_boulernm ( MMG5_pMesh  mesh,
int  start,
int  ip,
int *  ng,
int *  nr 
)
Parameters
meshpointer toward the mesh structure.
startindex of the starting tetrahedra.
iplocal index of the point in the tetrahedra start.
ngpointer toward the number of ridges.
nrpointer toward the number of reference edges.
Returns
ns the number of special edges passing through ip, -1 if fail.

Count the numer of ridges and reference edges incident to the vertex ip when ip is non-manifold.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_boulesurfvolp()

int _MMG5_boulesurfvolp ( MMG5_pMesh  mesh,
int  start,
int  ip,
int  iface,
int *  listv,
int *  ilistv,
int *  lists,
int *  ilists,
int  isnm 
)
Parameters
meshpointer toward the mesh structure.
startindex of the starting tetra.
ipindex in start of the looked point.
ifaceindex in start of the starting face.
listvpointer toward the computed volumic ball.
ilistvpointer toward the computed volumic ball size.
listspointer toward the computed surfacic ball.
ilistspointer toward the computed surfacic ball size.
isnmis the looked point ip non-manifold?
Returns
-1 if fail, 1 otherwise.

Compute the volumic ball of a SURFACE point p, as well as its surfacic ball, starting from tetra start, with point ip, and face if in tetra volumic ball. listv[k] = 4*number of tet + index of point surfacic ball. lists[k] = 4*number of tet + index of face.

Warning
Don't work for a non-manifold point if start has an adjacent through iface (for example : a non-manifold subdomain). Thus, if ip is non-manifold, must be called only if start has no adjacent through iface.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_bouletrid()

int _MMG5_bouletrid ( MMG5_pMesh  mesh,
int  start,
int  iface,
int  ip,
int *  il1,
int *  l1,
int *  il2,
int *  l2,
int *  ip0,
int *  ip1 
)
Parameters
meshpointer toward the mesh structure.
startindex of the starting tetrahedron.
ipindex of the looked ridge point.
ifaceindex in start of the starting face.
il1pointer toward the first ball size.
l1pointer toward the first computed ball (associated to n_1's side).
il2pointer toward the second ball size.
l2pointer toward the second computed ball (associated to n_2's side).
ip0index of the first extremity of the ridge.
ip1index of the second extremity of the ridge.
Returns
0 if fail, 1 otherwise.

Computation of the two surface balls of a ridge point: the list l1 is associated to the normal of face iface. ip0 and ip1 are the indices of the 2 ending point of the ridge. Both lists are returned enumerated in direct order.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_boulevolp()

int _MMG5_boulevolp ( MMG5_pMesh  mesh,
int  start,
int  ip,
int *  list 
)
Parameters
meshpointer toward the mesh structure.
startindex of the starting tetrahedra.
iplocal index of the point in the tetrahedra start.
listpointer toward the list of the tetra in the volumic ball of ip.
Returns
0 if fail and the number of the tetra in the ball otherwise.

Fill the volumic ball (i.e. filled with tetrahedra) of point ip in tetra start. Results are stored under the form $4*kel + jel$, kel = number of the tetra, jel = local index of p within kel.

Here is the caller graph for this function:

◆ _MMG5_caltet33_ani()

double _MMG5_caltet33_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTetra  pt 
)
inline
Parameters
meshpointer toward the mesh structure.
metpointer toward the meric structure.
ptpointer toward a tetrahedra.
Returns
The anisotropic quality of the tet or 0.0 if fail.

Compute the quality of the tet pt with respect to the anisotropic metric met. $ Q = V_met(K) / (sum(len(edge_K)^2)^(3/2) $ and for a calssic storage of metrics at ridges.

Todo:
test with the square of this measure
Here is the caller graph for this function:

◆ _MMG5_cavity_ani()

int _MMG5_cavity_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  iel,
int  ip,
int *  list,
int  lon,
double  volmin 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the sol structure.
ieltetra index.
ippoint local index in iel.
listpointer toward the list of tetra in the shell of edge where ip will be inserted.
lonnumber of tetra in the list.
Returns
ilist number of tetra inside the cavity or -ilist if one of the tet of the cavity is required.

Mark elements in cavity and update the list of tetra in the cavity.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_cavity_iso()

int _MMG5_cavity_iso ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
int  iel,
int  ip,
int *  list,
int  lon,
double  volmin 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the sol structure.
ieltetra index.
ippoint local index in iel.
listpointer toward the list of tetra in the shell of edge where ip will be inserted.
lonnumber of tetra in the list.
Returns
ilist number of tetra inside the cavity or -ilist if one of the tet of the cavity is required.

Mark elements in cavity and update the list of tetra in the cavity.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_cenrad_ani()

int _MMG5_cenrad_ani ( MMG5_pMesh  mesh,
double *  ct,
double *  m,
double *  c,
double *  rad 
)
Parameters
meshpointer toward the mesh structure.
ctcoordinates of vertices of the element.
mmetric at the point for which we compute the cavity.
ccenter of circumscribing circle to the element.
radsquared radius of circumscribing circle to the element.
Returns
0 if failed, 1 otherwise.

Compute radius (squared) and center of circumscribing circle to the element for an anisotropic metric m.

Here is the caller graph for this function:

◆ _MMG5_cenrad_iso()

int _MMG5_cenrad_iso ( MMG5_pMesh  mesh,
double *  ct,
double *  c,
double *  rad 
)
Parameters
meshpointer toward the mesh structure.
ctcoordinates of vertices of the element.
ccenter of circumscribing circle to the element.
radsquared radius of circumscribing circle to the element.
Returns
0 if failed, 1 otherwise.

Compute radius (squared) and center of circumscribing circle to the element.

Here is the caller graph for this function:

◆ _MMG5_chkBdryTria()

int _MMG5_chkBdryTria ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
1 if success, 0 otherwise.

Check the matching between actual and given number of faces in the mesh: Count the number of faces in mesh and compare this number to the number of given triangles. If the founded number exceed the given one, add the missing boundary triangles. Do nothing otherwise.

Step 1: scan the mesh and count the boundaries

Step 2: detect the extra boundaries (that will be ignored) provided by the user
Step 3: add the missing boundary triangles or, if the mesh contains prisms, set to required the triangles at interface betwen prisms and tet

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkcol_bdy()

int _MMG5_chkcol_bdy ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  iface,
char  iedg,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
char  typchk 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element in which we collapse.
ifaceface through wich we perform the collapse
iedgedge to collapse (in local face num)
listvpointer toward the list of the tetra in the ball of p0.
ilistvnumber of tetra in the ball of p0.
listspointer toward the surfacic ball of p0.
ilistsnumber of tetra in the surfacic ball of p0.
typchktypchk type of checking permformed for edge length (hmax or _MMG3D_LLONG criterion).
Returns
1 if success, 0 if the point cannot be collapsed, -1 if fail.

Check whether collapse ip -> iq could be performed, ip boundary point ; 'mechanical' tests (positive jacobian) are not performed here ; iface = boundary face on which lie edge iedg - in local face num. (pq, or ia in local tet notation).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkcol_int()

int _MMG5_chkcol_int ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  iface,
char  iedg,
int *  list,
int  ilist,
char  typchk 
)

Check whether collapse ip -> iq could be performed, ip internal ; 'mechanical' tests (positive jacobian) are not performed here

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkedg() [1/2]

char _MMG5_chkedg ( MMG5_pMesh  mesh,
MMG5_pTria  pt,
char  ori,
double  ,
double  ,
int   
)
Here is the caller graph for this function:

◆ _MMG5_chkedg() [2/2]

char _MMG5_chkedg ( MMG5_pMesh  mesh,
MMG5_Tria pt,
char  ori,
double  hmax,
double  hausd,
int  locPar 
)
Parameters
meshpointer toward the mesh structure.
ptpointer toward the triangle.
oriorientation of the triangle (1 for direct orientation, 0 otherwise).
hmaxmaximal edge length.
hausdmaximal hausdorff distance.
locPar1 if hmax and hausd are locals parameters.
Returns
0 if error.
edges of the triangle pt that need to be split.

Find edges of (virtual) triangle pt that need to be split with respect to the Hausdorff criterion.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkfemtopo()

int _MMG5_chkfemtopo ( MMG5_pMesh  mesh)

Count the number of tetras that have several boundary faces, as well as the number of internal edges connecting points of the boundary

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkmani()

int _MMG5_chkmani ( MMG5_pMesh  mesh)

Check whether implicit surface enclosed in volume is orientable

First test : check whether a tetra has 4 boundary faces
Second test : Check whether configuration is manifold in each ball

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkmanicoll()

int _MMG5_chkmanicoll ( MMG5_pMesh  mesh,
int  k,
int  iface,
int  iedg,
int  ndepmin,
int  ndepplus,
char  isminp,
char  isplp 
)
Parameters
meshpointer toward the mesh structure.
kindex of element in which we collapse.
ifaceface through wich we perform the collapse
iedgedge to collapse
ndepminindex of an elt with ref MG_MINUS and outside the shell of edge.
ndepplusndex of an elt with ref MG_PLUS and outside the shell of edge.
isminp1 if we have found a tetra with ref MG_MINUS
isplp1 if we have found a tetra with ref MG_PLUS
Returns
0 if we create a non manifold situation, 1 otherwise

Check whether collapse of point np to nq does not create a non manifold situation at nq ndepmin, ndepplus = tetra of ref minus, plus in ball of np, not in shell of (np,nq).

First step : pile up tetras of future ball of nq, crossing through the shell of (np,nq), as long as they have same ref as ndepmin list[l] <= 0 if element of ball of np, >= 0, if element of ball of nq
Second step : same process, starting with a tetra of different reference, in the ball of np

Here is the caller graph for this function:

◆ _MMG5_chkptonbdy()

int _MMG5_chkptonbdy ( MMG5_pMesh  mesh,
int  np 
)

Search boundary faces containing point np.

Returns
0 if fail, 1 otherwise
Warning
Not used.
Here is the call graph for this function:

◆ _MMG5_chkswpbdy()

int _MMG5_chkswpbdy ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  list,
int  ilist,
int  it1,
int  it2,
char  typchk 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
listpointer toward the shell of the edge.
ilistpointer toward the size of the shell of the edge.
it1first element of the open shell.
it2last element of the open shell.
typchktype of checking permformed for edge length (hmin or LSHORT criterion).
Returns
-1 if fail, 0 if we can not swap the edge, 1 otherwise.

Check whether edge whose shell is provided should be swapped for geometric approximation purposes (the 2 surface triangles are also provided).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_chkswpgen()

int _MMG5_chkswpgen ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  start,
int  ia,
int *  ilist,
int *  list,
double  crit,
char  typchk 
)
Parameters
meshpointer toward the mesh structure
metpointer toward the metric structure.
starttetrahedra in which the swap should be performed
iaedge that we want to swap
ilistpointer to store the size of the shell of the edge
listpointer to store the shell of the edge
critimprovment coefficient
Returns
0 if fail, the index of point corresponding to the swapped configuration otherwise ( $4*k+i$).
Parameters
typchktype of checking permformed for edge length (hmin or LSHORT criterion).

Check whether swap of edge ia in start should be performed, and return $4*k+i$ the index of point corresponding to the swapped configuration. The shell of edge is built during the process.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_cntbdypt()

int _MMG5_cntbdypt ( MMG5_pMesh  mesh,
int  nump 
)

Count how many boundary faces share point nump.

Warning
Not used.
Here is the call graph for this function:

◆ _MMG5_colver()

int _MMG5_colver ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  list,
int  ilist,
char  indq,
char  typchk 
)
Parameters
meshpointer toward the mesh
metpointer toward the metric
listpointer toward the ball of the point
ilistnumber of elements in the ball of the point
indqlocal index of the point on which we collapse
typchktype of check performed depending on the remeshing step
Returns
np the index of the collpased point if success, 0 if we cannot collapse, -1 if we fail.

Collapse vertex p = list[0]%4 of tetra list[0]/4 over vertex indq of tetra list[0]/4. Only physical tests (positive jacobian) are done (i.e. approximation of the surface, etc... must be performed outside).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_coquil()

int _MMG5_coquil ( MMG5_pMesh  mesh,
int  start,
int  ia,
int *  list 
)
Parameters
meshpointer toward the mesh structure
startindex of the starting tetra
iaindex of the edge
listlist of tetra sharing the edge ia
Returns
2*ilist if shell is closed, 2*ilist +1 otherwise, 0 if one of the tet of the shell is required

Find all tets sharing edge ia of tetra start.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_coquilface()

int _MMG5_coquilface ( MMG5_pMesh  mesh,
int  start,
char  iface,
int  ia,
int *  list,
int *  it1,
int *  it2,
int  silent 
)
Parameters
meshpointer toward the mesh structure.
startindex of the starting tetrahedron.
ifaceindex of the boundary face from which we come.
iaindex of edge whose shell is computed (in tetra).
listpointer toward the list of tetra in the shell (to fill).
it1pointer toward the index of the first boundary face sharing ia (to fill).
it2pointer toward the index of the second boundary face sharing ia (to fill).
silentif 1, print error message for more than 2 boundary triangles in the shell
Returns
-1 if fail, $2*ilist$ if shell is closed, $2*ilist+1$ otherwise.

Find all tets sharing edge ia of tetra start, and stores boundary faces when met. $ it1 $ and $ it2 = 6*iel + iface$, iel = index of tetra, iface = index of face in tetra.

Warning
Don't work if ia has only one boundary face in its shell.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_coquilFaceErrorMessage()

void _MMG5_coquilFaceErrorMessage ( MMG5_pMesh  mesh,
int  k1,
int  k2 
)
Parameters
meshpointer toward the mesh structure.
k1should contain a tetra index.
k2should contain a tetra index different from k2.

Print an error message if _MMG5_coquilFace detect a boundary topology problem.

Here is the caller graph for this function:

◆ _MMG5_coquilTravel()

int16_t _MMG5_coquilTravel ( MMG5_pMesh  mesh,
int  na,
int  nb,
int *  adj,
int *  piv,
char *  iface,
int *  i 
)
Parameters
meshpointer toward the mesh structure.
naglobal index of edge extremity.
nbglobal index of edge extremity.
adjstarting tetrahedron at the begining and finish tet at the end.
pivglobal index of the vertex opposite to the travelling face (updated for the finish tet at the end).
ifaceprevious traveling face of the tet (suspected to be boundary), updated.
ilocal index of the edge $[na,nb]$ in tet adj.
Returns
the tag of the face iface of the tetra adj or 0 if the tetra is not boundary.

Travel around the edge $[na,nb]$ from tetra adj and through the face piv.

Here is the caller graph for this function:

◆ _MMG5_countelt()

int _MMG5_countelt ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
double *  weightelt,
long *  npcible 
)

Approximation of the final number of vertex.

Warning
call MMG3D_hashTetra(mesh,1) or analysis before using
Todo:
Doxygen documentation
Here is the call graph for this function:

◆ _MMG5_defaultValues()

void _MMG5_defaultValues ( MMG5_pMesh  )

◆ _MMG5_delone()

int _MMG5_delone ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
int  ip,
int *  list,
int  ilist 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the solution structure.
ipindex of the point to insert.
listpointer toward the list of the tetra in the cavity (computed by _MMG5_cavity).
ilistnumber of tetra inside the cavity.
Returns
1 if sucess, 0 or -1 if fail.

Insertion of the vertex ip. The cavity of ip become its ball.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_deltag()

int _MMG5_deltag ( MMG5_pMesh  mesh,
int  start,
int  ia,
int16_t  tag 
)
inline
Parameters
meshpointer toward the mesh structure
startindex of the starting tetra
iaindex of the edge in tetra start that we want to modify
tagtag to remove
Returns
1 if success.

Remove the tag tag of edge ia in tetra start by travelling its shell.

Here is the caller graph for this function:

◆ _MMG5_denoisbdy()

int _MMG5_denoisbdy ( MMG5_pMesh  )

◆ _MMG5_dichodisp()

int _MMG5_dichodisp ( MMG5_pMesh  ,
double *   
)

◆ _MMG5_directsurfball()

int _MMG5_directsurfball ( MMG5_pMesh  mesh,
int  ip,
int *  list,
int  ilist,
double  n[3] 
)
inline

If need be, invert the travelling sense of surfacic ball so that it is travelled in the direct sense with respect to direction n anchored at point ip (ip = global num.): return 2 = orientation reversed, 1 otherwise

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_estavglen()

double _MMG5_estavglen ( MMG5_pMesh  )

◆ _MMG5_freeXPrisms()

void _MMG5_freeXPrisms ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.

Free xprism structure.

Here is the caller graph for this function:

◆ _MMG5_freeXTets()

void _MMG5_freeXTets ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.

Free xtetra structure.

Here is the caller graph for this function:

◆ _MMG5_gradsiz_ani()

int _MMG5_gradsiz_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
1

Enforces mesh gradation by truncating metric field.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_gradsiz_iso()

int _MMG5_gradsiz_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
0 if fail, 1 otherwise.

Enforce mesh gradation by truncating size map.

Here is the caller graph for this function:

◆ _MMG5_hashPop()

int _MMG5_hashPop ( _MMG5_Hash hash,
int  a,
int  b 
)

remove edge from hash table

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_hEdge()

int _MMG5_hEdge ( MMG5_pMesh  mesh,
MMG5_HGeom hash,
int  a,
int  b,
int  ref,
int16_t  tag 
)

store edge on geometry

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_hGeom()

int _MMG5_hGeom ( MMG5_pMesh  mesh)
Parameters
meshpointer toward he mesh structure.
Returns
0 if failed, 1 otherwise

Build hashtable for initial mesh edges.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_hGet()

int _MMG5_hGet ( MMG5_HGeom hash,
int  a,
int  b,
int *  ref,
int16_t *  tag 
)

get ref and tag to edge on geometry

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_hNew()

int _MMG5_hNew ( MMG5_pMesh  mesh,
MMG5_HGeom hash,
int  hsiz,
int  hmax 
)

to store edge on geometry

Here is the caller graph for this function:

◆ _MMG5_hPop()

int _MMG5_hPop ( MMG5_HGeom hash,
int  a,
int  b,
int *  ref,
int16_t *  tag 
)

remove edge from hash table

Here is the call graph for this function:

◆ _MMG5_hTag()

int _MMG5_hTag ( MMG5_HGeom hash,
int  a,
int  b,
int  ref,
int16_t  tag 
)

set tag to edge on geometry

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_Init_parameters()

void _MMG5_Init_parameters ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.

Initialization of the input parameters (stored in the Info structure).

MMG5_IPARAM_verbose = 1















MMG*_IPARAM_iso = 0














MMG5_IPARAM_mem = -1













MMG5_IPARAM_debug = 0












MMG5_IPARAM_npar = 0











MMG5_IPARAM_noinsert = 0










MMG5_IPARAM_noswap = 0









MMG5_IPARAM_nomove = 0








MMG5_IPARAM nmat = 0







MMG5_DPARAM_angleDetection = _MMG5_ANGEDG






MMG5_DPARAM_hmin = 0.001 $\times$ bounding box size;





MMG5_DPARAM_hmax = double of the bounding box size




MMG5_DPARAM_hsiz= -1.



MMG5_DPARAM_hausd = 0.01


MMG5_DPARAM_hgrad = 1.3

MMG5_PPARAM = NULL
MMG3D_IPARAM_lag = -1 used by mmg3d only but need to be negative in the scaleMesh function

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_interp4bar33_ani()

int _MMG5_interp4bar33_ani ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
int  ,
double *   
)

◆ _MMG5_interp4bar_ani()

int _MMG5_interp4bar_ani ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
int  ,
double *   
)

◆ _MMG5_interp4bar_iso()

int _MMG5_interp4bar_iso ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
int  ,
double *   
)

◆ _MMG5_intmet_ani()

int _MMG5_intmet_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  i,
int  ip,
double  s 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kelement index.
ilocal index of edge in k.
ipglobal index of the new point in which we want to compute the metric.
sinterpolation parameter (between 0 and 1).
Returns
0 if fail, 1 otherwise.

Interpolation of anisotropic sizemap at parameter s along edge i of elt k for a special storage of ridges metric (after defsiz call).

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_intmet_iso()

int _MMG5_intmet_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  i,
int  ip,
double  s 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kelement index.
ilocal index of edge in k.
ipglobal index of the new point in which we want to compute the metric.
sinterpolation parameter (between 0 and 1).
Returns
0 if fail, 1 otherwise.

Interpolation of anisotropic sizemap at parameter s along edge i of elt k.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_intregmet()

int _MMG5_intregmet ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
char  ,
double  ,
double *   
)

◆ _MMG5_intridmet()

int _MMG5_intridmet ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
int  ,
double  ,
double *  ,
double *   
)

◆ _MMG5_intvolmet()

int _MMG5_intvolmet ( MMG5_pMesh  ,
MMG5_pSol  ,
int  ,
char  ,
double  ,
double *   
)

◆ _MMG5_lapantilap()

int _MMG5_lapantilap ( MMG5_pMesh  ,
double *   
)

◆ _MMG5_lenedgCoor_iso()

double _MMG5_lenedgCoor_iso ( double *  ca,
double *  cb,
double *  ma,
double *  mb 
)
inline

Compute edge length from edge's coordinates.

Parameters
*capointer toward the coordinates of the first edge's extremity.
*cbpointer toward the coordinates of the second edge's extremity.
*mapointer toward the metric associated to the first edge's extremity.
*mbpointer toward the metric associated to the second edge's extremity.
Returns
edge length.

Compute length of edge $[ca,cb]$ (with ca and cb coordinates of edge extremities) according to the isotropic size prescription.

Here is the caller graph for this function:

◆ _MMG5_meancur()

int _MMG5_meancur ( MMG5_pMesh  mesh,
int  np,
double  c[3],
int  ilist,
int *  list,
double  h[3] 
)

◆ _MMG5_meansizreg_iso()

double _MMG5_meansizreg_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  nump,
int *  lists,
int  ilists,
double  hmin,
double  hmax 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
numpindex of point in which the size must be computed.
listspointer toward the surfacic ball of nump.
ilistssize of surfacic ball of nump.
hminminimal edge size.
hmaxmaximal edge size.
Returns
the isotropic size at the point.

For -nosurf option : define isotropic size at regular point nump, whose surfacic ball is provided. The size is computed as the mean of the length of the surface edges passing through nump.

Here is the caller graph for this function:

◆ _MMG5_memSize()

long long _MMG5_memSize ( void  )
Returns
the available memory size of the computer.

Compute the available memory size of the computer.

◆ _MMG5_mmg3d1_delone()

int _MMG5_mmg3d1_delone ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
0 if failed, 1 if success.

Main adaptation routine.

— stage 1: geometric mesh
— stage 2: computational mesh

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_mmg3d1_pattern()

int _MMG5_mmg3d1_pattern ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
Returns
0 if failed, 1 if success.

Main adaptation routine.

— stage 1: geometric mesh
— Stage 2: computational mesh

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_mmg3d3()

int _MMG5_mmg3d3 ( MMG5_pMesh  ,
MMG5_pSol  ,
MMG5_pSol   
)
Here is the caller graph for this function:

◆ _MMG5_mmg3dBezierCP()

int _MMG5_mmg3dBezierCP ( MMG5_pMesh  mesh,
MMG5_Tria pt,
_MMG5_pBezier  pb,
char  ori 
)
Parameters
meshpointer toward the mesh structure.
ptpointer toward the triangle structure.
pbpointer toward the computed Bezier structure.
oritriangle orientation.
Returns
1.

Compute Bezier control points on triangle pt (cf. Vlachos)

Todo:
merge with the _MMG5_mmg3dBeizerCP function and remove the pointer toward this functions.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_mmg3dChkmsh()

int _MMG5_mmg3dChkmsh ( MMG5_pMesh  mesh,
int  severe,
int  base 
)
Parameters
meshpointer toward the mesh structure.
severelevel of performed check
baseunused argument.
Returns
0 if fail, 1 if success.

Check the mesh validity

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdynompt_ani()

int _MMG5_movbdynompt_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if fail, 1 if success.

Move boundary non manifold point, whose volumic and (exterior) surfacic balls are passed

Remarks
we don't check if we break the hausdorff criterion.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdynompt_iso()

int _MMG5_movbdynompt_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if fail, 1 if success.

Move boundary non manifold point, whose volumic and (exterior) surfacic balls are passed

Remarks
the metric is not interpolated at the new position.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdyrefpt_ani()

int _MMG5_movbdyrefpt_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if fail, 1 if success.
Remarks
we don't check if we break the hausdorff criterion.

Move boundary reference point, whose volumic and surfacic balls are passed.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdyrefpt_iso()

int _MMG5_movbdyrefpt_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if fail, 1 if success.

Move boundary reference point, whose volumic and surfacic balls are passed.

Remarks
the metric is not interpolated at the new position.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdyregpt_ani()

int _MMG5_movbdyregpt_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improveSurf,
int  improveVol 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if we can't move the point, 1 if we can.
Remarks
we don't check if we break the hausdorff criterion.
the metric is not interpolated at the new position.

Move boundary regular point, whose volumic and surfacic balls are passed.

Step 1 : rotation matrix that sends normal n to the third coordinate vector of R^3



Step 2 : rotation of the oriented surfacic ball with r : lispoi[k] is the common edge between faces lists[k-1] and lists[k]


Step 3 : Compute gradient towards optimal position = centre of mass of the ball, projected to tangent plane

Step 4 : locate new point in the ball, and compute its barycentric coordinates
Step 5 : come back to original problem, and compute patch in triangle iel

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdyregpt_iso()

int _MMG5_movbdyregpt_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improveSurf,
int  improveVol 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if we can not move, 1 if success, -1 if fail.

Move boundary regular point, whose volumic and surfacic balls are passed.

Remarks
the metric is not interpolated at the new position.

Step 1 : rotation matrix that sends normal n to the third coordinate vector of R^3



Step 2 : rotation of the oriented surfacic ball with r : lispoi[k] is the common edge between faces lists[k-1] and lists[k]


Step 3 : Compute optimal position to make current triangle equilateral, and average of these positions

Step 4 : locate new point in the ball, and compute its barycentric coordinates
Step 5 : come back to original problem, and compute patch in triangle iel

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdyridpt_ani()

int _MMG5_movbdyridpt_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal
Returns
0 if fail, 1 if success.
Remarks
we don't check if we break the hausdorff criterion.

Move boundary ridge point, whose volumic and surfacic balls are passed.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movbdyridpt_iso()

int _MMG5_movbdyridpt_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  listv,
int  ilistv,
int *  lists,
int  ilists,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listvpointer toward the volumic ball of the point.
ilistvsize of the volumic ball.
listspointer toward the surfacic ball of the point.
ilistssize of the surfacic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if fail, 1 if success.

Move boundary ridge point, whose volumic and surfacic balls are passed.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movintpt_ani()

int _MMG5_movintpt_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  list,
int  ilist,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listpointer toward the volumic ball of the point.
ilistsize of the volumic ball.
improveforce the new minimum element quality to be greater or equal than 1.02 of the old minimum element quality.
Returns
0 if we can't move the point, 1 if we can.

Move internal point whose volumic is passed.

Remarks
the metric is not interpolated at the new position.
we don't check if we break the hausdorff criterion.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movintpt_iso()

int _MMG5_movintpt_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  list,
int  ilist,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listpointer toward the volumic ball of the point.
ilistsize of the volumic ball.
improveforce the new minimum element quality to be greater or equal than 0.9 of the old minimum element quality.
Returns
0 if we can't move the point, 1 if we can.

Move internal point whose volumic is passed.

Remarks
the metric is not interpolated at the new position.
we don't check if we break the hausdorff criterion.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_movintptLES_iso()

int _MMG5_movintptLES_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int *  list,
int  ilist,
int  improve 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
listpointer toward the volumic ball of the point.
ilistsize of the volumic ball.
improveforce the new minimum element quality to be greater or equal than 0.9 of the old minimum element quality.
Returns
0 if we can't move the point, 1 if we can.

Move internal point whose volumic ball is passed (for LES optimization). The optimal point position is computed as the barycenter of the optimal point position for each tetra. The optimal point position for a tetra is the point located over the normal of the face at the face barycenter and at the distance 1 of the face.

Remarks
the metric is not interpolated at the new position.
we don't check if we break the hausdorff criterion.
Here is the call graph for this function:

◆ _MMG5_movtet()

int _MMG5_movtet ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
double  clickSurf,
double  clickVol,
int  moveVol,
int  improveSurf,
int  improveVolSurf,
int  improveVol,
int  maxit 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
clickSurftriangle quality threshold under which we want to move
clickVoltetra quality threshold under which we want to move
moveVolinternal move
improveSurfforbid surface degradation during the move
improveVolSurfforbid volume degradation during the surfacic move
improveVolforbid volume degradation during the move
maxitmaximum number of iteration.
Returns
-1 if failed, number of moved points otherwise.

Analyze tetrahedra and move points so as to make mesh more uniform.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_moymet()

int _MMG5_moymet ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTetra  pt,
double *  m1 
)
inline
Parameters
meshpointer toward the mesh structure.
metpointer toward the sol structure.
ptpointer toward a tetra.
m1computed metric.
Returns
the number of vertices used for the mean computation, 0 if fail.

Compute mean metric over the internal tetra pt. Do not take into account the metric values at ridges points (because we don't know how to build it).

Here is the caller graph for this function:

◆ _MMG5_norface()

int _MMG5_norface ( MMG5_pMesh  mesh,
int  k,
int  iface,
double  n[3] 
)
inline

Compute normal to face iface of tetra k, exterior to tetra k

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_nsort()

void _MMG5_nsort ( int  n,
double *  val,
char *  perm 
)
inline

naive (increasing) sorting algorithm, for very small tabs ; permutation is stored in perm

◆ _MMG5_openCoquilTravel()

void _MMG5_openCoquilTravel ( MMG5_pMesh  mesh,
int  na,
int  nb,
int *  adj,
int *  piv,
char *  iface,
int *  i 
)
Parameters
meshpointer toward the mesh structure.
naglobal index of edge extremity.
nbglobal index of edge extremity.
adjstarting tetrahedron at the begining and finish tet at the end.
pivglobal index of the vertex opposite to the travelling face (updated for the finish tet at the end).
ifacetraveling face of the tet (suspected to be boundary), updated.
ilocal index of the edge $[na,nb]$ in tet adj.

Travel around the edge $[na,nb]$ from tetra adj and through the face piv. The shell of the edge is open and the tetra adj has no neighbour through the face iface.

Here is the caller graph for this function:

◆ _MMG5_orcal_poi()

double _MMG5_orcal_poi ( double  a[3],
double  b[3],
double  c[3],
double  d[3] 
)

◆ _MMG5_ppgdisp()

int _MMG5_ppgdisp ( MMG5_pMesh  ,
double *   
)

◆ _MMG5_printTetra()

void _MMG5_printTetra ( MMG5_pMesh  mesh,
char *  fileName 
)

Debug function (not use in clean code): print mesh->tetra structure

◆ _MMG5_saveDisp()

int _MMG5_saveDisp ( MMG5_pMesh  ,
MMG5_pSol   
)

◆ _MMG5_setNmTag()

int _MMG5_setNmTag ( MMG5_pMesh  mesh,
_MMG5_Hash hash 
)
Parameters
meshpointer towar the mesh structure.
hashedges hash table.
Returns
1 if success, 0 if failed.

Set tags to non-manifold edges and vertices. Not done before because we need the MMG5_xTetra table.

Warning
if fail, the edge hash table hash is not freed.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_settag()

int _MMG5_settag ( MMG5_pMesh  mesh,
int  start,
int  ia,
int16_t  tag,
int  edg 
)
inline

Set tag and edg of edge ia (if need be) in tetra start by travelling its shell

Here is the caller graph for this function:

◆ _MMG5_split1()

int _MMG5_split1 ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split 1 edge of tetra k.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split1b()

int _MMG5_split1b ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  list,
int  ret,
int  ip,
int  cas,
char  metRidTyp,
char  chkRidTet 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
listpointer toward the shell of edge.
retsize of the shell of edge.
ipidex of new point.
casflag to watch the length of the new edges.
metRidTypType of storage of ridges metrics: 0 for classic storage, 1 for special storage.
chkRidTetif 1, avoid the creation of a tet with 4 ridge vertices
Returns
-1 if we fail, 0 if we don't split the edge, 1 if success.

Split edge $list[0]\%6$, whose shell list is passed, introducing point ip Beware : shell has to be enumerated in ONLY ONE TRAVEL (always same sense).

2 different checks : 1) are we creating a too small edge (BUG_Split1b_SpereIso_0.125h_met) 2) in aniso and from the last wave of anatet(typchk=1): avoid the creation of a tetra with 4 ridge vertices.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split2()

int _MMG5_split2 ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split of two OPPOSITE edges

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split2sf()

int _MMG5_split2sf ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split of two edges that belong to a common face : 1 tetra becomes 3

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split3()

int _MMG5_split3 ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

1 face (3 edges) subdivided

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split3cone()

int _MMG5_split3cone ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split 3 edge in cone configuration

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split3op()

int _MMG5_split3op ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split 3 opposite edges in a tetra

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split4bar()

int _MMG5_split4bar ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ktetra index.
metRidTypmetric storage (classic or special)
Returns
0 if fail, index of created point otherwise (ib)

Split a tetra in 4 tetras by introducing its barycenter. FOR NOW : flags, that tell which edge should be split, are not updated (erased) : UPDATE NEEDED ?

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split4op()

int _MMG5_split4op ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split 4 edges in a configuration when no 3 edges lie on the same face

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split4sf()

int _MMG5_split4sf ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split 4 edges in a configuration when 3 lie on the same face

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split5()

int _MMG5_split5 ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

Split 5 edges

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_split6()

int _MMG5_split6 ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
int  vx[6],
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of element to split.
vx$vx[i]$ is the index of the point to add on the edge i.
metRidTypmetric storage (classic or special)
Returns
0 if fail, 1 otherwise

split all faces (6 edges)

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_splitedg()

int _MMG5_splitedg ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  iel,
int  iar,
double  crit 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ieltetra index
iaredge index of iel
critquality threshold.
Returns
-1 if lack of memory, 0 if we don't split the edge, ip if success.

Split edge iar of iel and verify that every new tet have a better quality than crit

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_srcbdy()

int _MMG5_srcbdy ( MMG5_pMesh  mesh,
int  start,
int  ia 
)

Identify whether edge ia in start is a boundary edge by unfolding its shell

Here is the caller graph for this function:

◆ _MMG5_srcface()

int _MMG5_srcface ( MMG5_pMesh  mesh,
int  n0,
int  n1,
int  n2 
)

◆ _MMG5_startedgsurfball()

int _MMG5_startedgsurfball ( MMG5_pMesh  mesh,
int  nump,
int  numq,
int *  list,
int  ilist 
)

If need be, reorder the surfacic ball of point ip, so that its first element has edge (p,q) (nump,q = global num) as edge _MMG5_iprv2[ip] of face iface. return 2 = orientation reversed, 1 otherwise

Here is the caller graph for this function:

◆ _MMG5_stiffelt()

int _MMG5_stiffelt ( MMG5_pMesh  ,
int  ,
double *  ,
double *   
)

◆ _MMG5_surftri()

double _MMG5_surftri ( MMG5_pMesh  ,
int  ,
int   
)

◆ _MMG5_swpbdy()

int _MMG5_swpbdy ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  list,
int  ret,
int  it1,
_MMG3D_pOctree  octree,
char  typchk 
)
Parameters
meshpointer toward the mesh structure
metpointer toward the solution structure
listpointer toward the shell of the edge
retdobble of the number of tetrahedra in the shell
it1boundary face carrying the beforehand tested terminal point for collapse
octreepointer toward the octree structure in Delaunay mode, NULL pointer in pattern mode.
typchktype of checking permformed for edge length (hmin or LSHORT criterion).
Returns
-1 if lack of memory, 0 if fail to swap, 1 otherwise

Swap boundary edge whose shell is provided.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_swpgen()

int _MMG5_swpgen ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  nconf,
int  ilist,
int *  list,
_MMG3D_pOctree  octree,
char  typchk 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the sol structure.
nconfconfiguration.
ilistnumber of tetrahedra in the shell of the edge that we want to swap.
listpointer toward the shell of the edge that we want to swap.
octreepointer toward the octree structure in Delaunay mode, NULL pointer in pattern mode.
typchktype of checking permformed for edge length (hmin or LSHORT criterion).
Returns
-1 if lack of memory, 0 if fail to swap, 1 otherwise.

Perform swap of edge whose shell is passed according to configuration nconf.

First step : split of edge (na,nb)
Second step : collapse of np towards enhancing configuration

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_swpmsh()

int _MMG5_swpmsh ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int  typchk 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure (only for delaunay).
typchktype of checking permformed for edge length (hmin or LSHORT criterion).
Returns
-1 if failed and swap number otherwise.

Search for boundary edges that could be swapped for geometric approximation.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_swptet()

int _MMG5_swptet ( MMG5_pMesh  mesh,
MMG5_pSol  met,
double  crit,
double  declic,
_MMG3D_pOctree  octree,
int  typchk 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
critcoefficient of quality improvment.
octreepointer toward the octree structure in delaunay mode and toward the NULL pointer otherwise
typchktype of checking permformed for edge length (hmin or LSHORT criterion).

Internal edge flipping.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _MMG5_tet2tri()

void _MMG5_tet2tri ( MMG5_pMesh  mesh,
int  k,
char  ie,
MMG5_Tria ptt 
)
Parameters
meshpointer toward the mesh structure.
ktetrahedron index.
ieface index of tetrahedron.
pttpointer toward the output triangle.

Set triangle corresponding to face ie of tetra k.

Here is the caller graph for this function:

◆ _MMG5_timestepMCF()

double _MMG5_timestepMCF ( MMG5_pMesh  ,
double   
)

◆ _MMG5_trydisp()

int _MMG5_trydisp ( MMG5_pMesh  ,
double *  ,
short   
)

◆ _MMG5_velextLS()

int _MMG5_velextLS ( MMG5_pMesh  ,
MMG5_pSol   
)

◆ _MMG5_volint()

double _MMG5_volint ( MMG5_pMesh  )

◆ _MMG5_warnOrientation()

static void _MMG5_warnOrientation ( MMG5_pMesh  mesh)
inlinestatic
Parameters
meshpointer toward the mesh structure.

Warn user that some tetrahedra of the mesh have been reoriented.

Here is the caller graph for this function:

◆ MMG3D_hashPrism()

int MMG3D_hashPrism ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
0 if failed, 1 otherwise.

Create table of adjacency for prisms.

Warning
check the hashtable efficiency
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG3D_movetetrapoints()

int MMG3D_movetetrapoints ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int  k 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
kindex of a tetra
Returns
1 if we move one of the vertices, 0 otherwise.

Try to move the vertices of the tetra k to improve its quality.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG3D_optbdry()

int MMG3D_optbdry ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree,
int  k 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
kindex of a tetra
Returns
1 if success, 0 if fail.

Try to optimize the tetra k. This tetra has a face on the boundary.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG3D_opttyp()

int MMG3D_opttyp ( MMG5_pMesh  mesh,
MMG5_pSol  met,
_MMG3D_pOctree  octree 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
octreepointer toward the octree structure.
Returns
0 if fail, number of improved elts otherwise.

Travel across the mesh to detect element with very bad quality (less than 0.2) and try to improve them by every means.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MMG3D_swap23()

int MMG3D_swap23 ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  k,
char  metRidTyp 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the sol structure.
kindex of the tetrahedron with multiple boundary faces (to be swapped).
metRidTypmetric storage (classic or special)
Returns
-1 if lack of memory, 0 if fail to swap, 1 otherwise.

Search an adjacent to the tetra k and perform swap 2->3 (the common face of the 2 tetra is destroyed and replaced by a common edge used by the three new elts).

Remarks
used in anatet4 to remove the tetra with multiple boundary faces.

Neighbouring element with which we will try to swap

Swap
Quality Update

Here is the call graph for this function:

Variable Documentation

◆ _MMG3D_octreein

int(* _MMG3D_octreein) (MMG5_pMesh,MMG5_pSol,_MMG3D_pOctree,int, double)

◆ _MMG5_arpt

const unsigned char _MMG5_arpt[4][3] = { {0,1,2}, {0,4,3}, {1,3,5}, {2,5,4} }
static

arpt[i]: edges passing through vertex i

◆ _MMG5_caltet

double(* _MMG5_caltet) (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)

◆ _MMG5_caltri

double(* _MMG5_caltri) (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)

◆ _MMG5_cavity

int(* _MMG5_cavity) (MMG5_pMesh,MMG5_pSol,int,int,int *, int,double)

◆ _MMG5_defsiz

int(* _MMG5_defsiz) (MMG5_pMesh,MMG5_pSol)

◆ _MMG5_gradsiz

int(* _MMG5_gradsiz) (MMG5_pMesh,MMG5_pSol)

◆ _MMG5_iare

const unsigned char _MMG5_iare[6][2] = { {0,1}, {0,2}, {0,3}, {1,2}, {1,3}, {2,3} }
static

vertices of extremities of the edges of the tetra

◆ _MMG5_iarf

const unsigned char _MMG5_iarf[4][3] = { {5,4,3}, {5,1,2}, {4,2,0}, {3,0,1} }
static

iarf[i]: edges of face opposite to vertex i

◆ _MMG5_iarf_pr

const unsigned char _MMG5_iarf_pr[5][5] = { {0,1,3,0}, {6,8,7,6}, {3,5,8,4}, {5,1,2,7},{0,4,6,2} }
static

iarf[i]: edges of face i for a prism

◆ _MMG5_iarfinv

const unsigned char _MMG5_iarfinv[4][6] = { {-1,-1,-1,2,1,0}, {-1,1,2,-1,-1,0},{2,-1,1,-1,0,-1},{1,2,-1,0,-1,-1}}
static

num of the j^th edge in the i^th face

◆ _MMG5_idir

const unsigned char _MMG5_idir[4][3] = { {1,2,3}, {0,3,2}, {0,1,3}, {0,2,1} }
static

idir[i]: vertices of face opposite to vertex i

◆ _MMG5_idir_pr

const unsigned char _MMG5_idir_pr[5][4] = { {0,1,2,0},{3,5,4,3},{1,4,5,2},{0,2,5,3},{0,3,4,1} }
static

idir[i]: vertices of face i for a prism

◆ _MMG5_idirinv

const char _MMG5_idirinv[4][4] = {{-1,0,1,2},{0,-1,2,1},{0,1,-1,2},{0,2,1,-1}}
static

◆ _MMG5_ifar

const unsigned char _MMG5_ifar[6][2] = { {2,3}, {1,3}, {1,2}, {0,3}, {0,2}, {0,1} }
static

ifar[i][]: faces sharing the ith edge of the tetra

◆ _MMG5_interp4bar

int(* _MMG5_interp4bar) (MMG5_pMesh, MMG5_pSol, int, int, double *)

◆ _MMG5_intmet

int(* _MMG5_intmet) (MMG5_pMesh, MMG5_pSol, int, char, int, double)

◆ _MMG5_inxt3

const unsigned char _MMG5_inxt3[7] = { 1,2,3,0,1,2,3 }
static

next vertex of tetra: {1,2,3,0,1,2,3}

◆ _MMG5_iprv3

const unsigned char _MMG5_iprv3[7] = { 3,0,1,2,3,0,1 }
static

previous vertex of tetra: {3,0,1,2,3,0,1}

◆ _MMG5_isar

const unsigned char _MMG5_isar[6][2] = { {2,3}, {3,1}, {1,2}, {0,3}, {2,0}, {0,1} }
static

isar[i][]: vertices of extremities of the edge opposite to the ith edge

◆ _MMG5_lenedg

double(* _MMG5_lenedg) (MMG5_pMesh,MMG5_pSol,int, MMG5_pTetra)

◆ _MMG5_lenedgspl

double(* _MMG5_lenedgspl) (MMG5_pMesh,MMG5_pSol,int, MMG5_pTetra)

◆ _MMG5_movbdynompt

int(* _MMG5_movbdynompt) (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree,int *, int, int *, int,int)

◆ _MMG5_movbdyrefpt

int(* _MMG5_movbdyrefpt) (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree,int *, int, int *, int,int)

◆ _MMG5_movbdyregpt

int(* _MMG5_movbdyregpt) (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree,int *, int, int *, int, int,int)

◆ _MMG5_movbdyridpt

int(* _MMG5_movbdyridpt) (MMG5_pMesh, MMG5_pSol, _MMG3D_pOctree,int *, int, int *, int,int)

◆ _MMG5_movintpt

int(* _MMG5_movintpt) (MMG5_pMesh,MMG5_pSol, _MMG3D_pOctree,int *, int, int)

◆ MMG5_permedge

const unsigned char MMG5_permedge[12][6]
static
Initial value:
= {
{0,1,2,3,4,5}, {1,2,0,5,3,4}, {2,0,1,4,5,3}, {0,4,3,2,1,5},
{3,0,4,1,5,2}, {4,3,0,5,2,1}, {1,3,5,0,2,4}, {3,5,1,4,0,2},
{5,1,3,2,4,0}, {2,5,4,1,0,3}, {4,2,5,0,3,1}, {5,4,2,3,1,0} }

Table that associates to each (even) permutation of the 4 vertices of a tetrahedron the corresponding permutation of its edges.
Labels : 0 : [0,1,2,3] 1 : [0,2,3,1] 2 : [0,3,1,2] 3 : [1,0,3,2] 4 : [1,2,0,3] 5 : [1,3,2,0] 6 : [2,0,1,3] 7 : [2,1,3,0] 8 : [2,3,0,1] 9 : [3,0,2,1] 10 : [3,1,0,2] 11 : [3,2,1,0] The edge 0 of the config 1 become the edge 1 of the reference config so permedge[1][0]=1 ...