This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/financial/financialFoam/financialFoam.C
2016-06-21 15:04:12 +02:00

89 lines
2.7 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
foam-extend is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
General Public License for more details.
You should have received a copy of the GNU General Public License
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
financialFoam
Description
Solves the Black-Scholes equation to price commodities.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "writeCellGraph.H"
#include "OSspecific.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createTime.H"
# include "createMesh.H"
# include "createFields.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< nl << "Calculating value(price of comodities)" << endl;
surfaceScalarField phi("phi", (sigmaSqr - r)*(Pf & mesh.Sf()));
volScalarField DV("DV", 0.5*sigmaSqr*sqr(P.component(Foam::vector::X)));
Info<< "Starting time loop\n" << endl;
while (runTime.loop())
{
delta == fvc::grad(V)().component(Foam::vector::X);
solve
(
fvm::ddt(V)
+ fvm::div(phi, V)
- fvm::Sp(fvc::div(phi), V)
- fvm::laplacian(DV, V)
==
- fvm::Sp(r, V)
);
runTime.write();
if (runTime.outputTime())
{
writeCellGraph(V, runTime.graphFormat());
writeCellGraph(delta, runTime.graphFormat());
}
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //