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