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-2018 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"
28 
29 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30 
31 template<>
33 (
34  const label patchi,
35  const label facei,
36  const direction,
37  const scalar value
38 )
39 {
40  if (psi_.needReference())
41  {
42  if (Pstream::master())
43  {
44  internalCoeffs_[patchi][facei] +=
45  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]];
46 
47  boundaryCoeffs_[patchi][facei] +=
48  diag()[psi_.mesh().boundary()[patchi].faceCells()[facei]]
49  *value;
50  }
51  }
52 }
53 
54 
55 template<>
58 (
59  const dictionary& solverControls
60 )
61 {
62  if (debug)
63  {
64  Info.masterStream(this->mesh().comm())
65  << "fvMatrix<scalar>::solver(const dictionary& solverControls) : "
66  "solver for fvMatrix<scalar>"
67  << endl;
68  }
69 
70  scalarField saveDiag(diag());
71  addBoundaryDiag(diag(), 0);
72 
74  (
76  (
77  *this,
79  (
80  psi_.name(),
81  *this,
82  boundaryCoeffs_,
83  internalCoeffs_,
84  psi_.boundaryField().scalarInterfaces(),
85  solverControls
86  )
87  )
88  );
89 
90  diag() = saveDiag;
91 
92  return solverPtr;
93 }
94 
95 
96 template<>
98 (
99  const dictionary& solverControls
100 )
101 {
104  (fvMat_.psi());
105 
106  scalarField saveDiag(fvMat_.diag());
107  fvMat_.addBoundaryDiag(fvMat_.diag(), 0);
108 
109  scalarField totalSource(fvMat_.source());
110  fvMat_.addBoundarySource(totalSource, false);
111 
112  // Assign new solver controls
113  solver_->read(solverControls);
114 
115  solverPerformance solverPerf = solver_->solve
116  (
117  psi.primitiveFieldRef(),
118  totalSource
119  );
120 
121  if (solverPerformance::debug)
122  {
123  solverPerf.print(Info.masterStream(fvMat_.mesh().comm()));
124  }
125 
126  fvMat_.diag() = saveDiag;
127 
129 
130  psi.mesh().setSolverPerformance(psi.name(), solverPerf);
131 
132  return solverPerf;
133 }
134 
135 
136 template<>
138 (
139  const dictionary& solverControls
140 )
141 {
142  if (debug)
143  {
144  Info.masterStream(this->mesh().comm())
145  << "fvMatrix<scalar>::solveSegregated"
146  "(const dictionary& solverControls) : "
147  "solving fvMatrix<scalar>"
148  << endl;
149  }
150 
153 
154  scalarField saveDiag(diag());
155  addBoundaryDiag(diag(), 0);
156 
157  scalarField totalSource(source_);
158  addBoundarySource(totalSource, false);
159 
160  // Solver call
162  (
163  psi.name(),
164  *this,
165  boundaryCoeffs_,
166  internalCoeffs_,
167  psi_.boundaryField().scalarInterfaces(),
168  solverControls
169  )->solve(psi.primitiveFieldRef(), totalSource);
170 
171  if (solverPerformance::debug)
172  {
173  solverPerf.print(Info.masterStream(mesh().comm()));
174  }
175 
176  diag() = saveDiag;
177 
179 
180  psi.mesh().setSolverPerformance(psi.name(), solverPerf);
181 
182  return solverPerf;
183 }
184 
185 
186 template<>
188 {
189  scalarField boundaryDiag(psi_.size(), 0.0);
190  addBoundaryDiag(boundaryDiag, 0);
191 
192  tmp<scalarField> tres
193  (
194  lduMatrix::residual
195  (
196  psi_.primitiveField(),
197  source_ - boundaryDiag*psi_.primitiveField(),
198  boundaryCoeffs_,
199  psi_.boundaryField().scalarInterfaces(),
200  0
201  )
202  );
203 
204  addBoundarySource(tres.ref());
205 
206  return tres;
207 }
208 
209 
210 template<>
212 {
213  tmp<volScalarField> tHphi
214  (
215  new volScalarField
216  (
217  IOobject
218  (
219  "H("+psi_.name()+')',
220  psi_.instance(),
221  psi_.mesh(),
222  IOobject::NO_READ,
223  IOobject::NO_WRITE
224  ),
225  psi_.mesh(),
226  dimensions_/dimVol,
227  extrapolatedCalculatedFvPatchScalarField::typeName
228  )
229  );
230  volScalarField& Hphi = tHphi.ref();
231 
232  Hphi.primitiveFieldRef() = (lduMatrix::H(psi_.primitiveField()) + source_);
233  addBoundarySource(Hphi.primitiveFieldRef());
234 
235  Hphi.primitiveFieldRef() /= psi_.mesh().V();
237 
238  return tHphi;
239 }
240 
241 
242 template<>
244 {
246  (
247  new volScalarField
248  (
249  IOobject
250  (
251  "H(1)",
252  psi_.instance(),
253  psi_.mesh(),
254  IOobject::NO_READ,
255  IOobject::NO_WRITE
256  ),
257  psi_.mesh(),
258  dimensions_/(dimVol*psi_.dimensions()),
259  extrapolatedCalculatedFvPatchScalarField::typeName
260  )
261  );
262  volScalarField& H1_ = tH1.ref();
263 
264  H1_.primitiveFieldRef() = lduMatrix::H1();
265  // addBoundarySource(Hphi.primitiveField());
266 
267  H1_.primitiveFieldRef() /= psi_.mesh().V();
268  H1_.correctBoundaryConditions();
269 
270  return tH1;
271 }
272 
273 
274 // ************************************************************************* //
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:297
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
uint8_t direction
Definition: direction.H:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
tmp< GeometricField< Type, fvPatchField, volMesh > > H() const
Return the H operation source.
Definition: fvMatrix.C:766
Solver class returned by the solver function.
Definition: fvMatrix.H:220
const dimensionSet dimVol(dimVolume)
Definition: dimensionSets.H:59
tmp< volScalarField > H1() const
Return H(1)
Definition: fvMatrix.C:828
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dynamicFvMesh & mesh
SolverPerformance< Type > solve()
Solve returning the solution statistics.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
rhoEqn solve()
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
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:33
SolverPerformance< Type > solveSegregated(const dictionary &)
Solve segregated returning the solution statistics.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
tmp< Field< Type > > residual() const
Return the matrix residual.