Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
Data Structures | Macros | Typedefs | Functions | Variables
mmgcommon_private.h File Reference
#include <assert.h>
#include <stdlib.h>
#include <stdio.h>
#include <limits.h>
#include <string.h>
#include <signal.h>
#include <ctype.h>
#include <float.h>
#include <math.h>
#include <complex.h>
#include "mmg/common/mmgcmakedefines.h"
#include "eigenv_private.h"
#include "libmmgcommon_private.h"
Include dependency graph for mmgcommon_private.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  MMG5_Bezier
 
struct  MMG5_iNode_s
 

Macros

#define MG_STR   "&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&"
 
#define MMG5_TRIA_LMAX   1024
 
#define MMG5_LMAX   10240
 
#define MMG5_LPARMAX   200
 
#define MG_SMSGN(a, b)   (((double)(a)*(double)(b) > (0.0)) ? (1) : (0))
 
#define MMG5_BOXSIZE   500
 
#define MMG5_MEMMAX   800
 
#define MMG5_BITWIZE_MB_TO_B   20
 
#define MMG5_MEMPERCENT   0.5
 
#define MMG5_UNSET   -1
 
#define MMG5_DISPREF   10
 
#define MMG5_MILLION   1048576
 
#define MMG5_ANGEDG   0.707106781186548 /*0.573576436351046 */
 
#define MMG5_ANGLIM   -0.999999
 
#define MMG5_ATHIRD   0.333333333333333
 
#define MMG5_EPSD   1.e-30
 
#define MMG5_EPSD2   1.0e-200
 
#define MMG5_EPS   1.e-06
 
#define MMG5_EPSOK   1.e-15
 
#define MMG5_NULKAL   1.e-30
 
#define MMG5_SQR32   0.866025403784439
 
#define A64TH   0.015625
 
#define A16TH   0.0625
 
#define A32TH   0.03125
 
#define MMG5_MEMMIN   38
 
#define MMG5_PATHSEP   '/'
 
#define MMG5_NONSET_MEM   -1
 
#define MMG5_NONSET_HMIN   -1
 
#define MMG5_NONSET_HMAX   -1
 
#define MMG5_NONSET_HSIZ   -1
 
#define MMG5_NONSET   -1
 
#define MMG5_HAUSD   0.01
 
#define MMG5_HGRAD   0.26236426446
 
#define MMG5_HGRADREQ   0.83290912294
 
#define MMG5_NOHGRAD   -1
 
#define MMG5_LAG   -1
 
#define MMG5_NR   -1
 
#define MMG5_LS   0.0
 
#define MMG5_PROCTREE   32
 
#define MMG5_OFF   0
 
#define MMG5_ON   1
 
#define MMG5_GAP   0.2
 
#define MMG5_HMINCOE   0.001
 
#define MMG5_HMAXCOE   2
 
#define MMG5_HMINMAXGAP   5
 
#define MMG5_FEM   1
 
#define MMG5_FILESTR_LGTH   128 /** Maximal length of a line in input file */
 
#define MG_MAX(a, b)   (((a) > (b)) ? (a) : (b))
 
#define MG_MIN(a, b)   (((a) < (b)) ? (a) : (b))
 
#define MG_NOTAG   (0)
 
#define MG_REF   (1 << 0)
 
#define MG_GEO   (1 << 1)
 
#define MG_REQ   (1 << 2)
 
#define MG_NOM   (1 << 3)
 
#define MG_BDY   (1 << 4)
 
#define MG_CRN   (1 << 5)
 
#define MG_NOSURF   (1 << 6)
 
#define MG_OPNBDY   (1 << 7)
 
#define MG_OLDPARBDY   (1 << 11)
 
#define MG_PARBDYBDY   (1 << 12)
 
#define MG_PARBDY   (1 << 13)
 
#define MG_NUL   (1 << 14)
 
#define MG_Vert   (1 << 0 )
 
#define MG_Tria   (1 << 1 )
 
#define MG_Tetra   (1 << 2 )
 
#define MG_Edge   (1 << 3 )
 
#define MG_VOK(ppt)   (ppt && ((ppt)->tag < MG_NUL))
 
#define MG_EOK(pt)   (pt && ((pt)->v[0] > 0))
 
#define MG_SIN(tag)   ((tag & MG_CRN) || (tag & MG_REQ))
 
#define MG_SIN_OR_NOM(tag)   ( MG_SIN(tag) || (tag & MG_NOM) )
 
#define MG_RID(tag)   ( ( !( MG_SIN_OR_NOM(tag)) ) && ( tag & MG_GEO ) )
 
#define MG_EDG(tag)   ((tag & MG_GEO) || (tag & MG_REF))
 
#define MG_GEO_OR_NOM(tag)   (( tag & MG_GEO ) || ( tag & MG_NOM ))
 
#define MG_EDG_OR_NOM(tag)   ( MG_EDG(tag) || (tag & MG_NOM ) )
 
#define MG_TRUE_BDY(tag)   ( (tag & MG_BDY) && !(tag & MG_PARBDY) )
 
#define MG_SET(flag, bit)   ((flag) |= (1 << (bit)))
 
#define MG_CLR(flag, bit)   ((flag) &= ~(1 << (bit)))
 
#define MG_GET(flag, bit)   ((flag) & (1 << (bit)))
 
#define MMG5_KA   7
 
#define MMG5_KB   11
 
#define MMG5_SW   4
 
#define MMG5_SD   8
 
#define _LIBMMG5_RETURN(mesh, sol, met, val)
 
#define MMG5_CHK_MEM(mesh, size, string, law)
 
#define MMG5_DEL_MEM(mesh, ptr)
 
#define MMG5_ADD_MEM(mesh, size, message, law)
 
#define MMG5_SAFE_FREE(ptr)
 
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
 
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
 
#define MMG5_SAFE_REALLOC(ptr, prevSize, newSize, type, message, law)
 
#define MMG5_CHK_INT32_OVERFLOW(wantedGap, oldSiz, coef, shift, law)
 
#define MMG5_SAFE_RECALLOC(ptr, prevSize, newSize, type, message, law)
 
#define MMG5_TAB_RECALLOC(mesh, ptr, initSize, wantedGap, type, message, law)
 
#define MMG5_INCREASE_MEM_MESSAGE()
 
#define MMG5_SAFELL2LCAST(longlongval)   (((longlongval) > (LONG_MAX)) ? 0 : ((long)(longlongval)))
 
#define MMG5_SAFELL2ICAST(longlongval)   (((longlongval) > (INT_MAX)) ? 0 : ((int)(longlongval)))
 
#define MMG_FREAD(ptr, size, count, stream)
 
#define CV_VA_NUM_ARGS_HELPER(_1, _2, _3, _4, _5, _6, _7, _8, _9, _10, N, ...)   N
 
#define CV_VA_NUM_ARGS(...)   CV_VA_NUM_ARGS_HELPER(__VA_ARGS__, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0)
 
#define MMG_FSCANF(stream, format, ...)
 
#define FUNCTION_POINTER(fproto)    MMG_EXTERN fproto MMG_ASSIGN_NULL
 
#define FORTRAN_NAME(nu, nl, pl, pc)
 Adds function definitions. More...
 
#define FORTRAN_VARIADIC(nu, nl, pl, body)
 Adds function definitions. More...
 

Typedefs

typedef MMG5_BezierMMG5_pBezier
 
typedef struct MMG5_iNode_s MMG5_iNode
 

Functions

static void * mycalloc (size_t c, size_t s)
 
static void * mymalloc (size_t s)
 
static void * myrealloc (void *ptr_in, size_t s, size_t oldsize)
 
static size_t myfree (void *ptr)
 
static void MMG5_warnScotch (MMG5_pMesh mesh)
 
static void MMG5_excfun (int sigid)
 
void MMG5_version (MMG5_pMesh, char *)
 Functions needed by libraries API. More...
 
void MMG5_nsort (int8_t, double *, int8_t *)
 
void MMG5_nperm (int8_t n, int8_t shift, int8_t stride, double *val, double *oldval, int8_t *perm)
 
double MMG5_det3pt1vec (double c0[3], double c1[3], double c2[3], double v[3])
 
double MMG5_det4pt (double c0[3], double c1[3], double c2[3], double c3[3])
 
int MMG5_devangle (double *n1, double *n2, double crit)
 
double MMG5_orvol (MMG5_pPoint point, MMG5_int *v)
 
int MMG5_Add_inode (MMG5_pMesh mesh, MMG5_iNode **liLi, int val)
 
int MMG5_eigenvmatsym2d (MMG5_pMesh mesh, double m[], double lambda[], double v[][2])
 
int MMG5_eigenvmatsym3d (MMG5_pMesh mesh, double m[], double lambda[], double v[][3])
 
int MMG5_eigenvmatnonsym2d (MMG5_pMesh mesh, double m[], double lambda[], double v[][2])
 
int MMG5_eigenvmatnonsym3d (MMG5_pMesh mesh, double m[], double lambda[], double v[][3])
 
void MMG5_bezierEdge (MMG5_pMesh, MMG5_int, MMG5_int, double *, double *, int8_t, double *)
 
int MMG5_buildridmet (MMG5_pMesh, MMG5_pSol, MMG5_int, double, double, double, double *, double[3][3])
 
int MMG5_buildridmetfic (MMG5_pMesh, double *, double *, double, double, double, double *)
 
int MMG5_buildridmetnor (MMG5_pMesh, MMG5_pSol, MMG5_int, double *, double *, double[3][3])
 
void MMG5_check_hminhmax (MMG5_pMesh mesh, int8_t sethmin, int8_t sethmax)
 
int MMG5_paratmet (double c0[3], double n0[3], double m[6], double c1[3], double n1[3], double mt[6])
 
void MMG5_transpose3d (double m[3][3])
 
void MMG5_dotprod (int8_t dim, double *a, double *b, double *result)
 
void MMG5_crossprod3d (double *a, double *b, double *result)
 
void MMG5_mn (double m[6], double n[6], double mn[9])
 
int MMG5_rmtr (double r[3][3], double m[6], double mr[6])
 
int MMG5_boundingBox (MMG5_pMesh mesh)
 
int MMG5_boulep (MMG5_pMesh mesh, MMG5_int start, int ip, MMG5_int *, MMG5_int *list, MMG5_int *tlist)
 
int MMG5_boulec (MMG5_pMesh, MMG5_int *, MMG5_int, int ip, double *tt)
 
int MMG5_boulen (MMG5_pMesh, MMG5_int *, MMG5_int, int ip, double *nn)
 
int MMG5_bouler (MMG5_pMesh, MMG5_int *, MMG5_int, int ip, MMG5_int *, MMG5_int *, int *, int *, int)
 
int MMG5_boulet (MMG5_pMesh mesh, MMG5_int start, int ip, MMG5_int *list, int8_t s, int8_t *opn)
 
double MMG5_caltri33_ani (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria pt)
 
