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-2021 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;
114  clearFieldGroups();
115 
116  // Check currently available fields
117  forAll(fields_, fieldi)
118  {
119  const word& fieldName = fields_[fieldi];
120 
121  if (mesh_.objectRegistry::found(fieldName))
122  {
123  nFields += appendFieldGroup
124  (
125  fieldName,
126  mesh_.find(fieldName)()->type()
127  );
128  }
129  }
130 
131  return nFields;
132 }
133 
134 
135 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:236
static const char *const typeName
Definition: Field.H:106
fieldGroup< symmTensor > surfaceSymmTensorFields_
Definition: probes.H:122
fieldGroup< tensor > tensorFields_
Definition: probes.H:116
void clearFieldGroups()
Clear old field groups.
label classifyFields()
Classify field types, returns the number of fields.
fieldGroup< tensor > surfaceTensorFields_
Definition: probes.H:123
fieldGroup< symmTensor > symmTensorFields_
Definition: probes.H:115
fieldGroup< sphericalTensor > surfaceSphericalTensorFields_
Definition: probes.H:121
fieldGroup< scalar > surfaceScalarFields_
Categorised scalar/vector/tensor surf fields.
Definition: probes.H:119
fieldGroup< vector > surfaceVectorFields_
Definition: probes.H:120
fieldGroup< scalar > scalarFields_
Categorised scalar/vector/tensor vol fields.
Definition: probes.H:112
label appendFieldGroup(const word &fieldName, const word &fieldType)
Append fieldName to the appropriate group.
fieldGroup< sphericalTensor > sphericalTensorFields_
Definition: probes.H:114
fieldGroup< vector > vectorFields_
Definition: probes.H:113
A class for handling words, derived from string.
Definition: word.H:62
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
Operations on lists of strings.
Foam::surfaceFields.