Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
|
#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"
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_XREG 0.4 |
#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_OVERLAP (1 << 14) |
#define | MG_NUL (1 << 15) |
#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. | |
#define | FORTRAN_VARIADIC(nu, nl, pl, body) |
Adds function definitions. | |
Typedefs | |
typedef MMG5_Bezier * | MMG5_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. | |
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. | |
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, uint16_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} |
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 194 of file mmgcommon_private.h.
#define A16TH 0.0625 |
Definition at line 108 of file mmgcommon_private.h.
#define A32TH 0.03125 |
Definition at line 109 of file mmgcommon_private.h.
#define A64TH 0.015625 |
Definition at line 107 of file mmgcommon_private.h.
#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 455 of file mmgcommon_private.h.
#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 452 of file mmgcommon_private.h.
#define FORTRAN_NAME | ( | nu, | |
nl, | |||
pl, | |||
pc | |||
) |
Adds function definitions.
nu | function name in upper case. |
nl | function name in lower case. |
pl | type of arguments. |
pc | name of arguments. |
Adds function definitions with upcase, underscore and double underscore to match any fortran compiler.
Definition at line 551 of file mmgcommon_private.h.
#define FORTRAN_VARIADIC | ( | nu, | |
nl, | |||
pl, | |||
body | |||
) |
Adds function definitions.
nu | function name in upper case. |
nl | function name in lower case. |
pl | type of arguments. |
body | body of the function. |
Adds function definitions with upcase, underscore and double underscore to match any fortran compiler.
Definition at line 573 of file mmgcommon_private.h.
#define FUNCTION_POINTER | ( | fproto | ) | MMG_EXTERN fproto MMG_ASSIGN_NULL |
fproto | function 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:
Definition at line 532 of file mmgcommon_private.h.
#define MG_BDY (1 << 4) |
16 boundary entity
Definition at line 149 of file mmgcommon_private.h.
#define MG_CLR | ( | flag, | |
bit | |||
) | ((flag) &= ~(1 << (bit))) |
bit number bit is set to 0
Definition at line 181 of file mmgcommon_private.h.
#define MG_CRN (1 << 5) |
32 corner
Definition at line 150 of file mmgcommon_private.h.
Edge or Ridge
Definition at line 173 of file mmgcommon_private.h.
Edge, ridge or non-manifold
Definition at line 175 of file mmgcommon_private.h.
#define MG_Edge (1 << 3 ) |
8 local parameter applied over edge
Definition at line 163 of file mmgcommon_private.h.
#define MG_EOK | ( | pt | ) | (pt && ((pt)->v[0] > 0)) |
Element OK
Definition at line 167 of file mmgcommon_private.h.
#define MG_GEO (1 << 1) |
2 geometric ridge
Definition at line 146 of file mmgcommon_private.h.
Ridge or non-manifold
Definition at line 174 of file mmgcommon_private.h.
#define MG_GET | ( | flag, | |
bit | |||
) | ((flag) & (1 << (bit))) |
return bit number bit value
Definition at line 182 of file mmgcommon_private.h.
#define MG_MAX | ( | a, | |
b | |||
) | (((a) > (b)) ? (a) : (b)) |
Definition at line 140 of file mmgcommon_private.h.
#define MG_MIN | ( | a, | |
b | |||
) | (((a) < (b)) ? (a) : (b)) |
Definition at line 141 of file mmgcommon_private.h.
#define MG_NOM (1 << 3) |
8 non manifold
Definition at line 148 of file mmgcommon_private.h.
#define MG_NOSURF (1 << 6) |
64 freezed boundary
Definition at line 151 of file mmgcommon_private.h.
#define MG_NOTAG (0) |
Definition at line 144 of file mmgcommon_private.h.
#define MG_NUL (1 << 15) |
32768 vertex removed
Definition at line 157 of file mmgcommon_private.h.
#define MG_OLDPARBDY (1 << 11) |
2048 old parallel boundary
Definition at line 153 of file mmgcommon_private.h.
#define MG_OPNBDY (1 << 7) |
128 open boundary
Definition at line 152 of file mmgcommon_private.h.
#define MG_OVERLAP (1 << 14) |
16384 elements on overlap
Definition at line 156 of file mmgcommon_private.h.
#define MG_PARBDY (1 << 13) |
8192 parallel boundary
Definition at line 155 of file mmgcommon_private.h.
#define MG_PARBDYBDY (1 << 12) |
4096 parallel boundary over a boundary
Definition at line 154 of file mmgcommon_private.h.
#define MG_REF (1 << 0) |
1 edge reference
Definition at line 145 of file mmgcommon_private.h.
#define MG_REQ (1 << 2) |
4 required entity
Definition at line 147 of file mmgcommon_private.h.
#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 171 of file mmgcommon_private.h.
#define MG_SET | ( | flag, | |
bit | |||
) | ((flag) |= (1 << (bit))) |
bit number bit is set to 1
Definition at line 180 of file mmgcommon_private.h.
Corner or Required
Definition at line 169 of file mmgcommon_private.h.
Corner, Required or non-manifold
Definition at line 170 of file mmgcommon_private.h.
#define MG_SMSGN | ( | a, | |
b | |||
) | (((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.
#define MG_STR "&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&" |
Definition at line 54 of file mmgcommon_private.h.
#define MG_Tetra (1 << 2 ) |
4 local parameter applied over tetrahedron
Definition at line 162 of file mmgcommon_private.h.
#define MG_Tria (1 << 1 ) |
2 local parameter applied over triangle
Definition at line 161 of file mmgcommon_private.h.
Non parbdy boundary point (true bdy)
Definition at line 176 of file mmgcommon_private.h.
#define MG_Vert (1 << 0 ) |
1 local parameter applied over vertex
Definition at line 160 of file mmgcommon_private.h.
#define MG_VOK | ( | ppt | ) | (ppt && ((ppt)->tag < MG_NUL)) |
Vertex OK
Definition at line 166 of file mmgcommon_private.h.
#define MMG5_ADD_MEM | ( | mesh, | |
size, | |||
message, | |||
law | |||
) |
Increment memory counter memCur and check if we don't overflow the maximum authorizied memory memMax.
Definition at line 302 of file mmgcommon_private.h.
#define MMG5_ANGEDG 0.707106781186548 /*0.573576436351046 */ |
Definition at line 90 of file mmgcommon_private.h.
#define MMG5_ANGLIM -0.999999 |
Definition at line 91 of file mmgcommon_private.h.
#define MMG5_ATHIRD 0.333333333333333 |
Definition at line 92 of file mmgcommon_private.h.
#define MMG5_BITWIZE_MB_TO_B 20 |
Bitwise convertion from Mo to O
Definition at line 77 of file mmgcommon_private.h.
#define MMG5_BOXSIZE 500 |
size of box for renumbering with scotch.
Definition at line 73 of file mmgcommon_private.h.
#define MMG5_CHK_INT32_OVERFLOW | ( | wantedGap, | |
oldSiz, | |||
coef, | |||
shift, | |||
law | |||
) |
Check for int32 overflow when trying to reallocate coef*oldSiz+shift array
Definition at line 353 of file mmgcommon_private.h.
#define MMG5_CHK_MEM | ( | mesh, | |
size, | |||
string, | |||
law | |||
) |
Check if used memory overflow maximal authorized memory. Execute the command law if lack of memory.
Definition at line 215 of file mmgcommon_private.h.
#define MMG5_DEL_MEM | ( | mesh, | |
ptr | |||
) |
Free pointer ptr of mesh structure and compute the new used memory.
Definition at line 293 of file mmgcommon_private.h.
#define MMG5_DISPREF 10 |
Definition at line 84 of file mmgcommon_private.h.
#define MMG5_EPS 1.e-06 |
Definition at line 96 of file mmgcommon_private.h.
#define MMG5_EPSD 1.e-30 |
Definition at line 94 of file mmgcommon_private.h.
#define MMG5_EPSD2 1.0e-200 |
Definition at line 95 of file mmgcommon_private.h.
#define MMG5_EPSOK 1.e-15 |
Definition at line 97 of file mmgcommon_private.h.
#define MMG5_FEM 1 |
defaut value for FEM mode
Definition at line 136 of file mmgcommon_private.h.
#define MMG5_FILESTR_LGTH 128 /** Maximal length of a line in input file */ |
Definition at line 137 of file mmgcommon_private.h.
#define MMG5_GAP 0.2 |
gap value for reallocation
Definition at line 132 of file mmgcommon_private.h.
#define MMG5_HAUSD 0.01 |
default value for hausdorff param
Definition at line 121 of file mmgcommon_private.h.
#define MMG5_HGRAD 0.26236426446 |
default value for gradation (1.3)
Definition at line 122 of file mmgcommon_private.h.
#define MMG5_HGRADREQ 0.83290912294 |
default value for required gradation (2.3)
Definition at line 123 of file mmgcommon_private.h.
#define MMG5_HMAXCOE 2 |
coefficient to compute the default hmax value
Definition at line 134 of file mmgcommon_private.h.
#define MMG5_HMINCOE 0.001 |
coefficient to compute the default hmin value
Definition at line 133 of file mmgcommon_private.h.
#define MMG5_HMINMAXGAP 5 |
imposed gap between hmin and hmax if hmax<hmin
Definition at line 135 of file mmgcommon_private.h.
#define MMG5_INCREASE_MEM_MESSAGE | ( | ) |
Error message when lack of memory
Definition at line 432 of file mmgcommon_private.h.
#define MMG5_KA 7 |
Key for hash tables.
Definition at line 184 of file mmgcommon_private.h.
#define MMG5_KB 11 |
Key for hash tables.
Definition at line 185 of file mmgcommon_private.h.
#define MMG5_LAG -1 |
default value for lagrangian option
Definition at line 125 of file mmgcommon_private.h.
#define MMG5_LMAX 10240 |
Maximum array size when calling boulep.
Definition at line 64 of file mmgcommon_private.h.
#define MMG5_LPARMAX 200 |
Maximal number of local parameters
Definition at line 67 of file mmgcommon_private.h.
#define MMG5_LS 0.0 |
default level-set to discretize
Definition at line 127 of file mmgcommon_private.h.
#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.
#define MMG5_MEMMIN 38 |
minimal memory needed to store the mesh/sol names
Definition at line 111 of file mmgcommon_private.h.
#define MMG5_MEMPERCENT 0.5 |
Percent of RAM used by default
Definition at line 78 of file mmgcommon_private.h.
#define MMG5_MILLION 1048576 |
Definition at line 87 of file mmgcommon_private.h.
#define MMG5_NOHGRAD -1 |
disable gradation
Definition at line 124 of file mmgcommon_private.h.
#define MMG5_NONSET -1 |
Definition at line 120 of file mmgcommon_private.h.
#define MMG5_NONSET_HMAX -1 |
hmax value for unspecified hmax size
Definition at line 118 of file mmgcommon_private.h.
#define MMG5_NONSET_HMIN -1 |
hmin value for unspecified hmin size
Definition at line 117 of file mmgcommon_private.h.
#define MMG5_NONSET_HSIZ -1 |
hsiz value for unspecified hsiz map
Definition at line 119 of file mmgcommon_private.h.
#define MMG5_NONSET_MEM -1 |
mem value for unspecified max memory
Definition at line 116 of file mmgcommon_private.h.
#define MMG5_NR -1 |
no ridge detection
Definition at line 126 of file mmgcommon_private.h.
#define MMG5_NULKAL 1.e-30 |
Definition at line 98 of file mmgcommon_private.h.
#define MMG5_OFF 0 |
0
Definition at line 130 of file mmgcommon_private.h.
#define MMG5_ON 1 |
1
Definition at line 131 of file mmgcommon_private.h.
#define MMG5_PATHSEP '/' |
Definition at line 113 of file mmgcommon_private.h.
#define MMG5_PROCTREE 32 |
default size of the PROctree
Definition at line 129 of file mmgcommon_private.h.
#define MMG5_SAFE_CALLOC | ( | ptr, | |
size, | |||
type, | |||
law | |||
) |
Safe allocation with calloc
Definition at line 316 of file mmgcommon_private.h.
#define MMG5_SAFE_FREE | ( | ptr | ) |
Safe deallocation
Definition at line 309 of file mmgcommon_private.h.
#define MMG5_SAFE_MALLOC | ( | ptr, | |
size, | |||
type, | |||
law | |||
) |
Safe allocation with malloc
Definition at line 326 of file mmgcommon_private.h.
#define MMG5_SAFE_REALLOC | ( | ptr, | |
prevSize, | |||
newSize, | |||
type, | |||
message, | |||
law | |||
) |
Safe reallocation
Definition at line 337 of file mmgcommon_private.h.
#define MMG5_SAFE_RECALLOC | ( | ptr, | |
prevSize, | |||
newSize, | |||
type, | |||
message, | |||
law | |||
) |
safe reallocation with memset at 0 for the new values of tab
Definition at line 383 of file mmgcommon_private.h.
#define MMG5_SAFELL2ICAST | ( | longlongval | ) | (((longlongval) > (INT_MAX)) ? 0 : ((int)(longlongval))) |
Definition at line 439 of file mmgcommon_private.h.
#define MMG5_SAFELL2LCAST | ( | longlongval | ) | (((longlongval) > (LONG_MAX)) ? 0 : ((long)(longlongval))) |
Definition at line 438 of file mmgcommon_private.h.
#define MMG5_SD 8 |
Definition at line 189 of file mmgcommon_private.h.
#define MMG5_SQR32 0.866025403784439 |
Definition at line 100 of file mmgcommon_private.h.
#define MMG5_SW 4 |
Definition at line 188 of file mmgcommon_private.h.
#define MMG5_TAB_RECALLOC | ( | mesh, | |
ptr, | |||
initSize, | |||
wantedGap, | |||
type, | |||
message, | |||
law | |||
) |
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 406 of file mmgcommon_private.h.
#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.
#define MMG5_UNSET -1 |
Definition at line 81 of file mmgcommon_private.h.
#define MMG5_XREG 0.4 |
default relaxation parameter for coordinate regularization
Definition at line 128 of file mmgcommon_private.h.
#define MMG_FREAD | ( | ptr, | |
size, | |||
count, | |||
stream | |||
) |
check the return value of fread
Definition at line 442 of file mmgcommon_private.h.
#define MMG_FSCANF | ( | stream, | |
format, | |||
... | |||
) |
check the return value of fscanf
Definition at line 474 of file mmgcommon_private.h.
typedef struct MMG5_iNode_s MMG5_iNode |
typedef MMG5_Bezier* MMG5_pBezier |
Definition at line 617 of file mmgcommon_private.h.
double MMG2D_quickarea | ( | double | a[2], |
double | b[2], | ||
double | c[2] | ||
) |
void MMG5_2d3dUsage | ( | void | ) |
Print help for common options between 2D and 3D.
Definition at line 294 of file libtools.c.
int MMG5_2dSolTruncature_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the solution structure. |
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.
int MMG5_3dSolTruncature_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the solution structure. |
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.
int MMG5_Add_inode | ( | MMG5_pMesh | mesh, |
MMG5_iNode ** | liLi, | ||
int | val | ||
) |
mesh | pointer to the mesh structure (for count of used memory). |
liLi | pointer to the address of the root of the linked list. |
val | value to add to the linked list. |
Add a node with value val to a sorted linked list with unique entries.
Definition at line 68 of file apptools.c.
void MMG5_advancedUsage | ( | void | ) |
Print help for advanced users of mmg.
Definition at line 305 of file libtools.c.
void MMG5_bezierEdge | ( | MMG5_pMesh | , |
MMG5_int | , | ||
MMG5_int | , | ||
double * | , | ||
double * | , | ||
int8_t | , | ||
double * | |||
) |
int MMG5_boulec | ( | MMG5_pMesh | mesh, |
MMG5_int * | adjt, | ||
MMG5_int | start, | ||
int | ip, | ||
double * | tt | ||
) |
mesh | pointer to the mesh structure. |
adjt | pointer to the table of triangle adjacency. |
start | index of triangle where we start to work. |
ip | index of vertex where the tangent is computed. |
tt | pointer to the computed tangent. |
Compute the tangent to the curve at point ip.
Definition at line 199 of file boulep.c.
int MMG5_boulen | ( | MMG5_pMesh | mesh, |
MMG5_int * | adjt, | ||
MMG5_int | start, | ||
int | ip, | ||
double * | nn | ||
) |
mesh | pointer to the mesh structure. |
adjt | pointer to the table of triangle adjacency. |
start | index of triangle where we start to work. |
ip | local index of vertex where the normal is computed. |
nn | pointer to the computed tangent. |
Compute average normal of triangles sharing P without crossing ridge.
Definition at line 119 of file boulep.c.
int MMG5_boulep | ( | MMG5_pMesh | mesh, |
MMG5_int | start, | ||
int | ip, | ||
MMG5_int * | adja, | ||
MMG5_int * | list, | ||
MMG5_int * | tlist | ||
) |
mesh | pointer to mesh structure. |
start | a triangle to which ip belongs. |
ip | local point index |
adja | pointer to the adjacency array. |
list | pointer to the list of points connected to ip. |
tlist | pointer to the list of triangles sharing ip. |
Return all vertices connected to ip (with list[0] = ip) and all triangles sharing ip.
Definition at line 52 of file boulep.c.
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 | ||
) |
mesh | pointer to the mesh structure. |
adjt | pointer to the table of triangle adjacency. |
start | index of triangle where we start to work. |
ip | index of vertex on which we work. |
list | pointer to the computed list of GEO vertices incident to ip. |
listref | pointer to the corresponding edge references |
ng | pointer to the number of ridges. |
nr | pointer to the number of reference edges. |
lmax | maxmum size for the ball of the point 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.
int MMG5_boulet | ( | MMG5_pMesh | mesh, |
MMG5_int | start, | ||
int | ip, | ||
MMG5_int * | list, | ||
int8_t | s, | ||
int8_t * | opn | ||
) |
mesh | pointer to the mesh structure. |
start | index of triangle to start. |
ip | index of point for wich we compute the ball. |
list | pointer to the computed ball of ip. |
s | 1 if called from mmgs, 0 if called from mmg2d. |
opn | 0 for a closed ball, 1 for an open ball. |
Find all triangles sharing ip, \(list[0] =\) start do not stop when crossing ridge.
Definition at line 363 of file boulep.c.
int MMG5_boundingBox | ( | MMG5_pMesh | mesh | ) |
void MMG5_build3DMetric | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
MMG5_int | ip, | ||
double | dbuf[6] | ||
) |
mesh | pointer to the mesh structure |
sol | pointer to the sol structure. |
index | of point in which we want to build the metric |
dbuf | builded metric |
Build the metric at point ip depending with its type (ridge/not ridge).
Definition at line 1524 of file inout.c.
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] | ||
) |
double MMG5_caltri33_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pTria | pt | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the meric structure. |
pt | pointer to the triangle structure. |
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.
double MMG5_caltri_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pTria | ptt | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the meric structure. |
ptt | pointer to the triangle structure. |
Compute the quality of the surface triangle ptt with respect to an anisotropic metric.
Definition at line 115 of file quality.c.
|
inline |
mesh | pointer to the mesh structure. |
met | pointer to the meric structure. |
ptt | pointer to the triangle structure. |
Compute the quality of the surface triangle ptt with respect to an isotropic metric.
Definition at line 199 of file quality.c.
void MMG5_check_hminhmax | ( | MMG5_pMesh | mesh, |
int8_t | sethmin, | ||
int8_t | sethmax | ||
) |
mesh | pointer to the mesh structure. |
sethmin | 1 if hmin is setted by the user. |
sethmax | 1 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.
int MMG5_check_readedMesh | ( | MMG5_pMesh | mesh, |
MMG5_int | nref | ||
) |
mesh | pointer to an Mmg mesh |
nref | pointer to the number of negative refs (replaced by abolute values). |
Check the tetra orientation, print warning it negative refs have been detected, mark points as used.
Definition at line 526 of file inout.c.
int MMG5_check_setted_hminhmax | ( | MMG5_pMesh | mesh | ) |
int MMG5_chkmaniball | ( | MMG5_pMesh | mesh, |
MMG5_int | start, | ||
int8_t | istart | ||
) |
mesh | pointer to the mesh structure. |
start | index of starting tria. |
istart | local index of point that we check (in tria start) |
Check whether the ball of vertex i in tria start is manifold;
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 1304 of file mmg2.c.
int MMG5_chkmanimesh | ( | MMG5_pMesh | mesh | ) |
mesh | pointer to the mesh. |
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 1424 of file mmg2.c.
int MMG5_chkMetricType | ( | MMG5_pMesh | mesh, |
int * | type, | ||
int * | entities, | ||
FILE * | inm | ||
) |
mesh | pointer to the mesh structure. |
type | type of the metric |
entities | entities on which the metric applies (should be MMG5_Vertex) |
inm | metric file |
Check metric data:
Definition at line 2663 of file inout.c.
int MMG5_compute_meanMetricAtMarkedPoints_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
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.
int MMG5_compute_meanMetricAtMarkedPoints_iso | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
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.
int MMG5_countLocalParamAtTri | ( | MMG5_pMesh | mesh, |
MMG5_iNode ** | bdryRefs | ||
) |
mesh | pointer to the mesh structure. |
bdryRefs | pointer to the list of the boundary references. |
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.
void MMG5_crossprod3d | ( | double * | a, |
double * | b, | ||
double * | result | ||
) |
int MMG5_defsiz_startingMessage | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
const char * | funcname | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
funcname | name of the calling function |
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.
void MMG5_defUninitSize | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
int8_t | ismet | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
ismet | 1 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.
|
inline |
|
inline |
int MMG5_devangle | ( | double * | n1, |
double * | n2, | ||
double | crit | ||
) |
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.
mesh | pointer to the mesh structure |
disp | pointer to the displacement field |
lastt | 0 if a movement is possible, pointer to the last tested fraction otherwise |
shortmax | maximal parameter t (MMG2D_SHORTMAX or MMG3D_SHORTMAX) |
chkmovmesh | function that has to be called to check motion validity |
Generic function to compute the the largest fraction t that makes the motion along disp valid.
Definition at line 49 of file mmg3.c.
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 | ||
) |
mesh | pointer to the mesh structure. |
ned | edges number. |
avlen | pointer to the average edges lengths. |
amin | index of first extremity of the smallest edge. |
bmin | index of second extremity of the smallest edge. |
lmin | smallest edge length. |
amax | index of first extremity of the largest edge. |
bmax | index of second extremity of the largest edge. |
lmax | largest edge length. |
nullEdge | number of edges for which we are unable to compute the length |
bd | pointer to the table of the quality span. |
hl | pointer to the table that store the number of edges for eac |
shift | value to shift the target lenght interval span of quality |
Display histogram of edge length.
Definition at line 252 of file quality.c.
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 | ||
) |
ned | edges number. |
amin | index of first extremity of the smallest edge. |
bmin | index of second extremity of the smallest edge. |
lmin | smallest edge length. |
amax | index of first extremity of the largest edge. |
bmax | index of second extremity of the largest edge. |
lmax | largest edge length. |
nullEdge | number of edges for which we are unable to compute the length |
bd | pointer to the table of the quality span. |
hl | pointer to the table that store the number of edges for eac |
shift | value to shift the target lenght interval span of quality |
imprim | verbosity level |
Display histogram of edge length without the histo header
Definition at line 294 of file quality.c.
void MMG5_dotprod | ( | int8_t | dim, |
double * | a, | ||
double * | b, | ||
double * | result | ||
) |
int MMG5_eigenvmatnonsym2d | ( | MMG5_pMesh | mesh, |
double | m[], | ||
double | lambda[], | ||
double | v[][2] | ||
) |
mesh | pointer to the mesh structure. |
m | matrix array. |
lambda | eigenvalues array. |
v | double array of eigenvectors. |
Recompose a 2x2 non-symmetric matrix from its eigenvalue decomposition.
Definition at line 184 of file mettools.c.
int MMG5_eigenvmatnonsym3d | ( | MMG5_pMesh | mesh, |
double | m[], | ||
double | lambda[], | ||
double | v[][3] | ||
) |
mesh | pointer to the mesh structure. |
m | matrix array. |
lambda | eigenvalues array. |
v | double array of eigenvectors. |
Recompose a 3x3 non-symmetric matrix from its eigenvalue decomposition.
Definition at line 209 of file mettools.c.
int MMG5_eigenvmatsym2d | ( | MMG5_pMesh | mesh, |
double | m[], | ||
double | lambda[], | ||
double | v[][2] | ||
) |
mesh | pointer to the mesh structure. |
m | matrix array. |
lambda | eigenvalues array. |
v | double array of eigenvectors. |
Recompose a 2x2 symmetric matrix from its eigenvalue decomposition.
Definition at line 144 of file mettools.c.
int MMG5_eigenvmatsym3d | ( | MMG5_pMesh | mesh, |
double | m[], | ||
double | lambda[], | ||
double | v[][3] | ||
) |
mesh | pointer to the mesh structure. |
m | matrix array. |
lambda | eigenvalues array. |
v | double array of eigenvectors. |
Recompose a 3x3 symmetric matrix from its eigenvalue decomposition.
Definition at line 164 of file mettools.c.
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] | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
pt | pointer to the tria on which we integrate. |
p0 | pointer to the point that we want to move. |
pb | bezier patch of the triangle. |
r | rotation matrix that sends the normal at point p0 to e_z. |
gv | centre of mass that we want to update using the computed element weight. |
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.
|
inlinestatic |
sigid | signal number. |
Signal handling: specify error messages depending from catched signal.
Definition at line 505 of file mmgcommon_private.h.
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 | ||
) |
mesh | pointer to the mesh structure (for count of used memory). |
liLi | pointer to 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.
int MMG5_getStartRef | ( | MMG5_pMesh | mesh, |
MMG5_int | ref, | ||
MMG5_int * | pref | ||
) |
mesh | pointer to the mesh |
ref | final reference for which we are searching the initial one |
pref | pointer to the reference of the parent material. |
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.
MMG5_int MMG5_grad2metSurf | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pTria | pt, | ||
MMG5_int | np1, | ||
MMG5_int | np2 | ||
) |
mesh | pointer to the mesh. |
met | pointer to the metric structure. |
pt | pointer to a triangle. |
np1 | global index of the first extremity of the edge. |
np2 | global index of the second extremity of the edge. |
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.
Definition at line 975 of file anisosiz.c.
int MMG5_grad2metSurfreq | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pTria | pt, | ||
MMG5_int | npmaster, | ||
MMG5_int | npslave | ||
) |
mesh | pointer to the mesh. |
met | pointer to the metric structure. |
pt | pointer to the processed triangle. |
npmaster | edge extremity that cannot be modified |
npslave | edge extremity to modify to respect the gradation. |
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.
Definition at line 1951 of file anisosiz.c.
void MMG5_gradation_info | ( | MMG5_pMesh | mesh | ) |
void MMG5_gradEigenvreq | ( | double * | dm, |
double * | dn, | ||
double | difsiz, | ||
int8_t | dir, | ||
int8_t * | ier | ||
) |
dm | eigenvalues of the first matrix (not modified) |
dn | eigenvalues of the second matrix (modified) |
difsiz | maximal size gap authorized by the gradation. |
dir | direction in which the sizes are graded. |
ier | 2 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.
MMG5_int MMG5_gradsiz_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
int * | it | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
it | number of performed iteration (to fill) |
Standard gradation procedure.
Mark the edges belonging to a required entity
Definition at line 2259 of file anisosiz.c.
int MMG5_gradsiz_iso | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh |
met | pointer to the metric |
Isotropic mesh gradation routine. The points belonging to a required edge are treated in gradsizreq_iso.
Definition at line 278 of file isosiz.c.
int MMG5_gradsizreq_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
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.
int MMG5_gradsizreq_iso | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh |
met | pointer to the metric |
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.
int MMG5_hashEdge | ( | MMG5_pMesh | mesh, |
MMG5_Hash * | hash, | ||
MMG5_int | a, | ||
MMG5_int | b, | ||
MMG5_int | k | ||
) |
mesh | pointer to the mesh structure. |
hash | pointer to the hash table of edges. |
a | index of the first extremity of the edge. |
b | index of the second extremity of the edge. |
k | index of point along the edge. |
Add edge \([a;b]\) to the hash table.
Definition at line 346 of file hash.c.
int MMG5_hashEdgeTag | ( | MMG5_pMesh | mesh, |
MMG5_Hash * | hash, | ||
MMG5_int | a, | ||
MMG5_int | b, | ||
uint16_t | tag | ||
) |
mesh | pointer to the mesh structure. |
hash | pointer to the hash table of edges. |
a | index of the first extremity of the edge. |
b | index of the second extremity of the edge. |
tag | edge tag |
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 441 of file hash.c.
MMG5_int MMG5_hashFace | ( | MMG5_pMesh | mesh, |
MMG5_Hash * | hash, | ||
MMG5_int | ia, | ||
MMG5_int | ib, | ||
MMG5_int | ic, | ||
MMG5_int | k | ||
) |
mesh | pointer to the mesh. |
hash | pointer to the hash table to fill. |
ia | first vertex of face to hash. |
ib | second vertex of face to hash. |
ic | third vertex of face to hash. |
k | index of face to hash. |
Definition at line 286 of file hash.c.
MMG5_int MMG5_hashGet | ( | MMG5_Hash * | hash, |
MMG5_int | a, | ||
MMG5_int | b | ||
) |
hash | pointer to the hash table of edges. |
a | index of the first extremity of the edge. |
b | index of the second extremity of the edge. |
Find the index of point stored along \([a;b]\).
Definition at line 501 of file hash.c.
int MMG5_hashNew | ( | MMG5_pMesh | mesh, |
MMG5_Hash * | hash, | ||
MMG5_int | hsiz, | ||
MMG5_int | hmax | ||
) |
int MMG5_hashUpdate | ( | MMG5_Hash * | hash, |
MMG5_int | a, | ||
MMG5_int | b, | ||
MMG5_int | k | ||
) |
mesh | pointer to the mesh structure. |
hash | pointer to the hash table of edges. |
a | index of the first extremity of the edge. |
b | index of the second extremity of the edge. |
k | new index of point along the edge. |
Update the index of the point stored along the edge \([a;b]\)
Definition at line 404 of file hash.c.
int MMG5_interp_iso | ( | double * | ma, |
double * | mb, | ||
double * | mp, | ||
double | t | ||
) |
int MMG5_interpreg_ani | ( | MMG5_pMesh | , |
MMG5_pSol | , | ||
MMG5_pTria | , | ||
int8_t | , | ||
double | , | ||
double * | mr | ||
) |
int MMG5_intersecmet22 | ( | MMG5_pMesh | mesh, |
double * | m, | ||
double * | n, | ||
double * | mr | ||
) |
mesh | pointer to the mesh structure. |
m | pointer to a \((2x2)\) metric. |
n | pointer to a \((2x2)\) metric. |
mr | computed \((2x2)\) metric. |
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.
int MMG5_intmetsavedir | ( | MMG5_pMesh | mesh, |
double * | m, | ||
double * | n, | ||
double * | mr | ||
) |
mesh | pointer to the mesh structure. |
m | pointer to the first metric to intersect. |
n | pointer to the second metric to intersect. |
mr | pointer to the computed intersected metric. |
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.
int MMG5_intridmet | ( | MMG5_pMesh | , |
MMG5_pSol | , | ||
MMG5_int | , | ||
MMG5_int | , | ||
double | , | ||
double * | , | ||
double * | |||
) |
int MMG5_invmat | ( | double * | m, |
double * | mi | ||
) |
int MMG5_invmat22 | ( | double | m[2][2], |
double | mi[2][2] | ||
) |
int MMG5_invmat33 | ( | double | m[3][3], |
double | mi[3][3] | ||
) |
int MMG5_invmatg | ( | double | m[9], |
double | mi[9] | ||
) |
int MMG5_isLevelSet | ( | MMG5_pMesh | mesh, |
MMG5_int | ref0, | ||
MMG5_int | ref1 | ||
) |
mesh | pointer to the mesh structure. |
ref0 | reference of the first tetrahedron sharing the face. |
ref1 | reference of the second tetrahedron sharing the face.. |
Identify whether a face is on the discrete level set or not.
Definition at line 472 of file mmg2.c.
int MMG5_ismaniball | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
MMG5_int | start, | ||
int8_t | istart | ||
) |
mesh | pointer to the mesh structure. |
sol | pointer to the level-set values. |
start | index of the starting tria |
istart | local index (inside the tria start) of the vertex that we check. |
Check whether snapping the value of vertex istart of start to 0 exactly leads to a non manifold situation.
Definition at line 623 of file mmg2.c.
int MMG5_isNotSplit | ( | MMG5_pMesh | mesh, |
MMG5_int | ref | ||
) |
mesh | pointer to the mesh structure. |
ref | initial reference. |
Identify whether an entity with reference ref should not be split.
Definition at line 448 of file mmg2.c.
int MMG5_isSplit | ( | MMG5_pMesh | mesh, |
MMG5_int | ref, | ||
MMG5_int * | refint, | ||
MMG5_int * | refext | ||
) |
mesh | pointer to the mesh structure. |
ref | initial reference. |
refint | internal reference after ls discretization. |
refext | external reference after ls discretization. |
Identify whether an entity with reference ref should be split, and the labels of the resulting entities.
Definition at line 412 of file mmg2.c.
void MMG5_keep_subdomainElts | ( | MMG5_pMesh | mesh, |
int | nsd, | ||
int(*)(MMG5_pMesh, MMG5_int) | delElt | ||
) |
void MMG5_lagUsage | ( | void | ) |
Print help for lagrangian motion option.
Definition at line 278 of file libtools.c.
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 | ||
) |
mesh | pointer to the mesh |
filename | pointer to the name of file |
inm | pointer to the file pointer |
posNodes | pointer to the position of nodes data in file |
posElts | pointer to the position of elts data in file |
posNodeData | pointer to the list of the positions of data in file |
bin | 1 if binary format |
nelts | number of elements in file |
nsol | number of data in file |
Begin to read mesh at MSH file format. Read the mesh size informations.
Definition at line 278 of file inout.c.
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 | ||
) |
mesh | pointer to the mesh |
sol | pointer to the solutions array |
inm | pointer to the file pointer |
posNodes | position of nodes data in file |
posElts | position of elts data in file |
posNodeData | position of solution data in file |
bin | 1 if binary format |
nelts | number of elements in file |
nsols | number of silutions in file |
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
Definition at line 664 of file inout.c.
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 | ||
) |
filename | name of file. |
meshDim | mesh dimenson. |
inm | allocatable pointer to the FILE structure |
ver | file version (1=simple precision, 2=double) |
bin | 1 if the file is a binary |
iswp | 1 or 0 depending on the endianness (binary only) |
np | number of solutions of each type |
dim | solution dimension |
nsols | number of solutions of different types in the file |
type | type of solutions |
posnp | pointer to the position of the point list in the file |
imprim | verbosity |
Open the "filename" solution file and read the file header.
Definition at line 2076 of file inout.c.
int MMG5_loadVtuMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
void MMG5_mark_pointsOnReqEdge_fromTria | ( | MMG5_pMesh | mesh | ) |
void MMG5_mark_usedVertices | ( | MMG5_pMesh | mesh, |
void(*)(MMG5_pMesh, MMG5_int) | delPt | ||
) |
void MMG5_mark_verticesAsUnused | ( | MMG5_pMesh | mesh | ) |
void MMG5_memOption_memSet | ( | MMG5_pMesh | mesh | ) |
size_t MMG5_memSize | ( | void | ) |
int MMG5_minQualCheck | ( | MMG5_int | iel, |
double | minqual, | ||
double | alpha | ||
) |
iel | index of the worst tetra of the mesh |
minqual | quality of the worst tetra of the mesh (will be normalized by alpha) |
alpha | normalisation parameter for the quality |
Print warning or error messages depending on the quality of the worst tetra of the mesh.
Definition at line 343 of file quality.c.
void MMG5_mmgDefaultValues | ( | MMG5_pMesh | mesh | ) |
mesh | pointer to the mesh structure. |
Print the default parameters values.
Definition at line 70 of file libtools.c.
int MMG5_mmgHashTria | ( | MMG5_pMesh | mesh, |
MMG5_int * | adjt, | ||
MMG5_Hash * | hash, | ||
int | chkISO | ||
) |
mesh | pointer to the mesh structure. |
adjt | pointer to the adjacency table of the surfacic mesh. |
hash | pointer to the edge hash table. |
chkISO | flag to say if we check ISO references (so if we come from mmg3d). |
Create surface adjacency
Definition at line 58 of file hash.c.
void MMG5_mmgInit_parameters | ( | MMG5_pMesh | mesh | ) |
int MMG5_mmgIntextmet | ( | MMG5_pMesh | , |
MMG5_pSol | , | ||
MMG5_int | , | ||
double * | , | ||
double * | |||
) |
int MMG5_mmgIntmet33_ani | ( | double * | m, |
double * | n, | ||
double * | mr, | ||
double | s | ||
) |
m | input metric. |
n | input metric. |
mr | computed output metric. |
s | parameter coordinate for the interpolation of metrics m and n. |
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.
void MMG5_mmgUsage | ( | char * | prog | ) |
*prog | pointer to the program name. |
Print help for common options of the 3 codes (first section).
Definition at line 212 of file libtools.c.
void MMG5_mn | ( | double | m[6], |
double | n[6], | ||
double | mn[9] | ||
) |
int MMG5_MultiMat_init | ( | MMG5_pMesh | mesh | ) |
mesh | pointer to the mesh structure. |
Initialize handling of multimaterial mode.
An indexed table of materials has been provided by the MMG5_Mat datatype in the form:
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
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):
... | ... 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?
Why child materials should be different?
Definition at line 326 of file mmg2.c.
|
inline |
|
inline |
mesh | pointer to the mesh stucture. |
ip1 | first point of face. |
ip2 | second point of face. |
ip3 | third point of face. |
n | pointer to store the computed normal. |
Compute non-normalized face normal given three points on the surface.
Definition at line 127 of file tools.c.
|
inline |
mesh | pointer to the mesh stucture. |
ip1 | first point of face. |
ip2 | second point of face. |
ip3 | third point of face. |
n | pointer to store the computed normal. |
Compute normalized face normal given three points on the surface.
Definition at line 183 of file tools.c.
|
inline |
|
inline |
n | array size |
shift | shift to apply when taking array value |
stride | stride to apply when taking array value |
val | array of double precision floating points |
oldval | array to store input values |
perm | permutation array |
Naively permute a small array. Use shift and stride to eventually permute matrix columns.
Definition at line 80 of file tools.c.
|
inline |
n | array size |
val | array of double precision floating points |
perm | permutation array |
naive (increasing) sorting algorithm, for very small tabs ; permutation is stored in perm
Definition at line 49 of file tools.c.
|
inline |
void MMG5_paramUsage1 | ( | void | ) |
Print help for common parameters options of the 3 codes (first section).
Definition at line 242 of file libtools.c.
void MMG5_paramUsage2 | ( | void | ) |
Print help for common options of the 3 codes (second section).
Definition at line 261 of file libtools.c.
int MMG5_paratmet | ( | double | c0[3], |
double | n0[3], | ||
double | m[6], | ||
double | c1[3], | ||
double | n1[3], | ||
double | mt[6] | ||
) |
c0 | table of the coordinates of the starting point. |
n0 | normal at the starting point. |
m | metric to be transported. |
c1 | table of the coordinates of the ending point. |
n1 | normal at the ending point. |
mt | computed metric. |
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.
void MMG5_printMetStats | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
void MMG5_printSolStats | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol | ||
) |
void MMG5_printTria | ( | MMG5_pMesh | mesh, |
char * | fileName | ||
) |
int MMG5_readDoubleSol3D | ( | MMG5_pSol | sol, |
FILE * | inm, | ||
int | bin, | ||
int | iswp, | ||
MMG5_int | pos | ||
) |
sol | pointer to an allocatable sol structure. |
inm | pointer to the solution file |
bin | 1 if binary file |
iswp | Endianess |
index | of the readed solution |
Read the solution value for vertex of index pos in double precision.
Definition at line 2271 of file inout.c.
int MMG5_readFloatSol3D | ( | MMG5_pSol | sol, |
FILE * | inm, | ||
int | bin, | ||
int | iswp, | ||
int | pos | ||
) |
sol | pointer to an allocatable sol structure. |
inm | pointer to the solution file |
bin | 1 if binary file |
iswp | Endianess |
index | of the readed solution |
Read the solution value for vertex of index pos in floating precision.
Definition at line 2222 of file inout.c.
int MMG5_regnor | ( | MMG5_pMesh | mesh | ) |
int MMG5_reset_metricAtReqEdges_surf | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
int8_t | ismet | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
ismet | 1 if user provided metric |
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.
int MMG5_resetRef_ls | ( | MMG5_pMesh | mesh | ) |
int MMG5_resetRef_lssurf | ( | MMG5_pMesh | mesh | ) |
double MMG5_ridSizeInNormalDir | ( | MMG5_pMesh | , |
int | , | ||
double * | , | ||
MMG5_pBezier | , | ||
double | , | ||
double | |||
) |
double MMG5_ridSizeInTangentDir | ( | MMG5_pMesh | mesh, |
MMG5_pPoint | p0, | ||
MMG5_int | idp, | ||
MMG5_int * | iprid, | ||
double | isqhmin, | ||
double | isqhmax | ||
) |
mesh | pointer to the mesh structure. |
p0 | pointer to the point at which we define the metric. |
idp | global index of the point at which we define the metric. |
iprid | pointer to the two extremities of the ridge. |
isqhmin | minimum edge size. |
isqhmax | maximum edge size. |
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.
int MMG5_rmc | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
mesh | pointer to the mesh |
sol | pointer to the level-set |
Removal of small parasitic components (bubbles of material, etc) with volume less than mesh->info.rmc * volume of the mesh.
Definition at line 912 of file mmg2.c.
|
inline |
|
inline |
int MMG5_saveDisp | ( | MMG5_pMesh | mesh, |
MMG5_pSol | disp | ||
) |
int MMG5_saveMshMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename, | ||
int | metricData | ||
) |
mesh | pointer to the mesh structure. |
sol | pointer to an array of solutions. |
filename | name of file. |
metricData | 1 if the data saved is a metric (if only 1 data) |
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 1583 of file inout.c.
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 | ||
) |
mesh | pointer to the mesh structure. |
inm | pointer to the opened file unit. |
ver | file version (1=simple precision, 2=double). |
bin | 1 if the file is a binary. |
bpos | cumulative field position for binary Medit format. |
nsols | number of solutions of different types in the file. |
nsolsAtTetra | number of solutions at tetra in the file. |
entities | kind of entity on which the solution applies. |
type | type of solutions. |
size | size of solutions. |
Save the number and type of solutions at Tetrahedron (not used by Mmg).
Definition at line 2595 of file inout.c.
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 | ||
) |
mesh | pointer to the mesh structure. |
inm | pointer to the opened file unit. |
ver | file version (1=simple precision, 2=double). |
bin | 1 if the file is a binary. |
bpos | cumulative field position for binary Medit format. |
nsols | number of solutions of different types in the file. |
nsolsAtTriangles | number of solutions at triangles in the file. |
entities | kind of entity on which the solution applies. |
type | type of solutions. |
size | size of solutions. |
Save the number and type of solutions at Triangles (not used by Mmg).
Definition at line 2526 of file inout.c.
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 | ||
) |
mesh | pointer to the mesh structure. |
filename | name of file. |
inm | allocatable pointer to the FILE structure. |
ver | file version (1=simple precision, 2=double). |
bin | 1 if the file is a binary. |
bpos | cumulative field position for binary Medit format. |
np | number of solutions of each type. |
dim | solution dimension. |
nsols | number of solutions of different types in the file. |
entities | kind of entity on which the solution applies (Vertex or Tetra) |
type | type of solutions. |
size | size of solutions. |
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 2389 of file inout.c.
int MMG5_scale_meshAndSol | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pSol | sol, | ||
double * | dd | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to a metric |
sol | pointer to a solution structure (level-set or displacement). |
dd | pointer to the scaling value (to fill) |
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.
int MMG5_scale_scalarMetric | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
double | dd | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
dd | scaling value. |
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.
int MMG5_scotchCall | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pSol | fields, | ||
MMG5_int * | permNodGlob | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the solution structure. |
fields | pointer to an array of solution fields (non mandatory) |
permNodGlob | store the global permutation of nodes (if provided). |
Call scotch renumbering.
Definition at line 231 of file librnbg.c.
void MMG5_Set_commonFunc | ( | void | ) |
int MMG5_setref_ls | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
int MMG5_setref_lssurf | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
int MMG5_simred2d | ( | MMG5_pMesh | mesh, |
double * | m, | ||
double * | n, | ||
double | dm[2], | ||
double | dn[2], | ||
double | vp[2][2] | ||
) |
mesh | pointer to the mesh |
m | first matrix |
n | second matrix |
dm | eigenvalues of m in the coreduction basis (to fill) |
dn | eigenvalues of n in the coreduction basis (to fill) |
vp | coreduction basis (to fill) |
Perform simultaneous reduction of matrices m and n.
Definition at line 1326 of file anisosiz.c.
int MMG5_simred3d | ( | MMG5_pMesh | mesh, |
double * | m, | ||
double * | n, | ||
double | dm[3], | ||
double | dn[3], | ||
double | vp[3][3] | ||
) |
mesh | pointer to the mesh |
m | first matrix |
n | second matrix |
dm | eigenvalues of m in the coreduction basis (to fill) |
dn | eigenvalues of n in the coreduction basis (to fill) |
vp | coreduction basis (to fill) |
Perform simultaneous reduction of matrices m and n.
Definition at line 1416 of file anisosiz.c.
int MMG5_snpval_ls | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
mesh | pointer to the mesh structure. |
sol | pointer to the level-set function. |
Snap values of the level set function very close to 0 to exactly 0, and prevent nonmanifold patterns from being generated.
Definition at line 503 of file mmg2.c.
int MMG5_snpval_lssurf | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
mesh | pointer to the mesh structure. |
sol | pointer to the level-set |
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.
int MMG5_solTruncature_iso | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the solution structure. |
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.
int MMG5_solveDefmetrefSys | ( | MMG5_pMesh | , |
MMG5_pPoint | , | ||
MMG5_int * | , | ||
double | r[3][3], | ||
double * | , | ||
double * | , | ||
double * | , | ||
double * | , | ||
double | , | ||
double | , | ||
double | |||
) |
int MMG5_solveDefmetregSys | ( | MMG5_pMesh | , |
double | r[3][3], | ||
double * | , | ||
double * | , | ||
double * | , | ||
double * | , | ||
double | , | ||
double | , | ||
double | |||
) |
int MMG5_sum_reqEdgeLengthsAtPoint | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_int | ip0, | ||
MMG5_int | ip1 | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
ip0 | index of the first edge extremity |
ip1 | index of the second edge extremity |
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.
double MMG5_surftri33_ani | ( | MMG5_pMesh | , |
MMG5_pTria | , | ||
double * | , | ||
double * | , | ||
double * | |||
) |
double MMG5_surftri_ani | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pTria | ptt | ||
) |
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
ptt | pointer to the triangle structure. |
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.
double MMG5_surftri_iso | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pTria | ptt | ||
) |
int MMG5_swapbin | ( | int | sbin | ) |
MMG5_int MMG5_swapbin_int | ( | MMG5_int | sbin | ) |
double MMG5_swapd | ( | double | sbin | ) |
float MMG5_swapf | ( | float | sbin | ) |
|
inline |
int MMG5_test_crossprod3d | ( | ) |
int MMG5_test_dotprod | ( | ) |
int MMG5_test_eigenvmatnonsym2d | ( | MMG5_pMesh | mesh, |
double * | mex, | ||
double | lambdaex[], | ||
double | vpex[][2], | ||
double | ivpex[][2] | ||
) |
mesh | pointer to the mesh structure. |
mex | test matrix array. |
lambdaex | exact eigenvalues array. |
vpex | double array of exact right eigenvectors. |
ivpex | double array of exact right eigenvectors inverse. |
For a 2x2 non-symmetric matrix, Test:
Recompose matrix from its eigendecomposition
Compute eigendecomposition
Compute both eigendecomposition and recomposition, and check matrix
Definition at line 379 of file mettools.c.
int MMG5_test_eigenvmatnonsym3d | ( | MMG5_pMesh | mesh, |
double * | mex, | ||
double | lambdaex[], | ||
double | vpex[][3], | ||
double | ivpex[][3] | ||
) |
mesh | pointer to the mesh structure. |
mex | test matrix array. |
lambdaex | exact eigenvalues array. |
vpex | double array of exact right eigenvectors. |
ivpex | double array of exact right eigenvectors inverse. |
For a 3x3 non-symmetric matrix, Test:
Recompose matrix from its eigendecomposition
Compute eigendecomposition
Compute both eigendecomposition and recomposition, and check matrix
Definition at line 536 of file mettools.c.
int MMG5_test_eigenvmatsym2d | ( | MMG5_pMesh | mesh, |
double * | mex, | ||
double | lambdaex[], | ||
double | vpex[][2] | ||
) |
mesh | pointer to the mesh structure. |
mex | test matrix array. |
lambdaex | exact eigenvalues array. |
vpex | double array of exact eigenvectors. |
For a 2x2 symmetric matrix, Test:
Recompose matrix from its eigendecomposition
Compute eigendecomposition
Compute both eigendecomposition and recomposition, and check matrix
Definition at line 301 of file mettools.c.
int MMG5_test_eigenvmatsym3d | ( | MMG5_pMesh | mesh, |
double * | mex, | ||
double | lambdaex[], | ||
double | vpex[][3] | ||
) |
mesh | pointer to the mesh structure. |
mex | test matrix array. |
lambdaex | exact eigenvalues array. |
vpex | double array of exact eigenvectors. |
For a 3x3 symmetric matrix, Test:
Recompose matrix from its eigendecomposition
Compute eigendecomposition
Compute both eigendecomposition and recomposition, and check matrix
Definition at line 456 of file mettools.c.
int MMG5_test_intersecmet22 | ( | MMG5_pMesh | mesh | ) |
mesh | pointer to the mesh structure. |
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.
int MMG5_test_intersecmet33 | ( | MMG5_pMesh | mesh | ) |
mesh | pointer to the mesh structure. |
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.
int MMG5_test_invmat22 | ( | ) |
int MMG5_test_invmat33 | ( | ) |
|
inline |
int MMG5_test_mn | ( | ) |
|
inline |
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.
int MMG5_test_simred2d | ( | MMG5_pMesh | mesh, |
double * | mex, | ||
double * | nex, | ||
double * | dmex, | ||
double * | dnex, | ||
double | vpex[][2] | ||
) |
mesh | pointer to the mesh structure |
mex | first symmetric test matrix |
nex | second symmetric test matrix |
dm | diagonalization of the first matrix on the reduction basis |
dn | diagonalization of the second matrix on the reduction basis |
vp | simultaneous reduction basis (stored by columns) |
For a couple of 2x2 symmetric matrices, Test:
Compute simultaneous reduction
Recompose matrices from diagonal values
Definition at line 1585 of file anisosiz.c.
int MMG5_test_simred3d | ( | MMG5_pMesh | mesh, |
double * | mex, | ||
double * | nex, | ||
double * | dmex, | ||
double * | dnex, | ||
double | vpex[][3] | ||
) |
mesh | pointer to the mesh structure |
mex | first symmetric test matrix |
nex | second symmetric test matrix |
dm | diagonalization of the first matrix on the reduction basis |
dn | diagonalization of the second matrix on the reduction basis |
vp | simultaneous reduction basis (stored by columns) |
For a couple of 3x3 symmetric matrices, Test:
Compute simultaneous reduction
Recompose matrices from diagonal values
Definition at line 1663 of file anisosiz.c.
int MMG5_test_transpose3d | ( | ) |
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.
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.
void MMG5_transpose3d | ( | double | m[3][3] | ) |
int MMG5_truncate_met3d | ( | MMG5_pSol | met, |
MMG5_int | ip, | ||
double | isqhmin, | ||
double | isqhmax | ||
) |
met | pointer to metric. |
ip | pointer to global index of point on which metric has to be truncated |
isqhmin | inverse square of hmin (min edge size) |
isqhmax | inverse square of hmax (max edge size) |
Truncation of anisotropic metric at point ip by hmin and hmax.
Definition at line 148 of file scalem.c.
|
inline |
m | first matrix |
n | second matrix |
dm | eigenvalues of m in the coreduction basis |
dn | eigenvalues of n in the coreduction basis |
vp | coreduction basis |
ier | flag of the updated sizes: (ier & 1) if we dm has been modified, (ier & 2) if dn has been modified. |
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.
|
inline |
m | first matrix |
n | second matrix |
dm | eigenvalues of m in the coreduction basis |
dn | eigenvalues of n in the coreduction basis |
vp | coreduction basis |
ier | flag of the updated sizes: (ier & 1) if we dm has been modified, (ier & 2) if dn has been modified. |
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.
int MMG5_updatemetreq_ani | ( | double * | n, |
double | dn[2], | ||
double | vp[2][2] | ||
) |
n | matrix to update |
dn | eigenvalues of n in the coreduction basis |
vp | coreduction basis |
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.
void MMG5_version | ( | MMG5_pMesh | mesh, |
char * | dim | ||
) |
Functions needed by libraries API.
mesh | pointer to the mesh |
dim | string dontaining the dimension (3D,2D or S) |
Print MMG release and date
Definition at line 43 of file libtools.c.
|
inlinestatic |
Inlined functions for libraries and executables Warn user that we overflow asked memory during scotch call
Definition at line 489 of file mmgcommon_private.h.
void MMG5_writeDoubleSol3D | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
FILE * | inm, | ||
int | bin, | ||
MMG5_int | pos, | ||
int | metricData | ||
) |
mesh | pointer to the mesh structure |
sol | pointer to an allocatable sol structure. |
inm | pointer to the solution file |
bin | 1 if binary file |
pos | of the writted solution |
metricData | 1 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 2320 of file inout.c.
int MMG5_writeLocalParamAtTri | ( | MMG5_pMesh | mesh, |
MMG5_iNode * | bdryRefs, | ||
FILE * | out | ||
) |
mesh | pointer to the mesh structure. |
bdryRefs | pointer to the list of the boundary references. |
out | pointer to the file in which to write. |
Write the local default values at triangles in the parameter file.
Definition at line 186 of file apptools.c.
|
inlinestatic |
Definition at line 228 of file mmgcommon_private.h.
|
inlinestatic |
Definition at line 279 of file mmgcommon_private.h.
|
inlinestatic |
|
inlinestatic |
|
static |
next vertex of triangle: {1,2,0}
Definition at line 585 of file mmgcommon_private.h.
|
static |
previous vertex of triangle: {2,0,1}
Definition at line 586 of file mmgcommon_private.h.