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