diff --git a/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.C b/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.C
index aa0f3829a..4608a67fb 100644
--- a/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.C
+++ b/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.C
@@ -25,6 +25,8 @@ License
 
 #include "wideBandAbsorptionEmission.H"
 #include "addToRunTimeSelectionTable.H"
+#include "basicSpecieMixture.H"
+#include "unitConversion.H"
 
 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
 
@@ -56,12 +58,7 @@ Foam::radiation::wideBandAbsorptionEmission::wideBandAbsorptionEmission
     coeffsDict_((dict.optionalSubDict(typeName + "Coeffs"))),
     speciesNames_(0),
     specieIndex_(label(0)),
-    lookUpTable_
-    (
-        fileName(coeffsDict_.lookup("lookUpTableFileName")),
-        mesh.time().constant(),
-        mesh
-    ),
+    lookUpTablePtr_(),
     thermo_(mesh.lookupObject<fluidThermo>(basicThermo::dictName)),
     Yj_(nSpecies_),
     totalWaveLength_(0)
@@ -95,7 +92,7 @@ Foam::radiation::wideBandAbsorptionEmission::wideBandAbsorptionEmission
                 if (!speciesNames_.found(key))
                 {
                     FatalErrorInFunction
-                        << "specie: " << key << "is not in all the bands"
+                        << "specie: " << key << " is not in all the bands"
                         << nl << exit(FatalError);
                 }
             }
@@ -106,36 +103,82 @@ Foam::radiation::wideBandAbsorptionEmission::wideBandAbsorptionEmission
     }
     nBands_ = nBand;
 
+    if (coeffsDict_.found("lookUpTableFileName"))
+    {
+        const word name = coeffsDict_.lookup("lookUpTableFileName");
+        if (name != "none")
+        {
+            lookUpTablePtr_.set
+            (
+                new interpolationLookUpTable<scalar>
+                (
+                    fileName(coeffsDict_.lookup("lookUpTableFileName")),
+                    mesh.time().constant(),
+                    mesh
+                )
+            );
+
+            if (!mesh.foundObject<volScalarField>("ft"))
+            {
+                FatalErrorInFunction
+                    << "specie ft is not present to use with "
+                    << "lookUpTableFileName " << nl
+                    << exit(FatalError);
+            }
+        }
+    }
+
     // Check that all the species on the dictionary are present in the
-    // look-up table  and save the corresponding indices of the look-up table
+    // look-up table and save the corresponding indices of the look-up table
 
     label j = 0;
     forAllConstIter(HashTable<label>, speciesNames_, iter)
     {
-        if (lookUpTable_.found(iter.key()))
+        if (!lookUpTablePtr_.empty())
         {
-            label index = lookUpTable_.findFieldIndex(iter.key());
-            Info<< "specie: " << iter.key() << " found in look-up table "
-                << " with index: " << index << endl;
-            specieIndex_[iter()] = index;
+            if (lookUpTablePtr_().found(iter.key()))
+            {
+                label index = lookUpTablePtr_().findFieldIndex(iter.key());
+
+                Info<< "specie: " << iter.key() << " found on look-up table "
+                    << " with index: " << index << endl;
+
+                specieIndex_[iter()] = index;
+            }
+            else if (mesh.foundObject<volScalarField>(iter.key()))
+            {
+                Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
+                specieIndex_[iter()] = 0;
+                j++;
+                Info<< "specie: " << iter.key() << " is being solved" << endl;
+            }
+            else
+            {
+                FatalErrorInFunction
+                    << "specie: " << iter.key()
+                    << " is neither in look-up table: "
+                    << lookUpTablePtr_().tableName()
+                    << " nor is being solved" << nl
+                    << exit(FatalError);
+            }
         }
         else if (mesh.foundObject<volScalarField>(iter.key()))
         {
             Yj_.set(j, &mesh.lookupObjectRef<volScalarField>(iter.key()));
-
-            specieIndex_[iter()] = 0.0;
+            specieIndex_[iter()] = 0;
             j++;
-            Info<< "species: " << iter.key() << " is being solved" << endl;
         }
         else
         {
             FatalErrorInFunction
-                << "specie: " << iter.key()
-                << " is neither in look-up table : "
-                << lookUpTable_.tableName() << " nor is being solved"
+                << " there is no lookup table and the specie" << nl
+                << iter.key() << nl
+                << " is not found " << nl
                 << exit(FatalError);
+
         }
     }
