ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / sampling / sampledSurface / sampledPlane / sampledPlane.C
blob3622dcb5853da8da9ab2739687b790fdc952e960
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 \*---------------------------------------------------------------------------*/
26 #include "sampledPlane.H"
27 #include "dictionary.H"
28 #include "polyMesh.H"
29 #include "volFields.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 namespace Foam
37     defineTypeNameAndDebug(sampledPlane, 0);
38     addNamedToRunTimeSelectionTable(sampledSurface, sampledPlane, word, plane);
41 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
43 Foam::sampledPlane::sampledPlane
45     const word& name,
46     const polyMesh& mesh,
47     const plane& planeDesc,
48     const keyType& zoneKey
51     sampledSurface(name, mesh),
52     cuttingPlane(planeDesc),
53     zoneKey_(zoneKey),
54     needsUpdate_(true)
56     if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
57     {
58         Info<< "cellZone " << zoneKey_
59             << " not found - using entire mesh" << endl;
60     }
64 Foam::sampledPlane::sampledPlane
66     const word& name,
67     const polyMesh& mesh,
68     const dictionary& dict
71     sampledSurface(name, mesh, dict),
72     cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
73     zoneKey_(keyType::null),
74     needsUpdate_(true)
76     // make plane relative to the coordinateSystem (Cartesian)
77     // allow lookup from global coordinate systems
78     if (dict.found("coordinateSystem"))
79     {
80         coordinateSystem cs(dict, mesh);
82         point  base = cs.globalPosition(planeDesc().refPoint());
83         vector norm = cs.globalVector(planeDesc().normal());
85         // assign the plane description
86         static_cast<plane&>(*this) = plane(base, norm);
87     }
89     dict.readIfPresent("zone", zoneKey_);
91     if (debug && zoneKey_.size() && mesh.cellZones().findIndex(zoneKey_) < 0)
92     {
93         Info<< "cellZone " << zoneKey_
94             << " not found - using entire mesh" << endl;
95     }
99 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
101 Foam::sampledPlane::~sampledPlane()
105 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
107 bool Foam::sampledPlane::needsUpdate() const
109     return needsUpdate_;
113 bool Foam::sampledPlane::expire()
115     // already marked as expired
116     if (needsUpdate_)
117     {
118         return false;
119     }
121     sampledSurface::clearGeom();
123     needsUpdate_ = true;
124     return true;
128 bool Foam::sampledPlane::update()
130     if (!needsUpdate_)
131     {
132         return false;
133     }
135     sampledSurface::clearGeom();
137     labelList selectedCells = mesh().cellZones().findMatching(zoneKey_).used();
139     if (selectedCells.empty())
140     {
141         reCut(mesh(), true);    // always triangulate. Note:Make option?
142     }
143     else
144     {
145         reCut(mesh(), true, selectedCells);
146     }
148     if (debug)
149     {
150         print(Pout);
151         Pout<< endl;
152     }
154     needsUpdate_ = false;
155     return true;
159 Foam::tmp<Foam::scalarField> Foam::sampledPlane::sample
161     const volScalarField& vField
162 ) const
164     return sampleField(vField);
168 Foam::tmp<Foam::vectorField> Foam::sampledPlane::sample
170     const volVectorField& vField
171 ) const
173     return sampleField(vField);
177 Foam::tmp<Foam::sphericalTensorField> Foam::sampledPlane::sample
179     const volSphericalTensorField& vField
180 ) const
182     return sampleField(vField);
186 Foam::tmp<Foam::symmTensorField> Foam::sampledPlane::sample
188     const volSymmTensorField& vField
189 ) const
191     return sampleField(vField);
195 Foam::tmp<Foam::tensorField> Foam::sampledPlane::sample
197     const volTensorField& vField
198 ) const
200     return sampleField(vField);
204 Foam::tmp<Foam::scalarField> Foam::sampledPlane::interpolate
206     const interpolation<scalar>& interpolator
207 ) const
209     return interpolateField(interpolator);
213 Foam::tmp<Foam::vectorField> Foam::sampledPlane::interpolate
215     const interpolation<vector>& interpolator
216 ) const
218     return interpolateField(interpolator);
221 Foam::tmp<Foam::sphericalTensorField> Foam::sampledPlane::interpolate
223     const interpolation<sphericalTensor>& interpolator
224 ) const
226     return interpolateField(interpolator);
230 Foam::tmp<Foam::symmTensorField> Foam::sampledPlane::interpolate
232     const interpolation<symmTensor>& interpolator
233 ) const
235     return interpolateField(interpolator);
239 Foam::tmp<Foam::tensorField> Foam::sampledPlane::interpolate
241     const interpolation<tensor>& interpolator
242 ) const
244     return interpolateField(interpolator);
248 void Foam::sampledPlane::print(Ostream& os) const
250     os  << "sampledPlane: " << name() << " :"
251         << "  base:" << refPoint()
252         << "  normal:" << normal()
253         << "  faces:" << faces().size()
254         << "  points:" << points().size();
258 // ************************************************************************* //