Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / OpenFOAM / dimensionSet / dimensionSet.C
blob8a440a74eb8497cf7725f89f4f69c556a64129d5
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 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
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 \*---------------------------------------------------------------------------*/
26 #include "dimensionSet.H"
27 #include "dimensionedScalar.H"
28 #include "OStringStream.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     for (int Dimension=0; Dimension<nDimensions; ++Dimension)
83     {
84         // ie, mag(exponents_[Dimension]) > smallExponent
85         if
86         (
87             exponents_[Dimension] > smallExponent
88          || exponents_[Dimension] < -smallExponent
89         )
90         {
91             return false;
92         }
93     }
95     return true;
99 void Foam::dimensionSet::reset(const dimensionSet& ds)
101     for (int Dimension=0; Dimension<nDimensions; ++Dimension)
102     {
103         exponents_[Dimension] = ds.exponents_[Dimension];
104     }
108 Foam::string Foam::dimensionSet::asText() const
110     OStringStream buf;
112     bool Dimensionless = true;
114     for (int Dimension=0; Dimension < dimensionSet::nDimensions-1; ++Dimension)
115     {
116         const scalar& expt = exponents_[Dimension];
118         if (expt < smallExponent && expt > -smallExponent)
119         {
120             continue;
121         }
123         if (Dimensionless)
124         {
125             Dimensionless = false;
126         }
127         else
128         {
129             buf << ' ';
130         }
132         // note: currently only handle SI
133         switch (Dimension)
134         {
135             case MASS:
136                 buf << "kg";
137                 break;
139             case LENGTH:
140                 buf << "m";
141                 break;
143             case TIME:
144                 buf << "s";
145                 break;
147             case TEMPERATURE:
148                 buf << "K";
149                 break;
151             case MOLES:
152                 buf << "mol";
153                 break;
155             case CURRENT:
156                 buf << "A";
157                 break;
159             case LUMINOUS_INTENSITY:
160                 buf << "Cd";
161                 break;
163             default:
164                 buf << "??";  // this shouldn't be - flag as being weird
165                 break;
166         }
168         if (expt != 1)
169         {
170             buf << '^' << expt;
171         }
172     }
174     if (Dimensionless)
175     {
176         return "none";
177     }
178     else
179     {
180         return buf.str();
181     }
185 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
187 Foam::scalar Foam::dimensionSet::operator[](const dimensionType type) const
189     return exponents_[type];
193 Foam::scalar& Foam::dimensionSet::operator[](const dimensionType type)
195     return exponents_[type];
199 bool Foam::dimensionSet::operator==(const dimensionSet& ds) const
201     for (int Dimension=0; Dimension < nDimensions; ++Dimension)
202     {
203         if
204         (
205             mag(exponents_[Dimension] - ds.exponents_[Dimension])
206           > smallExponent
207         )
208         {
209             return false;
210         }
211     }
213     return true;
217 bool Foam::dimensionSet::operator!=(const dimensionSet& ds) const
219     return !(operator==(ds));
223 bool Foam::dimensionSet::operator=(const dimensionSet& ds) const
225     if (dimensionSet::debug && *this != ds)
226     {
227         FatalErrorIn("dimensionSet::operator=(const dimensionSet&) const")
228             << "Different dimensions for =" << endl
229             << "     dimensions : " << *this << " = " << ds << endl
230             << abort(FatalError);
231     }
233     return true;
237 bool Foam::dimensionSet::operator+=(const dimensionSet& ds) const
239     if (dimensionSet::debug && *this != ds)
240     {
241         FatalErrorIn("dimensionSet::operator+=(const dimensionSet&) const")
242             << "Different dimensions for +=" << endl
243             << "     dimensions : " << *this << " = " << ds << endl
244             << abort(FatalError);
245     }
247     return true;
251 bool Foam::dimensionSet::operator-=(const dimensionSet& ds) const
253     if (dimensionSet::debug && *this != ds)
254     {
255         FatalErrorIn("dimensionSet::operator-=(const dimensionSet&) const")
256             << "Different dimensions for -=" << endl
257             << "     dimensions : " << *this << " = " << ds << endl
258             << abort(FatalError);
259     }
261     return true;
265 bool Foam::dimensionSet::operator*=(const dimensionSet& ds)
267     reset((*this)*ds);
269     return true;
273 bool Foam::dimensionSet::operator/=(const dimensionSet& ds)
275     reset((*this)/ds);
277     return true;
281 // * * * * * * * * * * * * * * * Friend functions * * * * * * * * * * * * * * //
283 Foam::dimensionSet Foam::max(const dimensionSet& ds1, const dimensionSet& ds2)
285     if (dimensionSet::debug && ds1 != ds2)
286     {
287         FatalErrorIn("max(const dimensionSet&, const dimensionSet&)")
288             << "Arguments of max have different dimensions" << endl
289             << "     dimensions : " << ds1 << " and " << ds2 << endl
290             << abort(FatalError);
291     }
293     return ds1;
297 Foam::dimensionSet Foam::min(const dimensionSet& ds1, const dimensionSet& ds2)
299     if (dimensionSet::debug && ds1 != ds2)
300     {
301         FatalErrorIn("min(const dimensionSet&, const dimensionSet&)")
302             << "Arguments of min have different dimensions" << endl
303             << "     dimensions : " << ds1 << " and " << ds2 << endl
304             << abort(FatalError);
305     }
307     return ds1;
311 Foam::dimensionSet Foam::cmptMultiply
313     const dimensionSet& ds1,
314     const dimensionSet& ds2
317     return ds1*ds2;
321 Foam::dimensionSet Foam::cmptDivide
323     const dimensionSet& ds1,
324     const dimensionSet& ds2
327     return ds1/ds2;
331 Foam::dimensionSet Foam::pow(const dimensionSet& ds, const scalar p)
333     dimensionSet dimPow
334     (
335         ds[dimensionSet::MASS]*p,
336         ds[dimensionSet::LENGTH]*p,
337         ds[dimensionSet::TIME]*p,
338         ds[dimensionSet::TEMPERATURE]*p,
339         ds[dimensionSet::MOLES]*p,
340         ds[dimensionSet::CURRENT]*p,
341         ds[dimensionSet::LUMINOUS_INTENSITY]*p
342     );
344     return dimPow;
348 Foam::dimensionSet Foam::pow
350     const dimensionSet& ds,
351     const dimensionedScalar& dS
354     if (dimensionSet::debug && !dS.dimensions().dimensionless())
355     {
356         FatalErrorIn("pow(const dimensionSet&, const dimensionedScalar&)")
357             << "Exponent of pow is not dimensionless"
358             << abort(FatalError);
359     }
361     dimensionSet dimPow
362     (
363         ds[dimensionSet::MASS]*dS.value(),
364         ds[dimensionSet::LENGTH]*dS.value(),
365         ds[dimensionSet::TIME]*dS.value(),
366         ds[dimensionSet::TEMPERATURE]*dS.value(),
367         ds[dimensionSet::MOLES]*dS.value(),
368         ds[dimensionSet::CURRENT]*dS.value(),
369         ds[dimensionSet::LUMINOUS_INTENSITY]*dS.value()
370     );
372     return dimPow;
376 Foam::dimensionSet Foam::pow
378     const dimensionedScalar& dS,
379     const dimensionSet& ds
382     if
383     (
384         dimensionSet::debug
385      && !dS.dimensions().dimensionless()
386      && !ds.dimensionless()
387     )
388     {
389         FatalErrorIn("pow(const dimensionedScalar&, const dimensionSet&)")
390             << "Argument or exponent of pow not dimensionless" << endl
391             << abort(FatalError);
392     }
394     return ds;
398 Foam::dimensionSet Foam::sqr(const dimensionSet& ds)
400     return pow(ds, 2);
404 Foam::dimensionSet Foam::pow3(const dimensionSet& ds)
406     return pow(ds, 3);
410 Foam::dimensionSet Foam::pow4(const dimensionSet& ds)
412     return pow(ds, 4);
416 Foam::dimensionSet Foam::pow5(const dimensionSet& ds)
418     return pow(ds, 5);
422 Foam::dimensionSet Foam::pow6(const dimensionSet& ds)
424     return pow(ds, 6);
428 Foam::dimensionSet Foam::pow025(const dimensionSet& ds)
430     return sqrt(sqrt(ds));
434 Foam::dimensionSet Foam::sqrt(const dimensionSet& ds)
436     return pow(ds, 0.5);
440 Foam::dimensionSet Foam::magSqr(const dimensionSet& ds)
442     return pow(ds, 2);
446 Foam::dimensionSet Foam::mag(const dimensionSet& ds)
448     return ds;
452 Foam::dimensionSet Foam::sign(const dimensionSet&)
454     return dimless;
458 Foam::dimensionSet Foam::pos(const dimensionSet&)
460     return dimless;
464 Foam::dimensionSet Foam::neg(const dimensionSet&)
466     return dimless;
470 Foam::dimensionSet Foam::inv(const dimensionSet& ds)
472     return dimless/ds;
476 Foam::dimensionSet Foam::trans(const dimensionSet& ds)
478     if (dimensionSet::debug && !ds.dimensionless())
479     {
480         FatalErrorIn("trans(const dimensionSet&)")
481             << "Argument of trancendental function not dimensionless"
482             << abort(FatalError);
483     }
485     return ds;
489 Foam::dimensionSet Foam::transform(const dimensionSet& ds)
491     return ds;
495 // * * * * * * * * * * * * * * * Friend Operators  * * * * * * * * * * * * * //
497 Foam::dimensionSet Foam::operator-(const dimensionSet& ds)
499     return ds;
503 Foam::dimensionSet Foam::operator+
505     const dimensionSet& ds1,
506     const dimensionSet& ds2
509     dimensionSet dimSum(ds1);
511     if (dimensionSet::debug && ds1 != ds2)
512     {
513         FatalErrorIn
514             ("operator+(const dimensionSet&, const dimensionSet&)")
515             << "LHS and RHS of + have different dimensions" << endl
516             << "     dimensions : " << ds1 << " + " << ds2 << endl
517             << abort(FatalError);
518     }
520     return dimSum;
524 Foam::dimensionSet Foam::operator-
526     const dimensionSet& ds1,
527     const dimensionSet& ds2
530     dimensionSet dimDifference(ds1);
532     if (dimensionSet::debug && ds1 != ds2)
533     {
534         FatalErrorIn
535             ("operator-(const dimensionSet&, const dimensionSet&)")
536             << "LHS and RHS of - have different dimensions" << endl
537             << "     dimensions : " << ds1 << " - " << ds2 << endl
538             << abort(FatalError);
539     }
541     return dimDifference;
545 Foam::dimensionSet Foam::operator*
547     const dimensionSet& ds1,
548     const dimensionSet& ds2
551     dimensionSet dimProduct(ds1);
553     for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
554     {
555         dimProduct.exponents_[Dimension] += ds2.exponents_[Dimension];
556     }
558     return dimProduct;
562 Foam::dimensionSet Foam::operator/
564     const dimensionSet& ds1,
565     const dimensionSet& ds2
568     dimensionSet dimQuotient(ds1);
570     for (int Dimension=0; Dimension<dimensionSet::nDimensions; Dimension++)
571     {
572         dimQuotient.exponents_[Dimension] -= ds2.exponents_[Dimension];
573     }
575     return dimQuotient;
579 Foam::dimensionSet Foam::operator&
581     const dimensionSet& ds1,
582     const dimensionSet& ds2
585     return ds1*ds2;
589 Foam::dimensionSet Foam::operator^
591     const dimensionSet& ds1,
592     const dimensionSet& ds2
595     return ds1*ds2;
599 Foam::dimensionSet Foam::operator&&
601     const dimensionSet& ds1,
602     const dimensionSet& ds2
605     return ds1*ds2;
609 // ************************************************************************* //