Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / fields / pointPatchFields / pointPatchField / pointPatchField.C
blob8e7ca0e25041158383225e24d96cdc19ca54e11d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 OpenCFD Ltd.
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 "pointPatchField.H"
27 #include "pointMesh.H"
28 #include "dictionary.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
35 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
37 template<class Type>
38 pointPatchField<Type>::pointPatchField
40     const pointPatch& p,
41     const DimensionedField<Type, pointMesh>& iF
44     patch_(p),
45     internalField_(iF),
46     updated_(false),
47     patchType_(word::null)
51 template<class Type>
52 pointPatchField<Type>::pointPatchField
54     const pointPatch& p,
55     const DimensionedField<Type, pointMesh>& iF,
56     const dictionary& dict
59     patch_(p),
60     internalField_(iF),
61     updated_(false),
62     patchType_(dict.lookupOrDefault<word>("patchType", word::null))
66 template<class Type>
67 pointPatchField<Type>::pointPatchField
69     const pointPatchField<Type>& ptf
72     patch_(ptf.patch_),
73     internalField_(ptf.internalField_),
74     updated_(false),
75     patchType_(ptf.patchType_)
79 template<class Type>
80 pointPatchField<Type>::pointPatchField
82     const pointPatchField<Type>& ptf,
83     const DimensionedField<Type, pointMesh>& iF
86     patch_(ptf.patch_),
87     internalField_(iF),
88     updated_(false),
89     patchType_(ptf.patchType_)
93 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
95 template<class Type>
96 const objectRegistry& pointPatchField<Type>::db() const
98     return patch_.boundaryMesh().mesh()();
102 template<class Type>
103 void pointPatchField<Type>::write(Ostream& os) const
105     os.writeKeyword("type") << type() << token::END_STATEMENT << nl;
107     if (patchType_.size())
108     {
109         os.writeKeyword("patchType") << patchType_
110             << token::END_STATEMENT << nl;
111     }
115 template<class Type>
116 tmp<Field<Type> > pointPatchField<Type>::patchInternalField() const
118     return patchInternalField(internalField());
122 template<class Type>
123 template<class Type1>
124 tmp<Field<Type1> > pointPatchField<Type>::patchInternalField
126     const Field<Type1>& iF,
127     const labelList& meshPoints
128 ) const
130     // Check size
131     if (iF.size() != internalField().size())
132     {
133         FatalErrorIn
134         (
135             "tmp<Field<Type1> > pointPatchField<"
136             "Type>::"
137             "patchInternalField(const Field<Type1>& iF) const"
138         )   << "given internal field does not correspond to the mesh. "
139             << "Field size: " << iF.size()
140             << " mesh size: " << internalField().size()
141             << abort(FatalError);
142     }
144     return tmp<Field<Type1> >(new Field<Type1>(iF, meshPoints));
148 template<class Type>
149 template<class Type1>
150 tmp<Field<Type1> > pointPatchField<Type>::patchInternalField
152     const Field<Type1>& iF
153 ) const
155     return patchInternalField(iF, patch().meshPoints());
159 template<class Type>
160 template<class Type1>
161 void pointPatchField<Type>::addToInternalField
163     Field<Type1>& iF,
164     const Field<Type1>& pF
165 ) const
167     // Check size
168     if (iF.size() != internalField().size())
169     {
170         FatalErrorIn
171         (
172             "void pointPatchField<Type>::"
173             "addToInternalField("
174             "Field<Type1>& iF, const Field<Type1>& iF) const"
175         )   << "given internal field does not correspond to the mesh. "
176             << "Field size: " << iF.size()
177             << " mesh size: " << internalField().size()
178             << abort(FatalError);
179     }
181     if (pF.size() != size())
182     {
183         FatalErrorIn
184         (
185             "void pointPatchField<Type>::"
186             "addToInternalField("
187             "Field<Type1>& iF, const Field<Type1>& iF) const"
188         )   << "given patch field does not correspond to the mesh. "
189             << "Field size: " << pF.size()
190             << " mesh size: " << size()
191             << abort(FatalError);
192     }
194     // Get the addressing
195     const labelList& mp = patch().meshPoints();
197     forAll(mp, pointI)
198     {
199         iF[mp[pointI]] += pF[pointI];
200     }
204 template<class Type>
205 template<class Type1>
206 void pointPatchField<Type>::addToInternalField
208     Field<Type1>& iF,
209     const Field<Type1>& pF,
210     const labelList& points
211 ) const
213     // Check size
214     if (iF.size() != internalField().size())
215     {
216         FatalErrorIn
217         (
218             "void pointPatchField<Type>::"
219             "addToInternalField("
220             "Field<Type1>& iF, const Field<Type1>& iF, const labelList&) const"
221         )   << "given internal field does not correspond to the mesh. "
222             << "Field size: " << iF.size()
223             << " mesh size: " << internalField().size()
224             << abort(FatalError);
225     }
227     if (pF.size() != size())
228     {
229         FatalErrorIn
230         (
231             "void pointPatchField<Type>::"
232             "addToInternalField("
233             "Field<Type1>& iF, const Field<Type1>& iF, const labelList&) const"
234         )   << "given patch field does not correspond to the mesh. "
235             << "Field size: " << pF.size()
236             << " mesh size: " << size()
237             << abort(FatalError);
238     }
240     // Get the addressing
241     const labelList& mp = patch().meshPoints();
243     forAll(points, i)
244     {
245         label pointI = points[i];
246         iF[mp[pointI]] += pF[pointI];
247     }
251 template<class Type>
252 template<class Type1>
253 void pointPatchField<Type>::setInInternalField
255     Field<Type1>& iF,
256     const Field<Type1>& pF,
257     const labelList& meshPoints
258 ) const
260     // Check size
261     if (iF.size() != internalField().size())
262     {
263         FatalErrorIn
264         (
265             "void pointPatchField<Type>::"
266             "setInInternalField("
267             "Field<Type1>& iF, const Field<Type1>& iF) const"
268         )   << "given internal field does not correspond to the mesh. "
269             << "Field size: " << iF.size()
270             << " mesh size: " << internalField().size()
271             << abort(FatalError);
272     }
274     if (pF.size() != meshPoints.size())
275     {
276         FatalErrorIn
277         (
278             "void pointPatchField<Type>::"
279             "setInInternalField("
280             "Field<Type1>& iF, const Field<Type1>& iF) const"
281         )   << "given patch field does not correspond to the meshPoints. "
282             << "Field size: " << pF.size()
283             << " meshPoints size: " << size()
284             << abort(FatalError);
285     }
287     forAll(meshPoints, pointI)
288     {
289         iF[meshPoints[pointI]] = pF[pointI];
290     }
294 template<class Type>
295 template<class Type1>
296 void pointPatchField<Type>::setInInternalField
298     Field<Type1>& iF,
299     const Field<Type1>& pF
300 ) const
302     setInInternalField(iF, pF, patch().meshPoints());
306 template<class Type>
307 void pointPatchField<Type>::evaluate(const Pstream::commsTypes)
309     if (!updated_)
310     {
311         updateCoeffs();
312     }
314     updated_ = false;
318 // * * * * * * * * * * * * * * * IOstream Operators  * * * * * * * * * * * * //
320 template<class Type>
321 Ostream& operator<<
323     Ostream& os,
324     const pointPatchField<Type>& ptf
327     ptf.write(os);
329     os.check("Ostream& operator<<(Ostream&, const pointPatchField<Type>&)");
331     return os;
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
337 } // End namespace Foam
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 #include "pointPatchFieldNew.C"
343 // ************************************************************************* //