Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
main.c
Go to the documentation of this file.
1/* =============================================================================
2** This file is part of the mmg software package for the
3** mesh modification.
4** Copyright (c) Bx INP/Inria/UBordeaux/UPMC, 2004- .
5**
6** mmg is free software: you can redistribute it and/or modify it
7** under the terms of the GNU Lesser General Public License as published
8** by the Free Software Foundation, either version 3 of the License, or
9** (at your option) any later version.
10**
11** mmg is distributed in the hope that it will be useful, but WITHOUT
12** ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13** FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
14** License for more details.
15**
16** You should have received a copy of the GNU Lesser General Public
17** License and of the GNU General Public License along with mmg (in
18** files COPYING.LESSER and COPYING). If not, see
19** <http://www.gnu.org/licenses/>. Please read their terms carefully and
20** use this copy of the mmg distribution only if you accept them.
21** =============================================================================
22*/
23
34
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 "libmmgs.h"
47// if the header file is in "include/mmg/mmgs"
48#include "mmg/mmgs/libmmgs.h"
49
50#define MAX0(a,b) (((a) > (b)) ? (a) : (b))
51#define MAX3(a,b,c) (((MAX0(a,b)) > c) ? (MAX0(a,b)) : c)
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, nt, na, nc, nr, nreq,ref,Tria[3], Edge[2];
61 int typEntity, typSol;
62 int *corner, *required, *ridge;
63 double Point[3],Sol;
64 char *fileout,*solout;
65
66 fprintf(stdout," -- TEST MMGSLIB \n");
67
68 if ( argc != 2 ) {
69 printf(" Usage: %s fileout\n",argv[0]);
70 return(1);
71 }
72
73 /* Name and path of the mesh file */
74 fileout = (char *) calloc(strlen(argv[1]) + 6, sizeof(char));
75 if ( fileout == NULL ) {
76 perror(" ## Memory problem: calloc");
77 exit(EXIT_FAILURE);
78 }
79 strcpy(fileout,argv[1]);
80 strcat(fileout,".mesh");
81
82 solout = (char *) calloc(strlen(argv[1]) + 5, sizeof(char));
83 if ( solout == NULL ) {
84 perror(" ## Memory problem: calloc");
85 exit(EXIT_FAILURE);
86 }
87 strcpy(solout,argv[1]);
88 strcat(solout,".sol");
89
92 /* args of InitMesh:
93 * MMG5_ARG_start: we start to give the args of a variadic func
94 * MMG5_ARG_ppMesh: next arg will be a pointer over a MMG5_pMesh
95 * &mmgMesh: pointer toward your MMG5_pMesh (that store your mesh)
96 * MMG5_ARG_ppMet: next arg will be a pointer over a MMG5_pSol storing a metric
97 * &mmgSol: pointer toward your MMG5_pSol (that store your metric) */
98 mmgMesh = NULL;
99 mmgSol = NULL;
100
102 MMG5_ARG_ppMesh,&mmgMesh,MMG5_ARG_ppMet,&mmgSol,
104
108
111 if ( MMGS_Set_meshSize(mmgMesh,12,20,0) != 1 ) exit(EXIT_FAILURE);
112
115 if ( MMGS_Set_vertex(mmgMesh,0 ,0 ,0 ,0, 1) != 1 ) exit(EXIT_FAILURE);
116 if ( MMGS_Set_vertex(mmgMesh,0.5,0 ,0 ,0, 2) != 1 ) exit(EXIT_FAILURE);
117 if ( MMGS_Set_vertex(mmgMesh,0.5,0 ,1 ,0, 3) != 1 ) exit(EXIT_FAILURE);
118 if ( MMGS_Set_vertex(mmgMesh,0 ,0 ,1 ,0, 4) != 1 ) exit(EXIT_FAILURE);
119 if ( MMGS_Set_vertex(mmgMesh,0 ,1 ,0 ,0, 5) != 1 ) exit(EXIT_FAILURE);
120 if ( MMGS_Set_vertex(mmgMesh,0.5,1 ,0 ,0, 6) != 1 ) exit(EXIT_FAILURE);
121 if ( MMGS_Set_vertex(mmgMesh,0.5,1 ,1 ,0, 7) != 1 ) exit(EXIT_FAILURE);
122 if ( MMGS_Set_vertex(mmgMesh,0 ,1 ,1 ,0, 8) != 1 ) exit(EXIT_FAILURE);
123 if ( MMGS_Set_vertex(mmgMesh,1 ,0 ,0 ,0, 9) != 1 ) exit(EXIT_FAILURE);
124 if ( MMGS_Set_vertex(mmgMesh,1 ,1 ,0 ,0, 10) != 1 ) exit(EXIT_FAILURE);
125 if ( MMGS_Set_vertex(mmgMesh,1 ,0 ,1 ,0, 11) != 1 ) exit(EXIT_FAILURE);
126 if ( MMGS_Set_vertex(mmgMesh,1 ,1 ,1 ,0, 12) != 1 ) exit(EXIT_FAILURE);
127
128
131 if ( MMGS_Set_triangle(mmgMesh, 1, 4, 8, 3, 1) != 1 ) exit(EXIT_FAILURE);
132 if ( MMGS_Set_triangle(mmgMesh, 1, 2, 4, 3, 2) != 1 ) exit(EXIT_FAILURE);
133 if ( MMGS_Set_triangle(mmgMesh, 8, 3, 7, 3, 3) != 1 ) exit(EXIT_FAILURE);
134 if ( MMGS_Set_triangle(mmgMesh, 5, 8, 6, 3, 4) != 1 ) exit(EXIT_FAILURE);
135 if ( MMGS_Set_triangle(mmgMesh, 5, 6, 2, 3, 5) != 1 ) exit(EXIT_FAILURE);
136 if ( MMGS_Set_triangle(mmgMesh, 5, 2, 1, 3, 6) != 1 ) exit(EXIT_FAILURE);
137 if ( MMGS_Set_triangle(mmgMesh, 5, 1, 8, 3, 7) != 1 ) exit(EXIT_FAILURE);
138 if ( MMGS_Set_triangle(mmgMesh, 7, 6, 8, 3, 8) != 1 ) exit(EXIT_FAILURE);
139 if ( MMGS_Set_triangle(mmgMesh, 4, 3, 8, 3, 9) != 1 ) exit(EXIT_FAILURE);
140 if ( MMGS_Set_triangle(mmgMesh, 2, 3, 4, 3,10) != 1 ) exit(EXIT_FAILURE);
141 if ( MMGS_Set_triangle(mmgMesh, 9, 3, 2, 4,11) != 1 ) exit(EXIT_FAILURE);
142 if ( MMGS_Set_triangle(mmgMesh, 11, 9, 12, 4,12) != 1 ) exit(EXIT_FAILURE);
143 if ( MMGS_Set_triangle(mmgMesh, 7, 11, 12, 4,13) != 1 ) exit(EXIT_FAILURE);
144 if ( MMGS_Set_triangle(mmgMesh, 6, 7, 10, 4,14) != 1 ) exit(EXIT_FAILURE);
145 if ( MMGS_Set_triangle(mmgMesh, 6, 10, 9, 4,15) != 1 ) exit(EXIT_FAILURE);
146 if ( MMGS_Set_triangle(mmgMesh, 6, 9, 2, 4,16) != 1 ) exit(EXIT_FAILURE);
147 if ( MMGS_Set_triangle(mmgMesh, 12, 10, 7, 4,17) != 1 ) exit(EXIT_FAILURE);
148 if ( MMGS_Set_triangle(mmgMesh, 12, 9, 10, 4,18) != 1 ) exit(EXIT_FAILURE);
149 if ( MMGS_Set_triangle(mmgMesh, 3, 11, 7, 4,19) != 1 ) exit(EXIT_FAILURE);
150 if ( MMGS_Set_triangle(mmgMesh, 9, 11, 3, 4,20) != 1 ) exit(EXIT_FAILURE);
151
152
156
160 if ( MMGS_Set_solSize(mmgMesh,mmgSol,MMG5_Vertex,12,MMG5_Scalar) != 1 )
161 exit(EXIT_FAILURE);
162
164 for(k=1 ; k<=12 ; k++) {
165 if ( MMGS_Set_scalarSol(mmgSol,0.5,k) != 1 ) exit(EXIT_FAILURE);
166 }
167
169 if ( MMGS_Chk_meshData(mmgMesh,mmgSol) != 1 ) exit(EXIT_FAILURE);
170
173 ier = MMGS_mmgslib(mmgMesh,mmgSol);
174
175 if ( ier == MMG5_STRONGFAILURE ) {
176 fprintf(stdout,"BAD ENDING OF MMGSLIB: UNABLE TO SAVE MESH\n");
177 return(ier);
178 } else if ( ier == MMG5_LOWFAILURE )
179 fprintf(stdout,"BAD ENDING OF MMGSLIB\n");
180
181
187
189 if( !(inm = fopen(fileout,"w")) ) {
190 fprintf(stderr," ** UNABLE TO OPEN OUTPUT MESH FILE.\n");
191 exit(EXIT_FAILURE);
192 }
193 fprintf(inm,"MeshVersionFormatted 2\n");
194 fprintf(inm,"\nDimension 3\n");
195
197 if ( MMGS_Get_meshSize(mmgMesh,&np,&nt,&na) !=1 ) exit(EXIT_FAILURE);
198
199 /* Table to know if a vertex is corner */
200 corner = (int*)calloc(np+1,sizeof(int));
201 if ( !corner ) {
202 perror(" ## Memory problem: calloc");
203 exit(EXIT_FAILURE);
204 }
205 /* Table to know if a vertex/tetra/tria/edge is required */
206 required = (int*)calloc(MAX3(np,nt,na)+1 ,sizeof(int));
207 if ( !required ) {
208 perror(" ## Memory problem: calloc");
209 exit(EXIT_FAILURE);
210 }
211 /* Table to know if a coponant is corner and/or required */
212 ridge = (int*)calloc(na+1 ,sizeof(int));
213 if ( !ridge ) {
214 perror(" ## Memory problem: calloc");
215 exit(EXIT_FAILURE);
216 }
217
218 nreq = 0; nc = 0;
219 fprintf(inm,"\nVertices\n%"MMG5_PRId"\n",np);
220 for(k=1; k<=np; k++) {
222 if ( MMGS_Get_vertex(mmgMesh,&(Point[0]),&(Point[1]),&(Point[2]),
223 &ref,&(corner[k]),&(required[k])) != 1 )
224 exit(EXIT_FAILURE);
225 fprintf(inm,"%.15lg %.15lg %.15lg %"MMG5_PRId" \n",Point[0],Point[1],Point[2],ref);
226 if ( corner[k] ) nc++;
227 if ( required[k] ) nreq++;
228 }
229 fprintf(inm,"\nCorners\n%"MMG5_PRId"\n",nc);
230 for(k=1; k<=np; k++) {
231 if ( corner[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
232 }
233 fprintf(inm,"\nRequiredVertices\n%"MMG5_PRId"\n",nreq);
234 for(k=1; k<=np; k++) {
235 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
236 }
237 free(corner);
238 corner = NULL;
239
240 nreq = 0;
241 fprintf(inm,"\nTriangles\n%"MMG5_PRId"\n",nt);
242 for(k=1; k<=nt; k++) {
244 if ( MMGS_Get_triangle(mmgMesh,&(Tria[0]),&(Tria[1]),&(Tria[2]),
245 &ref,&(required[k])) != 1 )
246 exit(EXIT_FAILURE);
247 fprintf(inm,"%"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" \n",
248 Tria[0],Tria[1],Tria[2],ref);
249 if ( required[k] ) nreq++;
250 }
251 fprintf(inm,"\nRequiredTriangles\n%"MMG5_PRId"\n",nreq);
252 for(k=1; k<=nt; k++) {
253 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
254 }
255
256 nreq = 0;nr = 0;
257 fprintf(inm,"\nEdges\n%"MMG5_PRId"\n",na);
258 for(k=1; k<=na; k++) {
260 if ( MMGS_Get_edge(mmgMesh,&(Edge[0]),&(Edge[1]),&ref,
261 &(ridge[k]),&(required[k])) != 1 ) exit(EXIT_FAILURE);
262 fprintf(inm,"%"MMG5_PRId" %"MMG5_PRId" %"MMG5_PRId" \n",Edge[0],Edge[1],ref);
263 if ( ridge[k] ) nr++;
264 if ( required[k] ) nreq++;
265 }
266 fprintf(inm,"\nRequiredEdges\n%"MMG5_PRId"\n",nreq);
267 for(k=1; k<=na; k++) {
268 if ( required[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
269 }
270 fprintf(inm,"\nRidges\n%"MMG5_PRId"\n",nr);
271 for(k=1; k<=na; k++) {
272 if ( ridge[k] ) fprintf(inm,"%"MMG5_PRId" \n",k);
273 }
274
275 fprintf(inm,"\nEnd\n");
276 fclose(inm);
277
278 free(required);
279 required = NULL;
280 free(ridge);
281 ridge = NULL;
282
284 if( !(inm = fopen(solout,"w")) ) {
285 fprintf(stderr," ** UNABLE TO OPEN OUTPUT SOL FILE.\n");
286 exit(EXIT_FAILURE);
287 }
288 fprintf(inm,"MeshVersionFormatted 2\n");
289 fprintf(inm,"\nDimension 3\n");
290
293 if ( MMGS_Get_solSize(mmgMesh,mmgSol,&typEntity,&np,&typSol) != 1 )
294 exit(EXIT_FAILURE);
295
296 if ( ( typEntity != MMG5_Vertex ) || ( typSol != MMG5_Scalar ) )
297 exit(EXIT_FAILURE);
298
299 fprintf(inm,"\nSolAtVertices\n%"MMG5_PRId"\n",np);
300 fprintf(inm,"1 1 \n\n");
301 for(k=1; k<=np; k++) {
303 if ( MMGS_Get_scalarSol(mmgSol,&Sol) != 1 ) exit(EXIT_FAILURE);
304 fprintf(inm,"%.15lg \n",Sol);
305 }
306 fprintf(inm,"\nEnd\n");
307 fclose(inm);
308
311 MMG5_ARG_ppMesh,&mmgMesh,MMG5_ARG_ppMet,&mmgSol,
313
314 free(fileout);
315 fileout = NULL;
316
317 free(solout);
318 solout = NULL;
319
320 return(ier);
321}
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.
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.
int MMGS_Init_mesh(const int starter,...)
Initialize a mesh structure and optionally the associated solution and metric structures.
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.
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.
int MMGS_Free_all(const int starter,...)
Deallocations before return.
int MMGS_Get_scalarSol(MMG5_pSol met, double *s)
Get the next element of a scalar solution structure.
int MMGS_Set_scalarSol(MMG5_pSol met, double s, MMG5_int pos)
Set a single element of a scalar solution structure.
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.
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.
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.
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.
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.
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.
int ier
#define MAX3(a, b, c)
Definition main.c:51
program main
Example for using mmglib (basic use)
Definition main.F90:6
int MMGS_mmgslib(MMG5_pMesh mesh, MMG5_pSol met)
Main "program" for mesh adaptation.
Definition libmmgs.c:545
#define MMG5_ARG_ppMesh
#define MMG5_ARG_end
#define MMG5_STRONGFAILURE
Definition libmmgtypes.h:65
#define MMG5_LOWFAILURE
Definition libmmgtypes.h:57
@ MMG5_Scalar
MMG5_Sol * MMG5_pSol
MMG5_Mesh * MMG5_pMesh
#define MMG5_ARG_start
Definition libmmgtypes.h:93
@ MMG5_Vertex
#define MMG5_ARG_ppMet