double MMG5_caltri_ani (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
 
double MMG5_caltri_iso (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
 
void MMG5_defUninitSize (MMG5_pMesh mesh, MMG5_pSol met, int8_t ismet)
 
void MMG5_displayLengthHisto (MMG5_pMesh, MMG5_int, double *, MMG5_int, MMG5_int, double, MMG5_int, MMG5_int, double, int, double *, MMG5_int *, int8_t)
 
void MMG5_displayLengthHisto_internal (MMG5_int, MMG5_int, MMG5_int, double, MMG5_int, MMG5_int, double, MMG5_int, double *, MMG5_int *, int8_t, int)
 
short MMG5_dikmov (MMG5_pMesh, MMG5_pSol, short *, short, MMG5_int chkmovmesh(MMG5_pMesh, MMG5_pSol, short, MMG5_int *))
 common functions for lagrangian meshing. More...
 
int MMG5_minQualCheck (MMG5_int iel, double minqual, double alpha)
 
int MMG5_elementWeight (MMG5_pMesh, MMG5_pSol, MMG5_pTria, MMG5_pPoint, MMG5_Bezier *, double r[3][3], double gv[2])
 
void MMG5_fillDefmetregSys (MMG5_int, MMG5_pPoint, int, MMG5_Bezier, double r[3][3], double *, double *, double *, double *)
 
void MMG5_Free_ilinkedList (MMG5_pMesh mesh, MMG5_iNode *liLi)
 
MMG5_int MMG5_grad2metSurf (MMG5_pMesh, MMG5_pSol, MMG5_pTria, MMG5_int, MMG5_int)
 
int MMG5_grad2metSurfreq (MMG5_pMesh, MMG5_pSol, MMG5_pTria, MMG5_int, MMG5_int)
 
MMG5_int MMG5_hashFace (MMG5_pMesh, MMG5_Hash *, MMG5_int, MMG5_int, MMG5_int, MMG5_int)
 
int MMG5_hashEdge (MMG5_pMesh mesh, MMG5_Hash *hash, MMG5_int a, MMG5_int b, MMG5_int k)
 
int MMG5_hashUpdate (MMG5_Hash *hash, MMG5_int a, MMG5_int b, MMG5_int k)
 
int MMG5_hashEdgeTag (MMG5_pMesh mesh, MMG5_Hash *hash, MMG5_int a, MMG5_int b, int16_t k)
 
MMG5_int MMG5_hashGet (MMG5_Hash *hash, MMG5_int a, MMG5_int b)
 
int MMG5_hashNew (MMG5_pMesh mesh, MMG5_Hash *hash, MMG5_int hsiz, MMG5_int hmax)
 
int MMG5_intmetsavedir (MMG5_pMesh mesh, double *m, double *n, double *mr)
 
int MMG5_intridmet (MMG5_pMesh, MMG5_pSol, MMG5_int, MMG5_int, double, double *, double *)
 
int MMG5_mmgIntmet33_ani (double *, double *, double *, double)
 
int MMG5_ismaniball (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int start, int8_t istart)
 
int MMG5_mmgIntextmet (MMG5_pMesh, MMG5_pSol, MMG5_int, double *, double *)
 
size_t MMG5_memSize (void)
 
void MMG5_memOption_memSet (MMG5_pMesh mesh)
 
void MMG5_mmgDefaultValues (MMG5_pMesh mesh)
 
int MMG5_mmgHashTria (MMG5_pMesh mesh, MMG5_int *adja, MMG5_Hash *, int chkISO)
 
void MMG5_mmgInit_parameters (MMG5_pMesh mesh)
 
void MMG5_mmgUsage (char *prog)
 
void MMG5_paramUsage1 (void)
 
void MMG5_paramUsage2 (void)
 
void MMG5_2d3dUsage (void)
 
void MMG5_lagUsage (void)
 
void MMG5_advancedUsage (void)
 
int MMG5_nonUnitNorPts (MMG5_pMesh, MMG5_int, MMG5_int, MMG5_int, double *)
 
double MMG5_nonorsurf (MMG5_pMesh mesh, MMG5_pTria pt)
 
int MMG5_norpts (MMG5_pMesh, MMG5_int, MMG5_int, MMG5_int, double *)
 
int MMG5_nortri (MMG5_pMesh mesh, MMG5_pTria pt, double *n)
 
void MMG5_printTria (MMG5_pMesh mesh, char *fileName)
 
int MMG5_rotmatrix (double n[3], double r[3][3])
 
int MMG5_invmat (double *m, double *mi)
 
int MMG5_invmatg (double m[9], double mi[9])
 
int MMG5_invmat33 (double m[3][3], double mi[3][3])
 
int MMG5_invmat22 (double m[2][2], double mi[2][2])
 
int MMG5_regnor (MMG5_pMesh mesh)
 
double MMG5_ridSizeInNormalDir (MMG5_pMesh, int, double *, MMG5_pBezier, double, double)
 
double MMG5_ridSizeInTangentDir (MMG5_pMesh, MMG5_pPoint, MMG5_int, MMG5_int *, double, double)
 
int MMG5_scale_meshAndSol (MMG5_pMesh, MMG5_pSol, MMG5_pSol, double *)
 
int MMG5_scale_scalarMetric (MMG5_pMesh, MMG5_pSol, double)
 
int MMG5_scotchCall (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pSol fields, MMG5_int *)
 
int MMG5_check_setted_hminhmax (MMG5_pMesh mesh)
 
int MMG5_solTruncature_iso (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_2dSolTruncature_ani (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_3dSolTruncature_ani (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_truncate_met3d (MMG5_pSol met, MMG5_int ip, double isqhmin, double isqhmax)
 
int MMG5_solveDefmetregSys (MMG5_pMesh, double r[3][3], double *, double *, double *, double *, double, double, double)
 
int MMG5_solveDefmetrefSys (MMG5_pMesh, MMG5_pPoint, MMG5_int *, double r[3][3], double *, double *, double *, double *, double, double, double)
 
double MMG5_surftri_ani (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
 
double MMG5_surftri33_ani (MMG5_pMesh, MMG5_pTria, double *, double *, double *)
 
double MMG5_surftri_iso (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
 
int MMG5_sys33sym (double a[6], double b[3], double r[3])
 
int MMG5_interpreg_ani (MMG5_pMesh, MMG5_pSol, MMG5_pTria, int8_t, double, double *mr)
 
int MMG5_interp_iso (double *ma, double *mb, double *mp, double t)
 
int MMG5_intersecmet22 (MMG5_pMesh mesh, double *m, double *n, double *mr)
 
int MMG5_countLocalParamAtTri (MMG5_pMesh, MMG5_iNode **)
 
int MMG5_writeLocalParamAtTri (MMG5_pMesh, MMG5_iNode *, FILE *)
 
double MMG2D_quickarea (double a[2], double b[2], double c[2])
 
void MMG5_build3DMetric (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ip, double dbuf[6])
 
int MMG5_loadVtuMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
 
int MMG5_loadMshMesh_part1 (MMG5_pMesh mesh, const char *filename, FILE **inm, long *posNodes, long *posElts, long **posNodeData, int *bin, int *iswp, MMG5_int *nelts, int *nsols)
 
int MMG5_check_readedMesh (MMG5_pMesh mesh, MMG5_int nref)
 
int MMG5_loadMshMesh_part2 (MMG5_pMesh mesh, MMG5_pSol *sol, FILE **inm, const long posNodes, const long posElts, const long *posNodeData, const int bin, const int iswp, const MMG5_int nelts, const int nsols)
 
int MMG5_saveMshMesh (MMG5_pMesh, MMG5_pSol *, const char *, int)
 
int MMG5_saveDisp (MMG5_pMesh, MMG5_pSol)
 
int MMG5_loadSolHeader (const char *, int, FILE **, int *, int *, int *, MMG5_int *, int *, int *, int **, long *, int)
 
int MMG5_chkMetricType (MMG5_pMesh mesh, int *type, int *, FILE *inm)
 
int MMG5_readFloatSol3D (MMG5_pSol, FILE *, int, int, int)
 
int MMG5_readDoubleSol3D (MMG5_pSol, FILE *, int, int, MMG5_int)
 
int MMG5_saveSolHeader (MMG5_pMesh, const char *, FILE **, int, int *, MMG5_int *, MMG5_int, int, int, int *, int *, int *)
 
int MMG5_saveSolAtTrianglesHeader (MMG5_pMesh, FILE *, int, int, MMG5_int *, int, int, int *, int *, int *)
 
int MMG5_saveSolAtTetrahedraHeader (MMG5_pMesh, FILE *, int, int, MMG5_int *, int, int, int *, int *, int *)
 
void MMG5_writeDoubleSol3D (MMG5_pMesh, MMG5_pSol, FILE *, int, MMG5_int, int)
 
void MMG5_printMetStats (MMG5_pMesh mesh, MMG5_pSol met)
 
void MMG5_printSolStats (MMG5_pMesh mesh, MMG5_pSol *sol)
 
int MMG5_defsiz_startingMessage (MMG5_pMesh, MMG5_pSol, const char *funcname)
 
void MMG5_gradation_info (MMG5_pMesh)
 
int MMG5_sum_reqEdgeLengthsAtPoint (MMG5_pMesh, MMG5_pSol, MMG5_int ip0, MMG5_int ip1)
 
int MMG5_compute_meanMetricAtMarkedPoints_iso (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_compute_meanMetricAtMarkedPoints_ani (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_reset_metricAtReqEdges_surf (MMG5_pMesh, MMG5_pSol, int8_t)
 
void MMG5_mark_pointsOnReqEdge_fromTria (MMG5_pMesh mesh)
 
int MMG5_gradsiz_iso (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_gradsizreq_iso (MMG5_pMesh, MMG5_pSol)
 
MMG5_int MMG5_gradsiz_ani (MMG5_pMesh mesh, MMG5_pSol met, int *it)
 
int MMG5_gradsizreq_ani (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_simred2d (MMG5_pMesh, double *, double *, double dm[2], double dn[2], double vp[2][2])
 
int MMG5_simred3d (MMG5_pMesh mesh, double *m, double *n, double dm[3], double dn[3], double vp[3][3])
 
int MMG5_updatemet2d_ani (double *m, double *n, double dm[2], double dn[2], double vp[2][2], int8_t ier)
 
int MMG5_updatemet3d_ani (double *m, double *n, double dm[3], double dn[3], double vp[3][3], int8_t ier)
 
void MMG5_gradEigenvreq (double *dm, double *dn, double, int8_t, int8_t *)
 
int MMG5_updatemetreq_ani (double *n, double dn[2], double vp[2][2])
 
int MMG5_swapbin (int sbin)
 
MMG5_int MMG5_swapbin_int (MMG5_int sbin)
 
float MMG5_swapf (float sbin)
 
double MMG5_swapd (double sbin)
 
int MMG5_MultiMat_init (MMG5_pMesh)
 
int MMG5_isLevelSet (MMG5_pMesh, MMG5_int, MMG5_int)
 
int MMG5_isSplit (MMG5_pMesh, MMG5_int, MMG5_int *, MMG5_int *)
 
int MMG5_isNotSplit (MMG5_pMesh, MMG5_int)
 
int MMG5_getStartRef (MMG5_pMesh, MMG5_int, MMG5_int *)
 
int MMG5_snpval_ls (MMG5_pMesh mesh, MMG5_pSol sol)
 
int MMG5_snpval_lssurf (MMG5_pMesh mesh, MMG5_pSol sol)
 
int MMG5_rmc (MMG5_pMesh, MMG5_pSol)
 
int MMG5_resetRef_ls (MMG5_pMesh)
 
int MMG5_resetRef_lssurf (MMG5_pMesh)
 
int MMG5_setref_ls (MMG5_pMesh mesh, MMG5_pSol sol)
 
int MMG5_setref_lssurf (MMG5_pMesh mesh, MMG5_pSol sol)
 
int MMG5_chkmaniball (MMG5_pMesh mesh, MMG5_int start, int8_t istart)
 
int MMG5_chkmanimesh (MMG5_pMesh mesh)
 
double MMG5_test_mat_error (int8_t nelem, double m1[], double m2[])
 
int MMG5_test_invmat22 ()
 
int MMG5_test_invmat33 ()
 
int MMG5_test_eigenvmatsym2d (MMG5_pMesh mesh, double *mex, double lambdaex[], double vpex[][2])
 
int MMG5_test_eigenvmatnonsym2d (MMG5_pMesh mesh, double *mex, double lambdaex[], double vpex[][2], double ivpex[][2])
 
int MMG5_test_eigenvmatsym3d (MMG5_pMesh mesh, double *mex, double lambdaex[], double vpex[][3])
 
int MMG5_test_eigenvmatnonsym3d (MMG5_pMesh mesh, double *mex, double lambdaex[], double vpex[][3], double ivpex[][3])
 
int MMG5_test_transpose3d ()
 
int MMG5_test_dotprod ()
 
int MMG5_test_crossprod3d ()
 
int MMG5_test_mn ()
 
int MMG5_test_rmtr ()
 
int MMG5_test_rotmatrix ()
 
int MMG5_test_simred2d (MMG5_pMesh mesh, double *mex, double *nex, double *dmex, double *dnex, double vpex[][2])
 
int MMG5_test_simred3d (MMG5_pMesh mesh, double *mex, double *nex, double *dmex, double *dnex, double vpex[][3])
 
int MMG5_test_updatemet2d_ani ()
 
int MMG5_test_updatemet3d_ani ()
 
int MMG5_test_intersecmet22 (MMG5_pMesh mesh)
 
int MMG5_test_intersecmet33 (MMG5_pMesh mesh)
 
void MMG5_mark_verticesAsUnused (MMG5_pMesh mesh)
 
void MMG5_mark_usedVertices (MMG5_pMesh mesh, void(*delPt)(MMG5_pMesh, MMG5_int))
 
void MMG5_keep_subdomainElts (MMG5_pMesh, int, int(*delElt)(MMG5_pMesh, MMG5_int))
 
void MMG5_Set_commonFunc (void)
 

Variables

static const uint8_t MMG5_inxt2 [6] = {1,2,0,1,2}
 
static const uint8_t MMG5_iprv2 [3] = {2,0,1}
 

Macro Definition Documentation

◆ _LIBMMG5_RETURN

#define _LIBMMG5_RETURN (   mesh,
  sol,
  met,
  val 
)
Value:
do \
{ \
signal(SIGABRT,SIG_DFL); \
signal(SIGFPE,SIG_DFL); \
signal(SIGILL,SIG_DFL); \
signal(SIGSEGV,SIG_DFL); \
signal(SIGTERM,SIG_DFL); \
signal(SIGINT,SIG_DFL); \
mesh->npi = mesh->np; \
mesh->nti = mesh->nt; \
mesh->nai = mesh->na; \
mesh->nei = mesh->ne; \
mesh->xt = 0; \
if ( sol ) { sol->npi = sol->np; } \
if ( met ) { met->npi = met->np; } \
return val; \
}while(0)
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh * mesh
MMG5_int xt
Definition: libmmgtypes.h:620
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_int nei
Definition: libmmgtypes.h:612
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_int nti
Definition: libmmgtypes.h:612
MMG5_int npi
Definition: libmmgtypes.h:612
MMG5_int nai
Definition: libmmgtypes.h:612
MMG5_int na
Definition: libmmgtypes.h:612

Reset the customized signals and set the internal counters of points, edges, tria and tetra to the suitable value (needed by users to recover their mesh using the API)

Definition at line 192 of file mmgcommon_private.h.

◆ A16TH

#define A16TH   0.0625

Definition at line 108 of file mmgcommon_private.h.

◆ A32TH

#define A32TH   0.03125

Definition at line 109 of file mmgcommon_private.h.

◆ A64TH

#define A64TH   0.015625

Definition at line 107 of file mmgcommon_private.h.

◆ CV_VA_NUM_ARGS

#define CV_VA_NUM_ARGS (   ...)    CV_VA_NUM_ARGS_HELPER(__VA_ARGS__, 10, 9, 8, 7, 6, 5, 4, 3, 2, 1, 0)

count the number of variadic arguments provided to the macro

Definition at line 453 of file mmgcommon_private.h.

◆ CV_VA_NUM_ARGS_HELPER

#define CV_VA_NUM_ARGS_HELPER (   _1,
  _2,
  _3,
  _4,
  _5,
  _6,
  _7,
  _8,
  _9,
  _10,
  N,
  ... 
)    N

macro to help to count the number of variadic arguments

Definition at line 450 of file mmgcommon_private.h.

◆ FORTRAN_NAME

#define FORTRAN_NAME (   nu,
  nl,
  pl,
  pc 
)
Value:
void nu pl; \
void nl pl \
{ nu pc; } \
void nl##_ pl \
{ nu pc; } \
void nl##__ pl \
{ nu pc; } \
void nu pl

Adds function definitions.

Parameters
nufunction name in upper case.
nlfunction name in lower case.
pltype of arguments.
pcname of arguments.
Note
Macro coming from Scotch library.

Adds function definitions with upcase, underscore and double underscore to match any fortran compiler.

Definition at line 549 of file mmgcommon_private.h.

◆ FORTRAN_VARIADIC

#define FORTRAN_VARIADIC (   nu,
  nl,
  pl,
  body 
)
Value:
void nu pl \
{ body } \
void nl pl \
{ body } \
void nl##_ pl \
{ body } \
void nl##__ pl \
{ body } \

Adds function definitions.

Parameters
nufunction name in upper case.
nlfunction name in lower case.
pltype of arguments.
bodybody of the function.

Adds function definitions with upcase, underscore and double underscore to match any fortran compiler.

Definition at line 571 of file mmgcommon_private.h.

◆ FUNCTION_POINTER

#define FUNCTION_POINTER (   fproto)     MMG_EXTERN fproto MMG_ASSIGN_NULL
Parameters
fprotofunction prototype

Expand automatically prototype of function pointer in .h/.c files depending on the definition of the MMG_EXTERN and MMG_ASSIGN_NULL preprocessor variables:

  • MMG_EXTERN is setted to "extern" in the .h file and empty in the .c one;
  • MMG_ASSIGN_NULL is empty in .h file and setted to =NULL in .c one.

Definition at line 530 of file mmgcommon_private.h.

◆ MG_BDY

#define MG_BDY   (1 << 4)

16 boundary entity

Definition at line 148 of file mmgcommon_private.h.

◆ MG_CLR

#define MG_CLR (   flag,
  bit 
)    ((flag) &= ~(1 << (bit)))

bit number bit is set to 0

Definition at line 179 of file mmgcommon_private.h.

◆ MG_CRN

#define MG_CRN   (1 << 5)

32 corner

Definition at line 149 of file mmgcommon_private.h.

◆ MG_EDG

#define MG_EDG (   tag)    ((tag & MG_GEO) || (tag & MG_REF))

Edge or Ridge

Definition at line 171 of file mmgcommon_private.h.

◆ MG_EDG_OR_NOM

#define MG_EDG_OR_NOM (   tag)    ( MG_EDG(tag) || (tag & MG_NOM ) )

Edge, ridge or non-manifold

Definition at line 173 of file mmgcommon_private.h.

◆ MG_Edge

#define MG_Edge   (1 << 3 )

8 local parameter applied over edge

Definition at line 161 of file mmgcommon_private.h.

◆ MG_EOK

#define MG_EOK (   pt)    (pt && ((pt)->v[0] > 0))

Element OK

Definition at line 165 of file mmgcommon_private.h.

◆ MG_GEO

#define MG_GEO   (1 << 1)

2 geometric ridge

Definition at line 145 of file mmgcommon_private.h.

◆ MG_GEO_OR_NOM

#define MG_GEO_OR_NOM (   tag)    (( tag & MG_GEO ) || ( tag & MG_NOM ))

Ridge or non-manifold

Definition at line 172 of file mmgcommon_private.h.

◆ MG_GET

#define MG_GET (   flag,
  bit 
)    ((flag) & (1 << (bit)))

return bit number bit value

Definition at line 180 of file mmgcommon_private.h.

◆ MG_MAX

#define MG_MAX (   a,
 
)    (((a) > (b)) ? (a) : (b))

Definition at line 139 of file mmgcommon_private.h.

◆ MG_MIN

#define MG_MIN (   a,
 
)    (((a) < (b)) ? (a) : (b))

Definition at line 140 of file mmgcommon_private.h.

◆ MG_NOM

#define MG_NOM   (1 << 3)

8 non manifold

Definition at line 147 of file mmgcommon_private.h.

◆ MG_NOSURF

#define MG_NOSURF   (1 << 6)

64 freezed boundary

Definition at line 150 of file mmgcommon_private.h.

◆ MG_NOTAG

#define MG_NOTAG   (0)

Definition at line 143 of file mmgcommon_private.h.

◆ MG_NUL

#define MG_NUL   (1 << 14)

16384 vertex removed

Definition at line 155 of file mmgcommon_private.h.

◆ MG_OLDPARBDY

#define MG_OLDPARBDY   (1 << 11)

2048 old parallel boundary

Definition at line 152 of file mmgcommon_private.h.

◆ MG_OPNBDY

#define MG_OPNBDY   (1 << 7)

128 open boundary

Definition at line 151 of file mmgcommon_private.h.

◆ MG_PARBDY

#define MG_PARBDY   (1 << 13)

8192 parallel boundary

Definition at line 154 of file mmgcommon_private.h.

◆ MG_PARBDYBDY

#define MG_PARBDYBDY   (1 << 12)

4096 parallel boundary over a boundary

Definition at line 153 of file mmgcommon_private.h.

◆ MG_REF

#define MG_REF   (1 << 0)

1 edge reference

Definition at line 144 of file mmgcommon_private.h.

◆ MG_REQ

#define MG_REQ   (1 << 2)

4 required entity

Definition at line 146 of file mmgcommon_private.h.

◆ MG_RID

#define MG_RID (   tag)    ( ( !( MG_SIN_OR_NOM(tag)) ) && ( tag & MG_GEO ) )

Non-singular ridge point (so ridge metric in aniso mode)

Definition at line 169 of file mmgcommon_private.h.

◆ MG_SET

#define MG_SET (   flag,
  bit 
)    ((flag) |= (1 << (bit)))

bit number bit is set to 1

Definition at line 178 of file mmgcommon_private.h.

◆ MG_SIN

#define MG_SIN (   tag)    ((tag & MG_CRN) || (tag & MG_REQ))

Corner or Required

Definition at line 167 of file mmgcommon_private.h.

◆ MG_SIN_OR_NOM

#define MG_SIN_OR_NOM (   tag)    ( MG_SIN(tag) || (tag & MG_NOM) )

Corner, Required or non-manifold

Definition at line 168 of file mmgcommon_private.h.

◆ MG_SMSGN

#define MG_SMSGN (   a,
 
)    (((double)(a)*(double)(b) > (0.0)) ? (1) : (0))

Check if a and b have the same sign

Definition at line 70 of file mmgcommon_private.h.

◆ MG_STR

#define MG_STR   "&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&"

Definition at line 54 of file mmgcommon_private.h.

◆ MG_Tetra

#define MG_Tetra   (1 << 2 )

4 local parameter applied over tetrahedron

Definition at line 160 of file mmgcommon_private.h.

◆ MG_Tria

#define MG_Tria   (1 << 1 )

2 local parameter applied over triangle

Definition at line 159 of file mmgcommon_private.h.

◆ MG_TRUE_BDY

#define MG_TRUE_BDY (   tag)    ( (tag & MG_BDY) && !(tag & MG_PARBDY) )

Non parbdy boundary point (true bdy)

Definition at line 174 of file mmgcommon_private.h.

◆ MG_Vert

#define MG_Vert   (1 << 0 )

1 local parameter applied over vertex

Definition at line 158 of file mmgcommon_private.h.

◆ MG_VOK

#define MG_VOK (   ppt)    (ppt && ((ppt)->tag < MG_NUL))

Vertex OK

Definition at line 164 of file mmgcommon_private.h.

◆ MMG5_ADD_MEM

#define MMG5_ADD_MEM (   mesh,
  size,
  message,
  law 
)
Value:
do \
{ \
(mesh)->memCur += (size); \
MMG5_CHK_MEM(mesh,size,message,law); \
}while(0)

