Skip to content
Snippets Groups Projects
Commit 700e17bd authored by Andrew Heather's avatar Andrew Heather
Browse files

INT: Initial check-in of IH Cantabria streamFunction wave generation

model and test case

Code supplied by Gabriel BARAJAS OJEDA
parent ac20f63c
Branches
Tags
No related merge requests found
Showing
with 1214 additions and 0 deletions
......@@ -8,6 +8,7 @@ waveGenerationModels/derived/Boussinesq/BoussinesqWaveModel.C
waveGenerationModels/derived/cnoidal/cnoidalWaveModel.C
waveGenerationModels/derived/Grimshaw/GrimshawWaveModel.C
waveGenerationModels/derived/McCowan/McCowanWaveModel.C
waveGenerationModels/derived/streamFunction/streamFunctionWaveModel.C
waveGenerationModels/derived/StokesII/StokesIIWaveModel.C
waveGenerationModels/derived/StokesI/StokesIWaveModel.C
waveGenerationModels/derived/StokesV/StokesVWaveModel.C
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016-2017 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2015 IH-Cantabria
-------------------------------------------------------------------------------
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 "streamFunctionWaveModel.H"
#include "mathematicalConstants.H"
#include "addToRunTimeSelectionTable.H"
using namespace Foam::constant;
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace waveModels
{
defineTypeNameAndDebug(streamFunction, 0);
addToRunTimeSelectionTable
(
waveModel,
streamFunction,
patch
);
}
}
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
Foam::scalar Foam::waveModels::streamFunction::eta
(
const scalar h,
const scalar kx,
const scalar ky,
const scalar T,
const scalar x,
const scalar y,
const scalar omega,
const scalar t,
const scalar phase
) const
{
const scalar k = sqrt(kx*kx + ky*ky);
scalar strfnAux = 0.0;
forAll(Ejs_, iterSF)
{
strfnAux += Ejs_[iterSF]*cos((iterSF + 1)
*(kx*x + ky*y - omega*t + phase));
}
return (1/k)*strfnAux;
}
Foam::vector Foam::waveModels::streamFunction::Uf
(
const scalar h,
const scalar kx,
const scalar ky,
const scalar T,
const scalar x,
const scalar y,
const scalar omega,
const scalar t,
const scalar phase,
const scalar z
) const
{
const scalar k = sqrt(kx*kx + ky*ky);
const scalar phaseTot = kx*x + ky*y - omega*t + phase;
scalar u = 0.0;
scalar w = 0.0;
forAll(Bjs_, iterSF2)
{
u += (iterSF2 + 1)*Bjs_[iterSF2] *cosh((iterSF2 + 1)*k*z)
/cosh((iterSF2 + 1)*k*h)*cos((iterSF2 + 1)*phaseTot);
w += (iterSF2 + 1)*Bjs_[iterSF2]*sinh((iterSF2 + 1)*k*z)
/cosh((iterSF2 + 1)*k*h)*sin((iterSF2 + 1)*phaseTot);
}
u = waveLengthSF_/T - uMean_ + sqrt(mag(g_)/k)*u;
w = sqrt(mag(g_)/k)*w;
scalar v = u*sin(waveAngle_);
u *= cos(waveAngle_);
return vector(u, v, w);
}
// * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
void Foam::waveModels::streamFunction::setLevel
(
const scalar t,
const scalar tCoeff,
scalarField& level
) const
{
const scalar waveOmega = mathematical::twoPi/wavePeriod_;
const scalar waveK = mathematical::twoPi/waveLengthSF_;
const scalar waveKx = waveK*cos(waveAngle_);
const scalar waveKy = waveK*sin(waveAngle_);
forAll(level, paddlei)
{
const scalar eta =
this->eta
(
waterDepthRef_,
waveKx,
waveKy,
wavePeriod_,
xPaddle_[paddlei],
yPaddle_[paddlei],
waveOmega,
t,
wavePhase_
);
level[paddlei] = waterDepthRef_ + tCoeff*eta;
}
}
void Foam::waveModels::streamFunction::setVelocity
(
const scalar t,
const scalar tCoeff,
const scalarField& level
)
{
const scalar waveOmega = mathematical::twoPi/wavePeriod_;
const scalar waveK = mathematical::twoPi/waveLengthSF_;
const scalar waveKx = waveK*cos(waveAngle_);
const scalar waveKy = waveK*sin(waveAngle_);
forAll(U_, facei)
{
// Fraction of geometry represented by paddle - to be set
scalar fraction = 1;
// Height - to be set
scalar z = 0;
setPaddlePropeties(level, facei, fraction, z);
if (fraction > 0)
{
const label paddlei = faceToPaddle_[facei];
const vector Uf = this->Uf
(
waterDepthRef_,
waveKx,
waveKy,
wavePeriod_,
xPaddle_[paddlei],
yPaddle_[paddlei],
waveOmega,
t,
wavePhase_,
z
);
U_[facei] = fraction*Uf*tCoeff;
}
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::waveModels::streamFunction::streamFunction
(
const dictionary& dict,
const fvMesh& mesh,
const polyPatch& patch,
const bool readFields
)
:
regularWaveModel(dict, mesh, patch, false),
uMean_(0),
waveLengthSF_(0),
Bjs_( List<scalar> (1, -1.0) ),
Ejs_( List<scalar> (1, -1.0) )
{
if (readFields)
{
readDict(dict);
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::waveModels::streamFunction::~streamFunction()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::waveModels::streamFunction::readDict(const dictionary& overrideDict)
{
if (regularWaveModel::readDict(overrideDict))
{
lookup("uMean") >> uMean_;
lookup("waveLengthSF") >> waveLengthSF_;
lookup("Bjs") >> Bjs_;
lookup("Ejs") >> Ejs_;
return true;
}
return false;
}
void Foam::waveModels::streamFunction::info(Ostream& os) const
{
regularWaveModel::info(os);
os << " uMean : " << uMean_ << nl
<< " wave Length streamFunction : " << waveLengthSF_ << nl
<< " Bj coefficients : " << Bjs_ << nl
<< " Ej coefficients : " << Ejs_ << nl;
}
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2016-2017 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2015 IH-Cantabria
-------------------------------------------------------------------------------
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::waveModels::streamFunction
Description
streamFunction wave model
\*---------------------------------------------------------------------------*/
#ifndef waveModels_streamFunction_H
#define waveModels_streamFunction_H
#include "regularWaveModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace waveModels
{
/*---------------------------------------------------------------------------*\
Class streamFunction Declaration
\*---------------------------------------------------------------------------*/
class streamFunction
:
public regularWaveModel
{
private:
// Private Member Functions
//- Wave height
virtual scalar eta
(
const scalar h,
const scalar kx,
const scalar ky,
const scalar T,
const scalar x,
const scalar y,
const scalar omega,
const scalar t,
const scalar phase
) const;
//- Wave velocity
virtual vector Uf
(
const scalar d,
const scalar kx,
const scalar ky,
const scalar T,
const scalar x,
const scalar y,
const scalar omega,
const scalar t,
const scalar phase,
const scalar z
) const;
protected:
// Protected data
//- Mean fluid speed in frame of reference (stream function)
scalar uMean_;
//-Stream waveLength
scalar waveLengthSF_;
//- Stream Function Bj coefficients
scalarList Bjs_;
//- Stream Function Ej coefficients
scalarList Ejs_;
// Protected Member Functions
//- Set the water level
virtual void setLevel
(
const scalar t,
const scalar tCoeff,
scalarField& level
) const;
//- Calculate the wave model velocity
virtual void setVelocity
(
const scalar t,
const scalar tCoeff,
const scalarField& level
);
public:
//- Runtime type information
TypeName("streamFunction");
//- Constructor
streamFunction
(
const dictionary& dict,
const fvMesh& mesh,
const polyPatch& patch,
const bool readFields = true
);
//- Destructor
virtual ~streamFunction();
// Public Member Functions
//- Read from dictionary
virtual bool readDict(const dictionary& overrideDict);
//- Info
virtual void info(Ostream& os) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace waveModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volVectorField;
location "0";
object U;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 1 -1 0 0 0 0];
internalField uniform (0 0 0);
boundaryField
{
inlet
{
type waveVelocity;
value uniform (0 0 0);
}
outlet
{
type waveVelocity;
value uniform (0 0 0);
}
sides
{
type empty;
}
ground
{
type fixedValue;
value uniform (0 0 0);
}
top
{
type pressureInletOutletVelocity;
value uniform (0 0 0);
}
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volScalarField;
object alpha.water;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 0 0 0 0 0 0];
internalField uniform 0;
boundaryField
{
inlet
{
type waveAlpha;
value uniform 0;
}
outlet
{
type zeroGradient;
}
ground
{
type zeroGradient;
}
sides
{
type empty;
}
top
{
type inletOutlet;
inletValue uniform 0;
value uniform 0;
}
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class volScalarField;
object p_rgh;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [1 -1 -2 0 0 0 0];
internalField uniform 0;
boundaryField
{
inlet
{
type fixedFluxPressure;
value uniform 0;
}
outlet
{
type fixedFluxPressure;
value uniform 0;
}
ground
{
type fixedFluxPressure;
value uniform 0;
}
sides
{
type empty;
}
top
{
type totalPressure;
p0 uniform 0;
}
}
// ************************************************************************* //
#!/bin/sh
cd ${0%/*} || exit 1 # run from this directory
# Source tutorial run functions
. $WM_PROJECT_DIR/bin/tools/CleanFunctions
rm -rf 0
cleanCase
#!/bin/sh
cd ${0%/*} || exit 1 # run from this directory
# Source tutorial run functions
. $WM_PROJECT_DIR/bin/tools/RunFunctions
restore0Dir
runApplication blockMesh
runApplication decomposePar
runParallel setFields
runParallel $(getApplication)
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class uniformDimensionedVectorField;
location "constant";
object g;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
dimensions [0 1 -2 0 0 0 0];
value ( 0 0 -9.81 );
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "constant";
object transportProperties;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
phases (water air);
water
{
transportModel Newtonian;
nu 1e-06;
rho 1000;
}
air
{
transportModel Newtonian;
nu 1.48e-05;
rho 1;
}
sigma 0.07;
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "constant";
object turbulenceProperties;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
simulationType laminar;
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "constant";
object wavesProperties;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
inlet
{
alpha alpha.water;
waveModel streamFunction;
nPaddle 1;
waveHeight 0.1517;
waveAngle 0.0;
rampTime 6.034;
activeAbsorption yes;
wavePeriod 3.017;
uMean 2.0825;
waveLengthSF 6.2832;
Bjs
10
(
8.6669014e-002
2.4849799e-002
7.7446850e-003
2.3355420e-003
6.4497731e-004
1.5205114e-004
2.5433769e-005
-2.2045436e-007
-2.8711504e-006
-1.2287334e-006
);
Ejs
10
(
5.6009609e-002
3.1638171e-002
1.5375952e-002
7.1743178e-003
3.3737077e-003
1.6324880e-003
8.2331980e-004
4.4403497e-004
2.7580059e-004
2.2810557e-004
);
}
outlet
{
alpha alpha.water;
waveModel shallowWaterAbsorption;
nPaddle 1;
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
object blockMeshDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
convertToMeters 1;
vertices
(
( 0.0 0.0 0.0)
( 40.0 0.0 0.0)
( 40.0 0.01 0.0)
( 0.0 0.01 0.0)
( 0.0 0.0 0.8)
( 40.0 0.0 0.8)
( 40.0 0.01 0.8)
( 0.0 0.01 0.8)
);
blocks
(
hex (0 1 2 3 4 5 6 7) (2000 1 80) simpleGrading (1 1 1)
);
edges
(
);
boundary
(
inlet
{
type patch;
faces
(
(0 4 7 3)
);
}
outlet
{
type patch;
faces
(
(1 5 6 2)
);
}
ground
{
type wall;
faces
(
(0 1 2 3)
);
}
top
{
type patch;
faces
(
(4 5 6 7)
);
}
sides
{
type empty;
faces
(
(0 1 5 4)
(3 2 6 7)
);
}
);
mergePatchPairs
(
);
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object controlDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
application interFoam;
startFrom latestTime;
startTime 0;
stopAt endTime;
endTime 30.0;
deltaT 0.01;
writeControl adjustableRunTime;
writeInterval 0.033;
purgeWrite 0;
writeFormat ascii;
writePrecision 6;
writeCompression off;
timeFormat general;
timePrecision 6;
runTimeModifiable yes;
adjustTimeStep yes;
maxCo 0.65;
maxAlphaCo 0.65;
maxDeltaT 0.05;
functions
{
line
{
type sets;
libs ("libsampling.so");
enabled true;
writeControl writeTime;
writeInterval 1;
interpolationScheme cellPoint;
setFormat raw;
sets
(
line1
{
type uniform;
axis distance;
start ( 1.0 0.005 0.0 );
end ( 1.0 0.005 1.0 );
nPoints 1001;
}
line2
{
type uniform;
axis distance;
start ( 2.0 0.005 0.0 );
end ( 2.0 0.005 1.0 );
nPoints 1001;
}
line3
{
type uniform;
axis distance;
start ( 3.0 0.005 0.0 );
end ( 3.0 0.005 1.0 );
nPoints 1001;
}
line4
{
type uniform;
axis distance;
start ( 5.0 0.005 0.0 );
end ( 5.0 0.005 1.0 );
nPoints 1001;
}
line5
{
type uniform;
axis distance;
start ( 7.5 0.005 0.0 );
end ( 7.5 0.005 1.0 );
nPoints 1001;
}
line6
{
type uniform;
axis distance;
start ( 10.0 0.005 0.0 );
end ( 10.0 0.005 1.0 );
nPoints 1001;
}
line7
{
type uniform;
axis distance;
start ( 12.0 0.005 0.0 );
end ( 12.0 0.005 1.0 );
nPoints 1001;
}
line8
{
type uniform;
axis distance;
start ( 15.0 0.005 0.0 );
end ( 15.0 0.005 1.0 );
nPoints 1001;
}
line9
{
type uniform;
axis distance;
start ( 20.0 0.005 0.0 );
end ( 20.0 0.005 1.0 );
nPoints 1001;
}
line10
{
type uniform;
axis distance;
start ( 28.0 0.005 0.0 );
end ( 28.0 0.005 1.0 );
nPoints 1001;
}
);
fixedLocations false;
fields
(
U alpha.water
);
}
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
object decomposeParDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
numberOfSubdomains 2;
method hierarchical;
hierarchicalCoeffs
{
n (2 1 1);
delta 0.001;
order xyz;
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object fvSchemes;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
ddtSchemes
{
default Euler;
}
gradSchemes
{
default Gauss linear;
}
divSchemes
{
div(rhoPhi,U) Gauss linearUpwind grad(U);
div(phi,alpha) Gauss vanLeer;
div(phirb,alpha) Gauss linear;
div(((rho*nuEff)*dev2(T(grad(U))))) Gauss linear;
}
laplacianSchemes
{
default Gauss linear orthogonal;
}
interpolationSchemes
{
default linear;
}
snGradSchemes
{
default orthogonal;
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object fvSolution;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
solvers
{
"alpha.water.*"
{
nAlphaCorr 1;
nAlphaSubCycles 3;
cAlpha 1;
}
"pcorr.*"
{
solver PCG;
preconditioner DIC;
tolerance 1e-6;
relTol 0;
}
p_rgh
{
solver PCG;
preconditioner DIC;
tolerance 1e-6;
relTol 0.1;
}
p_rghFinal
{
solver GAMG;
smoother DIC;
tolerance 1e-7;
relTol 0;
}
U
{
solver PBiCG;
preconditioner DILU;
tolerance 1e-6;
relTol 0.1;
}
UFinal
{
solver PBiCG;
preconditioner DILU;
tolerance 1e-6;
relTol 0;
}
}
PIMPLE
{
momentumPredictor no;
nCorrectors 2;
nNonOrthogonalCorrectors 0;
}
// ************************************************************************* //
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: plus |
| \\ / A nd | Web: www.OpenFOAM.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object setFieldsDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
defaultFieldValues
(
volScalarFieldValue alpha.water 0
);
regions
(
boxToCell
{
box (0 0 0) (40.0 1.0 0.4046);
fieldValues
(
volScalarFieldValue alpha.water 1
);
}
);
// ************************************************************************* //
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment