NEMoSys  0.63.0
A modular, extensible resource with robust automated mesh generation, mesh quality analysis, adaptive mesh refinement, and data transfer between arbitrary meshes.
dataSetRegionBoundaryFilter.H
Go to the documentation of this file.
1 /*******************************************************************************
2 * Promesh *
3 * Copyright (C) 2022, IllinoisRocstar LLC. All rights reserved. *
4 * *
5 * Promesh is the property of IllinoisRocstar LLC. *
6 * *
7 * IllinoisRocstar LLC *
8 * Champaign, IL *
9 * www.illinoisrocstar.com *
10 * promesh@illinoisrocstar.com *
11 *******************************************************************************/
12 /*******************************************************************************
13 * This file is part of Promesh *
14 * *
15 * This version of Promesh is free software: you can redistribute it and/or *
16 * modify it under the terms of the GNU Lesser General Public License as *
17 * published by the Free Software Foundation, either version 3 of the License, *
18 * or (at your option) any later version. *
19 * *
20 * Promesh is distributed in the hope that it will be useful, but WITHOUT ANY *
21 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS *
22 * FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more *
23 * details. *
24 * *
25 * You should have received a copy of the GNU Lesser General Public License *
26 * along with this program. If not, see <https://www.gnu.org/licenses/>. *
27 * *
28 *******************************************************************************/
29 /*
30  * Inspired by vtkDataSetRegionSurfaceFilter, but with an algorithm using STL
31  * instead of vtkFastGeomQuadStruct.
32  */
33 
34 #ifndef NEMOSYS_DATASETREGIONBOUNDARYFILTER_H_
35 #define NEMOSYS_DATASETREGIONBOUNDARYFILTER_H_
36 
37 #include "nemosys_export.h"
38 #include <vtkPolyDataAlgorithm.h>
39 
40 #include <string>
41 
42 namespace NEM {
43 namespace MSH {
44 
45 /**
46  * @brief Like vtkDataSetRegionSurfaceFilter, this filter extracts boundaries of
47  * materials, including both interfaces and external boundaries. Unlike
48  * vtkDataSetRegionSurfaceFilter, can extract boundary edges for 2d data sets.
49  * @details Output edges/faces are represented once for interfaces, with the
50  * orientation chosen by the input cell with lower material number. If input is
51  * a vtkPointSet, output points are shallow copy; otherwise, input points are
52  * deep copy - in either case, all points are copied, even those unused in
53  * output. Output point data is shallow copy of input point data; input cell
54  * data is not copied. Unlike vtkDataSetRegionSurfaceFilter, point ordering
55  * should match the edge/face described in OrigCellArrayName and
56  * CellSideArrayName arrays.
57  */
58 class NEMOSYS_EXPORT dataSetRegionBoundaryFilter : public vtkPolyDataAlgorithm {
59  public:
61  vtkTypeMacro(dataSetRegionBoundaryFilter, vtkPolyDataAlgorithm)
62 
63  //@{
64  /**
65  * @brief If 2-D, finds edges; if 3-D, finds faces. Default is 3.
66  */
67  vtkSetClampMacro(Dimension, int, 2, 3);
68  vtkGetMacro(Dimension, int);
69  //@}
70 
71  //@{
72  /**
73  * @brief Input material/region array name. Default is "Material".
74  * @details Assumed to be vtkIntArray in input cell data. Interfaces between
75  * materials are included in output. Should not use the value -1.
76  */
77  vtkSetMacro(MaterialArrayName, const std::string &);
78  vtkGetMacro(MaterialArrayName, const std::string &);
79  //@}
80 
81  //@{
82  /**
83  * @brief Output region array name. Default is "Region".
84  * @details If not empty, array will be added as vtkIntArray in output cell
85  * data. To recover relation to material, see RegionToMaterialArrayName
86  */
87  vtkSetMacro(BoundaryRegionArrayName, const std::string &);
88  vtkGetMacro(BoundaryRegionArrayName, const std::string &);
89  //@}
90 
91  //@{
92  /**
93  * @brief Output region to material array name. Default is "RegionToMaterial".
94  * @details If not empty, array will be added as vtkIntArray with two
95  * components in output field data. Indexed by region (value in
96  * BoundaryRegionArrayName array); value corresponds to material(s), with -1
97  * corresponding to no material (for regions bounding only one material).
98  */
99  vtkSetMacro(RegionToMaterialArrayName, const std::string &);
100  vtkGetMacro(RegionToMaterialArrayName, const std::string &);
101  //@}
102 
103  //@{
104  /**
105  * @brief Output original cell id array name. Default is "OrigCellIds".
106  * @details If not empty, array will be added as vtkIdTypeArray with two
107  * components in output cell data. For each output edge/face, contains the
108  * index of the input cell. Second component is -1 if external edge/face.
109  */
110  vtkSetMacro(OrigCellArrayName, const std::string &);
111  vtkGetMacro(OrigCellArrayName, const std::string &);
112  //@}
113 
114  //@{
115  /**
116  * @brief Output original cell edge/face id array name. Default is
117  * "CellFaceIds".
118  * @details If not empty, array will be added as vtkIntArray with two
119  * components in output cell data. For each output edge/face, contains the
120  * index of the corresponding edge/face of the cell given by the
121  * OrigCellArrayName array. Second component is -1 if external edge/face.
122  */
123  vtkSetMacro(CellSideArrayName, const std::string &);
124  vtkGetMacro(CellSideArrayName, const std::string &);
125  //@}
126 
127  protected:
128  dataSetRegionBoundaryFilter() = default;
129 
130  int RequestData(vtkInformation *request, vtkInformationVector **inputVector,
131  vtkInformationVector *outputVector) override;
132  int FillInputPortInformation(int port, vtkInformation *info) override;
133 
134  int Dimension{3};
135  std::string MaterialArrayName{"Material"};
136  std::string RegionToMaterialArrayName{"RegionToMaterial"};
137  std::string OrigCellArrayName{"OrigCellIds"};
138  std::string CellSideArrayName{"CellFaceIds"};
139  std::string BoundaryRegionArrayName{"Region"};
140 };
141 
142 } // namespace MSH
143 } // namespace NEM
144 
145 #endif // NEMOSYS_DATASETREGIONBOUNDARYFILTER_H_
geoMeshBase * New(MeshType meshType)
Create a new mesh object.
STL namespace.
Like vtkDataSetRegionSurfaceFilter, this filter extracts boundaries of materials, including both inte...