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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "sampledIsoSurfaceCell.H"
27 #include "dictionary.H"
28 #include "volFields.H"
29 #include "volPointInterpolation.H"
30 #include "addToRunTimeSelectionTable.H"
32 #include "isoSurfaceCell.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(sampledIsoSurfaceCell, 0);
39 addNamedToRunTimeSelectionTable
42 sampledIsoSurfaceCell,
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 bool Foam::sampledIsoSurfaceCell::updateGeometry() const
52 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
55 if (fvm.time().timeIndex() == prevTimeIndex_)
60 prevTimeIndex_ = fvm.time().timeIndex();
62 // Clear any stored topo
66 sampledSurface::clearGeom();
68 // Optionally read volScalarField
69 autoPtr<volScalarField> readFieldPtr_;
71 // 1. see if field in database
72 // 2. see if field can be read
73 const volScalarField* cellFldPtr = NULL;
74 if (fvm.foundObject<volScalarField>(isoField_))
78 Info<< "sampledIsoSurfaceCell::updateGeometry() : lookup "
82 cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
86 // Bit of a hack. Read field and store.
90 Info<< "sampledIsoSurfaceCell::updateGeometry() : reading "
91 << isoField_ << " from time " <<fvm.time().timeName()
102 fvm.time().timeName(),
112 cellFldPtr = readFieldPtr_.operator->();
114 const volScalarField& cellFld = *cellFldPtr;
116 tmp<pointScalarField> pointFld
118 volPointInterpolation::New(fvm).interpolate(cellFld)
123 //- From point field and interpolated cell.
124 scalarField cellAvg(fvm.nCells(), scalar(0.0));
125 labelField nPointCells(fvm.nCells(), 0);
127 for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
129 const labelList& pCells = fvm.pointCells(pointI);
133 label cellI = pCells[i];
135 cellAvg[cellI] += pointFld().internalField()[pointI];
136 nPointCells[cellI]++;
140 forAll(cellAvg, cellI)
142 cellAvg[cellI] /= nPointCells[cellI];
145 const isoSurfaceCell iso
149 pointFld().internalField(),
154 const_cast<sampledIsoSurfaceCell&>
157 ).triSurface::operator=(iso);
158 meshCells_ = iso.meshCells();
162 //- Direct from cell field and point field. Gives bad continuity.
163 const isoSurfaceCell iso
166 cellFld.internalField(),
167 pointFld().internalField(),
172 const_cast<sampledIsoSurfaceCell&>
175 ).triSurface::operator=(iso);
176 meshCells_ = iso.meshCells();
182 Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
184 << " regularise : " << regularise_ << nl
185 << " average : " << average_ << nl
186 << " isoField : " << isoField_ << nl
187 << " isoValue : " << isoVal_ << nl
188 << " points : " << points().size() << nl
189 << " tris : " << triSurface::size() << nl
190 << " cut cells : " << meshCells_.size() << endl;
197 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
199 Foam::sampledIsoSurfaceCell::sampledIsoSurfaceCell
202 const polyMesh& mesh,
203 const dictionary& dict
206 sampledSurface(name, mesh, dict),
207 isoField_(dict.lookup("isoField")),
208 isoVal_(readScalar(dict.lookup("isoValue"))),
209 regularise_(dict.lookupOrDefault("regularise", true)),
210 average_(dict.lookupOrDefault("average", true)),
211 zoneKey_(keyType::null),
216 // dict.readIfPresent("zone", zoneKey_);
218 // if (debug && zoneKey_.size() && mesh.cellZones().findZoneID(zoneKey_) < 0)
220 // Info<< "cellZone " << zoneKey_
221 // << " not found - using entire mesh" << endl;
226 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
228 Foam::sampledIsoSurfaceCell::~sampledIsoSurfaceCell()
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 bool Foam::sampledIsoSurfaceCell::needsUpdate() const
236 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
238 return fvm.time().timeIndex() != prevTimeIndex_;
242 bool Foam::sampledIsoSurfaceCell::expire()
246 // Clear derived data
247 sampledSurface::clearGeom();
249 // already marked as expired
250 if (prevTimeIndex_ == -1)
261 bool Foam::sampledIsoSurfaceCell::update()
263 return updateGeometry();
267 Foam::tmp<Foam::scalarField>
268 Foam::sampledIsoSurfaceCell::sample
270 const volScalarField& vField
273 return sampleField(vField);
277 Foam::tmp<Foam::vectorField>
278 Foam::sampledIsoSurfaceCell::sample
280 const volVectorField& vField
283 return sampleField(vField);
287 Foam::tmp<Foam::sphericalTensorField>
288 Foam::sampledIsoSurfaceCell::sample
290 const volSphericalTensorField& vField
293 return sampleField(vField);
297 Foam::tmp<Foam::symmTensorField>
298 Foam::sampledIsoSurfaceCell::sample
300 const volSymmTensorField& vField
303 return sampleField(vField);
307 Foam::tmp<Foam::tensorField>
308 Foam::sampledIsoSurfaceCell::sample
310 const volTensorField& vField
313 return sampleField(vField);
317 Foam::tmp<Foam::scalarField>
318 Foam::sampledIsoSurfaceCell::interpolate
320 const interpolation<scalar>& interpolator
323 return interpolateField(interpolator);
327 Foam::tmp<Foam::vectorField>
328 Foam::sampledIsoSurfaceCell::interpolate
330 const interpolation<vector>& interpolator
333 return interpolateField(interpolator);
336 Foam::tmp<Foam::sphericalTensorField>
337 Foam::sampledIsoSurfaceCell::interpolate
339 const interpolation<sphericalTensor>& interpolator
342 return interpolateField(interpolator);
346 Foam::tmp<Foam::symmTensorField>
347 Foam::sampledIsoSurfaceCell::interpolate
349 const interpolation<symmTensor>& interpolator
352 return interpolateField(interpolator);
356 Foam::tmp<Foam::tensorField>
357 Foam::sampledIsoSurfaceCell::interpolate
359 const interpolation<tensor>& interpolator
362 return interpolateField(interpolator);
366 void Foam::sampledIsoSurfaceCell::print(Ostream& os) const
368 os << "sampledIsoSurfaceCell: " << name() << " :"
369 << " field:" << isoField_
370 << " value:" << isoVal_;
371 //<< " faces:" << faces().size() // possibly no geom yet
372 //<< " points:" << points().size();
376 // ************************************************************************* //