Loading...
Searching...
No Matches
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-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33template<class Type>
36 scalarSquareMatrix(mSize),
37 source_(mSize)
38{}
39
40
41template<class Type>
43(
44 const label mSize,
45 const scalar coeffVal,
46 const Type& sourceVal
47)
49 scalarSquareMatrix(mSize, coeffVal),
50 source_(mSize, sourceVal)
51{}
52
53
54template<class Type>
58 source_(n, Foam::zero{})
59{}
60
61
62template<class Type>
64(
65 const scalarSquareMatrix& matrix,
66 const Field<Type>& source
67)
70 source_(source)
71{}
72
73
74template<class Type>
76:
78 source_(is)
79{}
80
81
82// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83
84template<class Type>
86{
87 scalarSquareMatrix tmpMatrix = *this;
88 Field<Type> sourceSol = source_;
89
90 Foam::solve(tmpMatrix, sourceSol);
91
92 return sourceSol;
93}
94
95
96template<class Type>
98{
99 scalarSquareMatrix luMatrix = *this;
100 Field<Type> sourceSol = source_;
101
102 Foam::LUsolve(luMatrix, sourceSol);
103
104 return sourceSol;
105}
106
107
108// * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
109
110template<class Type>
112{
113 if (this == &m)
115 return; // Self-assignment is a no-op
116 }
117
118 if (m() != m.m())
119 {
121 << "Different size matrices"
122 << abort(FatalError);
123 }
124
125 if (source_.size() != m.source_.size())
126 {
128 << "Different size source vectors"
129 << abort(FatalError);
130 }
131
133 source_ = m.source_;
134}
135
136
137// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
138
139template<class Type>
140Foam::simpleMatrix<Type> Foam::operator+
141(
142 const simpleMatrix<Type>& m1,
143 const simpleMatrix<Type>& m2
144)
145{
146 return simpleMatrix<Type>
147 (
148 static_cast<const scalarSquareMatrix&>(m1)
149 + static_cast<const scalarSquareMatrix&>(m2),
150 m1.source_ + m2.source_
151 );
152}
153
154
155template<class Type>
156Foam::simpleMatrix<Type> Foam::operator-
157(
158 const simpleMatrix<Type>& m1,
159 const simpleMatrix<Type>& m2
160)
161{
162 return simpleMatrix<Type>
163 (
164 static_cast<const scalarSquareMatrix&>(m1)
165 - static_cast<const scalarSquareMatrix&>(m2),
166 m1.source_ - m2.source_
167 );
168}
169
170
171template<class Type>
172Foam::simpleMatrix<Type> Foam::operator*
173(
174 const scalar s,
175 const simpleMatrix<Type>& m
176)
177{
178 return simpleMatrix<Type>(s*m.matrix_, s*m.source_);
179}
180
181
182// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
183
184template<class Type>
185Foam::Ostream& Foam::operator<<
186(
187 Ostream& os,
188 const simpleMatrix<Type>& m
189)
190{
192 return os;
193}
194
195
196// ************************************************************************* //
label n
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
label m() const noexcept
The number of rows.
Definition Matrix.H:261
label n() const noexcept
The number of columns.
Definition Matrix.H:271
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
SquareMatrix & operator=(const SquareMatrix &)=default
A simple square matrix solver with scalar coefficients.
void operator=(const simpleMatrix< Type > &)
Copy assignment.
simpleMatrix(const label n)
Construct given size.
Field< Type > solve() const
Solve the matrix using Gaussian elimination with pivoting and return the solution.
Field< Type > LUsolve() const
Solve the matrix using LU decomposition with pivoting and return the solution.
Field< Type > & source() noexcept
Return access to the source.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
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))
Namespace for OpenFOAM.
errorManip< error > abort(error &err)
Definition errorManip.H:139
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
SquareMatrix< scalar > scalarSquareMatrix
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting returning the LU form of the matrix and the sol...
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50