1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
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
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"
29 #include "demandDrivenData.H"
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
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_);
52 void Foam::sampledSurface::makeSf() const
54 // It is an error to recalculate if the pointer is already set
57 FatalErrorIn("Foam::sampledSurface::makeSf()")
58 << "face area vectors already exist"
62 const faceList& theFaces = faces();
63 SfPtr_ = new vectorField(theFaces.size());
65 vectorField& values = *SfPtr_;
66 forAll(theFaces, faceI)
68 values[faceI] = theFaces[faceI].normal(points());
73 void Foam::sampledSurface::makeMagSf() const
75 // It is an error to recalculate if the pointer is already set
78 FatalErrorIn("Foam::sampledSurface::makeMagSf()")
79 << "mag face areas already exist"
83 const faceList& theFaces = faces();
84 magSfPtr_ = new scalarField(theFaces.size());
86 scalarField& values = *magSfPtr_;
87 forAll(theFaces, faceI)
89 values[faceI] = theFaces[faceI].mag(points());
94 void Foam::sampledSurface::makeCf() const
96 // It is an error to recalculate if the pointer is already set
99 FatalErrorIn("Foam::sampledSurface::makeCf()")
100 << "face centres already exist"
101 << abort(FatalError);
104 const faceList& theFaces = faces();
105 CfPtr_ = new vectorField(theFaces.size());
107 vectorField& values = *CfPtr_;
108 forAll(theFaces, faceI)
110 values[faceI] = theFaces[faceI].centre(points());
115 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
118 Foam::autoPtr<Foam::sampledSurface>
119 Foam::sampledSurface::New
122 const polyMesh& mesh,
123 const dictionary& dict
126 word sampleType(dict.lookup("type"));
129 Info<< "Selecting sampledType " << sampleType << endl;
132 wordConstructorTable::iterator cstrIter =
133 wordConstructorTablePtr_->find(sampleType);
135 if (cstrIter == wordConstructorTablePtr_->end())
139 "sampledSurface::New"
140 "(const word&, const polyMesh&, const dictionary&)"
141 ) << "Unknown sample type " << sampleType
143 << "Valid sample types : " << endl
144 << wordConstructorTablePtr_->toc()
148 return autoPtr<sampledSurface>(cstrIter()(name, mesh, dict));
151 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
153 Foam::sampledSurface::sampledSurface
169 // Construct from dictionary
170 Foam::sampledSurface::sampledSurface
173 const polyMesh& mesh,
174 const dictionary& dict
179 interpolate_(dict.lookupOrDefault("interpolate", false)),
185 dict.readIfPresent("name", name_);
189 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
191 Foam::sampledSurface::~sampledSurface()
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
198 const Foam::vectorField& Foam::sampledSurface::Sf() const
209 const Foam::scalarField& Foam::sampledSurface::magSf() const
220 const Foam::vectorField& Foam::sampledSurface::Cf() const
231 Foam::scalar Foam::sampledSurface::area() const
235 area_ = sum(magSf());
236 reduce(area_, sumOp<scalar>());
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)
252 res[faceI] = field[faceI];
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);
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);
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);
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);
295 void Foam::sampledSurface::print(Ostream& os) const
300 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
302 Foam::Ostream& Foam::operator<<(Ostream &os, const sampledSurface& s)
305 os.check("Ostream& operator<<(Ostream&, const sampledSurface&");
309 // ************************************************************************* //