fvcMagSqrGradGrad.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 -------------------------------------------------------------------------------
10 License
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 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fvcMagSqrGradGrad.H"
29 #include "fvcGrad.H"
30 #include "fvMesh.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fvc
40 {
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 template<class Type>
45 tmp<volScalarField> magSqrGradGrad
46 (
48 )
49 {
50  tmp<volScalarField> tMagSqrGradGrad
51  (
53  );
54 
55  // Loop over other vector field components
56  for (direction cmpt = 1; cmpt < pTraits<Type>::nComponents; cmpt++)
57  {
58  tMagSqrGradGrad.ref() +=
59  magSqr(fvc::grad(fvc::grad(vf.component(cmpt))))();
60  }
61 
62  return tMagSqrGradGrad;
63 }
64 
65 
66 template<class Type>
69 (
71 )
72 {
73  tmp<volScalarField> tMagSqrGradGrad(fvc::magSqrGradGrad(tvf()));
74  tvf.clear();
75  return tMagSqrGradGrad;
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 } // End namespace fvc
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 } // End namespace Foam
86 
87 // ************************************************************************* //
Foam::GeometricField::component
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
fvcGrad.H
Calculate the gradient of the given field.
Foam::fvc::magSqrGradGrad
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcMagSqrGradGrad.C:46
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::GeometricField< Type, fvPatchField, volMesh >
fvcMagSqrGradGrad.H
Calculate the magnitiude of the square of the gradient of the gradient of the given volField.