All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
fieldToCell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 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  defineTypeNameAndDebug(fieldToCell, 0);
40  addToRunTimeSelectionTable(topoSetSource, fieldToCell, word);
41  addToRunTimeSelectionTable(topoSetSource, fieldToCell, istream);
42 }
43 
44 
45 Foam::topoSetSource::addToUsageTable Foam::fieldToCell::usage_
46 (
47  fieldToCell::typeName,
48  "\n Usage: fieldToCell field min max\n\n"
49  " Select all cells with field value >= min and <= max\n\n"
50 );
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::fieldToCell::applyToSet
56 (
57  const topoSetSource::setAction action,
58  const scalarField& field,
59  topoSet& set
60 ) const
61 {
62  Info<< " Field min:" << min(field)
63  << " max:" << max(field) << endl;
64 
65  if ((action == topoSetSource::NEW) || (action == topoSetSource::ADD))
66  {
67  Info<< " Adding all cells with value of field " << fieldName_
68  << " within range " << min_ << ".." << max_ << endl;
69 
70  forAll(field, celli)
71  {
72  if (field[celli] >= min_ && field[celli] <= max_)
73  {
74  set.insert(celli);
75  }
76  }
77  }
78  else if (action == topoSetSource::DELETE)
79  {
80  Info<< " Removing all cells with value of field " << fieldName_
81  << " within range " << min_ << ".." << max_ << endl;
82 
83  forAll(field, celli)
84  {
85  if (field[celli] >= min_ && field[celli] <= max_)
86  {
87  set.erase(celli);
88  }
89  }
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const polyMesh& mesh,
99  const word& fieldName,
100  const scalar min,
101  const scalar max
102 )
103 :
104  topoSetSource(mesh),
105  fieldName_(fieldName),
106  min_(min),
107  max_(max)
108 {}
109 
110 
112 (
113  const polyMesh& mesh,
114  const dictionary& dict
115 )
116 :
117  topoSetSource(mesh),
118  fieldName_(dict.lookup("field")),
119  min_(dict.lookup<scalar>("min")),
120  max_(dict.lookup<scalar>("max"))
121 {}
122 
123 
125 (
126  const polyMesh& mesh,
127  Istream& is
128 )
129 :
130  topoSetSource(mesh),
131  fieldName_(checkIs(is)),
132  min_(readScalar(checkIs(is))),
133  max_(readScalar(checkIs(is)))
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
138 
140 {}
141 
142 
143 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144 
145 void Foam::fieldToCell::applyToSet
146 (
147  const topoSetSource::setAction action,
148  topoSet& set
149 ) const
150 {
151  // Try to load field
153  (
154  fieldName_,
155  mesh().time().timeName(),
156  mesh(),
159  false
160  );
161 
162  // Note: should check for volScalarField but that introduces dependency
163  // on volMesh so just use another type with processor-local scope
164  if (!fieldObject.typeHeaderOk<labelIOList>(false))
165  {
167  << "Cannot read field " << fieldName_
168  << " from time " << mesh().time().timeName() << endl;
169  }
170  else if (fieldObject.headerClassName() == "volScalarField")
171  {
172  IFstream str(typeFilePath<labelIOList>(fieldObject));
173 
174  // Read dictionary
175  fieldDictionary fieldDict(fieldObject, fieldObject.headerClassName());
176 
177  scalarField internalVals("internalField", fieldDict, mesh().nCells());
178 
179  applyToSet(action, internalVals, set);
180  }
181  else if (fieldObject.headerClassName() == "volVectorField")
182  {
183  IFstream str(typeFilePath<labelIOList>(fieldObject));
184 
185  // Read dictionary
186  fieldDictionary fieldDict(fieldObject, fieldObject.headerClassName());
187 
188  vectorField internalVals("internalField", fieldDict, mesh().nCells());
189 
190  applyToSet(action, mag(internalVals), set);
191  }
192  else
193  {
195  << "Cannot handle fields of type " << fieldObject.headerClassName()
196  << endl;
197  }
198 }
199 
200 
201 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual ~fieldToCell()
Destructor.
Definition: fieldToCell.C:139
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
IOobject fieldObject(fieldNames[var2field[nVar]], runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE)
Macros for easy insertion into run-time selection tables.
Base class of a source for a topoSet.
Definition: topoSetSource.H:63
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
setAction
Enumeration defining the valid actions.
Definition: topoSetSource.H:82
word timeName
Definition: getTimeIndex.H:3
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
Read field as dictionary (without mesh).
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
Input from file stream.
Definition: IFstream.H:81
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
#define WarningInFunction
Report a warning using Foam::Warning.
Class with constructor to add usage string to table.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
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:97
const word & headerClassName() const
Return name of the class name read from header.
Definition: IOobject.H:309
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812