Test-fvc.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-2013 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  test
26 
27 Description
28  Finite volume method test code.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 int main(int argc, char *argv[])
37 {
38  #include "setRootCase.H"
39 
40  #include "createTime.H"
41  #include "createMesh.H"
42 
43  volScalarField fx(pow(mesh.C().component(vector::X), 1));
44  fx.write();
45  volScalarField gradx4(fvc::grad(fx)().component(vector::X));
46  gradx4.write();
47 
48  volVectorField curlC(fvc::curl(1.0*mesh.C()));
49  curlC.write();
50 
51  surfaceScalarField xf(mesh.Cf().component(vector::X));
52  surfaceScalarField xf4(pow(xf, 4));
53 
54  for (int i=1; i<xf4.size()-1; i++)
55  {
56  scalar gradx4a = (xf4[i] - xf4[i-1])/(xf[i] - xf[i-1]);
57  Info<< (gradx4a - gradx4[i])/gradx4a << endl;
58  }
59 
60  Info<< "end" << endl;
61 }
62 
63 
64 // ************************************************************************* //
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:380
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
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
Foam::fvc::curl
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcCurl.C:45
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
main
int main(int argc, char *argv[])
Definition: Test-fvc.C:33
Foam::Info
messageStream Info
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
setRootCase.H
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
createMesh.H
createTime.H
fvCFD.H
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52