Increment memory counter memCur and check if we don't overflow the maximum authorizied memory memMax.

Definition at line 300 of file mmgcommon_private.h.

◆ MMG5_ANGEDG

#define MMG5_ANGEDG   0.707106781186548 /*0.573576436351046 */

Definition at line 90 of file mmgcommon_private.h.

◆ MMG5_ANGLIM

#define MMG5_ANGLIM   -0.999999

Definition at line 91 of file mmgcommon_private.h.

◆ MMG5_ATHIRD

#define MMG5_ATHIRD   0.333333333333333

Definition at line 92 of file mmgcommon_private.h.

◆ MMG5_BITWIZE_MB_TO_B

#define MMG5_BITWIZE_MB_TO_B   20

Bitwise convertion from Mo to O

Definition at line 77 of file mmgcommon_private.h.

◆ MMG5_BOXSIZE

#define MMG5_BOXSIZE   500

size of box for renumbering with scotch.

Definition at line 73 of file mmgcommon_private.h.

◆ MMG5_CHK_INT32_OVERFLOW

#define MMG5_CHK_INT32_OVERFLOW (   wantedGap,
  oldSiz,
  coef,
  shift,
  law 
)
Value:
do \
{ \
/* Check for int32 overflow */ \
if ( sizeof(MMG5_int) == sizeof(int32_t) ) { \
MMG5_int gap_loc = (MMG5_int)((wantedGap) * (oldSiz)); \
if ( !gap_loc ) gap_loc = 1; \
\
int32_t max_ne = (INT32_MAX-(shift))/(coef); \
if ( max_ne < (oldSiz)+gap_loc ) { \
/* Detected overflow, target maximal possible size */ \
gap_loc = max_ne-(oldSiz); \
if ( gap_loc <=0 ) { \
/* No possibe realloc without int overflow */ \
fprintf(stderr," ## Error: %s: %d: Unable to reallocate adja array" \
" without int overflow.\n",__func__,__LINE__); \
gap_loc = 0; \
law; \
} \
else { \
wantedGap = (float)gap_loc/(float)oldSiz; \
printf("wantGap has been modified %15f\n",wantedGap); \
wantedGap = (double)gap_loc/(double)oldSiz; \
printf("DwantGap has been modified %15fl\n",wantedGap); \
} \
} \
} \
}while(0)

Check for int32 overflow when trying to reallocate coef*oldSiz+shift array

Definition at line 351 of file mmgcommon_private.h.

◆ MMG5_CHK_MEM

#define MMG5_CHK_MEM (   mesh,
  size,
  string,
  law 
)
Value:
do \
{ \
if ( (mesh)->memCur > (mesh)->memMax ) { \
fprintf(stderr," ## Error:"); \
fprintf(stderr," unable to allocate %s.\n",string); \
fprintf(stderr," ## Check the mesh size or "); \
fprintf(stderr,"increase maximal authorized memory with the -m option.\n"); \
(mesh)->memCur -= (size); \
law; \
} \
}while(0)

Check if used memory overflow maximal authorized memory. Execute the command law if lack of memory.

Definition at line 213 of file mmgcommon_private.h.

◆ MMG5_DEL_MEM

#define MMG5_DEL_MEM (   mesh,
  ptr 
)
Value:
do \
{ \
size_t size_to_free = myfree(ptr); \
(mesh)->memCur -= size_to_free; \
ptr = NULL; \
}while(0)
static size_t myfree(void *ptr)

Free pointer ptr of mesh structure and compute the new used memory.

Definition at line 291 of file mmgcommon_private.h.

◆ MMG5_DISPREF

#define MMG5_DISPREF   10

Definition at line 84 of file mmgcommon_private.h.

◆ MMG5_EPS

#define MMG5_EPS   1.e-06

Definition at line 96 of file mmgcommon_private.h.

◆ MMG5_EPSD

#define MMG5_EPSD   1.e-30

Definition at line 94 of file mmgcommon_private.h.

◆ MMG5_EPSD2

#define MMG5_EPSD2   1.0e-200

Definition at line 95 of file mmgcommon_private.h.

◆ MMG5_EPSOK

#define MMG5_EPSOK   1.e-15

Definition at line 97 of file mmgcommon_private.h.

◆ MMG5_FEM

#define MMG5_FEM   1

defaut value for FEM mode

Definition at line 135 of file mmgcommon_private.h.

◆ MMG5_FILESTR_LGTH

#define MMG5_FILESTR_LGTH   128 /** Maximal length of a line in input file */

Definition at line 136 of file mmgcommon_private.h.

◆ MMG5_GAP

#define MMG5_GAP   0.2

gap value for reallocation

Definition at line 131 of file mmgcommon_private.h.

◆ MMG5_HAUSD

#define MMG5_HAUSD   0.01

default value for hausdorff param

Definition at line 121 of file mmgcommon_private.h.

◆ MMG5_HGRAD

#define MMG5_HGRAD   0.26236426446

default value for gradation (1.3)

Definition at line 122 of file mmgcommon_private.h.

◆ MMG5_HGRADREQ

#define MMG5_HGRADREQ   0.83290912294

default value for required gradation (2.3)

Definition at line 123 of file mmgcommon_private.h.

◆ MMG5_HMAXCOE

#define MMG5_HMAXCOE   2

coefficient to compute the default hmax value

Definition at line 133 of file mmgcommon_private.h.

◆ MMG5_HMINCOE

#define MMG5_HMINCOE   0.001

coefficient to compute the default hmin value

Definition at line 132 of file mmgcommon_private.h.

◆ MMG5_HMINMAXGAP

#define MMG5_HMINMAXGAP   5

imposed gap between hmin and hmax if hmax<hmin

Definition at line 134 of file mmgcommon_private.h.

◆ MMG5_INCREASE_MEM_MESSAGE

#define MMG5_INCREASE_MEM_MESSAGE ( )
Value:
do \
{ \
printf(" ## Check the mesh size or increase maximal"); \
printf(" authorized memory with the -m option.\n"); \
} while(0)

Error message when lack of memory

Definition at line 430 of file mmgcommon_private.h.

◆ MMG5_KA

#define MMG5_KA   7

Key for hash tables.

Definition at line 182 of file mmgcommon_private.h.

◆ MMG5_KB

#define MMG5_KB   11

Key for hash tables.

Definition at line 183 of file mmgcommon_private.h.

◆ MMG5_LAG

#define MMG5_LAG   -1

default value for lagrangian option

Definition at line 125 of file mmgcommon_private.h.

◆ MMG5_LMAX

#define MMG5_LMAX   10240

Maximum array size when calling boulep.

Definition at line 64 of file mmgcommon_private.h.

◆ MMG5_LPARMAX

#define MMG5_LPARMAX   200

Maximal number of local parameters

Definition at line 67 of file mmgcommon_private.h.

◆ MMG5_LS

#define MMG5_LS   0.0

default level-set to discretize

Definition at line 127 of file mmgcommon_private.h.

◆ MMG5_MEMMAX

#define MMG5_MEMMAX   800

Maximal memory used if available memory compitation fail. Default mem if unable to compute memMax

Definition at line 76 of file mmgcommon_private.h.

◆ MMG5_MEMMIN

#define MMG5_MEMMIN   38

minimal memory needed to store the mesh/sol names

Definition at line 111 of file mmgcommon_private.h.

◆ MMG5_MEMPERCENT

#define MMG5_MEMPERCENT   0.5

Percent of RAM used by default

Definition at line 78 of file mmgcommon_private.h.

◆ MMG5_MILLION

#define MMG5_MILLION   1048576

Definition at line 87 of file mmgcommon_private.h.

◆ MMG5_NOHGRAD

#define MMG5_NOHGRAD   -1

disable gradation

Definition at line 124 of file mmgcommon_private.h.

◆ MMG5_NONSET

#define MMG5_NONSET   -1

Definition at line 120 of file mmgcommon_private.h.

◆ MMG5_NONSET_HMAX

#define MMG5_NONSET_HMAX   -1

hmax value for unspecified hmax size

Definition at line 118 of file mmgcommon_private.h.

◆ MMG5_NONSET_HMIN

#define MMG5_NONSET_HMIN   -1

hmin value for unspecified hmin size

Definition at line 117 of file mmgcommon_private.h.

◆ MMG5_NONSET_HSIZ

#define MMG5_NONSET_HSIZ   -1

hsiz value for unspecified hsiz map

Definition at line 119 of file mmgcommon_private.h.

◆ MMG5_NONSET_MEM

#define MMG5_NONSET_MEM   -1

mem value for unspecified max memory

Definition at line 116 of file mmgcommon_private.h.

◆ MMG5_NR

#define MMG5_NR   -1

no ridge detection

Definition at line 126 of file mmgcommon_private.h.

◆ MMG5_NULKAL

#define MMG5_NULKAL   1.e-30

Definition at line 98 of file mmgcommon_private.h.

◆ MMG5_OFF

#define MMG5_OFF   0

0

Definition at line 129 of file mmgcommon_private.h.

◆ MMG5_ON

#define MMG5_ON   1

1

Definition at line 130 of file mmgcommon_private.h.

◆ MMG5_PATHSEP

#define MMG5_PATHSEP   '/'

Definition at line 113 of file mmgcommon_private.h.

◆ MMG5_PROCTREE

#define MMG5_PROCTREE   32

default size of the PROctree

Definition at line 128 of file mmgcommon_private.h.

◆ MMG5_SAFE_CALLOC

#define MMG5_SAFE_CALLOC (   ptr,
  size,
  type,
  law 
)
Value:
do \
{ \
ptr = (type*)mycalloc(size,sizeof(type)); \
if ( !ptr ) { \
perror(" ## Memory problem: calloc"); \
law; \
} \
}while(0)
static void * mycalloc(size_t c, size_t s)

Safe allocation with calloc

Definition at line 314 of file mmgcommon_private.h.

◆ MMG5_SAFE_FREE

#define MMG5_SAFE_FREE (   ptr)
Value:
do \
{ \
myfree(ptr); \
ptr = NULL; \
}while(0)

Safe deallocation

Definition at line 307 of file mmgcommon_private.h.

◆ MMG5_SAFE_MALLOC

#define MMG5_SAFE_MALLOC (   ptr,
  size,
  type,
  law 
)
Value:
do \
{ \
size_t size_to_allocate = (size)*sizeof(type); \
ptr = (type*)mymalloc(size_to_allocate); \
if ( !ptr ) { \
perror(" ## Memory problem: malloc"); \
law; \
} \
}while(0)
static void * mymalloc(size_t s)

Safe allocation with malloc

Definition at line 324 of file mmgcommon_private.h.

◆ MMG5_SAFE_REALLOC

#define MMG5_SAFE_REALLOC (   ptr,
  prevSize,
  newSize,
  type,
  message,
  law 
)
Value:
do \
{ \
type* tmp; \
size_t size_to_allocate = (newSize)*sizeof(type); \
\
tmp = (type *)myrealloc((ptr),size_to_allocate,(prevSize)*sizeof(type)); \
if ( !tmp ) { \
MMG5_SAFE_FREE(ptr); \
perror(" ## Memory problem: realloc"); \
law; \
} \
\
(ptr) = tmp; \
}while(0)
tmp[*strlen0]
static void * myrealloc(void *ptr_in, size_t s, size_t oldsize)

Safe reallocation

Definition at line 335 of file mmgcommon_private.h.

◆ MMG5_SAFE_RECALLOC

#define MMG5_SAFE_RECALLOC (   ptr,
  prevSize,
  newSize,
  type,
  message,
  law 
)
Value:
do \
{ \
type* tmp; \
size_t size_to_allocate = (newSize)*sizeof(type); \
\
tmp = (type *)myrealloc((ptr),size_to_allocate,(prevSize)*sizeof(type)); \
if ( !tmp ) { \
MMG5_SAFE_FREE(ptr); \
perror(" ## Memory problem: realloc"); \
law; \
} \
else { \
(ptr) = tmp; \
assert(ptr); \
if ( newSize > prevSize ) { \
memset(&((ptr)[prevSize]),0,((newSize)-(prevSize))*sizeof(type)); \
} \
} \
}while(0)

safe reallocation with memset at 0 for the new values of tab

Definition at line 381 of file mmgcommon_private.h.

◆ MMG5_SAFELL2ICAST

#define MMG5_SAFELL2ICAST (   longlongval)    (((longlongval) > (INT_MAX)) ? 0 : ((int)(longlongval)))

Definition at line 437 of file mmgcommon_private.h.

◆ MMG5_SAFELL2LCAST

#define MMG5_SAFELL2LCAST (   longlongval)    (((longlongval) > (LONG_MAX)) ? 0 : ((long)(longlongval)))

Definition at line 436 of file mmgcommon_private.h.

◆ MMG5_SD

#define MMG5_SD   8

Definition at line 187 of file mmgcommon_private.h.

◆ MMG5_SQR32

#define MMG5_SQR32   0.866025403784439

Definition at line 100 of file mmgcommon_private.h.

◆ MMG5_SW

#define MMG5_SW   4

Definition at line 186 of file mmgcommon_private.h.

