fvMatrix.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-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 "volFields.H"
30#include "surfaceFields.H"
34#include "IndirectList.H"
35#include "UniformList.H"
36#include "demandDrivenData.H"
37
38#include "cyclicFvPatchField.H"
41
43
44// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45
46template<class Type>
47template<class Type2>
49(
50 const labelUList& addr,
51 const Field<Type2>& pf,
52 Field<Type2>& intf
53) const
54{
55 if (addr.size() != pf.size())
56 {
58 << "addressing (" << addr.size()
59 << ") and field (" << pf.size() << ") are different sizes" << endl
60 << abort(FatalError);
61 }
62
63 forAll(addr, facei)
64 {
65 intf[addr[facei]] += pf[facei];
66 }
67}
68
69
70template<class Type>
71template<class Type2>
73(
74 const labelUList& addr,
75 const tmp<Field<Type2>>& tpf,
76 Field<Type2>& intf
77) const
78{
79 addToInternalField(addr, tpf(), intf);
80 tpf.clear();
81}
82
83
84template<class Type>
85template<class Type2>
87(
88 const labelUList& addr,
89 const Field<Type2>& pf,
90 Field<Type2>& intf
91) const
92{
93 if (addr.size() != pf.size())
94 {
96 << "addressing (" << addr.size()
97 << ") and field (" << pf.size() << ") are different sizes" << endl
98 << abort(FatalError);
99 }
100
101 forAll(addr, facei)
102 {
103 intf[addr[facei]] -= pf[facei];
104 }
105}
106
107
108template<class Type>
109template<class Type2>
111(
112 const labelUList& addr,
113 const tmp<Field<Type2>>& tpf,
114 Field<Type2>& intf
115) const
116{
117 subtractFromInternalField(addr, tpf(), intf);
118 tpf.clear();
119}
120
121
122template<class Type>
124(
125 scalarField& diag,
126 const direction solveCmpt
127) const
128{
129 for (label fieldi = 0; fieldi < nMatrices(); ++fieldi)
130 {
131 const auto& bpsi = this->psi(fieldi).boundaryField();
132
133 forAll(bpsi, ptfi)
134 {
135 const label patchi = globalPatchID(fieldi, ptfi);
136
137 if (patchi != -1)
138 {
139 addToInternalField
140 (
141 lduAddr().patchAddr(patchi),
142 internalCoeffs_[patchi].component(solveCmpt),
143 diag
144 );
145 }
146 }
147 }
148}
149
150
151template<class Type>
153{
154 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
155 {
156 const auto& bpsi = this->psi(fieldi).boundaryField();
157
158 forAll(bpsi, ptfi)
159 {
160 const label patchi = globalPatchID(fieldi, ptfi);
161 if (patchi != -1)
162 {
163 addToInternalField
164 (
165 lduAddr().patchAddr(patchi),
166 cmptAv(internalCoeffs_[patchi]),
167 diag
168 );
169 }
170 }
171 }
172}
173
174
175template<class Type>
177(
178 Field<Type>& source,
179 const bool couples
180) const
181{
182 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
183 {
184 const auto& bpsi = this->psi(fieldi).boundaryField();
185
186 forAll(bpsi, ptfi)
187 {
188 const fvPatchField<Type>& ptf = bpsi[ptfi];
189
190 const label patchi = globalPatchID(fieldi, ptfi);
191
192 if (patchi != -1)
193 {
194 const Field<Type>& pbc = boundaryCoeffs_[patchi];
195
196 if (!ptf.coupled())
197 {
198 addToInternalField
199 (
200 lduAddr().patchAddr(patchi),
201 pbc,
202 source
203 );
204 }
205 else if (couples)
206 {
207 const tmp<Field<Type>> tpnf = ptf.patchNeighbourField();
208 const Field<Type>& pnf = tpnf();
209
210 const labelUList& addr = lduAddr().patchAddr(patchi);
211
212 forAll(addr, facei)
213 {
214 source[addr[facei]] +=
215 cmptMultiply(pbc[facei], pnf[facei]);
216 }
217 }
218 }
219 }
220 }
221}
222
223
224template<class Type>
225template<template<class> class ListType>
227(
228 const labelUList& cellLabels,
229 const ListType<Type>& values
230)
231{
232 const fvMesh& mesh = psi_.mesh();
233
234 const cellList& cells = mesh.cells();
235 const labelUList& own = mesh.owner();
236 const labelUList& nei = mesh.neighbour();
237
238 scalarField& Diag = diag();
240 const_cast
241 <
243 >(psi_).primitiveFieldRef();
244
245 forAll(cellLabels, i)
246 {
247 const label celli = cellLabels[i];
248 const Type& value = values[i];
249
250 psi[celli] = value;
251 source_[celli] = value*Diag[celli];
252
253 if (symmetric() || asymmetric())
254 {
255 for (const label facei : cells[celli])
256 {
257 if (mesh.isInternalFace(facei))
258 {
259 if (symmetric())
260 {
261 if (celli == own[facei])
262 {
263 source_[nei[facei]] -= upper()[facei]*value;
264 }
265 else
266 {
267 source_[own[facei]] -= upper()[facei]*value;
268 }
269
270 upper()[facei] = 0.0;
271 }
272 else
273 {
274 if (celli == own[facei])
275 {
276 source_[nei[facei]] -= lower()[facei]*value;
277 }
278 else
279 {
280 source_[own[facei]] -= upper()[facei]*value;
281 }
282
283 upper()[facei] = 0.0;
284 lower()[facei] = 0.0;
285 }
286 }
287 else
288 {
289 const label patchi = mesh.boundaryMesh().whichPatch(facei);
290
291 if (internalCoeffs_[patchi].size())
292 {
293 const label patchFacei =
294 mesh.boundaryMesh()[patchi].whichFace(facei);
295
296 internalCoeffs_[patchi][patchFacei] = Zero;
297 boundaryCoeffs_[patchi][patchFacei] = Zero;
298 }
299 }
300 }
301 }
302 }
303}
304
305
306template<class Type>
308{
309 const auto& bpsi = this->psi(fieldi).boundaryField();
310
311 word idName;
312 forAll(bpsi, patchi)
313 {
314 if (bpsi[patchi].useImplicit())
315 {
316 if (debug)
317 {
318 Pout<< "fvMatrix<Type>::checkImplicit "
319 << " field:" << this->psi(fieldi).name()
320 << " on mesh:"
321 << this->psi(fieldi).mesh().name()
322 << " patch:" << bpsi[patchi].patch().name()
323 << endl;
324 }
325
326 idName += Foam::name(patchi);
327 useImplicit_ = true;
328 }
329 }
330
331 if (useImplicit_)
332 {
333 lduAssemblyName_ = word("lduAssembly") + idName;
334 }
335
336 return !idName.empty();
337}
338
339
340// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341
342template<class Type>
344(
346 const dimensionSet& ds
347)
348:
349 lduMatrix(psi.mesh()),
350 psi_(psi),
351 useImplicit_(false),
352 lduAssemblyName_(),
353 nMatrix_(0),
354 dimensions_(ds),
355 source_(psi.size(), Zero),
356 internalCoeffs_(psi.mesh().boundary().size()),
357 boundaryCoeffs_(psi.mesh().boundary().size()),
358 faceFluxCorrectionPtr_(nullptr)
359{
361 << "Constructing fvMatrix<Type> for field " << psi_.name() << endl;
362
364
365 forAll(psi.mesh().boundary(), patchi)
366 {
367 internalCoeffs_.set
368 (
369 patchi,
370 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
371 );
372
373 boundaryCoeffs_.set
374 (
375 patchi,
376 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
377 );
378 }
379
380 auto& psiRef = this->psi(0);
381 const label currentStatePsi = psiRef.eventNo();
382 psiRef.boundaryFieldRef().updateCoeffs();
383 psiRef.eventNo() = currentStatePsi;
384}
385
386
387template<class Type>
389:
390 lduMatrix(fvm),
391 psi_(fvm.psi_),
392 useImplicit_(fvm.useImplicit_),
393 lduAssemblyName_(fvm.lduAssemblyName_),
394 nMatrix_(fvm.nMatrix_),
395 dimensions_(fvm.dimensions_),
396 source_(fvm.source_),
397 internalCoeffs_(fvm.internalCoeffs_),
398 boundaryCoeffs_(fvm.boundaryCoeffs_),
399 faceFluxCorrectionPtr_(nullptr)
400{
402 << "Copying fvMatrix<Type> for field " << psi_.name() << endl;
403
404 if (fvm.faceFluxCorrectionPtr_)
405 {
406 faceFluxCorrectionPtr_ =
408 (
409 *(fvm.faceFluxCorrectionPtr_)
410 );
411 }
412}
413
414
415template<class Type>
417:
418 lduMatrix(tmat.constCast(), tmat.movable()),
419 psi_(tmat().psi_),
420 useImplicit_(tmat().useImplicit_),
421 lduAssemblyName_(tmat().lduAssemblyName_),
422 nMatrix_(tmat().nMatrix_),
423 dimensions_(tmat().dimensions_),
424 source_(tmat.constCast().source_, tmat.movable()),
425 internalCoeffs_(tmat.constCast().internalCoeffs_, tmat.movable()),
426 boundaryCoeffs_(tmat.constCast().boundaryCoeffs_, tmat.movable()),
427 faceFluxCorrectionPtr_(nullptr)
428{
430 << "Copy/move fvMatrix<Type> for field " << psi_.name() << endl;
431
432 if (tmat().faceFluxCorrectionPtr_)
433 {
434 if (tmat.movable())
435 {
436 faceFluxCorrectionPtr_ = tmat().faceFluxCorrectionPtr_;
437 tmat().faceFluxCorrectionPtr_ = nullptr;
438 }
439 else
440 {
441 faceFluxCorrectionPtr_ =
443 (
444 *(tmat().faceFluxCorrectionPtr_)
445 );
446 }
447 }
448
449 tmat.clear();
450}
451
452
453// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * * //
454
455template<class Type>
457{
459 << "Destroying fvMatrix<Type> for field " << psi_.name() << endl;
460
461 deleteDemandDrivenData(faceFluxCorrectionPtr_);
462 subMatrices_.clear();
463}
464
465
466// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
467
468template<class Type>
470(
471 lduInterfaceFieldPtrsList& interfaces,
472 PtrDynList<lduInterfaceField>& newInterfaces
473)
474{
475 interfaces.setSize(internalCoeffs_.size());
476 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
477 {
478 const auto& bpsi = this->psi(fieldi).boundaryField();
479 lduInterfaceFieldPtrsList fieldInterfaces(bpsi.scalarInterfaces());
480
481 forAll (fieldInterfaces, patchi)
482 {
483 label globalPatchID = lduMeshPtr()->patchMap()[fieldi][patchi];
484
485 if (globalPatchID != -1)
486 {
487 if (fieldInterfaces.set(patchi))
488 {
489 if (isA<cyclicLduInterfaceField>(bpsi[patchi]))
490 {
491 newInterfaces.append
492 (
494 (
495 refCast<const fvPatch>
496 (
497 lduMeshPtr()->interfaces()[globalPatchID]
498 ),
499 bpsi[patchi].internalField()
500 )
501 );
502 interfaces.set(globalPatchID, &newInterfaces.last());
503
504 }
505 else if (isA<cyclicAMILduInterfaceField>(bpsi[patchi]))
506 {
507 newInterfaces.append
508 (
510 (
511 refCast<const fvPatch>
512 (
513 lduMeshPtr()->interfaces()[globalPatchID]
514 ),
515 bpsi[patchi].internalField()
516 )
517 );
518 interfaces.set(globalPatchID, &newInterfaces.last());
519 }
520 else if (isA<cyclicACMILduInterfaceField>(bpsi[patchi]))
521 {
522 newInterfaces.append
523 (
525 (
526 refCast<const fvPatch>
527 (
528 lduMeshPtr()->interfaces()[globalPatchID]
529 ),
530 bpsi[patchi].internalField()
531 )
532 );
533 interfaces.set(globalPatchID, &newInterfaces.last());
534 }
535 else
536 {
537 interfaces.set(globalPatchID, &fieldInterfaces[patchi]);
538 }
539 }
540 }
541 }
542 }
543}
544
545
546template<class Type>
548(
549 label fieldi,
550 const FieldField<Field, Type>& fluxContrib,
552 bool internal
553) const
554{
555 const lduPrimitiveMeshAssembly* ptr = lduMeshPtr();
556
557 const labelList& patchMap = ptr->patchMap()[fieldi];
558
559 forAll(contrib, patchi)
560 {
561 const label globalPtchId = patchMap[patchi];
562
563 if (globalPtchId != -1)
564 {
565 // Cache contrib before overwriting
566 const Field<Type> saveContrib(fluxContrib[globalPtchId]);
567 contrib[patchi].setSize(psi_.boundaryField()[patchi].size()),
568 contrib[patchi] = pTraits<Type>::zero;
569
570 if (internal)
571 {
572 contrib[patchi] =
574 (
575 saveContrib,
576 psi_.boundaryField()[patchi].patchInternalField()
577 );
578 }
579 else
580 {
581 if (this->psi(fieldi).boundaryField()[patchi].coupled())
582 {
583 contrib[patchi] =
585 (
586 saveContrib,
587 psi_.boundaryField()[patchi].patchNeighbourField()
588 );
589 }
590 }
591 }
592 else if (globalPtchId == -1)
593 {
594 const polyPatch& pp =
595 this->psi(fieldi).mesh().boundaryMesh()[patchi];
596
597 if (pp.masterImplicit())
598 {
599 label virtualPatch =
600 ptr->patchLocalToGlobalMap()[fieldi][patchi];
601
602 const label nbrPatchId = pp.neighbPolyPatchID();
603
604 // Copy contrib before overwriting
605 const Field<Type> saveContrib(fluxContrib[virtualPatch]);
606
607 Field<Type>& coeffs = contrib[patchi];
608 Field<Type>& nbrCoeffs = contrib[nbrPatchId];
609
610 coeffs.setSize(psi_.boundaryField()[patchi].size());
611 nbrCoeffs.setSize(psi_.boundaryField()[nbrPatchId].size());
612
613 coeffs = pTraits<Type>::zero;
614 nbrCoeffs = pTraits<Type>::zero;
615
616 // nrb cells
617 const labelList& nbrCellIds =
618 ptr->cellBoundMap()[fieldi][patchi];
619
620 const labelList& cellIds =
621 ptr->cellBoundMap()[fieldi][nbrPatchId];
622
624 this->psi(fieldi);
625
626 forAll(saveContrib, subFaceI)
627 {
628 const label faceId =
629 ptr->facePatchFaceMap()[fieldi][patchi][subFaceI];
630 const label nbrFaceId =
631 ptr->facePatchFaceMap()[fieldi][nbrPatchId][subFaceI];
632
633 const label nbrCellId = nbrCellIds[subFaceI];
634 const label cellId = cellIds[subFaceI];
635
636 if (internal)
637 {
638 coeffs[faceId] +=
639 cmptMultiply(saveContrib[subFaceI], psi[cellId]);
640
641 nbrCoeffs[nbrFaceId] +=
642 cmptMultiply(saveContrib[subFaceI], psi[nbrCellId]);
643 }
644 else //boundary
645 {
646 coeffs[faceId] +=
647 cmptMultiply(saveContrib[subFaceI], psi[nbrCellId]);
648
649 nbrCoeffs[nbrFaceId] +=
650 cmptMultiply(saveContrib[subFaceI], psi[cellId]);
651 }
652 }
653 }
654 }
655 }
656}
657
658
659template<class Type>
661{
662 // If it is a multi-fvMatrix needs correct internalCoeffs and
663 // boundaryCoeffs size
664 if (nMatrix_ > 0)
665 {
666 label interfaceI(0);
667 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
668 {
669 const auto& psi = this->psi(fieldi);
670
671 forAll(psi.mesh().boundary(), patchi)
672 {
673 interfaceI++;
674 }
675 }
676 internalCoeffs_.setSize(interfaceI);
677 boundaryCoeffs_.setSize(interfaceI);
678
679 interfaceI = 0;
680 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
681 {
682 const auto& psi = this->psi(fieldi);
683
684 forAll(psi.mesh().boundary(), patchi)
685 {
686 internalCoeffs_.set
687 (
688 interfaceI,
689 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
690 );
691
692 boundaryCoeffs_.set
693 (
694 interfaceI,
695 new Field<Type>(psi.mesh().boundary()[patchi].size(), Zero)
696 );
697 interfaceI++;
698 }
699 }
700 }
701
702 for (label i=0; i < nMatrices(); ++i)
703 {
704 const auto& bpsi = this->psi(i).boundaryField();
705
706 // Cache to-be implicit internal/boundary
708 FieldField<Field, Type> internal(bpsi.size());
709
710 label implicit = 0;
711 forAll(bpsi, patchI)
712 {
713 label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
714 if (globalPatchId == -1)
715 {
717 (
718 implicit,
719 matrix(i).boundaryCoeffs()[patchI].clone()
720 );
721 internal.set
722 (
723 implicit,
724 matrix(i).internalCoeffs()[patchI].clone()
725 );
726 implicit++;
727 }
728 }
729
730 // Update non-implicit patches (re-order)
731 forAll(bpsi, patchI)
732 {
733 label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
734 if (globalPatchId != -1)
735 {
736 if (matrix(i).internalCoeffs().set(patchI))
737 {
738 internalCoeffs_.set
739 (
740 globalPatchId,
741 matrix(i).internalCoeffs()[patchI].clone()
742 );
743 }
744
745 if (matrix(i).boundaryCoeffs().set(patchI))
746 {
747 boundaryCoeffs_.set
748 (
749 globalPatchId,
750 matrix(i).boundaryCoeffs()[patchI].clone()
751 );
752 }
753 }
754 }
755
756 // Store implicit patches at the end of the list
757 implicit = 0;
758 forAll(bpsi, patchI)
759 {
760 label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
761 if (globalPatchId == -1)
762 {
763 const label implicitPatchId =
764 lduMeshPtr()->patchLocalToGlobalMap()[i][patchI];
765
766 internalCoeffs_.set
767 (
768 implicitPatchId, internal[implicit].clone()
769 );
770 boundaryCoeffs_.set
771 (
772 implicitPatchId, boundary[implicit].clone()
773 );
774
775 implicit++;
776 }
777 }
778 }
779
780// forAll(internalCoeffs_, patchI)
781// {
782// DebugVar(patchI)
783// DebugVar(internalCoeffs_[patchI])
784// DebugVar(boundaryCoeffs_[patchI])
785// }
786}
787
788
789template<class Type>
791{
792 for (label i=0; i < nMatrices(); ++i)
793 {
794 forAll(psi(i).boundaryField(), patchI)
795 {
796 label globalPatchId = lduMeshPtr()->patchMap()[i][patchI];
797
798 if (globalPatchId == -1)
799 {
800 psi(i).boundaryFieldRef()[patchI].manipulateMatrix
801 (
802 *this,
803 i,
804 cmp
805 );
806 }
807 }
808 }
809}
810
811
812template<class Type>
814{
815 const labelListList& faceMap = lduMeshPtr()->faceMap();
816 const labelList& cellMap = lduMeshPtr()->cellOffsets();
817
818 label newFaces = lduMeshPtr()->lduAddr().upperAddr().size();
819 label newCells = lduMeshPtr()->lduAddr().size();
820
821 scalarField lowerAssemb(newFaces, Zero);
822 scalarField upperAssemb(newFaces, Zero);
823 scalarField diagAssemb(newCells, Zero);
824 Field<Type> sourceAssemb(newCells, Zero);
825
826 bool asymmetricAssemby = false;
827 for (label i=0; i < nMatrices(); ++i)
828 {
829 if (matrix(i).asymmetric())
830 {
831 asymmetricAssemby = true;
832 }
833 }
834 // Move append contents into intermediate list
835 for (label i=0; i < nMatrices(); ++i)
836 {
837 if (asymmetricAssemby)
838 {
839 const scalarField lowerSub(matrix(i).lower());
840 forAll(lowerSub, facei)
841 {
842 lowerAssemb[faceMap[i][facei]] = lowerSub[facei];
843 }
844 }
845
846 const scalarField upperSub(matrix(i).upper());
847 const scalarField diagSub(matrix(i).diag());
848 const Field<Type> sourceSub(matrix(i).source());
849
850 forAll(upperSub, facei)
851 {
852 upperAssemb[faceMap[i][facei]] = upperSub[facei];
853 }
854
855 forAll(diagSub, celli)
856 {
857 const label globalCelli = cellMap[i] + celli;
858 diagAssemb[globalCelli] = diagSub[celli];
859 sourceAssemb[globalCelli] = sourceSub[celli];
860 }
861 }
862
863 if (asymmetricAssemby)
864 {
865 lower().setSize(newFaces, Zero);
866 lower() = lowerAssemb;
867 }
868 upper().setSize(newFaces, Zero);
869 upper() = upperAssemb;
870
871 diag().setSize(newCells, Zero);
872 diag() = diagAssemb;
873
874 source().setSize(newCells, Zero);
875 source() = sourceAssemb;
876}
877
878
879template<class Type>
881{
882 const lduPrimitiveMeshAssembly* lduAssemMeshPtr =
883 psi_.mesh().thisDb().objectRegistry::template findObject
884 <
886 > (lduAssemblyName_);
887
888 return const_cast<lduPrimitiveMeshAssembly*>(lduAssemMeshPtr);
889}
890
891
892template<class Type>
894{
895 return
896 (
897 psi_.mesh().thisDb().objectRegistry::template cfindObject
898 <
900 > (lduAssemblyName_)
901 );
902}
903
904
905template<class Type>
907{
908 lduPrimitiveMeshAssembly* ptr = lduMeshPtr();
909
911 (
912 lduAssemblyName_,
913 psi_.mesh().time().timeName(),
914 psi_.mesh(),
917 );
918
919 UPtrList<lduMesh> uMeshPtr(nMatrices());
920
922 uFieldPtr(nMatrices());
923
924 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
925 {
926 const fvMesh& meshi = this->psi(fieldi).mesh();
927 uMeshPtr.set
928 (
929 fieldi,
930 &const_cast<fvMesh&>(meshi)
931 );
932 uFieldPtr.set(fieldi, &this->psi(fieldi));
933 }
934
935 if (!ptr)
936 {
937 lduPrimitiveMeshAssembly* lduAssemMeshPtr =
938 new lduPrimitiveMeshAssembly(io, uMeshPtr);
939
940 lduAssemMeshPtr->store();
941 lduAssemMeshPtr->update(uFieldPtr);
942
943 Info
944 << "Creating lduPrimitiveAssembly: " << lduAssemblyName_ << endl;
945 }
946 else if
947 (
948 psi_.mesh().changing() && !psi_.mesh().upToDatePoints(*ptr)
949 )
950 {
951 // Clear losortPtr_, ownerStartPtr_, losortStartPtr_
952 ptr->lduAddr().clearOut();
953 ptr->update(uFieldPtr);
954 psi_.mesh().setUpToDatePoints(*ptr);
955
956 Info
957 << "Updating lduPrimitiveAssembly: " << lduAssemblyName_ << endl;
958 }
959 else
960 {
961 Info
962 << "Using lduPrimitiveAssembly: " << lduAssemblyName_ << endl;
963 }
964}
965
966
967template<class Type>
969(
970 const labelUList& cellLabels,
971 const Type& value
972)
973{
974 this->setValuesFromList(cellLabels, UniformList<Type>(value));
975}
976
977
978template<class Type>
980(
981 const labelUList& cellLabels,
982 const UList<Type>& values
983)
984{
985 this->setValuesFromList(cellLabels, values);
986}
987
988
989template<class Type>
991(
992 const labelUList& cellLabels,
993 const UIndirectList<Type>& values
994)
995{
996 this->setValuesFromList(cellLabels, values);
997}
998
999
1000template<class Type>
1002(
1003 const label celli,
1004 const Type& value,
1005 const bool forceReference
1006)
1007{
1008 if ((forceReference || psi_.needReference()) && celli >= 0)
1009 {
1010 source()[celli] += diag()[celli]*value;
1011 diag()[celli] += diag()[celli];
1012 }
1013}
1014
1015
1016template<class Type>
1018(
1019 const labelUList& cellLabels,
1020 const Type& value,
1021 const bool forceReference
1022)
1023{
1024 if (forceReference || psi_.needReference())
1025 {
1026 forAll(cellLabels, celli)
1027 {
1028 const label cellId = cellLabels[celli];
1029 if (cellId >= 0)
1030 {
1031 source()[cellId] += diag()[cellId]*value;
1032 diag()[cellId] += diag()[cellId];
1033 }
1034 }
1035 }
1036}
1037
1038
1039template<class Type>
1041(
1042 const labelUList& cellLabels,
1043 const UList<Type>& values,
1044 const bool forceReference
1045)
1046{
1047 if (forceReference || psi_.needReference())
1048 {
1049 forAll(cellLabels, celli)
1050 {
1051 const label cellId = cellLabels[celli];
1052 if (cellId >= 0)
1053 {
1054 source()[cellId] += diag()[cellId]*values[celli];
1055 diag()[cellId] += diag()[cellId];
1056 }
1057 }
1058 }
1059}
1060
1061
1062template<class Type>
1064{
1065 subMatrices_.append(matrix.clone());
1066 ++nMatrix_;
1067
1068 if (dimensions_ != matrix.dimensions())
1069 {
1071 << "incompatible dimensions for matrix addition "
1072 << endl << " "
1073 << "[" << dimensions_ << " ] "
1074 << " [" << matrix.dimensions() << " ]"
1075 << abort(FatalError);
1076 }
1077
1078 for (label fieldi = 0; fieldi < nMatrices(); fieldi++)
1079 {
1080 if (checkImplicit(fieldi))
1081 {
1082 break;
1083 }
1084 }
1085
1086 internalCoeffs_.clear();
1087 boundaryCoeffs_.clear();
1088}
1089
1090
1091template<class Type>
1093{
1094 if (alpha <= 0)
1095 {
1096 return;
1097 }
1098
1100 << "Relaxing " << psi_.name() << " by " << alpha << endl;
1101
1102 Field<Type>& S = source();
1103 scalarField& D = diag();
1104
1105 // Store the current unrelaxed diagonal for use in updating the source
1106 scalarField D0(D);
1107
1108 // Calculate the sum-mag off-diagonal from the interior faces
1109 scalarField sumOff(D.size(), Zero);
1110 sumMagOffDiag(sumOff);
1111
1112 // Handle the boundary contributions to the diagonal
1113 forAll(psi_.boundaryField(), patchi)
1114 {
1115 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1116
1117 if (ptf.size())
1118 {
1119 const labelUList& pa = lduAddr().patchAddr(patchi);
1120 Field<Type>& iCoeffs = internalCoeffs_[patchi];
1121
1122 if (ptf.coupled())
1123 {
1124 const Field<Type>& pCoeffs = boundaryCoeffs_[patchi];
1125
1126 // For coupled boundaries add the diagonal and
1127 // off-diagonal contributions
1128 forAll(pa, face)
1129 {
1130 D[pa[face]] += component(iCoeffs[face], 0);
1131 sumOff[pa[face]] += mag(component(pCoeffs[face], 0));
1132 }
1133 }
1134 else
1135 {
1136 // For non-coupled boundaries add the maximum magnitude diagonal
1137 // contribution to ensure stability
1138 forAll(pa, face)
1139 {
1140 D[pa[face]] += cmptMax(cmptMag(iCoeffs[face]));
1141 }
1142 }
1143 }
1144 }
1145
1146
1147 if (debug)
1148 {
1149 // Calculate amount of non-dominance.
1150 label nNon = 0;
1151 scalar maxNon = 0.0;
1152 scalar sumNon = 0.0;
1153 forAll(D, celli)
1154 {
1155 scalar d = (sumOff[celli] - D[celli])/mag(D[celli]);
1156
1157 if (d > 0)
1158 {
1159 nNon++;
1160 maxNon = max(maxNon, d);
1161 sumNon += d;
1162 }
1163 }
1164
1165 reduce(nNon, sumOp<label>(), UPstream::msgType(), psi_.mesh().comm());
1166 reduce
1167 (
1168 maxNon,
1169 maxOp<scalar>(),
1171 psi_.mesh().comm()
1172 );
1173 reduce
1174 (
1175 sumNon,
1176 sumOp<scalar>(),
1178 psi_.mesh().comm()
1179 );
1180 sumNon /= returnReduce
1181 (
1182 D.size(),
1183 sumOp<label>(),
1185 psi_.mesh().comm()
1186 );
1187
1189 << "Matrix dominance test for " << psi_.name() << nl
1190 << " number of non-dominant cells : " << nNon << nl
1191 << " maximum relative non-dominance : " << maxNon << nl
1192 << " average relative non-dominance : " << sumNon << nl
1193 << endl;
1194 }
1195
1196
1197 // Ensure the matrix is diagonally dominant...
1198 // Assumes that the central coefficient is positive and ensures it is
1199 forAll(D, celli)
1200 {
1201 D[celli] = max(mag(D[celli]), sumOff[celli]);
1202 }
1203
1204 // ... then relax
1205 D /= alpha;
1206
1207 // Now remove the diagonal contribution from coupled boundaries
1208 forAll(psi_.boundaryField(), patchi)
1209 {
1210 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1211
1212 if (ptf.size())
1213 {
1214 const labelUList& pa = lduAddr().patchAddr(patchi);
1215 Field<Type>& iCoeffs = internalCoeffs_[patchi];
1216
1217 if (ptf.coupled())
1218 {
1219 forAll(pa, face)
1220 {
1221 D[pa[face]] -= component(iCoeffs[face], 0);
1222 }
1223 }
1224 else
1225 {
1226 forAll(pa, face)
1227 {
1228 D[pa[face]] -= cmptMin(iCoeffs[face]);
1229 }
1230 }
1231 }
1232 }
1233
1234 // Finally add the relaxation contribution to the source.
1235 S += (D - D0)*psi_.primitiveField();
1236}
1237
1238
1239template<class Type>
1241{
1242 word name = psi_.select
1243 (
1244 psi_.mesh().data::template getOrDefault<bool>
1245 ("finalIteration", false)
1246 );
1247
1248 if (psi_.mesh().relaxEquation(name))
1249 {
1250 relax(psi_.mesh().equationRelaxationFactor(name));
1251 }
1252}
1253
1254
1255template<class Type>
1257(
1259 Boundary& bFields
1260)
1261{
1262 forAll(bFields, patchi)
1263 {
1264 bFields[patchi].manipulateMatrix(*this);
1265 }
1266}
1267
1268
1269template<class Type>
1271{
1272 tmp<scalarField> tdiag(new scalarField(diag()));
1273 addCmptAvBoundaryDiag(tdiag.ref());
1274 return tdiag;
1275}
1276
1277
1278template<class Type>
1280{
1282
1283 forAll(psi_.boundaryField(), patchi)
1284 {
1285 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1286
1287 if (!ptf.coupled() && ptf.size())
1288 {
1289 addToInternalField
1290 (
1291 lduAddr().patchAddr(patchi),
1292 internalCoeffs_[patchi],
1293 tdiag.ref()
1294 );
1295 }
1296 }
1297
1298 return tdiag;
1299}
1300
1301
1302template<class Type>
1304{
1306 (
1307 new volScalarField
1308 (
1309 IOobject
1310 (
1311 "A("+psi_.name()+')',
1312 psi_.instance(),
1313 psi_.mesh(),
1316 ),
1317 psi_.mesh(),
1318 dimensions_/psi_.dimensions()/dimVol,
1319 extrapolatedCalculatedFvPatchScalarField::typeName
1320 )
1321 );
1322
1323 tAphi.ref().primitiveFieldRef() = D()/psi_.mesh().V();
1324 tAphi.ref().correctBoundaryConditions();
1325
1326 return tAphi;
1327}
1328
1329
1330template<class Type>
1333{
1335 (
1337 (
1338 IOobject
1339 (
1340 "H("+psi_.name()+')',
1341 psi_.instance(),
1342 psi_.mesh(),
1345 ),
1346 psi_.mesh(),
1347 dimensions_/dimVol,
1348 extrapolatedCalculatedFvPatchScalarField::typeName
1349 )
1350 );
1352
1353 // Loop over field components
1354 for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
1355 {
1356 scalarField psiCmpt(psi_.primitiveField().component(cmpt));
1357
1358 scalarField boundaryDiagCmpt(psi_.size(), Zero);
1359 addBoundaryDiag(boundaryDiagCmpt, cmpt);
1360 boundaryDiagCmpt.negate();
1361 addCmptAvBoundaryDiag(boundaryDiagCmpt);
1362
1363 Hphi.primitiveFieldRef().replace(cmpt, boundaryDiagCmpt*psiCmpt);
1364 }
1365
1366 Hphi.primitiveFieldRef() += lduMatrix::H(psi_.primitiveField()) + source_;
1367 addBoundarySource(Hphi.primitiveFieldRef());
1368
1369 Hphi.primitiveFieldRef() /= psi_.mesh().V();
1371
1372 typename Type::labelType validComponents
1373 (
1374 psi_.mesh().template validComponents<Type>()
1375 );
1376
1377 for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
1378 {
1379 if (validComponents[cmpt] == -1)
1380 {
1381 Hphi.replace
1382 (
1383 cmpt,
1385 );
1386 }
1387 }
1388
1389 return tHphi;
1390}
1391
1392
1393template<class Type>
1395{
1397 (
1398 new volScalarField
1399 (
1400 IOobject
1401 (
1402 "H(1)",
1403 psi_.instance(),
1404 psi_.mesh(),
1407 ),
1408 psi_.mesh(),
1409 dimensions_/(dimVol*psi_.dimensions()),
1410 extrapolatedCalculatedFvPatchScalarField::typeName
1411 )
1412 );
1413 volScalarField& H1_ = tH1.ref();
1414
1416
1417 forAll(psi_.boundaryField(), patchi)
1418 {
1419 const fvPatchField<Type>& ptf = psi_.boundaryField()[patchi];
1420
1421 if (ptf.coupled() && ptf.size())
1422 {
1423 addToInternalField
1424 (
1425 lduAddr().patchAddr(patchi),
1426 boundaryCoeffs_[patchi].component(0),
1427 H1_
1428 );
1429 }
1430 }
1431
1432 H1_.primitiveFieldRef() /= psi_.mesh().V();
1434
1435 return tH1;
1436}
1437
1438
1439
1440template<class Type>
1443flux() const
1444{
1445 if (!psi_.mesh().fluxRequired(psi_.name()))
1446 {
1448 << "flux requested but " << psi_.name()
1449 << " not specified in the fluxRequired sub-dictionary"
1450 " of fvSchemes."
1451 << abort(FatalError);
1452 }
1453
1454 if (nMatrices() > 1)
1455 {
1457 << "Flux requested but " << psi_.name()
1458 << " can't handle multiple fvMatrix."
1459 << abort(FatalError);
1460 }
1461
1462 // construct GeometricField<Type, fvsPatchField, surfaceMesh>
1464 (
1466 (
1467 IOobject
1468 (
1469 "flux("+psi_.name()+')',
1470 psi_.instance(),
1471 psi_.mesh(),
1474 ),
1475 psi_.mesh(),
1476 dimensions()
1477 )
1478 );
1480 tfieldFlux.ref();
1481
1482 fieldFlux.setOriented();
1483
1484 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
1485 {
1486 fieldFlux.primitiveFieldRef().replace
1487 (
1488 cmpt,
1489 lduMatrix::faceH(psi_.primitiveField().component(cmpt))
1490 );
1491 }
1492
1493 FieldField<Field, Type> InternalContrib = internalCoeffs_;
1494
1495 label fieldi = 0;
1496 if (!useImplicit_)
1497 {
1498 forAll(InternalContrib, patchi)
1499 {
1500 InternalContrib[patchi] =
1502 (
1503 InternalContrib[patchi],
1504 psi_.boundaryField()[patchi].patchInternalField()
1505 );
1506 }
1507 }
1508 else
1509 {
1510 FieldField<Field, Type> fluxInternalContrib(internalCoeffs_);
1511
1512 mapContributions(fieldi, fluxInternalContrib, InternalContrib, true);
1513 }
1514
1515 FieldField<Field, Type> NeighbourContrib = boundaryCoeffs_;
1516
1517 if (!useImplicit_)
1518 {
1519 forAll(NeighbourContrib, patchi)
1520 {
1521 if (psi_.boundaryField()[patchi].coupled())
1522 {
1523 NeighbourContrib[patchi] =
1525 (
1526 NeighbourContrib[patchi],
1527 psi_.boundaryField()[patchi].patchNeighbourField()
1528 );
1529 }
1530 }
1531 }
1532 else
1533 {
1534 FieldField<Field, Type> fluxBoundaryContrib(boundaryCoeffs_);
1535
1536 mapContributions(fieldi, fluxBoundaryContrib, NeighbourContrib, false);
1537 }
1538
1540 Boundary& ffbf = fieldFlux.boundaryFieldRef();
1541
1542 forAll(ffbf, patchi)
1543 {
1544 ffbf[patchi] = InternalContrib[patchi] - NeighbourContrib[patchi];
1545 //DebugVar(gSum(ffbf[patchi]))
1546 }
1547
1548 if (faceFluxCorrectionPtr_)
1549 {
1550 fieldFlux += *faceFluxCorrectionPtr_;
1551 }
1552
1553 return tfieldFlux;
1554}
1555
1556
1557template<class Type>
1559{
1560 return psi_.mesh().solverDict
1561 (
1562 psi_.select
1563 (
1564 psi_.mesh().data::template getOrDefault<bool>
1565 ("finalIteration", false)
1566 )
1567 );
1568}
1569
1570
1571// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1572
1573template<class Type>
1575{
1576 if (this == &fvmv)
1577 {
1578 return; // Self-assignment is a no-op
1579 }
1580
1581 if (&psi_ != &(fvmv.psi_))
1582 {
1584 << "different fields"
1585 << abort(FatalError);
1586 }
1587
1588 dimensions_ = fvmv.dimensions_;
1590 source_ = fvmv.source_;
1591 internalCoeffs_ = fvmv.internalCoeffs_;
1592 boundaryCoeffs_ = fvmv.boundaryCoeffs_;
1593
1594 if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1595 {
1596 *faceFluxCorrectionPtr_ = *fvmv.faceFluxCorrectionPtr_;
1597 }
1598 else if (fvmv.faceFluxCorrectionPtr_)
1599 {
1600 faceFluxCorrectionPtr_ =
1602 (*fvmv.faceFluxCorrectionPtr_);
1603 }
1604
1605 useImplicit_ = fvmv.useImplicit_;
1606 lduAssemblyName_ = fvmv.lduAssemblyName_;
1607}
1608
1609
1610template<class Type>
1612{
1613 operator=(tfvmv());
1614 tfvmv.clear();
1615}
1616
1617
1618template<class Type>
1620{
1622 source_.negate();
1623 internalCoeffs_.negate();
1624 boundaryCoeffs_.negate();
1625
1626 if (faceFluxCorrectionPtr_)
1627 {
1628 faceFluxCorrectionPtr_->negate();
1629 }
1630}
1631
1632
1633template<class Type>
1635{
1636 checkMethod(*this, fvmv, "+=");
1637
1638 dimensions_ += fvmv.dimensions_;
1640 source_ += fvmv.source_;
1641 internalCoeffs_ += fvmv.internalCoeffs_;
1642 boundaryCoeffs_ += fvmv.boundaryCoeffs_;
1643
1644 useImplicit_ = fvmv.useImplicit_;
1645 lduAssemblyName_ = fvmv.lduAssemblyName_;
1646 nMatrix_ = fvmv.nMatrix_;
1647
1648 if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1649 {
1650 *faceFluxCorrectionPtr_ += *fvmv.faceFluxCorrectionPtr_;
1651 }
1652 else if (fvmv.faceFluxCorrectionPtr_)
1653 {
1654 faceFluxCorrectionPtr_ = new
1656 (
1657 *fvmv.faceFluxCorrectionPtr_
1658 );
1659 }
1660}
1661
1662
1663template<class Type>
1665{
1666 operator+=(tfvmv());
1667 tfvmv.clear();
1668}
1669
1670
1671template<class Type>
1673{
1674 checkMethod(*this, fvmv, "-=");
1675
1676 dimensions_ -= fvmv.dimensions_;
1678 source_ -= fvmv.source_;
1679 internalCoeffs_ -= fvmv.internalCoeffs_;
1680 boundaryCoeffs_ -= fvmv.boundaryCoeffs_;
1681
1682 useImplicit_ = fvmv.useImplicit_;
1683 lduAssemblyName_ = fvmv.lduAssemblyName_;
1684 nMatrix_ = fvmv.nMatrix_;
1685
1686 if (faceFluxCorrectionPtr_ && fvmv.faceFluxCorrectionPtr_)
1687 {
1688 *faceFluxCorrectionPtr_ -= *fvmv.faceFluxCorrectionPtr_;
1689 }
1690 else if (fvmv.faceFluxCorrectionPtr_)
1691 {
1692 faceFluxCorrectionPtr_ =
1694 (-*fvmv.faceFluxCorrectionPtr_);
1695 }
1696}
1697
1698
1699template<class Type>
1701{
1702 operator-=(tfvmv());
1703 tfvmv.clear();
1704}
1705
1706
1707template<class Type>
1709(
1711)
1712{
1713 checkMethod(*this, su, "+=");
1714 source() -= su.mesh().V()*su.field();
1715}
1716
1717
1718template<class Type>
1720(
1722)
1723{
1724 operator+=(tsu());
1725 tsu.clear();
1726}
1727
1728
1729template<class Type>
1731(
1733)
1734{
1735 operator+=(tsu());
1736 tsu.clear();
1737}
1738
1739
1740template<class Type>
1742(
1744)
1745{
1746 checkMethod(*this, su, "-=");
1747 source() += su.mesh().V()*su.field();
1748}
1749
1750
1751template<class Type>
1753(
1755)
1756{
1757 operator-=(tsu());
1758 tsu.clear();
1759}
1760
1761
1762template<class Type>
1764(
1766)
1767{
1768 operator-=(tsu());
1769 tsu.clear();
1770}
1771
1772
1773template<class Type>
1775(
1776 const dimensioned<Type>& su
1777)
1778{
1779 source() -= psi().mesh().V()*su;
1780}
1781
1782
1783template<class Type>
1785(
1786 const dimensioned<Type>& su
1787)
1788{
1789 source() += psi().mesh().V()*su;
1790}
1791
1792
1793template<class Type>
1795{}
1796
1797
1798template<class Type>
1800{}
1801
1802
1803template<class Type>
1805(
1806 const volScalarField::Internal& dsf
1807)
1808{
1809 dimensions_ *= dsf.dimensions();
1810 lduMatrix::operator*=(dsf.field());
1811 source_ *= dsf.field();
1812
1813 forAll(boundaryCoeffs_, patchi)
1814 {
1815 scalarField pisf
1816 (
1817 dsf.mesh().boundary()[patchi].patchInternalField(dsf.field())
1818 );
1819
1820 internalCoeffs_[patchi] *= pisf;
1821 boundaryCoeffs_[patchi] *= pisf;
1822 }
1823
1824 if (faceFluxCorrectionPtr_)
1825 {
1827 << "cannot scale a matrix containing a faceFluxCorrection"
1828 << abort(FatalError);
1829 }
1830}
1831
1832
1833template<class Type>
1835(
1837)
1838{
1839 operator*=(tfld());
1840 tfld.clear();
1841}
1842
1843
1844template<class Type>
1846(
1847 const tmp<volScalarField>& tfld
1848)
1849{
1850 operator*=(tfld());
1851 tfld.clear();
1852}
1853
1854
1855template<class Type>
1857(
1858 const dimensioned<scalar>& ds
1859)
1860{
1861 dimensions_ *= ds.dimensions();
1862 lduMatrix::operator*=(ds.value());
1863 source_ *= ds.value();
1864 internalCoeffs_ *= ds.value();
1865 boundaryCoeffs_ *= ds.value();
1866
1867 if (faceFluxCorrectionPtr_)
1868 {
1869 *faceFluxCorrectionPtr_ *= ds.value();
1870 }
1871}
1872
1873
1874// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
1875
1876template<class Type>
1878(
1879 const fvMatrix<Type>& mat1,
1880 const fvMatrix<Type>& mat2,
1881 const char* op
1882)
1883{
1884 if (&mat1.psi() != &mat2.psi())
1885 {
1887 << "Incompatible fields for operation\n "
1888 << "[" << mat1.psi().name() << "] "
1889 << op
1890 << " [" << mat2.psi().name() << "]"
1891 << abort(FatalError);
1892 }
1893
1894 if
1895 (
1897 && mat1.dimensions() != mat2.dimensions()
1898 )
1899 {
1901 << "Incompatible dimensions for operation\n "
1902 << "[" << mat1.psi().name() << mat1.dimensions()/dimVolume << " ] "
1903 << op
1904 << " [" << mat2.psi().name() << mat2.dimensions()/dimVolume << " ]"
1905 << abort(FatalError);
1906 }
1907}
1908
1909
1910template<class Type>
1912(
1913 const fvMatrix<Type>& mat,
1915 const char* op
1916)
1917{
1918 if
1919 (
1921 && mat.dimensions()/dimVolume != fld.dimensions()
1922 )
1923 {
1925 << "Incompatible dimensions for operation\n "
1926 << "[" << mat.psi().name() << mat.dimensions()/dimVolume << " ] "
1927 << op
1928 << " [" << fld.name() << fld.dimensions() << " ]"
1929 << abort(FatalError);
1930 }
1931}
1932
1933
1934template<class Type>
1936(
1937 const fvMatrix<Type>& mat,
1938 const dimensioned<Type>& dt,
1939 const char* op
1940)
1941{
1942 if
1943 (
1945 && mat.dimensions()/dimVolume != dt.dimensions()
1946 )
1947 {
1949 << "Incompatible dimensions for operation\n "
1950 << "[" << mat.psi().name() << mat.dimensions()/dimVolume << " ] "
1951 << op
1952 << " [" << dt.name() << dt.dimensions() << " ]"
1953 << abort(FatalError);
1954 }
1955}
1956
1957
1958template<class Type>
1960(
1961 fvMatrix<Type>& fvm,
1962 const dictionary& solverControls
1963)
1964{
1965 return fvm.solve(solverControls);
1966}
1967
1968template<class Type>
1970(
1971 const tmp<fvMatrix<Type>>& tmat,
1972 const dictionary& solverControls
1973)
1974{
1975 SolverPerformance<Type> solverPerf(tmat.constCast().solve(solverControls));
1976
1977 tmat.clear();
1978
1979 return solverPerf;
1980}
1981
1982
1983template<class Type>
1984Foam::SolverPerformance<Type> Foam::solve(fvMatrix<Type>& fvm)
1985{
1986 return fvm.solve();
1987}
1988
1989template<class Type>
1990Foam::SolverPerformance<Type> Foam::solve(const tmp<fvMatrix<Type>>& tmat)
1991{
1992 SolverPerformance<Type> solverPerf(tmat.constCast().solve());
1993
1994 tmat.clear();
1995
1996 return solverPerf;
1997}
1998
1999
2000template<class Type>
2002(
2003 const fvMatrix<Type>& A
2004)
2005{
2006 tmp<Foam::fvMatrix<Type>> tAcorr = A - (A & A.psi());
2007
2008 // Delete the faceFluxCorrection from the correction matrix
2009 // as it does not have a clear meaning or purpose
2010 deleteDemandDrivenData(tAcorr.ref().faceFluxCorrectionPtr());
2011
2012 return tAcorr;
2013}
2014
2015
2016template<class Type>
2018(
2019 const tmp<fvMatrix<Type>>& tA
2020)
2021{
2022 tmp<Foam::fvMatrix<Type>> tAcorr = tA - (tA() & tA().psi());
2023
2024 // Delete the faceFluxCorrection from the correction matrix
2025 // as it does not have a clear meaning or purpose
2026 deleteDemandDrivenData(tAcorr.ref().faceFluxCorrectionPtr());
2027
2028 return tAcorr;
2029}
2030
2031
2032// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
2033
2034template<class Type>
2035Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2036(
2037 const fvMatrix<Type>& A,
2038 const fvMatrix<Type>& B
2039)
2040{
2041 checkMethod(A, B, "==");
2042 return (A - B);
2043}
2044
2045template<class Type>
2046Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2047(
2048 const tmp<fvMatrix<Type>>& tA,
2049 const fvMatrix<Type>& B
2050)
2051{
2052 checkMethod(tA(), B, "==");
2053 return (tA - B);
2054}
2055
2056template<class Type>
2057Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2058(
2059 const fvMatrix<Type>& A,
2060 const tmp<fvMatrix<Type>>& tB
2061)
2062{
2063 checkMethod(A, tB(), "==");
2064 return (A - tB);
2065}
2066
2067template<class Type>
2068Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2069(
2070 const tmp<fvMatrix<Type>>& tA,
2071 const tmp<fvMatrix<Type>>& tB
2072)
2073{
2074 checkMethod(tA(), tB(), "==");
2075 return (tA - tB);
2076}
2077
2078template<class Type>
2079Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2080(
2081 const fvMatrix<Type>& A,
2082 const DimensionedField<Type, volMesh>& su
2083)
2084{
2085 checkMethod(A, su, "==");
2086 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2087 tC.ref().source() += su.mesh().V()*su.field();
2088 return tC;
2089}
2090
2091template<class Type>
2092Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2093(
2094 const fvMatrix<Type>& A,
2095 const tmp<DimensionedField<Type, volMesh>>& tsu
2096)
2097{
2098 checkMethod(A, tsu(), "==");
2099 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2100 tC.ref().source() += tsu().mesh().V()*tsu().field();
2101 tsu.clear();
2102 return tC;
2103}
2104
2105template<class Type>
2106Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2107(
2108 const fvMatrix<Type>& A,
2109 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2110)
2111{
2112 checkMethod(A, tsu(), "==");
2113 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2114 tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2115 tsu.clear();
2116 return tC;
2117}
2118
2119template<class Type>
2120Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2121(
2122 const tmp<fvMatrix<Type>>& tA,
2123 const DimensionedField<Type, volMesh>& su
2124)
2125{
2126 checkMethod(tA(), su, "==");
2127 tmp<fvMatrix<Type>> tC(tA.ptr());
2128 tC.ref().source() += su.mesh().V()*su.field();
2129 return tC;
2130}
2131
2132template<class Type>
2133Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2134(
2135 const tmp<fvMatrix<Type>>& tA,
2136 const tmp<DimensionedField<Type, volMesh>>& tsu
2137)
2138{
2139 checkMethod(tA(), tsu(), "==");
2140 tmp<fvMatrix<Type>> tC(tA.ptr());
2141 tC.ref().source() += tsu().mesh().V()*tsu().field();
2142 tsu.clear();
2143 return tC;
2144}
2145
2146template<class Type>
2147Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2148(
2149 const tmp<fvMatrix<Type>>& tA,
2150 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2151)
2152{
2153 checkMethod(tA(), tsu(), "==");
2154 tmp<fvMatrix<Type>> tC(tA.ptr());
2155 tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2156 tsu.clear();
2157 return tC;
2158}
2159
2160template<class Type>
2161Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2162(
2163 const fvMatrix<Type>& A,
2164 const dimensioned<Type>& su
2165)
2166{
2167 checkMethod(A, su, "==");
2168 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2169 tC.ref().source() += A.psi().mesh().V()*su.value();
2170 return tC;
2171}
2172
2173template<class Type>
2174Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2175(
2176 const tmp<fvMatrix<Type>>& tA,
2177 const dimensioned<Type>& su
2178)
2179{
2180 checkMethod(tA(), su, "==");
2181 tmp<fvMatrix<Type>> tC(tA.ptr());
2182 tC.ref().source() += tC().psi().mesh().V()*su.value();
2183 return tC;
2184}
2185
2186template<class Type>
2187Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2188(
2189 const fvMatrix<Type>& A,
2190 const Foam::zero
2191)
2192{
2193 return A;
2194}
2195
2196
2197template<class Type>
2198Foam::tmp<Foam::fvMatrix<Type>> Foam::operator==
2199(
2200 const tmp<fvMatrix<Type>>& tA,
2201 const Foam::zero
2202)
2203{
2204 return tA;
2205}
2206
2207
2208template<class Type>
2209Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2210(
2211 const fvMatrix<Type>& A
2212)
2213{
2214 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2215 tC.ref().negate();
2216 return tC;
2217}
2218
2219template<class Type>
2220Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2221(
2222 const tmp<fvMatrix<Type>>& tA
2223)
2224{
2225 tmp<fvMatrix<Type>> tC(tA.ptr());
2226 tC.ref().negate();
2227 return tC;
2228}
2229
2230
2231template<class Type>
2232Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2233(
2234 const fvMatrix<Type>& A,
2235 const fvMatrix<Type>& B
2236)
2237{
2238 checkMethod(A, B, "+");
2239 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2240 tC.ref() += B;
2241 return tC;
2242}
2243
2244template<class Type>
2245Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2246(
2247 const tmp<fvMatrix<Type>>& tA,
2248 const fvMatrix<Type>& B
2249)
2250{
2251 checkMethod(tA(), B, "+");
2252 tmp<fvMatrix<Type>> tC(tA.ptr());
2253 tC.ref() += B;
2254 return tC;
2255}
2256
2257template<class Type>
2258Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2259(
2260 const fvMatrix<Type>& A,
2261 const tmp<fvMatrix<Type>>& tB
2262)
2263{
2264 checkMethod(A, tB(), "+");
2265 tmp<fvMatrix<Type>> tC(tB.ptr());
2266 tC.ref() += A;
2267 return tC;
2268}
2269
2270template<class Type>
2271Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2272(
2273 const tmp<fvMatrix<Type>>& tA,
2274 const tmp<fvMatrix<Type>>& tB
2275)
2276{
2277 checkMethod(tA(), tB(), "+");
2278 tmp<fvMatrix<Type>> tC(tA.ptr());
2279 tC.ref() += tB();
2280 tB.clear();
2281 return tC;
2282}
2283
2284template<class Type>
2285Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2286(
2287 const fvMatrix<Type>& A,
2288 const DimensionedField<Type, volMesh>& su
2289)
2290{
2291 checkMethod(A, su, "+");
2292 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2293 tC.ref().source() -= su.mesh().V()*su.field();
2294 return tC;
2295}
2296
2297template<class Type>
2298Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2299(
2300 const fvMatrix<Type>& A,
2301 const tmp<DimensionedField<Type, volMesh>>& tsu
2302)
2303{
2304 checkMethod(A, tsu(), "+");
2305 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2306 tC.ref().source() -= tsu().mesh().V()*tsu().field();
2307 tsu.clear();
2308 return tC;
2309}
2310
2311template<class Type>
2312Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2313(
2314 const fvMatrix<Type>& A,
2315 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2316)
2317{
2318 checkMethod(A, tsu(), "+");
2319 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2320 tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2321 tsu.clear();
2322 return tC;
2323}
2324
2325template<class Type>
2326Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2327(
2328 const tmp<fvMatrix<Type>>& tA,
2329 const DimensionedField<Type, volMesh>& su
2330)
2331{
2332 checkMethod(tA(), su, "+");
2333 tmp<fvMatrix<Type>> tC(tA.ptr());
2334 tC.ref().source() -= su.mesh().V()*su.field();
2335 return tC;
2336}
2337
2338template<class Type>
2339Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2340(
2341 const tmp<fvMatrix<Type>>& tA,
2342 const tmp<DimensionedField<Type, volMesh>>& tsu
2343)
2344{
2345 checkMethod(tA(), tsu(), "+");
2346 tmp<fvMatrix<Type>> tC(tA.ptr());
2347 tC.ref().source() -= tsu().mesh().V()*tsu().field();
2348 tsu.clear();
2349 return tC;
2350}
2351
2352template<class Type>
2353Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2354(
2355 const tmp<fvMatrix<Type>>& tA,
2356 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2357)
2358{
2359 checkMethod(tA(), tsu(), "+");
2360 tmp<fvMatrix<Type>> tC(tA.ptr());
2361 tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2362 tsu.clear();
2363 return tC;
2364}
2365
2366template<class Type>
2367Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2368(
2369 const DimensionedField<Type, volMesh>& su,
2370 const fvMatrix<Type>& A
2371)
2372{
2373 checkMethod(A, su, "+");
2374 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2375 tC.ref().source() -= su.mesh().V()*su.field();
2376 return tC;
2377}
2378
2379template<class Type>
2380Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2381(
2382 const tmp<DimensionedField<Type, volMesh>>& tsu,
2383 const fvMatrix<Type>& A
2384)
2385{
2386 checkMethod(A, tsu(), "+");
2387 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2388 tC.ref().source() -= tsu().mesh().V()*tsu().field();
2389 tsu.clear();
2390 return tC;
2391}
2392
2393template<class Type>
2394Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2395(
2396 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2397 const fvMatrix<Type>& A
2398)
2399{
2400 checkMethod(A, tsu(), "+");
2401 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2402 tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2403 tsu.clear();
2404 return tC;
2405}
2406
2407template<class Type>
2408Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2409(
2410 const DimensionedField<Type, volMesh>& su,
2411 const tmp<fvMatrix<Type>>& tA
2412)
2413{
2414 checkMethod(tA(), su, "+");
2415 tmp<fvMatrix<Type>> tC(tA.ptr());
2416 tC.ref().source() -= su.mesh().V()*su.field();
2417 return tC;
2418}
2419
2420template<class Type>
2421Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2422(
2423 const tmp<DimensionedField<Type, volMesh>>& tsu,
2424 const tmp<fvMatrix<Type>>& tA
2425)
2426{
2427 checkMethod(tA(), tsu(), "+");
2428 tmp<fvMatrix<Type>> tC(tA.ptr());
2429 tC.ref().source() -= tsu().mesh().V()*tsu().field();
2430 tsu.clear();
2431 return tC;
2432}
2433
2434template<class Type>
2435Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2436(
2437 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2438 const tmp<fvMatrix<Type>>& tA
2439)
2440{
2441 checkMethod(tA(), tsu(), "+");
2442 tmp<fvMatrix<Type>> tC(tA.ptr());
2443 tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2444 tsu.clear();
2445 return tC;
2446}
2447
2448
2449template<class Type>
2450Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2451(
2452 const fvMatrix<Type>& A,
2453 const fvMatrix<Type>& B
2454)
2455{
2456 checkMethod(A, B, "-");
2457 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2458 tC.ref() -= B;
2459 return tC;
2460}
2461
2462template<class Type>
2463Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2464(
2465 const tmp<fvMatrix<Type>>& tA,
2466 const fvMatrix<Type>& B
2467)
2468{
2469 checkMethod(tA(), B, "-");
2470 tmp<fvMatrix<Type>> tC(tA.ptr());
2471 tC.ref() -= B;
2472 return tC;
2473}
2474
2475template<class Type>
2476Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2477(
2478 const fvMatrix<Type>& A,
2479 const tmp<fvMatrix<Type>>& tB
2480)
2481{
2482 checkMethod(A, tB(), "-");
2483 tmp<fvMatrix<Type>> tC(tB.ptr());
2484 tC.ref() -= A;
2485 tC.ref().negate();
2486 return tC;
2487}
2488
2489template<class Type>
2490Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2491(
2492 const tmp<fvMatrix<Type>>& tA,
2493 const tmp<fvMatrix<Type>>& tB
2494)
2495{
2496 checkMethod(tA(), tB(), "-");
2497 tmp<fvMatrix<Type>> tC(tA.ptr());
2498 tC.ref() -= tB();
2499 tB.clear();
2500 return tC;
2501}
2502
2503template<class Type>
2504Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2505(
2506 const fvMatrix<Type>& A,
2507 const DimensionedField<Type, volMesh>& su
2508)
2509{
2510 checkMethod(A, su, "-");
2511 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2512 tC.ref().source() += su.mesh().V()*su.field();
2513 return tC;
2514}
2515
2516template<class Type>
2517Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2518(
2519 const fvMatrix<Type>& A,
2520 const tmp<DimensionedField<Type, volMesh>>& tsu
2521)
2522{
2523 checkMethod(A, tsu(), "-");
2524 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2525 tC.ref().source() += tsu().mesh().V()*tsu().field();
2526 tsu.clear();
2527 return tC;
2528}
2529
2530template<class Type>
2531Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2532(
2533 const fvMatrix<Type>& A,
2534 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2535)
2536{
2537 checkMethod(A, tsu(), "-");
2538 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2539 tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2540 tsu.clear();
2541 return tC;
2542}
2543
2544template<class Type>
2545Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2546(
2547 const tmp<fvMatrix<Type>>& tA,
2548 const DimensionedField<Type, volMesh>& su
2549)
2550{
2551 checkMethod(tA(), su, "-");
2552 tmp<fvMatrix<Type>> tC(tA.ptr());
2553 tC.ref().source() += su.mesh().V()*su.field();
2554 return tC;
2555}
2556
2557template<class Type>
2558Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2559(
2560 const tmp<fvMatrix<Type>>& tA,
2561 const tmp<DimensionedField<Type, volMesh>>& tsu
2562)
2563{
2564 checkMethod(tA(), tsu(), "-");
2565 tmp<fvMatrix<Type>> tC(tA.ptr());
2566 tC.ref().source() += tsu().mesh().V()*tsu().field();
2567 tsu.clear();
2568 return tC;
2569}
2570
2571template<class Type>
2572Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2573(
2574 const tmp<fvMatrix<Type>>& tA,
2575 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu
2576)
2577{
2578 checkMethod(tA(), tsu(), "-");
2579 tmp<fvMatrix<Type>> tC(tA.ptr());
2580 tC.ref().source() += tsu().mesh().V()*tsu().primitiveField();
2581 tsu.clear();
2582 return tC;
2583}
2584
2585template<class Type>
2586Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2587(
2588 const DimensionedField<Type, volMesh>& su,
2589 const fvMatrix<Type>& A
2590)
2591{
2592 checkMethod(A, su, "-");
2593 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2594 tC.ref().negate();
2595 tC.ref().source() -= su.mesh().V()*su.field();
2596 return tC;
2597}
2598
2599template<class Type>
2600Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2601(
2602 const tmp<DimensionedField<Type, volMesh>>& tsu,
2603 const fvMatrix<Type>& A
2604)
2605{
2606 checkMethod(A, tsu(), "-");
2607 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2608 tC.ref().negate();
2609 tC.ref().source() -= tsu().mesh().V()*tsu().field();
2610 tsu.clear();
2611 return tC;
2612}
2613
2614template<class Type>
2615Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2616(
2617 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2618 const fvMatrix<Type>& A
2619)
2620{
2621 checkMethod(A, tsu(), "-");
2622 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2623 tC.ref().negate();
2624 tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2625 tsu.clear();
2626 return tC;
2627}
2628
2629template<class Type>
2630Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2631(
2632 const DimensionedField<Type, volMesh>& su,
2633 const tmp<fvMatrix<Type>>& tA
2634)
2635{
2636 checkMethod(tA(), su, "-");
2637 tmp<fvMatrix<Type>> tC(tA.ptr());
2638 tC.ref().negate();
2639 tC.ref().source() -= su.mesh().V()*su.field();
2640 return tC;
2641}
2642
2643template<class Type>
2644Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2645(
2646 const tmp<DimensionedField<Type, volMesh>>& tsu,
2647 const tmp<fvMatrix<Type>>& tA
2648)
2649{
2650 checkMethod(tA(), tsu(), "-");
2651 tmp<fvMatrix<Type>> tC(tA.ptr());
2652 tC.ref().negate();
2653 tC.ref().source() -= tsu().mesh().V()*tsu().field();
2654 tsu.clear();
2655 return tC;
2656}
2657
2658template<class Type>
2659Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2660(
2661 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tsu,
2662 const tmp<fvMatrix<Type>>& tA
2663)
2664{
2665 checkMethod(tA(), tsu(), "-");
2666 tmp<fvMatrix<Type>> tC(tA.ptr());
2667 tC.ref().negate();
2668 tC.ref().source() -= tsu().mesh().V()*tsu().primitiveField();
2669 tsu.clear();
2670 return tC;
2671}
2672
2673template<class Type>
2674Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2675(
2676 const fvMatrix<Type>& A,
2677 const dimensioned<Type>& su
2678)
2679{
2680 checkMethod(A, su, "+");
2681 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2682 tC.ref().source() -= su.value()*A.psi().mesh().V();
2683 return tC;
2684}
2685
2686template<class Type>
2687Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2688(
2689 const tmp<fvMatrix<Type>>& tA,
2690 const dimensioned<Type>& su
2691)
2692{
2693 checkMethod(tA(), su, "+");
2694 tmp<fvMatrix<Type>> tC(tA.ptr());
2695 tC.ref().source() -= su.value()*tC().psi().mesh().V();
2696 return tC;
2697}
2698
2699template<class Type>
2700Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2701(
2702 const dimensioned<Type>& su,
2703 const fvMatrix<Type>& A
2704)
2705{
2706 checkMethod(A, su, "+");
2707 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2708 tC.ref().source() -= su.value()*A.psi().mesh().V();
2709 return tC;
2710}
2711
2712template<class Type>
2713Foam::tmp<Foam::fvMatrix<Type>> Foam::operator+
2714(
2715 const dimensioned<Type>& su,
2716 const tmp<fvMatrix<Type>>& tA
2717)
2718{
2719 checkMethod(tA(), su, "+");
2720 tmp<fvMatrix<Type>> tC(tA.ptr());
2721 tC.ref().source() -= su.value()*tC().psi().mesh().V();
2722 return tC;
2723}
2724
2725template<class Type>
2726Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2727(
2728 const fvMatrix<Type>& A,
2729 const dimensioned<Type>& su
2730)
2731{
2732 checkMethod(A, su, "-");
2733 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2734 tC.ref().source() += su.value()*tC().psi().mesh().V();
2735 return tC;
2736}
2737
2738template<class Type>
2739Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2740(
2741 const tmp<fvMatrix<Type>>& tA,
2742 const dimensioned<Type>& su
2743)
2744{
2745 checkMethod(tA(), su, "-");
2746 tmp<fvMatrix<Type>> tC(tA.ptr());
2747 tC.ref().source() += su.value()*tC().psi().mesh().V();
2748 return tC;
2749}
2750
2751template<class Type>
2752Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2753(
2754 const dimensioned<Type>& su,
2755 const fvMatrix<Type>& A
2756)
2757{
2758 checkMethod(A, su, "-");
2759 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2760 tC.ref().negate();
2761 tC.ref().source() -= su.value()*A.psi().mesh().V();
2762 return tC;
2763}
2764
2765template<class Type>
2766Foam::tmp<Foam::fvMatrix<Type>> Foam::operator-
2767(
2768 const dimensioned<Type>& su,
2769 const tmp<fvMatrix<Type>>& tA
2770)
2771{
2772 checkMethod(tA(), su, "-");
2773 tmp<fvMatrix<Type>> tC(tA.ptr());
2774 tC.ref().negate();
2775 tC.ref().source() -= su.value()*tC().psi().mesh().V();
2776 return tC;
2777}
2778
2779
2780template<class Type>
2781Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2782(
2783 const volScalarField::Internal& dsf,
2784 const fvMatrix<Type>& A
2785)
2786{
2787 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2788 tC.ref() *= dsf;
2789 return tC;
2790}
2791
2792template<class Type>
2793Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2794(
2795 const tmp<volScalarField::Internal>& tdsf,
2796 const fvMatrix<Type>& A
2797)
2798{
2799 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2800 tC.ref() *= tdsf;
2801 return tC;
2802}
2803
2804template<class Type>
2805Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2806(
2807 const tmp<volScalarField>& tvsf,
2808 const fvMatrix<Type>& A
2809)
2810{
2811 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2812 tC.ref() *= tvsf;
2813 return tC;
2814}
2815
2816template<class Type>
2817Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2818(
2819 const volScalarField::Internal& dsf,
2820 const tmp<fvMatrix<Type>>& tA
2821)
2822{
2823 tmp<fvMatrix<Type>> tC(tA.ptr());
2824 tC.ref() *= dsf;
2825 return tC;
2826}
2827
2828template<class Type>
2829Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2830(
2831 const tmp<volScalarField::Internal>& tdsf,
2832 const tmp<fvMatrix<Type>>& tA
2833)
2834{
2835 tmp<fvMatrix<Type>> tC(tA.ptr());
2836 tC.ref() *= tdsf;
2837 return tC;
2838}
2839
2840template<class Type>
2841Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2842(
2843 const tmp<volScalarField>& tvsf,
2844 const tmp<fvMatrix<Type>>& tA
2845)
2846{
2847 tmp<fvMatrix<Type>> tC(tA.ptr());
2848 tC.ref() *= tvsf;
2849 return tC;
2850}
2851
2852template<class Type>
2853Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2854(
2855 const dimensioned<scalar>& ds,
2856 const fvMatrix<Type>& A
2857)
2858{
2859 tmp<fvMatrix<Type>> tC(new fvMatrix<Type>(A));
2860 tC.ref() *= ds;
2861 return tC;
2862}
2863
2864template<class Type>
2865Foam::tmp<Foam::fvMatrix<Type>> Foam::operator*
2866(
2867 const dimensioned<scalar>& ds,
2868 const tmp<fvMatrix<Type>>& tA
2869)
2870{
2871 tmp<fvMatrix<Type>> tC(tA.ptr());
2872 tC.ref() *= ds;
2873 return tC;
2874}
2875
2876
2877template<class Type>
2879Foam::operator&
2880(
2881 const fvMatrix<Type>& M,
2882 const DimensionedField<Type, volMesh>& psi
2883)
2884{
2885 tmp<GeometricField<Type, fvPatchField, volMesh>> tMphi
2886 (
2887 new GeometricField<Type, fvPatchField, volMesh>
2888 (
2889 IOobject
2890 (
2891 "M&" + psi.name(),
2892 psi.instance(),
2893 psi.mesh(),
2896 ),
2897 psi.mesh(),
2898 M.dimensions()/dimVol,
2899 extrapolatedCalculatedFvPatchScalarField::typeName
2900 )
2901 );
2902 GeometricField<Type, fvPatchField, volMesh>& Mphi = tMphi.ref();
2903
2904 // Loop over field components
2905 if (M.hasDiag())
2906 {
2907 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
2908 {
2909 scalarField psiCmpt(psi.field().component(cmpt));
2910 scalarField boundaryDiagCmpt(M.diag());
2911 M.addBoundaryDiag(boundaryDiagCmpt, cmpt);
2912 Mphi.primitiveFieldRef().replace(cmpt, -boundaryDiagCmpt*psiCmpt);
2913 }
2914 }
2915 else
2916 {
2917 Mphi.primitiveFieldRef() = Zero;
2918 }
2919
2920 Mphi.primitiveFieldRef() += M.lduMatrix::H(psi.field()) + M.source();
2921 M.addBoundarySource(Mphi.primitiveFieldRef());
2922
2923 Mphi.primitiveFieldRef() /= -psi.mesh().V();
2924 Mphi.correctBoundaryConditions();
2925
2926 return tMphi;
2927}
2928
2929template<class Type>
2931Foam::operator&
2932(
2933 const fvMatrix<Type>& M,
2934 const tmp<DimensionedField<Type, volMesh>>& tpsi
2935)
2936{
2937 tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = M & tpsi();
2938 tpsi.clear();
2939 return tMpsi;
2940}
2941
2942template<class Type>
2944Foam::operator&
2945(
2946 const fvMatrix<Type>& M,
2947 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tpsi
2948)
2949{
2950 tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = M & tpsi();
2951 tpsi.clear();
2952 return tMpsi;
2953}
2954
2955template<class Type>
2957Foam::operator&
2958(
2959 const tmp<fvMatrix<Type>>& tM,
2960 const DimensionedField<Type, volMesh>& psi
2961)
2962{
2963 tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & psi;
2964 tM.clear();
2965 return tMpsi;
2966}
2967
2968template<class Type>
2970Foam::operator&
2971(
2972 const tmp<fvMatrix<Type>>& tM,
2973 const tmp<DimensionedField<Type, volMesh>>& tpsi
2974)
2975{
2976 tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2977 tM.clear();
2978 tpsi.clear();
2979 return tMpsi;
2980}
2981
2982template<class Type>
2984Foam::operator&
2985(
2986 const tmp<fvMatrix<Type>>& tM,
2987 const tmp<GeometricField<Type, fvPatchField, volMesh>>& tpsi
2988)
2989{
2990 tmp<GeometricField<Type, fvPatchField, volMesh>> tMpsi = tM() & tpsi();
2991 tM.clear();
2992 tpsi.clear();
2993 return tMpsi;
2994}
2995
2996
2997// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2998
2999template<class Type>
3000Foam::Ostream& Foam::operator<<(Ostream& os, const fvMatrix<Type>& fvm)
3001{
3002 os << static_cast<const lduMatrix&>(fvm) << nl
3003 << fvm.dimensions_ << nl
3004 << fvm.source_ << nl
3005 << fvm.internalCoeffs_ << nl
3006 << fvm.boundaryCoeffs_ << endl;
3007
3009
3010 return os;
3011}
3012
3013
3014// * * * * * * * * * * * * * * * * Solvers * * * * * * * * * * * * * * * * * //
3015
3016#include "fvMatrixSolve.C"
3017
3018// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
#define M(I)
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 dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
void setOriented(const bool oriented=true) noexcept
Set the oriented flag.
const Field< Type > & field() const
Return field.
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, volMesh > Internal
The internal field type from which this GeometricField is derived.
void replace(const direction d, const GeometricField< cmptType, PatchField, GeoMesh > &gcf)
Replace specified field component with content from another field.
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
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:330
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A dynamically resizable PtrList with allocation management.
Definition: PtrDynList.H:64
void append(T *ptr)
Append an element to the end of the list.
Definition: PtrDynListI.H:235
void setSize(const label newLen)
Same as resize()
Definition: PtrList.H:151
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
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
const T * set(const label i) const
Definition: UPtrList.H:248
void setSize(const label n)
Alias for resize()
Definition: UPtrList.H:261
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
T & last()
Return reference to the last element of the list.
Definition: UPtrListI.H:176
A single value that is represented as a list with an operator[] to access the value....
Definition: UniformList.H:55
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
This boundary condition enforces a cyclic condition between a pair of boundaries, whereby communicati...
This boundary condition enforces a cyclic condition between a pair of boundaries.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
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 Type & value() const
Return const reference to value.
const word & name() const
Return const reference to name.
tmp< GeometricField< Type, faPatchField, areaMesh > > H() const
Return the H operation source.
Definition: faMatrix.C:633
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:38
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void operator-=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1672
tmp< Field< Type > > DD() const
Return the matrix Type diagonal.
Definition: fvMatrix.C:1279
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:1443
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:1303
void manipulateMatrix(direction cmp)
Manipulate matrix.
Definition: fvMatrix.C:790
void addToInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Add patch contribution to internal field.
Definition: fvMatrix.C:49
void mapContributions(label fieldi, const FieldField< Field, Type > &fluxContrib, FieldField< Field, Type > &contrib, bool internal) const
Add internal and boundary contribution to local patches.
Definition: fvMatrix.C:548
void setReferences(const labelUList &cellLabels, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1018
void transferFvMatrixCoeffs()
Transfer lower, upper, diag and source to this fvMatrix.
Definition: fvMatrix.C:813
void relax()
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1240
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:1394
void operator+=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1634
void createOrUpdateLduPrimitiveAssembly()
Create or update ldu assembly.
Definition: fvMatrix.C:906
virtual ~fvMatrix()
Destructor.
Definition: fvMatrix.C:456
void setValues(const labelUList &cellLabels, const Type &value)
Definition: fvMatrix.C:969
const dimensionSet & dimensions() const noexcept
Definition: fvMatrix.H:453
void addCmptAvBoundaryDiag(scalarField &diag) const
Definition: fvMatrix.C:152
void addBoundarySource(Field< Type > &source, const bool couples=true) const
Definition: fvMatrix.C:177
void setReference(const label celli, const Type &value, const bool forceReference=false)
Set reference level for solution.
Definition: fvMatrix.C:1002
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:1257
void setInterfaces(lduInterfaceFieldPtrsList &, PtrDynList< lduInterfaceField > &newInterfaces)
Set interfaces.
Definition: fvMatrix.C:470
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:1332
void addBoundaryDiag(scalarField &diag, const direction cmpt) const
Definition: fvMatrix.C:124
void negate()
Inplace negate.
Definition: fvMatrix.C:1619
tmp< fvMatrix< Type > > clone() const
Construct and return a clone.
Definition: fvMatrix.H:323
void addFvMatrix(fvMatrix< Type > &matrix)
Add fvMatrix.
Definition: fvMatrix.C:1063
tmp< scalarField > D() const
Return the matrix scalar diagonal.
Definition: fvMatrix.C:1270
void subtractFromInternalField(const labelUList &addr, const Field< Type2 > &pf, Field< Type2 > &intf) const
Subtract patch contribution from internal field.
Definition: fvMatrix.C:87
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:412
bool checkImplicit(const label fieldi=0)
Name the implicit assembly addressing.
Definition: fvMatrix.C:307
lduPrimitiveMeshAssembly * lduMeshPtr()
Access to lduPrimitiveMeshAssembly.
Definition: fvMatrix.C:880
void operator=(const fvMatrix< Type > &)
Definition: fvMatrix.C:1574
void setValuesFromList(const labelUList &cellLabels, const ListType< Type > &values)
Set solution in given cells to the specified values.
Definition: fvMatrix.C:227
const dictionary & solverDict() const
Return the solver dictionary taking into account finalIteration.
Definition: fvMatrix.C:1558
void setBounAndInterCoeffs()
Manipulate boundary/internal coeffs for coupling.
Definition: fvMatrix.C:660
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:453
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:350
void clearOut()
Clear additional addressing.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition: lduMatrix.H:84
tmp< scalarField > H1() const
void operator=(const lduMatrix &)
tmp< Field< Type > > faceH(const Field< Type > &) const
void operator*=(const scalarField &)
void operator+=(const lduMatrix &)
void operator-=(const lduMatrix &)
virtual const objectRegistry & thisDb() const
Return the object registry.
Definition: lduMesh.C:42
An assembly of lduMatrix that is specific inter-region coupling through mapped patches.
const labelListListList & cellBoundMap() const
Return patch local sub-face to nbrCellId map.
const labelListList & patchMap() const
Return patchMap.
const labelListListList & facePatchFaceMap() const
Return patch local sub-face to local patch face map.
const labelListList & patchLocalToGlobalMap() const
Return patchLocalToGlobalMap.
void update(UPtrList< GeometricField< Type, fvPatchField, volMesh > > &psis)
Update mappings.
static const lduMesh & mesh(const lduMesh &mesh0, const PtrList< lduPrimitiveMesh > &otherMeshes, const label meshI)
Select either mesh0 (meshI is 0) or otherMeshes[meshI-1].
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
virtual label neighbPolyPatchID() const
Return nbr patchID.
Definition: polyPatch.H:332
virtual bool masterImplicit() const
Return implicit master.
Definition: polyPatch.H:346
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 for handling words, derived from Foam::string.
Definition: word.H:68
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
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
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
OBJstream os(runTime.globalPath()/outputName)
const cellShapeList & cells
label cellId
label faceId(-1)
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
#define FUNCTION_NAME
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
void checkMethod(const faMatrix< Type > &, const faMatrix< Type > &, const char *)
Definition: faMatrix.C:1050
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
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)
void cmptMin(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
messageStream Info
Information stream (stdout output on master, null elsewhere)
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)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
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)
const dimensionSet dimVol(dimVolume)
Older spelling for dimVolume.
Definition: dimensionSets.H:64
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.
volScalarField & alpha
const dimensionedScalar & D
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59
Foam::surfaceFields.
conserve primitiveFieldRef()+