All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
extendedFaceStencilScalarGradOpt.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::extendedFaceStencilScalarGradOpt
28 Description
29  Methods for optimal calculating of directional derivative.
30  With parallel realisation.
31 \*---------------------------------------------------------------------------*/
32 
33 
34 #include "leastSquaresStencilOpt.H"
35 #include "polyMesh.H"
36 #include "fvMesh.H"
37 #include "word.H"
38 #include "IOstream.H"
39 #include "Ostream.H"
40 #include <HashTable.H>
41 
42 #include "emptyFvPatch.H"
43 #include "coupledFvPatch.H"
44 #include "wedgeFvPatch.H"
45 #include "symmetryFvPatch.H"
46 #include "symmetryPlaneFvPatch.H"
47 
48 //- Calculate gradient of volume scalar function on the faces
49 //
50 // \param iF Internal scalar field.
51 // Allowable values: constant reference to the volScalarField.
52 //
53 // \return Gradient of iF (vector field) which was computed on the faces of mesh.
54 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::leastSquaresOpt::Grad(const volScalarField& iF)
55 {
56  surfaceScalarField sF = linearInterpolate(iF);
57  return Grad(iF,sF);
58 };
59 
60 Foam::tmp<Foam::surfaceVectorField>
61 Foam::fvsc::leastSquaresOpt::Grad(const volScalarField& iF, const surfaceScalarField& sF)
62 {
63 
64  tmp<surfaceVectorField> tgradIF(0.0*nf_*fvc::snGrad(iF));
65  surfaceVectorField& gradIF = tgradIF.ref();
66  //scalarField tField = sF;
67  surfaceScalarField tField = sF*0;
68 
69  faceScalarDer(iF.primitiveField(),sF.primitiveField(),0,tField);
70  gradIF.primitiveFieldRef().replace(0, tField);
71  faceScalarDer(iF.primitiveField(),sF.primitiveField(),1,tField);
72  gradIF.primitiveFieldRef().replace(1, tField);
73  faceScalarDer(iF.primitiveField(),sF.primitiveField(),2,tField);
74  gradIF.primitiveFieldRef().replace(2, tField);
75 
76  //update boundary field
77  forAll(mesh_.boundaryMesh(), ipatch)
78  {
79  bool notConstrain = true;
80  const fvPatch& fvp = mesh_.boundary()[ipatch];
81  if
82  (
83  isA<emptyFvPatch>(fvp) ||
84  isA<wedgeFvPatch>(fvp) ||
85  isA<coupledFvPatch>(fvp) ||
86  isA<symmetryFvPatch>(fvp) ||
87  isA<symmetryPlaneFvPatch>(fvp)
88  )
89  {
90  notConstrain = false;
91  }
92 
93  if (notConstrain)
94  {
95  gradIF.boundaryFieldRef()[ipatch] =
96  nf_.boundaryField()[ipatch] *
97  iF.boundaryField()[ipatch].snGrad();
98  }
99  }
100 
101  if(!Pstream::parRun())
102  {
103  return tgradIF;
104  }
105 
106  /*
107  *
108  * Update processor patches for parallel case
109  *
110  */
111 
112  //allocate storage for near-patch field
113 
114  List3<scalar> procVfValues(nProcPatches_); //array of values from neighb. processors
115  formVfValues(iF,procVfValues);
116 
117  faceScalarDer(procVfValues,sF,0,tField);
118  gradIF.boundaryFieldRef().replace(0, tField.boundaryFieldRef());
119  faceScalarDer(procVfValues,sF,1,tField);
120  gradIF.boundaryFieldRef().replace(1, tField.boundaryFieldRef());
121  faceScalarDer(procVfValues,sF,2,tField);
122  gradIF.boundaryFieldRef().replace(2, tField.boundaryFieldRef());
123 
124  return tgradIF;
125 
126 
127 };
128 
129 //
130 //END-OF-FILE
131 //
132 
void formVfValues(const volScalarField &iF, List3< scalar > &procVfValues)
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.
forAll(Y, i)
Definition: QGDYEqn.H:36