◆ MMG5_TAB_RECALLOC

#define MMG5_TAB_RECALLOC (   mesh,
  ptr,
  initSize,
  wantedGap,
  type,
  message,
  law 
)
Value:
do \
{ \
MMG5_int gap; \
\
assert ( mesh->memCur < mesh->memMax ); \
\
gap = (MMG5_int)(floor(wantedGap * initSize)); \
if ( !gap ) gap = 1; \
\
if ( mesh->memMax < mesh->memCur + gap*sizeof(type) ) { \
gap = (MMG5_int)((mesh->memMax-mesh->memCur)/sizeof(type)); \
if(gap<1) { \
fprintf(stderr," ## Error:"); \
fprintf(stderr," unable to allocate %s.\n",message); \
fprintf(stderr," ## Check the mesh size or "); \
fprintf(stderr,"increase maximal authorized memory with the -m option.\n"); \
law; \
} \
} \
\
MMG5_ADD_MEM(mesh,gap*sizeof(type),message,law); \
MMG5_SAFE_RECALLOC((ptr),initSize+1,initSize+gap+1,type,message,law); \
initSize = initSize+gap; \
}while(0);
size_t memCur
Definition: libmmgtypes.h:607
size_t memMax
Definition: libmmgtypes.h:606

Reallocation of ptr of type type at size (initSize+wantedGap*initSize) if possible or at maximum available size if not. Execute the command law if reallocation failed. Memset to 0 for the new values of table.

Definition at line 404 of file mmgcommon_private.h.

◆ MMG5_TRIA_LMAX

#define MMG5_TRIA_LMAX   1024

Maximum array size when storing list of tria containing a vertex.

Definition at line 59 of file mmgcommon_private.h.

◆ MMG5_UNSET

#define MMG5_UNSET   -1

Definition at line 81 of file mmgcommon_private.h.

◆ MMG_FREAD

#define MMG_FREAD (   ptr,
  size,
  count,
  stream 
)
Value:
do \
{ \
\
if ( count != fread(ptr,size,count,stream) ) { \
fputs ( "Reading error", stderr ); \
return -1; \
} \
} while(0);

check the return value of fread

Definition at line 440 of file mmgcommon_private.h.

◆ MMG_FSCANF

#define MMG_FSCANF (   stream,
  format,
  ... 
)
Value:
do \
{ \
int io_count = fscanf(stream,format,__VA_ARGS__); \
int args_count = CV_VA_NUM_ARGS(__VA_ARGS__); \
if ( 0 > io_count ) { \
fprintf (stderr, "Reading error: fscanf counts %d args\n",io_count); \
return -1; \
} \
} while(0);
#define CV_VA_NUM_ARGS(...)

check the return value of fscanf

Remarks
don't work without any variadic arg
don't work on MSVC because variadic args are not expanded

Definition at line 472 of file mmgcommon_private.h.

Typedef Documentation

◆ MMG5_iNode

typedef struct MMG5_iNode_s MMG5_iNode

◆ MMG5_pBezier

Definition at line 615 of file mmgcommon_private.h.

Function Documentation

◆ MMG2D_quickarea()

double MMG2D_quickarea ( double  a[2],
double  b[2],
double  c[2] 
)
Parameters
apoint coordinates
bpoint coor
cpoint coor

Compute tria area.

Definition at line 971 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_2d3dUsage()

void MMG5_2d3dUsage ( void  )

Print help for common options between 2D and 3D.

Definition at line 291 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_2dSolTruncature_ani()

int MMG5_2dSolTruncature_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the solution structure.
Returns
0 if fail, 1 if succeed.

Compute hmin and hmax values from unflagged points (if not setted by the user), check the values and truncate the 2D metric.

Definition at line 365 of file scalem.c.

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

◆ MMG5_3dSolTruncature_ani()

int MMG5_3dSolTruncature_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the solution structure.
Returns
0 if fail, 1 if succeed.

Compute hmin and hmax values from unflagged points (if not setted by the user), check the values and truncate the 3D metric.

Definition at line 453 of file scalem.c.

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

◆ MMG5_Add_inode()

int MMG5_Add_inode ( MMG5_pMesh  mesh,
MMG5_iNode **  liLi,
int  val 
)
Parameters
meshpointer toward the mesh structure (for count of used memory).
liLipointer toward the address of the root of the linked list.
valvalue to add to the linked list.
Returns
1 if the node is inserted, 0 if the node is not inserted, -1 if fail.

Add a node with value val to a sorted linked list with unique entries.

Remarks
as the linked list had unique entries, we don't insert a node if it exists.

Definition at line 68 of file apptools.c.

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

◆ MMG5_advancedUsage()

void MMG5_advancedUsage ( void  )

Print help for advanced users of mmg.

Definition at line 302 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_bezierEdge()

void MMG5_bezierEdge ( MMG5_pMesh  ,
MMG5_int  ,
MMG5_int  ,
double *  ,
double *  ,
int8_t  ,
double *   
)

◆ MMG5_boulec()

int MMG5_boulec ( MMG5_pMesh  mesh,
MMG5_int *  adjt,
MMG5_int  start,
int  ip,
double *  tt 
)
Parameters
meshpointer toward the mesh structure.
adjtpointer toward the table of triangle adjacency.
startindex of triangle where we start to work.
ipindex of vertex where the tangent is computed.
ttpointer toward the computed tangent.
Returns
0 if fail, 1 otherwise.

Compute the tangent to the curve at point ip.

Definition at line 199 of file boulep.c.

Here is the caller graph for this function:

◆ MMG5_boulen()

int MMG5_boulen ( MMG5_pMesh  mesh,
MMG5_int *  adjt,
MMG5_int  start,
int  ip,
double *  nn 
)
Parameters
meshpointer toward the mesh structure.
adjtpointer toward the table of triangle adjacency.
startindex of triangle where we start to work.
iplocal index of vertex where the normal is computed.
nnpointer toward the computed tangent.
Returns
0 if fail, 1 otherwise.

Compute average normal of triangles sharing P without crossing ridge.

Definition at line 119 of file boulep.c.

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

◆ MMG5_boulep()

int MMG5_boulep ( MMG5_pMesh  mesh,
MMG5_int  start,
int  ip,
MMG5_int *  adja,
MMG5_int *  list,
MMG5_int *  tlist 
)
Parameters
meshpointer toward mesh structure.
starta triangle to which ip belongs.
iplocal point index
adjapointer toward the adjacency array.
listpointer toward the list of points connected to ip.
tlistpointer toward the list of triangles sharing ip.
Returns
-ilist if buffer overflow, ilist otherwise.

Return all vertices connected to ip (with list[0] = ip) and all triangles sharing ip.

Definition at line 52 of file boulep.c.

Here is the caller graph for this function:

◆ MMG5_bouler()

int MMG5_bouler ( MMG5_pMesh  mesh,
MMG5_int *  adjt,
MMG5_int  start,
int  ip,
MMG5_int *  list,
MMG5_int *  listref,
int *  ng,
int *  nr,
int  lmax 
)
Parameters
meshpointer toward the mesh structure.
adjtpointer toward the table of triangle adjacency.
startindex of triangle where we start to work.
ipindex of vertex on which we work.
listpointer toward the computed list of GEO vertices incident to ip.
listrefpointer toward the corresponding edge references
ngpointer toward the number of ridges.
nrpointer toward the number of reference edges.
lmaxmaxmum size for the ball of the point ip.
Returns
The number of edges incident to the vertex ip.

Store edges and count the number of ridges and reference edges incident to the vertex ip.

Definition at line 287 of file boulep.c.

Here is the caller graph for this function:

◆ MMG5_boulet()

int MMG5_boulet ( MMG5_pMesh  mesh,
MMG5_int  start,
int  ip,
MMG5_int *  list,
int8_t  s,
int8_t *  opn 
)
Parameters
meshpointer toward the mesh structure.
startindex of triangle to start.
ipindex of point for wich we compute the ball.
listpointer toward the computed ball of ip.
s1 if called from mmgs, 0 if called from mmg2d.
opn0 for a closed ball, 1 for an open ball.
Returns
the size of the computed ball or 0 if fail.

Find all triangles sharing ip, $list[0] =$ start do not stop when crossing ridge.

Definition at line 363 of file boulep.c.

Here is the caller graph for this function:

◆ MMG5_boundingBox()

int MMG5_boundingBox ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
1 if success, 0 if fail (computed bounding box too small).

Compute the mesh bounding box and fill the min, max and delta fields of the MMG5_info structure.

Definition at line 46 of file scalem.c.

Here is the caller graph for this function:

◆ MMG5_build3DMetric()

void MMG5_build3DMetric ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
MMG5_int  ip,
double  dbuf[6] 
)
Parameters
meshpointer toward the mesh structure
solpointer toward the sol structure.
indexof point in which we want to build the metric
dbufbuilded metric

Build the metric at point ip depending with its type (ridge/not ridge).

Definition at line 1528 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_buildridmet()

int MMG5_buildridmet ( MMG5_pMesh  ,
MMG5_pSol  ,
MMG5_int  ,
double  ,
double  ,
double  ,
double *  ,
double  [3][3] 
)

◆ MMG5_buildridmetfic()

int MMG5_buildridmetfic ( MMG5_pMesh  ,
double *  ,
double *  ,
double  ,
double  ,
double  ,
double *   
)

◆ MMG5_buildridmetnor()

int MMG5_buildridmetnor ( MMG5_pMesh  ,
MMG5_pSol  ,
MMG5_int  ,
double *  ,
double *  ,
double  [3][3] 
)

◆ MMG5_caltri33_ani()

double MMG5_caltri33_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  pt 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the meric structure.
ptpointer toward the triangle structure.
Returns
The computed quality.

Compute the quality of the surface triangle ptt with respect to an anisotropic metric and a classic storage of the ridges metrics.

Definition at line 47 of file quality.c.

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

◆ MMG5_caltri_ani()

double MMG5_caltri_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  ptt 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the meric structure.
pttpointer toward the triangle structure.
Returns
The computed quality.

Compute the quality of the surface triangle ptt with respect to an anisotropic metric.

Warning
The quality is computed as if the triangle is a "straight" triangle.

Definition at line 115 of file quality.c.

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

◆ MMG5_caltri_iso()

double MMG5_caltri_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  ptt 
)
inline
Parameters
meshpointer toward the mesh structure.
metpointer toward the meric structure.
pttpointer toward the triangle structure.
Returns
The computed quality.

Compute the quality of the surface triangle ptt with respect to an isotropic metric.

Definition at line 199 of file quality.c.

Here is the caller graph for this function:

◆ MMG5_check_hminhmax()

void MMG5_check_hminhmax ( MMG5_pMesh  mesh,
int8_t  sethmin,
int8_t  sethmax 
)
Parameters
meshpointer toward the mesh structure.
sethmin1 if hmin is setted by the user.
sethmax1 if hmax is setted by the user.

Check the compatibility between the automatically computed hmin/hmax values and the user settings.

Definition at line 89 of file scalem.c.

Here is the caller graph for this function:

◆ MMG5_check_readedMesh()

int MMG5_check_readedMesh ( MMG5_pMesh  mesh,
MMG5_int  nref 
)
Parameters
meshpointer toward an Mmg mesh
nrefpointer toward the number of negative refs (replaced by abolute values).
Returns
1 if success, 0 otherwise

Check the tetra orientation, print warning it negative refs have been detected, mark points as used and remove elt refs in iso mode.

Definition at line 526 of file inout.c.

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

◆ MMG5_check_setted_hminhmax()

int MMG5_check_setted_hminhmax ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
0 if fail, 1 if succeed

Check that hmin (resp. hmax) is not user setted if it is negative.

Definition at line 116 of file scalem.c.

Here is the caller graph for this function:

◆ MMG5_chkmaniball()

int MMG5_chkmaniball ( MMG5_pMesh  mesh,
MMG5_int  start,
int8_t  istart 
)
Parameters
meshpointer toward the mesh structure.
startindex of starting tria.
istartlocal index of point that we check (in tria start)
Returns
1 if the ball is manifold, 0 otherwise.

Check whether the ball of vertex i in tria start is manifold;

Warning
inxt[i] is an edge belonging to the implicit boundary.

Step 1: Travel while another part of the implicit boundary is not met

Normal or multi-material mode: check for change in triangle references

Input reference preservation mode (mmgs –keep-ref option): Check if we cross an isoref edge

Step 2: Check why the loop has ended

a./ Case where a boundary is hit: travel in the other sense from start, and make sure that a boundary is hit too

Ball is manifold if tested point is connected to two external edges

b./ General case: go on travelling until another implicit boundary is met

Ball is non-manifold if at least 3 boundary segments meeting at p

Definition at line 1303 of file mmg2.c.

Here is the caller graph for this function:

◆ MMG5_chkmanimesh()

int MMG5_chkmanimesh ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh.
Returns
1 if the mesh is manifold, 0 otherwise.

Check whether the resulting two subdomains occupying mesh are manifold.

First check: check whether one triangle in the mesh has 3 boundary faces

Second check: check whether the configuration is manifold in the ball of each point; each vertex on the implicit boundary is caught in such a way that i1 inxt[i1] is one edge of the implicit boundary

Definition at line 1423 of file mmg2.c.

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

◆ MMG5_chkMetricType()

int MMG5_chkMetricType ( MMG5_pMesh  mesh,
int *  type,
int *  entities,
FILE *  inm 
)
Parameters
meshpointer toward the mesh structure.
typetype of the metric
entitiesentities on which the metric applies (should be MMG5_Vertex)
inmmetric file
Returns
1 if success, -1 if fail

Check metric data:

  1. check that the metric applies on vertices;
  2. check that the type of the metric is compatible with the remeshing mode. If not, close the metric file (note that if type is an allocated array, you must unallocate it outside).

Definition at line 2661 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_compute_meanMetricAtMarkedPoints_ani()

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

Compute the mean metric at mesh points with a non-nul s field. At the beginning, for a given point ip, $ met->m[met->size * ip] $ contains the sum of n metrics and the s field of ip contains the number of metrics summed in the point. Set the flag of the processed points to 3.

Definition at line 2205 of file anisosiz.c.

Here is the caller graph for this function:

◆ MMG5_compute_meanMetricAtMarkedPoints_iso()

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

Compute the mean metric at mesh points with a non-nul s field. At the beginning, the metric of a given point contains the sum of n metrics and the s field of the point the number of metrics summed in the point. Set the flag of the processed points to 3.

Definition at line 167 of file isosiz.c.

Here is the caller graph for this function:

◆ MMG5_countLocalParamAtTri()

int MMG5_countLocalParamAtTri ( MMG5_pMesh  mesh,
MMG5_iNode **  bdryRefs 
)
Parameters
meshpointer toward the mesh structure.
bdryRefspointer toward the list of the boundary references.
Returns
npar, the number of local parameters at triangles if success, 0 otherwise.

Count the local default values at triangles and fill the list of the boundary references.

Count the number of different boundary references and list it

Definition at line 142 of file apptools.c.

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

◆ MMG5_crossprod3d()

void MMG5_crossprod3d ( double *  a,
double *  b,
double *  result 
)
Parameters
afirst array
bsecond array
resultcross product of the two arrays

Compute cross product of two double precision arrays in 3D.

Definition at line 300 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_defsiz_startingMessage()

