ENH: liquidFilm subModels: add force models

This commit is contained in:
sergio
2020-10-26 10:18:55 -07:00
committed by Sergio Ferraris
parent e30cbfc6db
commit f98633e427
11 changed files with 1082 additions and 5 deletions

View File

@ -13,9 +13,6 @@ KirchhoffShell/KirchhoffShell.C
derivedFvPatchFields/thermalShell/thermalShellFvPatchScalarField.C
derivedFvPatchFields/vibrationShell/vibrationShellFvPatchScalarField.C
liquidFilm/liquidFilmBase.C
liquidFilm/liquidFilmBaseNew.C
/* Sub-Model */
liquidFilm/subModels/kinematic/filmTurbulenceModel/filmTurbulenceModel/filmTurbulenceModel.C
@ -23,12 +20,17 @@ liquidFilm/subModels/kinematic/filmTurbulenceModel/filmTurbulenceModel/filmTurbu
liquidFilm/subModels/kinematic/filmTurbulenceModel/laminar/laminar.C
liquidFilm/subModels/kinematic/injectionModel/injectionModelList/injectionModelList.C
liquidFilm/subModels/kinematic/injectionModel/injectionModel/injectionModel.C
liquidFilm/subModels/kinematic/injectionModel/injectionModel/injectionModelNew.C
liquidFilm/subModels/kinematic/injectionModel/curvatureSeparation/curvatureSeparation.C
liquidFilm/subModels/kinematic/force/forceList/forceList.C
liquidFilm/subModels/kinematic/force/force/force.C
liquidFilm/subModels/kinematic/force/force/forceNew.C
liquidFilm/subModels/kinematic/force/contactAngleForces/contactAngleForce/contactAngleForce.C
liquidFilm/subModels/kinematic/force/contactAngleForces/perturbedTemperatureDependent/perturbedTemperatureDependentContactAngleForce.C
liquidFilm/subModels/filmSubModelBase.C
liquidFilm/liquidFilmBase.C

View File

@ -5,6 +5,7 @@ EXE_INC = \
-I$(LIB_SRC)/meshTools/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/thermophysicalProperties/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/specie/lnInclude \
-I$(LIB_SRC)/lagrangian/distributionModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/incompressible/lnInclude \
@ -18,4 +19,5 @@ LIB_LIBS = \
-lmeshTools \
-lthermophysicalProperties \
-lspecie \
-lfaOptions
-lfaOptions \
-ldistributionModels

View File

