1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "sampledIsoSurfaceCell.H"
28 #include "dictionary.H"
29 #include "volFields.H"
30 #include "volPointInterpolation.H"
31 #include "addToRunTimeSelectionTable.H"
33 #include "isoSurfaceCell.H"
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 defineTypeNameAndDebug(sampledIsoSurfaceCell, 0);
40 addNamedToRunTimeSelectionTable(sampledSurface, sampledIsoSurfaceCell, word, isoSurfaceCell);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 bool Foam::sampledIsoSurfaceCell::updateGeometry() const
47 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
50 if (fvm.time().timeIndex() == prevTimeIndex_)
55 prevTimeIndex_ = fvm.time().timeIndex();
57 // Clear any stored topo
60 // Optionally read volScalarField
61 autoPtr<volScalarField> readFieldPtr_;
63 // 1. see if field in database
64 // 2. see if field can be read
65 const volScalarField* cellFldPtr = NULL;
66 if (fvm.foundObject<volScalarField>(isoField_))
70 Info<< "sampledIsoSurfaceCell::updateGeometry() : lookup "
74 cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
78 // Bit of a hack. Read field and store.
82 Info<< "sampledIsoSurfaceCell::updateGeometry() : reading "
83 << isoField_ << " from time " <<fvm.time().timeName()
94 fvm.time().timeName(),
104 cellFldPtr = readFieldPtr_.operator->();
106 const volScalarField& cellFld = *cellFldPtr;
108 tmp<pointScalarField> pointFld
110 volPointInterpolation::New(fvm).interpolate(cellFld)
115 //- From point field and interpolated cell.
116 scalarField cellAvg(fvm.nCells(), scalar(0.0));
117 labelField nPointCells(fvm.nCells(), 0);
119 for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
121 const labelList& pCells = fvm.pointCells(pointI);
125 label cellI = pCells[i];
127 cellAvg[cellI] += pointFld().internalField()[pointI];
128 nPointCells[cellI]++;
132 forAll(cellAvg, cellI)
134 cellAvg[cellI] /= nPointCells[cellI];
137 const isoSurfaceCell iso
141 pointFld().internalField(),
146 const_cast<sampledIsoSurfaceCell&>
149 ).triSurface::operator=(iso);
150 meshCells_ = iso.meshCells();
154 //- Direct from cell field and point field. Gives bad continuity.
155 const isoSurfaceCell iso
158 cellFld.internalField(),
159 pointFld().internalField(),
164 const_cast<sampledIsoSurfaceCell&>
167 ).triSurface::operator=(iso);
168 meshCells_ = iso.meshCells();
174 Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
176 << " regularise : " << regularise_ << nl
177 << " average : " << average_ << nl
178 << " isoField : " << isoField_ << nl
179 << " isoValue : " << isoVal_ << nl
180 << " points : " << points().size() << nl
181 << " tris : " << triSurface::size() << nl
182 << " cut cells : " << meshCells_.size() << endl;
189 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191 Foam::sampledIsoSurfaceCell::sampledIsoSurfaceCell
194 const polyMesh& mesh,
195 const dictionary& dict
198 sampledSurface(name, mesh, dict),
199 isoField_(dict.lookup("isoField")),
200 isoVal_(readScalar(dict.lookup("isoValue"))),
201 regularise_(dict.lookupOrDefault("regularise", true)),
202 average_(dict.lookupOrDefault("average", true)),
203 zoneName_(word::null),
208 // dict.readIfPresent("zone", zoneName_);
210 // if (debug && zoneName_.size())
212 // if (mesh.cellZones().findZoneID(zoneName_) < 0)
214 // Info<< "cellZone \"" << zoneName_
215 // << "\" not found - using entire mesh" << endl;
221 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
223 Foam::sampledIsoSurfaceCell::~sampledIsoSurfaceCell()
227 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
229 bool Foam::sampledIsoSurfaceCell::needsUpdate() const
231 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
233 return fvm.time().timeIndex() != prevTimeIndex_;
237 bool Foam::sampledIsoSurfaceCell::expire()
241 // already marked as expired
242 if (prevTimeIndex_ == -1)
253 bool Foam::sampledIsoSurfaceCell::update()
255 return updateGeometry();
259 Foam::tmp<Foam::scalarField>
260 Foam::sampledIsoSurfaceCell::sample
262 const volScalarField& vField
265 return sampleField(vField);
269 Foam::tmp<Foam::vectorField>
270 Foam::sampledIsoSurfaceCell::sample
272 const volVectorField& vField
275 return sampleField(vField);
279 Foam::tmp<Foam::sphericalTensorField>
280 Foam::sampledIsoSurfaceCell::sample
282 const volSphericalTensorField& vField
285 return sampleField(vField);
289 Foam::tmp<Foam::symmTensorField>
290 Foam::sampledIsoSurfaceCell::sample
292 const volSymmTensorField& vField
295 return sampleField(vField);
299 Foam::tmp<Foam::tensorField>
300 Foam::sampledIsoSurfaceCell::sample
302 const volTensorField& vField
305 return sampleField(vField);
309 Foam::tmp<Foam::scalarField>
310 Foam::sampledIsoSurfaceCell::interpolate
312 const interpolation<scalar>& interpolator
315 return interpolateField(interpolator);
319 Foam::tmp<Foam::vectorField>
320 Foam::sampledIsoSurfaceCell::interpolate
322 const interpolation<vector>& interpolator
325 return interpolateField(interpolator);
328 Foam::tmp<Foam::sphericalTensorField>
329 Foam::sampledIsoSurfaceCell::interpolate
331 const interpolation<sphericalTensor>& interpolator
334 return interpolateField(interpolator);
338 Foam::tmp<Foam::symmTensorField>
339 Foam::sampledIsoSurfaceCell::interpolate
341 const interpolation<symmTensor>& interpolator
344 return interpolateField(interpolator);
348 Foam::tmp<Foam::tensorField>
349 Foam::sampledIsoSurfaceCell::interpolate
351 const interpolation<tensor>& interpolator
354 return interpolateField(interpolator);
358 void Foam::sampledIsoSurfaceCell::print(Ostream& os) const
360 os << "sampledIsoSurfaceCell: " << name() << " :"
361 << " field:" << isoField_
362 << " value:" << isoVal_;
363 //<< " faces:" << faces().size() // possibly no geom yet
364 //<< " points:" << points().size();
368 // ************************************************************************* //