simpleMatrix.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 "simpleMatrix.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class Type>
32 :
33  scalarSquareMatrix(mSize),
34  source_(mSize)
35 {}
36 
37 
38 template<class Type>
40 (
41  const label mSize,
42  const scalar coeffVal,
43  const Type& sourceVal
44 )
45 :
46  scalarSquareMatrix(mSize, coeffVal),
47  source_(mSize, sourceVal)
48 {}
49 
50 
51 template<class Type>
53 (
54  const scalarSquareMatrix& matrix,
55  const Field<Type>& source
56 )
57 :
58  scalarSquareMatrix(matrix),
59  source_(source)
60 {}
61 
62 
63 template<class Type>
65 :
67  source_(is)
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 template<class Type>
75 {
76  scalarSquareMatrix tmpMatrix = *this;
77  Field<Type> sourceSol = source_;
78 
79  Foam::solve(tmpMatrix, sourceSol);
80 
81  return sourceSol;
82 }
83 
84 
85 template<class Type>
87 {
88  scalarSquareMatrix luMatrix = *this;
89  Field<Type> sourceSol = source_;
90 
91  Foam::LUsolve(luMatrix, sourceSol);
92 
93  return sourceSol;
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
98 
99 template<class Type>
101 {
102  if (this == &m)
103  {
105  << "Attempted assignment to self"
106  << abort(FatalError);
107  }
108 
109  if (m() != m.m())
110  {
112  << "Different size matrices"
113  << abort(FatalError);
114  }
115 
116  if (source_.size() != m.source_.size())
117  {
119  << "Different size source vectors"
121  }
122 
124  source_ = m.source_;
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
129 
130 template<class Type>
131 Foam::simpleMatrix<Type> Foam::operator+
132 (
133  const simpleMatrix<Type>& m1,
134  const simpleMatrix<Type>& m2
135 )
136 {
137  return simpleMatrix<Type>
138  (
139  static_cast<const scalarSquareMatrix&>(m1)
140  + static_cast<const scalarSquareMatrix&>(m2),
141  m1.source_ + m2.source_
142  );
143 }
144 
145 
146 template<class Type>
147 Foam::simpleMatrix<Type> Foam::operator-
148 (
149  const simpleMatrix<Type>& m1,
150  const simpleMatrix<Type>& m2
151 )
152 {
153  return simpleMatrix<Type>
154  (
155  static_cast<const scalarSquareMatrix&>(m1)
156  - static_cast<const scalarSquareMatrix&>(m2),
157  m1.source_ - m2.source_
158  );
159 }
160 
161 
162 template<class Type>
163 Foam::simpleMatrix<Type> Foam::operator*
164 (
165  const scalar s,
166  const simpleMatrix<Type>& m
167 )
168 {
169  return simpleMatrix<Type>(s*m.matrix_, s*m.source_);
170 }
171 
172 
173 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
174 
175 template<class Type>
176 Foam::Ostream& Foam::operator<<
177 (
178  Ostream& os,
179  const simpleMatrix<Type>& m
180 )
181 {
182  os << static_cast<const scalarSquareMatrix&>(m) << nl << m.source_;
183  return os;
184 }
185 
186 
187 // ************************************************************************* //
Pre-declare SubField and related Field type.
Definition: Field.H:83
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
label m() const
Return the number of rows.
Definition: MatrixI.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
void operator=(const zero)
Assignment of all elements to zero.
A simple square matrix solver with scalar coefficients.
Definition: simpleMatrix.H:65
void operator=(const simpleMatrix< Type > &)
Definition: simpleMatrix.C:100
Field< Type > solve() const
Solve the matrix using Gaussian elimination with pivoting.
Definition: simpleMatrix.C:74
simpleMatrix(const label)
Construct given size.
Definition: simpleMatrix.C:31
Field< Type > LUsolve() const
Solve the matrix using LU decomposition with pivoting.
Definition: simpleMatrix.C:86
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
SquareMatrix< scalar > scalarSquareMatrix
error FatalError
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
static const char nl
Definition: Ostream.H:266
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.