gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
GModelIO_GEO.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #include <sstream>
7 #include <stdlib.h>
8 #include "GmshConfig.h"
9 #include "GmshMessage.h"
10 #include "GModel.h"
11 #include "GModelIO_GEO.h"
12 #include "Geo.h"
13 #include "GeoInterpolation.h"
14 #include "OS.h"
15 #include "OpenFile.h"
16 #include "Numeric.h"
17 #include "ListUtils.h"
18 #include "gmshVertex.h"
19 #include "gmshFace.h"
20 #include "gmshEdge.h"
21 #include "gmshRegion.h"
22 #include "Context.h"
23 #include "Parser.h"
24 
25 #if defined(HAVE_MESH)
26 #include "Field.h"
27 #endif
28 
30 {
33 
34  Points = Tree_Create(sizeof(Vertex *), CompareVertex);
35  Curves = Tree_Create(sizeof(Curve *), CompareCurve);
40 
41  PhysicalGroups = List_Create(5, 5, sizeof(PhysicalGroup *));
42  DelPhysicalGroups = List_Create(5, 5, sizeof(PhysicalGroup *));
43 
48 
49  _changed = true;
50 }
51 
53 {
56 
69 
78 
83 
84  _changed = true;
85 }
86 
87 void GEO_Internals::setMaxTag(int dim, int val)
88 {
89  switch(dim) {
90  case 0: _maxPointNum = val; break;
91  case 1: _maxLineNum = val; break;
92  case -1: _maxLineLoopNum = val; break;
93  case 2: _maxSurfaceNum = val; break;
94  case -2: _maxSurfaceLoopNum = val; break;
95  case 3: _maxVolumeNum = val; break;
96  }
97 }
98 
99 int GEO_Internals::getMaxTag(int dim) const
100 {
101  switch(dim) {
102  case 0: return _maxPointNum;
103  case 1: return _maxLineNum;
104  case -1: return _maxLineLoopNum;
105  case 2: return _maxSurfaceNum;
106  case -2: return _maxSurfaceLoopNum;
107  case 3: return _maxVolumeNum;
108  default: return 0;
109  }
110 }
111 
112 bool GEO_Internals::addVertex(int &tag, double x, double y, double z, double lc)
113 {
114  if(tag >= 0 && FindPoint(tag)) {
115  Msg::Error("GEO point with tag %d already exists", tag);
116  return false;
117  }
118  if(tag < 0) tag = getMaxTag(0) + 1;
119  if(!lc) lc = MAX_LC;
120  Vertex *v = CreateVertex(tag, x, y, z, lc, 1.0);
121  Tree_Add(Points, &v);
122  _changed = true;
123  return true;
124 }
125 
126 bool GEO_Internals::addVertex(int &tag, double x, double y,
127  gmshSurface *surface, double lc)
128 {
129  if(tag >= 0 && FindPoint(tag)) {
130  Msg::Error("GEO point with tag %d already exists", tag);
131  return false;
132  }
133  if(tag < 0) tag = getMaxTag(0) + 1;
134  if(!lc) lc = MAX_LC;
135  Vertex *v = CreateVertex(tag, x, y, surface, lc);
136  Tree_Add(Points, &v);
137  _changed = true;
138  return true;
139 }
140 
141 bool GEO_Internals::addLine(int &tag, int startTag, int endTag)
142 {
143  std::vector<int> points;
144  points.push_back(startTag);
145  points.push_back(endTag);
146  return addLine(tag, points);
147 }
148 
149 bool GEO_Internals::addLine(int &tag, const std::vector<int> &pointTags)
150 {
151  if(tag >= 0 && FindCurve(tag)) {
152  Msg::Error("GEO curve with tag %d already exists", tag);
153  return false;
154  }
155  if(pointTags.size() < 2) {
156  Msg::Error("Line requires 2 points");
157  return false;
158  }
159  if(tag < 0) tag = getMaxTag(1) + 1;
160  List_T *tmp = List_Create(2, 2, sizeof(int));
161  for(std::size_t i = 0; i < pointTags.size(); i++) {
162  int t = pointTags[i];
163  List_Add(tmp, &t);
164  }
165  bool ok = true;
166  Curve *c =
167  CreateCurve(tag, MSH_SEGM_LINE, 1, tmp, nullptr, -1, -1, 0., 1., ok);
168  Tree_Add(Curves, &c);
170  List_Delete(tmp);
171  _changed = true;
172  return ok;
173 }
174 
175 bool GEO_Internals::addCircleArc(int &tag, int startTag, int centerTag,
176  int endTag, double nx, double ny, double nz)
177 {
178  if(tag >= 0 && FindCurve(tag)) {
179  Msg::Error("GEO curve with tag %d already exists", tag);
180  return false;
181  }
182  if(tag < 0) tag = getMaxTag(1) + 1;
183  List_T *tmp = List_Create(3, 2, sizeof(int));
184  List_Add(tmp, &startTag);
185  List_Add(tmp, &centerTag);
186  List_Add(tmp, &endTag);
187  bool ok = true;
188  Curve *c =
189  CreateCurve(tag, MSH_SEGM_CIRC, 2, tmp, nullptr, -1, -1, 0., 1., ok);
190  if(nx || ny || nz) {
191  c->Circle.n[0] = nx;
192  c->Circle.n[1] = ny;
193  c->Circle.n[2] = nz;
194  EndCurve(c);
195  }
196  Tree_Add(Curves, &c);
197  Curve *rc = CreateReversedCurve(c);
198  if(nx || ny || nz) {
199  rc->Circle.n[0] = nx;
200  rc->Circle.n[1] = ny;
201  rc->Circle.n[2] = nz;
202  EndCurve(rc);
203  }
204  List_Delete(tmp);
205  _changed = true;
206  return ok;
207 }
208 
209 bool GEO_Internals::addEllipseArc(int &tag, int startTag, int centerTag,
210  int majorTag, int endTag, double nx,
211  double ny, double nz)
212 {
213  if(tag >= 0 && FindCurve(tag)) {
214  Msg::Error("GEO curve with tag %d already exists", tag);
215  return false;
216  }
217  if(tag < 0) tag = getMaxTag(1) + 1;
218  List_T *tmp = List_Create(3, 2, sizeof(int));
219  List_Add(tmp, &startTag);
220  List_Add(tmp, &centerTag);
221  List_Add(tmp, &majorTag);
222  List_Add(tmp, &endTag);
223  bool ok = true;
224  Curve *c =
225  CreateCurve(tag, MSH_SEGM_ELLI, 2, tmp, nullptr, -1, -1, 0., 1., ok);
226  if(nx || ny || nz) {
227  c->Circle.n[0] = nx;
228  c->Circle.n[1] = ny;
229  c->Circle.n[2] = nz;
230  EndCurve(c);
231  }
232  Tree_Add(Curves, &c);
233  Curve *rc = CreateReversedCurve(c);
234  if(nx || ny || nz) {
235  rc->Circle.n[0] = nx;
236  rc->Circle.n[1] = ny;
237  rc->Circle.n[2] = nz;
238  EndCurve(rc);
239  }
240  List_Delete(tmp);
241  _changed = true;
242  return ok;
243 }
244 
245 bool GEO_Internals::addSpline(int &tag, const std::vector<int> &pointTags)
246 {
247  if(tag >= 0 && FindCurve(tag)) {
248  Msg::Error("GEO curve with tag %d already exists", tag);
249  return false;
250  }
251  if(pointTags.size() < 2) {
252  Msg::Error("Spline curve requires at least 2 control points");
253  return false;
254  }
255  if(tag < 0) tag = getMaxTag(1) + 1;
256  List_T *tmp = List_Create(2, 2, sizeof(int));
257  for(std::size_t i = 0; i < pointTags.size(); i++) {
258  int t = pointTags[i];
259  List_Add(tmp, &t);
260  }
261  bool ok = true;
262  Curve *c =
263  CreateCurve(tag, MSH_SEGM_SPLN, 3, tmp, nullptr, -1, -1, 0., 1., ok);
264  Tree_Add(Curves, &c);
266  List_Delete(tmp);
267  _changed = true;
268  return ok;
269 }
270 
271 bool GEO_Internals::addBezier(int &tag, const std::vector<int> &pointTags)
272 {
273  if(tag >= 0 && FindCurve(tag)) {
274  Msg::Error("GEO curve with tag %d already exists", tag);
275  return false;
276  }
277  if(tag < 0) tag = getMaxTag(1) + 1;
278  if(pointTags.size() < 2) {
279  Msg::Error("Bezier curve requires at least 2 control points");
280  return false;
281  }
282  List_T *tmp = List_Create(2, 2, sizeof(int));
283  for(std::size_t i = 0; i < pointTags.size(); i++) {
284  int t = pointTags[i];
285  List_Add(tmp, &t);
286  }
287  bool ok = true;
288  Curve *c =
289  CreateCurve(tag, MSH_SEGM_BEZIER, 2, tmp, nullptr, -1, -1, 0., 1., ok);
290  Tree_Add(Curves, &c);
292  List_Delete(tmp);
293  _changed = true;
294  return ok;
295 }
296 
297 bool GEO_Internals::addBSpline(int &tag, const std::vector<int> &pointTags,
298  const std::vector<double> &seqknots)
299 {
300  if(tag >= 0 && FindCurve(tag)) {
301  Msg::Error("GEO curve with tag %d already exists", tag);
302  return false;
303  }
304  if(pointTags.size() < 2) {
305  Msg::Error("BSpline curve requires at least 2 control points");
306  return false;
307  }
308  if(tag < 0) tag = getMaxTag(1) + 1;
309  List_T *tmp = List_Create(2, 2, sizeof(int));
310  for(std::size_t i = 0; i < pointTags.size(); i++) {
311  int t = pointTags[i];
312  List_Add(tmp, &t);
313  }
314  Curve *c = nullptr;
315  bool ok = true;
316  if(seqknots.empty()) {
317  c = CreateCurve(tag, MSH_SEGM_BSPLN, 2, tmp, nullptr, -1, -1, 0., 1., ok);
318  }
319  else {
320  int order = seqknots.size() - pointTags.size() - 1;
321  List_T *knotsList = List_Create(2, 2, sizeof(double));
322  for(std::size_t i = 0; i < seqknots.size(); i++) {
323  double d = seqknots[i];
324  List_Add(knotsList, &d);
325  }
326  c = CreateCurve(tag, MSH_SEGM_NURBS, order, tmp, knotsList, -1, -1, 0., 1.,
327  ok);
328  }
329  Tree_Add(Curves, &c);
331  List_Delete(tmp);
332  _changed = true;
333  return ok;
334 }
335 
337  const std::vector<int> &curveTags,
338  int numIntervals, bool bspline)
339 {
340  if(tag >= 0 && FindCurve(tag)) {
341  Msg::Error("GEO curve with tag %d already exists", tag);
342  return false;
343  }
344  if(curveTags.empty()) {
345  Msg::Error("Compound spline curve requires at least 1 input curve");
346  return false;
347  }
348  if(numIntervals < 0) {
349  Msg::Error("Negative number of intervals in compound spline");
350  return false;
351  }
352  if(tag < 0) tag = getMaxTag(1) + 1;
353 
354  List_T *tmp =
355  List_Create((numIntervals + 1) * curveTags.size(), 2, sizeof(int));
356  for(std::size_t i = 0; i < curveTags.size(); i++) {
357  Curve *c = FindCurve(curveTags[i]);
358  if(!c) {
359  Msg::Error("Unknown GEO curve with tag %d", curveTags[i]);
360  return false;
361  }
362  if(i == 0 && c->beg) List_Add(tmp, &c->beg->Num);
363  for(int j = 1; j < numIntervals; j++) {
364  double u = (double)j / (double)(numIntervals);
365  Vertex V = InterpolateCurve(c, u, 0);
366  double lc = (1 - u) * c->beg->lc + u * c->end->lc;
367  int tag = getMaxTag(0) + 1;
368  Vertex *v = CreateVertex(tag, V.Pos.X, V.Pos.Y, V.Pos.Z, lc, 1.0);
369  Tree_Add(Points, &v);
370  List_Add(tmp, &v->Num);
371  }
372  if(c->end) List_Add(tmp, &c->end->Num);
373  }
374  bool ok = true;
375  Curve *c;
376  if(bspline)
377  c = CreateCurve(tag, MSH_SEGM_BSPLN, 2, tmp, nullptr, -1, -1, 0., 1., ok);
378  else // often too oscillatory for non-uniform distribution of control points
379  c = CreateCurve(tag, MSH_SEGM_SPLN, 3, tmp, nullptr, -1, -1, 0., 1., ok);
380  Tree_Add(Curves, &c);
382  List_Delete(tmp);
383  _changed = true;
384  return ok;
385 }
386 
388  const std::vector<int> &curveTags,
389  int numIntervals)
390 {
391  return _addCompoundSpline(tag, curveTags, numIntervals, false);
392 }
393 
395  const std::vector<int> &curveTags,
396  int numIntervals)
397 {
398  return _addCompoundSpline(tag, curveTags, numIntervals, true);
399 }
400 
401 bool GEO_Internals::addCurveLoop(int &tag, const std::vector<int> &curveTags,
402  bool reorient)
403 {
404  if(tag >= 0 && FindEdgeLoop(tag)) {
405  Msg::Error("GEO curve loop with tag %d already exists", tag);
406  return false;
407  }
408  if(tag < 0) tag = getMaxTag(-1) + 1;
409  List_T *tmp = List_Create(2, 2, sizeof(int));
410  for(std::size_t i = 0; i < curveTags.size(); i++) {
411  int t = curveTags[i];
412  List_Add(tmp, &t);
413  }
414  bool ok = SortEdgesInLoop(tag, tmp, reorient);
415  EdgeLoop *l = CreateEdgeLoop(tag, tmp);
416  Tree_Add(EdgeLoops, &l);
417  List_Delete(tmp);
418  _changed = true;
419  return ok;
420 }
421 
423  bool operator()(const Vertex *v1, const Vertex *v2) const
424  {
425  return v1->Num < v2->Num;
426  }
427 };
428 
429 static bool SortCurvesConsecutive(const std::vector<Curve *> &e,
430  std::vector<std::vector<Vertex *> > &vs)
431 {
432  vs.clear();
433  if(e.empty()) return true;
434  std::map<Vertex *, std::pair<Vertex *, Vertex *>, VertexPtrLessThan> c;
435 
436  for(size_t i = 0; i < e.size(); i++) {
437  Vertex *v0 = e[i]->beg;
438  Vertex *v1 = e[i]->end;
439 
440  if(!v0 || !v1) {
441  Msg::Warning("Skipping GEO curve %d without begin or end point in curve "
442  "loop detection",
443  e[i]->Num);
444  continue;
445  }
446 
447  if(v0 == v1) { // periodic curve
448  std::vector<Vertex *> v = {v0, v1};
449  vs.push_back(v);
450  continue;
451  }
452 
453  auto it0 = c.find(v0), it1 = c.find(v1);
454  if(it0 == c.end())
455  c[v0] = std::make_pair(v1, (Vertex *)nullptr);
456  else {
457  if(it0->second.second == nullptr) { it0->second.second = v1; }
458  else {
459  Msg::Debug("A list of curves has points that are adjacent to 3 curves");
460  return false;
461  }
462  }
463  if(it1 == c.end())
464  c[v1] = std::make_pair(v0, (Vertex *)nullptr);
465  else {
466  if(it1->second.second == nullptr) { it1->second.second = v0; }
467  else {
468  Msg::Debug("Wrong topology for a list of curves");
469  Msg::Debug("Point %d is adjacent to more than 2 points %d %d", v1->Num,
470  it1->second.first->Num, it1->second.second->Num);
471  return false;
472  }
473  }
474  }
475 
476  while(!c.empty()) {
477  std::vector<Vertex *> v;
478  Vertex *start = nullptr;
479  {
480  auto it = c.begin();
481  start = it->first;
482  for(; it != c.end(); ++it) {
483  if(it->second.second == nullptr) {
484  start = it->first;
485  break;
486  }
487  }
488  }
489 
490  auto its = c.find(start);
491 
492  Vertex *prev =
493  (its->second.second == start) ? its->second.first : its->second.second;
494  Vertex *current = start;
495 
496  do {
497  if(c.size() == 0) {
498  Msg::Warning("Wrong topology in a curve loop");
499  return false;
500  }
501  v.push_back(current);
502  auto it = c.find(current);
503  if(it == c.end() || it->first == nullptr) {
504  Msg::Error("Impossible to find point %d", current->Num);
505  return false;
506  }
507  Vertex *v1 = it->second.first;
508  Vertex *v2 = it->second.second;
509  c.erase(it);
510  Vertex *temp = current;
511  if(v1 == prev)
512  current = v2;
513  else if(v2 == prev)
514  current = v1;
515  else {
516  break;
517  }
518  prev = temp;
519  if(current == start) { v.push_back(current); }
520  } while(current != start && current != nullptr);
521  if(v.size() > 2 && v[v.size() - 2] == v[v.size() - 1]) {
522  v.erase(v.begin() + v.size() - 1);
523  }
524  vs.push_back(v);
525  }
526  return true;
527 }
528 
529 bool GEO_Internals::addCurveLoops(const std::vector<int> &curveTags,
530  std::vector<int> &curveLoopTags)
531 {
532  curveLoopTags.clear();
533  std::vector<Curve *> curves;
534  std::multimap<std::pair<Vertex *, Vertex *>, Curve *> pairs;
535  for(auto j : curveTags) {
536  Curve *c = FindCurve(j);
537  if(!c) {
538  Msg::Error("Unknown GEO curve %d", j);
539  return false;
540  }
541  if(!c->beg || !c->end) {
542  Msg::Error("Cannot create curve loops using curve %d without begin or "
543  "end point",
544  c->Num);
545  return false;
546  }
547  pairs.insert(std::make_pair(std::make_pair(c->beg, c->end), c));
548  curves.push_back(c);
549  }
550  std::vector<std::vector<Vertex *> > vs;
551  if(!SortCurvesConsecutive(curves, vs)) {
552  Msg::Error("Could not sort curves while creating curve loops");
553  return false;
554  }
555 
556  for(std::size_t i = 0; i < vs.size(); i++) {
557  if(vs[i].size() < 2 || vs[i][0] != vs[i][vs[i].size() - 1]) {
558  Msg::Warning("Skipping invalid loop with %lu points", vs[i].size());
559  continue;
560  }
561  else {
562  List_T *tmp = List_Create(2, 2, sizeof(int));
563  for(std::size_t j = 0; j < vs[i].size() - 1; j++) {
564  std::pair<Vertex *, Vertex *> p(vs[i][j], vs[i][j + 1]);
565  int num = 0;
566  auto it = pairs.find(p);
567  if(it != pairs.end()) {
568  num = it->second->Num;
569  pairs.erase(it);
570  }
571  else {
572  std::pair<Vertex *, Vertex *> p2(vs[i][j + 1], vs[i][j]);
573  it = pairs.find(p2);
574  if(it != pairs.end()) {
575  num = -it->second->Num;
576  pairs.erase(it);
577  }
578  }
579  if(num) List_Add(tmp, &num);
580  }
581  int tag = getMaxTag(-1) + 1;
582  EdgeLoop *l = CreateEdgeLoop(tag, tmp);
583  Tree_Add(EdgeLoops, &l);
584  curveLoopTags.push_back(tag);
585  List_Delete(tmp);
586  }
587  }
588  _changed = true;
589  return curveLoopTags.empty() ? false : true;
590 }
591 
592 bool GEO_Internals::addPlaneSurface(int &tag, const std::vector<int> &wireTags)
593 {
594  if(tag >= 0 && FindSurface(tag)) {
595  Msg::Error("GEO surface with tag %d already exists", tag);
596  return false;
597  }
598  if(tag < 0) tag = getMaxTag(2) + 1;
599  if(wireTags.empty()) {
600  Msg::Error("Plane surface requires at least one curve loop");
601  return false;
602  }
603  List_T *tmp = List_Create(2, 2, sizeof(int));
604  for(std::size_t i = 0; i < wireTags.size(); i++) {
605  int t = wireTags[i];
606  List_Add(tmp, &t);
607  }
609  bool ok = SetSurfaceGeneratrices(s, tmp);
610  List_Delete(tmp);
611  EndSurface(s);
612  Tree_Add(Surfaces, &s);
613  _changed = true;
614  return ok;
615 }
616 
618 {
619  if(tag >= 0 && FindSurface(tag)) {
620  Msg::Error("GEO surface with tag %d already exists", tag);
621  return false;
622  }
623  if(tag < 0) tag = getMaxTag(2) + 1;
625  Tree_Add(Surfaces, &s);
626  _changed = true;
627  return true;
628 }
629 
631  const std::vector<int> &wireTags,
632  int sphereCenterTag)
633 {
634  if(tag >= 0 && FindSurface(tag)) {
635  Msg::Error("GEO surface with tag %d already exists", tag);
636  return false;
637  }
638  if(tag < 0) tag = getMaxTag(2) + 1;
639  if(wireTags.empty()) {
640  Msg::Error("Surface requires at least one curve loop");
641  return false;
642  }
643  int ll = (int)std::abs(wireTags[0]);
644  EdgeLoop *el = FindEdgeLoop(ll);
645  if(!el) {
646  Msg::Error("Unknown curve loop %d", ll);
647  return false;
648  }
649  int j = List_Nbr(el->Curves), type = MSH_SURF_PLAN;
650  if(j == 4)
651  type = MSH_SURF_REGL;
652  else if(j == 3)
653  type = MSH_SURF_TRIC;
654  else {
655  Msg::Error("Wrong definition of surface %d: %d borders instead of 3 or 4",
656  tag, j);
657  return false;
658  }
659  List_T *tmp = List_Create(2, 2, sizeof(int));
660  for(std::size_t i = 0; i < wireTags.size(); i++) {
661  int t = wireTags[i];
662  List_Add(tmp, &t);
663  }
664  Surface *s = CreateSurface(tag, type);
665  bool ok = SetSurfaceGeneratrices(s, tmp);
666  List_Delete(tmp);
667  EndSurface(s);
668  if(sphereCenterTag >= 0) {
669  s->InSphereCenter = FindPoint(sphereCenterTag);
670  if(!s->InSphereCenter) {
671  Msg::Error("Unknown sphere center point %d", sphereCenterTag);
672  ok = false;
673  }
674  }
675  Tree_Add(Surfaces, &s);
676  _changed = true;
677  return ok;
678 }
679 
681  const std::vector<int> &surfaceTags)
682 {
683  if(tag >= 0 && FindSurfaceLoop(tag)) {
684  Msg::Error("GEO surface loop with tag %d already exists", tag);
685  return false;
686  }
687  if(tag < 0) tag = getMaxTag(-2) + 1;
688 
689  List_T *tmp = List_Create(2, 2, sizeof(int));
690  for(std::size_t i = 0; i < surfaceTags.size(); i++) {
691  int t = surfaceTags[i];
692  List_Add(tmp, &t);
693  }
694  SurfaceLoop *l = CreateSurfaceLoop(tag, tmp);
695  Tree_Add(SurfaceLoops, &l);
696  List_Delete(tmp);
697  _changed = true;
698  return true;
699 }
700 
701 bool GEO_Internals::addVolume(int &tag, const std::vector<int> &shellTags)
702 {
703  if(tag >= 0 && FindVolume(tag)) {
704  Msg::Error("GEO volume with tag %d already exists", tag);
705  return false;
706  }
707  if(tag < 0) tag = getMaxTag(3) + 1;
708 
709  List_T *tmp = List_Create(2, 2, sizeof(int));
710  for(std::size_t i = 0; i < shellTags.size(); i++) {
711  int t = shellTags[i];
712  List_Add(tmp, &t);
713  }
714  Volume *v = CreateVolume(tag, MSH_VOLUME);
715  bool ok = SetVolumeSurfaces(v, tmp);
716  List_Delete(tmp);
717  Tree_Add(Volumes, &v);
718  _changed = true;
719  return ok;
720 }
721 
723  const std::vector<std::pair<int, int> > &inDimTags,
724  double x, double y, double z, double dx, double dy,
725  double dz, double ax, double ay, double az,
726  double angle,
727  std::vector<std::pair<int, int> > &outDimTags,
728  ExtrudeParams *e)
729 {
730  List_T *in = List_Create(inDimTags.size() + 1, 10, sizeof(Shape));
731  List_T *out = List_Create(3 * inDimTags.size() + 1, 10, sizeof(Shape));
732 
733  for(std::size_t i = 0; i < inDimTags.size(); i++) {
734  int dim = inDimTags[i].first;
735  int tag = inDimTags[i].second;
736  Shape s;
737  s.Type = (dim == 3) ? MSH_VOLUME :
738  (dim == 2) ? MSH_SURF_PLAN :
739  (dim == 1) ? MSH_SEGM_LINE :
740  MSH_POINT;
741  s.Num = tag;
742  List_Add(in, &s);
743  }
744 
745  if(mode == 0) { // extrude
746  ExtrudeShapes(TRANSLATE, in, dx, dy, dz, 0., 0., 0., 0., 0., 0., 0., e,
747  out);
748  }
749  else if(mode == 1) { // revolve
750  ExtrudeShapes(ROTATE, in, 0., 0., 0., ax, ay, az, x, y, z, angle, e, out);
751  }
752  else if(mode == 2) { // extrude+revolve
753  ExtrudeShapes(TRANSLATE_ROTATE, in, dx, dy, dz, ax, ay, az, x, y, z, angle,
754  e, out);
755  }
756  else if(mode == 3) { // boundary layer
757  ExtrudeShapes(BOUNDARY_LAYER, in, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., e,
758  out);
759  }
760 
761  for(int i = 0; i < List_Nbr(out); i++) {
762  Shape s;
763  List_Read(out, i, &s);
764  int dim = s.Type / 100 - 1;
765  if(dim >= 0 && dim <= 3)
766  outDimTags.push_back(std::make_pair(dim, s.Num));
767  }
768  List_Delete(in);
769  List_Delete(out);
770  _changed = true;
771  return true;
772 }
773 
774 bool GEO_Internals::extrude(const std::vector<std::pair<int, int> > &inDimTags,
775  double dx, double dy, double dz,
776  std::vector<std::pair<int, int> > &outDimTags,
777  ExtrudeParams *e)
778 {
779  return _extrude(0, inDimTags, 0., 0., 0., dx, dy, dz, 0., 0., 0., 0.,
780  outDimTags, e);
781 }
782 
783 bool GEO_Internals::revolve(const std::vector<std::pair<int, int> > &inDimTags,
784  double x, double y, double z, double ax, double ay,
785  double az, double angle,
786  std::vector<std::pair<int, int> > &outDimTags,
787  ExtrudeParams *e)
788 {
789  return _extrude(1, inDimTags, x, y, z, 0., 0., 0., ax, ay, az, angle,
790  outDimTags, e);
791 }
792 
793 bool GEO_Internals::twist(const std::vector<std::pair<int, int> > &inDimTags,
794  double x, double y, double z, double dx, double dy,
795  double dz, double ax, double ay, double az,
796  double angle,
797  std::vector<std::pair<int, int> > &outDimTags,
798  ExtrudeParams *e)
799 {
800  return _extrude(2, inDimTags, x, y, z, dx, dy, dz, ax, ay, az, angle,
801  outDimTags, e);
802 }
803 
805  const std::vector<std::pair<int, int> > &inDimTags,
806  std::vector<std::pair<int, int> > &outDimTags, ExtrudeParams *e)
807 {
808  return _extrude(3, inDimTags, 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.,
809  outDimTags, e);
810 }
811 
813  const std::vector<std::pair<int, int> > &dimTags,
814  double x, double y, double z, double dx,
815  double dy, double dz, double a, double b,
816  double c, double d)
817 {
818  List_T *list = List_Create(dimTags.size() + 1, 10, sizeof(Shape));
819  for(std::size_t i = 0; i < dimTags.size(); i++) {
820  int dim = dimTags[i].first;
821  int tag = dimTags[i].second;
822  Shape s;
823  s.Type = (dim == 3) ? MSH_VOLUME :
824  (dim == 2) ? MSH_SURF_PLAN :
825  (dim == 1) ? MSH_SEGM_LINE :
826  MSH_POINT;
827  s.Num = tag;
828  List_Add(list, &s);
829  }
830  bool ok = true;
831  switch(mode) {
832  case 0: ok = TranslateShapes(dx, dy, dz, list); break;
833  case 1: ok = RotateShapes(dx, dy, dz, x, y, z, a, list); break;
834  case 2: ok = DilatShapes(x, y, z, a, b, c, list); break;
835  case 3: ok = SymmetryShapes(a, b, c, d, list); break;
836  }
837  List_Delete(list);
838  _changed = true;
839  return ok;
840 }
841 
842 bool GEO_Internals::translate(const std::vector<std::pair<int, int> > &dimTags,
843  double dx, double dy, double dz)
844 {
845  return _transform(0, dimTags, 0, 0, 0, dx, dy, dz, 0, 0, 0, 0);
846 }
847 
848 bool GEO_Internals::rotate(const std::vector<std::pair<int, int> > &dimTags,
849  double x, double y, double z, double ax, double ay,
850  double az, double angle)
851 {
852  return _transform(1, dimTags, x, y, z, ax, ay, az, angle, 0, 0, 0);
853 }
854 
855 bool GEO_Internals::dilate(const std::vector<std::pair<int, int> > &dimTags,
856  double x, double y, double z, double a, double b,
857  double c)
858 {
859  return _transform(2, dimTags, x, y, z, 0, 0, 0, a, b, c, 0);
860 }
861 
862 bool GEO_Internals::symmetry(const std::vector<std::pair<int, int> > &dimTags,
863  double a, double b, double c, double d)
864 {
865  return _transform(3, dimTags, 0, 0, 0, 0, 0, 0, a, b, c, d);
866 }
867 
868 bool GEO_Internals::splitCurve(int tag, const std::vector<int> &pointTags,
869  std::vector<int> &curveTags)
870 {
871  List_T *tmp = List_Create(10, 10, sizeof(int));
872  for(std::size_t i = 0; i < pointTags.size(); i++) {
873  int t = pointTags[i];
874  List_Add(tmp, &t);
875  }
876  List_T *curves = List_Create(10, 10, sizeof(Curve *));
877  bool ok = SplitCurve(tag, tmp, curves);
878  for(int i = 0; i < List_Nbr(curves); i++) {
879  Curve *c;
880  List_Read(curves, i, &c);
881  curveTags.push_back(c->Num);
882  }
883  List_Delete(tmp);
884  List_Delete(curves);
885  _changed = true;
886  return ok;
887 }
888 
890  const std::vector<int> &curveTags, int surfaceTag,
891  std::vector<int> &pointTags)
892 {
893  List_T *curves = List_Create(10, 10, sizeof(double));
894  List_T *shapes = List_Create(10, 10, sizeof(Shape));
895  for(std::size_t i = 0; i < curveTags.size(); i++) {
896  double d = curveTags[i];
897  List_Add(curves, &d);
898  }
899  bool ok = IntersectCurvesWithSurface(curves, surfaceTag, shapes);
900  for(int i = 0; i < List_Nbr(shapes); i++) {
901  Shape s;
902  List_Read(shapes, i, &s);
903  if(s.Type == MSH_POINT) { pointTags.push_back(s.Num); }
904  else {
905  Msg::Error("Degenerated curve-surface intersection not implemented");
906  return false;
907  }
908  }
909  _changed = true;
910  return ok;
911 }
912 
913 bool GEO_Internals::copy(const std::vector<std::pair<int, int> > &inDimTags,
914  std::vector<std::pair<int, int> > &outDimTags)
915 {
916  bool ok = true;
917  for(std::size_t i = 0; i < inDimTags.size(); i++) {
918  int dim = inDimTags[i].first;
919  int tag = inDimTags[i].second;
920  if(dim == 0) {
921  Vertex *v = FindPoint(tag);
922  if(!v) {
923  Msg::Error("Unknown GEO point %d", tag);
924  ok = false;
925  }
926  else {
927  Vertex *newv = DuplicateVertex(v);
928  outDimTags.push_back(std::make_pair(0, newv->Num));
929  }
930  }
931  else if(dim == 1) {
932  Curve *c = FindCurve(tag);
933  if(!c) {
934  Msg::Error("Unknown GEO curve %d", tag);
935  ok = false;
936  }
937  else {
938  Curve *newc = DuplicateCurve(c);
939  outDimTags.push_back(std::make_pair(1, newc->Num));
940  }
941  }
942  else if(dim == 2) {
943  Surface *s = FindSurface(tag);
944  if(!s) {
945  Msg::Error("Unknown GEO surface %d", tag);
946  ok = false;
947  }
948  else {
949  Surface *news = DuplicateSurface(s);
950  outDimTags.push_back(std::make_pair(2, news->Num));
951  }
952  }
953  else if(dim == 3) {
954  Volume *v = FindVolume(tag);
955  if(!v) {
956  Msg::Error("Unknown GEO volume %d", tag);
957  ok = false;
958  }
959  else {
960  Volume *newv = DuplicateVolume(v);
961  outDimTags.push_back(std::make_pair(3, newv->Num));
962  }
963  }
964  }
965  _changed = true;
966  return ok;
967 }
968 
969 bool GEO_Internals::remove(int dim, int tag, bool recursive)
970 {
971  switch(dim) {
972  case 0: DeletePoint(tag, recursive); break;
973  case 1:
974  DeleteCurve(tag, recursive);
975  DeleteCurve(-tag, recursive);
976  break;
977  case 2: DeleteSurface(tag, recursive); break;
978  case 3: DeleteVolume(tag, recursive); break;
979  }
980  _changed = true;
981  return true;
982 }
983 
984 bool GEO_Internals::remove(const std::vector<std::pair<int, int> > &dimTags,
985  bool recursive)
986 {
987  for(std::size_t i = 0; i < dimTags.size(); i++)
988  remove(dimTags[i].first, dimTags[i].second, recursive);
989  return true;
990 }
991 
993 {
997  _changed = true;
998 }
999 
1000 bool GEO_Internals::modifyPhysicalGroup(int dim, int tag, int op,
1001  const std::vector<int> &tags)
1002 {
1003  int type;
1004  std::string str;
1005  switch(dim) {
1006  case 0:
1007  type = MSH_PHYSICAL_POINT;
1008  str = "point";
1009  break;
1010  case 1:
1011  type = MSH_PHYSICAL_LINE;
1012  str = "curve";
1013  break;
1014  case 2:
1015  type = MSH_PHYSICAL_SURFACE;
1016  str = "surface";
1017  break;
1018  case 3:
1019  type = MSH_PHYSICAL_VOLUME;
1020  str = "volume";
1021  break;
1022  default: return false;
1023  }
1024 
1025  PhysicalGroup *p = FindPhysicalGroup(tag, type);
1026  if(p && op == 0) {
1027  Msg::Error("Physical %s %d already exists", str.c_str(), tag);
1028  return false;
1029  }
1030  else if(!p && op == 2) {
1031  // we call this in gmsh::model::removePhysicalGroup(), so it's not an error
1032  // if the group does not exist
1033  return true;
1034  }
1035  else if(!p && op > 0) {
1036  Msg::Error("Physical %s %d does not exist", str.c_str(), tag);
1037  return false;
1038  }
1039  else if(op == 0) {
1040  List_T *tmp = List_Create(10, 10, sizeof(int));
1041  for(std::size_t i = 0; i < tags.size(); i++) {
1042  int t = tags[i];
1043  List_Add(tmp, &t);
1044  }
1045  p = CreatePhysicalGroup(tag, type, tmp);
1046  List_Delete(tmp);
1047  List_Add(PhysicalGroups, &p);
1048  }
1049  else if(op == 1) {
1050  for(std::size_t i = 0; i < tags.size(); i++) {
1051  int t = tags[i];
1052  List_Add(p->Entities, &t);
1053  }
1054  }
1055  else if(op == 2) {
1056  for(std::size_t i = 0; i < tags.size(); i++) {
1057  int t = tags[i];
1058  List_Suppress(p->Entities, &t, fcmp_int);
1059  }
1060  if(!List_Nbr(p->Entities) || tags.empty()) {
1061  switch(dim) {
1062  case 0: DeletePhysicalPoint(tag); break;
1063  case 1: DeletePhysicalLine(tag); break;
1064  case 2: DeletePhysicalSurface(tag); break;
1065  case 3: DeletePhysicalVolume(tag); break;
1066  }
1067  }
1068  }
1069  else {
1070  Msg::Error("Unsupported operation on physical %s %d", str.c_str(), tag);
1071  return false;
1072  }
1073  _changed = true;
1074  return true;
1075 }
1076 
1078 {
1080  _changed = true;
1081 }
1082 
1083 bool GEO_Internals::mergeVertices(const std::vector<int> &tags)
1084 {
1085  if(tags.size() < 2) return true;
1086  Vertex *target = FindPoint(tags[0]);
1087  if(!target) {
1088  Msg::Error("Unknown GEO point %d", tags[0]);
1089  return false;
1090  }
1091 
1092  double x = target->Pos.X, y = target->Pos.Y, z = target->Pos.Z;
1093  for(std::size_t i = 1; i < tags.size(); i++) {
1094  Vertex *source = FindPoint(tags[i]);
1095  if(!source) {
1096  Msg::Error("Unknown GEO point %d", tags[i]);
1097  return false;
1098  }
1099  source->Typ = target->Typ;
1100  source->Pos.X = x;
1101  source->Pos.Y = y;
1102  source->Pos.Z = z;
1103  source->boundaryLayerIndex = target->boundaryLayerIndex;
1104  }
1105  ExtrudeParams::normalsCoherence.push_back(SPoint3(x, y, z));
1107  _changed = true;
1108  return true;
1109 }
1110 
1111 void GEO_Internals::setCompoundMesh(int dim, const std::vector<int> &tags)
1112 {
1113  _meshCompounds.insert(std::make_pair(dim, tags));
1114  _changed = true;
1115 }
1116 
1117 void GEO_Internals::setMeshSize(int dim, int tag, double size)
1118 {
1119  if(dim != 0) {
1120  Msg::Error("Setting mesh size only available on GEO points");
1121  return;
1122  }
1123  Vertex *v = FindPoint(tag);
1124  if(v) v->lc = size;
1125  _changed = true;
1126 }
1127 
1128 void GEO_Internals::setDegenerated(int dim, int tag)
1129 {
1130  if(dim != 1) return;
1131  Curve *c = FindCurve(tag);
1132  if(c) c->degenerated = true;
1133  _changed = true;
1134 }
1135 
1136 void GEO_Internals::setTransfiniteLine(int tag, int nPoints, int type,
1137  double coef)
1138 {
1139  if(!tag) {
1140  List_T *tmp = Tree2List(Curves);
1141  for(int i = 0; i < List_Nbr(tmp); i++) {
1142  Curve *c;
1143  List_Read(tmp, i, &c);
1144  c->Method = MESH_TRANSFINITE;
1145  c->nbPointsTransfinite = (nPoints > 2) ? nPoints : 2;
1146  c->typeTransfinite = type;
1147  c->coeffTransfinite = coef;
1148  }
1149  List_Delete(tmp);
1150  }
1151  else {
1152  Curve *c = FindCurve(tag);
1153  if(c) {
1154  c->Method = MESH_TRANSFINITE;
1155  c->nbPointsTransfinite = (nPoints > 2) ? nPoints : 2;
1156  c->typeTransfinite = type;
1157  c->coeffTransfinite = coef;
1158  }
1159  }
1160  _changed = true;
1161 }
1162 
1163 void GEO_Internals::setTransfiniteSurface(int tag, int arrangement,
1164  const std::vector<int> &cornerTags)
1165 {
1166  if(!tag) {
1167  List_T *tmp = Tree2List(Surfaces);
1168  for(int i = 0; i < List_Nbr(tmp); i++) {
1169  Surface *s;
1170  List_Read(tmp, i, &s);
1171  s->Method = MESH_TRANSFINITE;
1172  s->Recombine_Dir = arrangement;
1173  List_Reset(s->TrsfPoints);
1174  }
1175  List_Delete(tmp);
1176  }
1177  else {
1178  Surface *s = FindSurface(tag);
1179  if(s) {
1180  s->Method = MESH_TRANSFINITE;
1181  s->Recombine_Dir = arrangement;
1182  List_Reset(s->TrsfPoints);
1183  if(cornerTags.empty() || cornerTags.size() == 3 ||
1184  cornerTags.size() == 4) {
1185  for(std::size_t j = 0; j < cornerTags.size(); j++) {
1186  Vertex *v = FindPoint(std::abs(cornerTags[j]));
1187  if(v)
1188  List_Add(s->TrsfPoints, &v);
1189  else
1190  Msg::Error("Unknown GEO point %d", cornerTags[j]);
1191  }
1192  }
1193  else {
1194  Msg::Error("Transfinite surface requires 3 or 4 corner points");
1195  }
1196  }
1197  }
1198  _changed = true;
1199 }
1200 
1202  const std::vector<int> &cornerTags)
1203 {
1204  if(!tag) {
1205  List_T *tmp = Tree2List(Volumes);
1206  for(int i = 0; i < List_Nbr(tmp); i++) {
1207  Volume *v;
1208  List_Read(tmp, i, &v);
1209  v->Method = MESH_TRANSFINITE;
1210  List_Reset(v->TrsfPoints);
1211  }
1212  List_Delete(tmp);
1213  }
1214  else {
1215  Volume *v = FindVolume(tag);
1216  if(v) {
1217  v->Method = MESH_TRANSFINITE;
1218  List_Reset(v->TrsfPoints);
1219  if(cornerTags.empty() || cornerTags.size() == 6 ||
1220  cornerTags.size() == 8) {
1221  for(std::size_t i = 0; i < cornerTags.size(); i++) {
1222  Vertex *vert = FindPoint(std::abs(cornerTags[i]));
1223  if(vert)
1224  List_Add(v->TrsfPoints, &vert);
1225  else
1226  Msg::Error("Unknown GEO point %d", cornerTags[i]);
1227  }
1228  }
1229  }
1230  }
1231  _changed = true;
1232 }
1233 
1235 {
1236  if(!tag) {
1237  List_T *tmp = Tree2List(Volumes);
1238  for(int i = 0; i < List_Nbr(tmp); i++) {
1239  Volume *v;
1240  List_Read(tmp, i, &v);
1242  }
1243  List_Delete(tmp);
1244  }
1245  else {
1246  Volume *v = FindVolume(tag);
1247  if(v) v->QuadTri = TRANSFINITE_QUADTRI_1;
1248  }
1249  _changed = true;
1250 }
1251 
1252 void GEO_Internals::setRecombine(int dim, int tag, double angle)
1253 {
1254  if(dim == 2) {
1255  if(!tag) {
1256  List_T *tmp = Tree2List(Surfaces);
1257  for(int i = 0; i < List_Nbr(tmp); i++) {
1258  Surface *s;
1259  List_Read(tmp, i, &s);
1260  s->Recombine = 1;
1261  s->RecombineAngle = angle;
1262  }
1263  List_Delete(tmp);
1264  }
1265  else {
1266  Surface *s = FindSurface(tag);
1267  if(s) {
1268  s->Recombine = 1;
1269  s->RecombineAngle = angle;
1270  }
1271  }
1272  }
1273  else if(dim == 3) {
1274  if(!tag) {
1275  List_T *tmp = Tree2List(Volumes);
1276  for(int i = 0; i < List_Nbr(tmp); i++) {
1277  Volume *v;
1278  List_Read(tmp, i, &v);
1279  v->Recombine3D = 1;
1280  }
1281  List_Delete(tmp);
1282  }
1283  else {
1284  Volume *v = FindVolume(tag);
1285  if(v) { v->Recombine3D = 1; }
1286  }
1287  }
1288  _changed = true;
1289 }
1290 
1291 void GEO_Internals::setSmoothing(int tag, int val)
1292 {
1293  if(!tag) {
1294  List_T *tmp = Tree2List(Surfaces);
1295  for(int i = 0; i < List_Nbr(tmp); i++) {
1296  Surface *s;
1297  List_Read(tmp, i, &s);
1298  s->TransfiniteSmoothing = val;
1299  }
1300  List_Delete(tmp);
1301  }
1302  else {
1303  Surface *s = FindSurface(tag);
1304  if(s) s->TransfiniteSmoothing = val;
1305  }
1306  _changed = true;
1307 }
1308 
1309 void GEO_Internals::setReverseMesh(int dim, int tag, bool val)
1310 {
1311  if(dim == 1) {
1312  if(!tag) {
1313  List_T *tmp = Tree2List(Curves);
1314  for(int i = 0; i < List_Nbr(tmp); i++) {
1315  Curve *c;
1316  List_Read(tmp, i, &c);
1317  c->ReverseMesh = val ? 1 : 0;
1318  }
1319  List_Delete(tmp);
1320  }
1321  else {
1322  Curve *c = FindCurve(tag);
1323  if(c) c->ReverseMesh = val ? 1 : 0;
1324  }
1325  }
1326  else if(dim == 2) {
1327  if(!tag) {
1328  List_T *tmp = Tree2List(Surfaces);
1329  for(int i = 0; i < List_Nbr(tmp); i++) {
1330  Surface *s;
1331  List_Read(tmp, i, &s);
1332  s->ReverseMesh = val ? 1 : 0;
1333  }
1334  List_Delete(tmp);
1335  }
1336  else {
1337  Surface *s = FindSurface(tag);
1338  if(s) s->ReverseMesh = val ? 1 : 0;
1339  }
1340  }
1341  _changed = true;
1342 }
1343 
1344 void GEO_Internals::setMeshAlgorithm(int dim, int tag, int val)
1345 {
1346  if(dim == 2) {
1347  Surface *s = FindSurface(tag);
1348  if(s) s->MeshAlgorithm = val;
1349  }
1350  _changed = true;
1351 }
1352 
1353 void GEO_Internals::setMeshSizeFromBoundary(int dim, int tag, int val)
1354 {
1355  if(dim == 2) {
1356  Surface *s = FindSurface(tag);
1357  if(s) s->MeshSizeFromBoundary = val;
1358  }
1359  _changed = true;
1360 }
1361 
1362 bool sortEntities(const std::pair<int, int> &a,
1363  const std::pair<int, int> &b)
1364 {
1365  if(a.first != b.first)
1366  return a.first > b.first;
1367  return a.second < b.second;
1368 }
1369 
1370 void GEO_Internals::synchronize(GModel *model, bool resetMeshAttributes)
1371 {
1372  Msg::Debug("Syncing GEO_Internals with GModel");
1373 
1374  // if the entities do not exist in GModel, we create them; if they exist as
1375  // GEO entities in GModel but don't exist (anymore) in the internal CAD data,
1376  // we remove them. And if they exist in both the internal CAD data and in the
1377  // GModel, we update the pointer and the underlying dependencies (e.g. surface
1378  // boundaries): this is necessary because a GEO entity can change (while
1379  // keeping the same tag), due e.g. to ReplaceDuplicates.
1380  //
1381  // We also remove any entities of type "UnknownModel": these are discrete
1382  // entities, which are also stored in GEO_Internals so that they can be
1383  // combined with GEO entities; but they are not GmshModel entities.
1384  std::vector<std::pair<int, int> > toRemove;
1385  for(auto it = model->firstVertex(); it != model->lastVertex(); ++it) {
1386  GVertex *gv = *it;
1387  if(gv->getNativeType() == GEntity::GmshModel ||
1389  if(!FindPoint(gv->tag()))
1390  toRemove.push_back(std::make_pair(0, gv->tag()));
1391  }
1392  }
1393  for(auto it = model->firstEdge(); it != model->lastEdge(); ++it) {
1394  GEdge *ge = *it;
1395  if(ge->getNativeType() == GEntity::GmshModel ||
1397  if(!FindCurve(ge->tag()))
1398  toRemove.push_back(std::make_pair(1, ge->tag()));
1399  }
1400  }
1401  for(auto it = model->firstFace(); it != model->lastFace(); ++it) {
1402  GFace *gf = *it;
1403  if(gf->getNativeType() == GEntity::GmshModel ||
1405  if(!FindSurface(gf->tag()))
1406  toRemove.push_back(std::make_pair(2, gf->tag()));
1407  }
1408  }
1409  for(auto it = model->firstRegion(); it != model->lastRegion(); ++it) {
1410  GRegion *gr = *it;
1411  if(gr->getNativeType() == GEntity::GmshModel ||
1413  if(!FindVolume(gr->tag()))
1414  toRemove.push_back(std::make_pair(3, gr->tag()));
1415  }
1416  }
1417  std::sort(toRemove.begin(), toRemove.end(), sortEntities);
1418  Msg::Debug("Sync will try to remove %d model entities", toRemove.size());
1419 
1420  std::vector<GEntity*> removed;
1421  model->remove(toRemove, removed);
1422  Msg::Debug("Destroying %lu model entities during first pass", removed.size());
1423  for(std::size_t i = 0; i < removed.size(); i++) delete removed[i];
1424 
1425  if(Tree_Nbr(Points)) {
1426  List_T *points = Tree2List(Points);
1427  for(int i = 0; i < List_Nbr(points); i++) {
1428  Vertex *p;
1429  List_Read(points, i, &p);
1430  GVertex *v = model->getVertexByTag(p->Num);
1431  if(!v) {
1432  v = new gmshVertex(model, p);
1433  model->add(v);
1434  }
1435  else {
1436  if(v->getNativeType() == GEntity::GmshModel)
1437  ((gmshVertex *)v)->resetNativePtr(p);
1438  if(resetMeshAttributes) v->resetMeshAttributes();
1439  }
1440  }
1441  List_Delete(points);
1442  }
1443 
1444  if(Tree_Nbr(Curves)) {
1445  List_T *curves = Tree2List(Curves);
1446  for(int i = 0; i < List_Nbr(curves); i++) {
1447  Curve *c;
1448  List_Read(curves, i, &c);
1449  if(c->Num >= 0) {
1450  GEdge *e = model->getEdgeByTag(c->Num);
1451  GVertex *beg = nullptr, *end = nullptr;
1452  if(c->begByTag) beg = model->getVertexByTag(c->begByTag);
1453  if(!beg && c->beg) beg = model->getVertexByTag(c->beg->Num);
1454  if(c->endByTag) end = model->getVertexByTag(c->endByTag);
1455  if(!end && c->end) end = model->getVertexByTag(c->end->Num);
1456  if(!e && beg && end) {
1457  e = new gmshEdge(model, c, beg, end);
1458  model->add(e);
1459  }
1460  else if(!e) {
1461  e = new gmshEdge(model, c, nullptr, nullptr);
1462  model->add(e);
1463  }
1464  else {
1465  if(e->getNativeType() == GEntity::GmshModel) {
1466  if(beg && end)
1467  ((gmshEdge *)e)->resetNativePtr(c, beg, end);
1468  else
1469  ((gmshEdge *)e)->resetNativePtr(c, nullptr, nullptr);
1470  }
1471  if(resetMeshAttributes) e->resetMeshAttributes();
1472  }
1473  if(c->degenerated) e->setTooSmall(true);
1474  }
1475  }
1476  List_Delete(curves);
1477  }
1478 
1479  if(Tree_Nbr(Surfaces)) {
1480  List_T *surfaces = Tree2List(Surfaces);
1481  for(int i = 0; i < List_Nbr(surfaces); i++) {
1482  Surface *s;
1483  List_Read(surfaces, i, &s);
1484  GFace *f = model->getFaceByTag(s->Num);
1485  if(!f) {
1486  f = new gmshFace(model, s);
1487  model->add(f);
1488  }
1489  else {
1490  if(f->getNativeType() == GEntity::GmshModel)
1491  ((gmshFace *)f)->resetNativePtr(s);
1492  if(resetMeshAttributes) { f->resetMeshAttributes(); }
1493  }
1494  }
1495  List_Delete(surfaces);
1496  }
1497 
1498  if(Tree_Nbr(Volumes)) {
1499  List_T *volumes = Tree2List(Volumes);
1500  for(int i = 0; i < List_Nbr(volumes); i++) {
1501  Volume *v;
1502  List_Read(volumes, i, &v);
1503  GRegion *r = model->getRegionByTag(v->Num);
1504  if(!r) {
1505  r = new gmshRegion(model, v);
1506  model->add(r);
1507  }
1508  else {
1509  if(r->getNativeType() == GEntity::GmshModel)
1510  ((gmshRegion *)r)->resetNativePtr(v);
1511  if(resetMeshAttributes) r->resetMeshAttributes();
1512  }
1513  }
1514  List_Delete(volumes);
1515  }
1516 
1517  // In the first removal pass, entities which were on the boundary of a
1518  // higher-dimensional entity that is not removed, will not have been removed
1519  // due to the dependency that still existed in the topological GModel. After
1520  // recreation of the topology of these higher dimensional entities
1521  // (cf. resetNativePtr() above), we can try to re-remove the lower dimension
1522  // entities
1523  removed.clear();
1524  model->remove(toRemove, removed);
1525  Msg::Debug("Destroying %lu model entities during second pass", removed.size());
1526  for(std::size_t i = 0; i < removed.size(); i++) delete removed[i];
1527 
1528  // delete all physical groups before sync only if there is no mesh (if there
1529  // is a mesh, it could have been loaded from a file with physical groups - we
1530  // don't want to remove those)
1531  if(!model->getNumMeshElements()) model->removePhysicalGroups();
1532  // we might want to store physical groups directly in GModel; but I guess this
1533  // is OK for now:
1534  for(int i = 0; i < List_Nbr(PhysicalGroups); i++) {
1535  PhysicalGroup *p;
1536  List_Read(PhysicalGroups, i, &p);
1537  for(int j = 0; j < List_Nbr(p->Entities); j++) {
1538  int num;
1539  List_Read(p->Entities, j, &num);
1540  GEntity *ge = nullptr;
1541  const char *name = "";
1542  int tag = CTX::instance()->geom.orientedPhysicals ? abs(num) : num;
1543  switch(p->Typ) {
1544  case MSH_PHYSICAL_POINT:
1545  ge = model->getVertexByTag(tag);
1546  name = "point";
1547  break;
1548  case MSH_PHYSICAL_LINE:
1549  ge = model->getEdgeByTag(tag);
1550  name = "curve";
1551  break;
1552  case MSH_PHYSICAL_SURFACE:
1553  ge = model->getFaceByTag(tag);
1554  name = "surface";
1555  break;
1556  case MSH_PHYSICAL_VOLUME:
1557  ge = model->getRegionByTag(tag);
1558  name = "volume";
1559  break;
1560  }
1561  if(ge) {
1562  int pnum = CTX::instance()->geom.orientedPhysicals ?
1563  (gmsh_sign(num) * p->Num) : p->Num;
1564  if(std::find(ge->physicals.begin(), ge->physicals.end(), pnum) ==
1565  ge->physicals.end()) {
1566  ge->physicals.push_back(pnum);
1567  }
1568  }
1569  else {
1570  Msg::Warning("Skipping unknown %s %d in physical %s %d",
1571  name, tag, name, p->Num);
1572  }
1573  }
1574  }
1575 
1576  // we might want to store mesh compounds directly in GModel; but this is OK
1577  // for now.
1578  for(auto it = _meshCompounds.begin(); it != _meshCompounds.end(); ++it) {
1579  int dim = it->first;
1580  std::vector<int> compound = it->second;
1581  std::vector<GEntity *> ents;
1582  for(std::size_t i = 0; i < compound.size(); i++) {
1583  int tag = compound[i];
1584  GEntity *ent = nullptr;
1585  const char *name = "";
1586  switch(dim) {
1587  case 1: ent = model->getEdgeByTag(tag); name = "curve"; break;
1588  case 2: ent = model->getFaceByTag(tag); name = "surface"; break;
1589  case 3: ent = model->getRegionByTag(tag); name = "volume"; break;
1590  default:
1591  Msg::Error("Compound mesh constraint with dimension %d", dim);
1592  continue;
1593  }
1594  if(ent) {
1595  ents.push_back(ent);
1596  }
1597  else {
1598  Msg::Warning("Skipping unknown %d %d in compound", name, tag);
1599  }
1600  }
1601  for(std::size_t i = 0; i < ents.size(); i++) { ents[i]->compound = ents; }
1602  }
1603 
1604  // recompute global boundind box in CTX
1605  SetBoundingBox();
1606 
1607  Msg::Debug("GModel imported:");
1608  Msg::Debug("%d points", model->getNumVertices());
1609  Msg::Debug("%d curves", model->getNumEdges());
1610  Msg::Debug("%d surfaces", model->getNumFaces());
1611  Msg::Debug("%d volumes", model->getNumRegions());
1612 
1613  _changed = false;
1614 }
1615 
1616 bool GEO_Internals::getVertex(int tag, double &x, double &y, double &z)
1617 {
1618  Vertex *v = FindPoint(tag);
1619  if(v) {
1620  x = v->Pos.X;
1621  y = v->Pos.Y;
1622  z = v->Pos.Z;
1623  return true;
1624  }
1625  return false;
1626 }
1627 
1629  int pointTag)
1630 {
1631  Vertex *v1 = FindPoint(centerTag);
1632  if(!v1) {
1633  Msg::Error("Unknown sphere center point %d", centerTag);
1634  return nullptr;
1635  }
1636  Vertex *v2 = FindPoint(pointTag);
1637  if(!v2) {
1638  Msg::Error("Unknown sphere point %d", pointTag);
1639  return nullptr;
1640  }
1641  return gmshSphere::NewSphere(
1642  tag, v1->Pos.X, v1->Pos.Y, v1->Pos.Z,
1643  sqrt((v2->Pos.X - v1->Pos.X) * (v2->Pos.X - v1->Pos.X) +
1644  (v2->Pos.Y - v1->Pos.Y) * (v2->Pos.Y - v1->Pos.Y) +
1645  (v2->Pos.Z - v1->Pos.Z) * (v2->Pos.Z - v1->Pos.Z)));
1646 }
1647 
1649  int pointTag)
1650 {
1651  Vertex *v1 = FindPoint(centerTag);
1652  if(!v1) {
1653  Msg::Error("Unknown polar sphere center point %d", centerTag);
1654  return nullptr;
1655  }
1656  Vertex *v2 = FindPoint(pointTag);
1657  if(!v2) {
1658  Msg::Error("Unknown polar sphere point %d", pointTag);
1659  return nullptr;
1660  }
1662  tag, v1->Pos.X, v1->Pos.Y, v1->Pos.Z,
1663  sqrt((v2->Pos.X - v1->Pos.X) * (v2->Pos.X - v1->Pos.X) +
1664  (v2->Pos.Y - v1->Pos.Y) * (v2->Pos.Y - v1->Pos.Y) +
1665  (v2->Pos.Z - v1->Pos.Z) * (v2->Pos.Z - v1->Pos.Z)));
1666 }
1667 
1668 // GModel interface
1669 
1671 
1673 {
1674  if(_geo_internals) delete _geo_internals;
1675  _geo_internals = nullptr;
1676 }
1677 
1678 #if defined(HAVE_MESH)
1679 
1680 class writeFieldOptionGEO {
1681 private:
1682  FILE *geo;
1683  Field *field;
1684 
1685 public:
1686  writeFieldOptionGEO(FILE *fp, Field *_field)
1687  {
1688  geo = fp ? fp : stdout;
1689  field = _field;
1690  }
1691  void operator()(std::pair<std::string, FieldOption *> it)
1692  {
1693  std::string v;
1694  it.second->getTextRepresentation(v);
1695  fprintf(geo, "Field[%i].%s = %s;\n", field->id, it.first.c_str(),
1696  v.c_str());
1697  }
1698 };
1699 
1700 class writeFieldGEO {
1701 private:
1702  FILE *geo;
1703 
1704 public:
1705  writeFieldGEO(FILE *fp) { geo = fp ? fp : stdout; }
1706  void operator()(std::pair<const int, Field *> it)
1707  {
1708  fprintf(geo, "Field[%i] = %s;\n", it.first, it.second->getName());
1709  std::for_each(it.second->options.begin(), it.second->options.end(),
1710  writeFieldOptionGEO(geo, it.second));
1711  }
1712 };
1713 
1714 #endif
1715 
1717 private:
1718  FILE *geo;
1719  int dim;
1721  std::map<int, std::string> &oldLabels;
1722  std::map<std::pair<int, int>, std::string> &newLabels;
1723 
1724 public:
1725  writePhysicalGroupGEO(FILE *fp, int i, bool labels,
1726  std::map<int, std::string> &o,
1727  std::map<std::pair<int, int>, std::string> &n)
1728  : dim(i), printLabels(labels), oldLabels(o), newLabels(n)
1729  {
1730  geo = fp ? fp : stdout;
1731  }
1732  void operator()(std::pair<const int, std::vector<GEntity *> > &g)
1733  {
1734  std::string oldName, newName;
1735  if(printLabels) {
1736  if(newLabels.count(std::make_pair(dim, g.first))) {
1737  newName = newLabels[std::make_pair(dim, g.first)];
1738  }
1739  else if(oldLabels.count(g.first)) {
1740  oldName = oldLabels[g.first];
1741  fprintf(geo, "%s = %d;\n", oldName.c_str(), g.first);
1742  }
1743  }
1744 
1745  switch(dim) {
1746  case 0: fprintf(geo, "Physical Point"); break;
1747  case 1: fprintf(geo, "Physical Curve"); break;
1748  case 2: fprintf(geo, "Physical Surface"); break;
1749  case 3: fprintf(geo, "Physical Volume"); break;
1750  }
1751 
1752  if(oldName.size())
1753  fprintf(geo, "(%s) = {", oldName.c_str());
1754  else if(newName.size())
1755  fprintf(geo, "(\"%s\") = {", newName.c_str());
1756  else
1757  fprintf(geo, "(%d) = {", g.first);
1758  for(std::size_t i = 0; i < g.second.size(); i++) {
1759  if(i) fprintf(geo, ", ");
1760  fprintf(geo, "%d", g.second[i]->tag());
1761  }
1762  fprintf(geo, "};\n");
1763  }
1764 };
1765 
1766 static bool skipRegion(GRegion *gr)
1767 {
1768  if(gr->physicals.size()) return false;
1769  return true;
1770 }
1771 
1772 static bool skipFace(GFace *gf)
1773 {
1774  if(gf->physicals.size()) return false;
1775  std::list<GRegion *> regions(gf->regions());
1776  for(auto itr = regions.begin(); itr != regions.end(); itr++) {
1777  if(!skipRegion(*itr)) return false;
1778  }
1779  return true;
1780 }
1781 
1782 static bool skipEdge(GEdge *ge)
1783 {
1784  if(ge->physicals.size()) return false;
1785  std::vector<GFace *> faces(ge->faces());
1786  for(auto itf = faces.begin(); itf != faces.end(); itf++) {
1787  if(!skipFace(*itf)) return false;
1788  }
1789  return true;
1790 }
1791 
1792 static bool skipVertex(GVertex *gv)
1793 {
1794  if(gv->physicals.size()) return false;
1795  std::vector<GEdge *> const &edges = gv->edges();
1796  for(auto ite = edges.begin(); ite != edges.end(); ite++) {
1797  if(!skipEdge(*ite)) return false;
1798  }
1799  return true;
1800 }
1801 
1802 int GModel::writeGEO(const std::string &name, bool printLabels,
1803  bool onlyPhysicals)
1804 {
1805  FILE *fp = Fopen(name.c_str(), "w");
1806  if(!fp) {
1807  Msg::Error("Could not open file '%s'", name.c_str());
1808  return 0;
1809  }
1810 
1811  std::map<double, std::string> meshSizeParameters;
1812  int cpt = 0;
1813  for(auto it = firstVertex(); it != lastVertex(); it++) {
1814  double val = (*it)->prescribedMeshSizeAtVertex();
1815  if(meshSizeParameters.find(val) == meshSizeParameters.end()) {
1816  std::ostringstream paramName;
1817  paramName << "cl__" << ++cpt;
1818  fprintf(fp, "%s = %.16g;\n", paramName.str().c_str(), val);
1819  meshSizeParameters.insert(std::make_pair(val, paramName.str()));
1820  }
1821  }
1822 
1823  for(auto it = firstVertex(); it != lastVertex(); it++) {
1824  double val = (*it)->prescribedMeshSizeAtVertex();
1825  if(!onlyPhysicals || !skipVertex(*it))
1826  (*it)->writeGEO(fp, meshSizeParameters[val]);
1827  }
1828  for(auto it = firstEdge(); it != lastEdge(); it++) {
1829  if(!onlyPhysicals || !skipEdge(*it)) (*it)->writeGEO(fp);
1830  }
1831  for(auto it = firstFace(); it != lastFace(); it++) {
1832  if(!onlyPhysicals || !skipFace(*it)) (*it)->writeGEO(fp);
1833  }
1834  for(auto it = firstRegion(); it != lastRegion(); it++) {
1835  if(!onlyPhysicals || !skipRegion(*it)) (*it)->writeGEO(fp);
1836  }
1837 
1838  std::map<int, std::string> labels;
1839 #if defined(HAVE_PARSER)
1840  // get "old-style" labels from parser
1841  for(auto it = gmsh_yysymbols.begin(); it != gmsh_yysymbols.end(); ++it)
1842  for(std::size_t i = 0; i < it->second.value.size(); i++)
1843  labels[(int)it->second.value[i]] = it->first;
1844 #endif
1845 
1846  std::map<int, std::vector<GEntity *> > groups[4];
1847  getPhysicalGroups(groups);
1848  for(int i = 0; i < 4; i++)
1849  std::for_each(
1850  groups[i].begin(), groups[i].end(),
1851  writePhysicalGroupGEO(fp, i, printLabels, labels, _physicalNames));
1852 
1853 #if defined(HAVE_MESH)
1854  std::for_each(getFields()->begin(), getFields()->end(), writeFieldGEO(fp));
1855  if(getFields()->getBackgroundField() > 0)
1856  fprintf(fp, "Background Field = %i;\n", getFields()->getBackgroundField());
1857 #endif
1858 
1859  fclose(fp);
1860  return 1;
1861 }
1862 
1863 int GModel::writePY(const std::string &name, bool printLabels,
1864  bool onlyPhysicals)
1865 {
1866  FILE *fp = Fopen(name.c_str(), "w");
1867  if(!fp) {
1868  Msg::Error("Could not open file '%s'", name.c_str());
1869  return 0;
1870  }
1871 
1872  fprintf(fp, "import gmsh\n");
1873  fprintf(fp, "gmsh.initialize()\n");
1874 
1875  std::map<double, std::string> meshSizeParameters;
1876  int cpt = 0;
1877  for(auto it = firstVertex(); it != lastVertex(); it++) {
1878  double val = (*it)->prescribedMeshSizeAtVertex();
1879  if(meshSizeParameters.find(val) == meshSizeParameters.end()) {
1880  std::ostringstream paramName;
1881  paramName << "cl__" << ++cpt;
1882  fprintf(fp, "%s = %.16g\n", paramName.str().c_str(), val);
1883  meshSizeParameters.insert(std::make_pair(val, paramName.str()));
1884  }
1885  }
1886 
1887  for(auto it = firstVertex(); it != lastVertex(); it++) {
1888  double val = (*it)->prescribedMeshSizeAtVertex();
1889  if(!onlyPhysicals || !skipVertex(*it))
1890  (*it)->writePY(fp, meshSizeParameters[val]);
1891  }
1892  for(auto it = firstEdge(); it != lastEdge(); it++) {
1893  if(!onlyPhysicals || !skipEdge(*it)) (*it)->writePY(fp);
1894  }
1895  for(auto it = firstFace(); it != lastFace(); it++) {
1896  if(!onlyPhysicals || !skipFace(*it)) (*it)->writePY(fp);
1897  }
1898  for(auto it = firstRegion(); it != lastRegion(); it++) {
1899  if(!onlyPhysicals || !skipRegion(*it)) (*it)->writePY(fp);
1900  }
1901 
1902  fprintf(fp, "gmsh.model.geo.synchronize()\n");
1903  fprintf(fp, "gmsh.fltk.run()\n");
1904  fprintf(fp, "gmsh.finalize()\n");
1905 
1906  fclose(fp);
1907  return 1;
1908 }
1909 
1911 {
1912  if(_geo_internals) delete _geo_internals;
1914 
1915  for(auto it = firstVertex(); it != lastVertex(); it++) {
1916  Vertex *v = CreateVertex((*it)->tag(), (*it)->x(), (*it)->y(), (*it)->z(),
1917  (*it)->prescribedMeshSizeAtVertex(), 1.0);
1919  }
1920 
1921  for(auto it = firstEdge(); it != lastEdge(); it++) {
1922  if((*it)->geomType() == GEntity::DiscreteCurve) {
1923  bool ok = true;
1924  Curve *c = CreateCurve((*it)->tag(), MSH_SEGM_DISCRETE, 1, nullptr,
1925  nullptr, -1, -1, 0., 1., ok);
1926  c->Control_Points = List_Create(2, 1, sizeof(Vertex *));
1927  GVertex *gvb = (*it)->getBeginVertex();
1928  if(gvb) {
1929  Vertex *v = FindPoint(gvb->tag());
1930  if(v) {
1931  List_Add(c->Control_Points, &v);
1932  c->beg = v;
1933  }
1934  else {
1935  Msg::Error("Unknown GEO point %d", gvb->tag());
1936  }
1937  }
1938  else {
1939  Msg::Warning("Discrete curve %d has no begin point", (*it)->tag());
1940  }
1941  GVertex *gve = (*it)->getEndVertex();
1942  if(gve) {
1943  Vertex *v = FindPoint(gve->tag());
1944  if(v) {
1945  List_Add(c->Control_Points, &v);
1946  c->end = v;
1947  }
1948  else {
1949  Msg::Error("Unknown GEO point %d", gve->tag());
1950  }
1951  }
1952  else {
1953  Msg::Warning("Discrete curve %d has no end point", (*it)->tag());
1954  }
1955  EndCurve(c);
1958  }
1959  }
1960 
1961  for(auto it = firstFace(); it != lastFace(); it++) {
1962  if((*it)->geomType() == GEntity::DiscreteSurface) {
1963  Surface *s = CreateSurface((*it)->tag(), MSH_SURF_DISCRETE);
1964  std::vector<GEdge *> const &edges = (*it)->edges();
1965  s->Generatrices = List_Create(edges.size() + 1, 1, sizeof(Curve *));
1966  for(auto ite = edges.begin(); ite != edges.end(); ite++) {
1967  Curve *c = FindCurve((*ite)->tag());
1968  if(c) { List_Add(s->Generatrices, &c); }
1969  else {
1970  Msg::Error("Unknown GEO curve %d", (*ite)->tag());
1971  }
1972  }
1973  EndSurface(s);
1975  }
1976  }
1977 
1978  for(auto it = firstRegion(); it != lastRegion(); it++) {
1979  if((*it)->geomType() == GEntity::DiscreteVolume) {
1980  Volume *v = CreateVolume((*it)->tag(), MSH_VOLUME_DISCRETE);
1981  std::vector<GFace *> faces = (*it)->faces();
1982  v->Surfaces = List_Create(faces.size() + 1, 1, sizeof(Surface *));
1983  for(auto itf = faces.begin(); itf != faces.end(); itf++) {
1984  Surface *s = FindSurface((*itf)->tag());
1985  if(s) { List_Add(v->Surfaces, &s); }
1986  else {
1987  Msg::Error("Unknown GEO surface %d", (*itf)->tag());
1988  }
1989  }
1991  }
1992  }
1993 
1994  Msg::Debug("Geo internal model has:");
1995  Msg::Debug("%d points", Tree_Nbr(_geo_internals->Points));
1996  Msg::Debug("%d curves", Tree_Nbr(_geo_internals->Curves));
1997  Msg::Debug("%d surfaces", Tree_Nbr(_geo_internals->Surfaces));
1998  Msg::Debug("%d volumes", Tree_Nbr(_geo_internals->Volumes));
1999 
2000  return 1;
2001 }
Volume
Definition: Geo.h:140
GEO_Internals::Curves
Tree_T * Curves
Definition: GModelIO_GEO.h:18
MSH_SEGM_DISCRETE
#define MSH_SEGM_DISCRETE
Definition: GeoDefines.h:31
Surface::TrsfPoints
List_T * TrsfPoints
Definition: Geo.h:122
GEO_Internals::removeAllDuplicates
void removeAllDuplicates()
Definition: GModelIO_GEO.cpp:1077
Surface::Num
int Num
Definition: Geo.h:113
Geo.h
GeoInterpolation.h
GEO_Internals::PhysicalGroups
List_T * PhysicalGroups
Definition: GModelIO_GEO.h:20
tags
static std::map< SPoint2, unsigned int > tags
Definition: drawGraph2d.cpp:400
MSH_SEGM_BEZIER
#define MSH_SEGM_BEZIER
Definition: GeoDefines.h:29
GEO_Internals::addBSpline
bool addBSpline(int &tag, const std::vector< int > &pointTags, const std::vector< double > &seqknots=std::vector< double >())
Definition: GModelIO_GEO.cpp:297
GEO_Internals::Volumes
Tree_T * Volumes
Definition: GModelIO_GEO.h:18
GModel::firstEdge
eiter firstEdge()
Definition: GModel.h:356
Volume::QuadTri
int QuadTri
Definition: Geo.h:146
Field.h
CreateVolume
Volume * CreateVolume(int Num, int Typ)
Definition: Geo.cpp:617
GEntity::UnknownModel
@ UnknownModel
Definition: GEntity.h:80
GEO_Internals::DelCurves
Tree_T * DelCurves
Definition: GModelIO_GEO.h:19
gmsh_sign
int gmsh_sign(T const &value)
Definition: Numeric.h:28
GEO_Internals::getMaxTag
int getMaxTag(int dim) const
Definition: GModelIO_GEO.cpp:99
SortCurvesConsecutive
static bool SortCurvesConsecutive(const std::vector< Curve * > &e, std::vector< std::vector< Vertex * > > &vs)
Definition: GModelIO_GEO.cpp:429
List_Action
void List_Action(List_T *liste, void(*action)(void *data, void *dummy))
Definition: ListUtils.cpp:275
PhysicalGroup::Num
int Num
Definition: Geo.h:156
SymmetryShapes
bool SymmetryShapes(double A, double B, double C, double D, List_T *shapes)
Definition: Geo.cpp:1608
GEO_Internals::boundaryLayer
bool boundaryLayer(const std::vector< std::pair< int, int > > &inDimTags, std::vector< std::pair< int, int > > &outDimTags, ExtrudeParams *e=0)
Definition: GModelIO_GEO.cpp:804
GEO_Internals::extrude
bool extrude(const std::vector< std::pair< int, int > > &inDimTags, double dx, double dy, double dz, std::vector< std::pair< int, int > > &outDimTags, ExtrudeParams *e=0)
Definition: GModelIO_GEO.cpp:774
CreateReversedCurve
Curve * CreateReversedCurve(Curve *c)
Definition: Geo.cpp:1147
DeleteVolume
void DeleteVolume(int iv, bool recursive)
Definition: Geo.cpp:1064
sortEntities
bool sortEntities(const std::pair< int, int > &a, const std::pair< int, int > &b)
Definition: GModelIO_GEO.cpp:1362
MSH_PHYSICAL_SURFACE
#define MSH_PHYSICAL_SURFACE
Definition: GeoDefines.h:45
gmshEdge.h
Curve
Definition: Geo.h:74
GEO_Internals::_maxPointNum
int _maxPointNum
Definition: GModelIO_GEO.h:24
FindSurfaceLoop
SurfaceLoop * FindSurfaceLoop(int inum)
Definition: Geo.cpp:744
gmshVertex.h
GEO_Internals::SurfaceLoops
Tree_T * SurfaceLoops
Definition: GModelIO_GEO.h:18
GEdge::resetMeshAttributes
virtual void resetMeshAttributes()
Definition: GEdge.cpp:187
GEO_Internals::rotate
bool rotate(const std::vector< std::pair< int, int > > &dimTags, double x, double y, double z, double ax, double ay, double az, double angle)
Definition: GModelIO_GEO.cpp:848
GFace
Definition: GFace.h:33
GModel::writeGEO
int writeGEO(const std::string &name, bool printLabels=true, bool onlyPhysicals=false)
Definition: GModelIO_GEO.cpp:1802
GEO_Internals::_maxPhysicalNum
int _maxPhysicalNum
Definition: GModelIO_GEO.h:25
GEO_Internals::setMaxTag
void setMaxTag(int dim, int val)
Definition: GModelIO_GEO.cpp:87
angle
double angle(const SVector3 &a, const SVector3 &b)
Definition: SVector3.h:157
VertexPtrLessThan
Definition: GModelIO_GEO.cpp:422
GVertex::resetMeshAttributes
virtual void resetMeshAttributes()
Definition: GVertex.cpp:30
GEO_Internals::setTransfiniteVolumeQuadTri
void setTransfiniteVolumeQuadTri(int tag)
Definition: GModelIO_GEO.cpp:1234
GEO_Internals::setTransfiniteLine
void setTransfiniteLine(int tag, int nPoints, int type, double coef)
Definition: GModelIO_GEO.cpp:1136
OS.h
MSH_SURF_PLAN
#define MSH_SURF_PLAN
Definition: GeoDefines.h:33
GEO_Internals::setMeshSizeFromBoundary
void setMeshSizeFromBoundary(int dim, int tag, int val)
Definition: GModelIO_GEO.cpp:1353
GEO_Internals::modifyPhysicalGroup
bool modifyPhysicalGroup(int dim, int tag, int op, const std::vector< int > &tags)
Definition: GModelIO_GEO.cpp:1000
MSH_PHYSICAL_LINE
#define MSH_PHYSICAL_LINE
Definition: GeoDefines.h:44
Msg::Debug
static void Debug(const char *fmt,...)
Definition: GmshMessage.cpp:752
GEO_Internals::setTransfiniteSurface
void setTransfiniteSurface(int tag, int arrangement, const std::vector< int > &cornerTags)
Definition: GModelIO_GEO.cpp:1163
c
static double c(int i, int j, fullMatrix< double > &CA, const std::vector< SPoint3 > &P, const std::vector< SPoint3 > &Q)
Definition: discreteFrechetDistance.cpp:15
MESH_TRANSFINITE
#define MESH_TRANSFINITE
Definition: GmshDefines.h:259
GEO_Internals::_freeAll
void _freeAll()
Definition: GModelIO_GEO.cpp:52
CompareEdgeLoop
int CompareEdgeLoop(const void *a, const void *b)
Definition: Geo.cpp:67
FreePhysicalGroup
void FreePhysicalGroup(void *a, void *b)
Definition: Geo.cpp:158
GEO_Internals::addCompoundBSpline
bool addCompoundBSpline(int &tag, const std::vector< int > &curveTags, int numPoints)
Definition: GModelIO_GEO.cpp:394
DuplicateCurve
Curve * DuplicateCurve(Curve *c)
Definition: Geo.cpp:830
GModel::remove
bool remove(GRegion *r)
Definition: GModel.cpp:435
GEO_Internals::setRecombine
void setRecombine(int dim, int tag, double angle)
Definition: GModelIO_GEO.cpp:1252
SplitCurve
bool SplitCurve(int line_id, List_T *vertices_id, List_T *curves)
Definition: Geo.cpp:3433
List_T
Definition: ListUtils.h:9
Msg::Warning
static void Warning(const char *fmt,...)
Definition: GmshMessage.cpp:543
List_Suppress
int List_Suppress(List_T *liste, void *data, int(*fcmp)(const void *a, const void *b))
Definition: ListUtils.cpp:236
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
ListUtils.h
fcmp_int
int fcmp_int(const void *a, const void *b)
Definition: ListUtils.cpp:26
GEntity::physicals
std::vector< int > physicals
Definition: GEntity.h:65
List_Reset
void List_Reset(List_T *liste)
Definition: ListUtils.cpp:269
SPoint3
Definition: SPoint3.h:14
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
FreeVolume
void FreeVolume(void *a, void *b)
Definition: Geo.cpp:635
List_Nbr
int List_Nbr(List_T *liste)
Definition: ListUtils.cpp:106
GModel::getFaceByTag
GFace * getFaceByTag(int n) const
Definition: GModel.cpp:326
Vertex
Definition: Geo.h:29
PhysicalGroup::Entities
List_T * Entities
Definition: Geo.h:158
MSH_POINT
#define MSH_POINT
Definition: GeoDefines.h:16
CompareVertex
int CompareVertex(const void *a, const void *b)
Definition: Geo.cpp:30
Tree2List
List_T * Tree2List(Tree_T *pTree)
Definition: TreeUtils.cpp:110
FreeVertex
void FreeVertex(void *a, void *b)
Definition: Geo.cpp:133
List_Create
List_T * List_Create(int n, int incr, int size)
Definition: ListUtils.cpp:46
GEO_Internals::DelPhysicalGroups
List_T * DelPhysicalGroups
Definition: GModelIO_GEO.h:20
GEO_Internals::addSurfaceFilling
bool addSurfaceFilling(int &tag, const std::vector< int > &wireTags, int sphereCenterTag=-1)
Definition: GModelIO_GEO.cpp:630
SortEdgesInLoop
bool SortEdgesInLoop(int num, List_T *edges, bool reorient)
Definition: Geo.cpp:3597
GVertex::edges
virtual std::vector< GEdge * > const & edges() const
Definition: GVertex.h:54
gmshRegion
Definition: gmshRegion.h:13
GModelIO_GEO.h
Surface::MeshSizeFromBoundary
int MeshSizeFromBoundary
Definition: Geo.h:130
gmshEdge
Definition: gmshEdge.h:13
GEntity::getNativeType
virtual ModelType getNativeType() const
Definition: GEntity.h:268
GEO_Internals::DelSurfaces
Tree_T * DelSurfaces
Definition: GModelIO_GEO.h:19
GModel::firstVertex
viter firstVertex()
Definition: GModel.h:357
CreateVertex
Vertex * CreateVertex(int Num, double X, double Y, double Z, double lc, double u)
Definition: Geo.cpp:106
GModel::deleteGEOInternals
void deleteGEOInternals()
Definition: GModelIO_GEO.cpp:1672
GModel::getEdgeByTag
GEdge * getEdgeByTag(int n) const
Definition: GModel.cpp:336
CreatePhysicalGroup
PhysicalGroup * CreatePhysicalGroup(int Num, int typ, List_T *intlist)
Definition: Geo.cpp:142
GmshMessage.h
edges
static int edges[6][2]
Definition: meshGRegionLocalMeshMod.cpp:23
GEntity::GmshModel
@ GmshModel
Definition: GEntity.h:81
MAX_LC
#define MAX_LC
Definition: GEntity.h:19
GModel::getNumMeshElements
std::size_t getNumMeshElements(int dim=-1) const
Definition: GModel.cpp:1540
GEO_Internals::addEllipseArc
bool addEllipseArc(int &tag, int startTag, int centerTag, int majorTag, int endTag, double nx=0., double ny=0., double nz=0.)
Definition: GModelIO_GEO.cpp:209
GEO_Internals::addLine
bool addLine(int &tag, int startTag, int endTag)
Definition: GModelIO_GEO.cpp:141
GModel::removePhysicalGroups
void removePhysicalGroups()
Definition: GModel.cpp:893
GEntity
Definition: GEntity.h:31
FindVolume
Volume * FindVolume(int inum)
Definition: Geo.cpp:722
SetBoundingBox
void SetBoundingBox(double xmin, double xmax, double ymin, double ymax, double zmin, double zmax)
Definition: OpenFile.cpp:115
GFace::regions
std::list< GRegion * > regions() const
Definition: GFace.h:92
MSH_SEGM_SPLN
#define MSH_SEGM_SPLN
Definition: GeoDefines.h:21
MSH_VOLUME_DISCRETE
#define MSH_VOLUME_DISCRETE
Definition: GeoDefines.h:41
ReplaceAllDuplicates
static void ReplaceAllDuplicates(std::vector< std::map< int, int > > &report)
Definition: Geo.cpp:2507
writePhysicalGroupGEO::writePhysicalGroupGEO
writePhysicalGroupGEO(FILE *fp, int i, bool labels, std::map< int, std::string > &o, std::map< std::pair< int, int >, std::string > &n)
Definition: GModelIO_GEO.cpp:1725
skipVertex
static bool skipVertex(GVertex *gv)
Definition: GModelIO_GEO.cpp:1792
TRANSLATE_ROTATE
#define TRANSLATE_ROTATE
Definition: ExtrudeParams.h:23
Vertex::boundaryLayerIndex
int boundaryLayerIndex
Definition: Geo.h:41
gmshFace
Definition: gmshFace.h:13
InterpolateCurve
Vertex InterpolateCurve(Curve *c, double u, int const derivee)
Definition: GeoInterpolation.cpp:441
GEO_Internals::addVolume
bool addVolume(int &tag, const std::vector< int > &shellTags)
Definition: GModelIO_GEO.cpp:701
GEntity::DiscreteCurve
@ DiscreteCurve
Definition: GEntity.h:104
GModel::_physicalNames
std::map< std::pair< int, int >, std::string > _physicalNames
Definition: GModel.h:150
Fopen
FILE * Fopen(const char *f, const char *mode)
Definition: OS.cpp:273
CreateSurface
Surface * CreateSurface(int Num, int Typ)
Definition: Geo.cpp:579
Volume::Surfaces
List_T * Surfaces
Definition: Geo.h:149
Vertex::Num
int Num
Definition: Geo.h:31
Shape::Type
int Type
Definition: GeoDefines.h:11
List_Add
void List_Add(List_T *liste, void *data)
Definition: ListUtils.cpp:90
GEdge::setTooSmall
void setTooSmall(bool const b)
Definition: GEdge.h:183
GModel::writePY
int writePY(const std::string &name, bool printLabels=true, bool onlyPhysicals=false)
Definition: GModelIO_GEO.cpp:1863
CTX::instance
static CTX * instance()
Definition: Context.cpp:122
GEO_Internals::Points
Tree_T * Points
Definition: GModelIO_GEO.h:18
GModel::getPhysicalGroups
void getPhysicalGroups(std::map< int, std::vector< GEntity * > > groups[4]) const
Definition: GModel.cpp:837
GModel::lastFace
fiter lastFace()
Definition: GModel.h:359
GEO_Internals::_maxLineNum
int _maxLineNum
Definition: GModelIO_GEO.h:24
ExtrudeParams::normalsCoherence
static std::vector< SPoint3 > normalsCoherence
Definition: ExtrudeParams.h:65
GEO_Internals::addSpline
bool addSpline(int &tag, const std::vector< int > &pointTags)
Definition: GModelIO_GEO.cpp:245
DeleteSurface
void DeleteSurface(int is, bool recursive)
Definition: Geo.cpp:1018
Volume::Method
int Method
Definition: Geo.h:144
gmshFace.h
gmshSurface
Definition: gmshSurface.h:22
DeletePhysicalSurface
void DeletePhysicalSurface(int num)
Definition: Geo.cpp:1125
writePhysicalGroupGEO::oldLabels
std::map< int, std::string > & oldLabels
Definition: GModelIO_GEO.cpp:1721
FindCurve
static Curve * FindCurve(int inum, Tree_T *t)
Definition: Geo.cpp:694
skipRegion
static bool skipRegion(GRegion *gr)
Definition: GModelIO_GEO.cpp:1766
GEO_Internals::Surfaces
Tree_T * Surfaces
Definition: GModelIO_GEO.h:18
GEO_Internals::addDiscreteSurface
bool addDiscreteSurface(int &tag)
Definition: GModelIO_GEO.cpp:617
CompareSurfaceLoop
int CompareSurfaceLoop(const void *a, const void *b)
Definition: Geo.cpp:60
GModel::lastVertex
viter lastVertex()
Definition: GModel.h:361
gmsh_yysymbols
std::map< std::string, gmsh_yysymbol > gmsh_yysymbols
Definition: Gmsh.tab.cpp:643
FindPhysicalGroup
PhysicalGroup * FindPhysicalGroup(int num, int type)
Definition: Geo.cpp:755
GVertex
Definition: GVertex.h:23
Tree_Action
void Tree_Action(Tree_T *tree, void(*action)(void *data, void *dummy))
Definition: TreeUtils.cpp:100
CompareSurface
int CompareSurface(const void *a, const void *b)
Definition: Geo.cpp:81
TranslateShapes
bool TranslateShapes(double X, double Y, double Z, List_T *shapes)
Definition: Geo.cpp:1549
Tree_Add
void * Tree_Add(Tree_T *tree, void *data)
Definition: TreeUtils.cpp:37
GEntity::DiscreteVolume
@ DiscreteVolume
Definition: GEntity.h:120
Numeric.h
GEO_Internals::DelPoints
Tree_T * DelPoints
Definition: GModelIO_GEO.h:19
GModel::getNumVertices
std::size_t getNumVertices() const
Definition: GModel.h:347
DeleteCurve
void DeleteCurve(int ip, bool recursive)
Definition: Geo.cpp:981
FindEdgeLoop
EdgeLoop * FindEdgeLoop(int inum)
Definition: Geo.cpp:733
MSH_PHYSICAL_VOLUME
#define MSH_PHYSICAL_VOLUME
Definition: GeoDefines.h:46
GModel
Definition: GModel.h:44
GEO_Internals::DelVolumes
Tree_T * DelVolumes
Definition: GModelIO_GEO.h:19
GEO_Internals::_transform
bool _transform(int mode, const std::vector< std::pair< int, int > > &dimTags, double x, double y, double z, double dx, double dy, double dz, double a, double b, double c, double d)
Definition: GModelIO_GEO.cpp:812
GEO_Internals::symmetry
bool symmetry(const std::vector< std::pair< int, int > > &dimTags, double a, double b, double c, double d)
Definition: GModelIO_GEO.cpp:862
gmshSphere::NewSphere
static gmshSurface * NewSphere(int _iSphere, double _x, double _y, double _z, double _r)
Definition: gmshSurface.cpp:38
GEO_Internals::_changed
bool _changed
Definition: GModelIO_GEO.h:28
MSH_SEGM_NURBS
#define MSH_SEGM_NURBS
Definition: GeoDefines.h:28
Vertex::Typ
int Typ
Definition: Geo.h:32
writePhysicalGroupGEO::printLabels
bool printLabels
Definition: GModelIO_GEO.cpp:1720
GEO_Internals::splitCurve
bool splitCurve(int tag, const std::vector< int > &pointTags, std::vector< int > &curveTags)
Definition: GModelIO_GEO.cpp:868
contextGeometryOptions::orientedPhysicals
int orientedPhysicals
Definition: Context.h:111
GEO_Internals::setCompoundMesh
void setCompoundMesh(int dim, const std::vector< int > &tags)
Definition: GModelIO_GEO.cpp:1111
GEO_Internals::setDegenerated
void setDegenerated(int dim, int tag)
Definition: GModelIO_GEO.cpp:1128
GEO_Internals::setMeshAlgorithm
void setMeshAlgorithm(int dim, int tag, int val)
Definition: GModelIO_GEO.cpp:1344
Surface
Definition: Geo.h:111
GModel::firstFace
fiter firstFace()
Definition: GModel.h:355
MSH_SURF_DISCRETE
#define MSH_SURF_DISCRETE
Definition: GeoDefines.h:38
writePhysicalGroupGEO::newLabels
std::map< std::pair< int, int >, std::string > & newLabels
Definition: GModelIO_GEO.cpp:1722
GEO_Internals::getVertex
bool getVertex(int tag, double &x, double &y, double &z)
Definition: GModelIO_GEO.cpp:1616
List_Delete
void List_Delete(List_T *liste)
Definition: ListUtils.cpp:66
FreeCurve
void FreeCurve(void *a, void *b)
Definition: Geo.cpp:567
GModel::createGEOInternals
void createGEOInternals()
Definition: GModelIO_GEO.cpp:1670
TRANSFINITE_QUADTRI_1
#define TRANSFINITE_QUADTRI_1
Definition: GmshDefines.h:268
Shape
Definition: GeoDefines.h:9
GModel::add
bool add(GRegion *r)
Definition: GModel.h:394
Parser.h
GModel::getNumRegions
std::size_t getNumRegions() const
Definition: GModel.h:344
Vertex::lc
double lc
Definition: Geo.h:33
MSH_SURF_REGL
#define MSH_SURF_REGL
Definition: GeoDefines.h:34
writePhysicalGroupGEO::geo
FILE * geo
Definition: GModelIO_GEO.cpp:1718
GEO_Internals::addVertex
bool addVertex(int &tag, double x, double y, double z, double lc)
Definition: GModelIO_GEO.cpp:112
GModel::getFields
FieldManager * getFields()
Definition: GModel.h:325
PhysicalGroup
Definition: Geo.h:154
GEO_Internals::_maxSurfaceLoopNum
int _maxSurfaceLoopNum
Definition: GModelIO_GEO.h:25
ROTATE
#define ROTATE
Definition: ExtrudeParams.h:22
GEntity::tag
int tag() const
Definition: GEntity.h:280
GModel::getNumFaces
std::size_t getNumFaces() const
Definition: GModel.h:345
EndCurve
bool EndCurve(Curve *c)
Definition: Geo.cpp:225
EdgeLoop
Definition: Geo.h:105
skipEdge
static bool skipEdge(GEdge *ge)
Definition: GModelIO_GEO.cpp:1782
Field
Definition: Field.h:103
GEO_Internals::dilate
bool dilate(const std::vector< std::pair< int, int > > &dimTags, double x, double y, double z, double a, double b, double c)
Definition: GModelIO_GEO.cpp:855
GEO_Internals::_allocateAll
void _allocateAll()
Definition: GModelIO_GEO.cpp:29
RotateShapes
bool RotateShapes(double Ax, double Ay, double Az, double Px, double Py, double Pz, double alpha, List_T *shapes)
Definition: Geo.cpp:1580
GEO_Internals::newGeometrySphere
gmshSurface * newGeometrySphere(int tag, int centerTag, int pointTag)
Definition: GModelIO_GEO.cpp:1628
GModel::exportDiscreteGEOInternals
int exportDiscreteGEOInternals()
Definition: GModelIO_GEO.cpp:1910
CTX::geom
contextGeometryOptions geom
Definition: Context.h:311
DeletePhysicalPoint
void DeletePhysicalPoint(int num)
Definition: Geo.cpp:1103
GRegion
Definition: GRegion.h:28
Tree_Delete
void Tree_Delete(Tree_T *tree)
Definition: TreeUtils.cpp:23
PhysicalGroup::Typ
int Typ
Definition: Geo.h:157
VertexPtrLessThan::operator()
bool operator()(const Vertex *v1, const Vertex *v2) const
Definition: GModelIO_GEO.cpp:423
GEO_Internals::addCurveLoops
bool addCurveLoops(const std::vector< int > &curveTags, std::vector< int > &curveLoopTags)
Definition: GModelIO_GEO.cpp:529
GModel::getNumEdges
std::size_t getNumEdges() const
Definition: GModel.h:346
GEdge::faces
virtual std::vector< GFace * > faces() const
Definition: GEdge.h:113
GEntity::DiscreteSurface
@ DiscreteSurface
Definition: GEntity.h:117
MSH_SEGM_CIRC
#define MSH_SEGM_CIRC
Definition: GeoDefines.h:22
MSH_SURF_TRIC
#define MSH_SURF_TRIC
Definition: GeoDefines.h:35
IntersectCurvesWithSurface
bool IntersectCurvesWithSurface(List_T *curve_ids, int surface_id, List_T *shapes)
Definition: Geo.cpp:3560
writePhysicalGroupGEO
Definition: GModelIO_GEO.cpp:1716
MSH_PHYSICAL_POINT
#define MSH_PHYSICAL_POINT
Definition: GeoDefines.h:43
gmshRegion.h
EndSurface
void EndSurface(Surface *s)
Definition: Geo.cpp:426
Surface::MeshAlgorithm
int MeshAlgorithm
Definition: Geo.h:130
GModel::firstRegion
riter firstRegion()
Definition: GModel.h:354
Shape::Num
int Num
Definition: GeoDefines.h:12
GEO_Internals::_addCompoundSpline
bool _addCompoundSpline(int &tag, const std::vector< int > &curveTags, int numPoints, bool bspline)
Definition: GModelIO_GEO.cpp:336
Surface::TransfiniteSmoothing
int TransfiniteSmoothing
Definition: Geo.h:119
CreateEdgeLoop
EdgeLoop * CreateEdgeLoop(int Num, List_T *intlist)
Definition: Geo.cpp:168
GModel::lastRegion
riter lastRegion()
Definition: GModel.h:358
GEO_Internals::translate
bool translate(const std::vector< std::pair< int, int > > &dimTags, double dx, double dy, double dz)
Definition: GModelIO_GEO.cpp:842
GEO_Internals::addCircleArc
bool addCircleArc(int &tag, int startTag, int centerTag, int endTag, double nx=0., double ny=0., double nz=0.)
Definition: GModelIO_GEO.cpp:175
CreateSurfaceLoop
SurfaceLoop * CreateSurfaceLoop(int Num, List_T *intlist)
Definition: Geo.cpp:193
Coord::X
double X
Definition: Geo.h:26
EdgeLoop::Curves
List_T * Curves
Definition: Geo.h:108
Surface::RecombineAngle
double RecombineAngle
Definition: Geo.h:118
Surface::Recombine
int Recombine
Definition: Geo.h:116
Context.h
GEO_Internals::intersectCurvesWithSurface
bool intersectCurvesWithSurface(const std::vector< int > &curveTags, int surfaceTag, std::vector< int > &pointTags)
Definition: GModelIO_GEO.cpp:889
GEO_Internals::mergeVertices
bool mergeVertices(const std::vector< int > &tags)
Definition: GModelIO_GEO.cpp:1083
DilatShapes
bool DilatShapes(double X, double Y, double Z, double A, double B, double C, List_T *shapes)
Definition: Geo.cpp:1564
DeletePhysicalVolume
void DeletePhysicalVolume(int num)
Definition: Geo.cpp:1136
GEO_Internals::setReverseMesh
void setReverseMesh(int dim, int tag, bool val=1)
Definition: GModelIO_GEO.cpp:1309
GEO_Internals::revolve
bool revolve(const std::vector< std::pair< int, int > > &inDimTags, double x, double y, double z, double ax, double ay, double az, double angle, std::vector< std::pair< int, int > > &outDimTags, ExtrudeParams *e=0)
Definition: GModelIO_GEO.cpp:783
z
const double z
Definition: GaussQuadratureQuad.cpp:56
Field::id
int id
Definition: Field.h:111
Curve::Num
int Num
Definition: Geo.h:76
GEO_Internals::copy
bool copy(const std::vector< std::pair< int, int > > &inDimTags, std::vector< std::pair< int, int > > &outDimTags)
Definition: GModelIO_GEO.cpp:913
GEO_Internals::twist
bool twist(const std::vector< std::pair< int, int > > &inDimTags, double x, double y, double z, double dx, double dy, double dz, double ax, double ay, double az, double angle, std::vector< std::pair< int, int > > &outDimTags, ExtrudeParams *e=0)
Definition: GModelIO_GEO.cpp:793
GModel::faces
std::set< GFace *, GEntityPtrLessThan > faces
Definition: GModel.h:144
Tree_Nbr
int Tree_Nbr(Tree_T *tree)
Definition: TreeUtils.cpp:46
FreeEdgeLoop
void FreeEdgeLoop(void *a, void *b)
Definition: Geo.cpp:183
GRegion::resetMeshAttributes
virtual void resetMeshAttributes()
Definition: GRegion.cpp:170
GEO_Internals::addCompoundSpline
bool addCompoundSpline(int &tag, const std::vector< int > &curveTags, int numPoints)
Definition: GModelIO_GEO.cpp:387
GEO_Internals::addPlaneSurface
bool addPlaneSurface(int &tag, const std::vector< int > &wireTags)
Definition: GModelIO_GEO.cpp:592
GEO_Internals::remove
bool remove(int dim, int tag, bool recursive=false)
Definition: GModelIO_GEO.cpp:969
Coord::Z
double Z
Definition: Geo.h:26
GEdge
Definition: GEdge.h:26
MSH_SEGM_BSPLN
#define MSH_SEGM_BSPLN
Definition: GeoDefines.h:27
GEO_Internals::setSmoothing
void setSmoothing(int tag, int val)
Definition: GModelIO_GEO.cpp:1291
BOUNDARY_LAYER
#define BOUNDARY_LAYER
Definition: ExtrudeParams.h:24
DeletePhysicalLine
void DeletePhysicalLine(int num)
Definition: Geo.cpp:1114
writePhysicalGroupGEO::dim
int dim
Definition: GModelIO_GEO.cpp:1719
GEO_Internals::addSurfaceLoop
bool addSurfaceLoop(int &tag, const std::vector< int > &surfaceTags)
Definition: GModelIO_GEO.cpp:680
FreeSurface
void FreeSurface(void *a, void *b)
Definition: Geo.cpp:604
DuplicateVertex
Vertex * DuplicateVertex(Vertex *v)
Definition: Geo.cpp:778
Coord::Y
double Y
Definition: Geo.h:26
GEO_Internals::setMeshSize
void setMeshSize(int dim, int tag, double size)
Definition: GModelIO_GEO.cpp:1117
MSH_VOLUME
#define MSH_VOLUME
Definition: GeoDefines.h:40
GEO_Internals::_extrude
bool _extrude(int mode, const std::vector< std::pair< int, int > > &inDimTags, double x, double y, double z, double dx, double dy, double dz, double ax, double ay, double az, double angle, std::vector< std::pair< int, int > > &outDimTags, ExtrudeParams *e=0)
Definition: GModelIO_GEO.cpp:722
Vertex::Pos
Coord Pos
Definition: Geo.h:34
ExtrudeShapes
void ExtrudeShapes(int type, List_T *list_in, double T0, double T1, double T2, double A0, double A1, double A2, double X0, double X1, double X2, double alpha, ExtrudeParams *e, List_T *list_out)
Definition: Geo.cpp:3221
SurfaceLoop
Definition: Geo.h:134
Volume::TrsfPoints
List_T * TrsfPoints
Definition: Geo.h:148
CreateCurve
Curve * CreateCurve(int Num, int Typ, int Order, List_T *Liste, List_T *Knots, int p1, int p2, double u1, double u2, bool &ok)
Definition: Geo.cpp:445
GModel::lastEdge
eiter lastEdge()
Definition: GModel.h:360
Surface::Generatrices
List_T * Generatrices
Definition: Geo.h:120
GModel.h
DuplicateSurface
Surface * DuplicateSurface(Surface *s)
Definition: Geo.cpp:895
Volume::Recombine3D
int Recombine3D
Definition: Geo.h:145
SetVolumeSurfaces
bool SetVolumeSurfaces(Volume *v, List_T *loops)
Definition: Geo.cpp:3734
DeletePoint
void DeletePoint(int ip, bool recursive)
Definition: Geo.cpp:956
GModel::getRegionByTag
GRegion * getRegionByTag(int n) const
Definition: GModel.cpp:316
FreeSurfaceLoop
void FreeSurfaceLoop(void *a, void *b)
Definition: Geo.cpp:208
MSH_SEGM_ELLI
#define MSH_SEGM_ELLI
Definition: GeoDefines.h:24
Surface::Recombine_Dir
int Recombine_Dir
Definition: Geo.h:117
FindPoint
static Vertex * FindPoint(int inum, Tree_T *t)
Definition: Geo.cpp:680
GEO_Internals::resetPhysicalGroups
void resetPhysicalGroups()
Definition: GModelIO_GEO.cpp:992
Surface::ReverseMesh
int ReverseMesh
Definition: Geo.h:130
Tree_Create
Tree_T * Tree_Create(int size, int(*fcmp)(const void *a, const void *b))
Definition: TreeUtils.cpp:15
TRANSLATE
#define TRANSLATE
Definition: ExtrudeParams.h:21
CircParam::n
double n[3]
Definition: Geo.h:71
GEO_Internals::synchronize
void synchronize(GModel *model, bool resetMeshAttributes=true)
Definition: GModelIO_GEO.cpp:1370
CompareVolume
int CompareVolume(const void *a, const void *b)
Definition: Geo.cpp:88
Surface::InSphereCenter
Vertex * InSphereCenter
Definition: Geo.h:123
GEO_Internals::EdgeLoops
Tree_T * EdgeLoops
Definition: GModelIO_GEO.h:18
skipFace
static bool skipFace(GFace *gf)
Definition: GModelIO_GEO.cpp:1772
GModel::getVertexByTag
GVertex * getVertexByTag(int n) const
Definition: GModel.cpp:346
Volume::Num
int Num
Definition: Geo.h:142
GModel::_geo_internals
GEO_Internals * _geo_internals
Definition: GModel.h:119
GEO_Internals
Definition: GModelIO_GEO.h:15
ExtrudeParams
Definition: ExtrudeParams.h:26
GEO_Internals::_maxSurfaceNum
int _maxSurfaceNum
Definition: GModelIO_GEO.h:24
GModel::edges
std::set< GEdge *, GEntityPtrLessThan > edges
Definition: GModel.h:145
GEO_Internals::setTransfiniteVolume
void setTransfiniteVolume(int tag, const std::vector< int > &cornerTags)
Definition: GModelIO_GEO.cpp:1201
writePhysicalGroupGEO::operator()
void operator()(std::pair< const int, std::vector< GEntity * > > &g)
Definition: GModelIO_GEO.cpp:1732
List_Read
void List_Read(List_T *liste, int index, void *data)
Definition: ListUtils.cpp:111
OpenFile.h
FindSurface
static Surface * FindSurface(int inum, Tree_T *t)
Definition: Geo.cpp:708
CompareCurve
int CompareCurve(const void *a, const void *b)
Definition: Geo.cpp:74
MSH_SEGM_LINE
#define MSH_SEGM_LINE
Definition: GeoDefines.h:20
GEO_Internals::_maxVolumeNum
int _maxVolumeNum
Definition: GModelIO_GEO.h:25
GEO_Internals::newGeometryPolarSphere
gmshSurface * newGeometryPolarSphere(int tag, int centerTag, int pointTag)
Definition: GModelIO_GEO.cpp:1648
gmshVertex
Definition: gmshVertex.h:13
GEO_Internals::addCurveLoop
bool addCurveLoop(int &tag, const std::vector< int > &curveTags, bool reorient=false)
Definition: GModelIO_GEO.cpp:401
GEO_Internals::_meshCompounds
std::multimap< int, std::vector< int > > _meshCompounds
Definition: GModelIO_GEO.h:23
GEO_Internals::addBezier
bool addBezier(int &tag, const std::vector< int > &pointTags)
Definition: GModelIO_GEO.cpp:271
Surface::Method
int Method
Definition: Geo.h:115
DuplicateVolume
Volume * DuplicateVolume(Volume *v)
Definition: Geo.cpp:942
faces
static int faces[4][3]
Definition: meshGRegionDelaunayInsertion.cpp:165
GEO_Internals::_maxLineLoopNum
int _maxLineLoopNum
Definition: GModelIO_GEO.h:24
Curve::Circle
CircParam Circle
Definition: Geo.h:92
gmshPolarSphere::NewPolarSphere
static gmshSurface * NewPolarSphere(int _iSphere, double _x, double _y, double _z, double _r)
Definition: gmshSurface.cpp:70
SetSurfaceGeneratrices
bool SetSurfaceGeneratrices(Surface *s, List_T *loops)
Definition: Geo.cpp:3674