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() = default;
86 
87  //- Reset format and field list
88  void clear()
89  {
91  formatter.clear();
92  }
93 
94  void setFormatter(const word& writeFormat, const dictionary& dict)
95  {
96  formatter = writer<Type>::New(writeFormat, dict);
97  }
98  };
99 
100 
101  //- Class used for sampling volFields
102  template<class Type>
103  class volFieldSampler
104  :
105  public List<Field<Type>>
106  {
107  //- Name of this collection of values
108  const word name_;
109 
110  public:
111 
112  //- Construct interpolating field to the sampleSets
113  volFieldSampler
114  (
115  const word& interpolationScheme,
117  const PtrList<sampledSet>&
118  );
119 
120  //- Construct mapping field to the sampleSets
121  volFieldSampler
122  (
124  const PtrList<sampledSet>&
125  );
126 
127  //- Construct from components
128  volFieldSampler
129  (
130  const List<Field<Type>>& values,
131  const word& name
132  );
133 
134  //- Return the field name
135  const word& name() const
136  {
137  return name_;
138  }
139  };
140 
141 
142  // Static Data Members
143 
144  //- Output verbosity
145  static bool verbose_;
146 
147 
148  // Private Data
149 
150  //- Const reference to fvMesh
151  const fvMesh& mesh_;
152 
153  //- Keep the dictionary to recreate sets for moving mesh cases
154  dictionary dict_;
155 
156  //- Load fields from files (not from objectRegistry)
157  bool loadFromFiles_;
158 
159  //- Output path
160  fileName outputPath_;
161 
162  //- Mesh search engine
163  meshSearch searchEngine_;
164 
165 
166  // Read from dictionary
167 
168  //- Names of fields to sample
169  wordRes fieldSelection_;
170 
171  //- Interpolation scheme to use
172  word interpolationScheme_;
173 
174  //- Output format to use
175  word writeFormat_;
176 
177  //- Dictionary containing writer options
178  dictionary writeFormatOptions_;
179 
180 
181  // Categorized scalar/vector/tensor fields
182 
183  fieldGroup<scalar> scalarFields_;
184  fieldGroup<vector> vectorFields_;
185  fieldGroup<sphericalTensor> sphericalTensorFields_;
186  fieldGroup<symmTensor> symmTensorFields_;
187  fieldGroup<tensor> tensorFields_;
188 
189 
190  // Merging structures
191 
192  PtrList<coordSet> masterSampledSets_;
193  labelListList indexSets_;
194 
195 
196  // Private Member Functions
197 
198  //- Clear old field groups
199  void clearFieldGroups();
200 
201  //- Classify field types, returns the number of fields
202  label classifyFields();
203 
204  //- Combine points from all processors. Sort by curveDist and produce
205  //- index list. Valid result only on master processor.
206  void combineSampledSets
207  (
208  PtrList<coordSet>& masterSampledSets,
209  labelListList& indexSets
210  );
211 
212  //- Combine values from all processors.
213  // Valid result only on master processor.
214  template<class T>
215  void combineSampledValues
216  (
217  const PtrList<volFieldSampler<T>>& sampledFields,
218  const labelListList& indexSets,
219  PtrList<volFieldSampler<T>>& masterFields
220  );
221 
222  //- Write set on master, return fileName
223  template<class Type>
224  fileName writeSampleFile
225  (
226  const coordSet& masterSampleSet,
227  const PtrList<volFieldSampler<Type>>& masterFields,
228  const label setI,
229  const fileName& timeDir,
230  const writer<Type>& formatter
231  );
232 
233  template<class Type>
234  void sampleAndWrite(fieldGroup<Type>& fields);
235 
236 
237  //- No copy construct
238  sampledSets(const sampledSets&) = delete;
239 
240  //- No copy assignment
241  void operator=(const sampledSets&) = delete;
242 
243 
244 public:
245 
246  //- Runtime type information
247  TypeName("sets");
248 
249 
250  // Constructors
251 
252  //- Construct from Time and dictionary
254  (
255  const word& name,
256  const Time& time,
257  const dictionary& dict
258  );
259 
260  //- Construct for given objectRegistry and dictionary
261  // allow the possibility to load fields from files
263  (
264  const word& name,
265  const objectRegistry&,
266  const dictionary&,
267  const bool loadFromFiles = false
268  );
269 
270 
271  //- Destructor
272  virtual ~sampledSets() = default;
273 
274 
275  // Member Functions
276 
277  //- Enable/disable verbose output
278  // \return old value
279  bool verbose(const bool on);
280 
281  //- Read the sampledSets
282  virtual bool read(const dictionary&);
283 
284  //- Execute, currently does nothing
285  virtual bool execute();
286 
287  //- Sample and write
288  virtual bool write();
289 
290  //- Correct for mesh changes
291  void correct();
292 
293  //- Update for changes of mesh
294  virtual void updateMesh(const mapPolyMesh&);
295 
296  //- Update for mesh point-motion
297  virtual void movePoints(const polyMesh&);
298 
299  //- Update for changes of mesh due to readUpdate
300  virtual void readUpdate(const polyMesh::readUpdateState state);
301 };
302 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 } // End namespace Foam
307 
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
309 
310 #ifdef NoRepository
311  #include "sampledSetsTemplates.C"
312 #endif
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 #endif
317 
318 // ************************************************************************* //
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:300
Foam::sampledSets::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:282
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:65
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:73
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:169
interpolation.H
Foam::DynamicList::clear
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:391
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:155
field
rDeltaTY field()
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
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:123
Foam::sampledSets::~sampledSets
virtual ~sampledSets()=default
Destructor.
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:81
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:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sampledSets::correct
void correct()
Correct for mesh changes.
Definition: sampledSets.C:265
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:163
Foam::functionObject::name
const word & name() const noexcept
Return the name of this functionObject.
Definition: functionObject.C:143
meshSearch.H
Foam::functionObjects::timeFunctionObject::time
const Time & time() const
Return time database.
Definition: timeFunctionObject.H:96
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:291
Foam::autoPtr::clear
void clear() noexcept
Same as reset(nullptr)
Definition: autoPtr.H:176
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:218
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:38