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