1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
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
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/>.
25 Foam::sampledTriSurfaceMesh
28 A sampledSurface from a triSurfaceMesh. It samples on the points/triangles
31 - it either samples cells or (non-coupled) boundary faces
34 - source=cells, interpolate=false:
35 finds per triangle centre the nearest cell centre and uses its value
36 - source=cells, interpolate=true
37 finds per triangle centre the nearest cell centre.
38 Per surface point checks if this nearest cell is the one containing
39 point; otherwise projects the point onto the nearest point on
40 the boundary of the cell (to make sure interpolateCellPoint
41 gets a valid location)
43 - source=boundaryFaces, interpolate=false:
44 finds per triangle centre the nearest point on the boundary
45 (uncoupled faces only) and uses the value (or 0 if the nearest
46 is on an empty boundary)
47 - source=boundaryFaces, interpolate=true:
48 finds per triangle centre the nearest point on the boundary
49 (uncoupled faces only).
50 Per surface point projects the point onto this boundary face
51 (to make sure interpolateCellPoint gets a valid location)
53 - since it finds a nearest per triangle each triangle is guaranteed
54 to be on one processor only. So after stitching (by sampledSurfaces)
55 the original surface should be complete.
58 sampledTriSurfaceMesh.C
60 \*---------------------------------------------------------------------------*/
62 #ifndef sampledTriSurfaceMesh_H
63 #define sampledTriSurfaceMesh_H
65 #include "sampledSurface.H"
66 #include "triSurfaceMesh.H"
67 #include "MeshedSurface.H"
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 /*---------------------------------------------------------------------------*\
77 Class sampledTriSurfaceMesh Declaration
78 \*---------------------------------------------------------------------------*/
80 class sampledTriSurfaceMesh
82 public sampledSurface,
83 public MeshedSurface<face>
86 //- Types of communications
95 //- Private typedefs for convenience
96 typedef MeshedSurface<face> MeshStorage;
101 static const NamedEnum<samplingSource, 2> samplingSourceNames_;
103 //- Surface to sample on
104 const triSurfaceMesh surface_;
106 //- Whether to sample internal cell values or boundary values
107 const samplingSource sampleSource_;
109 //- Track if the surface needs an update
110 mutable bool needsUpdate_;
112 //- Search tree for all non-coupled boundary faces
113 mutable autoPtr<indexedOctree<treeDataFace> > boundaryTreePtr_;
115 //- From local surface triangle to mesh cell/face.
116 labelList sampleElements_;
118 //- Local points to sample per point
119 pointField samplePoints_;
122 // Private Member Functions
124 //- Get tree of all non-coupled boundary faces
125 const indexedOctree<treeDataFace>& nonCoupledboundaryTree() const;
127 //- sample field on faces
128 template <class Type>
129 tmp<Field<Type> > sampleField
131 const GeometricField<Type, fvPatchField, volMesh>& vField
135 template <class Type>
137 interpolateField(const interpolation<Type>&) const;
141 //- Runtime type information
142 TypeName("sampledTriSurfaceMesh");
147 //- Construct from components
148 sampledTriSurfaceMesh
151 const polyMesh& mesh,
152 const word& surfaceName,
153 const samplingSource sampleSource
156 //- Construct from dictionary
157 sampledTriSurfaceMesh
160 const polyMesh& mesh,
161 const dictionary& dict
166 virtual ~sampledTriSurfaceMesh();
171 //- Does the surface need an update?
172 virtual bool needsUpdate() const;
174 //- Mark the surface as needing an update.
175 // May also free up unneeded data.
176 // Return false if surface was already marked as expired.
177 virtual bool expire();
179 //- Update the surface as required.
180 // Do nothing (and return false) if no update was needed
181 virtual bool update();
184 //- Points of surface
185 virtual const pointField& points() const
187 return MeshStorage::points();
191 virtual const faceList& faces() const
193 return MeshStorage::faces();
197 //- sample field on surface
198 virtual tmp<scalarField> sample
200 const volScalarField&
203 //- sample field on surface
204 virtual tmp<vectorField> sample
206 const volVectorField&
209 //- sample field on surface
210 virtual tmp<sphericalTensorField> sample
212 const volSphericalTensorField&
215 //- sample field on surface
216 virtual tmp<symmTensorField> sample
218 const volSymmTensorField&
221 //- sample field on surface
222 virtual tmp<tensorField> sample
224 const volTensorField&
228 //- interpolate field on surface
229 virtual tmp<scalarField> interpolate
231 const interpolation<scalar>&
235 //- interpolate field on surface
236 virtual tmp<vectorField> interpolate
238 const interpolation<vector>&
241 //- interpolate field on surface
242 virtual tmp<sphericalTensorField> interpolate
244 const interpolation<sphericalTensor>&
247 //- interpolate field on surface
248 virtual tmp<symmTensorField> interpolate
250 const interpolation<symmTensor>&
253 //- interpolate field on surface
254 virtual tmp<tensorField> interpolate
256 const interpolation<tensor>&
260 virtual void print(Ostream&) const;
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 } // End namespace Foam
268 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 # include "sampledTriSurfaceMeshTemplates.C"
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
278 // ************************************************************************* //