pointPatchField.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) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
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 
26 #include "pointPatchField.H"
27 #include "pointMesh.H"
28 #include "dictionary.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
37 template<class Type>
39 (
40  const pointPatch& p,
41  const DimensionedField<Type, pointMesh>& iF
42 )
43 :
44  patch_(p),
45  internalField_(iF),
46  updated_(false),
47  patchType_(word::null)
48 {}
49 
50 
51 template<class Type>
53 (
54  const pointPatch& p,
55  const DimensionedField<Type, pointMesh>& iF,
56  const dictionary& dict
57 )
58 :
59  patch_(p),
60  internalField_(iF),
61  updated_(false),
62  patchType_(dict.lookupOrDefault<word>("patchType", word::null))
63 {}
64 
65 
66 template<class Type>
68 (
69  const pointPatchField<Type>& ptf,
70  const pointPatch& p,
71  const DimensionedField<Type, pointMesh>& iF,
72  const pointPatchFieldMapper&
73 )
74 :
75  patch_(p),
76  internalField_(iF),
77  updated_(false),
78  patchType_(ptf.patchType_)
79 {}
80 
81 
82 template<class Type>
84 (
85  const pointPatchField<Type>& ptf
86 )
87 :
88  patch_(ptf.patch_),
89  internalField_(ptf.internalField_),
90  updated_(false),
91  patchType_(ptf.patchType_)
92 {}
93 
94 
95 template<class Type>
97 (
98  const pointPatchField<Type>& ptf,
99  const DimensionedField<Type, pointMesh>& iF
100 )
101 :
102  patch_(ptf.patch_),
103  internalField_(iF),
104  updated_(false),
105  patchType_(ptf.patchType_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
111 template<class Type>
113 {
114  return patch_.boundaryMesh().mesh()();
115 }
116 
117 
118 template<class Type>
120 {
121  os.writeKeyword("type") << type() << token::END_STATEMENT << nl;
122 
123  if (patchType_.size())
124  {
125  os.writeKeyword("patchType") << patchType_
126  << token::END_STATEMENT << nl;
127  }
128 }
129 
130 
131 template<class Type>
133 {
134  return patchInternalField(internalField());
135 }
136 
137 
138 template<class Type>
139 template<class Type1>
141 (
142  const Field<Type1>& iF,
143  const labelList& meshPoints
144 ) const
145 {
146  // Check size
147  if (iF.size() != internalField().size())
148  {
150  << "given internal field does not correspond to the mesh. "
151  << "Field size: " << iF.size()
152  << " mesh size: " << internalField().size()
153  << abort(FatalError);
154  }
155 
156  return tmp<Field<Type1> >(new Field<Type1>(iF, meshPoints));
157 }
158 
159 
160 template<class Type>
161 template<class Type1>
163 (
164  const Field<Type1>& iF
165 ) const
166 {
167  return patchInternalField(iF, patch().meshPoints());
168 }
169 
170 
171 template<class Type>
172 template<class Type1>
174 (
175  Field<Type1>& iF,
176  const Field<Type1>& pF
177 ) const
178 {
179  // Check size
180  if (iF.size() != internalField().size())
181  {
183  << "given internal field does not correspond to the mesh. "
184  << "Field size: " << iF.size()
185  << " mesh size: " << internalField().size()
186  << abort(FatalError);
187  }
188 
189  if (pF.size() != size())
190  {
192  << "given patch field does not correspond to the mesh. "
193  << "Field size: " << pF.size()
194  << " mesh size: " << size()
195  << abort(FatalError);
196  }
197 
198  // Get the addressing
199  const labelList& mp = patch().meshPoints();
200 
201  forAll(mp, pointI)
202  {
203  iF[mp[pointI]] += pF[pointI];
204  }
205 }
206 
207 
208 template<class Type>
209 template<class Type1>
211 (
212  Field<Type1>& iF,
213  const Field<Type1>& pF,
214  const labelList& points
215 ) const
216 {
217  // Check size
218  if (iF.size() != internalField().size())
219  {
221  << "given internal field does not correspond to the mesh. "
222  << "Field size: " << iF.size()
223  << " mesh size: " << internalField().size()
224  << abort(FatalError);
225  }
226 
227  if (pF.size() != size())
228  {
230  << "given patch field does not correspond to the mesh. "
231  << "Field size: " << pF.size()
232  << " mesh size: " << size()
233  << abort(FatalError);
234  }
235 
236  // Get the addressing
237  const labelList& mp = patch().meshPoints();
238 
239  forAll(points, i)
240  {
241  label pointI = points[i];
242  iF[mp[pointI]] += pF[pointI];
243  }
244 }
245 
246 
247 template<class Type>
248 template<class Type1>
250 (
251  Field<Type1>& iF,
252  const Field<Type1>& pF,
253  const labelList& meshPoints
254 ) const
255 {
256  // Check size
257  if (iF.size() != internalField().size())
258  {
260  << "given internal field does not correspond to the mesh. "
261  << "Field size: " << iF.size()
262  << " mesh size: " << internalField().size()
263  << abort(FatalError);
264  }
265 
266  if (pF.size() != meshPoints.size())
267  {
269  << "given patch field does not correspond to the meshPoints. "
270  << "Field size: " << pF.size()
271  << " meshPoints size: " << size()
272  << abort(FatalError);
273  }
274 
275  forAll(meshPoints, pointI)
276  {
277  iF[meshPoints[pointI]] = pF[pointI];
278  }
279 }
280 
281 
282 template<class Type>
283 template<class Type1>
285 (
286  Field<Type1>& iF,
287  const Field<Type1>& pF
288 ) const
289 {
290  setInInternalField(iF, pF, patch().meshPoints());
291 }
292 
293 
294 template<class Type>
296 {
297  if (!updated_)
298  {
299  updateCoeffs();
300  }
301 
302  updated_ = false;
303 }
304 
305 
306 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
307 
308 template<class Type>
309 Ostream& operator<<
310 (
311  Ostream& os,
312  const pointPatchField<Type>& ptf
313 )
314 {
315  ptf.write(os);
316 
317  os.check("Ostream& operator<<(Ostream&, const pointPatchField<Type>&)");
318 
319  return os;
320 }
321 
322 
323 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
324 
325 } // End namespace Foam
326 
327 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
328 
329 #include "pointPatchFieldNew.C"
330 
331 // ************************************************************************* //
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
p
p
Definition: pEqn.H:62
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
pointPatchField.H
Foam::pointPatchField::pointPatchField
pointPatchField(const pointPatch &, const DimensionedField< Type, pointMesh > &)
Construct from patch and internal field.
pointPatchFieldNew.C
Foam::pointPatchField
Abstract base class for point-mesh patch fields.
Definition: pointMVCWeight.H:58
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Foam::FatalError
error FatalError
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
internalField
conserve internalField()+
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
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
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
points
const pointField & points
Definition: gmvOutputHeader.H:1
dictionary.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
write
Tcoeff write()
pointMesh.H