int MMG5_defsiz_startingMessage ( MMG5_pMesh  mesh,
MMG5_pSol  met,
const char *  funcname 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
funcnamename of the calling function
Returns
1 if success, 0 if fail.

Print that we enter in the defsiz function in high verbosity level and check the hmax value.

Definition at line 77 of file isosiz.c.

Here is the caller graph for this function:

◆ MMG5_defUninitSize()

void MMG5_defUninitSize ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int8_t  ismet 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ismet1 if user provided metric.

Search for points with unintialized metric and define anisotropic size at this points.

Definition at line 228 of file anisosiz.c.

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

◆ MMG5_det3pt1vec()

double MMG5_det3pt1vec ( double  c0[3],
double  c1[3],
double  c2[3],
double  v[3] 
)
inline

Compute 3 * 3 determinant : det(c1-c0,c2-c0,v)

Definition at line 920 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_det4pt()

double MMG5_det4pt ( double  c0[3],
double  c1[3],
double  c2[3],
double  c3[3] 
)
inline

Compute 3 * 3 determinant : det(c1-c0,c2-c0,c3-c0)

Definition at line 932 of file tools.c.

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

◆ MMG5_devangle()

int MMG5_devangle ( double *  n1,
double *  n2,
double  crit 
)
Parameters
n1first normal
n2second normal
critridge threshold
Returns
1 if success, 0 if fail

Check if the angle between n1 and n2 is larger than the ridge criterion. If yes, return 1, 0 otherwise (ridge creation).

Definition at line 103 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_dikmov()

short MMG5_dikmov ( MMG5_pMesh  mesh,
MMG5_pSol  disp,
short *  lastt,
short  shortmax,
MMG5_int   chkmovmeshMMG5_pMesh, MMG5_pSol, short, MMG5_int * 
)

common functions for lagrangian meshing.

Author
Charles Dapogny (UPMC)
Cécile Dobrzynski (Bx INP/Inria/UBordeaux)
Pascal Frey (UPMC)
Algiane Froehly (Inria/UBordeaux)
Version
5
Parameters
meshpointer toward the mesh structure
disppointer toward the displacement field
lastt0 if a movement is possible, pointer toward the last tested fraction otherwise
shortmaxmaximal parameter t (MMG2D_SHORTMAX or MMG3D_SHORTMAX)
chkmovmeshfunction that has to be called to check motion validity
Returns
the largest fraction t allowing a valid motion

Generic function to compute the the largest fraction t that makes the motion along disp valid.

Definition at line 49 of file mmg3.c.

Here is the caller graph for this function:

◆ MMG5_displayLengthHisto()

void MMG5_displayLengthHisto ( MMG5_pMesh  mesh,
MMG5_int  ned,
double *  avlen,
MMG5_int  amin,
MMG5_int  bmin,
double  lmin,
MMG5_int  amax,
MMG5_int  bmax,
double  lmax,
int  nullEdge,
double *  bd,
MMG5_int *  hl,
int8_t  shift 
)
Parameters
meshpointer toward the mesh structure.
nededges number.
avlenpointer toward the average edges lengths.
aminindex of first extremity of the smallest edge.
bminindex of second extremity of the smallest edge.
lminsmallest edge length.
amaxindex of first extremity of the largest edge.
bmaxindex of second extremity of the largest edge.
lmaxlargest edge length.
nullEdgenumber of edges for which we are unable to compute the length
bdpointer toward the table of the quality span.
hlpointer toward the table that store the number of edges for eac
shiftvalue to shift the target lenght interval span of quality

Display histogram of edge length.

Definition at line 252 of file quality.c.

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

◆ MMG5_displayLengthHisto_internal()

void MMG5_displayLengthHisto_internal ( MMG5_int  ned,
MMG5_int  amin,
MMG5_int  bmin,
double  lmin,
MMG5_int  amax,
MMG5_int  bmax,
double  lmax,
MMG5_int  nullEdge,
double *  bd,
MMG5_int *  hl,
int8_t  shift,
int  imprim 
)
Parameters
nededges number.
aminindex of first extremity of the smallest edge.
bminindex of second extremity of the smallest edge.
lminsmallest edge length.
amaxindex of first extremity of the largest edge.
bmaxindex of second extremity of the largest edge.
lmaxlargest edge length.
nullEdgenumber of edges for which we are unable to compute the length
bdpointer toward the table of the quality span.
hlpointer toward the table that store the number of edges for eac
shiftvalue to shift the target lenght interval span of quality
imprimverbosity level

Display histogram of edge length without the histo header

Definition at line 294 of file quality.c.

Here is the caller graph for this function:

◆ MMG5_dotprod()

void MMG5_dotprod ( int8_t  dim,
double *  a,
double *  b,
double *  result 
)
Parameters
dimsize of the array
afirst array
bsecond array
resultscalar product of the two arrays

Compute scalar product of two double precision arrays.

Definition at line 265 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_eigenvmatnonsym2d()

int MMG5_eigenvmatnonsym2d ( MMG5_pMesh  mesh,
double  m[],
double  lambda[],
double  v[][2] 
)
Parameters
meshpointer toward the mesh structure.
mmatrix array.
lambdaeigenvalues array.
vdouble array of eigenvectors.
Returns
1 if success, 0 if failure.

Recompose a 2x2 non-symmetric matrix from its eigenvalue decomposition.

Warning
Eigenvectors in Mmg are stored as matrix rows (the first dimension of the double array spans the number of eigenvectors, the second dimension spans the number of entries of each eigenvector).

Definition at line 184 of file mettools.c.

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

◆ MMG5_eigenvmatnonsym3d()

int MMG5_eigenvmatnonsym3d ( MMG5_pMesh  mesh,
double  m[],
double  lambda[],
double  v[][3] 
)
Parameters
meshpointer toward the mesh structure.
mmatrix array.
lambdaeigenvalues array.
vdouble array of eigenvectors.
Returns
1 if success, 0 if failure.

Recompose a 3x3 non-symmetric matrix from its eigenvalue decomposition.

Warning
Eigenvectors in Mmg are stored as matrix rows (the first dimension of the double array spans the number of eigenvectors, the second dimension spans the number of entries of each eigenvector).

Definition at line 209 of file mettools.c.

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

◆ MMG5_eigenvmatsym2d()

int MMG5_eigenvmatsym2d ( MMG5_pMesh  mesh,
double  m[],
double  lambda[],
double  v[][2] 
)
Parameters
meshpointer toward the mesh structure.
mmatrix array.
lambdaeigenvalues array.
vdouble array of eigenvectors.
Returns
1 if success, 0 if failure.

Recompose a 2x2 symmetric matrix from its eigenvalue decomposition.

Warning
Eigenvectors in Mmg are stored as matrix rows (the first dimension of the double array spans the number of eigenvectors, the second dimension spans the number of entries of each eigenvector).

Definition at line 144 of file mettools.c.

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

◆ MMG5_eigenvmatsym3d()

int MMG5_eigenvmatsym3d ( MMG5_pMesh  mesh,
double  m[],
double  lambda[],
double  v[][3] 
)
Parameters
meshpointer toward the mesh structure.
mmatrix array.
lambdaeigenvalues array.
vdouble array of eigenvectors.
Returns
1 if success, 0 if failure.

Recompose a 3x3 symmetric matrix from its eigenvalue decomposition.

Warning
Eigenvectors in Mmg are stored as matrix rows (the first dimension of the double array spans the number of eigenvectors, the second dimension spans the number of entries of each eigenvector).

Definition at line 164 of file mettools.c.

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

◆ MMG5_elementWeight()

int MMG5_elementWeight ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  pt,
MMG5_pPoint  p0,
MMG5_Bezier pb,
double  r[3][3],
double  gv[2] 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ptpointer toward the tria on which we integrate.
p0pointer toward the point that we want to move.
pbbezier patch of the triangle.
rrotation matrix that sends the normal at point p0 to e_z.
gvcentre of mass that we want to update using the computed element weight.
Returns
0 if fail, 1 otherwise.

Compute integral of sqrt(T^J(xi) M(P(xi)) J(xi)) * P(xi) over the triangle.

Definition at line 53 of file anisomovpt.c.

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

◆ MMG5_excfun()

static void MMG5_excfun ( int  sigid)
inlinestatic
Parameters
sigidsignal number.

Signal handling: specify error messages depending from catched signal.

Definition at line 503 of file mmgcommon_private.h.

Here is the caller graph for this function:

◆ MMG5_fillDefmetregSys()

void MMG5_fillDefmetregSys ( MMG5_int  ,
MMG5_pPoint  ,
int  ,
MMG5_Bezier  ,
double  r[3][3],
double *  ,
double *  ,
double *  ,
double *   
)

◆ MMG5_Free_ilinkedList()

void MMG5_Free_ilinkedList ( MMG5_pMesh  mesh,
MMG5_iNode liLi 
)
Parameters
meshpointer toward the mesh structure (for count of used memory).
liLipointer toward the root of the linked list.

Free the memory used by the linked list whose root is liLi.

Definition at line 120 of file apptools.c.

Here is the caller graph for this function:

◆ MMG5_getStartRef()

int MMG5_getStartRef ( MMG5_pMesh  mesh,
MMG5_int  ref,
MMG5_int *  pref 
)
Parameters
meshpointer toward the mesh
reffinal reference for which we are searching the initial one
prefpointer to the reference of the parent material.
Returns
1 if found, 0 otherwise.

Retrieve the starting domain reference (parent material) associated to the split reference ref. Allow the call in non-multimaterial mode.

Definition at line 213 of file mmg2.c.

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

◆ MMG5_grad2metSurf()

MMG5_int MMG5_grad2metSurf ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  pt,
MMG5_int  np1,
MMG5_int  np2 
)
Parameters
meshpointer toward the mesh.
metpointer toward the metric structure.
ptpointer toward a triangle.
np1global index of the first extremity of the edge.
np2global index of the second extremity of the edge.
Returns
-1 if no gradation is needed, else index of graded point.

Enforces gradation of metric in one extremity of edge $f[ np1; np2]$f in tria pt with respect to the other, along the direction of the associated support curve first, then along the normal direction.

Warning
The gradation along the direction normal to the surface is made in an "isotropic way".
Remarks
ALGIANE: a mettre à plat : dans le cas d'une métrique très aniso avec la grande taille quasiment dans la direction de l'arête on se retrouve à modifier la grande taille uniquement (car proche de l'arête) sauf que cette modification n'a quasi pas d'influence sur le calcul de la longueur d'arête.

Definition at line 975 of file anisosiz.c.

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

◆ MMG5_grad2metSurfreq()

int MMG5_grad2metSurfreq ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  pt,
MMG5_int  npmaster,
MMG5_int  npslave 
)
Parameters
meshpointer toward the mesh.
metpointer toward the metric structure.
ptpointer toward the processed triangle.
npmasteredge extremity that cannot be modified
npslaveedge extremity to modify to respect the gradation.
Returns
0 if no gradation is needed, 1 otherwise.

Enforces gradation of metric of the extremity ±a npslave of edge $f[ npmaster; npslave]$f in tria pt with respect to the other, along the direction of the associated support curve first, then along the normal direction.

Warning
The gradation along the direction normal to the surface is made in an "isotropic way".

Definition at line 1951 of file anisosiz.c.

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

◆ MMG5_gradation_info()

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

Print gradation values (depending on the verbosity).

Definition at line 96 of file isosiz.c.

Here is the caller graph for this function:

◆ MMG5_gradEigenvreq()

void MMG5_gradEigenvreq ( double *  dm,
double *  dn,
double  difsiz,
int8_t  dir,
int8_t *  ier 
)
Parameters
dmeigenvalues of the first matrix (not modified)
dneigenvalues of the second matrix (modified)
difsizmaximal size gap authorized by the gradation.
dirdirection in which the sizes are graded.
ier2 if dn has been updated, 0 otherwise.

Gradation of size dn = 1/sqrt(eigenv of the tensor) for required points in the idir direction.

Definition at line 1883 of file anisosiz.c.

Here is the caller graph for this function:

◆ MMG5_gradsiz_ani()

MMG5_int MMG5_gradsiz_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int *  it 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
itnumber of performed iteration (to fill)
Returns
nup, the number of points updated.

Standard gradation procedure.

Mark the edges belonging to a required entity

Definition at line 2259 of file anisosiz.c.

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
metpointer toward the metric
Returns
0 if fail, 1 otherwise

Isotropic mesh gradation routine. The points belonging to a required edge are treated in gradsizreq_iso.

Definition at line 278 of file isosiz.c.

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

◆ MMG5_gradsizreq_ani()

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

Enforces mesh gradation (on required entities) by truncating metric field.

Mark the edges belonging to a required entity (already done if the classic gradation is enabled)

Definition at line 2322 of file anisosiz.c.

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

◆ MMG5_gradsizreq_iso()

int MMG5_gradsizreq_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh
metpointer toward the metric
Returns
the number of updated metrics.

Isotropic mesh gradation routine. The points belonging to a required entity are treated in gradsizreq_iso.

Mark the edges belonging to a required entity

Update the sizes and mark the treated points

Definition at line 372 of file isosiz.c.

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

◆ MMG5_hashEdge()

int MMG5_hashEdge ( MMG5_pMesh  mesh,
MMG5_Hash hash,
MMG5_int  a,
MMG5_int  b,
MMG5_int  k 
)
Parameters
meshpointer toward the mesh structure.
hashpointer toward the hash table of edges.
aindex of the first extremity of the edge.
bindex of the second extremity of the edge.
kindex of point along the edge.
Returns
1 if success, 0 if fail.

Add edge $[a;b]$ to the hash table.

Definition at line 345 of file hash.c.

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

◆ MMG5_hashEdgeTag()

int MMG5_hashEdgeTag ( MMG5_pMesh  mesh,
MMG5_Hash hash,
MMG5_int  a,
MMG5_int  b,
int16_t  tag 
)
Parameters
meshpointer toward the mesh structure.
hashpointer toward the hash table of edges.
aindex of the first extremity of the edge.
bindex of the second extremity of the edge.
tagedge tag
Returns
the edge tag if success, 0 if fail.

Add edge $[a;b]$ to the hash table if it doesn't exist and store the edge tag. If the edge exist, add the new tag to the already stored tags.

Definition at line 437 of file hash.c.

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

◆ MMG5_hashFace()

MMG5_int MMG5_hashFace ( MMG5_pMesh  mesh,
MMG5_Hash hash,
MMG5_int  ia,
MMG5_int  ib,
MMG5_int  ic,
MMG5_int  k 
)
Parameters
meshpointer toward the mesh.
hashpointer toward the hash table to fill.
iafirst vertex of face to hash.
ibsecond vertex of face to hash.
icthird vertex of face to hash.
kindex of face to hash.
Returns
0 if fail, -1 if the face is newly hashed, index of the first face hashed if another face with same vertices exist.

Definition at line 286 of file hash.c.

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

◆ MMG5_hashGet()

MMG5_int MMG5_hashGet ( MMG5_Hash hash,
MMG5_int  a,
MMG5_int  b 
)
Parameters
hashpointer toward the hash table of edges.
aindex of the first extremity of the edge.
bindex of the second extremity of the edge.
Returns
the index of point stored along $[a;b]$.

Find the index of point stored along $[a;b]$.

Definition at line 495 of file hash.c.

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

◆ MMG5_hashNew()

int MMG5_hashNew ( MMG5_pMesh  mesh,
MMG5_Hash hash,
MMG5_int  hsiz,
MMG5_int  hmax 
)
Parameters
meshpointer toward the mesh structure.
hashpointer toward the hash table of edges.
hsizinitial size of hash table.
hmaxmaximal size of hash table.
Returns
1 if success, 0 if fail.

Hash edges or faces.

Definition at line 526 of file hash.c.

Here is the caller graph for this function:

◆ MMG5_hashUpdate()

int MMG5_hashUpdate ( MMG5_Hash hash,
MMG5_int  a,
MMG5_int  b,
MMG5_int  k 
)
Parameters
meshpointer toward the mesh structure.
hashpointer toward the hash table of edges.
aindex of the first extremity of the edge.
bindex of the second extremity of the edge.
knew index of point along the edge.
Returns
1 if success, 0 if fail (edge is not found).

Update the index of the point stored along the edge $[a;b]$

Definition at line 400 of file hash.c.

Here is the caller graph for this function:

◆ MMG5_interp_iso()

int MMG5_interp_iso ( double *  ma,
double *  mb,
double *  mp,
double  t 
)
Parameters
mapointer on a metric
mbpointer on a metric
mppointer on the computed interpolated metric
tinterpolation parameter (comprise between 0 and 1)

Linear interpolation of isotropic sizemap along an edge

Definition at line 484 of file intmet.c.

Here is the caller graph for this function:

◆ MMG5_interpreg_ani()

