simpleMatrix.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 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 \*---------------------------------------------------------------------------*/
28 
29 #include "simpleMatrix.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 :
36  scalarSquareMatrix(mSize),
37  source_(mSize)
38 {}
39 
40 
41 template<class Type>
43 (
44  const label mSize,
45  const scalar coeffVal,
46  const Type& sourceVal
47 )
48 :
49  scalarSquareMatrix(mSize, coeffVal),
50  source_(mSize, sourceVal)
51 {}
52 
53 
54 template<class Type>
56 (
57  const scalarSquareMatrix& matrix,
58  const Field<Type>& source
59 )
60 :
61  scalarSquareMatrix(matrix),
62  source_(source)
63 {}
64 
65 
66 template<class Type>
68 :
70  source_(is)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
76 template<class Type>
78 {
79  scalarSquareMatrix tmpMatrix = *this;
80  Field<Type> sourceSol = source_;
81 
82  Foam::solve(tmpMatrix, sourceSol);
83 
84  return sourceSol;
85 }
86 
87 
88 template<class Type>
90 {
91  scalarSquareMatrix luMatrix = *this;
92  Field<Type> sourceSol = source_;
93 
94  Foam::LUsolve(luMatrix, sourceSol);
95 
96  return sourceSol;
97 }
98 
99 
100 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
101 
102 template<class Type>
104 {
105  if (this == &m)
106  {
107  return; // Self-assignment is a no-op
108  }
109 
110  if (m() != m.m())
111  {
113  << "Different size matrices"
114  << abort(FatalError);
115  }
116 
117  if (source_.size() != m.source_.size())
118  {
120  << "Different size source vectors"
121  << abort(FatalError);
122  }
123 
124  scalarSquareMatrix::operator=(m);
125  source_ = m.source_;
126 }
127 
128 
129 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
130 
131 template<class Type>
132 Foam::simpleMatrix<Type> Foam::operator+
133 (
134  const simpleMatrix<Type>& m1,
135  const simpleMatrix<Type>& m2
136 )
137 {
138  return simpleMatrix<Type>
139  (
140  static_cast<const scalarSquareMatrix&>(m1)
141  + static_cast<const scalarSquareMatrix&>(m2),
142  m1.source_ + m2.source_
143  );
144 }
145 
146 
147 template<class Type>
148 Foam::simpleMatrix<Type> Foam::operator-
149 (
150  const simpleMatrix<Type>& m1,
151  const simpleMatrix<Type>& m2
152 )
153 {
154  return simpleMatrix<Type>
155  (
156  static_cast<const scalarSquareMatrix&>(m1)
157  - static_cast<const scalarSquareMatrix&>(m2),
158  m1.source_ - m2.source_
159  );
160 }
161 
162 
163 template<class Type>
164 Foam::simpleMatrix<Type> Foam::operator*
165 (
166  const scalar s,
167  const simpleMatrix<Type>& m
168 )
169 {
170  return simpleMatrix<Type>(s*m.matrix_, s*m.source_);
171 }
172 
173 
174 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
175 
176 template<class Type>
177 Foam::Ostream& Foam::operator<<
178 (
179  Ostream& os,
180  const simpleMatrix<Type>& m
181 )
182 {
183  os << static_cast<const scalarSquareMatrix&>(m) << nl << m.source_;
184  return os;
185 }
186 
187 
188 // ************************************************************************* //
Foam::simpleMatrix::operator=
void operator=(const simpleMatrix< Type > &)
Definition: simpleMatrix.C:103
Foam::LUsolve
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Definition: scalarMatricesTemplates.C:215
Foam::simpleMatrix
A simple square matrix solver with scalar coefficients.
Definition: simpleMatrix.H:49
s
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))
Definition: gmvOutputSpray.H:25
Foam::simpleMatrix::LUsolve
Field< Type > LUsolve() const
Solve the matrix using LU decomposition with pivoting.
Definition: simpleMatrix.C:89
Foam::Field
Generic templated field type.
Definition: Field.H:63
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::simpleMatrix::simpleMatrix
simpleMatrix(const label)
Construct given size.
Definition: simpleMatrix.C:34
Foam::scalarSquareMatrix
SquareMatrix< scalar > scalarSquareMatrix
Definition: scalarMatrices.H:57
Foam::FatalError
error FatalError
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::SquareMatrix< scalar >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::nl
constexpr char nl
Definition: Ostream.H:385
simpleMatrix.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::simpleMatrix::solve
Field< Type > solve() const
Solve the matrix using Gaussian elimination with pivoting.
Definition: simpleMatrix.C:77