1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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 "emptyPolyPatch.H"
27 #include "commSchedule.H"
28 #include "globalMeshData.H"
29 #include "cyclicPolyPatch.H"
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 template<class Type, template<class> class PatchField, class GeoMesh>
34 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
35 GeometricBoundaryField
37 const BoundaryMesh& bmesh,
38 const DimensionedField<Type, GeoMesh>& field,
39 const word& patchFieldType
42 FieldField<PatchField, Type>(bmesh.size()),
47 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
48 "GeometricBoundaryField::"
49 "GeometricBoundaryField(const BoundaryMesh&, "
50 "const Field<Type>&, const word&)"
54 forAll(bmesh_, patchi)
70 template<class Type, template<class> class PatchField, class GeoMesh>
71 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
72 GeometricBoundaryField
74 const BoundaryMesh& bmesh,
75 const DimensionedField<Type, GeoMesh>& field,
76 const wordList& patchFieldTypes,
77 const wordList& constraintTypes
80 FieldField<PatchField, Type>(bmesh.size()),
85 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
86 "GeometricBoundaryField::"
87 "GeometricBoundaryField(const BoundaryMesh&, "
88 "const Field<Type>&, const wordList&, const wordList&)"
94 patchFieldTypes.size() != this->size()
95 || (constraintTypes.size() && (constraintTypes.size() != this->size()))
100 "GeometricField<Type, PatchField, GeoMesh>::"
101 "GeometricBoundaryField::"
102 "GeometricBoundaryField(const BoundaryMesh&, "
103 "const Field<Type>&, const wordList&, const wordList&)"
104 ) << "Incorrect number of patch type specifications given" << nl
105 << " Number of patches in mesh = " << bmesh.size()
106 << " number of patch type specifications = "
107 << patchFieldTypes.size()
108 << abort(FatalError);
111 if (constraintTypes.size())
113 forAll(bmesh_, patchi)
118 PatchField<Type>::New
120 patchFieldTypes[patchi],
121 constraintTypes[patchi],
130 forAll(bmesh_, patchi)
135 PatchField<Type>::New
137 patchFieldTypes[patchi],
147 template<class Type, template<class> class PatchField, class GeoMesh>
148 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
149 GeometricBoundaryField
151 const BoundaryMesh& bmesh,
152 const DimensionedField<Type, GeoMesh>& field,
153 const PtrList<PatchField<Type> >& ptfl
156 FieldField<PatchField, Type>(bmesh.size()),
161 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
162 "GeometricBoundaryField::"
163 "GeometricBoundaryField(const BoundaryMesh&, "
164 "const Field<Type>&, const PatchField<Type>List&)"
168 forAll(bmesh_, patchi)
170 this->set(patchi, ptfl[patchi].clone(field));
175 template<class Type, template<class> class PatchField, class GeoMesh>
176 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
177 GeometricBoundaryField
179 const DimensionedField<Type, GeoMesh>& field,
180 const typename GeometricField<Type, PatchField, GeoMesh>::
181 GeometricBoundaryField& btf
184 FieldField<PatchField, Type>(btf.size()),
189 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
190 "GeometricBoundaryField::"
191 "GeometricBoundaryField(const GeometricBoundaryField<Type, "
192 "PatchField, BoundaryMesh>&)"
196 forAll(bmesh_, patchi)
198 this->set(patchi, btf[patchi].clone(field));
204 // Dangerous because Field may be set to a field which gets deleted.
205 // Need new type of GeometricBoundaryField, one which IS part of a geometric
206 // field for which snGrad etc. may be called and a free standing
207 // GeometricBoundaryField for which such operations are unavailable.
208 template<class Type, template<class> class PatchField, class GeoMesh>
209 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
210 GeometricBoundaryField
212 const typename GeometricField<Type, PatchField, GeoMesh>::
213 GeometricBoundaryField& btf
216 FieldField<PatchField, Type>(btf),
221 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
222 "GeometricBoundaryField::"
223 "GeometricBoundaryField(const GeometricBoundaryField<Type, "
224 "PatchField, BoundaryMesh>&)"
230 template<class Type, template<class> class PatchField, class GeoMesh>
231 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
232 GeometricBoundaryField
234 const BoundaryMesh& bmesh,
235 const DimensionedField<Type, GeoMesh>& field,
236 const dictionary& dict
239 FieldField<PatchField, Type>(bmesh.size()),
244 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
245 "GeometricBoundaryField::"
246 "GeometricBoundaryField"
247 "(const BoundaryMesh&, const Field<Type>&, const dictionary&)"
251 forAll(bmesh_, patchi)
253 if (bmesh_[patchi].type() != emptyPolyPatch::typeName)
257 bmesh_[patchi].type() == cyclicPolyPatch::typeName
258 && !dict.found(bmesh_[patchi].name())
263 "GeometricField<Type, PatchField, GeoMesh>::\n"
264 "GeometricBoundaryField::GeometricBoundaryField\n"
266 " const BoundaryMesh&,\n"
267 " const DimensionedField<Type, GeoMesh>&,\n"
268 " const dictionary&\n"
271 ) << "Cannot find patchField entry for cyclic "
272 << bmesh_[patchi].name() << endl
273 << "Is your field uptodate with split cyclics?" << endl
274 << "Run foamUpgradeCyclics to convert mesh and fields"
275 << " to split cyclics." << exit(FatalIOError);
281 PatchField<Type>::New
285 dict.subDict(bmesh_[patchi].name())
294 PatchField<Type>::New
296 emptyPolyPatch::typeName,
306 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
308 template<class Type, template<class> class PatchField, class GeoMesh>
309 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
314 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
315 "GeometricBoundaryField::"
316 "updateCoeffs()" << endl;
319 forAll(*this, patchi)
321 this->operator[](patchi).updateCoeffs();
326 template<class Type, template<class> class PatchField, class GeoMesh>
327 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
332 Info<< "GeometricField<Type, PatchField, GeoMesh>::"
333 "GeometricBoundaryField::"
334 "evaluate()" << endl;
339 Pstream::defaultCommsType == Pstream::blocking
340 || Pstream::defaultCommsType == Pstream::nonBlocking
343 forAll(*this, patchi)
345 this->operator[](patchi).initEvaluate(Pstream::defaultCommsType);
348 // Block for any outstanding requests
352 && Pstream::defaultCommsType == Pstream::nonBlocking
355 Pstream::waitRequests();
358 forAll(*this, patchi)
360 this->operator[](patchi).evaluate(Pstream::defaultCommsType);
363 else if (Pstream::defaultCommsType == Pstream::scheduled)
365 const lduSchedule& patchSchedule =
366 bmesh_.mesh().globalData().patchSchedule();
368 forAll(patchSchedule, patchEvali)
370 if (patchSchedule[patchEvali].init)
372 this->operator[](patchSchedule[patchEvali].patch)
373 .initEvaluate(Pstream::scheduled);
377 this->operator[](patchSchedule[patchEvali].patch)
378 .evaluate(Pstream::scheduled);
384 FatalErrorIn("GeometricBoundaryField::evaluate()")
385 << "Unsuported communications type "
386 << Pstream::commsTypeNames[Pstream::defaultCommsType]
392 template<class Type, template<class> class PatchField, class GeoMesh>
394 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
397 const FieldField<PatchField, Type>& pff = *this;
399 wordList Types(pff.size());
403 Types[patchi] = pff[patchi].type();
410 template<class Type, template<class> class PatchField, class GeoMesh>
411 typename Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField
412 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
413 boundaryInternalField() const
415 typename GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField
416 BoundaryInternalField(*this);
418 forAll(BoundaryInternalField, patchi)
420 BoundaryInternalField[patchi] ==
421 this->operator[](patchi).patchInternalField();
424 return BoundaryInternalField;
428 template<class Type, template<class> class PatchField, class GeoMesh>
429 Foam::lduInterfaceFieldPtrsList
430 Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
433 lduInterfaceFieldPtrsList interfaces(this->size());
435 forAll(interfaces, patchi)
437 if (isA<lduInterfaceField>(this->operator[](patchi)))
442 &refCast<const lduInterfaceField>(this->operator[](patchi))
451 template<class Type, template<class> class PatchField, class GeoMesh>
452 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
453 writeEntry(const word& keyword, Ostream& os) const
455 os << keyword << nl << token::BEGIN_BLOCK << incrIndent << nl;
457 forAll(*this, patchi)
459 os << indent << this->operator[](patchi).patch().name() << nl
460 << indent << token::BEGIN_BLOCK << nl
461 << incrIndent << this->operator[](patchi) << decrIndent
462 << indent << token::END_BLOCK << endl;
465 os << decrIndent << token::END_BLOCK << endl;
467 // Check state of IOstream
470 "GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::"
471 "writeEntry(const word& keyword, Ostream& os) const"
476 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
478 template<class Type, template<class> class PatchField, class GeoMesh>
479 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
482 const typename GeometricField<Type, PatchField, GeoMesh>::
483 GeometricBoundaryField& bf
486 FieldField<PatchField, Type>::operator=(bf);
490 template<class Type, template<class> class PatchField, class GeoMesh>
491 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
494 const FieldField<PatchField, Type>& ptff
497 FieldField<PatchField, Type>::operator=(ptff);
501 template<class Type, template<class> class PatchField, class GeoMesh>
502 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
508 FieldField<PatchField, Type>::operator=(t);
512 // Forced assignments
513 template<class Type, template<class> class PatchField, class GeoMesh>
514 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
517 const typename GeometricField<Type, PatchField, GeoMesh>::
518 GeometricBoundaryField& bf
521 forAll((*this), patchI)
523 this->operator[](patchI) == bf[patchI];
528 template<class Type, template<class> class PatchField, class GeoMesh>
529 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
532 const FieldField<PatchField, Type>& ptff
535 forAll((*this), patchI)
537 this->operator[](patchI) == ptff[patchI];
542 template<class Type, template<class> class PatchField, class GeoMesh>
543 void Foam::GeometricField<Type, PatchField, GeoMesh>::GeometricBoundaryField::
549 forAll((*this), patchI)
551 this->operator[](patchI) == t;
556 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
558 template<class Type, template<class> class PatchField, class GeoMesh>
559 Foam::Ostream& Foam::operator<<
562 const typename GeometricField<Type, PatchField, GeoMesh>::
563 GeometricBoundaryField& bf
566 os << static_cast<const FieldField<PatchField, Type>&>(bf);
571 // ************************************************************************* //