probesGrouping.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-2018 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 "probes.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "IOobjectList.H"
30 #include "stringListOps.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
35 {
37  vectorFields_.clear();
38  sphericalTensorFields_.clear();
39  symmTensorFields_.clear();
40  tensorFields_.clear();
41 
43  surfaceVectorFields_.clear();
46  surfaceTensorFields_.clear();
47 }
48 
49 
51 (
52  const word& fieldName,
53  const word& fieldType
54 )
55 {
56  if (fieldType == volScalarField::typeName)
57  {
58  scalarFields_.append(fieldName);
59  return 1;
60  }
61  else if (fieldType == volVectorField::typeName)
62  {
63  vectorFields_.append(fieldName);
64  return 1;
65  }
66  else if (fieldType == volSphericalTensorField::typeName)
67  {
68  sphericalTensorFields_.append(fieldName);
69  return 1;
70  }
71  else if (fieldType == volSymmTensorField::typeName)
72  {
73  symmTensorFields_.append(fieldName);
74  return 1;
75  }
76  else if (fieldType == volTensorField::typeName)
77  {
78  tensorFields_.append(fieldName);
79  return 1;
80  }
81  else if (fieldType == surfaceScalarField::typeName)
82  {
83  surfaceScalarFields_.append(fieldName);
84  return 1;
85  }
86  else if (fieldType == surfaceVectorField::typeName)
87  {
88  surfaceVectorFields_.append(fieldName);
89  return 1;
90  }
91  else if (fieldType == surfaceSphericalTensorField::typeName)
92  {
93  surfaceSphericalTensorFields_.append(fieldName);
94  return 1;
95  }
96  else if (fieldType == surfaceSymmTensorField::typeName)
97  {
98  surfaceSymmTensorFields_.append(fieldName);
99  return 1;
100  }
101  else if (fieldType == surfaceTensorField::typeName)
102  {
103  surfaceTensorFields_.append(fieldName);
104  return 1;
105  }
106 
107  return 0;
108 }
109 
110 
112 {
113  label nFields = 0;
115 
116  if (loadFromFiles_)
117  {
118  // check files for a particular time
120  wordList allFields = objects.sortedNames();
121 
122  labelList indices = findStrings(fieldSelection_, allFields);
123 
124  forAll(indices, fieldi)
125  {
126  const word& fieldName = allFields[indices[fieldi]];
127 
128  nFields += appendFieldGroup
129  (
130  fieldName,
131  objects.find(fieldName)()->headerClassName()
132  );
133  }
134  }
135  else
136  {
137  // check currently available fields
138  wordList allFields = mesh_.sortedNames();
139  labelList indices = findStrings(fieldSelection_, allFields);
140 
141  forAll(indices, fieldi)
142  {
143  const word& fieldName = allFields[indices[fieldi]];
144 
145  nFields += appendFieldGroup
146  (
147  fieldName,
148  mesh_.find(fieldName)()->type()
149  );
150  }
151  }
152 
153  return nFields;
154 }
155 
156 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
fieldGroup< vector > surfaceVectorFields_
Definition: probes.H:124
void clearFieldGroups()
Clear old field groups.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
fieldGroup< scalar > surfaceScalarFields_
Categorized scalar/vector/tensor surf fields.
Definition: probes.H:123
fieldGroup< sphericalTensor > surfaceSphericalTensorFields_
Definition: probes.H:125
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
static const char *const typeName
Definition: Field.H:105
fieldGroup< symmTensor > surfaceSymmTensorFields_
Definition: probes.H:126
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Operations on lists of strings.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
fieldGroup< tensor > surfaceTensorFields_
Definition: probes.H:127
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
A class for handling words, derived from string.
Definition: word.H:59
wordList sortedNames() const
Return the sorted list of names of the IOobjects.
virtual const word & type() const =0
Runtime type information.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
fieldGroup< scalar > scalarFields_
Categorized scalar/vector/tensor vol fields.
Definition: probes.H:116
const fvMesh & mesh_
Const reference to fvMesh.
Definition: probes.H:92
label appendFieldGroup(const word &fieldName, const word &fieldType)
Append fieldName to the appropriate group.
fieldGroup< tensor > tensorFields_
Definition: probes.H:120
fieldGroup< symmTensor > symmTensorFields_
Definition: probes.H:119
objects
wordReList fieldSelection_
Names of fields to probe.
Definition: probes.H:101
fieldGroup< vector > vectorFields_
Definition: probes.H:117
bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: probes.H:95
label classifyFields()
Classify field types, returns the number of fields.
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
fieldGroup< sphericalTensor > sphericalTensorFields_
Definition: probes.H:118