lduMatrixOperations.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-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2024 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 Description
28  lduMatrix member operations.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "lduMatrix.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
37 {
38  const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
39  const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
40  scalarField& Diag = diag();
41 
42  const labelUList& l = lduAddr().lowerAddr();
43  const labelUList& u = lduAddr().upperAddr();
44 
45  for (label face=0; face<l.size(); face++)
46  {
47  Diag[l[face]] += Lower[face];
48  Diag[u[face]] += Upper[face];
49  }
50 }
51 
52 
54 {
55  const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
56  const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
57  scalarField& Diag = diag();
58 
59  const labelUList& l = lduAddr().lowerAddr();
60  const labelUList& u = lduAddr().upperAddr();
61 
62  for (label face=0; face<l.size(); face++)
63  {
64  Diag[l[face]] -= Lower[face];
65  Diag[u[face]] -= Upper[face];
66  }
67 }
68 
69 
71 (
72  scalarField& sumOff
73 ) const
74 {
75  const scalarField& Lower = const_cast<const lduMatrix&>(*this).lower();
76  const scalarField& Upper = const_cast<const lduMatrix&>(*this).upper();
77 
78  const labelUList& l = lduAddr().lowerAddr();
79  const labelUList& u = lduAddr().upperAddr();
80 
81  for (label face = 0; face < l.size(); face++)
82  {
83  sumOff[u[face]] += mag(Lower[face]);
84  sumOff[l[face]] += mag(Upper[face]);
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
90 
91 void Foam::lduMatrix::operator=(const lduMatrix& A)
92 {
93  if (this == &A)
94  {
95  return; // Self-assignment is a no-op
96  }
97 
98  if (A.hasLower())
99  {
100  lower() = A.lower();
101  }
102  else
103  {
104  lowerPtr_.reset(nullptr);
105  }
106 
107  if (A.hasUpper())
108  {
109  upper() = A.upper();
110  }
111  else
112  {
113  upperPtr_.reset(nullptr);
114  }
115 
116  if (A.hasDiag())
117  {
118  diag() = A.diag();
119  }
120 }
121 
122 
123 void Foam::lduMatrix::operator=(lduMatrix&& A)
124 {
125  if (this == &A)
126  {
127  return; // Self-assignment is a no-op
128  }
129 
130  diagPtr_ = std::move(A.diagPtr_);
131  upperPtr_ = std::move(A.upperPtr_);
132  lowerPtr_ = std::move(A.lowerPtr_);
133 }
134 
135 
137 {
138  if (diagPtr_)
139  {
140  diagPtr_->negate();
141  }
142 
143  if (upperPtr_)
144  {
145  upperPtr_->negate();
146  }
147 
148  if (lowerPtr_)
149  {
150  lowerPtr_->negate();
151  }
152 }
153 
154 
155 void Foam::lduMatrix::operator+=(const lduMatrix& A)
156 {
157  if (A.hasDiag())
158  {
159  diag() += A.diag();
160  }
161 
162  if (symmetric() && A.symmetric())
163  {
164  upper() += A.upper();
165  }
166  else if (symmetric() && A.asymmetric())
167  {
168  if (upperPtr_)
169  {
170  lower();
171  }
172  else
173  {
174  upper();
175  }
176 
177  upper() += A.upper();
178  lower() += A.lower();
179  }
180  else if (asymmetric() && A.symmetric())
181  {
182  if (A.hasUpper())
183  {
184  lower() += A.upper();
185  upper() += A.upper();
186  }
187  else
188  {
189  lower() += A.lower();
190  upper() += A.lower();
191  }
192 
193  }
194  else if (asymmetric() && A.asymmetric())
195  {
196  lower() += A.lower();
197  upper() += A.upper();
198  }
199  else if (diagonal())
200  {
201  if (A.hasUpper())
202  {
203  upper() = A.upper();
204  }
205 
206  if (A.hasLower())
207  {
208  lower() = A.lower();
209  }
210  }
211  else if (A.diagonal())
212  {
213  }
214  else
215  {
216  if (debug > 1)
217  {
219  << "Unknown matrix type combination" << nl
220  << " this : " << this->matrixTypeName()
221  << " A : " << A.matrixTypeName() << endl;
222  }
223  }
224 }
225 
226 
227 void Foam::lduMatrix::operator-=(const lduMatrix& A)
228 {
229  if (A.diagPtr_)
230  {
231  diag() -= A.diag();
232  }
233 
234  if (symmetric() && A.symmetric())
235  {
236  upper() -= A.upper();
237  }
238  else if (symmetric() && A.asymmetric())
239  {
240  if (upperPtr_)
241  {
242  lower();
243  }
244  else
245  {
246  upper();
247  }
248 
249  upper() -= A.upper();
250  lower() -= A.lower();
251  }
252  else if (asymmetric() && A.symmetric())
253  {
254  if (A.hasUpper())
255  {
256  lower() -= A.upper();
257  upper() -= A.upper();
258  }
259  else
260  {
261  lower() -= A.lower();
262  upper() -= A.lower();
263  }
264 
265  }
266  else if (asymmetric() && A.asymmetric())
267  {
268  lower() -= A.lower();
269  upper() -= A.upper();
270  }
271  else if (diagonal())
272  {
273  if (A.hasUpper())
274  {
275  upper() = -A.upper();
276  }
277 
278  if (A.hasLower())
279  {
280  lower() = -A.lower();
281  }
282  }
283  else if (A.diagonal())
284  {
285  }
286  else
287  {
288  if (debug > 1)
289  {
291  << "Unknown matrix type combination" << nl
292  << " this : " << this->matrixTypeName()
293  << " A : " << A.matrixTypeName() << endl;
294  }
295  }
296 }
297 
298 
300 {
301  if (diagPtr_)
302  {
303  *diagPtr_ *= sf;
304  }
305 
306  // Non-uniform scaling causes a symmetric matrix
307  // to become asymmetric
308  if (symmetric() || asymmetric())
309  {
310  scalarField& upper = this->upper();
311  scalarField& lower = this->lower();
312 
313  const labelUList& l = lduAddr().lowerAddr();
314  const labelUList& u = lduAddr().upperAddr();
315 
316  for (label face=0; face<upper.size(); face++)
317  {
318  upper[face] *= sf[l[face]];
319  }
320 
321  for (label face=0; face<lower.size(); face++)
322  {
323  lower[face] *= sf[u[face]];
324  }
325  }
326 }
327 
328 
329 void Foam::lduMatrix::operator*=(scalar s)
330 {
331  if (diagPtr_)
332  {
333  *diagPtr_ *= s;
334  }
335 
336  if (upperPtr_)
337  {
338  *upperPtr_ *= s;
339  }
340 
341  if (lowerPtr_)
342  {
343  *lowerPtr_ *= s;
344  }
345 }
346 
347 
348 // ************************************************************************* //
const scalarField & diag() const
Definition: lduMatrix.C:163
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1187
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
void sumMagOffDiag(scalarField &sumOff) const
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
const scalarField & lower() const
Definition: lduMatrix.C:268
void operator=(const lduMatrix &)
Copy assignment.
void operator*=(const scalarField &)
int debug
Static debugging option.
lduMatrix is a general matrix class in which the coefficients are stored as three arrays...
Definition: lduMatrix.H:80
const scalarField & upper() const
Definition: lduMatrix.C:203
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1171
#define WarningInFunction
Report a warning using Foam::Warning.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:755
void operator-=(const lduMatrix &)
void operator+=(const lduMatrix &)
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)