00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "sampledCuttingPlane.H"
00027 #include <OpenFOAM/dictionary.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <finiteVolume/volPointInterpolation.H>
00030 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00031 #include <finiteVolume/fvMesh.H>
00032 #include <sampling/isoSurface.H>
00033
00034
00035
00036 namespace Foam
00037 {
00038 defineTypeNameAndDebug(sampledCuttingPlane, 0);
00039 addNamedToRunTimeSelectionTable
00040 (
00041 sampledSurface,
00042 sampledCuttingPlane,
00043 word,
00044 cuttingPlane
00045 );
00046 }
00047
00048
00049
00050 void Foam::sampledCuttingPlane::createGeometry()
00051 {
00052 if (debug)
00053 {
00054 Pout<< "sampledCuttingPlane::createGeometry :updating geometry."
00055 << endl;
00056 }
00057
00058
00059 facesPtr_.clear();
00060 isoSurfPtr_.ptr();
00061 pointDistance_.clear();
00062 cellDistancePtr_.clear();
00063
00064
00065 clearGeom();
00066
00067
00068 if (zoneID_.index() != -1 && !subMeshPtr_.valid())
00069 {
00070 const polyBoundaryMesh& patches = mesh().boundaryMesh();
00071
00072
00073 label exposedPatchI = patches.findPatchID(exposedPatchName_);
00074
00075 if (debug)
00076 {
00077 Info<< "Allocating subset of size "
00078 << mesh().cellZones()[zoneID_.index()].size()
00079 << " with exposed faces into patch "
00080 << patches[exposedPatchI].name() << endl;
00081 }
00082
00083 subMeshPtr_.reset
00084 (
00085 new fvMeshSubset(static_cast<const fvMesh&>(mesh()))
00086 );
00087 subMeshPtr_().setLargeCellSubset
00088 (
00089 labelHashSet(mesh().cellZones()[zoneID_.index()]),
00090 exposedPatchI
00091 );
00092 }
00093
00094
00095
00096 const fvMesh& fvm =
00097 (
00098 subMeshPtr_.valid()
00099 ? subMeshPtr_().subMesh()
00100 : static_cast<const fvMesh&>(mesh())
00101 );
00102
00103
00104
00105
00106
00107 cellDistancePtr_.reset
00108 (
00109 new volScalarField
00110 (
00111 IOobject
00112 (
00113 "cellDistance",
00114 fvm.time().timeName(),
00115 fvm.time(),
00116 IOobject::NO_READ,
00117 IOobject::NO_WRITE,
00118 false
00119 ),
00120 fvm,
00121 dimensionedScalar("zero", dimLength, 0)
00122 )
00123 );
00124 volScalarField& cellDistance = cellDistancePtr_();
00125
00126
00127 {
00128 const pointField& cc = fvm.cellCentres();
00129 scalarField& fld = cellDistance.internalField();
00130
00131 forAll(cc, i)
00132 {
00133
00134 fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
00135 }
00136 }
00137
00138
00139 {
00140 forAll(cellDistance.boundaryField(), patchI)
00141 {
00142 if
00143 (
00144 isA<emptyFvPatchScalarField>
00145 (
00146 cellDistance.boundaryField()[patchI]
00147 )
00148 )
00149 {
00150 cellDistance.boundaryField().set
00151 (
00152 patchI,
00153 new calculatedFvPatchScalarField
00154 (
00155 fvm.boundary()[patchI],
00156 cellDistance
00157 )
00158 );
00159
00160 const polyPatch& pp = fvm.boundary()[patchI].patch();
00161 pointField::subField cc = pp.patchSlice(fvm.faceCentres());
00162
00163 fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
00164 fld.setSize(pp.size());
00165 forAll(fld, i)
00166 {
00167 fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
00168 }
00169 }
00170 else
00171 {
00172 const pointField& cc = fvm.C().boundaryField()[patchI];
00173 fvPatchScalarField& fld = cellDistance.boundaryField()[patchI];
00174
00175 forAll(fld, i)
00176 {
00177 fld[i] = (cc[i] - plane_.refPoint()) & plane_.normal();
00178 }
00179 }
00180 }
00181 }
00182
00183
00184
00185
00186
00187
00188
00189 pointDistance_.setSize(fvm.nPoints());
00190 {
00191 const pointField& pts = fvm.points();
00192
00193 forAll(pointDistance_, i)
00194 {
00195 pointDistance_[i] = (pts[i] - plane_.refPoint()) & plane_.normal();
00196 }
00197 }
00198
00199
00200 if (debug)
00201 {
00202 Pout<< "Writing cell distance:" << cellDistance.objectPath() << endl;
00203 cellDistance.write();
00204 pointScalarField pDist
00205 (
00206 IOobject
00207 (
00208 "pointDistance",
00209 fvm.time().timeName(),
00210 fvm.time(),
00211 IOobject::NO_READ,
00212 IOobject::NO_WRITE,
00213 false
00214 ),
00215 pointMesh::New(fvm),
00216 dimensionedScalar("zero", dimLength, 0)
00217 );
00218 pDist.internalField() = pointDistance_;
00219
00220 Pout<< "Writing point distance:" << pDist.objectPath() << endl;
00221 pDist.write();
00222 }
00223
00224
00225
00226 isoSurfPtr_.reset
00227 (
00228 new isoSurface
00229 (
00230 cellDistance,
00231 pointDistance_,
00232 0.0,
00233 regularise_
00234 )
00235 );
00236
00237 if (debug)
00238 {
00239 print(Pout);
00240 Pout<< endl;
00241 }
00242 }
00243
00244
00245
00246
00247 Foam::sampledCuttingPlane::sampledCuttingPlane
00248 (
00249 const word& name,
00250 const polyMesh& mesh,
00251 const dictionary& dict
00252 )
00253 :
00254 sampledSurface(name, mesh, dict),
00255 plane_(dict),
00256 mergeTol_(dict.lookupOrDefault("mergeTol", 1E-6)),
00257 regularise_(dict.lookupOrDefault("regularise", true)),
00258 zoneID_(dict.lookupOrDefault("zone", word::null), mesh.cellZones()),
00259 exposedPatchName_(word::null),
00260 needsUpdate_(true),
00261 subMeshPtr_(NULL),
00262 cellDistancePtr_(NULL),
00263 isoSurfPtr_(NULL),
00264 facesPtr_(NULL)
00265 {
00266 if (zoneID_.index() != -1)
00267 {
00268 dict.lookup("exposedPatchName") >> exposedPatchName_;
00269
00270 if (mesh.boundaryMesh().findPatchID(exposedPatchName_) == -1)
00271 {
00272 FatalErrorIn
00273 (
00274 "sampledCuttingPlane::sampledCuttingPlane"
00275 "(const word&, const polyMesh&, const dictionary&)"
00276 ) << "Cannot find patch " << exposedPatchName_
00277 << " in which to put exposed faces." << endl
00278 << "Valid patches are " << mesh.boundaryMesh().names()
00279 << exit(FatalError);
00280 }
00281
00282 if (debug && zoneID_.index() != -1)
00283 {
00284 Info<< "Restricting to cellZone " << zoneID_.name()
00285 << " with exposed internal faces into patch "
00286 << exposedPatchName_ << endl;
00287 }
00288 }
00289 }
00290
00291
00292
00293
00294 Foam::sampledCuttingPlane::~sampledCuttingPlane()
00295 {}
00296
00297
00298
00299
00300 bool Foam::sampledCuttingPlane::needsUpdate() const
00301 {
00302 return needsUpdate_;
00303 }
00304
00305
00306 bool Foam::sampledCuttingPlane::expire()
00307 {
00308 if (debug)
00309 {
00310 Pout<< "sampledCuttingPlane::expire :"
00311 << " have-facesPtr_:" << facesPtr_.valid()
00312 << " needsUpdate_:" << needsUpdate_ << endl;
00313 }
00314
00315
00316 facesPtr_.clear();
00317
00318
00319 clearGeom();
00320
00321
00322 if (needsUpdate_)
00323 {
00324 return false;
00325 }
00326
00327 needsUpdate_ = true;
00328 return true;
00329 }
00330
00331
00332 bool Foam::sampledCuttingPlane::update()
00333 {
00334 if (debug)
00335 {
00336 Pout<< "sampledCuttingPlane::update :"
00337 << " have-facesPtr_:" << facesPtr_.valid()
00338 << " needsUpdate_:" << needsUpdate_ << endl;
00339 }
00340
00341 if (!needsUpdate_)
00342 {
00343 return false;
00344 }
00345
00346 createGeometry();
00347
00348 needsUpdate_ = false;
00349 return true;
00350 }
00351
00352
00353 Foam::tmp<Foam::scalarField>
00354 Foam::sampledCuttingPlane::sample
00355 (
00356 const volScalarField& vField
00357 ) const
00358 {
00359 return sampleField(vField);
00360 }
00361
00362
00363 Foam::tmp<Foam::vectorField>
00364 Foam::sampledCuttingPlane::sample
00365 (
00366 const volVectorField& vField
00367 ) const
00368 {
00369 return sampleField(vField);
00370 }
00371
00372
00373 Foam::tmp<Foam::sphericalTensorField>
00374 Foam::sampledCuttingPlane::sample
00375 (
00376 const volSphericalTensorField& vField
00377 ) const
00378 {
00379 return sampleField(vField);
00380 }
00381
00382
00383 Foam::tmp<Foam::symmTensorField>
00384 Foam::sampledCuttingPlane::sample
00385 (
00386 const volSymmTensorField& vField
00387 ) const
00388 {
00389 return sampleField(vField);
00390 }
00391
00392
00393 Foam::tmp<Foam::tensorField>
00394 Foam::sampledCuttingPlane::sample
00395 (
00396 const volTensorField& vField
00397 ) const
00398 {
00399 return sampleField(vField);
00400 }
00401
00402
00403 Foam::tmp<Foam::scalarField>
00404 Foam::sampledCuttingPlane::interpolate
00405 (
00406 const interpolation<scalar>& interpolator
00407 ) const
00408 {
00409 return interpolateField(interpolator);
00410 }
00411
00412
00413 Foam::tmp<Foam::vectorField>
00414 Foam::sampledCuttingPlane::interpolate
00415 (
00416 const interpolation<vector>& interpolator
00417 ) const
00418 {
00419 return interpolateField(interpolator);
00420 }
00421
00422 Foam::tmp<Foam::sphericalTensorField>
00423 Foam::sampledCuttingPlane::interpolate
00424 (
00425 const interpolation<sphericalTensor>& interpolator
00426 ) const
00427 {
00428 return interpolateField(interpolator);
00429 }
00430
00431
00432 Foam::tmp<Foam::symmTensorField>
00433 Foam::sampledCuttingPlane::interpolate
00434 (
00435 const interpolation<symmTensor>& interpolator
00436 ) const
00437 {
00438 return interpolateField(interpolator);
00439 }
00440
00441
00442 Foam::tmp<Foam::tensorField>
00443 Foam::sampledCuttingPlane::interpolate
00444 (
00445 const interpolation<tensor>& interpolator
00446 ) const
00447 {
00448 return interpolateField(interpolator);
00449 }
00450
00451
00452 void Foam::sampledCuttingPlane::print(Ostream& os) const
00453 {
00454 os << "sampledCuttingPlane: " << name() << " :"
00455 << " plane:" << plane_
00456 << " faces:" << faces().size()
00457 << " points:" << points().size();
00458 }
00459
00460
00461