From 2f8e51ae7c9415450f38e891fba508e1aa87799a Mon Sep 17 00:00:00 2001 From: Henry <Henry> Date: Fri, 8 May 2015 16:40:03 +0100 Subject: [PATCH] surfaceFilmModels: Replaced removeInjection with patchInjection patchInjection accumulates mass exchanges on a per-patch basis which are reported by injectionModelList --- .../primitives/subModelBase/subModelBase.C | 8 +- .../primitives/subModelBase/subModelBase.H | 12 +- .../subModelBase/subModelBaseTemplates.C | 34 +--- src/regionModels/surfaceFilmModels/Make/files | 2 +- .../injectionModel/injectionModel.H | 32 +-- .../injectionModelList/injectionModelList.C | 37 ++-- .../patchInjection/patchInjection.C | 191 ++++++++++++++++++ .../patchInjection.H} | 59 +++--- .../removeInjection/removeInjection.C | 121 ----------- 9 files changed, 284 insertions(+), 212 deletions(-) create mode 100644 src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.C rename src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/{removeInjection/removeInjection.H => patchInjection/patchInjection.H} (66%) delete mode 100644 src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.C diff --git a/src/OpenFOAM/primitives/subModelBase/subModelBase.C b/src/OpenFOAM/primitives/subModelBase/subModelBase.C index 41b8af6a29e..4a7c7017eb2 100644 --- a/src/OpenFOAM/primitives/subModelBase/subModelBase.C +++ b/src/OpenFOAM/primitives/subModelBase/subModelBase.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -158,9 +158,7 @@ bool Foam::subModelBase::active() const void Foam::subModelBase::cacheFields(const bool) -{ - // do nothing -} +{} bool Foam::subModelBase::outputTime() const @@ -171,8 +169,6 @@ bool Foam::subModelBase::outputTime() const void Foam::subModelBase::write(Ostream& os) const { - // not writing complete cloud dictionary, only coeffs -// os << dict_; os << coeffDict_; } diff --git a/src/OpenFOAM/primitives/subModelBase/subModelBase.H b/src/OpenFOAM/primitives/subModelBase/subModelBase.H index a6fbc65615a..92627782e9d 100644 --- a/src/OpenFOAM/primitives/subModelBase/subModelBase.H +++ b/src/OpenFOAM/primitives/subModelBase/subModelBase.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -51,8 +51,6 @@ namespace Foam class subModelBase { -private: - // Private Member Functions //- Disallow default bitwise assignment @@ -181,6 +179,10 @@ public: // Model properties + //- Retrieve generic property from the sub-model + template<class Type> + void getModelProperty(const word& entryName, Type& value) const; + //- Retrieve generic property from the sub-model template<class Type> Type getModelProperty @@ -189,10 +191,6 @@ public: const Type& defaultValue = pTraits<Type>::zero ) const; - //- Retrieve generic property from the sub-model - template<class Type> - void getModelProperty(const word& entryName, Type& value) const; - //- Add generic property to the sub-model template<class Type> void setModelProperty(const word& entryName, const Type& value); diff --git a/src/OpenFOAM/primitives/subModelBase/subModelBaseTemplates.C b/src/OpenFOAM/primitives/subModelBase/subModelBaseTemplates.C index df881986e8e..4264d4fb136 100644 --- a/src/OpenFOAM/primitives/subModelBase/subModelBaseTemplates.C +++ b/src/OpenFOAM/primitives/subModelBase/subModelBaseTemplates.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2013 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2013-2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -80,52 +80,38 @@ void Foam::subModelBase::setBaseProperty template<class Type> -Type Foam::subModelBase::getModelProperty +void Foam::subModelBase::getModelProperty ( const word& entryName, - const Type& defaultValue + Type& value ) const { - Type result = defaultValue; - if (properties_.found(baseName_)) { const dictionary& baseDict = properties_.subDict(baseName_); if (inLine() && baseDict.found(modelName_)) { - baseDict.subDict(modelName_).readIfPresent(entryName, result); + baseDict.subDict(modelName_).readIfPresent(entryName, value); } else if (baseDict.found(modelType_)) { - baseDict.subDict(modelType_).readIfPresent(entryName, result); + baseDict.subDict(modelType_).readIfPresent(entryName, value); } } - - return result; } template<class Type> -void Foam::subModelBase::getModelProperty +Type Foam::subModelBase::getModelProperty ( const word& entryName, - Type& value + const Type& defaultValue ) const { - if (properties_.found(baseName_)) - { - const dictionary& baseDict = properties_.subDict(baseName_); - - if (inLine() && baseDict.found(modelName_)) - { - baseDict.subDict(modelName_).readIfPresent(entryName, value); - } - else if (baseDict.found(modelType_)) - { - baseDict.subDict(modelType_).readIfPresent(entryName, value); - } - } + Type result = defaultValue; + getModelProperty(entryName, result); + return result; } diff --git a/src/regionModels/surfaceFilmModels/Make/files b/src/regionModels/surfaceFilmModels/Make/files index 5c535dbf6d2..15accdfa79c 100644 --- a/src/regionModels/surfaceFilmModels/Make/files +++ b/src/regionModels/surfaceFilmModels/Make/files @@ -20,7 +20,7 @@ $(KINEMATICMODELS)/injectionModel/injectionModel/injectionModel.C $(KINEMATICMODELS)/injectionModel/injectionModel/injectionModelNew.C $(KINEMATICMODELS)/injectionModel/injectionModelList/injectionModelList.C $(KINEMATICMODELS)/injectionModel/drippingInjection/drippingInjection.C -$(KINEMATICMODELS)/injectionModel/removeInjection/removeInjection.C +$(KINEMATICMODELS)/injectionModel/patchInjection/patchInjection.C $(KINEMATICMODELS)/injectionModel/curvatureSeparation/curvatureSeparation.C $(KINEMATICMODELS)/filmThermoModel/filmThermoModel/filmThermoModel.C diff --git a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModel/injectionModel.H b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModel/injectionModel.H index 70ec873ce9c..b365ff737c3 100644 --- a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModel/injectionModel.H +++ b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModel/injectionModel.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -58,8 +58,6 @@ class injectionModel : public filmSubModelBase { -private: - // Private data //- Injected mass @@ -106,6 +104,7 @@ public: (owner, dict) ); + // Constructors //- Construct null @@ -137,18 +136,21 @@ public: // Member Functions - // Evolution - - //- Correct - virtual void correct - ( - scalarField& availableMass, - scalarField& massToInject, - scalarField& diameterToInject - ) = 0; - - //- Return the total mass injected - scalar injectedMassTotal() const; + //- Correct + virtual void correct + ( + scalarField& availableMass, + scalarField& massToInject, + scalarField& diameterToInject + ) = 0; + + //- Return the total mass injected + virtual scalar injectedMassTotal() const; + + //- Accumulate the total mass injected for the patches into the + // scalarField provided + virtual void patchInjectedMassTotals(scalarField& patchMasses) const + {} }; diff --git a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModelList/injectionModelList.C b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModelList/injectionModelList.C index 5cdebac7c0c..132106a3a54 100644 --- a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModelList/injectionModelList.C +++ b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/injectionModelList/injectionModelList.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -118,44 +118,57 @@ void injectionModelList::correct forAll(patchIDs, i) { - label patchI = patchIDs[i]; + label patchi = patchIDs[i]; massInjected_[i] = - massInjected_[i] + sum(massToInject.boundaryField()[patchI]); + massInjected_[i] + sum(massToInject.boundaryField()[patchi]); } } void injectionModelList::info(Ostream& os) { - scalar injectedMass = 0.0; + const polyBoundaryMesh& pbm = owner().regionMesh().boundaryMesh(); + + scalar injectedMass = 0; + scalarField patchInjectedMasses(pbm.size(), 0); + forAll(*this, i) { const injectionModel& im = operator[](i); injectedMass += im.injectedMassTotal(); + im.patchInjectedMassTotals(patchInjectedMasses); } os << indent << "injected mass = " << injectedMass << nl; - scalarList mass0(massInjected_.size(), 0.0); - this->getModelProperty("massInjected", mass0); + forAll(pbm, patchi) + { + if (mag(patchInjectedMasses[patchi]) > VSMALL) + { + os << indent << indent << "from patch " << pbm[patchi].name() + << " = " << patchInjectedMasses[patchi] << nl; + } + } - scalarList mass(massInjected_); + scalarField mass0(massInjected_.size(), 0.0); + this->getBaseProperty("massInjected", mass0); + + scalarField mass(massInjected_); Pstream::listCombineGather(mass, plusEqOp<scalar>()); - mass = mass + mass0; + mass += mass0; - const polyBoundaryMesh& pbm = owner().regionMesh().boundaryMesh(); const labelList& patchIDs = owner().intCoupledPatchIDs(); forAll(patchIDs, i) { - label patchI = patchIDs[i]; - Info<< indent << " - patch: " << pbm[patchI].name() << ": " + label patchi = patchIDs[i]; + Info<< indent << " - patch: " << pbm[patchi].name() << ": " << mass[i] << endl; } if (owner().time().outputTime()) { - setModelProperty("massInjected", mass); + setBaseProperty("massInjected", mass); massInjected_ = 0.0; } } diff --git a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.C b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.C new file mode 100644 index 00000000000..218c727a86e --- /dev/null +++ b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.C @@ -0,0 +1,191 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation + \\/ M anipulation | +------------------------------------------------------------------------------- +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 "patchInjection.H" +#include "addToRunTimeSelectionTable.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +namespace Foam +{ +namespace regionModels +{ +namespace surfaceFilmModels +{ + +// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // + +defineTypeNameAndDebug(patchInjection, 0); +addToRunTimeSelectionTable(injectionModel, patchInjection, dictionary); + +// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // + +patchInjection::patchInjection +( + surfaceFilmModel& owner, + const dictionary& dict +) +: + injectionModel(type(), owner, dict), + deltaStable_(coeffDict_.lookupOrDefault<scalar>("deltaStable", 0.0)) +{ + const polyBoundaryMesh& pbm = owner.regionMesh().boundaryMesh(); + patchIDs_.setSize(pbm.size()); + + if (coeffDict_.found("patches")) + { + wordReList patchNames(coeffDict_.lookup("patches")); + const labelHashSet patchSet = pbm.patchSet(patchNames); + + Info<< " applying to patches:" << nl; + + label pidi = 0; + forAllConstIter(labelHashSet, patchSet, iter) + { + label patchi = iter.key(); + patchIDs_[pidi++] = patchi; + Info<< " " << pbm[patchi].name() << endl; + } + patchIDs_.setSize(pidi); + patchInjectedMasses_.setSize(pidi, 0); + } + else + { + Info<< " applying to all patches" << endl; + + forAll(patchIDs_, patchi) + { + patchIDs_[patchi] = patchi; + } + + patchInjectedMasses_.setSize(pbm.size(), 0); + } +} + + +// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // + +patchInjection::~patchInjection() +{} + + +// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // + +void patchInjection::correct +( + scalarField& availableMass, + scalarField& massToInject, + scalarField& diameterToInject +) +{ + const scalarField& delta = owner().delta(); + const scalarField& rho = owner().rho(); + const scalarField& magSf = owner().magSf(); + + const polyBoundaryMesh& pbm = owner().regionMesh().boundaryMesh(); + + forAll(patchIDs_, pidi) + { + label patchi = patchIDs_[pidi]; + const polyPatch& pp = pbm[patchi]; + const labelList& faceCells = pp.faceCells(); + + // Accumulate the total mass removed from patch + scalar dMassPatch = 0; + + forAll(faceCells, fci) + { + label celli = faceCells[fci]; + + scalar ddelta = max(0.0, delta[celli] - deltaStable_); + scalar dMass = ddelta*rho[celli]*magSf[celli]; + massToInject[celli] += dMass; + availableMass[celli] -= dMass; + dMassPatch += dMass; + } + + patchInjectedMasses_[pidi] += dMassPatch; + addToInjectedMass(dMassPatch); + } + + injectionModel::correct(); + + if (outputTime()) + { + scalarField patchInjectedMasses0 + ( + getModelProperty<scalarField> + ( + "patchInjectedMasses", + scalarField(patchInjectedMasses_.size(), 0) + ) + ); + + scalarField patchInjectedMassTotals(patchInjectedMasses_); + Pstream::listCombineGather(patchInjectedMassTotals, plusEqOp<scalar>()); + patchInjectedMasses0 += patchInjectedMassTotals; + + setModelProperty<scalarField> + ( + "patchInjectedMasses", + patchInjectedMasses0 + ); + + patchInjectedMasses_ = 0; + } +} + + +void patchInjection::patchInjectedMassTotals(scalarField& patchMasses) const +{ + scalarField patchInjectedMasses + ( + getModelProperty<scalarField> + ( + "patchInjectedMasses", + scalarField(patchInjectedMasses_.size(), 0) + ) + ); + + scalarField patchInjectedMassTotals(patchInjectedMasses_); + //combineReduce(patchInjectedMassTotals, plusEqOp<scalarField>()); + Pstream::listCombineGather(patchInjectedMassTotals, plusEqOp<scalar>()); + + forAll(patchIDs_, pidi) + { + label patchi = patchIDs_[pidi]; + patchMasses[patchi] += + patchInjectedMasses[pidi] + patchInjectedMassTotals[pidi]; + } +} + + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +} // End namespace surfaceFilmModels +} // End namespace regionModels +} // End namespace Foam + +// ************************************************************************* // diff --git a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.H b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.H similarity index 66% rename from src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.H rename to src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.H index d9927e35cf8..0e49627d13e 100644 --- a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.H +++ b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/patchInjection/patchInjection.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -22,18 +22,19 @@ License along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. Class - Foam::removeInjection + Foam::patchInjection Description - All mass available to be removed from the system is removed. + Remove and inject the mass in the film as it passes over the selected + patches. SourceFiles - removeInjection.C + patchInjection.C \*---------------------------------------------------------------------------*/ -#ifndef removeInjection_H -#define removeInjection_H +#ifndef patchInjection_H +#define patchInjection_H #include "injectionModel.H" @@ -47,22 +48,20 @@ namespace surfaceFilmModels { /*---------------------------------------------------------------------------*\ - Class removeInjection Declaration + Class patchInjection Declaration \*---------------------------------------------------------------------------*/ -class removeInjection +class patchInjection : public injectionModel { -private: - // Private member functions //- Disallow default bitwise copy construct - removeInjection(const removeInjection&); + patchInjection(const patchInjection&); //- Disallow default bitwise assignment - void operator=(const removeInjection&); + void operator=(const patchInjection&); protected: @@ -71,37 +70,45 @@ protected: // this threhold value scalar deltaStable_; - //- Mask per cell to indicate whether mass can be removed - scalarField mask_; + //- List of patch IDs at which the film is removed + labelList patchIDs_; + + //- Injected mass for each patch at which the film is removed + scalarField patchInjectedMasses_; public: //- Runtime type information - TypeName("removeInjection"); + TypeName("patchInjection"); // Constructors //- Construct from surface film model - removeInjection(surfaceFilmModel& owner, const dictionary& dict); + patchInjection(surfaceFilmModel& owner, const dictionary& dict); //- Destructor - virtual ~removeInjection(); + virtual ~patchInjection(); // Member Functions - // Evolution - - //- Correct - virtual void correct - ( - scalarField& availableMass, - scalarField& massToInject, - scalarField& diameterToInject - ); + //- Correct + virtual void correct + ( + scalarField& availableMass, + scalarField& massToInject, + scalarField& diameterToInject + ); + + //- Accumulate the total mass injected for the patches into the + // scalarField provided + virtual void patchInjectedMassTotals + ( + scalarField& patchMasses + ) const; }; diff --git a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.C b/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.C deleted file mode 100644 index fe6b5730b34..00000000000 --- a/src/regionModels/surfaceFilmModels/submodels/kinematic/injectionModel/removeInjection/removeInjection.C +++ /dev/null @@ -1,121 +0,0 @@ -/*---------------------------------------------------------------------------*\ - ========= | - \\ / F ield | OpenFOAM: The Open Source CFD Toolbox - \\ / O peration | - \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation - \\/ M anipulation | -------------------------------------------------------------------------------- -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 "removeInjection.H" -#include "addToRunTimeSelectionTable.H" - -// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // - -namespace Foam -{ -namespace regionModels -{ -namespace surfaceFilmModels -{ - -// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // - -defineTypeNameAndDebug(removeInjection, 0); -addToRunTimeSelectionTable(injectionModel, removeInjection, dictionary); - -// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * // - -removeInjection::removeInjection -( - surfaceFilmModel& owner, - const dictionary& dict -) -: - injectionModel(type(), owner, dict), - deltaStable_(coeffDict_.lookupOrDefault<scalar>("deltaStable", 0.0)), - mask_(owner.regionMesh().nCells(), -1) -{ - wordReList patches; - if (coeffDict_.readIfPresent("patches", patches)) - { - Info<< " applying to patches:" << nl; - const polyBoundaryMesh& pbm = owner.regionMesh().boundaryMesh(); - const labelHashSet patchSet = pbm.patchSet(patches); - - forAllConstIter(labelHashSet, patchSet, iter) - { - label patchI = iter.key(); - const polyPatch& pp = pbm[patchI]; - UIndirectList<scalar>(mask_, pp.faceCells()) = 1.0; - - Info<< " " << pp.name() << endl; - } - } - else - { - Info<< " applying to all patches" << endl; - mask_ = 1.0; - } -} - - -// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * // - -removeInjection::~removeInjection() -{} - - -// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * // - -void removeInjection::correct -( - scalarField& availableMass, - scalarField& massToInject, - scalarField& diameterToInject -) -{ - const scalarField& delta = owner().delta(); - const scalarField& rho = owner().rho(); - const scalarField& magSf = owner().magSf(); - - forAll(delta, cellI) - { - if (mask_[cellI] > 0) - { - scalar ddelta = max(0.0, delta[cellI] - deltaStable_); - scalar dMass = ddelta*rho[cellI]*magSf[cellI]; - massToInject[cellI] += dMass; - availableMass[cellI] -= dMass; - - addToInjectedMass(dMass); - } - } - - injectionModel::correct(); -} - - -// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // - -} // End namespace surfaceFilmModels -} // End namespace regionModels -} // End namespace Foam - -// ************************************************************************* // -- GitLab