FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

liquidMixture.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "liquidMixture.H"
00027 #include <OpenFOAM/dictionary.H>
00028 #include <specie/specie.H>
00029 
00030 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00031 
00032 const Foam::scalar Foam::liquidMixture::TrMax = 0.999;
00033 
00034 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00035 
00036 Foam::liquidMixture::liquidMixture
00037 (
00038     const dictionary& thermophysicalProperties
00039 )
00040 :
00041     components_(thermophysicalProperties.lookup("liquidComponents")),
00042     properties_(components_.size())
00043 {
00044     // use sub-dictionary "liquidProperties" if possible to avoid
00045     // collisions with identically named gas-phase entries
00046     // (eg, H2O liquid vs. gas)
00047     forAll(components_, i)
00048     {
00049         const dictionary* subDictPtr = thermophysicalProperties.subDictPtr
00050         (
00051             "liquidProperties"
00052         );
00053 
00054         if (subDictPtr)
00055         {
00056             properties_.set
00057             (
00058                 i,
00059                 liquid::New(subDictPtr->lookup(components_[i]))
00060             );
00061         }
00062         else
00063         {
00064             properties_.set
00065             (
00066                 i,
00067                 liquid::New(thermophysicalProperties.lookup(components_[i]))
00068             );
00069         }
00070     }
00071 }
00072 
00073 
00074 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
00075 
00076 Foam::autoPtr<Foam::liquidMixture> Foam::liquidMixture::New
00077 (
00078     const dictionary& thermophysicalProperties
00079 )
00080 {
00081     return autoPtr<liquidMixture>(new liquidMixture(thermophysicalProperties));
00082 }
00083 
00084 
00085 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00086 
00087 Foam::scalar Foam::liquidMixture::Tc
00088 (
00089     const scalarField& x
00090 ) const
00091 {
00092 
00093     scalar vTc = 0.0;
00094     scalar vc = 0.0;
00095 
00096     forAll(properties_, i)
00097     {
00098         scalar x1 = x[i]*properties_[i].Vc();
00099         vc += x1;
00100         vTc += x1*properties_[i].Tc();
00101     }
00102 
00103     return vTc/vc;
00104 }
00105 
00106 
00107 Foam::scalar Foam::liquidMixture::Tpc
00108 (
00109     const scalarField& x
00110 ) const
00111 {
00112     scalar Tpc = 0.0;
00113     forAll(properties_, i)
00114     {
00115         Tpc += x[i]*properties_[i].Tc();
00116     }
00117 
00118     return Tpc;
00119 }
00120 
00121 
00122 Foam::scalar Foam::liquidMixture::Ppc
00123 (
00124     const scalarField& x
00125 ) const
00126 {
00127     scalar Vc = 0.0;
00128     scalar Zc = 0.0;
00129     forAll(properties_, i)
00130     {
00131         Vc += x[i]*properties_[i].Vc();
00132         Zc += x[i]*properties_[i].Zc();
00133     }
00134 
00135     return specie::RR*Zc*Tpc(x)/Vc;
00136 }
00137 
00138 
00139 Foam::scalar Foam::liquidMixture::omega
00140 (
00141     const scalarField& x
00142 ) const
00143 {
00144     scalar omega = 0.0;
00145     forAll(properties_, i)
00146     {
00147         omega += x[i]*properties_[i].omega();
00148     }
00149 
00150     return omega;
00151 }
00152 
00153 
00154 Foam::scalarField Foam::liquidMixture::Xs
00155 (
00156     const scalar p,
00157     const scalar Tg,
00158     const scalar Tl,
00159     const scalarField& xg,
00160     const scalarField& xl
00161 ) const
00162 {
00163     scalarField xs(xl.size(), 0.0);
00164 
00165     // Raoult's Law
00166     forAll(xs, i)
00167     {
00168         scalar Ti = min(TrMax*properties_[i].Tc(), Tl);
00169         xs[i] = properties_[i].pv(p, Ti)*xl[i]/p;
00170     }
00171     return xs;
00172 }
00173 
00174 
00175 Foam::scalar Foam::liquidMixture::W
00176 (
00177     const scalarField& x
00178 ) const
00179 {
00180     scalar W = 0.0;
00181     forAll(properties_, i)
00182     {
00183         W += x[i]*properties_[i].W();
00184     }
00185 
00186     return W;
00187 }
00188 
00189 
00190 Foam::scalarField Foam::liquidMixture::Y
00191 (
00192     const scalarField& X
00193 ) const
00194 {
00195     scalarField Y(X/W(X));
00196 
00197     forAll(Y, i)
00198     {
00199         Y[i] *= properties_[i].W();
00200     }
00201 
00202     return Y;
00203 }
00204 
00205 
00206 Foam::scalarField Foam::liquidMixture::X
00207 (
00208     const scalarField& Y
00209 ) const
00210 {
00211     scalarField X(Y.size());
00212     scalar Winv = 0.0;
00213     forAll(X, i)
00214     {
00215         Winv += Y[i]/properties_[i].W();
00216         X[i] = Y[i]/properties_[i].W();
00217     }
00218 
00219     return X/Winv;
00220 }
00221 
00222 
00223 Foam::scalar Foam::liquidMixture::rho
00224 (
00225     const scalar p,
00226     const scalar T,
00227     const scalarField& x
00228 ) const
00229 {
00230     scalar v = 0.0;
00231 
00232     forAll(properties_, i)
00233     {
00234         if (x[i] > SMALL)
00235         {
00236             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00237             scalar rho = SMALL + properties_[i].rho(p, Ti);
00238             v += x[i]*properties_[i].W()/rho;
00239         }
00240     }
00241 
00242     return W(x)/v;
00243 }
00244 
00245 
00246 Foam::scalar Foam::liquidMixture::pv
00247 (
00248     const scalar p,
00249     const scalar T,
00250     const scalarField& x
00251 ) const
00252 {
00253     scalar pv = 0.0;
00254 
00255     forAll(properties_, i)
00256     {
00257         if (x[i] > SMALL)
00258         {
00259             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00260             pv += x[i]*properties_[i].pv(p, Ti)*properties_[i].W();
00261         }
00262     }
00263 
00264     return pv/W(x);
00265 }
00266 
00267 
00268 Foam::scalar Foam::liquidMixture::hl
00269 (
00270     const scalar p,
00271     const scalar T,
00272     const scalarField& x
00273 ) const
00274 {
00275     scalar hl = 0.0;
00276 
00277     forAll(properties_, i)
00278     {
00279         if (x[i] > SMALL)
00280         {
00281             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00282             hl += x[i]*properties_[i].hl(p, Ti)*properties_[i].W();
00283         }
00284     }
00285 
00286     return hl/W(x);
00287 }
00288 
00289 
00290 Foam::scalar Foam::liquidMixture::cp
00291 (
00292     const scalar p,
00293     const scalar T,
00294     const scalarField& x
00295 ) const
00296 {
00297     scalar cp = 0.0;
00298 
00299     forAll(properties_, i)
00300     {
00301         if (x[i] > SMALL)
00302         {
00303             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00304             cp += x[i]*properties_[i].cp(p, Ti)*properties_[i].W();
00305         }
00306     }
00307 
00308     return cp/W(x);
00309 }
00310 
00311 
00312 Foam::scalar Foam::liquidMixture::sigma
00313 (
00314     const scalar p,
00315     const scalar T,
00316     const scalarField& x
00317 ) const
00318 {
00319     // sigma is based on surface mole fractions
00320     // which is estimated from Raoult's Law
00321     scalar sigma = 0.0;
00322     scalarField Xs(x.size(), 0.0);
00323     scalar XsSum = 0.0;
00324     forAll(properties_, i)
00325     {
00326         scalar Ti = min(TrMax*properties_[i].Tc(), T);
00327         scalar Pvs = properties_[i].pv(p, Ti);
00328         scalar xs = x[i]*Pvs/p;
00329         XsSum += xs;
00330         Xs[i] = xs;
00331     }
00332 
00333     forAll(properties_, i)
00334     {
00335         if (Xs[i] > SMALL)
00336         {
00337             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00338             sigma += (Xs[i]/XsSum)*properties_[i].sigma(p, Ti);
00339         }
00340     }
00341 
00342     return sigma;
00343 }
00344 
00345 
00346 Foam::scalar Foam::liquidMixture::mu
00347 (
00348     const scalar p,
00349     const scalar T,
00350     const scalarField& x
00351 ) const
00352 {
00353     scalar mu = 0.0;
00354 
00355     forAll(properties_, i)
00356     {
00357         if (x[i] > SMALL)
00358         {
00359             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00360             mu += x[i]*log(properties_[i].mu(p, Ti));
00361         }
00362     }
00363 
00364     return exp(mu);
00365 }
00366 
00367 
00368 Foam::scalar Foam::liquidMixture::K
00369 (
00370     const scalar p,
00371     const scalar T,
00372     const scalarField& x
00373 ) const
00374 {
00375     // calculate superficial volume fractions phii
00376     scalarField phii(x.size(), 0.0);
00377     scalar pSum = 0.0;
00378 
00379     forAll(properties_, i)
00380     {
00381         scalar Ti = min(TrMax*properties_[i].Tc(), T);
00382 
00383         scalar Vi = properties_[i].W()/properties_[i].rho(p, Ti);
00384         phii[i] = x[i]*Vi;
00385         pSum += phii[i];
00386     }
00387 
00388     forAll(phii, i)
00389     {
00390         phii[i] /= pSum;
00391     }
00392 
00393     scalar K = 0.0;
00394 
00395     forAll(properties_, i)
00396     {
00397         scalar Ti = min(TrMax*properties_[i].Tc(), T);
00398 
00399         forAll(properties_, j)
00400         {
00401             scalar Tj = min(TrMax*properties_[j].Tc(), T);
00402 
00403             scalar Kij =
00404                 2.0
00405                /(
00406                     1.0/properties_[i].K(p, Ti)
00407                   + 1.0/properties_[j].K(p, Tj)
00408                 );
00409             K += phii[i]*phii[j]*Kij;
00410         }
00411     }
00412 
00413     return K;
00414 }
00415 
00416 
00417 Foam::scalar Foam::liquidMixture::D
00418 (
00419     const scalar p,
00420     const scalar T,
00421     const scalarField& x
00422 ) const
00423 {
00424     // Blanc's law
00425     scalar Dinv = 0.0;
00426 
00427     forAll(properties_, i)
00428     {
00429         if (x[i] > SMALL)
00430         {
00431             scalar Ti = min(TrMax*properties_[i].Tc(), T);
00432             Dinv += x[i]/properties_[i].D(p, Ti);
00433         }
00434     }
00435 
00436     return 1.0/Dinv;
00437 }
00438 
00439 
00440 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines