financialFoam.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 financialFoam
28
29Group
30 grpFinancialSolvers
31
32Description
33 Solves the Black-Scholes equation to price commodities.
34
35\*---------------------------------------------------------------------------*/
36
37#include "fvCFD.H"
38#include "writeCellGraph.H"
39#include "OSspecific.H"
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43int main(int argc, char *argv[])
44{
45 argList::addNote
46 (
47 "Solves the Black-Scholes equation to price commodities."
48 );
49
50 #define NO_CONTROL
51 #include "postProcess.H"
52
53 #include "addCheckCaseOptions.H"
54 #include "setRootCaseLists.H"
55 #include "createTime.H"
56 #include "createMesh.H"
57 #include "createFields.H"
58
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61 Info<< nl << "Calculating value(price of comodities)" << endl;
62
63 surfaceScalarField phi("phi", (sigmaSqr - r)*(Pf & mesh.Sf()));
64
65 volScalarField DV("DV", 0.5*sigmaSqr*sqr(P.component(Foam::vector::X)));
66
67 Info<< "Starting time loop\n" << endl;
68
69 while (runTime.loop())
70 {
71 delta == fvc::grad(V)().component(Foam::vector::X);
72
73 solve
74 (
75 fvm::ddt(V)
76 + fvm::div(phi, V)
77 - fvm::Sp(fvc::div(phi), V)
78 - fvm::laplacian(DV, V)
79 ==
80 - fvm::Sp(r, V)
81 );
82
83 runTime.write();
84
85 if (runTime.writeTime())
86 {
87 writeCellGraph(V, runTime.graphFormat());
88 writeCellGraph(delta, runTime.graphFormat());
89 }
90
91 runTime.printExecutionTime(Info);
92 }
93
94 Info<< "End\n" << endl;
95
96 return 0;
97}
98
99
100// ************************************************************************* //
scalar delta
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Required Classes.
surfaceScalarField & phi
dynamicFvMesh & mesh
engineTime & runTime
Required Variables.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
messageStream Info
Information stream (stdout output on master, null elsewhere)
void writeCellGraph(const volScalarField &vsf, const word &graphFormat)
Definition: writeCellGraph.C:9
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Execute application functionObjects to post-process existing results.
CEqn solve()