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
0f8d7e67
Commit
0f8d7e67
authored
Jun 25, 2015
by
Henry Weller
Browse files
reactingTwoPhaseEulerFoam: Construct MRF and fvOptions in phaseSystem
parent
550ba440
Changes
15
Hide whitespace changes
Inline
Side-by-side
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/createMRF.H
deleted
100644 → 0
View file @
550ba440
IOMRFZoneList
MRF
(
mesh
);
MRF
.
correctBoundaryVelocity
(
U1
);
MRF
.
correctBoundaryVelocity
(
U2
);
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/generatePairsAndSubModels
deleted
100644 → 0
View file @
550ba440
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseCompressibleTurbulenceModels/Make/options
View file @
0f8d7e67
...
...
@@ -8,4 +8,5 @@ EXE_INC = \
-I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
-I$(LIB_SRC)/TurbulenceModels/phaseCompressible/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \
-I$(LIB_SRC)/fvOptions/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/PhaseSystems/HeatAndMassTransferPhaseSystem/HeatAndMassTransferPhaseSystem.C
View file @
0f8d7e67
...
...
@@ -171,13 +171,10 @@ Foam::HeatAndMassTransferPhaseSystem<BasePhaseSystem>::dmdt
template
<
class
BasePhaseSystem
>
Foam
::
autoPtr
<
Foam
::
phaseSystem
::
momentumTransferTable
>
Foam
::
HeatAndMassTransferPhaseSystem
<
BasePhaseSystem
>::
momentumTransfer
(
IOMRFZoneList
&
MRF
)
const
Foam
::
HeatAndMassTransferPhaseSystem
<
BasePhaseSystem
>::
momentumTransfer
()
const
{
autoPtr
<
phaseSystem
::
momentumTransferTable
>
eqnsPtr
(
BasePhaseSystem
::
momentumTransfer
(
MRF
));
eqnsPtr
(
BasePhaseSystem
::
momentumTransfer
());
phaseSystem
::
momentumTransferTable
&
eqns
=
eqnsPtr
();
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/PhaseSystems/HeatAndMassTransferPhaseSystem/HeatAndMassTransferPhaseSystem.H
View file @
0f8d7e67
...
...
@@ -142,10 +142,8 @@ public:
)
const
;
//- Return the momentum transfer matrices
virtual
autoPtr
<
phaseSystem
::
momentumTransferTable
>
momentumTransfer
(
IOMRFZoneList
&
mrfZones
)
const
;
virtual
autoPtr
<
phaseSystem
::
momentumTransferTable
>
momentumTransfer
()
const
;
//- Return the heat transfer matrices
virtual
autoPtr
<
phaseSystem
::
heatTransferTable
>
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/PhaseSystems/MomentumTransferPhaseSystem/MomentumTransferPhaseSystem.C
View file @
0f8d7e67
...
...
@@ -32,7 +32,6 @@ License
#include
"wallLubricationModel.H"
#include
"turbulentDispersionModel.H"
#include
"IOMRFZoneList.H"
#include
"HashPtrTable.H"
#include
"fvmDdt.H"
...
...
@@ -348,10 +347,7 @@ Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::D
template
<
class
BasePhaseSystem
>
Foam
::
autoPtr
<
Foam
::
phaseSystem
::
momentumTransferTable
>
Foam
::
MomentumTransferPhaseSystem
<
BasePhaseSystem
>::
momentumTransfer
(
IOMRFZoneList
&
MRF
)
const
Foam
::
MomentumTransferPhaseSystem
<
BasePhaseSystem
>::
momentumTransfer
()
const
{
// Create a momentum transfer matrix for each phase
autoPtr
<
phaseSystem
::
momentumTransferTable
>
eqnsPtr
...
...
@@ -452,7 +448,7 @@ Foam::MomentumTransferPhaseSystem<BasePhaseSystem>::momentumTransfer
-
fvm
::
Sp
(
fvc
::
div
(
phi
),
U
)
-
otherPhase
->
DUDt
()
)
-
MRF
.
DDt
(
Vm
,
U
-
otherPhase
->
U
());
-
this
->
MRF
_
.
DDt
(
Vm
,
U
-
otherPhase
->
U
());
Swap
(
phase
,
otherPhase
);
}
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/PhaseSystems/MomentumTransferPhaseSystem/MomentumTransferPhaseSystem.H
View file @
0f8d7e67
...
...
@@ -171,10 +171,8 @@ public:
virtual
tmp
<
volScalarField
>
D
(
const
phasePairKey
&
key
)
const
;
//- Return the momentum transfer matrices
virtual
autoPtr
<
phaseSystem
::
momentumTransferTable
>
momentumTransfer
(
IOMRFZoneList
&
mrfZones
)
const
;
virtual
autoPtr
<
phaseSystem
::
momentumTransferTable
>
momentumTransfer
()
const
;
//- Read base phaseProperties dictionary
virtual
bool
read
();
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/phaseModel/MovingPhaseModel/MovingPhaseModel.C
View file @
0f8d7e67
...
...
@@ -42,12 +42,10 @@ License
template
<
class
BasePhaseModel
>
Foam
::
tmp
<
Foam
::
surfaceScalarField
>
Foam
::
MovingPhaseModel
<
BasePhaseModel
>::
generatePhi
(
const
word
&
phiName
,
const
volVectorField
&
U
)
const
Foam
::
MovingPhaseModel
<
BasePhaseModel
>::
phi
(
const
volVectorField
&
U
)
const
{
word
phiName
(
IOobject
::
groupName
(
"phi"
,
this
->
name
()));
IOobject
phiHeader
(
phiName
,
...
...
@@ -60,22 +58,21 @@ Foam::MovingPhaseModel<BasePhaseModel>::generatePhi
{
Info
<<
"Reading face flux field "
<<
phiName
<<
endl
;
return
tmp
<
surfaceScalarField
>
return
tmp
<
surfaceScalarField
>
(
new
surfaceScalarField
(
new
surfaceScalarField
IOobject
(
IOobject
(
phiName
,
U
.
mesh
().
time
().
timeName
(),
U
.
mesh
(),
IOobject
::
MUST_READ
,
IOobject
::
AUTO_WRITE
),
U
.
mesh
()
)
);
phiName
,
U
.
mesh
().
time
().
timeName
(),
U
.
mesh
(),
IOobject
::
MUST_READ
,
IOobject
::
AUTO_WRITE
),
U
.
mesh
()
)
);
}
else
{
...
...
@@ -142,14 +139,7 @@ Foam::MovingPhaseModel<BasePhaseModel>::MovingPhaseModel
),
fluid
.
mesh
()
),
phi_
(
generatePhi
(
IOobject
::
groupName
(
"phi"
,
this
->
name
()),
U_
)
),
phi_
(
phi
(
U_
)),
alphaPhi_
(
IOobject
...
...
@@ -225,6 +215,8 @@ void Foam::MovingPhaseModel<BasePhaseModel>::correct()
{
BasePhaseModel
::
correct
();
this
->
fluid
().
MRF
().
correctBoundaryVelocity
(
U_
);
continuityError_
=
fvc
::
ddt
(
*
this
,
this
->
thermo
().
rho
())
+
fvc
::
div
(
alphaRhoPhi_
);
}
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/phaseModel/MovingPhaseModel/MovingPhaseModel.H
View file @
0f8d7e67
...
...
@@ -91,12 +91,8 @@ class MovingPhaseModel
// Private static member functions
//- Generate the flux field
tmp
<
surfaceScalarField
>
generatePhi
(
const
word
&
phiName
,
const
volVectorField
&
U
)
const
;
//- Calculate and return the flux field
tmp
<
surfaceScalarField
>
phi
(
const
volVectorField
&
U
)
const
;
public:
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/phaseSystem/phaseSystem.C
View file @
0f8d7e67
...
...
@@ -183,7 +183,10 @@ Foam::phaseSystem::phaseSystem
),
mesh
,
dimensionedScalar
(
"dpdt"
,
dimPressure
/
dimTime
,
0
)
)
),
MRF_
(
mesh_
),
fvOptions_
(
mesh_
)
{
// Blending methods
forAllConstIter
(
dictionary
,
subDict
(
"blending"
),
iter
)
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/phaseSystem/phaseSystem.H
View file @
0f8d7e67
...
...
@@ -41,6 +41,8 @@ SourceFiles
#include
"phaseModel.H"
#include
"phasePair.H"
#include
"orderedPhasePair.H"
#include
"IOMRFZoneList.H"
#include
"fvIOoptionList.H"
#include
"volFields.H"
#include
"surfaceFields.H"
...
...
@@ -55,7 +57,6 @@ class blendingMethod;
template
<
class
modelType
>
class
BlendedInterfacialModel
;
class
surfaceTensionModel
;
class
aspectRatioModel
;
class
IOMRFZoneList
;
/*---------------------------------------------------------------------------*\
Class phaseSystem Declaration
...
...
@@ -174,6 +175,12 @@ protected:
//- Rate of change of pressure
volScalarField
dpdt_
;
//- Optional MRF zones
IOMRFZoneList
MRF_
;
//- Optional FV-options
fv
::
IOoptionList
fvOptions_
;
//- Blending methods
blendingMethodTable
blendingMethods_
;
...
...
@@ -316,10 +323,7 @@ public:
virtual
tmp
<
volScalarField
>
dmdt
(
const
phasePairKey
&
key
)
const
=
0
;
//- Return the momentum transfer matrices
virtual
autoPtr
<
momentumTransferTable
>
momentumTransfer
(
IOMRFZoneList
&
MRF
)
const
=
0
;
virtual
autoPtr
<
momentumTransferTable
>
momentumTransfer
()
const
=
0
;
//- Return the heat transfer matrices
virtual
autoPtr
<
heatTransferTable
>
heatTransfer
()
const
=
0
;
...
...
@@ -363,6 +367,12 @@ public:
//- Access the rate of change of the pressure
inline
volScalarField
&
dpdt
();
//- Return MRF zones
inline
const
IOMRFZoneList
&
MRF
()
const
;
//- Optional FV-options
inline
fv
::
IOoptionList
&
fvOptions
();
//- Access a sub model between a phase pair
template
<
class
modelType
>
const
modelType
&
lookupSubModel
(
const
phasePair
&
key
)
const
;
...
...
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/phaseSystems/phaseSystem/phaseSystemI.H
View file @
0f8d7e67
...
...
@@ -55,4 +55,16 @@ inline Foam::volScalarField& Foam::phaseSystem::dpdt()
}
inline
const
Foam
::
IOMRFZoneList
&
Foam
::
phaseSystem
::
MRF
()
const
{
return
MRF_
;
}
inline
Foam
::
fv
::
IOoptionList
&
Foam
::
phaseSystem
::
fvOptions
()
{
return
fvOptions_
;
}
// ************************************************************************* //
applications/solvers/multiphase/reactingTwoPhaseEulerFoam/reactingTwoPhaseEulerFoam.C
View file @
0f8d7e67
...
...
@@ -37,7 +37,6 @@ Description
#include
"twoPhaseSystem.H"
#include
"PhaseCompressibleTurbulenceModel.H"
#include
"pimpleControl.H"
#include
"fvIOoptionList.H"
#include
"fixedFluxPressureFvPatchScalarField.H"
#include
"HashPtrTable.H"
...
...
@@ -53,8 +52,6 @@ int main(int argc, char *argv[])
pimpleControl
pimple
(
mesh
);
#include
"createFields.H"
#include
"createMRF.H"
#include
"createFvOptions.H"
#include
"initContinuityErrs.H"
#include
"readTimeControls.H"
#include
"CourantNos.H"
...
...
@@ -75,6 +72,9 @@ int main(int argc, char *argv[])
#include
"pUf/createDDtU.H"
const
IOMRFZoneList
&
MRF
=
fluid
.
MRF
();
fv
::
IOoptionList
&
fvOptions
=
fluid
.
fvOptions
();
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info
<<
"
\n
Starting time loop
\n
"
<<
endl
;
...
...
src/finiteVolume/cfdTools/general/MRF/MRFZoneList.C
View file @
0f8d7e67
...
...
@@ -164,7 +164,7 @@ void Foam::MRFZoneList::addAcceleration
Foam
::
tmp
<
Foam
::
volVectorField
>
Foam
::
MRFZoneList
::
DDt
(
const
volVectorField
&
U
)
)
const
{
tmp
<
volVectorField
>
tacceleration
(
...
...
@@ -195,7 +195,7 @@ Foam::tmp<Foam::volVectorField> Foam::MRFZoneList::DDt
(
const
volScalarField
&
rho
,
const
volVectorField
&
U
)
)
const
{
return
rho
*
DDt
(
U
);
}
...
...
src/finiteVolume/cfdTools/general/MRF/MRFZoneList.H
View file @
0f8d7e67
...
...
@@ -114,14 +114,14 @@ public:
tmp
<
volVectorField
>
DDt
(
const
volVectorField
&
U
);
)
const
;
//- Return the frame acceleration
tmp
<
volVectorField
>
DDt
(
const
volScalarField
&
rho
,
const
volVectorField
&
U
);
)
const
;
//- Make the given absolute velocity relative within the MRF region
void
makeRelative
(
volVectorField
&
U
)
const
;
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment