diff --git a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H
index 6a4a61b6479b324804e7336e8eeaf32ce34fbbce..3ece75e1780725d6b08c69841097291f61fd3017 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcel.H
@@ -6,7 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
-    Copyright (C) 2016-2019 OpenCFD Ltd.
+    Copyright (C) 2016-2020 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -339,6 +339,18 @@ public:
             //- Return const access to mass fractions of mixture []
             inline const scalarField& Y() const;
 
+            //- Return const access to mass fractions of gases
+            //  Note: for compatibilty only - returns Y()
+            inline const scalarField& YGas() const;
+
+            //- Return const access to mass fractions of liquids
+            //  Note: for compatibilty only - returns Y()
+            inline const scalarField& YLiquid() const;
+
+            //- Return const access to mass fractions of solids
+            //  Note: for compatibilty only - returns Y()
+            inline const scalarField& YSolid() const;
+
 
         // Edit
 
diff --git a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H
index 45a0ef0900bdd4befbb4ff500358ecd57e6d5cb6..b8beb3ed16c5225015f5448c6dac5e62f11f9e98 100644
--- a/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H
+++ b/src/lagrangian/intermediate/parcels/Templates/ReactingParcel/ReactingParcelI.H
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2017 OpenFOAM Foundation
+    Copyright (C) 2020 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -170,6 +171,29 @@ inline const Foam::scalarField& Foam::ReactingParcel<ParcelType>::Y() const
 }
 
 
