BUG: PointEdgeWave : n cyclics bool instead of label
[OpenFOAM-1.6.x.git] / src / sampling / sampledSet / sampledSets / sampledSets.H
blobda38f9e1c50ef5d43adf9426fb0517dbadd94a9a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Class
26     Foam::sampledSets
28 Description
29     Set of sets to sample.
30     Call sampledSets.write() to sample&write files.
32 SourceFiles
33     sampledSets.C
35 \*---------------------------------------------------------------------------*/
37 #ifndef sampledSets_H
38 #define sampledSets_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"
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 namespace Foam
52 class objectRegistry;
53 class dictionary;
54 class fvMesh;
56 /*---------------------------------------------------------------------------*\
57                       Class sampledSets Declaration
58 \*---------------------------------------------------------------------------*/
60 class sampledSets
62     public PtrList<sampledSet>
64     // Private classes
66         //- Class used for grouping field types
67         template<class Type>
68         class fieldGroup
69         :
70             public wordList
71         {
72         public:
74             //- Set formatter
75             autoPtr<writer<Type> > formatter;
77             //- Construct null
78             fieldGroup()
79             :
80                 wordList(0),
81                 formatter(NULL)
82             {}
84             void clear()
85             {
86                 wordList::clear();
87                 formatter.clear();
88             }
89         };
92         //- Class used for sampling volFields
93         template <class Type>
94         class volFieldSampler
95         :
96             public List<Field<Type> >
97         {
98             //- Name of this collection of values
99             const word name_;
101         public:
103             //- Construct interpolating field to the sampleSets
104             volFieldSampler
105             (
106                 const word& interpolationScheme,
107                 const GeometricField<Type, fvPatchField, volMesh>& field,
108                 const PtrList<sampledSet>&
109             );
111             //- Construct mapping field to the sampleSets
112             volFieldSampler
113             (
114                 const GeometricField<Type, fvPatchField, volMesh>& field,
115                 const PtrList<sampledSet>&
116             );
118             //- Construct from components
119             volFieldSampler
120             (
121                 const List<Field<Type> >& values,
122                 const word& name
123             );
125             //- Return the field name
126             const word& name() const
127             {
128                 return name_;
129             }
130         };
133     // Static data members
135         //- output verbosity
136         static bool verbose_;
139     // Private data
141         //- Name of this set of sets,
142         //  Also used as the name of the sampledSets directory.
143         word name_;
145         //- Const reference to fvMesh
146         const fvMesh& mesh_;
148         //- Keep the dictionary to recreate sets for moving mesh cases
149         dictionary dict_;
151         //- Load fields from files (not from objectRegistry)
152         bool loadFromFiles_;
154         //- Output path
155         fileName outputPath_;
157         //- Mesh search engine
158         meshSearch searchEngine_;
161         // Read from dictonary
163             //- Names of fields to sample
164             wordList fieldNames_;
166             //- Interpolation scheme to use
167             word interpolationScheme_;
169             //- Output format to use
170             word writeFormat_;
173         // Categorized scalar/vector/tensor fields
175             fieldGroup<scalar> scalarFields_;
176             fieldGroup<vector> vectorFields_;
177             fieldGroup<sphericalTensor> sphericalTensorFields_;
178             fieldGroup<symmTensor> symmTensorFields_;
179             fieldGroup<tensor> tensorFields_;
182         // Merging structures
184             PtrList<coordSet> masterSampledSets_;
185             labelListList indexSets_;
188     // Private Member Functions
190         //- Classify field types, return true if nFields > 0
191         bool checkFieldTypes();
193         //- Find the fields in the list of the given type, return count
194         template<class Type>
195         label grep
196         (
197             fieldGroup<Type>& fieldList,
198             const wordList& fieldTypes
199         ) const;
201         //- Combine points from all processors. Sort by curveDist and produce
202         //  index list. Valid result only on master processor.
203         void combineSampledSets
204         (
205             PtrList<coordSet>& masterSampledSets,
206             labelListList& indexSets
207         );
209         //- Combine values from all processors.
210         //  Valid result only on master processor.
211         template<class T>
212         void combineSampledValues
213         (
214             const PtrList<volFieldSampler<T> >& sampledFields,
215             const labelListList& indexSets,
216             PtrList<volFieldSampler<T> >& masterFields
217         );
219         template<class Type>
220         void writeSampleFile
221         (
222             const coordSet& masterSampleSet,
223             const PtrList<volFieldSampler<Type> >& masterFields,
224             const label setI,
225             const fileName& timeDir,
226             const writer<Type>& formatter
227         );
229         template<class Type>
230         void sampleAndWrite(fieldGroup<Type>& fields);
233         //- Disallow default bitwise copy construct and assignment
234         sampledSets(const sampledSets&);
235         void operator=(const sampledSets&);
238 public:
240     //- Runtime type information
241     TypeName("sets");
244     // Constructors
246         //- Construct for given objectRegistry and dictionary
247         //  allow the possibility to load fields from files
248         sampledSets
249         (
250             const word& name,
251             const objectRegistry&,
252             const dictionary&,
253             const bool loadFromFiles = false
254         );
257     // Destructor
259         virtual ~sampledSets();
262     // Member Functions
264         //- Return name of the set of probes
265         virtual const word& name() const
266         {
267             return name_;
268         }
270         //- set verbosity level
271         void verbose(const bool verbosity = true);
273         //- Execute, currently does nothing
274         virtual void execute();
276         //- Execute at the final time-loop, currently does nothing
277         virtual void end();
279         //- Sample and write
280         virtual void write();
282         //- Read the sampledSets
283         virtual void read(const dictionary&);
285         //- Correct for mesh changes
286         void correct();
288         //- Update for changes of mesh
289         virtual void updateMesh(const mapPolyMesh&);
291         //- Update for mesh point-motion
292         virtual void movePoints(const pointField&);
294         //- Update for changes of mesh due to readUpdate
295         virtual void readUpdate(const polyMesh::readUpdateState state);
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 } // End namespace Foam
303 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 #ifdef NoRepository
306 #   include "sampledSetsTemplates.C"
307 #endif
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 #endif
313 // ************************************************************************* //