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 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
combinePatchFaces
Description
Checks for multiple patch faces on same cell and combines them.
Multiple patch faces can result from e.g. removal of refined
neighbouring cells, leaving 4 exposed faces with same owner.
Rules for merging:
- only boundary faces (since multiple internal faces between two cells
not allowed anyway)
- faces have to have same owner
- faces have to be connected via edge which are not features (so angle
between them < feature angle)
- outside of faces has to be single loop
- outside of face should not be (or just slightly) concave (so angle
between consecutive edges < concaveangle
E.g. to allow all faces on same patch to be merged:
combinePatchFaces 180 -concaveAngle 90
\*---------------------------------------------------------------------------*/
#include "PstreamReduceOps.H"
#include "argList.H"
#include "Time.H"
#include "polyTopoChange.H"
#include "polyModifyFace.H"
#include "polyAddFace.H"
#include "combineFaces.H"
#include "removePoints.H"
#include "polyMesh.H"
#include "mapPolyMesh.H"
#include "unitConversion.H"
#include "motionSmoother.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Merge faces on the same patch (usually from exposing refinement)
// Can undo merges if these cause problems.
label mergePatchFaces
(
const scalar minCos,
const scalar concaveSin,
const autoPtr<IOdictionary>& qualDictPtr,
const Time& runTime,
polyMesh& mesh
)
{
// Patch face merging engine
combineFaces faceCombiner(mesh);
// Get all sets of faces that can be merged
labelListList allFaceSets(faceCombiner.getMergeSets(minCos, concaveSin));
label nFaceSets = returnReduce(allFaceSets.size(), sumOp<label>());
Info<< "Merging " << nFaceSets << " sets of faces." << endl;
if (nFaceSets > 0)
{
// Store the faces of the face sets
List<faceList> allFaceSetsFaces(allFaceSets.size());
forAll(allFaceSets, setI)
{
allFaceSetsFaces[setI] = UIndirectList<face>
(
mesh.faces(),
allFaceSets[setI]
);
}
autoPtr<mapPolyMesh> map;
{
// Topology changes container
polyTopoChange meshMod(mesh);
// Merge all faces of a set into the first face of the set.
faceCombiner.setRefinement(allFaceSets, meshMod);
// Change the mesh (no inflation)
map = meshMod.changeMesh(mesh, false, true);
// Update fields
mesh.updateMesh(map);
// Move mesh (since morphing does not do this)
if (map().hasMotionPoints())
{
mesh.movePoints(map().preMotionPoints());
}
else
{
// Delete mesh volumes. No other way to do this?
mesh.clearOut();
}
}
// Check for errors and undo
// ~~~~~~~~~~~~~~~~~~~~~~~~~
// Faces in error.
labelHashSet errorFaces;
if (qualDictPtr.valid())
{
motionSmoother::checkMesh(false, mesh, qualDictPtr(), errorFaces);
}
else
{
mesh.checkFacePyramids(false, -SMALL, &errorFaces);
}
// Sets where the master is in error
labelHashSet errorSets;
forAll(allFaceSets, setI)
{
label newMasterI = map().reverseFaceMap()[allFaceSets[setI][0]];
if (errorFaces.found(newMasterI))
{
errorSets.insert(setI);
}
}
label nErrorSets = returnReduce(errorSets.size(), sumOp<label>());
Info<< "Detected " << nErrorSets
<< " error faces on boundaries that have been merged."
<< " These will be restored to their original faces."
<< endl;
if (nErrorSets > 0)
{
// Renumber stored faces to new vertex numbering.
forAllConstIter(labelHashSet, errorSets, iter)
{
label setI = iter.key();
faceList& setFaceVerts = allFaceSetsFaces[setI];
forAll(setFaceVerts, i)
{
inplaceRenumber(map().reversePointMap(), setFaceVerts[i]);
// Debug: check that all points are still there.
forAll(setFaceVerts[i], j)
{
label newVertI = setFaceVerts[i][j];
if (newVertI < 0)
{
FatalErrorInFunction
<< "In set:" << setI << " old face labels:"
<< allFaceSets[setI] << " new face vertices:"
<< setFaceVerts[i] << " are unmapped vertices!"
<< abort(FatalError);
}
}
}
}
// Topology changes container
polyTopoChange meshMod(mesh);
// Restore faces
forAllConstIter(labelHashSet, errorSets, iter)
{
label setI = iter.key();
const labelList& setFaces = allFaceSets[setI];
const faceList& setFaceVerts = allFaceSetsFaces[setI];
label newMasterI = map().reverseFaceMap()[setFaces[0]];
// Restore. Get face properties.
label own = mesh.faceOwner()[newMasterI];
label zoneID = mesh.faceZones().whichZone(newMasterI);
bool zoneFlip = false;
if (zoneID >= 0)
{
const faceZone& fZone = mesh.faceZones()[zoneID];
zoneFlip = fZone.flipMap()[fZone.whichFace(newMasterI)];
}
label patchID = mesh.boundaryMesh().whichPatch(newMasterI);
Pout<< "Restoring new master face " << newMasterI
<< " to vertices " << setFaceVerts[0] << endl;
// Modify the master face.
meshMod.setAction
(
polyModifyFace
(
setFaceVerts[0], // original face
newMasterI, // label of face
own, // owner
-1, // neighbour
false, // face flip
patchID, // patch for face
false, // remove from zone
zoneID, // zone for face
zoneFlip // face flip in zone
)
);
// Add the previously removed faces
for (label i = 1; i < setFaces.size(); i++)
{
Pout<< "Restoring removed face " << setFaces[i]
<< " with vertices " << setFaceVerts[i] << endl;
meshMod.setAction
(
polyAddFace
(
setFaceVerts[i], // vertices
own, // owner,
-1, // neighbour,
-1, // masterPointID,
-1, // masterEdgeID,
newMasterI, // masterFaceID,
false, // flipFaceFlux,
patchID, // patchID,
zoneID, // zoneID,
zoneFlip // zoneFlip
)
);
}
}
// Change the mesh (no inflation)
map = meshMod.changeMesh(mesh, false, true);
// Update fields
mesh.updateMesh(map);
// Move mesh (since morphing does not do this)
if (map().hasMotionPoints())
{
mesh.movePoints(map().preMotionPoints());
}
else
{
// Delete mesh volumes. No other way to do this?
mesh.clearOut();
}
}
}
else
{
Info<< "No faces merged ..." << endl;
}
return nFaceSets;
}
// Remove points not used by any face or points used by only two faces where
// the edges are in line
label mergeEdges(const scalar minCos, polyMesh& mesh)
{
Info<< "Merging all points on surface that" << nl
<< "- are used by only two boundary faces and" << nl
<< "- make an angle with a cosine of more than " << minCos
<< "." << nl << endl;
// Point removal analysis engine
removePoints pointRemover(mesh);
// Count usage of points
boolList pointCanBeDeleted;
label nRemove = pointRemover.countPointUsage(minCos, pointCanBeDeleted);
if (nRemove > 0)
{
Info<< "Removing " << nRemove
<< " straight edge points ..." << endl;
// Topology changes container
polyTopoChange meshMod(mesh);
pointRemover.setRefinement(pointCanBeDeleted, meshMod);
// Change the mesh (no inflation)
autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false, true);
// Update fields
mesh.updateMesh(map);
// Move mesh (since morphing does not do this)
if (map().hasMotionPoints())
{
mesh.movePoints(map().preMotionPoints());
}
else
{
// Delete mesh volumes. No other way to do this?
mesh.clearOut();
}
}
else
{
Info<< "No straight edges simplified and no points removed ..." << endl;
}
return nRemove;
}
int main(int argc, char *argv[])
{
#include "addOverwriteOption.H"
argList::validArgs.append("featureAngle [0..180]");
argList::addOption
(
"concaveAngle",
"degrees",
"specify concave angle [0..180] (default: 30 degrees)"
);
argList::addBoolOption
(
"meshQuality",
"read user-defined mesh quality criterions from system/meshQualityDict"
);
#include "setRootCase.H"
#include "createTime.H"
runTime.functionObjects().off();
#include "createPolyMesh.H"
const word oldInstance = mesh.pointsInstance();
const scalar featureAngle = args.argRead<scalar>(1);
const scalar minCos = Foam::cos(degToRad(featureAngle));
// Sin of angle between two consecutive edges on a face.
// If sin(angle) larger than this the face will be considered concave.
scalar concaveAngle = args.optionLookupOrDefault("concaveAngle", 30.0);
scalar concaveSin = Foam::sin(degToRad(concaveAngle));
const bool overwrite = args.optionFound("overwrite");
const bool meshQuality = args.optionFound("meshQuality");
Info<< "Merging all faces of a cell" << nl
<< " - which are on the same patch" << nl
<< " - which make an angle < " << featureAngle << " degrees"
<< nl
<< " (cos:" << minCos << ')' << nl
<< " - even when resulting face becomes concave by more than "
<< concaveAngle << " degrees" << nl
<< " (sin:" << concaveSin << ')' << nl
<< endl;
autoPtr<IOdictionary> qualDict;
if (meshQuality)
{
Info<< "Enabling user-defined geometry checks." << nl << endl;
qualDict.reset
(
new IOdictionary
(
IOobject
(
"meshQualityDict",
mesh.time().system(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE
)
)
);
}
if (!overwrite)
{
runTime++;
}
// Merge faces on same patch
label nChanged = mergePatchFaces
(
minCos,
concaveSin,
qualDict,
runTime,
mesh
);
// Merge points on straight edges and remove unused points
if (qualDict.valid())
{
Info<< "Merging all 'loose' points on surface edges, "
<< "regardless of the angle they make." << endl;
// Surface bnound to be used to extrude. Merge all loose points.
nChanged += mergeEdges(-1, mesh);
}
else
{
nChanged += mergeEdges(minCos, mesh);
}
if (nChanged > 0)
{
if (overwrite)
{
mesh.setInstance(oldInstance);
}
Info<< "Writing morphed mesh to time " << runTime.timeName() << endl;
mesh.write();
}
else
{
Info<< "Mesh unchanged." << endl;
}
Info<< "\nEnd\n" << endl;
return 0;
}
// ************************************************************************* //
|