/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2013-2016 OpenFOAM Foundation ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM 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. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see . Application Test-fieldMapping Description Test app for mapping of fields. \*---------------------------------------------------------------------------*/ #include "argList.H" #include "timeSelector.H" #include "fvMesh.H" #include "volFields.H" #include "Time.H" #include "OFstream.H" #include "meshTools.H" #include "removeFaces.H" #include "mapPolyMesh.H" #include "polyTopoChange.H" #include "fvCFD.H" #include "Random.H" using namespace Foam; // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // Main program: int main(int argc, char *argv[]) { timeSelector::addOptions_singleTime(); // Single-time options argList::addBoolOption ( "inflate", "Use inflation/deflation for deleting cells" ); #include "setRootCase.H" #include "createTime.H" // Allow override of time from specified time options, or no-op timeSelector::setTimeIfPresent(runTime, args); #include "createMesh.H" const bool inflate = args.found("inflate"); if (inflate) { Info<< "Deleting cells using inflation/deflation" << nl << endl; } else { Info<< "Deleting cells, introducing points at new position" << nl << endl; } Random rndGen(0); // Test mapping // ------------ // Mapping is volume averaged // 1. uniform field stays uniform volScalarField one ( IOobject ( "one", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE ), mesh, dimensionedScalar("one", dimless, 1.0), fvPatchFieldBase::zeroGradientType() ); Info<< "Writing one field " << one.name() << " in " << runTime.timeName() << endl; one.write(); // 2. linear profile gets preserved volScalarField ccX ( IOobject ( "ccX", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE ), mesh.C().component(0) ); Info<< "Writing x component of cell centres to " << ccX.name() << " in " << runTime.timeName() << endl; ccX.write(); // Uniform surface field surfaceScalarField surfaceOne ( IOobject ( "surfaceOne", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE ), mesh, dimensionedScalar("one", dimless, 1.0) ); Info<< "Writing surface one field " << surfaceOne.name() << " in " << runTime.timeName() << endl; surfaceOne.write(); // Force allocation of V. Important for any mesh changes since otherwise // old time volumes are not stored const scalar totalVol = gSum(mesh.V()); // Face removal engine. No checking for not merging boundary faces. removeFaces faceRemover(mesh, GREAT); // Comparison for inequality const auto isNotEqual = notEqualOp(1e-10); while (runTime.loop()) { Info<< "Time = " << runTime.timeName() << nl << endl; if (!mesh.nInternalFaces()) { break; } // Remove face label candidateFacei = rndGen.position