snappyHexMesh: Write correct refinement files once only

The files relating to the hex refinement are written out explicitly both by
snappyHexMesh and dynamicRefineFvMesh and hence should be set "NO_WRITE" rather
than "AUTO_WRITE" to avoid writing them twice.  This change corrects the
handling of the "refinementHistory" file which should not be written by
snappyHexMesh.
This commit is contained in:
Henry Weller
2017-01-24 08:15:43 +00:00
parent dd0dddf46a
commit 9129f5afc3
3 changed files with 25 additions and 43 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -1938,7 +1938,7 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh, const bool readHistory)
polyMesh::meshSubDir,
mesh_,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
labelList(mesh_.nCells(), 0)
),
@ -1951,7 +1951,7 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh, const bool readHistory)
polyMesh::meshSubDir,
mesh_,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
labelList(mesh_.nPoints(), 0)
),
@ -1964,7 +1964,7 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh, const bool readHistory)
polyMesh::meshSubDir,
mesh_,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
dimensionedScalar("level0Edge", dimLength, getLevel0EdgeLength())
),
@ -1977,9 +1977,10 @@ Foam::hexRef8::hexRef8(const polyMesh& mesh, const bool readHistory)
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
(readHistory ? mesh_.nCells() : 0) // All cells visible if not be read
// All cells visible if not read or readHistory = false
(readHistory ? mesh_.nCells() : 0)
),
faceRemover_(mesh_, GREAT), // merge boundary faces wherever possible
savedPointLevel_(0),
@ -2063,7 +2064,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
cellLevel
),
@ -2076,7 +2077,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
pointLevel
),
@ -2089,7 +2090,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
dimensionedScalar
(
@ -2107,7 +2108,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
history
),
@ -2171,7 +2172,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
cellLevel
),
@ -2184,7 +2185,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
pointLevel
),
@ -2197,7 +2198,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
dimensionedScalar
(
@ -2215,7 +2216,7 @@ Foam::hexRef8::hexRef8
polyMesh::meshSubDir,
mesh_,
IOobject::NO_READ,
IOobject::AUTO_WRITE
IOobject::NO_WRITE
),
List<refinementHistory::splitCell8>(0),
labelList(0),
@ -3075,7 +3076,7 @@ Foam::labelList Foam::hexRef8::consistentSlowRefinement2
// fMesh.time().timeName(),
// fMesh,
// IOobject::NO_READ,
// IOobject::AUTO_WRITE,
// IOobject::NO_WRITE,
// false
// ),
// fMesh,

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -335,7 +335,6 @@ void Foam::refinementHistory::freeSplitCell(const label index)
}
// Mark entry in splitCells. Recursively mark its parent and subs.
void Foam::refinementHistory::markSplit
(
const label index,
@ -372,7 +371,6 @@ void Foam::refinementHistory::markSplit
}
// Mark index and all its descendants
void Foam::refinementHistory::mark
(
const label val,
@ -762,7 +760,6 @@ Foam::refinementHistory::refinementHistory
}
// Construct as copy
Foam::refinementHistory::refinementHistory
(
const IOobject& io,
@ -783,7 +780,6 @@ Foam::refinementHistory::refinementHistory
}
// Construct from multiple
Foam::refinementHistory::refinementHistory
(
const IOobject& io,
@ -903,7 +899,6 @@ Foam::refinementHistory::refinementHistory
}
// Construct from Istream
Foam::refinementHistory::refinementHistory(const IOobject& io, Istream& is)
:
regIOobject(io),
@ -1180,7 +1175,6 @@ void Foam::refinementHistory::updateMesh(const mapPolyMesh& map)
}
// Update numbering for subsetting
void Foam::refinementHistory::subset
(
const labelList& pointMap,

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -207,8 +207,6 @@ void Foam::meshRefinement::calcNeighbourData
}
// Find intersections of edges (given by their two endpoints) with surfaces.
// Returns first intersection if there are more than one.
void Foam::meshRefinement::updateIntersections(const labelList& changedFaces)
{
const pointField& cellCentres = mesh_.cellCentres();
@ -609,8 +607,6 @@ void Foam::meshRefinement::setInstance(const fileName& inst)
}
// Remove cells. Put exposedFaces (output of getExposedFaces(cellsToRemove))
// into exposedPatchIDs.
Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::doRemoveCells
(
const labelList& cellsToRemove,
@ -670,7 +666,6 @@ Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::doRemoveCells
}
// Split faces
Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::splitFaces
(
const labelList& splitFaces,
@ -684,7 +679,6 @@ Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::splitFaces
label facei = splitFaces[i];
const face& f = mesh_.faces()[facei];
// Split as start and end index in face
const labelPair& split = splits[i];
@ -737,9 +731,12 @@ Foam::autoPtr<Foam::mapPolyMesh> Foam::meshRefinement::splitFaces
}
Pout<< "face:" << facei << " verts:" << f
<< " split into f0:" << f0
<< " f1:" << f1 << endl;
if (debug)
{
Pout<< "face:" << facei << " verts:" << f
<< " split into f0:" << f0
<< " f1:" << f1 << endl;
}
// Change/add faces
meshMod.modifyFace
@ -753,6 +750,7 @@ Pout<< "face:" << facei << " verts:" << f
zoneI, // zone for face
zoneFlip // face flip in zone
);
meshMod.addFace
(
f1, // modified face
@ -1160,7 +1158,6 @@ Pout<< "face:" << facei << " verts:" << f
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
// Construct from components
Foam::meshRefinement::meshRefinement
(
fvMesh& mesh,
@ -1381,11 +1378,6 @@ Foam::autoPtr<Foam::mapDistributePolyMesh> Foam::meshRefinement::balance
if (Pstream::parRun())
{
//if (debug_)
//{
// const_cast<Time&>(mesh_.time())++;
//}
// Wanted distribution
labelList distribution;
@ -1648,7 +1640,6 @@ Foam::autoPtr<Foam::mapDistributePolyMesh> Foam::meshRefinement::balance
}
// Helper function to get intersected faces
Foam::labelList Foam::meshRefinement::intersectedFaces() const
{
label nBoundaryFaces = 0;
@ -1675,7 +1666,6 @@ Foam::labelList Foam::meshRefinement::intersectedFaces() const
}
// Helper function to get points used by faces
Foam::labelList Foam::meshRefinement::intersectedPoints() const
{
const faceList& faces = mesh_.faces();
@ -1788,7 +1778,6 @@ Foam::autoPtr<Foam::indirectPrimitivePatch> Foam::meshRefinement::makePatch
}
// Construct pointVectorField with correct boundary conditions
Foam::tmp<Foam::pointVectorField> Foam::meshRefinement::makeDisplacementField
(
const pointMesh& pMesh,
@ -2088,7 +2077,6 @@ Foam::label Foam::meshRefinement::appendPatch
}
// Adds patch if not yet there. Returns patchID.
Foam::label Foam::meshRefinement::addPatch
(
fvMesh& mesh,
@ -2451,7 +2439,6 @@ void Foam::meshRefinement::distribute(const mapDistributePolyMesh& map)
}
// Update local data for a mesh change
void Foam::meshRefinement::updateMesh
(
const mapPolyMesh& map,