FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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 Usage
31  - financialFoam [OPTION]
32 
33  @param -case <dir> \n
34  Specify the case directory
35 
36  @param -parallel \n
37  Run the case in parallel
38 
39  @param -help \n
40  Display short usage message
41 
42  @param -doc \n
43  Display Doxygen documentation page
44 
45  @param -srcDoc \n
46  Display source code
47 
48 \*---------------------------------------------------------------------------*/
49 
50 #include <finiteVolume/fvCFD.H>
52 #include <OpenFOAM/OSspecific.H>
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 int main(int argc, char *argv[])
57 {
58  #include <OpenFOAM/setRootCase.H>
59 
60  #include <OpenFOAM/createTime.H>
61  #include <OpenFOAM/createMesh.H>
62  #include "createFields.H"
63 
64  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66  Info<< nl << "Calculating value(price of comodities)" << endl;
67 
68  surfaceScalarField phi("phi", (sigmaSqr - r)*(Pf & mesh.Sf()));
69 
70  volScalarField DV("DV", 0.5*sigmaSqr*sqr(P.component(Foam::vector::X)));
71 
72  Info<< "Starting time loop\n" << endl;
73 
74  while (runTime.loop())
75  {
76  delta == fvc::grad(V)().component(Foam::vector::X);
77 
78  solve
79  (
80  fvm::ddt(V)
81  + fvm::div(phi, V)
82  - fvm::Sp(fvc::div(phi), V)
83  - fvm::laplacian(DV, V)
84  ==
85  - fvm::Sp(r, V)
86  );
87 
88  runTime.write();
89 
90  if (runTime.outputTime())
91  {
92  writeCellGraph(V, runTime.graphFormat());
93  writeCellGraph(delta, runTime.graphFormat());
94  }
95 
96  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
97  << " ClockTime = " << runTime.elapsedClockTime() << " s"
98  << nl << endl;
99  }
100 
101  Info<< "End\n" << endl;
102 
103  return 0;
104 }
105 
106 
107 // ************************ vim: set sw=4 sts=4 et: ************************ //