mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: singleCellMesh: do all times, write to different region
This commit is contained in:
@ -25,25 +25,32 @@ Application
|
|||||||
singleCellMesh
|
singleCellMesh
|
||||||
|
|
||||||
Description
|
Description
|
||||||
Removes all but one cells of the mesh. Used to generate mesh and fields
|
Reads all fields and maps them to a mesh with all internal faces removed
|
||||||
|
(singleCellFvMesh) which gets written to region "singleCell".
|
||||||
|
|
||||||
|
Used to generate mesh and fields
|
||||||
that can be used for boundary-only data.
|
that can be used for boundary-only data.
|
||||||
Might easily result in illegal mesh though so only look at boundaries
|
Might easily result in illegal mesh though so only look at boundaries
|
||||||
in paraview.
|
in paraview.
|
||||||
|
|
||||||
\*---------------------------------------------------------------------------*/
|
\*---------------------------------------------------------------------------*/
|
||||||
|
|
||||||
|
|
||||||
#include "argList.H"
|
#include "argList.H"
|
||||||
#include "fvMesh.H"
|
#include "fvMesh.H"
|
||||||
#include "volFields.H"
|
#include "volFields.H"
|
||||||
#include "Time.H"
|
#include "Time.H"
|
||||||
#include "ReadFields.H"
|
#include "ReadFields.H"
|
||||||
#include "singleCellFvMesh.H"
|
#include "singleCellFvMesh.H"
|
||||||
|
#include "timeSelector.H"
|
||||||
|
|
||||||
using namespace Foam;
|
using namespace Foam;
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
|
// Name of region to create
|
||||||
|
const string singleCellName = "singleCell";
|
||||||
|
|
||||||
|
|
||||||
template<class GeoField>
|
template<class GeoField>
|
||||||
void interpolateFields
|
void interpolateFields
|
||||||
(
|
(
|
||||||
@ -65,19 +72,82 @@ void interpolateFields
|
|||||||
|
|
||||||
int main(int argc, char *argv[])
|
int main(int argc, char *argv[])
|
||||||
{
|
{
|
||||||
# include "addOverwriteOption.H"
|
// constant, not false
|
||||||
# include "addTimeOptions.H"
|
timeSelector::addOptions(true, false);
|
||||||
|
|
||||||
# include "setRootCase.H"
|
# include "setRootCase.H"
|
||||||
# include "createTime.H"
|
# include "createTime.H"
|
||||||
// Get times list
|
|
||||||
instantList Times = runTime.times();
|
|
||||||
# include "checkTimeOptions.H"
|
|
||||||
runTime.setTime(Times[startTime], startTime);
|
|
||||||
# include "createMesh.H"
|
|
||||||
const word oldInstance = mesh.pointsInstance();
|
|
||||||
|
|
||||||
const bool overwrite = args.optionFound("overwrite");
|
instantList timeDirs = timeSelector::select0(runTime, args);
|
||||||
|
|
||||||
|
# include "createNamedMesh.H"
|
||||||
|
|
||||||
|
if (regionName == singleCellName)
|
||||||
|
{
|
||||||
|
FatalErrorIn(args.executable())
|
||||||
|
<< "Cannot convert region " << singleCellName
|
||||||
|
<< " since result would overwrite it. Please rename your region."
|
||||||
|
<< exit(FatalError);
|
||||||
|
}
|
||||||
|
|
||||||
|
// Create the mesh
|
||||||
|
Info<< "Creating singleCell mesh" << nl << endl;
|
||||||
|
autoPtr<singleCellFvMesh> scMesh
|
||||||
|
(
|
||||||
|
new singleCellFvMesh
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
singleCellName,
|
||||||
|
mesh.polyMesh::instance(),
|
||||||
|
runTime,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::AUTO_WRITE
|
||||||
|
),
|
||||||
|
mesh
|
||||||
|
)
|
||||||
|
);
|
||||||
|
// For convenience create any fv* files
|
||||||
|
if (!exists(scMesh().fvSolution::objectPath()))
|
||||||
|
{
|
||||||
|
mkDir(scMesh().fvSolution::path());
|
||||||
|
ln("../fvSolution", scMesh().fvSolution::objectPath());
|
||||||
|
}
|
||||||
|
if (!exists(scMesh().fvSchemes::objectPath()))
|
||||||
|
{
|
||||||
|
mkDir(scMesh().fvSolution::path());
|
||||||
|
ln("../fvSchemes", scMesh().fvSchemes::objectPath());
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
forAll(timeDirs, timeI)
|
||||||
|
{
|
||||||
|
runTime.setTime(timeDirs[timeI], timeI);
|
||||||
|
|
||||||
|
Info<< nl << "Time = " << runTime.timeName() << endl;
|
||||||
|
|
||||||
|
|
||||||
|
// Check for new mesh
|
||||||
|
if (mesh.readUpdate() != polyMesh::UNCHANGED)
|
||||||
|
{
|
||||||
|
Info<< "Detected changed mesh. Recreating singleCell mesh." << endl;
|
||||||
|
scMesh.clear(); // remove any registered objects
|
||||||
|
scMesh.reset
|
||||||
|
(
|
||||||
|
new singleCellFvMesh
|
||||||
|
(
|
||||||
|
IOobject
|
||||||
|
(
|
||||||
|
singleCellName,
|
||||||
|
mesh.polyMesh::instance(),
|
||||||
|
runTime,
|
||||||
|
IOobject::NO_READ,
|
||||||
|
IOobject::AUTO_WRITE
|
||||||
|
),
|
||||||
|
mesh
|
||||||
|
)
|
||||||
|
);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
// Read objects in time directory
|
// Read objects in time directory
|
||||||
@ -99,38 +169,18 @@ int main(int argc, char *argv[])
|
|||||||
PtrList<volTensorField> vtFlds;
|
PtrList<volTensorField> vtFlds;
|
||||||
ReadFields(mesh, objects, vtFlds);
|
ReadFields(mesh, objects, vtFlds);
|
||||||
|
|
||||||
|
|
||||||
if (!overwrite)
|
|
||||||
{
|
|
||||||
runTime++;
|
|
||||||
}
|
|
||||||
|
|
||||||
// Create the mesh
|
|
||||||
singleCellFvMesh scMesh
|
|
||||||
(
|
|
||||||
IOobject
|
|
||||||
(
|
|
||||||
mesh.name(),
|
|
||||||
mesh.polyMesh::instance(),
|
|
||||||
runTime,
|
|
||||||
IOobject::NO_READ,
|
|
||||||
IOobject::AUTO_WRITE
|
|
||||||
),
|
|
||||||
mesh
|
|
||||||
);
|
|
||||||
|
|
||||||
|
|
||||||
// Map and store the fields on the scMesh.
|
// Map and store the fields on the scMesh.
|
||||||
interpolateFields(scMesh, vsFlds);
|
interpolateFields(scMesh(), vsFlds);
|
||||||
interpolateFields(scMesh, vvFlds);
|
interpolateFields(scMesh(), vvFlds);
|
||||||
interpolateFields(scMesh, vstFlds);
|
interpolateFields(scMesh(), vstFlds);
|
||||||
interpolateFields(scMesh, vsymtFlds);
|
interpolateFields(scMesh(), vsymtFlds);
|
||||||
interpolateFields(scMesh, vtFlds);
|
interpolateFields(scMesh(), vtFlds);
|
||||||
|
|
||||||
|
|
||||||
// Write
|
// Write
|
||||||
Info<< "Writing mesh to time " << runTime.timeName() << endl;
|
Info<< "Writing mesh to time " << runTime.timeName() << endl;
|
||||||
scMesh.write();
|
scMesh().write();
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
Info<< "End\n" << endl;
|
Info<< "End\n" << endl;
|
||||||
|
|||||||
Reference in New Issue
Block a user