sampledSets.H
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-2016 OpenFOAM Foundation
9  Copyright (C) 2015-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 Class
28  Foam::sampledSets
29 
30 Description
31  Set of sets to sample.
32  Call sampledSets.write() to sample and write files.
33 
34 SourceFiles
35  sampledSets.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef sampledSets_H
40 #define sampledSets_H
41 
42 #include "regionFunctionObject.H"
43 #include "sampledSet.H"
44 #include "volFieldsFwd.H"
45 #include "meshSearch.H"
46 #include "interpolation.H"
47 #include "coordSet.H"
48 #include "writer.H"
49 #include "wordRes.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 namespace Foam
54 {
55 
56 // Forward declarations
57 class Time;
58 class objectRegistry;
59 class dictionary;
60 class fvMesh;
61 
62 /*---------------------------------------------------------------------------*\
63  Class sampledSets Declaration
64 \*---------------------------------------------------------------------------*/
65 
66 class sampledSets
67 :
69  public PtrList<sampledSet>
70 {
71  // Private classes
72 
73  //- Class used for grouping field types
74  template<class Type>
75  class fieldGroup
76  :
77  public DynamicList<word>
78  {
79  public:
80 
81  //- The set formatter
82  autoPtr<writer<Type>> formatter;
83 
84  //- Construct null
85  fieldGroup()
86  :
88  formatter(nullptr)
89  {}
90 
91  //- Construct for a particular format
92  fieldGroup(const word& writeFormat)
93  :
95  formatter(writer<Type>::New(writeFormat))
96  {}
97 
98  //- Reset format and field list
99  void clear()
100  {
102  formatter.clear();
103  }
104 
105  //- Assign a new formatter
106  void operator=(const word& writeFormat)
107  {
108  formatter = writer<Type>::New(writeFormat);
109  }
110  };
111 
112 
113  //- Class used for sampling volFields
114  template<class Type>
115  class volFieldSampler
116  :
117  public List<Field<Type>>
118  {
119  //- Name of this collection of values
120  const word name_;
121 
122  public:
123 
124  //- Construct interpolating field to the sampleSets
125  volFieldSampler
126  (
127  const word& interpolationScheme,
129  const PtrList<sampledSet>&
130  );
131 
132  //- Construct mapping field to the sampleSets
133  volFieldSampler
134  (
136  const PtrList<sampledSet>&
137  );
138 
139  //- Construct from components
140  volFieldSampler
141  (
142  const List<Field<Type>>& values,
143  const word& name
144  );
145 
146  //- Return the field name
147  const word& name() const
148  {
149  return name_;
150  }
151  };
152 
153 
154  // Static data members
155 
156  //- Output verbosity
157  static bool verbose_;
158 
159 
160  // Private data
161 
162  //- Const reference to fvMesh
163  const fvMesh& mesh_;
164 
165  //- Keep the dictionary to recreate sets for moving mesh cases
166  dictionary dict_;
167 
168  //- Load fields from files (not from objectRegistry)
169  bool loadFromFiles_;
170 
171  //- Output path
172  fileName outputPath_;
173 
174  //- Mesh search engine
175  meshSearch searchEngine_;
176 
177 
178  // Read from dictonary
179 
180  //- Names of fields to sample
181  wordRes fieldSelection_;
182 
183  //- Interpolation scheme to use
184  word interpolationScheme_;
185 
186  //- Output format to use
187  word writeFormat_;
188 
189 
190  // Categorized scalar/vector/tensor fields
191 
192  fieldGroup<scalar> scalarFields_;
193  fieldGroup<vector> vectorFields_;
194  fieldGroup<sphericalTensor> sphericalTensorFields_;
195  fieldGroup<symmTensor> symmTensorFields_;
196  fieldGroup<tensor> tensorFields_;
197 
198 
199  // Merging structures
200 
201  PtrList<coordSet> masterSampledSets_;
202  labelListList indexSets_;
203 
204 
205  // Private Member Functions
206 
207  //- Clear old field groups
208  void clearFieldGroups();
209 
210  //- Classify field types, returns the number of fields
211  label classifyFields();
212 
213  //- Combine points from all processors. Sort by curveDist and produce
214  // index list. Valid result only on master processor.
215  void combineSampledSets
216  (
217  PtrList<coordSet>& masterSampledSets,
218  labelListList& indexSets
219  );
220 
221  //- Combine values from all processors.
222  // Valid result only on master processor.
223  template<class T>
224  void combineSampledValues
225  (
226  const PtrList<volFieldSampler<T>>& sampledFields,
227  const labelListList& indexSets,
228  PtrList<volFieldSampler<T>>& masterFields
229  );
230 
231  //- Write set on master, return fileName
232  template<class Type>
233  fileName writeSampleFile
234  (
235  const coordSet& masterSampleSet,
236  const PtrList<volFieldSampler<Type>>& masterFields,
237  const label setI,
238  const fileName& timeDir,
239  const writer<Type>& formatter
240  );
241 
242  template<class Type>
243  void sampleAndWrite(fieldGroup<Type>& fields);
244 
245 
246  //- No copy construct
247  sampledSets(const sampledSets&) = delete;
248 
249  //- No copy assignment
250  void operator=(const sampledSets&) = delete;
251 
252 
253 public:
254 
255  //- Runtime type information
256  TypeName("sets");
257 
258 
259  // Constructors
260 
261  //- Construct from Time and dictionary
263  (
264  const word& name,
265  const Time& time,
266  const dictionary& dict
267  );
268 
269  //- Construct for given objectRegistry and dictionary
270  // allow the possibility to load fields from files
272  (
273  const word& name,
274  const objectRegistry&,
275  const dictionary&,
276  const bool loadFromFiles = false
277  );
278 
279 
280  //- Destructor
281  virtual ~sampledSets() = default;
282 
283 
284  // Member Functions
285 
286  //- Enable/disable verbose output
287  // \return old value
288  bool verbose(const bool on);
289 
290  //- Read the sampledSets
291  virtual bool read(const dictionary&);
292 
293  //- Execute, currently does nothing
294  virtual bool execute();
295 
296  //- Sample and write
297  virtual bool write();
298 
299  //- Correct for mesh changes
300  void correct();
301 
302  //- Update for changes of mesh
303  virtual void updateMesh(const mapPolyMesh&);
304 
305  //- Update for mesh point-motion
306  virtual void movePoints(const polyMesh&);
307 
308  //- Update for changes of mesh due to readUpdate
309  virtual void readUpdate(const polyMesh::readUpdateState state);
310 };
311 
312 
313 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314 
315 } // End namespace Foam
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 #ifdef NoRepository
320  #include "sampledSetsTemplates.C"
321 #endif
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 #endif
326 
327 // ************************************************************************* //
volFieldsFwd.H
wordRes.H
Foam::sampledSets::readUpdate
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:298
Foam::sampledSets::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:280
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::sampledSets
Set of sets to sample. Call sampledSets.write() to sample and write files.
Definition: sampledSets.H:65
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:55
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
coordSet.H
Foam::sampledSets::write
virtual bool write()
Sample and write.
Definition: sampledSets.C:167
interpolation.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
sampledSet.H
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::sampledSets::verbose
bool verbose(const bool on)
Enable/disable verbose output.
Definition: sampledSets.C:153
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:348
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
regionFunctionObject.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::sampledSets::~sampledSets
virtual ~sampledSets()=default
Destructor.
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:80
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam::PtrList< sampledSet >::clear
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:97
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sampledSets::correct
void correct()
Correct for mesh changes.
Definition: sampledSets.C:263
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::sampledSets::execute
virtual bool execute()
Execute, currently does nothing.
Definition: sampledSets.C:161
meshSearch.H
Foam::functionObjects::timeFunctionObject::time
const Time & time() const
Return time database.
Definition: timeFunctionObject.H:96
Foam::functionObject::name
const word & name() const
Return the name of this functionObject.
Definition: functionObject.C:131
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::sampledSets::TypeName
TypeName("sets")
Runtime type information.
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::sampledSets::movePoints
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:289
Foam::autoPtr::clear
void clear() noexcept
Same as reset(nullptr)
Definition: autoPtr.H:181
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::functionObjects::regionFunctionObject
Specialization of Foam::functionObject for a region and providing a reference to the region Foam::obj...
Definition: regionFunctionObject.H:90
sampledSetsTemplates.C
writer.H
Foam::GeometricField< Type, fvPatchField, volMesh >
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::sampledSets::read
virtual bool read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:216
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:38