+template<class ParcelType>
+inline const Foam::scalarField& Foam::ReactingParcel<ParcelType>::YGas() const
+{
+    return Y_;
+}
+
+
+template<class ParcelType>
+inline const Foam::scalarField&
+Foam::ReactingParcel<ParcelType>::YLiquid() const
+{
+    return Y_;
+}
+
+
+template<class ParcelType>
+inline const Foam::scalarField&
+Foam::ReactingParcel<ParcelType>::YSolid() const
+{
+    return Y_;
+}
+ 
+
 template<class ParcelType>
 inline Foam::scalar& Foam::ReactingParcel<ParcelType>::mass0()
 {
diff --git a/src/lagrangian/intermediate/parcels/derived/basicReactingMultiphaseParcel/makeBasicReactingMultiphaseParcelSubmodels.C b/src/lagrangian/intermediate/parcels/derived/basicReactingMultiphaseParcel/makeBasicReactingMultiphaseParcelSubmodels.C
index 9b5db0fc485bfdb5123317703ac21e6890d1f6e1..dc329d39278f16a6b1ed384a0ddf4594c8f37c23 100644
--- a/src/lagrangian/intermediate/parcels/derived/basicReactingMultiphaseParcel/makeBasicReactingMultiphaseParcelSubmodels.C
+++ b/src/lagrangian/intermediate/parcels/derived/basicReactingMultiphaseParcel/makeBasicReactingMultiphaseParcelSubmodels.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2015 OpenFOAM Foundation
+    Copyright (C) 2020 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -27,7 +28,7 @@ License
 
 #include "basicReactingMultiphaseCloud.H"
 
-#include "makeParcelCloudFunctionObjects.H"
+#include "makeReactingParcelCloudFunctionObjects.H" // Reacting variant
 
 // Kinematic
 #include "makeThermoParcelForces.H" // thermo variant
diff --git a/src/lagrangian/intermediate/parcels/derived/basicReactingParcel/makeBasicReactingParcelSubmodels.C b/src/lagrangian/intermediate/parcels/derived/basicReactingParcel/makeBasicReactingParcelSubmodels.C
index b0d8029d75c4a2b6b441f88d8e6fd2ad2251c02a..09f7f95d28849c60125016c1c82fe518a07c5072 100644
--- a/src/lagrangian/intermediate/parcels/derived/basicReactingParcel/makeBasicReactingParcelSubmodels.C
+++ b/src/lagrangian/intermediate/parcels/derived/basicReactingParcel/makeBasicReactingParcelSubmodels.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2015 OpenFOAM Foundation
+    Copyright (C) 2020 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -27,7 +28,7 @@ License
 
 #include "basicReactingCloud.H"
 
-#include "makeParcelCloudFunctionObjects.H"
+#include "makeReactingParcelCloudFunctionObjects.H" // Reacting variant
 
 // Kinematic
 #include "makeThermoParcelForces.H" // thermo variant
diff --git a/src/lagrangian/intermediate/parcels/include/makeReactingParcelCloudFunctionObjects.H b/src/lagrangian/intermediate/parcels/include/makeReactingParcelCloudFunctionObjects.H
new file mode 100644
index 0000000000000000000000000000000000000000..ca86c65d63191fb64f0a53e6ae571716f3a979c6
--- /dev/null
+++ b/src/lagrangian/intermediate/parcels/include/makeReactingParcelCloudFunctionObjects.H
@@ -0,0 +1,65 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2011-2018 OpenFOAM Foundation
+    Copyright (C) 2020 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef makeReactingParcelCloudFunctionObjects_H
+#define makeReactingParcelCloudFunctionObjects_H
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#include "FacePostProcessing.H"
+#include "ParticleCollector.H"
+#include "ParticleErosion.H"
+#include "ParticleTracks.H"
+#include "ParticleTrap.H"
+#include "PatchCollisionDensity.H"
+#include "PatchPostProcessing.H"
+#include "VoidFraction.H"
+#include "WeberNumberReacting.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#define makeParcelCloudFunctionObjects(CloudType)                              \
+                                                                               \
+    makeCloudFunctionObject(CloudType);                                        \
+                                                                               \
+    makeCloudFunctionObjectType(FacePostProcessing, CloudType);                \
+    makeCloudFunctionObjectType(ParticleCollector, CloudType);                 \
+    makeCloudFunctionObjectType(ParticleErosion, CloudType);                   \
+    makeCloudFunctionObjectType(ParticleTracks, CloudType);                    \
+    makeCloudFunctionObjectType(ParticleTrap, CloudType);                      \
+    makeCloudFunctionObjectType(PatchCollisionDensity, CloudType);             \
+    makeCloudFunctionObjectType(PatchPostProcessing, CloudType);               \
+    makeCloudFunctionObjectType(VoidFraction, CloudType);                      \
+    makeCloudFunctionObjectType(WeberNumberReacting, CloudType);
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/lagrangian/intermediate/submodels/CloudFunctionObjects/WeberNumber/WeberNumberReacting.C b/src/lagrangian/intermediate/submodels/CloudFunctionObjects/WeberNumber/WeberNumberReacting.C
new file mode 100644
index 0000000000000000000000000000000000000000..f5a72900e8e42bfe3a3c365b5f69b842fd7ca60a
--- /dev/null
+++ b/src/lagrangian/intermediate/submodels/CloudFunctionObjects/WeberNumber/WeberNumberReacting.C
@@ -0,0 +1,124 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2020 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+\*---------------------------------------------------------------------------*/
+
+#include "WeberNumberReacting.H"
+#include "SLGThermo.H"
+
+// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
+
+template<class CloudType>
+Foam::WeberNumberReacting<CloudType>::WeberNumberReacting
+(
+    const dictionary& dict,
+    CloudType& owner,
+    const word& modelName
+)
+:
+    CloudFunctionObject<CloudType>(dict, owner, modelName, typeName)
+{}
+
+
+template<class CloudType>
+Foam::WeberNumberReacting<CloudType>::WeberNumberReacting
+(
+    const WeberNumberReacting<CloudType>& we
+)
+:
+    CloudFunctionObject<CloudType>(we)
+{}
+
+
+// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
+
+template<class CloudType>
+void Foam::WeberNumberReacting<CloudType>::postEvolve
+(
+    const typename parcelType::trackingData& td
+)
+{
+    const auto& c = this->owner();
+
+    if (!c.template foundObject<IOField<scalar>>("We"))
+    {
+        IOField<scalar>* WePtr =
+            new IOField<scalar>
+            (
+                IOobject
+                (
+                    "We",
+                    c.time().timeName(),
+                    c,
+                    IOobject::NO_READ
+                )
+            );
+
+        WePtr->store();
+    }
+
+    auto& We = c.template lookupObjectRef<IOField<scalar>>("We");
+    We.setSize(c.size());
+
+    const auto& thermo = c.db().template lookupObject<SLGThermo>("SLGThermo");
+    const auto& liquids = thermo.liquids();
+
+    const auto& UInterp = td.UInterp();
+    const auto& pInterp = td.pInterp();
+    const auto& rhoInterp = td.rhoInterp();
+
+    label parceli = 0;
+    forAllConstIters(c, parcelIter)
+    {
+        const parcelType& p = parcelIter();
+
+        const auto& coords = p.coordinates();
+        const auto& tetIs = p.currentTetIndices();
+
+        const vector Uc(UInterp.interpolate(coords, tetIs));
+
+        const scalar pc =
+            max
+            (
+                pInterp.interpolate(coords, tetIs),
+                c.constProps().pMin()
+            );
+
+        const scalar rhoc(rhoInterp.interpolate(coords, tetIs));
+        const scalarField X(liquids.X(p.YLiquid()));
+        const scalar sigma = liquids.sigma(pc, p.T(), X);
+
+        We[parceli++] = rhoc*magSqr(p.U() - Uc)*p.d()/sigma;
+    }
+
+
+    if (c.size() && c.time().writeTime())
+    {
+        We.write();
+    }
+}
+
+
+// ************************************************************************* //
diff --git a/src/lagrangian/intermediate/submodels/CloudFunctionObjects/WeberNumber/WeberNumberReacting.H b/src/lagrangian/intermediate/submodels/CloudFunctionObjects/WeberNumber/WeberNumberReacting.H
new file mode 100644
index 0000000000000000000000000000000000000000..b686fbe04176f6a82497201699ca8926fe2d3a50
--- /dev/null
+++ b/src/lagrangian/intermediate/submodels/CloudFunctionObjects/WeberNumber/WeberNumberReacting.H
@@ -0,0 +1,121 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | www.openfoam.com
+     \\/     M anipulation  |
+-------------------------------------------------------------------------------
+    Copyright (C) 2020 OpenCFD Ltd.
+-------------------------------------------------------------------------------
+License
+    This file is part of OpenFOAM.
+
+    OpenFOAM is free software: you can redistribute it and/or modify it
+    under the terms of the GNU General Public License as published by
+    the Free Software Foundation, either version 3 of the License, or
+    (at your option) any later version.
+
+    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
+    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
+    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
+    for more details.
+
+    You should have received a copy of the GNU General Public License
+    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
+
+Class
+    Foam::WeberNumberReacting
+
+Group
+    grpLagrangianIntermediateFunctionObjects
+
+Description
+    Creates particle Weber number field on the cloud
+
+SourceFiles
+    WeberNumberReacting.C
+
+\*---------------------------------------------------------------------------*/
+
+#ifndef WeberNumberReacting_H
+#define WeberNumberReacting_H
+
+#include "CloudFunctionObject.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+namespace Foam
+{
+
+/*---------------------------------------------------------------------------*\
+                       Class WeberNumberReacting Declaration
+\*---------------------------------------------------------------------------*/
+
+template<class CloudType>
+class WeberNumberReacting
+:
+    public CloudFunctionObject<CloudType>
+{
+    // Private Data
+
+        // Typedefs
+
+            //- Convenience typedef for parcel type
+            typedef typename CloudType::parcelType parcelType;
+
+
+public:
+
+    //- Runtime type information
+    TypeName("WeberNumber");
+
+
+    // Constructors
+
+        //- Construct from dictionary
+        WeberNumberReacting
+        (
+            const dictionary& dict,
+            CloudType& owner,
+            const word& modelName
+        );
+
+        //- Construct copy
+        WeberNumberReacting(const WeberNumberReacting<CloudType>& vf);
+
+        //- Construct and return a clone
+        virtual autoPtr<CloudFunctionObject<CloudType>> clone() const
+        {
+            return autoPtr<CloudFunctionObject<CloudType>>
+            (
+                new WeberNumberReacting<CloudType>(*this)
+            );
+        }
+
+
+    //- Destructor
+    virtual ~WeberNumberReacting() = default;
+
+
+    // Member Functions
+
+        //- Post-evolve hook
+        virtual void postEvolve(const typename parcelType::trackingData& td);
+};
+
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+} // End namespace Foam
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#ifdef NoRepository
+    #include "WeberNumberReacting.C"
+#endif
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+#endif
+
+// ************************************************************************* //
diff --git a/src/lagrangian/spray/parcels/derived/basicSprayParcel/makeBasicSprayParcelSubmodels.C b/src/lagrangian/spray/parcels/derived/basicSprayParcel/makeBasicSprayParcelSubmodels.C
index 9ad6f70cff7829e3f4a888ec8408cf0a2fb1cd86..0cb36971ca1e24d29bdd16dc81d60d197e6d988b 100644
--- a/src/lagrangian/spray/parcels/derived/basicSprayParcel/makeBasicSprayParcelSubmodels.C
+++ b/src/lagrangian/spray/parcels/derived/basicSprayParcel/makeBasicSprayParcelSubmodels.C
@@ -6,6 +6,7 @@
      \\/     M anipulation  |
 -------------------------------------------------------------------------------
     Copyright (C) 2011-2015 OpenFOAM Foundation
+    Copyright (C) 2020 OpenCFD Ltd.
 -------------------------------------------------------------------------------
 License
     This file is part of OpenFOAM.
@@ -27,7 +28,7 @@ License
 
 #include "basicSprayCloud.H"
 
-#include "makeParcelCloudFunctionObjects.H"
+#include "makeReactingParcelCloudFunctionObjects.H" // Reacting variant
 
 // Kinematic
 #include "makeThermoParcelForces.H" // thermo variant
diff --git a/tutorials/lagrangian/sprayFoam/aachenBomb/constant/sprayCloudProperties b/tutorials/lagrangian/sprayFoam/aachenBomb/constant/sprayCloudProperties
index 412e822d05281bfc1ac7ad1a8a04ec938a322dab..a613d069b1de1a9f43620e2d73b443db36827f7e 100644
--- a/tutorials/lagrangian/sprayFoam/aachenBomb/constant/sprayCloudProperties
+++ b/tutorials/lagrangian/sprayFoam/aachenBomb/constant/sprayCloudProperties
@@ -228,7 +228,12 @@ subModels
 
 
 cloudFunctions
-{}
+{
+    WeberNumber1
+    {
+        type    WeberNumber;
+    }
+}
 
 
 // ************************************************************************* //