mmg2d
mmgcommon.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 "eigenv.h"
#include "libmmgcommon.h"
Include dependency graph for mmgcommon.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
 
struct  MMG5_dNode_s
 

Macros

#define MG_STR   "&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&"
 
#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 MG_PLUS   2
 
#define MG_MINUS   3
 
#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_LMAX   10240
 
#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_EDG(tag)   ((tag & MG_GEO) || (tag & MG_REF))
 
#define MG_SIN(tag)   ((tag & MG_CRN) || (tag & MG_REQ))
 
#define MG_RID(tag)   ((!( MG_SIN(tag)||(tag & MG_NOM))) && tag & MG_GEO )
 
#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_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 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
 
typedef struct MMG5_dNode_s MMG5_dNode
 

Enumerations

enum  MMG5_Format {
  MMG5_FMT_MeditASCII, MMG5_FMT_MeditBinary, MMG5_FMT_GmshASCII, MMG5_FMT_GmshBinary,
  MMG5_FMT_VtkPvtp, MMG5_FMT_VtkPvtu, MMG5_FMT_VtkVtu, MMG5_FMT_VtkVtp,
  MMG5_FMT_VtkVtk, MMG5_FMT_Tetgen, MMG5_FMT_Unknown
}
 Type of supported file format. More...
 

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_excfun (int sigid)
 
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, int *v)
 
int MMG5_Add_inode (MMG5_pMesh mesh, MMG5_iNode **liLi, int val)
 
int MMG5_Add_dnode (MMG5_pMesh mesh, MMG5_dNode **liLi, int, double)
 
void MMG5_bezierEdge (MMG5_pMesh, int, int, double *, double *, int8_t, double *)
 
int MMG5_buildridmet (MMG5_pMesh, MMG5_pSol, 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, 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_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, int start, int ip, int *, int *list)
 
int MMG5_boulec (MMG5_pMesh, int *, int, int i, double *tt)
 
int MMG5_boulen (MMG5_pMesh, int *, int, int i, double *nn)
 
int MMG5_bouler (MMG5_pMesh, int *, int, int i, int *, int *, int *, int *, int)
 
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, int, double *, int, int, double, int, int, double, int, double *, int *, int8_t)
 
void MMG5_displayLengthHisto_internal (int, int, int, double, int, int, double, int, double *, int *, int8_t, int)
 
