FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

sampledPlane.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "sampledPlane.H"
00027 #include <OpenFOAM/dictionary.H>
00028 #include <OpenFOAM/polyMesh.H>
00029 #include <finiteVolume/volFields.H>
00030 
00031 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00032 
00033 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00034 
00035 namespace Foam
00036 {
00037     defineTypeNameAndDebug(sampledPlane, 0);
00038     addNamedToRunTimeSelectionTable(sampledSurface, sampledPlane, word, plane);
00039 }
00040 
00041 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00042 
00043 Foam::sampledPlane::sampledPlane
00044 (
00045     const word& name,
00046     const polyMesh& mesh,
00047     const plane& planeDesc,
00048     const word& zoneName
00049 )
00050 :
00051     sampledSurface(name, mesh),
00052     cuttingPlane(planeDesc),
00053     zoneName_(zoneName),
00054     needsUpdate_(true)
00055 {
00056     if (debug && zoneName_.size())
00057     {
00058         if (mesh.cellZones().findZoneID(zoneName_) < 0)
00059         {
00060             Info<< "cellZone \"" << zoneName_
00061                 << "\" not found - using entire mesh" << endl;
00062         }
00063     }
00064 }
00065 
00066 
00067 Foam::sampledPlane::sampledPlane
00068 (
00069     const word& name,
00070     const polyMesh& mesh,
00071     const dictionary& dict
00072 )
00073 :
00074     sampledSurface(name, mesh, dict),
00075     cuttingPlane(plane(dict.lookup("basePoint"), dict.lookup("normalVector"))),
00076     zoneName_(word::null),
00077     needsUpdate_(true)
00078 {
00079     // make plane relative to the coordinateSystem (Cartesian)
00080     // allow lookup from global coordinate systems
00081     if (dict.found("coordinateSystem"))
00082     {
00083         coordinateSystem cs(dict, mesh);
00084 
00085         point  base = cs.globalPosition(planeDesc().refPoint());
00086         vector norm = cs.globalVector(planeDesc().normal());
00087 
00088         // assign the plane description
00089         static_cast<plane&>(*this) = plane(base, norm);
00090     }
00091 
00092     dict.readIfPresent("zone", zoneName_);
00093 
00094     if (debug && zoneName_.size())
00095     {
00096         if (mesh.cellZones().findZoneID(zoneName_) < 0)
00097         {
00098             Info<< "cellZone \"" << zoneName_
00099                 << "\" not found - using entire mesh" << endl;
00100         }
00101     }
00102 
00103 }
00104 
00105 
00106 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00107 
00108 Foam::sampledPlane::~sampledPlane()
00109 {}
00110 
00111 
00112 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00113 
00114 bool Foam::sampledPlane::needsUpdate() const
00115 {
00116     return needsUpdate_;
00117 }
00118 
00119 
00120 bool Foam::sampledPlane::expire()
00121 {
00122     // already marked as expired
00123     if (needsUpdate_)
00124     {
00125         return false;
00126     }
00127 
00128     sampledSurface::clearGeom();
00129 
00130     needsUpdate_ = true;
00131     return true;
00132 }
00133 
00134 
00135 bool Foam::sampledPlane::update()
00136 {
00137     if (!needsUpdate_)
00138     {
00139         return false;
00140     }
00141 
00142     sampledSurface::clearGeom();
00143 
00144     label zoneId = -1;
00145     if (zoneName_.size())
00146     {
00147         zoneId = mesh().cellZones().findZoneID(zoneName_);
00148     }
00149 
00150     if (zoneId < 0)
00151     {
00152         reCut(mesh());
00153     }
00154     else
00155     {
00156         reCut(mesh(), mesh().cellZones()[zoneId]);
00157     }
00158 
00159     if (debug)
00160     {
00161         print(Pout);
00162         Pout << endl;
00163     }
00164 
00165     needsUpdate_ = false;
00166     return true;
00167 }
00168 
00169 
00170 Foam::tmp<Foam::scalarField>
00171 Foam::sampledPlane::sample
00172 (
00173     const volScalarField& vField
00174 ) const
00175 {
00176     return sampleField(vField);
00177 }
00178 
00179 
00180 Foam::tmp<Foam::vectorField>
00181 Foam::sampledPlane::sample
00182 (
00183     const volVectorField& vField
00184 ) const
00185 {
00186     return sampleField(vField);
00187 }
00188 
00189 
00190 Foam::tmp<Foam::sphericalTensorField>
00191 Foam::sampledPlane::sample
00192 (
00193     const volSphericalTensorField& vField
00194 ) const
00195 {
00196     return sampleField(vField);
00197 }
00198 
00199 
00200 Foam::tmp<Foam::symmTensorField>
00201 Foam::sampledPlane::sample
00202 (
00203     const volSymmTensorField& vField
00204 ) const
00205 {
00206     return sampleField(vField);
00207 }
00208 
00209 
00210 Foam::tmp<Foam::tensorField>
00211 Foam::sampledPlane::sample
00212 (
00213     const volTensorField& vField
00214 ) const
00215 {
00216     return sampleField(vField);
00217 }
00218 
00219 
00220 Foam::tmp<Foam::scalarField>
00221 Foam::sampledPlane::interpolate
00222 (
00223     const interpolation<scalar>& interpolator
00224 ) const
00225 {
00226     return interpolateField(interpolator);
00227 }
00228 
00229 
00230 Foam::tmp<Foam::vectorField>
00231 Foam::sampledPlane::interpolate
00232 (
00233     const interpolation<vector>& interpolator
00234 ) const
00235 {
00236     return interpolateField(interpolator);
00237 }
00238 
00239 Foam::tmp<Foam::sphericalTensorField>
00240 Foam::sampledPlane::interpolate
00241 (
00242     const interpolation<sphericalTensor>& interpolator
00243 ) const
00244 {
00245     return interpolateField(interpolator);
00246 }
00247 
00248 
00249 Foam::tmp<Foam::symmTensorField>
00250 Foam::sampledPlane::interpolate
00251 (
00252     const interpolation<symmTensor>& interpolator
00253 ) const
00254 {
00255     return interpolateField(interpolator);
00256 }
00257 
00258 
00259 Foam::tmp<Foam::tensorField>
00260 Foam::sampledPlane::interpolate
00261 (
00262     const interpolation<tensor>& interpolator
00263 ) const
00264 {
00265     return interpolateField(interpolator);
00266 }
00267 
00268 
00269 void Foam::sampledPlane::print(Ostream& os) const
00270 {
00271     os  << "sampledPlane: " << name() << " :"
00272         << "  base:" << refPoint()
00273         << "  normal:" << normal()
00274         << "  faces:" << faces().size()
00275         << "  points:" << points().size();
00276 }
00277 
00278 
00279 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines