initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / OpenFOAM / dimensionSet / dimensionSet.C
blob9fba472d9801d92233fefc7e7fe941cf7d6e3f69
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 \*---------------------------------------------------------------------------*/
27 #include "dimensionSet.H"
28 #include "dimensionedScalar.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 defineTypeNameAndDebug(Foam::dimensionSet, 1);
33 const Foam::scalar Foam::dimensionSet::smallExponent = SMALL;
36 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
38 Foam::dimensionSet::dimensionSet
40     const scalar mass,
41     const scalar length,
42     const scalar time,
43     const scalar temperature,
44     const scalar moles,
45     const scalar current,
46     const scalar luminousIntensity
49     exponents_[MASS] = mass;
50     exponents_[LENGTH] = length;
51     exponents_[TIME] = time;
52     exponents_[TEMPERATURE] = temperature;
53     exponents_[MOLES] = moles;
54     exponents_[CURRENT] = current;
55     exponents_[LUMINOUS_INTENSITY] = luminousIntensity;
59 Foam::dimensionSet::dimensionSet
61     const scalar mass,
62     const scalar length,
63     const scalar time,
64     const scalar temperature,
65     const scalar moles
68     exponents_[MASS] = mass;
69     exponents_[LENGTH] = length;
70     exponents_[TIME] = time;
71     exponents_[TEMPERATURE] = temperature;
72     exponents_[MOLES] = moles;
73     exponents_[CURRENT] = 0;
74     exponents_[LUMINOUS_INTENSITY] = 0;
78 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
80 bool Foam::dimensionSet::dimensionless() const
82     bool Dimensionless = true;
84     for (int Dimension=0; Dimension<nDimensions; Dimension++)
85     {
86         Dimensionless = Dimensionless &&
87         (
88             exponents_[Dimension] < smallExponent
89          && exponents_[Dimension] > -smallExponent
90         );
91     }
93     return Dimensionless;
97 void Foam::dimensionSet::reset(const dimensionSet& ds)
99     for (int Dimension=0; Dimension<nDimensions; Dimension++)
100     {
101         exponents_[Dimension] = ds.exponents_[Dimension];
102     }
106 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
108 Foam::scalar Foam::dimensionSet::operator[](const dimensionType type) const
110     return exponents_[type];
113 Foam::scalar& Foam::dimensionSet::operator[](const dimensionType type)
115     return exponents_[type];
119 bool Foam::dimensionSet::operator==(const dimensionSet& ds) const
121     bool equall = true;
123     for (int Dimension=0; Dimension<nDimensions; Dimension++)
124     {
125         equall = equall &&
126             (mag(exponents_[Dimension] - ds.exponents_[Dimension])
127           < smallExponent);
128     }
130     return equall;
133 bool Foam::dimensionSet::operator!=(const dimensionSet& ds) const
135     return !(operator==(ds));
139 bool Foam::dimensionSet::operator=(const dimensionSet& ds) const
141     if (dimensionSet::debug && *this != ds)
142     {
143         FatalErrorIn("dimensionSet::operator=(const dimensionSet& ds) const")
144             << "Different dimensions for =" << endl
145             << "     dimensions : " << *this << " = " << ds << endl
146             << abort(FatalError);
147     }
149     return true;
153 bool Foam::dimensionSet::operator+=(const dimensionSet& ds) const
155     if (dimensionSet::debug && *this != ds)
156     {
157         FatalErrorIn("dimensionSet::operator+=(const dimensionSet& ds) const")
158             << "Different dimensions for +=" << endl
159             << "     dimensions : " << *this << " = " << ds << endl
160             << abort(FatalError);
161     }
163     return true;
166 bool Foam::dimensionSet::operator-=(const dimensionSet& ds) const
168     if (dimensionSet::debug && *this != ds)
169     {
170         FatalErrorIn("dimensionSet::operator-=(const dimensionSet& ds) const")
171             << "Different dimensions for -=" << endl
172             << "     dimensions : " << *this << " = " << ds << endl
173             << abort(FatalError);
174     }
176     return true;
179 bool Foam::dimensionSet::operator*=(const dimensionSet& ds)
181     reset((*this)*ds);
183     return true;
186 bool Foam::dimensionSet::operator/=(const dimensionSet& ds)
188     reset((*this)/ds);
190     return true;
194 // * * * * * * * * * * * * * * * Friend functions * * * * * * * * * * * * * * //
196 Foam::dimensionSet Foam::max(const dimensionSet& ds1, const dimensionSet& ds2)
198     if (dimensionSet::debug && ds1 != ds2)
199     {
200         FatalErrorIn("max(const dimensionSet& ds1, const dimensionSet& ds2)")
201             << "Arguments of max have different dimensions" << endl
202             << "     dimensions : " << ds1 << " and " << ds2 << endl
203             << abort(FatalError);
204     }
206     return ds1;
209 Foam::dimensionSet Foam::min(const dimensionSet& ds1, const dimensionSet& ds2)
211     if (dimensionSet::debug && ds1 != ds2)
212     {
213         FatalErrorIn("min(const dimensionSet& ds1, const dimensionSet& ds2)")
214             << "Arguments of min have different dimensions" << endl
215             << "     dimensions : " << ds1 << " and " << ds2 << endl
216             << abort(FatalError);
217     }
219     return ds1;
223 Foam::dimensionSet Foam::cmptMultiply
225     const dimensionSet& ds1,
226     const dimensionSet& ds2
229     return ds1*ds2;
233 Foam::dimensionSet Foam::cmptDivide
235     const dimensionSet& ds1,
236     const dimensionSet& ds2
239     return ds1/ds2;
243 Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
245     dimensionSet dimPow
246     (
247         ds[dimensionSet::MASS]*p,
248         ds[dimensionSet::LENGTH]*p,
249         ds[dimensionSet::TIME]*p,
250         ds[dimensionSet::TEMPERATURE]*p,
251         ds[dimensionSet::MOLES]*p,
252         ds[dimensionSet::CURRENT]*p,
253         ds[dimensionSet::LUMINOUS_INTENSITY]*p
254     );
256     return dimPow;
259 Foam::dimensionSet Foam::pow
261     const dimensionSet& ds,
262     const dimensionedScalar& dS
265     if (dimensionSet::debug && !dS.dimensions().dimensionless())
266     {
267         FatalErrorIn("pow(const dimensionSet& ds, const dimensionedScalar& dS)")
268             << "Exponent of pow are not dimensionless"
269             << abort(FatalError);
270     }
272     dimensionSet dimPow
273     (
274         ds[dimensionSet::MASS]*dS.value(),
275         ds[dimensionSet::LENGTH]*dS.value(),
276         ds[dimensionSet::TIME]*dS.value(),
277         ds[dimensionSet::TEMPERATURE]*dS.value(),
278         ds[dimensionSet::MOLES]*dS.value(),
279         ds[dimensionSet::CURRENT]*dS.value(),
280         ds[dimensionSet::LUMINOUS_INTENSITY]*dS.value()
281     );
283     return dimPow;
286 Foam::dimensionSet Foam::pow
288     const dimensionedScalar& dS,
289     const dimensionSet& ds
292     if
293     (
294         dimensionSet::debug
295      && !dS.dimensions().dimensionless()
296      && !ds.dimensionless())
297     {
298         FatalErrorIn("pow(const dimensionedScalar& dS, const dimensionSet& ds)")
299             << "Argument or exponent of pow not dimensionless" << endl
300             << abort(FatalError);
301     }
303     return ds;
307 Foam::dimensionSet Foam::sqr(const dimensionSet& ds)
309     return pow(ds, 2);
312 Foam::dimensionSet Foam::pow3(const dimensionSet& ds)
314     return pow(ds, 3);
317 Foam::dimensionSet Foam::pow4(const dimensionSet& ds)
319     return pow(ds, 4);
322 Foam::dimensionSet Foam::pow5(const dimensionSet& ds)
324     return pow(ds, 5);
327 Foam::dimensionSet Foam::pow6(const dimensionSet& ds)
329     return pow(ds, 6);
332 Foam::dimensionSet Foam::sqrt(const dimensionSet& ds)
334     return pow(ds, 0.5);
337 Foam::dimensionSet Foam::magSqr(const dimensionSet& ds)
339     return pow(ds, 2);
342 Foam::dimensionSet Foam::mag(const dimensionSet& ds)
344     return ds;
347 Foam::dimensionSet Foam::sign(const dimensionSet&)
349     return dimless;
352 Foam::dimensionSet Foam::pos(const dimensionSet&)
354     return dimless;
357 Foam::dimensionSet Foam::neg(const dimensionSet&)
359     return dimless;
362 Foam::dimensionSet Foam::inv(const dimensionSet& ds)
364     return dimless/ds;
367 Foam::dimensionSet Foam::trans(const dimensionSet& ds)
369     if (dimensionSet::debug && !ds.dimensionless())
370     {
371         FatalErrorIn("trans(const dimensionSet& ds)")
372             << "Argument of trancendental function not dimensionless"
373             << abort(FatalError);
374     }
376     return ds;
379 Foam::dimensionSet Foam::transform(const dimensionSet& ds)
381     return ds;
385 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
387 Foam::dimensionSet Foam::operator-(const dimensionSet& ds)
389     return ds;
392 Foam::dimensionSet Foam::operator+
394     const dimensionSet& ds1,
395     const dimensionSet& ds2
398     dimensionSet dimSum(ds1);
400     if (dimensionSet::debug && ds1 != ds2)
401     {
402         FatalErrorIn
403             ("operator+(const dimensionSet& ds1, const dimensionSet& ds2)")
404             << "LHS and RHS of + have different dimensions" << endl
405             << "     dimensions : " << ds1 << " + " << ds2 << endl
406             << abort(FatalError);
407     }
409     return dimSum;
412 Foam::dimensionSet Foam::operator-
414     const dimensionSet& ds1,
415     const dimensionSet& ds2
418     dimensionSet dimDifference(ds1);
420     if (dimensionSet::debug && ds1 != ds2)
421     {
422         FatalErrorIn
423             ("operator-(const dimensionSet& ds1, const dimensionSet& ds2)")
424             << "LHS and RHS of - have different dimensions" << endl
425             << "     dimensions : " << ds1 << " - " << ds2 << endl
426             << abort(FatalError);
427     }
429     return dimDifference;
432 Foam::dimensionSet Foam::operator*
434     const dimensionSet& ds1,
435     const dimensionSet& ds2
438     dimensionSet dimProduct(ds1);
440     for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
441     {
442         dimProduct.exponents_[Dimension] += ds2.exponents_[Dimension];
443     }
445     return dimProduct;
448 Foam::dimensionSet Foam::operator/
450     const dimensionSet& ds1,
451     const dimensionSet& ds2
454     dimensionSet dimQuotient(ds1);
456     for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
457     {
458         dimQuotient.exponents_[Dimension] -= ds2.exponents_[Dimension];
459     }
461     return dimQuotient;
465 Foam::dimensionSet Foam::operator&
467     const dimensionSet& ds1,
468     const dimensionSet& ds2
471     return ds1*ds2;
474 Foam::dimensionSet Foam::operator^
476     const dimensionSet& ds1,
477     const dimensionSet& ds2
480     return ds1*ds2;
483 Foam::dimensionSet Foam::operator&&
485     const dimensionSet& ds1,
486     const dimensionSet& ds2
489     return ds1*ds2;
493 // ************************************************************************* //