ENH: RASModel.C: clipping input to log
[OpenFOAM-1.7.x.git] / src / sampling / sampledSet / sampledSets / sampledSets.H
blobe7c28698b39955f26d25beda4c79712d73fcaaac
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Class
25     Foam::sampledSets
27 Description
28     Set of sets to sample.
29     Call sampledSets.write() to sample&write files.
31 SourceFiles
32     sampledSets.C
34 \*---------------------------------------------------------------------------*/
36 #ifndef sampledSets_H
37 #define sampledSets_H
39 #include "sampledSet.H"
40 #include "volFieldsFwd.H"
41 #include "meshSearch.H"
42 #include "interpolation.H"
43 #include "coordSet.H"
44 #include "writer.H"
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 namespace Foam
51 class objectRegistry;
52 class dictionary;
53 class fvMesh;
55 /*---------------------------------------------------------------------------*\
56                       Class sampledSets Declaration
57 \*---------------------------------------------------------------------------*/
59 class sampledSets
61     public PtrList<sampledSet>
63     // Private classes
65         //- Class used for grouping field types
66         template<class Type>
67         class fieldGroup
68         :
69             public wordList
70         {
71         public:
73             //- Set formatter
74             autoPtr<writer<Type> > formatter;
76             //- Construct null
77             fieldGroup()
78             :
79                 wordList(0),
80                 formatter(NULL)
81             {}
83             void clear()
84             {
85                 wordList::clear();
86                 formatter.clear();
87             }
88         };
91         //- Class used for sampling volFields
92         template <class Type>
93         class volFieldSampler
94         :
95             public List<Field<Type> >
96         {
97             //- Name of this collection of values
98             const word name_;
100         public:
102             //- Construct interpolating field to the sampleSets
103             volFieldSampler
104             (
105                 const word& interpolationScheme,
106                 const GeometricField<Type, fvPatchField, volMesh>& field,
107                 const PtrList<sampledSet>&
108             );
110             //- Construct mapping field to the sampleSets
111             volFieldSampler
112             (
113                 const GeometricField<Type, fvPatchField, volMesh>& field,
114                 const PtrList<sampledSet>&
115             );
117             //- Construct from components
118             volFieldSampler
119             (
120                 const List<Field<Type> >& values,
121                 const word& name
122             );
124             //- Return the field name
125             const word& name() const
126             {
127                 return name_;
128             }
129         };
132     // Static data members
134         //- output verbosity
135         static bool verbose_;
138     // Private data
140         //- Name of this set of sets,
141         //  Also used as the name of the sampledSets directory.
142         word name_;
144         //- Const reference to fvMesh
145         const fvMesh& mesh_;
147         //- Keep the dictionary to recreate sets for moving mesh cases
148         dictionary dict_;
150         //- Load fields from files (not from objectRegistry)
151         bool loadFromFiles_;
153         //- Output path
154         fileName outputPath_;
156         //- Mesh search engine
157         meshSearch searchEngine_;
160         // Read from dictonary
162             //- Names of fields to sample
163             wordList fieldNames_;
165             //- Interpolation scheme to use
166             word interpolationScheme_;
168             //- Output format to use
169             word writeFormat_;
172         // Categorized scalar/vector/tensor fields
174             fieldGroup<scalar> scalarFields_;
175             fieldGroup<vector> vectorFields_;
176             fieldGroup<sphericalTensor> sphericalTensorFields_;
177             fieldGroup<symmTensor> symmTensorFields_;
178             fieldGroup<tensor> tensorFields_;
181         // Merging structures
183             PtrList<coordSet> masterSampledSets_;
184             labelListList indexSets_;
187     // Private Member Functions
189         //- Classify field types, return true if nFields > 0
190         bool checkFieldTypes();
192         //- Find the fields in the list of the given type, return count
193         template<class Type>
194         label grep
195         (
196             fieldGroup<Type>& fieldList,
197             const wordList& fieldTypes
198         ) const;
200         //- Combine points from all processors. Sort by curveDist and produce
201         //  index list. Valid result only on master processor.
202         void combineSampledSets
203         (
204             PtrList<coordSet>& masterSampledSets,
205             labelListList& indexSets
206         );
208         //- Combine values from all processors.
209         //  Valid result only on master processor.
210         template<class T>
211         void combineSampledValues
212         (
213             const PtrList<volFieldSampler<T> >& sampledFields,
214             const labelListList& indexSets,
215             PtrList<volFieldSampler<T> >& masterFields
216         );
218         template<class Type>
219         void writeSampleFile
220         (
221             const coordSet& masterSampleSet,
222             const PtrList<volFieldSampler<Type> >& masterFields,
223             const label setI,
224             const fileName& timeDir,
225             const writer<Type>& formatter
226         );
228         template<class Type>
229         void sampleAndWrite(fieldGroup<Type>& fields);
232         //- Disallow default bitwise copy construct and assignment
233         sampledSets(const sampledSets&);
234         void operator=(const sampledSets&);
237 public:
239     //- Runtime type information
240     TypeName("sets");
243     // Constructors
245         //- Construct for given objectRegistry and dictionary
246         //  allow the possibility to load fields from files
247         sampledSets
248         (
249             const word& name,
250             const objectRegistry&,
251             const dictionary&,
252             const bool loadFromFiles = false
253         );
256     // Destructor
258         virtual ~sampledSets();
261     // Member Functions
263         //- Return name of the set of probes
264         virtual const word& name() const
265         {
266             return name_;
267         }
269         //- set verbosity level
270         void verbose(const bool verbosity = true);
272         //- Execute, currently does nothing
273         virtual void execute();
275         //- Execute at the final time-loop, currently does nothing
276         virtual void end();
278         //- Sample and write
279         virtual void write();
281         //- Read the sampledSets
282         virtual void read(const dictionary&);
284         //- Correct for mesh changes
285         void correct();
287         //- Update for changes of mesh
288         virtual void updateMesh(const mapPolyMesh&);
290         //- Update for mesh point-motion
291         virtual void movePoints(const pointField&);
293         //- Update for changes of mesh due to readUpdate
294         virtual void readUpdate(const polyMesh::readUpdateState state);
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 } // End namespace Foam
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
304 #ifdef NoRepository
305 #   include "sampledSetsTemplates.C"
306 #endif
308 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 #endif
312 // ************************************************************************* //