All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
QGDInterpolate.H
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  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Description
29  Creates interpolation instances templated for QGD solver
30 
31 SourceFiles
32  QGDInterpolate.H
33 
34 \*---------------------------------------------------------------------------*/
35 namespace Foam
36 {
37 
38 template<class T>
39 tmp<GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh> >
40 qgdInterpolate(const GeometricField<T, Foam::fvPatchField, Foam::volMesh > &psi)
41 {
42  if
43  (
44  psi.mesh().schemesDict().subDict("interpolationSchemes").found
45  ("interpolate("+psi.name()+")")
46  )
47  {
48  return fvc::interpolate(psi);
49  }
50  else if
51  (
52  psi.mesh().schemesDict().subDict("interpolationSchemes").found
53  ("default")
54  )
55  {
56  if
57  (
58  psi.mesh().schemesDict().subDict("interpolationSchemes").template get<word>
59  ("default") == "none"
60  )
61  {
62  return linearInterpolate(psi);
63  }
64  return fvc::interpolate(psi);
65  }
66  return linearInterpolate(psi);
67 }
68 
69 template<class T>
70 tmp<GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh> >
71 qgdInterpolate(const tmp<GeometricField<T, Foam::fvPatchField, Foam::volMesh >> &tpsi)
72 {
73  return qgdInterpolate(tpsi());
74 }
75 
76 template<class T>
77 tmp<GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh> >
78 qgdFlux
79 (
80  const GeometricField<scalar, Foam::fvsPatchField, Foam::surfaceMesh>& flux,
81  const GeometricField<T, Foam::fvPatchField, Foam::volMesh > &psi,
82  const GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh>& psif,
83  const word fluxName
84 )
85 {
86  if (psi.mesh().schemesDict().subDict("divSchemes").found(fluxName))
87  {
88  //if
89  //(
90  // psi.mesh().schemesDict().subDict("divSchemes").template get<word>(fluxName)
91  // ==
92  // "linear"
93  //)
94  //{
95  // return flux*psif;
96  //}
97  return fvc::flux
98  (
99  flux,
100  psi,
101  fluxName
102  );
103  }
104  return flux*psif;
105 }
107 template<class T>
108 tmp<GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh> >
109 qgdFlux
110 (
111  const GeometricField<scalar, Foam::fvsPatchField, Foam::surfaceMesh>& flux,
112  const GeometricField<T, Foam::fvPatchField, Foam::volMesh > &psi,
113  const GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh>& psif
114 )
115 {
116  word fluxName = "div("+flux.name()+","+psi.name()+")";
117  return qgdFlux(flux,psi,psif,fluxName);
118 }
119 
121 template<class T>
122 tmp<GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh> >
123 qgdFlux
124 (
125  const tmp<GeometricField<scalar, Foam::fvsPatchField, Foam::surfaceMesh>>& flux,
126  const tmp<GeometricField<T, Foam::fvPatchField, Foam::volMesh >> &psi,
127  const tmp<GeometricField<T, Foam::fvsPatchField, Foam::surfaceMesh>>& psif
128 )
129 {
130  return qgdFlux(flux(),psi(),psif());
131 }
132 
133 }
134 
135 //
136 //END-OF-FILE
137 //
138 
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &dir, const word &reconFieldName=word::null)
Interpolate field vf according to direction dir.
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdFlux(const GeometricField< scalar, Foam::fvsPatchField, Foam::surfaceMesh > &flux, const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi, const GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > &psif, const word fluxName)
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdInterpolate(const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi)
const volScalarField & psi
Definition: createFields.H:20