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