Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
475e8857
Commit
475e8857
authored
Jan 02, 2014
by
william
Browse files
ENH: Abstracted and made run-time selectable the lift models in twoPhaseEulerFoam
parent
d3539642
Changes
15
Hide whitespace changes
Inline
Side-by-side
applications/solvers/multiphase/twoPhaseEulerFoam/interfacialModels/Make/files
View file @
475e8857
...
...
@@ -12,4 +12,8 @@ heatTransferModels/heatTransferModel/heatTransferModel.C
heatTransferModels/heatTransferModel/newHeatTransferModel.C
heatTransferModels/RanzMarshall/RanzMarshall.C
liftModels/liftModel/liftModel.C
liftModels/liftModel/newLiftModel.C
liftModels/constantCoefficient/constantCoefficient.C
LIB = $(FOAM_LIBBIN)/libcompressibleEulerianInterfacialModels
applications/solvers/multiphase/twoPhaseEulerFoam/interfacialModels/liftModels/constantCoefficient/constantCoefficient.C
0 → 100644
View file @
475e8857
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 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
"constantCoefficient.H"
#include
"addToRunTimeSelectionTable.H"
#include
"fvc.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace
Foam
{
namespace
liftModels
{
defineTypeNameAndDebug
(
constantCoefficient
,
0
);
addToRunTimeSelectionTable
(
liftModel
,
constantCoefficient
,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam
::
liftModels
::
constantCoefficient
::
constantCoefficient
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
)
:
liftModel
(
dict
,
alpha1
,
phase1
,
phase2
),
Cl_
(
"Cl"
,
dimless
,
dict
.
lookup
(
"Cl"
))
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam
::
liftModels
::
constantCoefficient
::~
constantCoefficient
()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam
::
tmp
<
Foam
::
volVectorField
>
Foam
::
liftModels
::
constantCoefficient
::
F
(
const
volVectorField
&
U
)
const
{
return
Cl_
*
(
phase1_
*
phase1_
.
rho
()
+
phase2_
*
phase2_
.
rho
())
*
(
(
phase1_
.
U
()
-
phase2_
.
U
())
^
fvc
::
curl
(
U
)
);
}
// ************************************************************************* //
applications/solvers/multiphase/twoPhaseEulerFoam/interfacialModels/liftModels/constantCoefficient/constantCoefficient.H
0 → 100644
View file @
475e8857
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 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/>.
Class
Foam::liftModels::constantCoefficient
Description
SourceFiles
constantCoefficient.C
\*---------------------------------------------------------------------------*/
#ifndef constantCoefficient_H
#define constantCoefficient_H
#include
"liftModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace
Foam
{
namespace
liftModels
{
/*---------------------------------------------------------------------------*\
Class constantCoefficient Declaration
\*---------------------------------------------------------------------------*/
class
constantCoefficient
:
public
liftModel
{
// Private data
//- Constant lift coefficient
dimensionedScalar
Cl_
;
public:
//- Runtime type information
TypeName
(
"constantCoefficient"
);
// Constructors
//- Construct from components
constantCoefficient
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
);
//- Destructor
virtual
~
constantCoefficient
();
// Member Functions
//- Lift force
tmp
<
volVectorField
>
F
(
const
volVectorField
&
U
)
const
;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
}
// End namespace liftModels
}
// End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
applications/solvers/multiphase/twoPhaseEulerFoam/interfacialModels/liftModels/liftModel/liftModel.C
0 → 100644
View file @
475e8857
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 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
"liftModel.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace
Foam
{
defineTypeNameAndDebug
(
liftModel
,
0
);
defineRunTimeSelectionTable
(
liftModel
,
dictionary
);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam
::
liftModel
::
liftModel
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
)
:
dict_
(
dict
),
alpha1_
(
alpha1
),
phase1_
(
phase1
),
phase2_
(
phase2
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam
::
liftModel
::~
liftModel
()
{}
// ************************************************************************* //
applications/solvers/multiphase/twoPhaseEulerFoam/interfacialModels/liftModels/liftModel/liftModel.H
0 → 100644
View file @
475e8857
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 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/>.
Class
Foam::liftModel
Description
SourceFiles
liftModel.C
newLiftModel.C
\*---------------------------------------------------------------------------*/
#ifndef liftModel_H
#define liftModel_H
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#include
"dictionary.H"
#include
"phaseModel.H"
#include
"runTimeSelectionTables.H"
namespace
Foam
{
/*---------------------------------------------------------------------------*\
Class liftModel Declaration
\*---------------------------------------------------------------------------*/
class
liftModel
{
protected:
// Protected data
const
dictionary
&
dict_
;
const
volScalarField
&
alpha1_
;
const
phaseModel
&
phase1_
;
const
phaseModel
&
phase2_
;
public:
//- Runtime type information
TypeName
(
"liftModel"
);
// Declare runtime construction
declareRunTimeSelectionTable
(
autoPtr
,
liftModel
,
dictionary
,
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
),
(
dict
,
alpha1
,
phase1
,
phase2
)
);
// Constructors
liftModel
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
);
//- Destructor
virtual
~
liftModel
();
// Selectors
static
autoPtr
<
liftModel
>
New
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
);
// Member Functions
//- Lift force
virtual
tmp
<
volVectorField
>
F
(
const
volVectorField
&
U
)
const
=
0
;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
}
// End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
applications/solvers/multiphase/twoPhaseEulerFoam/interfacialModels/liftModels/liftModel/newLiftModel.C
0 → 100644
View file @
475e8857
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 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
"liftModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Foam
::
autoPtr
<
Foam
::
liftModel
>
Foam
::
liftModel
::
New
(
const
dictionary
&
dict
,
const
volScalarField
&
alpha1
,
const
phaseModel
&
phase1
,
const
phaseModel
&
phase2
)
{
word
liftModelType
(
dict
.
subDict
(
phase1
.
name
()).
lookup
(
"type"
)
);
Info
<<
"Selecting liftModel for phase "
<<
phase1
.
name
()
<<
": "
<<
liftModelType
<<
endl
;
dictionaryConstructorTable
::
iterator
cstrIter
=
dictionaryConstructorTablePtr_
->
find
(
liftModelType
);
if
(
cstrIter
==
dictionaryConstructorTablePtr_
->
end
())
{
FatalErrorIn
(
"liftModel::New"
)
<<
"Unknown liftModelType type "
<<
liftModelType
<<
endl
<<
endl
<<
"Valid liftModel types are : "
<<
endl
<<
dictionaryConstructorTablePtr_
->
sortedToc
()
<<
exit
(
FatalError
);
}
return
cstrIter
()
(
dict
.
subDict
(
phase1
.
name
()).
subDict
(
liftModelType
+
"Coeffs"
),
alpha1
,
phase1
,
phase2
);
}
// ************************************************************************* //
applications/solvers/multiphase/twoPhaseEulerFoam/twoPhaseSystem/twoPhaseSystem.C
View file @
475e8857
...
...
@@ -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
-2014
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -111,13 +111,6 @@ Foam::twoPhaseSystem::twoPhaseSystem
lookup
(
"Cvm"
)
),
Cl_
(
"Cl"
,
dimless
,
lookup
(
"Cl"
)
),
drag1_
(
dragModel
::
New
...
...
@@ -162,6 +155,28 @@ Foam::twoPhaseSystem::twoPhaseSystem
)
),
lift1_
(
liftModel
::
New
(
subDict
(
"lift"
),
phase1_
,
phase1_
,
phase2_
)
),
lift2_
(
liftModel
::
New
(
subDict
(
"lift"
),
phase2_
,
phase2_
,
phase1_
)
),
dispersedPhase_
(
lookup
(
"dispersedPhase"
)),
residualPhaseFraction_
...
...
@@ -311,11 +326,28 @@ Foam::tmp<Foam::volVectorField> Foam::twoPhaseSystem::liftForce
);
volVectorField
&
liftForce
=
tliftForce
();
volVectorField
Ur
(
phase1_
.
U
()
-
phase2_
.
U
());
liftForce
=
Cl_
*
(
phase1_
*
phase1_
.
rho
()
+
phase2_
*
phase2_
.
rho
())
*
(
Ur
^
fvc
::
curl
(
U
));
if
(
dispersedPhase_
==
phase1_
.
name
())
{
liftForce
=
lift1
().
F
(
U
);
}
else
if
(
dispersedPhase_
==
phase2_
.
name
())
{
liftForce
=
lift2
().
F
(
U
);
}
else
if
(
dispersedPhase_
==
"both"
)
{
liftForce
=
(
phase2_
*
lift1
().
F
(
U
)
+
phase1_
*
lift2
().
F
(
U
)
);
}
else
{
FatalErrorIn
(
"twoPhaseSystem::liftForce()"
)
<<
"dispersedPhase: "
<<
dispersedPhase_
<<
" is incorrect"
<<
exit
(
FatalError
);
}
// Remove lift at fixed-flux boundaries
forAll
(
phase1_
.
phi
().
boundaryField
(),
patchi
)
...
...
@@ -631,7 +663,6 @@ bool Foam::twoPhaseSystem::read()
lookup
(
"sigma"
)
>>
sigma_
;
lookup
(
"Cvm"
)
>>
Cvm_
;
lookup
(
"Cl"
)
>>
Cl_
;
// drag1_->read(*this);
// drag2_->read(*this);
...
...
applications/solvers/multiphase/twoPhaseEulerFoam/twoPhaseSystem/twoPhaseSystem.H
View file @
475e8857
...
...
@@ -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
-2014
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -45,6 +45,7 @@ SourceFiles
#include
"IOdictionary.H"
#include
"phaseModel.H"
#include
"dragModel.H"
#include
"liftModel.H"
#include
"heatTransferModel.H"
#include
"volFields.H"
#include
"surfaceFields.H"
...
...
@@ -57,7 +58,7 @@ namespace Foam
// Forward declarations
class
dragModel
;
class
heatTransferModel
;
class
liftModel
;
/*---------------------------------------------------------------------------*\
Class twoPhaseSystem Declaration
...
...
@@ -69,31 +70,57 @@ class twoPhaseSystem
{
// Private data
//- Reference to the mesh
const
fvMesh
&
mesh_
;