All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
extendedFaceStencilFormVfValues.H
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::extendedFaceStencilFormVfValues
28 Description
29  Function for filling of array which contain which will contain
30  the cells necessary for use by several processors.
31 Source file
32  extendedFaceStencilScalarDer.C
33  extendedFaceStencilScalarGradOpt.C
34  leastSquaresStencilOpt.C
35 \*---------------------------------------------------------------------------*/
36 
37 void formVfValues (const volScalarField& iF,List3<scalar>& procVfValues)
38 {
39  label cellId = -1;
40  forAll(procPairs_, patchI)
41  {
42  if (procPairs_[patchI] > -1)
43  {
44  procVfValues[patchI].resize(procWf2_[patchI].size());
45  forAll(procVfValues[patchI], faceI)
46  {
47  procVfValues[patchI][faceI].resize(procWf2_[patchI][faceI].size());
48  procVfValues[patchI][faceI] = 0.0; //make values zero
49  forAll(myProcPatchCells_[patchI][faceI], cellI)
50  {
51  cellId = myProcPatchCells_[patchI][faceI][cellI];
52  procVfValues[patchI][faceI][cellI] = iF.primitiveField()[cellId];
53  }
54  }
55  }
56  }
57 
58  //Step 1. Send field data to neighbouring processors (non-blocking mode)
59  PstreamBuffers pBuffers(Pstream::commsTypes::nonBlocking);
60  forAll(procPairs_, procI)
61  {
62  label procId = neigProcs_[procI];
63  DynamicList<scalar> locVf;
64 
65  if (procPairs_[procI] > -1) //patch proc pair
66  {
67  forAll(procVfValues[procI], faceI)
68  {
69  for(
70  label
71  cellI = 0;
72  cellI <= ownEnd_[procI][faceI];
73  cellI++
74  )
75  {
76  locVf.append(procVfValues[procI][faceI][cellI]);
77  }
78  }
79  }
80  else //corner connected process
81  {
82  label cellId = -1;
83  label addrId = corProcIds_[procId];
84  forAll(corCellIds_[addrId], iCellId)
85  {
86  cellId = corCellIds_[addrId][iCellId];
87  locVf.append(iF.primitiveField()[cellId]);
88  }
89  }
90 
91  UOPstream oProcStr(procId, pBuffers);
92  oProcStr << locVf;
93  }
94 
95  //Step 2. Recieve field data from neighbouring processors
96  pBuffers.finishedSends();
97  label iCorProc = 0;
98  forAll(procPairs_, procI)
99  {
100  //data size from processor, that is neighbouring through the patch
101  label procId = neigProcs_[procI];
102 
103  UIPstream iProcStr(procId, pBuffers);
104  List<scalar> locVf (iProcStr);
105 
106  if (procPairs_[procI] > -1)
107  {
108  label iVf = 0;
109  forAll(neiStart_[procI], iFace)
110  {
111  for(
112  label
113  iCell=neiStart_[procI][iFace];
114  iCell<=neiEnd_[procI][iFace];
115  iCell++
116  )
117  {
118  procVfValues[procI][iFace][iCell] =
119  locVf[iVf];
120  iVf++;
121  }
122  }
123  }
124  else
125  {
126  label patchNo = -1;
127  label faceNo = -1;
128  label cellNo = -1;
129  label offset = -1;
130 
131  const List<Triple<label> >& addr = corAddr_[iCorProc];
132 // Pout << "corAddr for scalar: " << "corAddt["<< iCorProc << "] = " << corAddr_[iCorProc] << endl;
133 
134  forAll(addr, iVal)
135  {
136  patchNo = addr[iVal][0];
137  faceNo = addr[iVal][1];
138  cellNo = addr[iVal][2];
139 
140  offset = corStart_[patchNo][faceNo];
141  procVfValues[patchNo][faceNo][cellNo+offset] = locVf[iVal];
142  }
143  iCorProc++;
144  }
145 
146  }
147 
148 // Info << "procVfValues = " << procVfValues << endl;
149 }
150 
151 
152 template<class FieldType>
153 void formVfValues (const GeometricField<FieldType, fvPatchField, volMesh>& iF,List<List3<scalar>>& procVfValues)
154 {
155  //set values from this domain
156  label cellId = -1;
157  label nComps = pTraits<FieldType>::nComponents;
158 
159  procVfValues.resize(nComps);
160  for(label compI = 0;compI <= nComps-1;compI++)
161  {
162  forAll(procPairs_, patchI)
163  {
164  if (procPairs_[patchI] > -1)
165  {
166  procVfValues[compI].resize(procVfValues[compI].size()+1);
167  procVfValues[compI][patchI].resize(procWf2_[patchI].size());
168 
169  forAll(procVfValues[compI][patchI], faceI)
170  {
171  procVfValues[compI][patchI][faceI].resize(procWf2_[patchI][faceI].size());
172 // procVfValues[patchI][faceI] = 0.0; //make values zero
173  forAll(procVfValues[compI][patchI][faceI], cellI)
174  {
175  procVfValues[compI][patchI][faceI][cellI] =0.0; //make values zero
176 
177 // cellId = myProcPatchCells_[patchI][faceI][cellI];
178  procVfValues[compI][patchI][faceI][cellI] = iF.primitiveField()[cellId].component(compI);
179 
180 
181  }
182  }
183  }
184  }
185  }
186 
187  //Step 1. Send field data to neighbouring processors (non-blocking mode)
188 
189  PstreamBuffers pBuffers(Pstream::commsTypes::nonBlocking);
190 
191  forAll(procPairs_, procI)
192  {
193  label procId = neigProcs_[procI];
194  List<DynamicList<scalar>> locVf;
195  locVf.resize(nComps);
196 
197  for(label compI = 0;compI <= nComps-1;compI++)
198  {
199  if (procPairs_[procI] > -1) //patch proc pair
200  {
201  forAll(procVfValues[compI][procI], faceI)
202  {
203  for(
204  label
205  cellI = 0;
206  cellI <= ownEnd_[procI][faceI];
207  cellI++
208  )
209  {
210 
211  locVf[compI].append(procVfValues[compI][procI][faceI][cellI]);
212  }
213  }
214 
215  }
216 
217  else //corner connected process
218  {
219  label cellId = -1;
220  label addrId = corProcIds_[procId];
221  forAll(corCellIds_[addrId], iCellId)
222  {
223  cellId = corCellIds_[addrId][iCellId];
224  locVf[compI].append(iF.primitiveField()[cellId].component(compI));
225  }
226  }
227  }
228  UOPstream oProcStr(procId, pBuffers);
229  oProcStr << locVf;
230  }
231  pBuffers.finishedSends();
232  //Step 2. Recieve field data from neighbouring processors
233  label iCorProc = 0;
234  forAll(procPairs_, procI)
235  {
236 // data size from processor, that is neighbouring through the patch
237 
238  label procId = neigProcs_[procI];
239 
240  UIPstream iProcStr(procId, pBuffers);
241 
242  List<List<scalar>> locVf(iProcStr);
243 
244  if (procPairs_[procI] > -1)
245  {
246  label iVf = 0;
247  forAll(neiStart_[procI], iFace)
248  {
249  for(
250  label
251  iCell=neiStart_[procI][iFace];
252  iCell<=neiEnd_[procI][iFace];
253  iCell++
254  )
255  {
256  for(label compI = 0;compI <= nComps-1;compI++)
257  procVfValues[compI][procI][iFace][iCell] =
258  locVf[compI][iVf];
259  iVf++;
260  }
261  }
262  }
263  else
264  {
265  label patchNo = -1;
266  label faceNo = -1;
267  label cellNo = -1;
268  label offset = -1;
269 
270  const List<Triple<label> >& addr = corAddr_[iCorProc];
271 
272 // Pout << "corAddr for template: " << "corAddt["<< iCorProc << "] = " << corAddr_[iCorProc] << endl;
273  forAll(addr, iVal)
274  {
275  patchNo = addr[iVal][0];
276  faceNo = addr[iVal][1];
277  cellNo = addr[iVal][2];
278  offset = corStart_[patchNo][faceNo];
279 
280  for(label compI = 0;compI <= nComps-1;compI++)
281  procVfValues[compI][patchNo][faceNo][cellNo+offset] = locVf[compI][iVal];
282  }
283  iCorProc++;
284 
285  }
286  }
287 
288 };
void formVfValues(const volScalarField &iF, List3< scalar > &procVfValues)
forAll(Y, i)
Definition: QGDYEqn.H:36