ENH: patchCloudSet: new sampledSet
[OpenFOAM-1.7.x.git] / src / sampling / sampledSurface / sampledPlane / sampledPlane.C
blob1efcf81e329da698beae5f0e509f1dc9a014e6b8
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 \*---------------------------------------------------------------------------*/
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 word& zoneName
51     sampledSurface(name, mesh),
52     cuttingPlane(planeDesc),
53     zoneName_(zoneName),
54     needsUpdate_(true)
56     if (debug && zoneName_.size())
57     {
58         if (mesh.cellZones().findZoneID(zoneName_) < 0)
59         {
60             Info<< "cellZone \"" << zoneName_
61                 << "\" not found - using entire mesh" << endl;
62         }
63     }
67 Foam::sampledPlane::sampledPlane
69     const word& name,
70     const polyMesh& mesh,
71     const dictionary& dict
74     sampledSurface(name, mesh, dict),
75     cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
76     zoneName_(word::null),
77     needsUpdate_(true)
79     // make plane relative to the coordinateSystem (Cartesian)
80     // allow lookup from global coordinate systems
81     if (dict.found("coordinateSystem"))
82     {
83         coordinateSystem cs(dict, mesh);
85         point  base = cs.globalPosition(planeDesc().refPoint());
86         vector norm = cs.globalVector(planeDesc().normal());
88         // assign the plane description
89         static_cast<plane&>(*this) = plane(base, norm);
90     }
92     dict.readIfPresent("zone", zoneName_);
94     if (debug && zoneName_.size())
95     {
96         if (mesh.cellZones().findZoneID(zoneName_) < 0)
97         {
98             Info<< "cellZone \"" << zoneName_
99                 << "\" not found - using entire mesh" << endl;
100         }
101     }
106 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
108 Foam::sampledPlane::~sampledPlane()
112 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
114 bool Foam::sampledPlane::needsUpdate() const
116     return needsUpdate_;
120 bool Foam::sampledPlane::expire()
122     // already marked as expired
123     if (needsUpdate_)
124     {
125         return false;
126     }
128     sampledSurface::clearGeom();
130     needsUpdate_ = true;
131     return true;
135 bool Foam::sampledPlane::update()
137     if (!needsUpdate_)
138     {
139         return false;
140     }
142     sampledSurface::clearGeom();
144     label zoneId = -1;
145     if (zoneName_.size())
146     {
147         zoneId = mesh().cellZones().findZoneID(zoneName_);
148     }
150     if (zoneId < 0)
151     {
152         reCut(mesh());
153     }
154     else
155     {
156         reCut(mesh(), mesh().cellZones()[zoneId]);
157     }
159     if (debug)
160     {
161         print(Pout);
162         Pout << endl;
163     }
165     needsUpdate_ = false;
166     return true;
170 Foam::tmp<Foam::scalarField>
171 Foam::sampledPlane::sample
173     const volScalarField& vField
174 ) const
176     return sampleField(vField);
180 Foam::tmp<Foam::vectorField>
181 Foam::sampledPlane::sample
183     const volVectorField& vField
184 ) const
186     return sampleField(vField);
190 Foam::tmp<Foam::sphericalTensorField>
191 Foam::sampledPlane::sample
193     const volSphericalTensorField& vField
194 ) const
196     return sampleField(vField);
200 Foam::tmp<Foam::symmTensorField>
201 Foam::sampledPlane::sample
203     const volSymmTensorField& vField
204 ) const
206     return sampleField(vField);
210 Foam::tmp<Foam::tensorField>
211 Foam::sampledPlane::sample
213     const volTensorField& vField
214 ) const
216     return sampleField(vField);
220 Foam::tmp<Foam::scalarField>
221 Foam::sampledPlane::interpolate
223     const interpolation<scalar>& interpolator
224 ) const
226     return interpolateField(interpolator);
230 Foam::tmp<Foam::vectorField>
231 Foam::sampledPlane::interpolate
233     const interpolation<vector>& interpolator
234 ) const
236     return interpolateField(interpolator);
239 Foam::tmp<Foam::sphericalTensorField>
240 Foam::sampledPlane::interpolate
242     const interpolation<sphericalTensor>& interpolator
243 ) const
245     return interpolateField(interpolator);
249 Foam::tmp<Foam::symmTensorField>
250 Foam::sampledPlane::interpolate
252     const interpolation<symmTensor>& interpolator
253 ) const
255     return interpolateField(interpolator);
259 Foam::tmp<Foam::tensorField>
260 Foam::sampledPlane::interpolate
262     const interpolation<tensor>& interpolator
263 ) const
265     return interpolateField(interpolator);
269 void Foam::sampledPlane::print(Ostream& os) const
271     os  << "sampledPlane: " << name() << " :"
272         << "  base:" << refPoint()
273         << "  normal:" << normal()
274         << "  faces:" << faces().size()
275         << "  points:" << points().size();
279 // ************************************************************************* //