BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / OpenFOAM / fields / Fields / complexFields / complexFields.C
blobf10ed55e861234182b7ca5a48ae9bfe4d2c77fb5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
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
19     for more details.
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 Description
25     Specialisation of Field\<T\> for complex and complexVector.
27 \*---------------------------------------------------------------------------*/
29 #include "complexFields.H"
30 #include "addToRunTimeSelectionTable.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 defineCompoundTypeName(List<complex>, complexList);
40 addCompoundToRunTimeSelectionTable(List<complex>, complexList);
42 complexField ComplexField(const UList<scalar>& re, const UList<scalar>& im)
44     complexField cf(re.size());
46     forAll(cf, i)
47     {
48         cf[i].Re() = re[i];
49         cf[i].Im() = im[i];
50     }
52     return cf;
56 complexField ReComplexField(const UList<scalar>& sf)
58     complexField cf(sf.size());
60     forAll(cf, i)
61     {
62         cf[i].Re() = sf[i];
63         cf[i].Im() = 0.0;
64     }
66     return cf;
70 complexField ImComplexField(const UList<scalar>& sf)
72     complexField cf(sf.size());
74     forAll(cf, i)
75     {
76         cf[i].Re() = 0.0;
77         cf[i].Im() = sf[i];
78     }
80     return cf;
84 scalarField ReImSum(const UList<complex>& cf)
86     scalarField sf(cf.size());
88     forAll(sf, i)
89     {
90         sf[i] = cf[i].Re() + cf[i].Im();
91     }
93     return sf;
97 scalarField Re(const UList<complex>& cf)
99     scalarField sf(cf.size());
101     forAll(sf, i)
102     {
103         sf[i] = cf[i].Re();
104     }
106     return sf;
110 scalarField Im(const UList<complex>& cf)
112     scalarField sf(cf.size());
114     forAll(sf, i)
115     {
116         sf[i] = cf[i].Im();
117     }
119     return sf;
123 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 defineCompoundTypeName(List<complexVector>, complexVectorList);
126 addCompoundToRunTimeSelectionTable(List<complexVector>, complexVectorList);
128 complexVectorField ComplexField
130     const UList<vector>& re,
131     const UList<vector>& im
134     complexVectorField cvf(re.size());
136     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
137     {
138         forAll(cvf, i)
139         {
140             cvf[i].component(cmpt).Re() = re[i].component(cmpt);
141             cvf[i].component(cmpt).Im() = im[i].component(cmpt);
142         }
143     }
145     return cvf;
149 complexVectorField ReComplexField(const UList<vector>& vf)
151     complexVectorField cvf(vf.size());
153     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
154     {
155         forAll(cvf, i)
156         {
157             cvf[i].component(cmpt).Re() = vf[i].component(cmpt);
158             cvf[i].component(cmpt).Im() = 0.0;
159         }
160     }
162     return cvf;
166 complexVectorField ImComplexField(const UList<vector>& vf)
168     complexVectorField cvf(vf.size());
170     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
171     {
172         forAll(cvf, i)
173         {
174             cvf[i].component(cmpt).Re() = 0.0;
175             cvf[i].component(cmpt).Im() = vf[i].component(cmpt);
176         }
177     }
179     return cvf;
183 vectorField ReImSum(const UList<complexVector>& cvf)
185     vectorField vf(cvf.size());
187     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
188     {
189         forAll(cvf, i)
190         {
191             vf[i].component(cmpt) =
192                 cvf[i].component(cmpt).Re() + cvf[i].component(cmpt).Im();
193         }
194     }
196     return vf;
200 vectorField Re(const UList<complexVector>& cvf)
202     vectorField vf(cvf.size());
204     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
205     {
206         forAll(cvf, i)
207         {
208             vf[i].component(cmpt) = cvf[i].component(cmpt).Re();
209         }
210     }
212     return vf;
216 vectorField Im(const UList<complexVector>& cvf)
218     vectorField vf(cvf.size());
220     for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
221     {
222         forAll(cvf, i)
223         {
224             vf[i].component(cmpt) = cvf[i].component(cmpt).Im();
225         }
226     }
228     return vf;
232 complexVectorField operator^
234     const UList<vector>& vf,
235     const UList<complexVector>& cvf
238     return ComplexField(vf^Re(cvf), vf^Im(cvf));
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244 } // End namespace Foam
246 // ************************************************************************* //