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-2019 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 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
1234 
1235 template<class Type>
1236 void Foam::checkMethod
1238  const fvMatrix<Type>& fvm1,
1239  const fvMatrix<Type>& fvm2,
1240  const char* op
1241 )
1242 {
1243  if (&fvm1.psi() != &fvm2.psi())
1244  {
1246  << "incompatible fields for operation "
1247  << endl << " "
1248  << "[" << fvm1.psi().name() << "] "
1249  << op
1250  << " [" << fvm2.psi().name() << "]"
1251  << abort(FatalError);
1252  }
1253 
1254  if (dimensionSet::debug && fvm1.dimensions() != fvm2.dimensions())
1255  {
1257  << "incompatible dimensions for operation "
1258  << endl << " "
1259  << "[" << fvm1.psi().name() << fvm1.dimensions()/dimVolume << " ] "
1260  << op
1261  << " [" << fvm2.psi().name() << fvm2.dimensions()/dimVolume << " ]"
1262  << abort(FatalError);
1263  }
1264 }
1265 
1266 
1267 template<class Type>
1268 void Foam::checkMethod
1270  const fvMatrix<Type>& fvm,
1272  const char* op
1273 )
1274 {
1275  if (dimensionSet::debug && fvm.dimensions()/dimVolume != df.dimensions())
1276  {
1278  << endl << " "
1279  << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
1280  << op
1281  << " [" << df.name() << df.dimensions() << " ]"
1282  << abort(FatalError);
1283  }
1284 }
1285 
1286 
1287 template<class Type>
1288 void Foam::checkMethod
1290  const fvMatrix<Type>& fvm,
1291  const dimensioned<Type>& dt,
1292  const char* op
1293 )
1294 {
1295  if (dimensionSet::debug && fvm.dimensions()/dimVolume != dt.dimensions())
1296  {
1298  << "incompatible dimensions for operation "
1299  << endl << " "
1300  << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
1301  << op
1302  << " [" << dt.name() << dt.dimensions() << " ]"
1303  << abort(FatalError);
1304  }
1305 }
1306 
1307 
1308 template<class Type>
1310 (
1311  const fvMatrix<Type>& A
1312 )
1313 {
1314  tmp<Foam::fvMatrix<Type>> tAcorr = A - (A & A.psi());
1315 
1316  // Delete the faceFluxCorrection from the correction matrix
1317  // as it does not have a clear meaning or purpose
1318  deleteDemandDrivenData(tAcorr.ref().faceFluxCorrectionPtr());
1319 
1320  return tAcorr;
1321 }
1322 
1323 
1324 template<class Type>
1326 (
1327  const tmp<fvMatrix<Type>>& tA
1328 )
1329 {
1330  tmp<Foam::fvMatrix<Type>> tAcorr = tA - (tA() & tA().psi());
1331 
1332  // Delete the faceFluxCorrection from the correction matrix
1333  // as it does not have a clear meaning or purpose
1334  deleteDemandDrivenData(tAcorr.ref().faceFluxCorrectionPtr());
1335 
1336  return tAcorr;
1337 }
1338 
1339 
1340 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
1341 
1342 template<class Type>
1343 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1344 (
1345  const fvMatrix<Type>& A,
1346  const fvMatrix<Type>& B
1347 )
1348 {
1349  checkMethod(A, B, "==");
1350  return (A - B);
1351 }
1352 
1353 template<class Type>
1354 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1355 (
1356  const tmp<fvMatrix<Type>>& tA,
1357  const fvMatrix<Type>& B
1358 )
1359 {
1360  checkMethod(tA(), B, "==");
1361  return (tA - B);
1362 }
1363 
1364 template<class Type>
1365 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1366 (
1367  const fvMatrix<Type>& A,
1368  const tmp<fvMatrix<Type>>& tB
1369 )
1370 {
1371  checkMethod(A, tB(), "==");
1372  return (A - tB);
1373 }
1374 
1375 template<class Type>
1376 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1377 (
1378  const tmp<fvMatrix<Type>>& tA,
1379  const tmp<fvMatrix<Type>>& tB
1380 )
1381 {
1382  checkMethod(tA(), tB(), "==");
1383  return (tA - tB);
1384 }
1385 
1386 template<class Type>
1387 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1388 (
1389  const fvMatrix<Type>& A,
1391 )
1392 {
1393  checkMethod(A, su, "==");
1394  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1395  tC.ref().source() += su.mesh().V()*su.field();
1396  return tC;
1397 }
1398 
1399 template<class Type>
1400 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1401 (
1402  const fvMatrix<Type>& A,
1404 )
1405 {
1406  checkMethod(A, tsu(), "==");
1407  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1408  tC.ref().source() += tsu().mesh().V()*tsu().field();
1409  tsu.clear();
1410  return tC;
1411 }
1412 
1413 template<class Type>
1414 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1415 (
1416  const fvMatrix<Type>& A,
1418 )
1419 {
1420  checkMethod(A, tsu(), "==");
1421  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1422  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1423  tsu.clear();
1424  return tC;
1425 }
1426 
1427 template<class Type>
1428 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1429 (
1430  const tmp<fvMatrix<Type>>& tA,
1432 )
1433 {
1434  checkMethod(tA(), su, "==");
1435  tmp<fvMatrix<Type>> tC(tA.ptr());
1436  tC.ref().source() += su.mesh().V()*su.field();
1437  return tC;
1438 }
1439 
1440 template<class Type>
1441 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1442 (
1443  const tmp<fvMatrix<Type>>& tA,
1445 )
1446 {
1447  checkMethod(tA(), tsu(), "==");
1448  tmp<fvMatrix<Type>> tC(tA.ptr());
1449  tC.ref().source() += tsu().mesh().V()*tsu().field();
1450  tsu.clear();
1451  return tC;
1452 }
1453 
1454 template<class Type>
1455 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1456 (
1457  const tmp<fvMatrix<Type>>& tA,
1459 )
1460 {
1461  checkMethod(tA(), tsu(), "==");
1462  tmp<fvMatrix<Type>> tC(tA.ptr());
1463  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1464  tsu.clear();
1465  return tC;
1466 }
1467 
1468 template<class Type>
1469 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1470 (
1471  const fvMatrix<Type>& A,
1472  const dimensioned<Type>& su
1473 )
1474 {
1475  checkMethod(A, su, "==");
1476  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1477  tC.ref().source() += A.psi().mesh().V()*su.value();
1478  return tC;
1479 }
1480 
1481 template<class Type>
1482 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1483 (
1484  const tmp<fvMatrix<Type>>& tA,
1485  const dimensioned<Type>& su
1486 )
1487 {
1488  checkMethod(tA(), su, "==");
1489  tmp<fvMatrix<Type>> tC(tA.ptr());
1490  tC.ref().source() += tC().psi().mesh().V()*su.value();
1491  return tC;
1492 }
1493 
1494 template<class Type>
1495 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1496 (
1497  const fvMatrix<Type>& A,
1498  const zero&
1499 )
1500 {
1501  return A;
1502 }
1503 
1504 
1505 template<class Type>
1506 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
1507 (
1508  const tmp<fvMatrix<Type>>& tA,
1509  const zero&
1510 )
1511 {
1512  return tA;
1513 }
1514 
1515 
1516 template<class Type>
1517 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1518 (
1519  const fvMatrix<Type>& A
1520 )
1521 {
1522  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1523  tC.ref().negate();
1524  return tC;
1525 }
1526 
1527 template<class Type>
1528 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1529 (
1530  const tmp<fvMatrix<Type>>& tA
1531 )
1532 {
1533  tmp<fvMatrix<Type>> tC(tA.ptr());
1534  tC.ref().negate();
1535  return tC;
1536 }
1537 
1538 
1539 template<class Type>
1540 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1541 (
1542  const fvMatrix<Type>& A,
1543  const fvMatrix<Type>& B
1544 )
1545 {
1546  checkMethod(A, B, "+");
1547  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1548  tC.ref() += B;
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 fvMatrix<Type>& B
1557 )
1558 {
1559  checkMethod(tA(), B, "+");
1560  tmp<fvMatrix<Type>> tC(tA.ptr());
1561  tC.ref() += B;
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 tmp<fvMatrix<Type>>& tB
1570 )
1571 {
1572  checkMethod(A, tB(), "+");
1573  tmp<fvMatrix<Type>> tC(tB.ptr());
1574  tC.ref() += A;
1575  return tC;
1576 }
1577 
1578 template<class Type>
1579 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1580 (
1581  const tmp<fvMatrix<Type>>& tA,
1582  const tmp<fvMatrix<Type>>& tB
1583 )
1584 {
1585  checkMethod(tA(), tB(), "+");
1586  tmp<fvMatrix<Type>> tC(tA.ptr());
1587  tC.ref() += tB();
1588  tB.clear();
1589  return tC;
1590 }
1591 
1592 template<class Type>
1593 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1594 (
1595  const fvMatrix<Type>& A,
1597 )
1598 {
1599  checkMethod(A, su, "+");
1600  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1601  tC.ref().source() -= su.mesh().V()*su.field();
1602  return tC;
1603 }
1604 
1605 template<class Type>
1606 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1607 (
1608  const fvMatrix<Type>& A,
1610 )
1611 {
1612  checkMethod(A, tsu(), "+");
1613  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1614  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1615  tsu.clear();
1616  return tC;
1617 }
1618 
1619 template<class Type>
1620 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1621 (
1622  const fvMatrix<Type>& A,
1624 )
1625 {
1626  checkMethod(A, tsu(), "+");
1627  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1628  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1629  tsu.clear();
1630  return tC;
1631 }
1632 
1633 template<class Type>
1634 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1635 (
1636  const tmp<fvMatrix<Type>>& tA,
1638 )
1639 {
1640  checkMethod(tA(), su, "+");
1641  tmp<fvMatrix<Type>> tC(tA.ptr());
1642  tC.ref().source() -= su.mesh().V()*su.field();
1643  return tC;
1644 }
1645 
1646 template<class Type>
1647 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1648 (
1649  const tmp<fvMatrix<Type>>& tA,
1651 )
1652 {
1653  checkMethod(tA(), tsu(), "+");
1654  tmp<fvMatrix<Type>> tC(tA.ptr());
1655  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1656  tsu.clear();
1657  return tC;
1658 }
1659 
1660 template<class Type>
1661 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1662 (
1663  const tmp<fvMatrix<Type>>& tA,
1665 )
1666 {
1667  checkMethod(tA(), tsu(), "+");
1668  tmp<fvMatrix<Type>> tC(tA.ptr());
1669  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1670  tsu.clear();
1671  return tC;
1672 }
1673 
1674 template<class Type>
1675 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1676 (
1678  const fvMatrix<Type>& A
1679 )
1680 {
1681  checkMethod(A, su, "+");
1682  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1683  tC.ref().source() -= su.mesh().V()*su.field();
1684  return tC;
1685 }
1686 
1687 template<class Type>
1688 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1689 (
1691  const fvMatrix<Type>& A
1692 )
1693 {
1694  checkMethod(A, tsu(), "+");
1695  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1696  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1697  tsu.clear();
1698  return tC;
1699 }
1700 
1701 template<class Type>
1702 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1703 (
1705  const fvMatrix<Type>& A
1706 )
1707 {
1708  checkMethod(A, tsu(), "+");
1709  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1710  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1711  tsu.clear();
1712  return tC;
1713 }
1714 
1715 template<class Type>
1716 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1717 (
1719  const tmp<fvMatrix<Type>>& tA
1720 )
1721 {
1722  checkMethod(tA(), su, "+");
1723  tmp<fvMatrix<Type>> tC(tA.ptr());
1724  tC.ref().source() -= su.mesh().V()*su.field();
1725  return tC;
1726 }
1727 
1728 template<class Type>
1729 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1730 (
1732  const tmp<fvMatrix<Type>>& tA
1733 )
1734 {
1735  checkMethod(tA(), tsu(), "+");
1736  tmp<fvMatrix<Type>> tC(tA.ptr());
1737  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1738  tsu.clear();
1739  return tC;
1740 }
1741 
1742 template<class Type>
1743 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1744 (
1746  const tmp<fvMatrix<Type>>& tA
1747 )
1748 {
1749  checkMethod(tA(), tsu(), "+");
1750  tmp<fvMatrix<Type>> tC(tA.ptr());
1751  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1752  tsu.clear();
1753  return tC;
1754 }
1755 
1756 
1757 template<class Type>
1758 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1759 (
1760  const fvMatrix<Type>& A,
1761  const fvMatrix<Type>& B
1762 )
1763 {
1764  checkMethod(A, B, "-");
1765  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1766  tC.ref() -= B;
1767  return tC;
1768 }
1769 
1770 template<class Type>
1771 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1772 (
1773  const tmp<fvMatrix<Type>>& tA,
1774  const fvMatrix<Type>& B
1775 )
1776 {
1777  checkMethod(tA(), B, "-");
1778  tmp<fvMatrix<Type>> tC(tA.ptr());
1779  tC.ref() -= B;
1780  return tC;
1781 }
1782 
1783 template<class Type>
1784 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1785 (
1786  const fvMatrix<Type>& A,
1787  const tmp<fvMatrix<Type>>& tB
1788 )
1789 {
1790  checkMethod(A, tB(), "-");
1791  tmp<fvMatrix<Type>> tC(tB.ptr());
1792  tC.ref() -= A;
1793  tC.ref().negate();
1794  return tC;
1795 }
1796 
1797 template<class Type>
1798 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1799 (
1800  const tmp<fvMatrix<Type>>& tA,
1801  const tmp<fvMatrix<Type>>& tB
1802 )
1803 {
1804  checkMethod(tA(), tB(), "-");
1805  tmp<fvMatrix<Type>> tC(tA.ptr());
1806  tC.ref() -= tB();
1807  tB.clear();
1808  return tC;
1809 }
1810 
1811 template<class Type>
1812 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1813 (
1814  const fvMatrix<Type>& A,
1816 )
1817 {
1818  checkMethod(A, su, "-");
1819  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1820  tC.ref().source() += su.mesh().V()*su.field();
1821  return tC;
1822 }
1823 
1824 template<class Type>
1825 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1826 (
1827  const fvMatrix<Type>& A,
1829 )
1830 {
1831  checkMethod(A, tsu(), "-");
1832  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1833  tC.ref().source() += tsu().mesh().V()*tsu().field();
1834  tsu.clear();
1835  return tC;
1836 }
1837 
1838 template<class Type>
1839 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1840 (
1841  const fvMatrix<Type>& A,
1843 )
1844 {
1845  checkMethod(A, tsu(), "-");
1846  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1847  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1848  tsu.clear();
1849  return tC;
1850 }
1851 
1852 template<class Type>
1853 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1854 (
1855  const tmp<fvMatrix<Type>>& tA,
1857 )
1858 {
1859  checkMethod(tA(), su, "-");
1860  tmp<fvMatrix<Type>> tC(tA.ptr());
1861  tC.ref().source() += su.mesh().V()*su.field();
1862  return tC;
1863 }
1864 
1865 template<class Type>
1866 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1867 (
1868  const tmp<fvMatrix<Type>>& tA,
1870 )
1871 {
1872  checkMethod(tA(), tsu(), "-");
1873  tmp<fvMatrix<Type>> tC(tA.ptr());
1874  tC.ref().source() += tsu().mesh().V()*tsu().field();
1875  tsu.clear();
1876  return tC;
1877 }
1878 
1879 template<class Type>
1880 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1881 (
1882  const tmp<fvMatrix<Type>>& tA,
1884 )
1885 {
1886  checkMethod(tA(), tsu(), "-");
1887  tmp<fvMatrix<Type>> tC(tA.ptr());
1888  tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
1889  tsu.clear();
1890  return tC;
1891 }
1892 
1893 template<class Type>
1894 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1895 (
1897  const fvMatrix<Type>& A
1898 )
1899 {
1900  checkMethod(A, su, "-");
1901  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1902  tC.ref().negate();
1903  tC.ref().source() -= su.mesh().V()*su.field();
1904  return tC;
1905 }
1906 
1907 template<class Type>
1908 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1909 (
1911  const fvMatrix<Type>& A
1912 )
1913 {
1914  checkMethod(A, tsu(), "-");
1915  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1916  tC.ref().negate();
1917  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1918  tsu.clear();
1919  return tC;
1920 }
1921 
1922 template<class Type>
1923 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1924 (
1926  const fvMatrix<Type>& A
1927 )
1928 {
1929  checkMethod(A, tsu(), "-");
1930  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1931  tC.ref().negate();
1932  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1933  tsu.clear();
1934  return tC;
1935 }
1936 
1937 template<class Type>
1938 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1939 (
1941  const tmp<fvMatrix<Type>>& tA
1942 )
1943 {
1944  checkMethod(tA(), su, "-");
1945  tmp<fvMatrix<Type>> tC(tA.ptr());
1946  tC.ref().negate();
1947  tC.ref().source() -= su.mesh().V()*su.field();
1948  return tC;
1949 }
1950 
1951 template<class Type>
1952 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1953 (
1955  const tmp<fvMatrix<Type>>& tA
1956 )
1957 {
1958  checkMethod(tA(), tsu(), "-");
1959  tmp<fvMatrix<Type>> tC(tA.ptr());
1960  tC.ref().negate();
1961  tC.ref().source() -= tsu().mesh().V()*tsu().field();
1962  tsu.clear();
1963  return tC;
1964 }
1965 
1966 template<class Type>
1967 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
1968 (
1970  const tmp<fvMatrix<Type>>& tA
1971 )
1972 {
1973  checkMethod(tA(), tsu(), "-");
1974  tmp<fvMatrix<Type>> tC(tA.ptr());
1975  tC.ref().negate();
1976  tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
1977  tsu.clear();
1978  return tC;
1979 }
1980 
1981 template<class Type>
1982 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1983 (
1984  const fvMatrix<Type>& A,
1985  const dimensioned<Type>& su
1986 )
1987 {
1988  checkMethod(A, su, "+");
1989  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
1990  tC.ref().source() -= su.value()*A.psi().mesh().V();
1991  return tC;
1992 }
1993 
1994 template<class Type>
1995 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
1996 (
1997  const tmp<fvMatrix<Type>>& tA,
1998  const dimensioned<Type>& su
1999 )
2000 {
2001  checkMethod(tA(), su, "+");
2002  tmp<fvMatrix<Type>> tC(tA.ptr());
2003  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2004  return tC;
2005 }
2006 
2007 template<class Type>
2008 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2009 (
2010  const dimensioned<Type>& su,
2011  const fvMatrix<Type>& A
2012 )
2013 {
2014  checkMethod(A, su, "+");
2015  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2016  tC.ref().source() -= su.value()*A.psi().mesh().V();
2017  return tC;
2018 }
2019 
2020 template<class Type>
2021 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2022 (
2023  const dimensioned<Type>& su,
2024  const tmp<fvMatrix<Type>>& tA
2025 )
2026 {
2027  checkMethod(tA(), su, "+");
2028  tmp<fvMatrix<Type>> tC(tA.ptr());
2029  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2030  return tC;
2031 }
2032 
2033 template<class Type>
2034 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2035 (
2036  const fvMatrix<Type>& A,
2037  const dimensioned<Type>& su
2038 )
2039 {
2040  checkMethod(A, su, "-");
2041  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2042  tC.ref().source() += su.value()*tC().psi().mesh().V();
2043  return tC;
2044 }
2045 
2046 template<class Type>
2047 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2048 (
2049  const tmp<fvMatrix<Type>>& tA,
2050  const dimensioned<Type>& su
2051 )
2052 {
2053  checkMethod(tA(), su, "-");
2054  tmp<fvMatrix<Type>> tC(tA.ptr());
2055  tC.ref().source() += su.value()*tC().psi().mesh().V();
2056  return tC;
2057 }
2058 
2059 template<class Type>
2060 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2061 (
2062  const dimensioned<Type>& su,
2063  const fvMatrix<Type>& A
2064 )
2065 {
2066  checkMethod(A, su, "-");
2067  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2068  tC.ref().negate();
2069  tC.ref().source() -= su.value()*A.psi().mesh().V();
2070  return tC;
2071 }
2072 
2073 template<class Type>
2074 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2075 (
2076  const dimensioned<Type>& su,
2077  const tmp<fvMatrix<Type>>& tA
2078 )
2079 {
2080  checkMethod(tA(), su, "-");
2081  tmp<fvMatrix<Type>> tC(tA.ptr());
2082  tC.ref().negate();
2083  tC.ref().source() -= su.value()*tC().psi().mesh().V();
2084  return tC;
2085 }
2086 
2087 
2088 template<class Type>
2089 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2090 (
2091  const volScalarField::Internal& dsf,
2092  const fvMatrix<Type>& A
2093 )
2094 {
2095  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2096  tC.ref() *= dsf;
2097  return tC;
2098 }
2099 
2100 template<class Type>
2101 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2102 (
2103  const tmp<volScalarField::Internal>& tdsf,
2104  const fvMatrix<Type>& A
2105 )
2106 {
2107  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2108  tC.ref() *= tdsf;
2109  return tC;
2110 }
2111 
2112 template<class Type>
2113 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2114 (
2115  const tmp<volScalarField>& tvsf,
2116  const fvMatrix<Type>& A
2117 )
2118 {
2119  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2120  tC.ref() *= tvsf;
2121  return tC;
2122 }
2123 
2124 template<class Type>
2125 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2126 (
2127  const volScalarField::Internal& dsf,
2128  const tmp<fvMatrix<Type>>& tA
2129 )
2130 {
2131  tmp<fvMatrix<Type>> tC(tA.ptr());
2132  tC.ref() *= dsf;
2133  return tC;
2134 }
2135 
2136 template<class Type>
2137 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2138 (
2139  const tmp<volScalarField::Internal>& tdsf,
2140  const tmp<fvMatrix<Type>>& tA
2141 )
2142 {
2143  tmp<fvMatrix<Type>> tC(tA.ptr());
2144  tC.ref() *= tdsf;
2145  return tC;
2146 }
2147 
2148 template<class Type>
2149 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2150 (
2151  const tmp<volScalarField>& tvsf,
2152  const tmp<fvMatrix<Type>>& tA
2153 )
2154 {
2155  tmp<fvMatrix<Type>> tC(tA.ptr());
2156  tC.ref() *= tvsf;
2157  return tC;
2158 }
2159 
2160 template<class Type>
2161 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2162 (
2163  const dimensioned<scalar>& ds,
2164  const fvMatrix<Type>& A
2165 )
2166 {
2167  tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2168  tC.ref() *= ds;
2169  return tC;
2170 }
2171 
2172 template<class Type>
2173 Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2174 (
2175  const dimensioned<scalar>& ds,
2176  const tmp<fvMatrix<Type>>& tA
2177 )
2178 {
2179  tmp<fvMatrix<Type>> tC(tA.ptr());
2180  tC.ref() *= ds;
2181  return tC;
2182 }
2183 
2184 
2185 template<class Type>
2187 Foam::operator&
2188 (
2189  const fvMatrix<Type>& M,
2191 )
2192 {
2194  (
2196  (
2197  "M&" + psi.name(),
2198  psi.mesh(),
2199  M.dimensions()/dimVol,
2200  extrapolatedCalculatedFvPatchScalarField::typeName
2201  )
2202  );
2204 
2205  // Loop over field components
2206  if (M.hasDiag())
2207  {
2208  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2209  {
2210  scalarField psiCmpt(psi.field().component(cmpt));
2211  scalarField boundaryDiagCmpt(M.diag());
2212  M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2213  Mphi.primitiveFieldRef().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2214  }
2215  }
2216  else
2217  {
2218  Mphi.primitiveFieldRef() = Zero;
2219  }
2220 
2221  Mphi.primitiveFieldRef() += M.lduMatrix::H(psi.field()) + M.source();
2222  M.addBoundarySource(Mphi.primitiveFieldRef());
2223 
2224  Mphi.primitiveFieldRef() /= -psi.mesh().V();
2226 
2227  return tMphi;
2228 }
2229 
2230 template<class Type>
2232 Foam::operator&
2233 (
2234  const fvMatrix<Type>& M,
2236 )
2237 {
2239  tpsi.clear();
2240  return tMpsi;
2241 }
2242 
2243 template<class Type>
2245 Foam::operator&
2246 (
2247  const fvMatrix<Type>& M,
2249 )
2250 {
2252  tpsi.clear();
2253  return tMpsi;
2254 }
2255 
2256 template<class Type>
2258 Foam::operator&
2259 (
2260  const tmp<fvMatrix<Type>>& tM,
2262 )
2263 {
2265  tM.clear();
2266  return tMpsi;
2267 }
2268 
2269 template<class Type>
2271 Foam::operator&
2272 (
2273  const tmp<fvMatrix<Type>>& tM,
2275 )
2276 {
2277  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2278  tM.clear();
2279  tpsi.clear();
2280  return tMpsi;
2281 }
2282 
2283 template<class Type>
2285 Foam::operator&
2286 (
2287  const tmp<fvMatrix<Type>>& tM,
2289 )
2290 {
2291  tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2292  tM.clear();
2293  tpsi.clear();
2294  return tMpsi;
2295 }
2296 
2297 
2298 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2299 
2300 template<class Type>
2301 Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
2302 {
2303  os << static_cast<const lduMatrix&>(fvm) << nl
2304  << fvm.dimensions_ << nl
2305  << fvm.source_ << nl
2306  << fvm.internalCoeffs_ << nl
2307  << fvm.boundaryCoeffs_ << endl;
2308 
2309  os.check("Ostream& operator<<(Ostream&, fvMatrix<Type>&");
2310 
2311  return os;
2312 }
2313 
2314 
2315 // * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
2316 
2317 #include "fvMatrixSolve.C"
2318 
2319 // ************************************************************************* //
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:230
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:319
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:174
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
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 dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
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
const dimensionedScalar & c
Speed of light in a vacuum.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
tmp< scalarField > H1() const
void sumMagOffDiag(scalarField &sumOff) const
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 dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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
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:412
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:484
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:313
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
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.
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:1237
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 &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:198
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:299
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