@ -0,0 +1,195 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "contactAngleForce.H"
#include "addToRunTimeSelectionTable.H"
#include "faCFD.H"
#include "unitConversion.H"
#include "meshWavePatchDistMethod.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(contactAngleForce, 0);
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
void contactAngleForce::initialise()
{
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
contactAngleForce::contactAngleForce
(
const word& typeName,
liquidFilmBase& film,
const dictionary& dict
)
:
force(typeName, film, dict),
Ccf_(coeffDict_.get<scalar>("Ccf")),
mask_
(
IOobject
(
typeName + ":fContactForceMask",
film.primaryMesh().time().timeName(),
film.primaryMesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
film.regionMesh(),
dimensionedScalar("mask", dimless, 1.0)
)
{
initialise();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
contactAngleForce::~contactAngleForce()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
tmp<faVectorMatrix> contactAngleForce::correct(areaVectorField& U)
{
tmp<areaVectorField> tForce
(
new areaVectorField
(
IOobject
(
typeName + ":contactForce",
film().primaryMesh().time().timeName(),
film().primaryMesh(),
IOobject::NO_READ,
IOobject::NO_WRITE
),
film().regionMesh(),
dimensionedVector(dimForce/dimDensity/dimArea, Zero)
)
);
vectorField& force = tForce.ref().primitiveFieldRef();
const labelUList& own = film().regionMesh().owner();
const labelUList& nbr = film().regionMesh().neighbour();
const DimensionedField<scalar, areaMesh>& magSf = film().regionMesh().S();
tmp<areaScalarField> talpha = film().alpha();
const areaScalarField& sigma = film().sigma();
const areaScalarField& rhof = film().rho();
tmp<areaScalarField> ttheta = theta();
const areaScalarField& theta = ttheta();
const areaVectorField gradAlpha(fac::grad(talpha()));
forAll(nbr, edgei)
{
const label faceO = own[edgei];
const label faceN = nbr[edgei];
label facei = -1;
if ((talpha()[faceO] > 0.5) && (talpha()[faceN] < 0.5))
{
facei = faceO;
}
else if ((talpha()[faceO] < 0.5) && (talpha()[faceN] > 0.5))
{
facei = faceN;
}
if (facei != -1 && mask_[facei] > 0.5)
{
const scalar invDx = film().regionMesh().deltaCoeffs()[edgei];
const vector n
(
gradAlpha[facei]/(mag(gradAlpha[facei]) + ROOTVSMALL)
);
const scalar cosTheta = cos(degToRad(theta[facei]));
force[facei] +=
Ccf_*n*sigma[facei]*(1 - cosTheta)/invDx/rhof[facei];
}
}
forAll(sigma.boundaryField(), patchi)
{
const faPatchField<scalar>& alphaPf = sigma.boundaryField()[patchi];
const faPatchField<scalar>& sigmaPf = sigma.boundaryField()[patchi];
const labelUList& faces = alphaPf.patch().edgeFaces();
forAll(sigmaPf, edgei)
{
label face0 = faces[edgei];
force[face0] = vector::zero;
}
}
force /= magSf.field();
if (film().regionMesh().time().writeTime())
{
tForce().write();
gradAlpha.write();
}
tmp<faVectorMatrix> tfvm
(
new faVectorMatrix(U, dimForce/dimDensity)
);
tfvm.ref() += tForce;
return tfvm;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,127 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::regionModels::areaSurfaceFilmModels::contactAngleForce
Description
Base-class for film contact angle force models.
The effect of the contact angle force can be ignored over a specified
distance from patches.
SourceFiles
contactAngleForce.C
\*---------------------------------------------------------------------------*/
#ifndef contactAngleForce_H
#define contactAngleForce_H
#include "force.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
/*---------------------------------------------------------------------------*\
Class contactAngleForce Declaration
\*---------------------------------------------------------------------------*/
class contactAngleForce
:
public force
{
// Private Data
//- Coefficient applied to the contact angle force
scalar Ccf_;
//- Mask over which force is applied
areaScalarField mask_;
// Private Member Functions
//- Initialise
void initialise();
//- No copy construct
contactAngleForce(const contactAngleForce&) = delete;
//- No copy assignment
void operator=(const contactAngleForce&) = delete;
protected:
//- Return the contact angle field
virtual tmp<areaScalarField> theta() const = 0;
public:
//- Runtime type information
TypeName("contactAngle");
// Constructors
//- Construct from surface film model
contactAngleForce
(
const word& typeName,
liquidFilmBase& film,
const dictionary& dict
);
//- Destructor
virtual ~contactAngleForce();
// Member Functions
//- Correct
virtual tmp<faVectorMatrix> correct(areaVectorField& U);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -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 "perturbedTemperatureDependentContactAngleForce.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(perturbedTemperatureDependentContactAngleForce, 0);
addToRunTimeSelectionTable
(
force,
perturbedTemperatureDependentContactAngleForce,
dictionary
);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
perturbedTemperatureDependentContactAngleForce::
perturbedTemperatureDependentContactAngleForce
(
liquidFilmBase& film,
const dictionary& dict
)
:
contactAngleForce(typeName, film, dict),
thetaPtr_(Function1<scalar>::New("theta", coeffDict_)),
rndGen_(label(0)),
distribution_
(
distributionModel::New
(
coeffDict_.subDict("distribution"),
rndGen_
)
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
perturbedTemperatureDependentContactAngleForce::
~perturbedTemperatureDependentContactAngleForce()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
tmp<areaScalarField> perturbedTemperatureDependentContactAngleForce::
theta() const
{
tmp<areaScalarField> ttheta
(
new areaScalarField
(
IOobject
(
typeName + ":theta",
film().primaryMesh().time().timeName(),
film().primaryMesh()
),
film().regionMesh(),
dimensionedScalar(dimless, Zero)
)
);
areaScalarField& theta = ttheta.ref();
scalarField& thetai = theta.ref();
const areaScalarField& T = film().Tf();
// Initialize with the function of temperature
thetai = thetaPtr_->value(T());
// Add the stochastic perturbation
forAll(thetai, facei)
{
thetai[facei] += distribution_->sample();
}
return ttheta;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace surfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,139 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::regionModels::areaSurfaceFilmModels::
perturbedTemperatureDependentContactAngleForce
Description
Temperature dependent contact angle force with a stochastic perturbation.
The contact angle in degrees is specified as a \c Foam::Function1 type,
to enable the use of, e.g. \c constant, \c polynomial, \c table values
and the stochastic perturbation obtained from a
\c Foam::distributionModels::distributionModel.
See also
- Foam::regionModels::areaSurfaceFilmModels::contactAngleForce
- areaSurfaceFilmModels::temperatureDependentContactAngleForce
- Foam::regionModels::areaSurfaceFilmModels::distributionContactAngleForce
- Foam::Function1Types
- Foam::distributionModel
SourceFiles
perturbedTemperatureDependentContactAngleForce.C
\*---------------------------------------------------------------------------*/
#ifndef perturbedTemperatureDependentContactAngleForce_H
#define perturbedTemperatureDependentContactAngleForce_H
#include "contactAngleForce.H"
#include "Function1.H"
#include "distributionModel.H"
#include "Random.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
/*---------------------------------------------------------------------------*\
Class perturbedTemperatureDependentContactAngleForce Declaration
\*---------------------------------------------------------------------------*/
class perturbedTemperatureDependentContactAngleForce
:
public contactAngleForce
{
// Private Data
//- Contact angle function
autoPtr<Function1<scalar>> thetaPtr_;
//- Random number generator
Random rndGen_;
//- Parcel size PDF model
const autoPtr<distributionModel> distribution_;
// Private Member Functions
//- No copy construct
perturbedTemperatureDependentContactAngleForce
(
const perturbedTemperatureDependentContactAngleForce&
) = delete;
//- No copy assignment
void operator=
(
const perturbedTemperatureDependentContactAngleForce&
) = delete;
protected:
//- Return the contact angle field
virtual tmp<areaScalarField> theta() const;
public:
//- Runtime type information
TypeName("perturbedTemperatureDependentContactAngle");
// Constructors
//- Construct from surface film model
perturbedTemperatureDependentContactAngleForce
(
liquidFilmBase& film,
const dictionary& dict
);
//- Destructor
virtual ~perturbedTemperatureDependentContactAngleForce();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,75 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "force.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(force, 0);
defineRunTimeSelectionTable(force, dictionary);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
force::force(liquidFilmBase& film)
:
filmSubModelBase(film)
{}
force::force
(
const word& modelType,
liquidFilmBase& film,
const dictionary& dict
)
:
filmSubModelBase(film, dict, typeName, modelType)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
force::~force()
{}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,138 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::regionModels::areaSurfaceFilmModels::force
Description
Base class for film (stress-based) force models
SourceFiles
force.C
forceNew.C
\*---------------------------------------------------------------------------*/
#ifndef force_H
#define force_H
#include "filmSubModelBase.H"
#include "runTimeSelectionTables.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
/*---------------------------------------------------------------------------*\
Class force Declaration
\*---------------------------------------------------------------------------*/
class force
:
public filmSubModelBase
{
// Private Member Functions
//- No copy construct
force(const force&) = delete;
//- No copy assignment
void operator=(const force&) = delete;
public:
//- Runtime type information
TypeName("force");
// Declare runtime constructor selection table
declareRunTimeSelectionTable
(
autoPtr,
force,
dictionary,
(
liquidFilmBase& film,
const dictionary& dict
),
(film, dict)
);
// Constructors
//- Construct null
force(liquidFilmBase& film);
//- Construct from type name, dictionary and surface film model
force
(
const word& modelType,
liquidFilmBase& film,
const dictionary& dict
);
// Selectors
//- Return a reference to the selected force model
static autoPtr<force> New
(
liquidFilmBase& film,
const dictionary& dict,
const word& modelType
);
//- Destructor
virtual ~force();
// Member Functions
// Evolution
//- Correct
virtual tmp<faVectorMatrix> correct(areaVectorField& U) = 0;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,73 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "force.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
autoPtr<force> force::New
(
liquidFilmBase& model,
const dictionary& dict,
const word& modelType
)
{
Info<< " " << modelType << endl;
auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
if (!cstrIter.found())
{
FatalIOErrorInLookup
(
dict,
"force",
modelType,
*dictionaryConstructorTablePtr_
) << exit(FatalIOError);
}
return autoPtr<force>(cstrIter()(model, dict));
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,104 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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 "forceList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
forceList::forceList(liquidFilmBase& film)
:
PtrList<force>()
{}
forceList::forceList
(
liquidFilmBase& film,
const dictionary& dict
)
:
PtrList<force>()
{
const wordList models(dict.lookup("forces"));
Info<< " Selecting film force models" << endl;
if (models.size() > 0)
{
this->setSize(models.size());
forAll(models, i)
{
set(i, force::New(film, dict, models[i]));
}
}
else
{
Info<< " none" << endl;
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
forceList::~forceList()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
tmp<faVectorMatrix> forceList::correct(areaVectorField& U)
{
tmp<faVectorMatrix> tResult
(
new faVectorMatrix(U, dimForce/dimDensity)
);
faVectorMatrix& result = tResult.ref();
forAll(*this, i)
{
result += this->operator[](i).correct(U);
}
return tResult;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// ************************************************************************* //

View File

@ -0,0 +1,98 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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::regionModels::surfaceFilmModels::forceList
Description
List container for film sources
SourceFiles
forceList.C
\*---------------------------------------------------------------------------*/
#ifndef forceList_H
#define forceList_H
#include "PtrList.H"
#include "force.H"
#include "faCFD.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace areaSurfaceFilmModels
{
/*---------------------------------------------------------------------------*\
Class forceList Declaration
\*---------------------------------------------------------------------------*/
class forceList
:
public PtrList<force>
{
public:
// Constructors
//- Construct null
forceList(liquidFilmBase& film);
//- Construct from type name, dictionary and surface film model
forceList
(
liquidFilmBase& film,
const dictionary& dict
);
//- Destructor
virtual ~forceList();
// Member Functions
//- Return (net) force system
tmp<faVectorMatrix> correct(areaVectorField& U);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace areaSurfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //