All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
extendedFaceStencilScalarDer.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  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  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
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  You should have received a copy of the GNU General Public License
23  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 Group grpleastSquaresOpt
25  This group contains common part of QGD solvers.
26 Class
27  Foam::fvsc::leastSquaresOpt::extendedFaceStencilScalarDer
28 Description
29  Methods for optimal calculating of directional derivative.
30  With parallel realisation.
31 \*---------------------------------------------------------------------------*/
32 
33 #include "leastSquaresStencilOpt.H"
34 #include "polyMesh.H"
35 #include "fvMesh.H"
36 #include "word.H"
37 #include "IOstream.H"
38 #include "Ostream.H"
39 #include <HashTable.H>
40 
41 void Foam::fvsc::leastSquaresOpt::faceScalarDer(const Field<scalar>& iF,const Field<scalar>& sF,int com, surfaceScalarField& rField)
42 {
43  forAll(sF, facei)
44  {
45  rField[facei] = 0.0;
46  forAll(GdfAll_[facei], i)
47  {
48  rField[facei] += wf2All_[facei][i]*GdfAll_[facei][i].component(com)*(iF[neighbourCells_[facei][i]] - sF[facei]);
49  }
50  }
51 
52 };
53 
54 void Foam::fvsc::leastSquaresOpt::faceScalarDer(const tmp<Field<scalar>>& tiF,const tmp<Field<scalar>>& tsF, int com, tmp<surfaceScalarField>& trField )
55 {
56 };
57 
58 
59 void Foam::fvsc::leastSquaresOpt::faceScalarDer(const List3<scalar>& procVfValues,const surfaceScalarField& sF,int derComp, surfaceScalarField& rField)
60 {
61  scalar gf = 0.0;
62  forAll(procPairs_, patchI)
63  {
64  label procPatchId = procPairs_[patchI];
65  if (procPatchId > -1)
66  {
67  fvsPatchScalarField& pgradf = rField.boundaryFieldRef()[procPatchId];
68  const List2<scalar> & pvf = procVfValues[patchI];
69  const List2<scalar> & pwf2= procWf2_[patchI];
70  const List2<vector> & pgdf= procGdf_[patchI];
71 
72  forAll(pgradf, iFace)
73  {
74  gf = 0.0;
75  forAll(procGdf_[patchI][iFace], i)
76  {
77  gf += pwf2[iFace][i]*pgdf[iFace][i].component(derComp)*
78  (pvf[iFace][i] - sF.boundaryField()[procPatchId][iFace]);
79  }
80  rField[iFace] = gf;
81  }
82  }
83  }
84 };
85 
86 void Foam::fvsc::leastSquaresOpt::faceScalarDer(const tmp<List3<scalar>>& tprocVfValues,const tmp<surfaceScalarField>&, int derComp, tmp<surfaceScalarField>& trField )
87 {
88 };
89 
void faceScalarDer(const Field< scalar > &iF, const Field< scalar > &sF, int com, surfaceScalarField &rField)
forAll(Y, i)
Definition: QGDYEqn.H:36