fvMatrix.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 // * * * * * * * * * * * * * Private 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  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_).internalField();
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] =
253 
254  boundaryCoeffs_[patchi][patchFacei] =
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(), pTraits<Type>::zero),
277  internalCoeffs_(psi.mesh().boundary().size()),
278  boundaryCoeffs_(psi.mesh().boundary().size()),
279  faceFluxCorrectionPtr_(NULL)
280 {
281  if (debug)
282  {
283  Info<< "fvMatrix<Type>(GeometricField<Type, fvPatchField, volMesh>&,"
284  " const dimensionSet&) : "
285  "constructing fvMatrix<Type> for field " << psi_.name()
286  << endl;
287  }
288 
289  // Initialise coupling coefficients
290  forAll(psi.mesh().boundary(), patchI)
291  {
292  internalCoeffs_.set
293  (
294  patchI,
295  new Field<Type>
296  (
297  psi.mesh().boundary()[patchI].size(),
298  pTraits<Type>::zero
299  )
300  );
301 
302  boundaryCoeffs_.set
303  (
304  patchI,
305  new Field<Type>
306  (
307  psi.mesh().boundary()[patchI].size(),
308  pTraits<Type>::zero
309  )
310  );
311  }
312 
313  // Update the boundary coefficients of psi without changing its event No.
314  GeometricField<Type, fvPatchField, volMesh>& psiRef =
315  const_cast<GeometricField<Type, fvPatchField, volMesh>&>(psi_);
316 
317  label currentStatePsi = psiRef.eventNo();
318  psiRef.boundaryField().updateCoeffs();
319  psiRef.eventNo() = currentStatePsi;
320 }
321 
322 
323 template<class Type>
324 Foam::fvMatrix<Type>::fvMatrix(const fvMatrix<Type>& fvm)
325 :
326  refCount(),
327  lduMatrix(fvm),
328  psi_(fvm.psi_),
329  dimensions_(fvm.dimensions_),
330  source_(fvm.source_),
331  internalCoeffs_(fvm.internalCoeffs_),
332  boundaryCoeffs_(fvm.boundaryCoeffs_),
333  faceFluxCorrectionPtr_(NULL)
334 {
335  if (debug)
336  {
337  Info<< "fvMatrix<Type>::fvMatrix(const fvMatrix<Type>&) : "
338  << "copying fvMatrix<Type> for field " << psi_.name()
339  << endl;
340  }
341 
342  if (fvm.faceFluxCorrectionPtr_)
343  {
344  faceFluxCorrectionPtr_ = new
345  GeometricField<Type, fvsPatchField, surfaceMesh>
346  (
347  *(fvm.faceFluxCorrectionPtr_)
348  );
349  }
350 }
351 
352 
353 #ifndef NoConstructFromTmp
354 template<class Type>
355 Foam::fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type> >& tfvm)
356 :
357  refCount(),
358  lduMatrix
359  (
360  const_cast<fvMatrix<Type>&>(tfvm()),
361  tfvm.isTmp()
362  ),
363  psi_(tfvm().psi_),
364  dimensions_(tfvm().dimensions_),
365  source_
366  (
367  const_cast<fvMatrix<Type>&>(tfvm()).source_,
368  tfvm.isTmp()
369  ),
370  internalCoeffs_
371  (
372  const_cast<fvMatrix<Type>&>(tfvm()).internalCoeffs_,
373  tfvm.isTmp()
374  ),
375  boundaryCoeffs_
376  (
377  const_cast<fvMatrix<Type>&>(tfvm()).boundaryCoeffs_,
378  tfvm.isTmp()
379  ),
380  faceFluxCorrectionPtr_(NULL)
381 {
382  if (debug)
383  {
384  Info<< "fvMatrix<Type>::fvMatrix(const tmp<fvMatrix<Type> >&) : "
385  << "copying fvMatrix<Type> for field " << psi_.name()
386  << endl;
387  }
388 
389  if (tfvm().faceFluxCorrectionPtr_)
390  {
391  if (tfvm.isTmp())
392  {
393  faceFluxCorrectionPtr_ = tfvm().faceFluxCorrectionPtr_;
394  tfvm().faceFluxCorrectionPtr_ = NULL;
395  }
396  else
397  {
398  faceFluxCorrectionPtr_ = new
399  GeometricField<Type, fvsPatchField, surfaceMesh>
400  (
401  *(tfvm().faceFluxCorrectionPtr_)
402  );
403  }
404  }
405 
406  tfvm.clear();
407 }
408 #endif
409 
410 
411 template<class Type>
413 (
414  const GeometricField<Type, fvPatchField, volMesh>& psi,
415  Istream& is
416 )
417 :
418  lduMatrix(psi.mesh()),
419  psi_(psi),
420  dimensions_(is),
421  source_(is),
422  internalCoeffs_(psi.mesh().boundary().size()),
423  boundaryCoeffs_(psi.mesh().boundary().size()),
424  faceFluxCorrectionPtr_(NULL)
425 {
426  if (debug)
427  {
428  Info<< "fvMatrix<Type>"
429  "(GeometricField<Type, fvPatchField, volMesh>&, Istream&) : "
430  "constructing fvMatrix<Type> for field " << psi_.name()
431  << endl;
432  }
433 
434  // Initialise coupling coefficients
435  forAll(psi.mesh().boundary(), patchI)
436  {
437  internalCoeffs_.set
438  (
439  patchI,
440  new Field<Type>
441  (
442  psi.mesh().boundary()[patchI].size(),
443  pTraits<Type>::zero
444  )
445  );
446 
447  boundaryCoeffs_.set
448  (
449  patchI,
450  new Field<Type>
451  (
452  psi.mesh().boundary()[patchI].size(),
453  pTraits<Type>::zero
454  )
455  );
456  }
457 
458 }
459 
460 
461 template<class Type>
463 {
464  if (debug)
465  {
466  Info<< "fvMatrix<Type>::~fvMatrix<Type>() : "
467  << "destroying fvMatrix<Type> for field " << psi_.name()
468  << endl;
469  }
470 
471  if (faceFluxCorrectionPtr_)
472  {
473  delete faceFluxCorrectionPtr_;
474  }
475 }
476 
477 
478 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
479 
480 template<class Type>
482 (
483  const labelUList& cellLabels,
484  const UList<Type>& values
485 )
486 {
487  this->setValuesFromList(cellLabels, values);
488 }
489 
490 
491 template<class Type>
493 (
494  const labelUList& cellLabels,
495  const UIndirectList<Type>& values
496 )
497 {
498  this->setValuesFromList(cellLabels, values);
499 }
500 
501 
502 template<class Type>
504 (
505  const label celli,
506  const Type& value,
507  const bool forceReference
508 )
509 {
510  if ((forceReference || psi_.needReference()) && celli >= 0)
511  {
512  source()[celli] += diag()[celli]*value;
513  diag()[celli] += diag()[celli];
514  }
515 }
516 
517 
518 template<class Type>
520 {
521  if (alpha <= 0)
522  {
523  return;
524  }
525 
526  if (debug)
527  {
529  << "Relaxing " << psi_.name() << " by " << alpha
530  << endl;
531  }
532 
533  Field<Type>& S = source();
534  scalarField& D = diag();
535 
536  // Store the current unrelaxed diagonal for use in updating the source
537  scalarField D0(D);
538 
539  // Calculate the sum-mag off-diagonal from the interior faces
540  scalarField sumOff(D.size(), 0.0);
541  sumMagOffDiag(sumOff);
542 
543  // Handle the boundary contributions to the diagonal
544  forAll(psi_.boundaryField(), patchI)
545  {
546  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
547 
548  if (ptf.size())
549  {
550  const labelUList& pa = lduAddr().patchAddr(patchI);
551  Field<Type>& iCoeffs = internalCoeffs_[patchI];
552 
553  if (ptf.coupled())
554  {
555  const Field<Type>& pCoeffs = boundaryCoeffs_[patchI];
556 
557  // For coupled boundaries add the diagonal and
558  // off-diagonal contributions
559  forAll(pa, face)
560  {
561  D[pa[face]] += component(iCoeffs[face], 0);
562  sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
563  }
564  }
565  else
566  {
567  // For non-coupled boundaries add the maximum magnitude diagonal
568  // contribution to ensure stability
569  forAll(pa, face)
570  {
571  D[pa[face]] += cmptMax(cmptMag(iCoeffs[face]));
572  }
573  }
574  }
575  }
576 
577 
578  if (debug)
579  {
580  // Calculate amount of non-dominance.
581  label nNon = 0;
582  scalar maxNon = 0.0;
583  scalar sumNon = 0.0;
584  forAll(D, celli)
585  {
586  scalar d = (sumOff[celli] - D[celli])/mag(D[celli]);
587 
588  if (d > 0)
589  {
590  nNon++;
591  maxNon = max(maxNon, d);
592  sumNon += d;
593  }
594  }
595 
596  reduce(nNon, sumOp<label>(), UPstream::msgType(), psi_.mesh().comm());
597  reduce
598  (
599  maxNon,
600  maxOp<scalar>(),
601  UPstream::msgType(),
602  psi_.mesh().comm()
603  );
604  reduce
605  (
606  sumNon,
607  sumOp<scalar>(),
608  UPstream::msgType(),
609  psi_.mesh().comm()
610  );
611  sumNon /= returnReduce
612  (
613  D.size(),
614  sumOp<label>(),
615  UPstream::msgType(),
616  psi_.mesh().comm()
617  );
618 
620  << "Matrix dominance test for " << psi_.name() << nl
621  << " number of non-dominant cells : " << nNon << nl
622  << " maximum relative non-dominance : " << maxNon << nl
623  << " average relative non-dominance : " << sumNon << nl
624  << endl;
625  }
626 
627 
628  // Ensure the matrix is diagonally dominant...
629  // Assumes that the central coefficient is positive and ensures it is
630  forAll(D, celli)
631  {
632  D[celli] = max(mag(D[celli]), sumOff[celli]);
633  }
634 
635  // ... then relax
636  D /= alpha;
637 
638  // Now remove the diagonal contribution from coupled boundaries
639  forAll(psi_.boundaryField(), patchI)
640  {
641  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
642 
643  if (ptf.size())
644  {
645  const labelUList& pa = lduAddr().patchAddr(patchI);
646  Field<Type>& iCoeffs = internalCoeffs_[patchI];
647 
648  if (ptf.coupled())
649  {
650  forAll(pa, face)
651  {
652  D[pa[face]] -= component(iCoeffs[face], 0);
653  }
654  }
655  else
656  {
657  forAll(pa, face)
658  {
659  D[pa[face]] -= cmptMin(iCoeffs[face]);
660  }
661  }
662  }
663  }
664 
665  // Finally add the relaxation contribution to the source.
666  S += (D - D0)*psi_.internalField();
667 }
668 
669 
670 template<class Type>
672 {
673  word name = psi_.select
674  (
675  psi_.mesh().data::template lookupOrDefault<bool>
676  ("finalIteration", false)
677  );
678 
679  if (psi_.mesh().relaxEquation(name))
680  {
681  relax(psi_.mesh().equationRelaxationFactor(name));
682  }
683 }
684 
685 
686 template<class Type>
688 (
690  GeometricBoundaryField& bFields
691 )
692 {
693  forAll(bFields, patchI)
694  {
695  bFields[patchI].manipulateMatrix(*this);
696  }
697 }
698 
699 
700 template<class Type>
702 {
703  tmp<scalarField> tdiag(new scalarField(diag()));
704  addCmptAvBoundaryDiag(tdiag());
705  return tdiag;
706 }
707 
708 
709 template<class Type>
711 {
713 
714  forAll(psi_.boundaryField(), patchI)
715  {
716  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
717 
718  if (!ptf.coupled() && ptf.size())
719  {
720  addToInternalField
721  (
722  lduAddr().patchAddr(patchI),
723  internalCoeffs_[patchI],
724  tdiag()
725  );
726  }
727  }
728 
729  return tdiag;
730 }
731 
732 
733 template<class Type>
735 {
736  tmp<volScalarField> tAphi
737  (
738  new volScalarField
739  (
740  IOobject
741  (
742  "A("+psi_.name()+')',
743  psi_.instance(),
744  psi_.mesh(),
745  IOobject::NO_READ,
746  IOobject::NO_WRITE
747  ),
748  psi_.mesh(),
749  dimensions_/psi_.dimensions()/dimVol,
750  zeroGradientFvPatchScalarField::typeName
751  )
752  );
753 
754  tAphi().internalField() = D()/psi_.mesh().V();
755  tAphi().correctBoundaryConditions();
756 
757  return tAphi;
758 }
759 
760 
761 template<class Type>
764 {
766  (
768  (
769  IOobject
770  (
771  "H("+psi_.name()+')',
772  psi_.instance(),
773  psi_.mesh(),
774  IOobject::NO_READ,
775  IOobject::NO_WRITE
776  ),
777  psi_.mesh(),
778  dimensions_/dimVol,
779  zeroGradientFvPatchScalarField::typeName
780  )
781  );
783 
784  // Loop over field components
785  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
786  {
787  scalarField psiCmpt(psi_.internalField().component(cmpt));
788 
789  scalarField boundaryDiagCmpt(psi_.size(), 0.0);
790  addBoundaryDiag(boundaryDiagCmpt, cmpt);
791  boundaryDiagCmpt.negate();
792  addCmptAvBoundaryDiag(boundaryDiagCmpt);
793 
794  Hphi.internalField().replace(cmpt, boundaryDiagCmpt*psiCmpt);
795  }
796 
797  Hphi.internalField() += lduMatrix::H(psi_.internalField()) + source_;
798  addBoundarySource(Hphi.internalField());
799 
800  Hphi.internalField() /= psi_.mesh().V();
802 
803  typename Type::labelType validComponents
804  (
805  psi_.mesh().template validComponents<Type>()
806  );
807 
808  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
809  {
810  if (validComponents[cmpt] == -1)
811  {
812  Hphi.replace
813  (
814  cmpt,
815  dimensionedScalar("0", Hphi.dimensions(), 0.0)
816  );
817  }
818  }
819 
820  return tHphi;
821 }
822 
823 
824 template<class Type>
826 {
828  (
829  new volScalarField
830  (
831  IOobject
832  (
833  "H(1)",
834  psi_.instance(),
835  psi_.mesh(),
836  IOobject::NO_READ,
837  IOobject::NO_WRITE
838  ),
839  psi_.mesh(),
840  dimensions_/(dimVol*psi_.dimensions()),
841  zeroGradientFvPatchScalarField::typeName
842  )
843  );
844  volScalarField& H1_ = tH1();
845 
846  H1_.internalField() = lduMatrix::H1();
847 
848  forAll(psi_.boundaryField(), patchI)
849  {
850  const fvPatchField<Type>& ptf = psi_.boundaryField()[patchI];
851 
852  if (ptf.coupled() && ptf.size())
853  {
854  addToInternalField
855  (
856  lduAddr().patchAddr(patchI),
857  boundaryCoeffs_[patchI].component(0),
858  H1_
859  );
860  }
861  }
862 
863  H1_.internalField() /= psi_.mesh().V();
865 
866  return tH1;
867 }
868 
869 
870 template<class Type>
873 flux() const
874 {
875  if (!psi_.mesh().fluxRequired(psi_.name()))
876  {
878  << "flux requested but " << psi_.name()
879  << " not specified in the fluxRequired sub-dictionary"
880  " of fvSchemes."
881  << abort(FatalError);
882  }
883 
884  // construct GeometricField<Type, fvsPatchField, surfaceMesh>
886  (
888  (
889  IOobject
890  (
891  "flux("+psi_.name()+')',
892  psi_.instance(),
893  psi_.mesh(),
894  IOobject::NO_READ,
895  IOobject::NO_WRITE
896  ),
897  psi_.mesh(),
898  dimensions()
899  )
900  );
901  GeometricField<Type, fvsPatchField, surfaceMesh>& fieldFlux = tfieldFlux();
902 
903  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
904  {
905  fieldFlux.internalField().replace
906  (
907  cmpt,
908  lduMatrix::faceH(psi_.internalField().component(cmpt))
909  );
910  }
911 
912  FieldField<Field, Type> InternalContrib = internalCoeffs_;
913 
914  forAll(InternalContrib, patchI)
915  {
916  InternalContrib[patchI] =
918  (
919  InternalContrib[patchI],
920  psi_.boundaryField()[patchI].patchInternalField()
921  );
922  }
923 
924  FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
925 
926  forAll(NeighbourContrib, patchI)
927  {
928  if (psi_.boundaryField()[patchI].coupled())
929  {
930  NeighbourContrib[patchI] =
932  (
933  NeighbourContrib[patchI],
934  psi_.boundaryField()[patchI].patchNeighbourField()
935  );
936  }
937  }
938 
939  forAll(fieldFlux.boundaryField(), patchI)
940  {
941  fieldFlux.boundaryField()[patchI] =
942  InternalContrib[patchI] - NeighbourContrib[patchI];
943  }
944 
945  if (faceFluxCorrectionPtr_)
946  {
947  fieldFlux += *faceFluxCorrectionPtr_;
948  }
949 
950  return tfieldFlux;
951 }
952 
953 
954 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
955 
956 template<class Type>
958 {
959  if (this == &fvmv)
960  {
962  << "attempted assignment to self"
963  << abort(FatalError);
964  }
965 
966  if (&psi_ != &(fvmv.psi_))
967  {
969  << "different fields"
970  << abort(FatalError);
971  }
972 
973  dimensions_ = fvmv.dimensions_;
974  lduMatrix::operator=(fvmv);
975  source_ = fvmv.source_;
976  internalCoeffs_ = fvmv.internalCoeffs_;
977  boundaryCoeffs_ = fvmv.boundaryCoeffs_;
978 
979  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
980  {
981  *faceFluxCorrectionPtr_ = *fvmv.faceFluxCorrectionPtr_;
982  }
983  else if (fvmv.faceFluxCorrectionPtr_)
984  {
985  faceFluxCorrectionPtr_ =
986  new GeometricField<Type, fvsPatchField, surfaceMesh>
987  (*fvmv.faceFluxCorrectionPtr_);
988  }
989 }
990 
991 
992 template<class Type>
993 void Foam::fvMatrix<Type>::operator=(const tmp<fvMatrix<Type> >& tfvmv)
994 {
995  operator=(tfvmv());
996  tfvmv.clear();
997 }
998 
999 
1000 template<class Type>
1002 {
1004  source_.negate();
1005  internalCoeffs_.negate();
1006  boundaryCoeffs_.negate();
1007 
1008  if (faceFluxCorrectionPtr_)
1009  {
1010  faceFluxCorrectionPtr_->negate();
1011  }
1012 }
1013 
1014 
1015 template<class Type>
1017 {
1018  checkMethod(*this, fvmv, "+=");
1019 
1020  dimensions_ += fvmv.dimensions_;
1021  lduMatrix::operator+=(fvmv);
1022  source_ += fvmv.source_;
1023  internalCoeffs_ += fvmv.internalCoeffs_;
1024  boundaryCoeffs_ += fvmv.boundaryCoeffs_;
1025 
1026  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1027  {
1028  *faceFluxCorrectionPtr_ += *fvmv.faceFluxCorrectionPtr_;
1029  }
1030  else if (fvmv.faceFluxCorrectionPtr_)
1031  {
1032  faceFluxCorrectionPtr_ = new
1033  GeometricField<Type, fvsPatchField, surfaceMesh>
1034  (
1036  );
1037  }
1038 }
1039 
1040 
1041 template<class Type>
1042 void Foam::fvMatrix<Type>::operator+=(const tmp<fvMatrix<Type> >& tfvmv)
1043 {
1044  operator+=(tfvmv());
1045  tfvmv.clear();
1046 }
1047 
1048 
1049 template<class Type>
1050 void Foam::fvMatrix<Type>::operator-=(const fvMatrix<Type>& fvmv)
1051 {
1052  checkMethod(*this, fvmv, "-=");
1053 
1054  dimensions_ -= fvmv.dimensions_;
1055  lduMatrix::operator-=(fvmv);
1056  source_ -= fvmv.source_;
1057  internalCoeffs_ -= fvmv.internalCoeffs_;
1058  boundaryCoeffs_ -= fvmv.boundaryCoeffs_;
1059 
1060  if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1061  {
1062  *faceFluxCorrectionPtr_ -= *fvmv.faceFluxCorrectionPtr_;
1063  }
1064  else if (fvmv.faceFluxCorrectionPtr_)
1065  {
1066  faceFluxCorrectionPtr_ =
1067  new GeometricField<Type, fvsPatchField, surfaceMesh>
1068  (-*fvmv.faceFluxCorrectionPtr_);
1069  }
1070 }
1071 
1072 
1073 template<class Type>
1074 void Foam::fvMatrix<Type>::operator-=(const tmp<fvMatrix<Type> >& tfvmv)
1075 {
1076  operator-=(tfvmv());
1077  tfvmv.clear();
1078 }
1079 
1080 
1081 template<class Type>
1083 (
1084  const DimensionedField<Type, volMesh>& su
1085 )
1086 {
1087  checkMethod(*this, su, "+=");
1088  source() -= su.mesh().V()*su.field();
1089 }
1090 
1091 
1092 template<class Type>
1094 (
1095  const tmp<DimensionedField<Type, volMesh> >& tsu
1096 )
1097 {
1098  operator+=(tsu());
1099  tsu.clear();
1100 }
1101 
1102 
1103 template<class Type>
1105 (
1106  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1107 )
1108 {
1109  operator+=(tsu());
1110  tsu.clear();
1111 }
1112 
1113 
1114 template<class Type>
1116 (
1117  const DimensionedField<Type, volMesh>& su
1118 )
1119 {
1120  checkMethod(*this, su, "-=");
1121  source() += su.mesh().V()*su.field();
1122 }
1123 
1124 
1125 template<class Type>
1127 (
1128  const tmp<DimensionedField<Type, volMesh> >& tsu
1129 )
1130 {
1131  operator-=(tsu());
1132  tsu.clear();
1133 }
1134 
1135 
1136 template<class Type>
1138 (
1139  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1140 )
1141 {
1142  operator-=(tsu());
1143  tsu.clear();
1144 }
1145 
1146 
1147 template<class Type>
1149 (
1150  const dimensioned<Type>& su
1151 )
1152 {
1153  source() -= psi().mesh().V()*su;
1154 }
1155 
1156 
1157 template<class Type>
1159 (
1160  const dimensioned<Type>& su
1161 )
1162 {
1163  source() += psi().mesh().V()*su;
1164 }
1165 
1166 
1167 template<class Type>
1170  const zero&
1171 )
1172 {}
1173 
1174 
1175 template<class Type>
1178  const zero&
1179 )
1180 {}
1181 
1182 
1183 template<class Type>
1185 (
1187 )
1188 {
1189  dimensions_ *= dsf.dimensions();
1190  lduMatrix::operator*=(dsf.field());
1191  source_ *= dsf.field();
1192 
1193  forAll(boundaryCoeffs_, patchI)
1194  {
1195  scalarField pisf
1196  (
1197  dsf.mesh().boundary()[patchI].patchInternalField(dsf.field())
1198  );
1199 
1200  internalCoeffs_[patchI] *= pisf;
1201  boundaryCoeffs_[patchI] *= pisf;
1202  }
1203 
1204  if (faceFluxCorrectionPtr_)
1205  {
1207  << "cannot scale a matrix containing a faceFluxCorrection"
1208  << abort(FatalError);
1209  }
1210 }
1211 
1212 
1213 template<class Type>
1215 (
1216  const tmp<DimensionedField<scalar, volMesh> >& tdsf
1217 )
1218 {
1219  operator*=(tdsf());
1220  tdsf.clear();
1221 }
1222 
1223 
1224 template<class Type>
1226 (
1227  const tmp<volScalarField>& tvsf
1228 )
1229 {
1230  operator*=(tvsf());
1231  tvsf.clear();
1232 }
1233 
1234 
1235 template<class Type>
1237 (
1238  const dimensioned<scalar>& ds
1239 )
1240 {
1241  dimensions_ *= ds.dimensions();
1242  lduMatrix::operator*=(ds.value());
1243  source_ *= ds.value();
1244  internalCoeffs_ *= ds.value();
1245  boundaryCoeffs_ *= ds.value();
1246 
1247  if (faceFluxCorrectionPtr_)
1248  {
1249  *faceFluxCorrectionPtr_ *= ds.value();
1250  }
1251 }
1252 
1253 
1254 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
1255 
1256 template<class Type>
1257 void Foam::checkMethod
1259  const fvMatrix<Type>& fvm1,
1260  const fvMatrix<Type>& fvm2,
1261  const char* op
1262 )
1263 {
1264  if (&fvm1.psi() != &fvm2.psi())
1265  {
1267  << "incompatible fields for operation "
1268  << endl << " "
1269  << "[" << fvm1.psi().name() << "] "
1270  << op
1271  << " [" << fvm2.psi().name() << "]"
1272  << abort(FatalError);
1273  }
1274 
1275  if (dimensionSet::debug && fvm1.dimensions() != fvm2.dimensions())
1276  {
1278  << "incompatible dimensions for operation "
1279  << endl << " "
1280  << "[" << fvm1.psi().name() << fvm1.dimensions()/dimVolume << " ] "
1281  << op
1282  << " [" << fvm2.psi().name() << fvm2.dimensions()/dimVolume << " ]"
1283  << abort(FatalError);
1284  }
1285 }
1286 
1287 
1288 template<class Type>
1289 void Foam::checkMethod
1291  const fvMatrix<Type>& fvm,
1293  const char* op
1294 )
1295 {
1296  if (dimensionSet::debug && fvm.dimensions()/dimVolume != df.dimensions())
1297  {
1299  << endl << " "
1300  << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
1301  << op
1302  << " [" << df.name() << df.dimensions() << " ]"
1303  << abort(FatalError);
1304  }
1305 }
1306 
1307 
1308 template<class Type>
1309 void Foam::checkMethod
1311  const fvMatrix<Type>& fvm,
1312  const dimensioned<Type>& dt,
1313  const char* op
1314 )
1315 {
1316  if (dimensionSet::debug && fvm.dimensions()/dimVolume != dt.dimensions())
1317  {
1319  << "incompatible dimensions for operation "
1320  << endl << " "
1321  << "[" << fvm.psi().name() << fvm.dimensions()/dimVolume << " ] "
1322  << op
1323  << " [" << dt.name() << dt.dimensions() << " ]"
1324  << abort(FatalError);
1325  }
1326 }
1327 
1328 
1329 template<class Type>
1331 (
1332  fvMatrix<Type>& fvm,
1333  const dictionary& solverControls
1334 )
1335 {
1336  return fvm.solve(solverControls);
1337 }
1338 
1339 template<class Type>
1341 (
1342  const tmp<fvMatrix<Type> >& tfvm,
1343  const dictionary& solverControls
1344 )
1345 {
1346  SolverPerformance<Type> solverPerf =
1347  const_cast<fvMatrix<Type>&>(tfvm()).solve(solverControls);
1348 
1349  tfvm.clear();
1350 
1351  return solverPerf;
1352 }
1353 
1354 
1355 template<class Type>
1356 Foam::SolverPerformance<Type> Foam::solve(fvMatrix<Type>& fvm)
1357 {
1358  return fvm.solve();
1359 }
1360 
1361 template<class Type>
1362 Foam::SolverPerformance<Type> Foam::solve(const tmp<fvMatrix<Type> >& tfvm)
1363 {
1364  SolverPerformance<Type> solverPerf =
1365  const_cast<fvMatrix<Type>&>(tfvm()).solve();
1366 
1367  tfvm.clear();
1368 
1369  return solverPerf;
1370 }
1371 
1372 
1373 template<class Type>
1375 (
1376  const fvMatrix<Type>& A
1377 )
1378 {
1379  tmp<Foam::fvMatrix<Type> > tAcorr = A - (A & A.psi());
1380 
1381  if
1382  (
1383  (A.hasUpper() || A.hasLower())
1384  && A.psi().mesh().fluxRequired(A.psi().name())
1385  )
1386  {
1387  tAcorr().faceFluxCorrectionPtr() = (-A.flux()).ptr();
1388  }
1389 
1390  return tAcorr;
1391 }
1392 
1393 
1394 template<class Type>
1396 (
1397  const tmp<fvMatrix<Type> >& tA
1398 )
1399 {
1400  tmp<Foam::fvMatrix<Type> > tAcorr = tA - (tA() & tA().psi());
1401 
1402  // Note the matrix coefficients are still that of matrix A
1403  const fvMatrix<Type>& A = tAcorr();
1404 
1405  if
1406  (
1407  (A.hasUpper() || A.hasLower())
1408  && A.psi().mesh().fluxRequired(A.psi().name())
1409  )
1410  {
1411  tAcorr().faceFluxCorrectionPtr() = (-A.flux()).ptr();
1412  }
1413 
1414  return tAcorr;
1415 }
1416 
1417 
1418 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
1419 
1420 template<class Type>
1421 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1422 (
1423  const fvMatrix<Type>& A,
1424  const fvMatrix<Type>& B
1425 )
1426 {
1427  checkMethod(A, B, "==");
1428  return (A - B);
1429 }
1430 
1431 template<class Type>
1432 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1433 (
1434  const tmp<fvMatrix<Type> >& tA,
1435  const fvMatrix<Type>& B
1436 )
1437 {
1438  checkMethod(tA(), B, "==");
1439  return (tA - B);
1440 }
1441 
1442 template<class Type>
1443 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1444 (
1445  const fvMatrix<Type>& A,
1446  const tmp<fvMatrix<Type> >& tB
1447 )
1448 {
1449  checkMethod(A, tB(), "==");
1450  return (A - tB);
1451 }
1452 
1453 template<class Type>
1454 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1455 (
1456  const tmp<fvMatrix<Type> >& tA,
1457  const tmp<fvMatrix<Type> >& tB
1458 )
1459 {
1460  checkMethod(tA(), tB(), "==");
1461  return (tA - tB);
1462 }
1463 
1464 template<class Type>
1465 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1466 (
1467  const fvMatrix<Type>& A,
1468  const DimensionedField<Type, volMesh>& su
1469 )
1470 {
1471  checkMethod(A, su, "==");
1472  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1473  tC().source() += su.mesh().V()*su.field();
1474  return tC;
1475 }
1476 
1477 template<class Type>
1478 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1479 (
1480  const fvMatrix<Type>& A,
1481  const tmp<DimensionedField<Type, volMesh> >& tsu
1482 )
1483 {
1484  checkMethod(A, tsu(), "==");
1485  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1486  tC().source() += tsu().mesh().V()*tsu().field();
1487  tsu.clear();
1488  return tC;
1489 }
1490 
1491 template<class Type>
1492 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1493 (
1494  const fvMatrix<Type>& A,
1495  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1496 )
1497 {
1498  checkMethod(A, tsu(), "==");
1499  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1500  tC().source() += tsu().mesh().V()*tsu().internalField();
1501  tsu.clear();
1502  return tC;
1503 }
1504 
1505 template<class Type>
1506 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1507 (
1508  const tmp<fvMatrix<Type> >& tA,
1509  const DimensionedField<Type, volMesh>& su
1510 )
1511 {
1512  checkMethod(tA(), su, "==");
1513  tmp<fvMatrix<Type> > tC(tA.ptr());
1514  tC().source() += su.mesh().V()*su.field();
1515  return tC;
1516 }
1517 
1518 template<class Type>
1519 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1520 (
1521  const tmp<fvMatrix<Type> >& tA,
1522  const tmp<DimensionedField<Type, volMesh> >& tsu
1523 )
1524 {
1525  checkMethod(tA(), tsu(), "==");
1526  tmp<fvMatrix<Type> > tC(tA.ptr());
1527  tC().source() += tsu().mesh().V()*tsu().field();
1528  tsu.clear();
1529  return tC;
1530 }
1531 
1532 template<class Type>
1533 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1534 (
1535  const tmp<fvMatrix<Type> >& tA,
1536  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1537 )
1538 {
1539  checkMethod(tA(), tsu(), "==");
1540  tmp<fvMatrix<Type> > tC(tA.ptr());
1541  tC().source() += tsu().mesh().V()*tsu().internalField();
1542  tsu.clear();
1543  return tC;
1544 }
1545 
1546 template<class Type>
1547 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1548 (
1549  const fvMatrix<Type>& A,
1550  const dimensioned<Type>& su
1551 )
1552 {
1553  checkMethod(A, su, "==");
1554  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1555  tC().source() += A.psi().mesh().V()*su.value();
1556  return tC;
1557 }
1558 
1559 template<class Type>
1560 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1561 (
1562  const tmp<fvMatrix<Type> >& tA,
1563  const dimensioned<Type>& su
1564 )
1565 {
1566  checkMethod(tA(), su, "==");
1567  tmp<fvMatrix<Type> > tC(tA.ptr());
1568  tC().source() += tC().psi().mesh().V()*su.value();
1569  return tC;
1570 }
1571 
1572 template<class Type>
1573 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1574 (
1575  const fvMatrix<Type>& A,
1576  const zero&
1577 )
1578 {
1579  return A;
1580 }
1581 
1582 
1583 template<class Type>
1584 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator==
1585 (
1586  const tmp<fvMatrix<Type> >& tA,
1587  const zero&
1588 )
1589 {
1590  return tA;
1591 }
1592 
1593 
1594 template<class Type>
1595 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1596 (
1597  const fvMatrix<Type>& A
1598 )
1599 {
1600  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1601  tC().negate();
1602  return tC;
1603 }
1604 
1605 template<class Type>
1606 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1607 (
1608  const tmp<fvMatrix<Type> >& tA
1609 )
1610 {
1611  tmp<fvMatrix<Type> > tC(tA.ptr());
1612  tC().negate();
1613  return tC;
1614 }
1615 
1616 
1617 template<class Type>
1618 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1619 (
1620  const fvMatrix<Type>& A,
1621  const fvMatrix<Type>& B
1622 )
1623 {
1624  checkMethod(A, B, "+");
1625  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1626  tC() += B;
1627  return tC;
1628 }
1629 
1630 template<class Type>
1631 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1632 (
1633  const tmp<fvMatrix<Type> >& tA,
1634  const fvMatrix<Type>& B
1635 )
1636 {
1637  checkMethod(tA(), B, "+");
1638  tmp<fvMatrix<Type> > tC(tA.ptr());
1639  tC() += B;
1640  return tC;
1641 }
1642 
1643 template<class Type>
1644 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1645 (
1646  const fvMatrix<Type>& A,
1647  const tmp<fvMatrix<Type> >& tB
1648 )
1649 {
1650  checkMethod(A, tB(), "+");
1651  tmp<fvMatrix<Type> > tC(tB.ptr());
1652  tC() += A;
1653  return tC;
1654 }
1655 
1656 template<class Type>
1657 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1658 (
1659  const tmp<fvMatrix<Type> >& tA,
1660  const tmp<fvMatrix<Type> >& tB
1661 )
1662 {
1663  checkMethod(tA(), tB(), "+");
1664  tmp<fvMatrix<Type> > tC(tA.ptr());
1665  tC() += tB();
1666  tB.clear();
1667  return tC;
1668 }
1669 
1670 template<class Type>
1671 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1672 (
1673  const fvMatrix<Type>& A,
1674  const DimensionedField<Type, volMesh>& su
1675 )
1676 {
1677  checkMethod(A, su, "+");
1678  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1679  tC().source() -= su.mesh().V()*su.field();
1680  return tC;
1681 }
1682 
1683 template<class Type>
1684 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1685 (
1686  const fvMatrix<Type>& A,
1687  const tmp<DimensionedField<Type, volMesh> >& tsu
1688 )
1689 {
1690  checkMethod(A, tsu(), "+");
1691  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1692  tC().source() -= tsu().mesh().V()*tsu().field();
1693  tsu.clear();
1694  return tC;
1695 }
1696 
1697 template<class Type>
1698 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1699 (
1700  const fvMatrix<Type>& A,
1701  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1702 )
1703 {
1704  checkMethod(A, tsu(), "+");
1705  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1706  tC().source() -= tsu().mesh().V()*tsu().internalField();
1707  tsu.clear();
1708  return tC;
1709 }
1710 
1711 template<class Type>
1712 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1713 (
1714  const tmp<fvMatrix<Type> >& tA,
1715  const DimensionedField<Type, volMesh>& su
1716 )
1717 {
1718  checkMethod(tA(), su, "+");
1719  tmp<fvMatrix<Type> > tC(tA.ptr());
1720  tC().source() -= su.mesh().V()*su.field();
1721  return tC;
1722 }
1723 
1724 template<class Type>
1725 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1726 (
1727  const tmp<fvMatrix<Type> >& tA,
1728  const tmp<DimensionedField<Type, volMesh> >& tsu
1729 )
1730 {
1731  checkMethod(tA(), tsu(), "+");
1732  tmp<fvMatrix<Type> > tC(tA.ptr());
1733  tC().source() -= tsu().mesh().V()*tsu().field();
1734  tsu.clear();
1735  return tC;
1736 }
1737 
1738 template<class Type>
1739 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1740 (
1741  const tmp<fvMatrix<Type> >& tA,
1742  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1743 )
1744 {
1745  checkMethod(tA(), tsu(), "+");
1746  tmp<fvMatrix<Type> > tC(tA.ptr());
1747  tC().source() -= tsu().mesh().V()*tsu().internalField();
1748  tsu.clear();
1749  return tC;
1750 }
1751 
1752 template<class Type>
1753 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1754 (
1755  const DimensionedField<Type, volMesh>& su,
1756  const fvMatrix<Type>& A
1757 )
1758 {
1759  checkMethod(A, su, "+");
1760  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1761  tC().source() -= su.mesh().V()*su.field();
1762  return tC;
1763 }
1764 
1765 template<class Type>
1766 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1767 (
1768  const tmp<DimensionedField<Type, volMesh> >& tsu,
1769  const fvMatrix<Type>& A
1770 )
1771 {
1772  checkMethod(A, tsu(), "+");
1773  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1774  tC().source() -= tsu().mesh().V()*tsu().field();
1775  tsu.clear();
1776  return tC;
1777 }
1778 
1779 template<class Type>
1780 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1781 (
1782  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
1783  const fvMatrix<Type>& A
1784 )
1785 {
1786  checkMethod(A, tsu(), "+");
1787  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1788  tC().source() -= tsu().mesh().V()*tsu().internalField();
1789  tsu.clear();
1790  return tC;
1791 }
1792 
1793 template<class Type>
1794 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1795 (
1796  const DimensionedField<Type, volMesh>& su,
1797  const tmp<fvMatrix<Type> >& tA
1798 )
1799 {
1800  checkMethod(tA(), su, "+");
1801  tmp<fvMatrix<Type> > tC(tA.ptr());
1802  tC().source() -= su.mesh().V()*su.field();
1803  return tC;
1804 }
1805 
1806 template<class Type>
1807 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1808 (
1809  const tmp<DimensionedField<Type, volMesh> >& tsu,
1810  const tmp<fvMatrix<Type> >& tA
1811 )
1812 {
1813  checkMethod(tA(), tsu(), "+");
1814  tmp<fvMatrix<Type> > tC(tA.ptr());
1815  tC().source() -= tsu().mesh().V()*tsu().field();
1816  tsu.clear();
1817  return tC;
1818 }
1819 
1820 template<class Type>
1821 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
1822 (
1823  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
1824  const tmp<fvMatrix<Type> >& tA
1825 )
1826 {
1827  checkMethod(tA(), tsu(), "+");
1828  tmp<fvMatrix<Type> > tC(tA.ptr());
1829  tC().source() -= tsu().mesh().V()*tsu().internalField();
1830  tsu.clear();
1831  return tC;
1832 }
1833 
1834 
1835 template<class Type>
1836 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1837 (
1838  const fvMatrix<Type>& A,
1839  const fvMatrix<Type>& B
1840 )
1841 {
1842  checkMethod(A, B, "-");
1843  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1844  tC() -= B;
1845  return tC;
1846 }
1847 
1848 template<class Type>
1849 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1850 (
1851  const tmp<fvMatrix<Type> >& tA,
1852  const fvMatrix<Type>& B
1853 )
1854 {
1855  checkMethod(tA(), B, "-");
1856  tmp<fvMatrix<Type> > tC(tA.ptr());
1857  tC() -= B;
1858  return tC;
1859 }
1860 
1861 template<class Type>
1862 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1863 (
1864  const fvMatrix<Type>& A,
1865  const tmp<fvMatrix<Type> >& tB
1866 )
1867 {
1868  checkMethod(A, tB(), "-");
1869  tmp<fvMatrix<Type> > tC(tB.ptr());
1870  tC() -= A;
1871  tC().negate();
1872  return tC;
1873 }
1874 
1875 template<class Type>
1876 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1877 (
1878  const tmp<fvMatrix<Type> >& tA,
1879  const tmp<fvMatrix<Type> >& tB
1880 )
1881 {
1882  checkMethod(tA(), tB(), "-");
1883  tmp<fvMatrix<Type> > tC(tA.ptr());
1884  tC() -= tB();
1885  tB.clear();
1886  return tC;
1887 }
1888 
1889 template<class Type>
1890 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1891 (
1892  const fvMatrix<Type>& A,
1893  const DimensionedField<Type, volMesh>& su
1894 )
1895 {
1896  checkMethod(A, su, "-");
1897  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1898  tC().source() += su.mesh().V()*su.field();
1899  return tC;
1900 }
1901 
1902 template<class Type>
1903 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1904 (
1905  const fvMatrix<Type>& A,
1906  const tmp<DimensionedField<Type, volMesh> >& tsu
1907 )
1908 {
1909  checkMethod(A, tsu(), "-");
1910  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1911  tC().source() += tsu().mesh().V()*tsu().field();
1912  tsu.clear();
1913  return tC;
1914 }
1915 
1916 template<class Type>
1917 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1918 (
1919  const fvMatrix<Type>& A,
1920  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1921 )
1922 {
1923  checkMethod(A, tsu(), "-");
1924  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1925  tC().source() += tsu().mesh().V()*tsu().internalField();
1926  tsu.clear();
1927  return tC;
1928 }
1929 
1930 template<class Type>
1931 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1932 (
1933  const tmp<fvMatrix<Type> >& tA,
1934  const DimensionedField<Type, volMesh>& su
1935 )
1936 {
1937  checkMethod(tA(), su, "-");
1938  tmp<fvMatrix<Type> > tC(tA.ptr());
1939  tC().source() += su.mesh().V()*su.field();
1940  return tC;
1941 }
1942 
1943 template<class Type>
1944 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1945 (
1946  const tmp<fvMatrix<Type> >& tA,
1947  const tmp<DimensionedField<Type, volMesh> >& tsu
1948 )
1949 {
1950  checkMethod(tA(), tsu(), "-");
1951  tmp<fvMatrix<Type> > tC(tA.ptr());
1952  tC().source() += tsu().mesh().V()*tsu().field();
1953  tsu.clear();
1954  return tC;
1955 }
1956 
1957 template<class Type>
1958 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1959 (
1960  const tmp<fvMatrix<Type> >& tA,
1961  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu
1962 )
1963 {
1964  checkMethod(tA(), tsu(), "-");
1965  tmp<fvMatrix<Type> > tC(tA.ptr());
1966  tC().source() += tsu().mesh().V()*tsu().internalField();
1967  tsu.clear();
1968  return tC;
1969 }
1970 
1971 template<class Type>
1972 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1973 (
1974  const DimensionedField<Type, volMesh>& su,
1975  const fvMatrix<Type>& A
1976 )
1977 {
1978  checkMethod(A, su, "-");
1979  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1980  tC().negate();
1981  tC().source() -= su.mesh().V()*su.field();
1982  return tC;
1983 }
1984 
1985 template<class Type>
1986 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
1987 (
1988  const tmp<DimensionedField<Type, volMesh> >& tsu,
1989  const fvMatrix<Type>& A
1990 )
1991 {
1992  checkMethod(A, tsu(), "-");
1993  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
1994  tC().negate();
1995  tC().source() -= tsu().mesh().V()*tsu().field();
1996  tsu.clear();
1997  return tC;
1998 }
1999 
2000 template<class Type>
2001 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2002 (
2003  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
2004  const fvMatrix<Type>& A
2005 )
2006 {
2007  checkMethod(A, tsu(), "-");
2008  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2009  tC().negate();
2010  tC().source() -= tsu().mesh().V()*tsu().internalField();
2011  tsu.clear();
2012  return tC;
2013 }
2014 
2015 template<class Type>
2016 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2017 (
2018  const DimensionedField<Type, volMesh>& su,
2019  const tmp<fvMatrix<Type> >& tA
2020 )
2021 {
2022  checkMethod(tA(), su, "-");
2023  tmp<fvMatrix<Type> > tC(tA.ptr());
2024  tC().negate();
2025  tC().source() -= su.mesh().V()*su.field();
2026  return tC;
2027 }
2028 
2029 template<class Type>
2030 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2031 (
2032  const tmp<DimensionedField<Type, volMesh> >& tsu,
2033  const tmp<fvMatrix<Type> >& tA
2034 )
2035 {
2036  checkMethod(tA(), tsu(), "-");
2037  tmp<fvMatrix<Type> > tC(tA.ptr());
2038  tC().negate();
2039  tC().source() -= tsu().mesh().V()*tsu().field();
2040  tsu.clear();
2041  return tC;
2042 }
2043 
2044 template<class Type>
2045 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2046 (
2047  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tsu,
2048  const tmp<fvMatrix<Type> >& tA
2049 )
2050 {
2051  checkMethod(tA(), tsu(), "-");
2052  tmp<fvMatrix<Type> > tC(tA.ptr());
2053  tC().negate();
2054  tC().source() -= tsu().mesh().V()*tsu().internalField();
2055  tsu.clear();
2056  return tC;
2057 }
2058 
2059 template<class Type>
2060 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
2061 (
2062  const fvMatrix<Type>& A,
2063  const dimensioned<Type>& su
2064 )
2065 {
2066  checkMethod(A, su, "+");
2067  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2068  tC().source() -= su.value()*A.psi().mesh().V();
2069  return tC;
2070 }
2071 
2072 template<class Type>
2073 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
2074 (
2075  const tmp<fvMatrix<Type> >& tA,
2076  const dimensioned<Type>& su
2077 )
2078 {
2079  checkMethod(tA(), su, "+");
2080  tmp<fvMatrix<Type> > tC(tA.ptr());
2081  tC().source() -= su.value()*tC().psi().mesh().V();
2082  return tC;
2083 }
2084 
2085 template<class Type>
2086 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
2087 (
2088  const dimensioned<Type>& su,
2089  const fvMatrix<Type>& A
2090 )
2091 {
2092  checkMethod(A, su, "+");
2093  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2094  tC().source() -= su.value()*A.psi().mesh().V();
2095  return tC;
2096 }
2097 
2098 template<class Type>
2099 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator+
2100 (
2101  const dimensioned<Type>& su,
2102  const tmp<fvMatrix<Type> >& tA
2103 )
2104 {
2105  checkMethod(tA(), su, "+");
2106  tmp<fvMatrix<Type> > tC(tA.ptr());
2107  tC().source() -= su.value()*tC().psi().mesh().V();
2108  return tC;
2109 }
2110 
2111 template<class Type>
2112 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2113 (
2114  const fvMatrix<Type>& A,
2115  const dimensioned<Type>& su
2116 )
2117 {
2118  checkMethod(A, su, "-");
2119  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2120  tC().source() += su.value()*tC().psi().mesh().V();
2121  return tC;
2122 }
2123 
2124 template<class Type>
2125 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2126 (
2127  const tmp<fvMatrix<Type> >& tA,
2128  const dimensioned<Type>& su
2129 )
2130 {
2131  checkMethod(tA(), su, "-");
2132  tmp<fvMatrix<Type> > tC(tA.ptr());
2133  tC().source() += su.value()*tC().psi().mesh().V();
2134  return tC;
2135 }
2136 
2137 template<class Type>
2138 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2139 (
2140  const dimensioned<Type>& su,
2141  const fvMatrix<Type>& A
2142 )
2143 {
2144  checkMethod(A, su, "-");
2145  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2146  tC().negate();
2147  tC().source() -= su.value()*A.psi().mesh().V();
2148  return tC;
2149 }
2150 
2151 template<class Type>
2152 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator-
2153 (
2154  const dimensioned<Type>& su,
2155  const tmp<fvMatrix<Type> >& tA
2156 )
2157 {
2158  checkMethod(tA(), su, "-");
2159  tmp<fvMatrix<Type> > tC(tA.ptr());
2160  tC().negate();
2161  tC().source() -= su.value()*tC().psi().mesh().V();
2162  return tC;
2163 }
2164 
2165 
2166 template<class Type>
2167 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2168 (
2169  const DimensionedField<scalar, volMesh>& dsf,
2170  const fvMatrix<Type>& A
2171 )
2172 {
2173  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2174  tC() *= dsf;
2175  return tC;
2176 }
2177 
2178 template<class Type>
2179 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2180 (
2181  const tmp< DimensionedField<scalar, volMesh> >& tdsf,
2182  const fvMatrix<Type>& A
2183 )
2184 {
2185  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2186  tC() *= tdsf;
2187  return tC;
2188 }
2189 
2190 template<class Type>
2191 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2192 (
2193  const tmp<volScalarField>& tvsf,
2194  const fvMatrix<Type>& A
2195 )
2196 {
2197  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2198  tC() *= tvsf;
2199  return tC;
2200 }
2201 
2202 template<class Type>
2203 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2204 (
2205  const DimensionedField<scalar, volMesh>& dsf,
2206  const tmp<fvMatrix<Type> >& tA
2207 )
2208 {
2209  tmp<fvMatrix<Type> > tC(tA.ptr());
2210  tC() *= dsf;
2211  return tC;
2212 }
2213 
2214 template<class Type>
2215 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2216 (
2217  const tmp<DimensionedField<scalar, volMesh> >& tdsf,
2218  const tmp<fvMatrix<Type> >& tA
2219 )
2220 {
2221  tmp<fvMatrix<Type> > tC(tA.ptr());
2222  tC() *= tdsf;
2223  return tC;
2224 }
2225 
2226 template<class Type>
2227 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2228 (
2229  const tmp<volScalarField>& tvsf,
2230  const tmp<fvMatrix<Type> >& tA
2231 )
2232 {
2233  tmp<fvMatrix<Type> > tC(tA.ptr());
2234  tC() *= tvsf;
2235  return tC;
2236 }
2237 
2238 template<class Type>
2239 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2240 (
2241  const dimensioned<scalar>& ds,
2242  const fvMatrix<Type>& A
2243 )
2244 {
2245  tmp<fvMatrix<Type> > tC(new fvMatrix<Type>(A));
2246  tC() *= ds;
2247  return tC;
2248 }
2249 
2250 template<class Type>
2251 Foam::tmp<Foam::fvMatrix<Type> > Foam::operator*
2252 (
2253  const dimensioned<scalar>& ds,
2254  const tmp<fvMatrix<Type> >& tA
2255 )
2256 {
2257  tmp<fvMatrix<Type> > tC(tA.ptr());
2258  tC() *= ds;
2259  return tC;
2260 }
2261 
2262 
2263 template<class Type>
2265 Foam::operator&
2266 (
2267  const fvMatrix<Type>& M,
2268  const DimensionedField<Type, volMesh>& psi
2269 )
2270 {
2271  tmp<GeometricField<Type, fvPatchField, volMesh> > tMphi
2272  (
2273  new GeometricField<Type, fvPatchField, volMesh>
2274  (
2275  IOobject
2276  (
2277  "M&" + psi.name(),
2278  psi.instance(),
2279  psi.mesh(),
2280  IOobject::NO_READ,
2281  IOobject::NO_WRITE
2282  ),
2283  psi.mesh(),
2284  M.dimensions()/dimVol,
2285  zeroGradientFvPatchScalarField::typeName
2286  )
2287  );
2288  GeometricField<Type, fvPatchField, volMesh>& Mphi = tMphi();
2289 
2290  // Loop over field components
2291  if (M.hasDiag())
2292  {
2293  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2294  {
2295  scalarField psiCmpt(psi.field().component(cmpt));
2296  scalarField boundaryDiagCmpt(M.diag());
2297  M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2298  Mphi.internalField().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2299  }
2300  }
2301  else
2302  {
2303  Mphi.internalField() = pTraits<Type>::zero;
2304  }
2305 
2306  Mphi.internalField() += M.lduMatrix::H(psi.field()) + M.source();
2307  M.addBoundarySource(Mphi.internalField());
2308 
2309  Mphi.internalField() /= -psi.mesh().V();
2310  Mphi.correctBoundaryConditions();
2311 
2312  return tMphi;
2313 }
2314 
2315 template<class Type>
2317 Foam::operator&
2318 (
2319  const fvMatrix<Type>& M,
2320  const tmp<DimensionedField<Type, volMesh> >& tpsi
2321 )
2322 {
2323  tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = M & tpsi();
2324  tpsi.clear();
2325  return tMpsi;
2326 }
2327 
2328 template<class Type>
2330 Foam::operator&
2331 (
2332  const fvMatrix<Type>& M,
2333  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tpsi
2334 )
2335 {
2336  tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = M & tpsi();
2337  tpsi.clear();
2338  return tMpsi;
2339 }
2340 
2341 template<class Type>
2343 Foam::operator&
2344 (
2345  const tmp<fvMatrix<Type> >& tM,
2346  const DimensionedField<Type, volMesh>& psi
2347 )
2348 {
2349  tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & psi;
2350  tM.clear();
2351  return tMpsi;
2352 }
2353 
2354 template<class Type>
2356 Foam::operator&
2357 (
2358  const tmp<fvMatrix<Type> >& tM,
2359  const tmp<DimensionedField<Type, volMesh> >& tpsi
2360 )
2361 {
2362  tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & tpsi();
2363  tM.clear();
2364  tpsi.clear();
2365  return tMpsi;
2366 }
2367 
2368 template<class Type>
2370 Foam::operator&
2371 (
2372  const tmp<fvMatrix<Type> >& tM,
2373  const tmp<GeometricField<Type, fvPatchField, volMesh> >& tpsi
2374 )
2375 {
2376  tmp<GeometricField<Type, fvPatchField, volMesh> > tMpsi = tM() & tpsi();
2377  tM.clear();
2378  tpsi.clear();
2379  return tMpsi;
2380 }
2381 
2382 
2383 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2384 
2385 template<class Type>
2386 Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
2387 {
2388  os << static_cast<const lduMatrix&>(fvm) << nl
2389  << fvm.dimensions_ << nl
2390  << fvm.source_ << nl
2391  << fvm.internalCoeffs_ << nl
2392  << fvm.boundaryCoeffs_ << endl;
2393 
2394  os.check("Ostream& operator<<(Ostream&, fvMatrix<Type>&");
2395 
2396  return os;
2397 }
2398 
2399 
2400 // * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
2401 
2402 #include "fvMatrixSolve.C"
2403 
2404 // ************************************************************************* //
Foam::fvPatchField< Type >
volFields.H
Foam::fvMatrix::H
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:763
Foam::fvMatrix::addBoundarySource
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:145
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
fvMatrixSolve.C
InfoInFunction
#define InfoInFunction
Report a information message using Foam::Info.
Definition: messageStream.H:281
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:41
relax
p relax()
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvMatrix::dimensions_
dimensionSet dimensions_
Dimension set.
Definition: fvMatrix.H:126
UIndirectList.H
Foam::fvMatrix::addCmptAvBoundaryDiag
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:129
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::checkMethod
void checkMethod(const fvMatrix< Type > &, const fvMatrix< Type > &, const char *)
Definition: fvMatrix.C:1258
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:262
Foam::fvMatrix::D
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition: fvMatrix.C:701
Foam::fvMatrix::subtractFromInternalField
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:75
calculatedFvPatchFields.H
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::fvMatrix::operator=
void operator=(const fvMatrix< Type > &)
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:77
Foam::fvMatrix::dimensions
const dimensionSet & dimensions() const
Definition: fvMatrix.H:286
boundary
faceListList boundary(nPatches)
Foam::fvMatrix::relax
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:671
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:235
Foam::fvMatrix::negate
void negate()
Definition: fvMatrix.C:1001
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::fvMatrix::addBoundaryDiag
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:111
Foam::fvMatrix::setReference
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:504
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:116
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
A
simpleMatrix< scalar > A(Nc)
Foam::fvMatrix::faceFluxCorrectionPtr_
GeometricField< Type, fvsPatchField, surfaceMesh > * faceFluxCorrectionPtr_
Face flux field for non-orthogonal correction.
Definition: fvMatrix.H:142
Foam::cmptMin
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:269
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:359
coupledFvPatchFields.H
Foam::cmptMax
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:224
solve
rhoEqn solve()
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:342
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::fvMatrix::setValues
void setValues(const labelUList &cells, const UList< Type > &values)
Set solution in given cells to the specified values.
Foam::label
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
Foam::fvMatrix::operator-=
void operator-=(const fvMatrix< Type > &)
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvMatrix::internalCoeffs_
FieldField< Field, Type > internalCoeffs_
Boundary scalar field containing pseudo-matrix coeffs.
Definition: fvMatrix.H:133
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvMatrix::H1
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:825
Foam::fvMatrix::~fvMatrix
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:462
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::Field::replace
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
Foam::fvMatrix::boundaryManipulate
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::GeometricBoundaryField &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:688
Foam::negate
void negate(FieldField< Field, Type > &res, const FieldField< Field, Type > &f)
Foam::fvMatrix::boundaryCoeffs_
FieldField< Field, Type > boundaryCoeffs_
Boundary scalar field containing pseudo-matrix coeffs.
Definition: fvMatrix.H:137
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
Foam::FatalError
error FatalError
Foam::Field::negate
void negate()
Negate this field.
Definition: Field.C:635
Foam::fvMatrix::psi_
const GeometricField< Type, fvPatchField, volMesh > & psi_
Const reference to GeometricField<Type, fvPatchField, volMesh>
Definition: fvMatrix.H:123
Foam::fvMatrix::DD
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:710
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned< Type >
Foam::fvMatrix::A
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:734
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::tmp::ptr
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:146
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::DimensionedField::dimensions
const dimensionSet & dimensions() const
Return dimensions.
Definition: DimensionedFieldI.H:46
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:408
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
internalField
conserve internalField()+
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvMatrix::setValuesFromList
void setValuesFromList(const labelUList &cells, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:178
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::cmptAv
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:252
Foam::sumOp
Definition: ops.H:162
Foam::fvMatrix::operator+=
void operator+=(const fvMatrix< Type > &)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
scalarField
volScalarField scalarField(fieldObject, mesh)
M
#define M(I)
Foam::UList
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
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::fvMatrix< Type >
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::GeometricField::replace
void replace(const direction, const GeometricField< cmptType, PatchField, GeoMesh > &)
Foam::dimVol
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::labelUList
UList< label > labelUList
Definition: UList.H:63
Foam::UList::size
label size() const
Return the number of elements in the UList.
Definition: UListI.H:299
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:248
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvMatrix::addToInternalField
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:38
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:49
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
zeroGradientFvPatchFields.H
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::fvMatrix::flux
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:873
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::fvMatrix::fvMatrix
fvMatrix(const GeometricField< Type, fvPatchField, volMesh > &, const dimensionSet &)
Construct given a field to solve for.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::fvMatrix::source_
Field< Type > source_
Source term.
Definition: fvMatrix.H:129
Foam::zero
A class representing the concept of 0 used to avoid unnecessary manipulations for objects that are kn...
Definition: zero.H:47