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 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670
|
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2017 OpenCFD Ltd.
-------------------------------------------------------------------------------
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
stitchMesh
Group
grpMeshManipulationUtilities
Description
'Stitches' a mesh.
Takes a mesh and two patches and merges the faces on the two patches
(if geometrically possible) so the faces become internal.
Can do
- 'perfect' match: faces and points on patches align exactly. Order might
be different though.
- 'integral' match: where the surfaces on both patches exactly
match but the individual faces not
- 'partial' match: where the non-overlapping part of the surface remains
in the respective patch.
Note : Is just a front-end to perfectInterface/slidingInterface.
Comparable to running a meshModifier of the form
(if masterPatch is called "M" and slavePatch "S"):
\verbatim
couple
{
type slidingInterface;
masterFaceZoneName MSMasterZone
slaveFaceZoneName MSSlaveZone
cutPointZoneName MSCutPointZone
cutFaceZoneName MSCutFaceZone
masterPatchName M;
slavePatchName S;
typeOfMatch partial or integral
}
\endverbatim
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "polyTopoChanger.H"
#include "mapPolyMesh.H"
#include "slidingInterface.H"
#include "perfectInterface.H"
#include "IOobjectList.H"
#include "ReadFields.H"
#include <numeric>
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Checks whether patch present and non-zero
bool checkPatch(const polyBoundaryMesh& bMesh, const word& name)
{
const label patchi = bMesh.findPatchID(name);
if (patchi == -1)
{
Info<< "No patch " << name << " in mesh" << nl
<< "Known patches: " << bMesh.names() << endl;
return false;
}
if (bMesh[patchi].empty())
{
Info<< "Patch " << name << " has zero size" << nl;
return false;
}
return true;
}
int main(int argc, char *argv[])
{
argList::addNote
(
"Merge the faces on specified patches (if geometrically possible)"
" so that the\n"
"faces become internal.\n"
"This utility can be called without arguments (uses stitchMeshDict)"
" or with\n"
"two arguments (master/slave patch names)."
);
argList::noParallel();
argList::noFunctionObjects(); // Never use function objects
#include "addOverwriteOption.H"
#include "addRegionOption.H"
argList::addOption("dict", "file", "Use alternative stitchMeshDict");
argList::addBoolOption
(
"integral",
"Couple integral master/slave patches (2 argument mode: default)"
);
argList::addBoolOption
(
"partial",
"Couple partially overlapping master/slave patches (2 argument mode)"
);
argList::addBoolOption
(
"perfect",
"Couple perfectly aligned master/slave patches (2 argument mode)"
);
argList::addBoolOption
(
"intermediate",
"Write intermediate stages, not just the final result"
);
argList::addOption
(
"toleranceDict",
"file",
"Dictionary file with tolerances"
);
// The arguments are optional (non-mandatory) when using dictionary mode
argList::noMandatoryArgs();
argList::addArgument
(
"master",
"The master patch name (non-dictionary mode)"
);
argList::addArgument
(
"slave",
"The slave patch name (non-dictionary mode)"
);
#include "setRootCase.H"
// We now handle checking args and general sanity etc.
const bool useCommandArgs = (args.size() > 1);
if (useCommandArgs)
{
if (args.found("dict"))
{
FatalErrorInFunction
<< "Cannot specify both dictionary and command-line arguments"
<< nl
<< endl;
}
// If we have arguments - we require all arguments!
if (!args.check(true, false))
{
FatalError.exit();
}
}
else
{
// Carp about inapplicable options
if (args.found("integral"))
{
FatalErrorInFunction
<< "Only specify -integral with command-line arguments"
<< endl;
}
if (args.found("partial"))
{
FatalErrorInFunction
<< "Only specify -partial with command-line arguments"
<< endl;
}
if (args.found("perfect"))
{
FatalErrorInFunction
<< "Only specify -perfect with command-line arguments"
<< endl;
}
}
#include "createTime.H"
#include "createNamedMesh.H"
const word oldInstance = mesh.pointsInstance();
const bool intermediate = args.found("intermediate");
const bool overwrite = args.found("overwrite");
const word dictName("stitchMeshDict");
// A validated input dictionary
dictionary validatedDict;
if (useCommandArgs)
{
// Command argument driven:
const int integralCover = args.found("integral");
const int partialCover = args.found("partial");
const int perfectCover = args.found("perfect");
if ((integralCover + partialCover + perfectCover) > 1)
{
FatalErrorInFunction
<< "Can only specify one of -integral | -partial | -perfect."
<< nl
<< "Use perfect match option if the patches perfectly align"
<< " (both vertex positions and face centres)" << endl
<< exit(FatalError);
}
// Patch names
const word masterPatchName(args[1]);
const word slavePatchName(args[2]);
// Patch names
Info<< " " << masterPatchName
<< " / " << slavePatchName << nl;
// Bail out if either patch has problems
if
(
!checkPatch(mesh.boundaryMesh(), masterPatchName)
|| !checkPatch(mesh.boundaryMesh(), slavePatchName)
)
{
FatalErrorInFunction
<< "Cannot continue"
<< exit(FatalError);
return 1;
}
// Input was validated
dictionary dict;
if (perfectCover)
{
dict.add("match", word("perfect"));
}
else if (partialCover)
{
dict.add
(
"match",
slidingInterface::typeOfMatchNames[slidingInterface::PARTIAL]
);
}
else
{
dict.add
(
"match",
slidingInterface::typeOfMatchNames[slidingInterface::INTEGRAL]
);
}
// Patch names
dict.add("master", masterPatchName);
dict.add("slave", slavePatchName);
validatedDict.add("stitchMesh", dict);
}
else
{
// dictionary-driven:
#include "setSystemRunTimeDictionaryIO.H"
Info<< "Reading " << dictName;
IOdictionary stitchDict(dictIO);
Info<< " with " << stitchDict.size() << " entries" << nl;
// Suppress duplicate names
wordHashSet requestedPatches;
for (const entry& dEntry : stitchDict)
{
if (!dEntry.isDict())
{
Info<< "Ignoring non-dictionary entry: "
<< dEntry.keyword() << nl;
continue;
}
const keyType& key = dEntry.keyword();
const dictionary& dict = dEntry.dict();
// Match type
word matchName;
if (dict.readIfPresent("match", matchName))
{
if
(
matchName != "perfect"
&& !slidingInterface::typeOfMatchNames.found(matchName)
)
{
Info<< "Error: unknown match type - " << matchName
<< " should be one of "
<< slidingInterface::typeOfMatchNames.toc() << nl;
continue;
}
}
// Patch names
const word masterPatchName(dict.get<word>("master"));
const word slavePatchName(dict.get<word>("slave"));
// Patch names
Info<< " " << masterPatchName
<< " / " << slavePatchName << nl;
if (!requestedPatches.insert(masterPatchName))
{
Info<< "Error: patch specified multiple times - "
<< masterPatchName << nl;
continue;
}
if (!requestedPatches.insert(slavePatchName))
{
Info<< "Error: patch specified multiple times - "
<< slavePatchName << nl;
requestedPatches.erase(masterPatchName);
continue;
}
// Bail out if either patch has problems
if
(
!checkPatch(mesh.boundaryMesh(), masterPatchName)
|| !checkPatch(mesh.boundaryMesh(), slavePatchName)
)
{
requestedPatches.erase(masterPatchName);
requestedPatches.erase(slavePatchName);
continue;
}
// Input was validated
validatedDict.add(key, dict);
}
}
const label nActions = validatedDict.size();
Info<< nl << nActions << " validated actions" << endl;
if (!nActions)
{
Info<<"\nStopping" << nl << endl;
return 1;
}
// ------------------------------------------
// This is where the real work begins
// set up the tolerances for the sliding mesh
dictionary slidingTolerances;
if (args.found("toleranceDict"))
{
IOdictionary toleranceFile
(
IOobject
(
args["toleranceDict"],
runTime.constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::NO_WRITE
)
);
slidingTolerances += toleranceFile;
}
// Search for list of objects for this time
IOobjectList objects(mesh, runTime.timeName());
// Read all current fvFields so they will get mapped
Info<< "Reading all current volfields" << endl;
PtrList<volScalarField> volScalarFields;
ReadFields(mesh, objects, volScalarFields);
PtrList<volVectorField> volVectorFields;
ReadFields(mesh, objects, volVectorFields);
PtrList<volSphericalTensorField> volSphericalTensorFields;
ReadFields(mesh, objects, volSphericalTensorFields);
PtrList<volSymmTensorField> volSymmTensorFields;
ReadFields(mesh, objects, volSymmTensorFields);
PtrList<volTensorField> volTensorFields;
ReadFields(mesh, objects, volTensorFields);
//- Uncomment if you want to interpolate surface fields (usually a bad idea)
//Info<< "Reading all current surfaceFields" << endl;
//PtrList<surfaceScalarField> surfaceScalarFields;
//ReadFields(mesh, objects, surfaceScalarFields);
//
//PtrList<surfaceVectorField> surfaceVectorFields;
//ReadFields(mesh, objects, surfaceVectorFields);
//
//PtrList<surfaceTensorField> surfaceTensorFields;
//ReadFields(mesh, objects, surfaceTensorFields);
// Increase precision for output mesh points
IOstream::defaultPrecision(max(10u, IOstream::defaultPrecision()));
polyTopoChanger stitcher(mesh, IOobject::NO_READ);
// Step through the topology changes
label actioni = 0;
for (const entry& dEntry : validatedDict)
{
const dictionary& dict = dEntry.dict();
// Match type
bool perfect = false;
slidingInterface::typeOfMatch matchType = slidingInterface::PARTIAL;
word matchName;
if (dict.readIfPresent("match", matchName))
{
if (matchName == "perfect")
{
perfect = true;
}
else
{
matchType = slidingInterface::typeOfMatchNames[matchName];
}
}
// Patch names
const word masterPatchName(dict.get<word>("master"));
const word slavePatchName(dict.get<word>("slave"));
// Zone names
const word mergePatchName(masterPatchName + slavePatchName);
const word cutZoneName(mergePatchName + "CutFaceZone");
Info<< nl << "========================================" << nl;
// Information messages
if (perfect)
{
Info<< "Coupling PERFECTLY aligned patches "
<< masterPatchName << " / " << slavePatchName << nl << nl
<< "Resulting (internal) faces in faceZone "
<< cutZoneName << nl << nl
<< "The patch vertices and face centres must align within a"
<< " tolerance relative to the minimum edge length on the patch"
<< nl << endl;
}
else if (matchType == slidingInterface::INTEGRAL)
{
Info<< "Coupling INTEGRALLY matching of patches "
<< masterPatchName << " / " << slavePatchName << nl << nl
<< "Resulting (internal) faces in faceZone "
<< cutZoneName << nl << nl
<< "The overall area covered by both patches should be"
<< " identical!" << endl
<< "If this is not the case use partial"
<< nl << endl;
}
else
{
Info<< "Coupling PARTIALLY overlapping patches "
<< masterPatchName << " / " << slavePatchName << nl
<< "Resulting internal faces in faceZone "
<< cutZoneName << nl
<< "Uncovered faces remain in their patch"
<< nl << endl;
}
// Master/slave patches
const polyPatch& masterPatch = mesh.boundaryMesh()[masterPatchName];
const polyPatch& slavePatch = mesh.boundaryMesh()[slavePatchName];
mesh.pointZones().clearAddressing();
mesh.faceZones().clearAddressing();
mesh.cellZones().clearAddressing();
// Lists of master and slave faces:
labelList faceIds;
// Markup master face ids
faceIds.setSize(masterPatch.size());
std::iota(faceIds.begin(), faceIds.end(), masterPatch.start());
stitcher.clear();
stitcher.setSize(1);
if (perfect)
{
// Add new (empty) zone for resulting internal faces
mesh.faceZones()
(
cutZoneName,
true // verbose
).resetAddressing(std::move(faceIds), false);
// Add the perfect interface mesh modifier
stitcher.set
(
0,
new perfectInterface
(
"couple" + Foam::name(actioni),
0,
stitcher,
cutZoneName,
masterPatchName,
slavePatchName
)
);
}
else
{
mesh.pointZones()
(
mergePatchName + "CutPointZone",
true // verbose
) = labelList();
mesh.faceZones()
(
mergePatchName + "MasterZone",
true // verbose
).resetAddressing(std::move(faceIds), false);
// Markup slave face ids
faceIds.setSize(slavePatch.size());
std::iota(faceIds.begin(), faceIds.end(), slavePatch.start());
mesh.faceZones()
(
mergePatchName + "SlaveZone",
true // verbose
).resetAddressing(std::move(faceIds), false);
// Add empty zone for cut faces
mesh.faceZones()
(
cutZoneName,
true // verbose
).resetAddressing(labelList(), false);
// Add the sliding interface mesh modifier
stitcher.set
(
0,
new slidingInterface
(
"couple" + Foam::name(actioni),
0,
stitcher,
mergePatchName + "MasterZone",
mergePatchName + "SlaveZone",
mergePatchName + "CutPointZone",
cutZoneName,
masterPatchName,
slavePatchName,
matchType, // integral or partial
true // couple/decouple mode
)
);
static_cast<slidingInterface&>(stitcher[0]).setTolerances
(
slidingTolerances,
true
);
}
++actioni;
// Advance time for intermediate results or only on final
if (!overwrite && (intermediate || actioni == nActions))
{
++runTime;
}
// Execute all polyMeshModifiers
autoPtr<mapPolyMesh> morphMap = stitcher.changeMesh(true);
mesh.movePoints(morphMap->preMotionPoints());
// Write mesh
if (overwrite)
{
mesh.setInstance(oldInstance);
stitcher.instance() = oldInstance;
}
if (intermediate || actioni == nActions)
{
Info<< nl << "Writing polyMesh to time "
<< runTime.timeName() << endl;
// Bypass runTime write (since only writes at writeTime)
if
(
!runTime.objectRegistry::writeObject
(
runTime.writeFormat(),
IOstream::currentVersion,
runTime.writeCompression(),
true
)
)
{
FatalErrorInFunction
<< "Failed writing polyMesh."
<< exit(FatalError);
}
mesh.faceZones().write();
mesh.pointZones().write();
mesh.cellZones().write();
// Write fields
runTime.write();
}
}
Info<< "\nEnd\n" << endl;
return 0;
}
// ************************************************************************* //
|