createFields.H
Go to the documentation of this file.
1 Info<< "Reading financial properties\n" << endl;
2 
3 IOdictionary financialProperties
4 (
5  IOobject
6  (
7  "financialProperties",
8  runTime.constant(),
9  mesh,
10  IOobject::MUST_READ_IF_MODIFIED,
11  IOobject::NO_WRITE
12  )
13 );
14 
15 dimensionedScalar strike
16 (
17  "strike",
18  dimLength,
19  financialProperties
20 );
21 
23 (
24  "r",
26  financialProperties
27 );
28 
30 (
31  "sigma",
32  dimensionSet(0, 0, -0.5, 0, 0),
33  financialProperties
34 );
35 
36 dimensionedScalar sigmaSqr = sqr(sigma);
37 
38 
39 Info<< nl << "Reading field V" << endl;
40 
42 (
43  IOobject
44  (
45  "V",
46  runTime.timeName(),
47  mesh,
48  IOobject::MUST_READ,
49  IOobject::AUTO_WRITE
50  ),
51  mesh
52 );
53 
54 
56 (
57  IOobject
58  (
59  "Pf",
60  runTime.timeName(),
61  mesh,
62  IOobject::NO_READ,
63  IOobject::NO_WRITE
64  ),
65  mesh.Cf()
66 );
67 
68 
70 (
71  IOobject
72  (
73  "P",
74  runTime.timeName(),
75  mesh,
76  IOobject::NO_READ,
77  IOobject::NO_WRITE
78  ),
79  mesh.C()
80 );
81 
82 V == max
83 (
84  P.component(Foam::vector::X) - strike,
85  dimensionedScalar(V.dimensions(), Zero)
86 );
87 
89 (
90  IOobject
91  (
92  "delta",
93  runTime.timeName(),
94  mesh,
95  IOobject::NO_READ,
96  IOobject::AUTO_WRITE
97  ),
98  fvc::grad(V)().component(Foam::vector::X)
99 );
const vector delta(L.x()/scalar(N.x()), L.y()/scalar(N.y()), L.z()/scalar(N.z()))
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127