initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / fields / fvPatchFields / constraint / jumpCyclic / jumpCyclicFvPatchField.C
blob75dc758026903bbce67229f7c5a11d61d6510fd8
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 "jumpCyclicFvPatchField.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
34 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
36 template<class Type>
37 jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
39     const fvPatch& p,
40     const DimensionedField<Type, volMesh>& iF
43     cyclicFvPatchField<Type>(p, iF)
47 template<class Type>
48 jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
50     const jumpCyclicFvPatchField<Type>& ptf,
51     const fvPatch& p,
52     const DimensionedField<Type, volMesh>& iF,
53     const fvPatchFieldMapper& mapper
56     cyclicFvPatchField<Type>(ptf, p, iF, mapper)
60 template<class Type>
61 jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
63     const fvPatch& p,
64     const DimensionedField<Type, volMesh>& iF,
65     const dictionary& dict
68     cyclicFvPatchField<Type>(p, iF, dict)
70     // Call this evaluation in derived classes
71     //this->evaluate(Pstream::blocking);
75 template<class Type>
76 jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
78     const jumpCyclicFvPatchField<Type>& ptf
81     cyclicLduInterfaceField(),
82     cyclicFvPatchField<Type>(ptf)
86 template<class Type>
87 jumpCyclicFvPatchField<Type>::jumpCyclicFvPatchField
89     const jumpCyclicFvPatchField<Type>& ptf,
90     const DimensionedField<Type, volMesh>& iF
93     cyclicFvPatchField<Type>(ptf, iF)
97 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
99 template<class Type>
100 tmp<Field<Type> > jumpCyclicFvPatchField<Type>::patchNeighbourField() const
102     const Field<Type>& iField = this->internalField();
103     const unallocLabelList& faceCells = this->cyclicPatch().faceCells();
105     tmp<Field<Type> > tpnf(new Field<Type>(this->size()));
106     Field<Type>& pnf = tpnf();
108     tmp<Field<scalar> > tjf = jump();
109     const Field<scalar>& jf = tjf();
111     label sizeby2 = this->size()/2;
113     if (this->doTransform())
114     {
115         for (label facei=0; facei<sizeby2; facei++)
116         {
117             pnf[facei] = transform
118             (
119                 this->forwardT()[0], iField[faceCells[facei + sizeby2]]
120             ) - jf[facei];
122             pnf[facei + sizeby2] = transform
123             (
124                 this->reverseT()[0], iField[faceCells[facei]] + jf[facei]
125             );
126         }
127     }
128     else
129     {
130         for (label facei=0; facei<sizeby2; facei++)
131         {
132             pnf[facei] = iField[faceCells[facei + sizeby2]] - jf[facei];
133             pnf[facei + sizeby2] = iField[faceCells[facei]] + jf[facei];
134         }
135     }
137     return tpnf;
141 template<class Type>
142 void jumpCyclicFvPatchField<Type>::updateInterfaceMatrix
144     const scalarField& psiInternal,
145     scalarField& result,
146     const lduMatrix&,
147     const scalarField& coeffs,
148     const direction cmpt,
149     const Pstream::commsTypes
150 ) const
152     scalarField pnf(this->size());
154     label sizeby2 = this->size()/2;
155     const unallocLabelList& faceCells = this->cyclicPatch().faceCells();
157     if (long(&psiInternal) == long(&this->internalField()))
158     {
159         tmp<Field<scalar> > tjf = jump();
160         const Field<scalar>& jf = tjf();
162         for (label facei=0; facei<sizeby2; facei++)
163         {
164             pnf[facei] = psiInternal[faceCells[facei + sizeby2]] - jf[facei];
165             pnf[facei + sizeby2] = psiInternal[faceCells[facei]] + jf[facei];
166         }
167     }
168     else
169     {
170         for (label facei=0; facei<sizeby2; facei++)
171         {
172             pnf[facei] = psiInternal[faceCells[facei + sizeby2]];
173             pnf[facei + sizeby2] = psiInternal[faceCells[facei]];
174         }
175     }
177     // Transform according to the transformation tensors
178     this->transformCoupleField(pnf, cmpt);
180     // Multiply the field by coefficients and add into the result
181     forAll(faceCells, elemI)
182     {
183         result[faceCells[elemI]] -= coeffs[elemI]*pnf[elemI];
184     }
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 } // End namespace Foam
192 // ************************************************************************* //