faMatrix.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 Wikki Ltd
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "areaFields.H"
30 #include "edgeFields.H"
34 #include "IndirectList.H"
35 #include "UniformList.H"
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 template<class Type>
41 template<class Type2>
43 (
44  const labelUList& addr,
45  const Field<Type2>& pf,
46  Field<Type2>& intf
47 ) const
48 {
49  if (addr.size() != pf.size())
50  {
52  << "addressing (" << addr.size()
53  << ") and field (" << pf.size() << ") are different sizes" << endl
54  << abort(FatalError);
55  }
56 
57  forAll(addr, faceI)
58  {
59  intf[addr[faceI]] += pf[faceI];
60  }
61 }
62 
63 
64 template<class Type>
65 template<class Type2>
67 (
68  const labelUList& addr,
69  const tmp<Field<Type2>>& tpf,
70  Field<Type2>& intf
71 ) const
72 {
73  addToInternalField(addr, tpf(), intf);
74  tpf.clear();
75 }
76 
77 
78 template<class Type>
79 template<class Type2>
81 (
82  const labelUList& addr,
83  const Field<Type2>& pf,
84  Field<Type2>& intf
85 ) const
86 {
87  if (addr.size() != pf.size())
88  {
90  << "addressing (" << addr.size()
91  << ") and field (" << pf.size() << ") are different sizes" << endl
92  << abort(FatalError);
93  }
94 
95  forAll(addr, faceI)
96  {
97  intf[addr[faceI]] -= pf[faceI];
98  }
99 }
100 
101 
102 template<class Type>
103 template<class Type2>
105 (
106  const labelUList& addr,
107  const tmp<Field<Type2>>& tpf,
108  Field<Type2>& intf
109 ) const
110 {
111  subtractFromInternalField(addr, tpf(), intf);
112  tpf.clear();
113 }
114 
115 
116 template<class Type>
118 (
119  scalarField& diag,
120  const direction solveCmpt
121 ) const
122 {
123  forAll(internalCoeffs_, patchI)
124  {
125  addToInternalField
126  (
127  lduAddr().patchAddr(patchI),
128  internalCoeffs_[patchI].component(solveCmpt),
130  );
131  }
132 }
133 
134 
135 template<class Type>
137 {
138  forAll(internalCoeffs_, patchI)
139  {
140  addToInternalField
141  (
142  lduAddr().patchAddr(patchI),
143  cmptAv(internalCoeffs_[patchI]),
144  diag
145  );
146  }
147 }
148 
149 
150 template<class Type>
152 (
153  Field<Type>& source,
154  const bool couples
155 ) const
156 {
157  forAll(psi_.boundaryField(), patchI)
158  {
159  const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
160  const Field<Type>& pbc = boundaryCoeffs_[patchI];
161 
162  if (!ptf.coupled())
163  {
164  addToInternalField(lduAddr().patchAddr(patchI), pbc, source);
165  }
166  else if (couples)
167  {
168  tmp<Field<Type>> tpnf = ptf.patchNeighbourField();
169  const Field<Type>& pnf = tpnf();
170 
171  const labelUList& addr = lduAddr().patchAddr(patchI);
172 
173  forAll(addr, facei)
174  {
175  source[addr[facei]] += cmptMultiply(pbc[facei], pnf[facei]);
176  }
177  }
178  }
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
183 
184 template<class Type>
186 (
188  const dimensionSet& ds
189 )
190 :
191  lduMatrix(psi.mesh()),
192  psi_(psi),
193  dimensions_(ds),
194  source_(psi.size(), Zero),
195  internalCoeffs_(psi.mesh().boundary().size()),
196  boundaryCoeffs_(psi.mesh().boundary().size()),
197  faceFluxCorrectionPtr_(nullptr)
198 {
200  << "constructing faMatrix<Type> for field " << psi_.name()
201  << endl;
202 
203  // Initialise coupling coefficients
204  forAll(psi.mesh().boundary(), patchi)
205  {
206  internalCoeffs_.set
207  (
208  patchi,
209  new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
210  );
211 
212  boundaryCoeffs_.set
213  (
214  patchi,
215  new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
216  );
217  }
218 
219  // Update the boundary coefficients of psi without changing its event No.
220  auto& psiRef =
221  const_cast<GeometricField<Type, faPatchField, areaMesh>&>(psi_);
222 
223  const label currentStatePsi = psiRef.eventNo();
224  psiRef.boundaryFieldRef().updateCoeffs();
225  psiRef.eventNo() = currentStatePsi;
226 }
227 
228 
229 template<class Type>
231 :
232  lduMatrix(fam),
233  psi_(fam.psi_),
234  dimensions_(fam.dimensions_),
235  source_(fam.source_),
236  internalCoeffs_(fam.internalCoeffs_),
237  boundaryCoeffs_(fam.boundaryCoeffs_),
238  faceFluxCorrectionPtr_(nullptr)
239 {
241  << "Copying faMatrix<Type> for field " << psi_.name() << endl;
242 
243  if (fam.faceFluxCorrectionPtr_)
244  {
245  faceFluxCorrectionPtr_ =
247  (
248  *(fam.faceFluxCorrectionPtr_)
249  );
250  }
251 }
252 
253 
254 template<class Type>
255 Foam::faMatrix<Type>::faMatrix(const tmp<faMatrix<Type>>& tmat)
256 :
257  lduMatrix(tmat.constCast(), tmat.movable()),
258  psi_(tmat().psi_),
259  dimensions_(tmat().dimensions_),
260  source_(tmat.constCast().source_, tmat.movable()),
261  internalCoeffs_(tmat.constCast().internalCoeffs_, tmat.movable()),
262  boundaryCoeffs_(tmat.constCast().boundaryCoeffs_, tmat.movable()),
263  faceFluxCorrectionPtr_(nullptr)
264 {
266  << "Copy/Move faMatrix<Type> for field " << psi_.name() << endl;
267 
268  if (tmat().faceFluxCorrectionPtr_)
269  {
270  if (tmat.movable())
271  {
272  faceFluxCorrectionPtr_ = tmat().faceFluxCorrectionPtr_;
273  tmat().faceFluxCorrectionPtr_ = nullptr;
274  }
275  else
276  {
277  faceFluxCorrectionPtr_ =
278  new GeometricField<Type, faePatchField, edgeMesh>
279  (
280  *(tmat().faceFluxCorrectionPtr_)
281  );
282  }
283  }
284 
285  tmat.clear();
286 }
287 
288 
289 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
290 
291 template<class Type>
293 {
295  << "Destroying faMatrix<Type> for field " << psi_.name() << endl;
296 
297  deleteDemandDrivenData(faceFluxCorrectionPtr_);
298 }
300 
301 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
302 
303 template<class Type>
304 template<template<class> class ListType>
306 (
307  const labelUList& faceLabels,
308  const ListType<Type>& values
309 )
310 {
311 #if 0 /* Specific to foam-extend */
312  // Record face labels of eliminated equations
313  for (const label i : faceLabels)
314  {
315  this->eliminatedEqns().insert(i);
316  }
317 #endif
318 
319  const faMesh& mesh = psi_.mesh();
320 
321  const labelListList& edges = mesh.patch().faceEdges();
322  const labelUList& own = mesh.owner();
323  const labelUList& nei = mesh.neighbour();
324 
325  scalarField& Diag = diag();
326  Field<Type>& psi =
327  const_cast
328  <
330  >(psi_).primitiveFieldRef();
331 
332  // Following actions:
333  // - adjust local field psi
334  // - set local matrix to be diagonal (so adjust source)
335  // - cut connections to neighbours
336  // - make (on non-adjusted cells) contribution explicit
337 
338  if (symmetric() || asymmetric())
339  {
340  forAll(faceLabels, i)
341  {
342  const label facei = faceLabels[i];
343  const Type& value = values[i];
344 
345  for (const label edgei : edges[facei])
346  {
347  if (mesh.isInternalEdge(edgei))
348  {
349  if (symmetric())
350  {
351  if (facei == own[edgei])
352  {
353  source_[nei[edgei]] -= upper()[edgei]*value;
354  }
355  else
356  {
357  source_[own[edgei]] -= upper()[edgei]*value;
358  }
359 
360  upper()[edgei] = 0.0;
361  }
362  else
363  {
364  if (facei == own[edgei])
365  {
366  source_[nei[edgei]] -= lower()[edgei]*value;
367  }
368  else
369  {
370  source_[own[edgei]] -= upper()[edgei]*value;
371  }
372 
373  upper()[edgei] = 0.0;
374  lower()[edgei] = 0.0;
375  }
376  }
377  else
378  {
379  const label patchi = mesh.boundary().whichPatch(edgei);
380 
381  if (internalCoeffs_[patchi].size())
382  {
383  const label patchEdgei =
384  mesh.boundary()[patchi].whichEdge(edgei);
385 
386  internalCoeffs_[patchi][patchEdgei] = Zero;
387  boundaryCoeffs_[patchi][patchEdgei] = Zero;
388  }
389  }
390  }
391  }
392  }
393 
394  // Note: above loop might have affected source terms on adjusted cells
395  // so make sure to adjust them afterwards
396  forAll(faceLabels, i)
397  {
398  const label facei = faceLabels[i];
399  const Type& value = values[i];
400 
401  psi[facei] = value;
402  source_[facei] = value*Diag[facei];
403  }
404 }
405 
406 
407 template<class Type>
409 (
410  const labelUList& faceLabels,
411  const Type& value
412 )
413 {
414  this->setValuesFromList(faceLabels, UniformList<Type>(value));
415 }
416 
417 
418 template<class Type>
420 (
421  const labelUList& faceLabels,
422  const UList<Type>& values
423 )
424 {
425  this->setValuesFromList(faceLabels, values);
426 }
427 
428 
429 template<class Type>
431 (
432  const labelUList& faceLabels,
434 )
435 {
436  this->setValuesFromList(faceLabels, values);
437 }
438 
439 
440 template<class Type>
442 (
443  const label facei,
444  const Type& value,
445  const bool forceReference
446 )
447 {
448  if ((forceReference || psi_.needReference()) && facei >= 0)
449  {
450  if (Pstream::master())
451  {
452  source()[facei] += diag()[facei]*value;
453  diag()[facei] += diag()[facei];
454  }
455  }
456 }
457 
458 
459 template<class Type>
461 (
462  const labelUList& faceLabels,
463  const Type& value,
464  const bool forceReference
465 )
466 {
467  if (forceReference || psi_.needReference())
468  {
469  forAll(faceLabels, facei)
470  {
471  const label faceId = faceLabels[facei];
472  if (faceId >= 0)
473  {
474  source()[faceId] += diag()[faceId]*value;
475  diag()[faceId] += diag()[faceId];
476  }
477  }
478  }
479 }
480 
481 
482 template<class Type>
484 (
485  const labelUList& faceLabels,
486  const UList<Type>& values,
487  const bool forceReference
488 )
489 {
490  if (forceReference || psi_.needReference())
491  {
492  forAll(faceLabels, facei)
493  {
494  const label faceId = faceLabels[facei];
495  if (faceId >= 0)
496  {
497  source()[faceId] += diag()[faceId]*values[facei];
498  diag()[faceId] += diag()[faceId];
499  }
500  }
501  }
502 }
503 
504 
505 template<class Type>
506 void Foam::faMatrix<Type>::relax(const scalar alpha)
507 {
508  if (alpha <= 0)
509  {
510  return;
511  }
512 
513  Field<Type>& S = source();
514  scalarField& D = diag();
515 
516  // Store the current unrelaxed diagonal for use in updating the source
517  scalarField D0(D);
518 
519  // Calculate the sum-mag off-diagonal from the interior faces
520  scalarField sumOff(D.size(), Zero);
521  sumMagOffDiag(sumOff);
522 
523  // Handle the boundary contributions to the diagonal
524  forAll(psi_.boundaryField(), patchI)
525  {
526  const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
527 
528  if (ptf.size())
529  {
530  const labelUList& pa = lduAddr().patchAddr(patchI);
531  Field<Type>& iCoeffs = internalCoeffs_[patchI];
532 
533  if (ptf.coupled())
534  {
535  const Field<Type>& pCoeffs = boundaryCoeffs_[patchI];
536 
537  // For coupled boundaries add the diagonal and
538  // off-diagonal contributions
539  forAll(pa, face)
540  {
541  D[pa[face]] += component(iCoeffs[face], 0);
542  sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
543  }
544  }
545  else
546  {
547  // For non-coupled boundaries subtract the diagonal
548  // contribution off-diagonal sum which avoids having to remove
549  // it from the diagonal later.
550  // Also add the source contribution from the relaxation
551  forAll(pa, face)
552  {
553  Type iCoeff0 = iCoeffs[face];
554  iCoeffs[face] = cmptMag(iCoeffs[face]);
555  sumOff[pa[face]] -= cmptMin(iCoeffs[face]);
556  iCoeffs[face] /= alpha;
557  S[pa[face]] +=
558  cmptMultiply(iCoeffs[face] - iCoeff0, psi_[pa[face]]);
559  }
560  }
561  }
562  }
563 
564  // Ensure the matrix is diagonally dominant...
565  max(D, D, sumOff);
566 
567  // ... then relax
568  D /= alpha;
569 
570  // Now remove the diagonal contribution from coupled boundaries
571  forAll(psi_.boundaryField(), patchI)
572  {
573  const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
574 
575  if (ptf.size())
576  {
577  const labelUList& pa = lduAddr().patchAddr(patchI);
578  Field<Type>& iCoeffs = internalCoeffs_[patchI];
579 
580  if (ptf.coupled())
581  {
582  forAll(pa, face)
583  {
584  D[pa[face]] -= component(iCoeffs[face], 0);
585  }
586  }
587  }
588  }
590  // Finally add the relaxation contribution to the source.
591  S += (D - D0)*psi_.primitiveField();
592 }
593 
594 
595 template<class Type>
597 {
598  if (psi_.mesh().relaxEquation(psi_.name()))
599  {
600  relax(psi_.mesh().equationRelaxationFactor(psi_.name()));
601  }
602  else
603  {
605  << "Relaxation factor for field " << psi_.name()
606  << " not found. Relaxation will not be used." << endl;
607  }
608 }
609 
610 
611 template<class Type>
613 {
614  auto tdiag = tmp<scalarField>::New(diag());
615  addCmptAvBoundaryDiag(tdiag.ref());
616  return tdiag;
617 }
618 
619 
620 template<class Type>
622 {
623  auto tAphi = areaScalarField::New
624  (
625  "A(" + psi_.name() + ')',
626  psi_.mesh(),
627  dimensions_/psi_.dimensions()/dimArea,
629  );
630 
631  tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().S();
632  tAphi.ref().correctBoundaryConditions();
634  return tAphi;
635 }
636 
637 
638 template<class Type>
641 {
643  (
644  "H(" + psi_.name() + ')',
645  psi_.mesh(),
646  dimensions_/dimArea,
648  );
649  auto& Hphi = tHphi.ref();
650 
651  // Loop over field components
652  for (direction cmpt=0; cmpt<Type::nComponents; ++cmpt)
653  {
654  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
655 
656  scalarField boundaryDiagCmpt(psi_.size(), Zero);
657  addBoundaryDiag(boundaryDiagCmpt, cmpt);
658  boundaryDiagCmpt.negate();
659  addCmptAvBoundaryDiag(boundaryDiagCmpt);
660 
661  Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
662  }
663 
664  Hphi.primitiveFieldRef() += lduMatrix::H(psi_.primitiveField()) + source_;
665  addBoundarySource(Hphi.primitiveFieldRef());
666 
667  Hphi.primitiveFieldRef() /= psi_.mesh().S();
670  return tHphi;
671 }
672 
673 
674 template<class Type>
677 {
678  if (!psi_.mesh().fluxRequired(psi_.name()))
679  {
681  << "flux requested but " << psi_.name()
682  << " not specified in the fluxRequired sub-dictionary of faSchemes"
683  << abort(FatalError);
684  }
685 
687  (
688  "flux(" + psi_.name() + ')',
689  psi_.mesh(),
690  dimensions()
691  );
692  auto& fieldFlux = tfieldFlux.ref();
693 
694  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
695  {
696  fieldFlux.primitiveFieldRef().replace
697  (
698  cmpt,
699  lduMatrix::faceH(psi_.primitiveField().component(cmpt))
700  );
701  }
702 
703  FieldField<Field, Type> InternalContrib = internalCoeffs_;
704 
705  forAll(InternalContrib, patchI)
706  {
707  InternalContrib[patchI] =
709  (
710  InternalContrib[patchI],
711  psi_.boundaryField()[patchI].patchInternalField()
712  );
713  }
714 
715  FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
716 
717  forAll(NeighbourContrib, patchI)
718  {
719  if (psi_.boundaryField()[patchI].coupled())
720  {
721  NeighbourContrib[patchI] =
723  (
724  NeighbourContrib[patchI],
725  psi_.boundaryField()[patchI].patchNeighbourField()
726  );
727  }
728  }
729 
730  forAll(fieldFlux.boundaryField(), patchI)
731  {
732  fieldFlux.boundaryFieldRef()[patchI] =
733  InternalContrib[patchI] - NeighbourContrib[patchI];
734  }
735 
736  if (faceFluxCorrectionPtr_)
737  {
738  fieldFlux += *faceFluxCorrectionPtr_;
739  }
740 
741  return tfieldFlux;
742 }
743 
744 
745 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
746 
747 template<class Type>
749 {
750  if (this == &famv)
751  {
752  return; // Self-assignment is a no-op
753  }
754 
755  if (&psi_ != &(famv.psi_))
756  {
758  << "different fields"
759  << abort(FatalError);
760  }
761 
762  lduMatrix::operator=(famv);
763  source_ = famv.source_;
764  internalCoeffs_ = famv.internalCoeffs_;
765  boundaryCoeffs_ = famv.boundaryCoeffs_;
766 
767  if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
768  {
769  *faceFluxCorrectionPtr_ = *famv.faceFluxCorrectionPtr_;
770  }
771  else if (famv.faceFluxCorrectionPtr_)
772  {
773  faceFluxCorrectionPtr_ =
775  (*famv.faceFluxCorrectionPtr_);
776  }
777 }
778 
779 
780 template<class Type>
781 void Foam::faMatrix<Type>::operator=(const tmp<faMatrix<Type>>& tfamv)
782 {
783  operator=(tfamv());
784  tfamv.clear();
785 }
786 
787 
788 template<class Type>
790 {
792  source_.negate();
793  internalCoeffs_.negate();
794  boundaryCoeffs_.negate();
795 
796  if (faceFluxCorrectionPtr_)
797  {
798  faceFluxCorrectionPtr_->negate();
799  }
800 }
801 
802 
803 template<class Type>
804 void Foam::faMatrix<Type>::operator+=(const faMatrix<Type>& famv)
805 {
806  checkMethod(*this, famv, "+=");
807 
808  dimensions_ += famv.dimensions_;
809  lduMatrix::operator+=(famv);
810  source_ += famv.source_;
811  internalCoeffs_ += famv.internalCoeffs_;
812  boundaryCoeffs_ += famv.boundaryCoeffs_;
813 
814  if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
815  {
816  *faceFluxCorrectionPtr_ += *famv.faceFluxCorrectionPtr_;
817  }
818  else if (famv.faceFluxCorrectionPtr_)
819  {
820  faceFluxCorrectionPtr_ = new
821  GeometricField<Type, faePatchField, edgeMesh>
822  (
823  *famv.faceFluxCorrectionPtr_
824  );
825  }
826 }
827 
828 
829 template<class Type>
830 void Foam::faMatrix<Type>::operator+=(const tmp<faMatrix<Type>>& tfamv)
831 {
832  operator+=(tfamv());
833  tfamv.clear();
834 }
835 
836 
837 template<class Type>
839 {
840  checkMethod(*this, famv, "+=");
841 
842  dimensions_ -= famv.dimensions_;
843  lduMatrix::operator-=(famv);
844  source_ -= famv.source_;
845  internalCoeffs_ -= famv.internalCoeffs_;
846  boundaryCoeffs_ -= famv.boundaryCoeffs_;
847 
848  if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
849  {
850  *faceFluxCorrectionPtr_ -= *famv.faceFluxCorrectionPtr_;
851  }
852  else if (famv.faceFluxCorrectionPtr_)
853  {
854  faceFluxCorrectionPtr_ =
856  (-*famv.faceFluxCorrectionPtr_);
857  }
858 }
859 
860 
861 template<class Type>
862 void Foam::faMatrix<Type>::operator-=(const tmp<faMatrix<Type>>& tfamv)
863 {
864  operator-=(tfamv());
865  tfamv.clear();
866 }
867 
868 
869 template<class Type>
871 (
873 )
874 {
875  checkMethod(*this, su, "+=");
876  source() -= su.mesh().S()*su.field();
877 }
878 
879 
880 template<class Type>
882 (
884 )
885 {
886  operator+=(tsu());
887  tsu.clear();
888 }
889 
890 
891 template<class Type>
893 (
895 )
896 {
897  operator+=(tsu());
898  tsu.clear();
899 }
900 
901 
902 template<class Type>
904 (
906 )
907 {
908  checkMethod(*this, su, "-=");
909  source() += su.mesh().S()*su.field();
910 }
911 
912 
913 template<class Type>
915 (
917 )
918 {
919  operator-=(tsu());
920  tsu.clear();
921 }
922 
923 
924 template<class Type>
926 (
928 )
929 {
930  operator-=(tsu());
931  tsu.clear();
932 }
933 
934 
935 template<class Type>
937 (
938  const dimensioned<Type>& su
939 )
940 {
941  source() -= psi().mesh().S()*su;
942 }
943 
944 
945 template<class Type>
947 (
948  const dimensioned<Type>& su
949 )
950 {
951  source() += psi().mesh().S()*su;
952 }
953 
955 template<class Type>
957 {}
958 
959 
960 template<class Type>
962 {}
963 
964 
965 template<class Type>
967 (
968  const areaScalarField::Internal& dsf
969 )
970 {
971  dimensions_ *= dsf.dimensions();
972  lduMatrix::operator*=(dsf.field());
973  source_ *= dsf.field();
974 
975  forAll(boundaryCoeffs_, patchi)
976  {
977  const scalarField pisf
978  (
979  dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
980  );
981  internalCoeffs_[patchi] *= pisf;
982  boundaryCoeffs_[patchi] *= pisf;
983  }
984 
985  if (faceFluxCorrectionPtr_)
986  {
988  << "cannot scale a matrix containing a faceFluxCorrection"
990  }
991 }
992 
993 
994 template<class Type>
996 (
997  const tmp<areaScalarField::Internal>& tfld
998 )
999 {
1000  operator*=(tfld());
1001  tfld.clear();
1002 }
1003 
1004 
1005 template<class Type>
1007 (
1008  const tmp<areaScalarField>& tfld
1009 )
1010 {
1011  operator*=(tfld());
1012  tfld.clear();
1013 }
1014 
1015 
1016 template<class Type>
1018 (
1019  const dimensioned<scalar>& ds
1020 )
1021 {
1022  dimensions_ *= ds.dimensions();
1023  lduMatrix::operator*=(ds.value());
1024  source_ *= ds.value();
1025  internalCoeffs_ *= ds.value();
1026  boundaryCoeffs_ *= ds.value();
1027 
1028  if (faceFluxCorrectionPtr_)
1029  {
1030  *faceFluxCorrectionPtr_ *= ds.value();
1031  }
1033 
1034 
1035 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
1036 
1037 template<class Type>
1038 void Foam::checkMethod
1039 (
1040  const faMatrix<Type>& mat1,
1041  const faMatrix<Type>& mat2,
1042  const char* op
1043 )
1044 {
1045  if (&mat1.psi() != &mat2.psi())
1046  {
1048  << "Incompatible fields for operation\n "
1049  << "[" << mat1.psi().name() << "] "
1050  << op
1051  << " [" << mat2.psi().name() << "]"
1052  << abort(FatalError);
1053  }
1054 
1055  if
1056  (
1058  && mat1.dimensions() != mat2.dimensions()
1059  )
1060  {
1062  << "Incompatible dimensions for operation\n "
1063  << "[" << mat1.psi().name() << mat1.dimensions()/dimArea << " ] "
1064  << op
1065  << " [" << mat2.psi().name() << mat2.dimensions()/dimArea << " ]"
1067  }
1068 }
1069 
1070 
1071 template<class Type>
1072 void Foam::checkMethod
1073 (
1074  const faMatrix<Type>& mat,
1075  const DimensionedField<Type, areaMesh>& fld,
1076  const char* op
1077 )
1078 {
1079  if
1080  (
1082  && mat.dimensions()/dimArea != fld.dimensions()
1083  )
1084  {
1086  << "Incompatible dimensions for operation\n "
1087  << "[" << mat.psi().name() << mat.dimensions()/dimArea << " ] "
1088  << op
1089  << " [" << fld.name() << fld.dimensions() << " ]"
1091  }
1092 }
1093 
1094 
1095 template<class Type>
1096 void Foam::checkMethod
1097 (
1098  const faMatrix<Type>& mat,
1099  const dimensioned<Type>& dt,
1100  const char* op
1101 )
1102 {
1103  if
1104  (
1106  && mat.dimensions()/dimArea != dt.dimensions()
1107  )
1108  {
1110  << "Incompatible dimensions for operation\n "
1111  << "[" << mat.psi().name() << mat.dimensions()/dimArea << " ] "
1112  << op
1113  << " [" << dt.name() << dt.dimensions() << " ]"
1114  << abort(FatalError);
1115  }
1116 }
1117 
1118 
1119 template<class Type>
1121 (
1122  faMatrix<Type>& fam,
1123  const dictionary& solverControls
1124 )
1125 {
1126  return fam.solve(solverControls);
1127 }
1128 
1129 
1130 template<class Type>
1132 (
1133  const tmp<faMatrix<Type>>& tmat,
1134  const dictionary& solverControls
1135 )
1136 {
1137  SolverPerformance<Type> solverPerf(tmat.constCast().solve(solverControls));
1138 
1139  tmat.clear();
1140 
1141  return solverPerf;
1142 }
1143 
1144 
1145 template<class Type>
1147 {
1148  return fam.solve();
1149 }
1150 
1151 
1152 template<class Type>
1153 Foam::SolverPerformance<Type> Foam::solve(const tmp<faMatrix<Type>>& tmat)
1154 {
1155  SolverPerformance<Type> solverPerf(tmat.constCast().solve());
1156 
1157  tmat.clear();
1158 
1159  return solverPerf;
1160 }
1161 
1162 
1163 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
1164 
1165 template<class Type>
1166 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1167 (
1168  const faMatrix<Type>& A,
1169  const faMatrix<Type>& B
1170 )
1171 {
1172  checkMethod(A, B, "+");
1173  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1174  tC.ref() += B;
1175  return tC;
1176 }
1177 
1178 
1179 template<class Type>
1180 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1181 (
1182  const tmp<faMatrix<Type>>& tA,
1183  const faMatrix<Type>& B
1184 )
1185 {
1186  checkMethod(tA(), B, "+");
1187  tmp<faMatrix<Type>> tC(tA.ptr());
1188  tC.ref() += B;
1189  return tC;
1190 }
1191 
1192 
1193 template<class Type>
1194 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1195 (
1196  const faMatrix<Type>& A,
1197  const tmp<faMatrix<Type>>& tB
1198 )
1199 {
1200  checkMethod(A, tB(), "+");
1201  tmp<faMatrix<Type>> tC(tB.ptr());
1202  tC.ref() += A;
1203  return tC;
1204 }
1205 
1206 
1207 template<class Type>
1208 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1209 (
1210  const tmp<faMatrix<Type>>& tA,
1211  const tmp<faMatrix<Type>>& tB
1212 )
1213 {
1214  checkMethod(tA(), tB(), "+");
1215  tmp<faMatrix<Type>> tC(tA.ptr());
1216  tC.ref() += tB();
1217  tB.clear();
1218  return tC;
1219 }
1220 
1221 
1222 template<class Type>
1223 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1224 (
1225  const faMatrix<Type>& A
1226 )
1227 {
1228  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1229  tC.ref().negate();
1230  return tC;
1231 }
1232 
1233 
1234 template<class Type>
1235 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1236 (
1237  const tmp<faMatrix<Type>>& tA
1238 )
1239 {
1240  tmp<faMatrix<Type>> tC(tA.ptr());
1241  tC.ref().negate();
1242  return tC;
1243 }
1244 
1245 
1246 template<class Type>
1247 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1248 (
1249  const faMatrix<Type>& A,
1250  const faMatrix<Type>& B
1251 )
1252 {
1253  checkMethod(A, B, "-");
1254  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1255  tC.ref() -= B;
1256  return tC;
1257 }
1258 
1259 
1260 template<class Type>
1261 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1262 (
1263  const tmp<faMatrix<Type>>& tA,
1264  const faMatrix<Type>& B
1265 )
1266 {
1267  checkMethod(tA(), B, "-");
1268  tmp<faMatrix<Type>> tC(tA.ptr());
1269  tC.ref() -= B;
1270  return tC;
1271 }
1272 
1273 
1274 template<class Type>
1275 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1276 (
1277  const faMatrix<Type>& A,
1278  const tmp<faMatrix<Type>>& tB
1279 )
1280 {
1281  checkMethod(A, tB(), "-");
1282  tmp<faMatrix<Type>> tC(tB.ptr());
1283  tC.ref() -= A;
1284  tC.ref().negate();
1285  return tC;
1286 }
1287 
1288 
1289 template<class Type>
1290 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1291 (
1292  const tmp<faMatrix<Type>>& tA,
1293  const tmp<faMatrix<Type>>& tB
1294 )
1295 {
1296  checkMethod(tA(), tB(), "-");
1297  tmp<faMatrix<Type>> tC(tA.ptr());
1298  tC.ref() -= tB();
1299  tB.clear();
1300  return tC;
1301 }
1302 
1303 
1304 template<class Type>
1305 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1306 (
1307  const faMatrix<Type>& A,
1308  const faMatrix<Type>& B
1309 )
1310 {
1311  checkMethod(A, B, "==");
1312  return (A - B);
1313 }
1314 
1315 
1316 template<class Type>
1317 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1318 (
1319  const tmp<faMatrix<Type>>& tA,
1320  const faMatrix<Type>& B
1321 )
1322 {
1323  checkMethod(tA(), B, "==");
1324  return (tA - B);
1325 }
1326 
1327 
1328 template<class Type>
1329 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1330 (
1331  const faMatrix<Type>& A,
1332  const tmp<faMatrix<Type>>& tB
1333 )
1334 {
1335  checkMethod(A, tB(), "==");
1336  return (A - tB);
1337 }
1338 
1339 
1340 template<class Type>
1341 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1342 (
1343  const tmp<faMatrix<Type>>& tA,
1344  const tmp<faMatrix<Type>>& tB
1345 )
1346 {
1347  checkMethod(tA(), tB(), "==");
1348  return (tA - tB);
1349 }
1350 
1351 
1352 template<class Type>
1353 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1354 (
1355  const faMatrix<Type>& A,
1356  const DimensionedField<Type, areaMesh>& su
1357 )
1358 {
1359  checkMethod(A, su, "+");
1360  auto tC(A.clone());
1361  tC.ref().source() -= su.mesh().S()*su.field();
1362  return tC;
1363 }
1364 
1365 
1366 template<class Type>
1367 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1368 (
1369  const faMatrix<Type>& A,
1370  const tmp<DimensionedField<Type, areaMesh>>& tsu
1371 )
1372 {
1373  checkMethod(A, tsu(), "+");
1374  auto tC(A.clone());
1375  tC.ref().source() -= tsu().mesh().S()*tsu().field();
1376  tsu.clear();
1377  return tC;
1378 }
1379 
1380 
1381 template<class Type>
1382 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1383 (
1384  const faMatrix<Type>& A,
1385  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1386 )
1387 {
1388  checkMethod(A, tsu(), "+");
1389  auto tC(A.clone());
1390  tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1391  tsu.clear();
1392  return tC;
1393 }
1394 
1395 
1396 template<class Type>
1397 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1398 (
1399  const tmp<faMatrix<Type>>& tA,
1400  const DimensionedField<Type, areaMesh>& su
1401 )
1402 {
1403  checkMethod(tA(), su, "+");
1404  tmp<faMatrix<Type>> tC(tA.ptr());
1405  tC.ref().source() -= su.mesh().S()*su.field();
1406  return tC;
1407 }
1408 
1409 
1410 template<class Type>
1411 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1412 (
1413  const tmp<faMatrix<Type>>& tA,
1414  const tmp<DimensionedField<Type, areaMesh>>& tsu
1415 )
1416 {
1417  checkMethod(tA(), tsu(), "+");
1418  tmp<faMatrix<Type>> tC(tA.ptr());
1419  tC.ref().source() -= tsu().mesh().S()*tsu().field();
1420  tsu.clear();
1421  return tC;
1422 }
1423 
1424 
1425 template<class Type>
1426 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1427 (
1428  const tmp<faMatrix<Type>>& tA,
1429  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1430 )
1431 {
1432  checkMethod(tA(), tsu(), "+");
1433  tmp<faMatrix<Type>> tC(tA.ptr());
1434  tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1435  tsu.clear();
1436  return tC;
1437 }
1438 
1439 
1440 template<class Type>
1441 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1442 (
1443  const DimensionedField<Type, areaMesh>& su,
1444  const faMatrix<Type>& A
1445 )
1446 {
1447  checkMethod(A, su, "+");
1448  auto tC(A.clone());
1449  tC.ref().source() -= su.mesh().S()*su.field();
1450  return tC;
1451 }
1452 
1453 
1454 template<class Type>
1455 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1456 (
1457  const tmp<DimensionedField<Type, areaMesh>>& tsu,
1458  const faMatrix<Type>& A
1459 )
1460 {
1461  checkMethod(A, tsu(), "+");
1462  auto tC(A.clone());
1463  tC.ref().source() -= tsu().mesh().S()*tsu().field();
1464  tsu.clear();
1465  return tC;
1466 }
1467 
1468 
1469 template<class Type>
1470 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1471 (
1472  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1473  const faMatrix<Type>& A
1474 )
1475 {
1476  checkMethod(A, tsu(), "+");
1477  auto tC(A.clone());
1478  tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1479  tsu.clear();
1480  return tC;
1481 }
1482 
1483 
1484 template<class Type>
1485 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1486 (
1487  const DimensionedField<Type, areaMesh>& su,
1488  const tmp<faMatrix<Type>>& tA
1489 )
1490 {
1491  checkMethod(tA(), su, "+");
1492  tmp<faMatrix<Type>> tC(tA.ptr());
1493  tC.ref().source() -= su.mesh().S()*su.field();
1494  return tC;
1495 }
1496 
1497 
1498 template<class Type>
1499 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1500 (
1501  const tmp<DimensionedField<Type, areaMesh>>& tsu,
1502  const tmp<faMatrix<Type>>& tA
1503 )
1504 {
1505  checkMethod(tA(), tsu(), "+");
1506  tmp<faMatrix<Type>> tC(tA.ptr());
1507  tC.ref().source() -= tsu().mesh().S()*tsu().field();
1508  tsu.clear();
1509  return tC;
1510 }
1511 
1512 
1513 template<class Type>
1514 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1515 (
1516  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1517  const tmp<faMatrix<Type>>& tA
1518 )
1519 {
1520  checkMethod(tA(), tsu(), "+");
1521  tmp<faMatrix<Type>> tC(tA.ptr());
1522  tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1523  tsu.clear();
1524  return tC;
1525 }
1526 
1527 
1528 template<class Type>
1529 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1530 (
1531  const faMatrix<Type>& A,
1532  const DimensionedField<Type, areaMesh>& su
1533 )
1534 {
1535  checkMethod(A, su, "-");
1536  auto tC(A.clone());
1537  tC.ref().source() += su.mesh().S()*su.field();
1538  return tC;
1539 }
1540 
1541 
1542 template<class Type>
1543 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1544 (
1545  const faMatrix<Type>& A,
1546  const tmp<DimensionedField<Type, areaMesh>>& tsu
1547 )
1548 {
1549  checkMethod(A, tsu(), "-");
1550  auto tC(A.clone());
1551  tC.ref().source() += tsu().mesh().S()*tsu().field();
1552  tsu.clear();
1553  return tC;
1554 }
1555 
1556 
1557 template<class Type>
1558 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1559 (
1560  const faMatrix<Type>& A,
1561  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1562 )
1563 {
1564  checkMethod(A, tsu(), "-");
1565  auto tC(A.clone());
1566  tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1567  tsu.clear();
1568  return tC;
1569 }
1570 
1571 
1572 template<class Type>
1573 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1574 (
1575  const tmp<faMatrix<Type>>& tA,
1576  const DimensionedField<Type, areaMesh>& su
1577 )
1578 {
1579  checkMethod(tA(), su, "-");
1580  tmp<faMatrix<Type>> tC(tA.ptr());
1581  tC.ref().source() += su.mesh().S()*su.field();
1582  return tC;
1583 }
1584 
1585 
1586 template<class Type>
1587 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1588 (
1589  const tmp<faMatrix<Type>>& tA,
1590  const tmp<DimensionedField<Type, areaMesh>>& tsu
1591 )
1592 {
1593  checkMethod(tA(), tsu(), "-");
1594  tmp<faMatrix<Type>> tC(tA.ptr());
1595  tC.ref().source() += tsu().mesh().S()*tsu().field();
1596  tsu.clear();
1597  return tC;
1598 }
1599 
1600 
1601 template<class Type>
1602 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1603 (
1604  const tmp<faMatrix<Type>>& tA,
1605  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1606 )
1607 {
1608  checkMethod(tA(), tsu(), "-");
1609  tmp<faMatrix<Type>> tC(tA.ptr());
1610  tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1611  tsu.clear();
1612  return tC;
1613 }
1614 
1615 
1616 template<class Type>
1617 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1618 (
1619  const DimensionedField<Type, areaMesh>& su,
1620  const faMatrix<Type>& A
1621 )
1622 {
1623  checkMethod(A, su, "-");
1624  auto tC(A.clone());
1625  tC.ref().negate();
1626  tC.ref().source() -= su.mesh().S()*su.field();
1627  return tC;
1628 }
1629 
1630 
1631 template<class Type>
1632 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1633 (
1634  const tmp<DimensionedField<Type, areaMesh>>& tsu,
1635  const faMatrix<Type>& A
1636 )
1637 {
1638  checkMethod(A, tsu(), "-");
1639  auto tC(A.clone());
1640  tC.ref().negate();
1641  tC.ref().source() -= tsu().mesh().S()*tsu().field();
1642  tsu.clear();
1643  return tC;
1644 }
1645 
1646 
1647 template<class Type>
1648 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1649 (
1650  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1651  const faMatrix<Type>& A
1652 )
1653 {
1654  checkMethod(A, tsu(), "-");
1655  auto tC(A.clone());
1656  tC.ref().negate();
1657  tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1658  tsu.clear();
1659  return tC;
1660 }
1661 
1662 
1663 template<class Type>
1664 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1665 (
1666  const DimensionedField<Type, areaMesh>& su,
1667  const tmp<faMatrix<Type>>& tA
1668 )
1669 {
1670  checkMethod(tA(), su, "-");
1671  tmp<faMatrix<Type>> tC(tA.ptr());
1672  tC.ref().negate();
1673  tC.ref().source() -= su.mesh().S()*su.field();
1674  return tC;
1675 }
1676 
1677 
1678 template<class Type>
1679 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1680 (
1681  const tmp<DimensionedField<Type, areaMesh>>& tsu,
1682  const tmp<faMatrix<Type>>& tA
1683 )
1684 {
1685  checkMethod(tA(), tsu(), "-");
1686  tmp<faMatrix<Type>> tC(tA.ptr());
1687  tC.ref().negate();
1688  tC.ref().source() -= tsu().mesh().S()*tsu().field();
1689  tsu.clear();
1690  return tC;
1691 }
1692 
1693 
1694 template<class Type>
1695 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1696 (
1697  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1698  const tmp<faMatrix<Type>>& tA
1699 )
1700 {
1701  checkMethod(tA(), tsu(), "-");
1702  tmp<faMatrix<Type>> tC(tA.ptr());
1703  tC.ref().negate();
1704  tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1705  tsu.clear();
1706  return tC;
1707 }
1708 
1709 
1710 template<class Type>
1711 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1712 (
1713  const faMatrix<Type>& A,
1714  const dimensioned<Type>& su
1715 )
1716 {
1717  checkMethod(A, su, "+");
1718  auto tC(A.clone());
1719  tC.ref().source() -= su.value()*A.psi().mesh().S();
1720  return tC;
1721 }
1722 
1723 
1724 template<class Type>
1725 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1726 (
1727  const tmp<faMatrix<Type>>& tA,
1728  const dimensioned<Type>& su
1729 )
1730 {
1731  checkMethod(tA(), su, "+");
1732  tmp<faMatrix<Type>> tC(tA.ptr());
1733  tC.ref().source() -= su.value()*tC().psi().mesh().S();
1734  return tC;
1735 }
1736 
1737 
1738 template<class Type>
1739 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1740 (
1741  const dimensioned<Type>& su,
1742  const faMatrix<Type>& A
1743 )
1744 {
1745  checkMethod(A, su, "+");
1746  auto tC(A.clone());
1747  tC.ref().source() -= su.value()*A.psi().mesh().S();
1748  return tC;
1749 }
1750 
1751 
1752 template<class Type>
1753 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1754 (
1755  const dimensioned<Type>& su,
1756  const tmp<faMatrix<Type>>& tA
1757 )
1758 {
1759  checkMethod(tA(), su, "+");
1760  tmp<faMatrix<Type>> tC(tA.ptr());
1761  tC.ref().source() -= su.value()*tC().psi().mesh().S();
1762  return tC;
1763 }
1764 
1765 
1766 template<class Type>
1767 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1768 (
1769  const faMatrix<Type>& A,
1770  const dimensioned<Type>& su
1771 )
1772 {
1773  checkMethod(A, su, "-");
1774  auto tC(A.clone());
1775  tC.ref().source() += su.value()*tC().psi().mesh().S();
1776  return tC;
1777 }
1778 
1779 
1780 template<class Type>
1781 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1782 (
1783  const tmp<faMatrix<Type>>& tA,
1784  const dimensioned<Type>& su
1785 )
1786 {
1787  checkMethod(tA(), su, "-");
1788  tmp<faMatrix<Type>> tC(tA.ptr());
1789  tC.ref().source() += su.value()*tC().psi().mesh().S();
1790  return tC;
1791 }
1792 
1793 
1794 template<class Type>
1795 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1796 (
1797  const dimensioned<Type>& su,
1798  const faMatrix<Type>& A
1799 )
1800 {
1801  checkMethod(A, su, "-");
1802  auto tC(A.clone());
1803  tC.ref().negate();
1804  tC.ref().source() -= su.value()*A.psi().mesh().S();
1805  return tC;
1806 }
1807 
1808 
1809 template<class Type>
1810 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1811 (
1812  const dimensioned<Type>& su,
1813  const tmp<faMatrix<Type>>& tA
1814 )
1815 {
1816  checkMethod(tA(), su, "-");
1817  tmp<faMatrix<Type>> tC(tA.ptr());
1818  tC.ref().negate();
1819  tC.ref().source() -= su.value()*tC().psi().mesh().S();
1820  return tC;
1821 }
1822 
1823 
1824 template<class Type>
1825 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1826 (
1827  const faMatrix<Type>& A,
1828  const DimensionedField<Type, areaMesh>& su
1829 )
1830 {
1831  checkMethod(A, su, "==");
1832  auto tC(A.clone());
1833  tC.ref().source() += su.mesh().S()*su.field();
1834  return tC;
1835 }
1836 
1837 
1838 template<class Type>
1839 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1840 (
1841  const faMatrix<Type>& A,
1842  const tmp<DimensionedField<Type, areaMesh>>& tsu
1843 )
1844 {
1845  checkMethod(A, tsu(), "==");
1846  auto tC(A.clone());
1847  tC.ref().source() += tsu().mesh().S()*tsu().field();
1848  tsu.clear();
1849  return tC;
1850 }
1851 
1852 
1853 template<class Type>
1854 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1855 (
1856  const faMatrix<Type>& A,
1857  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1858 )
1859 {
1860  checkMethod(A, tsu(), "==");
1861  auto tC(A.clone());
1862  tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1863  tsu.clear();
1864  return tC;
1865 }
1866 
1867 
1868 template<class Type>
1869 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1870 (
1871  const tmp<faMatrix<Type>>& tA,
1872  const DimensionedField<Type, areaMesh>& su
1873 )
1874 {
1875  checkMethod(tA(), su, "==");
1876  tmp<faMatrix<Type>> tC(tA.ptr());
1877  tC.ref().source() += su.mesh().S()*su.field();
1878  return tC;
1879 }
1880 
1881 
1882 template<class Type>
1883 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1884 (
1885  const tmp<faMatrix<Type>>& tA,
1886  const tmp<DimensionedField<Type, areaMesh>>& tsu
1887 )
1888 {
1889  checkMethod(tA(), tsu(), "==");
1890  tmp<faMatrix<Type>> tC(tA.ptr());
1891  tC.ref().source() += tsu().mesh().S()*tsu().field();
1892  tsu.clear();
1893  return tC;
1894 }
1895 
1896 
1897 template<class Type>
1898 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1899 (
1900  const tmp<faMatrix<Type>>& tA,
1901  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1902 )
1903 {
1904  checkMethod(tA(), tsu(), "==");
1905  tmp<faMatrix<Type>> tC(tA.ptr());
1906  tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1907  tsu.clear();
1908  return tC;
1909 }
1910 
1911 
1912 template<class Type>
1913 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1914 (
1915  const faMatrix<Type>& A,
1916  const dimensioned<Type>& su
1917 )
1918 {
1919  checkMethod(A, su, "==");
1920  auto tC(A.clone());
1921  tC.ref().source() += A.psi().mesh().S()*su.value();
1922  return tC;
1923 }
1924 
1925 
1926 template<class Type>
1927 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1928 (
1929  const tmp<faMatrix<Type>>& tA,
1930  const dimensioned<Type>& su
1931 )
1932 {
1933  checkMethod(tA(), su, "==");
1934  tmp<faMatrix<Type>> tC(tA.ptr());
1935  tC.ref().source() += tC().psi().mesh().S()*su.value();
1936  return tC;
1937 }
1938 
1939 
1940 template<class Type>
1941 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1942 (
1943  const faMatrix<Type>& A,
1944  const Foam::zero
1945 )
1946 {
1947  return A;
1948 }
1949 
1950 
1951 template<class Type>
1952 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1953 (
1954  const tmp<faMatrix<Type>>& tA,
1955  const Foam::zero
1956 )
1957 {
1958  return tA;
1959 }
1960 
1961 
1962 template<class Type>
1963 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1964 (
1965  const areaScalarField::Internal& dsf,
1966  const faMatrix<Type>& A
1967 )
1968 {
1969  auto tC(A.clone());
1970  tC.ref() *= dsf;
1971  return tC;
1972 }
1973 
1974 
1975 template<class Type>
1976 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1977 (
1978  const tmp<areaScalarField::Internal>& tdsf,
1979  const faMatrix<Type>& A
1980 )
1981 {
1982  auto tC(A.clone());
1983  tC.ref() *= tdsf;
1984  return tC;
1985 }
1986 
1987 
1988 template<class Type>
1989 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1990 (
1991  const tmp<areaScalarField>& tvsf,
1992  const faMatrix<Type>& A
1993 )
1994 {
1995  auto tC(A.clone());
1996  tC.ref() *= tvsf;
1997  return tC;
1998 }
1999 
2000 
2001 template<class Type>
2002 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2003 (
2004  const areaScalarField::Internal& dsf,
2005  const tmp<faMatrix<Type>>& tA
2006 )
2007 {
2008  tmp<faMatrix<Type>> tC(tA.ptr());
2009  tC.ref() *= dsf;
2010  return tC;
2011 }
2012 
2013 
2014 template<class Type>
2015 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2016 (
2017  const tmp<areaScalarField::Internal>& tdsf,
2018  const tmp<faMatrix<Type>>& tA
2019 )
2020 {
2021  tmp<faMatrix<Type>> tC(tA.ptr());
2022  tC.ref() *= tdsf;
2023  return tC;
2024 }
2025 
2026 
2027 template<class Type>
2028 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2029 (
2030  const tmp<areaScalarField>& tvsf,
2031  const tmp<faMatrix<Type>>& tA
2032 )
2033 {
2034  tmp<faMatrix<Type>> tC(tA.ptr());
2035  tC.ref() *= tvsf;
2036  return tC;
2037 }
2038 
2039 
2040 template<class Type>
2041 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2042 (
2043  const dimensioned<scalar>& ds,
2044  const faMatrix<Type>& A
2045 )
2046 {
2047  auto tC(A.clone());
2048  tC.ref() *= ds;
2049  return tC;
2050 }
2051 
2052 
2053 template<class Type>
2054 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2055 (
2056  const dimensioned<scalar>& ds,
2057  const tmp<faMatrix<Type>>& tA
2058 )
2059 {
2060  tmp<faMatrix<Type>> tC(tA.ptr());
2061  tC.ref() *= ds;
2062  return tC;
2063 }
2064 
2065 
2066 template<class Type>
2068 Foam::operator&
2069 (
2070  const faMatrix<Type>& M,
2071  const DimensionedField<Type, areaMesh>& psi
2072 )
2073 {
2074  auto tMphi = tmp<GeometricField<Type, faPatchField, areaMesh>>::New
2075  (
2076  IOobject
2077  (
2078  "M&" + psi.name(),
2079  psi.instance(),
2080  psi.mesh().mesh(),
2083  ),
2084  psi.mesh(),
2085  M.dimensions()/dimArea
2086  );
2087  auto& Mphi = tMphi.ref();
2088 
2089  // Loop over field components
2090  if (M.hasDiag())
2091  {
2092  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2093  {
2094  scalarField psiCmpt(psi.field().component(cmpt));
2095  scalarField boundaryDiagCmpt(M.diag());
2096  M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2097  Mphi.primitiveFieldRef().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2098  }
2099  }
2100  else
2101  {
2102  Mphi.primitiveFieldRef() = Zero;
2103  }
2104 
2105  Mphi.primitiveFieldRef() += M.lduMatrix::H(psi.field()) + M.source();
2106  M.addBoundarySource(Mphi.primitiveFieldRef());
2107 
2108  Mphi.primitiveFieldRef() /= -psi.mesh().S();
2109  Mphi.correctBoundaryConditions();
2110 
2111  return tMphi;
2112 }
2113 
2114 
2115 template<class Type>
2117 Foam::operator&
2118 (
2119  const faMatrix<Type>& M,
2120  const tmp<DimensionedField<Type, areaMesh>>& tpsi
2121 )
2122 {
2123  tmp<GeometricField<Type, faPatchField, areaMesh>> tMpsi = M & tpsi();
2124  tpsi.clear();
2125  return tMpsi;
2126 }
2127 
2128 
2129 template<class Type>
2131 Foam::operator&
2132 (
2133  const faMatrix<Type>& M,
2134  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tpsi
2135 )
2136 {
2137  tmp<GeometricField<Type, faPatchField, areaMesh>> tMpsi = M & tpsi();
2138  tpsi.clear();
2139  return tMpsi;
2140 }
2141 
2142 
2143 template<class Type>
2145 Foam::operator&
2146 (
2147  const tmp<faMatrix<Type>>& tM,
2148  const DimensionedField<Type, areaMesh>& psi
2149 )
2150 {
2151  tmp<GeometricField<Type, faPatchField, areaMesh>> tMpsi = tM() & psi;
2152  tM.clear();
2153  return tMpsi;
2154 }
2155 
2156 
2157 template<class Type>
2159 Foam::operator&
2160 (
2161  const tmp<faMatrix<Type>>& tM,
2162  const tmp<DimensionedField<Type, areaMesh>>& tpsi
2163 )
2164 {
2165  tmp<GeometricField<Type, faPatchField, areaMesh>> tMpsi = tM() & tpsi();
2166  tM.clear();
2167  tpsi.clear();
2168  return tMpsi;
2169 }
2170 
2171 
2172 template<class Type>
2174 Foam::operator&
2175 (
2176  const tmp<faMatrix<Type>>& tM,
2177  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tpsi
2178 )
2179 {
2180  tmp<GeometricField<Type, faPatchField, areaMesh>> tMpsi = tM() & tpsi();
2181  tM.clear();
2182  tpsi.clear();
2183  return tMpsi;
2184 }
2185 
2186 
2187 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2188 
2189 template<class Type>
2190 Foam::Ostream& Foam::operator<<(Ostream& os, const faMatrix<Type>& fam)
2191 {
2192  os << static_cast<const lduMatrix&>(fam) << nl
2193  << fam.dimensions_ << nl
2194  << fam.source_ << nl
2195  << fam.internalCoeffs_ << nl
2196  << fam.boundaryCoeffs_ << endl;
2197 
2199 
2200  return os;
2201 }
2202 
2203 
2204 // * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
2205 
2206 #include "faMatrixSolve.C"
2207 
2208 // ************************************************************************* //
void negate()
Inplace negate.
Definition: faMatrix.C:782
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:88
faceListList boundary
void setValues(const labelUList &faceLabels, const Type &value)
Set solution in given faces to the specified value and eliminate the corresponding equations from the...
Definition: faMatrix.C:402
void operator-=(const faMatrix< Type > &)
Definition: faMatrix.C:831
tmp< Field< Type > > faceH(const Field< Type > &) const
tmp< scalarField > D() const
Return the matrix diagonal.
Definition: faMatrix.C:605
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:48
label faceId(-1)
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &f1)
void operator+=(const faMatrix< Type > &)
Definition: faMatrix.C:797
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UEqn relax()
faMatrix(const GeometricField< Type, faPatchField, areaMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition: faMatrix.C:179
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
const labelListList & faceEdges() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:45
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: faMatrix.C:36
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:162
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:244
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: faMatrix.C:111
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: faPatchField.H:197
void setReferences(const labelUList &faceLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: faMatrix.C:454
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1186
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Generic dimensioned Type class.
Ignore writing from objectRegistry::writeObject()
labelList faceLabels(nFaceLabels)
DimensionedField< scalar, areaMesh > Internal
The internal field type from which this GeometricField is derived.
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: faMatrix.C:145
conserve primitiveFieldRef()+
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:549
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
tmp< GeometricField< Type, faePatchField, edgeMesh > > flux() const
Return the face-flux field from the matrix.
Definition: faMatrix.C:669
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: faMatrix.C:74
dynamicFvMesh & mesh
void relax()
Relax matrix (for steady-state solution).
Definition: faMatrix.C:589
void operator=(const faMatrix< Type > &)
Definition: faMatrix.C:741
Generic templated field type.
Definition: Field.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
#define DebugInFunction
Report an information message using Foam::Info.
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1032
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:212
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void setValuesFromList(const labelUList &faceLabels, const ListType< Type > &values)
Set solution in given faces to the specified values.
Definition: faMatrix.C:299
void operator=(const lduMatrix &)
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
const Mesh & mesh() const noexcept
Return mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
void operator*=(const scalarField &)
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
tmp< areaScalarField > A() const
Return the central coefficient.
Definition: faMatrix.C:614
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:233
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:541
void setReference(const label facei, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: faMatrix.C:435
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const Field< Type > & field() const noexcept
Return const-reference to the field values.
Template functions to aid in the implementation of demand driven data.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Finite-Area matrix basic solvers.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1170
static bool checking() noexcept
True if dimension checking is enabled (the usual default)
Definition: dimensionSet.H:241
A List with indirect addressing. Like IndirectList but does not store addressing. ...
Definition: faMatrix.H:56
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1037
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
A single value that is represented as a list with an operator[] to access the value. This can be useful for templated operations expecting a list accessor.
Definition: UniformList.H:48
A special matrix type and solver, designed for finite area solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: faMatricesFwd.H:37
void correctBoundaryConditions()
Correct boundary field.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:397
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:352
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition: tmpI.H:298
const volScalarField & psi
tmp< Field< Type > > H(const Field< Type > &) const
const dimensionedScalar & D
void operator-=(const lduMatrix &)
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: faMatrix.C:129
void operator+=(const lduMatrix &)
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:265
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
void deleteDemandDrivenData(DataPtr &dataPtr)
Calculate the matrix for the second temporal derivative.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
#define M(I)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133
virtual ~faMatrix()
Destructor.
Definition: faMatrix.C:285