Files
OpenFOAM-12/applications/utilities/preProcessing/mapFields/mapMeshes.C
Henry Weller ed7e703040 Time::timeName(): no longer needed, calls replaced by name()
The timeName() function simply returns the dimensionedScalar::name() which holds
the user-time name of the current time and now that timeName() is no longer
virtual the dimensionedScalar::name() can be called directly.  The timeName()
function implementation is maintained for backward-compatibility.
2022-11-30 15:53:51 +00:00

278 lines
7.2 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "mapMeshes.H"
#include "surfaceMesh.H"
#include "processorFvPatch.H"
#include "mapLagrangian.H"
#include "MapVolFields.H"
#include "MapConsistentVolFields.H"
#include "UnMapped.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
void Foam::mapConsistentMesh
(
const fvMesh& meshSource,
const fvMesh& meshTarget,
const meshToMesh0::order& mapOrder
)
{
// Create the interpolation scheme
meshToMesh0 meshToMesh0Interp(meshSource, meshTarget);
Info<< nl
<< "Consistently creating and mapping fields for time "
<< meshSource.time().name() << nl << endl;
{
// Search for list of objects for this time
IOobjectList objects(meshSource, meshSource.time().name());
// Map volFields
// ~~~~~~~~~~~~~
MapConsistentVolFields<scalar>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapConsistentVolFields<vector>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapConsistentVolFields<sphericalTensor>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapConsistentVolFields<symmTensor>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapConsistentVolFields<tensor>
(
objects,
meshToMesh0Interp,
mapOrder
);
}
{
// Search for list of target objects for this time
IOobjectList objects(meshTarget, meshTarget.time().name());
// Mark surfaceFields as unmapped
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
UnMapped<surfaceScalarField>(objects);
UnMapped<surfaceVectorField>(objects);
UnMapped<surfaceSphericalTensorField>(objects);
UnMapped<surfaceSymmTensorField>(objects);
UnMapped<surfaceTensorField>(objects);
// Mark pointFields as unmapped
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
UnMapped<pointScalarField>(objects);
UnMapped<pointVectorField>(objects);
UnMapped<pointSphericalTensorField>(objects);
UnMapped<pointSymmTensorField>(objects);
UnMapped<pointTensorField>(objects);
}
mapLagrangian(meshToMesh0Interp);
}
void Foam::mapSubMesh
(
const fvMesh& meshSource,
const fvMesh& meshTarget,
const HashTable<word>& patchMap,
const wordList& cuttingPatches,
const meshToMesh0::order& mapOrder
)
{
// Create the interpolation scheme
meshToMesh0 meshToMesh0Interp
(
meshSource,
meshTarget,
patchMap,
cuttingPatches
);
Info<< nl
<< "Mapping fields for time " << meshSource.time().name()
<< nl << endl;
{
// Search for list of source objects for this time
IOobjectList objects(meshSource, meshSource.time().name());
// Map volFields
// ~~~~~~~~~~~~~
MapVolFields<scalar>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapVolFields<vector>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapVolFields<sphericalTensor>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapVolFields<symmTensor>
(
objects,
meshToMesh0Interp,
mapOrder
);
MapVolFields<tensor>
(
objects,
meshToMesh0Interp,
mapOrder
);
}
{
// Search for list of target objects for this time
IOobjectList objects(meshTarget, meshTarget.time().name());
// Mark surfaceFields as unmapped
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
UnMapped<surfaceScalarField>(objects);
UnMapped<surfaceVectorField>(objects);
UnMapped<surfaceSphericalTensorField>(objects);
UnMapped<surfaceSymmTensorField>(objects);
UnMapped<surfaceTensorField>(objects);
// Mark pointFields as unmapped
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
UnMapped<pointScalarField>(objects);
UnMapped<pointVectorField>(objects);
UnMapped<pointSphericalTensorField>(objects);
UnMapped<pointSymmTensorField>(objects);
UnMapped<pointTensorField>(objects);
}
mapLagrangian(meshToMesh0Interp);
}
void Foam::mapConsistentSubMesh
(
const fvMesh& meshSource,
const fvMesh& meshTarget,
const meshToMesh0::order& mapOrder
)
{
HashTable<word> patchMap;
HashTable<label> cuttingPatchTable;
forAll(meshTarget.boundary(), patchi)
{
if (!isA<processorFvPatch>(meshTarget.boundary()[patchi]))
{
patchMap.insert
(
meshTarget.boundary()[patchi].name(),
meshTarget.boundary()[patchi].name()
);
}
else
{
cuttingPatchTable.insert
(
meshTarget.boundaryMesh()[patchi].name(),
-1
);
}
}
mapSubMesh
(
meshSource,
meshTarget,
patchMap,
cuttingPatchTable.toc(),
mapOrder
);
}
Foam::wordList Foam::addProcessorPatches
(
const fvMesh& meshTarget,
const wordList& cuttingPatches
)
{
// Add the processor patches to the cutting list
HashTable<label> cuttingPatchTable;
forAll(cuttingPatches, i)
{
cuttingPatchTable.insert(cuttingPatches[i], i);
}
forAll(meshTarget.boundary(), patchi)
{
if (isA<processorFvPatch>(meshTarget.boundary()[patchi]))
{
if
(
!cuttingPatchTable.found
(
meshTarget.boundaryMesh()[patchi].name()
)
)
{
cuttingPatchTable.insert
(
meshTarget.boundaryMesh()[patchi].name(),
-1
);
}
}
}
return cuttingPatchTable.toc();
}
// ************************************************************************* //