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