int MMG5_interpreg_ani ( MMG5_pMesh  ,
MMG5_pSol  ,
MMG5_pTria  ,
int8_t  ,
double  ,
double *  mr 
)

◆ MMG5_intersecmet22()

int MMG5_intersecmet22 ( MMG5_pMesh  mesh,
double *  m,
double *  n,
double *  mr 
)
Parameters
meshpointer toward the mesh structure.
mpointer toward a $(2x2)$ metric.
npointer toward a $(2x2)$ metric.
mrcomputed $(2x2)$ metric.
Returns
0 if fail, 1 otherwise.

Compute the intersected (2 x 2) metric from metrics m and n : take simultaneous reduction, and proceed to truncation in sizes.

Definition at line 814 of file mettools.c.

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

◆ MMG5_intmetsavedir()

int MMG5_intmetsavedir ( MMG5_pMesh  mesh,
double *  m,
double *  n,
double *  mr 
)
Parameters
meshpointer toward the mesh structure.
mpointer toward the first metric to intersect.
npointer toward the second metric to intersect.
mrpointer toward the computed intersected metric.
Returns
1.

Compute the intersected (2 x 2) metric between metrics m and n, PRESERVING the directions of m. Result is stored in mr.

Definition at line 635 of file mettools.c.

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

◆ MMG5_intridmet()

int MMG5_intridmet ( MMG5_pMesh  ,
MMG5_pSol  ,
MMG5_int  ,
MMG5_int  ,
double  ,
double *  ,
double *   
)

◆ MMG5_invmat()

int MMG5_invmat ( double *  m,
double *  mi 
)
Parameters
mpointer toward a 3x3 symetric matrix
mipointer toward the computed 3x3 matrix.

Invert m (3x3 symetric matrix) and store the result on mi

Definition at line 562 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_invmat22()

int MMG5_invmat22 ( double  m[2][2],
double  mi[2][2] 
)
Parameters
minitial matrix.
miinverted matrix.

Invert 2x2 non-symmetric matrix stored in 2 dimensions

Definition at line 745 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_invmat33()

int MMG5_invmat33 ( double  m[3][3],
double  mi[3][3] 
)
Parameters
minitial matrix.
miinverted matrix.

Invert 3x3 non-symmetric matrix stored in 2 dimensions

Definition at line 653 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_invmatg()

int MMG5_invmatg ( double  m[9],
double  mi[9] 
)
Parameters
minitial matrix.
miinverted matrix.

Invert 3x3 non-symmetric matrix.

Definition at line 613 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_isLevelSet()

int MMG5_isLevelSet ( MMG5_pMesh  mesh,
MMG5_int  ref0,
MMG5_int  ref1 
)
Parameters
meshpointer toward the mesh structure.
ref0reference of the first tetrahedron sharing the face.
ref1reference of the second tetrahedron sharing the face..
Returns
1 if face is on the discrete level set, 0 if not.

Identify whether a face is on the discrete level set or not.

Definition at line 471 of file mmg2.c.

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

◆ MMG5_ismaniball()

int MMG5_ismaniball ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
MMG5_int  start,
int8_t  istart 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the level-set values.
startindex of the starting tria
istartlocal index (inside the tria start) of the vertex that we check.
Returns
1 if success, 0 if fail

Check whether snapping the value of vertex istart of start to 0 exactly leads to a non manifold situation.

Warning
: we assume that the triangle start has vertex istart with value 0 and the other two with changing values.

Definition at line 622 of file mmg2.c.

Here is the caller graph for this function:

◆ MMG5_isNotSplit()

int MMG5_isNotSplit ( MMG5_pMesh  mesh,
MMG5_int  ref 
)
Parameters
meshpointer toward the mesh structure.
refinitial reference.
Returns
1 if entity cannot be split, 0 if can be split.

Identify whether an entity with reference ref should not be split.

Definition at line 447 of file mmg2.c.

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

◆ MMG5_isSplit()

int MMG5_isSplit ( MMG5_pMesh  mesh,
MMG5_int  ref,
MMG5_int *  refint,
MMG5_int *  refext 
)
Parameters
meshpointer toward the mesh structure.
refinitial reference.
refintinternal reference after ls discretization.
refextexternal reference after ls discretization.
Returns
1 if entity can be splitted, 0 if cannot be splitted.

Identify whether an entity with reference ref should be split, and the labels of the resulting entities.

Definition at line 411 of file mmg2.c.

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

◆ MMG5_keep_subdomainElts()

void MMG5_keep_subdomainElts ( MMG5_pMesh  mesh,
int  nsd,
int(*)(MMG5_pMesh, MMG5_int)  delElt 
)
Parameters
meshpointer toward the mesh structure.
nsdsubdomain index.
delEltfunction to call to delete elt.

Remove triangles that do not belong to subdomain of index nsd

Definition at line 1072 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_lagUsage()

void MMG5_lagUsage ( void  )

Print help for lagrangian motion option.

Definition at line 275 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_loadMshMesh_part1()

int MMG5_loadMshMesh_part1 ( MMG5_pMesh  mesh,
const char *  filename,
FILE **  inm,
long *  posNodes,
long *  posElts,
long **  posNodeData,
int *  bin,
int *  iswp,
MMG5_int *  nelts,
int *  nsols 
)
Parameters
meshpointer toward the mesh
filenamepointer toward the name of file
inmpointer toward the file pointer
posNodespointer toward the position of nodes data in file
posEltspointer toward the position of elts data in file
posNodeDatapointer toward the list of the positions of data in file
bin1 if binary format
neltsnumber of elements in file
nsolnumber of data in file
Returns
1 if success, 0 if file is not found, -1 if fail.

Begin to read mesh at MSH file format. Read the mesh size informations.

Remarks
For now only intput files containing int32 integers are supported.

Definition at line 278 of file inout.c.

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

◆ MMG5_loadMshMesh_part2()

int MMG5_loadMshMesh_part2 ( MMG5_pMesh  mesh,
MMG5_pSol sol,
FILE **  inm,
const long  posNodes,
const long  posElts,
const long *  posNodeData,
const int  bin,
const int  iswp,
const MMG5_int  nelts,
const int  nsols 
)
Parameters
meshpointer toward the mesh
solpointer toward the solutions array
inmpointer toward the file pointer
posNodesposition of nodes data in file
posEltsposition of elts data in file
posNodeDataposition of solution data in file
bin1 if binary format
neltsnumber of elements in file
nsolsnumber of silutions in file
Returns
1 if success, -1 if fail.

End to read mesh and solution array at MSH file format after the mesh/solution array alloc.

Remarks
For now only intput files containing int32 integers are supported.

Second step: read the nodes and elements

Read the solution at nodes

Definition at line 668 of file inout.c.

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

◆ MMG5_loadSolHeader()

int MMG5_loadSolHeader ( const char *  filename,
int  meshDim,
FILE **  inm,
int *  ver,
int *  bin,
int *  iswp,
MMG5_int *  np,
int *  dim,
int *  nsols,
int **  type,
long *  posnp,
int  imprim 
)
Parameters
filenamename of file.
meshDimmesh dimenson.
inmallocatable pointer toward the FILE structure
verfile version (1=simple precision, 2=double)
bin1 if the file is a binary
iswp1 or 0 depending on the endianness (binary only)
npnumber of solutions of each type
dimsolution dimension
nsolsnumber of solutions of different types in the file
typetype of solutions
posnppointer toward the position of the point list in the file
imprimverbosity
Returns
-1 data invalid or we fail, 0 no file, 1 ok.

Open the "filename" solution file and read the file header.

Definition at line 2080 of file inout.c.

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

◆ MMG5_loadVtuMesh()

int MMG5_loadVtuMesh ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
const char *  filename 
)

◆ MMG5_mark_pointsOnReqEdge_fromTria()

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

Set the s field of the points that belongs to a required edge to 1, set it to 0 otherwise.

Definition at line 243 of file isosiz.c.

Here is the caller graph for this function:

◆ MMG5_mark_usedVertices()

void MMG5_mark_usedVertices ( MMG5_pMesh  mesh,
void(*)(MMG5_pMesh, MMG5_int)  delPt 
)
Parameters
meshpointer toward the mesh structure.
delPtfunction to call to delete point.

Mark the mesh vertices that belong to triangles or quadrangles as used (for Mmgs or Mmg2d).

Definition at line 1018 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_mark_verticesAsUnused()

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

Mark all mesh vertices as unused.

Definition at line 994 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_memOption_memSet()

void MMG5_memOption_memSet ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure

Set the memMax value to its "true" value if memory asked by user. Here the MMG5_MEMPERCENT coef is already applied on memMax.

Definition at line 891 of file tools.c.

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

◆ MMG5_memSize()

size_t MMG5_memSize ( void  )
Returns
the available memory size of the computer.

Compute the available memory size of the computer.

Definition at line 852 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_minQualCheck()

int MMG5_minQualCheck ( MMG5_int  iel,
double  minqual,
double  alpha 
)
Parameters
ielindex of the worst tetra of the mesh
minqualquality of the worst tetra of the mesh (will be normalized by alpha)
alphanormalisation parameter for the quality
Returns
1 if success, 0 if fail (the quality is lower than MMG5_NULKAL).

Print warning or error messages depending on the quality of the worst tetra of the mesh.

Definition at line 343 of file quality.c.

Here is the caller graph for this function:

◆ MMG5_mmgDefaultValues()

void MMG5_mmgDefaultValues ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh structure.
Returns
0 if fail, 1 if success.

Print the default parameters values.

Definition at line 70 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_mmgHashTria()

int MMG5_mmgHashTria ( MMG5_pMesh  mesh,
MMG5_int *  adjt,
MMG5_Hash hash,
int  chkISO 
)
Parameters
meshpointer toward the mesh structure.
adjtpointer toward the adjacency table of the surfacic mesh.
hashpointer toward the edge hash table.
chkISOflag to say if we check ISO references (so if we come from mmg3d).
Returns
1 if success, 0 otherwise.

Create surface adjacency

Remarks
the ph->s field computation is useless in mmgs.
: as all triangles are mesh boundaries, we do not need to mark their adges as MG_BDY so the MG_BDY tag may be used inside geometrical triangles (external non-parallel, or internal parallel) to tag edges on the intersection with purely parallel (non-geometrical) triangles. The MG_PARBDYBDY tag is also added, as it does not have a supporting triangle to inherit this tag from.

Definition at line 58 of file hash.c.

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

◆ MMG5_mmgInit_parameters()

void MMG5_mmgInit_parameters ( MMG5_pMesh  mesh)

◆ MMG5_mmgIntextmet()

int MMG5_mmgIntextmet ( MMG5_pMesh  ,
MMG5_pSol  ,
MMG5_int  ,
double *  ,
double *   
)

◆ MMG5_mmgIntmet33_ani()

int MMG5_mmgIntmet33_ani ( double *  m,
double *  n,
double *  mr,
double  s 
)
Parameters
minput metric.
ninput metric.
mrcomputed output metric.
sparameter coordinate for the interpolation of metrics m and n.
Returns
0 if fail, 1 otherwise.

Compute the interpolated $(3 x 3)$ metric from metrics m and n, at parameter s : $ mr = (1-s)*m +s*n $, both metrics being expressed in the simultaneous reduction basis: linear interpolation of sizes.

Definition at line 50 of file intmet.c.

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

◆ MMG5_mmgUsage()

void MMG5_mmgUsage ( char *  prog)
Parameters
*progpointer toward the program name.

Print help for common options of the 3 codes (first section).

Definition at line 210 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_mn()

void MMG5_mn ( double  m[6],
double  n[6],
double  mn[9] 
)
Parameters
msymmetric matrix
nsymmetric matrix
mnresult

Compute product m*n (mn stored by rows for consistency with MMG5_eigenv3d).

Definition at line 337 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_MultiMat_init()

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

Initialize handling of multimaterial mode.

An indexed table of materials has been provided by the MMG5_Mat datatype in the form:

index | dospl | ref | rin | rex

0 | dospl_0 | ref_0 | rin_0 | rex_0 ... | ... | ... | ... | k | dospl_k | ref_k | rin_k | rex_k ... | ... | ... | ... | n-1 | dospl_{n-1} | ref_{n-1} | rin_{n-1} | rex_{n-1}

where dospl is the split/preserve attribute of the material, and rin,rex are its child materials (if dospl). Viceversa, ref is the parent material for rin,rex.

Here a lookup table for material references is built through trivial hashing for all references (both parent and child materials) with the key:

key = ref - ref_min, ref_min = min_{k = 0,n-1} (ref_k, rin_k, rex_k)

For all references, it is important to store

  • the index k of the row in the original table,
  • the characteristic attribute (parent-split, parent-preserve, child-interior, child-exterior) of the material.

Since dospl = 0 or 1, MG_PLUS = 2, and MG_MINUS = 3

we can store the attribute as dospl (for the parent material) or MG_MINUS/ MG_PLUS (for the child materials), with its value ranging between 0 and 3.

A convenient entry to store both the index and the attribute in the lookup table is thus:

entry = 4*(index+1) + attribute

leading to a lookup table in the form (the key ordering is only an example):

key | entry

... | ... ref_k | 4*(k+1)+dospl_k ... | ... rin_k | 4*(k+1)+MG_MINUS ... | ... rex_k | 4*(k+1)+MG_PLUS ... | ...

where the index and the attribute of the material can be retrieved as

index = entry / 4 -1 attribute = entry % 4

What if two materials have the same reference?

  • child references should be distinct (the entry in the lookup table would be overwritten),
  • a parent material can have itself as child (a positive attribute would say it should be split, the attribute value would say if it is interior or exterior). Thus, each of the maps parent->child_in and parent->child_ext is injective, but a fixed point is allowed.

Why child materials should be different?

  • because on failure it is necessary to recover the parent references and apply them on the tetrahedra to restore the starting materials distribution.

Definition at line 326 of file mmg2.c.

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

◆ MMG5_nonorsurf()

double MMG5_nonorsurf ( MMG5_pMesh  mesh,
MMG5_pTria  pt 
)
inline
Parameters
meshpointer toward the mesh stucture.
pttriangle for which we compute the surface.
Returns
the computed surface

Compute non-oriented surface area of a triangle.

Definition at line 160 of file tools.c.

Here is the call graph for this function:

◆ MMG5_nonUnitNorPts()

int MMG5_nonUnitNorPts ( MMG5_pMesh  mesh,
MMG5_int  ip1,
MMG5_int  ip2,
MMG5_int  ip3,
double *  n 
)
inline
Parameters
meshpointer toward the mesh stucture.
ip1first point of face.
ip2second point of face.
ip3third point of face.
npointer to store the computed normal.
Returns
1

Compute non-normalized face normal given three points on the surface.

Definition at line 127 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_norpts()

int MMG5_norpts ( MMG5_pMesh  mesh,
MMG5_int  ip1,
MMG5_int  ip2,
MMG5_int  ip3,
double *  n 
)
inline
Parameters
meshpointer toward the mesh stucture.
ip1first point of face.
ip2second point of face.
ip3third point of face.
npointer to store the computed normal.
Returns
1 if success, 0 otherwise.

Compute normalized face normal given three points on the surface.

Definition at line 183 of file tools.c.

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

◆ MMG5_nortri()

int MMG5_nortri ( MMG5_pMesh  mesh,
MMG5_pTria  pt,
double *  n 
)
inline
Parameters
meshpointer toward the mesh stucture.
ptpointer toward the triangle structure.
npointer to store the computed normal.
Returns
1

Compute triangle normal.

Definition at line 209 of file tools.c.

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

◆ MMG5_nperm()

void MMG5_nperm ( int8_t  n,
int8_t  shift,
int8_t  stride,
double *  val,
double *  oldval,
int8_t *  perm 
)
inline
Parameters
narray size
shiftshift to apply when taking array value
stridestride to apply when taking array value
valarray of double precision floating points
oldvalarray to store input values
permpermutation array

Naively permute a small array. Use shift and stride to eventually permute matrix columns.

Remarks
to use only on very small arrays such as metric tensors (size lower than int8_max).

Definition at line 80 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_nsort()

