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