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 "sampledPatchInternalField.H"
27 #include "dictionary.H"
29 #include "polyPatch.H"
30 #include "volFields.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 defineTypeNameAndDebug(sampledPatchInternalField, 0);
39 addNamedToRunTimeSelectionTable
42 sampledPatchInternalField,
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 Foam::sampledPatchInternalField::sampledPatchInternalField
55 const dictionary& dict
58 sampledPatch(name, mesh, dict),
59 mappers_(patchIDs().size())
61 const scalar distance = readScalar(dict.lookup("distance"));
65 label patchI = patchIDs()[i];
69 new directMappedPatchBase
71 mesh.boundaryMesh()[patchI],
72 mesh.name(), // sampleRegion
73 directMappedPatchBase::NEARESTCELL, // sampleMode
74 word::null, // samplePatch
75 -distance // sample inside my domain
82 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
84 Foam::sampledPatchInternalField::~sampledPatchInternalField()
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 Foam::tmp<Foam::scalarField> Foam::sampledPatchInternalField::sample
92 const volScalarField& vField
95 return sampleField(vField);
99 Foam::tmp<Foam::vectorField> Foam::sampledPatchInternalField::sample
101 const volVectorField& vField
104 return sampleField(vField);
107 Foam::tmp<Foam::sphericalTensorField> Foam::sampledPatchInternalField::sample
109 const volSphericalTensorField& vField
112 return sampleField(vField);
116 Foam::tmp<Foam::symmTensorField> Foam::sampledPatchInternalField::sample
118 const volSymmTensorField& vField
121 return sampleField(vField);
125 Foam::tmp<Foam::tensorField> Foam::sampledPatchInternalField::sample
127 const volTensorField& vField
130 return sampleField(vField);
134 Foam::tmp<Foam::scalarField> Foam::sampledPatchInternalField::interpolate
136 const interpolation<scalar>& interpolator
139 return interpolateField(interpolator);
143 Foam::tmp<Foam::vectorField> Foam::sampledPatchInternalField::interpolate
145 const interpolation<vector>& interpolator
148 return interpolateField(interpolator);
151 Foam::tmp<Foam::sphericalTensorField>
152 Foam::sampledPatchInternalField::interpolate
154 const interpolation<sphericalTensor>& interpolator
157 return interpolateField(interpolator);
161 Foam::tmp<Foam::symmTensorField> Foam::sampledPatchInternalField::interpolate
163 const interpolation<symmTensor>& interpolator
166 return interpolateField(interpolator);
170 Foam::tmp<Foam::tensorField> Foam::sampledPatchInternalField::interpolate
172 const interpolation<tensor>& interpolator
175 return interpolateField(interpolator);
179 void Foam::sampledPatchInternalField::print(Ostream& os) const
181 os << "sampledPatchInternalField: " << name() << " :"
182 << " patches:" << patchNames()
183 << " faces:" << faces().size()
184 << " points:" << points().size();
188 // ************************************************************************* //