patchCorrectedInterpolation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(patchCorrectedInterpolation, 0);
34 
36  (
37  motionInterpolation,
38  patchCorrectedInterpolation,
39  Istream
40  );
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 (
48  Istream& entry
49 ) const
50 {
51  List<List<word> > patchGroupNames(entry);
52 
53  labelListList patchGroups(patchGroupNames.size());
54 
55  forAll(patchGroupNames, patchI)
56  {
57  patchGroups[patchI].resize(patchGroupNames[patchI].size());
58 
59  forAll(patchGroupNames[patchI], patchJ)
60  {
61  patchGroups[patchI][patchJ] =
62  mesh().boundaryMesh().findPatchID
63  (
64  patchGroupNames[patchI][patchJ]
65  );
66 
67  if (patchGroups[patchI][patchJ] == -1)
68  {
70  << "patch \"" << patchGroupNames[patchI][patchJ]
71  << "\" not found" << exit(FatalError);
72  }
73  }
74  }
75 
76  return patchGroups;
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const fvMesh& mesh,
85  Istream& entry
86 )
87 :
89  patchGroups_(getPatchGroups(entry))
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
100 
102 (
103  const volScalarField& cellDisplacement,
104  pointScalarField& pointDisplacement
105 ) const
106 {
107  interpolateType(cellDisplacement, pointDisplacement);
108 }
109 
110 
112 (
113  const volVectorField& cellDisplacement,
114  pointVectorField& pointDisplacement
115 ) const
116 {
117  interpolateType(cellDisplacement, pointDisplacement);
118 }
119 
120 
121 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:65
Foam::patchCorrectedInterpolation::~patchCorrectedInterpolation
virtual ~patchCorrectedInterpolation()
Destructor.
Definition: patchCorrectedInterpolation.C:95
Foam::patchCorrectedInterpolation::getPatchGroups
labelListList getPatchGroups(Istream &entry) const
Get patch groups from the input stream.
Definition: patchCorrectedInterpolation.C:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::patchCorrectedInterpolation::interpolate
virtual void interpolate(const volScalarField &, pointScalarField &) const
Interpolate the given scalar cell displacement.
Definition: patchCorrectedInterpolation.C:102
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::patchCorrectedInterpolation::patchCorrectedInterpolation
patchCorrectedInterpolation(const fvMesh &mesh, Istream &entry)
Construct from an fvMesh and an Istream.
Definition: patchCorrectedInterpolation.C:83
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::List::resize
void resize(const label)
Alias for setSize(const label)
Foam::FatalError
error FatalError
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::motionInterpolation
Base class for interpolation of cell displacement fields, generated by fvMotionSolvers,...
Definition: motionInterpolation.H:50
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
patchCorrectedInterpolation.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)