fvMatrix.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
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 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFields.H"
30 #include "surfaceFields.H"
33 #include "coupledFvPatchFields.H"
34 #include "IndirectList.H"
35 #include "UniformList.H"
36 
37 #include "cyclicFvPatchField.H"
38 #include "cyclicAMIFvPatchField.H"
39 #include "cyclicACMIFvPatchField.H"
40 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
45 template<class Type>
46 template<class Type2>
48 (
49  const labelUList& addr,
50  const Field<Type2>& pf,
51  Field<Type2>& intf
52 ) const
53 {
54  if (addr.size() != pf.size())
55  {
57  << "addressing (" << addr.size()
58  << ") and field (" << pf.size() << ") are different sizes" << endl
59  << abort(FatalError);
60  }
61 
62  forAll(addr, facei)
63  {
64  intf[addr[facei]] += pf[facei];
65  }
66 }
67 
68 
69 template<class Type>
70 template<class Type2>
72 (
73  const labelUList& addr,
74  const tmp<Field<Type2>>& tpf,
75  Field<Type2>& intf
76 ) const
77 {
78  addToInternalField(addr, tpf(), intf);
79  tpf.clear();
80 }
81 
82 
83 template<class Type>
84 template<class Type2>
86 (
87  const labelUList& addr,
88  const Field<Type2>& pf,
89  Field<Type2>& intf
90 ) const
91 {
92  if (addr.size() != pf.size())
93  {
95  << "addressing (" << addr.size()
96  << ") and field (" << pf.size() << ") are different sizes" << endl
97  << abort(FatalError);
98  }
99 
100  forAll(addr, facei)
101  {
102  intf[addr[facei]] -= pf[facei];
103  }
104 }
105 
106 
107 template<class Type>
108 template<class Type2>
110 (
111  const labelUList& addr,
112  const tmp<Field<Type2>>& tpf,
113  Field<Type2>& intf
114 ) const
115 {
116  subtractFromInternalField(addr, tpf(), intf);
117  tpf.clear();
118 }
119 
120 
121 template<class Type>
123 (
124  scalarField& diag,
125  const direction solveCmpt
126 ) const
127 {
128  for (label fieldi = 0; fieldi < nMatrices(); ++fieldi)
129  {
130  const auto& bpsi = this->psi(fieldi).boundaryField();
131 
132  forAll(bpsi, ptfi)
133  {
134  const label patchi = globalPatchID(fieldi, ptfi);
135 
136  if (patchi != -1)
137  {
138  addToInternalField
139  (
140  lduAddr().patchAddr(patchi),
141  internalCoeffs_[patchi].component(solveCmpt),
142  diag
143  );
144  }
145  }
146  }
147 }
148 
149 
150 template<class Type>
152 {
153  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
154  {
155  const auto& bpsi = this->psi(fieldi).boundaryField();
156 
157  forAll(bpsi, ptfi)
158  {
159  const label patchi = globalPatchID(fieldi, ptfi);
160  if (patchi != -1)
161  {
162  addToInternalField
163  (
164  lduAddr().patchAddr(patchi),
165  cmptAv(internalCoeffs_[patchi]),
166  diag
167  );
168  }
169  }
170  }
171 }
172 
173 
174 template<class Type>
176 (
177  Field<Type>& source,
178  const bool couples
179 ) const
180 {
181  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
182  {
183  const auto& bpsi = this->psi(fieldi).boundaryField();
184 
185  forAll(bpsi, ptfi)
186  {
187  const fvPatchField<Type>& ptf = bpsi[ptfi];
188 
189  const label patchi = globalPatchID(fieldi, ptfi);
190 
191  if (patchi != -1)
192  {
193  const Field<Type>& pbc = boundaryCoeffs_[patchi];
194 
195  if (!ptf.coupled())
196  {
197  addToInternalField
198  (
199  lduAddr().patchAddr(patchi),
200  pbc,
201  source
202  );
203  }
204  else if (couples)
205  {
206  const tmp<Field<Type>> tpnf = ptf.patchNeighbourField();
207  const Field<Type>& pnf = tpnf();
208 
209  const labelUList& addr = lduAddr().patchAddr(patchi);
210 
211  forAll(addr, facei)
212  {
213  source[addr[facei]] +=
214  cmptMultiply(pbc[facei], pnf[facei]);
215  }
216  }
217  }
218  }
219  }
220 }
221 
222 
223 template<class Type>
224 template<template<class> class ListType>
226 (
227  const labelUList& cellLabels,
228  const ListType<Type>& values
229 )
230 {
231  const fvMesh& mesh = psi_.mesh();
232 
233  const cellList& cells = mesh.cells();
234  const labelUList& own = mesh.owner();
235  const labelUList& nei = mesh.neighbour();
236 
237  scalarField& Diag = diag();
238  Field<Type>& psi =
239  const_cast
240  <
241  GeometricField<Type, fvPatchField, volMesh>&
242  >(psi_).primitiveFieldRef();
243 
244 
245  // Following actions:
246  // - adjust local field psi
247  // - set local matrix to be diagonal (so adjust source)
248  // - cut connections to neighbours
249  // - make (on non-adjusted cells) contribution explicit
250 
251  if (symmetric() || asymmetric())
252  {
253  forAll(cellLabels, i)
254  {
255  const label celli = cellLabels[i];
256  const Type& value = values[i];
257 
258  for (const label facei : cells[celli])
259  {
260  const label patchi = mesh.boundaryMesh().patchID(facei);
261 
262  if (patchi == -1)
263  {
264  if (symmetric())
265  {
266  if (celli == own[facei])
267  {
268  source_[nei[facei]] -= upper()[facei]*value;
269  }
270  else
271  {
272  source_[own[facei]] -= upper()[facei]*value;
273  }
274 
275  upper()[facei] = 0.0;
276  }
277  else
278  {
279  if (celli == own[facei])
280  {
281  source_[nei[facei]] -= lower()[facei]*value;
282  }
283  else
284  {
285  source_[own[facei]] -= upper()[facei]*value;
286  }
287 
288  upper()[facei] = 0.0;
289  lower()[facei] = 0.0;
290  }
291  }
292  else
293  {
294  if (internalCoeffs_[patchi].size())
295  {
296  const label patchFacei =
297  mesh.boundaryMesh()[patchi].whichFace(facei);
298 
299  internalCoeffs_[patchi][patchFacei] = Zero;
300  boundaryCoeffs_[patchi][patchFacei] = Zero;
301  }
302  }
303  }
304  }
305  }
306 
307  // Note: above loop might have affected source terms on adjusted cells
308  // so make sure to adjust them afterwards
309  forAll(cellLabels, i)
310  {
311  const label celli = cellLabels[i];
312  const Type& value = values[i];
313 
314  psi[celli] = value;
315  source_[celli] = value*Diag[celli];
316  }
317 }
318 
319 
320 template<class Type>
321 bool Foam::fvMatrix<Type>::checkImplicit(const label fieldi)
322 {
323  const auto& bpsi = this->psi(fieldi).boundaryField();
324 
325  word idName;
326  forAll(bpsi, patchi)
327  {
328  if (bpsi[patchi].useImplicit())
329  {
330  if (debug)
331  {
332  Pout<< "fvMatrix<Type>::checkImplicit "
333  << " field:" << this->psi(fieldi).name()
334  << " on mesh:"
335  << this->psi(fieldi).mesh().name()
336  << " patch:" << bpsi[patchi].patch().name()
337  << endl;
338  }
339 
340  idName += Foam::name(patchi);
341  useImplicit_ = true;
342  }
343  }
344 
345  if (useImplicit_)
346  {
347  lduAssemblyName_ = word("lduAssembly") + idName;
348  }
349 
350  return !idName.empty();
351 }
352 
353 
354 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
355 
356 template<class Type>
358 (
360  const dimensionSet& ds
361 )
362 :
363  lduMatrix(psi.mesh()),
364  psi_(psi),
365  useImplicit_(false),
366  lduAssemblyName_(),
367  nMatrix_(0),
368  dimensions_(ds),
369  source_(psi.size(), Zero),
370  internalCoeffs_(psi.mesh().boundary().size()),
371  boundaryCoeffs_(psi.mesh().boundary().size())
372 {
374  << "Constructing fvMatrix<Type> for field " << psi_.name() << endl;
375 
376  checkImplicit();
377 
378  forAll(psi.mesh().boundary(), patchi)
379  {
380  internalCoeffs_.set
381  (
382  patchi,
383  new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
384  );
385 
386  boundaryCoeffs_.set
387  (
388  patchi,
389  new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
390  );
391  }
392 
393  auto& psiRef = this->psi(0);
394  const label currentStatePsi = psiRef.eventNo();
395  psiRef.boundaryFieldRef().updateCoeffs();
396  psiRef.eventNo() = currentStatePsi;
397 }
398 
399 
400 template<class Type>
402 :
403  lduMatrix(fvm),
404  psi_(fvm.psi_),
405  useImplicit_(fvm.useImplicit_),
406  lduAssemblyName_(fvm.lduAssemblyName_),
407  nMatrix_(fvm.nMatrix_),
408  dimensions_(fvm.dimensions_),
409  source_(fvm.source_),
410  internalCoeffs_(fvm.internalCoeffs_),
411  boundaryCoeffs_(fvm.boundaryCoeffs_)
412 {
414  << "Copying fvMatrix<Type> for field " << psi_.name() << endl;
415 
416  if (fvm.faceFluxCorrectionPtr_)
417  {
418  faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
419  (
420  *(fvm.faceFluxCorrectionPtr_)
421  );
422  }
423 }
424 
425 
426 template<class Type>
427 Foam::fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type>>& tmat)
428 :
429  lduMatrix(tmat.constCast(), tmat.movable()),
430  psi_(tmat().psi_),
431  useImplicit_(tmat().useImplicit_),
432  lduAssemblyName_(tmat().lduAssemblyName_),
433  nMatrix_(tmat().nMatrix_),
434  dimensions_(tmat().dimensions_),
435  source_(tmat.constCast().source_, tmat.movable()),
436  internalCoeffs_(tmat.constCast().internalCoeffs_, tmat.movable()),
437  boundaryCoeffs_(tmat.constCast().boundaryCoeffs_, tmat.movable())
438 {
440  << "Copy/move fvMatrix<Type> for field " << psi_.name() << endl;
441 
442  if (tmat().faceFluxCorrectionPtr_)
443  {
444  if (tmat.movable())
445  {
446  faceFluxCorrectionPtr_ =
447  std::move(tmat.constCast().faceFluxCorrectionPtr_);
448  }
449  else if (tmat().faceFluxCorrectionPtr_)
450  {
451  faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
452  (
453  *(tmat().faceFluxCorrectionPtr_)
454  );
455  }
456  }
457 
458  tmat.clear();
459 }
460 
461 
462 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
463 
464 template<class Type>
466 {
468  << "Destroying fvMatrix<Type> for field " << psi_.name() << endl;
469 
470  subMatrices_.clear();
471 }
472 
473 
474 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
475 
476 template<class Type>
478 (
479  lduInterfaceFieldPtrsList& interfaces,
480  PtrDynList<lduInterfaceField>& newInterfaces
481 )
482 {
483  interfaces.setSize(internalCoeffs_.size());
484  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
485  {
486  const auto& bpsi = this->psi(fieldi).boundaryField();
487  lduInterfaceFieldPtrsList fieldInterfaces(bpsi.scalarInterfaces());
488 
489  forAll (fieldInterfaces, patchi)
490  {
491  label globalPatchID = lduMeshPtr()->patchMap()[fieldi][patchi];
492 
493  if (globalPatchID != -1)
494  {
495  if (fieldInterfaces.set(patchi))
496  {
497  if (isA<cyclicLduInterfaceField>(bpsi[patchi]))
498  {
499  newInterfaces.append
500  (
502  (
503  refCast<const fvPatch>
504  (
505  lduMeshPtr()->interfaces()[globalPatchID]
506  ),
507  bpsi[patchi].internalField()
508  )
509  );
510  interfaces.set(globalPatchID, &newInterfaces.last());
511 
512  }
513  else if (isA<cyclicAMILduInterfaceField>(bpsi[patchi]))
514  {
515  newInterfaces.append
516  (
518  (
519  refCast<const fvPatch>
520  (
521  lduMeshPtr()->interfaces()[globalPatchID]
522  ),
523  bpsi[patchi].internalField()
524  )
525  );
526  interfaces.set(globalPatchID, &newInterfaces.last());
527  }
528  else if (isA<cyclicACMILduInterfaceField>(bpsi[patchi]))
529  {
530  newInterfaces.append
531  (
533  (
534  refCast<const fvPatch>
535  (
536  lduMeshPtr()->interfaces()[globalPatchID]
537  ),
538  bpsi[patchi].internalField()
539  )
540  );
541  interfaces.set(globalPatchID, &newInterfaces.last());
542  }
543  else
544  {
545  interfaces.set(globalPatchID, &fieldInterfaces[patchi]);
546  }
547  }
548  }
549  }
550  }
551 }
552 
553 
554 template<class Type>
556 (
557  label fieldi,
558  const FieldField<Field, Type>& fluxContrib,
559  FieldField<Field, Type>& contrib,
560  bool internal
561 ) const
562 {
563  const lduPrimitiveMeshAssembly* ptr = lduMeshPtr();
564 
565  const labelList& patchMap = ptr->patchMap()[fieldi];
566 
567  forAll(contrib, patchi)
568  {
569  const label globalPtchId = patchMap[patchi];
570 
571  if (globalPtchId != -1)
572  {
573  // Cache contrib before overwriting
574  const Field<Type> saveContrib(fluxContrib[globalPtchId]);
575  contrib[patchi].setSize(psi_.boundaryField()[patchi].size()),
576  contrib[patchi] = pTraits<Type>::zero;
577 
578  if (internal)
579  {
580  contrib[patchi] =
582  (
583  saveContrib,
584  psi_.boundaryField()[patchi].patchInternalField()
585  );
586  }
587  else
588  {
589  if (this->psi(fieldi).boundaryField()[patchi].coupled())
590  {
591  contrib[patchi] =
593  (
594  saveContrib,
595  psi_.boundaryField()[patchi].patchNeighbourField()
596  );
597  }
598  }
599  }
600  else if (globalPtchId == -1)
601  {
602  const polyPatch& pp =
603  this->psi(fieldi).mesh().boundaryMesh()[patchi];
604 
605  if (pp.masterImplicit())
606  {
607  label virtualPatch =
608  ptr->patchLocalToGlobalMap()[fieldi][patchi];
609 
610  const label nbrPatchId = pp.neighbPolyPatchID();
611 
612  // Copy contrib before overwriting
613  const Field<Type> saveContrib(fluxContrib[virtualPatch]);
614 
615  Field<Type>& coeffs = contrib[patchi];
616  Field<Type>& nbrCoeffs = contrib[nbrPatchId];
617 
618  coeffs.setSize(psi_.boundaryField()[patchi].size());
619  nbrCoeffs.setSize(psi_.boundaryField()[nbrPatchId].size());
620 
621  coeffs = pTraits<Type>::zero;
622  nbrCoeffs = pTraits<Type>::zero;
623 
624  // nrb cells
625  const labelList& nbrCellIds =
626  ptr->cellBoundMap()[fieldi][patchi];
627 
628  const labelList& cellIds =
629  ptr->cellBoundMap()[fieldi][nbrPatchId];
630 
631  const GeometricField<Type, fvPatchField, volMesh>& psi =
632  this->psi(fieldi);
633 
634  forAll(saveContrib, subFaceI)
635  {
636  const label faceId =
637  ptr->facePatchFaceMap()[fieldi][patchi][subFaceI];
638  const label nbrFaceId =
639  ptr->facePatchFaceMap()[fieldi][nbrPatchId][subFaceI];
640 
641  const label nbrCellId = nbrCellIds[subFaceI];
642  const label cellId = cellIds[subFaceI];
643 
644  if (internal)
645  {
646  coeffs[faceId] +=
647  cmptMultiply(saveContrib[subFaceI], psi[cellId]);
648 
649  nbrCoeffs[nbrFaceId] +=
650  cmptMultiply(saveContrib[subFaceI], psi[nbrCellId]);
651  }
652  else //boundary
653  {
654  coeffs[faceId] +=
655  cmptMultiply(saveContrib[subFaceI], psi[nbrCellId]);
656 
657  nbrCoeffs[nbrFaceId] +=
658  cmptMultiply(saveContrib[subFaceI], psi[cellId]);
659  }
660  }
661  }
662  }
663  }
664 }
665 
666 
667 template<class Type>
669 {
670  // If it is a multi-fvMatrix needs correct internalCoeffs and
671  // boundaryCoeffs size
672  if (nMatrix_ > 0)
673  {
674  label interfaceI(0);
675  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
676  {
677  const auto& psi = this->psi(fieldi);
678 
679  forAll(psi.mesh().boundary(), patchi)
680  {
681  interfaceI++;
682  }
683  }
684  internalCoeffs_.setSize(interfaceI);
685  boundaryCoeffs_.setSize(interfaceI);
686 
687  interfaceI = 0;
688  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
689  {
690  const auto& psi = this->psi(fieldi);
691 
692  forAll(psi.mesh().boundary(), patchi)
693  {
694  internalCoeffs_.set
695  (
696  interfaceI,
697  new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
698  );
699 
700  boundaryCoeffs_.set
701  (
702  interfaceI,
703  new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
704  );
705  interfaceI++;
706  }
707  }
708  }
709 
710  for (label i=0; i < nMatrices(); ++i)
711  {
712  const auto& bpsi = this->psi(i).boundaryField();
713 
714  // Cache to-be implicit internal/boundary
715  FieldField<Field, Type> boundary(bpsi.size());
716  FieldField<Field, Type> internal(bpsi.size());
717 
718  label implicit = 0;
719  forAll(bpsi, patchI)
720  {
721  label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
722  if (globalPatchId == -1)
723  {
724  boundary.set
725  (
726  implicit,
727  matrix(i).boundaryCoeffs()[patchI].clone()
728  );
729  internal.set
730  (
731  implicit,
732  matrix(i).internalCoeffs()[patchI].clone()
733  );
734  implicit++;
735  }
736  }
737 
738  // Update non-implicit patches (re-order)
739  forAll(bpsi, patchI)
740  {
741  label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
742  if (globalPatchId != -1)
743  {
744  if (matrix(i).internalCoeffs().set(patchI))
745  {
746  internalCoeffs_.set
747  (
748  globalPatchId,
749  matrix(i).internalCoeffs()[patchI].clone()
750  );
751  }
752 
753  if (matrix(i).boundaryCoeffs().set(patchI))
754  {
755  boundaryCoeffs_.set
756  (
757  globalPatchId,
758  matrix(i).boundaryCoeffs()[patchI].clone()
759  );
760  }
761  }
762  }
763 
764  // Store implicit patches at the end of the list
765  implicit = 0;
766  forAll(bpsi, patchI)
767  {
768  label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
769  if (globalPatchId == -1)
770  {
771  const label implicitPatchId =
772  lduMeshPtr()->patchLocalToGlobalMap()[i][patchI];
773 
774  internalCoeffs_.set
775  (
776  implicitPatchId, internal[implicit].clone()
777  );
778  boundaryCoeffs_.set
779  (
780  implicitPatchId, boundary[implicit].clone()
781  );
782 
783  implicit++;
784  }
785  }
786  }
787 
788 // forAll(internalCoeffs_, patchI)
789 // {
790 // DebugVar(patchI)
791 // DebugVar(internalCoeffs_[patchI])
792 // DebugVar(boundaryCoeffs_[patchI])
793 // }
794 }
795 
796 
797 template<class Type>
799 {
800  for (label i=0; i < nMatrices(); ++i)
801  {
802  forAll(psi(i).boundaryField(), patchI)
803  {
804  label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
805 
806  if (globalPatchId == -1)
807  {
808  psi(i).boundaryFieldRef()[patchI].manipulateMatrix
809  (
810  *this,
811  i,
812  cmp
813  );
814  }
815  }
816  }
817 }
818 
819 
820 template<class Type>
822 {
823  const labelListList& faceMap = lduMeshPtr()->faceMap();
824  const labelList& cellMap = lduMeshPtr()->cellOffsets();
825 
826  label newFaces = lduMeshPtr()->lduAddr().upperAddr().size();
827  label newCells = lduMeshPtr()->lduAddr().size();
828 
829  scalarField lowerAssemb(newFaces, Zero);
830  scalarField upperAssemb(newFaces, Zero);
831  scalarField diagAssemb(newCells, Zero);
832  Field<Type> sourceAssemb(newCells, Zero);
833 
834  bool asymmetricAssemby = false;
835  for (label i=0; i < nMatrices(); ++i)
836  {
837  if (matrix(i).asymmetric())
838  {
839  asymmetricAssemby = true;
840  }
841  }
842  // Move append contents into intermediate list
843  for (label i=0; i < nMatrices(); ++i)
844  {
845  if (asymmetricAssemby)
846  {
847  const scalarField lowerSub(matrix(i).lower());
848  forAll(lowerSub, facei)
849  {
850  lowerAssemb[faceMap[i][facei]] = lowerSub[facei];
851  }
852  }
853 
854  const scalarField upperSub(matrix(i).upper());
855  const scalarField diagSub(matrix(i).diag());
856  const Field<Type> sourceSub(matrix(i).source());
857 
858  forAll(upperSub, facei)
859  {
860  upperAssemb[faceMap[i][facei]] = upperSub[facei];
861  }
862 
863  forAll(diagSub, celli)
864  {
865  const label globalCelli = cellMap[i] + celli;
866  diagAssemb[globalCelli] = diagSub[celli];
867  sourceAssemb[globalCelli] = sourceSub[celli];
868  }
869  }
870 
871  if (asymmetricAssemby)
872  {
873  lower().setSize(newFaces, Zero);
874  lower() = lowerAssemb;
875  }
876  upper().setSize(newFaces, Zero);
877  upper() = upperAssemb;
878 
879  diag().setSize(newCells, Zero);
880  diag() = diagAssemb;
882  source().setSize(newCells, Zero);
883  source() = sourceAssemb;
884 }
885 
886 
887 template<class Type>
889 {
890  return
891  (
892  psi_.mesh().thisDb().objectRegistry::template getObjectPtr
893  <
895  > (lduAssemblyName_)
896  );
897 }
898 
899 
900 template<class Type>
902 {
903  return
904  (
905  psi_.mesh().thisDb().objectRegistry::template cfindObject
906  <
908  > (lduAssemblyName_)
909  );
910 }
911 
912 
913 template<class Type>
915 {
916  lduPrimitiveMeshAssembly* ptr = lduMeshPtr();
917 
918  IOobject io
919  (
920  lduAssemblyName_,
921  psi_.mesh().time().timeName(),
922  psi_.mesh().thisDb(),
926  );
927 
928  UPtrList<lduMesh> uMeshPtr(nMatrices());
929 
931  uFieldPtr(nMatrices());
932 
933  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
934  {
935  const fvMesh& meshi = this->psi(fieldi).mesh();
936  uMeshPtr.set
937  (
938  fieldi,
939  &const_cast<fvMesh&>(meshi)
940  );
941  uFieldPtr.set(fieldi, &this->psi(fieldi));
942  }
943 
944  if (!ptr)
945  {
946  lduPrimitiveMeshAssembly* lduAssemMeshPtr =
947  new lduPrimitiveMeshAssembly(io, uMeshPtr);
948 
949  lduAssemMeshPtr->store();
950  lduAssemMeshPtr->update(uFieldPtr);
951 
952  Info
953  << "Creating lduPrimitiveAssembly: " << lduAssemblyName_ << endl;
954  }
955  else if
956  (
957  psi_.mesh().changing() && !psi_.mesh().upToDatePoints(*ptr)
958  )
959  {
960  // Clear losortPtr_, ownerStartPtr_, losortStartPtr_
961  ptr->lduAddr().clearOut();
962  ptr->update(uFieldPtr);
963  psi_.mesh().setUpToDatePoints(*ptr);
964 
965  Info
966  << "Updating lduPrimitiveAssembly: " << lduAssemblyName_ << endl;
967  }
968  else
969  {
970  Info
971  << "Using lduPrimitiveAssembly: " << lduAssemblyName_ << endl;
972  }
973 }
974 
975 
976 template<class Type>
978 (
979  const labelUList& cellLabels,
980  const Type& value
981 )
982 {
983  this->setValuesFromList(cellLabels, UniformList<Type>(value));
984 }
985 
986 
987 template<class Type>
989 (
990  const labelUList& cellLabels,
991  const UList<Type>& values
992 )
993 {
994  this->setValuesFromList(cellLabels, values);
995 }
996 
997 
998 template<class Type>
1000 (
1001  const labelUList& cellLabels,
1003 )
1005  this->setValuesFromList(cellLabels, values);
1006 }
1007 
1008 
1009 template<class Type>
1011 (
1012  const label celli,
1013  const Type& value,
1014  const bool forceReference
1015 )
1016 {
1017  if ((forceReference || psi_.needReference()) && celli >= 0)
1018  {
1019  source()[celli] += diag()[celli]*value;
1020  diag()[celli] += diag()[celli];
1021  }
1022 }
1023 
1024 
1025 template<class Type>
1027 (
1028  const labelUList& cellLabels,
1029  const Type& value,
1030  const bool forceReference
1031 )
1032 {
1033  if (forceReference || psi_.needReference())
1034  {
1035  forAll(cellLabels, celli)
1036  {
1037  const label cellId = cellLabels[celli];
1038  if (cellId >= 0)
1039  {
1040  source()[cellId] += diag()[cellId]*value;
1041  diag()[cellId] += diag()[cellId];
1042  }
1043  }
1044  }
1045 }
1046 
1047 
1048 template<class Type>
1050 (
1051  const labelUList& cellLabels,
1052  const UList<Type>& values,
1053  const bool forceReference
1054 )
1055 {
1056  if (forceReference || psi_.needReference())
1057  {
1058  forAll(cellLabels, celli)
1059  {
1060  const label cellId = cellLabels[celli];
1061  if (cellId >= 0)
1062  {
1063  source()[cellId] += diag()[cellId]*values[celli];
1064  diag()[cellId] += diag()[cellId];
1065  }
1066  }
1067  }
1068 }
1069 
1070 
1071 template<class Type>
1072 void Foam::fvMatrix<Type>::addFvMatrix(fvMatrix& matrix)
1073 {
1074  subMatrices_.append(matrix.clone());
1075  ++nMatrix_;
1076 
1077  if (dimensions_ != matrix.dimensions())
1078  {
1080  << "incompatible dimensions for matrix addition "
1081  << endl << " "
1082  << "[" << dimensions_ << " ] "
1083  << " [" << matrix.dimensions() << " ]"
1084  << abort(FatalError);
1085  }
1086 
1087  for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
1088  {
1089  if (checkImplicit(fieldi))
1090  {
1091  break;
1092  }
1093  }
1095  internalCoeffs_.clear();
1096  boundaryCoeffs_.clear();
1097 }
1098 
1099 
1100 template<class Type>
1101 void Foam::fvMatrix<Type>::relax(const scalar alpha)
1102 {
1103  if (alpha <= 0)
1104  {
1105  return;
1106  }
1107 
1109  << "Relaxing " << psi_.name() << " by " << alpha << endl;
1110 
1111  Field<Type>& S = source();
1112  scalarField& D = diag();
1113 
1114  // Store the current unrelaxed diagonal for use in updating the source
1115  scalarField D0(D);
1116 
1117  // Calculate the sum-mag off-diagonal from the interior faces
1118  scalarField sumOff(D.size(), Zero);
1119  sumMagOffDiag(sumOff);
1120 
1121  // Handle the boundary contributions to the diagonal
1122  forAll(psi_.boundaryField(), patchi)
1123  {
1124  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1125 
1126  if (ptf.size())
1127  {
1128  const labelUList& pa = lduAddr().patchAddr(patchi);
1129  Field<Type>& iCoeffs = internalCoeffs_[patchi];
1130 
1131  if (ptf.coupled())
1132  {
1133  const Field<Type>& pCoeffs = boundaryCoeffs_[patchi];
1134 
1135  // For coupled boundaries add the diagonal and
1136  // off-diagonal contributions
1137  forAll(pa, face)
1138  {
1139  D[pa[face]] += component(iCoeffs[face], 0);
1140  sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
1141  }
1142  }
1143  else
1144  {
1145  // For non-coupled boundaries add the maximum magnitude diagonal
1146  // contribution to ensure stability
1147  forAll(pa, face)
1148  {
1149  D[pa[face]] += cmptMax(cmptMag(iCoeffs[face]));
1150  }
1151  }
1152  }
1153  }
1154 
1155 
1156  if (debug)
1157  {
1158  // Calculate amount of non-dominance.
1159  label nNon = 0;
1160  scalar maxNon = 0.0;
1161  scalar sumNon = 0.0;
1162  forAll(D, celli)
1163  {
1164  scalar d = (sumOff[celli] - D[celli])/mag(D[celli]);
1165 
1166  if (d > 0)
1167  {
1168  nNon++;
1169  maxNon = max(maxNon, d);
1170  sumNon += d;
1171  }
1172  }
1173 
1174  reduce(nNon, sumOp<label>(), UPstream::msgType(), psi_.mesh().comm());
1175  reduce
1176  (
1177  maxNon,
1178  maxOp<scalar>(),
1180  psi_.mesh().comm()
1181  );
1182  reduce
1183  (
1184  sumNon,
1185  sumOp<scalar>(),
1187  psi_.mesh().comm()
1188  );
1189  sumNon /= returnReduce
1190  (
1191  D.size(),
1192  sumOp<label>(),
1194  psi_.mesh().comm()
1195  );
1196 
1198  << "Matrix dominance test for " << psi_.name() << nl
1199  << " number of non-dominant cells : " << nNon << nl
1200  << " maximum relative non-dominance : " << maxNon << nl
1201  << " average relative non-dominance : " << sumNon << nl
1202  << endl;
1203  }
1204 
1205 
1206  // Ensure the matrix is diagonally dominant...
1207  // Assumes that the central coefficient is positive and ensures it is
1208  forAll(D, celli)
1209  {
1210  D[celli] = max(mag(D[celli]), sumOff[celli]);
1211  }
1212 
1213  // ... then relax
1214  D /= alpha;
1215 
1216  // Now remove the diagonal contribution from coupled boundaries
1217  forAll(psi_.boundaryField(), patchi)
1218  {
1219  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1220 
1221  if (ptf.size())
1222  {
1223  const labelUList& pa = lduAddr().patchAddr(patchi);
1224  Field<Type>& iCoeffs = internalCoeffs_[patchi];
1225 
1226  if (ptf.coupled())
1227  {
1228  forAll(pa, face)
1229  {
1230  D[pa[face]] -= component(iCoeffs[face], 0);
1231  }
1232  }
1233  else
1234  {
1235  forAll(pa, face)
1236  {
1237  D[pa[face]] -= cmptMin(iCoeffs[face]);
1238  }
1239  }
1240  }
1241  }
1243  // Finally add the relaxation contribution to the source.
1244  S += (D - D0)*psi_.primitiveField();
1245 }
1246 
1247 
1248 template<class Type>
1250 {
1251  word name = psi_.select
1252  (
1253  psi_.mesh().data().isFinalIteration()
1254  );
1255 
1256  scalar relaxCoeff = 0;
1257 
1258  if (psi_.mesh().relaxEquation(name, relaxCoeff))
1259  {
1260  relax(relaxCoeff);
1261  }
1262 }
1263 
1264 
1265 template<class Type>
1267 (
1268  typename GeometricField<Type, fvPatchField, volMesh>::
1269  Boundary& bFields
1270 )
1271 {
1272  forAll(bFields, patchi)
1273  {
1274  bFields[patchi].manipulateMatrix(*this);
1275  }
1276 }
1277 
1278 
1279 template<class Type>
1281 {
1282  auto tdiag = tmp<scalarField>::New(diag());
1283  addCmptAvBoundaryDiag(tdiag.ref());
1284  return tdiag;
1285 }
1286 
1287 
1288 template<class Type>
1290 {
1292 
1293  forAll(psi_.boundaryField(), patchi)
1294  {
1295  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1296 
1297  if (!ptf.coupled() && ptf.size())
1298  {
1299  addToInternalField
1300  (
1301  lduAddr().patchAddr(patchi),
1302  internalCoeffs_[patchi],
1303  tdiag.ref()
1304  );
1305  }
1306  }
1307 
1308  return tdiag;
1309 }
1310 
1311 
1312 template<class Type>
1314 {
1315  auto tAphi = volScalarField::New
1316  (
1317  "A(" + psi_.name() + ')',
1318  psi_.mesh(),
1319  dimensions_/psi_.dimensions()/dimVol,
1321  );
1322 
1323  tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().V();
1324  tAphi.ref().correctBoundaryConditions();
1326  return tAphi;
1327 }
1328 
1329 
1330 template<class Type>
1333 {
1335  (
1336  "H(" + psi_.name() + ')',
1337  psi_.mesh(),
1338  dimensions_/dimVol,
1340  );
1341  auto& Hphi = tHphi.ref();
1342 
1343  // Loop over field components
1344  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
1345  {
1346  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
1347 
1348  scalarField boundaryDiagCmpt(psi_.size(), Zero);
1349  addBoundaryDiag(boundaryDiagCmpt, cmpt);
1350  boundaryDiagCmpt.negate();
1351  addCmptAvBoundaryDiag(boundaryDiagCmpt);
1352 
1353  Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
1354  }
1355 
1356  Hphi.primitiveFieldRef() += lduMatrix::H(psi_.primitiveField()) + source_;
1357  addBoundarySource(Hphi.primitiveFieldRef());
1358 
1359  Hphi.primitiveFieldRef() /= psi_.mesh().V();
1361 
1362  typename Type::labelType validComponents
1363  (
1364  psi_.mesh().template validComponents<Type>()
1365  );
1366 
1367  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
1368  {
1369  if (validComponents[cmpt] == -1)
1370  {
1371  Hphi.replace
1372  (
1373  cmpt,
1375  );
1376  }
1377  }
1378 
1379  return tHphi;
1380 }
1381 
1382 
1383 template<class Type>
1385 {
1386  auto tH1 = volScalarField::New
1387  (
1388  "H(1)",
1389  psi_.mesh(),
1390  dimensions_/(dimVol*psi_.dimensions()),
1392  );
1393  auto& H1_ = tH1.ref();
1394 
1395  H1_.primitiveFieldRef() = lduMatrix::H1();
1396 
1397  forAll(psi_.boundaryField(), patchi)
1398  {
1399  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1400 
1401  if (ptf.coupled() && ptf.size())
1402  {
1403  addToInternalField
1404  (
1405  lduAddr().patchAddr(patchi),
1406  boundaryCoeffs_[patchi].component(0),
1407  H1_
1408  );
1409  }
1410  }
1411 
1412  H1_.primitiveFieldRef() /= psi_.mesh().V();
1413  H1_.correctBoundaryConditions();
1414 
1415  return tH1;
1416 }
1417 
1418 
1419 template<class Type>
1422 flux() const
1423 {
1424  if (!psi_.mesh().fluxRequired(psi_.name()))
1425  {
1427  << "flux requested but " << psi_.name()
1428  << " not specified in the fluxRequired sub-dictionary"
1429  " of fvSchemes."
1430  << abort(FatalError);
1431  }
1432 
1433  if (nMatrices() > 1)
1434  {
1436  << "Flux requested but " << psi_.name()
1437  << " can't handle multiple fvMatrix."
1438  << abort(FatalError);
1439  }
1440 
1442  (
1443  "flux(" + psi_.name() + ')',
1444  psi_.mesh(),
1445  dimensions(),
1446  lduMatrix::faceH(psi_.primitiveField())
1447  );
1448  auto& fieldFlux = tfieldFlux.ref();
1449  fieldFlux.setOriented();
1450 
1451 
1452  FieldField<Field, Type> InternalContrib = internalCoeffs_;
1453 
1454  label fieldi = 0;
1455  if (!useImplicit_)
1456  {
1457  forAll(InternalContrib, patchi)
1458  {
1459  InternalContrib[patchi] =
1460  cmptMultiply
1461  (
1462  InternalContrib[patchi],
1463  psi_.boundaryField()[patchi].patchInternalField()
1464  );
1465  }
1466  }
1467  else
1468  {
1469  FieldField<Field, Type> fluxInternalContrib(internalCoeffs_);
1470 
1471  mapContributions(fieldi, fluxInternalContrib, InternalContrib, true);
1472  }
1473 
1474  FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
1475 
1476  if (!useImplicit_)
1477  {
1478  forAll(NeighbourContrib, patchi)
1479  {
1480  if (psi_.boundaryField()[patchi].coupled())
1481  {
1482  NeighbourContrib[patchi] =
1483  cmptMultiply
1484  (
1485  NeighbourContrib[patchi],
1486  psi_.boundaryField()[patchi].patchNeighbourField()
1487  );
1488  }
1489  }
1490  }
1491  else
1492  {
1493  FieldField<Field, Type> fluxBoundaryContrib(boundaryCoeffs_);
1494 
1495  mapContributions(fieldi, fluxBoundaryContrib, NeighbourContrib, false);
1496  }
1497 
1498 
1499  {
1500  auto& ffbf = fieldFlux.boundaryFieldRef();
1501 
1502  forAll(ffbf, patchi)
1503  {
1504  ffbf[patchi] = InternalContrib[patchi] - NeighbourContrib[patchi];
1505  //DebugVar(gSum(ffbf[patchi]))
1506  }
1507  }
1508 
1509  if (faceFluxCorrectionPtr_)
1510  {
1511  fieldFlux += *faceFluxCorrectionPtr_;
1512  }
1514  return tfieldFlux;
1515 }
1516 
1517 
1518 template<class Type>
1520 (
1521  const word& name
1522 ) const
1523 {
1524  return psi_.mesh().solverDict(name);
1525 }
1526 
1527 
1528 template<class Type>
1530 {
1531  return psi_.mesh().solverDict
1532  (
1533  psi_.select(psi_.mesh().data().isFinalIteration())
1534  );
1535 }
1536 
1537 
1538 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1539 
1540 template<class Type>
1542 {
1543  if (this == &fvmv)
1544  {
1545  return; // Self-assignment is a no-op
1546  }
1547 
1548  if (&psi_ != &(fvmv.psi_))
1549  {
1551  << "different fields"
1552  << abort(FatalError);
1553  }
1554 
1555  dimensions_ = fvmv.dimensions_;
1556  lduMatrix::operator=(fvmv);
1557  source_ = fvmv.source_;
1558  internalCoeffs_ = fvmv.internalCoeffs_;
1559  boundaryCoeffs_ = fvmv.boundaryCoeffs_;
1560 
1561  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1562  {
1563  *faceFluxCorrectionPtr_ = *fvmv.faceFluxCorrectionPtr_;
1564  }
1565  else if (fvmv.faceFluxCorrectionPtr_)
1566  {
1567  faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
1568  (
1569  *fvmv.faceFluxCorrectionPtr_
1570  );
1571  }
1573  useImplicit_ = fvmv.useImplicit_;
1574  lduAssemblyName_ = fvmv.lduAssemblyName_;
1575 }
1576 
1577 
1578 template<class Type>
1581  operator=(tfvmv());
1582  tfvmv.clear();
1583 }
1584 
1585 
1586 template<class Type>
1588 {
1590  source_.negate();
1591  internalCoeffs_.negate();
1592  boundaryCoeffs_.negate();
1593 
1594  if (faceFluxCorrectionPtr_)
1595  {
1596  faceFluxCorrectionPtr_->negate();
1597  }
1598 }
1599 
1600 
1601 template<class Type>
1602 void Foam::fvMatrix<Type>::operator+=(const fvMatrix<Type>& fvmv)
1603 {
1604  checkMethod(*this, fvmv, "+=");
1605 
1606  dimensions_ += fvmv.dimensions_;
1607  lduMatrix::operator+=(fvmv);
1608  source_ += fvmv.source_;
1609  internalCoeffs_ += fvmv.internalCoeffs_;
1610  boundaryCoeffs_ += fvmv.boundaryCoeffs_;
1611 
1612  useImplicit_ = fvmv.useImplicit_;
1613  lduAssemblyName_ = fvmv.lduAssemblyName_;
1614  nMatrix_ = fvmv.nMatrix_;
1615 
1616  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1617  {
1618  *faceFluxCorrectionPtr_ += *fvmv.faceFluxCorrectionPtr_;
1619  }
1620  else if (fvmv.faceFluxCorrectionPtr_)
1621  {
1622  faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
1623  (
1624  *fvmv.faceFluxCorrectionPtr_
1625  );
1626  }
1627 }
1628 
1629 
1630 template<class Type>
1631 void Foam::fvMatrix<Type>::operator+=(const tmp<fvMatrix<Type>>& tfvmv)
1633  operator+=(tfvmv());
1634  tfvmv.clear();
1635 }
1636 
1637 
1638 template<class Type>
1640 {
1641  checkMethod(*this, fvmv, "-=");
1642 
1643  dimensions_ -= fvmv.dimensions_;
1644  lduMatrix::operator-=(fvmv);
1645  source_ -= fvmv.source_;
1646  internalCoeffs_ -= fvmv.internalCoeffs_;
1647  boundaryCoeffs_ -= fvmv.boundaryCoeffs_;
1648 
1649  useImplicit_ = fvmv.useImplicit_;
1650  lduAssemblyName_ = fvmv.lduAssemblyName_;
1651  nMatrix_ = fvmv.nMatrix_;
1652 
1653  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1654  {
1655  *faceFluxCorrectionPtr_ -= *fvmv.faceFluxCorrectionPtr_;
1656  }
1657  else if (fvmv.faceFluxCorrectionPtr_)
1658  {
1659  faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
1660  (
1661  -*fvmv.faceFluxCorrectionPtr_
1662  );
1663  }
1664 }
1665 
1666 
1667 template<class Type>
1668 void Foam::fvMatrix<Type>::operator-=(const tmp<fvMatrix<Type>>& tfvmv)
1669 {
1670  operator-=(tfvmv());
1671  tfvmv.clear();
1672 }
1673 
1674 
1675 template<class Type>
1677 (
1679 )
1680 {
1681  checkMethod(*this, su, "+=");
1682  source() -= su.mesh().V()*su.field();
1683 }
1684 
1685 
1686 template<class Type>
1688 (
1690 )
1691 {
1692  operator+=(tsu());
1693  tsu.clear();
1694 }
1695 
1696 
1697 template<class Type>
1699 (
1701 )
1702 {
1703  operator+=(tsu());
1704  tsu.clear();
1705 }
1706 
1707 
1708 template<class Type>
1710 (
1712 )
1713 {
1714  checkMethod(*this, su, "-=");
1715  source() += su.mesh().V()*su.field();
1716 }
1717 
1718 
1719 template<class Type>
1721 (
1723 )
1724 {
1725  operator-=(tsu());
1726  tsu.clear();
1727 }
1728 
1729 
1730 template<class Type>
1732 (
1734 )
1735 {
1736  operator-=(tsu());
1737  tsu.clear();
1738 }
1739 
1740 
1741 template<class Type>
1743 (
1744  const dimensioned<Type>& su
1745 )
1747  source() -= psi().mesh().V()*su;
1748 }
1749 
1750 
1751 template<class Type>
1753 (
1754  const dimensioned<Type>& su
1755 )
1757  source() += psi().mesh().V()*su;
1758 }
1759 
1760 
1761 template<class Type>
1763 (
1764  const volScalarField::Internal& dsf
1765 )
1766 {
1767  dimensions_ *= dsf.dimensions();
1768  lduMatrix::operator*=(dsf.field());
1769  source_ *= dsf.field();
1770 
1771  forAll(boundaryCoeffs_, patchi)
1772  {
1773  scalarField pisf
1774  (
1775  dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
1776  );
1777 
1778  internalCoeffs_[patchi] *= pisf;
1779  boundaryCoeffs_[patchi] *= pisf;
1780  }
1781 
1782  if (faceFluxCorrectionPtr_)
1783  {
1785  << "cannot scale a matrix containing a faceFluxCorrection"
1787  }
1788 }
1789 
1790 
1791 template<class Type>
1793 (
1794  const tmp<volScalarField::Internal>& tfld
1795 )
1796 {
1797  operator*=(tfld());
1798  tfld.clear();
1799 }
1800 
1801 
1802 template<class Type>
1804 (
1805  const tmp<volScalarField>& tfld
1806 )
1807 {
1808  operator*=(tfld());
1809  tfld.clear();
1810 }
1811 
1812 
1813 template<class Type>
1815 (
1816  const dimensioned<scalar>& ds
1817 )
1818 {
1819  dimensions_ *= ds.dimensions();
1820  lduMatrix::operator*=(ds.value());
1821  source_ *= ds.value();
1822  internalCoeffs_ *= ds.value();
1823  boundaryCoeffs_ *= ds.value();
1824 
1825  if (faceFluxCorrectionPtr_)
1826  {
1827  *faceFluxCorrectionPtr_ *= ds.value();
1828  }
1830 
1831 
1832 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
1833 
1834 template<class Type>
1835 void Foam::checkMethod
1836 (
1837  const fvMatrix<Type>& mat1,
1838  const fvMatrix<Type>& mat2,
1839  const char* op
1840 )
1841 {
1842  if (&mat1.psi() != &mat2.psi())
1843  {
1845  << "Incompatible fields for operation\n "
1846  << "[" << mat1.psi().name() << "] "
1847  << op
1848  << " [" << mat2.psi().name() << "]"
1849  << abort(FatalError);
1850  }
1851 
1852  if
1853  (
1855  && mat1.dimensions() != mat2.dimensions()
1856  )
1857  {
1859  << "Incompatible dimensions for operation\n "
1860  << "[" << mat1.psi().name() << mat1.dimensions()/dimVolume << " ] "
1861  << op
1862  << " [" << mat2.psi().name() << mat2.dimensions()/dimVolume << " ]"
1864  }
1865 }
1866 
1867 
1868 template<class Type>
1869 void Foam::checkMethod
1870 (
1871  const fvMatrix<Type>& mat,
1872  const DimensionedField<Type, volMesh>& fld,
1873  const char* op
1874 )
1875 {
1876  if
1877  (
1879  && mat.dimensions()/dimVolume != fld.dimensions()
1880  )
1881  {
1883  << "Incompatible dimensions for operation\n "
1884  << "[" << mat.psi().name() << mat.dimensions()/dimVolume << " ] "
1885  << op
1886  << " [" << fld.name() << fld.dimensions() << " ]"
1888  }
1889 }
1890 
1891 
1892 template<class Type>
1893 void Foam::checkMethod
1894 (
1895  const fvMatrix<Type>& mat,
1896  const dimensioned<Type>& dt,
1897  const char* op
1898 )
1899 {
1900  if
1901  (
1903  && mat.dimensions()/dimVolume != dt.dimensions()
1904  )
1905  {
1907  << "Incompatible dimensions for operation\n "
1908  << "[" << mat.psi().name() << mat.dimensions()/dimVolume << " ] "
1909  << op
1910  << " [" << dt.name() << dt.dimensions() << " ]"
1911  << abort(FatalError);
1912  }
1913 }
1914 
1915 
1916 template<class Type>
1918 (
1919  fvMatrix<Type>& mat,
1920  const dictionary& solverControls
1921 )
1922 {
1923  return mat.solve(solverControls);
1924 }
1925 
1926 template<class Type>
1928 (
1929  const tmp<fvMatrix<Type>>& tmat,
1930  const dictionary& solverControls
1931 )
1932 {
1933  SolverPerformance<Type> solverPerf(tmat.constCast().solve(solverControls));
1934 
1935  tmat.clear();
1936 
1937  return solverPerf;
1938 }
1939 
1940 
1941 template<class Type>
1943 (
1944  fvMatrix<Type>& mat,
1945  const word& name
1946 )
1947 {
1948  return mat.solve(name);
1949 }
1950 
1951 template<class Type>
1953 (
1954  const tmp<fvMatrix<Type>>& tmat,
1955  const word& name
1956 )
1957 {
1958  SolverPerformance<Type> solverPerf(tmat.constCast().solve(name));
1959 
1960  tmat.clear();
1961 
1962  return solverPerf;
1963 }
1964 
1965 
1966 template<class Type>
1967 Foam::SolverPerformance<Type> Foam::solve(fvMatrix<Type>& mat)
1968 {
1969  return mat.solve();
1970 }
1971 
1972 template<class Type>
1973 Foam::SolverPerformance<Type> Foam::solve(const tmp<fvMatrix<Type>>& tmat)
1974 {
1975  SolverPerformance<Type> solverPerf(tmat.constCast().solve());
1976 
1977  tmat.clear();
1978 
1979  return solverPerf;
1980 }
1981 
1982 
1983 template<class Type>
1985 (
1986  const fvMatrix<Type>& A
1987 )
1988 {
1989  tmp<Foam::fvMatrix<Type>> tAcorr = A - (A & A.psi());
1990 
1991  // Delete the faceFluxCorrection from the correction matrix
1992  // as it does not have a clear meaning or purpose
1993  tAcorr.ref().faceFluxCorrectionPtr(nullptr);
1994 
1995  return tAcorr;
1996 }
1997 
1998 
1999 template<class Type>
2001 (
2002  const tmp<fvMatrix<Type>>& tA
2003 )
2004 {
2005  tmp<Foam::fvMatrix<Type>> tAcorr = tA - (tA() & tA().psi());
2006 
2007  // Delete the faceFluxCorrection from the correction matrix
2008  // as it does not have a clear meaning or purpose
2009  tAcorr.ref().faceFluxCorrectionPtr(nullptr);
2010 
2011  return tAcorr;
2012 }
2013 
2014 
2015 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
2016 
2017 template<class Type>
2018 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2019 (
2020  const fvMatrix<Type>& A,
2021  const fvMatrix<Type>& B
2022 )
2023 {
2024  checkMethod(A, B, "==");
2025  return (A - B);
2026 }
2027 
2028 template<class Type>
2029 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2030 (
2031  const tmp<fvMatrix<Type>>& tA,
2032  const fvMatrix<Type>& B
2033 )
2034 {
2035  checkMethod(tA(), B, "==");
2036  return (tA - B);
2037 }
2038 
2039 template<class Type>
2040 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2041 (
2042  const fvMatrix<Type>& A,
2043  const tmp<fvMatrix<Type>>& tB
2044 )
2045 {
2046  checkMethod(A, tB(), "==");
2047  return (A - tB);
2048 }
2049 
2050 template<class Type>
2051 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2052 (
2053  const tmp<fvMatrix<Type>>& tA,
2054  const tmp<fvMatrix<Type>>& tB
2055 )
2056 {
2057  checkMethod(tA(), tB(), "==");
2058  return (tA - tB);
2059 }
2060 
2061 template<class Type>
2062 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2063 (
2064  const fvMatrix<Type>& A,
2065  const DimensionedField<Type, volMesh>& su
2066 )
2067 {
2068  checkMethod(A, su, "==");
2069  auto tC = tmp<fvMatrix<Type>>::New(A);
2070  tC.ref().source() += su.mesh().V()*su.field();
2071  return tC;
2072 }
2073 
2074 template<class Type>
2075 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2076 (
2077  const fvMatrix<Type>& A,
2078  const tmp<DimensionedField<Type, volMesh>>& tsu
2079 )
2080 {
2081  checkMethod(A, tsu(), "==");
2082  auto tC = tmp<fvMatrix<Type>>::New(A);
2083  tC.ref().source() += tsu().mesh().V()*tsu().field();
2084  tsu.clear();
2085  return tC;
2086 }
2087 
2088 template<class Type>
2089 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2090 (
2091  const fvMatrix<Type>& A,
2092  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2093 )
2094 {
2095  checkMethod(A, tsu(), "==");
2096  auto tC = tmp<fvMatrix<Type>>::New(A);
2097  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2098  tsu.clear();
2099  return tC;
2100 }
2101 
2102 template<class Type>
2103 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2104 (
2105  const tmp<fvMatrix<Type>>& tA,
2106  const DimensionedField<Type, volMesh>& su
2107 )
2108 {
2109  checkMethod(tA(), su, "==");
2110  tmp<fvMatrix<Type>> tC(tA.ptr());
2111  tC.ref().source() += su.mesh().V()*su.field();
2112  return tC;
2113 }
2114 
2115 template<class Type>
2116 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2117 (
2118  const tmp<fvMatrix<Type>>& tA,
2119  const tmp<DimensionedField<Type, volMesh>>& tsu
2120 )
2121 {
2122  checkMethod(tA(), tsu(), "==");
2123  tmp<fvMatrix<Type>> tC(tA.ptr());
2124  tC.ref().source() += tsu().mesh().V()*tsu().field();
2125  tsu.clear();
2126  return tC;
2127 }
2128 
2129 template<class Type>
2130 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2131 (
2132  const tmp<fvMatrix<Type>>& tA,
2133  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2134 )
2135 {
2136  checkMethod(tA(), tsu(), "==");
2137  tmp<fvMatrix<Type>> tC(tA.ptr());
2138  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2139  tsu.clear();
2140  return tC;
2141 }
2142 
2143 template<class Type>
2144 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2145 (
2146  const fvMatrix<Type>& A,
2147  const dimensioned<Type>& su
2148 )
2149 {
2150  checkMethod(A, su, "==");
2151  auto tC = tmp<fvMatrix<Type>>::New(A);
2152  tC.ref().source() += A.psi().mesh().V()*su.value();
2153  return tC;
2154 }
2155 
2156 template<class Type>
2157 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2158 (
2159  const tmp<fvMatrix<Type>>& tA,
2160  const dimensioned<Type>& su
2161 )
2162 {
2163  checkMethod(tA(), su, "==");
2164  tmp<fvMatrix<Type>> tC(tA.ptr());
2165  tC.ref().source() += tC().psi().mesh().V()*su.value();
2166  return tC;
2167 }
2168 
2169 template<class Type>
2170 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2171 (
2172  const fvMatrix<Type>& A,
2173  const Foam::zero
2174 )
2175 {
2176  return A;
2177 }
2178 
2179 
2180 template<class Type>
2181 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2182 (
2183  const tmp<fvMatrix<Type>>& tA,
2184  const Foam::zero
2185 )
2186 {
2187  return tA;
2188 }
2189 
2190 
2191 template<class Type>
2192 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2193 (
2194  const fvMatrix<Type>& A
2195 )
2196 {
2197  auto tC = tmp<fvMatrix<Type>>::New(A);
2198  tC.ref().negate();
2199  return tC;
2200 }
2201 
2202 template<class Type>
2203 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2204 (
2205  const tmp<fvMatrix<Type>>& tA
2206 )
2207 {
2208  tmp<fvMatrix<Type>> tC(tA.ptr());
2209  tC.ref().negate();
2210  return tC;
2211 }
2212 
2213 
2214 template<class Type>
2215 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2216 (
2217  const fvMatrix<Type>& A,
2218  const fvMatrix<Type>& B
2219 )
2220 {
2221  checkMethod(A, B, "+");
2222  auto tC = tmp<fvMatrix<Type>>::New(A);
2223  tC.ref() += B;
2224  return tC;
2225 }
2226 
2227 template<class Type>
2228 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2229 (
2230  const tmp<fvMatrix<Type>>& tA,
2231  const fvMatrix<Type>& B
2232 )
2233 {
2234  checkMethod(tA(), B, "+");
2235  tmp<fvMatrix<Type>> tC(tA.ptr());
2236  tC.ref() += B;
2237  return tC;
2238 }
2239 
2240 template<class Type>
2241 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2242 (
2243  const fvMatrix<Type>& A,
2244  const tmp<fvMatrix<Type>>& tB
2245 )
2246 {
2247  checkMethod(A, tB(), "+");
2248  tmp<fvMatrix<Type>> tC(tB.ptr());
2249  tC.ref() += A;
2250  return tC;
2251 }
2252 
2253 template<class Type>
2254 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2255 (
2256  const tmp<fvMatrix<Type>>& tA,
2257  const tmp<fvMatrix<Type>>& tB
2258 )
2259 {
2260  checkMethod(tA(), tB(), "+");
2261  tmp<fvMatrix<Type>> tC(tA.ptr());
2262  tC.ref() += tB();
2263  tB.clear();
2264  return tC;
2265 }
2266 
2267 template<class Type>
2268 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2269 (
2270  const fvMatrix<Type>& A,
2271  const DimensionedField<Type, volMesh>& su
2272 )
2273 {
2274  checkMethod(A, su, "+");
2275  auto tC = tmp<fvMatrix<Type>>::New(A);
2276  tC.ref().source() -= su.mesh().V()*su.field();
2277  return tC;
2278 }
2279 
2280 template<class Type>
2281 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2282 (
2283  const fvMatrix<Type>& A,
2284  const tmp<DimensionedField<Type, volMesh>>& tsu
2285 )
2286 {
2287  checkMethod(A, tsu(), "+");
2288  auto tC = tmp<fvMatrix<Type>>::New(A);
2289  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2290  tsu.clear();
2291  return tC;
2292 }
2293 
2294 template<class Type>
2295 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2296 (
2297  const fvMatrix<Type>& A,
2298  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2299 )
2300 {
2301  checkMethod(A, tsu(), "+");
2302  auto tC = tmp<fvMatrix<Type>>::New(A);
2303  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2304  tsu.clear();
2305  return tC;
2306 }
2307 
2308 template<class Type>
2309 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2310 (
2311  const tmp<fvMatrix<Type>>& tA,
2312  const DimensionedField<Type, volMesh>& su
2313 )
2314 {
2315  checkMethod(tA(), su, "+");
2316  tmp<fvMatrix<Type>> tC(tA.ptr());
2317  tC.ref().source() -= su.mesh().V()*su.field();
2318  return tC;
2319 }
2320 
2321 template<class Type>
2322 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2323 (
2324  const tmp<fvMatrix<Type>>& tA,
2325  const tmp<DimensionedField<Type, volMesh>>& tsu
2326 )
2327 {
2328  checkMethod(tA(), tsu(), "+");
2329  tmp<fvMatrix<Type>> tC(tA.ptr());
2330  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2331  tsu.clear();
2332  return tC;
2333 }
2334 
2335 template<class Type>
2336 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2337 (
2338  const tmp<fvMatrix<Type>>& tA,
2339  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2340 )
2341 {
2342  checkMethod(tA(), tsu(), "+");
2343  tmp<fvMatrix<Type>> tC(tA.ptr());
2344  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2345  tsu.clear();
2346  return tC;
2347 }
2348 
2349 template<class Type>
2350 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2351 (
2352  const DimensionedField<Type, volMesh>& su,
2353  const fvMatrix<Type>& A
2354 )
2355 {
2356  checkMethod(A, su, "+");
2357  auto tC = tmp<fvMatrix<Type>>::New(A);
2358  tC.ref().source() -= su.mesh().V()*su.field();
2359  return tC;
2360 }
2361 
2362 template<class Type>
2363 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2364 (
2365  const tmp<DimensionedField<Type, volMesh>>& tsu,
2366  const fvMatrix<Type>& A
2367 )
2368 {
2369  checkMethod(A, tsu(), "+");
2370  auto tC = tmp<fvMatrix<Type>>::New(A);
2371  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2372  tsu.clear();
2373  return tC;
2374 }
2375 
2376 template<class Type>
2377 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2378 (
2379  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2380  const fvMatrix<Type>& A
2381 )
2382 {
2383  checkMethod(A, tsu(), "+");
2384  auto tC = tmp<fvMatrix<Type>>::New(A);
2385  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2386  tsu.clear();
2387  return tC;
2388 }
2389 
2390 template<class Type>
2391 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2392 (
2393  const DimensionedField<Type, volMesh>& su,
2394  const tmp<fvMatrix<Type>>& tA
2395 )
2396 {
2397  checkMethod(tA(), su, "+");
2398  tmp<fvMatrix<Type>> tC(tA.ptr());
2399  tC.ref().source() -= su.mesh().V()*su.field();
2400  return tC;
2401 }
2402 
2403 template<class Type>
2404 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2405 (
2406  const tmp<DimensionedField<Type, volMesh>>& tsu,
2407  const tmp<fvMatrix<Type>>& tA
2408 )
2409 {
2410  checkMethod(tA(), tsu(), "+");
2411  tmp<fvMatrix<Type>> tC(tA.ptr());
2412  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2413  tsu.clear();
2414  return tC;
2415 }
2416 
2417 template<class Type>
2418 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2419 (
2420  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2421  const tmp<fvMatrix<Type>>& tA
2422 )
2423 {
2424  checkMethod(tA(), tsu(), "+");
2425  tmp<fvMatrix<Type>> tC(tA.ptr());
2426  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2427  tsu.clear();
2428  return tC;
2429 }
2430 
2431 
2432 template<class Type>
2433 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2434 (
2435  const fvMatrix<Type>& A,
2436  const fvMatrix<Type>& B
2437 )
2438 {
2439  checkMethod(A, B, "-");
2440  auto tC = tmp<fvMatrix<Type>>::New(A);
2441  tC.ref() -= B;
2442  return tC;
2443 }
2444 
2445 template<class Type>
2446 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2447 (
2448  const tmp<fvMatrix<Type>>& tA,
2449  const fvMatrix<Type>& B
2450 )
2451 {
2452  checkMethod(tA(), B, "-");
2453  tmp<fvMatrix<Type>> tC(tA.ptr());
2454  tC.ref() -= B;
2455  return tC;
2456 }
2457 
2458 template<class Type>
2459 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2460 (
2461  const fvMatrix<Type>& A,
2462  const tmp<fvMatrix<Type>>& tB
2463 )
2464 {
2465  checkMethod(A, tB(), "-");
2466  tmp<fvMatrix<Type>> tC(tB.ptr());
2467  tC.ref() -= A;
2468  tC.ref().negate();
2469  return tC;
2470 }
2471 
2472 template<class Type>
2473 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2474 (
2475  const tmp<fvMatrix<Type>>& tA,
2476  const tmp<fvMatrix<Type>>& tB
2477 )
2478 {
2479  checkMethod(tA(), tB(), "-");
2480  tmp<fvMatrix<Type>> tC(tA.ptr());
2481  tC.ref() -= tB();
2482  tB.clear();
2483  return tC;
2484 }
2485 
2486 template<class Type>
2487 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2488 (
2489  const fvMatrix<Type>& A,
2490  const DimensionedField<Type, volMesh>& su
2491 )
2492 {
2493  checkMethod(A, su, "-");
2494  auto tC = tmp<fvMatrix<Type>>::New(A);
2495  tC.ref().source() += su.mesh().V()*su.field();
2496  return tC;
2497 }
2498 
2499 template<class Type>
2500 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2501 (
2502  const fvMatrix<Type>& A,
2503  const tmp<DimensionedField<Type, volMesh>>& tsu
2504 )
2505 {
2506  checkMethod(A, tsu(), "-");
2507  auto tC = tmp<fvMatrix<Type>>::New(A);
2508  tC.ref().source() += tsu().mesh().V()*tsu().field();
2509  tsu.clear();
2510  return tC;
2511 }
2512 
2513 template<class Type>
2514 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2515 (
2516  const fvMatrix<Type>& A,
2517  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2518 )
2519 {
2520  checkMethod(A, tsu(), "-");
2521  auto tC = tmp<fvMatrix<Type>>::New(A);
2522  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2523  tsu.clear();
2524  return tC;
2525 }
2526 
2527 template<class Type>
2528 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2529 (
2530  const tmp<fvMatrix<Type>>& tA,
2531  const DimensionedField<Type, volMesh>& su
2532 )
2533 {
2534  checkMethod(tA(), su, "-");
2535  tmp<fvMatrix<Type>> tC(tA.ptr());
2536  tC.ref().source() += su.mesh().V()*su.field();
2537  return tC;
2538 }
2539 
2540 template<class Type>
2541 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2542 (
2543  const tmp<fvMatrix<Type>>& tA,
2544  const tmp<DimensionedField<Type, volMesh>>& tsu
2545 )
2546 {
2547  checkMethod(tA(), tsu(), "-");
2548  tmp<fvMatrix<Type>> tC(tA.ptr());
2549  tC.ref().source() += tsu().mesh().V()*tsu().field();
2550  tsu.clear();
2551  return tC;
2552 }
2553 
2554 template<class Type>
2555 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2556 (
2557  const tmp<fvMatrix<Type>>& tA,
2558  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2559 )
2560 {
2561  checkMethod(tA(), tsu(), "-");
2562  tmp<fvMatrix<Type>> tC(tA.ptr());
2563  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2564  tsu.clear();
2565  return tC;
2566 }
2567 
2568 template<class Type>
2569 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2570 (
2571  const DimensionedField<Type, volMesh>& su,
2572  const fvMatrix<Type>& A
2573 )
2574 {
2575  checkMethod(A, su, "-");
2576  auto tC = tmp<fvMatrix<Type>>::New(A);
2577  tC.ref().negate();
2578  tC.ref().source() -= su.mesh().V()*su.field();
2579  return tC;
2580 }
2581 
2582 template<class Type>
2583 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2584 (
2585  const tmp<DimensionedField<Type, volMesh>>& tsu,
2586  const fvMatrix<Type>& A
2587 )
2588 {
2589  checkMethod(A, tsu(), "-");
2590  auto tC = tmp<fvMatrix<Type>>::New(A);
2591  tC.ref().negate();
2592  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2593  tsu.clear();
2594  return tC;
2595 }
2596 
2597 template<class Type>
2598 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2599 (
2600  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2601  const fvMatrix<Type>& A
2602 )
2603 {
2604  checkMethod(A, tsu(), "-");
2605  auto tC = tmp<fvMatrix<Type>>::New(A);
2606  tC.ref().negate();
2607  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2608  tsu.clear();
2609  return tC;
2610 }
2611 
2612 template<class Type>
2613 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2614 (
2615  const DimensionedField<Type, volMesh>& su,
2616  const tmp<fvMatrix<Type>>& tA
2617 )
2618 {
2619  checkMethod(tA(), su, "-");
2620  tmp<fvMatrix<Type>> tC(tA.ptr());
2621  tC.ref().negate();
2622  tC.ref().source() -= su.mesh().V()*su.field();
2623  return tC;
2624 }
2625 
2626 template<class Type>
2627 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2628 (
2629  const tmp<DimensionedField<Type, volMesh>>& tsu,
2630  const tmp<fvMatrix<Type>>& tA
2631 )
2632 {
2633  checkMethod(tA(), tsu(), "-");
2634  tmp<fvMatrix<Type>> tC(tA.ptr());
2635  tC.ref().negate();
2636  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2637  tsu.clear();
2638  return tC;
2639 }
2640 
2641 template<class Type>
2642 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2643 (
2644  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2645  const tmp<fvMatrix<Type>>& tA
2646 )
2647 {
2648  checkMethod(tA(), tsu(), "-");
2649  tmp<fvMatrix<Type>> tC(tA.ptr());
2650  tC.ref().negate();
2651  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2652  tsu.clear();
2653  return tC;
2654 }
2655 
2656 template<class Type>
2657 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2658 (
2659  const fvMatrix<Type>& A,
2660  const dimensioned<Type>& su
2661 )
2662 {
2663  checkMethod(A, su, "+");
2664  auto tC = tmp<fvMatrix<Type>>::New(A);
2665  tC.ref().source() -= su.value()*A.psi().mesh().V();
2666  return tC;
2667 }
2668 
2669 template<class Type>
2670 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2671 (
2672  const tmp<fvMatrix<Type>>& tA,
2673  const dimensioned<Type>& su
2674 )
2675 {
2676  checkMethod(tA(), su, "+");
2677  tmp<fvMatrix<Type>> tC(tA.ptr());
2678  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2679  return tC;
2680 }
2681 
2682 template<class Type>
2683 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2684 (
2685  const dimensioned<Type>& su,
2686  const fvMatrix<Type>& A
2687 )
2688 {
2689  checkMethod(A, su, "+");
2690  auto tC = tmp<fvMatrix<Type>>::New(A);
2691  tC.ref().source() -= su.value()*A.psi().mesh().V();
2692  return tC;
2693 }
2694 
2695 template<class Type>
2696 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2697 (
2698  const dimensioned<Type>& su,
2699  const tmp<fvMatrix<Type>>& tA
2700 )
2701 {
2702  checkMethod(tA(), su, "+");
2703  tmp<fvMatrix<Type>> tC(tA.ptr());
2704  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2705  return tC;
2706 }
2707 
2708 template<class Type>
2709 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2710 (
2711  const fvMatrix<Type>& A,
2712  const dimensioned<Type>& su
2713 )
2714 {
2715  checkMethod(A, su, "-");
2716  auto tC = tmp<fvMatrix<Type>>::New(A);
2717  tC.ref().source() += su.value()*tC().psi().mesh().V();
2718  return tC;
2719 }
2720 
2721 template<class Type>
2722 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2723 (
2724  const tmp<fvMatrix<Type>>& tA,
2725  const dimensioned<Type>& su
2726 )
2727 {
2728  checkMethod(tA(), su, "-");
2729  tmp<fvMatrix<Type>> tC(tA.ptr());
2730  tC.ref().source() += su.value()*tC().psi().mesh().V();
2731  return tC;
2732 }
2733 
2734 template<class Type>
2735 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2736 (
2737  const dimensioned<Type>& su,
2738  const fvMatrix<Type>& A
2739 )
2740 {
2741  checkMethod(A, su, "-");
2742  auto tC = tmp<fvMatrix<Type>>::New(A);
2743  tC.ref().negate();
2744  tC.ref().source() -= su.value()*A.psi().mesh().V();
2745  return tC;
2746 }
2747 
2748 template<class Type>
2749 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2750 (
2751  const dimensioned<Type>& su,
2752  const tmp<fvMatrix<Type>>& tA
2753 )
2754 {
2755  checkMethod(tA(), su, "-");
2756  tmp<fvMatrix<Type>> tC(tA.ptr());
2757  tC.ref().negate();
2758  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2759  return tC;
2760 }
2761 
2762 
2763 template<class Type>
2764 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2765 (
2766  const volScalarField::Internal& dsf,
2767  const fvMatrix<Type>& A
2768 )
2769 {
2770  auto tC = tmp<fvMatrix<Type>>::New(A);
2771  tC.ref() *= dsf;
2772  return tC;
2773 }
2774 
2775 template<class Type>
2776 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2777 (
2778  const tmp<volScalarField::Internal>& tdsf,
2779  const fvMatrix<Type>& A
2780 )
2781 {
2782  auto tC = tmp<fvMatrix<Type>>::New(A);
2783  tC.ref() *= tdsf;
2784  return tC;
2785 }
2786 
2787 template<class Type>
2788 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2789 (
2790  const tmp<volScalarField>& tvsf,
2791  const fvMatrix<Type>& A
2792 )
2793 {
2794  auto tC = tmp<fvMatrix<Type>>::New(A);
2795  tC.ref() *= tvsf;
2796  return tC;
2797 }
2798 
2799 template<class Type>
2800 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2801 (
2802  const volScalarField::Internal& dsf,
2803  const tmp<fvMatrix<Type>>& tA
2804 )
2805 {
2806  tmp<fvMatrix<Type>> tC(tA.ptr());
2807  tC.ref() *= dsf;
2808  return tC;
2809 }
2810 
2811 template<class Type>
2812 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2813 (
2814  const tmp<volScalarField::Internal>& tdsf,
2815  const tmp<fvMatrix<Type>>& tA
2816 )
2817 {
2818  tmp<fvMatrix<Type>> tC(tA.ptr());
2819  tC.ref() *= tdsf;
2820  return tC;
2821 }
2822 
2823 template<class Type>
2824 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2825 (
2826  const tmp<volScalarField>& tvsf,
2827  const tmp<fvMatrix<Type>>& tA
2828 )
2829 {
2830  tmp<fvMatrix<Type>> tC(tA.ptr());
2831  tC.ref() *= tvsf;
2832  return tC;
2833 }
2834 
2835 template<class Type>
2836 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2837 (
2838  const dimensioned<scalar>& ds,
2839  const fvMatrix<Type>& A
2840 )
2841 {
2842  auto tC = tmp<fvMatrix<Type>>::New(A);
2843  tC.ref() *= ds;
2844  return tC;
2845 }
2846 
2847 template<class Type>
2848 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2849 (
2850  const dimensioned<scalar>& ds,
2851  const tmp<fvMatrix<Type>>& tA
2852 )
2853 {
2854  tmp<fvMatrix<Type>> tC(tA.ptr());
2855  tC.ref() *= ds;
2856  return tC;
2857 }
2858 
2859 
2860 template<class Type>
2862 Foam::operator&
2863 (
2864  const fvMatrix<Type>& M,
2865  const DimensionedField<Type, volMesh>& psi
2866 )
2867 {
2869  (
2870  "M&" + psi.name(),
2871  psi.mesh(),
2872  M.dimensions()/dimVol,
2874  );
2875  auto& Mphi = tMphi.ref();
2876 
2877  // Loop over field components
2878  if (M.hasDiag())
2879  {
2880  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2881  {
2882  scalarField psiCmpt(psi.field().component(cmpt));
2883  scalarField boundaryDiagCmpt(M.diag());
2884  M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2885  Mphi.primitiveFieldRef().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2886  }
2887  }
2888  else
2889  {
2890  Mphi.primitiveFieldRef() = Zero;
2891  }
2892 
2893  Mphi.primitiveFieldRef() += M.lduMatrix::H(psi.field()) + M.source();
2894  M.addBoundarySource(Mphi.primitiveFieldRef());
2895 
2896  Mphi.primitiveFieldRef() /= -psi.mesh().V();
2897  Mphi.correctBoundaryConditions();
2898 
2899  return tMphi;
2900 }
2901 
2902 template<class Type>
2904 Foam::operator&
2905 (
2906  const fvMatrix<Type>& M,
2907  const tmp<DimensionedField<Type, volMesh>>& tpsi
2908 )
2909 {
2910  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = M & tpsi();
2911  tpsi.clear();
2912  return tMpsi;
2913 }
2914 
2915 template<class Type>
2917 Foam::operator&
2918 (
2919  const fvMatrix<Type>& M,
2920  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tpsi
2921 )
2922 {
2923  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = M & tpsi();
2924  tpsi.clear();
2925  return tMpsi;
2926 }
2927 
2928 template<class Type>
2930 Foam::operator&
2931 (
2932  const tmp<fvMatrix<Type>>& tM,
2933  const DimensionedField<Type, volMesh>& psi
2934 )
2935 {
2936  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & psi;
2937  tM.clear();
2938  return tMpsi;
2939 }
2940 
2941 template<class Type>
2943 Foam::operator&
2944 (
2945  const tmp<fvMatrix<Type>>& tM,
2946  const tmp<DimensionedField<Type, volMesh>>& tpsi
2947 )
2948 {
2949  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2950  tM.clear();
2951  tpsi.clear();
2952  return tMpsi;
2953 }
2954 
2955 template<class Type>
2957 Foam::operator&
2958 (
2959  const tmp<fvMatrix<Type>>& tM,
2960  const tmp<GeometricField<Type, fvPatchField, volMesh>>& tpsi
2961 )
2962 {
2963  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2964  tM.clear();
2965  tpsi.clear();
2966  return tMpsi;
2967 }
2968 
2969 
2970 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2971 
2972 template<class Type>
2973 Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
2974 {
2975  os << static_cast<const lduMatrix&>(fvm) << nl
2976  << fvm.dimensions_ << nl
2977  << fvm.source_ << nl
2978  << fvm.internalCoeffs_ << nl
2979  << fvm.boundaryCoeffs_ << endl;
2980 
2982 
2983  return os;
2984 }
2986 
2987 // * * * * * * * * * * * * * Expression Templates * * * * * * * * * * * * * //
2988 
2989 template<class Type>
2990 template<typename E>
2992 (
2995  <
2996  E,
2997  typename E::DiagExpr,
2998  typename E::UpperExpr,
2999  typename E::LowerExpr,
3000  typename E::FaceFluxExpr,
3001  typename E::SourceExpr
3002  >& expr
3003 )
3004 :
3005  lduMatrix(psi.mesh()),
3006  psi_(psi),
3007  useImplicit_(false),
3008  lduAssemblyName_(),
3009  nMatrix_(0),
3010  dimensions_(expr.dimensions()),
3011  source_(psi.size(), Zero),
3012  internalCoeffs_(psi.mesh().boundary().size()),
3013  boundaryCoeffs_(psi.mesh().boundary().size())
3014 {
3016  << "Constructing fvMatrix<Type> from expression for field "
3017  << psi_.name() << endl;
3018 
3019  checkImplicit();
3020 
3021  // Fill in diag,upper,lower etc.
3022  expr.evaluate(*this);
3023 
3024  auto& psiRef = this->psi(0);
3025  const label currentStatePsi = psiRef.eventNo();
3026  psiRef.boundaryFieldRef().updateCoeffs();
3027  psiRef.eventNo() = currentStatePsi;
3028 }
3029 
3030 
3031 template<class Type>
3034 {
3036 }
3037 
3038 
3039 template<class Type>
3040 template<typename E>
3042 (
3044  <
3045  E,
3046  typename E::DiagExpr,
3047  typename E::UpperExpr,
3048  typename E::LowerExpr,
3049  typename E::FaceFluxExpr,
3050  typename E::SourceExpr
3051  >& expr
3052 )
3053 {
3054  expr.evaluate(*this);
3055 }
3056 
3057 
3058 // * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
3059 
3060 #include "fvMatrixSolve.C"
3061 
3062 // ************************************************************************* //
void mapContributions(label fieldi, const FieldField< Field, Type > &fluxContrib, FieldField< Field, Type > &contrib, bool internal) const
Add internal and boundary contribution to local patches.
Definition: fvMatrix.C:549
Foam::surfaceFields.
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:41
faceListList boundary
void update(UPtrList< GeometricField< Type, fvPatchField, volMesh >> &psis)
Update mappings.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
tmp< Field< Type > > faceH(const Field< Type > &) const
virtual bool coupled() const
True if the patch field is coupled.
Definition: fvPatchField.H:242
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
void setValuesFromList(const labelUList &cellLabels, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:219
uint8_t direction
Definition: direction.H:46
label faceId(-1)
T & last()
Return reference to the last element of the list.
Definition: UPtrList.H:863
Expression wrap of const reference to fvMatrix.
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:116
void createOrUpdateLduPrimitiveAssembly()
Create or update ldu assembly.
Definition: fvMatrix.C:907
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UEqn relax()
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1242
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:130
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:600
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
void negate()
Inplace negate.
Definition: fvMatrix.C:1580
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:205
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:144
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1004
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition: fvMatrix.C:791
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:519
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1325
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:61
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1632
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:65
bool checkImplicit(const label fieldi=0)
Name the implicit assembly addressing.
Definition: fvMatrix.C:314
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1926
Generic GeometricField class.
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Generic dimensioned Type class.
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1377
Ignore writing from objectRegistry::writeObject()
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:79
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
tmp< scalarField > H1() const
conserve primitiveFieldRef()+
UList< label > labelUList
A UList of labels.
Definition: UList.H:75
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition: fvMatrix.C:1273
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1534
void reduce(T &value, [[maybe_unused]] BinaryOp bop, [[maybe_unused]] const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:299
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:487
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
void addFvMatrix(fvMatrix< Type > &matrix)
Add fvMatrix.
Definition: fvMatrix.C:1065
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition: Field.H:69
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition: fvMatrix.C:471
ListConstRefWrap< typename Expr::value_type > expr
Fully self-contained constant field wrapper. Not needed?
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
void setReferences(const labelUList &cellLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1020
#define DebugInFunction
Report an information message using Foam::Info.
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1026
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
const auto & io
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:458
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:215
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
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
void transferFvMatrixCoeffs()
Transfer lower, upper, diag and source to this fvMatrix.
Definition: fvMatrix.C:814
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void operator=(const lduMatrix &)
Copy assignment.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
This boundary condition enforces a cyclic condition between a pair of boundaries. ...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:607
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition: fvMatrix.C:351
void operator*=(const scalarField &)
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
Expression::fvMatrixConstRefWrap< fvMatrix< Type > > expr() const
Wrap value as expression.
Definition: fvMatrix.C:3026
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
const dimensionSet & dimensions() const noexcept
Return dimensions.
A dynamically resizable PtrList with allocation management.
Definition: PtrDynList.H:48
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition: fvMatrix.C:881
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1595
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:1282
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:80
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1306
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
bool set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:469
static bool checking() noexcept
True if dimension checking is enabled (the usual default)
Definition: dimensionSet.H:241
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const DynamicField< Type > & field() const noexcept
Return const-reference to the primitive field values.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:52
label cellId
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Nothing to be read.
const Mesh & mesh() const noexcept
Return const reference to mesh.
A single value that is represented as a list with an operator[] to access the value. This can be useful for templated operations expecting a list accessor.
Definition: UniformList.H:48
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:169
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:289
const volScalarField & psi
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1415
void setValues(const labelUList &cellLabels, const Type &value)
Set solution in given cells to the specified value and eliminate the corresponding equations from the...
Definition: fvMatrix.C:971
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1260
tmp< Field< Type > > H(const Field< Type > &) const
const dimensionedScalar & D
void operator-=(const lduMatrix &)
void operator+=(const lduMatrix &)
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:256
List< label > labelList
A List of labels.
Definition: List.H:62
void clearOut()
Clear additional addressing.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:842
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:194
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:188
Request registration (bool: true)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
void setBounAndInterCoeffs()
Manipulate boundary/internal coeffs for coupling.
Definition: fvMatrix.C:661
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
void setSize(label n)
Alias for resize()
Definition: List.H:536
#define M(I)
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const dictionary & solverDict() const
Return the solver dictionary for psi, taking into account finalIteration.
Definition: fvMatrix.C:1522
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
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
void append(autoPtr< T > &ptr)
Move append an element to the end of the list.
Definition: PtrDynList.H:398
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.