possibility to use old and new surface projection strategy
[engrid.git] / src / libengrid / surfacealgorithm.cpp
blobf88d2e78bdaf791623b5aa98c89158eedee9ed08
1 //
2 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 // + +
4 // + This file is part of enGrid. +
5 // + +
6 // + Copyright 2008-2012 enGits GmbH +
7 // + +
8 // + enGrid is free software: you can redistribute it and/or modify +
9 // + it under the terms of the GNU General Public License as published by +
10 // + the Free Software Foundation, either version 3 of the License, or +
11 // + (at your option) any later version. +
12 // + +
13 // + enGrid is distributed in the hope that it will be useful, +
14 // + but WITHOUT ANY WARRANTY; without even the implied warranty of +
15 // + MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +
16 // + GNU General Public License for more details. +
17 // + +
18 // + You should have received a copy of the GNU General Public License +
19 // + along with enGrid. If not, see <http://www.gnu.org/licenses/>. +
20 // + +
21 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
22 //
23 #include "surfacealgorithm.h"
25 #include "insertpoints.h"
26 #include "removepoints.h"
27 #include "updatedesiredmeshdensity.h"
28 #include "smoothingutilities.h"
29 #include "swaptriangles.h"
30 #include "laplacesmoother.h"
31 #include "guimainwindow.h"
34 SurfaceAlgorithm::SurfaceAlgorithm()
36 EG_TYPENAME;
37 getSet("surface meshing", "maximal number of iterations", 5, m_NumMaxIter);
38 getSet("surface meshing", "number of smoothing steps" , 2, m_NumSmoothSteps);
39 getSet("surface meshing", "number of Delaunay sweeps" , 1, m_NumDelaunaySweeps);
40 m_NodesPerQuarterCircle = 0;
41 m_RespectFeatureEdgesForDeleteNodes = false;
42 m_FeatureAngleForDeleteNodes = deg2rad(45);
43 m_PerformGeometricTests = true;
44 m_UseProjectionForSmoothing = true;
45 m_UseNormalCorrectionForSmoothing = false;
46 m_AllowFeatureEdgeSwapping = true;
47 m_AllowSmallAreaSwapping = false;
48 m_GrowthFactor = 1.5;
49 m_FeatureResolution2D = 0;
50 m_FeatureResolution3D = 0;
51 setDeleteNodesOn();
52 setInsertNodesOn();
55 void SurfaceAlgorithm::readSettings()
57 QString buffer = GuiMainWindow::pointer()->getXmlSection("engrid/surface/settings").replace("\n", " ");
58 QTextStream in(&buffer, QIODevice::ReadOnly);
59 in >> m_MaxEdgeLength;
60 in >> m_MinEdgeLength;
61 in >> m_GrowthFactor;
62 in >> m_NodesPerQuarterCircle;
63 int num_bcs;
64 in >> num_bcs;
65 QVector<int> tmp_bcs;
66 GuiMainWindow::pointer()->getAllBoundaryCodes(tmp_bcs);
67 m_BoundaryCodes.clear();
68 if (num_bcs == tmp_bcs.size()) {
69 foreach (int bc, tmp_bcs) {
70 int check_state;
71 in >> check_state;
72 if (check_state == 1) {
73 m_BoundaryCodes.insert(bc);
77 if (!in.atEnd()) {
78 in >> m_FeatureResolution2D;
79 in >> m_FeatureResolution3D;
83 void SurfaceAlgorithm::prepare()
85 setAllCells();
86 readSettings();
88 EG_VTKDCN(vtkCharArray, node_type, m_Grid, "node_type");//node type
90 updateNodeInfo();
94 void SurfaceAlgorithm::computeMeshDensity()
96 ///\todo Optimize by using only one loop through nodes!
97 UpdateDesiredMeshDensity update_desired_mesh_density;
98 update_desired_mesh_density.setGrid(m_Grid);
99 update_desired_mesh_density.setVertexMeshDensityVector(m_VMDvector);
100 update_desired_mesh_density.setMaxEdgeLength(m_MaxEdgeLength);
101 update_desired_mesh_density.setMinEdgeLength(m_MinEdgeLength);
102 update_desired_mesh_density.setNodesPerQuarterCircle(m_NodesPerQuarterCircle);
103 update_desired_mesh_density.setCellGrowthFactor(m_GrowthFactor);
104 update_desired_mesh_density.setBoundaryCodes(m_BoundaryCodes);
105 update_desired_mesh_density.setFeatureResolution2D(m_FeatureResolution2D);
106 update_desired_mesh_density.setFeatureResolution3D(m_FeatureResolution3D);
107 update_desired_mesh_density();
110 void SurfaceAlgorithm::swap()
112 SwapTriangles swap;
113 swap.setGrid(m_Grid);
114 swap.setRespectBC(true);
115 swap.setFeatureSwap(m_AllowFeatureEdgeSwapping);
116 swap.setFeatureAngle(m_FeatureAngle);
117 swap.setMaxNumLoops(m_NumDelaunaySweeps);
118 swap.setSmallAreaSwap(m_AllowSmallAreaSwapping);
119 swap.setBCodesFeatureDefinition(m_BCodeFeatureDefinition);
120 QSet<int> rest_bcs = GuiMainWindow::pointer()->getAllBoundaryCodes();
121 rest_bcs -= m_BoundaryCodes;
122 swap.setBoundaryCodes(rest_bcs);
123 swap();
126 void SurfaceAlgorithm::smooth(int N_iter, bool correct_curveture)
128 LaplaceSmoother lap;
129 lap.setGrid(m_Grid);
130 QVector<vtkIdType> cls;
131 getSurfaceCells(m_BoundaryCodes, cls, m_Grid);
132 lap.setCells(cls);
133 lap.setNumberOfIterations(N_iter);
134 lap.setBoundaryCodes(m_BoundaryCodes);//IMPORTANT: so that unselected nodes become fixed when node types are updated!
135 lap.setCorrectCurvature(correct_curveture);
136 lap.setBCodesFeatureDefinition(m_BCodeFeatureDefinition);
137 if (m_UseProjectionForSmoothing) {
138 lap.setProjectionOn();
139 } else {
140 lap.setProjectionOff();
142 if (m_UseNormalCorrectionForSmoothing) {
143 lap.setNormalCorrectionOn();
144 } else {
145 lap.setNormalCorrectionOff();
147 lap();
148 m_SmoothSuccess = lap.succeeded();
151 int SurfaceAlgorithm::insertNodes()
153 if (m_InsertNodes) {
154 InsertPoints insert_points;
155 insert_points.setGrid(m_Grid);
156 insert_points.setBoundaryCodes(m_BoundaryCodes);
157 insert_points();
158 return insert_points.getNumInserted();
160 return 0;
163 int SurfaceAlgorithm::deleteNodes()
165 if (m_DeleteNodes) {
166 RemovePoints remove_points;
167 remove_points.setGrid(m_Grid);
168 remove_points.setBoundaryCodes(m_BoundaryCodes);
169 remove_points.setStretchingFactor(m_StretchingFactor);
170 remove_points.setFeatureAngle(m_FeatureAngle);
171 if (m_RespectFeatureEdgesForDeleteNodes) {
172 remove_points.setProtectFeatureEdgesOn();
173 } else {
174 remove_points.setProtectFeatureEdgesOff();
176 //remove_points.setFeatureAngle(m_FeatureAngleForDeleteNodes);
177 if (m_PerformGeometricTests) {
178 remove_points.setPerformGeometricChecksOn();
179 } else {
180 remove_points.setPerformGeometricChecksOff();
182 remove_points();
183 return remove_points.getNumRemoved();
185 return 0;