ENH: sampling: add offsetMode to patchInternalField sampleSurface
[OpenFOAM-2.0.x.git] / src / sampling / sampledSurface / sampledPatchInternalField / sampledPatchInternalField.C
blob970b1522e8a2d22bb4e72803770815b8b002fcef
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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"
28 #include "polyMesh.H"
29 #include "polyPatch.H"
30 #include "volFields.H"
32 #include "addToRunTimeSelectionTable.H"
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 namespace Foam
38     defineTypeNameAndDebug(sampledPatchInternalField, 0);
39     addNamedToRunTimeSelectionTable
40     (
41         sampledSurface,
42         sampledPatchInternalField,
43         word,
44         patchInternalField
45     );
49 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
51 Foam::sampledPatchInternalField::sampledPatchInternalField
53     const word& name,
54     const polyMesh& mesh,
55     const dictionary& dict
58     sampledPatch(name, mesh, dict),
59     mappers_(patchIDs().size())
61     directMappedPatchBase::offsetMode mode = directMappedPatchBase::NORMAL;
62     if (dict.found("offsetMode"))
63     {
64         mode = directMappedPatchBase::offsetModeNames_.read
65         (
66             dict.lookup("offsetMode")
67         );
68     }
70     switch (mode)
71     {
72         case directMappedPatchBase::NORMAL:
73         {
74             const scalar distance = readScalar(dict.lookup("distance"));
75             forAll(patchIDs(), i)
76             {
77                 mappers_.set
78                 (
79                     i,
80                     new directMappedPatchBase
81                     (
82                         mesh.boundaryMesh()[patchIDs()[i]],
83                         mesh.name(),                        // sampleRegion
84                         directMappedPatchBase::NEARESTCELL, // sampleMode
85                         word::null,                         // samplePatch
86                         -distance                  // sample inside my domain
87                     )
88                 );
89             }
90         }
91         break;
93         case directMappedPatchBase::UNIFORM:
94         {
95             const point offset(dict.lookup("offset"));
96             forAll(patchIDs(), i)
97             {
98                 mappers_.set
99                 (
100                     i,
101                     new directMappedPatchBase
102                     (
103                         mesh.boundaryMesh()[patchIDs()[i]],
104                         mesh.name(),                        // sampleRegion
105                         directMappedPatchBase::NEARESTCELL, // sampleMode
106                         word::null,                         // samplePatch
107                         offset                  // sample inside my domain
108                     )
109                 );
110             }
111         }
112         break;
114         case directMappedPatchBase::NONUNIFORM:
115         {
116             const pointField offsets(dict.lookup("offsets"));
117             forAll(patchIDs(), i)
118             {
119                 mappers_.set
120                 (
121                     i,
122                     new directMappedPatchBase
123                     (
124                         mesh.boundaryMesh()[patchIDs()[i]],
125                         mesh.name(),                        // sampleRegion
126                         directMappedPatchBase::NEARESTCELL, // sampleMode
127                         word::null,                         // samplePatch
128                         offsets                  // sample inside my domain
129                     )
130                 );
131             }
132         }
133         break;
134     }
138 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
140 Foam::sampledPatchInternalField::~sampledPatchInternalField()
144 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
146 Foam::tmp<Foam::scalarField> Foam::sampledPatchInternalField::sample
148     const volScalarField& vField
149 ) const
151     return sampleField(vField);
155 Foam::tmp<Foam::vectorField> Foam::sampledPatchInternalField::sample
157     const volVectorField& vField
158 ) const
160     return sampleField(vField);
163 Foam::tmp<Foam::sphericalTensorField> Foam::sampledPatchInternalField::sample
165     const volSphericalTensorField& vField
166 ) const
168     return sampleField(vField);
172 Foam::tmp<Foam::symmTensorField> Foam::sampledPatchInternalField::sample
174     const volSymmTensorField& vField
175 ) const
177     return sampleField(vField);
181 Foam::tmp<Foam::tensorField> Foam::sampledPatchInternalField::sample
183     const volTensorField& vField
184 ) const
186     return sampleField(vField);
190 Foam::tmp<Foam::scalarField> Foam::sampledPatchInternalField::interpolate
192     const interpolation<scalar>& interpolator
193 ) const
195     return interpolateField(interpolator);
199 Foam::tmp<Foam::vectorField> Foam::sampledPatchInternalField::interpolate
201     const interpolation<vector>& interpolator
202 ) const
204     return interpolateField(interpolator);
207 Foam::tmp<Foam::sphericalTensorField>
208 Foam::sampledPatchInternalField::interpolate
210     const interpolation<sphericalTensor>& interpolator
211 ) const
213     return interpolateField(interpolator);
217 Foam::tmp<Foam::symmTensorField> Foam::sampledPatchInternalField::interpolate
219     const interpolation<symmTensor>& interpolator
220 ) const
222     return interpolateField(interpolator);
226 Foam::tmp<Foam::tensorField> Foam::sampledPatchInternalField::interpolate
228     const interpolation<tensor>& interpolator
229 ) const
231     return interpolateField(interpolator);
235 void Foam::sampledPatchInternalField::print(Ostream& os) const
237     os  << "sampledPatchInternalField: " << name() << " :"
238         << "  patches:" << patchNames()
239         << "  faces:" << faces().size()
240         << "  points:" << points().size();
244 // ************************************************************************* //