File: write_MEDIT.h

package info (click to toggle)
cgal 6.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 144,912 kB
  • sloc: cpp: 810,858; ansic: 208,477; sh: 493; python: 411; makefile: 286; javascript: 174
file content (109 lines) | stat: -rw-r--r-- 4,259 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
// Copyright (c) 2019-2024  GeometryFactory Sarl (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL: https://github.com/CGAL/cgal/blob/v6.1/Constrained_triangulation_3/include/CGAL/IO/write_MEDIT.h $
// $Id: include/CGAL/IO/write_MEDIT.h b26b07a1242 $
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
// Author(s)     : Jane Tournois

#include "CGAL/unordered_flat_map.h"
#include <CGAL/Conforming_constrained_Delaunay_triangulation_3.h>

#include <CGAL/Named_function_parameters.h>
#include <CGAL/boost/graph/named_params_helper.h>

#include <CGAL/IO/File_medit.h>

#include <stack>
#include <ostream>

namespace CGAL
{
namespace IO
{
/*!
 * @ingroup PkgCDT3IOFunctions
 * @brief outputs a conforming constrained Delaunay triangulation to
 * the MEDIT (`.mesh`) file format.
 *        See \cgalCite{frey:inria-00069921} for a comprehensive description of this
 *        file format.
 * @param os the output stream
 * @param ccdt the conforming constrained Delaunay triangulation to be written
 * \param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
 *
 * \cgalNamedParamsBegin
 *    \cgalParamNBegin{with_plc_face_id}
 *    \cgalParamDescription{a Boolean activating the numbering of PLC face identifiers in the output.
 *                          If `ccdt` was constructed with the `plc_face_id` property map given as a named parameter,
 *                          and this parameter is set to `true`,
 *                          the output will contain the corresponding patch identifier for each facet of the triangulation.
 *                          If this parameter is set to `false`, the output will not contain any patch identifier.
 *                          If `ccdt` was not constructed with the `plc_face_id` property map, and this parameter is
 *                          set to `true`, the output will contain a patch identifier for each facet of the triangulation.}
 *    \cgalParamType{Boolean}
 *    \cgalParamDefault{`false`}
 *   \cgalParamNEnd
 * \cgalNamedParamsEnd

 * \see \ref IOStreamMedit
 */
template <typename Traits,
          typename Tr,
          typename NamedParameters = parameters::Default_named_parameters>
void write_MEDIT(std::ostream& os,
                 const Conforming_constrained_Delaunay_triangulation_3<Traits, Tr>& ccdt,
                 const NamedParameters& np = parameters::default_values())
{
  const auto& tr = ccdt.triangulation();

  using Tr_ = typename cpp20::remove_cvref_t<decltype(tr)>;

  using Vertex_handle = typename Tr_::Vertex_handle;
  using Facet = typename Tr_::Facet;
  using Cell_handle = typename Tr_::Cell_handle;
  CGAL::unordered_flat_map<Cell_handle, bool> cells_in_domain;
  for(Cell_handle c : tr.all_cell_handles())
    cells_in_domain[c] = true;

  std::stack<Cell_handle> stack;
  stack.push(tr.infinite_cell());
  while(!stack.empty())
  {
    auto ch = stack.top();
    stack.pop();
    cells_in_domain[ch] = false;
    for(int i = 0; i < 4; ++i)
    {
      if(ccdt.is_facet_constrained(ch, i))
        continue;
      auto n = ch->neighbor(i);
      if(cells_in_domain[n])
        stack.push(n);
    }
  }

  using parameters::choose_parameter;
  using parameters::get_parameter;
  const bool has_plc_face_id
    = choose_parameter(get_parameter(np, internal_np::with_plc_face_id), false);

  auto plc_patch_map = boost::make_function_property_map<Facet>([&](const Facet& f)
    { return has_plc_face_id ? f.first->ccdt_3_data().face_constraint_index(f.second) + 1 : 1; });

  return SMDS_3::output_to_medit(os,
                                 tr,
                                 tr.finite_vertex_handles(),
                                 ccdt.constrained_facets(),
                                 tr.finite_cell_handles(),
                                 boost::make_function_property_map<Vertex_handle>([](Vertex_handle) { return 0; }),
                                 plc_patch_map,
                                 boost::make_function_property_map<Cell_handle>([&](Cell_handle ch) {
                                   return cells_in_domain[ch] ? 1 : 0;
                                 }));
}

}// end namespace IO
}// end namespace CGAL