mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
125 lines
3.4 KiB
C
125 lines
3.4 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd.
|
|
\\/ 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 "sampledSurfaces.H"
|
|
#include "volFields.H"
|
|
#include "IOobjectList.H"
|
|
#include "stringListOps.H"
|
|
|
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
|
|
|
void Foam::sampledSurfaces::clearFieldGroups()
|
|
{
|
|
scalarFields_.clear();
|
|
vectorFields_.clear();
|
|
sphericalTensorFields_.clear();
|
|
symmTensorFields_.clear();
|
|
tensorFields_.clear();
|
|
}
|
|
|
|
|
|
Foam::label Foam::sampledSurfaces::appendFieldGroup
|
|
(
|
|
const word& fieldName,
|
|
const word& fieldType
|
|
)
|
|
{
|
|
if (fieldType == volScalarField::typeName)
|
|
{
|
|
scalarFields_.append(fieldName);
|
|
return 1;
|
|
}
|
|
else if (fieldType == volVectorField::typeName)
|
|
{
|
|
vectorFields_.append(fieldName);
|
|
return 1;
|
|
}
|
|
else if (fieldType == volSphericalTensorField::typeName)
|
|
{
|
|
sphericalTensorFields_.append(fieldName);
|
|
return 1;
|
|
}
|
|
else if (fieldType == volSymmTensorField::typeName)
|
|
{
|
|
symmTensorFields_.append(fieldName);
|
|
return 1;
|
|
}
|
|
else if (fieldType == volTensorField::typeName)
|
|
{
|
|
tensorFields_.append(fieldName);
|
|
return 1;
|
|
}
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
Foam::label Foam::sampledSurfaces::classifyFields()
|
|
{
|
|
label nFields = 0;
|
|
clearFieldGroups();
|
|
|
|
// check files for a particular time
|
|
if (loadFromFiles_)
|
|
{
|
|
IOobjectList objects(mesh_, mesh_.time().timeName());
|
|
wordList allFields = objects.sortedNames();
|
|
|
|
labelList indices = findStrings(fieldSelection_, allFields);
|
|
|
|
forAll(indices, fieldI)
|
|
{
|
|
const word& fieldName = allFields[indices[fieldI]];
|
|
|
|
nFields += appendFieldGroup
|
|
(
|
|
fieldName,
|
|
objects.find(fieldName)()->headerClassName()
|
|
);
|
|
}
|
|
}
|
|
else
|
|
{
|
|
wordList allFields = mesh_.sortedNames();
|
|
labelList indices = findStrings(fieldSelection_, allFields);
|
|
|
|
forAll(indices, fieldI)
|
|
{
|
|
const word& fieldName = allFields[indices[fieldI]];
|
|
|
|
nFields += appendFieldGroup
|
|
(
|
|
fieldName,
|
|
mesh_.find(fieldName)()->type()
|
|
);
|
|
}
|
|
}
|
|
|
|
return nFields;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|