1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. 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(sampledSurface, sampledIsoSurfaceCell, word, isoSurfaceCell);
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 bool Foam::sampledIsoSurfaceCell::updateGeometry() const
46 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
49 if (fvm.time().timeIndex() == prevTimeIndex_)
54 prevTimeIndex_ = fvm.time().timeIndex();
56 // Clear any stored topo
59 // Optionally read volScalarField
60 autoPtr<volScalarField> readFieldPtr_;
62 // 1. see if field in database
63 // 2. see if field can be read
64 const volScalarField* cellFldPtr = NULL;
65 if (fvm.foundObject<volScalarField>(isoField_))
69 Info<< "sampledIsoSurfaceCell::updateGeometry() : lookup "
73 cellFldPtr = &fvm.lookupObject<volScalarField>(isoField_);
77 // Bit of a hack. Read field and store.
81 Info<< "sampledIsoSurfaceCell::updateGeometry() : reading "
82 << isoField_ << " from time " <<fvm.time().timeName()
93 fvm.time().timeName(),
103 cellFldPtr = readFieldPtr_.operator->();
105 const volScalarField& cellFld = *cellFldPtr;
107 tmp<pointScalarField> pointFld
109 volPointInterpolation::New(fvm).interpolate(cellFld)
114 //- From point field and interpolated cell.
115 scalarField cellAvg(fvm.nCells(), scalar(0.0));
116 labelField nPointCells(fvm.nCells(), 0);
118 for (label pointI = 0; pointI < fvm.nPoints(); pointI++)
120 const labelList& pCells = fvm.pointCells(pointI);
124 label cellI = pCells[i];
126 cellAvg[cellI] += pointFld().internalField()[pointI];
127 nPointCells[cellI]++;
131 forAll(cellAvg, cellI)
133 cellAvg[cellI] /= nPointCells[cellI];
136 const isoSurfaceCell iso
140 pointFld().internalField(),
145 const_cast<sampledIsoSurfaceCell&>
148 ).triSurface::operator=(iso);
149 meshCells_ = iso.meshCells();
153 //- Direct from cell field and point field. Gives bad continuity.
154 const isoSurfaceCell iso
157 cellFld.internalField(),
158 pointFld().internalField(),
163 const_cast<sampledIsoSurfaceCell&>
166 ).triSurface::operator=(iso);
167 meshCells_ = iso.meshCells();
173 Pout<< "sampledIsoSurfaceCell::updateGeometry() : constructed iso:"
175 << " regularise : " << regularise_ << nl
176 << " average : " << average_ << nl
177 << " isoField : " << isoField_ << nl
178 << " isoValue : " << isoVal_ << nl
179 << " points : " << points().size() << nl
180 << " tris : " << triSurface::size() << nl
181 << " cut cells : " << meshCells_.size() << endl;
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
190 Foam::sampledIsoSurfaceCell::sampledIsoSurfaceCell
193 const polyMesh& mesh,
194 const dictionary& dict
197 sampledSurface(name, mesh, dict),
198 isoField_(dict.lookup("isoField")),
199 isoVal_(readScalar(dict.lookup("isoValue"))),
200 regularise_(dict.lookupOrDefault("regularise", true)),
201 average_(dict.lookupOrDefault("average", true)),
202 zoneName_(word::null),
207 // dict.readIfPresent("zone", zoneName_);
209 // if (debug && zoneName_.size())
211 // if (mesh.cellZones().findZoneID(zoneName_) < 0)
213 // Info<< "cellZone \"" << zoneName_
214 // << "\" not found - using entire mesh" << endl;
220 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
222 Foam::sampledIsoSurfaceCell::~sampledIsoSurfaceCell()
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
228 bool Foam::sampledIsoSurfaceCell::needsUpdate() const
230 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
232 return fvm.time().timeIndex() != prevTimeIndex_;
236 bool Foam::sampledIsoSurfaceCell::expire()
240 // already marked as expired
241 if (prevTimeIndex_ == -1)
252 bool Foam::sampledIsoSurfaceCell::update()
254 return updateGeometry();
258 Foam::tmp<Foam::scalarField>
259 Foam::sampledIsoSurfaceCell::sample
261 const volScalarField& vField
264 return sampleField(vField);
268 Foam::tmp<Foam::vectorField>
269 Foam::sampledIsoSurfaceCell::sample
271 const volVectorField& vField
274 return sampleField(vField);
278 Foam::tmp<Foam::sphericalTensorField>
279 Foam::sampledIsoSurfaceCell::sample
281 const volSphericalTensorField& vField
284 return sampleField(vField);
288 Foam::tmp<Foam::symmTensorField>
289 Foam::sampledIsoSurfaceCell::sample
291 const volSymmTensorField& vField
294 return sampleField(vField);
298 Foam::tmp<Foam::tensorField>
299 Foam::sampledIsoSurfaceCell::sample
301 const volTensorField& vField
304 return sampleField(vField);
308 Foam::tmp<Foam::scalarField>
309 Foam::sampledIsoSurfaceCell::interpolate
311 const interpolation<scalar>& interpolator
314 return interpolateField(interpolator);
318 Foam::tmp<Foam::vectorField>
319 Foam::sampledIsoSurfaceCell::interpolate
321 const interpolation<vector>& interpolator
324 return interpolateField(interpolator);
327 Foam::tmp<Foam::sphericalTensorField>
328 Foam::sampledIsoSurfaceCell::interpolate
330 const interpolation<sphericalTensor>& interpolator
333 return interpolateField(interpolator);
337 Foam::tmp<Foam::symmTensorField>
338 Foam::sampledIsoSurfaceCell::interpolate
340 const interpolation<symmTensor>& interpolator
343 return interpolateField(interpolator);
347 Foam::tmp<Foam::tensorField>
348 Foam::sampledIsoSurfaceCell::interpolate
350 const interpolation<tensor>& interpolator
353 return interpolateField(interpolator);
357 void Foam::sampledIsoSurfaceCell::print(Ostream& os) const
359 os << "sampledIsoSurfaceCell: " << name() << " :"
360 << " field:" << isoField_
361 << " value:" << isoVal_;
362 //<< " faces:" << faces().size() // possibly no geom yet
363 //<< " points:" << points().size();
367 // ************************************************************************* //