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-2020 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"
33 #include "UIndirectList.H"
34 #include "UniformList.H"
35 #include "demandDrivenData.H"
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 template<class Type>
40 template<class Type2>
42 (
43  const labelUList& addr,
44  const Field<Type2>& pf,
45  Field<Type2>& intf
46 ) const
47 {
48  if (addr.size() != pf.size())
49  {
51  << "addressing (" << addr.size()
52  << ") and field (" << pf.size() << ") are different sizes" << endl
53  << abort(FatalError);
54  }
55 
56  forAll(addr, faceI)
57  {
58  intf[addr[faceI]] += pf[faceI];
59  }
60 }
61 
62 
63 template<class Type>
64 template<class Type2>
66 (
67  const labelUList& addr,
68  const tmp<Field<Type2>>& tpf,
69  Field<Type2>& intf
70 ) const
71 {
72  addToInternalField(addr, tpf(), intf);
73  tpf.clear();
74 }
75 
76 
77 template<class Type>
78 template<class Type2>
80 (
81  const labelUList& addr,
82  const Field<Type2>& pf,
83  Field<Type2>& intf
84 ) const
85 {
86  if (addr.size() != pf.size())
87  {
89  << "addressing (" << addr.size()
90  << ") and field (" << pf.size() << ") are different sizes" << endl
91  << abort(FatalError);
92  }
93 
94  forAll(addr, faceI)
95  {
96  intf[addr[faceI]] -= pf[faceI];
97  }
98 }
99 
100 
101 template<class Type>
102 template<class Type2>
104 (
105  const labelUList& addr,
106  const tmp<Field<Type2>>& tpf,
107  Field<Type2>& intf
108 ) const
109 {
110  subtractFromInternalField(addr, tpf(), intf);
111  tpf.clear();
112 }
113 
114 
115 template<class Type>
117 (
118  scalarField& diag,
119  const direction solveCmpt
120 ) const
121 {
122  forAll(internalCoeffs_, patchI)
123  {
124  addToInternalField
125  (
126  lduAddr().patchAddr(patchI),
127  internalCoeffs_[patchI].component(solveCmpt),
128  diag
129  );
130  }
131 }
132 
133 
134 template<class Type>
136 {
137  forAll(internalCoeffs_, patchI)
138  {
139  addToInternalField
140  (
141  lduAddr().patchAddr(patchI),
142  cmptAv(internalCoeffs_[patchI]),
143  diag
144  );
145  }
146 }
147 
148 
149 template<class Type>
151 (
152  Field<Type>& source,
153  const bool couples
154 ) const
155 {
156  forAll(psi_.boundaryField(), patchI)
157  {
158  const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
159  const Field<Type>& pbc = boundaryCoeffs_[patchI];
160 
161  if (!ptf.coupled())
162  {
163  addToInternalField(lduAddr().patchAddr(patchI), pbc, source);
164  }
165  else if (couples)
166  {
167  tmp<Field<Type>> tpnf = ptf.patchNeighbourField();
168  const Field<Type>& pnf = tpnf();
169 
170  const labelUList& addr = lduAddr().patchAddr(patchI);
171 
172  forAll(addr, facei)
173  {
174  source[addr[facei]] += cmptMultiply(pbc[facei], pnf[facei]);
175  }
176  }
177  }
178 }
179 
180 
181 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
182 
183 template<class Type>
185 (
187  const dimensionSet& ds
188 )
189 :
190  lduMatrix(psi.mesh()),
191  psi_(psi),
192  dimensions_(ds),
193  source_(psi.size(), Zero),
194  internalCoeffs_(psi.mesh().boundary().size()),
195  boundaryCoeffs_(psi.mesh().boundary().size()),
196  faceFluxCorrectionPtr_(nullptr)
197 {
199  << "constructing faMatrix<Type> for field " << psi_.name()
200  << endl;
201 
202  // Initialise coupling coefficients
203  forAll(psi.mesh().boundary(), patchI)
204  {
205  internalCoeffs_.set
206  (
207  patchI,
208  new Field<Type>
209  (
210  psi.mesh().boundary()[patchI].size(),
211  Zero
212  )
213  );
214 
215  boundaryCoeffs_.set
216  (
217  patchI,
218  new Field<Type>
219  (
220  psi.mesh().boundary()[patchI].size(),
221  Zero
222  )
223  );
224  }
225 
226  // Update the boundary coefficients of psi without changing its event No.
227  auto& psiRef =
229 
230  const label currentStatePsi = psiRef.eventNo();
231  psiRef.boundaryFieldRef().updateCoeffs();
232  psiRef.eventNo() = currentStatePsi;
233 }
234 
235 
236 template<class Type>
238 :
239  refCount(),
240  lduMatrix(fam),
241  psi_(fam.psi_),
242  dimensions_(fam.dimensions_),
243  source_(fam.source_),
244  internalCoeffs_(fam.internalCoeffs_),
245  boundaryCoeffs_(fam.boundaryCoeffs_),
246  faceFluxCorrectionPtr_(nullptr)
247 {
249  << "copying faMatrix<Type> for field " << psi_.name()
250  << endl;
251 
252  if (fam.faceFluxCorrectionPtr_)
253  {
254  faceFluxCorrectionPtr_ = new
256  (
257  *(fam.faceFluxCorrectionPtr_)
258  );
259  }
260 }
261 
262 
263 template<class Type>
265 (
267  Istream& is
268 )
269 :
270  lduMatrix(psi.mesh()),
271  psi_(psi),
272  dimensions_(is),
273  source_(is),
274  internalCoeffs_(psi.mesh().boundary().size()),
275  boundaryCoeffs_(psi.mesh().boundary().size()),
276  faceFluxCorrectionPtr_(nullptr)
277 {
279  << "constructing faMatrix<Type> for field " << psi_.name()
280  << endl;
281 
282  // Initialise coupling coefficients
283  forAll(psi.mesh().boundary(), patchI)
284  {
285  internalCoeffs_.set
286  (
287  patchI,
288  new Field<Type>
289  (
290  psi.mesh().boundary()[patchI].size(),
291  Zero
292  )
293  );
294 
295  boundaryCoeffs_.set
296  (
297  patchI,
298  new Field<Type>
299  (
300  psi.mesh().boundary()[patchI].size(),
301  Zero
302  )
303  );
304  }
305 }
306 
307 
308 template<class Type>
310 {
311  return tmp<faMatrix<Type>>
312  (
313  new faMatrix<Type>(*this)
314  );
315 }
316 
317 
318 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
319 
320 
321 template<class Type>
323 {
325  << "Destroying faMatrix<Type> for field " << psi_.name() << endl;
326 
327  deleteDemandDrivenData(faceFluxCorrectionPtr_);
328 }
329 
330 
331 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
332 
333 template<class Type>
334 template<template<class> class ListType>
336 (
337  const labelUList& faceLabels,
338  const ListType<Type>& values
339 )
340 {
341 #if 0 /* Specific to foam-extend */
342  // Record face labels of eliminated equations
343  for (const label i : faceLabels)
344  {
345  this->eliminatedEqns().insert(i);
346  }
347 #endif
348 
349  const faMesh& mesh = psi_.mesh();
350 
351  const labelListList& edges = mesh.patch().faceEdges();
352  const labelUList& own = mesh.owner();
353  const labelUList& nei = mesh.neighbour();
354 
355  scalarField& Diag = diag();
356  Field<Type>& psi =
357  const_cast
358  <
360  >(psi_).primitiveFieldRef();
361 
362  forAll(faceLabels, i)
363  {
364  const label facei = faceLabels[i];
365  const Type& value = values[i];
366 
367  psi[facei] = value;
368  source_[facei] = value*Diag[facei];
369 
370  if (symmetric() || asymmetric())
371  {
372  for (const label edgei : edges[facei])
373  {
374  if (mesh.isInternalEdge(edgei))
375  {
376  if (symmetric())
377  {
378  if (facei == own[edgei])
379  {
380  source_[nei[edgei]] -= upper()[edgei]*value;
381  }
382  else
383  {
384  source_[own[edgei]] -= upper()[edgei]*value;
385  }
386 
387  upper()[edgei] = 0.0;
388  }
389  else
390  {
391  if (facei == own[edgei])
392  {
393  source_[nei[edgei]] -= lower()[edgei]*value;
394  }
395  else
396  {
397  source_[own[edgei]] -= upper()[edgei]*value;
398  }
399 
400  upper()[edgei] = 0.0;
401  lower()[edgei] = 0.0;
402  }
403  }
404  else
405  {
406  const label patchi = mesh.boundary().whichPatch(edgei);
407 
408  if (internalCoeffs_[patchi].size())
409  {
410  const label patchEdgei =
411  mesh.boundary()[patchi].whichEdge(edgei);
412 
413  internalCoeffs_[patchi][patchEdgei] = Zero;
414  boundaryCoeffs_[patchi][patchEdgei] = Zero;
415  }
416  }
417  }
418  }
419  }
420 }
421 
422 
423 template<class Type>
425 (
426  const labelUList& faceLabels,
427  const Type& value
428 )
429 {
430  this->setValuesFromList(faceLabels, UniformList<Type>(value));
431 }
432 
433 
434 template<class Type>
436 (
437  const labelUList& faceLabels,
438  const UList<Type>& values
439 )
440 {
441  this->setValuesFromList(faceLabels, values);
442 }
443 
444 
445 template<class Type>
447 (
448  const labelUList& faceLabels,
450 )
451 {
452  this->setValuesFromList(faceLabels, values);
453 }
454 
455 
456 template<class Type>
458 (
459  const label facei,
460  const Type& value,
461  const bool forceReference
462 )
463 {
464  if ((forceReference || psi_.needReference()) && facei >= 0)
465  {
466  if (Pstream::master())
467  {
468  source()[facei] += diag()[facei]*value;
469  diag()[facei] += diag()[facei];
470  }
471  }
472 }
473 
474 
475 template<class Type>
477 (
478  const labelUList& faceLabels,
479  const Type& value,
480  const bool forceReference
481 )
482 {
483  if (forceReference || psi_.needReference())
484  {
485  forAll(faceLabels, facei)
486  {
487  const label faceId = faceLabels[facei];
488  if (faceId >= 0)
489  {
490  source()[faceId] += diag()[faceId]*value;
491  diag()[faceId] += diag()[faceId];
492  }
493  }
494  }
495 }
496 
497 
498 template<class Type>
500 (
501  const labelUList& faceLabels,
502  const UList<Type>& values,
503  const bool forceReference
504 )
505 {
506  if (forceReference || psi_.needReference())
507  {
508  forAll(faceLabels, facei)
509  {
510  const label faceId = faceLabels[facei];
511  if (faceId >= 0)
512  {
513  source()[faceId] += diag()[faceId]*values[facei];
514  diag()[faceId] += diag()[faceId];
515  }
516  }
517  }
518 }
519 
520 
521 template<class Type>
523 {
524  if (alpha <= 0)
525  {
526  return;
527  }
528 
529  Field<Type>& S = source();
530  scalarField& D = diag();
531 
532  // Store the current unrelaxed diagonal for use in updating the source
533  scalarField D0(D);
534 
535  // Calculate the sum-mag off-diagonal from the interior faces
536  scalarField sumOff(D.size(), Zero);
537  sumMagOffDiag(sumOff);
538 
539  // Handle the boundary contributions to the diagonal
540  forAll(psi_.boundaryField(), patchI)
541  {
542  const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
543 
544  if (ptf.size())
545  {
546  const labelUList& pa = lduAddr().patchAddr(patchI);
547  Field<Type>& iCoeffs = internalCoeffs_[patchI];
548 
549  if (ptf.coupled())
550  {
551  const Field<Type>& pCoeffs = boundaryCoeffs_[patchI];
552 
553  // For coupled boundaries add the diagonal and
554  // off-diagonal contributions
555  forAll(pa, face)
556  {
557  D[pa[face]] += component(iCoeffs[face], 0);
558  sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
559  }
560  }
561  else
562  {
563  // For non-coupled boundaries subtract the diagonal
564  // contribution off-diagonal sum which avoids having to remove
565  // it from the diagonal later.
566  // Also add the source contribution from the relaxation
567  forAll(pa, face)
568  {
569  Type iCoeff0 = iCoeffs[face];
570  iCoeffs[face] = cmptMag(iCoeffs[face]);
571  sumOff[pa[face]] -= cmptMin(iCoeffs[face]);
572  iCoeffs[face] /= alpha;
573  S[pa[face]] +=
574  cmptMultiply(iCoeffs[face] - iCoeff0, psi_[pa[face]]);
575  }
576  }
577  }
578  }
579 
580  // Ensure the matrix is diagonally dominant...
581  max(D, D, sumOff);
582 
583  // ... then relax
584  D /= alpha;
585 
586  // Now remove the diagonal contribution from coupled boundaries
587  forAll(psi_.boundaryField(), patchI)
588  {
589  const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
590 
591  if (ptf.size())
592  {
593  const labelUList& pa = lduAddr().patchAddr(patchI);
594  Field<Type>& iCoeffs = internalCoeffs_[patchI];
595 
596  if (ptf.coupled())
597  {
598  forAll(pa, face)
599  {
600  D[pa[face]] -= component(iCoeffs[face], 0);
601  }
602  }
603  }
604  }
605 
606  // Finally add the relaxation contribution to the source.
607  S += (D - D0)*psi_.internalField();
608 }
609 
610 
611 template<class Type>
613 {
614  if (psi_.mesh().relaxEquation(psi_.name()))
615  {
616  relax(psi_.mesh().equationRelaxationFactor(psi_.name()));
617  }
618  else
619  {
621  << "Relaxation factor for field " << psi_.name()
622  << " not found. Relaxation will not be used." << endl;
623  }
624 }
625 
626 
627 template<class Type>
629 {
630  tmp<scalarField> tdiag(new scalarField(diag()));
631  addCmptAvBoundaryDiag(tdiag.ref());
632  return tdiag;
633 }
634 
635 
636 template<class Type>
638 {
640  (
641  new areaScalarField
642  (
643  IOobject
644  (
645  "A("+psi_.name()+')',
646  psi_.instance(),
647  psi_.db()
648  ),
649  psi_.mesh(),
650  dimensions_/psi_.dimensions()/dimArea,
651  zeroGradientFaPatchScalarField::typeName
652  )
653  );
654 
655  tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().S();
656  tAphi.ref().correctBoundaryConditions();
657 
658  return tAphi;
659 }
660 
661 
662 template<class Type>
665 {
667  (
669  (
670  IOobject
671  (
672  "H("+psi_.name()+')',
673  psi_.instance(),
674  psi_.db()
675  ),
676  psi_.mesh(),
677  dimensions_/dimArea,
678  zeroGradientFaPatchScalarField::typeName
679  )
680  );
682 
683  // Loop over field components
684  for (direction cmpt=0; cmpt<Type::nComponents; ++cmpt)
685  {
686  scalarField psiCmpt(psi_.primitiveField().component(cmpt));
687 
688  scalarField boundaryDiagCmpt(psi_.size(), Zero);
689  addBoundaryDiag(boundaryDiagCmpt, cmpt);
690  boundaryDiagCmpt.negate();
691  addCmptAvBoundaryDiag(boundaryDiagCmpt);
692 
693  Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
694  }
695 
696  Hphi.primitiveFieldRef() += lduMatrix::H(psi_.internalField()) + source_;
697  addBoundarySource(Hphi.primitiveFieldRef());
698 
699  Hphi.primitiveFieldRef() /= psi_.mesh().S();
701 
702  return tHphi;
703 }
704 
705 
706 template<class Type>
709 {
710  if (!psi_.mesh().fluxRequired(psi_.name()))
711  {
713  << "flux requested but " << psi_.name()
714  << " not specified in the fluxRequired sub-dictionary of faSchemes"
715  << abort(FatalError);
716  }
717 
718  // construct GeometricField<Type, faePatchField, edgeMesh>
720  (
722  (
723  IOobject
724  (
725  "flux("+psi_.name()+')',
726  psi_.instance(),
727  psi_.db()
728  ),
729  psi_.mesh(),
730  dimensions()
731  )
732  );
734  tfieldFlux.ref();
735 
736  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
737  {
738  fieldFlux.primitiveFieldRef().replace
739  (
740  cmpt,
741  lduMatrix::faceH(psi_.primitiveField().component(cmpt))
742  );
743  }
744 
745  FieldField<Field, Type> InternalContrib = internalCoeffs_;
746 
747  forAll(InternalContrib, patchI)
748  {
749  InternalContrib[patchI] =
751  (
752  InternalContrib[patchI],
753  psi_.boundaryField()[patchI].patchInternalField()
754  );
755  }
756 
757  FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
758 
759  forAll(NeighbourContrib, patchI)
760  {
761  if (psi_.boundaryField()[patchI].coupled())
762  {
763  NeighbourContrib[patchI] =
765  (
766  NeighbourContrib[patchI],
767  psi_.boundaryField()[patchI].patchNeighbourField()
768  );
769  }
770  }
771 
772  forAll(fieldFlux.boundaryField(), patchI)
773  {
774  fieldFlux.boundaryFieldRef()[patchI] =
775  InternalContrib[patchI] - NeighbourContrib[patchI];
776  }
777 
778  if (faceFluxCorrectionPtr_)
779  {
780  fieldFlux += *faceFluxCorrectionPtr_;
781  }
782 
783  return tfieldFlux;
784 }
785 
786 
787 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
788 
789 template<class Type>
791 {
792  if (this == &famv)
793  {
794  return; // Self-assignment is a no-op
795  }
796 
797  if (&psi_ != &(famv.psi_))
798  {
800  << "different fields"
801  << abort(FatalError);
802  }
803 
804  lduMatrix::operator=(famv);
805  source_ = famv.source_;
806  internalCoeffs_ = famv.internalCoeffs_;
807  boundaryCoeffs_ = famv.boundaryCoeffs_;
808 
809  if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
810  {
811  *faceFluxCorrectionPtr_ = *famv.faceFluxCorrectionPtr_;
812  }
813  else if (famv.faceFluxCorrectionPtr_)
814  {
815  faceFluxCorrectionPtr_ =
817  (*famv.faceFluxCorrectionPtr_);
818  }
819 }
820 
821 
822 template<class Type>
824 {
825  operator=(tfamv());
826  tfamv.clear();
827 }
828 
829 
830 template<class Type>
832 {
834  source_.negate();
835  internalCoeffs_.negate();
836  boundaryCoeffs_.negate();
837 
838  if (faceFluxCorrectionPtr_)
839  {
840  faceFluxCorrectionPtr_->negate();
841  }
842 }
843 
844 
845 template<class Type>
847 {
848  checkMethod(*this, famv, "+=");
849 
850  dimensions_ += famv.dimensions_;
851  lduMatrix::operator+=(famv);
852  source_ += famv.source_;
853  internalCoeffs_ += famv.internalCoeffs_;
854  boundaryCoeffs_ += famv.boundaryCoeffs_;
855 
856  if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
857  {
858  *faceFluxCorrectionPtr_ += *famv.faceFluxCorrectionPtr_;
859  }
860  else if (famv.faceFluxCorrectionPtr_)
861  {
862  faceFluxCorrectionPtr_ = new
864  (
865  *famv.faceFluxCorrectionPtr_
866  );
867  }
868 }
869 
870 
871 template<class Type>
873 {
874  operator+=(tfamv());
875  tfamv.clear();
876 }
877 
878 
879 template<class Type>
881 {
882  checkMethod(*this, famv, "+=");
883 
884  dimensions_ -= famv.dimensions_;
885  lduMatrix::operator-=(famv);
886  source_ -= famv.source_;
887  internalCoeffs_ -= famv.internalCoeffs_;
888  boundaryCoeffs_ -= famv.boundaryCoeffs_;
889 
890  if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
891  {
892  *faceFluxCorrectionPtr_ -= *famv.faceFluxCorrectionPtr_;
893  }
894  else if (famv.faceFluxCorrectionPtr_)
895  {
896  faceFluxCorrectionPtr_ =
898  (-*famv.faceFluxCorrectionPtr_);
899  }
900 }
901 
902 
903 template<class Type>
905 {
906  operator-=(tfamv());
907  tfamv.clear();
908 }
909 
910 
911 template<class Type>
913 (
915 )
916 {
917  checkMethod(*this, su, "+=");
918  source() -= su.mesh().S()*su.internalField();
919 }
920 
921 
922 template<class Type>
924 (
926 )
927 {
928  operator+=(tsu());
929  tsu.clear();
930 }
931 
932 
933 template<class Type>
935 (
937 )
938 {
939  checkMethod(*this, su, "-=");
940  source() += su.mesh().S()*su.internalField();
941 }
942 
943 
944 template<class Type>
946 (
948 )
949 {
950  operator-=(tsu());
951  tsu.clear();
952 }
953 
954 
955 template<class Type>
957 (
958  const dimensioned<Type>& su
959 )
960 {
961  source() -= su.mesh().S()*su;
962 }
963 
964 
965 template<class Type>
967 (
968  const dimensioned<Type>& su
969 )
970 {
971  source() += su.mesh().S()*su;
972 }
973 
974 
975 template<class Type>
977 (
978  const areaScalarField& vsf
979 )
980 {
981  dimensions_ *= vsf.dimensions();
982  lduMatrix::operator*=(vsf.internalField());
983  source_ *= vsf.internalField();
984 
985  forAll(boundaryCoeffs_, patchI)
986  {
987  const scalarField psf(vsf.boundaryField()[patchI].patchInternalField());
988  internalCoeffs_[patchI] *= psf;
989  boundaryCoeffs_[patchI] *= psf;
990  }
991 
992  if (faceFluxCorrectionPtr_)
993  {
995  << "cannot scale a matrix containing a faceFluxCorrection"
996  << abort(FatalError);
997  }
998 }
999 
1000 
1001 template<class Type>
1004  const tmp<areaScalarField>& tvsf
1005 )
1006 {
1007  operator*=(tvsf());
1008  tvsf.clear();
1009 }
1010 
1011 
1012 template<class Type>
1015  const dimensioned<scalar>& ds
1016 )
1017 {
1018  dimensions_ *= ds.dimensions();
1019  lduMatrix::operator*=(ds.value());
1020  source_ *= ds.value();
1021  internalCoeffs_ *= ds.value();
1022  boundaryCoeffs_ *= ds.value();
1023 
1024  if (faceFluxCorrectionPtr_)
1025  {
1026  *faceFluxCorrectionPtr_ *= ds.value();
1027  }
1028 }
1029 
1030 
1031 // * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
1032 
1033 template<class Type>
1034 void Foam::checkMethod
1036  const faMatrix<Type>& fam1,
1037  const faMatrix<Type>& fam2,
1038  const char* op
1039 )
1040 {
1041  if (&fam1.psi() != &fam2.psi())
1042  {
1044  << "incompatible fields for operation "
1045  << endl << " "
1046  << "[" << fam1.psi().name() << "] "
1047  << op
1048  << " [" << fam2.psi().name() << "]"
1049  << abort(FatalError);
1050  }
1051 
1052  if (dimensionSet::debug && fam1.dimensions() != fam2.dimensions())
1053  {
1055  << "incompatible dimensions for operation "
1056  << endl << " "
1057  << "[" << fam1.psi().name() << fam1.dimensions()/dimArea << " ] "
1058  << op
1059  << " [" << fam2.psi().name() << fam2.dimensions()/dimArea << " ]"
1060  << abort(FatalError);
1061  }
1062 }
1063 
1064 
1065 template<class Type>
1066 void Foam::checkMethod
1068  const faMatrix<Type>& fam,
1070  const char* op
1071 )
1072 {
1073  if (dimensionSet::debug && fam.dimensions()/dimArea != vf.dimensions())
1074  {
1076  << "incompatible dimensions for operation "
1077  << endl << " "
1078  << "[" << fam.psi().name() << fam.dimensions()/dimArea << " ] "
1079  << op
1080  << " [" << vf.name() << vf.dimensions() << " ]"
1081  << abort(FatalError);
1082  }
1083 }
1084 
1085 
1086 template<class Type>
1087 void Foam::checkMethod
1089  const faMatrix<Type>& fam,
1090  const dimensioned<Type>& dt,
1091  const char* op
1092 )
1093 {
1094  if (dimensionSet::debug && fam.dimensions()/dimArea != dt.dimensions())
1095  {
1097  << "incompatible dimensions for operation "
1098  << endl << " "
1099  << "[" << fam.psi().name() << fam.dimensions()/dimArea << " ] "
1100  << op
1101  << " [" << dt.name() << dt.dimensions() << " ]"
1102  << abort(FatalError);
1103  }
1104 }
1105 
1106 
1107 template<class Type>
1109 (
1110  faMatrix<Type>& fam,
1111  Istream& solverControls
1112 )
1113 {
1114  return fam.solve(solverControls);
1115 }
1116 
1117 
1118 template<class Type>
1120 (
1121  const tmp<faMatrix<Type>>& tfam,
1122  Istream& solverControls
1123 )
1124 {
1125  SolverPerformance<Type> solverPerf =
1126  const_cast<faMatrix<Type>&>(tfam()).solve(solverControls);
1127 
1128  tfam.clear();
1129  return solverPerf;
1130 }
1131 
1132 
1133 template<class Type>
1135 {
1136  return fam.solve();
1137 }
1138 
1139 
1140 template<class Type>
1141 Foam::SolverPerformance<Type> Foam::solve(const tmp<faMatrix<Type>>& tfam)
1142 {
1143  SolverPerformance<Type> solverPerf =
1144  const_cast<faMatrix<Type>&>(tfam()).solve();
1145 
1146  tfam.clear();
1147  return solverPerf;
1148 }
1149 
1150 
1151 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
1152 
1153 template<class Type>
1154 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1155 (
1156  const faMatrix<Type>& A,
1157  const faMatrix<Type>& B
1158 )
1159 {
1160  checkMethod(A, B, "+");
1161  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1162  tC.ref() += B;
1163  return tC;
1164 }
1165 
1166 
1167 template<class Type>
1168 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1169 (
1170  const tmp<faMatrix<Type>>& tA,
1171  const faMatrix<Type>& B
1172 )
1173 {
1174  checkMethod(tA(), B, "+");
1175  tmp<faMatrix<Type>> tC(tA.ptr());
1176  tC.ref() += B;
1177  return tC;
1178 }
1179 
1180 
1181 template<class Type>
1182 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1183 (
1184  const faMatrix<Type>& A,
1185  const tmp<faMatrix<Type>>& tB
1186 )
1187 {
1188  checkMethod(A, tB(), "+");
1189  tmp<faMatrix<Type>> tC(tB.ptr());
1190  tC.ref() += A;
1191  return tC;
1192 }
1193 
1194 
1195 template<class Type>
1196 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1197 (
1198  const tmp<faMatrix<Type>>& tA,
1199  const tmp<faMatrix<Type>>& tB
1200 )
1201 {
1202  checkMethod(tA(), tB(), "+");
1203  tmp<faMatrix<Type>> tC(tA.ptr());
1204  tC.ref() += tB();
1205  tB.clear();
1206  return tC;
1207 }
1208 
1209 
1210 template<class Type>
1211 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1212 (
1213  const faMatrix<Type>& A
1214 )
1215 {
1216  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1217  tC.ref().negate();
1218  return tC;
1219 }
1220 
1221 
1222 template<class Type>
1223 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1224 (
1225  const tmp<faMatrix<Type>>& tA
1226 )
1227 {
1228  tmp<faMatrix<Type>> tC(tA.ptr());
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 faMatrix<Type>& A,
1238  const faMatrix<Type>& B
1239 )
1240 {
1241  checkMethod(A, B, "-");
1242  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1243  tC.ref() -= B;
1244  return tC;
1245 }
1246 
1247 
1248 template<class Type>
1249 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1250 (
1251  const tmp<faMatrix<Type>>& tA,
1252  const faMatrix<Type>& B
1253 )
1254 {
1255  checkMethod(tA(), B, "-");
1256  tmp<faMatrix<Type>> tC(tA.ptr());
1257  tC.ref() -= B;
1258  return tC;
1259 }
1260 
1261 
1262 template<class Type>
1263 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1264 (
1265  const faMatrix<Type>& A,
1266  const tmp<faMatrix<Type>>& tB
1267 )
1268 {
1269  checkMethod(A, tB(), "-");
1270  tmp<faMatrix<Type>> tC(tB.ptr());
1271  tC.ref() -= A;
1272  tC.ref().negate();
1273  return tC;
1274 }
1275 
1276 
1277 template<class Type>
1278 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1279 (
1280  const tmp<faMatrix<Type>>& tA,
1281  const tmp<faMatrix<Type>>& tB
1282 )
1283 {
1284  checkMethod(tA(), tB(), "-");
1285  tmp<faMatrix<Type>> tC(tA.ptr());
1286  tC.ref() -= tB();
1287  tB.clear();
1288  return tC;
1289 }
1290 
1291 
1292 template<class Type>
1293 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1294 (
1295  const faMatrix<Type>& A,
1296  const faMatrix<Type>& B
1297 )
1298 {
1299  checkMethod(A, B, "==");
1300  return (A - B);
1301 }
1302 
1303 
1304 template<class Type>
1305 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1306 (
1307  const tmp<faMatrix<Type>>& tA,
1308  const faMatrix<Type>& B
1309 )
1310 {
1311  checkMethod(tA(), B, "==");
1312  return (tA - B);
1313 }
1314 
1315 
1316 template<class Type>
1317 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1318 (
1319  const faMatrix<Type>& A,
1320  const tmp<faMatrix<Type>>& tB
1321 )
1322 {
1323  checkMethod(A, tB(), "==");
1324  return (A - tB);
1325 }
1326 
1327 
1328 template<class Type>
1329 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1330 (
1331  const tmp<faMatrix<Type>>& tA,
1332  const tmp<faMatrix<Type>>& tB
1333 )
1334 {
1335  checkMethod(tA(), tB(), "==");
1336  return (tA - tB);
1337 }
1338 
1339 
1340 template<class Type>
1341 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1342 (
1343  const faMatrix<Type>& A,
1344  const GeometricField<Type, faPatchField, areaMesh>& su
1345 )
1346 {
1347  checkMethod(A, su, "+");
1348  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1349  tC.ref().source() -= su.mesh().S()*su.internalField();
1350  return tC;
1351 }
1352 
1353 
1354 template<class Type>
1355 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1356 (
1357  const tmp<faMatrix<Type>>& tA,
1358  const GeometricField<Type, faPatchField, areaMesh>& su
1359 )
1360 {
1361  checkMethod(tA(), su, "+");
1362  tmp<faMatrix<Type>> tC(tA.ptr());
1363  tC.ref().source() -= su.mesh().S()*su.internalField();
1364  return tC;
1365 }
1366 
1367 
1368 template<class Type>
1369 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1370 (
1371  const faMatrix<Type>& A,
1372  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1373 )
1374 {
1375  checkMethod(A, tsu(), "+");
1376  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1377  tC.ref().source() -= tsu().mesh().S()*tsu().internalField();
1378  tsu.clear();
1379  return tC;
1380 }
1381 
1382 
1383 template<class Type>
1384 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1385 (
1386  const tmp<faMatrix<Type>>& tA,
1387  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1388 )
1389 {
1390  checkMethod(tA(), tsu(), "+");
1391  tmp<faMatrix<Type>> tC(tA.ptr());
1392  tC.ref().source() -= tsu().mesh().S()*tsu().internalField();
1393  tsu.clear();
1394  return tC;
1395 }
1396 
1397 
1398 template<class Type>
1399 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1400 (
1401  const GeometricField<Type, faPatchField, areaMesh>& su,
1402  const faMatrix<Type>& A
1403 )
1404 {
1405  checkMethod(A, su, "+");
1406  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1407  tC.ref().source() -= su.mesh().S()*su.internalField();
1408  return tC;
1409 }
1410 
1411 
1412 template<class Type>
1413 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1414 (
1415  const GeometricField<Type, faPatchField, areaMesh>& su,
1416  const tmp<faMatrix<Type>>& tA
1417 )
1418 {
1419  checkMethod(tA(), su, "+");
1420  tmp<faMatrix<Type>> tC(tA.ptr());
1421  tC.ref().source() -= su.mesh().S()*su.internalField();
1422  return tC;
1423 }
1424 
1425 
1426 template<class Type>
1427 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1428 (
1429  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1430  const faMatrix<Type>& A
1431 )
1432 {
1433  checkMethod(A, tsu(), "+");
1434  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1435  tC.ref().source() -= tsu().mesh().S()*tsu().internalField();
1436  tsu.clear();
1437  return tC;
1438 }
1439 
1440 
1441 template<class Type>
1442 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1443 (
1444  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1445  const tmp<faMatrix<Type>>& tA
1446 )
1447 {
1448  checkMethod(tA(), tsu(), "+");
1449  tmp<faMatrix<Type>> tC(tA.ptr());
1450  tC.ref().source() -= tsu().mesh().S()*tsu().internalField();
1451  tsu.clear();
1452  return tC;
1453 }
1454 
1455 
1456 template<class Type>
1457 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1458 (
1459  const faMatrix<Type>& A,
1460  const GeometricField<Type, faPatchField, areaMesh>& su
1461 )
1462 {
1463  checkMethod(A, su, "-");
1464  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1465  tC.ref().source() += su.mesh().S()*su.internalField();
1466  return tC;
1467 }
1468 
1469 
1470 template<class Type>
1471 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1472 (
1473  const tmp<faMatrix<Type>>& tA,
1474  const GeometricField<Type, faPatchField, areaMesh>& su
1475 )
1476 {
1477  checkMethod(tA(), su, "-");
1478  tmp<faMatrix<Type>> tC(tA.ptr());
1479  tC.ref().source() += su.mesh().S()*su.internalField();
1480  return tC;
1481 }
1482 
1483 
1484 template<class Type>
1485 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1486 (
1487  const faMatrix<Type>& A,
1488  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1489 )
1490 {
1491  checkMethod(A, tsu(), "-");
1492  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1493  tC.ref().source() += tsu().mesh().S()*tsu().internalField();
1494  tsu.clear();
1495  return tC;
1496 }
1497 
1498 
1499 template<class Type>
1500 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1501 (
1502  const tmp<faMatrix<Type>>& tA,
1503  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1504 )
1505 {
1506  checkMethod(tA(), tsu(), "-");
1507  tmp<faMatrix<Type>> tC(tA.ptr());
1508  tC.ref().source() += tsu().mesh().S()*tsu().internalField();
1509  tsu.clear();
1510  return tC;
1511 }
1512 
1513 
1514 template<class Type>
1515 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1516 (
1517  const GeometricField<Type, faPatchField, areaMesh>& su,
1518  const faMatrix<Type>& A
1519 )
1520 {
1521  checkMethod(A, su, "-");
1522  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1523  tC.ref().negate();
1524  tC.ref().source() -= su.mesh().S()*su.internalField();
1525  return tC;
1526 }
1527 
1528 
1529 template<class Type>
1530 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1531 (
1532  const GeometricField<Type, faPatchField, areaMesh>& su,
1533  const tmp<faMatrix<Type>>& tA
1534 )
1535 {
1536  checkMethod(tA(), su, "-");
1537  tmp<faMatrix<Type>> tC(tA.ptr());
1538  tC.ref().negate();
1539  tC.ref().source() -= su.mesh().S()*su.internalField();
1540  return tC;
1541 }
1542 
1543 
1544 template<class Type>
1545 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1546 (
1547  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1548  const faMatrix<Type>& A
1549 )
1550 {
1551  checkMethod(A, tsu(), "-");
1552  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1553  tC.ref().negate();
1554  tC.ref().source() -= tsu().mesh().S()*tsu().internalField();
1555  tsu.clear();
1556  return tC;
1557 }
1558 
1559 
1560 template<class Type>
1561 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1562 (
1563  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1564  const tmp<faMatrix<Type>>& tA
1565 )
1566 {
1567  checkMethod(tA(), tsu(), "-");
1568  tmp<faMatrix<Type>> tC(tA.ptr());
1569  tC.ref().negate();
1570  tC.ref().source() -= tsu().mesh().S()*tsu().internalField();
1571  tsu.clear();
1572  return tC;
1573 }
1574 
1575 
1576 template<class Type>
1577 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1578 (
1579  const faMatrix<Type>& A,
1580  const dimensioned<Type>& su
1581 )
1582 {
1583  checkMethod(A, su, "+");
1584  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1585  tC.ref().source() -= su.value()*A.psi().mesh().S();
1586  return tC;
1587 }
1588 
1589 
1590 template<class Type>
1591 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1592 (
1593  const tmp<faMatrix<Type>>& tA,
1594  const dimensioned<Type>& su
1595 )
1596 {
1597  checkMethod(tA(), su, "+");
1598  tmp<faMatrix<Type>> tC(tA.ptr());
1599  tC.ref().source() -= su.value()*tC().psi().mesh().S();
1600  return tC;
1601 }
1602 
1603 
1604 template<class Type>
1605 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1606 (
1607  const dimensioned<Type>& su,
1608  const faMatrix<Type>& A
1609 )
1610 {
1611  checkMethod(A, su, "+");
1612  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1613  tC.ref().source() -= su.value()*A.psi().mesh().S();
1614  return tC;
1615 }
1616 
1617 
1618 template<class Type>
1619 Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1620 (
1621  const dimensioned<Type>& su,
1622  const tmp<faMatrix<Type>>& tA
1623 )
1624 {
1625  checkMethod(tA(), su, "+");
1626  tmp<faMatrix<Type>> tC(tA.ptr());
1627  tC.ref().source() -= su.value()*tC().psi().mesh().S();
1628  return tC;
1629 }
1630 
1631 
1632 template<class Type>
1633 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1634 (
1635  const faMatrix<Type>& A,
1636  const dimensioned<Type>& su
1637 )
1638 {
1639  checkMethod(A, su, "-");
1640  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1641  tC.ref().source() += su.value()*tC().psi().mesh().S();
1642  return tC;
1643 }
1644 
1645 
1646 template<class Type>
1647 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1648 (
1649  const tmp<faMatrix<Type>>& tA,
1650  const dimensioned<Type>& su
1651 )
1652 {
1653  checkMethod(tA(), su, "-");
1654  tmp<faMatrix<Type>> tC(tA.ptr());
1655  tC.ref().source() += su.value()*tC().psi().mesh().S();
1656  return tC;
1657 }
1658 
1659 
1660 template<class Type>
1661 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1662 (
1663  const dimensioned<Type>& su,
1664  const faMatrix<Type>& A
1665 )
1666 {
1667  checkMethod(A, su, "-");
1668  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1669  tC.ref().negate();
1670  tC.ref().source() -= su.value()*A.psi().mesh().S();
1671  return tC;
1672 }
1673 
1674 
1675 template<class Type>
1676 Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1677 (
1678  const dimensioned<Type>& su,
1679  const tmp<faMatrix<Type>>& tA
1680 )
1681 {
1682  checkMethod(tA(), su, "-");
1683  tmp<faMatrix<Type>> tC(tA.ptr());
1684  tC.ref().negate();
1685  tC.ref().source() -= su.value()*tC().psi().mesh().S();
1686  return tC;
1687 }
1688 
1689 
1690 template<class Type>
1691 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1692 (
1693  const faMatrix<Type>& A,
1694  const GeometricField<Type, faPatchField, areaMesh>& su
1695 )
1696 {
1697  checkMethod(A, su, "==");
1698  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1699  tC.ref().source() += su.mesh().S()*su.internalField();
1700  return tC;
1701 }
1702 
1703 
1704 template<class Type>
1705 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1706 (
1707  const tmp<faMatrix<Type>>& tA,
1708  const GeometricField<Type, faPatchField, areaMesh>& su
1709 )
1710 {
1711  checkMethod(tA(), su, "==");
1712  tmp<faMatrix<Type>> tC(tA.ptr());
1713  tC.ref().source() += su.mesh().S()*su.internalField();
1714  return tC;
1715 }
1716 
1717 
1718 template<class Type>
1719 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1720 (
1721  const faMatrix<Type>& A,
1722  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1723 )
1724 {
1725  checkMethod(A, tsu(), "==");
1726  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1727  tC.ref().source() += tsu().mesh().S()*tsu().internalField();
1728  tsu.clear();
1729  return tC;
1730 }
1731 
1732 
1733 template<class Type>
1734 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1735 (
1736  const tmp<faMatrix<Type>>& tA,
1737  const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1738 )
1739 {
1740  checkMethod(tA(), tsu(), "==");
1741  tmp<faMatrix<Type>> tC(tA.ptr());
1742  tC.ref().source() += tsu().mesh().S()*tsu().internalField();
1743  tsu.clear();
1744  return tC;
1745 }
1746 
1747 
1748 template<class Type>
1749 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1750 (
1751  const faMatrix<Type>& A,
1752  const dimensioned<Type>& su
1753 )
1754 {
1755  checkMethod(A, su, "==");
1756  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1757  tC.ref().source() += A.psi().mesh().S()*su.value();
1758  return tC;
1759 }
1760 
1761 
1762 template<class Type>
1763 Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1764 (
1765  const tmp<faMatrix<Type>>& tA,
1766  const dimensioned<Type>& su
1767 )
1768 {
1769  checkMethod(tA(), su, "==");
1770  tmp<faMatrix<Type>> tC(tA.ptr());
1771  tC.ref().source() += tC().psi().mesh().S()*su.value();
1772  return tC;
1773 }
1774 
1775 
1776 template<class Type>
1777 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1778 (
1779  const areaScalarField& vsf,
1780  const faMatrix<Type>& A
1781 )
1782 {
1783  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1784  tC.ref() *= vsf;
1785  return tC;
1786 }
1787 
1788 
1789 template<class Type>
1790 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1791 (
1792  const tmp<areaScalarField>& tvsf,
1793  const faMatrix<Type>& A
1794 )
1795 {
1796  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1797  tC.ref() *= tvsf;
1798  return tC;
1799 }
1800 
1801 
1802 template<class Type>
1803 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1804 (
1805  const areaScalarField& vsf,
1806  const tmp<faMatrix<Type>>& tA
1807 )
1808 {
1809  tmp<faMatrix<Type>> tC(tA.ptr());
1810  tC.ref() *= vsf;
1811  return tC;
1812 }
1813 
1814 
1815 template<class Type>
1816 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1817 (
1818  const tmp<areaScalarField>& tvsf,
1819  const tmp<faMatrix<Type>>& tA
1820 )
1821 {
1822  tmp<faMatrix<Type>> tC(tA.ptr());
1823  tC.ref() *= tvsf;
1824  return tC;
1825 }
1826 
1827 
1828 template<class Type>
1829 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1830 (
1831  const dimensioned<scalar>& ds,
1832  const faMatrix<Type>& A
1833 )
1834 {
1835  tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1836  tC.ref() *= ds;
1837  return tC;
1838 }
1839 
1840 
1841 template<class Type>
1842 Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1843 (
1844  const dimensioned<scalar>& ds,
1845  const tmp<faMatrix<Type>>& tA
1846 )
1847 {
1848  tmp<faMatrix<Type>> tC(tA.ptr());
1849  tC.ref() *= ds;
1850  return tC;
1851 }
1852 
1853 
1854 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1855 
1856 template<class Type>
1857 Foam::Ostream& Foam::operator<<(Ostream& os, const faMatrix<Type>& fam)
1858 {
1859  os << static_cast<const lduMatrix&>(fam) << nl
1860  << fam.dimensions_ << nl
1861  << fam.source_ << nl
1862  << fam.internalCoeffs_ << nl
1863  << fam.boundaryCoeffs_ << endl;
1864 
1865  os.check(FUNCTION_NAME);
1866 
1867  return os;
1868 }
1869 
1870 
1871 // * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
1872 
1873 #include "faMatrixSolve.C"
1874 
1875 // ************************************************************************* //
Foam::checkMethod
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1035
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::faMatrix::operator+=
void operator+=(const faMatrix< Type > &)
Definition: faMatrix.C:846
Foam::lduMatrix::operator-=
void operator-=(const lduMatrix &)
Definition: lduMatrixOperations.C:223
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::faMatrix::addToInternalField
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: faMatrix.C:42
Foam::faMatrix::dimensions
const dimensionSet & dimensions() const
Definition: faMatrix.H:244
Foam::faPatchField
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: areaFieldsFwd.H:50
Foam::component
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Definition: FieldFieldFunctions.C:44
Foam::faMatrix
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:59
UIndirectList.H
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:291
Foam::faMatrix::negate
void negate()
Definition: faMatrix.C:831
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::faMatrix::setValues
void setValues(const labelUList &faceLabels, const Type &value)
Definition: faMatrix.C:425
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::lduMatrix::operator+=
void operator+=(const lduMatrix &)
Definition: lduMatrixOperations.C:144
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::faMatrix::operator-=
void operator-=(const faMatrix< Type > &)
Definition: faMatrix.C:880
Foam::refCount
Reference counter for various OpenFOAM components.
Definition: refCount.H:50
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::lduMatrix
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:83
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
Foam::faMatrix::flux
tmp< GeometricField< Type, faePatchField, edgeMesh > > flux() const
Return the face-flux field from the matrix.
Definition: faMatrix.C:708
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::faMatrix::faMatrix
faMatrix(const GeometricField< Type, faPatchField, areaMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition: faMatrix.C:185
Foam::dimensionSet
Dimension set for the base types.
Definition: dimensionSet.H:65
primitiveFieldRef
conserve primitiveFieldRef()+
Foam::faMatrix::subtractFromInternalField
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: faMatrix.C:80
Foam::GeometricField::Boundary::updateCoeffs
void updateCoeffs()
Update the boundary condition coefficients.
Definition: GeometricBoundaryField.C:395
faMatrixSolve.C
Finite-Area matrix basic solvers.
Foam::faMatrix::setReferences
void setReferences(const labelUList &faceLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: faMatrix.C:477
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::UniformList< Type >
Foam::cmptMin
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:302
Foam::cmptMag
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:400
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
Foam::operator<<
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Foam::faMatrix::addBoundaryDiag
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: faMatrix.C:117
Foam::faMatrix::operator=
void operator=(const faMatrix< Type > &)
Definition: faMatrix.C:790
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::faPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: faPatchField.H:320
Foam::stringOps::lower
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1186
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:365
edgeFields.H
calculatedFaPatchFields.H
faceId
label faceId(-1)
Foam::primitiveMesh::faceEdges
const labelListList & faceEdges() const
Definition: primitiveMeshEdges.C:528
Foam::cmptAv
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh >> cmptAv(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:246
Foam::faMatrix::setValuesFromList
void setValuesFromList(const labelUList &faceLabels, const ListType< Type > &values)
Set solution in given faces to the specified values.
Definition: faMatrix.C:336
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::faMatrix::addBoundarySource
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: faMatrix.C:151
Foam::faMatrix::relax
void relax()
Relax matrix (for steady-state solution).
Definition: faMatrix.C:612
UniformList.H
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
areaFields.H
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:51
Foam::FatalError
error FatalError
Foam::fvMesh::neighbour
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:411
Foam::faMatrix::H
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:664
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::lduMatrix::negate
void negate()
Definition: lduMatrixOperations.C:125
Foam::tmp::ptr
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:259
Foam::areaScalarField
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:53
Foam::lduMatrix::H
tmp< Field< Type > > H(const Field< Type > &) const
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::fvMesh::owner
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:405
Foam::faMatrix::setReference
void setReference(const label facei, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: faMatrix.C:458
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
zeroGradientFaPatchFields.H
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::faPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: faPatchField.H:364
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:679
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::lduMatrix::operator=
void operator=(const lduMatrix &)
Definition: lduMatrixOperations.C:91
Foam::faMatrix::psi
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:239
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::List< labelList >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::lduMatrix::operator*=
void operator*=(const scalarField &)
Definition: lduMatrixOperations.C:302
Foam::faMatrix::addCmptAvBoundaryDiag
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: faMatrix.C:135
Foam::lduMatrix::faceH
tmp< Field< Type > > faceH(const Field< Type > &) const
Foam::UList< label >
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::faMatrix::~faMatrix
virtual ~faMatrix()
Destructor.
Definition: faMatrix.C:322
FUNCTION_NAME
#define FUNCTION_NAME
Definition: messageStream.H:270
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::faMesh
Finite area mesh. Used for 2-D non-Euclidian finite area method.
Definition: faMesh.H:77
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::faMatrix::D
tmp< scalarField > D() const
Return the matrix diagonal.
Definition: faMatrix.C:628
Foam::faMatrix::A
tmp< areaScalarField > A() const
Return the central coefficient.
Definition: faMatrix.C:637
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
fam
Calculate the matrix for the second temporal derivative.
Foam::stringOps::upper
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1202
Foam::faMatrix::clone
tmp< faMatrix< Type > > clone() const
Clone.
Definition: faMatrix.C:309
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:420
Foam::GeometricField
Generic GeometricField class.
Definition: areaFieldsFwd.H:53
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::SolverPerformance
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
Definition: SolverPerformance.H:51
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
relax
UEqn relax()