00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #include "autoHexMeshDriver.H"
00028 #include "fvMesh.H"
00029 #include "Time.H"
00030 #include "boundBox.H"
00031 #include "globalIndex.H"
00032 #include "wallPolyPatch.H"
00033 #include "mapDistributePolyMesh.H"
00034 #include "cellSet.H"
00035 #include "syncTools.H"
00036 #include "motionSmoother.H"
00037 #include "refinementParameters.H"
00038 #include "snapParameters.H"
00039 #include "layerParameters.H"
00040 #include "autoRefineDriver.H"
00041 #include "autoSnapDriver.H"
00042 #include "autoLayerDriver.H"
00043 #include "triSurfaceMesh.H"
00044
00045
00046
00047 namespace Foam
00048 {
00049 defineTypeNameAndDebug(autoHexMeshDriver, 0);
00050 }
00051
00052
00053
00054
00055
00056 Foam::scalar Foam::autoHexMeshDriver::getMergeDistance(const scalar mergeTol)
00057 const
00058 {
00059 const boundBox& meshBb = mesh_.bounds();
00060 scalar mergeDist = mergeTol * meshBb.mag();
00061 scalar writeTol = std::pow
00062 (
00063 scalar(10.0),
00064 -scalar(IOstream::defaultPrecision())
00065 );
00066
00067 Info<< nl
00068 << "Overall mesh bounding box : " << meshBb << nl
00069 << "Relative tolerance : " << mergeTol << nl
00070 << "Absolute matching distance : " << mergeDist << nl
00071 << endl;
00072
00073 if (mesh_.time().writeFormat() == IOstream::ASCII && mergeTol < writeTol)
00074 {
00075 FatalErrorIn("autoHexMeshDriver::getMergeDistance(const scalar) const")
00076 << "Your current settings specify ASCII writing with "
00077 << IOstream::defaultPrecision() << " digits precision." << endl
00078 << "Your merging tolerance (" << mergeTol << ") is finer than this."
00079 << endl
00080 << "Please change your writeFormat to binary"
00081 << " or increase the writePrecision" << endl
00082 << "or adjust the merge tolerance (-mergeTol)."
00083 << exit(FatalError);
00084 }
00085
00086 return mergeDist;
00087 }
00088
00089
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148 Foam::autoHexMeshDriver::autoHexMeshDriver
00149 (
00150 fvMesh& mesh,
00151 const bool overwrite,
00152 const dictionary& dict,
00153 const dictionary& decomposeDict
00154 )
00155 :
00156 mesh_(mesh),
00157 dict_(dict),
00158 debug_(readLabel(dict_.lookup("debug"))),
00159 mergeDist_(getMergeDistance(readScalar(dict_.lookup("mergeTolerance"))))
00160 {
00161 if (debug_ > 0)
00162 {
00163 meshRefinement::debug = debug_;
00164 autoHexMeshDriver::debug = debug_;
00165 autoRefineDriver::debug = debug;
00166 autoSnapDriver::debug = debug;
00167 autoLayerDriver::debug = debug;
00168 }
00169
00170 refinementParameters refineParams(dict, 1);
00171
00172 Info<< "Overall cell limit : "
00173 << refineParams.maxGlobalCells() << endl;
00174 Info<< "Per processor cell limit : "
00175 << refineParams.maxLocalCells() << endl;
00176 Info<< "Minimum number of cells to refine : "
00177 << refineParams.minRefineCells() << endl;
00178 Info<< "Curvature : "
00179 << refineParams.curvature() << nl << endl;
00180 Info<< "Layers between different refinement levels : "
00181 << refineParams.nBufferLayers() << endl;
00182
00183 PtrList<dictionary> shellDicts(dict_.lookup("refinementShells"));
00184
00185 PtrList<dictionary> surfaceDicts(dict_.lookup("surfaces"));
00186
00187
00188
00189
00190
00191 {
00192 Info<< "Reading all geometry." << endl;
00193
00194
00195 dictionary geometryDict;
00196
00197 forAll(shellDicts, shellI)
00198 {
00199 dictionary shellDict = shellDicts[shellI];
00200 const word name(shellDict.lookup("name"));
00201 shellDict.remove("name");
00202 shellDict.remove("level");
00203 shellDict.remove("refineInside");
00204 geometryDict.add(name, shellDict);
00205 }
00206
00207 forAll(surfaceDicts, surfI)
00208 {
00209 dictionary surfDict = surfaceDicts[surfI];
00210 const word name(string::validate<word>(surfDict.lookup("file")));
00211 surfDict.remove("file");
00212 surfDict.remove("regions");
00213 if (!surfDict.found("name"))
00214 {
00215 surfDict.add("name", name);
00216 }
00217 surfDict.add("type", triSurfaceMesh::typeName);
00218 geometryDict.add(name, surfDict);
00219 }
00220
00221 allGeometryPtr_.reset
00222 (
00223 new searchableSurfaces
00224 (
00225 IOobject
00226 (
00227 "abc",
00228
00229
00230 mesh_.time().constant(),
00231 "triSurface",
00232 mesh_.time(),
00233 IOobject::MUST_READ,
00234 IOobject::NO_WRITE
00235 ),
00236 geometryDict
00237 )
00238 );
00239
00240 Info<< "Read geometry in = "
00241 << mesh_.time().cpuTimeIncrement() << " s" << endl;
00242 }
00243
00244
00245
00246
00247
00248 {
00249 Info<< "Reading surfaces and constructing search trees." << endl;
00250
00251 surfacesPtr_.reset
00252 (
00253 new refinementSurfaces
00254 (
00255 allGeometryPtr_(),
00256 surfaceDicts
00257 )
00258 );
00259 Info<< "Read surfaces in = "
00260 << mesh_.time().cpuTimeIncrement() << " s" << endl;
00261 }
00262
00263
00264
00265
00266 {
00267 Info<< "Reading refinement shells." << endl;
00268 shellsPtr_.reset
00269 (
00270 new shellSurfaces
00271 (
00272 allGeometryPtr_(),
00273 shellDicts
00274 )
00275 );
00276 Info<< "Read refinement shells in = "
00277 << mesh_.time().cpuTimeIncrement() << " s" << endl;
00278
00280
00281
00282
00283
00284
00285
00286 Info<< "Setting refinement level of surface to be consistent"
00287 << " with shells." << endl;
00288 surfacesPtr_().setMinLevelFields(shells());
00289 Info<< "Checked shell refinement in = "
00290 << mesh_.time().cpuTimeIncrement() << " s" << endl;
00291 }
00292
00293
00294
00295 meshRefinement::checkCoupledFaceZones(mesh_);
00296
00297
00298
00299
00300
00301 {
00302 Info<< nl
00303 << "Determining initial surface intersections" << nl
00304 << "-----------------------------------------" << nl
00305 << endl;
00306
00307
00308 meshRefinerPtr_.reset
00309 (
00310 new meshRefinement
00311 (
00312 mesh,
00313 mergeDist_,
00314 overwrite,
00315 surfaces(),
00316 shells()
00317 )
00318 );
00319 Info<< "Calculated surface intersections in = "
00320 << mesh_.time().cpuTimeIncrement() << " s" << endl;
00321
00322
00323 meshRefinerPtr_().printMeshInfo(debug_, "Initial mesh");
00324
00325 meshRefinerPtr_().write
00326 (
00327 debug_&meshRefinement::OBJINTERSECTIONS,
00328 mesh_.time().path()/meshRefinerPtr_().timeName()
00329 );
00330 }
00331
00332
00333
00334
00335
00336 {
00337 Info<< nl
00338 << "Adding patches for surface regions" << nl
00339 << "----------------------------------" << nl
00340 << endl;
00341
00342
00343 globalToPatch_.setSize(surfaces().nRegions(), -1);
00344
00345 Info<< "Patch\tRegion" << nl
00346 << "-----\t------"
00347 << endl;
00348
00349 const labelList& surfaceGeometry = surfaces().surfaces();
00350 forAll(surfaceGeometry, surfI)
00351 {
00352 label geomI = surfaceGeometry[surfI];
00353
00354 const wordList& regNames = allGeometryPtr_().regionNames()[geomI];
00355
00356 Info<< surfaces().names()[surfI] << ':' << nl << nl;
00357
00358 forAll(regNames, i)
00359 {
00360 label patchI = meshRefinerPtr_().addMeshedPatch
00361 (
00362 regNames[i],
00363 wallPolyPatch::typeName
00364 );
00365
00366 Info<< patchI << '\t' << regNames[i] << nl;
00367
00368 globalToPatch_[surfaces().globalRegion(surfI, i)] = patchI;
00369 }
00370
00371 Info<< nl;
00372 }
00373 Info<< "Added patches in = "
00374 << mesh_.time().cpuTimeIncrement() << " s" << nl << endl;
00375 }
00376
00377
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420 {
00421
00422 decomposerPtr_ = decompositionMethod::New
00423 (
00424 decomposeDict,
00425 mesh_
00426 );
00427 decompositionMethod& decomposer = decomposerPtr_();
00428
00429
00430 if (Pstream::parRun() && !decomposer.parallelAware())
00431 {
00432 FatalErrorIn("autoHexMeshDriver::autoHexMeshDriver"
00433 "(const IOobject&, fvMesh&)")
00434 << "You have selected decomposition method "
00435 << decomposer.typeName
00436 << " which is not parallel aware." << endl
00437 << "Please select one that is (parMetis, hierarchical)"
00438 << exit(FatalError);
00439 }
00440
00441
00442 distributorPtr_.reset(new fvMeshDistribute(mesh_, mergeDist_));
00443 }
00444 }
00445
00446
00447
00448
00449 void Foam::autoHexMeshDriver::writeMesh(const string& msg) const
00450 {
00451 const meshRefinement& meshRefiner = meshRefinerPtr_();
00452
00453 meshRefiner.printMeshInfo(debug_, msg);
00454 Info<< "Writing mesh to time " << meshRefiner.timeName() << endl;
00455
00456 meshRefiner.write(meshRefinement::MESH|meshRefinement::SCALARLEVELS, "");
00457 if (debug_ & meshRefinement::OBJINTERSECTIONS)
00458 {
00459 meshRefiner.write
00460 (
00461 meshRefinement::OBJINTERSECTIONS,
00462 mesh_.time().path()/meshRefiner.timeName()
00463 );
00464 }
00465 Info<< "Written mesh in = "
00466 << mesh_.time().cpuTimeIncrement() << " s." << endl;
00467 }
00468
00469
00470 void Foam::autoHexMeshDriver::doMesh()
00471 {
00472 Switch wantRefine(dict_.lookup("doRefine"));
00473 Switch wantSnap(dict_.lookup("doSnap"));
00474 Switch wantLayers(dict_.lookup("doLayers"));
00475
00476 Info<< "Do refinement : " << wantRefine << nl
00477 << "Do snapping : " << wantSnap << nl
00478 << "Do layers : " << wantLayers << nl
00479 << endl;
00480
00481 if (wantRefine)
00482 {
00483 const dictionary& motionDict = dict_.subDict("motionDict");
00484
00485 autoRefineDriver refineDriver
00486 (
00487 meshRefinerPtr_(),
00488 decomposerPtr_(),
00489 distributorPtr_(),
00490 globalToPatch_
00491 );
00492
00493
00494 refinementParameters refineParams(dict_, 1);
00495
00496 refineDriver.doRefine(dict_, refineParams, wantSnap, motionDict);
00497
00498
00499 writeMesh("Refined mesh");
00500 }
00501
00502 if (wantSnap)
00503 {
00504 const dictionary& snapDict = dict_.subDict("snapDict");
00505 const dictionary& motionDict = dict_.subDict("motionDict");
00506
00507 autoSnapDriver snapDriver
00508 (
00509 meshRefinerPtr_(),
00510 globalToPatch_
00511 );
00512
00513
00514 snapParameters snapParams(snapDict, 1);
00515
00516 snapDriver.doSnap(snapDict, motionDict, snapParams);
00517
00518
00519 writeMesh("Snapped mesh");
00520 }
00521
00522 if (wantLayers)
00523 {
00524 const dictionary& motionDict = dict_.subDict("motionDict");
00525 const dictionary& shrinkDict = dict_.subDict("shrinkDict");
00526 PtrList<dictionary> surfaceDicts(dict_.lookup("surfaces"));
00527
00528 autoLayerDriver layerDriver(meshRefinerPtr_());
00529
00530
00531 layerParameters layerParams
00532 (
00533 surfaceDicts,
00534 surfacesPtr_(),
00535 globalToPatch_,
00536 shrinkDict,
00537 mesh_.boundaryMesh()
00538 );
00539
00540 layerDriver.doLayers
00541 (
00542 shrinkDict,
00543 motionDict,
00544 layerParams,
00545 decomposerPtr_(),
00546 distributorPtr_()
00547 );
00548
00549
00550 writeMesh("Layer mesh");
00551 }
00552 }
00553
00554
00555