int MMG5_minQualCheck (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 (int, MMG5_pPoint, int, MMG5_Bezier, double r[3][3], double *, double *, double *, double *)
 
void MMG5_Free_ilinkedList (MMG5_pMesh mesh, MMG5_iNode *liLi)
 
void MMG5_Free_dlinkedList (MMG5_pMesh mesh, MMG5_dNode *liLi)
 
int MMG5_grad2metSurf (MMG5_pMesh, MMG5_pSol, MMG5_pTria, int, int)
 
int MMG5_grad2metSurfreq (MMG5_pMesh, MMG5_pSol, MMG5_pTria, int, int)
 
char * MMG5_Get_filenameExt (char *filename)
 
char * MMG5_Get_basename (char *path)
 
char * MMG5_Get_path (char *path)
 
char * MMG5_Remove_ext (char *path, char *)
 
const char * MMG5_Get_formatName (enum MMG5_Format fmt)
 
int MMG5_Get_format (char *ptr, int)
 
int MMG5_hashEdge (MMG5_pMesh mesh, MMG5_Hash *hash, int a, int b, int k)
 
int MMG5_hashUpdate (MMG5_Hash *hash, int a, int b, int k)
 
int MMG5_hashEdgeTag (MMG5_pMesh mesh, MMG5_Hash *hash, int a, int b, int16_t k)
 
int MMG5_hashGet (MMG5_Hash *hash, int a, int b)
 
int MMG5_hashNew (MMG5_pMesh mesh, MMG5_Hash *hash, int hsiz, int hmax)
 
int MMG5_intmetsavedir (MMG5_pMesh mesh, double *m, double *n, double *mr)
 
int MMG5_intridmet (MMG5_pMesh, MMG5_pSol, int, int, double, double *, double *)
 
int MMG5_mmgIntmet33_ani (double *, double *, double *, double)
 
int MMG5_mmgIntextmet (MMG5_pMesh, MMG5_pSol, 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, 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, int, int, int, double *)
 
double MMG5_nonorsurf (MMG5_pMesh mesh, MMG5_pTria pt)
 
int MMG5_norpts (MMG5_pMesh, int, int, 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_regnor (MMG5_pMesh mesh)
 
double MMG5_ridSizeInNormalDir (MMG5_pMesh, int, double *, MMG5_pBezier, double, double)
 
double MMG5_ridSizeInTangentDir (MMG5_pMesh, MMG5_pPoint, int, int *, double, double)
 
int MMG5_scale_meshAndSol (MMG5_pMesh, MMG5_pSol, MMG5_pSol, double *, int8_t *, int8_t *)
 
int MMG5_scale_scalarMetric (MMG5_pMesh, MMG5_pSol, double, int8_t, int8_t)
 
int MMG5_scaleMesh (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol ls)
 
int MMG5_scotchCall (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pSol fields, int *)
 
void MMG5_solTruncatureForOptim (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_solveDefmetregSys (MMG5_pMesh, double r[3][3], double *, double *, double *, double *, double, double, double)
 
int MMG5_solveDefmetrefSys (MMG5_pMesh, MMG5_pPoint, 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_unscaleMesh (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol ls)
 
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, 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, int *nelts, int *nsols)
 
int MMG5_check_readedMesh (MMG5_pMesh mesh, 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 int nelts, const int nsols)
 
int MMG5_saveMshMesh (MMG5_pMesh, MMG5_pSol *, const char *, int)
 
int MMG5_loadSolHeader (const char *, int, FILE **, int *, int *, int *, int *, int *, int *, int **, long *, int)
 
int MMG5_chkMetricType (MMG5_pMesh mesh, int *type, FILE *inm)
 
int MMG5_readFloatSol3D (MMG5_pSol, FILE *, int, int, int)
 
int MMG5_readDoubleSol3D (MMG5_pSol, FILE *, int, int, int)
 
int MMG5_saveSolHeader (MMG5_pMesh, const char *, FILE **, int, int *, int, int, int, int *, int *)
 
void MMG5_writeDoubleSol3D (MMG5_pMesh, MMG5_pSol, FILE *, int, 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, int ip0, 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)
 
int MMG5_gradsiz_ani (MMG5_pMesh mesh, MMG5_pSol met, int *it)
 
int MMG5_gradsizreq_ani (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMG5_simred (MMG5_pMesh, double *, double *, double dm[2], double dn[2], double vp[2][2])
 
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)
 
float MMG5_swapf (float sbin)
 
double MMG5_swapd (double sbin)
 
int MMG5_isSplit (MMG5_pMesh, int, int *, int *)
 
int MMG5_getIniRef (MMG5_pMesh, int)
 
void MMG5_mark_verticesAsUnused (MMG5_pMesh mesh)
 
void MMG5_mark_usedVertices (MMG5_pMesh mesh, void(*delPt)(MMG5_pMesh, int))
 
void MMG5_keep_subdomainElts (MMG5_pMesh, int, int(*delElt)(MMG5_pMesh, 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}
 
int(* MMG5_chkmsh )(MMG5_pMesh, int, int)
 
int(* MMG5_bezierCP )(MMG5_pMesh, MMG5_Tria *, MMG5_pBezier, int8_t)
 
double(* MMG5_lenSurfEdg )(MMG5_pMesh mesh, MMG5_pSol sol, int, int, int8_t)
 
int(* MMG5_grad2met_ani )(MMG5_pMesh, MMG5_pSol, MMG5_pTria, int, int)
 
int(* MMG5_grad2metreq_ani )(MMG5_pMesh, MMG5_pSol, MMG5_pTria, int, int)
 
int(* MMG5_compute_meanMetricAtMarkedPoints )(MMG5_pMesh, MMG5_pSol)
 
int(* MMG5_indElt )(MMG5_pMesh mesh, int kel)
 
int(* MMG5_indPt )(MMG5_pMesh mesh, int kp)
 

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)

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)

◆ A16TH

#define A16TH   0.0625

◆ A32TH

#define A32TH   0.03125

◆ A64TH

#define A64TH   0.015625

◆ 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

◆ 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

◆ 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.

◆ 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.

◆ MG_BDY

#define MG_BDY   (1 << 4)

16 boundary entity

◆ MG_CLR

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

bit number bit is set to 0

◆ MG_CRN

#define MG_CRN   (1 << 5)

32 corner

◆ MG_EDG

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

Edge or Ridge

◆ MG_Edge

#define MG_Edge   (1 << 3 )

8 local parameter applied over edge

◆ MG_EOK

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

Element OK

◆ MG_GEO

#define MG_GEO   (1 << 1)

2 geometric ridge

◆ MG_GET

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

return bit number bit value

◆ MG_MAX

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

◆ MG_MIN

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

◆ MG_MINUS

#define MG_MINUS   3

◆ MG_NOM

#define MG_NOM   (1 << 3)

8 non manifold

◆ MG_NOSURF

#define MG_NOSURF   (1 << 6)

64 freezed boundary

◆ MG_NOTAG

#define MG_NOTAG   (0)

◆ MG_NUL

#define MG_NUL   (1 << 14)

16384 vertex removed

◆ MG_OLDPARBDY

#define MG_OLDPARBDY   (1 << 11)

2048 old parallel boundary

◆ MG_OPNBDY

#define MG_OPNBDY   (1 << 7)

128 open boundary

◆ MG_PARBDY

#define MG_PARBDY   (1 << 13)

8192 parallel boundary

◆ MG_PARBDYBDY

#define MG_PARBDYBDY   (1 << 12)

4096 parallel boundary over a boundary

◆ MG_PLUS

#define MG_PLUS   2

◆ MG_REF

#define MG_REF   (1 << 0)

1 edge reference

◆ MG_REQ

#define MG_REQ   (1 << 2)

4 required entity

◆ MG_RID

#define MG_RID (   tag)    ((!( MG_SIN(tag)||(tag & MG_NOM))) && tag & MG_GEO )

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

◆ MG_SET

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

bit number bit is set to 1

◆ MG_SIN

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

Corner or Required

◆ MG_SMSGN

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

Check if a and b have the same sign

◆ MG_STR

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

◆ MG_Tetra

#define MG_Tetra   (1 << 2 )

4 local parameter applied over tetrahedron

◆ MG_Tria

#define MG_Tria   (1 << 1 )

2 local parameter applied over triangle

◆ MG_Vert

#define MG_Vert   (1 << 0 )

1 local parameter applied over vertex

◆ MG_VOK

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

Vertex OK

◆ 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.

◆ MMG5_ANGEDG

#define MMG5_ANGEDG   0.707106781186548 /*0.573576436351046 */

◆ MMG5_ANGLIM

#define MMG5_ANGLIM   -0.999999

◆ MMG5_ATHIRD

#define MMG5_ATHIRD   0.333333333333333

◆ MMG5_BITWIZE_MB_TO_B

#define MMG5_BITWIZE_MB_TO_B   20

Bitwise convertion from Mo to O

◆ MMG5_BOXSIZE

#define MMG5_BOXSIZE   500

size of box for renumbering with scotch.

◆ 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.

◆ 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)

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

◆ MMG5_DISPREF

#define MMG5_DISPREF   10

◆ MMG5_EPS

#define MMG5_EPS   1.e-06

◆ MMG5_EPSD

#define MMG5_EPSD   1.e-30

◆ MMG5_EPSD2

#define MMG5_EPSD2   1.0e-200

◆ MMG5_EPSOK

#define MMG5_EPSOK   1.e-15

◆ MMG5_FEM

#define MMG5_FEM   1

defaut value for FEM mode

◆ MMG5_FILESTR_LGTH

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

◆ MMG5_GAP

#define MMG5_GAP   0.2

gap value for reallocation

◆ MMG5_HAUSD

#define MMG5_HAUSD   0.01

default value for hausdorff param

◆ MMG5_HGRAD

#define MMG5_HGRAD   0.26236426446

default value for gradation (1.3)

◆ MMG5_HGRADREQ

#define MMG5_HGRADREQ   0.83290912294

default value for required gradation (2.3)

◆ MMG5_HMAXCOE

#define MMG5_HMAXCOE   2

coefficient to compute the default hmax value

◆ MMG5_HMINCOE

#define MMG5_HMINCOE   0.001

coefficient to compute the default hmin value

◆ MMG5_HMINMAXGAP

#define MMG5_HMINMAXGAP   5

imposed gap between hmin and hmax if hmax<hmin

◆ 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

◆ MMG5_KA

#define MMG5_KA   7

Key for hash tables.

◆ MMG5_KB

#define MMG5_KB   11

Key for hash tables.

◆ MMG5_LAG

#define MMG5_LAG   -1

default value for lagrangian option

◆ MMG5_LMAX

#define MMG5_LMAX   10240

◆ MMG5_LS

#define MMG5_LS   0.0

default level-set to discretize

◆ MMG5_MEMMAX

#define MMG5_MEMMAX   800

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

◆ MMG5_MEMMIN

#define MMG5_MEMMIN   38

minimal memory needed to store the mesh/sol names

◆ MMG5_MEMPERCENT

#define MMG5_MEMPERCENT   0.5

Percent of RAM used by default

◆ MMG5_MILLION

#define MMG5_MILLION   1048576

◆ MMG5_NOHGRAD

#define MMG5_NOHGRAD   -1

disable gradation

◆ MMG5_NONSET

#define MMG5_NONSET   -1

◆ MMG5_NONSET_HMAX

#define MMG5_NONSET_HMAX   -1

hmax value for unspecified hmax size

◆ MMG5_NONSET_HMIN

#define MMG5_NONSET_HMIN   -1

hmin value for unspecified hmin size

◆ MMG5_NONSET_HSIZ

#define MMG5_NONSET_HSIZ   -1

hsiz value for unspecified hsiz map

◆ MMG5_NONSET_MEM

#define MMG5_NONSET_MEM   -1

mem value for unspecified max memory

◆ MMG5_NR

#define MMG5_NR   -1

no ridge detection

◆ MMG5_NULKAL

#define MMG5_NULKAL   1.e-30

◆ MMG5_OFF

#define MMG5_OFF   0

0

◆ MMG5_ON

#define MMG5_ON   1

1

◆ MMG5_PATHSEP

#define MMG5_PATHSEP   '/'

◆ MMG5_PROCTREE

#define MMG5_PROCTREE   32

default size of the PROctree

◆ 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)

