All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
extendedFaceStencilScalarGrad.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::leastSquares::extendedFaceStencilScalarGrad
26 \*---------------------------------------------------------------------------*/
27 
28 #include "leastSquaresStencil.H"
29 #include "polyMesh.H"
30 #include "fvMesh.H"
31 #include "word.H"
32 #include "IOstream.H"
33 #include "Ostream.H"
34 #include <HashTable.H>
35 
36 #include "emptyFvPatch.H"
37 #include "coupledFvPatch.H"
38 #include "wedgeFvPatch.H"
39 #include "symmetryFvPatch.H"
40 #include "symmetryPlaneFvPatch.H"
41 
42 //- Calculate gradient of volume scalar function on the faces
43 //
44 // \param iF Internal scalar field.
45 // Allowable values: constant reference to the volScalarField.
46 //
47 // \return Gradient of iF (vector field) which was computed on the faces of mesh.
48 Foam::tmp<Foam::surfaceVectorField> Foam::fvsc::leastSquares::Grad(const volScalarField& iF)
49 {
50  surfaceScalarField sF = linearInterpolate(iF);
51  surfaceScalarField sngF (fvc::snGrad(iF));
52 
53  tmp<surfaceVectorField> tgradIF(0.0 * nf_ * sngF);
54  surfaceVectorField& gradIF = tgradIF.ref();
55 
56  // List of faces
57  const faceList& faces = mesh_.faces();
58 
59  vector gf = vector::zero;
60  forAll(faces, facei)
61  {
62  if (mesh_.isInternalFace(facei))
63  {
64  gf = vector::zero;
65  forAll(GdfAll_[facei], i)
66  {
67  gf = gf + wf2All_[facei][i]*GdfAll_[facei][i]*(iF[neighbourCells_[facei][i]] - sF[facei]);
68  }
69 
70  gradIF[facei] = gf;
71  }
72  }
73 
74  //process faces with degenerate stencil
75  label dFaceId = -1;
77  {
78  dFaceId = internalDegFaces_[facei];
79  //gradIF[dFaceId] = nf_[dFaceId] * sngF[dFaceId];
80  gradIF[dFaceId] = sngF[dFaceId] * nf_[dFaceId];
81  }
82 
83  //update boundary field
84  forAll(mesh_.boundaryMesh(), ipatch)
85  {
86  bool notConstrain = true;
87  const fvPatch& fvp = mesh_.boundary()[ipatch];
88  if
89  (
90  isA<emptyFvPatch>(fvp) ||
91  isA<wedgeFvPatch>(fvp) ||
92  isA<coupledFvPatch>(fvp) ||
93  isA<symmetryFvPatch>(fvp) ||
94  isA<symmetryPlaneFvPatch>(fvp)
95 // fvp.coupled()
96  )
97  {
98  notConstrain = false;
99  }
100 
101  if (notConstrain)
102  {
103  gradIF.boundaryFieldRef()[ipatch] = nf_.boundaryField()[ipatch] *
104  iF.boundaryField()[ipatch].snGrad();
105 
106  }
107  }
108 
109  if(!Pstream::parRun())
110  {
111  return tgradIF;
112  }
113 
114  /*
115  *
116  * Update processor patches for parallel case
117  *
118  */
119  //allocate storage for near-patch field
120  List3<scalar> procVfValues(nProcPatches_); //array of values from neighb. processors
121  //set values from this domain
122  label cellId = -1;
123  forAll(procPairs_, patchI)
124  {
125  if (procPairs_[patchI] > -1)
126  {
127  procVfValues[patchI].resize(procWf2_[patchI].size());
128  forAll(procVfValues[patchI], faceI)
129  {
130  procVfValues[patchI][faceI].resize(procWf2_[patchI][faceI].size());
131  procVfValues[patchI][faceI] = 0.0; //make values zero
132  forAll(myProcPatchCells_[patchI][faceI], cellI)
133  {
134  cellId = myProcPatchCells_[patchI][faceI][cellI];
135  procVfValues[patchI][faceI][cellI] = iF.primitiveField()[cellId];
136  }
137  }
138  }
139  }
140 
141  //Step 1. Send field data to neighbouring processors (non-blocking mode)
142 
143  PstreamBuffers pBuffers(Pstream::commsTypes::nonBlocking);
144  forAll(procPairs_, procI)
145  {
146  label procId = neigProcs_[procI];
147 // label dataSz = 0;
148 
149  DynamicList<scalar> locVf;
150 
151  if (procPairs_[procI] > -1) //patch proc pair
152  {
153  forAll(procVfValues[procI], faceI)
154  {
155  for(
156  label
157  cellI = 0;
158  cellI <= ownEnd_[procI][faceI];
159  cellI++
160  )
161  {
162  locVf.append(procVfValues[procI][faceI][cellI]);
163  }
164 
165  }
166  }
167  else //corner connected process
168  {
169  label cellId = -1;
170  label addrId = corProcIds_[procId];
171  //label addrId = corProcIds_.capacity();
172  forAll(corCellIds_[addrId], iCellId)
173  {
174  cellId = corCellIds_[addrId][iCellId];
175  locVf.append(iF.primitiveField()[cellId]);
176  }
177  }
178 
179  UOPstream oProcStr(procId, pBuffers);
180  oProcStr << locVf;
181  }
182 
183  //Step 2. Recieve field data from neighbouring processors
184  pBuffers.finishedSends();
185  label iCorProc = 0;
186  forAll(procPairs_, procI)
187  {
188  label procId = neigProcs_[procI];
189 
190  UIPstream iProcStr(procId, pBuffers);
191  List<scalar> locVf (iProcStr);
192 
193  if (procPairs_[procI] > -1)
194  {
195  label iVf = 0;
196  forAll(neiStart_[procI], iFace)
197  {
198  for(
199  label
200  iCell=neiStart_[procI][iFace];
201  iCell<=neiEnd_[procI][iFace];
202  iCell++
203  )
204  {
205  procVfValues[procI][iFace][iCell] =
206  locVf[iVf];
207  iVf++;
208  }
209  }
210  }
211  else
212  {
213  label patchNo = -1;
214  label faceNo = -1;
215  label cellNo = -1;
216  label offset = -1;
217 
218  const List<Triple<label> >& addr = corAddr_[iCorProc];
219 
220  forAll(addr, iVal)
221  {
222  patchNo = addr[iVal][0];
223  faceNo = addr[iVal][1];
224  cellNo = addr[iVal][2];
225 
226  offset = corStart_[patchNo][faceNo];
227  procVfValues[patchNo][faceNo][cellNo+offset] = locVf[iVal];
228  }
229  iCorProc++;
230  }
231  }
232 
233  //Step 3. Calculate gradient at faces on processor patches
234  forAll(procPairs_, patchI)
235  {
236  label procPatchId = procPairs_[patchI];
237  if (procPatchId > -1)
238  {
239  fvsPatchVectorField& pgradf = gradIF.boundaryFieldRef()[procPatchId];
240  const List2<scalar> & pvf = procVfValues[patchI];
241  const List2<scalar> & pwf2= procWf2_[patchI];
242  const List2<vector> & pgdf= procGdf_[patchI];
243 
244  forAll(pgradf, iFace)
245  {
246  gf = vector::zero;
247  forAll(procGdf_[patchI][iFace], i)
248  {
249  gf = gf + pwf2[iFace][i]*pgdf[iFace][i]*
250  (pvf[iFace][i] - sF.boundaryField()[procPatchId][iFace]);
251  }
252  pgradf[iFace] = gf;
253  }
254 
255  //Update processor degenerate faces
256  const labelList& degProcFaces = procDegFaces_[patchI];
257  label degId = -1;
258  forAll(degProcFaces, iFace)
259  {
260  degId = degProcFaces[iFace];
261 
262  pgradf[degId] = nf_.boundaryField()[procPatchId][degId]*
263  sngF.boundaryField()[procPatchId][degId];
264  }
265  }
266  }
267 
268  return tgradIF;
269 };
270 
271 //
272 //END-OF-FILE
273 //
274 
275 
List< DynamicList< label > > procDegFaces_
tmp< surfaceVectorField > Grad(const volScalarField &iF)
Calculate gradient of volume scalar function on the faces.
labelHashTable< label > corProcIds_
List< face > faceList
surfaceVectorField nf_
Definition: fvscStencil.H:54
const fvMesh & mesh_
Definition: fvscStencil.H:44
List2< Triple< label > > corAddr_
forAll(Y, i)
Definition: QGDYEqn.H:36
DynamicList< label > internalDegFaces_