2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2013-12-11 16:09:41 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2016-06-20 15:00:40 +00:00
|
|
|
\\ / O peration | Version: 4.0
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2010-05-12 13:27:55 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2010-05-12 13:27:55 +00:00
|
|
|
under the terms of the GNU General Public License as published by the
|
2013-12-11 16:09:41 +00:00
|
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
2010-05-12 13:27:55 +00:00
|
|
|
option) any later version.
|
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend 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.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
2013-12-11 16:09:41 +00:00
|
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Application
|
|
|
|
reconstructParMesh
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
Author
|
|
|
|
Hrvoje Jasak, Wikki Ltd. All rights reserved
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
Description
|
|
|
|
Reconstructs a mesh using geometrical matching and catenation.
|
|
|
|
Use following topological changes in parallel to create global mesh
|
|
|
|
and xxxxProcAddressing files in the processor meshes.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "processorMeshesReconstructor.H"
|
|
|
|
#include "fvFieldReconstructor.H"
|
|
|
|
#include "pointFieldReconstructor.H"
|
|
|
|
#include "tetPointFieldReconstructor.H"
|
|
|
|
#include "reconstructLagrangian.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
#include "faCFD.H"
|
|
|
|
#include "faMesh.H"
|
|
|
|
#include "processorFaMeshes.H"
|
|
|
|
#include "faFieldReconstructor.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
// enable -constant ... if someone really wants it
|
|
|
|
// enable -zeroTime to prevent accidentally trashing the initial fields
|
2011-02-25 12:45:44 +00:00
|
|
|
timeSelector::addOptions(false, true);
|
2011-02-19 01:05:24 +00:00
|
|
|
argList::noParallel();
|
|
|
|
# include "addRegionOption.H"
|
|
|
|
argList::validOptions.insert("cellDist", "");
|
|
|
|
argList::validOptions.insert("fields", "\"(list of fields)\"");
|
|
|
|
argList::validOptions.insert("noLagrangian", "");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
# include "setRootCase.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
bool writeCellDist = args.optionFound("cellDist");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
# include "createTime.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
HashSet<word> selectedFields;
|
|
|
|
if (args.optionFound("fields"))
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
args.optionLookup("fields")() >> selectedFields;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
bool noLagrangian = args.optionFound("noLagrangian");
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Determine the processor count directly
|
|
|
|
label nProcs = 0;
|
|
|
|
while (isDir(args.path()/(word("processor") + name(nProcs))))
|
|
|
|
{
|
|
|
|
++nProcs;
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if (!nProcs)
|
|
|
|
{
|
|
|
|
FatalErrorIn(args.executable())
|
|
|
|
<< "No processor* directories found"
|
|
|
|
<< exit(FatalError);
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Create the processor databases
|
|
|
|
PtrList<Time> databases(nProcs);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
forAll (databases, procI)
|
|
|
|
{
|
|
|
|
databases.set
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
procI,
|
|
|
|
new Time
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
Time::controlDictName,
|
|
|
|
args.rootPath(),
|
|
|
|
args.caseName()/fileName(word("processor") + name(procI))
|
2010-05-12 13:27:55 +00:00
|
|
|
)
|
|
|
|
);
|
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// use the times list from the master processor
|
|
|
|
// and select a subset based on the command-line options
|
|
|
|
instantList timeDirs = timeSelector::select
|
2010-08-26 14:22:03 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
databases[0].times(),
|
|
|
|
args
|
2010-08-26 14:22:03 +00:00
|
|
|
);
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if (timeDirs.empty())
|
2010-08-26 14:22:03 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
FatalErrorIn(args.executable())
|
|
|
|
<< "No times selected"
|
|
|
|
<< exit(FatalError);
|
2010-08-26 14:22:03 +00:00
|
|
|
}
|
|
|
|
|
2017-01-11 12:03:23 +00:00
|
|
|
word regionName = polyMesh::defaultRegion;
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if (args.optionReadIfPresent("region", regionName))
|
2010-08-26 14:22:03 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
Info<< "Selecting region " << regionName << " for time = "
|
2017-01-11 12:03:23 +00:00
|
|
|
<< runTime.timeName() << nl << endl;
|
2011-02-19 01:05:24 +00:00
|
|
|
}
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Set all times on processor meshes equal to reconstructed mesh
|
|
|
|
forAll (databases, procI)
|
|
|
|
{
|
|
|
|
Info<< "Reading database for processor " << procI << endl;
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
databases[procI].setTime(runTime.timeName(), runTime.timeIndex());
|
2010-08-26 14:22:03 +00:00
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Read all meshes and addressing to reconstructed mesh
|
2018-02-28 18:39:10 +00:00
|
|
|
processorMeshesReconstructor procMeshes
|
|
|
|
(
|
|
|
|
databases,
|
2018-03-19 07:44:16 +00:00
|
|
|
regionName
|
2018-02-28 18:39:10 +00:00
|
|
|
);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2018-04-24 18:04:02 +00:00
|
|
|
// Get reconstructed mesh
|
2011-02-19 01:05:24 +00:00
|
|
|
autoPtr<fvMesh> meshPtr = procMeshes.reconstructMesh(runTime);
|
2010-08-26 14:22:03 +00:00
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Mesh write will be controlled by hand
|
|
|
|
meshPtr->write();
|
2011-03-22 12:29:57 +00:00
|
|
|
procMeshes.writeAddressing();
|
2011-02-19 01:05:24 +00:00
|
|
|
meshPtr->setMotionWriteOpt(IOobject::NO_WRITE);
|
|
|
|
meshPtr->setTopoWriteOpt(IOobject::NO_WRITE);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-03-22 12:29:57 +00:00
|
|
|
// Write cell decomposition
|
|
|
|
if (writeCellDist)
|
|
|
|
{
|
|
|
|
// Write as volScalarField for post-processing
|
|
|
|
Info<< "Writing cellDist to time " << runTime.timeName()
|
|
|
|
<< endl;
|
|
|
|
volScalarField cellDist
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"cellDist",
|
|
|
|
runTime.timeName(),
|
|
|
|
meshPtr(),
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
meshPtr(),
|
|
|
|
dimensionedScalar("cellDist", dimless, 0),
|
|
|
|
zeroGradientFvPatchScalarField::typeName
|
|
|
|
);
|
|
|
|
scalarField& cellDistIn = cellDist.internalField();
|
|
|
|
|
|
|
|
label cellI = 0;
|
|
|
|
|
|
|
|
forAll (procMeshes.meshes(), procI)
|
|
|
|
{
|
|
|
|
for
|
|
|
|
(
|
|
|
|
label i = 0;
|
|
|
|
i < procMeshes.meshes()[procI].nCells();
|
|
|
|
i++
|
|
|
|
)
|
|
|
|
{
|
|
|
|
cellDistIn[cellI] = procI;
|
|
|
|
cellI++;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
cellDist.write();
|
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Get region prefix for lagrangian
|
2010-05-12 13:27:55 +00:00
|
|
|
fileName regionPrefix = "";
|
2011-02-19 01:05:24 +00:00
|
|
|
if (regionName != fvMesh::defaultRegion)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
regionPrefix = regionName;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Loop over all times
|
|
|
|
forAll (timeDirs, timeI)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
// Set time for global database
|
|
|
|
runTime.setTime(timeDirs[timeI], timeI);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
Info << "Time = " << runTime.timeName() << endl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Set time for all databases
|
|
|
|
forAll (databases, procI)
|
|
|
|
{
|
|
|
|
databases[procI].setTime(timeDirs[timeI], timeI);
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
polyMesh::readUpdateState procStat = procMeshes.readUpdate();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if (procStat == polyMesh::UNCHANGED)
|
|
|
|
{
|
|
|
|
Info<< "Mesh unchanged" << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
meshPtr->setMotionWriteOpt(IOobject::NO_WRITE);
|
|
|
|
meshPtr->setTopoWriteOpt(IOobject::NO_WRITE);
|
|
|
|
}
|
|
|
|
else if (procStat == polyMesh::POINTS_MOVED)
|
|
|
|
{
|
|
|
|
Info<< "Mesh motion detected. Reconstruct motion points"
|
|
|
|
<< endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Reconstruct the points for moving mesh cases and write them out
|
|
|
|
procMeshes.reconstructPoints(meshPtr());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Set write options
|
|
|
|
meshPtr->setMotionWriteOpt(IOobject::AUTO_WRITE);
|
|
|
|
meshPtr->setTopoWriteOpt(IOobject::NO_WRITE);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Global mesh write
|
|
|
|
meshPtr->write();
|
|
|
|
}
|
|
|
|
else if
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
procStat == polyMesh::TOPO_CHANGE
|
|
|
|
|| procStat == polyMesh::TOPO_PATCH_CHANGE
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Info<< "Topological change detected. Reconstructing mesh"
|
|
|
|
<< endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Reconstruct mesh
|
|
|
|
meshPtr = procMeshes.reconstructMesh(runTime);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Set write options
|
|
|
|
meshPtr->setMotionWriteOpt(IOobject::AUTO_WRITE);
|
|
|
|
meshPtr->setTopoWriteOpt(IOobject::AUTO_WRITE);
|
|
|
|
procMeshes.writeAddressing();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Global mesh write
|
|
|
|
meshPtr->write();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Write out mapping in processor directories
|
|
|
|
forAll (databases, procI)
|
|
|
|
{
|
|
|
|
databases[procI].write();
|
|
|
|
}
|
|
|
|
|
|
|
|
if (writeCellDist)
|
|
|
|
{
|
2011-03-11 17:46:22 +00:00
|
|
|
// Write as volScalarField for post-processing
|
2011-03-22 12:29:57 +00:00
|
|
|
Info<< "Writing cellDist to time " << runTime.timeName()
|
|
|
|
<< endl;
|
2011-02-19 01:05:24 +00:00
|
|
|
volScalarField cellDist
|
|
|
|
(
|
|
|
|
IOobject
|
|
|
|
(
|
|
|
|
"cellDist",
|
|
|
|
runTime.timeName(),
|
|
|
|
meshPtr(),
|
|
|
|
IOobject::NO_READ,
|
|
|
|
IOobject::NO_WRITE
|
|
|
|
),
|
|
|
|
meshPtr(),
|
|
|
|
dimensionedScalar("cellDist", dimless, 0),
|
|
|
|
zeroGradientFvPatchScalarField::typeName
|
|
|
|
);
|
|
|
|
scalarField& cellDistIn = cellDist.internalField();
|
|
|
|
|
|
|
|
label cellI = 0;
|
|
|
|
|
|
|
|
forAll (procMeshes.meshes(), procI)
|
|
|
|
{
|
|
|
|
for
|
|
|
|
(
|
|
|
|
label i = 0;
|
|
|
|
i < procMeshes.meshes()[procI].nCells();
|
|
|
|
i++
|
|
|
|
)
|
|
|
|
{
|
|
|
|
cellDistIn[cellI] = procI;
|
|
|
|
cellI++;
|
|
|
|
}
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
cellDist.write();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
FatalErrorIn(args.executable())
|
2011-02-19 01:05:24 +00:00
|
|
|
<< "Unknown readUpdate state"
|
|
|
|
<< abort(FatalError);
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
fvMesh& mesh = meshPtr();
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Get list of objects from processor0 database
|
|
|
|
IOobjectList objects(procMeshes.meshes()[0], databases[0].timeName());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// If there are any FV fields, reconstruct them
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if
|
2010-12-02 14:33:56 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
objects.lookupClass(volScalarField::typeName).size()
|
|
|
|
|| objects.lookupClass(volVectorField::typeName).size()
|
|
|
|
|| objects.lookupClass(volSphericalTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(volSymmTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(volTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(surfaceScalarField::typeName).size()
|
|
|
|
|| objects.lookupClass(surfaceVectorField::typeName).size()
|
|
|
|
|| objects.lookupClass(surfaceSphericalTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(surfaceSymmTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(surfaceTensorField::typeName).size()
|
|
|
|
)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
Info << "Reconstructing FV fields" << nl << endl;
|
2010-12-02 14:33:56 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
fvFieldReconstructor fvReconstructor
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
procMeshes.faceProcAddressing(),
|
|
|
|
procMeshes.cellProcAddressing(),
|
|
|
|
procMeshes.boundaryProcAddressing()
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
fvReconstructor.reconstructFvVolumeFields<scalar>
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
fvReconstructor.reconstructFvVolumeFields<vector>
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
objects,
|
|
|
|
selectedFields
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
2011-02-19 01:05:24 +00:00
|
|
|
fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
|
2010-08-26 14:22:03 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
fvReconstructor.reconstructFvVolumeFields<symmTensor>
|
|
|
|
(
|
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
fvReconstructor.reconstructFvVolumeFields<tensor>
|
|
|
|
(
|
|
|
|
objects,
|
|
|
|
selectedFields
|
2010-08-26 14:22:03 +00:00
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
fvReconstructor.reconstructFvSurfaceFields<scalar>
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
objects,
|
|
|
|
selectedFields
|
2010-05-12 13:27:55 +00:00
|
|
|
);
|
2011-02-19 01:05:24 +00:00
|
|
|
fvReconstructor.reconstructFvSurfaceFields<vector>
|
|
|
|
(
|
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
|
|
|
|
(
|
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
fvReconstructor.reconstructFvSurfaceFields<symmTensor>
|
|
|
|
(
|
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
fvReconstructor.reconstructFvSurfaceFields<tensor>
|
|
|
|
(
|
|
|
|
objects,
|
|
|
|
selectedFields
|
|
|
|
);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
Info << "No FV fields" << nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// If there are any point fields, reconstruct them
|
|
|
|
if
|
|
|
|
(
|
|
|
|
objects.lookupClass(pointScalarField::typeName).size()
|
|
|
|
|| objects.lookupClass(pointVectorField::typeName).size()
|
|
|
|
|| objects.lookupClass(pointSphericalTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(pointSymmTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(pointTensorField::typeName).size()
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Info << "Reconstructing point fields" << nl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
pointMesh pMesh(mesh);
|
|
|
|
PtrList<pointMesh> pMeshes(procMeshes.meshes().size());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
forAll (pMeshes, procI)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
pMeshes.set(procI, new pointMesh(procMeshes.meshes()[procI]));
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
pointFieldReconstructor pointReconstructor
|
|
|
|
(
|
|
|
|
pMesh,
|
|
|
|
pMeshes,
|
|
|
|
procMeshes.pointProcAddressing(),
|
|
|
|
procMeshes.boundaryProcAddressing()
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
pointReconstructor.reconstructFields<scalar>(objects);
|
|
|
|
pointReconstructor.reconstructFields<vector>(objects);
|
|
|
|
pointReconstructor.reconstructFields<sphericalTensor>(objects);
|
|
|
|
pointReconstructor.reconstructFields<symmTensor>(objects);
|
|
|
|
pointReconstructor.reconstructFields<tensor>(objects);
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
2011-02-19 01:05:24 +00:00
|
|
|
else
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
Info << "No point fields" << nl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// If there are any tetFem fields, reconstruct them
|
|
|
|
if
|
|
|
|
(
|
|
|
|
objects.lookupClass(tetPointScalarField::typeName).size()
|
|
|
|
|| objects.lookupClass(tetPointVectorField::typeName).size()
|
|
|
|
|| objects.lookupClass(tetPointSphericalTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(tetPointSymmTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(tetPointTensorField::typeName).size()
|
|
|
|
|
|
|
|
|| objects.lookupClass(elementScalarField::typeName).size()
|
|
|
|
|| objects.lookupClass(elementVectorField::typeName).size()
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Info << "Reconstructing tet point fields" << nl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
tetPolyMesh tetMesh(mesh);
|
|
|
|
PtrList<tetPolyMesh> tetMeshes(procMeshes.meshes().size());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
forAll (tetMeshes, procI)
|
|
|
|
{
|
|
|
|
tetMeshes.set
|
|
|
|
(
|
|
|
|
procI,
|
|
|
|
new tetPolyMesh(procMeshes.meshes()[procI])
|
|
|
|
);
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
tetPointFieldReconstructor tetPointReconstructor
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
tetMesh,
|
|
|
|
tetMeshes,
|
|
|
|
procMeshes.pointProcAddressing(),
|
|
|
|
procMeshes.faceProcAddressing(),
|
|
|
|
procMeshes.cellProcAddressing(),
|
|
|
|
procMeshes.boundaryProcAddressing()
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// Reconstruct tet point fields
|
|
|
|
tetPointReconstructor.reconstructTetPointFields<scalar>(objects);
|
|
|
|
tetPointReconstructor.reconstructTetPointFields<vector>(objects);
|
|
|
|
tetPointReconstructor.
|
|
|
|
reconstructTetPointFields<sphericalTensor>(objects);
|
|
|
|
tetPointReconstructor.
|
|
|
|
reconstructTetPointFields<symmTensor>(objects);
|
|
|
|
tetPointReconstructor.reconstructTetPointFields<tensor>(objects);
|
|
|
|
|
|
|
|
tetPointReconstructor.reconstructElementFields<scalar>(objects);
|
|
|
|
tetPointReconstructor.reconstructElementFields<vector>(objects);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
Info << "No tetFem fields" << nl << endl;
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// If there are any clouds, reconstruct them.
|
|
|
|
// The problem is that a cloud of size zero will not get written so
|
|
|
|
// in pass 1 we determine the cloud names and per cloud name the
|
|
|
|
// fields. Note that the fields are stored as IOobjectList from
|
|
|
|
// the first processor that has them. They are in pass2 only used
|
|
|
|
// for name and type (scalar, vector etc).
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if (!noLagrangian)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
HashTable<IOobjectList> cloudObjects;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
forAll (databases, procI)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
fileNameList cloudDirs
|
|
|
|
(
|
|
|
|
readDir
|
|
|
|
(
|
|
|
|
databases[procI].timePath()/regionPrefix/cloud::prefix,
|
|
|
|
fileName::DIRECTORY
|
|
|
|
)
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
forAll (cloudDirs, i)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
// Check if we already have cloud objects for
|
|
|
|
// this cloudname
|
|
|
|
HashTable<IOobjectList>::const_iterator iter =
|
|
|
|
cloudObjects.find(cloudDirs[i]);
|
|
|
|
|
|
|
|
if (iter == cloudObjects.end())
|
|
|
|
{
|
|
|
|
// Do local scan for valid cloud objects
|
|
|
|
IOobjectList sprayObjs
|
|
|
|
(
|
|
|
|
procMeshes.meshes()[procI],
|
|
|
|
databases[procI].timeName(),
|
|
|
|
cloud::prefix/cloudDirs[i]
|
|
|
|
);
|
|
|
|
|
|
|
|
IOobject* positionsPtr = sprayObjs.lookup("positions");
|
|
|
|
|
|
|
|
if (positionsPtr)
|
|
|
|
{
|
|
|
|
cloudObjects.insert(cloudDirs[i], sprayObjs);
|
|
|
|
}
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
2011-02-19 01:05:24 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
if (cloudObjects.size())
|
|
|
|
{
|
|
|
|
// Pass2: reconstruct the cloud
|
|
|
|
forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
const word cloudName = string::validate<word>(iter.key());
|
|
|
|
|
|
|
|
// Objects (on arbitrary processor)
|
|
|
|
const IOobjectList& sprayObjs = iter();
|
|
|
|
|
|
|
|
Info<< "Reconstructing lagrangian fields for cloud "
|
|
|
|
<< cloudName << nl << endl;
|
|
|
|
|
|
|
|
reconstructLagrangianPositions
|
|
|
|
(
|
|
|
|
mesh,
|
|
|
|
cloudName,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
procMeshes.faceProcAddressing(),
|
|
|
|
procMeshes.cellProcAddressing()
|
|
|
|
);
|
|
|
|
reconstructLagrangianFields<label>
|
|
|
|
(
|
|
|
|
cloudName,
|
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
sprayObjs
|
|
|
|
);
|
|
|
|
reconstructLagrangianFields<scalar>
|
|
|
|
(
|
|
|
|
cloudName,
|
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
sprayObjs
|
|
|
|
);
|
|
|
|
reconstructLagrangianFields<vector>
|
|
|
|
(
|
|
|
|
cloudName,
|
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
sprayObjs
|
|
|
|
);
|
|
|
|
reconstructLagrangianFields<sphericalTensor>
|
|
|
|
(
|
|
|
|
cloudName,
|
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
sprayObjs
|
|
|
|
);
|
|
|
|
reconstructLagrangianFields<symmTensor>
|
|
|
|
(
|
|
|
|
cloudName,
|
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
sprayObjs
|
|
|
|
);
|
|
|
|
reconstructLagrangianFields<tensor>
|
|
|
|
(
|
|
|
|
cloudName,
|
|
|
|
mesh,
|
|
|
|
procMeshes.meshes(),
|
|
|
|
sprayObjs
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
2011-02-19 01:05:24 +00:00
|
|
|
Info << "No lagrangian fields" << nl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// If there are any FA fields, reconstruct them
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
if
|
|
|
|
(
|
|
|
|
objects.lookupClass(areaScalarField::typeName).size()
|
|
|
|
|| objects.lookupClass(areaVectorField::typeName).size()
|
|
|
|
|| objects.lookupClass(areaSphericalTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(areaSymmTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(areaTensorField::typeName).size()
|
|
|
|
|| objects.lookupClass(edgeScalarField::typeName).size()
|
|
|
|
)
|
|
|
|
{
|
|
|
|
Info << "Reconstructing FA fields" << nl << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
faMesh aMesh(mesh);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2018-02-28 18:39:10 +00:00
|
|
|
processorFaMeshes procFaMeshes(procMeshes.meshes(), true);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
faFieldReconstructor faReconstructor
|
2010-05-12 13:27:55 +00:00
|
|
|
(
|
2011-02-19 01:05:24 +00:00
|
|
|
aMesh,
|
|
|
|
procFaMeshes.meshes(),
|
|
|
|
procFaMeshes.edgeProcAddressing(),
|
|
|
|
procFaMeshes.faceProcAddressing(),
|
|
|
|
procFaMeshes.boundaryProcAddressing()
|
|
|
|
);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
faReconstructor.reconstructFaAreaFields<scalar>(objects);
|
|
|
|
faReconstructor.reconstructFaAreaFields<vector>(objects);
|
|
|
|
faReconstructor
|
|
|
|
.reconstructFaAreaFields<sphericalTensor>(objects);
|
|
|
|
faReconstructor.reconstructFaAreaFields<symmTensor>(objects);
|
|
|
|
faReconstructor.reconstructFaAreaFields<tensor>(objects);
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
faReconstructor.reconstructFaEdgeFields<scalar>(objects);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
Info << "No FA fields" << nl << endl;
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
// If there are any "uniform" directories copy them from
|
|
|
|
// the master processor
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-19 01:05:24 +00:00
|
|
|
fileName uniformDir0 = databases[0].timePath()/"uniform";
|
|
|
|
if (isDir(uniformDir0))
|
|
|
|
{
|
|
|
|
cp(uniformDir0, runTime.timePath());
|
|
|
|
}
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
2010-08-26 14:22:03 +00:00
|
|
|
Info<< "End.\n" << endl;
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
return 0;
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|