Loading...
Searching...
No Matches
faMatrix.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2017 Wikki Ltd
9 Copyright (C) 2019-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
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 "UniformList.H"
34
35// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36
37template<class Type>
38template<class Type2>
40(
41 const labelUList& addr,
42 const Field<Type2>& pf,
43 Field<Type2>& intf
44) const
45{
46 if (addr.size() != pf.size())
47 {
49 << "addressing (" << addr.size()
50 << ") and field (" << pf.size() << ") are different sizes"
51 << abort(FatalError);
52 }
53
54 forAll(addr, faceI)
55 {
56 intf[addr[faceI]] += pf[faceI];
57 }
58}
59
60
61template<class Type>
62template<class Type2>
64(
65 const labelUList& addr,
66 const tmp<Field<Type2>>& tpf,
67 Field<Type2>& intf
68) const
69{
70 addToInternalField(addr, tpf(), intf);
71 tpf.clear();
72}
73
74
75template<class Type>
76template<class Type2>
78(
79 const labelUList& addr,
80 const Field<Type2>& pf,
81 Field<Type2>& intf
82) const
83{
84 if (addr.size() != pf.size())
85 {
87 << "addressing (" << addr.size()
88 << ") and field (" << pf.size() << ") are different sizes"
89 << abort(FatalError);
90 }
91
92 forAll(addr, faceI)
93 {
94 intf[addr[faceI]] -= pf[faceI];
95 }
96}
97
98
99template<class Type>
100template<class Type2>
102(
103 const labelUList& addr,
104 const tmp<Field<Type2>>& tpf,
105 Field<Type2>& intf
106) const
108 subtractFromInternalField(addr, tpf(), intf);
109 tpf.clear();
110}
111
112
113template<class Type>
115(
117 const direction solveCmpt
118) const
119{
120 forAll(internalCoeffs_, patchI)
121 {
122 addToInternalField
123 (
124 lduAddr().patchAddr(patchI),
125 internalCoeffs_[patchI].component(solveCmpt),
127 );
128 }
129}
130
131
132template<class Type>
133void Foam::faMatrix<Type>::addCmptAvBoundaryDiag(scalarField& diag) const
134{
135 forAll(internalCoeffs_, patchI)
136 {
137 addToInternalField
138 (
139 lduAddr().patchAddr(patchI),
140 cmptAv(internalCoeffs_[patchI]),
142 );
143 }
144}
145
146
147template<class Type>
149(
150 Field<Type>& source,
151 const bool couples
152) const
153{
154 forAll(psi_.boundaryField(), patchI)
155 {
156 const faPatchField<Type>& ptf = psi_.boundaryField()[patchI];
157 const Field<Type>& pbc = boundaryCoeffs_[patchI];
158
159 if (!ptf.coupled())
160 {
161 addToInternalField(lduAddr().patchAddr(patchI), pbc, source);
162 }
163 else if (couples)
164 {
165 tmp<Field<Type>> tpnf = ptf.patchNeighbourField();
166 const Field<Type>& pnf = tpnf();
167
168 const labelUList& addr = lduAddr().patchAddr(patchI);
169
170 forAll(addr, facei)
171 {
172 source[addr[facei]] += cmptMultiply(pbc[facei], pnf[facei]);
173 }
174 }
176}
177
178
179// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
180
181template<class Type>
183(
184 const GeometricField<Type, faPatchField, areaMesh>& psi,
185 const dimensionSet& ds
186)
187:
188 lduMatrix(psi.mesh()),
189 psi_(psi),
190 dimensions_(ds),
191 source_(psi.size(), Zero),
192 internalCoeffs_(psi.mesh().boundary().size()),
193 boundaryCoeffs_(psi.mesh().boundary().size())
194{
196 << "constructing faMatrix<Type> for field " << psi_.name()
197 << endl;
198
199 // Initialise coupling coefficients
200 forAll(psi.mesh().boundary(), patchi)
201 {
202 internalCoeffs_.set
203 (
204 patchi,
205 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
206 );
207
208 boundaryCoeffs_.set
209 (
210 patchi,
211 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
212 );
213 }
214
215 // Update the boundary coefficients of psi without changing its event No.
216 auto& psiRef = psi_.constCast();
217
218 const label currentStatePsi = psiRef.eventNo();
219 psiRef.boundaryFieldRef().updateCoeffs();
220 psiRef.eventNo() = currentStatePsi;
221}
222
223
224template<class Type>
226:
227 lduMatrix(fam),
228 psi_(fam.psi_),
229 dimensions_(fam.dimensions_),
230 source_(fam.source_),
231 internalCoeffs_(fam.internalCoeffs_),
232 boundaryCoeffs_(fam.boundaryCoeffs_)
233{
235 << "Copying faMatrix<Type> for field " << psi_.name() << endl;
236
237 if (fam.faceFluxCorrectionPtr_)
238 {
239 faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
240 (
241 *(fam.faceFluxCorrectionPtr_)
242 );
243 }
244}
245
246
247template<class Type>
249:
250 lduMatrix(tmat.constCast(), tmat.movable()),
251 psi_(tmat().psi_),
252 dimensions_(tmat().dimensions_),
253 source_(tmat.constCast().source_, tmat.movable()),
254 internalCoeffs_(tmat.constCast().internalCoeffs_, tmat.movable()),
255 boundaryCoeffs_(tmat.constCast().boundaryCoeffs_, tmat.movable())
256{
258 << "Copy/Move faMatrix<Type> for field " << psi_.name() << endl;
259
260 if (tmat().faceFluxCorrectionPtr_)
261 {
262 if (tmat.movable())
263 {
264 faceFluxCorrectionPtr_ =
265 std::move(tmat.constCast().faceFluxCorrectionPtr_);
266 }
267 else if (tmat().faceFluxCorrectionPtr_)
268 {
269 faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
270 (
271 *(tmat().faceFluxCorrectionPtr_)
272 );
273 }
274 }
275
276 tmat.clear();
277}
278
279
280// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
281
282template<class Type>
284{
286 << "Destroying faMatrix<Type> for field " << psi_.name() << endl;
288
289
290// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
291
292template<class Type>
293template<template<class> class ListType>
295(
296 const labelUList& faceLabels,
297 const ListType<Type>& values
298)
299{
300#if 0 /* Specific to foam-extend */
301 // Record face labels of eliminated equations
302 for (const label i : faceLabels)
303 {
304 this->eliminatedEqns().insert(i);
305 }
306#endif
307
308 const faMesh& mesh = psi_.mesh();
309
310 const labelListList& edges = mesh.patch().faceEdges();
311 const labelUList& own = mesh.owner();
312 const labelUList& nei = mesh.neighbour();
313
314 scalarField& Diag = diag();
315 Field<Type>& psi = psi_.constCast().primitiveFieldRef();
316
317 // Following actions:
318 // - adjust local field psi
319 // - set local matrix to be diagonal (so adjust source)
320 // - cut connections to neighbours
321 // - make (on non-adjusted cells) contribution explicit
322
323 if (symmetric() || asymmetric())
324 {
326 {
327 const label facei = faceLabels[i];
328 const Type& value = values[i];
329
330 for (const label edgei : edges[facei])
331 {
332 if (mesh.isInternalEdge(edgei))
333 {
334 if (symmetric())
335 {
336 if (facei == own[edgei])
337 {
338 source_[nei[edgei]] -= upper()[edgei]*value;
339 }
340 else
341 {
342 source_[own[edgei]] -= upper()[edgei]*value;
343 }
344
345 upper()[edgei] = 0.0;
346 }
347 else
348 {
349 if (facei == own[edgei])
350 {
351 source_[nei[edgei]] -= lower()[edgei]*value;
352 }
353 else
354 {
355 source_[own[edgei]] -= upper()[edgei]*value;
356 }
357
358 upper()[edgei] = 0.0;
359 lower()[edgei] = 0.0;
360 }
361 }
362 else
363 {
364 const label patchi = mesh.boundary().whichPatch(edgei);
365
366 if (internalCoeffs_[patchi].size())
367 {
368 const label patchEdgei =
369 mesh.boundary()[patchi].whichEdge(edgei);
370
371 internalCoeffs_[patchi][patchEdgei] = Zero;
372 boundaryCoeffs_[patchi][patchEdgei] = Zero;
373 }
374 }
375 }
376 }
377 }
378
379 // Note: above loop might have affected source terms on adjusted cells
380 // so make sure to adjust them afterwards
382 {
383 const label facei = faceLabels[i];
384 const Type& value = values[i];
385
386 psi[facei] = value;
387 source_[facei] = value*Diag[facei];
388 }
389}
390
391
392template<class Type>
394(
395 const labelUList& faceLabels,
396 const Type& value
398{
400}
401
402
403template<class Type>
405(
406 const labelUList& faceLabels,
407 const UList<Type>& values
409{
410 this->setValuesFromList(faceLabels, values);
411}
412
413
414template<class Type>
416(
417 const labelUList& faceLabels,
418 const UIndirectList<Type>& values
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];
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];
485 }
486 }
487}
488
489
490template<class Type>
491void Foam::faMatrix<Type>::relax(const scalar alpha)
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 }
575 // Finally add the relaxation contribution to the source.
576 S += (D - D0)*psi_.primitiveField();
577}
578
579
580template<class Type>
582{
583 scalar relaxCoeff = 0;
584
585 if (psi_.mesh().relaxEquation(psi_.name(), relaxCoeff))
586 {
587 relax(relaxCoeff);
588 }
589 else
590 {
592 << "No relaxation specified for field " << psi_.name() << nl;
593 }
594}
595
596
597template<class Type>
599{
601 addCmptAvBoundaryDiag(tdiag.ref());
602 return tdiag;
603}
604
605
606template<class Type>
608{
609 auto tAphi = areaScalarField::New
610 (
611 "A(" + psi_.name() + ')',
612 psi_.mesh(),
613 dimensions_/psi_.dimensions()/dimArea,
615 );
616
617 tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().S();
618 tAphi.ref().correctBoundaryConditions();
620 return tAphi;
621}
622
623
624template<class Type>
627{
629 (
630 "H(" + psi_.name() + ')',
631 psi_.mesh(),
632 dimensions_/dimArea,
634 );
635 auto& Hphi = tHphi.ref();
636
637 // Loop over field components
638 for (direction cmpt=0; cmpt<Type::nComponents; ++cmpt)
639 {
640 scalarField psiCmpt(psi_.primitiveField().component(cmpt));
641
642 scalarField boundaryDiagCmpt(psi_.size(), Zero);
643 addBoundaryDiag(boundaryDiagCmpt, cmpt);
644 boundaryDiagCmpt.negate();
645 addCmptAvBoundaryDiag(boundaryDiagCmpt);
646
647 Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
648 }
649
650 Hphi.primitiveFieldRef() += lduMatrix::H(psi_.primitiveField()) + source_;
651 addBoundarySource(Hphi.primitiveFieldRef());
652
653 Hphi.primitiveFieldRef() /= psi_.mesh().S();
654 Hphi.correctBoundaryConditions();
656 return tHphi;
657}
658
659
660template<class Type>
663{
664 if (!psi_.mesh().fluxRequired(psi_.name()))
665 {
667 << "flux requested but " << psi_.name()
668 << " not specified in the fluxRequired sub-dictionary of faSchemes"
669 << abort(FatalError);
670 }
671
673 (
674 "flux(" + psi_.name() + ')',
675 psi_.mesh(),
676 dimensions(),
677 lduMatrix::faceH(psi_.primitiveField())
678 );
679 auto& fieldFlux = tfieldFlux.ref();
680 // not yet: fieldFlux.setOriented();
681
682
683 FieldField<Field, Type> InternalContrib = internalCoeffs_;
684
685 forAll(InternalContrib, patchI)
686 {
687 InternalContrib[patchI] =
689 (
690 InternalContrib[patchI],
691 psi_.boundaryField()[patchI].patchInternalField()
692 );
693 }
694
695 FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
696
697 forAll(NeighbourContrib, patchI)
698 {
699 if (psi_.boundaryField()[patchI].coupled())
700 {
701 NeighbourContrib[patchI] =
703 (
704 NeighbourContrib[patchI],
705 psi_.boundaryField()[patchI].patchNeighbourField()
706 );
707 }
708 }
709
710 {
711 auto& ffbf = fieldFlux.boundaryFieldRef();
712
713 forAll(ffbf, patchi)
714 {
715 ffbf[patchi] = InternalContrib[patchi] - NeighbourContrib[patchi];
716 }
717 }
718
719 if (faceFluxCorrectionPtr_)
720 {
721 fieldFlux += *faceFluxCorrectionPtr_;
723
724 return tfieldFlux;
725}
726
727
728template<class Type>
730(
731 const word& name
732) const
733{
734 return psi_.mesh().solverDict(name);
735}
736
737
738template<class Type>
740{
741 return psi_.mesh().solverDict
742 (
743 psi_.name()
744 );
745}
746
747
748// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
749
750template<class Type>
752{
753 if (this == &famv)
754 {
755 return; // Self-assignment is a no-op
756 }
757
758 if (&psi_ != &(famv.psi_))
759 {
761 << "different fields"
762 << abort(FatalError);
763 }
764
766 source_ = famv.source_;
767 internalCoeffs_ = famv.internalCoeffs_;
768 boundaryCoeffs_ = famv.boundaryCoeffs_;
769
770 if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
771 {
772 *faceFluxCorrectionPtr_ = *famv.faceFluxCorrectionPtr_;
773 }
774 else if (famv.faceFluxCorrectionPtr_)
775 {
776 faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
777 (
778 *famv.faceFluxCorrectionPtr_
779 );
780 }
781}
782
783
784template<class Type>
785void Foam::faMatrix<Type>::operator=(const tmp<faMatrix<Type>>& tfamv)
787 operator=(tfamv());
788 tfamv.clear();
789}
790
791
792template<class Type>
794{
796 source_.negate();
797 internalCoeffs_.negate();
798 boundaryCoeffs_.negate();
799
800 if (faceFluxCorrectionPtr_)
802 faceFluxCorrectionPtr_->negate();
803 }
804}
805
806
807template<class Type>
808void Foam::faMatrix<Type>::operator+=(const faMatrix<Type>& famv)
809{
810 checkMethod(*this, famv, "+=");
811
812 dimensions_ += famv.dimensions_;
814 source_ += famv.source_;
815 internalCoeffs_ += famv.internalCoeffs_;
816 boundaryCoeffs_ += famv.boundaryCoeffs_;
817
818 if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
819 {
820 *faceFluxCorrectionPtr_ += *famv.faceFluxCorrectionPtr_;
821 }
822 else if (famv.faceFluxCorrectionPtr_)
823 {
824 faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
825 (
826 *famv.faceFluxCorrectionPtr_
827 );
828 }
829}
830
831
832template<class Type>
833void Foam::faMatrix<Type>::operator+=(const tmp<faMatrix<Type>>& tfamv)
835 operator+=(tfamv());
836 tfamv.clear();
837}
838
839
840template<class Type>
842{
843 checkMethod(*this, famv, "+=");
844
845 dimensions_ -= famv.dimensions_;
847 source_ -= famv.source_;
848 internalCoeffs_ -= famv.internalCoeffs_;
849 boundaryCoeffs_ -= famv.boundaryCoeffs_;
850
851 if (faceFluxCorrectionPtr_ && famv.faceFluxCorrectionPtr_)
852 {
853 *faceFluxCorrectionPtr_ -= *famv.faceFluxCorrectionPtr_;
854 }
855 else if (famv.faceFluxCorrectionPtr_)
856 {
857 faceFluxCorrectionPtr_ = std::make_unique<faceFluxFieldType>
858 (
859 -*famv.faceFluxCorrectionPtr_
860 );
861 }
862}
863
864
865template<class Type>
866void Foam::faMatrix<Type>::operator-=(const tmp<faMatrix<Type>>& tfamv)
868 operator-=(tfamv());
869 tfamv.clear();
870}
871
872
873template<class Type>
875(
877)
879 checkMethod(*this, su, "+=");
880 source() -= su.mesh().S()*su.field();
881}
882
883
884template<class Type>
886(
888)
890 operator+=(tsu());
891 tsu.clear();
892}
893
894
895template<class Type>
897(
899)
901 operator+=(tsu());
902 tsu.clear();
903}
904
905
906template<class Type>
908(
910)
912 checkMethod(*this, su, "-=");
913 source() += su.mesh().S()*su.field();
914}
915
916
917template<class Type>
919(
921)
923 operator-=(tsu());
924 tsu.clear();
925}
926
927
928template<class Type>
930(
932)
934 operator-=(tsu());
935 tsu.clear();
936}
937
938
939template<class Type>
941(
942 const dimensioned<Type>& su
944{
945 source() -= psi().mesh().S()*su;
946}
947
948
949template<class Type>
951(
952 const dimensioned<Type>& su
954{
955 source() += psi().mesh().S()*su;
956}
957
958
959template<class Type>
961(
963)
964{
965 dimensions_ *= dsf.dimensions();
966 lduMatrix::operator*=(dsf.field());
967 source_ *= dsf.field();
968
969 forAll(boundaryCoeffs_, patchi)
970 {
971 const scalarField pisf
972 (
973 dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
974 );
975 internalCoeffs_[patchi] *= pisf;
976 boundaryCoeffs_[patchi] *= pisf;
977 }
978
979 if (faceFluxCorrectionPtr_)
980 {
982 << "cannot scale a matrix containing a faceFluxCorrection"
983 << abort(FatalError);
984 }
985}
986
987
988template<class Type>
990(
991 const tmp<areaScalarField::Internal>& tfld
992)
994 operator*=(tfld());
995 tfld.clear();
996}
997
998
999template<class Type>
1001(
1002 const tmp<areaScalarField>& tfld
1003)
1005 operator*=(tfld());
1006 tfld.clear();
1007}
1008
1009
1010template<class Type>
1012(
1013 const dimensioned<scalar>& ds
1014)
1015{
1016 dimensions_ *= ds.dimensions();
1017 lduMatrix::operator*=(ds.value());
1018 source_ *= ds.value();
1019 internalCoeffs_ *= ds.value();
1020 boundaryCoeffs_ *= ds.value();
1021
1022 if (faceFluxCorrectionPtr_)
1023 {
1024 *faceFluxCorrectionPtr_ *= ds.value();
1026}
1027
1028
1029// * * * * * * * * * * * * * * * Friend Functions * * * * * * * * * * * * * //
1030
1031template<class Type>
1033(
1034 const faMatrix<Type>& mat1,
1035 const faMatrix<Type>& mat2,
1036 const char* op
1037)
1038{
1039 if (&mat1.psi() != &mat2.psi())
1040 {
1042 << "Incompatible fields for operation\n "
1043 << "[" << mat1.psi().name() << "] "
1044 << op
1045 << " [" << mat2.psi().name() << "]"
1046 << abort(FatalError);
1047 }
1048
1049 if
1050 (
1052 && mat1.dimensions() != mat2.dimensions()
1053 )
1054 {
1056 << "Incompatible dimensions for operation\n "
1057 << "[" << mat1.psi().name() << mat1.dimensions()/dimArea << " ] "
1058 << op
1059 << " [" << mat2.psi().name() << mat2.dimensions()/dimArea << " ]"
1060 << abort(FatalError);
1061 }
1062}
1063
1064
1065template<class Type>
1067(
1068 const faMatrix<Type>& mat,
1070 const char* op
1071)
1072{
1073 if
1074 (
1076 && mat.dimensions()/dimArea != fld.dimensions()
1077 )
1078 {
1080 << "Incompatible dimensions for operation\n "
1081 << "[" << mat.psi().name() << mat.dimensions()/dimArea << " ] "
1082 << op
1083 << " [" << fld.name() << fld.dimensions() << " ]"
1084 << abort(FatalError);
1085 }
1086}
1087
1088
1089template<class Type>
1091(
1092 const faMatrix<Type>& mat,
1093 const dimensioned<Type>& dt,
1094 const char* op
1095)
1096{
1097 if
1098 (
1100 && mat.dimensions()/dimArea != dt.dimensions()
1101 )
1102 {
1104 << "Incompatible dimensions for operation\n "
1105 << "[" << mat.psi().name() << mat.dimensions()/dimArea << " ] "
1106 << op
1107 << " [" << dt.name() << dt.dimensions() << " ]"
1108 << abort(FatalError);
1109 }
1110}
1111
1112
1113template<class Type>
1114Foam::SolverPerformance<Type> Foam::solve
1115(
1116 faMatrix<Type>& mat,
1117 const dictionary& solverControls
1118)
1119{
1120 return mat.solve(solverControls);
1121}
1122
1123
1124template<class Type>
1125Foam::SolverPerformance<Type> Foam::solve
1126(
1127 const tmp<faMatrix<Type>>& tmat,
1128 const dictionary& solverControls
1129)
1130{
1131 SolverPerformance<Type> solverPerf(tmat.constCast().solve(solverControls));
1132
1133 tmat.clear();
1134
1135 return solverPerf;
1136}
1137
1138
1139template<class Type>
1140Foam::SolverPerformance<Type> Foam::solve
1141(
1142 faMatrix<Type>& mat,
1143 const word& name
1144)
1145{
1146 return mat.solve(name);
1147}
1148
1149
1150template<class Type>
1151Foam::SolverPerformance<Type> Foam::solve
1152(
1153 const tmp<faMatrix<Type>>& tmat,
1154 const word& name
1155)
1156{
1157 SolverPerformance<Type> solverPerf(tmat.constCast().solve(name));
1158
1159 tmat.clear();
1160
1161 return solverPerf;
1162}
1163
1164
1165template<class Type>
1166Foam::SolverPerformance<Type> Foam::solve(faMatrix<Type>& mat)
1167{
1168 return mat.solve();
1169}
1170
1171
1172template<class Type>
1173Foam::SolverPerformance<Type> Foam::solve(const tmp<faMatrix<Type>>& tmat)
1174{
1175 SolverPerformance<Type> solverPerf(tmat.constCast().solve());
1176
1177 tmat.clear();
1178
1179 return solverPerf;
1180}
1181
1182
1183// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
1184
1185template<class Type>
1186Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1187(
1188 const faMatrix<Type>& A,
1189 const faMatrix<Type>& B
1190)
1191{
1192 checkMethod(A, B, "+");
1193 auto tC = tmp<faMatrix<Type>>::New(A);
1194 tC.ref() += B;
1195 return tC;
1196}
1197
1198
1199template<class Type>
1200Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1201(
1202 const tmp<faMatrix<Type>>& tA,
1203 const faMatrix<Type>& B
1204)
1205{
1206 checkMethod(tA(), B, "+");
1207 tmp<faMatrix<Type>> tC(tA.ptr());
1208 tC.ref() += B;
1209 return tC;
1210}
1211
1212
1213template<class Type>
1214Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1215(
1216 const faMatrix<Type>& A,
1217 const tmp<faMatrix<Type>>& tB
1218)
1219{
1220 checkMethod(A, tB(), "+");
1221 tmp<faMatrix<Type>> tC(tB.ptr());
1222 tC.ref() += A;
1223 return tC;
1224}
1225
1226
1227template<class Type>
1228Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1229(
1230 const tmp<faMatrix<Type>>& tA,
1231 const tmp<faMatrix<Type>>& tB
1232)
1233{
1234 checkMethod(tA(), tB(), "+");
1235 tmp<faMatrix<Type>> tC(tA.ptr());
1236 tC.ref() += tB();
1237 tB.clear();
1238 return tC;
1239}
1240
1241
1242template<class Type>
1243Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1244(
1245 const faMatrix<Type>& A
1246)
1247{
1248 auto tC = tmp<faMatrix<Type>>::New(A);
1249 tC.ref().negate();
1250 return tC;
1251}
1252
1253
1254template<class Type>
1255Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1256(
1257 const tmp<faMatrix<Type>>& tA
1258)
1259{
1260 tmp<faMatrix<Type>> tC(tA.ptr());
1261 tC.ref().negate();
1262 return tC;
1263}
1264
1265
1266template<class Type>
1267Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1268(
1269 const faMatrix<Type>& A,
1270 const faMatrix<Type>& B
1271)
1272{
1273 checkMethod(A, B, "-");
1274 auto tC = tmp<faMatrix<Type>>::New(A);
1275 tC.ref() -= B;
1276 return tC;
1277}
1278
1279
1280template<class Type>
1281Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1282(
1283 const tmp<faMatrix<Type>>& tA,
1284 const faMatrix<Type>& B
1285)
1286{
1287 checkMethod(tA(), B, "-");
1288 tmp<faMatrix<Type>> tC(tA.ptr());
1289 tC.ref() -= B;
1290 return tC;
1291}
1292
1293
1294template<class Type>
1295Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1296(
1297 const faMatrix<Type>& A,
1298 const tmp<faMatrix<Type>>& tB
1299)
1300{
1301 checkMethod(A, tB(), "-");
1302 tmp<faMatrix<Type>> tC(tB.ptr());
1303 tC.ref() -= A;
1304 tC.ref().negate();
1305 return tC;
1306}
1307
1308
1309template<class Type>
1310Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1311(
1312 const tmp<faMatrix<Type>>& tA,
1313 const tmp<faMatrix<Type>>& tB
1314)
1315{
1316 checkMethod(tA(), tB(), "-");
1317 tmp<faMatrix<Type>> tC(tA.ptr());
1318 tC.ref() -= tB();
1319 tB.clear();
1320 return tC;
1321}
1322
1323
1324template<class Type>
1325Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1326(
1327 const faMatrix<Type>& A,
1328 const faMatrix<Type>& B
1329)
1330{
1331 checkMethod(A, B, "==");
1332 return (A - B);
1333}
1334
1335
1336template<class Type>
1337Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1338(
1339 const tmp<faMatrix<Type>>& tA,
1340 const faMatrix<Type>& B
1341)
1342{
1343 checkMethod(tA(), B, "==");
1344 return (tA - B);
1345}
1346
1347
1348template<class Type>
1349Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1350(
1351 const faMatrix<Type>& A,
1352 const tmp<faMatrix<Type>>& tB
1353)
1354{
1355 checkMethod(A, tB(), "==");
1356 return (A - tB);
1357}
1358
1359
1360template<class Type>
1361Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1362(
1363 const tmp<faMatrix<Type>>& tA,
1364 const tmp<faMatrix<Type>>& tB
1365)
1366{
1367 checkMethod(tA(), tB(), "==");
1368 return (tA - tB);
1369}
1370
1371
1372template<class Type>
1373Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1374(
1375 const faMatrix<Type>& A,
1377)
1378{
1379 checkMethod(A, su, "+");
1380 auto tC(A.clone());
1381 tC.ref().source() -= su.mesh().S()*su.field();
1382 return tC;
1383}
1384
1385
1386template<class Type>
1387Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1388(
1389 const faMatrix<Type>& A,
1391)
1392{
1393 checkMethod(A, tsu(), "+");
1394 auto tC(A.clone());
1395 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1396 tsu.clear();
1397 return tC;
1398}
1399
1400
1401template<class Type>
1402Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1403(
1404 const faMatrix<Type>& A,
1406)
1407{
1408 checkMethod(A, tsu(), "+");
1409 auto tC(A.clone());
1410 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1411 tsu.clear();
1412 return tC;
1413}
1414
1415
1416template<class Type>
1417Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1418(
1419 const tmp<faMatrix<Type>>& tA,
1421)
1422{
1423 checkMethod(tA(), su, "+");
1424 tmp<faMatrix<Type>> tC(tA.ptr());
1425 tC.ref().source() -= su.mesh().S()*su.field();
1426 return tC;
1427}
1428
1429
1430template<class Type>
1431Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1432(
1433 const tmp<faMatrix<Type>>& tA,
1435)
1436{
1437 checkMethod(tA(), tsu(), "+");
1438 tmp<faMatrix<Type>> tC(tA.ptr());
1439 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1440 tsu.clear();
1441 return tC;
1442}
1443
1444
1445template<class Type>
1446Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1447(
1448 const tmp<faMatrix<Type>>& tA,
1450)
1451{
1452 checkMethod(tA(), tsu(), "+");
1453 tmp<faMatrix<Type>> tC(tA.ptr());
1454 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1455 tsu.clear();
1456 return tC;
1457}
1458
1459
1460template<class Type>
1461Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1462(
1464 const faMatrix<Type>& A
1465)
1466{
1467 checkMethod(A, su, "+");
1468 auto tC(A.clone());
1469 tC.ref().source() -= su.mesh().S()*su.field();
1470 return tC;
1471}
1472
1473
1474template<class Type>
1475Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1476(
1478 const faMatrix<Type>& A
1479)
1480{
1481 checkMethod(A, tsu(), "+");
1482 auto tC(A.clone());
1483 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1484 tsu.clear();
1485 return tC;
1486}
1487
1488
1489template<class Type>
1490Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1491(
1493 const faMatrix<Type>& A
1494)
1495{
1496 checkMethod(A, tsu(), "+");
1497 auto tC(A.clone());
1498 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1499 tsu.clear();
1500 return tC;
1501}
1502
1503
1504template<class Type>
1505Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1506(
1508 const tmp<faMatrix<Type>>& tA
1509)
1510{
1511 checkMethod(tA(), su, "+");
1512 tmp<faMatrix<Type>> tC(tA.ptr());
1513 tC.ref().source() -= su.mesh().S()*su.field();
1514 return tC;
1515}
1516
1517
1518template<class Type>
1519Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1520(
1522 const tmp<faMatrix<Type>>& tA
1523)
1524{
1525 checkMethod(tA(), tsu(), "+");
1526 tmp<faMatrix<Type>> tC(tA.ptr());
1527 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1528 tsu.clear();
1529 return tC;
1530}
1531
1532
1533template<class Type>
1534Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1535(
1537 const tmp<faMatrix<Type>>& tA
1538)
1539{
1540 checkMethod(tA(), tsu(), "+");
1541 tmp<faMatrix<Type>> tC(tA.ptr());
1542 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1543 tsu.clear();
1544 return tC;
1545}
1546
1547
1548template<class Type>
1549Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1550(
1551 const faMatrix<Type>& A,
1553)
1554{
1555 checkMethod(A, su, "-");
1556 auto tC(A.clone());
1557 tC.ref().source() += su.mesh().S()*su.field();
1558 return tC;
1559}
1560
1561
1562template<class Type>
1563Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1564(
1565 const faMatrix<Type>& A,
1567)
1568{
1569 checkMethod(A, tsu(), "-");
1570 auto tC(A.clone());
1571 tC.ref().source() += tsu().mesh().S()*tsu().field();
1572 tsu.clear();
1573 return tC;
1574}
1575
1576
1577template<class Type>
1578Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1579(
1580 const faMatrix<Type>& A,
1582)
1583{
1584 checkMethod(A, tsu(), "-");
1585 auto tC(A.clone());
1586 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1587 tsu.clear();
1588 return tC;
1589}
1590
1591
1592template<class Type>
1593Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1594(
1595 const tmp<faMatrix<Type>>& tA,
1597)
1598{
1599 checkMethod(tA(), su, "-");
1600 tmp<faMatrix<Type>> tC(tA.ptr());
1601 tC.ref().source() += su.mesh().S()*su.field();
1602 return tC;
1603}
1604
1605
1606template<class Type>
1607Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1608(
1609 const tmp<faMatrix<Type>>& tA,
1611)
1612{
1613 checkMethod(tA(), tsu(), "-");
1614 tmp<faMatrix<Type>> tC(tA.ptr());
1615 tC.ref().source() += tsu().mesh().S()*tsu().field();
1616 tsu.clear();
1617 return tC;
1618}
1619
1620
1621template<class Type>
1622Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1623(
1624 const tmp<faMatrix<Type>>& tA,
1626)
1627{
1628 checkMethod(tA(), tsu(), "-");
1629 tmp<faMatrix<Type>> tC(tA.ptr());
1630 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1631 tsu.clear();
1632 return tC;
1633}
1634
1635
1636template<class Type>
1637Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1638(
1640 const faMatrix<Type>& A
1641)
1642{
1643 checkMethod(A, su, "-");
1644 auto tC(A.clone());
1645 tC.ref().negate();
1646 tC.ref().source() -= su.mesh().S()*su.field();
1647 return tC;
1648}
1649
1650
1651template<class Type>
1652Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1653(
1655 const faMatrix<Type>& A
1656)
1657{
1658 checkMethod(A, tsu(), "-");
1659 auto tC(A.clone());
1660 tC.ref().negate();
1661 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1662 tsu.clear();
1663 return tC;
1664}
1665
1666
1667template<class Type>
1668Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1669(
1671 const faMatrix<Type>& A
1672)
1673{
1674 checkMethod(A, tsu(), "-");
1675 auto tC(A.clone());
1676 tC.ref().negate();
1677 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1678 tsu.clear();
1679 return tC;
1680}
1681
1682
1683template<class Type>
1684Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1685(
1687 const tmp<faMatrix<Type>>& tA
1688)
1689{
1690 checkMethod(tA(), su, "-");
1691 tmp<faMatrix<Type>> tC(tA.ptr());
1692 tC.ref().negate();
1693 tC.ref().source() -= su.mesh().S()*su.field();
1694 return tC;
1695}
1696
1697
1698template<class Type>
1699Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1700(
1702 const tmp<faMatrix<Type>>& tA
1703)
1704{
1705 checkMethod(tA(), tsu(), "-");
1706 tmp<faMatrix<Type>> tC(tA.ptr());
1707 tC.ref().negate();
1708 tC.ref().source() -= tsu().mesh().S()*tsu().field();
1709 tsu.clear();
1710 return tC;
1711}
1712
1713
1714template<class Type>
1715Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1716(
1718 const tmp<faMatrix<Type>>& tA
1719)
1720{
1721 checkMethod(tA(), tsu(), "-");
1722 tmp<faMatrix<Type>> tC(tA.ptr());
1723 tC.ref().negate();
1724 tC.ref().source() -= tsu().mesh().S()*tsu().primitiveField();
1725 tsu.clear();
1726 return tC;
1727}
1728
1729
1730template<class Type>
1731Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1732(
1733 const faMatrix<Type>& A,
1734 const dimensioned<Type>& su
1735)
1736{
1737 checkMethod(A, su, "+");
1738 auto tC(A.clone());
1739 tC.ref().source() -= su.value()*A.psi().mesh().S();
1740 return tC;
1741}
1742
1743
1744template<class Type>
1745Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1746(
1747 const tmp<faMatrix<Type>>& tA,
1748 const dimensioned<Type>& su
1749)
1750{
1751 checkMethod(tA(), su, "+");
1752 tmp<faMatrix<Type>> tC(tA.ptr());
1753 tC.ref().source() -= su.value()*tC().psi().mesh().S();
1754 return tC;
1755}
1756
1757
1758template<class Type>
1759Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1760(
1761 const dimensioned<Type>& su,
1762 const faMatrix<Type>& A
1763)
1764{
1765 checkMethod(A, su, "+");
1766 auto tC(A.clone());
1767 tC.ref().source() -= su.value()*A.psi().mesh().S();
1768 return tC;
1769}
1770
1771
1772template<class Type>
1773Foam::tmp<Foam::faMatrix<Type>> Foam::operator+
1774(
1775 const dimensioned<Type>& su,
1776 const tmp<faMatrix<Type>>& tA
1777)
1778{
1779 checkMethod(tA(), su, "+");
1780 tmp<faMatrix<Type>> tC(tA.ptr());
1781 tC.ref().source() -= su.value()*tC().psi().mesh().S();
1782 return tC;
1783}
1784
1785
1786template<class Type>
1787Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1788(
1789 const faMatrix<Type>& A,
1790 const dimensioned<Type>& su
1791)
1792{
1793 checkMethod(A, su, "-");
1794 auto tC(A.clone());
1795 tC.ref().source() += su.value()*tC().psi().mesh().S();
1796 return tC;
1797}
1798
1799
1800template<class Type>
1801Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1802(
1803 const tmp<faMatrix<Type>>& tA,
1804 const dimensioned<Type>& su
1805)
1806{
1807 checkMethod(tA(), su, "-");
1808 tmp<faMatrix<Type>> tC(tA.ptr());
1809 tC.ref().source() += su.value()*tC().psi().mesh().S();
1810 return tC;
1811}
1812
1813
1814template<class Type>
1815Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1816(
1817 const dimensioned<Type>& su,
1818 const faMatrix<Type>& A
1819)
1820{
1821 checkMethod(A, su, "-");
1822 auto tC(A.clone());
1823 tC.ref().negate();
1824 tC.ref().source() -= su.value()*A.psi().mesh().S();
1825 return tC;
1826}
1827
1828
1829template<class Type>
1830Foam::tmp<Foam::faMatrix<Type>> Foam::operator-
1831(
1832 const dimensioned<Type>& su,
1833 const tmp<faMatrix<Type>>& tA
1834)
1835{
1836 checkMethod(tA(), su, "-");
1837 tmp<faMatrix<Type>> tC(tA.ptr());
1838 tC.ref().negate();
1839 tC.ref().source() -= su.value()*tC().psi().mesh().S();
1840 return tC;
1841}
1842
1843
1844template<class Type>
1845Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1846(
1847 const faMatrix<Type>& A,
1849)
1850{
1851 checkMethod(A, su, "==");
1852 auto tC(A.clone());
1853 tC.ref().source() += su.mesh().S()*su.field();
1854 return tC;
1855}
1856
1857
1858template<class Type>
1859Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1860(
1861 const faMatrix<Type>& A,
1863)
1864{
1865 checkMethod(A, tsu(), "==");
1866 auto tC(A.clone());
1867 tC.ref().source() += tsu().mesh().S()*tsu().field();
1868 tsu.clear();
1869 return tC;
1870}
1871
1872
1873template<class Type>
1874Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1875(
1876 const faMatrix<Type>& A,
1878)
1879{
1880 checkMethod(A, tsu(), "==");
1881 auto tC(A.clone());
1882 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1883 tsu.clear();
1884 return tC;
1885}
1886
1887
1888template<class Type>
1889Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1890(
1891 const tmp<faMatrix<Type>>& tA,
1893)
1894{
1895 checkMethod(tA(), su, "==");
1896 tmp<faMatrix<Type>> tC(tA.ptr());
1897 tC.ref().source() += su.mesh().S()*su.field();
1898 return tC;
1899}
1900
1901
1902template<class Type>
1903Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1904(
1905 const tmp<faMatrix<Type>>& tA,
1907)
1908{
1909 checkMethod(tA(), tsu(), "==");
1910 tmp<faMatrix<Type>> tC(tA.ptr());
1911 tC.ref().source() += tsu().mesh().S()*tsu().field();
1912 tsu.clear();
1913 return tC;
1914}
1915
1916
1917template<class Type>
1918Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1919(
1920 const tmp<faMatrix<Type>>& tA,
1922)
1923{
1924 checkMethod(tA(), tsu(), "==");
1925 tmp<faMatrix<Type>> tC(tA.ptr());
1926 tC.ref().source() += tsu().mesh().S()*tsu().primitiveField();
1927 tsu.clear();
1928 return tC;
1929}
1930
1931
1932template<class Type>
1933Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1934(
1935 const faMatrix<Type>& A,
1936 const dimensioned<Type>& su
1937)
1938{
1939 checkMethod(A, su, "==");
1940 auto tC(A.clone());
1941 tC.ref().source() += A.psi().mesh().S()*su.value();
1942 return tC;
1943}
1944
1945
1946template<class Type>
1947Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1948(
1949 const tmp<faMatrix<Type>>& tA,
1950 const dimensioned<Type>& su
1951)
1952{
1953 checkMethod(tA(), su, "==");
1954 tmp<faMatrix<Type>> tC(tA.ptr());
1955 tC.ref().source() += tC().psi().mesh().S()*su.value();
1956 return tC;
1957}
1958
1959
1960template<class Type>
1961Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1962(
1963 const faMatrix<Type>& A,
1964 const Foam::zero
1965)
1966{
1967 return A;
1968}
1969
1970
1971template<class Type>
1972Foam::tmp<Foam::faMatrix<Type>> Foam::operator==
1973(
1974 const tmp<faMatrix<Type>>& tA,
1975 const Foam::zero
1976)
1977{
1978 return tA;
1979}
1980
1981
1982template<class Type>
1983Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1984(
1985 const areaScalarField::Internal& dsf,
1986 const faMatrix<Type>& A
1987)
1988{
1989 auto tC(A.clone());
1990 tC.ref() *= dsf;
1991 return tC;
1992}
1993
1994
1995template<class Type>
1996Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
1997(
1999 const faMatrix<Type>& A
2000)
2001{
2002 auto tC(A.clone());
2003 tC.ref() *= tdsf;
2004 return tC;
2005}
2006
2007
2008template<class Type>
2009Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2010(
2011 const tmp<areaScalarField>& tvsf,
2012 const faMatrix<Type>& A
2013)
2014{
2015 auto tC(A.clone());
2016 tC.ref() *= tvsf;
2017 return tC;
2018}
2019
2020
2021template<class Type>
2022Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2023(
2024 const areaScalarField::Internal& dsf,
2025 const tmp<faMatrix<Type>>& tA
2026)
2027{
2028 tmp<faMatrix<Type>> tC(tA.ptr());
2029 tC.ref() *= dsf;
2030 return tC;
2031}
2032
2033
2034template<class Type>
2035Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2036(
2038 const tmp<faMatrix<Type>>& tA
2039)
2040{
2041 tmp<faMatrix<Type>> tC(tA.ptr());
2042 tC.ref() *= tdsf;
2043 return tC;
2044}
2045
2046
2047template<class Type>
2048Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2049(
2050 const tmp<areaScalarField>& tvsf,
2051 const tmp<faMatrix<Type>>& tA
2052)
2053{
2054 tmp<faMatrix<Type>> tC(tA.ptr());
2055 tC.ref() *= tvsf;
2056 return tC;
2057}
2058
2059
2060template<class Type>
2061Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2062(
2063 const dimensioned<scalar>& ds,
2064 const faMatrix<Type>& A
2065)
2066{
2067 auto tC(A.clone());
2068 tC.ref() *= ds;
2069 return tC;
2070}
2071
2072
2073template<class Type>
2074Foam::tmp<Foam::faMatrix<Type>> Foam::operator*
2075(
2076 const dimensioned<scalar>& ds,
2077 const tmp<faMatrix<Type>>& tA
2078)
2079{
2080 tmp<faMatrix<Type>> tC(tA.ptr());
2081 tC.ref() *= ds;
2082 return tC;
2083}
2084
2085
2086template<class Type>
2087Foam::tmp<Foam::GeometricField<Type, Foam::faPatchField, Foam::areaMesh>>
2088Foam::operator&
2089(
2090 const faMatrix<Type>& M,
2092)
2093{
2095 (
2096 "M&" + psi.name(),
2097 psi.mesh(),
2098 M.dimensions()/dimArea,
2100 );
2101 auto& Mphi = tMphi.ref();
2102
2103 // Loop over field components
2104 if (M.hasDiag())
2105 {
2106 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2107 {
2108 scalarField psiCmpt(psi.field().component(cmpt));
2109 scalarField boundaryDiagCmpt(M.diag());
2110 M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2111 Mphi.primitiveFieldRef().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2112 }
2113 }
2114 else
2115 {
2116 Mphi.primitiveFieldRef() = Zero;
2117 }
2118
2119 Mphi.primitiveFieldRef() += M.lduMatrix::H(psi.field()) + M.source();
2120 M.addBoundarySource(Mphi.primitiveFieldRef());
2121
2122 Mphi.primitiveFieldRef() /= -psi.mesh().S();
2123 Mphi.correctBoundaryConditions();
2124
2125 return tMphi;
2126}
2127
2128
2129template<class Type>
2130Foam::tmp<Foam::GeometricField<Type, Foam::faPatchField, Foam::areaMesh>>
2131Foam::operator&
2132(
2133 const faMatrix<Type>& M,
2135)
2136{
2138 tpsi.clear();
2139 return tMpsi;
2140}
2141
2142
2143template<class Type>
2144Foam::tmp<Foam::GeometricField<Type, Foam::faPatchField, Foam::areaMesh>>
2145Foam::operator&
2146(
2147 const faMatrix<Type>& M,
2149)
2150{
2152 tpsi.clear();
2153 return tMpsi;
2154}
2155
2156
2157template<class Type>
2158Foam::tmp<Foam::GeometricField<Type, Foam::faPatchField, Foam::areaMesh>>
2159Foam::operator&
2160(
2161 const tmp<faMatrix<Type>>& tM,
2163)
2164{
2166 tM.clear();
2167 return tMpsi;
2168}
2169
2170
2171template<class Type>
2172Foam::tmp<Foam::GeometricField<Type, Foam::faPatchField, Foam::areaMesh>>
2173Foam::operator&
2174(
2175 const tmp<faMatrix<Type>>& tM,
2177)
2178{
2180 tM.clear();
2181 tpsi.clear();
2182 return tMpsi;
2183}
2184
2185
2186template<class Type>
2187Foam::tmp<Foam::GeometricField<Type, Foam::faPatchField, Foam::areaMesh>>
2188Foam::operator&
2189(
2190 const tmp<faMatrix<Type>>& tM,
2192)
2193{
2195 tM.clear();
2196 tpsi.clear();
2197 return tMpsi;
2198}
2199
2200
2201// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2202
2203template<class Type>
2204Foam::Ostream& Foam::operator<<(Ostream& os, const faMatrix<Type>& fam)
2205{
2207 << fam.dimensions_ << nl
2208 << fam.source_ << nl
2209 << fam.internalCoeffs_ << nl
2210 << fam.boundaryCoeffs_ << endl;
2211
2213
2214 return os;
2215}
2216
2217
2218// * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
2219
2220#include "faMatrixSolve.C"
2221
2222// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
#define M(I)
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
labelList faceLabels(nFaceLabels)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DynamicField< Type > & field() const noexcept
Return const-reference to the primitive field values.
const Mesh & mesh() const noexcept
Return const reference to mesh.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void negate()
Inplace negate this field (negative).
Definition Field.C:592
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
Generic GeometricField class.
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
DimensionedField< scalar, areaMesh > Internal
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition IOstream.C:45
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
A List with indirect addressing. Like IndirectList but does not store addressing.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
A single value that is represented as a list with an operator[] to access the value....
Definition UniformList.H:49
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
static bool checking(bool on) noexcept
Turn dimension checking on/off.
Generic dimensioned Type class.
const dimensionSet & dimensions() const noexcept
Return const reference to dimensions.
const word & name() const noexcept
Return const reference to name.
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition faMatrix.H:108
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition faMatrix.C:33
void operator=(const faMatrix< Type > &)
Definition faMatrix.C:744
void relax()
Relax matrix (for steady-state solution).
Definition faMatrix.C:574
tmp< GeometricField< Type, faePatchField, edgeMesh > > flux() const
Return the face-flux field from the matrix.
Definition faMatrix.C:655
void setValues(const labelUList &faceLabels, const Type &value)
Set solution in given faces to the specified value and eliminate the corresponding equations from the...
Definition faMatrix.C:387
void setReference(const label facei, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition faMatrix.C:420
void setReferences(const labelUList &faceLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition faMatrix.C:439
void operator+=(const faMatrix< Type > &)
Definition faMatrix.C:801
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition faMatrix.H:348
const dimensionSet & dimensions() const noexcept
Definition faMatrix.H:353
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition faMatrix.C:619
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition faMatrix.C:126
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition faMatrix.C:142
virtual ~faMatrix()
Destructor.
Definition faMatrix.C:276
tmp< areaScalarField > A() const
Return the central coefficient.
Definition faMatrix.C:600
void operator-=(const faMatrix< Type > &)
Definition faMatrix.C:834
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition faMatrix.C:108
void negate()
Inplace negate.
Definition faMatrix.C:786
faMatrix(const GeometricField< Type, faPatchField, areaMesh > &psi, const dimensionSet &ds)
Construct given a field to solve for.
Definition faMatrix.C:176
void setValuesFromList(const labelUList &faceLabels, const ListType< Type > &values)
Set solution in given faces to the specified values.
Definition faMatrix.C:288
tmp< scalarField > D() const
Return the matrix diagonal.
Definition faMatrix.C:591
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition faMatrix.C:71
Field< Type > & source() noexcept
Definition faMatrix.H:358
const dictionary & solverDict() const
Return the solver dictionary for psi.
Definition faMatrix.C:732
void operator*=(const areaScalarField::Internal &)
Definition faMatrix.C:954
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition faMesh.H:140
virtual bool coupled() const
True if the patch field is coupled.
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
faPatchField<Type> abstract base class. This class gives a fat-interface to all derived classes cover...
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
lduMatrix(const lduMesh &mesh)
Construct (without coefficients) for an LDU addressed mesh.
Definition lduMatrix.C:54
void operator=(const lduMatrix &)
Copy assignment.
tmp< Field< Type > > faceH(const Field< Type > &) const
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition lduMatrix.H:769
bool asymmetric() const noexcept
Matrix is asymmetric (ie, full).
Definition lduMatrix.H:850
const scalarField & diag() const
Definition lduMatrix.C:195
bool symmetric() const noexcept
Matrix is symmetric.
Definition lduMatrix.H:842
const scalarField & upper() const
Definition lduMatrix.C:235
void operator*=(const scalarField &)
tmp< Field< Type > > H(const Field< Type > &) const
void operator+=(const lduMatrix &)
const scalarField & lower() const
Definition lduMatrix.C:306
void sumMagOffDiag(scalarField &sumOff) const
const lduMesh & mesh() const noexcept
Return the LDU mesh from which the addressing is obtained.
Definition lduMatrix.H:753
void operator-=(const lduMatrix &)
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T * ptr() const
Return managed pointer for reuse, or clone() the object reference.
Definition tmpI.H:256
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & psi
UEqn relax()
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Finite-Area matrix basic solvers.
OBJstream os(runTime.globalPath()/outputName)
auto & name
label faceId(-1)
#define FUNCTION_NAME
#define DebugInFunction
Report an information message using Foam::Info.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition faMatrix.C:1026
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
const dimensionSet dimArea(sqr(dimLength))
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).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
tmp< DimensionedField< typename DimensionedField< Type, GeoMesh >::cmptType, GeoMesh > > cmptAv(const DimensionedField< Type, GeoMesh > &f1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
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:299