Fix tutorials: coupled/conjugateHeatFoam/conjugateCavity: fix Allrun file
[OpenFOAM-1.6-ext.git] / src / finiteVolume / fields / fvPatchFields / constraint / ggi / ggiFvPatchField.C
blob712f7abca8b6e10543e634c0f8d0eab8e74a185d
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 Author
26     Hrvoje Jasak, Wikki Ltd.  All rights reserved
28 Contributor
29     Martin Beaudoin, Hydro-Quebec, (2008)
31 Note on parallelisation
32     In order to handle parallelisation correctly, I need to rely on the fact
33     that all patches that require a global gather-scatter come before
34     processor patches.  In that case, the communication pattern
35     will be correct without intervention.  HJ, 6/Aug/2009
37 \*---------------------------------------------------------------------------*/
39 #include "ggiFvPatchField.H"
40 #include "symmTransformField.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 namespace Foam
47 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
49 template<class Type>
50 ggiFvPatchField<Type>::ggiFvPatchField
52     const fvPatch& p,
53     const DimensionedField<Type, volMesh>& iF
56     coupledFvPatchField<Type>(p, iF),
57     ggiPatch_(refCast<const ggiFvPatch>(p))
61 template<class Type>
62 ggiFvPatchField<Type>::ggiFvPatchField
64     const fvPatch& p,
65     const DimensionedField<Type, volMesh>& iF,
66     const dictionary& dict
69     coupledFvPatchField<Type>(p, iF, dict, false),
70     ggiPatch_(refCast<const ggiFvPatch>(p))
72     if (!isType<ggiFvPatch>(p))
73     {
74         FatalIOErrorIn
75         (
76             "ggiFvPatchField<Type>::ggiFvPatchField\n"
77             "(\n"
78             "    const fvPatch& p,\n"
79             "    const DimensionedField<Type, volMesh>& iF,\n"
80             "    const dictionary& dict\n"
81             ")\n",
82             dict
83         )   << "patch " << this->patch().index() << " not ggi type. "
84             << "Patch type = " << p.type()
85             << exit(FatalIOError);
86     }
88     if (!dict.found("value"))
89     {
90         // Grab the internal value for initialisation. (?) HJ, 27/Feb/2009
91         fvPatchField<Type>::operator=(this->patchInternalField()());
92     }
96 template<class Type>
97 ggiFvPatchField<Type>::ggiFvPatchField
99     const ggiFvPatchField<Type>& ptf,
100     const fvPatch& p,
101     const DimensionedField<Type, volMesh>& iF,
102     const fvPatchFieldMapper& mapper
105     coupledFvPatchField<Type>(ptf, p, iF, mapper),
106     ggiPatch_(refCast<const ggiFvPatch>(p))
108     if (!isType<ggiFvPatch>(this->patch()))
109     {
110         FatalErrorIn
111         (
112             "ggiFvPatchField<Type>::ggiFvPatchField\n"
113             "(\n"
114             "    const ggiFvPatchField<Type>& ptf,\n"
115             "    const fvPatch& p,\n"
116             "    const DimensionedField<Type, volMesh>& iF,\n"
117             "    const fvPatchFieldMapper& mapper\n"
118             ")\n"
119         )   << "Field type does not correspond to patch type for patch "
120             << this->patch().index() << "." << endl
121             << "Field type: " << typeName << endl
122             << "Patch type: " << this->patch().type()
123             << exit(FatalError);
124     }
128 template<class Type>
129 ggiFvPatchField<Type>::ggiFvPatchField
131     const ggiFvPatchField<Type>& ptf,
132     const DimensionedField<Type, volMesh>& iF
135     ggiLduInterfaceField(),
136     coupledFvPatchField<Type>(ptf, iF),
137     ggiPatch_(refCast<const ggiFvPatch>(ptf.patch()))
141 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
143 // Return neighbour field
144 template<class Type>
145 tmp<Field<Type> > ggiFvPatchField<Type>::patchNeighbourField() const
147     const Field<Type>& iField = this->internalField();
149     // Get shadow face-cells and assemble shadow field
150     const unallocLabelList& sfc = ggiPatch_.shadow().faceCells();
152     Field<Type> sField(sfc.size());
154     forAll (sField, i)
155     {
156         sField[i] = iField[sfc[i]];
157     }
159     tmp<Field<Type> > tpnf(ggiPatch_.interpolate(sField));
160     Field<Type>& pnf = tpnf();
162     if (ggiPatch_.bridgeOverlap())
163     {
164         // Symmetry treatment used for overlap
165         vectorField nHat = this->patch().nf();
167         // Use mirrored neighbour field for interpolation
168         // HJ, 21/Jan/2009
169         Field<Type> bridgeField =
170             transform(I - 2.0*sqr(nHat), this->patchInternalField());
172         ggiPatch_.bridge(bridgeField, pnf);
173     }
175     return tpnf;
179 template<class Type>
180 void ggiFvPatchField<Type>::initEvaluate
182     const Pstream::commsTypes commsType
185     if (!this->updated())
186     {
187         this->updateCoeffs();
188     }
190     Field<Type> pf
191     (
192         this->patch().weights()*this->patchInternalField()
193       + (1.0 - this->patch().weights())*this->patchNeighbourField()
194     );
196     if (ggiPatch_.bridgeOverlap())
197     {
198         // Symmetry treatment used for overlap
199         vectorField nHat = this->patch().nf();
201         Field<Type> bridgeField =
202         (
203             this->patchInternalField()
204           + transform(I - 2.0*sqr(nHat), this->patchInternalField())
205         )/2.0;
207         ggiPatch_.bridge(bridgeField, pf);
208     }
210     Field<Type>::operator=(pf);
214 template<class Type>
215 void ggiFvPatchField<Type>::evaluate
217     const Pstream::commsTypes
220     if (!this->updated())
221     {
222         this->updateCoeffs();
223     }
227 template<class Type>
228 void ggiFvPatchField<Type>::initInterfaceMatrixUpdate
230     const scalarField& psiInternal,
231     scalarField& result,
232     const lduMatrix&,
233     const scalarField& coeffs,
234     const direction cmpt,
235     const Pstream::commsTypes commsType
236 ) const
238     // Get shadow face-cells and assemble shadow field
239     const unallocLabelList& sfc = ggiPatch_.shadow().faceCells();
241     scalarField sField(sfc.size());
243     forAll (sField, i)
244     {
245         sField[i] = psiInternal[sfc[i]];
246     }
248     scalarField pnf = ggiPatch_.interpolate(sField);
250     // Multiply the field by coefficients and add into the result
251     const unallocLabelList& fc = ggiPatch_.faceCells();
253     forAll(fc, elemI)
254     {
255         result[fc[elemI]] -= coeffs[elemI]*pnf[elemI];
256     }
260 template<class Type>
261 void ggiFvPatchField<Type>::updateInterfaceMatrix
263     const scalarField& psiInternal,
264     scalarField& result,
265     const lduMatrix&,
266     const scalarField& coeffs,
267     const direction cmpt,
268     const Pstream::commsTypes
269 ) const
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 } // End namespace Foam
277 // ************************************************************************* //