void MMG5_nsort ( int8_t  n,
double *  val,
int8_t *  perm 
)
inline
Parameters
narray size
valarray of double precision floating points
permpermutation array

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

Remarks
to use only on very small arrays such as metric tensors (size lower than int8_max).

Definition at line 49 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_orvol()

double MMG5_orvol ( MMG5_pPoint  point,
MMG5_int *  v 
)
inline
Parameters
pointPointer toward the points array
vpointer toward the point indices
Returns
the oriented volume of tetra

Compute oriented volume of a tetrahedron (x6)

Definition at line 951 of file tools.c.

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

◆ MMG5_paramUsage1()

void MMG5_paramUsage1 ( void  )

Print help for common parameters options of the 3 codes (first section).

Definition at line 239 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_paramUsage2()

void MMG5_paramUsage2 ( void  )

Print help for common options of the 3 codes (second section).

Definition at line 258 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_paratmet()

int MMG5_paratmet ( double  c0[3],
double  n0[3],
double  m[6],
double  c1[3],
double  n1[3],
double  mt[6] 
)
Parameters
c0table of the coordinates of the starting point.
n0normal at the starting point.
mmetric to be transported.
c1table of the coordinates of the ending point.
n1normal at the ending point.
mtcomputed metric.
Returns
0 if fail, 1 otherwise.

Parallel transport of a metric tensor field, attached to point c0, with normal n0, to point c1, with normal n1.

Definition at line 1390 of file mettools.c.

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

◆ MMG5_printMetStats()

void MMG5_printMetStats ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.

print metric statistics

Definition at line 2705 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_printSolStats()

