initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / sampling / sampledSurface / sampledSurface / sampledSurface.C
blob5cbe630486d1c91f828fbd1583df80b9f5d15741
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 \*---------------------------------------------------------------------------*/
27 #include "sampledSurface.H"
28 #include "polyMesh.H"
29 #include "demandDrivenData.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 namespace Foam
36     defineTypeNameAndDebug(sampledSurface, 0);
37     defineRunTimeSelectionTable(sampledSurface, word);
41 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
43 void Foam::sampledSurface::clearGeom() const
45     deleteDemandDrivenData(SfPtr_);
46     deleteDemandDrivenData(magSfPtr_);
47     deleteDemandDrivenData(CfPtr_);
48     area_ = -1;
52 void Foam::sampledSurface::makeSf() const
54     // It is an error to recalculate if the pointer is already set
55     if (SfPtr_)
56     {
57         FatalErrorIn("Foam::sampledSurface::makeSf()")
58             << "face area vectors already exist"
59             << abort(FatalError);
60     }
62     const faceList& theFaces = faces();
63     SfPtr_ = new vectorField(theFaces.size());
65     vectorField& values = *SfPtr_;
66     forAll(theFaces, faceI)
67     {
68         values[faceI] = theFaces[faceI].normal(points());
69     }
73 void Foam::sampledSurface::makeMagSf() const
75     // It is an error to recalculate if the pointer is already set
76     if (magSfPtr_)
77     {
78         FatalErrorIn("Foam::sampledSurface::makeMagSf()")
79             << "mag face areas already exist"
80             << abort(FatalError);
81     }
83     const faceList& theFaces = faces();
84     magSfPtr_ = new scalarField(theFaces.size());
86     scalarField& values = *magSfPtr_;
87     forAll(theFaces, faceI)
88     {
89         values[faceI] = theFaces[faceI].mag(points());
90     }
94 void Foam::sampledSurface::makeCf() const
96     // It is an error to recalculate if the pointer is already set
97     if (CfPtr_)
98     {
99         FatalErrorIn("Foam::sampledSurface::makeCf()")
100             << "face centres already exist"
101             << abort(FatalError);
102     }
104     const faceList& theFaces = faces();
105     CfPtr_ = new vectorField(theFaces.size());
107     vectorField& values = *CfPtr_;
108     forAll(theFaces, faceI)
109     {
110         values[faceI] = theFaces[faceI].centre(points());
111     }
115 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
118 Foam::autoPtr<Foam::sampledSurface>
119 Foam::sampledSurface::New
121     const word& name,
122     const polyMesh& mesh,
123     const dictionary& dict
126     word sampleType(dict.lookup("type"));
127     if (debug)
128     {
129         Info<< "Selecting sampledType " << sampleType << endl;
130     }
132     wordConstructorTable::iterator cstrIter =
133         wordConstructorTablePtr_->find(sampleType);
135     if (cstrIter == wordConstructorTablePtr_->end())
136     {
137         FatalErrorIn
138         (
139             "sampledSurface::New"
140             "(const word&, const polyMesh&, const dictionary&)"
141         )   << "Unknown sample type " << sampleType
142             << endl << endl
143             << "Valid sample types : " << endl
144             << wordConstructorTablePtr_->toc()
145             << exit(FatalError);
146     }
148     return autoPtr<sampledSurface>(cstrIter()(name, mesh, dict));
151 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
153 Foam::sampledSurface::sampledSurface
155     const word& name,
156     const polyMesh& mesh
159     name_(name),
160     mesh_(mesh),
161     interpolate_(false),
162     SfPtr_(NULL),
163     magSfPtr_(NULL),
164     CfPtr_(NULL),
165     area_(-1)
169 // Construct from dictionary
170 Foam::sampledSurface::sampledSurface
172     const word& name,
173     const polyMesh& mesh,
174     const dictionary& dict
177     name_(name),
178     mesh_(mesh),
179     interpolate_(dict.lookupOrDefault("interpolate", false)),
180     SfPtr_(NULL),
181     magSfPtr_(NULL),
182     CfPtr_(NULL),
183     area_(-1)
185     dict.readIfPresent("name", name_);
189 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
191 Foam::sampledSurface::~sampledSurface()
193     clearGeom();
196 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
198 const Foam::vectorField& Foam::sampledSurface::Sf() const
200     if (!SfPtr_)
201     {
202         makeSf();
203     }
205     return *SfPtr_;
209 const Foam::scalarField& Foam::sampledSurface::magSf() const
211     if (!magSfPtr_)
212     {
213         makeMagSf();
214     }
216     return *magSfPtr_;
220 const Foam::vectorField& Foam::sampledSurface::Cf() const
222     if (!CfPtr_)
223     {
224         makeCf();
225     }
227     return *CfPtr_;
231 Foam::scalar Foam::sampledSurface::area() const
233     if (area_ < 0)
234     {
235         area_ = sum(magSf());
236         reduce(area_, sumOp<scalar>());
237     }
239     return area_;
243 // do not project scalar - just copy values
244 Foam::tmp<Foam::Field<Foam::scalar> >
245 Foam::sampledSurface::project(const Field<scalar>& field) const
247     tmp<Field<scalar> > tRes(new Field<scalar>(faces().size()));
248     Field<scalar>& res = tRes();
250     forAll(faces(), faceI)
251     {
252         res[faceI] = field[faceI];
253     }
255     return tRes;
259 Foam::tmp<Foam::Field<Foam::scalar> >
260 Foam::sampledSurface::project(const Field<vector>& field) const
262     tmp<Field<scalar> > tRes(new Field<scalar>(faces().size()));
263     project(tRes(), field);
264     return tRes;
268 Foam::tmp<Foam::Field<Foam::vector> >
269 Foam::sampledSurface::project(const Field<sphericalTensor>& field) const
271     tmp<Field<vector> > tRes(new Field<vector>(faces().size()));
272     project(tRes(), field);
273     return tRes;
277 Foam::tmp<Foam::Field<Foam::vector> >
278 Foam::sampledSurface::project(const Field<symmTensor>& field) const
280     tmp<Field<vector> > tRes(new Field<vector>(faces().size()));
281     project(tRes(), field);
282     return tRes;
286 Foam::tmp<Foam::Field<Foam::vector> >
287 Foam::sampledSurface::project(const Field<tensor>& field) const
289     tmp<Field<vector> > tRes(new Field<vector>(faces().size()));
290     project(tRes(), field);
291     return tRes;
295 void Foam::sampledSurface::print(Ostream& os) const
297     os << type();
300 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
302 Foam::Ostream& Foam::operator<<(Ostream &os, const sampledSurface& s)
304     s.print(os);
305     os.check("Ostream& operator<<(Ostream&, const sampledSurface&");
306     return os;
309 // ************************************************************************* //