Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
Functions
libmmgs_tools.c File Reference

Tools functions for the mmgs library. More...

#include "libmmgs.h"
#include "libmmgs_private.h"
#include "inlined_functions_private.h"
#include "mmgsexterns_private.h"
#include "mmgexterns_private.h"
Include dependency graph for libmmgs_tools.c:

Go to the source code of this file.

Functions

void MMGS_setfunc (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMGS_usage (char *prog)
 
int MMGS_defaultValues (MMG5_pMesh mesh)
 
int MMGS_parsar (int argc, char *argv[], MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol)
 
int MMGS_freeLocalPar (MMG5_pMesh mesh)
 
int MMGS_stockOptions (MMG5_pMesh mesh, MMG5_Info *info)
 
void MMGS_destockOptions (MMG5_pMesh mesh, MMG5_Info *info)
 
int MMGS_Get_numberOfNonBdyEdges (MMG5_pMesh mesh, MMG5_int *nb_edges)
 
int MMGS_Get_nonBdyEdge (MMG5_pMesh mesh, MMG5_int *e0, MMG5_int *e1, MMG5_int *ref, MMG5_int idx)
 
int MMGS_Get_adjaTri (MMG5_pMesh mesh, MMG5_int kel, MMG5_int listri[3])
 Return adjacent elements of a triangle. More...
 
int MMGS_Get_adjaVerticesFast (MMG5_pMesh mesh, MMG5_int ip, MMG5_int start, MMG5_int lispoi[MMGS_LMAX])
 Return adjacent elements of a triangle. More...
 
static int MMGS_solTruncatureForOptim (MMG5_pMesh mesh, MMG5_pSol met, int ani)
 
int MMGS_doSol_iso (MMG5_pMesh mesh, MMG5_pSol met)
 
static int MMGS_unitTensor_3D (MMG5_pMesh mesh, MMG5_int k, int i, MMG5_pPoint p1, double *m)
 
static int MMGS_surfopenballRotation (MMG5_pMesh mesh, MMG5_pPoint p0, MMG5_int k, int i, int ilist, double r[3][3], double *lispoi, double n[3])
 
static int MMGS_unitTensor_2D (MMG5_pMesh mesh, MMG5_int k, int i, MMG5_pPoint p1, double *m, double isqhmax)
 
int MMGS_doSol_ani (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMGS_Set_constantSize (MMG5_pMesh mesh, MMG5_pSol met)
 
int MMGS_Compute_eigenv (double m[6], double lambda[3], double vp[3][3])
 
void MMGS_Free_solutions (MMG5_pMesh mesh, MMG5_pSol sol)
 
int MMGS_Clean_isoSurf (MMG5_pMesh mesh)
 

Detailed Description

Tools functions for the mmgs library.

Author
Charles Dapogny (UPMC)
Cécile Dobrzynski (Bx INP/Inria/UBordeaux)
Pascal Frey (UPMC)
Algiane Froehly (Inria/UBordeaux)
Version
5
Todo:
Doxygen documentation

Definition in file libmmgs_tools.c.

Function Documentation

◆ MMGS_Clean_isoSurf()

int MMGS_Clean_isoSurf ( MMG5_pMesh  mesh)
Parameters
meshpointer toward mesh sructure
Returns
1 if successful, 0 otherwise.

Clean data (triangles and edges) linked to isosurface.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_CLEAN_ISOSURF(mesh,retval)
MMG5_DATA_PTR_T, INTENT(IN) :: mesh
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 1586 of file libmmgs_tools.c.

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

◆ MMGS_Compute_eigenv()

int MMGS_Compute_eigenv ( double  m[6],
double  lambda[3],
double  vp[3][3] 
)
Parameters
mupper part of a symetric matric diagonalizable in |R
lambdaarray of the metric eigenvalues
vparray of the metric eigenvectors
Returns
the order of the eigenvalues

Compute the real eigenvalues and eigenvectors of a symetric matrice m whose upper part is provided (m11, m12, m13, m22, m23, m33 in this order). lambda[0] is the eigenvalue associated to the eigenvector ( v[0][0], v[0,1], v[0,2] ) in C and to the eigenvector v(1,:) in fortran lambda[1] is the eigenvalue associated to the eigenvector ( v[1][0], v[1,1], v[1,2] ) in C and to the eigenvector v(2,:) in fortran lambda[2] is the eigenvalue associated to the eigenvector ( v[2][0], v[2,1], v[2,2] ) in C and to the eigenvector v(3,:) in fortran

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_COMPUTE_EIGENV(m,lambda,vp,retval)
REAL(KIND=8), INTENT(IN) :: m(*)
REAL(KIND=8), INTENT(OUT) :: lambda(*),vp(*)
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 1553 of file libmmgs_tools.c.

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

◆ MMGS_defaultValues()

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

Print the default parameters values.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_DEFAULTVALUES(mesh,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 111 of file libmmgs_tools.c.

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

◆ MMGS_destockOptions()

void MMGS_destockOptions ( MMG5_pMesh  mesh,
MMG5_Info info 
)
Parameters
meshpointer toward the mesh structure.
infopointer toward the info structure.

Recover the info structure stored in the mesh structure.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_DESTOCKOPTIONS(mesh,info)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,info
END SUBROUTINE

Definition at line 511 of file libmmgs_tools.c.

Here is the caller graph for this function:

◆ MMGS_doSol_ani()

int MMGS_doSol_ani ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh
metpointer toward the metric
Returns
1 if succeed, 0 if fail
Remarks
need the normal at vertices.
hmax/hmin may be not setted here so we can't use it.

Compute the unit metric tensor at mesh vertices (from edges passing through the vertices).

Increment base marker to detect points already processed

Step 1: treat non-manifold points: travel triangle edges and add edge contribution to extremities (we have to do that because ball of non-manifold points can't be computed)

Corner point (no normal at vertex): if the corner defines an angle we can compute the 3D unit tensor (non singular), in the other cases (corner along a flat surface, for example because it is at the intersection of 3 specific edges or because it is provided by the user) we will project the edges onto the tangent plane and compute the 2D unit tensor (as for required and regular points)

Required point (no normal at vertex): projection of edges onto the tangent plane and computation of the 2D unit tensor

Ridge point (2 normals): normally we can compute the 3D unit tensor. If the angle is too flat, the computation fails and we try to compute the 2D unit tensor on the tangent plane of the point.

Regular or reference point: projection of edges onto the tangent plane and computation of the 2D unit metric tensor

Step 3: check computed metric

Step 4: computation of hmin/hmax (if needed) and metric truncature

Definition at line 1280 of file libmmgs_tools.c.

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

◆ MMGS_doSol_iso()

int MMGS_doSol_iso ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh
metpointer toward the metric
Returns
1 if succeed, 0 if fail
Remarks
need the normal at vertices.

Compute isotropic size map according to the mean of the length of the edges passing through a point.

Definition at line 805 of file libmmgs_tools.c.

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

◆ MMGS_Free_solutions()

void MMGS_Free_solutions ( MMG5_pMesh  mesh,
MMG5_pSol  sol 
)
Parameters
meshpointer toward the mesh structure
solpointer toward the solution structure

Free the solution.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_FREE_SOLUTIONS(mesh,sol)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
END SUBROUTINE

Definition at line 1559 of file libmmgs_tools.c.

Here is the caller graph for this function:

◆ MMGS_freeLocalPar()

int MMGS_freeLocalPar ( MMG5_pMesh  mesh)

Definition at line 490 of file libmmgs_tools.c.

◆ MMGS_Get_adjaTri()

int MMGS_Get_adjaTri ( MMG5_pMesh  mesh,
MMG5_int  kel,
MMG5_int  listri[3] 
)

Return adjacent elements of a triangle.

Parameters
meshpointer toward the mesh structure.
keltriangle index.
listripointer toward the table of the indices of the three adjacent triangles of the elt kel (the index is 0 if there is no adjacent).
Returns
1.

Find the indices of the 3 adjacent elements of triangle kel. $v_i = 0$ if the $i^{th}$ face has no adjacent element (so we are on a boundary face).

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_GET_ADJATRI(mesh,kel,listri,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN) :: kel
INTEGER(MMG5F_INT), DIMENSION(3), INTENT(OUT) :: listri
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 655 of file libmmgs_tools.c.

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

◆ MMGS_Get_adjaVerticesFast()

int MMGS_Get_adjaVerticesFast ( MMG5_pMesh  mesh,
MMG5_int  ip,
MMG5_int  start,
MMG5_int  lispoi[MMGS_LMAX] 
)

Return adjacent elements of a triangle.

Parameters
meshpointer toward the mesh structure.
ipvertex index.
startindex of a triangle holding ip.
lispoipointer toward an array of size MMGS_LMAX that will contain the indices of adjacent vertices to the vertex ip.
Returns
nbpoi the number of adjacent points if success, 0 if fail.

Find the indices of the adjacent vertices of the vertex ip of the triangle start.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_GET_ADJAVERTICESFAST(mesh,ip,start,lispoi,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN) :: ip,start
INTEGER(MMG5F_INT), DIMENSION(MMGS_LMAX), INTENT(OUT) :: lispoi
INTEGER(MMG5F_INT), INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 669 of file libmmgs_tools.c.

Here is the caller graph for this function:

◆ MMGS_Get_nonBdyEdge()

int MMGS_Get_nonBdyEdge ( MMG5_pMesh  mesh,
MMG5_int *  e0,
MMG5_int *  e1,
MMG5_int *  ref,
MMG5_int  idx 
)
Parameters
meshpointer toward the mesh structure.
e0pointer toward the first extremity of the edge.
e1pointer toward the second extremity of the edge.
refpointer toward the edge reference.
idxindex of the non boundary edge to get (between 1 and nb_edges)
Returns
0 if failed, 1 otherwise.

Get extremities e0, e1 and reference ref of the idx^th non boundary edge (for DG methods for example). An edge is boundary if it is located at the interface of 2 domains witch different references, if it belongs to one triangle only or if it is a singular edge (ridge or required).

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_GET_NONBDYEDGE(mesh,e0,e1,ref,idx,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(OUT):: e0,e1
INTEGER(MMG5F_INT) :: ref
INTEGER(MMG5F_INT), INTENT(IN) :: idx
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 613 of file libmmgs_tools.c.

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

◆ MMGS_Get_numberOfNonBdyEdges()

int MMGS_Get_numberOfNonBdyEdges ( MMG5_pMesh  mesh,
MMG5_int *  nb_edges 
)
Parameters
meshpointer toward the mesh structure.
nb_edgespointer toward the number of non boundary edges.
Returns
0 if failed, 1 otherwise.

Get the number of non boundary edges (for DG methods for example). An edge is boundary if it is located at the interface of 2 domains with different references, if it belongs to one triangle only or if it is a singular edge (ridge or required). Append these edges to the list of edge.

Warning
reallocate the edge array and append the internal edges. This may modify the behaviour of other functions.
Remarks
Fortran interface:

‍ SUBROUTINE MMGS_GET_NUMBEROFNONBDYEDGES(mesh,nb_edges,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(OUT):: nb_edges
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 517 of file libmmgs_tools.c.

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

◆ MMGS_parsar()

int MMGS_parsar ( int  argc,
char *  argv[],
MMG5_pMesh  mesh,
MMG5_pSol  met,
MMG5_pSol  sol 
)
Parameters
argcnumber of command line arguments.
argvcommand line arguments.
meshpointer toward the mesh structure.
metpointer toward the sol structure.
solpointer toward a level-set or displacement
Returns
1.

Store command line arguments.

Remarks
no matching fortran function.

Definition at line 126 of file libmmgs_tools.c.

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

◆ MMGS_Set_constantSize()

int MMGS_Set_constantSize ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)
Parameters
meshpointer toward the mesh structure
metpointer toward the sol structure
Returns
1 if success

Compute constant size map according to mesh->info.hsiz, mesh->info.hmin and mesh->info.hmax. Update this 3 value if not compatible.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_SET_CONSTANTSIZE(mesh,met,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 1525 of file libmmgs_tools.c.

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

◆ MMGS_setfunc()

void MMGS_setfunc ( MMG5_pMesh  mesh,
MMG5_pSol  met 
)

To associate function pointers without calling MMGS_mmgslib

Parameters
meshpointer toward the mesh structure (unused).
metpointer toward the sol structure (unused).

Set function pointers for caltet, lenedg, defsiz and gradsiz.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_SETFUNC(mesh,met)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
END SUBROUTINE

Definition at line 43 of file libmmgs_tools.c.

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

◆ MMGS_solTruncatureForOptim()

static int MMGS_solTruncatureForOptim ( MMG5_pMesh  mesh,
MMG5_pSol  met,
int  ani 
)
inlinestatic
Parameters
meshpointer toward the mesh structure.
metpointer toward the solution structure.
ani1 for aniso metric, 0 for iso one
Returns
0 if fail, 1 if succeed.

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.

Definition at line 756 of file libmmgs_tools.c.

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

◆ MMGS_stockOptions()

int MMGS_stockOptions ( MMG5_pMesh  mesh,
MMG5_Info info 
)
Parameters
meshpointer toward the mesh structure.
infopointer toward the info structure.
Returns
1.

Store the info structure in the mesh structure.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_STOCKOPTIONS(mesh,info,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,info
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 498 of file libmmgs_tools.c.

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

◆ MMGS_surfopenballRotation()

static int MMGS_surfopenballRotation ( MMG5_pMesh  mesh,
MMG5_pPoint  p0,
MMG5_int  k,
int  i,
int  ilist,
double  r[3][3],
double *  lispoi,
double  n[3] 
)
inlinestatic
Parameters
meshpointer toward the mesh structure.
p0starting point
kindex of starting element
ilocal index of p0 in k
ilistcomputed number of tria in the ball of p0
rrotation that send the normal at p0 onto the z vector
lipointrotated ball of point p0
nnormal at point p0
Returns
1 if success, 0 otherwise.

Compute the rotation matrix that sends the tangent plane at p0 onto z=0 and apply this rotation to the opened ball of p0.

Definition at line 996 of file libmmgs_tools.c.

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

◆ MMGS_unitTensor_2D()

static int MMGS_unitTensor_2D ( MMG5_pMesh  mesh,
MMG5_int  k,
int  i,
MMG5_pPoint  p1,
double *  m,
double  isqhmax 
)
inlinestatic
Parameters
meshpointer toward the mesh
kindex of starting triangle
ilocal index of point p1 in k
p1point on which we want to compute the 3D unit tensor
mpointer to store computed metric
isqhmaxsquared inverse of MMG5_HMAXCOE
Returns
1 if succeed, 0 if fail.

Compute the 2D unit tensor at point p1, the vertex number i of tetra k : edges of the ball of point are projected onto the tangent plane.

Step 1: compute ball of point

Step 2: Store or compute the suitable normal depending on the type of point we are processing.

Step 3: Rotation of the ball of p0 so lispoi will contain all the points of the ball of p0, rotated so that t_{p_0}S = [z = 0]

Step 4: computation of unit tensor

Definition at line 1094 of file libmmgs_tools.c.

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

◆ MMGS_unitTensor_3D()

static int MMGS_unitTensor_3D ( MMG5_pMesh  mesh,
MMG5_int  k,
int  i,
MMG5_pPoint  p1,
double *  m 
)
inlinestatic
Parameters
meshpointer toward the mesh
kindex of starting triangle
ilocal index of point p1 in k
p1point on which we want to compute the 3D unit tensor
mpointer to store computed metric
Returns
1 if succeed, 0 if the tensor of covariance is non invertible.

Compute the 3D unit tensor at point p1, the vertex number i of tetra k.

Remarks
If all edges starting from p1 belong to the same plane, the tensor of covariance is not invertible and m can't be computed.

Step 1: compute ball of point

Step 2: compute unit tensor

Definition at line 885 of file libmmgs_tools.c.

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

◆ MMGS_usage()

int MMGS_usage ( char *  prog)
Parameters
progpointer toward the program name.

Print help for mmgs options.

Remarks
Fortran interface:

‍ SUBROUTINE MMGS_USAGE(prog,strlen0,retval)
CHARACTER(LEN=*), INTENT(IN) :: prog
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE

Definition at line 83 of file libmmgs_tools.c.

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