initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / OpenFOAM / fields / Fields / complexFields / complexFields.C
blob019505aecfd330084a512b14821fd702c908bc6c
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 Description
26     Specialisation of Field\<T\> for complex and complexVector.
28 \*---------------------------------------------------------------------------*/
30 #include "complexFields.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 defineCompoundTypeName(List<complex>, complexList);
41 addCompoundToRunTimeSelectionTable(List<complex>, complexList);
43 complexField ComplexField(const UList<scalar>& re, const UList<scalar>& im)
45     complexField cf(re.size());
47     forAll(cf, i)
48     {
49         cf[i].Re() = re[i];
50         cf[i].Im() = im[i];
51     }
53     return cf;
57 complexField ReComplexField(const UList<scalar>& sf)
59     complexField cf(sf.size());
61     forAll(cf, i)
62     {
63         cf[i].Re() = sf[i];
64         cf[i].Im() = 0.0;
65     }
67     return cf;
71 complexField ImComplexField(const UList<scalar>& sf)
73     complexField cf(sf.size());
75     forAll(cf, i)
76     {
77         cf[i].Re() = 0.0;
78         cf[i].Im() = sf[i];
79     }
81     return cf;
85 scalarField ReImSum(const UList<complex>& cf)
87     scalarField sf(cf.size());
89     forAll(sf, i)
90     {
91         sf[i] = cf[i].Re() + cf[i].Im();
92     }
94     return sf;
98 scalarField Re(const UList<complex>& cf)
100     scalarField sf(cf.size());
102     forAll(sf, i)
103     {
104         sf[i] = cf[i].Re();
105     }
107     return sf;
111 scalarField Im(const UList<complex>& cf)
113     scalarField sf(cf.size());
115     forAll(sf, i)
116     {
117         sf[i] = cf[i].Im();
118     }
120     return sf;
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 defineCompoundTypeName(List<complexVector>, complexVectorList);
127 addCompoundToRunTimeSelectionTable(List<complexVector>, complexVectorList);
129 complexVectorField ComplexField(const UList<vector>& re, const UList<vector>& im)
131     complexVectorField cvf(re.size());
133     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
134     {
135         forAll(cvf, i)
136         {
137             cvf[i].component(cmpt).Re() = re[i].component(cmpt);
138             cvf[i].component(cmpt).Im() = im[i].component(cmpt);
139         }
140     }
142     return cvf;
146 complexVectorField ReComplexField(const UList<vector>& vf)
148     complexVectorField cvf(vf.size());
150     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
151     {
152         forAll(cvf, i)
153         {
154             cvf[i].component(cmpt).Re() = vf[i].component(cmpt);
155             cvf[i].component(cmpt).Im() = 0.0;
156         }
157     }
159     return cvf;
163 complexVectorField ImComplexField(const UList<vector>& vf)
165     complexVectorField cvf(vf.size());
167     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
168     {
169         forAll(cvf, i)
170         {
171             cvf[i].component(cmpt).Re() = 0.0;
172             cvf[i].component(cmpt).Im() = vf[i].component(cmpt);
173         }
174     }
176     return cvf;
180 vectorField ReImSum(const UList<complexVector>& cvf)
182     vectorField vf(cvf.size());
184     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
185     {
186         forAll(cvf, i)
187         {
188             vf[i].component(cmpt) =
189                 cvf[i].component(cmpt).Re() + cvf[i].component(cmpt).Im();
190         }
191     }
193     return vf;
197 vectorField Re(const UList<complexVector>& cvf)
199     vectorField vf(cvf.size());
201     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
202     {
203         forAll(cvf, i)
204         {
205             vf[i].component(cmpt) = cvf[i].component(cmpt).Re();
206         }
207     }
209     return vf;
213 vectorField Im(const UList<complexVector>& cvf)
215     vectorField vf(cvf.size());
217     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
218     {
219         forAll(cvf, i)
220         {
221             vf[i].component(cmpt) = cvf[i].component(cmpt).Im();
222         }
223     }
225     return vf;
229 complexVectorField operator^
231     const UList<vector>& vf,
232     const UList<complexVector>& cvf
235     return ComplexField(vf^Re(cvf), vf^Im(cvf));
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 } // End namespace Foam
243 // ************************************************************************* //