Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
77b7beb6
Commit
77b7beb6
authored
Jul 13, 2017
by
Henry Weller
Committed by
Andrew Heather
Jul 13, 2017
Browse files
swirlFlowRateInletVelocity: Added support for specifying the origin and axis of rotation
parent
065bfa26
Changes
2
Hide whitespace changes
Inline
Side-by-side
src/finiteVolume/fields/fvPatchFields/derived/swirlFlowRateInletVelocity/swirlFlowRateInletVelocityFvPatchVectorField.C
View file @
77b7beb6
...
...
@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-201
6
OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-201
7
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -42,6 +42,8 @@ swirlFlowRateInletVelocityFvPatchVectorField
fixedValueFvPatchField
<
vector
>
(
p
,
iF
),
phiName_
(
"phi"
),
rhoName_
(
"rho"
),
origin_
(),
axis_
(
Zero
),
flowRate_
(),
rpm_
()
{}
...
...
@@ -50,33 +52,51 @@ swirlFlowRateInletVelocityFvPatchVectorField
Foam
::
swirlFlowRateInletVelocityFvPatchVectorField
::
swirlFlowRateInletVelocityFvPatchVectorField
(
const
swirlFlowRateInletVelocityFvPatchVectorField
&
ptf
,
const
fvPatch
&
p
,
const
DimensionedField
<
vector
,
volMesh
>&
iF
,
const
fvPatchFieldMapper
&
mapper
const
dictionary
&
dict
)
:
fixedValueFvPatchField
<
vector
>
(
ptf
,
p
,
iF
,
mapper
),
phiName_
(
ptf
.
phiName_
),
rhoName_
(
ptf
.
rhoName_
),
flowRate_
(
ptf
.
flowRate_
,
false
),
rpm_
(
ptf
.
rpm_
,
false
)
fixedValueFvPatchField
<
vector
>
(
p
,
iF
,
dict
),
phiName_
(
dict
.
lookupOrDefault
<
word
>
(
"phi"
,
"phi"
)),
rhoName_
(
dict
.
lookupOrDefault
<
word
>
(
"rho"
,
"rho"
)),
origin_
(
dict
.
lookupOrDefault
(
"origin"
,
gSum
(
patch
().
Cf
()
*
patch
().
magSf
())
/
gSum
(
patch
().
magSf
())
)
),
axis_
(
dict
.
lookupOrDefault
(
"axis"
,
-
gSum
(
patch
().
Sf
())
/
gSum
(
patch
().
magSf
())
)
),
flowRate_
(
Function1
<
scalar
>::
New
(
"flowRate"
,
dict
)),
rpm_
(
Function1
<
scalar
>::
New
(
"rpm"
,
dict
))
{}
Foam
::
swirlFlowRateInletVelocityFvPatchVectorField
::
swirlFlowRateInletVelocityFvPatchVectorField
(
const
swirlFlowRateInletVelocityFvPatchVectorField
&
ptf
,
const
fvPatch
&
p
,
const
DimensionedField
<
vector
,
volMesh
>&
iF
,
const
dictionary
&
dict
const
fvPatchFieldMapper
&
mapper
)
:
fixedValueFvPatchField
<
vector
>
(
p
,
iF
,
dict
),
phiName_
(
dict
.
lookupOrDefault
<
word
>
(
"phi"
,
"phi"
)),
rhoName_
(
dict
.
lookupOrDefault
<
word
>
(
"rho"
,
"rho"
)),
flowRate_
(
Function1
<
scalar
>::
New
(
"flowRate"
,
dict
)),
rpm_
(
Function1
<
scalar
>::
New
(
"rpm"
,
dict
))
fixedValueFvPatchField
<
vector
>
(
ptf
,
p
,
iF
,
mapper
),
phiName_
(
ptf
.
phiName_
),
rhoName_
(
ptf
.
rhoName_
),
origin_
(
ptf
.
origin_
),
axis_
(
ptf
.
axis_
),
flowRate_
(
ptf
.
flowRate_
,
false
),
rpm_
(
ptf
.
rpm_
,
false
)
{}
...
...
@@ -89,6 +109,8 @@ swirlFlowRateInletVelocityFvPatchVectorField
fixedValueFvPatchField
<
vector
>
(
ptf
),
phiName_
(
ptf
.
phiName_
),
rhoName_
(
ptf
.
rhoName_
),
origin_
(
ptf
.
origin_
),
axis_
(
ptf
.
axis_
),
flowRate_
(
ptf
.
flowRate_
,
false
),
rpm_
(
ptf
.
rpm_
,
false
)
{}
...
...
@@ -104,6 +126,8 @@ swirlFlowRateInletVelocityFvPatchVectorField
fixedValueFvPatchField
<
vector
>
(
ptf
,
iF
),
phiName_
(
ptf
.
phiName_
),
rhoName_
(
ptf
.
rhoName_
),
origin_
(
ptf
.
origin_
),
axis_
(
ptf
.
axis_
),
flowRate_
(
ptf
.
flowRate_
,
false
),
rpm_
(
ptf
.
rpm_
,
false
)
{}
...
...
@@ -122,18 +146,16 @@ void Foam::swirlFlowRateInletVelocityFvPatchVectorField::updateCoeffs()
const
scalar
flowRate
=
flowRate_
->
value
(
t
);
const
scalar
rpm
=
rpm_
->
value
(
t
);
const
scalar
totArea
=
gSum
(
patch
().
magSf
());
const
scalar
totArea
=
gSum
(
patch
().
magSf
());
const
scalar
avgU
=
-
flowRate
/
totArea
;
const
vector
avgCenter
=
gSum
(
patch
().
Cf
()
*
patch
().
magSf
())
/
totArea
;
const
vector
avgNormal
=
gSum
(
patch
().
Sf
())
/
totArea
;
const
vector
axisHat
=
axis_
/
mag
(
axis_
);
// Update angular velocity - convert [rpm] to [rad/s]
tmp
<
vectorField
>
tangentialVelocity
(
(
rpm
*
constant
::
mathematical
::
pi
/
30
.
0
)
*
(
patch
().
Cf
()
-
avgCenter
)
^
avgNormal
);
(
axisHat
^
(
rpm
*
constant
::
mathematical
::
pi
/
30
.
0
)
*
(
patch
().
Cf
()
-
origin_
)
);
tmp
<
vectorField
>
n
=
patch
().
nf
();
...
...
@@ -175,6 +197,8 @@ void Foam::swirlFlowRateInletVelocityFvPatchVectorField::write
fvPatchField
<
vector
>::
write
(
os
);
writeEntryIfDifferent
<
word
>
(
os
,
"phi"
,
"phi"
,
phiName_
);
writeEntryIfDifferent
<
word
>
(
os
,
"rho"
,
"rho"
,
rhoName_
);
os
.
writeKeyword
(
"origin"
)
<<
origin_
<<
token
::
END_STATEMENT
<<
nl
;
os
.
writeKeyword
(
"axis"
)
<<
axis_
<<
token
::
END_STATEMENT
<<
nl
;
flowRate_
->
writeData
(
os
);
rpm_
->
writeData
(
os
);
writeEntry
(
"value"
,
os
);
...
...
src/finiteVolume/fields/fvPatchFields/derived/swirlFlowRateInletVelocity/swirlFlowRateInletVelocityFvPatchVectorField.H
View file @
77b7beb6
...
...
@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-201
6
OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-201
7
OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
...
...
@@ -42,6 +42,8 @@ Usage
Property | Description | Required | Default value
phi | flux field name | no | phi
rho | density field name | no | rho
origin | origin of rotation | no | patch centre
axis | axis of rotation | no | -(patch normal)
flowRate | flow rate profile | yes |
rpm | rotational speed profile | yes |
\endtable
...
...
@@ -97,6 +99,12 @@ class swirlFlowRateInletVelocityFvPatchVectorField
//- Name of the density field used to normalize the mass flux
const
word
rhoName_
;
//- Origin of the rotation
const
vector
origin_
;
//- Axis of the rotation
const
vector
axis_
;
//- Inlet integral flow rate
autoPtr
<
Function1
<
scalar
>>
flowRate_
;
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a 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