Safe allocation with calloc

◆ MMG5_SAFE_FREE

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

Safe deallocation

◆ 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)

Safe allocation with malloc

◆ 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)

Safe reallocation

◆ 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

◆ MMG5_SAFELL2ICAST

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

◆ MMG5_SAFELL2LCAST

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

◆ MMG5_SD

#define MMG5_SD   8

◆ MMG5_SQR32

#define MMG5_SQR32   0.866025403784439

◆ MMG5_SW

#define MMG5_SW   4

◆ MMG5_TAB_RECALLOC

#define MMG5_TAB_RECALLOC (   mesh,
  ptr,
  initSize,
  wantedGap,
  type,
  message,
  law 
)
Value:
do \
{ \
int gap; \
\
assert ( mesh->memCur < mesh->memMax ); \
\
gap = (int)(wantedGap * initSize); \
if ( !gap ) gap = 1; \
if ( mesh->memMax < mesh->memCur + gap*sizeof(type) ) { \
gap = (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);

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.

◆ MMG5_UNSET

#define MMG5_UNSET   -1

◆ 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

◆ 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);

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

Typedef Documentation

◆ MMG5_dNode

typedef struct MMG5_dNode_s MMG5_dNode

◆ MMG5_iNode

typedef struct MMG5_iNode_s MMG5_iNode

◆ MMG5_pBezier

Enumeration Type Documentation

◆ MMG5_Format

Type of supported file format.

Enumerator
MMG5_FMT_MeditASCII 

ASCII Medit (.mesh)

MMG5_FMT_MeditBinary 

Binary Medit (.meshb)

MMG5_FMT_GmshASCII 

ASCII Gmsh

MMG5_FMT_GmshBinary 

Binary Gmsh

MMG5_FMT_VtkPvtp 

VTK pvtp

MMG5_FMT_VtkPvtu 

VTK pvtu

MMG5_FMT_VtkVtu 

VTK vtu

MMG5_FMT_VtkVtp 

VTK vtp

MMG5_FMT_VtkVtk 

VTK vtk

MMG5_FMT_Tetgen 

Tetgen or Triangle

MMG5_FMT_Unknown 

Unrecognized

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.

Here is the caller graph for this function:

◆ MMG5_2d3dUsage()

void MMG5_2d3dUsage ( void  )

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

Here is the caller graph for this function:

◆ MMG5_Add_dnode()

int MMG5_Add_dnode ( MMG5_pMesh  mesh,
MMG5_dNode **  liLi,
int  k,
double  val 
)
Parameters
meshpointer toward the mesh structure (for count of used memory).
liLipointer toward the address of the root of the linked list.
kinteger value to add to the linked list.
valreal value 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 integer value k and real 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.
Here is the call 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.
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.

Here is the caller graph for this function:

◆ MMG5_bezierEdge()

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

◆ MMG5_boulec()

int MMG5_boulec ( MMG5_pMesh  mesh,
int *  adjt,
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.

◆ MMG5_boulen()

int MMG5_boulen ( MMG5_pMesh  mesh,
int *  adjt,
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.
ipindex 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.

Here is the call graph for this function:

◆ MMG5_boulep()

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

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

Here is the caller graph for this function:

◆ MMG5_bouler()

int MMG5_bouler ( MMG5_pMesh  mesh,
int *  adjt,
int  start,
int  ip,
int *  list,
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.

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.

Here is the caller graph for this function:

◆ MMG5_build3DMetric()

void MMG5_build3DMetric ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
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).

Here is the caller graph for this function:

◆ MMG5_buildridmet()

int MMG5_buildridmet ( MMG5_pMesh  ,
MMG5_pSol  ,
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  ,
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.

Here is the call 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.
Here is the call 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.

◆ 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.

Here is the caller graph for this function:

◆ MMG5_check_readedMesh()

int MMG5_check_readedMesh ( MMG5_pMesh  mesh,
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.

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,
FILE *  inm 
)
Parameters
meshpointer toward the mesh structure.
typetype of the metric
inmmetric file
Returns
1 if success, -1 if fail

Check if the type of the metric is compatible with the remeshing mode. If not, deallocate the type array and close the metric file.

Here is the call graph for this function:
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.

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.

Here is the caller graph for this function:

◆ MMG5_countLocalParamAtTri()

int MMG5_countLocalParamAtTri ( MMG5_pMesh  mesh,
MMG5_iNode **  bdryRefs 
)
inline
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

Here is the call graph for this function:
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.

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.

Here is the call 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)

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)

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).

◆ MMG5_displayLengthHisto()

void MMG5_displayLengthHisto ( MMG5_pMesh  mesh,
int  ned,
double *  avlen,
int  amin,
int  bmin,
double  lmin,
int  amax,
int  bmax,
double  lmax,
int  nullEdge,
double *  bd,
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.

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

◆ MMG5_displayLengthHisto_internal()

void MMG5_displayLengthHisto_internal ( int  ned,
int  amin,
int  bmin,
double  lmin,
int  amax,
int  bmax,
double  lmax,
int  nullEdge,
double *  bd,
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

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.

Here is the call graph for this function:

◆ MMG5_excfun()

static void MMG5_excfun ( int  sigid)
inlinestatic

Inlined functions for libraries and executables

Parameters
sigidsignal number.

Signal handling: specify error messages depending from catched signal.

◆ MMG5_fillDefmetregSys()

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

◆ MMG5_Free_dlinkedList()

void MMG5_Free_dlinkedList ( MMG5_pMesh  mesh,
MMG5_dNode 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.

◆ 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.

Here is the caller graph for this function:

◆ MMG5_Get_basename()

char* MMG5_Get_basename ( char *  path)
Parameters
pathstring containing a filename and its path
Returns
a pointer toward the allocated string that contains the file basename.

Extract basename from a path (allocate a string to store it).

Here is the caller graph for this function:

◆ MMG5_Get_filenameExt()

char* MMG5_Get_filenameExt ( char *  filename)
Parameters
filenamestring containing a filename
Returns
pointer toward the filename extension or toward the end of the string if no extension have been founded

Get the extension of the filename string. Do not consider '.o' as an extension.

Here is the caller graph for this function:

◆ MMG5_Get_format()

int MMG5_Get_format ( char *  ptr,
int  fmt 
)
Parameters
ptrpointer toward the file extension (dot included)
fmtdefault file format.
Returns
and index associated to the file format detected from the extension.

Get the wanted file format from the mesh extension. If fmt is provided, it is used as default file format (ptr==NULL), otherwise, the default file format is the medit one.

Here is the caller graph for this function:

◆ MMG5_Get_formatName()

const char* MMG5_Get_formatName ( enum MMG5_Format  fmt)
Parameters
fmtfile format.
Returns
The name of the file format in a string.

Print the name of the file format associated to fmt.

Here is the caller graph for this function:

◆ MMG5_Get_path()

char* MMG5_Get_path ( char *  path)
Parameters
pathstring containing a filename and its path
Returns
a pointer toward the path allocated here

Remove filename from a path and return the path in a newly allocated string.

Here is the call graph for this function:

◆ MMG5_getIniRef()

int MMG5_getIniRef ( MMG5_pMesh  mesh,
int  ref 
)
Parameters
meshpointer toward the mesh
reffinal reference for which we are searching the initial one
Returns
initial reference associated to ref if founded, ref if not founded.

Retrieve the initial domain reference associated to the (split) reference ref.

Here is the caller graph for this function:

◆ MMG5_grad2metSurf()

int MMG5_grad2metSurf ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  pt,
int  np1,
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.
Here is the call graph for this function:

◆ MMG5_grad2metSurfreq()

int MMG5_grad2metSurfreq ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pTria  pt,
int  npmaster,
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".
Here is the call 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).

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.

Here is the caller graph for this function:

◆ MMG5_gradsiz_ani()

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

Here is the call 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.

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)

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

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,
int  a,
int  b,
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.

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,
int  a,
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.

Here is the call graph for this function:

◆ MMG5_hashGet()

int MMG5_hashGet ( MMG5_Hash hash,
int  a,
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]$.

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,
int  hsiz,
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.

Here is the caller graph for this function:

◆ MMG5_hashUpdate()

int MMG5_hashUpdate ( MMG5_Hash hash,
int  a,
int  b,
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]$

◆ 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

◆ 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.

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.

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  ,
int  ,
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

◆ 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

◆ MMG5_invmatg()

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

Invert 3x3 non-symmetric matrix.

Here is the caller graph for this function:

◆ MMG5_isSplit()

int MMG5_isSplit ( MMG5_pMesh  mesh,
int  ref,
int *  refint,
int *  refext 
)
Parameters
meshpointer toward the mesh structure.
refinitial reference.
refintinternal reference after ls discretization.
refintinternal 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.

Here is the caller graph for this function:

◆ MMG5_keep_subdomainElts()

void MMG5_keep_subdomainElts ( MMG5_pMesh  mesh,
int  nsd,
int(*)(MMG5_pMesh, 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

Here is the caller graph for this function:

◆ MMG5_lagUsage()

void MMG5_lagUsage ( void  )

Print help for lagrangian motion option.

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,
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.

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 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, 0 if fail.

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

Second step: read the nodes and elements

Read the solution at nodes

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,
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.

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.

Here is the caller graph for this function:

◆ MMG5_mark_usedVertices()

void MMG5_mark_usedVertices ( MMG5_pMesh  mesh,
void(*)(MMG5_pMesh, 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).

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.

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.

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.

Here is the caller graph for this function:

◆ MMG5_minQualCheck()

int MMG5_minQualCheck ( 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.

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.

Here is the caller graph for this function:

◆ MMG5_mmgHashTria()

int MMG5_mmgHashTria ( MMG5_pMesh  mesh,
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.
Here is the call graph for this function:

◆ MMG5_mmgInit_parameters()

void MMG5_mmgInit_parameters ( MMG5_pMesh  mesh)

◆ MMG5_mmgIntextmet()

int MMG5_mmgIntextmet ( MMG5_pMesh  ,
MMG5_pSol  ,
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.

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).

Here is the caller graph for this function:

◆ MMG5_mn()

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

Compute product m*n (mn stored in columns: mn[1] = mn[1][0]).

◆ 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.

Here is the call graph for this function:

◆ MMG5_nonUnitNorPts()

int MMG5_nonUnitNorPts ( MMG5_pMesh  mesh,
int  ip1,
int  ip2,
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.

Here is the caller graph for this function:

◆ MMG5_norpts()

int MMG5_norpts ( MMG5_pMesh  mesh,
int  ip1,
int  ip2,
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.

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.

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

◆ MMG5_orvol()

double MMG5_orvol ( MMG5_pPoint  point,
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)

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).

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).

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.

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

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

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.

◆ MMG5_readDoubleSol3D()

int MMG5_readDoubleSol3D ( 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 double precision.

Here is the call 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.

Here is the call 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

Here is the call graph for this function:

◆ MMG5_Remove_ext()

char* MMG5_Remove_ext ( char *  path,
char *  ext 
)
Parameters
pathpath from which we want to remove the extension.
Returns
allocated string or NULL if the allocation fail.

Allocate a new string and copy path without extension in it.

Here is the call 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.

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,
int  idp,
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 of a ridge in the direction of the tangent of the ridge.

Here is the call 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

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.

Here is the caller graph for this function:

◆ 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"

Here is the call graph for this function:
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,
int  np,
int  dim,
int  nsols,
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.
npnumber of solutions of each type.
dimsolution dimension.
nsolsnumber of solutions of different types in the file.
typetype of solutions.
sizesize of solutions.
Returns
0 if unable to open the file, 1 if success.

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

Here is the call graph for this function:
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,
int8_t *  sethmin,
int8_t *  sethmax 
)
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)
sethminsetted to 1 if hmin must not be computed from the metric.
sethmaxsetted to 1 if hmax must not be computed from the metric.
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.

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,
int8_t  sethmin,
int8_t  sethmax 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
ddscaling value.
sethmin1 if hmin must not be automatically computed
sethmax1 if hmin must not be automatically computed
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.

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

◆ MMG5_scaleMesh()

int MMG5_scaleMesh ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward the metric structure.
solpointer toward a solution structure (level-set or displacement).
Returns
1 if success, 0 if fail (computed bounding box too small or one af the anisotropic input metric is not valid).

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

Here is the call graph for this function:

◆ MMG5_scotchCall()

int MMG5_scotchCall ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  fields,
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.

◆ MMG5_Set_commonFunc()

void MMG5_Set_commonFunc ( void  )

◆ MMG5_simred()

int MMG5_simred ( 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.

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

◆ MMG5_solTruncatureForOptim()

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

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

Warning
works only for a metric computed by the DoSol function because we suppose that we have a diagonal tensor in aniso.
Here is the call graph for this function:

◆ MMG5_solveDefmetrefSys()

int MMG5_solveDefmetrefSys ( MMG5_pMesh  ,
MMG5_pPoint  ,
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,
int  ip0,
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.

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).

Here is the call 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.

◆ MMG5_swapbin()

int MMG5_swapbin ( int  sbin)
Here is the caller graph for this function:

◆ MMG5_swapd()

double MMG5_swapd ( double  sbin)
Here is the caller graph for this function:

◆ MMG5_swapf()

float MMG5_swapf ( float  sbin)
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 $.

Here is the caller graph for this function:

◆ MMG5_unscaleMesh()

int MMG5_unscaleMesh ( MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure.
metpointer toward a metric.
solpointer toward a solution structure (level-set or displacement).
Returns
1.

Unscale the mesh and the size informations to their initial sizes.

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

Here is the caller graph for this function:

◆ MMG5_writeDoubleSol3D()

void MMG5_writeDoubleSol3D ( MMG5_pMesh  mesh,
MMG5_pSol  sol,
FILE *  inm,
int  bin,
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.

Here is the call graph for this function:

◆ MMG5_writeLocalParamAtTri()

int MMG5_writeLocalParamAtTri ( MMG5_pMesh  mesh,
MMG5_iNode bdryRefs,
FILE *  out 
)
inline
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.

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

◆ myfree()

static size_t myfree ( void *  ptr)
inlinestatic

◆ mymalloc()

static void* mymalloc ( size_t  s)
inlinestatic
Here is the caller graph for this function:

◆ myrealloc()

static void* myrealloc ( void *  ptr_in,
size_t  s,
size_t  oldsize 
)
inlinestatic
Here is the call graph for this function:

Variable Documentation

◆ MMG5_bezierCP

int(* MMG5_bezierCP) (MMG5_pMesh,MMG5_Tria *, MMG5_pBezier,int8_t)

◆ MMG5_chkmsh

int(* MMG5_chkmsh) (MMG5_pMesh, int, int)

◆ MMG5_compute_meanMetricAtMarkedPoints

int(* MMG5_compute_meanMetricAtMarkedPoints) (MMG5_pMesh, MMG5_pSol)

◆ MMG5_grad2met_ani

int(* MMG5_grad2met_ani) (MMG5_pMesh, MMG5_pSol, MMG5_pTria, int, int)

◆ MMG5_grad2metreq_ani

int(* MMG5_grad2metreq_ani) (MMG5_pMesh, MMG5_pSol, MMG5_pTria, int, int)

◆ MMG5_indElt

int(* MMG5_indElt) (MMG5_pMesh mesh, int kel)

◆ MMG5_indPt

int(* MMG5_indPt) (MMG5_pMesh mesh, int kp)

◆ MMG5_inxt2

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

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

◆ MMG5_iprv2

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

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

◆ MMG5_lenSurfEdg

double(* MMG5_lenSurfEdg) (MMG5_pMesh mesh, MMG5_pSol sol,int,int, int8_t)
CV_VA_NUM_ARGS
#define CV_VA_NUM_ARGS(...)
Definition: mmgcommon.h:413
myfree
static size_t myfree(void *ptr)
Definition: mmgcommon.h:267
MMG5_Mesh::memCur
size_t memCur
Definition: libmmgtypes.h:554
MMG5_Mesh::memMax
size_t memMax
Definition: libmmgtypes.h:553
sol
MMG5_pMesh MMG5_pSol * sol
Definition: API_functionsf_2d.c:63
MMG5_Mesh::np
int np
Definition: libmmgtypes.h:559
mymalloc
static void * mymalloc(size_t s)
Definition: mmgcommon.h:229
MMG5_Mesh::ne
int ne
Definition: libmmgtypes.h:559
mesh
MMG5_pMesh * mesh
Definition: API_functionsf_2d.c:63
MMG5_ADD_MEM
#define MMG5_ADD_MEM(mesh, size, message, law)
Definition: mmgcommon.h:290
MMG5_Mesh::npi
int npi
Definition: libmmgtypes.h:559
mycalloc
static void * mycalloc(size_t c, size_t s)
Definition: mmgcommon.h:216
MMG5_Mesh::na
int na
Definition: libmmgtypes.h:559
myrealloc
static void * myrealloc(void *ptr_in, size_t s, size_t oldsize)
Definition: mmgcommon.h:242
MMG5_Mesh::nt
int nt
Definition: libmmgtypes.h:559
MMG5_dNode_s::val
double val
Definition: mmgcommon.h:580
tmp
tmp[*strlen0]
Definition: API_functionsf_2d.c:771
if
if(!ier) exit(EXIT_FAILURE)