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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fieldToCell.H"
30 #include "polyMesh.H"
31 #include "cellSet.H"
32 #include "Time.H"
33 #include "IFstream.H"
34 #include "fieldDictionary.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41  defineTypeNameAndDebug(fieldToCell, 0);
42  addToRunTimeSelectionTable(topoSetSource, fieldToCell, word);
43  addToRunTimeSelectionTable(topoSetSource, fieldToCell, istream);
44  addToRunTimeSelectionTable(topoSetCellSource, fieldToCell, word);
45  addToRunTimeSelectionTable(topoSetCellSource, fieldToCell, istream);
47  (
48  topoSetCellSource,
49  fieldToCell,
50  word,
51  field
52  );
54  (
55  topoSetCellSource,
56  fieldToCell,
57  istream,
58  field
59  );
60 }
61 
62 
63 Foam::topoSetSource::addToUsageTable Foam::fieldToCell::usage_
64 (
65  fieldToCell::typeName,
66  "\n Usage: fieldToCell field min max\n\n"
67  " Select all cells with field value >= min and <= max\n\n"
68 );
69 
70 
71 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
72 
73 void Foam::fieldToCell::applyToSet
74 (
75  const topoSetSource::setAction action,
76  const scalarField& field,
77  topoSet& set
78 ) const
79 {
80  if (verbose_)
81  {
82  Info << " Field min:" << min(field) << " max:" << max(field) << nl;
83  }
84 
85  if (action == topoSetSource::ADD || action == topoSetSource::NEW)
86  {
87  if (verbose_)
88  {
89  Info<< " Adding all cells with value of field " << fieldName_
90  << " within range " << min_ << ".." << max_ << endl;
91  }
92 
93  forAll(field, celli)
94  {
95  if (field[celli] >= min_ && field[celli] <= max_)
96  {
97  set.set(celli);
98  }
99  }
100  }
101  else if (action == topoSetSource::SUBTRACT)
102  {
103  if (verbose_)
104  {
105  Info<< " Removing all cells with value of field " << fieldName_
106  << " within range " << min_ << ".." << max_ << endl;
107  }
108 
109  forAll(field, celli)
110  {
111  if (field[celli] >= min_ && field[celli] <= max_)
112  {
113  set.unset(celli);
114  }
115  }
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const polyMesh& mesh,
125  const word& fieldName,
126  const scalar min,
127  const scalar max
128 )
129 :
130  topoSetCellSource(mesh),
131  fieldName_(fieldName),
132  min_(min),
133  max_(max)
134 {
135  if (min_ > max_)
136  {
138  << "Input min value = " << min_ << " is larger than "
139  << "input max value = " << max_ << " for field = " << fieldName_
140  << endl;
141  }
142 }
143 
144 
146 (
147  const polyMesh& mesh,
148  const dictionary& dict
149 )
150 :
151  fieldToCell
152  (
153  mesh,
154  dict.get<word>("field"),
155  dict.get<scalar>("min"),
156  dict.get<scalar>("max")
157  )
158 {}
159 
160 
162 (
163  const polyMesh& mesh,
164  Istream& is
165 )
166 :
168  fieldName_(checkIs(is)),
169  min_(readScalar(checkIs(is))),
170  max_(readScalar(checkIs(is)))
171 {}
172 
173 
174 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175 
176 void Foam::fieldToCell::applyToSet
177 (
178  const topoSetSource::setAction action,
179  topoSet& set
180 ) const
181 {
182  // Try to load field
183  IOobject fieldObject
184  (
185  fieldName_,
186  mesh().time().timeName(),
187  mesh(),
191  );
192 
193  // Note: should check for volScalarField but that introduces dependency
194  // on volMesh so just use another type with processor-local scope
195  if (!fieldObject.typeHeaderOk<labelIOList>(false))
196  {
198  << "Cannot read field " << fieldName_
199  << " from time " << mesh().time().timeName() << endl;
200  }
201  else if ("volScalarField" == fieldObject.headerClassName())
202  {
203  // Note: cannot use volScalarField::typeName since that would
204  // introduce linkage problems (finiteVolume needs meshTools)
205 
206  IFstream str(fieldObject.typeFilePath<labelIOList>());
207 
208  // Read as dictionary
209  fieldDictionary fieldDict(fieldObject, fieldObject.headerClassName());
210 
211  scalarField internalVals("internalField", fieldDict, mesh().nCells());
212 
213  applyToSet(action, internalVals, set);
214  }
215  else if ("volVectorField" == fieldObject.headerClassName())
216  {
217  // Note: cannot use volVectorField::typeName since that would
218  // introduce linkage problems (finiteVolume needs meshTools)
219 
220  IFstream str(fieldObject.typeFilePath<labelIOList>());
221 
222  // Read as dictionary
223  fieldDictionary fieldDict(fieldObject, fieldObject.headerClassName());
224 
225  vectorField internalVals("internalField", fieldDict, mesh().nCells());
226 
227  applyToSet(action, mag(internalVals), set);
228  }
229  else
230  {
232  << "Cannot handle fields of type " << fieldObject.headerClassName()
233  << endl;
234  }
235 }
236 
237 
238 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
rDeltaTY field()
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
bool verbose_
Output verbosity (default: true)
Create a new set and ADD elements to it.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Add elements to current set.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
The topoSetCellSource is a intermediate class for handling topoSet sources for selecting cells...
Macros for easy insertion into run-time selection tables.
addNamedToRunTimeSelectionTable(topoSetCellSource, badQualityToCell, word, badQuality)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
word timeName
Definition: getTimeIndex.H:3
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
setAction
Enumeration defining various actions.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
defineTypeNameAndDebug(combustionModel, 0)
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:59
Subtract elements from current set.
#define WarningInFunction
Report a warning using Foam::Warning.
Class with constructor to add usage string to table.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
fieldToCell(const polyMesh &mesh, const word &fieldName, const scalar min, const scalar max)
Construct from components.
Definition: fieldToCell.C:116
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Do not request registration (bool: false)
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)