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 | Copyright (C) 2011 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 Group
28  grpFinancialSolvers
29 
30 Description
31  Solves the Black-Scholes equation to price commodities.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "writeCellGraph.H"
37 #include "OSspecific.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 int main(int argc, char *argv[])
42 {
43  #include "setRootCase.H"
44 
45  #include "createTime.H"
46  #include "createMesh.H"
47  #include "createFields.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< nl << "Calculating value(price of comodities)" << 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.outputTime())
76  {
77  writeCellGraph(V, runTime.graphFormat());
78  writeCellGraph(delta, runTime.graphFormat());
79  }
80 
81  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
82  << " ClockTime = " << runTime.elapsedClockTime() << " s"
83  << nl << endl;
84  }
85 
86  Info<< "End\n" << endl;
87 
88  return 0;
89 }
90 
91 
92 // ************************************************************************* //
OSspecific.H
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::writeCellGraph
void writeCellGraph(const volScalarField &vsf, const word &graphFormat)
Definition: writeCellGraph.C:14
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::Vector< scalar >::X
@ X
Definition: Vector.H:89
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
setRootCase.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
createMesh.H
createTime.H
fvCFD.H
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
writeCellGraph.H