oversetFvMeshBaseTemplates.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) 2014-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
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 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volFields.H"
29 #include "fvMatrix.H"
30 #include "oversetFvPatchField.H"
34 #include "syncTools.H"
35 
36 // * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
37 
38 template<class Type>
40 (
41  Field<Type>& coeffs,
42  const labelUList& types,
43  const scalarList& factor,
44  const bool setHoleCellValue,
45  const label celli,
46  const label facei
47 ) const
48 {
49  const label cType = types[celli];
50  const scalar f = factor[celli];
51 
52  if (cType == cellCellStencil::INTERPOLATED)
53  {
54  coeffs[facei] *= 1.0-f;
55  }
56  else if (cType == cellCellStencil::HOLE)
57  {
58  // Disconnect hole cell from influence of neighbour
59  coeffs[facei] = pTraits<Type>::zero;
60  }
61  else if (cType == cellCellStencil::SPECIAL)
62  {
63  if (setHoleCellValue)
64  {
65  // Behave like hole
66  coeffs[facei] = pTraits<Type>::zero;
67  }
68  else
69  {
70  // Behave like interpolated
71  coeffs[facei] *= 1.0-f;
72  }
73  }
74 }
75 
76 template<class GeoField, class PatchType>
78 (
79  typename GeoField::Boundary& bfld,
80  const bool typeOnly
81 )
82 {
83  // Alternative (C++14)
84  //
85  // bfld.evaluate_if
86  // (
87  // [typeOnly](const auto& pfld) -> bool
88  // {
89  // return (typeOnly == bool(isA<PatchType>(pfld)));
90  // },
91  // UPstream::defaultCommsType
92  // );
93 
95  const label startOfRequests = UPstream::nRequests();
96 
97  for (auto& pfld : bfld)
98  {
99  if (typeOnly == bool(isA<PatchType>(pfld)))
100  {
101  pfld.initEvaluate(commsType);
102  }
103  }
104 
105  // Wait for outstanding requests (non-blocking)
106  UPstream::waitRequests(startOfRequests);
107 
108  for (auto& pfld : bfld)
109  {
110  if (typeOnly == bool(isA<PatchType>(pfld)))
111  {
112  pfld.evaluate(commsType);
113  }
114  }
115 }
116 
117 
118 template<class Type>
120 (
121  const fvMatrix<Type>& m
122 ) const
123 {
124  // Determine normalisation. This is normally the original diagonal plus
125  // remote contributions. This needs to be stabilised for hole cells
126  // which can have a zero diagonal. Assume that if any component has
127  // a non-zero diagonal the cell does not need stabilisation.
128  tmp<scalarField> tnorm(tmp<scalarField>::New(m.diag()));
129  scalarField& norm = tnorm.ref();
130 
131  // Add remote coeffs to duplicate behaviour of fvMatrix::addBoundaryDiag
132  const FieldField<Field, Type>& internalCoeffs = m.internalCoeffs();
133  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
134  {
135  forAll(internalCoeffs, patchi)
136  {
137  const labelUList& fc = mesh_.lduAddr().patchAddr(patchi);
138  const Field<Type>& intCoeffs = internalCoeffs[patchi];
139  const scalarField cmptCoeffs(intCoeffs.component(cmpt));
140  forAll(fc, i)
141  {
142  norm[fc[i]] += cmptCoeffs[i];
143  }
144  }
145  }
146 
147  // Count number of problematic cells
148  label nZeroDiag = 0;
149  forAll(norm, celli)
150  {
151  const scalar& n = norm[celli];
152  if (magSqr(n) < sqr(SMALL))
153  {
154  //Pout<< "For field " << m.psi().name()
155  // << " have diagonal " << n << " for cell " << celli
156  // << " at:" << cellCentres()[celli] << endl;
157  nZeroDiag++;
158  }
159  }
160 
161  if (debug)
162  {
163  Pout<< "For field " << m.psi().name() << " have zero diagonals for "
164  << returnReduce(nZeroDiag, sumOp<label>()) << " cells" << endl;
165  }
166 
167  if (returnReduceOr(nZeroDiag))
168  {
169  // Walk out the norm across hole cells
170 
171  const labelList& own = mesh_.faceOwner();
172  const labelList& nei = mesh_.faceNeighbour();
173  const cellCellStencilObject& overlap = Stencil::New(mesh_);
174  const labelUList& types = overlap.cellTypes();
175 
176  // label nHoles = 0;
177  scalarField extrapolatedNorm(norm);
178  forAll(types, celli)
179  {
180  if (types[celli] == cellCellStencil::HOLE)
181  {
182  extrapolatedNorm[celli] = -GREAT;
183  // ++nHoles;
184  }
185  }
186 
187  bitSet isFront(mesh_.nFaces());
188  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
189  {
190  label ownType = types[own[facei]];
191  label neiType = types[nei[facei]];
192  if
193  (
194  (ownType == cellCellStencil::HOLE)
195  != (neiType == cellCellStencil::HOLE)
196  )
197  {
198  isFront.set(facei);
199  }
200  }
201  labelList nbrTypes;
202  syncTools::swapBoundaryCellList(mesh_, types, nbrTypes);
203  for
204  (
205  label facei = mesh_.nInternalFaces();
206  facei < mesh_.nFaces();
207  ++facei
208  )
209  {
210  const label ownType = types[own[facei]];
211  const label neiType = nbrTypes[facei-mesh_.nInternalFaces()];
212  if
213  (
214  (ownType == cellCellStencil::HOLE)
215  != (neiType == cellCellStencil::HOLE)
216  )
217  {
218  isFront.set(facei);
219  }
220  }
221 
222 
223  while (true)
224  {
225  scalarField nbrNorm;
226  syncTools::swapBoundaryCellList(mesh_, extrapolatedNorm, nbrNorm);
227 
228  bitSet newIsFront(mesh_.nFaces());
229  scalarField newNorm(extrapolatedNorm);
230 
231  label nChanged = 0;
232  for (const label facei : isFront)
233  {
234  if (extrapolatedNorm[own[facei]] == -GREAT)
235  {
236  // Average owner cell, add faces to newFront
237  newNorm[own[facei]] = cellAverage
238  (
239  types,
240  nbrTypes,
241  extrapolatedNorm,
242  nbrNorm,
243  own[facei],
244  newIsFront
245  );
246  nChanged++;
247  }
248  if
249  (
250  mesh_.isInternalFace(facei)
251  && extrapolatedNorm[nei[facei]] == -GREAT
252  )
253  {
254  // Average nei cell, add faces to newFront
255  newNorm[nei[facei]] = cellAverage
256  (
257  types,
258  nbrTypes,
259  extrapolatedNorm,
260  nbrNorm,
261  nei[facei],
262  newIsFront
263  );
264  nChanged++;
265  }
266  }
267 
268  if (!returnReduceOr(nChanged))
269  {
270  break;
271  }
272 
273  // Transfer new front
274  extrapolatedNorm.transfer(newNorm);
275  isFront.transfer(newIsFront);
276  syncTools::syncFaceList(mesh_, isFront, maxEqOp<unsigned int>());
277  }
278 
279 
280  forAll(norm, celli)
281  {
282  scalar& n = norm[celli];
283  if (magSqr(n) < sqr(SMALL))
284  {
285  //Pout<< "For field " << m.psi().name()
286  // << " for cell " << celli
287  // << " at:" << cellCentres()[celli]
288  // << " have norm " << n
289  // << " have extrapolated norm " << extrapolatedNorm[celli]
290  // << endl;
291  // Override the norm
292  n = extrapolatedNorm[celli];
293  }
294  }
295  }
296  return tnorm;
297 }
298 
299 
300 template<class Type>
302 (
303  fvMatrix<Type>& m,
304  const scalarField& normalisation,
305  const bool setHoleCellValue,
306  const Type& holeCellValue
307 ) const
308 {
309  const cellCellStencilObject& overlap = Stencil::New(mesh_);
310  const List<scalarList>& wghts = overlap.cellInterpolationWeights();
311  const labelListList& stencil = overlap.cellStencil();
312  const scalarList& factor = overlap.cellInterpolationWeight();
313  const labelUList& types = overlap.cellTypes();
314 
315 
316  // Force asymmetric matrix (if it wasn't already)
317  scalarField& lower = m.lower();
318  scalarField& upper = m.upper();
319  Field<Type>& source = m.source();
320  scalarField& diag = m.diag();
321 
322 
323  // Get the addressing. Note that the addressing is now extended with
324  // any interpolation faces.
325  const lduAddressing& addr = lduAddr();
326  const labelUList& upperAddr = addr.upperAddr();
327  const labelUList& lowerAddr = addr.lowerAddr();
328  const lduInterfacePtrsList& interfaces = allInterfaces_;
329 
330  if (!isA<fvMeshPrimitiveLduAddressing>(addr))
331  {
333  << "Problem : addressing is not fvMeshPrimitiveLduAddressing"
334  << exit(FatalError);
335  }
336 
337 
338  // 1. Adapt lduMatrix for additional faces and new ordering
339  upper.setSize(upperAddr.size(), 0.0);
340  inplaceReorder(reverseFaceMap_, upper);
341  lower.setSize(lowerAddr.size(), 0.0);
342  inplaceReorder(reverseFaceMap_, lower);
343 
344 
345  //const label nOldInterfaces = dynamicMotionSolverFvMesh::interfaces().size();
346  const label nOldInterfaces = mesh_.fvMesh::interfaces().size();
347 
348 
349  if (interfaces.size() > nOldInterfaces)
350  {
351  // Extend matrix coefficients
352  m.internalCoeffs().setSize(interfaces.size());
353  m.boundaryCoeffs().setSize(interfaces.size());
354 
355  // 1b. Adapt for additional interfaces
356  for
357  (
358  label patchi = nOldInterfaces;
359  patchi < interfaces.size();
360  patchi++
361  )
362  {
363  const labelUList& fc = interfaces[patchi].faceCells();
364  m.internalCoeffs().set(patchi, new Field<Type>(fc.size(), Zero));
365  m.boundaryCoeffs().set(patchi, new Field<Type>(fc.size(), Zero));
366  }
367 
368  // 1c. Adapt field for additional interfaceFields (note: solver uses
369  // GeometricField::scalarInterfaces() to get hold of interfaces)
370  typedef GeometricField<Type, fvPatchField, volMesh> GeoField;
371 
372  typename GeoField::Boundary& bfld =
373  const_cast<GeoField&>(m.psi()).boundaryFieldRef();
374 
375  bfld.setSize(interfaces.size());
376 
377 
378  // This gets quite interesting: we do not want to add additional
379  // fvPatches (since direct correspondence to polyMesh) so instead
380  // add a reference to an existing processor patch
381  label addPatchi = 0;
382  for (label patchi = 0; patchi < nOldInterfaces; patchi++)
383  {
384  if (isA<processorFvPatch>(bfld[patchi].patch()))
385  {
386  addPatchi = patchi;
387  break;
388  }
389  }
390 
391  for
392  (
393  label patchi = nOldInterfaces;
394  patchi < interfaces.size();
395  patchi++
396  )
397  {
398  bfld.set
399  (
400  patchi,
401  new calculatedProcessorFvPatchField<Type>
402  (
403  interfaces[patchi],
404  bfld[addPatchi].patch(), // dummy processorFvPatch
405  m.psi()
406  )
407  );
408  }
409  }
410 
411 
412  // 2. Adapt fvMatrix level: faceFluxCorrectionPtr
413  // Question: do we need to do this?
414  // This seems to be set/used only by the gaussLaplacianScheme and
415  // fvMatrix:correction, both of which are outside the linear solver.
416 
417 
418  // Clear out existing connections on cells to be interpolated
419  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
420  // Note: could avoid doing the zeroing of the new faces since these
421  // are set to zero anyway.
422 
423  forAll(upperAddr, facei)
424  {
425  const label l = lowerAddr[facei];
426  scaleConnection(upper, types, factor, setHoleCellValue, l, facei);
427  const label u = upperAddr[facei];
428  scaleConnection(lower, types, factor, setHoleCellValue, u, facei);
429  /*
430  if (types[upperAddr[facei]] == cellCellStencil::INTERPOLATED)
431  {
432  // Disconnect upper from lower
433  label celli = upperAddr[facei];
434  lower[facei] *= 1.0-factor[celli];
435  }
436  if (types[lowerAddr[facei]] == cellCellStencil::INTERPOLATED)
437  {
438  // Disconnect lower from upper
439  label celli = lowerAddr[facei];
440  upper[facei] *= 1.0-factor[celli];
441  }
442  */
443  }
444 
445  for (label patchi = 0; patchi < nOldInterfaces; ++patchi)
446  {
447  const labelUList& fc = addr.patchAddr(patchi);
448  Field<Type>& bCoeffs = m.boundaryCoeffs()[patchi];
449  Field<Type>& iCoeffs = m.internalCoeffs()[patchi];
450  forAll(fc, i)
451  {
452  scaleConnection(bCoeffs, types, factor, setHoleCellValue, fc[i], i);
453 
454  scaleConnection(iCoeffs, types, factor, setHoleCellValue, fc[i], i);
455  }
456  /*
457  const labelUList& fc = addr.patchAddr(patchi);
458  Field<Type>& intCoeffs = m.internalCoeffs()[patchi];
459  Field<Type>& bouCoeffs = m.boundaryCoeffs()[patchi];
460  forAll(fc, i)
461  {
462  label celli = fc[i];
463  {
464  if (types[celli] == cellCellStencil::INTERPOLATED)
465  {
466  scalar f = factor[celli];
467  intCoeffs[i] *= 1.0-f;
468  bouCoeffs[i] *= 1.0-f;
469  }
470  else if (types[celli] == cellCellStencil::HOLE)
471  {
472  intCoeffs[i] = pTraits<Type>::zero;
473  bouCoeffs[i] = pTraits<Type>::zero;
474  }
475  }
476  }
477  */
478  }
479 
480 
481 
482  // Modify matrix
483  // ~~~~~~~~~~~~~
484 
485  // Do hole cells. Note: maybe put into interpolationCells() loop above?
486  forAll(types, celli)
487  {
488  if
489  (
490  types[celli] == cellCellStencil::HOLE
491  || (setHoleCellValue && types[celli] == cellCellStencil::SPECIAL)
492  )
493  {
494  const Type wantedValue
495  (
496  setHoleCellValue
497  ? holeCellValue
498  : m.psi()[celli]
499  );
500  diag[celli] = normalisation[celli];
501  source[celli] = normalisation[celli]*wantedValue;
502  }
503  else if
504  (
505  types[celli] == cellCellStencil::INTERPOLATED
506  || (!setHoleCellValue && types[celli] == cellCellStencil::SPECIAL)
507  )
508  {
509  const scalar f = factor[celli];
510  const scalarList& w = wghts[celli];
511  const labelList& nbrs = stencil[celli];
512  const labelList& nbrFaces = stencilFaces_[celli];
513  const labelList& nbrPatches = stencilPatches_[celli];
514 
515  diag[celli] *= (1.0-f);
516  source[celli] *= (1.0-f);
517 
518  forAll(nbrs, nbri)
519  {
520  const label patchi = nbrPatches[nbri];
521  const label facei = nbrFaces[nbri];
522 
523  if (patchi == -1)
524  {
525  const label nbrCelli = nbrs[nbri];
526  // Add the coefficients
527  const scalar s = normalisation[celli]*f*w[nbri];
528 
529  scalar& u = upper[facei];
530  scalar& l = lower[facei];
531  if (celli < nbrCelli)
532  {
533  diag[celli] += s;
534  u += -s;
535  }
536  else
537  {
538  diag[celli] += s;
539  l += -s;
540  }
541  }
542  else
543  {
544  // Patch face. Store in boundaryCoeffs. Note sign change.
545  //const label globalCelli = globalCellIDs[nbrs[nbri]];
546  //const label proci =
547  // globalNumbering.whichProcID(globalCelli);
548  //const label remoteCelli =
549  // globalNumbering.toLocal(proci, globalCelli);
550  //
551  //Pout<< "for cell:" << celli
552  // << " need weight from remote slot:" << nbrs[nbri]
553  // << " proc:" << proci << " remote cell:" << remoteCelli
554  // << " patch:" << patchi
555  // << " patchFace:" << facei
556  // << " weight:" << w[nbri]
557  // << endl;
558 
559  const scalar s = normalisation[celli]*f*w[nbri];
560  m.boundaryCoeffs()[patchi][facei] += pTraits<Type>::one*s;
561  m.internalCoeffs()[patchi][facei] += pTraits<Type>::one*s;
562 
563  // Note: do NOT add to diagonal - this is in the
564  // internalCoeffs and gets added to the diagonal
565  // inside fvMatrix::solve
566  }
567  }
568  }
569  /*
570  label startLabel = ownerStartAddr[celli];
571  label endLabel = ownerStartAddr[celli + 1];
572 
573  for (label facei = startLabel; facei < endLabel; facei++)
574  {
575  upper[facei] = 0.0;
576  }
577 
578  startLabel = addr.losortStartAddr()[celli];
579  endLabel = addr.losortStartAddr()[celli + 1];
580 
581  for (label i = startLabel; i < endLabel; i++)
582  {
583  label facei = losortAddr[i];
584  lower[facei] = 0.0;
585  }
586 
587  diag[celli] = normalisation[celli];
588  source[celli] = normalisation[celli]*m.psi()[celli];
589  */
590  }
591 
592 
593  //const globalIndex globalNumbering(V().size());
594  //labelList globalCellIDs(overlap.cellInterpolationMap().constructSize());
595  //forAll(V(), cellI)
596  //{
597  // globalCellIDs[cellI] = globalNumbering.toGlobal(cellI);
598  //}
599  //overlap.cellInterpolationMap().distribute(globalCellIDs);
600 
601 /*
602  forAll(cellIDs, i)
603  {
604  label celli = cellIDs[i];
605 
606  const scalar f = factor[celli];
607  const scalarList& w = wghts[celli];
608  const labelList& nbrs = stencil[celli];
609  const labelList& nbrFaces = stencilFaces_[celli];
610  const labelList& nbrPatches = stencilPatches_[celli];
611 
612  if (types[celli] == cellCellStencil::HOLE)
613  {
614  FatalErrorInFunction << "Found HOLE cell " << celli
615  << " at:" << C()[celli]
616  << " . Should this be in interpolationCells()????"
617  << abort(FatalError);
618  }
619  else
620  {
621  // Create interpolation stencil
622 
623  diag[celli] *= (1.0-f);
624  source[celli] *= (1.0-f);
625 
626  forAll(nbrs, nbri)
627  {
628  label patchi = nbrPatches[nbri];
629  label facei = nbrFaces[nbri];
630 
631  if (patchi == -1)
632  {
633  label nbrCelli = nbrs[nbri];
634 
635  // Add the coefficients
636  const scalar s = normalisation[celli]*f*w[nbri];
637 
638  scalar& u = upper[facei];
639  scalar& l = lower[facei];
640  if (celli < nbrCelli)
641  {
642  diag[celli] += s;
643  u += -s;
644  }
645  else
646  {
647  diag[celli] += s;
648  l += -s;
649  }
650  }
651  else
652  {
653  // Patch face. Store in boundaryCoeffs. Note sign change.
654  //const label globalCelli = globalCellIDs[nbrs[nbri]];
655  //const label proci =
656  // globalNumbering.whichProcID(globalCelli);
657  //const label remoteCelli =
658  // globalNumbering.toLocal(proci, globalCelli);
659  //
660  //Pout<< "for cell:" << celli
661  // << " need weight from remote slot:" << nbrs[nbri]
662  // << " proc:" << proci << " remote cell:" << remoteCelli
663  // << " patch:" << patchi
664  // << " patchFace:" << facei
665  // << " weight:" << w[nbri]
666  // << endl;
667 
668  const scalar s = normalisation[celli]*f*w[nbri];
669  m.boundaryCoeffs()[patchi][facei] += pTraits<Type>::one*s;
670  m.internalCoeffs()[patchi][facei] += pTraits<Type>::one*s;
671 
672  // Note: do NOT add to diagonal - this is in the
673  // internalCoeffs and gets added to the diagonal
674  // inside fvMatrix::solve
675  }
676  }
677 
678  //if (mag(diag[celli]) < SMALL)
679  //{
680  // Pout<< "for cell:" << celli
681  // << " at:" << this->C()[celli]
682  // << " diag:" << diag[celli] << endl;
683  //
684  // forAll(nbrs, nbri)
685  // {
686  // label patchi = nbrPatches[nbri];
687  // label facei = nbrFaces[nbri];
688  //
689  // const label globalCelli = globalCellIDs[nbrs[nbri]];
690  // const label proci =
691  // globalNumbering.whichProcID(globalCelli);
692  // const label remoteCelli =
693  // globalNumbering.toLocal(proci, globalCelli);
694  //
695  // Pout<< " need weight from slot:" << nbrs[nbri]
696  // << " proc:" << proci << " remote cell:"
697  // << remoteCelli
698  // << " patch:" << patchi
699  // << " patchFace:" << facei
700  // << " weight:" << w[nbri]
701  // << endl;
702  // }
703  // Pout<< endl;
704  //}
705  }
706  }
707 */
708 }
709 
710 
711 template<class Type>
713 (
714  fvMatrix<Type>& m,
715  const dictionary& dict
716 ) const
717 {
719  // Check we're running with bcs that can handle implicit matrix manipulation
720  typename GeoField::Boundary& bpsi =
721  const_cast<GeoField&>(m.psi()).boundaryFieldRef();
722 
723  bool hasOverset = false;
724  forAll(bpsi, patchi)
725  {
726  if (isA<oversetFvPatchField<Type>>(bpsi[patchi]))
727  {
728  hasOverset = true;
729  break;
730  }
731  }
732 
733  if (!hasOverset)
734  {
735  if (debug)
736  {
737  Pout<< "oversetFvMeshBase::solveOverset() :"
738  << " bypassing matrix adjustment for field " << m.psi().name()
739  << endl;
740  }
741  //return dynamicMotionSolverFvMesh::solve(m, dict);
742  return mesh_.fvMesh::solve(m, dict);
743  }
744 
745  if (debug)
746  {
747  Pout<< "oversetFvMeshBase::solveOverset() :"
748  << " adjusting matrix for interpolation for field "
749  << m.psi().name() << endl;
750  }
751 
752  // Calculate stabilised diagonal as normalisation for interpolation
753  const scalarField norm(normalisation(m));
754 
755  if (debug && mesh_.time().writeTime())
756  {
757  volScalarField scale
758  (
759  IOobject
760  (
761  m.psi().name() + "_normalisation",
762  mesh_.time().timeName(),
763  mesh_,
767  ),
768  mesh_,
770  );
771  scale.ref().field() = norm;
773  <
775  oversetFvPatchField<scalar>
776  >(scale.boundaryFieldRef(), false);
777  scale.write();
778 
779  if (debug)
780  {
781  Pout<< "oversetFvMeshBase::solveOverset() :"
782  << " writing matrix normalisation for field " << m.psi().name()
783  << " to " << scale.name() << endl;
784  }
785  }
786 
787 
788  // Switch to extended addressing (requires mesh::update() having been
789  // called)
790  active(true);
791 
792  // Adapt matrix
793  scalarField oldUpper(m.upper());
794  scalarField oldLower(m.lower());
795  FieldField<Field, Type> oldInt(m.internalCoeffs());
796  FieldField<Field, Type> oldBou(m.boundaryCoeffs());
797 
798  Field<Type> oldSource(m.source());
799  scalarField oldDiag(m.diag());
800 
801  const label oldSize = bpsi.size();
802 
803  // Insert the interpolation into the matrix (done inside
804  // oversetFvPatchField<Type>::manipulateMatrix)
805  m.boundaryManipulate(bpsi);
806 
807  // Swap psi values so added patches have patchNeighbourField
808  correctBoundaryConditions<GeoField, calculatedProcessorFvPatchField<Type>>
809  (
810  bpsi,
811  true
812  );
813 
814  // Use lower level solver
815  //SolverPerformance<Type> s(dynamicMotionSolverFvMesh::solve(m, dict));
816  SolverPerformance<Type> s(mesh_.fvMesh::solve(m, dict));
817 
818  // Restore boundary field
819  bpsi.setSize(oldSize);
820 
821  // Restore matrix
822  m.upper().transfer(oldUpper);
823  m.lower().transfer(oldLower);
824 
825  m.source().transfer(oldSource);
826  m.diag().transfer(oldDiag);
827 
828  m.internalCoeffs().transfer(oldInt);
829  m.boundaryCoeffs().transfer(oldBou);
830 
831 
832  const cellCellStencilObject& overlap = Stencil::New(mesh_);
833  const labelUList& types = overlap.cellTypes();
834  const labelList& own = mesh_.faceOwner();
835  const labelList& nei = mesh_.faceNeighbour();
836 
837  auto& psi =
838  const_cast<GeometricField<Type, fvPatchField, volMesh>&>(m.psi());
839 
840  // Mirror hole cell values next to calculated cells
841  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
842  {
843  const label ownType = types[own[facei]];
844  const label neiType = types[nei[facei]];
845 
846  if
847  (
848  ownType == cellCellStencil::HOLE
849  && neiType == cellCellStencil::CALCULATED)
850  {
851  psi[own[facei]] = psi[nei[facei]];
852  }
853  else if
854  (
855  ownType == cellCellStencil::CALCULATED
856  && neiType == cellCellStencil::HOLE
857  )
858  {
859  psi[nei[facei]] = psi[own[facei]];
860  }
861  }
862 
863  // Switch to original addressing
864  active(false);
866  return s;
867 }
868 
869 
870 template<class Type>
872 (
873  Ostream& os,
874  const fvMatrix<Type>& m,
875  const lduAddressing& addr,
876  const lduInterfacePtrsList& interfaces
877 ) const
878 {
879  os << "** Matrix **" << endl;
880  const labelUList& upperAddr = addr.upperAddr();
881  const labelUList& lowerAddr = addr.lowerAddr();
882  const labelUList& ownerStartAddr = addr.ownerStartAddr();
883  const labelUList& losortAddr = addr.losortAddr();
884 
885  const scalarField& lower = m.lower();
886  const scalarField& upper = m.upper();
887  const Field<Type>& source = m.source();
888  const scalarField& diag = m.diag();
889 
890 
891  // Invert patch addressing
892  labelListList cellToPatch(addr.size());
893  labelListList cellToPatchFace(addr.size());
894  {
895  forAll(interfaces, patchi)
896  {
897  if (interfaces.set(patchi))
898  {
899  const labelUList& fc = interfaces[patchi].faceCells();
900 
901  forAll(fc, i)
902  {
903  cellToPatch[fc[i]].append(patchi);
904  cellToPatchFace[fc[i]].append(i);
905  }
906  }
907  }
908  }
909 
910  forAll(source, celli)
911  {
912  os << "cell:" << celli << " diag:" << diag[celli]
913  << " source:" << source[celli] << endl;
914 
915  label startLabel = ownerStartAddr[celli];
916  label endLabel = ownerStartAddr[celli + 1];
917 
918  for (label facei = startLabel; facei < endLabel; facei++)
919  {
920  os << " face:" << facei
921  << " nbr:" << upperAddr[facei] << " coeff:" << upper[facei]
922  << endl;
923  }
924 
925  startLabel = addr.losortStartAddr()[celli];
926  endLabel = addr.losortStartAddr()[celli + 1];
927 
928  for (label i = startLabel; i < endLabel; i++)
929  {
930  label facei = losortAddr[i];
931  os << " face:" << facei
932  << " nbr:" << lowerAddr[facei] << " coeff:" << lower[facei]
933  << endl;
934  }
935 
936  forAll(cellToPatch[celli], i)
937  {
938  label patchi = cellToPatch[celli][i];
939  label patchFacei = cellToPatchFace[celli][i];
940 
941  os << " patch:" << patchi
942  << " patchface:" << patchFacei
943  << " intcoeff:" << m.internalCoeffs()[patchi][patchFacei]
944  << " boucoeff:" << m.boundaryCoeffs()[patchi][patchFacei]
945  << endl;
946  }
947  }
948  forAll(m.internalCoeffs(), patchi)
949  {
950  if (m.internalCoeffs().set(patchi))
951  {
952  const labelUList& fc = addr.patchAddr(patchi);
953 
954  os << "patch:" << patchi
955  //<< " type:" << interfaces[patchi].type()
956  << " size:" << fc.size() << endl;
957  if
958  (
959  interfaces.set(patchi)
960  && isA<processorLduInterface>(interfaces[patchi])
961  )
962  {
963  const processorLduInterface& ppp =
964  refCast<const processorLduInterface>(interfaces[patchi]);
965  os << "(processor with my:" << ppp.myProcNo()
966  << " nbr:" << ppp.neighbProcNo()
967  << ")" << endl;
968  }
969 
970  forAll(fc, i)
971  {
972  os << " cell:" << fc[i]
973  << " int:" << m.internalCoeffs()[patchi][i]
974  << " boun:" << m.boundaryCoeffs()[patchi][i]
975  << endl;
976  }
977  }
978  }
979 
980 
981  lduInterfaceFieldPtrsList interfaceFields =
982  m.psi().boundaryField().scalarInterfaces();
983  forAll(interfaceFields, inti)
984  {
985  if (interfaceFields.set(inti))
986  {
987  os << "interface:" << inti
988  << " if type:" << interfaceFields[inti].interface().type()
989  << " fld type:" << interfaceFields[inti].type() << endl;
990  }
991  }
992 
993  os << "** End of Matrix **" << endl;
994 }
995 
996 
997 template<class GeoField>
999 {
1000  auto& bfld = fld.boundaryFieldRef();
1001 
1003 
1004  const label startOfRequests = UPstream::nRequests();
1005 
1006  for (auto& pfld : bfld)
1007  {
1008  if (pfld.coupled())
1009  {
1010  //Pout<< "initEval of " << pfld.patch().name() << endl;
1011  pfld.initEvaluate(commsType);
1012  }
1013  }
1014 
1015  // Wait for outstanding requests (non-blocking)
1016  UPstream::waitRequests(startOfRequests);
1017 
1018  for (auto& pfld : bfld)
1019  {
1020  if (pfld.coupled())
1021  {
1022  //Pout<< "eval of " << pfld.patch().name() << endl;
1023  pfld.evaluate(commsType);
1024  }
1025  }
1026 }
1027 
1028 
1029 template<class GeoField>
1030 void Foam::oversetFvMeshBase::checkCoupledBC(const GeoField& fld)
1031 {
1032  Pout<< "** starting checking of " << fld.name() << endl;
1033 
1034  GeoField fldCorr(fld.name()+"_correct", fld);
1035  correctCoupledBoundaryConditions(fldCorr);
1036 
1037  const typename GeoField::Boundary& bfld = fld.boundaryField();
1038  const typename GeoField::Boundary& bfldCorr = fldCorr.boundaryField();
1039 
1040  forAll(bfld, patchi)
1041  {
1042  const typename GeoField::Patch& pfld = bfld[patchi];
1043  const typename GeoField::Patch& pfldCorr = bfldCorr[patchi];
1044 
1045  Pout<< "Patch:" << pfld.patch().name() << endl;
1046 
1047  forAll(pfld, i)
1048  {
1049  if (pfld[i] != pfldCorr[i])
1050  {
1051  Pout<< " " << i << " orig:" << pfld[i]
1052  << " corrected:" << pfldCorr[i] << endl;
1053  }
1054  }
1055  }
1056  Pout<< "** end of " << fld.name() << endl;
1057 }
1058 
1059 
1060 // ************************************************************************* //
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
const scalarField & diag() const
Definition: lduMatrix.C:163
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition: syncTools.H:465
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
commsTypes
Communications types.
Definition: UPstream.H:77
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const labelUList & losortStartAddr() const
Return losort start addressing.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
static const cellCellStencilObject & New(const fvMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests)
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
cellMask correctBoundaryConditions()
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1187
void scaleConnection(Field< Type > &coeffs, const labelUList &types, const scalarList &factor, const bool setHoleCellValue, const label celli, const label facei) const
Freeze values at holes.
Generic GeometricField class.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
tmp< scalarField > normalisation(const fvMatrix< Type > &m) const
Determine normalisation for interpolation. This equals the original diagonal except stabilised for ze...
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:134
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
static void correctCoupledBoundaryConditions(GeoField &fld)
Debug: correct coupled bc.
SolverPerformance< Type > solveOverset(fvMatrix< Type > &m, const dictionary &) const
Solve given dictionary with settings.
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1561
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual const labelUList & cellTypes() const
Return the cell type list.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
label size() const noexcept
Return number of equations.
Generic templated field type.
Definition: Field.H:62
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
virtual const labelUList & upperAddr() const =0
Return upper addressing.
void addInterpolation(fvMatrix< Type > &m, const scalarField &normalisation, const bool setHoleCellValue, const Type &holeCellValue) const
Manipulate the matrix to add the interpolation/set hole.
const cellCellStencilObject & overlap
Definition: correctPhi.H:57
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
const scalarField & lower() const
Definition: lduMatrix.C:268
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
static void correctBoundaryConditions(typename GeoField::Boundary &bfld, const bool typeOnly)
Correct boundary conditions of certain type (typeOnly = true) or explicitly not of the type (typeOnly...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const Field< Type > & field() const noexcept
Return const-reference to the field values.
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:397
const scalarField & upper() const
Definition: lduMatrix.C:203
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1171
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
const std::string patch
OpenFOAM patch number as a std::string.
const labelUList & losortAddr() const
Return losort addressing.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
label n
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition: fvMatrix.H:547
const volScalarField & psi
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition: fvMatrix.H:565
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1261
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition: typeInfo.H:87
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void write(Ostream &, const fvMatrix< Type > &, const lduAddressing &, const lduInterfacePtrsList &) const
Debug: print matrix.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Do not request registration (bool: false)
const labelUList & ownerStartAddr() const
Return owner start addressing.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrList.H:366
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127