fieldToCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fieldToCell.H"
27 #include "polyMesh.H"
28 #include "cellSet.H"
29 #include "Time.H"
30 #include "IFstream.H"
31 #include "fieldDictionary.H"
32 
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 defineTypeNameAndDebug(fieldToCell, 0);
41 
42 addToRunTimeSelectionTable(topoSetSource, fieldToCell, word);
43 
44 addToRunTimeSelectionTable(topoSetSource, fieldToCell, istream);
45 
46 }
47 
48 
50 (
51  fieldToCell::typeName,
52  "\n Usage: fieldToCell field min max\n\n"
53  " Select all cells with field value >= min and <= max\n\n"
54 );
55 
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
60 (
61  const topoSetSource::setAction action,
62  const scalarField& field,
63  topoSet& set
64 ) const
65 {
66  Info<< " Field min:" << min(field)
67  << " max:" << max(field) << endl;
68 
69  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
70  {
71  Info<< " Adding all cells with value of field " << fieldName_
72  << " within range " << min_ << ".." << max_ << endl;
73 
74  forAll(field, cellI)
75  {
76  if (field[cellI] >= min_ && field[cellI] <= max_)
77  {
78  set.insert(cellI);
79  }
80  }
81  }
82  else if (action == topoSetSource::DELETE)
83  {
84  Info<< " Removing all cells with value of field " << fieldName_
85  << " within range " << min_ << ".." << max_ << endl;
86 
87  forAll(field, cellI)
88  {
89  if (field[cellI] >= min_ && field[cellI] <= max_)
90  {
91  set.erase(cellI);
92  }
93  }
94  }
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
100 // Construct from components
102 (
103  const polyMesh& mesh,
104  const word& fieldName,
105  const scalar min,
106  const scalar max
107 )
108 :
110  fieldName_(fieldName),
111  min_(min),
112  max_(max)
113 {}
114 
115 
116 // Construct from dictionary
118 (
119  const polyMesh& mesh,
120  const dictionary& dict
121 )
122 :
124  fieldName_(dict.lookup("fieldName")),
125  min_(readScalar(dict.lookup("min"))),
126  max_(readScalar(dict.lookup("max")))
127 {}
128 
129 
130 // Construct from Istream
132 (
133  const polyMesh& mesh,
134  Istream& is
135 )
136 :
138  fieldName_(checkIs(is)),
139  min_(readScalar(checkIs(is))),
140  max_(readScalar(checkIs(is)))
141 {}
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 (
154  const topoSetSource::setAction action,
155  topoSet& set
156 ) const
157 {
158 
159 // // Construct temporary fvMesh from polyMesh
160 // fvMesh fMesh
161 // (
162 // mesh(), // IOobject
163 // mesh().points(),
164 // mesh().faces(),
165 // mesh().cells()
166 // );
167 //
168 // const polyBoundaryMesh& patches = mesh().boundaryMesh();
169 //
170 // List<polyPatch*> newPatches(patches.size());
171 // forAll(patches, patchI)
172 // {
173 // const polyPatch& pp = patches[patchI];
174 //
175 // newPatches[patchI] =
176 // patches[patchI].clone
177 // (
178 // fMesh.boundaryMesh(),
179 // patchI,
180 // pp.size(),
181 // pp.start()
182 // ).ptr();
183 // }
184 // fMesh.addFvPatches(newPatches);
185 
186  // Try to load field
188  (
189  fieldName_,
190  mesh().time().timeName(),
191  mesh(),
194  false
195  );
196 
197  if (!fieldObject.headerOk())
198  {
200  << "Cannot read field " << fieldName_
201  << " from time " << mesh().time().timeName() << endl;
202  }
203  else if (fieldObject.headerClassName() == "volScalarField")
204  {
206 
207  // Read dictionary
209 
210  scalarField internalVals("internalField", fieldDict, mesh().nCells());
211 
212  applyToSet(action, internalVals, set);
213  }
214  else if (fieldObject.headerClassName() == "volVectorField")
215  {
217 
218  // Read dictionary
220 
221  vectorField internalVals("internalField", fieldDict, mesh().nCells());
222 
223  applyToSet(action, mag(internalVals), set);
224  }
225  else
226  {
228  << "Cannot handle fields of type " << fieldObject.headerClassName()
229  << endl;
230  }
231 }
232 
233 
234 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::filePath
fileName filePath() const
Return complete path + object name if the file exists.
Definition: IOobject.C:336
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
fieldToCell.H
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::fieldDictionary
Read field as dictionary (without mesh).
Definition: fieldDictionary.H:47
Foam::IFstream
Input from file stream.
Definition: IFstream.H:81
Foam::fieldToCell::fieldToCell
fieldToCell(const polyMesh &mesh, const word &fieldName, const scalar min, const scalar max)
Construct from components.
Definition: fieldToCell.C:102
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::fieldToCell::usage_
static addToUsageTable usage_
Add usage string.
Definition: fieldToCell.H:57
Foam::topoSetSource::addToUsageTable
Class with constructor to add usage string to table.
Definition: topoSetSource.H:100
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::topoSetSource::setAction
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
Foam::fieldToCell::applyToSet
void applyToSet(const topoSetSource::setAction action, const scalarField &field, topoSet &set) const
Depending on field values add to or delete from cellSet.
Definition: fieldToCell.C:60
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::IOobject::headerClassName
const word & headerClassName() const
Return name of the class name read from header.
Definition: IOobject.H:279
Foam::Info
messageStream Info
Foam::HashTable::erase
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Foam::topoSet
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
IFstream.H
Foam::fieldToCell::~fieldToCell
virtual ~fieldToCell()
Destructor.
Definition: fieldToCell.C:146
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::topoSetSource
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
fieldObject
IOobject fieldObject(fieldNames[var2field[nVar]], runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE)
cellSet.H
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
fieldDictionary.H