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