financialFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  financialFoam
26 
27 Description
28  Solves the Black-Scholes equation to price commodities.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "OSspecific.H"
34 #include "setWriter.H"
35 #include "writeFile.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 int main(int argc, char *argv[])
40 {
41  #define NO_CONTROL
42  #include "postProcess.H"
43 
44  #include "setRootCaseLists.H"
45  #include "createTime.H"
46  #include "createMesh.H"
47  #include "createFields.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< nl << "Calculating value(price of commodities)" << endl;
52 
53  surfaceScalarField phi("phi", (sigmaSqr - r)*(Pf & mesh.Sf()));
54 
55  volScalarField DV("DV", 0.5*sigmaSqr*sqr(P.component(Foam::vector::X)));
56 
57  Info<< "Starting time loop\n" << endl;
58 
59  while (runTime.loop())
60  {
62 
63  solve
64  (
65  fvm::ddt(V)
66  + fvm::div(phi, V)
67  - fvm::Sp(fvc::div(phi), V)
68  - fvm::laplacian(DV, V)
69  ==
70  - fvm::Sp(r, V)
71  );
72 
73  runTime.write();
74 
75  if (runTime.writeTime())
76  {
77  setWriter::New(runTime.graphFormat())->write
78  (
79  runTime.globalPath()
80  /functionObjects::writeFile::outputPrefix
81  /args.executable()
82  /runTime.timeName(),
83 
84  args.executable(),
85 
86  coordSet(true, word::null, mesh.C().primitiveField(), "x"),
87 
88  "V", V.primitiveField(),
89  "delta", delta.primitiveField()
90  );
91  }
92 
93  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
94  << " ClockTime = " << runTime.elapsedClockTime() << " s"
95  << nl << endl;
96  }
97 
98  Info<< "End\n" << endl;
99 
100  return 0;
101 }
102 
103 
104 // ************************************************************************* //
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
phi
Definition: correctPhi.H:3
static const char nl
Definition: Ostream.H:260
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
rhoEqn solve()
messageStream Info
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
Foam::argList args(argc, argv)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)