|
|
|
sampledCuttingPlane.HGo 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-2009 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 the 00013 Free Software Foundation; either version 2 of the License, or (at your 00014 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, write to the Free Software Foundation, 00023 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA 00024 00025 Class 00026 Foam::sampledCuttingPlane 00027 00028 Description 00029 A sampledSurface defined by a plane 00030 00031 SourceFiles 00032 sampledCuttingPlane.C 00033 00034 \*---------------------------------------------------------------------------*/ 00035 00036 #ifndef sampledCuttingPlane_H 00037 #define sampledCuttingPlane_H 00038 00039 #include "sampledSurface.H" 00040 #include "isoSurface.H" 00041 #include "plane.H" 00042 #include "ZoneIDs.H" 00043 #include "fvMeshSubset.H" 00044 00045 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00046 00047 namespace Foam 00048 { 00049 00050 /*---------------------------------------------------------------------------*\ 00051 Class sampledCuttingPlane Declaration 00052 \*---------------------------------------------------------------------------*/ 00053 00054 class sampledCuttingPlane 00055 : 00056 public sampledSurface 00057 { 00058 // Private data 00059 00060 //- Plane 00061 const plane plane_; 00062 00063 //- Merge tolerance 00064 const scalar mergeTol_; 00065 00066 //- Whether to coarsen 00067 const Switch regularise_; 00068 00069 //- zone name/index (if restricted to zones) 00070 mutable cellZoneID zoneID_; 00071 00072 //- for zones: patch to put exposed faces into 00073 mutable word exposedPatchName_; 00074 00075 //- Track if the surface needs an update 00076 mutable bool needsUpdate_; 00077 00078 00079 //- Optional subsetted mesh 00080 autoPtr<fvMeshSubset> subMeshPtr_; 00081 00082 //- Distance to cell centres 00083 autoPtr<volScalarField> cellDistancePtr_; 00084 00085 //- Distance to points 00086 scalarField pointDistance_; 00087 00088 //- Constructed iso surface 00089 autoPtr<isoSurface> isoSurfPtr_; 00090 00091 //- triangles converted to faceList 00092 mutable autoPtr<faceList> facesPtr_; 00093 00094 00095 // Private Member Functions 00096 00097 //- Create iso surface 00098 void createGeometry(); 00099 00100 //- sample field on faces 00101 template <class Type> 00102 tmp<Field<Type> > sampleField 00103 ( 00104 const GeometricField<Type, fvPatchField, volMesh>& vField 00105 ) const; 00106 00107 00108 template <class Type> 00109 tmp<Field<Type> > 00110 interpolateField(const interpolation<Type>&) const; 00111 00112 00113 public: 00114 00115 //- Runtime type information 00116 TypeName("sampledCuttingPlane"); 00117 00118 00119 // Constructors 00120 00121 //- Construct from dictionary 00122 sampledCuttingPlane 00123 ( 00124 const word& name, 00125 const polyMesh& mesh, 00126 const dictionary& dict 00127 ); 00128 00129 00130 // Destructor 00131 00132 virtual ~sampledCuttingPlane(); 00133 00134 00135 // Member Functions 00136 00137 //- Does the surface need an update? 00138 virtual bool needsUpdate() const; 00139 00140 //- Mark the surface as needing an update. 00141 // May also free up unneeded data. 00142 // Return false if surface was already marked as expired. 00143 virtual bool expire(); 00144 00145 //- Update the surface as required. 00146 // Do nothing (and return false) if no update was needed 00147 virtual bool update(); 00148 00149 //- Points of surface 00150 virtual const pointField& points() const 00151 { 00152 return surface().points(); 00153 } 00154 00155 //- Faces of surface 00156 virtual const faceList& faces() const 00157 { 00158 if (facesPtr_.empty()) 00159 { 00160 const triSurface& s = surface(); 00161 00162 facesPtr_.reset(new faceList(s.size())); 00163 00164 forAll(s, i) 00165 { 00166 facesPtr_()[i] = s[i].triFaceFace(); 00167 } 00168 } 00169 return facesPtr_; 00170 } 00171 00172 00173 const isoSurface& surface() const 00174 { 00175 return isoSurfPtr_(); 00176 } 00177 00178 //- sample field on surface 00179 virtual tmp<scalarField> sample 00180 ( 00181 const volScalarField& 00182 ) const; 00183 00184 //- sample field on surface 00185 virtual tmp<vectorField> sample 00186 ( 00187 const volVectorField& 00188 ) const; 00189 00190 //- sample field on surface 00191 virtual tmp<sphericalTensorField> sample 00192 ( 00193 const volSphericalTensorField& 00194 ) const; 00195 00196 //- sample field on surface 00197 virtual tmp<symmTensorField> sample 00198 ( 00199 const volSymmTensorField& 00200 ) const; 00201 00202 //- sample field on surface 00203 virtual tmp<tensorField> sample 00204 ( 00205 const volTensorField& 00206 ) const; 00207 00208 00209 //- interpolate field on surface 00210 virtual tmp<scalarField> interpolate 00211 ( 00212 const interpolation<scalar>& 00213 ) const; 00214 00215 //- interpolate field on surface 00216 virtual tmp<vectorField> interpolate 00217 ( 00218 const interpolation<vector>& 00219 ) const; 00220 00221 //- interpolate field on surface 00222 virtual tmp<sphericalTensorField> interpolate 00223 ( 00224 const interpolation<sphericalTensor>& 00225 ) const; 00226 00227 //- interpolate field on surface 00228 virtual tmp<symmTensorField> interpolate 00229 ( 00230 const interpolation<symmTensor>& 00231 ) const; 00232 00233 //- interpolate field on surface 00234 virtual tmp<tensorField> interpolate 00235 ( 00236 const interpolation<tensor>& 00237 ) const; 00238 00239 //- Write 00240 virtual void print(Ostream&) const; 00241 }; 00242 00243 00244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00245 00246 } // End namespace Foam 00247 00248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00249 00250 #ifdef NoRepository 00251 # include "sampledCuttingPlaneTemplates.C" 00252 #endif 00253 00254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // 00255 00256 #endif 00257 00258 // ************************************************************************* // |