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-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "functionObjectState.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 class objectRegistry;
54 class dictionary;
55 class fvMesh;
56 
57 /*---------------------------------------------------------------------------*\
58  Class sampledSets Declaration
59 \*---------------------------------------------------------------------------*/
60 
61 class sampledSets
62 :
63  public functionObjectState,
64  public PtrList<sampledSet>
65 {
66  // Private classes
67 
68  //- Class used for grouping field types
69  template<class Type>
70  class fieldGroup
71  :
72  public DynamicList<word>
73  {
74  public:
75 
76  //- The set formatter
78 
79  //- Construct null
80  fieldGroup()
81  :
82  DynamicList<word>(0),
83  formatter(NULL)
84  {}
85 
86  //- Construct for a particular format
87  fieldGroup(const word& writeFormat)
88  :
89  DynamicList<word>(0),
90  formatter(writer<Type>::New(writeFormat))
91  {}
92 
93  //- Reset format and field list
94  void clear()
95  {
97  formatter.clear();
98  }
99 
100  //- Assign a new formatter
101  void operator=(const word& writeFormat)
102  {
103  formatter = writer<Type>::New(writeFormat);
104  }
105 
106  };
107 
108 
109  //- Class used for sampling volFields
110  template<class Type>
111  class volFieldSampler
112  :
113  public List<Field<Type> >
114  {
115  //- Name of this collection of values
116  const word name_;
117 
118  public:
119 
120  //- Construct interpolating field to the sampleSets
122  (
123  const word& interpolationScheme,
125  const PtrList<sampledSet>&
126  );
127 
128  //- Construct mapping field to the sampleSets
130  (
132  const PtrList<sampledSet>&
133  );
134 
135  //- Construct from components
137  (
138  const List<Field<Type> >& values,
139  const word& name
140  );
141 
142  //- Return the field name
143  const word& name() const
144  {
145  return name_;
146  }
147  };
148 
149 
150  // Static data members
151 
152  //- Output verbosity
153  static bool verbose_;
154 
155 
156  // Private data
157 
158  //- Const reference to fvMesh
159  const fvMesh& mesh_;
160 
161  //- Keep the dictionary to recreate sets for moving mesh cases
163 
164  //- Load fields from files (not from objectRegistry)
165  bool loadFromFiles_;
166 
167  //- Output path
169 
170  //- Mesh search engine
172 
173 
174  // Read from dictonary
175 
176  //- Names of fields to sample
178 
179  //- Interpolation scheme to use
181 
182  //- Output format to use
184 
185 
186  // Categorized scalar/vector/tensor fields
187 
193 
194 
195  // Merging structures
196 
199 
200 
201  // Private Member Functions
202 
203  //- Clear old field groups
204  void clearFieldGroups();
205 
206  //- Append fieldName to the appropriate group
207  label appendFieldGroup(const word& fieldName, const word& fieldType);
208 
209  //- Classify field types, returns the number of fields
211 
212  //- Combine points from all processors. Sort by curveDist and produce
213  // index list. Valid result only on master processor.
214  void combineSampledSets
215  (
216  PtrList<coordSet>& masterSampledSets,
217  labelListList& indexSets
218  );
219 
220  //- Combine values from all processors.
221  // Valid result only on master processor.
222  template<class T>
224  (
225  const PtrList<volFieldSampler<T> >& sampledFields,
226  const labelListList& indexSets,
227  PtrList<volFieldSampler<T> >& masterFields
228  );
229 
230  template<class Type>
231  void writeSampleFile
232  (
233  const coordSet& masterSampleSet,
234  const PtrList<volFieldSampler<Type> >& masterFields,
235  const label setI,
236  const fileName& timeDir,
237  const writer<Type>& formatter
238  );
239 
240  template<class Type>
242 
243 
244  //- Disallow default bitwise copy construct and assignment
245  sampledSets(const sampledSets&);
246  void operator=(const sampledSets&);
247 
248 
249 public:
250 
251  //- Runtime type information
252  TypeName("sets");
253 
254 
255  // Constructors
256 
257  //- Construct for given objectRegistry and dictionary
258  // allow the possibility to load fields from files
260  (
261  const word& name,
262  const objectRegistry&,
263  const dictionary&,
264  const bool loadFromFiles = false
265  );
266 
267 
268  //- Destructor
269  virtual ~sampledSets();
270 
271 
272  // Member Functions
273 
274  //- Set verbosity level
275  void verbose(const bool verbosity = true);
276 
277  //- Execute, currently does nothing
278  virtual void execute();
279 
280  //- Execute at the final time-loop, currently does nothing
281  virtual void end();
282 
283  //- Called when time was set at the end of the Time::operator++
284  virtual void timeSet();
285 
286  //- Sample and write
287  virtual void write();
288 
289  //- Read the sampledSets
290  virtual void read(const dictionary&);
291 
292  //- Correct for mesh changes
293  void correct();
294 
295  //- Update for changes of mesh
296  virtual void updateMesh(const mapPolyMesh&);
297 
298  //- Update for mesh point-motion
299  virtual void movePoints(const polyMesh&);
300 
301  //- Update for changes of mesh due to readUpdate
302  virtual void readUpdate(const polyMesh::readUpdateState state);
303 };
304 
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 } // End namespace Foam
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
312 #ifdef NoRepository
313 # include "sampledSetsTemplates.C"
314 #endif
315 
316 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
317 
318 #endif
319 
320 // ************************************************************************* //
Foam::sampledSets::fieldGroup::clear
void clear()
Reset format and field list.
Definition: sampledSets.H:93
Foam::sampledSets::~sampledSets
virtual ~sampledSets()
Destructor.
Definition: sampledSets.C:167
Foam::sampledSets::mesh_
const fvMesh & mesh_
Const reference to fvMesh.
Definition: sampledSets.H:158
volFieldsFwd.H
Foam::sampledSets::fieldGroup::fieldGroup
fieldGroup(const word &writeFormat)
Construct for a particular format.
Definition: sampledSets.H:86
Foam::sampledSets::readUpdate
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate.
Definition: sampledSets.C:316
Foam::sampledSets::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition: sampledSets.C:304
Foam::sampledSets::execute
virtual void execute()
Execute, currently does nothing.
Definition: sampledSets.C:179
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:57
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::sampledSets::vectorFields_
fieldGroup< vector > vectorFields_
Definition: sampledSets.H:188
Foam::functionObjectState
Base class for function objects, adding functionality to read/write state information (data required ...
Definition: functionObjectState.H:54
Foam::sampledSets
Set of sets to sample. Call sampledSets.write() to sample&write files.
Definition: sampledSets.H:60
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::sampledSets::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: sampledSets.C:185
Foam::sampledSets::read
virtual void read(const dictionary &)
Read the sampledSets.
Definition: sampledSets.C:240
Foam::sampledSets::searchEngine_
meshSearch searchEngine_
Mesh search engine.
Definition: sampledSets.H:170
fields
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:127
coordSet.H
Foam::sampledSets::fieldGroup::formatter
autoPtr< writer< Type > > formatter
The set formatter.
Definition: sampledSets.H:76
Foam::sampledSets::verbose
void verbose(const bool verbosity=true)
Set verbosity level.
Definition: sampledSets.C:173
Foam::sampledSets::tensorFields_
fieldGroup< tensor > tensorFields_
Definition: sampledSets.H:191
Foam::sampledSets::volFieldSampler::name
const word & name() const
Return the field name.
Definition: sampledSets.H:142
Foam::sampledSets::write
virtual void write()
Sample and write.
Definition: sampledSets.C:197
functionObjectState.H
Foam::sampledSets::symmTensorFields_
fieldGroup< symmTensor > symmTensorFields_
Definition: sampledSets.H:190
Foam::sampledSets::clearFieldGroups
void clearFieldGroups()
Clear old field groups.
Definition: sampledSetsGrouping.C:33
interpolation.H
Foam::sampledSets::indexSets_
labelListList indexSets_
Definition: sampledSets.H:197
Foam::sampledSets::sampleAndWrite
void sampleAndWrite(fieldGroup< Type > &fields)
Definition: sampledSetsTemplates.C:228
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::sampledSets::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: sampledSets.C:191
Foam::sampledSets::outputPath_
fileName outputPath_
Output path.
Definition: sampledSets.H:167
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
sampledSet.H
Foam::label
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
Foam::Field< Type >
Foam::DynamicList::clear
void clear()
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:242
Foam::sampledSets::appendFieldGroup
label appendFieldGroup(const word &fieldName, const word &fieldType)
Append fieldName to the appropriate group.
Definition: sampledSetsGrouping.C:44
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::sampledSets::interpolationScheme_
word interpolationScheme_
Interpolation scheme to use.
Definition: sampledSets.H:179
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::sampledSets::classifyFields
label classifyFields()
Classify field types, returns the number of fields.
Definition: sampledSetsGrouping.C:79
Foam::sampledSets::volFieldSampler::volFieldSampler
volFieldSampler(const word &interpolationScheme, const GeometricField< Type, fvPatchField, volMesh > &field, const PtrList< sampledSet > &)
Construct interpolating field to the sampleSets.
Definition: sampledSetsTemplates.C:34
Foam::sampledSets::writeSampleFile
void writeSampleFile(const coordSet &masterSampleSet, const PtrList< volFieldSampler< Type > > &masterFields, const label setI, const fileName &timeDir, const writer< Type > &formatter)
Definition: sampledSetsTemplates.C:126
Foam::writer
Base class for graphics format writing. Entry points are.
Definition: writer.H:78
Foam::sampledSets::combineSampledSets
void combineSampledSets(PtrList< coordSet > &masterSampledSets, labelListList &indexSets)
Combine points from all processors. Sort by curveDist and produce.
Definition: sampledSets.C:46
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:49
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::sampledSets::correct
void correct()
Correct for mesh changes.
Definition: sampledSets.C:286
Foam::sampledSets::loadFromFiles_
bool loadFromFiles_
Load fields from files (not from objectRegistry)
Definition: sampledSets.H:164
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:88
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::sampledSets::operator=
void operator=(const sampledSets &)
Foam::sampledSets::sphericalTensorFields_
fieldGroup< sphericalTensor > sphericalTensorFields_
Definition: sampledSets.H:189
Foam::sampledSets::volFieldSampler::name_
const word name_
Name of this collection of values.
Definition: sampledSets.H:115
meshSearch.H
Foam::sampledSets::scalarFields_
fieldGroup< scalar > scalarFields_
Definition: sampledSets.H:187
Foam::sampledSets::combineSampledValues
void combineSampledValues(const PtrList< volFieldSampler< T > > &sampledFields, const labelListList &indexSets, PtrList< volFieldSampler< T > > &masterFields)
Combine values from all processors.
Definition: sampledSetsTemplates.C:178
Foam::sampledSets::fieldSelection_
wordReList fieldSelection_
Names of fields to sample.
Definition: sampledSets.H:176
Foam::sampledSets::dict_
dictionary dict_
Keep the dictionary to recreate sets for moving mesh cases.
Definition: sampledSets.H:161
Foam::sampledSets::volFieldSampler
Class used for sampling volFields.
Definition: sampledSets.H:110
wordReList.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sampledSets::TypeName
TypeName("sets")
Runtime type information.
Foam::sampledSets::fieldGroup::fieldGroup
fieldGroup()
Construct null.
Definition: sampledSets.H:79
Foam::sampledSets::writeFormat_
word writeFormat_
Output format to use.
Definition: sampledSets.H:182
Foam::sampledSets::fieldGroup
Class used for grouping field types.
Definition: sampledSets.H:69
Foam::sampledSets::movePoints
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition: sampledSets.C:310
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
sampledSetsTemplates.C
Foam::sampledSets::fieldGroup::operator=
void operator=(const word &writeFormat)
Assign a new formatter.
Definition: sampledSets.H:100
writer.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::sampledSets::verbose_
static bool verbose_
Output verbosity.
Definition: sampledSets.H:152
Foam::functionObjectState::name
const word & name() const
Return the name.
Definition: functionObjectState.C:58
Foam::sampledSets::masterSampledSets_
PtrList< coordSet > masterSampledSets_
Definition: sampledSets.H:196
Foam::sampledSets::sampledSets
sampledSets(const sampledSets &)
Disallow default bitwise copy construct and assignment.
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:35