+    
 }
 
 
@@ -151,11 +194,11 @@ Foam::radiation::wideBandAbsorptionEmission::~wideBandAbsorptionEmission()
 Foam::tmp<Foam::volScalarField>
 Foam::radiation::wideBandAbsorptionEmission::aCont(const label bandI) const
 {
+    const basicSpecieMixture& mixture =
+        dynamic_cast<const basicSpecieMixture&>(thermo_);
+
     const volScalarField& T = thermo_.T();
     const volScalarField& p = thermo_.p();
-    const volScalarField& ft = mesh_.lookupObject<volScalarField>("ft");
-
-    label nSpecies = speciesNames_.size();
 
     tmp<volScalarField> ta
     (
@@ -176,38 +219,47 @@ Foam::radiation::wideBandAbsorptionEmission::aCont(const label bandI) const
 
     scalarField& a = ta.ref().primitiveFieldRef();
 
-    forAll(a, i)
+    forAll(a, celli)
     {
-        const List<scalar>& species = lookUpTable_.lookUp(ft[i]);
-
-        for (label n=0; n<nSpecies; n++)
+        forAllConstIter(HashTable<label>, speciesNames_, iter)
         {
-            label l = 0;
-            scalar Yipi = 0.0;
+            label n = iter();
+            scalar Xipi = 0.0;
             if (specieIndex_[n] != 0)
             {
-                // moles x pressure [atm]
-                Yipi = species[specieIndex_[n]]*p[i]*9.869231e-6;
+	        const volScalarField& ft =
+		    mesh_.lookupObject<volScalarField>("ft");
+
+                const List<scalar>& Ynft = lookUpTablePtr_().lookUp(ft[celli]);
+                //moles x pressure [atm]
+                Xipi = Ynft[specieIndex_[n]]*paToAtm(p[celli]);
             }
             else
             {
-                // mass fraction from species being solved
-                Yipi = Yj_[l][i];
-                l++;
+                scalar invWt = 0.0;
+                forAll(mixture.Y(), s)
+                {
+                    invWt += mixture.Y(s)[celli]/mixture.W(s);
+                }
+
+                label index = mixture.species()[iter.key()];
+                scalar Xk = mixture.Y(index)[celli]/(mixture.W(index)*invWt);
+
+                Xipi = Xk*paToAtm(p[celli]);
             }
 
-            scalar Ti = T[i];
+            scalar Ti = T[celli];
 
             const absorptionCoeffs::coeffArray& b =
-                coeffs_[n][bandI].coeffs(T[i]);
+	        coeffs_[bandI][n].coeffs(T[celli]);
 
-            if (coeffs_[n][bandI].invTemp())
+            if (coeffs_[bandI][n].invTemp())
             {
-                Ti = 1.0/T[i];
+                Ti = 1.0/T[celli];
             }
 
-            a[i]+=
-                Yipi
+            a[celli]+=
+                Xipi
                *(
                     ((((b[5]*Ti + b[4])*Ti + b[3])*Ti + b[2])*Ti + b[1])*Ti
                   + b[0]
diff --git a/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.H b/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.H
index e1fb10cf0..880537b68 100644
--- a/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.H
+++ b/src/thermophysicalModels/radiation/submodels/absorptionEmissionModel/wideBandAbsorptionEmission/wideBandAbsorptionEmission.H
@@ -37,11 +37,8 @@ Description
 
     The emission constant proportionality is specified per band (EhrrCoeff).
 
-    The coefficients for the species in the lookup table have to be specified
-    for use in moles x P [atm].i.e. (k[i] = species[i]*p*9.869231e-6).
-
-    The coefficients for CO and soot or any other added are multiplied by the
-    respective mass fraction being solved.
+    The coefficients for the species have to be specified for use in 
+    moles x P [atm], i.e. (k[i] = species[i]*p*9.869231e-6).
 
     The look Up table file should be in the constant directory.
 
@@ -156,8 +153,8 @@ private:
         //- Proportion of the heat released rate emitted
         FixedList<scalar, maxBands_> iEhrrCoeffs_;
 
-        //- Lookup table of species related to ft
-        mutable interpolationLookUpTable<scalar> lookUpTable_;
+        //- Look-up table of species related to ft
+        mutable autoPtr<interpolationLookUpTable<scalar>> lookUpTablePtr_;
 
         //- Thermo package
         const fluidThermo& thermo_;
