Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
|
API headers and documentation for the mmgs library. More...
#include "mmg/mmgs/mmgs_export.h"
#include "mmg/common/libmmgtypes.h"
Go to the source code of this file.
Macros | |
#define | MMGS_LMAX 1024 |
Functions | |
LIBMMGS_EXPORT int | MMGS_Init_mesh (const int starter,...) |
Initialize a mesh structure and optionally the associated solution and metric structures. | |
LIBMMGS_EXPORT void | MMGS_Init_fileNames (MMG5_pMesh mesh, MMG5_pSol sol) |
Initialize file names to their default values. | |
LIBMMGS_EXPORT void | MMGS_Init_parameters (MMG5_pMesh mesh) |
Initialize the input parameters. | |
LIBMMGS_EXPORT int | MMGS_Set_inputMeshName (MMG5_pMesh mesh, const char *meshin) |
Set the name of the input mesh. | |
LIBMMGS_EXPORT int | MMGS_Set_outputMeshName (MMG5_pMesh mesh, const char *meshout) |
Set the name of the output mesh file. | |
LIBMMGS_EXPORT int | MMGS_Set_inputSolName (MMG5_pMesh mesh, MMG5_pSol sol, const char *solin) |
Set the name of the input solution file. | |
LIBMMGS_EXPORT int | MMGS_Set_outputSolName (MMG5_pMesh mesh, MMG5_pSol sol, const char *solout) |
Set the name of the output solution file. | |
LIBMMGS_EXPORT int | MMGS_Set_inputParamName (MMG5_pMesh mesh, const char *fparamin) |
Set the name of the input parameter file. | |
LIBMMGS_EXPORT int | MMGS_Set_solSize (MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol) |
Initialize an array of solution fields: set dimension, types and number of fields. | |
LIBMMGS_EXPORT int | MMGS_Set_solsAtVerticesSize (MMG5_pMesh mesh, MMG5_pSol *sol, int nsols, MMG5_int nentities, int *typSol) |
Initialize an array of solution fields defined at vertices: set dimension, types and number of fields. | |
LIBMMGS_EXPORT int | MMGS_Set_meshSize (MMG5_pMesh mesh, MMG5_int np, MMG5_int nt, MMG5_int na) |
Set the number of vertices, triangles and edges of the mesh and allocate the associated tables. | |
LIBMMGS_EXPORT int | MMGS_Set_vertex (MMG5_pMesh mesh, double c0, double c1, double c2, MMG5_int ref, MMG5_int pos) |
Set the coordinates of a single vertex. | |
LIBMMGS_EXPORT int | MMGS_Set_vertices (MMG5_pMesh mesh, double *vertices, MMG5_int *refs) |
Set the coordinates and references of all vertices in a mesh. | |
LIBMMGS_EXPORT int | MMGS_Set_triangle (MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int ref, MMG5_int pos) |
Set the coordinates and reference of a single triangle. | |
LIBMMGS_EXPORT int | MMGS_Set_triangles (MMG5_pMesh mesh, MMG5_int *tria, MMG5_int *refs) |
Set the vertices and references of all triangles in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Set_edge (MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int ref, MMG5_int pos) |
Set the vertices and reference of one edge in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Set_corner (MMG5_pMesh mesh, MMG5_int k) |
Assign the "corner" attribute to a vertex. | |
LIBMMGS_EXPORT int | MMGS_Unset_corner (MMG5_pMesh mesh, MMG5_int k) |
Remove the "corner" attribute from a vertex. | |
LIBMMGS_EXPORT int | MMGS_Set_requiredVertex (MMG5_pMesh mesh, MMG5_int k) |
Assign the "required" attribute to a vertex. | |
LIBMMGS_EXPORT int | MMGS_Unset_requiredVertex (MMG5_pMesh mesh, MMG5_int k) |
Remove the "required" attribute from a vertex. | |
LIBMMGS_EXPORT int | MMGS_Set_requiredTriangle (MMG5_pMesh mesh, MMG5_int k) |
Assign the "required" attribute to a triangle. | |
LIBMMGS_EXPORT int | MMGS_Unset_requiredTriangle (MMG5_pMesh mesh, MMG5_int k) |
Remove the "required" attribute from a vertex. | |
LIBMMGS_EXPORT int | MMGS_Set_ridge (MMG5_pMesh mesh, MMG5_int k) |
Assign the "ridge" attribute to an edge. | |
LIBMMGS_EXPORT int | MMGS_Unset_ridge (MMG5_pMesh mesh, MMG5_int k) |
Remove the "ridge" attribute from an edge. | |
LIBMMGS_EXPORT int | MMGS_Set_requiredEdge (MMG5_pMesh mesh, MMG5_int k) |
Assign the "required" attribute to an edge. | |
LIBMMGS_EXPORT int | MMGS_Unset_requiredEdge (MMG5_pMesh mesh, MMG5_int k) |
Remove the "required" attribute from an edge. | |
LIBMMGS_EXPORT int | MMGS_Set_edges (MMG5_pMesh mesh, MMG5_int *edges, MMG5_int *refs) |
Set the vertices and references of all edges in a mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_edges (MMG5_pMesh mesh, MMG5_int *edges, MMG5_int *refs, int *areRidges, int *areRequired) |
Get vertices, references and attributes of all edges in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Set_normalAtVertex (MMG5_pMesh mesh, MMG5_int k, double n0, double n1, double n2) |
Set the normal orientation at a single vertex. | |
double | MMGS_Get_triangleQuality (MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k) |
Get the quality measure of a triangle. | |
LIBMMGS_EXPORT int | MMGS_Set_scalarSol (MMG5_pSol met, double s, MMG5_int pos) |
Set a single element of a scalar solution structure. | |
LIBMMGS_EXPORT int | MMGS_Set_scalarSols (MMG5_pSol met, double *s) |
Set the values of all elements of a scalar solution structure. | |
LIBMMGS_EXPORT int | MMGS_Set_vectorSol (MMG5_pSol met, double vx, double vy, double vz, MMG5_int pos) |
Set a single element of a vector solution structure. | |
LIBMMGS_EXPORT int | MMGS_Set_vectorSols (MMG5_pSol met, double *sols) |
Set all elements of a vector solution structure. | |
LIBMMGS_EXPORT int | MMGS_Set_tensorSol (MMG5_pSol met, double m11, double m12, double m13, double m22, double m23, double m33, MMG5_int pos) |
Set a single element of a tensor solution structure. | |
LIBMMGS_EXPORT int | MMGS_Set_tensorSols (MMG5_pSol met, double *sols) |
Set all elements of a tensor solution structure. | |
LIBMMGS_EXPORT int | MMGS_Set_ithSol_inSolsAtVertices (MMG5_pSol sol, int i, double *s, MMG5_int pos) |
Set a single element of one out of multiple solution fields that are defined on vertices. | |
LIBMMGS_EXPORT int | MMGS_Set_ithSols_inSolsAtVertices (MMG5_pSol sol, int i, double *s) |
Set all elements of one out of multiple solution fields that are defined on vertices. | |
LIBMMGS_EXPORT int | MMGS_Chk_meshData (MMG5_pMesh mesh, MMG5_pSol met) |
Check if the numbers of given entities match with mesh and solution size and check mesh data. | |
LIBMMGS_EXPORT int | MMGS_Set_iparameter (MMG5_pMesh mesh, MMG5_pSol sol, int iparam, MMG5_int val) |
set an integer parameter of the remesher | |
LIBMMGS_EXPORT int | MMGS_Set_dparameter (MMG5_pMesh mesh, MMG5_pSol sol, int dparam, double val) |
set a real-valued parameter of the remesher | |
LIBMMGS_EXPORT int | MMGS_Set_localParameter (MMG5_pMesh mesh, MMG5_pSol sol, int typ, MMG5_int ref, double hmin, double hmax, double hausd) |
set a local parameter | |
LIBMMGS_EXPORT int | MMGS_Set_multiMat (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ref, int split, MMG5_int rmin, MMG5_int rplus) |
Set the reference mapping for the elements of reference ref in level-set discretization mode. | |
LIBMMGS_EXPORT int | MMGS_Set_lsBaseReference (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int br) |
Set a new level-set base reference. | |
LIBMMGS_EXPORT int | MMGS_Get_meshSize (MMG5_pMesh mesh, MMG5_int *np, MMG5_int *nt, MMG5_int *na) |
Get the number of vertices, triangles, and edges of the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_solSize (MMG5_pMesh mesh, MMG5_pSol sol, int *typEntity, MMG5_int *np, int *typSol) |
Get the number of elements, dimension, and type of a solution. | |
LIBMMGS_EXPORT int | MMGS_Get_solsAtVerticesSize (MMG5_pMesh mesh, MMG5_pSol *sol, int *nsols, MMG5_int *nentities, int *typSol) |
Get the number of elements, type, and dimensions of several solutions defined on vertices. | |
LIBMMGS_EXPORT int | MMGS_Get_vertex (MMG5_pMesh mesh, double *c0, double *c1, double *c2, MMG5_int *ref, int *isCorner, int *isRequired) |
Get the coordinates c0, c1,c2 and reference ref of the next vertex of mesh. | |
LIBMMGS_EXPORT int | MMGS_GetByIdx_vertex (MMG5_pMesh mesh, double *c0, double *c1, double *c2, MMG5_int *ref, int *isCorner, int *isRequired, MMG5_int idx) |
Get the coordinates and reference of a specific vertex in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_vertices (MMG5_pMesh mesh, double *vertices, MMG5_int *refs, int *areCorners, int *areRequired) |
Get the coordinates, references and attributes of all vertices in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_triangle (MMG5_pMesh mesh, MMG5_int *v0, MMG5_int *v1, MMG5_int *v2, MMG5_int *ref, int *isRequired) |
Get the vertices, reference, and required attribute of the next triangle in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_triangles (MMG5_pMesh mesh, MMG5_int *tria, MMG5_int *refs, int *areRequired) |
Get the vertices, references, and required attributes of all triangles in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_edge (MMG5_pMesh mesh, MMG5_int *e0, MMG5_int *e1, MMG5_int *ref, int *isRidge, int *isRequired) |
Get the vertices, reference, and attributes of the next edge in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_normalAtVertex (MMG5_pMesh mesh, MMG5_int k, double *n0, double *n1, double *n2) |
Get the normal orientation at an edge. | |
LIBMMGS_EXPORT int | MMGS_Get_scalarSol (MMG5_pSol met, double *s) |
Get the next element of a scalar solution structure. | |
LIBMMGS_EXPORT int | MMGS_Get_scalarSols (MMG5_pSol met, double *s) |
Get all elements of a scalar solution structure. | |
LIBMMGS_EXPORT int | MMGS_Get_vectorSol (MMG5_pSol met, double *vx, double *vy, double *vz) |
Get the next element of a vector solution structure. | |
LIBMMGS_EXPORT int | MMGS_Get_vectorSols (MMG5_pSol met, double *sols) |
Get all elements of a vector solution structure. | |
LIBMMGS_EXPORT int | MMGS_Get_tensorSol (MMG5_pSol met, double *m11, double *m12, double *m13, double *m22, double *m23, double *m33) |
Get the next element of a tensor solution structure. | |
LIBMMGS_EXPORT int | MMGS_Get_tensorSols (MMG5_pSol met, double *sols) |
Get all elements of a tensor solution field. | |
LIBMMGS_EXPORT int | MMGS_Get_ithSol_inSolsAtVertices (MMG5_pSol sol, int i, double *s, MMG5_int pos) |
Get one out of several solutions at a specific vertex. | |
LIBMMGS_EXPORT int | MMGS_Get_ithSols_inSolsAtVertices (MMG5_pSol sol, int i, double *s) |
Get one out of several solutions at all vertices in the mesh. | |
LIBMMGS_EXPORT int | MMGS_Get_iparameter (MMG5_pMesh mesh, MMG5_int iparam) |
Get the value of an integer parameter of the remesher. | |
LIBMMGS_EXPORT int | MMGS_loadMesh (MMG5_pMesh mesh, const char *filename) |
Load a mesh (in .mesh/.mesb format) from file. | |
LIBMMGS_EXPORT int | MMGS_loadVtpMesh (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename) |
Load a mesh and optionally a solution in VTP (VTK) format from file. | |
LIBMMGS_EXPORT int | MMGS_loadVtpMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Load a mesh and multiple solutions in VTP (VTK) format from file. | |
LIBMMGS_EXPORT int | MMGS_loadVtuMesh (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename) |
Load a mesh and possibly data in VTU (VTK) format from file. | |
LIBMMGS_EXPORT int | MMGS_loadVtuMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Load a mesh and multiple solutions in VTU (VTK) format from file. | |
LIBMMGS_EXPORT int | MMGS_loadVtkMesh (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename) |
Load a mesh and possibly data in VTK format from file. | |
LIBMMGS_EXPORT int | MMGS_loadVtkMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Load a mesh and multiple solutions in VTK format from file. | |
LIBMMGS_EXPORT int | MMGS_loadMshMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename) |
Load a mesh and possibly a solution in .msh format from file. | |
LIBMMGS_EXPORT int | MMGS_loadMshMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Load a mesh and all data from a file in MSH format. | |
LIBMMGS_EXPORT int | MMGS_loadGenericMesh (MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename) |
Load a mesh and all data from a file. The format will be guessed from the filename extension. | |
LIBMMGS_EXPORT int | MMGS_saveMesh (MMG5_pMesh mesh, const char *filename) |
Save a mesh in .mesh or .meshb format. | |
LIBMMGS_EXPORT int | MMGS_saveMshMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename) |
Write mesh and optionally one data field in MSH file format (.msh extension). | |
LIBMMGS_EXPORT int | MMGS_saveMshMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Save a mesh and multiple data fields in MSH format, ascii or binary depending on the filename extension. | |
LIBMMGS_EXPORT int | MMGS_saveVtkMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename) |
Write mesh and optionally one data field in Vtk file format (.vtk extension). | |
LIBMMGS_EXPORT int | MMGS_saveVtkMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Save a mesh and multiple data fields in VTK format. | |
LIBMMGS_EXPORT int | MMGS_saveVtuMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename) |
Write mesh and optionally one data field vtu Vtk file format (.vtu extension). | |
LIBMMGS_EXPORT int | MMGS_saveVtuMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Write a mesh and multiple data fields in vtu Vtk file format (.vtu extension). | |
LIBMMGS_EXPORT int | MMGS_saveVtpMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename) |
Save a mesh and optionally one data field in VTP format. | |
LIBMMGS_EXPORT int | MMGS_saveVtpMesh_and_allData (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Save a mesh and multiple data fields in VTP format. | |
LIBMMGS_EXPORT int | MMGS_saveGenericMesh (MMG5_pMesh mesh, MMG5_pSol sol, const char *filename) |
Save mesh data in a file whose format depends on the filename extension. | |
LIBMMGS_EXPORT int | MMGS_loadSol (MMG5_pMesh mesh, MMG5_pSol met, const char *filename) |
Load a metric field (or other solution) in medit's .sol format. | |
LIBMMGS_EXPORT int | MMGS_loadAllSols (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Load one or more solutions in a solution file in medit file format. | |
LIBMMGS_EXPORT int | MMGS_saveSol (MMG5_pMesh mesh, MMG5_pSol met, const char *filename) |
Write an isotropic or anisotropic metric in medit file format. | |
LIBMMGS_EXPORT int | MMGS_saveAllSols (MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename) |
Save one or more solutions in a solution file in medit file format. | |
LIBMMGS_EXPORT int | MMGS_Free_allSols (MMG5_pMesh mesh, MMG5_pSol *sol) |
Deallocate an array of solution fields. | |
LIBMMGS_EXPORT int | MMGS_Free_all (const int starter,...) |
Deallocations before return. | |
LIBMMGS_EXPORT int | MMGS_Free_structures (const int starter,...) |
Structure deallocations before return. | |
LIBMMGS_EXPORT int | MMGS_Free_names (const int starter,...) |
Structure deallocations before return. | |
LIBMMGS_EXPORT int | MMGS_mmgslib (MMG5_pMesh mesh, MMG5_pSol met) |
Main "program" for mesh adaptation. | |
LIBMMGS_EXPORT int | MMGS_mmgsls (MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pSol met) |
Main "program" for level-set discretization. | |
LIBMMGS_EXPORT void | MMGS_setfunc (MMG5_pMesh mesh, MMG5_pSol met) |
Set function pointers for caltet, lenedg, defsiz and gradsiz. | |
LIBMMGS_EXPORT int | MMGS_Get_numberOfNonBdyEdges (MMG5_pMesh mesh, MMG5_int *nb_edges) |
Get the number of non-boundary edges. | |
LIBMMGS_EXPORT int | MMGS_Get_nonBdyEdge (MMG5_pMesh mesh, MMG5_int *e0, MMG5_int *e1, MMG5_int *ref, MMG5_int idx) |
Get vertices and reference of a non-boundary edge. | |
LIBMMGS_EXPORT int | MMGS_Set_constantSize (MMG5_pMesh mesh, MMG5_pSol met) |
Compute a constant size map. | |
LIBMMGS_EXPORT int | MMGS_usage (char *prog) |
Print help for mmgs options. | |
LIBMMGS_EXPORT int | MMGS_parsar (int argc, char *argv[], MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol) |
Store command line arguments. | |
LIBMMGS_EXPORT int | MMGS_defaultValues (MMG5_pMesh mesh) |
Print the default parameter values. | |
LIBMMGS_EXPORT int | MMGS_stockOptions (MMG5_pMesh mesh, MMG5_Info *info) |
Store the info structure in the mesh structure. | |
LIBMMGS_EXPORT void | MMGS_destockOptions (MMG5_pMesh mesh, MMG5_Info *info) |
Recover the info structure stored in the mesh structure. | |
LIBMMGS_EXPORT int | MMGS_Get_adjaTri (MMG5_pMesh mesh, MMG5_int kel, MMG5_int listri[3]) |
Return adjacent triangles of a triangle. | |
LIBMMGS_EXPORT int | MMGS_Get_adjaVerticesFast (MMG5_pMesh mesh, MMG5_int ip, MMG5_int start, MMG5_int lispoi[MMGS_LMAX]) |
Find adjacent vertices of a triangle. | |
LIBMMGS_EXPORT int | MMGS_Compute_eigenv (double m[6], double lambda[3], double vp[3][3]) |
Compute the real eigenvalues and eigenvectors of a symmetric matrix. | |
LIBMMGS_EXPORT void | MMGS_Free_solutions (MMG5_pMesh mesh, MMG5_pSol sol) |
Free a solution. | |
LIBMMGS_EXPORT int | MMGS_Clean_isoSurf (MMG5_pMesh mesh) |
Clean data (triangles and edges) linked to isosurface. | |
LIBMMGS_EXPORT void | MMGS_Set_commonFunc (void) |
Set common function pointers between mmgs and mmg3d to the matching mmgs functions. | |
Variables | |
LIBMMGS_EXPORT int(* | MMGS_doSol )(MMG5_pMesh mesh, MMG5_pSol met) |
Compute an isotropic size map according to the mean of the length of the edges passing through a vertex. | |
API headers and documentation for the mmgs library.
These are the API functions for the mmgs library. These functions allow to load and save meshes and data defined on meshes; add, extract, or modify mesh data; and to call the library functions that perform remeshing and level-set discretization.
Meshes are here defined in terms of vertices and two-dimensional objects: triangles and quadrangles, which live in 3D space. Edges can also be represented. All of these entities can have a reference: an integer value that can serve as a group identifier. In addition mesh entities can have attributes such as "required" or "corner".
Data defined on meshes can be for example functions that are meant for level-set discretization and metric tensors that will govern edge lengths. These data can be scalar, vector, or (symmetric) tensor-valued; and there can be more than one data item associated with a mesh entity. These data are often referred to as solutions.
Two of the functions here are referred to as "programs", because they perform the tasks for which mmgs is meant: meshing and level-set discretization. The other functions merely serve to load and save data and to perform pre- and post-processing. These programs actually behave much like independent programs: they send diagnostic output to stdout and in rare cases they may call the exit() function.
Definition in file libmmgs.h.
#define MMGS_LMAX 1024 |
enum MMGS_Param |
Input parameters for mmg library.
Input parameters for mmg library. Options prefixed by MMGS_IPARAM asked for integers values ans options prefixed by MMGS_DPARAM asked for real values.
LIBMMGS_EXPORT int MMGS_Chk_meshData | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
Check if the numbers of given entities match with mesh and solution size and check mesh data.
mesh | pointer to the mesh structure. |
met | pointer to the solution structure. |
This function checks if the numbers of given entities match with the mesh and solution sizes and checks the mesh data. Use of this function is not mandatory.
SUBROUTINE MMGS_CHK_MESHDATA(mesh,met,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,met
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1260 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Clean_isoSurf | ( | MMG5_pMesh | mesh | ) |
Clean data (triangles and edges) linked to isosurface.
mesh | pointer to mesh structure |
SUBROUTINE MMGS_CLEAN_ISOSURF(mesh,retval)
MMG5_DATA_PTR_T, INTENT(IN) :: mesh
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1640 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Compute_eigenv | ( | double | m[6], |
double | lambda[3], | ||
double | vp[3][3] | ||
) |
Compute the real eigenvalues and eigenvectors of a symmetric matrix.
m | upper part of a symmetric matrix diagonalizable in |R |
lambda | array of eigenvalues |
vp | array of eigenvectors |
Compute the real eigenvalues and eigenvectors of a symmetric matrix 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
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 1607 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_defaultValues | ( | MMG5_pMesh | mesh | ) |
Print the default parameter values.
mesh | pointer to the mesh structure. |
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.
LIBMMGS_EXPORT void MMGS_destockOptions | ( | MMG5_pMesh | mesh, |
MMG5_Info * | info | ||
) |
Recover the info structure stored in the mesh structure.
mesh | pointer to the mesh structure. |
info | pointer to the info structure. |
SUBROUTINE MMGS_DESTOCKOPTIONS(mesh,info)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,info
END SUBROUTINE
Definition at line 565 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Free_all | ( | const int | starter, |
... | |||
) |
Deallocations before return.
starter | dummy argument used to initialize the variadic argument list. |
... | variadic arguments. |
For the MMGS_mmgslib function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppMet, &your_metric,MMG5_ARG_end).
For the MMGS_mmgsls function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppLs, &your_level_set,MMG5_ARG_end).
Here, your_mesh is a MMG5_pMesh, your_metric and your_level_set are MMG5_pSol.
Definition at line 1690 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Free_allSols | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol | ||
) |
Deallocate an array of solution fields.
mesh | pointer to the mesh structure. |
sol | pointer to an array of solution structures (that stores solution fields). |
SUBROUTINE MMGS_Free_allSols(mesh,sol,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1685 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Free_names | ( | const int | starter, |
... | |||
) |
Structure deallocations before return.
starter | dummy argument used to initialize the variadic argument list. |
... | variadic arguments. |
For the MMGS_mmgslib function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppMet, &your_metric,MMG5_ARG_end).
For the MMGS_mmgsls function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppLs, &your_level_set,MMG5_ARG_end).
Here, your_mesh is a MMG5_pMesh, your_metric and your_level_set are MMG5_pSol.
Definition at line 1718 of file API_functions_s.c.
LIBMMGS_EXPORT void MMGS_Free_solutions | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
Free a solution.
mesh | pointer to the mesh structure |
sol | pointer to the solution structure |
SUBROUTINE MMGS_FREE_SOLUTIONS(mesh,sol)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
END SUBROUTINE
Definition at line 1613 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Free_structures | ( | const int | starter, |
... | |||
) |
Structure deallocations before return.
starter | dummy argument used to initialize the variadic argument list. |
... | variadic arguments. |
For the MMGS_mmgslib function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppMet, &your_metric,MMG5_ARG_end).
For the MMGS_mmgsls function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppLs, &your_level_set,MMG5_ARG_end).
Here, your_mesh is a pointer to MMG5_pMesh and your_metric and your_level_set are pointers to MMG5_pSol.
Definition at line 1704 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_adjaTri | ( | MMG5_pMesh | mesh, |
MMG5_int | kel, | ||
MMG5_int | listri[3] | ||
) |
Return adjacent triangles of a triangle.
mesh | pointer to the mesh structure. |
kel | triangle index. |
listri | pointer to the array of indices of the three adjacent triangles of triangle kel (the index is 0 if there is no adjacent triangle). |
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).
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 709 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Get_adjaVerticesFast | ( | MMG5_pMesh | mesh, |
MMG5_int | ip, | ||
MMG5_int | start, | ||
MMG5_int | lispoi[MMGS_LMAX] | ||
) |
Find adjacent vertices of a triangle.
mesh | pointer to the mesh structure. |
ip | vertex index. |
start | index of a triangle holding ip. |
lispoi | pointer to an array of size MMGS_LMAX that will contain the indices of adjacent vertices to the vertex ip. |
Find the indices of the adjacent vertices of the vertex ip of the triangle start.
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 723 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Get_edge | ( | MMG5_pMesh | mesh, |
MMG5_int * | e0, | ||
MMG5_int * | e1, | ||
MMG5_int * | ref, | ||
int * | isRidge, | ||
int * | isRequired | ||
) |
Get the vertices, reference, and attributes of the next edge in the mesh.
mesh | pointer to the mesh structure. |
e0 | pointer to the index of the first vertex of the edge. |
e1 | pointer to the index of the second vertex of the edge. |
ref | pointer to the edge reference. |
isRidge | pointer to the flag saying if the edge is ridge. |
isRequired | pointer to the flag saying if the edge is required. |
This function retrieves the extremities e0, e1, reference ref, and attributes of the next edge of mesh. It is meant to be called in a loop over all edges. When it has been called as many times as there are edges in the mesh, the internal edge counter will be reset.
SUBROUTINE MMGS_GET_EDGE(mesh,e0,e1,ref,isRidge,isRequired,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(OUT):: e0,e1
INTEGER(MMG5F_INT) :: ref
INTEGER :: isRidge,isRequired
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 617 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_edges | ( | MMG5_pMesh | mesh, |
MMG5_int * | edges, | ||
MMG5_int * | refs, | ||
int * | areRidges, | ||
int * | areRequired | ||
) |
Get vertices, references and attributes of all edges in the mesh.
mesh | pointer to the mesh structure. |
edges | pointer to an array of edges. The vertices of edge i are stored in edges[(i-1)*2] and edges[(i-1)*2+1]. |
refs | edge references. refs[i-1] is the reference of edge i. |
areRidges | 1 if the edge is a ridge, 0 otherwise. |
areRequired | 1 if the edge is required, 0 otherwise. |
SUBROUTINE MMGS_GET_EDGES(mesh,edges,refs,areRidges,areRequired,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: refs(*),edges(*)
INTEGER, INTENT(OUT) :: areRequired(*),areRidges(*)
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 680 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_iparameter | ( | MMG5_pMesh | mesh, |
MMG5_int | iparam | ||
) |
Get the value of an integer parameter of the remesher.
mesh | pointer to the mesh structure. |
iparam | integer parameter to get (see MMGS_Param structure). |
This function retrieves the value of integer parameter iparam (see MMGS_Param for a list of parameters). It returns the value of the parameter, or zero if the value of iparam is not recognized.
SUBROUTINE MMGS_GET_IPARAMETER(mesh,iparam,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN) :: iparam
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1461 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_ithSol_inSolsAtVertices | ( | MMG5_pSol | sol, |
int | i, | ||
double * | s, | ||
MMG5_int | pos | ||
) |
Get one out of several solutions at a specific vertex.
sol | pointer to the array of solutions |
i | position of the solution field that we want to get. |
s | solution(s) at mesh vertex pos. The required size of this array depends on the type of solution. |
pos | index of the vertex on which we get the solution. |
This function retreives the value of field i of the solution array at vertex pos. (pos from 1 to the number of vertices included and i from 1 to the number of solutions). It works for any type of solution; the types are inferred from sol.
SUBROUTINE MMGS_GET_ITHSOL_INSOLSATVERTICES(sol,i,s,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: sol
INTEGER, INTENT(IN) :: i
INTEGER(MMG5F_INT), INTENT(IN) :: pos
REAL(KIND=8), DIMENSION(*),INTENT(OUT) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1173 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_ithSols_inSolsAtVertices | ( | MMG5_pSol | sol, |
int | i, | ||
double * | s | ||
) |
Get one out of several solutions at all vertices in the mesh.
sol | pointer to the array of solutions |
i | position of the solution field that we want to get. |
s | array of solutions at mesh vertices. The solution at vertex k is given by s[k-1] for a scalar sol, s[3*(k-1)]@3 for a vectorial solution and s[6*(k-1)]@6 for a tensor solution. |
This function retrieves the values of field i of the solution array sol (i from 1 to the number of solutions). It works for any type of solution; the type is inferred from sol.
SUBROUTINE MMGS_GET_ITHSOLS_INSOLSATVERTICES(sol,i,s,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: sol
INTEGER, INTENT(IN) :: i
REAL(KIND=8), DIMENSION(*),INTENT(OUT) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1231 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_meshSize | ( | MMG5_pMesh | mesh, |
MMG5_int * | np, | ||
MMG5_int * | nt, | ||
MMG5_int * | na | ||
) |
Get the number of vertices, triangles, and edges of the mesh.
recover data
mesh | pointer to the mesh structure. |
np | pointer to the number of vertices. |
nt | pointer to the number of triangles. |
na | pointer to the number of edges. |
SUBROUTINE MMGS_GET_MESHSIZE(mesh,np,nt,na,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT) :: np,nt,na
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 285 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_nonBdyEdge | ( | MMG5_pMesh | mesh, |
MMG5_int * | e0, | ||
MMG5_int * | e1, | ||
MMG5_int * | ref, | ||
MMG5_int | idx | ||
) |
Get vertices and reference of a non-boundary edge.
mesh | pointer to the mesh structure. |
e0 | pointer to the first extremity of the edge (a vertex number). |
e1 | pointer to the second extremity of the edge (a vertex number). |
ref | pointer to the edge reference. |
idx | index of the non boundary edge to get (between 1 and the number of edges) |
This function returns the vertices e0, e1 and reference ref of the non boundary edge idx. 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).
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 667 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Get_normalAtVertex | ( | MMG5_pMesh | mesh, |
MMG5_int | k, | ||
double * | n0, | ||
double * | n1, | ||
double * | n2 | ||
) |
Get the normal orientation at an edge.
mesh | pointer to the mesh structure. |
k | vertex number |
n0 | x componant of the normal at vertex k. |
n1 | y componant of the normal at vertex k. |
n2 | z componant of the normal at vertex k. |
This function retrieves the normal (n0,n1,n2) at vertex k.
SUBROUTINE MMGS_GET_NORMALATVERTEX(mesh,k,n0,n1,n2,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
REAL(KIND=8) :: n0,n1,n2
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 786 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_numberOfNonBdyEdges | ( | MMG5_pMesh | mesh, |
MMG5_int * | nb_edges | ||
) |
Get the number of non-boundary edges.
mesh | pointer to the mesh structure. |
the | number of edges pointer to the number of non boundary edges. |
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 edges.
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 571 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Get_scalarSol | ( | MMG5_pSol | met, |
double * | s | ||
) |
Get the next element of a scalar solution structure.
met | pointer to the sol structure. |
s | pointer to the scalar solution value. |
This function retrieves the next element s of the solution field met. It is meant to be called in a loop over all elements. When it has been called as many times as there are elements in the solution, the internal loop counter will be reset.
SUBROUTINE MMGS_GET_SCALARSOL(met,s,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), INTENT(OUT) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 858 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_scalarSols | ( | MMG5_pSol | met, |
double * | s | ||
) |
Get all elements of a scalar solution structure.
met | pointer to the solution structure. |
s | array of scalar solutions at mesh vertices. s[i-1] is the solution at vertex i. |
SUBROUTINE MMGS_GET_SCALARSOLS(met,s,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), DIMENSION(*),INTENT(OUT) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 904 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_solsAtVerticesSize | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
int * | nsols, | ||
MMG5_int * | nentities, | ||
int * | typSol | ||
) |
Get the number of elements, type, and dimensions of several solutions defined on vertices.
mesh | pointer to the mesh structure. |
sol | pointer to an array of sol structures. |
nsols | number of solutions per entity |
nentities | pointer to the number of entities. |
typSol | array of size MMG5_NSOL_MAX to store type of each solution (scalar, vectorial, ..., see MMG5_type for possible values). |
SUBROUTINE MMGS_GET_SOLSATVERTICESSIZE(mesh,sol,nsols,nentities,typSol,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER :: nsols
INTEGER(MMG5F_INT) :: nentities
INTEGER :: typSol(*)
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 256 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_solSize | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
int * | typEntity, | ||
MMG5_int * | np, | ||
int * | typSol | ||
) |
Get the number of elements, dimension, and type of a solution.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
typEntity | pointer to the type of entities to which solutions are applied. (see MMG5_entities for possible values) |
np | pointer to the number of elements in the solution. |
typSol | pointer to the type of the solution (MMG5_Scalar, MMG5_Vector, MMG5_Tensor, MMG5_Notype) |
SUBROUTINE MMGS_GET_SOLSIZE(mesh,sol,typEntity,np,typSol,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER :: typEntity,typSol
INTEGER(MMG5F_INT) :: np
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 230 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_tensorSol | ( | MMG5_pSol | met, |
double * | m11, | ||
double * | m12, | ||
double * | m13, | ||
double * | m22, | ||
double * | m23, | ||
double * | m33 | ||
) |
Get the next element of a tensor solution structure.
met | pointer to the sol structure. |
m11 | pointer to the position (1,1) in the solution tensor. |
m12 | pointer to the position (1,2) in the solution tensor. |
m13 | pointer to the position (1,3) in the solution tensor. |
m22 | pointer to the position (2,2) in the solution tensor. |
m23 | pointer to the position (2,3) in the solution tensor. |
m33 | pointer to the position (3,3) in the solution tensor. |
This function retrieves the next element \((m_{11},m_{12},m_{13},m_{22},m_{23},m_{33})\) of a tensor-valued solution field. It is meant to be called in a loop over all vertices. When it has been called as many times as there are elements in the solution, the internal loop counter will be reset.
SUBROUTINE MMGS_GET_TENSORSOL(met,m11,m12,m13,m22,m23,m33,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), INTENT(OUT) :: m11,m12,m13,m22,m23,m33
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1064 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_tensorSols | ( | MMG5_pSol | met, |
double * | sols | ||
) |
Get all elements of a tensor solution field.
met | pointer to the sol structure. |
sols | array of solution values. sols[6*(i-1)]@6 is the solution at vertex i. |
SUBROUTINE MMGS_GET_TENSORSOLS(met,sols,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), DIMENSION(*), INTENT(OUT) :: sols
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1127 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_triangle | ( | MMG5_pMesh | mesh, |
MMG5_int * | v0, | ||
MMG5_int * | v1, | ||
MMG5_int * | v2, | ||
MMG5_int * | ref, | ||
int * | isRequired | ||
) |
Get the vertices, reference, and required attribute of the next triangle in the mesh.
mesh | pointer to the mesh structure. |
v0 | pointer to the first vertex of the triangle. |
v1 | pointer to the second vertex of the triangle. |
v2 | pointer to the third vertex of the triangle. |
ref | pointer to the triangle reference. |
isRequired | pointer to the flag saying if the triangle is required. |
This function retrieves the vertices v0, v1, v2, reference ref, and required attribute isRequired of the next triangle of mesh. It is meant to be called in a loop over all triangles. When it has been called as many times as there are triangles, the internal loop counter will be reset.
SUBROUTINE MMGS_GET_TRIANGLE(mesh,v0,v1,v2,ref,isRequired,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(OUT):: v0,v1,v2,ref
INTEGER :: isRequired
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 494 of file API_functions_s.c.
double MMGS_Get_triangleQuality | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_int | k | ||
) |
Get the quality measure of a triangle.
mesh | pointer to the mesh structure. |
met | pointer to the metric structure. |
k | index of the triangle for which we want to get the quality. |
SUBROUTINE MMGS_GET_TRIANGLEQUALITY(mesh,met,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,met
INTEGER(MMG5F_INT), INTENT(IN):: k
REAL(KIND=8), INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 796 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_triangles | ( | MMG5_pMesh | mesh, |
MMG5_int * | tria, | ||
MMG5_int * | refs, | ||
int * | areRequired | ||
) |
Get the vertices, references, and required attributes of all triangles in the mesh.
mesh | pointer to the mesh structure. |
tria | pointer to an array of vertices Vertices of triangle i are stored in tria[(i-1)*3]@3. |
refs | pointer to the array of triangles references. refs[i-1] is the ref of triangle i. |
areRequired | pointer to an array of flags saying if triangles are required. areRequired[i-1]=1 if triangle i is required. |
! SUBROUTINE MMGS_GET_TRIANGLES(mesh,tria,refs,areRequired,retval)
! MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
! INTEGER(MMG5F_INT), DIMENSION(*),INTENT(OUT) :: tria
! INTEGER(MMG5F_INT), DIMENSION(*) :: refs
! INTEGER, DIMENSION(*) :: areRequired
! INTEGER, INTENT(OUT) :: retval
! END SUBROUTINE
Definition at line 559 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_vectorSol | ( | MMG5_pSol | met, |
double * | vx, | ||
double * | vy, | ||
double * | vz | ||
) |
Get the next element of a vector solution structure.
met | pointer to the sol structure. |
vx | x value of the vectorial solution. |
vy | y value of the vectorial solution. |
vz | z value of the vectorial solution. |
This function retrieves the next vector-valued element \((v_x,v_y,vz)\) of a solution field. It is meant to be called in a loop over all elements. When it has been called as many times as there are elements in the solution, the internal loop counter will be reset.
SUBROUTINE MMGS_GET_VECTORSOL(met,vx,vy,vz,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), INTENT(OUT) :: vx,vy,vz
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 952 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_vectorSols | ( | MMG5_pSol | met, |
double * | sols | ||
) |
Get all elements of a vector solution structure.
met | pointer to the sol structure. |
sols | array of solutions at mesh vertices. sols[3*(i-1)]@3 is the solution at vertex i. |
This function retrieves all elements of a vector-valued solution field.
SUBROUTINE MMGS_GET_VECTORSOLS(met,sols,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), DIMENSION(*),INTENT(OUT) :: sols
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1007 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_vertex | ( | MMG5_pMesh | mesh, |
double * | c0, | ||
double * | c1, | ||
double * | c2, | ||
MMG5_int * | ref, | ||
int * | isCorner, | ||
int * | isRequired | ||
) |
Get the coordinates c0, c1,c2 and reference ref of the next vertex of mesh.
mesh | pointer to the mesh structure. |
c0 | pointer to the coordinate of the vertex along the first dimension. |
c1 | pointer to the coordinate of the vertex along the second dimension. |
c2 | pointer to the coordinate of the vertex along the third dimension. |
ref | pointer to the vertex reference. |
isCorner | pointer to the flag saying if the vertex is corner. |
isRequired | pointer to the flag saying if the vertex is required. |
This function retrieves the coordinates c0, c1,c2, reference ref, and attributes of the next vertex of a mesh. It is meant to be used in a loop over all vertices. When this function has been called as many times as there are vertices, the internal loop counter will be reset. To obtain data for a specific vertex, the MMGS_GetByIdx_vertex function can be used instead.
SUBROUTINE MMGS_GET_VERTEX(mesh,c0,c1,c2,ref,isCorner,isRequired,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
REAL(KIND=8), INTENT(OUT) :: c0,c1,c2
INTEGER(MMG5F_INT) :: ref
INTEGER :: isCorner,isRequired
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 360 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Get_vertices | ( | MMG5_pMesh | mesh, |
double * | vertices, | ||
MMG5_int * | refs, | ||
int * | areCorners, | ||
int * | areRequired | ||
) |
Get the coordinates, references and attributes of all vertices in the mesh.
mesh | pointer to the mesh structure. |
vertices | pointer to the array of coordinates. The coordinates of vertex i are stored in vertices[(i-1)*3]@3. |
refs | pointer to the array of vertex references. The ref of vertex i is stored in refs[i-1]. |
areCorners | pointer to the array of flags saying if vertices are corners. areCorners[i-1]=1 if vertex i is corner. |
areRequired | pointer to the table of flags saying if vertices are required. areRequired[i-1]=1 if vertex i is required. |
! SUBROUTINE MMGS_GET_VERTICES(mesh,vertices,refs,areCorners,&
! areRequired,retval)
! MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
! REAL(KIND=8),DIMENSION(*), INTENT(OUT) :: vertices
! INTEGER(MMG5F_INT), DIMENSION(*) :: refs
! INTEGER, DIMENSION(*) :: areCorners,areRequired
! INTEGER, INTENT(OUT) :: retval
! END SUBROUTINE
Definition at line 420 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_GetByIdx_vertex | ( | MMG5_pMesh | mesh, |
double * | c0, | ||
double * | c1, | ||
double * | c2, | ||
MMG5_int * | ref, | ||
int * | isCorner, | ||
int * | isRequired, | ||
MMG5_int | idx | ||
) |
Get the coordinates and reference of a specific vertex in the mesh.
mesh | pointer to the mesh structure. |
c0 | pointer to the coordinate of the vertex along the first dimension. |
c1 | pointer to the coordinate of the vertex along the second dimension. |
c2 | pointer to the coordinate of the vertex along the third dimension. |
ref | pointer to the vertex reference. |
isCorner | pointer to the flag saying if the vertex is corner. |
isRequired | pointer to the flag saying if the vertex is required. |
idx | index of vertex to get. |
This function retrieves the coordinates c0, c1, c2 and reference ref of vertex idx of mesh, as well as its "corner" and "required" attributes.
SUBROUTINE MMGS_GETBYIDX_VERTEX(mesh,c0,c1,c2,ref,isCorner,isRequired,idx,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
REAL(KIND=8), INTENT(OUT) :: c0,c1,c2
INTEGER :: isCorner,isRequired
INTEGER(MMG5F_INT) :: ref,idx
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 387 of file API_functions_s.c.
LIBMMGS_EXPORT void MMGS_Init_fileNames | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol | ||
) |
Initialize file names to their default values.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
SUBROUTINE MMGS_INIT_FILENAMES(mesh,sol)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
END SUBROUTINE
Definition at line 56 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Init_mesh | ( | const int | starter, |
... | |||
) |
Initialize a mesh structure and optionally the associated solution and metric structures.
starter | dummy argument used to initialize the variadic argument list |
... | variadic arguments. |
For the MMGS_mmgslib function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppMet, &your_metric,MMG5_ARG_end).
For the MMGS_mmgsls function, you need to call the MMGS_Init_mesh function with the following arguments : MMGS_Init_mesh(MMG5_ARG_start,MMG5_ARG_ppMesh, &your_mesh, MMG5_ARG_ppLs, &your_level_set,MMG5_ARG_end).
Here, your_mesh is a MMG5_pMesh, your_metric and your_level_set are MMG5_pSol.
Definition at line 43 of file API_functions_s.c.
LIBMMGS_EXPORT void MMGS_Init_parameters | ( | MMG5_pMesh | mesh | ) |
Initialize the input parameters.
mesh | pointer to the mesh structure. |
Initialization of the input parameters (stored in the Info structure).
SUBROUTINE MMGS_INIT_PARAMETERS(mesh)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
END SUBROUTINE
Definition at line 84 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_loadAllSols | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Load one or more solutions in a solution file in medit file format.
mesh | pointer to the mesh structure. |
sol | pointer to the solutions array |
filename | name of the file to load. |
SUBROUTINE MMGS_LOADALLSOLS(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Read the file header
Sol tab allocation
Definition at line 1387 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_loadGenericMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Load a mesh and all data from a file. The format will be guessed from the filename extension.
mesh | pointer to the mesh structure. |
met | pointer to the metric structure or the NULL pointer. |
sol | pointer to the level-set structure or the NULL pointer. |
filename | name of the file to load. |
SUBROUTINE MMGS_LOADGENERICMESH(mesh,met,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 764 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_loadMesh | ( | MMG5_pMesh | mesh, |
const char * | filename | ||
) |
Load a mesh (in .mesh/.mesb format) from file.
mesh | pointer to the mesh structure. |
filename | name of the file to load. |
This function reads .mesh (ASCII) and .meshb (binary) files. If the name contains ".mesh" the file will be read as an ASCII file and if the name contains .meshb it be read as a binary file. If the file contains neither of these strings the function will first try to open "[filename].meshb" and if this fails it will try "[filename].mesh".
SUBROUTINE MMGS_LOADMESH(mesh,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 41 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_loadMshMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Load a mesh and possibly a solution in .msh format from file.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to load. |
Read a mesh and optionally one data field in MSH file format (.msh extension). We read only low-order vertices, edges, triangles, quadrangles, tetrahedra and prisms.
SUBROUTINE MMGS_LOADMSHMESH(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 645 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_loadMshMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Load a mesh and all data from a file in MSH format.
mesh | pointer to the mesh structure. |
sol | pointer to a list of solution structures. |
filename | name of the file to load. |
Read a mesh and multiple data in MSH file format (.msh extension). We read only low-order vertices, edges, triangles, quadrangles, tetrahedra and prisms.
SUBROUTINE MMGS_LOADMSHMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 706 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_loadSol | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
const char * | filename | ||
) |
Load a metric field (or other solution) in medit's .sol format.
mesh | pointer to the mesh structure. |
met | pointer to the sol structure. |
filename | name of the file to load. |
Load metric field. The solution file (in medit file format) must contain only 1 solution: the metric.
SUBROUTINE MMGS_LOADSOL(mesh,met,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Read the file header
Definition at line 1312 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_loadVtkMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Load a mesh and possibly data in VTK format from file.
mesh | pointer to the mesh structure. |
met | pointer to the metric structure or the NULL pointer. |
sol | pointer to the level-set structure or the NULL pointer. |
filename | name of the file to load. |
Read a mesh and optionally one data field in VTK vtk file format (.vtk extension). We read only low-order vertices, edges, triangles and quadrangles.
SUBROUTINE MMGS_LOADVTKMESH(mesh,met,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 155 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_loadVtkMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Load a mesh and multiple solutions in VTK format from file.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to load. |
Read a mesh and multiple data field in VTK vtk file format (.vtk extension). We read only low-order vertices, edges, triangles and quadrangles.
SUBROUTINE MMGS_LOADVTKMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 199 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_loadVtpMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Load a mesh and optionally a solution in VTP (VTK) format from file.
mesh | pointer to the mesh structure. |
met | pointer to the metric structure or the NULL pointer. |
sol | pointer to the level-set structure or the NULL pointer. |
filename | name of the file to load. |
This function reads a mesh and optionally one data field in VTK vtp file format (.vtp extension). We read only low-order vertices, edges, triangles and quadrangles.
SUBROUTINE MMGS_LOADVTPMESH(mesh,met,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 75 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_loadVtpMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Load a mesh and multiple solutions in VTP (VTK) format from file.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to load. |
Read a mesh and multiple data fields in VTK vtp file format (.vtp extension). We read only low-order vertices, edges, triangles and quadrangles.
SUBROUTINE MMGS_LOADVTPMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 119 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_loadVtuMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Load a mesh and possibly data in VTU (VTK) format from file.
mesh | pointer to the mesh structure. |
met | pointer to the metric structure or the NULL pointer. |
sol | pointer to the level-set structure or the NULL pointer. |
filename | name of the file to load. |
Read a mesh and optionally one data field in VTK vtu file format (.vtu extension). We read only low-order vertices, edges, triangles and quadrangles.
SUBROUTINE MMGS_LOADVTUMESH(mesh,met,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 235 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_loadVtuMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Load a mesh and multiple solutions in VTU (VTK) format from file.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to load. |
Read a mesh and multiple data field in VTK vtu file format (.vtu extension). We read only low-order vertices, edges, triangles and quadrangles.
SUBROUTINE MMGS_LOADVTUMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 279 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_mmgslib | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
Main "program" for mesh adaptation.
mesh | pointer to the mesh structure. |
met | pointer to the sol (metric) structure. |
Main program for the library.
SUBROUTINE MMGS_MMGSLIB(mesh,met,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
In debug mode, check that all structures are allocated
Free topologic tables (adja, xpoint, xtetra) resulting from a previous run
Definition at line 545 of file libmmgs.c.
LIBMMGS_EXPORT int MMGS_mmgsls | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
MMG5_pSol | met | ||
) |
Main "program" for level-set discretization.
mesh | pointer to the mesh structure. |
sol | pointer to the sol (level-set) structure. |
met | pointer to the sol (metric) structure (optionnal). |
Main program for level set discretization library. If a metric met is provided, use it to adapt the mesh.
SUBROUTINE MMGS_MMGSLS(mesh,sol,met,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
MMG5_DATA_PTR_T :: met
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
In debug mode, check that all structures are allocated
Free topologic tables (adja, xpoint, xtetra) resulting from a previous run
Definition at line 298 of file libmmgs.c.
LIBMMGS_EXPORT int MMGS_parsar | ( | int | argc, |
char * | argv[], | ||
MMG5_pMesh | mesh, | ||
MMG5_pSol | met, | ||
MMG5_pSol | sol | ||
) |
Store command line arguments.
argc | number of command line arguments. |
argv | command line arguments. |
mesh | pointer to the mesh structure. |
met | pointer to the sol structure. |
sol | pointer to a level-set or displacement |
Definition at line 126 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_saveAllSols | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Save one or more solutions in a solution file in medit file format.
mesh | pointer to the mesh structure. |
sol | pointer to the solutions array |
filename | name of the solution file. |
SUBROUTINE MMGS_SAVEALLSOLS(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1520 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_saveGenericMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Save mesh data in a file whose format depends on the filename extension.
mesh | pointer to the mesh structure. |
filename | name of the file to write. |
SUBROUTINE MMGS_SAVEGENERICMESH(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1609 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_saveMesh | ( | MMG5_pMesh | mesh, |
const char * | filename | ||
) |
Save a mesh in .mesh or .meshb format.
mesh | pointer to the mesh structure. |
filename | name of the file to load. |
This function saves a mesh in .mesh or .meshb format (depending on the filename extension).
SUBROUTINE MMGS_SAVEMESH(mesh,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 842 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_saveMshMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Write mesh and optionally one data field in MSH file format (.msh extension).
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to load. |
The file is saved in ASCII format for .msh extension, an in binary format for a .mshb extension.
SUBROUTINE MMGS_SAVEMSHMESH(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1302 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_saveMshMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Save a mesh and multiple data fields in MSH format, ascii or binary depending on the filename extension.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
This function saves a mesh and multiple data fields (that are considered as solutions and not metrics, thus, we do nothing over the ridge vertices) in MSH file format (.msh extension). The file is saved in ASCII format for .msh extension and in binary format for a .mshb extension.
SUBROUTINE MMGS_SAVEMSHMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1307 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_saveSol | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met, | ||
const char * | filename | ||
) |
Write an isotropic or anisotropic metric in medit file format.
mesh | pointer to the mesh structure. |
met | pointer to the sol structure. |
filename | name of the file to write. |
SUBROUTINE MMGS_SAVESOL(mesh,met,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1483 of file inout_s.c.
LIBMMGS_EXPORT int MMGS_saveVtkMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Write mesh and optionally one data field in Vtk file format (.vtk extension).
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
SUBROUTINE MMGS_SAVEVTKMESH(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 350 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_saveVtkMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Save a mesh and multiple data fields in VTK format.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
This function writes a mesh and a list of data fields in Vtk file format (.vtk extension).
SUBROUTINE MMGS_SAVEVTKMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 365 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_saveVtpMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Save a mesh and optionally one data field in VTP format.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
This function writes a mesh and optionally one data in polydata Vtk file format (.vtp extension).
SUBROUTINE MMGS_SAVEVTPMESH(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 385 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_saveVtpMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Save a mesh and multiple data fields in VTP format.
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
This function writes a mesh and multiple data fields in polydata Vtk file format (.vtp extension).
SUBROUTINE MMGS_SAVEVTPMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 400 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_saveVtuMesh | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | filename | ||
) |
Write mesh and optionally one data field vtu Vtk file format (.vtu extension).
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
SUBROUTINE MMGS_SAVEVTUMESH(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 315 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT int MMGS_saveVtuMesh_and_allData | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
const char * | filename | ||
) |
Write a mesh and multiple data fields in vtu Vtk file format (.vtu extension).
mesh | pointer to the mesh structure. |
sol | pointer to the solution structure. |
filename | name of the file to write. |
SUBROUTINE MMGS_SAVEVTUMESH_AND_ALLDATA(mesh,sol,filename,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: filename
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 330 of file inoutcpp_s.cpp.
LIBMMGS_EXPORT void MMGS_Set_commonFunc | ( | void | ) |
LIBMMGS_EXPORT int MMGS_Set_constantSize | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
Compute a constant size map.
mesh | pointer to the mesh structure |
met | pointer to the sol structure |
This function computes a constant size map according to mesh->info.hsiz, mesh->info.hmin and mesh->info.hmax. It updates these 3 values if they are not compatible.
SUBROUTINE MMGS_SET_CONSTANTSIZE(mesh,met,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1579 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Set_corner | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Assign the "corner" attribute to a vertex.
mesh | pointer to the mesh structure. |
k | vertex number. |
This function sets the corner attribute at vertex pos (pos from 1 to the number of vertices included).
SUBROUTINE MMGS_SET_CORNER(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 710 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_dparameter | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
int | dparam, | ||
double | val | ||
) |
set a real-valued parameter of the remesher
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure (unused). |
dparam | double parameter to set (see MMGS_Param for a list of parameters that can be set). |
val | value of the parameter. |
This function sets the double parameter dparam to value val.
SUBROUTINE MMGS_SET_DPARAMETER(mesh,sol,dparam,val,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
MMG5_DATA_PTR_T :: sol
INTEGER, INTENT(IN) :: dparam
REAL(KIND=8), INTENT(IN) :: val
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1511 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_edge | ( | MMG5_pMesh | mesh, |
MMG5_int | v0, | ||
MMG5_int | v1, | ||
MMG5_int | ref, | ||
MMG5_int | pos | ||
) |
Set the vertices and reference of one edge in the mesh.
mesh | pointer to the mesh structure. |
v0 | first extremity of the edge. |
v1 | second extremity of the edge. |
ref | edge reference. |
pos | edge position in the mesh. |
Assigns vertices v0, v1 and reference ref to the edge at position pos in the mesh structure (pos from 1 to the number of edges included).
SUBROUTINE MMGS_SET_EDGE(mesh,v0,v1,ref,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: v0,v1,ref,pos
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 585 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_edges | ( | MMG5_pMesh | mesh, |
MMG5_int * | edges, | ||
MMG5_int * | refs | ||
) |
Set the vertices and references of all edges in a mesh.
mesh | pointer to the mesh structure. |
edges | pointer to an array of edges. The vertices of edge i should be given in edges[(i-1)*2] and edges[(i-1)*2+1]. |
refs | pointer to an array of references. refs[i-1] is the reference of edge i. |
SUBROUTINE MMGS_SET_EDGES(mesh,edges,refs,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: edges(*),refs(*)
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 663 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_inputMeshName | ( | MMG5_pMesh | mesh, |
const char * | meshin | ||
) |
Set the name of the input mesh.
mesh | pointer to the mesh structure. |
meshin | input mesh name. |
SUBROUTINE MMGS_SET_INPUTMESHNAME(mesh,meshin,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
CHARACTER(LEN=*), INTENT(IN) :: meshin
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 63 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_inputParamName | ( | MMG5_pMesh | mesh, |
const char * | fparamin | ||
) |
Set the name of the input parameter file.
mesh | pointer to the mesh structure. |
fparamin | name of the input parameter file. |
SUBROUTINE MMGS_SET_INPUTPARAMNAME(mesh,fparamin,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
CHARACTER(LEN=*), INTENT(IN) :: fparamin
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 72 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_inputSolName | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | solin | ||
) |
Set the name of the input solution file.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
solin | name of the input solution file. |
SUBROUTINE MMGS_SET_INPUTSOLNAME(mesh,sol,solin,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: solin
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 68 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_iparameter | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
int | iparam, | ||
MMG5_int | val | ||
) |
set an integer parameter of the remesher
functions to set parameters
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure (unused). |
iparam | integer parameter to set (see MMGS_Param for a list of parameters that can be set). |
val | value for the parameter. |
This function sets integer parameter iparam to value val.
SUBROUTINE MMGS_SET_IPARAMETER(mesh,sol,iparam,val,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
MMG5_DATA_PTR_T :: sol
INTEGER, INTENT(IN) :: iparam
INTEGER(MMG5F_INT), INTENT(IN) :: val
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1305 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_ithSol_inSolsAtVertices | ( | MMG5_pSol | sol, |
int | i, | ||
double * | s, | ||
MMG5_int | pos | ||
) |
Set a single element of one out of multiple solution fields that are defined on vertices.
sol | pointer to the array of solutions |
i | position of the solution field that we want to set. |
s | solution(s) at mesh vertex pos. |
pos | index of the vertex on which we set the solution. |
Set values of the solution at field i of the solution array and at position \pos (pos from 1 to the number of vertices included and i from 1 to the number of solutions). The type of solution is determined from sol.
SUBROUTINE MMGS_SET_ITHSOL_INSOLSATVERTICES(sol,i,s,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: sol
INTEGER, INTENT(IN) :: i
INTEGER(MMG5F_INT), INTENT(IN) :: pos
REAL(KIND=8), DIMENSION(*),INTENT(OUT) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1146 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_ithSols_inSolsAtVertices | ( | MMG5_pSol | sol, |
int | i, | ||
double * | s | ||
) |
Set all elements of one out of multiple solution fields that are defined on vertices.
sol | pointer to the array of solutions |
i | position of the solution field that we want to set. |
s | array of solutions at mesh vertices. The solution at vertex k is given by s[k-1] for a scalar sol, s[3*(k-1)]@3 for a vectorial solution and s[6*(k-1)]@6 for a tensor solution. |
Set values of the solution at field i of the solution array (i from 1 to the number of solutions). The type of solution is determined from sol.
SUBROUTINE MMGS_SET_ITHSOLS_INSOLSATVERTICES(sol,i,s,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: sol
INTEGER, INTENT(IN) :: i
REAL(KIND=8), DIMENSION(*),INTENT(OUT) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1203 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_localParameter | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
int | typ, | ||
MMG5_int | ref, | ||
double | hmin, | ||
double | hmax, | ||
double | hausd | ||
) |
set a local parameter
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
typ | type of entity (triangle, edge,...). |
ref | reference of the entity. |
hmin | minimal edge length. |
hmax | maximal edge length. |
Hausdorff | distance. |
Set local parameters: set the hausdorff distance at hausd, the minmal edge length at hmin and the maximal edge length at hmax for all elements of type typ and reference ref.
SUBROUTINE MMGS_SET_LOCALPARAMETER(mesh,sol,typ,ref,hmin,hmax,hausd,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER, INTENT(IN) :: typ
INTEGER(MMG5F_INT), INTENT(IN) :: ref
REAL(KIND=8), INTENT(IN) :: hmin,hmax,hausd
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1591 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_lsBaseReference | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
MMG5_int | br | ||
) |
Set a new level-set base reference.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
br | new level-set base reference. |
Set a new level-set base reference of ref br in level-set discretization mode. Base references are boundary conditions to which an implicit domain can be attached. All implicit volumes that are not attached to listed base references are deleted as spurious volumes by the rmc option.
SUBROUTINE MMGS_SET_LSBASEREFERENCE(mesh,sol,br,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER(MMG5F_INT), INTENT(IN):: br
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1681 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_meshSize | ( | MMG5_pMesh | mesh, |
MMG5_int | np, | ||
MMG5_int | nt, | ||
MMG5_int | na | ||
) |
Set the number of vertices, triangles and edges of the mesh and allocate the associated tables.
mesh | pointer to the mesh structure. |
np | number of vertices. |
nt | number of triangles. |
na | number of edges. |
If called again, this function resets the whole mesh to reallocate it at the new size
SUBROUTINE MMGS_SET_MESHSIZE(mesh,np,nt,na,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT) :: np,nt,na
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 185 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_multiMat | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
MMG5_int | ref, | ||
int | split, | ||
MMG5_int | rmin, | ||
MMG5_int | rplus | ||
) |
Set the reference mapping for the elements of reference ref in level-set discretization mode.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
ref | input triangle reference. |
split | MMG5_MMAT_NoSplit if the entity must not be split, MMG5_MMAT_Split otherwise |
rmin | reference for the negative side after LS discretization |
rplus | reference for the positive side after LS discretization |
With this function you can determine which references will be given to the triangles on both sides of the level set, after discretization. Negative and positive here refer to areas where the function is smaller or larger, respectively, than the isovalue of the level set.
SUBROUTINE MMGS_SET_MULTIMAT(mesh,sol,ref,split,rmin,rplus,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER(MMG5F_INT), INTENT(IN):: ref,rmin,rplus
INTEGER, INTENT(IN) :: split
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1676 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_normalAtVertex | ( | MMG5_pMesh | mesh, |
MMG5_int | k, | ||
double | n0, | ||
double | n1, | ||
double | n2 | ||
) |
Set the normal orientation at a single vertex.
mesh | pointer to the mesh structure. |
k | vertex index |
n0 | x componant of the normal at vertex k. |
n1 | y componant of the normal at vertex k. |
n2 | z componant of the normal at vertex k. |
Set normal (n0,n1,n2) at vertex k.
SUBROUTINE MMGS_SET_NORMALATVERTEX(mesh,k,n0,n1,n2,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
REAL(KIND=8), INTENT(IN) :: n0,n1,n2
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 774 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_outputMeshName | ( | MMG5_pMesh | mesh, |
const char * | meshout | ||
) |
Set the name of the output mesh file.
mesh | pointer to the mesh structure. |
meshout | name of the output mesh file. |
SUBROUTINE MMGS_SET_OUTPUTMESHNAME(mesh,meshout,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh
CHARACTER(LEN=*), INTENT(IN) :: meshout
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 76 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_outputSolName | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
const char * | solout | ||
) |
Set the name of the output solution file.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
solout | name of the output solution file. |
SUBROUTINE MMGS_SET_OUTPUTSOLNAME(mesh,sol,solout,strlen0,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,sol
CHARACTER(LEN=*), INTENT(IN) :: solout
INTEGER, INTENT(IN) :: strlen0
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 81 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_requiredEdge | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Assign the "required" attribute to an edge.
mesh | pointer to the mesh structure. |
k | edge index. |
This function makes edge k a required edge. Required edges will not be modified by the remesher.
SUBROUTINE MMGS_SET_REQUIREDEDGE(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 762 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_requiredTriangle | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Assign the "required" attribute to a triangle.
mesh | pointer to the mesh structure. |
k | triangle index. |
This function sets the required attribute at triangle k.
SUBROUTINE MMGS_SET_REQUIREDTRIANGLE(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 735 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_requiredVertex | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Assign the "required" attribute to a vertex.
mesh | pointer to the mesh structure. |
k | vertex number. |
This function sets the required attribute at vertex k. Vertices with this attribute will not be modified by the remesher.
SUBROUTINE MMGS_SET_REQUIREDVERTEX(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 722 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_ridge | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Assign the "ridge" attribute to an edge.
mesh | pointer to the mesh structure. |
k | edge index. |
This function gives the ridge attribute to edge k. This influences how this edge is treated by the remesher.
SUBROUTINE MMGS_SET_RIDGE(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 751 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_scalarSol | ( | MMG5_pSol | met, |
double | s, | ||
MMG5_int | pos | ||
) |
Set a single element of a scalar solution structure.
met | pointer to the sol structure. |
s | solution scalar value. |
pos | position of the solution in the mesh. |
This function sets the scalar value s at position pos in the solution structure (pos from 1 to the number of vertices included).
SUBROUTINE MMGS_SET_SCALARSOL(met,s,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), INTENT(IN) :: s
INTEGER(MMG5F_INT), INTENT(IN):: pos
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 823 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_scalarSols | ( | MMG5_pSol | met, |
double * | s | ||
) |
Set the values of all elements of a scalar solution structure.
met | pointer to the sol structure. |
s | array of scalar solutions values. s[i-1] is the solution at vertex i. |
SUBROUTINE MMGS_SET_SCALARSOLS(met,s,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8),DIMENSION(*), INTENT(IN) :: s
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 887 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_solsAtVerticesSize | ( | MMG5_pMesh | mesh, |
MMG5_pSol * | sol, | ||
int | nsols, | ||
MMG5_int | nentities, | ||
int * | typSol | ||
) |
Initialize an array of solution fields defined at vertices: set dimension, types and number of fields.
mesh | pointer to the mesh structure. |
sol | pointer to an allocatable sol structure. |
nsols | number of solutions per entity |
nentities | number of entities |
typSol | Array of size nsol listing the type of the solutions (scalar, vectorial, ..., see MMG5_type for possible values) |
To use to initialize an array of solution fields (not used by Mmg itself).
SUBROUTINE MMGS_SET_SOLSATVERTICESSIZE(mesh,sol,nsols,nentities,typSol,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER, INTENT(IN) :: nsols
INTEGER(MMG5F_INT), INTENT(IN):: nentities
INTEGER, INTENT(IN) :: typSol(*)
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Sol tab allocation
Definition at line 140 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_solSize | ( | MMG5_pMesh | mesh, |
MMG5_pSol | sol, | ||
int | typEntity, | ||
MMG5_int | np, | ||
int | typSol | ||
) |
Initialize an array of solution fields: set dimension, types and number of fields.
mesh | pointer to the mesh structure. |
sol | pointer to the sol structure. |
typEntity | types of solution entities (vertices, triangles, ... see MMG5_entities for possible values). |
np | number of solutions. |
typSol | type of solution (scalar, vectorial, ..., see MMG5_type for possible values) |
To use to initialize an array of solution fields (not used by Mmg itself).
SUBROUTINE MMGS_SET_SOLSIZE(mesh,sol,typEntity,np,typSol,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh,sol
INTEGER, INTENT(IN) :: typEntity,typSol
INTEGER(MMG5F_INT), INTENT(IN):: np
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 93 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_tensorSol | ( | MMG5_pSol | met, |
double | m11, | ||
double | m12, | ||
double | m13, | ||
double | m22, | ||
double | m23, | ||
double | m33, | ||
MMG5_int | pos | ||
) |
Set a single element of a tensor solution structure.
met | pointer to the sol structure. |
m11 | value of the tensorial solution at position (1,1) in the tensor. |
m12 | value of the tensorial solution at position (1,2) in the tensor. |
m13 | value of the tensorial solution at position (1,3) in the tensor. |
m22 | value of the tensorial solution at position (2,2) in the tensor. |
m23 | value of the tensorial solution at position (2,3) in the tensor. |
m33 | value of the tensorial solution at position (3,3) in the tensor. |
pos | position of the solution in the mesh (begin to 1). |
This function sets a tensor value at position pos in solution structure (pos from 1 to the number of vertices included).
SUBROUTINE MMGS_SET_TENSORSOL(met,m11,m12,m13,m22,m23,m33,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), INTENT(IN) :: m11,m12,m13,m22,m23,m33
INTEGER(MMG5F_INT), INTENT(IN):: pos
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1022 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_tensorSols | ( | MMG5_pSol | met, |
double * | sols | ||
) |
Set all elements of a tensor solution structure.
met | pointer to the sol structure. |
sols | array of tensorial solutions. sols[6*(i-1)]@6 is the solution at vertex i |
SUBROUTINE MMGS_SET_TENSORSOLS(met,sols,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8),DIMENSION(*), INTENT(IN) :: sols
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 1101 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_triangle | ( | MMG5_pMesh | mesh, |
MMG5_int | v0, | ||
MMG5_int | v1, | ||
MMG5_int | v2, | ||
MMG5_int | ref, | ||
MMG5_int | pos | ||
) |
Set the coordinates and reference of a single triangle.
mesh | pointer to the mesh structure. |
v0 | first vertex of triangle. |
v1 | second vertex of triangle. |
v2 | third vertex of triangle. |
ref | triangle reference. |
pos | triangle position in the mesh. |
This function sets a triangle with vertices v0, v1, v2 and reference ref at position pos in the mesh structure (pos from 1 to the number of triangles included).
SUBROUTINE MMGS_SET_TRIANGLE(mesh,v0,v1,v2,ref,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: v0,v1,v2,ref,pos
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 456 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_triangles | ( | MMG5_pMesh | mesh, |
MMG5_int * | tria, | ||
MMG5_int * | refs | ||
) |
Set the vertices and references of all triangles in the mesh.
mesh | pointer to the mesh structure. |
tria | pointer to an array of vertex numbers Vertices of the \(i^{th}\) triangles are stored in tria[(i-1)*3]@3. |
refs | pointer to an array of triangle references. refs[i-1] is the reference of the \(i^{th}\) triangle. |
! SUBROUTINE MMGS_SET_TRIANGLES(mesh,tria,refs,retval)
! MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
! INTEGER(MMG5F_INT),DIMENSION(*), INTENT(IN) :: tria,refs
! INTEGER, INTENT(OUT) :: retval
! END SUBROUTINE
Definition at line 537 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_vectorSol | ( | MMG5_pSol | met, |
double | vx, | ||
double | vy, | ||
double | vz, | ||
MMG5_int | pos | ||
) |
Set a single element of a vector solution structure.
met | pointer to the sol structure. |
vx | x value of the vectorial solution. |
vy | y value of the vectorial solution. |
vz | z value of the vectorial solution. |
pos | position of the solution in the mesh (begin to 1). |
This function sets the vectorial value \((v_x,v_y,v_z)\) at position pos in the solution structure (pos from 1 to the number of vertices included).
SUBROUTINE MMGS_SET_VECTORSOL(met,vx,vy,vz,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8), INTENT(IN) :: vx,vy,vz
INTEGER(MMG5F_INT), INTENT(IN):: pos
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 914 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_vectorSols | ( | MMG5_pSol | met, |
double * | sols | ||
) |
Set all elements of a vector solution structure.
met | pointer to the sol structure. |
sols | array of vectorial solutions sols[3*(i-1)]@3 is the solution at vertex i |
This function sets a vector-valued solution at each element of solution structure.
SUBROUTINE MMGS_SET_VECTORSOLS(met,sols,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: met
REAL(KIND=8),DIMENSION(*), INTENT(IN) :: sols
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 984 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_vertex | ( | MMG5_pMesh | mesh, |
double | c0, | ||
double | c1, | ||
double | c2, | ||
MMG5_int | ref, | ||
MMG5_int | pos | ||
) |
Set the coordinates of a single vertex.
mesh | pointer to the mesh structure. |
c0 | coordinate of the vertex along the first dimension. |
c1 | coordinate of the vertex along the second dimension. |
c2 | coordinate of the vertex along the third dimension. |
ref | vertex reference. |
pos | position of the vertex in the mesh. |
Set vertex coordinates c0, c1,c2 and reference ref at position pos in the mesh structure (pos from 1 to the number of vertices included).
SUBROUTINE MMGS_SET_VERTEX(mesh,c0,c1,c2,ref,pos,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
REAL(KIND=8), INTENT(IN) :: c0,c1,c2
INTEGER(MMG5F_INT), INTENT(IN):: ref,pos
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 297 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Set_vertices | ( | MMG5_pMesh | mesh, |
double * | vertices, | ||
MMG5_int * | refs | ||
) |
Set the coordinates and references of all vertices in a mesh.
mesh | pointer to the mesh structure. |
vertices | array of vertex coordinates in the order \([x_1, y_1, z_1, x_2, \ldots, z_N]\) where \(N\) is the number of vertices. |
refs | array of references. The reference of vertex \(i\) is stored in refs[ \(i-1\)]. |
! SUBROUTINE MMGS_SET_VERTICES(mesh,vertices,refs,retval)
! MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
! REAL(KIND=8), DIMENSION(*),INTENT(IN) :: vertices
! INTEGER(MMG5F_INT),DIMENSION(*), INTENT(IN) :: refs
! INTEGER, INTENT(OUT) :: retval
! END SUBROUTINE
Definition at line 333 of file API_functions_s.c.
LIBMMGS_EXPORT void MMGS_setfunc | ( | MMG5_pMesh | mesh, |
MMG5_pSol | met | ||
) |
Set function pointers for caltet, lenedg, defsiz and gradsiz.
To associate function pointers without calling MMGS_mmgslib
mesh | pointer to the mesh structure (unused). |
met | pointer to the sol structure (unused). |
SUBROUTINE MMGS_SETFUNC(mesh,met)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
END SUBROUTINE
Definition at line 43 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_stockOptions | ( | MMG5_pMesh | mesh, |
MMG5_Info * | info | ||
) |
Store the info structure in the mesh structure.
mesh | pointer to the mesh structure. |
info | pointer to the info structure. |
SUBROUTINE MMGS_STOCKOPTIONS(mesh,info,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,info
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 552 of file libmmgs_tools.c.
LIBMMGS_EXPORT int MMGS_Unset_corner | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Remove the "corner" attribute from a vertex.
mesh | pointer to the mesh structure. |
k | vertex number. |
This function removes the corner attribute from vertex pos (from 1 to the number of vertices included).
SUBROUTINE MMGS_UNSET_CORNER(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 716 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Unset_requiredEdge | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Remove the "required" attribute from an edge.
mesh | pointer to the mesh structure. |
k | edge index. |
This function removes the "required" attribute from edge k.
SUBROUTINE MMGS_UNSET_REQUIREDEDGE(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 768 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Unset_requiredTriangle | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Remove the "required" attribute from a vertex.
mesh | pointer to the mesh structure. |
k | triangle index. |
This function removes the required attribute from triangle k.
SUBROUTINE MMGS_UNSET_REQUIREDTRIANGLE(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 743 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Unset_requiredVertex | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Remove the "required" attribute from a vertex.
mesh | pointer to the mesh structure. |
k | vertex number. |
This function removes the required attribute from vertex k.
SUBROUTINE MMGS_UNSET_REQUIREDVERTEX(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 729 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_Unset_ridge | ( | MMG5_pMesh | mesh, |
MMG5_int | k | ||
) |
Remove the "ridge" attribute from an edge.
mesh | pointer to the mesh structure. |
k | edge index. |
This function removes the ridge attribute from edge k.
SUBROUTINE MMGS_UNSET_RIDGE(mesh,k,retval)
MMG5_DATA_PTR_T,INTENT(INOUT) :: mesh
INTEGER(MMG5F_INT), INTENT(IN):: k
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 756 of file API_functions_s.c.
LIBMMGS_EXPORT int MMGS_usage | ( | char * | prog | ) |
Print help for mmgs options.
prog | pointer to the program name. |
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.
|
extern |
Compute an isotropic size map according to the mean of the length of the edges passing through a vertex.
mesh | pointer to the mesh structure |
met | pointer to the sol structure |
SUBROUTINE MMGS_DOSOL(mesh,met,retval)
MMG5_DATA_PTR_T, INTENT(INOUT) :: mesh,met
INTEGER, INTENT(OUT) :: retval
END SUBROUTINE
Definition at line 9 of file mmgsexterns.c.