fvMatrix.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "volFields.H"
27 #include "surfaceFields.H"
30 #include "coupledFvPatchFields.H"
31 #include "UIndirectList.H"
32 
33 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
36 template<class Type2>
38 (
39  const labelUList& addr,
40  const Field<Type2>& pf,
41  Field<Type2>& intf
42 ) const
43 {
44  if (addr.size() != pf.size())
45  {
47  << "sizes of addressing and field are different"
48  << abort(FatalError);
49  }
50 
51  forAll(addr, facei)
52  {
53  intf[addr[facei]] += pf[facei];
54  }
55 }
56 
57 
58 template<class Type>
59 template<class Type2>
61 (
62  const labelUList& addr,
63  const tmp<Field<Type2>>& tpf,
64  Field<Type2>& intf
65 ) const
66 {
67  addToInternalField(addr, tpf(), intf);
68  tpf.clear();
69 }
70 
71 
72 template<class Type>
73 template<class Type2>
75 (
76  const labelUList& addr,
77  const Field<Type2>& pf,
78  Field<Type2>& intf
79 ) const
80 {
81  if (addr.size() != pf.size())
82  {
84  << "sizes of addressing and field are different"
85  << abort(FatalError);
86  }
87 
88  forAll(addr, facei)
89  {
90  intf[addr[facei]] -= pf[facei];
91  }
92 }
93 
94 
95 template<class Type>
96 template<class Type2>
98 (
99  const labelUList& addr,
100  const tmp<Field<Type2>>& tpf,
101  Field<Type2>& intf
102 ) const
103 {
104  subtractFromInternalField(addr, tpf(), intf);
105  tpf.clear();
106 }
107 
108 
109 template<class Type>
111 (
112  scalarField& diag,
113  const direction solveCmpt
114 ) const
115 {
116  forAll(internalCoeffs_, patchi)
117  {
118  addToInternalField
119  (
120  lduAddr().patchAddr(patchi),
121  internalCoeffs_[patchi].component(solveCmpt),
122  diag
123  );
124  }
125 }
126 
127 
128 template<class Type>
130 {
131  forAll(internalCoeffs_, patchi)
132  {
133  addToInternalField
134  (
135  lduAddr().patchAddr(patchi),
136  cmptAv(internalCoeffs_[patchi]),
137  diag
138  );
139  }
140 }
141 
142 
143 template<class Type>
145 (
146  Field<Type>& source,
147  const bool couples
148 ) const
149 {
150  forAll(psi_.boundaryField(), patchi)
151  {
152  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
153  const Field<Type>& pbc = boundaryCoeffs_[patchi];
154 
155  if (!ptf.coupled())
156  {
157  addToInternalField(lduAddr().patchAddr(patchi), pbc, source);
158  }
159  else if (couples)
160  {
161  const tmp<Field<Type>> tpnf = ptf.patchNeighbourField();
162  const Field<Type>& pnf = tpnf();
163 
164  const labelUList& addr = lduAddr().patchAddr(patchi);
165 
166  forAll(addr, facei)
167  {
168  source[addr[facei]] += cmptMultiply(pbc[facei], pnf[facei]);
169  }
170  }
171  }
172 }
173 
174 
175 template<class Type>
176 template<template<class> class ListType>
178 (
179  const labelUList& cellLabels,
180  const ListType<Type>& values
181 )
182 {
183  const fvMesh& mesh = psi_.mesh();
184 
185  const cellList& cells = mesh.cells();
186  const labelUList& own = mesh.owner();
187  const labelUList& nei = mesh.neighbour();
188 
189  scalarField& Diag = diag();
190  Field<Type>& psi =
191  const_cast
192  <
194  >(psi_).primitiveFieldRef();
195 
196  forAll(cellLabels, i)
197  {
198  const label celli = cellLabels[i];
199  const Type& value = values[i];
200 
201  psi[celli] = value;
202  source_[celli] = value*Diag[celli];
203 
204  if (symmetric() || asymmetric())
205  {
206  const cell& c = cells[celli];
207 
208  forAll(c, j)
209  {
210  const label facei = c[j];
211 
212  if (mesh.isInternalFace(facei))
213  {
214  if (symmetric())
215  {
216  if (celli == own[facei])
217  {
218  source_[nei[facei]] -= upper()[facei]*value;
219  }
220  else
221  {
222  source_[own[facei]] -= upper()[facei]*value;
223  }
224 
225  upper()[facei] = 0.0;
226  }
227  else
228  {
229  if (celli == own[facei])
230  {
231  source_[nei[facei]] -= lower()[facei]*value;
232  }
233  else
234  {
235  source_[own[facei]] -= upper()[facei]*value;
236  }
237 
238  upper()[facei] = 0.0;
239  lower()[facei] = 0.0;
240  }
241  }
242  else
243  {
244  label patchi = mesh.boundaryMesh().whichPatch(facei);
245 
246  if (internalCoeffs_[patchi].size())
247  {
248  label patchFacei =
249  mesh.boundaryMesh()[patchi].whichFace(facei);
250 
251  internalCoeffs_[patchi][patchFacei] =
252  Zero;
253 
254  boundaryCoeffs_[patchi][patchFacei] =
255  Zero;
256  }
257  }
258  }
259  }
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
266 template<class Type>
268 (
270  const dimensionSet& ds
271 )
272 :
273  lduMatrix(psi.mesh()),
274  psi_(psi),
275  dimensions_(ds),
276  source_(psi.size(), Zero),
277  internalCoeffs_(psi.mesh().boundary().size()),
278  boundaryCoeffs_(psi.mesh().boundary().size()),
279  faceFluxCorrectionPtr_(nullptr)
280 {
281  if (debug)
282  {
284  << "Constructing fvMatrix<Type> for field " << psi_.name() << endl;
285  }
286 
287  // Initialise coupling coefficients
288  forAll(psi.mesh().boundary(), patchi)
289  {
290  internalCoeffs_.set
291  (
292  patchi,
293  new Field<Type>
294  (
295  psi.mesh().boundary()[patchi].size(),
296  Zero
297  )
298  );
299 
300  boundaryCoeffs_.set
301  (
302  patchi,
303  new Field<Type>
304  (
305  psi.mesh().boundary()[patchi].size(),
306  Zero
307  )
308  );
309  }
310 
311  // Update the boundary coefficients of psi without changing its event No.
314 
315  label currentStatePsi = psiRef.eventNo();
316  psiRef.boundaryFieldRef().updateCoeffs();
317  psiRef.eventNo() = currentStatePsi;
318 }
319 
320 
321 template<class Type>
323 :
324  tmp<fvMatrix<Type>>::refCount(),
325  lduMatrix(fvm),
326  psi_(fvm.psi_),
327  dimensions_(fvm.dimensions_),
328  source_(fvm.source_),
329  internalCoeffs_(fvm.internalCoeffs_),
330  boundaryCoeffs_(fvm.boundaryCoeffs_),
331  faceFluxCorrectionPtr_(nullptr)
332 {
333  if (debug)
334  {
336  << "Copying fvMatrix<Type> for field " << psi_.name() << endl;
337  }
338 
339  if (fvm.faceFluxCorrectionPtr_)
340  {
341  faceFluxCorrectionPtr_ = new
343  (
344  *(fvm.faceFluxCorrectionPtr_)
345  );
346  }
347 }
348 
349 
350 template<class Type>
352 :
353  lduMatrix
354  (
355  const_cast<fvMatrix<Type>&>(tfvm()),
356  tfvm.isTmp()
357  ),
358  psi_(tfvm().psi_),
359  dimensions_(tfvm().dimensions_),
360  source_
361  (
362  const_cast<fvMatrix<Type>&>(tfvm()).source_,
363  tfvm.isTmp()
364  ),
365  internalCoeffs_
366  (
367  const_cast<fvMatrix<Type>&>(tfvm()).internalCoeffs_,
368  tfvm.isTmp()
369  ),
370  boundaryCoeffs_
371  (
372  const_cast<fvMatrix<Type>&>(tfvm()).boundaryCoeffs_,
373  tfvm.isTmp()
374  ),
375  faceFluxCorrectionPtr_(nullptr)
376 {
377  if (debug)
378  {
380  << "Copying fvMatrix<Type> for field " << psi_.name() << endl;
381  }
382 
383  if (tfvm().faceFluxCorrectionPtr_)
384  {
385  if (tfvm.isTmp())
386  {
387  faceFluxCorrectionPtr_ = tfvm().faceFluxCorrectionPtr_;
388  tfvm().faceFluxCorrectionPtr_ = nullptr;
389  }
390  else
391  {
392  faceFluxCorrectionPtr_ = new
394  (
395  *(tfvm().faceFluxCorrectionPtr_)
396  );
397  }
398  }
399 
400  tfvm.clear();
401 }
402 
403 
404 template<class Type>
406 (
408  Istream& is
409 )
410 :
411  lduMatrix(psi.mesh()),
412  psi_(psi),
413  dimensions_(is),
414  source_(is),
415  internalCoeffs_(psi.mesh().boundary().size()),
416  boundaryCoeffs_(psi.mesh().boundary().size()),
417  faceFluxCorrectionPtr_(nullptr)
418 {
419  if (debug)
420  {
422  << "Constructing fvMatrix<Type> for field " << psi_.name() << endl;
423  }
424 
425  // Initialise coupling coefficients
426  forAll(psi.mesh().boundary(), patchi)
427  {
428  internalCoeffs_.set
429  (
430  patchi,
431  new Field<Type>
432  (
433  psi.mesh().boundary()[patchi].size(),
434  Zero
435  )
436  );
437 
438  boundaryCoeffs_.set
439  (
440  patchi,
441  new Field<Type>
442  (
443  psi.mesh().boundary()[patchi].size(),
444  Zero
445  )
446  );
447  }
448 
449 }
450 
451 
452 template<class Type>
454 {
455  return tmp<fvMatrix<Type>>
456  (
457  new fvMatrix<Type>(*this)
458  );
459 }
460 
461 
462 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
463 
464 template<class Type>
466 {
467  if (debug)
468  {
470  << "Destroying fvMatrix<Type> for field " << psi_.name() << endl;
471  }
472 
473  if (faceFluxCorrectionPtr_)
474  {
475  delete faceFluxCorrectionPtr_;
476  }
477 }
478 
479 
480 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
481 
482 template<class Type>
484 (
485  const labelUList& cellLabels,
486  const UList<Type>& values
487 )
488 {
489  this->setValuesFromList(cellLabels, values);
490 }
491 
492 
493 template<class Type>
495 (
496  const labelUList& cellLabels,
497  const UIndirectList<Type>& values
498 )
499 {
500  this->setValuesFromList(cellLabels, values);
501 }
502 
503 
504 template<class Type>
506 (
507  const label celli,
508  const Type& value,
509  const bool forceReference
510 )
511 {
512  if ((forceReference || psi_.needReference()) && celli >= 0)
513  {
514  source()[celli] += diag()[celli]*value;
515  diag()[celli] += diag()[celli];
516  }
517 }
518 
519 
520 template<class Type>
522 {
523  if (alpha <= 0)
524  {
525  return;
526  }
527 
528  if (debug)
529  {
531  << "Relaxing " << psi_.name() << " by " << alpha << endl;
532  }
533 
534  Field<Type>& S = source();
535  scalarField& D = diag();
536 
537  // Store the current unrelaxed diagonal for use in updating the source
538  scalarField D0(D);
539 
540  // Calculate the sum-mag off-diagonal from the interior faces
541  scalarField sumOff(D.size(), 0.0);
542  sumMagOffDiag(sumOff);
543 
544  // Handle the boundary contributions to the diagonal
545  forAll(psi_.boundaryField(), patchi)
546  {
547  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
548 
549  if (ptf.size())
550  {
551  const labelUList& pa = lduAddr().patchAddr(patchi);
552  Field<Type>& iCoeffs = internalCoeffs_[patchi];
553 
554  if (ptf.coupled())
555  {
556  const Field<Type>& pCoeffs = boundaryCoeffs_[patchi];
557 
558  // For coupled boundaries add the diagonal and
559  // off-diagonal contributions
560  forAll(pa, face)
561  {
562  D[pa[face]] += component(iCoeffs[face], 0);
563  sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
564  }
565  }
566  else
567  {
568  // For non-coupled boundaries add the maximum magnitude diagonal
569  // contribution to ensure stability
570  forAll(pa, face)
571  {
572  D[pa[face]] += cmptMax(cmptMag(iCoeffs[face]));
573  }
574  }
575  }
576  }
577 
578 
579  if (debug)
580  {
581  // Calculate amount of non-dominance.
582  label nNon = 0;
583  scalar maxNon = 0.0;
584  scalar sumNon = 0.0;
585  forAll(D, celli)
586  {
587  scalar d = (sumOff[celli] - D[celli])/mag(D[celli]);
588 
589  if (d > 0)
590  {
591  nNon++;
592  maxNon = max(maxNon, d);
593  sumNon += d;
594  }
595  }
596 
597  reduce(nNon, sumOp<label>(), UPstream::msgType(), psi_.mesh().comm());
598  reduce
599  (
600  maxNon,
601  maxOp<scalar>(),
603  psi_.mesh().comm()
604  );
605  reduce
606  (
607  sumNon,
608  sumOp<scalar>(),
610  psi_.mesh().comm()
611  );
612  sumNon /= returnReduce
613  (
614  D.size(),
615  sumOp<label>(),
617  psi_.mesh().comm()
618  );
619 
621  << "Matrix dominance test for " << psi_.name() << nl
622  << " number of non-dominant cells : " << nNon << nl
623  << " maximum relative non-dominance : " << maxNon << nl
624  << " average relative non-dominance : " << sumNon << nl
625  << endl;
626  }
627 
628 
629  // Ensure the matrix is diagonally dominant...
630  // Assumes that the central coefficient is positive and ensures it is
631  forAll(D, celli)
632  {
633  D[celli] = max(mag(D[celli]), sumOff[celli]);
634  }
635 
636  // ... then relax
637  D /= alpha;
638 
639  // Now remove the diagonal contribution from coupled boundaries
640  forAll(psi_.boundaryField(), patchi)
641  {
642  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
643 
644  if (ptf.size())
645  {
646  const labelUList& pa = lduAddr().patchAddr(patchi);
647  Field<Type>& iCoeffs = internalCoeffs_[patchi];
648 
649  if (ptf.coupled())
650  {
651  forAll(pa, face)
652  {
653  D[pa[face]] -= component(iCoeffs[face], 0);
654  }
655  }
656  else
657  {
658  forAll(pa, face)
659  {
660  D[pa[face]] -= cmptMin(iCoeffs[face]);
661  }
662  }
663  }
664  }
665 
666  // Finally add the relaxation contribution to the source.
667  S += (D - D0)*psi_.primitiveField();
668 }
669 
670 
671 template<class Type>
673 {
674  if
675  (
676  psi_.mesh().data::template lookupOrDefault<bool>
677  ("finalIteration", false)
678  && psi_.mesh().relaxEquation(psi_.name() + "Final")
679  )
680  {
681  relax(psi_.mesh().equationRelaxationFactor(psi_.name() + "Final"));
682  }
683  else if (psi_.mesh().relaxEquation(psi_.name()))
684  {
685  relax(psi_.mesh().equationRelaxationFactor(psi_.name()));
686  }
687 }
688 
689 
690 template<class Type>
692 (
694  Boundary& bFields
695 )
696 {
697  forAll(bFields, patchi)
698  {
699  bFields[patchi].manipulateMatrix(*this);
700  }
701 }
702 
703 
704 template<class Type>
706 {
707  tmp<scalarField> tdiag(new scalarField(diag()));
708  addCmptAvBoundaryDiag(tdiag.ref());
709  return tdiag;
710 }
711 
712 
713 template<class Type>
715 {
717 
718  forAll(psi_.boundaryField(), patchi)
719  {
720  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
721 
722  if (!ptf.coupled() && ptf.size())
723  {
725  (
726  lduAddr().patchAddr(patchi),
727  internalCoeffs_[patchi],
728  tdiag.ref()
729  );
730  }
731  }
732 
733  return tdiag;
734 }
735 
736 
737 template<class Type>
739 {
740  tmp<volScalarField> tAphi
741  (
743  (
744  "A("+psi_.name()+')',
745  psi_.mesh(),
746  dimensions_/psi_.dimensions()/dimVol,
747  extrapolatedCalculatedFvPatchScalarField::typeName
748  )
749  );
750 
751  tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().V();
752  tAphi.ref().correctBoundaryConditions();
753 
754  return tAphi;
755 }
756 
757 
758 template<class Type>
761 {
763  (
765  (
766  "H("+psi_.name()+')',
767  psi_.mesh(),
768  dimensions_/dimVol,
769  extrapolatedCalculatedFvPatchScalarField::typeName
770  )
771  );
773 
774  // Loop over field components
775  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
776  {
777  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
778 
779  scalarField boundaryDiagCmpt(psi_.size(), 0.0);
780  addBoundaryDiag(boundaryDiagCmpt, cmpt);
781  boundaryDiagCmpt.negate();
782  addCmptAvBoundaryDiag(boundaryDiagCmpt);
783 
784  Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
785  }
786 
787  Hphi.primitiveFieldRef() += lduMatrix::H(psi_.primitiveField()) + source_;
789 
790  Hphi.primitiveFieldRef() /= psi_.mesh().V();
792 
793  typename Type::labelType validComponents
794  (
795  psi_.mesh().template validComponents<Type>()
796  );
797 
798  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
799  {
800  if (validComponents[cmpt] == -1)
801  {
802  Hphi.replace
803  (
804  cmpt,
805  dimensionedScalar(Hphi.dimensions(), 0)
806  );
807  }
808  }
809 
810  return tHphi;
811 }
812 
813 
814 template<class Type>
816 {
818  (
820  (
821  "H(1)",
822  psi_.mesh(),
823  dimensions_/(dimVol*psi_.dimensions()),
824  extrapolatedCalculatedFvPatchScalarField::typeName
825  )
826  );
827  volScalarField& H1_ = tH1.ref();
828 
829  H1_.primitiveFieldRef() = lduMatrix::H1();
830 
831  forAll(psi_.boundaryField(), patchi)
832  {
833  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
834 
835  if (ptf.coupled() && ptf.size())
836  {
838  (
839  lduAddr().patchAddr(patchi),
840  boundaryCoeffs_[patchi].component(0),
841  H1_
842  );
843  }
844  }
845 
846  H1_.primitiveFieldRef() /= psi_.mesh().V();
847  H1_.correctBoundaryConditions();
848 
849  return tH1;
850 }
851 
852 
853 template<class Type>
856 flux() const
857 {
858  if (!psi_.mesh().fluxRequired(psi_.name()))
859  {
861  << "flux requested but " << psi_.name()
862  << " not specified in the fluxRequired sub-dictionary"
863  " of fvSchemes."
864  << abort(FatalError);
865  }
866 
867  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
869  (
871  (
872  "flux("+psi_.name()+')',
873  psi_.mesh(),
874  dimensions()
875  )
876  );
878  tfieldFlux.ref();
879 
880  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
881  {
882  fieldFlux.primitiveFieldRef().replace
883  (
884  cmpt,
885  lduMatrix::faceH(psi_.primitiveField().component(cmpt))
886  );
887  }
888 
889  FieldField<Field, Type> InternalContrib = internalCoeffs_;
890 
891  forAll(InternalContrib, patchi)
892  {
893  InternalContrib[patchi] =
895  (
896  InternalContrib[patchi],
897  psi_.boundaryField()[patchi].patchInternalField()
898  );
899  }
900 
901  FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
902 
903  forAll(NeighbourContrib, patchi)
904  {
905  if (psi_.boundaryField()[patchi].coupled())
906  {
907  NeighbourContrib[patchi] =
909  (
910  NeighbourContrib[patchi],
911  psi_.boundaryField()[patchi].patchNeighbourField()
912  );
913  }
914  }
915 
917  Boundary& ffbf = fieldFlux.boundaryFieldRef();
918 
919  forAll(ffbf, patchi)
920  {
921  ffbf[patchi] = InternalContrib[patchi] - NeighbourContrib[patchi];
922  }
923 
924  if (faceFluxCorrectionPtr_)
925  {
926  fieldFlux += *faceFluxCorrectionPtr_;
927  }
928 
929  return tfieldFlux;
930 }
931 
932 
933 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
934 
935 template<class Type>
937 {
938  if (this == &fvmv)
939  {
941  << "attempted assignment to self"
942  << abort(FatalError);
943  }
944 
945  if (&psi_ != &(fvmv.psi_))
946  {
948  << "different fields"
949  << abort(FatalError);
950  }
951 
952  dimensions_ = fvmv.dimensions_;
953  lduMatrix::operator=(fvmv);
954  source_ = fvmv.source_;
955  internalCoeffs_ = fvmv.internalCoeffs_;
956  boundaryCoeffs_ = fvmv.boundaryCoeffs_;
957 
958  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
959  {
960  *faceFluxCorrectionPtr_ = *fvmv.faceFluxCorrectionPtr_;
961  }
962  else if (fvmv.faceFluxCorrectionPtr_)
963  {
964  faceFluxCorrectionPtr_ =
966  (*fvmv.faceFluxCorrectionPtr_);
967  }
968 }
969 
970 
971 template<class Type>
973 {
974  operator=(tfvmv());
975  tfvmv.clear();
976 }
977 
978 
979 template<class Type>
981 {
983  source_.negate();
984  internalCoeffs_.negate();
985  boundaryCoeffs_.negate();
986 
987  if (faceFluxCorrectionPtr_)
988  {
989  faceFluxCorrectionPtr_->negate();
990  }
991 }
992 
993 
994 template<class Type>
996 {
997  checkMethod(*this, fvmv, "+=");
998 
999  dimensions_ += fvmv.dimensions_;
1000  lduMatrix::operator+=(fvmv);
1001  source_ += fvmv.source_;
1002  internalCoeffs_ += fvmv.internalCoeffs_;
1003  boundaryCoeffs_ += fvmv.boundaryCoeffs_;
1004 
1005  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1006  {
1007  *faceFluxCorrectionPtr_ += *fvmv.faceFluxCorrectionPtr_;
1008  }
1009  else if (fvmv.faceFluxCorrectionPtr_)
1010  {
1011  faceFluxCorrectionPtr_ = new
1013  (
1014  *fvmv.faceFluxCorrectionPtr_
1015  );
1016  }
1017 }
1018 
1019 
1020 template<class Type>
1022 {
1023  operator+=(tfvmv());
1024  tfvmv.clear();
1025 }
1026 
1027 
1028 template<class Type>
1030 {
1031  checkMethod(*this, fvmv, "-=");
1032 
1033  dimensions_ -= fvmv.dimensions_;
1034  lduMatrix::operator-=(fvmv);
1035  source_ -= fvmv.source_;
1036  internalCoeffs_ -= fvmv.internalCoeffs_;
1037  boundaryCoeffs_ -= fvmv.boundaryCoeffs_;
1038 
1039  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1040  {
1041  *faceFluxCorrectionPtr_ -= *fvmv.faceFluxCorrectionPtr_;
1042  }
1043  else if (fvmv.faceFluxCorrectionPtr_)
1044  {
1045  faceFluxCorrectionPtr_ =
1047  (-*fvmv.faceFluxCorrectionPtr_);
1048  }
1049 }
1050 
1051 
1052 template<class Type>
1054 {
1055  operator-=(tfvmv());
1056  tfvmv.clear();
1057 }
1058 
1059 
1060 template<class Type>
1061 void Foam::fvMatrix<Type>::operator+=
1064 )
1065 {
1066  checkMethod(*this, su, "+=");
1067  source() -= su.mesh().V()*su.field();
1068 }
1069 
1070 
1071 template<class Type>
1072 void Foam::fvMatrix<Type>::operator+=
1075 )
1076 {
1077  operator+=(tsu());
1078  tsu.clear();
1079 }
1080 
1081 
1082 template<class Type>
1083 void Foam::fvMatrix<Type>::operator+=
1086 )
1087 {
1088  operator+=(tsu());
1089  tsu.clear();
1090 }
1091 
1092 
1093 template<class Type>
1094 void Foam::fvMatrix<Type>::operator-=
1097 )
1098 {
1099  checkMethod(*this, su, "-=");
1100  source() += su.mesh().V()*su.field();
1101 }
1102 
1103 
1104 template<class Type>
1105 void Foam::fvMatrix<Type>::operator-=
1108 )
1109 {
1110  operator-=(tsu());
1111  tsu.clear();
1112 }
1113 
1114 
1115 template<class Type>
1116 void Foam::fvMatrix<Type>::operator-=
1119 )
1120 {
1121  operator-=(tsu());
1122  tsu.clear();
1123 }
1124 
1125 
1126 template<class Type>
1127 void Foam::fvMatrix<Type>::operator+=
1129  const dimensioned<Type>& su
1130 )
1131 {
1132  source() -= psi().mesh().V()*su;
1133 }
1134 
1135 
1136 template<class Type>
1137 void Foam::fvMatrix<Type>::operator-=
1139  const dimensioned<Type>& su
1140 )
1141 {
1142  source() += psi().mesh().V()*su;
1143 }
1144 
1145 
1146 template<class Type>
1147 void Foam::fvMatrix<Type>::operator+=
1149  const zero&
1150 )
1151 {}
1152 
1153 
1154 template<class Type>
1155 void Foam::fvMatrix<Type>::operator-=
1157  const zero&
1158 )
1159 {}
1160 
1161 
1162 template<class Type>
1163 void Foam::fvMatrix<Type>::operator*=
1165  const volScalarField::Internal& dsf
1166 )
1167 {
1168  dimensions_ *= dsf.dimensions();
1169  lduMatrix::operator*=(dsf.field());
1170  source_ *= dsf.field();
1171 
1172  forAll(boundaryCoeffs_, patchi)
1173  {
1174  scalarField pisf
1175  (
1176  dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
1177  );
1178 
1179  internalCoeffs_[patchi] *= pisf;
1180  boundaryCoeffs_[patchi] *= pisf;
1181  }
1182 
1183  if (faceFluxCorrectionPtr_)
1184  {
1186  << "cannot scale a matrix containing a faceFluxCorrection"
1187  << abort(FatalError);
1188  }
1189 }
1190 
1191 
1192 template<class Type>
1193 void Foam::fvMatrix<Type>::operator*=
1195  const tmp<volScalarField::Internal>& tdsf
1196 )
1197 {
1198  operator*=(tdsf());
1199  tdsf.clear();
1200 }
1201 
1202 
1203 template<class Type>
1204 void Foam::fvMatrix<Type>::operator*=
1206  const tmp<volScalarField>& tvsf
1207 )
1208 {
1209  operator*=(tvsf());
1210  tvsf.clear();
1211 }
1212 
1213 
1214 template<class Type>
1215 void Foam::fvMatrix<Type>::operator*=
1217  const dimensioned<scalar>& ds
1218 )
1219 {
1220  dimensions_ *= ds.dimensions();
1221  lduMatrix::operator*=(ds.value());
1222  source_ *= ds.value();
1223  internalCoeffs_ *= ds.value();
1224  boundaryCoeffs_ *= ds.value();
1225 
1226  if (faceFluxCorrectionPtr_)
1227  {
1228  *faceFluxCorrectionPtr_ *= ds.value();
1229  }
1230 }
1231 
1232 
1233 template<class Type>
1234 void Foam::fvMatrix<Type>::operator/=
1236  const volScalarField::Internal& dsf
1237 )
1238 {
1239  dimensions_ /= dsf.dimensions();
1240  lduMatrix::operator/=(dsf.field());
1241  source_ /= dsf.field();
1242 
1243  forAll(boundaryCoeffs_, patchi)
1244  {
1245  scalarField pisf
1246  (
1247  dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
1248  );
1249 
1250  internalCoeffs_[patchi] /= pisf;
1251  boundaryCoeffs_[patchi] /= pisf;
1252  }
1253 
1254  if (faceFluxCorrectionPtr_)
1255  {
1257  << "cannot scale a matrix containing a faceFluxCorrection"
1258  << abort(FatalError);
1259  }
1260 }
1261 
1262 
1263 template<class Type>
1264 void Foam::fvMatrix<Type>::operator/=
1266  const tmp<volScalarField::Internal>& tdsf
1267 )
1268 {
1269  operator/=(tdsf());
1270  tdsf.clear();
1271 }
1272 
1273 
1274 template<class Type>
1275 void Foam::fvMatrix<Type>::operator/=
1277  const tmp<volScalarField>& tvsf
1278 )
1279 {
1280  operator/=(tvsf());
1281  tvsf.clear();
1282 }
1283 
1284 
1285 template<class Type>
1286 void Foam::fvMatrix<Type>::operator/=
1288  const dimensioned<scalar>& ds
1289 )
1290 {
1291  dimensions_ /= ds.dimensions();
1292  lduMatrix::operator/=(ds.value());
1293  source_ /= ds.value();
1294  internalCoeffs_ /= ds.value();
1295  boundaryCoeffs_ /= ds.value();
1296 
1297  if (faceFluxCorrectionPtr_)
1298  {
1299  *faceFluxCorrectionPtr_ /= ds.value();
1300  }
1301 }
1302 
1303 
1304 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
1305 
1306 template<class Type>
1307 void Foam::checkMethod
1309  const fvMatrix<Type>& fvm1,
1310  const fvMatrix<Type>& fvm2,
1311  const char* op
1312 )
1313 {
1314  if (&fvm1.psi() != &fvm2.psi())
1315  {
1317  << "incompatible fields for operation "
1318  << endl << " "
1319  << "[" << fvm1.psi().name() << "] "
1320  << op
1321  << " [" << fvm2.psi().name() << "]"
1322  << abort(FatalError);
1323  }
1324 
1325  if (dimensionSet::debug && fvm1.dimensions() != fvm2.dimensions())
1326  {
1328  << "incompatible dimensions for operation "
1329  << endl << " "
1330  << "[" << fvm1.psi().name() << fvm1.dimensions()/dimVolume << " ] "
1331  << op
1332  << " [" << fvm2.psi().name() << fvm2.dimensions()/dimVolume << " ]"
1333  << abort(FatalError);
1334  }
1335 }
1336 
1337 
1338 template<class Type>
1339 void Foam::checkMethod
1341  const fvMatrix<Type>& fvm,
1343  const char* op
1344 )
1345 {
1346  if (dimensionSet::debug && fvm.dimensions()/dimVolume != df.dimensions())
1347  {
1349  << endl << " "
1350  << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
1351  << op
1352  << " [" << df.name() << df.dimensions() << " ]"
1353  << abort(FatalError);
1354  }
1355 }
1356 
1357 
1358 template<class Type>
1359 void Foam::checkMethod
1361  const fvMatrix<Type>& fvm,
1362  const dimensioned<Type>& dt,
1363  const char* op
1364 )
1365 {
1366  if (dimensionSet::debug && fvm.dimensions()/dimVolume != dt.dimensions())
1367  {
1369  << "incompatible dimensions for operation "
1370  << endl << " "
1371  << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
1372  << op
1373  << " [" << dt.name() << dt.dimensions() << " ]"
1374  << abort(FatalError);
1375  }
1376 }
1377 
1378 
1379 template<class Type>
1381 (
1382  const fvMatrix<Type>& A
1383 )
1384 {
1385  tmp<Foam::fvMatrix<Type>> tAcorr = A - (A & A.psi());
1386 
1387  // Delete the faceFluxCorrection from the correction matrix
1388  // as it does not have a clear meaning or purpose
1389  deleteDemandDrivenData(tAcorr.ref().faceFluxCorrectionPtr());
1390 
1391  return tAcorr;
1392 }
1393 
1394 
1395 template<class Type>
1397 (
1398  const tmp<fvMatrix<Type>>& tA
1399 )
1400 {
1401  tmp<Foam::fvMatrix<Type>> tAcorr = tA - (tA() & tA().psi());
1402 
1403  // Delete the faceFluxCorrection from the correction matrix
1404  // as it does not have a clear meaning or purpose
1405  deleteDemandDrivenData(tAcorr.ref().faceFluxCorrectionPtr());
1406 
1407  return tAcorr;
1408 }
1409 
1410 
1411 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
1412 
1413 template<class Type>
1414 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1415 (
1416  const fvMatrix<Type>& A,
1417  const fvMatrix<Type>& B
1418 )
1419 {
1420  checkMethod(A, B, "==");
1421  return (A - B);
1422 }
1423 
1424 template<class Type>
1425 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1426 (
1427  const tmp<fvMatrix<Type>>& tA,
1428  const fvMatrix<Type>& B
1429 )
1430 {
1431  checkMethod(tA(), B, "==");
1432  return (tA - B);
1433 }
1434 
1435 template<class Type>
1436 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1437 (
1438  const fvMatrix<Type>& A,
1439  const tmp<fvMatrix<Type>>& tB
1440 )
1441 {
1442  checkMethod(A, tB(), "==");
1443  return (A - tB);
1444 }
1445 
1446 template<class Type>
1447 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1448 (
1449  const tmp<fvMatrix<Type>>& tA,
1450  const tmp<fvMatrix<Type>>& tB
1451 )
1452 {
1453  checkMethod(tA(), tB(), "==");
1454  return (tA - tB);
1455 }
1456 
1457 template<class Type>
1458 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1459 (
1460  const fvMatrix<Type>& A,
1462 )
1463 {
1464  checkMethod(A, su, "==");
1465  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1466  tC.ref().source() += su.mesh().V()*su.field();
1467  return tC;
1468 }
1469 
1470 template<class Type>
1471 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1472 (
1473  const fvMatrix<Type>& A,
1475 )
1476 {
1477  checkMethod(A, tsu(), "==");
1478  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1479  tC.ref().source() += tsu().mesh().V()*tsu().field();
1480  tsu.clear();
1481  return tC;
1482 }
1483 
1484 template<class Type>
1485 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1486 (
1487  const fvMatrix<Type>& A,
1489 )
1490 {
1491  checkMethod(A, tsu(), "==");
1492  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1493  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1494  tsu.clear();
1495  return tC;
1496 }
1497 
1498 template<class Type>
1499 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1500 (
1501  const tmp<fvMatrix<Type>>& tA,
1503 )
1504 {
1505  checkMethod(tA(), su, "==");
1506  tmp<fvMatrix<Type>> tC(tA.ptr());
1507  tC.ref().source() += su.mesh().V()*su.field();
1508  return tC;
1509 }
1510 
1511 template<class Type>
1512 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1513 (
1514  const tmp<fvMatrix<Type>>& tA,
1516 )
1517 {
1518  checkMethod(tA(), tsu(), "==");
1519  tmp<fvMatrix<Type>> tC(tA.ptr());
1520  tC.ref().source() += tsu().mesh().V()*tsu().field();
1521  tsu.clear();
1522  return tC;
1523 }
1524 
1525 template<class Type>
1526 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1527 (
1528  const tmp<fvMatrix<Type>>& tA,
1530 )
1531 {
1532  checkMethod(tA(), tsu(), "==");
1533  tmp<fvMatrix<Type>> tC(tA.ptr());
1534  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1535  tsu.clear();
1536  return tC;
1537 }
1538 
1539 template<class Type>
1540 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1541 (
1542  const fvMatrix<Type>& A,
1543  const dimensioned<Type>& su
1544 )
1545 {
1546  checkMethod(A, su, "==");
1547  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1548  tC.ref().source() += A.psi().mesh().V()*su.value();
1549  return tC;
1550 }
1551 
1552 template<class Type>
1553 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1554 (
1555  const tmp<fvMatrix<Type>>& tA,
1556  const dimensioned<Type>& su
1557 )
1558 {
1559  checkMethod(tA(), su, "==");
1560  tmp<fvMatrix<Type>> tC(tA.ptr());
1561  tC.ref().source() += tC().psi().mesh().V()*su.value();
1562  return tC;
1563 }
1564 
1565 template<class Type>
1566 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1567 (
1568  const fvMatrix<Type>& A,
1569  const zero&
1570 )
1571 {
1572  return A;
1573 }
1574 
1575 
1576 template<class Type>
1577 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1578 (
1579  const tmp<fvMatrix<Type>>& tA,
1580  const zero&
1581 )
1582 {
1583  return tA;
1584 }
1585 
1586 
1587 template<class Type>
1588 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1589 (
1590  const fvMatrix<Type>& A
1591 )
1592 {
1593  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1594  tC.ref().negate();
1595  return tC;
1596 }
1597 
1598 template<class Type>
1599 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1600 (
1601  const tmp<fvMatrix<Type>>& tA
1602 )
1603 {
1604  tmp<fvMatrix<Type>> tC(tA.ptr());
1605  tC.ref().negate();
1606  return tC;
1607 }
1608 
1609 
1610 template<class Type>
1611 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1612 (
1613  const fvMatrix<Type>& A,
1614  const fvMatrix<Type>& B
1615 )
1616 {
1617  checkMethod(A, B, "+");
1618  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1619  tC.ref() += B;
1620  return tC;
1621 }
1622 
1623 template<class Type>
1624 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1625 (
1626  const tmp<fvMatrix<Type>>& tA,
1627  const fvMatrix<Type>& B
1628 )
1629 {
1630  checkMethod(tA(), B, "+");
1631  tmp<fvMatrix<Type>> tC(tA.ptr());
1632  tC.ref() += B;
1633  return tC;
1634 }
1635 
1636 template<class Type>
1637 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1638 (
1639  const fvMatrix<Type>& A,
1640  const tmp<fvMatrix<Type>>& tB
1641 )
1642 {
1643  checkMethod(A, tB(), "+");
1644  tmp<fvMatrix<Type>> tC(tB.ptr());
1645  tC.ref() += A;
1646  return tC;
1647 }
1648 
1649 template<class Type>
1650 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1651 (
1652  const tmp<fvMatrix<Type>>& tA,
1653  const tmp<fvMatrix<Type>>& tB
1654 )
1655 {
1656  checkMethod(tA(), tB(), "+");
1657  tmp<fvMatrix<Type>> tC(tA.ptr());
1658  tC.ref() += tB();
1659  tB.clear();
1660  return tC;
1661 }
1662 
1663 template<class Type>
1664 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1665 (
1666  const fvMatrix<Type>& A,
1668 )
1669 {
1670  checkMethod(A, su, "+");
1671  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1672  tC.ref().source() -= su.mesh().V()*su.field();
1673  return tC;
1674 }
1675 
1676 template<class Type>
1677 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1678 (
1679  const fvMatrix<Type>& A,
1681 )
1682 {
1683  checkMethod(A, tsu(), "+");
1684  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1685  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1686  tsu.clear();
1687  return tC;
1688 }
1689 
1690 template<class Type>
1691 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1692 (
1693  const fvMatrix<Type>& A,
1695 )
1696 {
1697  checkMethod(A, tsu(), "+");
1698  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1699  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1700  tsu.clear();
1701  return tC;
1702 }
1703 
1704 template<class Type>
1705 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1706 (
1707  const tmp<fvMatrix<Type>>& tA,
1709 )
1710 {
1711  checkMethod(tA(), su, "+");
1712  tmp<fvMatrix<Type>> tC(tA.ptr());
1713  tC.ref().source() -= su.mesh().V()*su.field();
1714  return tC;
1715 }
1716 
1717 template<class Type>
1718 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1719 (
1720  const tmp<fvMatrix<Type>>& tA,
1722 )
1723 {
1724  checkMethod(tA(), tsu(), "+");
1725  tmp<fvMatrix<Type>> tC(tA.ptr());
1726  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1727  tsu.clear();
1728  return tC;
1729 }
1730 
1731 template<class Type>
1732 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1733 (
1734  const tmp<fvMatrix<Type>>& tA,
1736 )
1737 {
1738  checkMethod(tA(), tsu(), "+");
1739  tmp<fvMatrix<Type>> tC(tA.ptr());
1740  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1741  tsu.clear();
1742  return tC;
1743 }
1744 
1745 template<class Type>
1746 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1747 (
1749  const fvMatrix<Type>& A
1750 )
1751 {
1752  checkMethod(A, su, "+");
1753  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1754  tC.ref().source() -= su.mesh().V()*su.field();
1755  return tC;
1756 }
1757 
1758 template<class Type>
1759 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1760 (
1762  const fvMatrix<Type>& A
1763 )
1764 {
1765  checkMethod(A, tsu(), "+");
1766  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1767  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1768  tsu.clear();
1769  return tC;
1770 }
1771 
1772 template<class Type>
1773 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1774 (
1776  const fvMatrix<Type>& A
1777 )
1778 {
1779  checkMethod(A, tsu(), "+");
1780  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1781  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1782  tsu.clear();
1783  return tC;
1784 }
1785 
1786 template<class Type>
1787 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1788 (
1790  const tmp<fvMatrix<Type>>& tA
1791 )
1792 {
1793  checkMethod(tA(), su, "+");
1794  tmp<fvMatrix<Type>> tC(tA.ptr());
1795  tC.ref().source() -= su.mesh().V()*su.field();
1796  return tC;
1797 }
1798 
1799 template<class Type>
1800 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1801 (
1803  const tmp<fvMatrix<Type>>& tA
1804 )
1805 {
1806  checkMethod(tA(), tsu(), "+");
1807  tmp<fvMatrix<Type>> tC(tA.ptr());
1808  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1809  tsu.clear();
1810  return tC;
1811 }
1812 
1813 template<class Type>
1814 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1815 (
1817  const tmp<fvMatrix<Type>>& tA
1818 )
1819 {
1820  checkMethod(tA(), tsu(), "+");
1821  tmp<fvMatrix<Type>> tC(tA.ptr());
1822  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1823  tsu.clear();
1824  return tC;
1825 }
1826 
1827 
1828 template<class Type>
1829 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1830 (
1831  const fvMatrix<Type>& A,
1832  const fvMatrix<Type>& B
1833 )
1834 {
1835  checkMethod(A, B, "-");
1836  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1837  tC.ref() -= B;
1838  return tC;
1839 }
1840 
1841 template<class Type>
1842 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1843 (
1844  const tmp<fvMatrix<Type>>& tA,
1845  const fvMatrix<Type>& B
1846 )
1847 {
1848  checkMethod(tA(), B, "-");
1849  tmp<fvMatrix<Type>> tC(tA.ptr());
1850  tC.ref() -= B;
1851  return tC;
1852 }
1853 
1854 template<class Type>
1855 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1856 (
1857  const fvMatrix<Type>& A,
1858  const tmp<fvMatrix<Type>>& tB
1859 )
1860 {
1861  checkMethod(A, tB(), "-");
1862  tmp<fvMatrix<Type>> tC(tB.ptr());
1863  tC.ref() -= A;
1864  tC.ref().negate();
1865  return tC;
1866 }
1867 
1868 template<class Type>
1869 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1870 (
1871  const tmp<fvMatrix<Type>>& tA,
1872  const tmp<fvMatrix<Type>>& tB
1873 )
1874 {
1875  checkMethod(tA(), tB(), "-");
1876  tmp<fvMatrix<Type>> tC(tA.ptr());
1877  tC.ref() -= tB();
1878  tB.clear();
1879  return tC;
1880 }
1881 
1882 template<class Type>
1883 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1884 (
1885  const fvMatrix<Type>& A,
1887 )
1888 {
1889  checkMethod(A, su, "-");
1890  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1891  tC.ref().source() += su.mesh().V()*su.field();
1892  return tC;
1893 }
1894 
1895 template<class Type>
1896 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1897 (
1898  const fvMatrix<Type>& A,
1900 )
1901 {
1902  checkMethod(A, tsu(), "-");
1903  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1904  tC.ref().source() += tsu().mesh().V()*tsu().field();
1905  tsu.clear();
1906  return tC;
1907 }
1908 
1909 template<class Type>
1910 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1911 (
1912  const fvMatrix<Type>& A,
1914 )
1915 {
1916  checkMethod(A, tsu(), "-");
1917  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1918  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1919  tsu.clear();
1920  return tC;
1921 }
1922 
1923 template<class Type>
1924 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1925 (
1926  const tmp<fvMatrix<Type>>& tA,
1928 )
1929 {
1930  checkMethod(tA(), su, "-");
1931  tmp<fvMatrix<Type>> tC(tA.ptr());
1932  tC.ref().source() += su.mesh().V()*su.field();
1933  return tC;
1934 }
1935 
1936 template<class Type>
1937 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1938 (
1939  const tmp<fvMatrix<Type>>& tA,
1941 )
1942 {
1943  checkMethod(tA(), tsu(), "-");
1944  tmp<fvMatrix<Type>> tC(tA.ptr());
1945  tC.ref().source() += tsu().mesh().V()*tsu().field();
1946  tsu.clear();
1947  return tC;
1948 }
1949 
1950 template<class Type>
1951 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1952 (
1953  const tmp<fvMatrix<Type>>& tA,
1955 )
1956 {
1957  checkMethod(tA(), tsu(), "-");
1958  tmp<fvMatrix<Type>> tC(tA.ptr());
1959  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1960  tsu.clear();
1961  return tC;
1962 }
1963 
1964 template<class Type>
1965 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1966 (
1968  const fvMatrix<Type>& A
1969 )
1970 {
1971  checkMethod(A, su, "-");
1972  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1973  tC.ref().negate();
1974  tC.ref().source() -= su.mesh().V()*su.field();
1975  return tC;
1976 }
1977 
1978 template<class Type>
1979 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1980 (
1982  const fvMatrix<Type>& A
1983 )
1984 {
1985  checkMethod(A, tsu(), "-");
1986  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1987  tC.ref().negate();
1988  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1989  tsu.clear();
1990  return tC;
1991 }
1992 
1993 template<class Type>
1994 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1995 (
1997  const fvMatrix<Type>& A
1998 )
1999 {
2000  checkMethod(A, tsu(), "-");
2001  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2002  tC.ref().negate();
2003  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2004  tsu.clear();
2005  return tC;
2006 }
2007 
2008 template<class Type>
2009 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2010 (
2012  const tmp<fvMatrix<Type>>& tA
2013 )
2014 {
2015  checkMethod(tA(), su, "-");
2016  tmp<fvMatrix<Type>> tC(tA.ptr());
2017  tC.ref().negate();
2018  tC.ref().source() -= su.mesh().V()*su.field();
2019  return tC;
2020 }
2021 
2022 template<class Type>
2023 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2024 (
2026  const tmp<fvMatrix<Type>>& tA
2027 )
2028 {
2029  checkMethod(tA(), tsu(), "-");
2030  tmp<fvMatrix<Type>> tC(tA.ptr());
2031  tC.ref().negate();
2032  tC.ref().source() -= tsu().mesh().V()*tsu().field();
2033  tsu.clear();
2034  return tC;
2035 }
2036 
2037 template<class Type>
2038 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2039 (
2041  const tmp<fvMatrix<Type>>& tA
2042 )
2043 {
2044  checkMethod(tA(), tsu(), "-");
2045  tmp<fvMatrix<Type>> tC(tA.ptr());
2046  tC.ref().negate();
2047  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2048  tsu.clear();
2049  return tC;
2050 }
2051 
2052 template<class Type>
2053 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2054 (
2055  const fvMatrix<Type>& A,
2056  const dimensioned<Type>& su
2057 )
2058 {
2059  checkMethod(A, su, "+");
2060  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2061  tC.ref().source() -= su.value()*A.psi().mesh().V();
2062  return tC;
2063 }
2064 
2065 template<class Type>
2066 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2067 (
2068  const tmp<fvMatrix<Type>>& tA,
2069  const dimensioned<Type>& su
2070 )
2071 {
2072  checkMethod(tA(), su, "+");
2073  tmp<fvMatrix<Type>> tC(tA.ptr());
2074  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2075  return tC;
2076 }
2077 
2078 template<class Type>
2079 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2080 (
2081  const dimensioned<Type>& su,
2082  const fvMatrix<Type>& A
2083 )
2084 {
2085  checkMethod(A, su, "+");
2086  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2087  tC.ref().source() -= su.value()*A.psi().mesh().V();
2088  return tC;
2089 }
2090 
2091 template<class Type>
2092 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2093 (
2094  const dimensioned<Type>& su,
2095  const tmp<fvMatrix<Type>>& tA
2096 )
2097 {
2098  checkMethod(tA(), su, "+");
2099  tmp<fvMatrix<Type>> tC(tA.ptr());
2100  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2101  return tC;
2102 }
2103 
2104 template<class Type>
2105 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2106 (
2107  const fvMatrix<Type>& A,
2108  const dimensioned<Type>& su
2109 )
2110 {
2111  checkMethod(A, su, "-");
2112  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2113  tC.ref().source() += su.value()*tC().psi().mesh().V();
2114  return tC;
2115 }
2116 
2117 template<class Type>
2118 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2119 (
2120  const tmp<fvMatrix<Type>>& tA,
2121  const dimensioned<Type>& su
2122 )
2123 {
2124  checkMethod(tA(), su, "-");
2125  tmp<fvMatrix<Type>> tC(tA.ptr());
2126  tC.ref().source() += su.value()*tC().psi().mesh().V();
2127  return tC;
2128 }
2129 
2130 template<class Type>
2131 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2132 (
2133  const dimensioned<Type>& su,
2134  const fvMatrix<Type>& A
2135 )
2136 {
2137  checkMethod(A, su, "-");
2138  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2139  tC.ref().negate();
2140  tC.ref().source() -= su.value()*A.psi().mesh().V();
2141  return tC;
2142 }
2143 
2144 template<class Type>
2145 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2146 (
2147  const dimensioned<Type>& su,
2148  const tmp<fvMatrix<Type>>& tA
2149 )
2150 {
2151  checkMethod(tA(), su, "-");
2152  tmp<fvMatrix<Type>> tC(tA.ptr());
2153  tC.ref().negate();
2154  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2155  return tC;
2156 }
2157 
2158 
2159 template<class Type>
2160 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2161 (
2162  const volScalarField::Internal& dsf,
2163  const fvMatrix<Type>& A
2164 )
2165 {
2166  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2167  tC.ref() *= dsf;
2168  return tC;
2169 }
2170 
2171 template<class Type>
2172 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2173 (
2174  const tmp<volScalarField::Internal>& tdsf,
2175  const fvMatrix<Type>& A
2176 )
2177 {
2178  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2179  tC.ref() *= tdsf;
2180  return tC;
2181 }
2182 
2183 template<class Type>
2184 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2185 (
2186  const tmp<volScalarField>& tvsf,
2187  const fvMatrix<Type>& A
2188 )
2189 {
2190  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2191  tC.ref() *= tvsf;
2192  return tC;
2193 }
2194 
2195 template<class Type>
2196 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2197 (
2198  const volScalarField::Internal& dsf,
2199  const tmp<fvMatrix<Type>>& tA
2200 )
2201 {
2202  tmp<fvMatrix<Type>> tC(tA.ptr());
2203  tC.ref() *= dsf;
2204  return tC;
2205 }
2206 
2207 template<class Type>
2208 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2209 (
2210  const tmp<volScalarField::Internal>& tdsf,
2211  const tmp<fvMatrix<Type>>& tA
2212 )
2213 {
2214  tmp<fvMatrix<Type>> tC(tA.ptr());
2215  tC.ref() *= tdsf;
2216  return tC;
2217 }
2218 
2219 template<class Type>
2220 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2221 (
2222  const tmp<volScalarField>& tvsf,
2223  const tmp<fvMatrix<Type>>& tA
2224 )
2225 {
2226  tmp<fvMatrix<Type>> tC(tA.ptr());
2227  tC.ref() *= tvsf;
2228  return tC;
2229 }
2230 
2231 template<class Type>
2232 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2233 (
2234  const dimensioned<scalar>& ds,
2235  const fvMatrix<Type>& A
2236 )
2237 {
2238  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2239  tC.ref() *= ds;
2240  return tC;
2241 }
2242 
2243 template<class Type>
2244 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2245 (
2246  const dimensioned<scalar>& ds,
2247  const tmp<fvMatrix<Type>>& tA
2248 )
2249 {
2250  tmp<fvMatrix<Type>> tC(tA.ptr());
2251  tC.ref() *= ds;
2252  return tC;
2253 }
2254 
2255 
2256 template<class Type>
2257 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2258 (
2259  const fvMatrix<Type>& A,
2260  const volScalarField::Internal& dsf
2261 )
2262 {
2263  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2264  tC.ref() /= dsf;
2265  return tC;
2266 }
2267 
2268 template<class Type>
2269 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2270 (
2271  const fvMatrix<Type>& A,
2272  const tmp<volScalarField::Internal>& tdsf
2273 )
2274 {
2275  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2276  tC.ref() /= tdsf;
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<volScalarField>& tvsf
2285 )
2286 {
2287  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2288  tC.ref() /= tvsf;
2289  return tC;
2290 }
2291 
2292 template<class Type>
2293 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2294 (
2295  const tmp<fvMatrix<Type>>& tA,
2296  const volScalarField::Internal& dsf
2297 )
2298 {
2299  tmp<fvMatrix<Type>> tC(tA.ptr());
2300  tC.ref() /= dsf;
2301  return tC;
2302 }
2303 
2304 template<class Type>
2305 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2306 (
2307  const tmp<fvMatrix<Type>>& tA,
2308  const tmp<volScalarField::Internal>& tdsf
2309 )
2310 {
2311  tmp<fvMatrix<Type>> tC(tA.ptr());
2312  tC.ref() /= tdsf;
2313  return tC;
2314 }
2315 
2316 template<class Type>
2317 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2318 (
2319  const tmp<fvMatrix<Type>>& tA,
2320  const tmp<volScalarField>& tvsf
2321 )
2322 {
2323  tmp<fvMatrix<Type>> tC(tA.ptr());
2324  tC.ref() /= tvsf;
2325  return tC;
2326 }
2327 
2328 template<class Type>
2329 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2330 (
2331  const fvMatrix<Type>& A,
2332  const dimensioned<scalar>& ds
2333 )
2334 {
2335  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2336  tC.ref() /= ds;
2337  return tC;
2338 }
2339 
2340 template<class Type>
2341 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator/
2342 (
2343  const tmp<fvMatrix<Type>>& tA,
2344  const dimensioned<scalar>& ds
2345 )
2346 {
2347  tmp<fvMatrix<Type>> tC(tA.ptr());
2348  tC.ref() /= ds;
2349  return tC;
2350 }
2351 
2352 
2353 template<class Type>
2355 Foam::operator&
2356 (
2357  const fvMatrix<Type>& M,
2359 )
2360 {
2362  (
2364  (
2365  "M&" + psi.name(),
2366  psi.mesh(),
2367  M.dimensions()/dimVol,
2368  extrapolatedCalculatedFvPatchScalarField::typeName
2369  )
2370  );
2372 
2373  // Loop over field components
2374  if (M.hasDiag())
2375  {
2376  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2377  {
2378  scalarField psiCmpt(psi.field().component(cmpt));
2379  scalarField boundaryDiagCmpt(M.diag());
2380  M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2381  Mphi.primitiveFieldRef().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2382  }
2383  }
2384  else
2385  {
2386  Mphi.primitiveFieldRef() = Zero;
2387  }
2388 
2389  Mphi.primitiveFieldRef() += M.lduMatrix::H(psi.field()) + M.source();
2390  M.addBoundarySource(Mphi.primitiveFieldRef());
2391 
2392  Mphi.primitiveFieldRef() /= -psi.mesh().V();
2394 
2395  return tMphi;
2396 }
2397 
2398 template<class Type>
2400 Foam::operator&
2401 (
2402  const fvMatrix<Type>& M,
2404 )
2405 {
2407  tpsi.clear();
2408  return tMpsi;
2409 }
2410 
2411 template<class Type>
2413 Foam::operator&
2414 (
2415  const fvMatrix<Type>& M,
2417 )
2418 {
2420  tpsi.clear();
2421  return tMpsi;
2422 }
2423 
2424 template<class Type>
2426 Foam::operator&
2427 (
2428  const tmp<fvMatrix<Type>>& tM,
2430 )
2431 {
2433  tM.clear();
2434  return tMpsi;
2435 }
2436 
2437 template<class Type>
2439 Foam::operator&
2440 (
2441  const tmp<fvMatrix<Type>>& tM,
2443 )
2444 {
2445  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2446  tM.clear();
2447  tpsi.clear();
2448  return tMpsi;
2449 }
2450 
2451 template<class Type>
2453 Foam::operator&
2454 (
2455  const tmp<fvMatrix<Type>>& tM,
2457 )
2458 {
2459  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2460  tM.clear();
2461  tpsi.clear();
2462  return tMpsi;
2463 }
2464 
2465 
2466 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2467 
2468 template<class Type>
2469 Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
2470 {
2471  os << static_cast<const lduMatrix&>(fvm) << nl
2472  << fvm.dimensions_ << nl
2473  << fvm.source_ << nl
2474  << fvm.internalCoeffs_ << nl
2475  << fvm.boundaryCoeffs_ << endl;
2476 
2477  os.check("Ostream& operator<<(Ostream&, fvMatrix<Type>&");
2478 
2479  return os;
2480 }
2481 
2482 
2483 // * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
2484 
2485 #include "fvMatrixSolve.C"
2486 
2487 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:38
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
tmp< Field< Type > > faceH(const Field< Type > &) const
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:111
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Reference counter for various OpenFOAM components.
Definition: refCount.H:49
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:672
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
void negate()
Definition: fvMatrix.C:980
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:43
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
uint8_t direction
Definition: direction.H:45
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:129
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:506
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:760
Traits class for primitives.
Definition: pTraits.H:50
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1029
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const cellList & cells() const
Generic GeometricField class.
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
Generic dimensioned Type class.
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:815
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:75
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
tmp< scalarField > H1() const
void sumMagOffDiag(scalarField &sumOff) const
const dimensionedScalar c
Speed of light in a vacuum.
conserve primitiveFieldRef()+
Generic field type.
Definition: FieldField.H:51
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition: fvMatrix.C:705
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:936
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:284
const dimensionSet & dimensions() const
Return dimensions.
Dimension set for the base types.
Definition: dimensionSet.H:120
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dynamicFvMesh & mesh
const cellShapeList & cells
void operator/=(const volScalarField::Internal &)
Definition: fvMatrix.C:1235
Pre-declare SubField and related Field type.
Definition: Field.H:56
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:411
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:484
const dimensionSet dimVol
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:465
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:309
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:72
label eventNo() const
Event number at last update.
Definition: regIOobjectI.H:83
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void operator=(const lduMatrix &)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > clone() const
Clone.
Definition: fvMatrix.C:453
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:260
void operator*=(const scalarField &)
Field< Type > & source()
Definition: fvMatrix.H:292
void operator/=(const scalarField &)
const Mesh & mesh() const
Return mesh.
void operator*=(const volScalarField::Internal &)
Definition: fvMatrix.C:1164
const word & name() const
Return const reference to name.
const Field< Type > & field() const
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
Internal & ref()
Return a reference to the dimensioned internal field.
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:995
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:714
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:738
label patchi
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:550
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return const reference to dimensions.
const dimensionSet dimVolume
void setValuesFromList(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:178
void correctBoundaryConditions()
Correct boundary field.
void checkMethod(const fvMatrix< Type > &, const fvMatrix< Type > &, const char *)
Definition: fvMatrix.C:1308
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:49
dimensioned< scalar > mag(const dimensioned< Type > &)
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:145
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:856
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:692
tmp< Field< Type > > H(const Field< Type > &) const
void operator-=(const lduMatrix &)
void operator+=(const lduMatrix &)
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:205
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField & diag()
Definition: lduMatrix.C:186
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
void deleteDemandDrivenData(DataPtr &dataPtr)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &, const dimensionSet &)
Construct given a field to solve for.
Definition: fvMatrix.C:268
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
#define M(I)
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionSet & dimensions() const
Definition: fvMatrix.H:287