Skip to content
Snippets Groups Projects
Commit cc825319 authored by mattijs's avatar mattijs
Browse files

ENH: singleCellMesh: do all times, write to different region

parent b723a6aa
Branches
Tags
No related merge requests found
...@@ -25,25 +25,32 @@ Application ...@@ -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,72 +72,115 @@ void interpolateFields ...@@ -65,72 +72,115 @@ 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");
// Read objects in time directory
IOobjectList objects(mesh, runTime.timeName());
// Read vol fields.
PtrList<volScalarField> vsFlds;
ReadFields(mesh, objects, vsFlds);
PtrList<volVectorField> vvFlds;
ReadFields(mesh, objects, vvFlds);
PtrList<volSphericalTensorField> vstFlds;
ReadFields(mesh, objects, vstFlds);
PtrList<volSymmTensorField> vsymtFlds; instantList timeDirs = timeSelector::select0(runTime, args);
ReadFields(mesh, objects, vsymtFlds);
PtrList<volTensorField> vtFlds; # include "createNamedMesh.H"
ReadFields(mesh, objects, vtFlds);
if (regionName == singleCellName)
if (!overwrite)
{ {
runTime++; FatalErrorIn(args.executable())
<< "Cannot convert region " << singleCellName
<< " since result would overwrite it. Please rename your region."
<< exit(FatalError);
} }
// Create the mesh // Create the mesh
singleCellFvMesh scMesh Info<< "Creating singleCell mesh" << nl << endl;
autoPtr<singleCellFvMesh> scMesh
( (
IOobject new singleCellFvMesh
( (
mesh.name(), IOobject
mesh.polyMesh::instance(), (
runTime, singleCellName,
IOobject::NO_READ, mesh.polyMesh::instance(),
IOobject::AUTO_WRITE runTime,
), IOobject::NO_READ,
mesh 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());
}
// Map and store the fields on the scMesh. forAll(timeDirs, timeI)
interpolateFields(scMesh, vsFlds); {
interpolateFields(scMesh, vvFlds); runTime.setTime(timeDirs[timeI], timeI);
interpolateFields(scMesh, vstFlds);
interpolateFields(scMesh, vsymtFlds); Info<< nl << "Time = " << runTime.timeName() << endl;
interpolateFields(scMesh, vtFlds);
// Check for new mesh
// Write if (mesh.readUpdate() != polyMesh::UNCHANGED)
Info<< "Writing mesh to time " << runTime.timeName() << endl; {
scMesh.write(); 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
IOobjectList objects(mesh, runTime.timeName());
// Read vol fields.
PtrList<volScalarField> vsFlds;
ReadFields(mesh, objects, vsFlds);
PtrList<volVectorField> vvFlds;
ReadFields(mesh, objects, vvFlds);
PtrList<volSphericalTensorField> vstFlds;
ReadFields(mesh, objects, vstFlds);
PtrList<volSymmTensorField> vsymtFlds;
ReadFields(mesh, objects, vsymtFlds);
PtrList<volTensorField> vtFlds;
ReadFields(mesh, objects, vtFlds);
// Map and store the fields on the scMesh.
interpolateFields(scMesh(), vsFlds);
interpolateFields(scMesh(), vvFlds);
interpolateFields(scMesh(), vstFlds);
interpolateFields(scMesh(), vsymtFlds);
interpolateFields(scMesh(), vtFlds);
// Write
Info<< "Writing mesh to time " << runTime.timeName() << endl;
scMesh().write();
}
Info<< "End\n" << endl; Info<< "End\n" << endl;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment