initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / OpenFOAM / fields / pointPatchFields / pointPatchField / pointPatchField.C
blobe697ccc20b44d9391f2d47528c4e27866481f6ee
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "pointPatchField.H"
28 #include "pointMesh.H"
29 #include "dictionary.H"
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace Foam
36 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
38 template<class Type>
39 pointPatchField<Type>::pointPatchField
41     const pointPatch& p,
42     const DimensionedField<Type, pointMesh>& iF
45     patch_(p),
46     internalField_(iF),
47     updated_(false),
48     patchType_(word::null)
52 template<class Type>
53 pointPatchField<Type>::pointPatchField
55     const pointPatch& p,
56     const DimensionedField<Type, pointMesh>& iF,
57     const dictionary& dict
60     patch_(p),
61     internalField_(iF),
62     updated_(false),
63     patchType_(dict.lookupOrDefault<word>("patchType", word::null))
67 template<class Type>
68 pointPatchField<Type>::pointPatchField
70     const pointPatchField<Type>& ptf
73     patch_(ptf.patch_),
74     internalField_(ptf.internalField_),
75     updated_(false),
76     patchType_(ptf.patchType_)
80 template<class Type>
81 pointPatchField<Type>::pointPatchField
83     const pointPatchField<Type>& ptf,
84     const DimensionedField<Type, pointMesh>& iF
87     patch_(ptf.patch_),
88     internalField_(iF),
89     updated_(false),
90     patchType_(ptf.patchType_)
94 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
96 template<class Type>
97 const objectRegistry& pointPatchField<Type>::db() const
99     return patch_.boundaryMesh().mesh()();
103 template<class Type>
104 void pointPatchField<Type>::write(Ostream& os) const
106     os.writeKeyword("type") << type() << token::END_STATEMENT << nl;
108     if (patchType_.size())
109     {
110         os.writeKeyword("patchType") << patchType_
111             << token::END_STATEMENT << nl;
112     }
116 template<class Type>
117 tmp<Field<Type> > pointPatchField<Type>::patchInternalField() const
119     return patchInternalField(internalField());
123 template<class Type>
124 template<class Type1>
125 tmp<Field<Type1> > pointPatchField<Type>::patchInternalField
127     const Field<Type1>& iF
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     // get addressing
145     const labelList& meshPoints = patch().meshPoints();
147     tmp<Field<Type1> > tvalues(new Field<Type1>(meshPoints.size()));
148     Field<Type1>& values = tvalues();
150     forAll (meshPoints, pointI)
151     {
152         values[pointI] = iF[meshPoints[pointI]];
153     }
155     return tvalues;
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>::setInInternalField
208     Field<Type1>& iF,
209     const Field<Type1>& pF,
210     const labelList& meshPoints
211 ) const
213     // Check size
214     if (iF.size() != internalField().size())
215     {
216         FatalErrorIn
217         (
218             "void pointPatchField<Type>::"
219             "setInInternalField("
220             "Field<Type1>& iF, const Field<Type1>& iF) 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() != meshPoints.size())
228     {
229         FatalErrorIn
230         (
231             "void pointPatchField<Type>::"
232             "setInInternalField("
233             "Field<Type1>& iF, const Field<Type1>& iF) const"
234         )   << "given patch field does not correspond to the meshPoints. "
235             << "Field size: " << pF.size()
236             << " meshPoints size: " << size()
237             << abort(FatalError);
238     }
240     forAll (meshPoints, pointI)
241     {
242         iF[meshPoints[pointI]] = pF[pointI];
243     }
247 template<class Type>
248 template<class Type1>
249 void pointPatchField<Type>::setInInternalField
251     Field<Type1>& iF,
252     const Field<Type1>& pF
253 ) const
255     setInInternalField(iF, pF, patch().meshPoints());
259 template<class Type>
260 void pointPatchField<Type>::evaluate(const Pstream::commsTypes)
262     if (!updated_)
263     {
264         updateCoeffs();
265     }
267     updated_ = false;
271 // * * * * * * * * * * * * * * * IOstream Operators  * * * * * * * * * * * * //
273 template<class Type>
274 Ostream& operator<<
276     Ostream& os,
277     const pointPatchField<Type>& ptf
280     ptf.write(os);
282     os.check("Ostream& operator<<(Ostream&, const pointPatchField<Type>&)");
284     return os;
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 } // End namespace Foam
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
294 #include "newPointPatchField.C"
296 // ************************************************************************* //