fvScalarMatrix.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "fvScalarMatrix.H"
27 #include "Residuals.H"
29 
30 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
31 
32 template<>
34 (
35  const label patchi,
36  const label facei,
37  const direction,
38  const scalar value
39 )
40 {
41  if (psi_.needReference())
42  {
43  if (Pstream::master())
44  {
45  internalCoeffs_[patchi][facei] +=
46  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
47 
48  boundaryCoeffs_[patchi][facei] +=
49  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
50  *value;
51  }
52  }
53 }
54 
55 
56 template<>
59 (
60  const dictionary& solverControls
61 )
62 {
63  if (debug)
64  {
65  Info.masterStream(this->mesh().comm())
66  << "fvMatrix<scalar>::solver(const dictionary& solverControls) : "
67  "solver for fvMatrix<scalar>"
68  << endl;
69  }
70 
71  scalarField saveDiag(diag());
72  addBoundaryDiag(diag(), 0);
73 
75  (
77  (
78  *this,
80  (
81  psi_.name(),
82  *this,
83  boundaryCoeffs_,
84  internalCoeffs_,
85  psi_.boundaryField().scalarInterfaces(),
86  solverControls
87  )
88  )
89  );
90 
91  diag() = saveDiag;
92 
93  return solverPtr;
94 }
95 
96 
97 template<>
99 (
100  const dictionary& solverControls
101 )
102 {
105  (fvMat_.psi());
106 
107  scalarField saveDiag(fvMat_.diag());
108  fvMat_.addBoundaryDiag(fvMat_.diag(), 0);
109 
110  scalarField totalSource(fvMat_.source());
111  fvMat_.addBoundarySource(totalSource, false);
112 
113  // Assign new solver controls
114  solver_->read(solverControls);
115 
116  solverPerformance solverPerf = solver_->solve
117  (
118  psi.primitiveFieldRef(),
119  totalSource
120  );
121 
122  if (solverPerformance::debug)
123  {
124  solverPerf.print(Info.masterStream(fvMat_.mesh().comm()));
125  }
126 
127  fvMat_.diag() = saveDiag;
128 
130 
131  Residuals<scalar>::append(psi.mesh(), solverPerf);
132 
133  return solverPerf;
134 }
135 
136 
137 template<>
139 (
140  const dictionary& solverControls
141 )
142 {
143  if (debug)
144  {
145  Info.masterStream(this->mesh().comm())
146  << "fvMatrix<scalar>::solveSegregated"
147  "(const dictionary& solverControls) : "
148  "solving fvMatrix<scalar>"
149  << endl;
150  }
151 
154 
155  scalarField saveDiag(diag());
156  addBoundaryDiag(diag(), 0);
157 
158  scalarField totalSource(source_);
159  addBoundarySource(totalSource, false);
160 
161  // Solver call
163  (
164  psi.name(),
165  *this,
166  boundaryCoeffs_,
167  internalCoeffs_,
168  psi_.boundaryField().scalarInterfaces(),
169  solverControls
170  )->solve(psi.primitiveFieldRef(), totalSource);
171 
172  if (solverPerformance::debug)
173  {
174  solverPerf.print(Info.masterStream(mesh().comm()));
175  }
176 
177  diag() = saveDiag;
178 
180 
181  Residuals<scalar>::append(psi.mesh(), solverPerf);
182 
183  return solverPerf;
184 }
185 
186 
187 template<>
189 {
190  scalarField boundaryDiag(psi_.size(), 0.0);
191  addBoundaryDiag(boundaryDiag, 0);
192 
193  tmp<scalarField> tres
194  (
195  lduMatrix::residual
196  (
197  psi_.primitiveField(),
198  source_ - boundaryDiag*psi_.primitiveField(),
199  boundaryCoeffs_,
200  psi_.boundaryField().scalarInterfaces(),
201  0
202  )
203  );
204 
205  addBoundarySource(tres.ref());
206 
207  return tres;
208 }
209 
210 
211 template<>
213 {
214  tmp<volScalarField> tHphi
215  (
217  (
218  "H("+psi_.name()+')',
219  psi_.mesh(),
220  dimensions_/dimVol,
221  extrapolatedCalculatedFvPatchScalarField::typeName
222  )
223  );
224  volScalarField& Hphi = tHphi.ref();
225 
226  Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
227  addBoundarySource(Hphi.primitiveFieldRef());
228 
229  Hphi.primitiveFieldRef() /= psi_.mesh().V();
231 
232  return tHphi;
233 }
234 
235 
236 template<>
238 {
240  (
242  (
243  "H(1)",
244  psi_.mesh(),
245  dimensions_/(dimVol*psi_.dimensions()),
246  extrapolatedCalculatedFvPatchScalarField::typeName
247  )
248  );
249  volScalarField& H1_ = tH1.ref();
250 
251  H1_.primitiveFieldRef() = lduMatrix::H1();
252  // addBoundarySource(Hphi.primitiveField());
253 
254  H1_.primitiveFieldRef() /= psi_.mesh().V();
255  H1_.correctBoundaryConditions();
256 
257  return tH1;
258 }
259 
260 
261 // ************************************************************************* //
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:760
Solver class returned by the solver function.
Definition: fvMatrix.H:220
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:815
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dynamicFvMesh & mesh
SolverPerformance< Type > solve()
Solve returning the solution statistics.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
const dimensionSet dimVol
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
Internal & ref()
Return a reference to the dimensioned internal field.
void print(Ostream &os) const
Print summary of solver performance to the given stream.
OSstream & masterStream(const label communicator)
Convert to OSstream.
Definition: messageStream.C:59
label patchi
rhoEqn solve()
A scalar instance of fvMatrix.
autoPtr< fvSolver > solver()
Construct and return the solver.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const volScalarField & psi
void setComponentReference(const label patchi, const label facei, const direction cmpt, const scalar value)
Set reference level for a component of the solution.
Definition: fvMatrixSolve.C:34
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< Field< Type > > residual() const
Return the matrix residual.
MeshObject to store the solver performance residuals of all the fields of the type it is instantiated...
Definition: Residuals.H:52