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