Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
main.c
Go to the documentation of this file.
1
35#include <assert.h>
36#include <stdio.h>
37#include <stdlib.h>
38#include <signal.h>
39#include <string.h>
40#include <ctype.h>
41#include <math.h>
42#include <float.h>
43
45// if the header file is in the "include" directory
46// #include "libmmg3d.h"
47// if the header file is in "include/mmg/mmg3d"
48#include "mmg/mmg3d/libmmg3d.h"
49
50#define MAX0(a,b) (((a) > (b)) ? (a) : (b))
51#define MAX4(a,b,c,d) (((MAX0(a,b)) > (MAX0(c,d))) ? (MAX0(a,b)) : (MAX0(c,d)))
52
53int main(int argc,char *argv[]) {
54 MMG5_pMesh mmgMesh;
55 MMG5_pSol mmgSol;
56 int ier;
57 /* To save final mesh in a file */
58 FILE* inm;
59 /* To manually recover the mesh */
60 MMG5_int k,np,ne,nt,na,nc,nr,nreq,ref,Tetra[4],Tria[3],Edge[2];
61 MMG5_int ktet[2];
62 int typEntity, typSol,iface[2];
63 int *corner, *required, *ridge;
64 double Point[3],Sol;
65 char *fileout, *solout;
66
67 fprintf(stdout," -- TEST MMG3DLIB \n");
68
69 if ( argc != 2 ) {
70 printf(" Usage: %s fileout\n",argv[0]);
71 return(1);
72 }
73
74 /* Name and path of the mesh file */
75 fileout = (char *) calloc(strlen(argv[1]) + 6, sizeof(char));
76 if ( fileout == NULL ) {
77 perror(" ## Memory problem: calloc");
78 exit(EXIT_FAILURE);
79 }
80 strcpy(fileout,argv[1]);
81 strcat(fileout,".mesh");
82
83 solout = (char *) calloc(strlen(argv[1]) + 5, sizeof(char));
84 if ( solout == NULL ) {
85 perror(" ## Memory problem: calloc");
86 exit(EXIT_FAILURE);
87 }
88 strcpy(solout,argv[1]);
89 strcat(solout,".sol");
90
91
92
95 /* args of InitMesh:
96 * MMG5_ARG_start: we start to give the args of a variadic func
97 * MMG5_ARG_ppMesh: next arg will be a pointer over a MMG5_pMesh
98 * &mmgMesh: pointer toward your MMG5_pMesh (that store your mesh)
99 * MMG5_ARG_ppMet: next arg will be a pointer over a MMG5_pSol storing a metric
100 * &mmgSol: pointer toward your MMG5_pSol (that store your metric) */
101 mmgMesh = NULL;
102 mmgSol = NULL;
103
105 MMG5_ARG_ppMesh,&mmgMesh,MMG5_ARG_ppMet,&mmgSol,
107
115 if ( MMG3D_Set_meshSize(mmgMesh,12,12,0,20,0,0) != 1 ) exit(EXIT_FAILURE);
116
119 if ( MMG3D_Set_vertex(mmgMesh,0 ,0 ,0 ,0, 1) != 1 ) exit(EXIT_FAILURE);
120 if ( MMG3D_Set_vertex(mmgMesh,0.5,0 ,0 ,0, 2) != 1 ) exit(EXIT_FAILURE);
121 if ( MMG3D_Set_vertex(mmgMesh,0.5,0 ,1 ,0, 3) != 1 ) exit(EXIT_FAILURE);
122 if ( MMG3D_Set_vertex(mmgMesh,0 ,0 ,1 ,0, 4) != 1 ) exit(EXIT_FAILURE);
123 if ( MMG3D_Set_vertex(mmgMesh,0 ,1 ,0 ,0, 5) != 1 ) exit(EXIT_FAILURE);
124 if ( MMG3D_Set_vertex(mmgMesh,0.5,1 ,0 ,0, 6) != 1 ) exit(EXIT_FAILURE);
125 if ( MMG3D_Set_vertex(mmgMesh,0.5,1 ,1 ,0, 7) != 1 ) exit(EXIT_FAILURE);
126 if ( MMG3D_Set_vertex(mmgMesh,0 ,1 ,1 ,0, 8) != 1 ) exit(EXIT_FAILURE);
127 if ( MMG3D_Set_vertex(mmgMesh,1 ,0 ,0 ,0, 9) != 1 ) exit(EXIT_FAILURE);
128 if ( MMG3D_Set_vertex(mmgMesh,1 ,1 ,0 ,0, 10) != 1 ) exit(EXIT_FAILURE);
129 if ( MMG3D_Set_vertex(mmgMesh,1 ,0 ,1 ,0, 11) != 1 ) exit(EXIT_FAILURE);
130 if ( MMG3D_Set_vertex(mmgMesh,1 ,1 ,1 ,0, 12) != 1 ) exit(EXIT_FAILURE);
131
134 if ( MMG3D_Set_tetrahedron(mmgMesh, 1, 4, 2, 8,1, 1) != 1 ) exit(EXIT_FAILURE);
135 if ( MMG3D_Set_tetrahedron(mmgMesh, 8, 3, 2, 7,1, 2) != 1 ) exit(EXIT_FAILURE);
136 if ( MMG3D_Set_tetrahedron(mmgMesh, 5, 2, 6, 8,1, 3) != 1 ) exit(EXIT_FAILURE);
137 if ( MMG3D_Set_tetrahedron(mmgMesh, 5, 8, 1, 2,1, 4) != 1 ) exit(EXIT_FAILURE);
138 if ( MMG3D_Set_tetrahedron(mmgMesh, 7, 2, 8, 6,1, 5) != 1 ) exit(EXIT_FAILURE);
139 if ( MMG3D_Set_tetrahedron(mmgMesh, 2, 4, 3, 8,1, 6) != 1 ) exit(EXIT_FAILURE);
140 if ( MMG3D_Set_tetrahedron(mmgMesh, 9, 2, 3, 7,2, 7) != 1 ) exit(EXIT_FAILURE);
141 if ( MMG3D_Set_tetrahedron(mmgMesh, 7, 11, 9, 12,2, 8) != 1 ) exit(EXIT_FAILURE);
142 if ( MMG3D_Set_tetrahedron(mmgMesh, 6, 9, 10, 7,2, 9) != 1 ) exit(EXIT_FAILURE);
143 if ( MMG3D_Set_tetrahedron(mmgMesh, 6, 7, 2, 9,2,10) != 1 ) exit(EXIT_FAILURE);
144 if ( MMG3D_Set_tetrahedron(mmgMesh, 12, 9, 7, 10,2,11) != 1 ) exit(EXIT_FAILURE);
145 if ( MMG3D_Set_tetrahedron(mmgMesh, 9, 3, 11, 7,2,12) != 1 ) exit(EXIT_FAILURE);
146
149 if ( MMG3D_Set_triangle(mmgMesh, 1, 4, 8, 3, 1) != 1 ) exit(EXIT_FAILURE);
150 if ( MMG3D_Set_triangle(mmgMesh, 1, 2, 4, 3, 2) != 1 ) exit(EXIT_FAILURE);
151 if ( MMG3D_Set_triangle(mmgMesh, 8, 3, 7, 3, 3) != 1 ) exit(EXIT_FAILURE);
152 if ( MMG3D_Set_triangle(mmgMesh, 5, 8, 6, 3, 4) != 1 ) exit(EXIT_FAILURE);
153 if ( MMG3D_Set_triangle(mmgMesh, 5, 6, 2, 3, 5) != 1 ) exit(EXIT_FAILURE);
154 if ( MMG3D_Set_triangle(mmgMesh, 5, 2, 1, 3, 6) != 1 ) exit(EXIT_FAILURE);
155 if ( MMG3D_Set_triangle(mmgMesh, 5, 1, 8, 3, 7) != 1 ) exit(EXIT_FAILURE);
156 if ( MMG3D_Set_triangle(mmgMesh, 7, 6, 8, 3, 8) != 1 ) exit(EXIT_FAILURE);
157 if ( MMG3D_Set_triangle(mmgMesh, 4, 3, 8, 3, 9) != 1 ) exit(EXIT_FAILURE);
158 if ( MMG3D_Set_triangle(mmgMesh, 2, 3, 4, 3,10) != 1 ) exit(EXIT_FAILURE);
159 if ( MMG3D_Set_triangle(mmgMesh, 9, 3, 2, 4,11) != 1 ) exit(EXIT_FAILURE);
160 if ( MMG3D_Set_triangle(mmgMesh, 11, 9, 12, 4,12) != 1 ) exit(EXIT_FAILURE);
161 if ( MMG3D_Set_triangle(mmgMesh, 7, 11, 12, 4,13) != 1 ) exit(EXIT_FAILURE);
162 if ( MMG3D_Set_triangle(mmgMesh, 6, 7, 10, 4,14) != 1 ) exit(EXIT_FAILURE);
163 if ( MMG3D_Set_triangle(mmgMesh, 6, 10, 9, 4,15) != 1 ) exit(EXIT_FAILURE);
164 if ( MMG3D_Set_triangle(mmgMesh, 6, 9, 2, 4,16) != 1 ) exit(EXIT_FAILURE);
165 if ( MMG3D_Set_triangle(mmgMesh, 12, 10, 7, 4,17) != 1 ) exit(EXIT_FAILURE);
166 if ( MMG3D_Set_triangle(mmgMesh, 12, 9, 10, 4,18) != 1 ) exit(EXIT_FAILURE);
167 if ( MMG3D_Set_triangle(mmgMesh, 3, 11, 7, 4,19) != 1 ) exit(EXIT_FAILURE);
168 if ( MMG3D_Set_triangle(mmgMesh, 9, 11, 3, 4,20) != 1 ) exit(EXIT_FAILURE);
169
170
178 if ( MMG3D_Set_solSize(mmgMesh,mmgSol,MMG5_Vertex,12,MMG5_Scalar) != 1 )
179 exit(EXIT_FAILURE);
180
182 for(k=1 ; k<=12 ; k++) {
183 if ( MMG3D_Set_scalarSol(mmgSol,0.5,k) != 1 ) exit(EXIT_FAILURE);
184 }
185
187 if ( MMG3D_Chk_meshData(mmgMesh,mmgSol) != 1 ) exit(EXIT_FAILURE);
188
191 ier = MMG3D_mmg3dlib(mmgMesh,mmgSol);
192
193 if ( ier == MMG5_STRONGFAILURE ) {
194 fprintf(stdout,"BAD ENDING OF MMG3DLIB: UNABLE TO SAVE MESH\n");
195 return(ier);
196 } else if ( ier == MMG5_LOWFAILURE )
197 fprintf(stdout,"BAD ENDING OF MMG3DLIB\n");
198
199
207 if( !(inm = fopen(fileout,"w")) ) {
208 fprintf(stderr," ** UNABLE TO OPEN OUTPUT MESH FILE.\n");
209 exit(EXIT_FAILURE);
210 }
211 fprintf(inm,"MeshVersionFormatted 2\n");
212 fprintf(inm,"\nDimension 3\n");
213
215 if ( MMG3D_Get_meshSize(mmgMesh,&np,&ne,NULL,&nt,NULL,&na) !=1 ) exit(EXIT_FAILURE);
216
217 /* Table to know if a vertex is corner */
218 corner = (int*)calloc(np+1,sizeof(int));
219 if ( !corner ) {
220 perror(" ## Memory problem: calloc");
221 exit(EXIT_FAILURE);
222 }
223 /* Table to know if a vertex/tetra/tria/edge is required */
224 required = (int*)calloc(MAX4(np,ne,nt,na)+1 ,sizeof(int));
225 if ( !required ) {
226 perror(" ## Memory problem: calloc");
227 exit(EXIT_FAILURE);
228 }
229 /* Table to know if an edge delimits a sharp angle */
230 ridge = (int*)calloc(na+1 ,sizeof(int));
231 if ( !ridge ) {
232 perror(" ## Memory problem: calloc");
233 exit(EXIT_FAILURE);
234 }
235
236 nreq = 0; nc = 0;
237 fprintf(inm,"\nVertices\n%"MMG5_PRId"\n",np);
238 for(k=1; k<=np; k++) {
240 if ( MMG3D_Get_vertex(mmgMesh,&(Point[0]),&(Point[1]),&(Point[2]),
241 &ref,&(corner[k]),&(required[k])) != 1 )
242 exit(EXIT_FAILURE);
243 fprintf(inm,"%.15lg %.15lg %.15lg %"MMG5_PRId" \n",Point[0],Point[1],Point[2],ref);
244 if ( corner[k] ) nc++;
245 if ( required[k] ) nreq++;
246 }
247 fprintf(inm,"\nCorners\n%"MMG5_PRId"\n",nc);
248 for(k=1; k<=np; k++) {
249 if ( corner[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
250 }
251 fprintf(inm,"\nRequiredVertices\n%"MMG5_PRId"\n",nreq);
252 for(k=1; k<=np; k++) {
253 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
254 }
255 free(corner);
256 corner = NULL;
257
258 nreq = 0;
259 fprintf(inm,"\nTriangles\n%"MMG5_PRId"\n",nt);
260 for(k=1; k<=nt; k++) {
262 if ( MMG3D_Get_triangle(mmgMesh,&(Tria[0]),&(Tria[1]),&(Tria[2]),
263 &ref,&(required[k])) != 1 )
264 exit(EXIT_FAILURE);
265 fprintf(inm,"%"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" \n",Tria[0],Tria[1],Tria[2],ref);
266 if ( required[k] ) nreq++;
267 }
268 fprintf(inm,"\nRequiredTriangles\n%"MMG5_PRId"\n",nreq);
269 for(k=1; k<=nt; k++) {
270 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
271 }
272
273 /* Facultative step : if you want to know with which tetrahedra a triangle is
274 * connected */
275 for(k=1; k<=nt; k++) {
276 ktet[0] = ktet[1] = 0;
277 iface[0] = iface[1] = 0;
278
279 if (! MMG3D_Get_tetFromTria(mmgMesh,k,ktet,iface) ) {
280 printf("Get tet from tria fail.\n");
281 return 0;
282 }
283 printf("Tria %"MMG5_PRId" is connected with tet %"MMG5_PRId" "
284 "(face %d) and %"MMG5_PRId" (face %d) \n",
285 k,ktet[0],iface[0],ktet[1],iface[1]);
286 }
287
288 nreq = 0;nr = 0;
289 fprintf(inm,"\nEdges\n%"MMG5_PRId"\n",na);
290 for(k=1; k<=na; k++) {
292 if ( MMG3D_Get_edge(mmgMesh,&(Edge[0]),&(Edge[1]),&ref,
293 &(ridge[k]),&(required[k])) != 1 ) exit(EXIT_FAILURE);
294 fprintf(inm,"%"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" \n",Edge[0],Edge[1],ref);
295 if ( ridge[k] ) nr++;
296 if ( required[k] ) nreq++;
297 }
298 fprintf(inm,"\nRequiredEdges\n%"MMG5_PRId"\n",nreq);
299 for(k=1; k<=na; k++) {
300 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
301 }
302 fprintf(inm,"\nRidges\n%"MMG5_PRId"\n",nr);
303 for(k=1; k<=na; k++) {
304 if ( ridge[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
305 }
306
307 nreq = 0;
308 fprintf(inm,"\nTetrahedra\n%"MMG5_PRId"\n",ne);
309 for(k=1; k<=ne; k++) {
311 if ( MMG3D_Get_tetrahedron(mmgMesh,&(Tetra[0]),&(Tetra[1]),&(Tetra[2]),&(Tetra[3]),
312 &ref,&(required[k])) != 1 ) exit(EXIT_FAILURE);
313 fprintf(inm,"%"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" \n",
314 Tetra[0],Tetra[1],Tetra[2],Tetra[3],ref);
315 if ( required[k] ) nreq++;
316 }
317 fprintf(inm,"\nRequiredTetrahedra\n%"MMG5_PRId"\n",nreq);
318 for(k=1; k<=ne; k++) {
319 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
320 }
321
322 fprintf(inm,"\nEnd\n");
323 fclose(inm);
324
325 free(required);
326 required = NULL;
327 free(ridge);
328 ridge = NULL;
329
331 if( !(inm = fopen(solout,"w")) ) {
332 fprintf(stderr," ** UNABLE TO OPEN OUTPUT FILE.\n");
333 exit(EXIT_FAILURE);
334 }
335 fprintf(inm,"MeshVersionFormatted 2\n");
336 fprintf(inm,"\nDimension 3\n");
337
340 if ( MMG3D_Get_solSize(mmgMesh,mmgSol,&typEntity,&np,&typSol) != 1 )
341 exit(EXIT_FAILURE);
342
343 if ( ( typEntity != MMG5_Vertex ) || ( typSol != MMG5_Scalar ) )
344 exit(EXIT_FAILURE);
345
346 fprintf(inm,"\nSolAtVertices\n%"MMG5_PRId"\n",np);
347 fprintf(inm,"1 1 \n\n");
348 for(k=1; k<=np; k++) {
350 if ( MMG3D_Get_scalarSol(mmgSol,&Sol) != 1 ) exit(EXIT_FAILURE);
351 fprintf(inm,"%.15lg \n",Sol);
352 }
353 fprintf(inm,"\nEnd\n");
354 fclose(inm);
355
358 MMG5_ARG_ppMesh,&mmgMesh,MMG5_ARG_ppMet,&mmgSol,
360
361 free(fileout);
362 fileout = NULL;
363
364 free(solout);
365 solout = NULL;
366
367
368 return(ier);
369}
int MMG3D_Init_mesh(const int starter,...)
Initialize a mesh structure and optionally the associated solution and metric structures.
int MMG3D_Get_edge(MMG5_pMesh mesh, MMG5_int *e0, MMG5_int *e1, MMG5_int *ref, int *isRidge, int *isRequired)
Get the vertices and reference of the next edge in the mesh.
int MMG3D_Get_triangle(MMG5_pMesh mesh, MMG5_int *v0, MMG5_int *v1, MMG5_int *v2, MMG5_int *ref, int *isRequired)
Get the vertices and reference of the next triangle in the mesh.
int MMG3D_Set_tetrahedron(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int v3, MMG5_int ref, MMG5_int pos)
set a single tetrahedron's vertices
int MMG3D_Set_triangle(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int ref, MMG5_int pos)
Set the vertices and reference of a single triangle in a mesh.
int MMG3D_Get_tetrahedron(MMG5_pMesh mesh, MMG5_int *v0, MMG5_int *v1, MMG5_int *v2, MMG5_int *v3, MMG5_int *ref, int *isRequired)
Get the vertices and reference of the next tetrahedron in the mesh.
int MMG3D_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 structure.
int MMG3D_Set_scalarSol(MMG5_pSol met, double s, MMG5_int pos)
Set a single element of a scalar solution structure.
int MMG3D_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.
int MMG3D_Chk_meshData(MMG5_pMesh mesh, MMG5_pSol met)
Check if the number of given entities match with mesh and sol size.
int MMG3D_Set_vertex(MMG5_pMesh mesh, double c0, double c1, double c2, MMG5_int ref, MMG5_int pos)
Set the coordinates of a single vertex.
int MMG3D_Get_scalarSol(MMG5_pSol met, double *s)
Get the next element of a scalar solution structure defined at vertices.
int MMG3D_Get_meshSize(MMG5_pMesh mesh, MMG5_int *np, MMG5_int *ne, MMG5_int *nprism, MMG5_int *nt, MMG5_int *nquad, MMG5_int *na)
Get the number of vertices, tetrahedra, prisms, triangles, quadrilaterals and edges of the mesh.
int MMG3D_Free_all(const int starter,...)
Deallocations before return.
int MMG3D_Set_meshSize(MMG5_pMesh mesh, MMG5_int np, MMG5_int ne, MMG5_int nprism, MMG5_int nt, MMG5_int nquad, MMG5_int na)
Set the number of vertices, tetrahedra, prisms, triangles, quadrilaterals, and edges of a mesh.
int MMG3D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
Initialize a solution field.
int ier
program main
Example for using mmglib (basic use)
Definition: main.F90:6
#define MAX4(a, b, c, d)
Definition: main.c:25
int MMG3D_mmg3dlib(MMG5_pMesh mesh, MMG5_pSol met)
Main "program" for the mesh adaptation library.
Definition: libmmg3d.c:975
LIBMMG3D_EXPORT int MMG3D_Get_tetFromTria(MMG5_pMesh mesh, MMG5_int ktri, MMG5_int *ktet, int *iface)
Get a tetrahedron given one of its triangles and the index by which it refers to this triangle (DEPRE...
#define MMG5_ARG_ppMesh
Definition: libmmgtypes.h:102
#define MMG5_ARG_end
Definition: libmmgtypes.h:179
#define MMG5_STRONGFAILURE
Definition: libmmgtypes.h:65
#define MMG5_LOWFAILURE
Definition: libmmgtypes.h:57
@ MMG5_Scalar
Definition: libmmgtypes.h:219
#define MMG5_ARG_start
Definition: libmmgtypes.h:93
@ MMG5_Vertex
Definition: libmmgtypes.h:230
#define MMG5_ARG_ppMet
Definition: libmmgtypes.h:122
MMG mesh structure.
Definition: libmmgtypes.h:613