Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
openfoam
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Deploy
Releases
Model registry
Analyze
Contributor analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Development
openfoam
Commits
9774bf2d
Commit
9774bf2d
authored
16 years ago
by
Andrew Heather
Browse files
Options
Downloads
Patches
Plain Diff
updaing to include compressible case
parent
776f31b9
Branches
Branches containing commit
Tags
Tags containing commit
No related merge requests found
Changes
1
Hide whitespace changes
Inline
Side-by-side
Showing
1 changed file
applications/utilities/preProcessing/applyWallFunctionBounaryConditions/applyWallFunctionBounaryConditions.C
+122
-19
122 additions, 19 deletions
...ionBounaryConditions/applyWallFunctionBounaryConditions.C
with
122 additions
and
19 deletions
applications/utilities/preProcessing/applyWallFunctionBounaryConditions/applyWallFunctionBounaryConditions.C
+
122
−
19
View file @
9774bf2d
...
@@ -26,10 +26,9 @@ Application
...
@@ -26,10 +26,9 @@ Application
applyWallFunctionBounaryConditions
applyWallFunctionBounaryConditions
Description
Description
Updates OpenFOAM incompressible RAS cases to use the new wall function
Updates OpenFOAM RAS cases to use the new wall function framework
framework
Attempts to determine whether case is compressible or incompressible, or
can be supplied with -compressible command line argument
NOTE: For incompressible RAS calculations ONLY
\*---------------------------------------------------------------------------*/
\*---------------------------------------------------------------------------*/
...
@@ -38,6 +37,7 @@ Description
...
@@ -38,6 +37,7 @@ Description
#include
"fvMesh.H"
#include
"fvMesh.H"
#include
"Time.H"
#include
"Time.H"
#include
"volFields.H"
#include
"volFields.H"
#include
"surfaceFields.H"
#include
"wallPolyPatch.H"
#include
"wallPolyPatch.H"
...
@@ -45,39 +45,124 @@ using namespace Foam;
...
@@ -45,39 +45,124 @@ using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
bool
caseIsCompressible
(
const
fvMesh
&
mesh
)
{
// Attempt flux field
IOobject
phiHeader
(
"phi"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
MUST_READ
,
IOobject
::
NO_WRITE
);
// Main program:
if
(
phiHeader
.
headerOk
())
{
surfaceScalarField
phi
(
phiHeader
,
mesh
);
if
(
phi
.
dimensions
()
==
dimDensity
*
dimVelocity
*
dimArea
)
{
return
true
;
}
}
void
createNut
(
const
fvMesh
&
mesh
)
// Attempt density field
IOobject
rhoHeader
(
"rho"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
MUST_READ
,
IOobject
::
NO_WRITE
);
if
(
rhoHeader
.
headerOk
())
{
volScalarField
rho
(
rhoHeader
,
mesh
);
if
(
rho
.
dimensions
()
==
dimDensity
)
{
return
true
;
}
}
// Attempt pressure field
IOobject
pHeader
(
"p"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
MUST_READ
,
IOobject
::
NO_WRITE
);
if
(
pHeader
.
headerOk
())
{
volScalarField
p
(
pHeader
,
mesh
);
if
(
p
.
dimensions
()
==
dimMass
/
sqr
(
dimTime
)
/
dimLength
)
{
return
true
;
}
}
// Attempt hydrostatic pressure field
IOobject
pdHeader
(
"pd"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
MUST_READ
,
IOobject
::
NO_WRITE
);
if
(
pdHeader
.
headerOk
())
{
volScalarField
pd
(
pdHeader
,
mesh
);
if
(
pd
.
dimensions
()
==
dimMass
/
sqr
(
dimTime
)
/
dimLength
)
{
return
true
;
}
}
// If none of the above are true, assume that the case is incompressible
return
false
;
}
void
createVolScalarField
(
const
fvMesh
&
mesh
,
const
word
&
fieldName
,
const
dimensionSet
&
dims
)
{
{
IOobject
nut
Header
IOobject
field
Header
(
(
"nut"
,
fieldName
,
mesh
.
time
().
timeName
(),
mesh
.
time
().
timeName
(),
mesh
,
mesh
,
IOobject
::
MUST_READ
,
IOobject
::
MUST_READ
,
IOobject
::
NO_WRITE
IOobject
::
NO_WRITE
);
);
if
(
!
nut
Header
.
headerOk
())
if
(
!
field
Header
.
headerOk
())
{
{
Info
<<
"Creating field
nut"
<<
nl
<<
endl
;
Info
<<
"Creating field
"
<<
fieldName
<<
nl
<<
endl
;
volScalarField
nut
volScalarField
field
(
(
IOobject
IOobject
(
(
"nut"
,
fieldName
,
mesh
.
time
().
timeName
(),
mesh
.
time
().
timeName
(),
mesh
,
mesh
,
IOobject
::
NO_READ
,
IOobject
::
NO_READ
,
IOobject
::
NO_WRITE
IOobject
::
NO_WRITE
),
),
mesh
,
mesh
,
dimensionedScalar
(
"zero"
,
dim
ensionSet
(
0
,
2
,
-
1
,
0
,
0
)
,
0
.
0
)
dimensionedScalar
(
"zero"
,
dim
s
,
0
.
0
)
);
);
nut
.
write
();
field
.
write
();
}
}
}
}
...
@@ -116,7 +201,7 @@ void replaceBoundaryType
...
@@ -116,7 +201,7 @@ void replaceBoundaryType
// Make a backup of the old field
// Make a backup of the old field
word
backupName
(
dict
.
name
()
+
".old"
);
word
backupName
(
dict
.
name
()
+
".old"
);
Info
<<
" copying
original
"
<<
dict
.
name
()
<<
" to "
Info
<<
" copying "
<<
dict
.
name
()
<<
" to "
<<
backupName
<<
endl
;
<<
backupName
<<
endl
;
IOdictionary
dictOld
=
dict
;
IOdictionary
dictOld
=
dict
;
dictOld
.
rename
(
backupName
);
dictOld
.
rename
(
backupName
);
...
@@ -149,18 +234,36 @@ int main(int argc, char *argv[])
...
@@ -149,18 +234,36 @@ int main(int argc, char *argv[])
{
{
# include "addTimeOptions.H"
# include "addTimeOptions.H"
argList
::
validOptions
.
insert
(
"compressible"
,
""
);
# include "setRootCase.H"
# include "setRootCase.H"
# include "createTime.H"
# include "createTime.H"
# include "createMesh.H"
# include "createMesh.H"
bool
compressible
=
args
.
options
().
found
(
"compressible"
);
Info
<<
"Updating turbulence fields to operate using new run time "
Info
<<
"Updating turbulence fields to operate using new run time "
<<
"selectable"
<<
nl
<<
"wall functions"
<<
nl
<<
nl
<<
"selectable"
<<
nl
<<
"wall functions"
<<
">>>>NOTE: only applicable to incompressible RAS models"
<<
nl
<<
endl
;
<<
nl
<<
endl
;
createNut
(
mesh
);
if
(
compressible
||
caseIsCompressible
(
mesh
))
{
Info
<<
"Case treated as compressible"
<<
nl
<<
endl
;
createVolScalarField
(
mesh
,
"mut"
,
dimArea
/
dimTime
*
dimDensity
);
replaceBoundaryType
(
mesh
,
"mut"
,
"mutWallFunction"
,
"0"
);
}
else
{
Info
<<
"Case treated as incompressible"
<<
nl
<<
endl
;
createVolScalarField
(
mesh
,
"nut"
,
dimArea
/
dimTime
);
replaceBoundaryType
(
mesh
,
"nut"
,
"nutWallFunction"
,
"0"
);
}
replaceBoundaryType
(
mesh
,
"nut"
,
"nutWallFunction"
,
"0"
);
replaceBoundaryType
(
mesh
,
"epsilon"
,
"epsilonWallFunction"
,
"0"
);
replaceBoundaryType
(
mesh
,
"epsilon"
,
"epsilonWallFunction"
,
"0"
);
replaceBoundaryType
(
mesh
,
"omega"
,
"omegaWallFunction"
,
"0"
);
replaceBoundaryType
(
mesh
,
"omega"
,
"omegaWallFunction"
,
"0"
);
replaceBoundaryType
(
mesh
,
"k"
,
"kQRWallFunction"
,
"0"
);
replaceBoundaryType
(
mesh
,
"k"
,
"kQRWallFunction"
,
"0"
);
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
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!
Save comment
Cancel
Please
register
or
sign in
to comment