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-plus
Commits
cec56de0
Commit
cec56de0
authored
Feb 16, 2015
by
Henry
Browse files
Compressible solvers: rhoMax and rhoMin now optional and may be supplied without dimensions
parent
cd5ade6c
Changes
6
Hide whitespace changes
Inline
Side-by-side
applications/solvers/compressible/rhoPimpleFoam/createFields.H
View file @
cec56de0
...
...
@@ -39,8 +39,27 @@
#include "compressibleCreatePhi.H"
dimensionedScalar
rhoMax
(
pimple
.
dict
().
lookup
(
"rhoMax"
));
dimensionedScalar
rhoMin
(
pimple
.
dict
().
lookup
(
"rhoMin"
));
dimensionedScalar
rhoMax
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
pimple
.
dict
(),
GREAT
,
dimDensity
)
);
dimensionedScalar
rhoMin
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
pimple
.
dict
(),
0
,
dimDensity
)
);
Info
<<
"Creating turbulence model
\n
"
<<
endl
;
autoPtr
<
compressible
::
turbulenceModel
>
turbulence
...
...
applications/solvers/compressible/rhoSimpleFoam/createFields.H
View file @
cec56de0
...
...
@@ -44,8 +44,27 @@
scalar
pRefValue
=
0
.
0
;
setRefCell
(
p
,
simple
.
dict
(),
pRefCell
,
pRefValue
);
dimensionedScalar
rhoMax
(
simple
.
dict
().
lookup
(
"rhoMax"
));
dimensionedScalar
rhoMin
(
simple
.
dict
().
lookup
(
"rhoMin"
));
dimensionedScalar
rhoMax
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
simple
.
dict
(),
GREAT
,
dimDensity
)
);
dimensionedScalar
rhoMin
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
simple
.
dict
(),
0
,
dimDensity
)
);
Info
<<
"Creating turbulence model
\n
"
<<
endl
;
autoPtr
<
compressible
::
RASModel
>
turbulence
...
...
applications/solvers/compressible/rhoSimpleFoam/rhoPorousSimpleFoam/createFields.H
View file @
cec56de0
...
...
@@ -43,8 +43,27 @@
scalar
pRefValue
=
0
.
0
;
setRefCell
(
p
,
simple
.
dict
(),
pRefCell
,
pRefValue
);
dimensionedScalar
rhoMax
(
simple
.
dict
().
lookup
(
"rhoMax"
));
dimensionedScalar
rhoMin
(
simple
.
dict
().
lookup
(
"rhoMin"
));
dimensionedScalar
rhoMax
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
simple
.
dict
(),
GREAT
,
dimDensity
)
);
dimensionedScalar
rhoMin
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
simple
.
dict
(),
0
,
dimDensity
)
);
Info
<<
"Creating turbulence model
\n
"
<<
endl
;
autoPtr
<
compressible
::
RASModel
>
turbulence
...
...
applications/solvers/lagrangian/reactingParcelFoam/createFields.H
View file @
cec56de0
...
...
@@ -58,18 +58,24 @@
dimensionedScalar
rhoMax
(
"rhoMax"
,
dimDensity
,
pimple
.
dict
()
.
lookupOrDefault
<
scalar
>
(
"rhoMax"
,
GREAT
)
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
pimple
.
dict
(),
GREAT
,
dimDensity
)
);
dimensionedScalar
rhoMin
(
"rhoMin"
,
dimDensity
,
pimple
.
dict
()
.
lookupOrDefault
<
scalar
>
(
"rhoMin"
,
-
GREAT
)
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
pimple
.
dict
(),
0
,
dimDensity
)
);
Info
<<
"Creating turbulence model
\n
"
<<
endl
;
...
...
applications/solvers/lagrangian/reactingParcelFoam/simpleReactingParcelFoam/createFields.H
View file @
cec56de0
...
...
@@ -56,8 +56,28 @@
#include "compressibleCreatePhi.H"
dimensionedScalar
rhoMax
(
simple
.
dict
().
lookup
(
"rhoMax"
));
dimensionedScalar
rhoMin
(
simple
.
dict
().
lookup
(
"rhoMin"
));
dimensionedScalar
rhoMax
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
simple
.
dict
(),
GREAT
,
dimDensity
)
);
dimensionedScalar
rhoMin
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
simple
.
dict
(),
0
,
dimDensity
)
);
Info
<<
"Creating turbulence model
\n
"
<<
endl
;
autoPtr
<
compressible
::
turbulenceModel
>
turbulence
...
...
applications/solvers/lagrangian/sprayFoam/createFields.H
View file @
cec56de0
...
...
@@ -34,7 +34,7 @@
"rho"
,
runTime
.
timeName
(),
mesh
,
IOobject
::
NO_
READ
,
IOobject
::
READ
_IF_PRESENT
,
IOobject
::
AUTO_WRITE
),
thermo
.
rho
()
...
...
@@ -56,6 +56,28 @@
#include "compressibleCreatePhi.H"
dimensionedScalar
rhoMax
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMax"
,
pimple
.
dict
(),
GREAT
,
dimDensity
)
);
dimensionedScalar
rhoMin
(
dimensionedScalar
::
lookupOrDefault
(
"rhoMin"
,
pimple
.
dict
(),
0
,
dimDensity
)
);
Info
<<
"Creating turbulence model
\n
"
<<
endl
;
autoPtr
<
compressible
::
turbulenceModel
>
turbulence
(
...
...
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