void MMG5_printSolStats ( MMG5_pMesh  mesh,
MMG5_pSol sol 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the solutions array.

print solutions statistics

Definition at line 2723 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_printTria()

void MMG5_printTria ( MMG5_pMesh  mesh,
char *  fileName 
)
Parameters
meshpointer toward the mesh structure.
fileNamepointer toward the file name.

Debug function (not use in clean code): write mesh->tria structure in file.

Definition at line 825 of file tools.c.

◆ MMG5_readDoubleSol3D()

int MMG5_readDoubleSol3D ( MMG5_pSol  sol,
FILE *  inm,
int  bin,
int  iswp,
MMG5_int  pos 
)
Parameters
solpointer toward an allocatable sol structure.
inmpointer toward the solution file
bin1 if binary file
iswpEndianess
indexof the readed solution
Returns
1 if success, -1 if fail

Read the solution value for vertex of index pos in double precision.

Definition at line 2269 of file inout.c.

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

◆ MMG5_readFloatSol3D()

int MMG5_readFloatSol3D ( MMG5_pSol  sol,
FILE *  inm,
int  bin,
int  iswp,
int  pos 
)
Parameters
solpointer toward an allocatable sol structure.
inmpointer toward the solution file
bin1 if binary file
iswpEndianess
indexof the readed solution
Returns
1 if success, -1 if fail

Read the solution value for vertex of index pos in floating precision.

Definition at line 2220 of file inout.c.

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

◆ MMG5_regnor()

int MMG5_regnor ( MMG5_pMesh  mesh)
Parameters
meshpointer toward a MMG5 mesh structure.
Returns
0 if fail, 1 otherwise.

Regularization procedure for derivatives, dual Laplacian

Definition at line 46 of file analys.c.

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

◆ MMG5_reset_metricAtReqEdges_surf()

int MMG5_reset_metricAtReqEdges_surf ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int8_t  ismet 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ismet1 if user provided metric
Returns
1 if success, 0 if fail.

For a triangle mesh, process the triangles and set to 0 the metrics at points that are at the extremities of a required edge.

Definition at line 204 of file isosiz.c.

Here is the caller graph for this function:

◆ MMG5_resetRef_ls()

int MMG5_resetRef_ls ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh

Reset mesh->info.isoref vertex and edge references to 0.

Definition at line 1188 of file mmg2.c.

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

◆ MMG5_resetRef_lssurf()

int MMG5_resetRef_lssurf ( MMG5_pMesh  mesh)
Parameters
meshpointer toward the mesh

Reset mesh->info.isoref vertex references to 0.

Definition at line 129 of file mmg2s.c.

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

◆ MMG5_ridSizeInNormalDir()

double MMG5_ridSizeInNormalDir ( MMG5_pMesh  ,
int  ,
double *  ,
MMG5_pBezier  ,
double  ,
double   
)

◆ MMG5_ridSizeInTangentDir()

double MMG5_ridSizeInTangentDir ( MMG5_pMesh  mesh,
MMG5_pPoint  p0,
MMG5_int  idp,
MMG5_int *  iprid,
double  isqhmin,
double  isqhmax 
)
Parameters
meshpointer toward the mesh structure.
p0pointer toward the point at which we define the metric.
idpglobal index of the point at which we define the metric.
ipridpointer toward the two extremities of the ridge.
isqhminminimum edge size.
isqhmaxmaximum edge size.
Returns
the computed ridge size in the tangent direction.

Compute the specific size that we want to apply to a ridge in the direction of the tangent of the ridge. This wanted size is computed as using the majoration of the haudorff distance between a triangle and its ideal curve approximation by the Hessian of the signed distance function to the ideal surface.

\[
d^H(\partial \Omega,S_T) \leq \displaystyle\frac{1}{2}\left(
\frac{d-1}{d} \right)^2 \max\limits_{T\in S_T} \max\limits_{x \in T
} \max\limits_{y,z \in T} \langle \left| H(d_\omega)(x) \right|
yz,yz\rangle.
\]

where $ d^H(\partial \Omega,S_T) $ is the distance between the triangle $S_T$ and the ideal boundary (reconstructed using cubic Bezier patches) $ \partial \Omega $, $ d$ is the mesh dimension and $
  H(d_\Omega) $ is the hessian matrix of the signed distance function to $
  \Omega $.

For all $ x \in \partial \Omega $, $H(d_\omega)(x)$ is the second fundamental form whose eigenvalues are the principal curvatures ( $\kappa_1
  $ and $ \kappa_2 $) of $ \partial \Omega $ at $ x $ so the previous formula can be rewritten as:

\[
d^H(\partial \Omega,S_T) \leq \displaystyle\frac{1}{2}\left(
\frac{d-1}{d} \right)^2 \max( \left|\kappa_1\right|,\left|\kappa_2\right|).
\]

As we want to respect the imposed the hausd threshold and we are interessed to get the wanted size along the tangent direction at the ridge point only, finally, m is computed as

\[
m = \kappa \frac{1}{2\textrm{hausd}}
\left(\frac{d-1}{d} \right)^2.
\]

See Theorem 1 of [3].

Definition at line 764 of file anisosiz.c.

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

◆ MMG5_rmc()

int MMG5_rmc ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh
solpointer toward the level-set
Returns
1 if success, 0 otherwise

Removal of small parasitic components (bubbles of material, etc) with volume less than mesh->info.rmc * volume of the mesh.

Definition at line 911 of file mmg2.c.

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

◆ MMG5_rmtr()

int MMG5_rmtr ( double  r[3][3],
double  m[6],
double  mr[6] 
)
inline
Parameters
r3x3 matrix
msymetric matrix
mrresult
Returns
1

Compute product R*M*tR when M is symmetric

Definition at line 405 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_rotmatrix()

int MMG5_rotmatrix ( double  n[3],
double  r[3][3] 
)
inline
Parameters
npointer toward the vector that we want to send on the third vector of canonical basis.
rcomputed rotation matrix.

Compute rotation matrix that sends vector n to the third vector of canonical basis.

Definition at line 467 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_saveDisp()

int MMG5_saveDisp ( MMG5_pMesh  mesh,
MMG5_pSol  disp 
)
Parameters
meshpointer toward the mesh structure
disppointer toward the displacement field
Returns
1 if success, 0 if fail.

For debugging purposes: save displacement field.

Definition at line 112 of file mmg3.c.

◆ MMG5_saveMshMesh()

int MMG5_saveMshMesh ( MMG5_pMesh  mesh,
MMG5_pSol sol,
const char *  filename,
int  metricData 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward an array of solutions.
filenamename of file.
metricData1 if the data saved is a metric (if only 1 data)
Returns
0 if failed, 1 otherwise.

Write mesh and a list of solutions at MSH file format (.msh extension). Write binary file for .mshb extension.and ASCII for .msh one.

First step: Count the number of elements of each type

Second step: save the elements at following format: "idx type tagNumber tag0 tag1... v0_elt v1_elt..."

Write solution

Save the solution at following format: "idx sol"

Definition at line 1587 of file inout.c.

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

◆ MMG5_saveSolAtTetrahedraHeader()

int MMG5_saveSolAtTetrahedraHeader ( MMG5_pMesh  mesh,
FILE *  inm,
int  ver,
int  bin,
MMG5_int *  bpos,
int  nsols,
int  nsolsAtTetra,
int *  entities,
int *  type,
int *  size 
)
Parameters
meshpointer toward the mesh structure.
inmpointer toward the opened file unit.
verfile version (1=simple precision, 2=double).
bin1 if the file is a binary.
bposcumulative field position for binary Medit format.
nsolsnumber of solutions of different types in the file.
nsolsAtTetranumber of solutions at tetra in the file.
entitieskind of entity on which the solution applies.
typetype of solutions.
sizesize of solutions.

Save the number and type of solutions at Tetrahedron (not used by Mmg).

Definition at line 2593 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_saveSolAtTrianglesHeader()

int MMG5_saveSolAtTrianglesHeader ( MMG5_pMesh  mesh,
FILE *  inm,
int  ver,
int  bin,
MMG5_int *  bpos,
int  nsols,
int  nsolsAtTriangles,
int *  entities,
int *  type,
int *  size 
)
Parameters
meshpointer toward the mesh structure.
inmpointer toward the opened file unit.
verfile version (1=simple precision, 2=double).
bin1 if the file is a binary.
bposcumulative field position for binary Medit format.
nsolsnumber of solutions of different types in the file.
nsolsAtTrianglesnumber of solutions at triangles in the file.
entitieskind of entity on which the solution applies.
typetype of solutions.
sizesize of solutions.

Save the number and type of solutions at Triangles (not used by Mmg).

Definition at line 2524 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_saveSolHeader()

int MMG5_saveSolHeader ( MMG5_pMesh  mesh,
const char *  filename,
FILE **  inm,
int  ver,
int *  bin,
MMG5_int *  bpos,
MMG5_int  np,
int  dim,
int  nsols,
int *  entities,
int *  type,
int *  size 
)
Parameters
meshpointer toward the mesh structure.
filenamename of file.
inmallocatable pointer toward the FILE structure.
verfile version (1=simple precision, 2=double).
bin1 if the file is a binary.
bposcumulative field position for binary Medit format.
npnumber of solutions of each type.
dimsolution dimension.
nsolsnumber of solutions of different types in the file.
entitieskind of entity on which the solution applies (Vertex or Tetra)
typetype of solutions.
sizesize of solutions.
Returns
0 if unable to open the file, 1 if success.

Open the "filename" solution file and save the file header and the number and type of solutions at vertices (Native solutions for Mmg).

Definition at line 2387 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_scale_meshAndSol()

int MMG5_scale_meshAndSol ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  sol,
double *  dd 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward a metric
solpointer toward a solution structure (level-set or displacement).
ddpointer toward the scaling value (to fill)
Returns
1 if success, 0 if fail.

Scale the mesh and the size informations between 0 and 1. Compute a default value for the hmin/hmax parameters if needed.

Definition at line 561 of file scalem.c.

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

◆ MMG5_scale_scalarMetric()

int MMG5_scale_scalarMetric ( MMG5_pMesh  mesh,
MMG5_pSol  met,
double  dd 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ddscaling value.
Returns
1 if success, 0 if fail.

Scale and truncate by hmin/hmax the scalar metric stored in met. If hmin/hmax are not provided by the user, it is automatically computed from the metric.

Definition at line 206 of file scalem.c.

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

◆ MMG5_scotchCall()

int MMG5_scotchCall ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  fields,
MMG5_int *  permNodGlob 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the solution structure.
fieldspointer toward an array of solution fields (non mandatory)
permNodGlobstore the global permutation of nodes (if provided).
Returns
0 if MMG5_renumbering fail (non conformal mesh), 1 otherwise (renumerotation success of renumerotation fail but the mesh is still conformal).

Call scotch renumbering.

Definition at line 231 of file librnbg.c.

Here is the caller graph for this function:

◆ MMG5_Set_commonFunc()

void MMG5_Set_commonFunc ( void  )

◆ MMG5_setref_ls()

int MMG5_setref_ls ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the level-set values.
Returns
1.

Set references to tris according to the sign of the level set function.

Definition at line 1224 of file mmg2.c.

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

◆ MMG5_setref_lssurf()

int MMG5_setref_lssurf ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)

Definition at line 160 of file mmg2s.c.

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

◆ MMG5_simred2d()

int MMG5_simred2d ( MMG5_pMesh  mesh,
double *  m,
double *  n,
double  dm[2],
double  dn[2],
double  vp[2][2] 
)
Parameters
meshpointer toward the mesh
mfirst matrix
nsecond matrix
dmeigenvalues of m in the coreduction basis (to fill)
dneigenvalues of n in the coreduction basis (to fill)
vpcoreduction basis (to fill)
Returns
0 if fail 1 otherwise.

Perform simultaneous reduction of matrices m and n.

Definition at line 1326 of file anisosiz.c.

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

◆ MMG5_simred3d()

int MMG5_simred3d ( MMG5_pMesh  mesh,
double *  m,
double *  n,
double  dm[3],
double  dn[3],
double  vp[3][3] 
)
Parameters
meshpointer toward the mesh
mfirst matrix
nsecond matrix
dmeigenvalues of m in the coreduction basis (to fill)
dneigenvalues of n in the coreduction basis (to fill)
vpcoreduction basis (to fill)
Returns
0 if fail 1 otherwise.

Perform simultaneous reduction of matrices m and n.

Definition at line 1416 of file anisosiz.c.

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

◆ MMG5_snpval_ls()

int MMG5_snpval_ls ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the level-set function.
Returns
1 if success, 0 if fail.

Snap values of the level set function very close to 0 to exactly 0, and prevent nonmanifold patterns from being generated.

Definition at line 502 of file mmg2.c.

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

◆ MMG5_snpval_lssurf()

int MMG5_snpval_lssurf ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
solpointer toward the level-set
Returns
1 if success, 0 if fail

Snap values of sol very close to 0 to 0 exactly in the case of surface ls splitting (to avoid very small triangles in cutting)

Definition at line 45 of file mmg2s.c.

Here is the caller graph for this function:

◆ MMG5_solTruncature_iso()

int MMG5_solTruncature_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the solution structure.
Returns
0 if fail, 1 if succeed.

Compute hmin and hmax values from unflagged points (if not setted by the user), check the values and truncate the metric.

Definition at line 299 of file scalem.c.

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

◆ MMG5_solveDefmetrefSys()

int MMG5_solveDefmetrefSys ( MMG5_pMesh  ,
MMG5_pPoint  ,
MMG5_int *  ,
double  r[3][3],
double *  ,
double *  ,
double *  ,
double *  ,
double  ,
double  ,
double   
)

◆ MMG5_solveDefmetregSys()

int MMG5_solveDefmetregSys ( MMG5_pMesh  ,
double  r[3][3],
double *  ,
double *  ,
double *  ,
double *  ,
double  ,
double  ,
double   
)

◆ MMG5_sum_reqEdgeLengthsAtPoint()

int MMG5_sum_reqEdgeLengthsAtPoint ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_int  ip0,
MMG5_int  ip1 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ip0index of the first edge extremity
ip1index of the second edge extremity
Returns
1 if success, 0 if fail.

Compute the euclidean length of the edge ip0 ip1, add this length to the metric of the edge extremities and increment the count of times we have processed this extremities.

Definition at line 129 of file isosiz.c.

Here is the caller graph for this function:

◆ MMG5_surftri33_ani()

double MMG5_surftri33_ani ( MMG5_pMesh  ,
MMG5_pTria  ,
double *  ,
double *  ,
double *   
)

◆ MMG5_surftri_ani()

double MMG5_surftri_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  ptt 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
pttpointer toward the triangle structure.
Returns
The double of the triangle area.

Compute the double of the area of the surface triangle ptt with respect to the anisotropic metric met (for special storage of ridges metrics).

Definition at line 124 of file anisosiz.c.

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

◆ MMG5_surftri_iso()

double MMG5_surftri_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  ptt 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the meric structure.
pttpointer toward the triangle structure.
Returns
The computed area.

Compute the area of the surface triangle ptt with respect to the isotropic metric met.

Definition at line 42 of file isosiz.c.

◆ MMG5_swapbin()

int MMG5_swapbin ( int  sbin)

swap bytes if needed (conversion from big/little endian toward little/big endian)

Definition at line 42 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_swapbin_int()

MMG5_int MMG5_swapbin_int ( MMG5_int  sbin)

swap bytes if needed (conversion from big/little endian toward little/big endian)

Definition at line 61 of file inout.c.

◆ MMG5_swapd()

double MMG5_swapd ( double  sbin)

swap bytes if needed (conversion from big/little endian toward little/big endian)

Definition at line 97 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_swapf()

float MMG5_swapf ( float  sbin)

swap bytes if needed (conversion from big/little endian toward little/big endian)

Definition at line 80 of file inout.c.

Here is the caller graph for this function:

◆ MMG5_sys33sym()

int MMG5_sys33sym ( double  a[6],
double  b[3],
double  r[3] 
)
inline
Parameters
amatrix to invert.
blast member.
rvector of unknowns.
Returns
0 if fail, 1 otherwise.

Solve $ 3\times 3$ symmetric system $ A . r = b $.

Definition at line 769 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_test_crossprod3d()

int MMG5_test_crossprod3d ( )
Returns
1 if success, 0 if fail.

Test vector cross product.

Definition at line 311 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_dotprod()

int MMG5_test_dotprod ( )
Returns
1 if success, 0 if fail.

Test vector scalar product.

Definition at line 276 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_eigenvmatnonsym2d()

int MMG5_test_eigenvmatnonsym2d ( MMG5_pMesh  mesh,
double *  mex,
double  lambdaex[],
double  vpex[][2],
double  ivpex[][2] 
)
Parameters
meshpointer toward the mesh structure.
mextest matrix array.
lambdaexexact eigenvalues array.
vpexdouble array of exact right eigenvectors.
ivpexdouble array of exact right eigenvectors inverse.
Returns
1 if success, 0 if failure.

For a 2x2 non-symmetric matrix, Test:

  • the recomposition of the matrix from its exact eigendecomposition;
  • the computation of the eigenvalues of the matrix;
  • the computation of the eigenvectors of the matrix;
  • the recomposition of the matrix from its numerical eigendecomposition.

Recompose matrix from its eigendecomposition

Compute eigendecomposition

Compute both eigendecomposition and recomposition, and check matrix

Definition at line 379 of file mettools.c.

Here is the call graph for this function:

◆ MMG5_test_eigenvmatnonsym3d()

int MMG5_test_eigenvmatnonsym3d ( MMG5_pMesh  mesh,
double *  mex,
double  lambdaex[],
double  vpex[][3],
double  ivpex[][3] 
)
Parameters
meshpointer toward the mesh structure.
mextest matrix array.
lambdaexexact eigenvalues array.
vpexdouble array of exact right eigenvectors.
ivpexdouble array of exact right eigenvectors inverse.
Returns
1 if success, 0 if failure.

For a 3x3 non-symmetric matrix, Test:

  • the recomposition of the matrix from its exact eigendecomposition;
  • the computation of the eigenvalues of the matrix;
  • the computation of the eigenvectors of the matrix;
  • the recomposition of the matrix from its numerical eigendecomposition.

Recompose matrix from its eigendecomposition

Compute eigendecomposition

Compute both eigendecomposition and recomposition, and check matrix

Definition at line 536 of file mettools.c.

Here is the call graph for this function:

◆ MMG5_test_eigenvmatsym2d()

int MMG5_test_eigenvmatsym2d ( MMG5_pMesh  mesh,
double *  mex,
double  lambdaex[],
double  vpex[][2] 
)
Parameters
meshpointer toward the mesh structure.
mextest matrix array.
lambdaexexact eigenvalues array.
vpexdouble array of exact eigenvectors.
Returns
1 if success, 0 if failure.

For a 2x2 symmetric matrix, Test:

  • the recomposition of the matrix from its exact eigendecomposition;
  • the computation of the eigenvalues of the matrix;
  • the computation of the eigenvectors of the matrix;
  • the recomposition of the matrix from its numerical eigendecomposition.

Recompose matrix from its eigendecomposition

Compute eigendecomposition

Compute both eigendecomposition and recomposition, and check matrix

Definition at line 301 of file mettools.c.

Here is the call graph for this function:

◆ MMG5_test_eigenvmatsym3d()

int MMG5_test_eigenvmatsym3d ( MMG5_pMesh  mesh,
double *  mex,
double  lambdaex[],
double  vpex[][3] 
)
Parameters
meshpointer toward the mesh structure.
mextest matrix array.
lambdaexexact eigenvalues array.
vpexdouble array of exact eigenvectors.
Returns
1 if success, 0 if failure.

For a 3x3 symmetric matrix, Test:

  • the recomposition of the matrix from its exact eigendecomposition;
  • the computation of the eigenvalues of the matrix;
  • the computation of the eigenvectors of the matrix; the recomposition of the matrix from its numerical eigendecomposition.

Recompose matrix from its eigendecomposition

Compute eigendecomposition

Compute both eigendecomposition and recomposition, and check matrix

Definition at line 456 of file mettools.c.

Here is the call graph for this function:

◆ MMG5_test_intersecmet22()

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

Test the intersection of (2 x 2) metrics.

Compute intersection m^{-1}n

Iteratively Recompute intersection of the result with (m,n), changing the matrix to be inverted.

Compute the intersection with n, invert n

Compute the intersection with n, invert intnum

Compute the intersection with n, invert n

Compute the intersection with m, invert intnum

Definition at line 962 of file mettools.c.

Here is the call graph for this function:

◆ MMG5_test_intersecmet33()

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

Test the intersection of (3 x 3) metrics.

Compute intersection m^{-1}n

Iteratively Recompute intersection of the result with (m,n), changing the matrix to be inverted.

Compute the intersection with n, invert n

Compute the intersection with n, invert intnum

Compute the intersection with n, invert n

Compute the intersection with m, invert intnum

Definition at line 1047 of file mettools.c.

Here is the call graph for this function:

◆ MMG5_test_invmat22()

int MMG5_test_invmat22 ( )

Test inversion of 2x2 non-symmetric matrix stored in 2 dimensions.

Definition at line 1140 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_invmat33()

int MMG5_test_invmat33 ( )

Test inversion of 3x3 non-symmetric matrix stored in 2 dimensions.

Definition at line 1165 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_mat_error()

double MMG5_test_mat_error ( int8_t  nelem,
double  m1[],
double  m2[] 
)
inline
Parameters
nelemnumber of matrix elements.
m1first matrix (single array).
m2second matrix (single array).

Compute maximum error between two matrices.

Definition at line 1123 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_test_mn()

int MMG5_test_mn ( )

Test product of 3x3 symmetric matrices.

Compute product m*n

Compute product n*m

Definition at line 357 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_rmtr()

int MMG5_test_rmtr ( )
inline

Test computation of product R*M*tR when M is symmetric

Compute transformation

Definition at line 435 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_rotmatrix()

int MMG5_test_rotmatrix ( )

Test computation of the rotation matrix that sends vector n to the third vector of canonical basis.

Rodrigues' rotation formula (transposed to give a map from n to [0,0,1]). Input vector must be a unit vector.

Approximate z-unit vector

Check orthonormality

Definition at line 512 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_simred2d()

int MMG5_test_simred2d ( MMG5_pMesh  mesh,
double *  mex,
double *  nex,
double *  dmex,
double *  dnex,
double  vpex[][2] 
)
Parameters
meshpointer toward the mesh structure
mexfirst symmetric test matrix
nexsecond symmetric test matrix
dmdiagonalization of the first matrix on the reduction basis
dndiagonalization of the second matrix on the reduction basis
vpsimultaneous reduction basis (stored by columns)
Returns
1 if success, 0 if fail

For a couple of 2x2 symmetric matrices, Test:

  • the computation of the simultaneous reduction values of the matrices;
  • the computation of the simultaneous reduction basis vectors..

Compute simultaneous reduction

Recompose matrices from diagonal values

Definition at line 1585 of file anisosiz.c.

Here is the call graph for this function:

◆ MMG5_test_simred3d()

int MMG5_test_simred3d ( MMG5_pMesh  mesh,
double *  mex,
double *  nex,
double *  dmex,
double *  dnex,
double  vpex[][3] 
)
Parameters
meshpointer toward the mesh structure
mexfirst symmetric test matrix
nexsecond symmetric test matrix
dmdiagonalization of the first matrix on the reduction basis
dndiagonalization of the second matrix on the reduction basis
vpsimultaneous reduction basis (stored by columns)
Returns
1 if success, 0 if fail

For a couple of 3x3 symmetric matrices, Test:

  • the computation of the simultaneous reduction values of the matrices;
  • the computation of the simultaneous reduction basis vectors..

Compute simultaneous reduction

Recompose matrices from diagonal values

Definition at line 1663 of file anisosiz.c.

Here is the call graph for this function:

◆ MMG5_test_transpose3d()

int MMG5_test_transpose3d ( )
Returns
1 if success, 0 if fail.

Test the transposition of a 3x3 matrix.

Definition at line 236 of file tools.c.

Here is the call graph for this function:

◆ MMG5_test_updatemet2d_ani()

int MMG5_test_updatemet2d_ani ( )

Test Update of the metrics = tP^-1 diag(d0,d1)P^-1, P = (vp[0], vp[1]) stored in columns in 2D.

Recompose matrices from exact simultaneous diagonalization

Definition at line 1805 of file anisosiz.c.

Here is the call graph for this function:

◆ MMG5_test_updatemet3d_ani()

int MMG5_test_updatemet3d_ani ( )

Test Update of the metrics = tP^-1 diag(d0,d1)P^-1, P = (vp[0], vp[1]) stored in columns in 3D.

Recompose matrices from exact simultaneous diagonalization

Definition at line 1841 of file anisosiz.c.

Here is the call graph for this function:

◆ MMG5_transpose3d()

void MMG5_transpose3d ( double  m[3][3])
Parameters
mmatrix (stored as double array)

Transpose a square matrix in 3D, stored as double array.

Definition at line 220 of file tools.c.

Here is the caller graph for this function:

◆ MMG5_truncate_met3d()

int MMG5_truncate_met3d ( MMG5_pSol  met,
MMG5_int  ip,
double  isqhmin,
double  isqhmax 
)
Parameters
metpointer toward metric.
ippointer toward global index of point on which metric has to be truncated
isqhmininverse square of hmin (min edge size)
isqhmaxinverse square of hmax (max edge size)
Returns
0 if fail, 1 if succeed

Truncation of anisotropic metric at point ip by hmin and hmax.

Definition at line 148 of file scalem.c.

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

◆ MMG5_updatemet2d_ani()

int MMG5_updatemet2d_ani ( double *  m,
double *  n,
double  dm[2],
double  dn[2],
double  vp[2][2],
int8_t  ier 
)
inline
Parameters
mfirst matrix
nsecond matrix
dmeigenvalues of m in the coreduction basis
dneigenvalues of n in the coreduction basis
vpcoreduction basis
ierflag of the updated sizes: (ier & 1) if we dm has been modified, (ier & 2) if dn has been modified.
Returns
0 if fail, 1 otherwise

Update of the metrics = tP^-1 diag(d0,d1)P^-1, P = (vp[0], vp[1]) stored in columns in 2D.

Definition at line 1742 of file anisosiz.c.

Here is the caller graph for this function:

◆ MMG5_updatemet3d_ani()

int MMG5_updatemet3d_ani ( double *  m,
double *  n,
double  dm[3],
double  dn[3],
double  vp[3][3],
int8_t  ier 
)
inline
Parameters
mfirst matrix
nsecond matrix
dmeigenvalues of m in the coreduction basis
dneigenvalues of n in the coreduction basis
vpcoreduction basis
ierflag of the updated sizes: (ier & 1) if we dm has been modified, (ier & 2) if dn has been modified.
Returns
0 if fail, 1 otherwise

Update of the metrics = tP^-1 diag(d0,d1)P^-1, P = (vp[0], vp[1]) stored in columns in 3D.

Definition at line 1783 of file anisosiz.c.

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

◆ MMG5_updatemetreq_ani()

int MMG5_updatemetreq_ani ( double *  n,
double  dn[2],
double  vp[2][2] 
)
Parameters
nmatrix to update
dneigenvalues of n in the coreduction basis
vpcoreduction basis
Returns
0 if fail, 1 otherwise

Update of the metric n = tP^-1 diag(dn0,dn1)P^-1, P = (vp[0], vp[1]) stored in columns

Definition at line 1914 of file anisosiz.c.

Here is the caller graph for this function:

◆ MMG5_version()

void MMG5_version ( MMG5_pMesh  mesh,
char *  dim 
)

Functions needed by libraries API.

Author
Charles Dapogny (UPMC)
Cécile Dobrzynski (Bx INP/Inria/UBordeaux)
Pascal Frey (UPMC)
Algiane Froehly (Inria/UBordeaux)
Version
5
Date
04 2015
Parameters
meshpointer toward the mesh
dimstring dontaining the dimension (3D,2D or S)

Print MMG release and date

Definition at line 43 of file libtools.c.

Here is the caller graph for this function:

◆ MMG5_warnScotch()

static void MMG5_warnScotch ( MMG5_pMesh  mesh)
inlinestatic

Inlined functions for libraries and executables Warn user that we overflow asked memory during scotch call

Definition at line 487 of file mmgcommon_private.h.

Here is the caller graph for this function:

◆ MMG5_writeDoubleSol3D()

void MMG5_writeDoubleSol3D ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
FILE *  inm,
int  bin,
MMG5_int  pos,
int  metricData 
)
Parameters
meshpointer toward the mesh structure
solpointer toward an allocatable sol structure.
inmpointer toward the solution file
bin1 if binary file
posof the writted solution
metricData1 if the data saved is a metric (if only 1 data)

Write the solution value for vertex of index pos in double precision.

Definition at line 2318 of file inout.c.

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

◆ MMG5_writeLocalParamAtTri()

int MMG5_writeLocalParamAtTri ( MMG5_pMesh  mesh,
MMG5_iNode bdryRefs,
FILE *  out 
)
Parameters
meshpointer toward the mesh structure.
bdryRefspointer toward the list of the boundary references.
outpointer toward the file in which to write.
Returns
1 if success, 0 otherwise.

Write the local default values at triangles in the parameter file.

Definition at line 186 of file apptools.c.

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

◆ mycalloc()

static void * mycalloc ( size_t  c,
size_t  s 
)
inlinestatic

Definition at line 226 of file mmgcommon_private.h.

◆ myfree()

static size_t myfree ( void *  ptr)
inlinestatic

Definition at line 277 of file mmgcommon_private.h.

◆ mymalloc()

static void * mymalloc ( size_t  s)
inlinestatic

Definition at line 239 of file mmgcommon_private.h.

Here is the caller graph for this function:

◆ myrealloc()

static void * myrealloc ( void *  ptr_in,
size_t  s,
size_t  oldsize 
)
inlinestatic

Definition at line 252 of file mmgcommon_private.h.

Here is the call graph for this function:

Variable Documentation

◆ MMG5_inxt2

const uint8_t MMG5_inxt2[6] = {1,2,0,1,2}
static

next vertex of triangle: {1,2,0}

Definition at line 583 of file mmgcommon_private.h.

◆ MMG5_iprv2

const uint8_t MMG5_iprv2[3] = {2,0,1}
static

previous vertex of triangle: {2,0,1}

Definition at line 584 of file mmgcommon_private.h.