faMatrix.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2017 Wikki Ltd
9 Copyright (C) 2019-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "IndirectList.H"
34#include "UniformList.H"
35#include "demandDrivenData.H"
36
37// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38
39template<class Type>
40template<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
63template<class Type>
64template<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
77template<class Type>
78template<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
101template<class Type>
102template<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
115template<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
134template<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
149template<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 {
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
183template<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>(psi.mesh().boundary()[patchi].size(), Zero)
209 );
210
211 boundaryCoeffs_.set
212 (
213 patchi,
214 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
215 );
216 }
217
218 // Update the boundary coefficients of psi without changing its event No.
219 auto& psiRef =
221
222 const label currentStatePsi = psiRef.eventNo();
223 psiRef.boundaryFieldRef().updateCoeffs();
224 psiRef.eventNo() = currentStatePsi;
225}
226
227
228template<class Type>
230:
231 lduMatrix(fam),
232 psi_(fam.psi_),
233 dimensions_(fam.dimensions_),
234 source_(fam.source_),
235 internalCoeffs_(fam.internalCoeffs_),
236 boundaryCoeffs_(fam.boundaryCoeffs_),
237 faceFluxCorrectionPtr_(nullptr)
238{
240 << "Copying faMatrix<Type> for field " << psi_.name() << endl;
241
242 if (fam.faceFluxCorrectionPtr_)
243 {
244 faceFluxCorrectionPtr_ =
246 (
247 *(fam.faceFluxCorrectionPtr_)
248 );
249 }
250}
251
252
253template<class Type>
255:
256 lduMatrix(tmat.constCast(), tmat.movable()),
257 psi_(tmat().psi_),
258 dimensions_(tmat().dimensions_),
259 source_(tmat.constCast().source_, tmat.movable()),
260 internalCoeffs_(tmat.constCast().internalCoeffs_, tmat.movable()),
261 boundaryCoeffs_(tmat.constCast().boundaryCoeffs_, tmat.movable()),
262 faceFluxCorrectionPtr_(nullptr)
263{
265 << "Copy/Move faMatrix<Type> for field " << psi_.name() << endl;
266
267 if (tmat().faceFluxCorrectionPtr_)
268 {
269 if (tmat.movable())
270 {
271 faceFluxCorrectionPtr_ = tmat().faceFluxCorrectionPtr_;
272 tmat().faceFluxCorrectionPtr_ = nullptr;
273 }
274 else
275 {
276 faceFluxCorrectionPtr_ =
278 (
279 *(tmat().faceFluxCorrectionPtr_)
280 );
281 }
282 }
283
284 tmat.clear();
285}
286
287
288// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
289
290template<class Type>
292{
294 << "Destroying faMatrix<Type> for field " << psi_.name() << endl;
295
296 deleteDemandDrivenData(faceFluxCorrectionPtr_);
297}
298
299
300// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
301
302template<class Type>
303template<template<class> class ListType>
305(
306 const labelUList& faceLabels,
307 const ListType<Type>& values
308)
309{
310#if 0 /* Specific to foam-extend */
311 // Record face labels of eliminated equations
312 for (const label i : faceLabels)
313 {
314 this->eliminatedEqns().insert(i);
315 }
316#endif
317
318 const faMesh& mesh = psi_.mesh();
319
320 const labelListList& edges = mesh.patch().faceEdges();
321 const labelUList& own = mesh.owner();
322 const labelUList& nei = mesh.neighbour();
323
324 scalarField& Diag = diag();
326 const_cast
327 <
329 >(psi_).primitiveFieldRef();
330
331 forAll(faceLabels, i)
332 {
333 const label facei = faceLabels[i];
334 const Type& value = values[i];
335
336 psi[facei] = value;
337 source_[facei] = value*Diag[facei];
338
339 if (symmetric() || asymmetric())
340 {
341 for (const label edgei : edges[facei])
342 {
343 if (mesh.isInternalEdge(edgei))
344 {
345 if (symmetric())
346 {
347 if (facei == own[edgei])
348 {
349 source_[nei[edgei]] -= upper()[edgei]*value;
350 }
351 else
352 {
353 source_[own[edgei]] -= upper()[edgei]*value;
354 }
355
356 upper()[edgei] = 0.0;
357 }
358 else
359 {
360 if (facei == own[edgei])
361 {
362 source_[nei[edgei]] -= lower()[edgei]*value;
363 }
364 else
365 {
366 source_[own[edgei]] -= upper()[edgei]*value;
367 }
368
369 upper()[edgei] = 0.0;
370 lower()[edgei] = 0.0;
371 }
372 }
373 else
374 {
375 const label patchi = mesh.boundary().whichPatch(edgei);
376
377 if (internalCoeffs_[patchi].size())
378 {
379 const label patchEdgei =
380 mesh.boundary()[patchi].whichEdge(edgei);
381
382 internalCoeffs_[patchi][patchEdgei] = Zero;
383 boundaryCoeffs_[patchi][patchEdgei] = Zero;
384 }
385 }
386 }
387 }
388 }
389}
390
391
392template<class Type>
394(
395 const labelUList& faceLabels,
396 const Type& value
397)
398{
399 this->setValuesFromList(faceLabels, UniformList<Type>(value));
400}
401
402
403template<class Type>
405(
406 const labelUList& faceLabels,
407 const UList<Type>& values
408)
409{
410 this->setValuesFromList(faceLabels, values);
411}
412
413
414template<class Type>
416(
417 const labelUList& faceLabels,
418 const UIndirectList<Type>& values
419)
420{
421 this->setValuesFromList(faceLabels, values);
422}
423
424
425template<class Type>
427(
428 const label facei,
429 const Type& value,
430 const bool forceReference
431)
432{
433 if ((forceReference || psi_.needReference()) && facei >= 0)
434 {
435 if (Pstream::master())
436 {
437 source()[facei] += diag()[facei]*value;
438 diag()[facei] += diag()[facei];
439 }
440 }
441}
442
443
444template<class Type>
446(
447 const labelUList& faceLabels,
448 const Type& value,
449 const bool forceReference
450)
451{
452 if (forceReference || psi_.needReference())
453 {
454 forAll(faceLabels, facei)
455 {
456 const label faceId = faceLabels[facei];
457 if (faceId >= 0)
458 {
459 source()[faceId] += diag()[faceId]*value;
460 diag()[faceId] += diag()[faceId];
461 }
462 }
463 }
464}
465
466
467template<class Type>
469(
470 const labelUList& faceLabels,
471 const UList<Type>& values,
472 const bool forceReference
473)
474{
475 if (forceReference || psi_.needReference())
476 {
477 forAll(faceLabels, facei)
478 {
479 const label faceId = faceLabels[facei];
480 if (faceId >= 0)
481 {
482 source()[faceId] += diag()[faceId]*values[facei];
483 diag()[faceId] += diag()[faceId];
484 }
485 }
486 }
487}
488
489
490template<class Type>
492{
493 if (alpha <= 0)
494 {
495 return;
496 }
497
498 Field<Type>& S = source();
499 scalarField& D = diag();
500
501 // Store the current unrelaxed diagonal for use in updating the source
502 scalarField D0(D);
503
504 // Calculate the sum-mag off-diagonal from the interior faces
505 scalarField sumOff(D.size(), Zero);
506 sumMagOffDiag(sumOff);
507
508 // Handle the boundary contributions to the diagonal
509 forAll(psi_.boundaryField(), patchI)
510 {
511 const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
512
513 if (ptf.size())
514 {
515 const labelUList& pa = lduAddr().patchAddr(patchI);
516 Field<Type>& iCoeffs = internalCoeffs_[patchI];
517
518 if (ptf.coupled())
519 {
520 const Field<Type>& pCoeffs = boundaryCoeffs_[patchI];
521
522 // For coupled boundaries add the diagonal and
523 // off-diagonal contributions
524 forAll(pa, face)
525 {
526 D[pa[face]] += component(iCoeffs[face], 0);
527 sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
528 }
529 }
530 else
531 {
532 // For non-coupled boundaries subtract the diagonal
533 // contribution off-diagonal sum which avoids having to remove
534 // it from the diagonal later.
535 // Also add the source contribution from the relaxation
536 forAll(pa, face)
537 {
538 Type iCoeff0 = iCoeffs[face];
539 iCoeffs[face] = cmptMag(iCoeffs[face]);
540 sumOff[pa[face]] -= cmptMin(iCoeffs[face]);
541 iCoeffs[face] /= alpha;
542 S[pa[face]] +=
543 cmptMultiply(iCoeffs[face] - iCoeff0, psi_[pa[face]]);
544 }
545 }
546 }
547 }
548
549 // Ensure the matrix is diagonally dominant...
550 max(D, D, sumOff);
551
552 // ... then relax
553 D /= alpha;
554
555 // Now remove the diagonal contribution from coupled boundaries
556 forAll(psi_.boundaryField(), patchI)
557 {
558 const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
559
560 if (ptf.size())
561 {
562 const labelUList& pa = lduAddr().patchAddr(patchI);
563 Field<Type>& iCoeffs = internalCoeffs_[patchI];
564
565 if (ptf.coupled())
566 {
567 forAll(pa, face)
568 {
569 D[pa[face]] -= component(iCoeffs[face], 0);
570 }
571 }
572 }
573 }
574
575 // Finally add the relaxation contribution to the source.
576 S += (D - D0)*psi_.primitiveField();
577}
578
579
580template<class Type>
582{
583 if (psi_.mesh().relaxEquation(psi_.name()))
584 {
585 relax(psi_.mesh().equationRelaxationFactor(psi_.name()));
586 }
587 else
588 {
590 << "Relaxation factor for field " << psi_.name()
591 << " not found. Relaxation will not be used." << endl;
592 }
593}
594
595
596template<class Type>
598{
599 tmp<scalarField> tdiag(new scalarField(diag()));
600 addCmptAvBoundaryDiag(tdiag.ref());
601 return tdiag;
602}
603
604
605template<class Type>
607{
609 (
611 (
613 (
614 "A("+psi_.name()+')',
615 psi_.instance(),
616 psi_.db()
617 ),
618 psi_.mesh(),
619 dimensions_/psi_.dimensions()/dimArea,
620 zeroGradientFaPatchScalarField::typeName
621 )
622 );
623
624 tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().S();
625 tAphi.ref().correctBoundaryConditions();
626
627 return tAphi;
628}
629
630
631template<class Type>
634{
636 (
638 (
640 (
641 "H("+psi_.name()+')',
642 psi_.instance(),
643 psi_.db()
644 ),
645 psi_.mesh(),
646 dimensions_/dimArea,
647 zeroGradientFaPatchScalarField::typeName
648 )
649 );
651
652 // Loop over field components
653 for (direction cmpt=0; cmpt<Type::nComponents; ++cmpt)
654 {
655 scalarField psiCmpt(psi_.primitiveField().component(cmpt));
656
657 scalarField boundaryDiagCmpt(psi_.size(), Zero);
658 addBoundaryDiag(boundaryDiagCmpt, cmpt);
659 boundaryDiagCmpt.negate();
660 addCmptAvBoundaryDiag(boundaryDiagCmpt);
661
662 Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
663 }
664
665 Hphi.primitiveFieldRef() += lduMatrix::H(psi_.primitiveField()) + source_;
666 addBoundarySource(Hphi.primitiveFieldRef());
667
668 Hphi.primitiveFieldRef() /= psi_.mesh().S();
670
671 return tHphi;
672}
673
674
675template<class Type>
678{
679 if (!psi_.mesh().fluxRequired(psi_.name()))
680 {
682 << "flux requested but " << psi_.name()
683 << " not specified in the fluxRequired sub-dictionary of faSchemes"
684 << abort(FatalError);
685 }
686
687 // construct GeometricField<Type, faePatchField, edgeMesh>
689 (
691 (
693 (
694 "flux("+psi_.name()+')',
695 psi_.instance(),
696 psi_.db()
697 ),
698 psi_.mesh(),
699 dimensions()
700 )
701 );
703 tfieldFlux.ref();
704
705 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
706 {
707 fieldFlux.primitiveFieldRef().replace
708 (
709 cmpt,
710 lduMatrix::faceH(psi_.primitiveField().component(cmpt))
711 );
712 }
713
714 FieldField<Field, Type> InternalContrib = internalCoeffs_;
715
716 forAll(InternalContrib, patchI)
717 {
718 InternalContrib[patchI] =
720 (
721 InternalContrib[patchI],
722 psi_.boundaryField()[patchI].patchInternalField()
723 );
724 }
725
726 FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
727
728 forAll(NeighbourContrib, patchI)
729 {
730 if (psi_.boundaryField()[patchI].coupled())
731 {
732 NeighbourContrib[patchI] =
734 (
735 NeighbourContrib[patchI],
736 psi_.boundaryField()[patchI].patchNeighbourField()
737 );
738 }
739 }
740
741 forAll(fieldFlux.boundaryField(), patchI)
742 {
743 fieldFlux.boundaryFieldRef()[patchI] =
744 InternalContrib[patchI] - NeighbourContrib[patchI];
745 }
746
747 if (faceFluxCorrectionPtr_)
748 {
749 fieldFlux += *faceFluxCorrectionPtr_;
750 }
751
752 return tfieldFlux;
753}
754
755
756// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
757
758template<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
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_ =
786 (*famv.faceFluxCorrectionPtr_);
787 }
788}
789
790
791template<class Type>
793{
794 operator=(tfamv());
795 tfamv.clear();
796}
797
798
799template<class Type>
801{
803 source_.negate();
804 internalCoeffs_.negate();
805 boundaryCoeffs_.negate();
806
807 if (faceFluxCorrectionPtr_)
808 {
809 faceFluxCorrectionPtr_->negate();
810 }
811}
812
813
814template<class Type>
816{
817 checkMethod(*this, famv, "+=");
818
819 dimensions_ += famv.dimensions_;
821 source_ += famv.source_;
822 internalCoeffs_ += famv.internalCoeffs_;
823 boundaryCoeffs_ += famv.boundaryCoeffs_;
824
825 if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
826 {
827 *faceFluxCorrectionPtr_ += *famv.faceFluxCorrectionPtr_;
828 }
829 else if (famv.faceFluxCorrectionPtr_)
830 {
831 faceFluxCorrectionPtr_ = new
833 (
834 *famv.faceFluxCorrectionPtr_
835 );
836 }
837}
838
839
840template<class Type>
842{
843 operator+=(tfamv());
844 tfamv.clear();
845}
846
847
848template<class Type>
850{
851 checkMethod(*this, famv, "+=");
852
853 dimensions_ -= famv.dimensions_;
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_ =
867 (-*famv.faceFluxCorrectionPtr_);
868 }
869}
870
871
872template<class Type>
874{
875 operator-=(tfamv());
876 tfamv.clear();
877}
878
879
880template<class Type>
882(
884)
885{
886 checkMethod(*this, su, "+=");
887 source() -= su.mesh().S()*su.field();
888}
889
890
891template<class Type>
893(
895)
896{
897 operator+=(tsu());
898 tsu.clear();
899}
900
901
902template<class Type>
904(
906)
907{
908 operator+=(tsu());
909 tsu.clear();
910}
911
912
913template<class Type>
915(
917)
918{
919 checkMethod(*this, su, "-=");
920 source() += su.mesh().S()*su.field();
921}
922
923
924template<class Type>
926(
928)
929{
930 operator-=(tsu());
931 tsu.clear();
932}
933
934
935template<class Type>
937(
939)
940{
941 operator-=(tsu());
942 tsu.clear();
943}
944
945
946template<class Type>
948(
949 const dimensioned<Type>& su
950)
951{
952 source() -= psi().mesh().S()*su;
953}
954
955
956template<class Type>
958(
959 const dimensioned<Type>& su
960)
961{
962 source() += psi().mesh().S()*su;
963}
964
965
966template<class Type>
968{}
969
970
971template<class Type>
973{}
974
975
976template<class Type>
978(
980)
981{
982 dimensions_ *= dsf.dimensions();
983 lduMatrix::operator*=(dsf.field());
984 source_ *= dsf.field();
985
986 forAll(boundaryCoeffs_, patchi)
987 {
988 const scalarField pisf
989 (
990 dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
991 );
992 internalCoeffs_[patchi] *= pisf;
993 boundaryCoeffs_[patchi] *= pisf;
994 }
995
996 if (faceFluxCorrectionPtr_)
997 {
999 << "cannot scale a matrix containing a faceFluxCorrection"
1000 << abort(FatalError);
1001 }
1002}
1003
1004
1005template<class Type>
1007(
1009)
1010{
1011 operator*=(tfld());
1012 tfld.clear();
1013}
1014
1015
1016template<class Type>
1018(
1019 const tmp<areaScalarField>& tfld
1020)
1021{
1022 operator*=(tfld());
1023 tfld.clear();
1024}
1025
1026
1027template<class Type>
1029(
1030 const dimensioned<scalar>& ds
1031)
1032{
1033 dimensions_ *= ds.dimensions();
1034 lduMatrix::operator*=(ds.value());
1035 source_ *= ds.value();
1036 internalCoeffs_ *= ds.value();
1037 boundaryCoeffs_ *= ds.value();
1038
1039 if (faceFluxCorrectionPtr_)
1040 {
1041 *faceFluxCorrectionPtr_ *= ds.value();
1042 }
1043}
1044
1045
1046// * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
1047
1048template<class Type>
1050(
1051 const faMatrix<Type>& mat1,
1052 const faMatrix<Type>& mat2,
1053 const char* op
1054)
1055{
1056 if (&mat1.psi() != &mat2.psi())
1057 {
1059 << "Incompatible fields for operation\n "
1060 << "[" << mat1.psi().name() << "] "
1061 << op
1062 << " [" << mat2.psi().name() << "]"
1063 << abort(FatalError);
1064 }
1065
1066 if
1067 (
1069 && mat1.dimensions() != mat2.dimensions()
1070 )
1071 {
1073 << "Incompatible dimensions for operation\n "
1074 << "[" << mat1.psi().name() << mat1.dimensions()/dimArea << " ] "
1075 << op
1076 << " [" << mat2.psi().name() << mat2.dimensions()/dimArea << " ]"
1077 << abort(FatalError);
1078 }
1079}
1080
1081
1082template<class Type>
1084(
1085 const faMatrix<Type>& mat,
1087 const char* op
1088)
1089{
1090 if
1091 (
1093 && mat.dimensions()/dimArea != fld.dimensions()
1094 )
1095 {
1097 << "Incompatible dimensions for operation\n "
1098 << "[" << mat.psi().name() << mat.dimensions()/dimArea << " ] "
1099 << op
1100 << " [" << fld.name() << fld.dimensions() << " ]"
1101 << abort(FatalError);
1102 }
1103}
1104
1105
1106template<class Type>
1108(
1109 const faMatrix<Type>& mat,
1110 const dimensioned<Type>& dt,
1111 const char* op
1112)
1113{
1114 if
1115 (
1117 && mat.dimensions()/dimArea != dt.dimensions()
1118 )
1119 {
1121 << "Incompatible dimensions for operation\n "
1122 << "[" << mat.psi().name() << mat.dimensions()/dimArea << " ] "
1123 << op
1124 << " [" << dt.name() << dt.dimensions() << " ]"
1125 << abort(FatalError);
1126 }
1127}
1128
1129
1130template<class Type>
1132(
1133 faMatrix<Type>& fam,
1134 const dictionary& solverControls
1135)
1136{
1137 return fam.solve(solverControls);
1138}
1139
1140
1141template<class Type>
1143(
1144 const tmp<faMatrix<Type>>& tmat,
1145 const dictionary& solverControls
1146)
1147{
1148 SolverPerformance<Type> solverPerf(tmat.constCast().solve(solverControls));
1149
1150 tmat.clear();
1151
1152 return solverPerf;
1153}
1154
1155
1156template<class Type>
1158{
1159 return fam.solve();
1160}
1161
1162
1163template<class Type>
1164Foam::SolverPerformance<Type> Foam::solve(const tmp<faMatrix<Type>>& tmat)
1165{
1166 SolverPerformance<Type> solverPerf(tmat.constCast().solve());
1167
1168 tmat.clear();
1169
1170 return solverPerf;
1171}
1172
1173
1174// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
1175
1176template<class Type>
1177Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1178(
1179 const faMatrix<Type>& A,
1180 const faMatrix<Type>& B
1181)
1182{
1183 checkMethod(A, B, "+");
1184 tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1185 tC.ref() += B;
1186 return tC;
1187}
1188
1189
1190template<class Type>
1191Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1192(
1193 const tmp<faMatrix<Type>>& tA,
1194 const faMatrix<Type>& B
1195)
1196{
1197 checkMethod(tA(), B, "+");
1198 tmp<faMatrix<Type>> tC(tA.ptr());
1199 tC.ref() += B;
1200 return tC;
1201}
1202
1203
1204template<class Type>
1205Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1206(
1207 const faMatrix<Type>& A,
1208 const tmp<faMatrix<Type>>& tB
1209)
1210{
1211 checkMethod(A, tB(), "+");
1212 tmp<faMatrix<Type>> tC(tB.ptr());
1213 tC.ref() += A;
1214 return tC;
1215}
1216
1217
1218template<class Type>
1219Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1220(
1221 const tmp<faMatrix<Type>>& tA,
1222 const tmp<faMatrix<Type>>& tB
1223)
1224{
1225 checkMethod(tA(), tB(), "+");
1226 tmp<faMatrix<Type>> tC(tA.ptr());
1227 tC.ref() += tB();
1228 tB.clear();
1229 return tC;
1230}
1231
1232
1233template<class Type>
1234Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1235(
1236 const faMatrix<Type>& A
1237)
1238{
1239 tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1240 tC.ref().negate();
1241 return tC;
1242}
1243
1244
1245template<class Type>
1246Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1247(
1248 const tmp<faMatrix<Type>>& tA
1249)
1250{
1251 tmp<faMatrix<Type>> tC(tA.ptr());
1252 tC.ref().negate();
1253 return tC;
1254}
1255
1256
1257template<class Type>
1258Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1259(
1260 const faMatrix<Type>& A,
1261 const faMatrix<Type>& B
1262)
1263{
1264 checkMethod(A, B, "-");
1265 tmp<faMatrix<Type>> tC(new faMatrix<Type>(A));
1266 tC.ref() -= B;
1267 return tC;
1268}
1269
1270
1271template<class Type>
1272Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1273(
1274 const tmp<faMatrix<Type>>& tA,
1275 const faMatrix<Type>& B
1276)
1277{
1278 checkMethod(tA(), B, "-");
1279 tmp<faMatrix<Type>> tC(tA.ptr());
1280 tC.ref() -= B;
1281 return tC;
1282}
1283
1284
1285template<class Type>
1286Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1287(
1288 const faMatrix<Type>& A,
1289 const tmp<faMatrix<Type>>& tB
1290)
1291{
1292 checkMethod(A, tB(), "-");
1293 tmp<faMatrix<Type>> tC(tB.ptr());
1294 tC.ref() -= A;
1295 tC.ref().negate();
1296 return tC;
1297}
1298
1299
1300template<class Type>
1301Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1302(
1303 const tmp<faMatrix<Type>>& tA,
1304 const tmp<faMatrix<Type>>& tB
1305)
1306{
1307 checkMethod(tA(), tB(), "-");
1308 tmp<faMatrix<Type>> tC(tA.ptr());
1309 tC.ref() -= tB();
1310 tB.clear();
1311 return tC;
1312}
1313
1314
1315template<class Type>
1316Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1317(
1318 const faMatrix<Type>& A,
1319 const faMatrix<Type>& B
1320)
1321{
1322 checkMethod(A, B, "==");
1323 return (A - B);
1324}
1325
1326
1327template<class Type>
1328Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1329(
1330 const tmp<faMatrix<Type>>& tA,
1331 const faMatrix<Type>& B
1332)
1333{
1334 checkMethod(tA(), B, "==");
1335 return (tA - B);
1336}
1337
1338
1339template<class Type>
1340Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1341(
1342 const faMatrix<Type>& A,
1343 const tmp<faMatrix<Type>>& tB
1344)
1345{
1346 checkMethod(A, tB(), "==");
1347 return (A - tB);
1348}
1349
1350
1351template<class Type>
1352Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1353(
1354 const tmp<faMatrix<Type>>& tA,
1355 const tmp<faMatrix<Type>>& tB
1356)
1357{
1358 checkMethod(tA(), tB(), "==");
1359 return (tA - tB);
1360}
1361
1362
1363template<class Type>
1364Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1365(
1366 const faMatrix<Type>& A,
1367 const DimensionedField<Type, areaMesh>& su
1368)
1369{
1370 checkMethod(A, su, "+");
1371 auto tC(A.clone());
1372 tC.ref().source() -= su.mesh().S()*su.field();
1373 return tC;
1374}
1375
1376
1377template<class Type>
1378Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1379(
1380 const faMatrix<Type>& A,
1381 const tmp<DimensionedField<Type, areaMesh>>& tsu
1382)
1383{
1384 checkMethod(A, tsu(), "+");
1385 auto tC(A.clone());
1386 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1387 tsu.clear();
1388 return tC;
1389}
1390
1391
1392template<class Type>
1393Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1394(
1395 const faMatrix<Type>& A,
1396 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1397)
1398{
1399 checkMethod(A, tsu(), "+");
1400 auto tC(A.clone());
1401 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1402 tsu.clear();
1403 return tC;
1404}
1405
1406
1407template<class Type>
1408Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1409(
1410 const tmp<faMatrix<Type>>& tA,
1411 const DimensionedField<Type, areaMesh>& su
1412)
1413{
1414 checkMethod(tA(), su, "+");
1415 tmp<faMatrix<Type>> tC(tA.ptr());
1416 tC.ref().source() -= su.mesh().S()*su.field();
1417 return tC;
1418}
1419
1420
1421template<class Type>
1422Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1423(
1424 const tmp<faMatrix<Type>>& tA,
1425 const tmp<DimensionedField<Type, areaMesh>>& tsu
1426)
1427{
1428 checkMethod(tA(), tsu(), "+");
1429 tmp<faMatrix<Type>> tC(tA.ptr());
1430 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1431 tsu.clear();
1432 return tC;
1433}
1434
1435
1436template<class Type>
1437Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1438(
1439 const tmp<faMatrix<Type>>& tA,
1440 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1441)
1442{
1443 checkMethod(tA(), tsu(), "+");
1444 tmp<faMatrix<Type>> tC(tA.ptr());
1445 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1446 tsu.clear();
1447 return tC;
1448}
1449
1450
1451template<class Type>
1452Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1453(
1454 const DimensionedField<Type, areaMesh>& su,
1455 const faMatrix<Type>& A
1456)
1457{
1458 checkMethod(A, su, "+");
1459 auto tC(A.clone());
1460 tC.ref().source() -= su.mesh().S()*su.field();
1461 return tC;
1462}
1463
1464
1465template<class Type>
1466Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1467(
1468 const tmp<DimensionedField<Type, areaMesh>>& tsu,
1469 const faMatrix<Type>& A
1470)
1471{
1472 checkMethod(A, tsu(), "+");
1473 auto tC(A.clone());
1474 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1475 tsu.clear();
1476 return tC;
1477}
1478
1479
1480template<class Type>
1481Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1482(
1483 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1484 const faMatrix<Type>& A
1485)
1486{
1487 checkMethod(A, tsu(), "+");
1488 auto tC(A.clone());
1489 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1490 tsu.clear();
1491 return tC;
1492}
1493
1494
1495template<class Type>
1496Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1497(
1498 const DimensionedField<Type, areaMesh>& su,
1499 const tmp<faMatrix<Type>>& tA
1500)
1501{
1502 checkMethod(tA(), su, "+");
1503 tmp<faMatrix<Type>> tC(tA.ptr());
1504 tC.ref().source() -= su.mesh().S()*su.field();
1505 return tC;
1506}
1507
1508
1509template<class Type>
1510Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1511(
1512 const tmp<DimensionedField<Type, areaMesh>>& tsu,
1513 const tmp<faMatrix<Type>>& tA
1514)
1515{
1516 checkMethod(tA(), tsu(), "+");
1517 tmp<faMatrix<Type>> tC(tA.ptr());
1518 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1519 tsu.clear();
1520 return tC;
1521}
1522
1523
1524template<class Type>
1525Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1526(
1527 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1528 const tmp<faMatrix<Type>>& tA
1529)
1530{
1531 checkMethod(tA(), tsu(), "+");
1532 tmp<faMatrix<Type>> tC(tA.ptr());
1533 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1534 tsu.clear();
1535 return tC;
1536}
1537
1538
1539template<class Type>
1540Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1541(
1542 const faMatrix<Type>& A,
1543 const DimensionedField<Type, areaMesh>& su
1544)
1545{
1546 checkMethod(A, su, "-");
1547 auto tC(A.clone());
1548 tC.ref().source() += su.mesh().S()*su.field();
1549 return tC;
1550}
1551
1552
1553template<class Type>
1554Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1555(
1556 const faMatrix<Type>& A,
1557 const tmp<DimensionedField<Type, areaMesh>>& tsu
1558)
1559{
1560 checkMethod(A, tsu(), "-");
1561 auto tC(A.clone());
1562 tC.ref().source() += tsu().mesh().S()*tsu().field();
1563 tsu.clear();
1564 return tC;
1565}
1566
1567
1568template<class Type>
1569Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1570(
1571 const faMatrix<Type>& A,
1572 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1573)
1574{
1575 checkMethod(A, tsu(), "-");
1576 auto tC(A.clone());
1577 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1578 tsu.clear();
1579 return tC;
1580}
1581
1582
1583template<class Type>
1584Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1585(
1586 const tmp<faMatrix<Type>>& tA,
1587 const DimensionedField<Type, areaMesh>& su
1588)
1589{
1590 checkMethod(tA(), su, "-");
1591 tmp<faMatrix<Type>> tC(tA.ptr());
1592 tC.ref().source() += su.mesh().S()*su.field();
1593 return tC;
1594}
1595
1596
1597template<class Type>
1598Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1599(
1600 const tmp<faMatrix<Type>>& tA,
1601 const tmp<DimensionedField<Type, areaMesh>>& tsu
1602)
1603{
1604 checkMethod(tA(), tsu(), "-");
1605 tmp<faMatrix<Type>> tC(tA.ptr());
1606 tC.ref().source() += tsu().mesh().S()*tsu().field();
1607 tsu.clear();
1608 return tC;
1609}
1610
1611
1612template<class Type>
1613Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1614(
1615 const tmp<faMatrix<Type>>& tA,
1616 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1617)
1618{
1619 checkMethod(tA(), tsu(), "-");
1620 tmp<faMatrix<Type>> tC(tA.ptr());
1621 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1622 tsu.clear();
1623 return tC;
1624}
1625
1626
1627template<class Type>
1628Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1629(
1630 const DimensionedField<Type, areaMesh>& su,
1631 const faMatrix<Type>& A
1632)
1633{
1634 checkMethod(A, su, "-");
1635 auto tC(A.clone());
1636 tC.ref().negate();
1637 tC.ref().source() -= su.mesh().S()*su.field();
1638 return tC;
1639}
1640
1641
1642template<class Type>
1643Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1644(
1645 const tmp<DimensionedField<Type, areaMesh>>& tsu,
1646 const faMatrix<Type>& A
1647)
1648{
1649 checkMethod(A, tsu(), "-");
1650 auto tC(A.clone());
1651 tC.ref().negate();
1652 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1653 tsu.clear();
1654 return tC;
1655}
1656
1657
1658template<class Type>
1659Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1660(
1661 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1662 const faMatrix<Type>& A
1663)
1664{
1665 checkMethod(A, tsu(), "-");
1666 auto tC(A.clone());
1667 tC.ref().negate();
1668 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1669 tsu.clear();
1670 return tC;
1671}
1672
1673
1674template<class Type>
1675Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1676(
1677 const DimensionedField<Type, areaMesh>& su,
1678 const tmp<faMatrix<Type>>& tA
1679)
1680{
1681 checkMethod(tA(), su, "-");
1682 tmp<faMatrix<Type>> tC(tA.ptr());
1683 tC.ref().negate();
1684 tC.ref().source() -= su.mesh().S()*su.field();
1685 return tC;
1686}
1687
1688
1689template<class Type>
1690Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1691(
1692 const tmp<DimensionedField<Type, areaMesh>>& tsu,
1693 const tmp<faMatrix<Type>>& tA
1694)
1695{
1696 checkMethod(tA(), tsu(), "-");
1697 tmp<faMatrix<Type>> tC(tA.ptr());
1698 tC.ref().negate();
1699 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1700 tsu.clear();
1701 return tC;
1702}
1703
1704
1705template<class Type>
1706Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1707(
1708 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu,
1709 const tmp<faMatrix<Type>>& tA
1710)
1711{
1712 checkMethod(tA(), tsu(), "-");
1713 tmp<faMatrix<Type>> tC(tA.ptr());
1714 tC.ref().negate();
1715 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1716 tsu.clear();
1717 return tC;
1718}
1719
1720
1721template<class Type>
1722Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1723(
1724 const faMatrix<Type>& A,
1725 const dimensioned<Type>& su
1726)
1727{
1728 checkMethod(A, su, "+");
1729 auto tC(A.clone());
1730 tC.ref().source() -= su.value()*A.psi().mesh().S();
1731 return tC;
1732}
1733
1734
1735template<class Type>
1736Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1737(
1738 const tmp<faMatrix<Type>>& tA,
1739 const dimensioned<Type>& su
1740)
1741{
1742 checkMethod(tA(), su, "+");
1743 tmp<faMatrix<Type>> tC(tA.ptr());
1744 tC.ref().source() -= su.value()*tC().psi().mesh().S();
1745 return tC;
1746}
1747
1748
1749template<class Type>
1750Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1751(
1752 const dimensioned<Type>& su,
1753 const faMatrix<Type>& A
1754)
1755{
1756 checkMethod(A, su, "+");
1757 auto tC(A.clone());
1758 tC.ref().source() -= su.value()*A.psi().mesh().S();
1759 return tC;
1760}
1761
1762
1763template<class Type>
1764Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1765(
1766 const dimensioned<Type>& su,
1767 const tmp<faMatrix<Type>>& tA
1768)
1769{
1770 checkMethod(tA(), su, "+");
1771 tmp<faMatrix<Type>> tC(tA.ptr());
1772 tC.ref().source() -= su.value()*tC().psi().mesh().S();
1773 return tC;
1774}
1775
1776
1777template<class Type>
1778Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1779(
1780 const faMatrix<Type>& A,
1781 const dimensioned<Type>& su
1782)
1783{
1784 checkMethod(A, su, "-");
1785 auto tC(A.clone());
1786 tC.ref().source() += su.value()*tC().psi().mesh().S();
1787 return tC;
1788}
1789
1790
1791template<class Type>
1792Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1793(
1794 const tmp<faMatrix<Type>>& tA,
1795 const dimensioned<Type>& su
1796)
1797{
1798 checkMethod(tA(), su, "-");
1799 tmp<faMatrix<Type>> tC(tA.ptr());
1800 tC.ref().source() += su.value()*tC().psi().mesh().S();
1801 return tC;
1802}
1803
1804
1805template<class Type>
1806Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1807(
1808 const dimensioned<Type>& su,
1809 const faMatrix<Type>& A
1810)
1811{
1812 checkMethod(A, su, "-");
1813 auto tC(A.clone());
1814 tC.ref().negate();
1815 tC.ref().source() -= su.value()*A.psi().mesh().S();
1816 return tC;
1817}
1818
1819
1820template<class Type>
1821Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1822(
1823 const dimensioned<Type>& su,
1824 const tmp<faMatrix<Type>>& tA
1825)
1826{
1827 checkMethod(tA(), su, "-");
1828 tmp<faMatrix<Type>> tC(tA.ptr());
1829 tC.ref().negate();
1830 tC.ref().source() -= su.value()*tC().psi().mesh().S();
1831 return tC;
1832}
1833
1834
1835template<class Type>
1836Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1837(
1838 const faMatrix<Type>& A,
1839 const DimensionedField<Type, areaMesh>& su
1840)
1841{
1842 checkMethod(A, su, "==");
1843 auto tC(A.clone());
1844 tC.ref().source() += su.mesh().S()*su.field();
1845 return tC;
1846}
1847
1848
1849template<class Type>
1850Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1851(
1852 const faMatrix<Type>& A,
1853 const tmp<DimensionedField<Type, areaMesh>>& tsu
1854)
1855{
1856 checkMethod(A, tsu(), "==");
1857 auto tC(A.clone());
1858 tC.ref().source() += tsu().mesh().S()*tsu().field();
1859 tsu.clear();
1860 return tC;
1861}
1862
1863
1864template<class Type>
1865Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1866(
1867 const faMatrix<Type>& A,
1868 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1869)
1870{
1871 checkMethod(A, tsu(), "==");
1872 auto tC(A.clone());
1873 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1874 tsu.clear();
1875 return tC;
1876}
1877
1878
1879template<class Type>
1880Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1881(
1882 const tmp<faMatrix<Type>>& tA,
1883 const DimensionedField<Type, areaMesh>& su
1884)
1885{
1886 checkMethod(tA(), su, "==");
1887 tmp<faMatrix<Type>> tC(tA.ptr());
1888 tC.ref().source() += su.mesh().S()*su.field();
1889 return tC;
1890}
1891
1892
1893template<class Type>
1894Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1895(
1896 const tmp<faMatrix<Type>>& tA,
1897 const tmp<DimensionedField<Type, areaMesh>>& tsu
1898)
1899{
1900 checkMethod(tA(), tsu(), "==");
1901 tmp<faMatrix<Type>> tC(tA.ptr());
1902 tC.ref().source() += tsu().mesh().S()*tsu().field();
1903 tsu.clear();
1904 return tC;
1905}
1906
1907
1908template<class Type>
1909Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1910(
1911 const tmp<faMatrix<Type>>& tA,
1912 const tmp<GeometricField<Type, faPatchField, areaMesh>>& tsu
1913)
1914{
1915 checkMethod(tA(), tsu(), "==");
1916 tmp<faMatrix<Type>> tC(tA.ptr());
1917 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1918 tsu.clear();
1919 return tC;
1920}
1921
1922
1923template<class Type>
1924Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1925(
1926 const faMatrix<Type>& A,
1927 const dimensioned<Type>& su
1928)
1929{
1930 checkMethod(A, su, "==");
1931 auto tC(A.clone());
1932 tC.ref().source() += A.psi().mesh().S()*su.value();
1933 return tC;
1934}
1935
1936
1937template<class Type>
1938Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1939(
1940 const tmp<faMatrix<Type>>& tA,
1941 const dimensioned<Type>& su
1942)
1943{
1944 checkMethod(tA(), su, "==");
1945 tmp<faMatrix<Type>> tC(tA.ptr());
1946 tC.ref().source() += tC().psi().mesh().S()*su.value();
1947 return tC;
1948}
1949
1950
1951template<class Type>
1952Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1953(
1954 const faMatrix<Type>& A,
1955 const Foam::zero
1956)
1957{
1958 return A;
1959}
1960
1961
1962template<class Type>
1963Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1964(
1965 const tmp<faMatrix<Type>>& tA,
1966 const Foam::zero
1967)
1968{
1969 return tA;
1970}
1971
1972
1973template<class Type>
1974Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1975(
1976 const areaScalarField::Internal& dsf,
1977 const faMatrix<Type>& A
1978)
1979{
1980 auto tC(A.clone());
1981 tC.ref() *= dsf;
1982 return tC;
1983}
1984
1985
1986template<class Type>
1987Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1988(
1989 const tmp<areaScalarField::Internal>& tdsf,
1990 const faMatrix<Type>& A
1991)
1992{
1993 auto tC(A.clone());
1994 tC.ref() *= tdsf;
1995 return tC;
1996}
1997
1998
1999template<class Type>
2000Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2001(
2002 const tmp<areaScalarField>& tvsf,
2003 const faMatrix<Type>& A
2004)
2005{
2006 auto tC(A.clone());
2007 tC.ref() *= tvsf;
2008 return tC;
2009}
2010
2011
2012template<class Type>
2013Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2014(
2015 const areaScalarField::Internal& dsf,
2016 const tmp<faMatrix<Type>>& tA
2017)
2018{
2019 tmp<faMatrix<Type>> tC(tA.ptr());
2020 tC.ref() *= dsf;
2021 return tC;
2022}
2023
2024
2025template<class Type>
2026Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2027(
2028 const tmp<areaScalarField::Internal>& tdsf,
2029 const tmp<faMatrix<Type>>& tA
2030)
2031{
2032 tmp<faMatrix<Type>> tC(tA.ptr());
2033 tC.ref() *= tdsf;
2034 return tC;
2035}
2036
2037
2038template<class Type>
2039Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2040(
2041 const tmp<areaScalarField>& tvsf,
2042 const tmp<faMatrix<Type>>& tA
2043)
2044{
2045 tmp<faMatrix<Type>> tC(tA.ptr());
2046 tC.ref() *= tvsf;
2047 return tC;
2048}
2049
2050
2051template<class Type>
2052Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2053(
2054 const dimensioned<scalar>& ds,
2055 const faMatrix<Type>& A
2056)
2057{
2058 auto tC(A.clone());
2059 tC.ref() *= ds;
2060 return tC;
2061}
2062
2063
2064template<class Type>
2065Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2066(
2067 const dimensioned<scalar>& ds,
2068 const tmp<faMatrix<Type>>& tA
2069)
2070{
2071 tmp<faMatrix<Type>> tC(tA.ptr());
2072 tC.ref() *= ds;
2073 return tC;
2074}
2075
2076
2077// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2078
2079template<class Type>
2080Foam::Ostream& Foam::operator<<(Ostream& os, const faMatrix<Type>& fam)
2081{
2082 os << static_cast<const lduMatrix&>(fam) << nl
2083 << fam.dimensions_ << nl
2084 << fam.source_ << nl
2085 << fam.internalCoeffs_ << nl
2086 << fam.boundaryCoeffs_ << endl;
2087
2089
2090 return os;
2091}
2092
2093
2094// * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
2095
2096#include "faMatrixSolve.C"
2097
2098// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Generic templated field type.
Definition: Field.H:82
void replace(const direction, const UList< cmptType > &)
Replace a component field of the field.
Definition: Field.C:557
void negate()
Inplace negate this field (negative).
Definition: Field.C:530
Generic GeometricField class.
DimensionedField< scalar, areaMesh > Internal
The internal field type from which this GeometricField is derived.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
A List with indirect addressing. Like IndirectList but does not store addressing.
Definition: IndirectList.H:79
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A single value that is represented as a list with an operator[] to access the value....
Definition: UniformList.H:55
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
static bool checking() noexcept
True if dimension checking is enabled (the usual default)
Definition: dimensionSet.H:229
Generic dimensioned Type class.
const dimensionSet & dimensions() const
Return const reference to dimensions.
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:76
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: faMatrix.C:42
void operator=(const faMatrix< Type > &)
Definition: faMatrix.C:759
void relax()
Relax matrix (for steady-state solution).
Definition: faMatrix.C:581
tmp< GeometricField< Type, faePatchField, edgeMesh > > flux() const
Return the face-flux field from the matrix.
Definition: faMatrix.C:677
void setValues(const labelUList &faceLabels, const Type &value)
Definition: faMatrix.C:394
void setReference(const label facei, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: faMatrix.C:427
void setReferences(const labelUList &faceLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: faMatrix.C:446
void operator+=(const faMatrix< Type > &)
Definition: faMatrix.C:815
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:270
const dimensionSet & dimensions() const noexcept
Definition: faMatrix.H:275
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: faMatrix.C:135
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: faMatrix.C:151
virtual ~faMatrix()
Destructor.
Definition: faMatrix.C:291
tmp< areaScalarField > A() const
Return the central coefficient.
Definition: faMatrix.C:606
void operator-=(const faMatrix< Type > &)
Definition: faMatrix.C:849
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: faMatrix.C:117
void negate()
Inplace negate.
Definition: faMatrix.C:800
void setValuesFromList(const labelUList &faceLabels, const ListType< Type > &values)
Set solution in given faces to the specified values.
Definition: faMatrix.C:305
tmp< scalarField > D() const
Return the matrix diagonal.
Definition: faMatrix.C:597
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: faMatrix.C:80
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
Definition: faPatchField.H:82
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: faPatchField.H:371
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: faPatchField.H:327
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:417
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:423
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
void operator=(const lduMatrix &)
tmp< Field< Type > > faceH(const Field< Type > &) const
void operator*=(const scalarField &)
void operator+=(const lduMatrix &)
void operator-=(const lduMatrix &)
const labelListList & faceEdges() const
label eventNo() const noexcept
Event number at last update.
Definition: regIOobjectI.H:191
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition: tmpI.H:255
T & ref() const
Definition: tmpI.H:227
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
const volScalarField & psi
UEqn relax()
faceListList boundary
dynamicFvMesh & mesh
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Finite-Area matrix basic solvers.
OBJstream os(runTime.globalPath()/outputName)
label faceId(-1)
#define FUNCTION_NAME
#define DebugInFunction
Report an information message using Foam::Info.
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1050
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &df)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
void deleteDemandDrivenData(DataPtr &dataPtr)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Calculate the matrix for the second temporal derivative.
volScalarField & alpha
const dimensionedScalar & D
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
conserve primitiveFieldRef()+