lduMatrix.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 "lduMatrix.H"
27 #include "IOstreams.H"
28 #include "Switch.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(lduMatrix, 1);
35 }
36 
37 
38 const Foam::label Foam::lduMatrix::solver::defaultMaxIter_ = 1000;
39 
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
44 :
45  lduMesh_(mesh),
46  lowerPtr_(nullptr),
47  diagPtr_(nullptr),
48  upperPtr_(nullptr)
49 {}
50 
51 
53 :
54  lduMesh_(A.lduMesh_),
55  lowerPtr_(nullptr),
56  diagPtr_(nullptr),
57  upperPtr_(nullptr)
58 {
59  if (A.lowerPtr_)
60  {
61  lowerPtr_ = new scalarField(*(A.lowerPtr_));
62  }
63 
64  if (A.diagPtr_)
65  {
66  diagPtr_ = new scalarField(*(A.diagPtr_));
67  }
68 
69  if (A.upperPtr_)
70  {
71  upperPtr_ = new scalarField(*(A.upperPtr_));
72  }
73 }
74 
75 
77 :
78  lduMesh_(A.lduMesh_),
79  lowerPtr_(nullptr),
80  diagPtr_(nullptr),
81  upperPtr_(nullptr)
82 {
83  if (reuse)
84  {
85  if (A.lowerPtr_)
86  {
87  lowerPtr_ = A.lowerPtr_;
88  A.lowerPtr_ = nullptr;
89  }
90 
91  if (A.diagPtr_)
92  {
93  diagPtr_ = A.diagPtr_;
94  A.diagPtr_ = nullptr;
95  }
96 
97  if (A.upperPtr_)
98  {
99  upperPtr_ = A.upperPtr_;
100  A.upperPtr_ = nullptr;
101  }
102  }
103  else
104  {
105  if (A.lowerPtr_)
106  {
107  lowerPtr_ = new scalarField(*(A.lowerPtr_));
108  }
109 
110  if (A.diagPtr_)
111  {
112  diagPtr_ = new scalarField(*(A.diagPtr_));
113  }
114 
115  if (A.upperPtr_)
116  {
117  upperPtr_ = new scalarField(*(A.upperPtr_));
118  }
119  }
120 }
121 
122 
124 :
125  lduMesh_(mesh),
126  lowerPtr_(nullptr),
127  diagPtr_(nullptr),
128  upperPtr_(nullptr)
129 {
130  Switch hasLow(is);
131  Switch hasDiag(is);
132  Switch hasUp(is);
133 
134  if (hasLow)
135  {
136  lowerPtr_ = new scalarField(is);
137  }
138  if (hasDiag)
139  {
140  diagPtr_ = new scalarField(is);
141  }
142  if (hasUp)
143  {
144  upperPtr_ = new scalarField(is);
145  }
146 }
147 
148 
150 {
151  if (lowerPtr_)
152  {
153  delete lowerPtr_;
154  }
155 
156  if (diagPtr_)
157  {
158  delete diagPtr_;
159  }
160 
161  if (upperPtr_)
162  {
163  delete upperPtr_;
164  }
165 }
166 
167 
169 {
170  if (!lowerPtr_)
171  {
172  if (upperPtr_)
173  {
174  lowerPtr_ = new scalarField(*upperPtr_);
175  }
176  else
177  {
178  lowerPtr_ = new scalarField(lduAddr().lowerAddr().size(), 0.0);
179  }
180  }
181 
182  return *lowerPtr_;
183 }
184 
185 
187 {
188  if (!diagPtr_)
189  {
190  diagPtr_ = new scalarField(lduAddr().size(), 0.0);
191  }
192 
193  return *diagPtr_;
194 }
195 
196 
198 {
199  if (!upperPtr_)
200  {
201  if (lowerPtr_)
202  {
203  upperPtr_ = new scalarField(*lowerPtr_);
204  }
205  else
206  {
207  upperPtr_ = new scalarField(lduAddr().lowerAddr().size(), 0.0);
208  }
209  }
210 
211  return *upperPtr_;
212 }
213 
214 
216 {
217  if (!lowerPtr_)
218  {
219  if (upperPtr_)
220  {
221  lowerPtr_ = new scalarField(*upperPtr_);
222  }
223  else
224  {
225  lowerPtr_ = new scalarField(nCoeffs, 0.0);
226  }
227  }
228 
229  return *lowerPtr_;
230 }
231 
232 
234 {
235  if (!diagPtr_)
236  {
237  diagPtr_ = new scalarField(size, 0.0);
238  }
239 
240  return *diagPtr_;
241 }
242 
243 
245 {
246  if (!upperPtr_)
247  {
248  if (lowerPtr_)
249  {
250  upperPtr_ = new scalarField(*lowerPtr_);
251  }
252  else
253  {
254  upperPtr_ = new scalarField(nCoeffs, 0.0);
255  }
256  }
257 
258  return *upperPtr_;
259 }
260 
261 
263 {
264  if (!lowerPtr_ && !upperPtr_)
265  {
267  << "lowerPtr_ or upperPtr_ unallocated"
268  << abort(FatalError);
269  }
270 
271  if (lowerPtr_)
272  {
273  return *lowerPtr_;
274  }
275  else
276  {
277  return *upperPtr_;
278  }
279 }
280 
281 
283 {
284  if (!diagPtr_)
285  {
287  << "diagPtr_ unallocated"
288  << abort(FatalError);
289  }
290 
291  return *diagPtr_;
292 }
293 
294 
296 {
297  if (!lowerPtr_ && !upperPtr_)
298  {
300  << "lowerPtr_ or upperPtr_ unallocated"
301  << abort(FatalError);
302  }
303 
304  if (upperPtr_)
305  {
306  return *upperPtr_;
307  }
308  else
309  {
310  return *lowerPtr_;
311  }
312 }
313 
314 
315 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
316 
318 {
319  Switch hasLow = ldum.hasLower();
320  Switch hasDiag = ldum.hasDiag();
321  Switch hasUp = ldum.hasUpper();
322 
323  os << hasLow << token::SPACE << hasDiag << token::SPACE
324  << hasUp << token::SPACE;
325 
326  if (hasLow)
327  {
328  os << ldum.lower();
329  }
330 
331  if (hasDiag)
332  {
333  os << ldum.diag();
334  }
335 
336  if (hasUp)
337  {
338  os << ldum.upper();
339  }
340 
341  os.check("Ostream& operator<<(Ostream&, const lduMatrix&");
342 
343  return os;
344 }
345 
346 
347 Foam::Ostream& Foam::operator<<(Ostream& os, const InfoProxy<lduMatrix>& ip)
348 {
349  const lduMatrix& ldum = ip.t_;
350 
351  Switch hasLow = ldum.hasLower();
352  Switch hasDiag = ldum.hasDiag();
353  Switch hasUp = ldum.hasUpper();
354 
355  os << "Lower:" << hasLow
356  << " Diag:" << hasDiag
357  << " Upper:" << hasUp << endl;
358 
359  if (hasLow)
360  {
361  os << "lower:" << ldum.lower().size() << endl;
362  }
363  if (hasDiag)
364  {
365  os << "diag :" << ldum.diag().size() << endl;
366  }
367  if (hasUp)
368  {
369  os << "upper:" << ldum.upper().size() << endl;
370  }
371 
372 
373  // if (hasLow)
374  //{
375  // os << "lower contents:" << endl;
376  // forAll(ldum.lower(), i)
377  // {
378  // os << "i:" << i << "\t" << ldum.lower()[i] << endl;
379  // }
380  // os << endl;
381  //}
382  // if (hasDiag)
383  //{
384  // os << "diag contents:" << endl;
385  // forAll(ldum.diag(), i)
386  // {
387  // os << "i:" << i << "\t" << ldum.diag()[i] << endl;
388  // }
389  // os << endl;
390  //}
391  // if (hasUp)
392  //{
393  // os << "upper contents:" << endl;
394  // forAll(ldum.upper(), i)
395  // {
396  // os << "i:" << i << "\t" << ldum.upper()[i] << endl;
397  // }
398  // os << endl;
399  //}
400 
401  os.check("Ostream& operator<<(Ostream&, const lduMatrix&");
402 
403  return os;
404 }
405 
406 
407 // ************************************************************************* //
bool hasUpper() const
Definition: lduMatrix.H:585
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
lduMatrix(const lduMesh &)
Construct given an LDU addressed mesh.
Definition: lduMatrix.C:43
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:59
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
bool hasLower() const
Definition: lduMatrix.H:590
scalarField & upper()
Definition: lduMatrix.C:197
fvMesh & mesh
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
bool hasDiag() const
Definition: lduMatrix.H:580
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
defineTypeNameAndDebug(combustionModel, 0)
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:79
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:550
scalarField & lower()
Definition: lduMatrix.C:168
Ostream & operator<<(Ostream &, const ensightPart &)
scalarField & diag()
Definition: lduMatrix.C:186
const lduMesh & mesh() const
Return the LDU mesh from which the addressing is obtained.
Definition: lduMatrix.H:544
Namespace for OpenFOAM.
~lduMatrix()
Destructor.
Definition: lduMatrix.C:149