G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsFlowSolverBase.h
Go to the documentation of this file.
1
12#pragma once
13
18
19namespace gismo
20{
21
25template<class T, int MatOrder>
27{
28
29protected: // *** Class members ***
30
31 gsFlowSolverParams<T> m_params;
34
35 gsMatrix<T> m_solution;
36 unsigned m_iterationNumber;
37 gsStopwatch m_clock;
38 T m_initAssembT, m_assembT, m_relNorm;
39
40 bool m_fileOutput, m_dispOutput;
41 std::ofstream m_outFile;
42 std::stringstream m_outStream;
43
44
45public: // *** Constructor/destructor ***
46
48 m_params(params)
49 {
50 m_assemblerPtr = NULL;
51 }
52
53 virtual ~gsFlowSolverBase()
54 {
55 if (m_assemblerPtr)
56 {
57 delete m_assemblerPtr;
58 m_assemblerPtr = NULL;
59 }
60
61 if (m_linSolverPtr)
62 {
63 delete m_linSolverPtr;
64 m_linSolverPtr = NULL;
65 }
66
67 if (m_outFile.is_open())
68 m_outFile.close();
69 }
70
71
72protected: // *** Member functions ***
73
75 virtual void initMembers();
76
78 virtual void updateSizes();
79
81 virtual void createOutputFile();
82
84 real_t stopwatchStart();
85
87 real_t stopwatchStop();
88
89
90public: // *** Member functions ***
91
93 virtual void initialize();
94
96 virtual void initIteration()
97 { getLinSolver()->setupSolver(getAssembler()->matrix()); }
98
102 { getLinSolver()->setupSolver(mat); }
103
106 virtual void applySolver(gsMatrix<T>& solution)
107 { getLinSolver()->applySolver(getAssembler()->matrix(), getAssembler()->rhs(), solution); }
108
111
114 virtual void nextIteration(const unsigned numberOfIterations);
115
120 void solve(const int maxIterations, const T epsilon = 1e-3, const int minIterations = 1);
121
124 virtual void updateAssembler(const gsMatrix<T>& sol, bool updateSol = true);
125
128 virtual void updateAssembler(bool updateSol = true)
129 { updateAssembler(m_solution, updateSol); }
130
132 T solutionChangeRelNorm() const;
133
137 T solutionChangeRelNorm(gsMatrix<T> solOld, gsMatrix<T> solNew) const;
138
142 virtual void writeSolChangeRelNorm(gsMatrix<T> solOld, gsMatrix<T> solNew);
143
145 virtual T residualRelNorm() const
146 { return residualRelNorm(m_solution); }
147
149 virtual T residualRelNorm(const gsMatrix<T>& solution) const
151
156 { return getAssembler()->constructSolution(m_solution, unk); }
157
161 virtual void markDofsAsEliminatedZeros(const std::vector< gsMatrix< index_t > > & boundaryDofs, const int unk);
162
168 int checkGeoJacobian(int npts = -1, T dist = -1, T tol = -1);
169
170public: // *** Getters/setters ***
171
173 virtual std::string getName() { return "gsFlowSolverBase"; }
174
176 virtual gsFlowSolverParams<T> getParams() { return m_params; }
177
179 virtual gsFlowAssemblerBase<T, MatOrder>* getAssembler() const { return m_assemblerPtr; }
180
182 gsFlowLinSystSolver<T, MatOrder>* getLinSolver() const { return m_linSolverPtr; }
183
185 const gsMatrix<T> & getSolution() const { return m_solution; }
186
189 virtual void setSolution(const gsMatrix<T> & solVector) { m_solution = solVector; }
190
192 unsigned getIterationNumber() const { return m_iterationNumber; }
193
195 int numDofs() const { return getAssembler()->numDofs(); }
196
198 virtual const T getInitAssemblyTime() const { return m_initAssembT; }
199
201 virtual const T getAssemblyTime() const { return m_initAssembT + m_assembT; }
202
204 virtual const T getSolverSetupTime() const { return getLinSolver()->getSolverSetupTime(); }
205
207 virtual const T getSolveTime() const { return getLinSolver()->getSolveTime(); }
208
209
210}; // gsFlowSolverBase
211
212} // namespace gismo
213
214#ifndef GISMO_BUILD_LIB
215#include GISMO_HPP_HEADER(gsFlowSolverBase.hpp)
216#endif
A scalar of vector field defined on a m_parametric geometry.
Definition gsField.h:55
A base class for all assemblers in gsIncompressibleFlow.
Definition gsFlowAssemblerBase.h:25
Interface for classes solving linear systems inside the incompressible flow solvers (classes derived ...
Definition gsFlowLinSystSolver.h:27
A base class for all flow solvers in gsIncompressibleFlow.
Definition gsFlowSolverBase.h:27
virtual void nextIteration()
Perform next iteration step.
Definition gsFlowSolverBase.h:110
virtual void setSolution(const gsMatrix< T > &solVector)
Set a given solution vector as current solution.
Definition gsFlowSolverBase.h:189
virtual gsFlowAssemblerBase< T, MatOrder > * getAssembler() const
Returns a pointer to the assembler.
Definition gsFlowSolverBase.h:179
virtual void initIteration(const gsSparseMatrix< T, MatOrder > &mat)
Prepare for the solution process.
Definition gsFlowSolverBase.h:101
virtual void initialize()
Initialize the solver.
Definition gsFlowSolverBase.hpp:96
T solutionChangeRelNorm() const
Compute and return the relative norm of the solution change.
Definition gsFlowSolverBase.hpp:158
virtual const T getInitAssemblyTime() const
Returns the time of the initial matrix assembly.
Definition gsFlowSolverBase.h:198
virtual void writeSolChangeRelNorm(gsMatrix< T > solOld, gsMatrix< T > solNew)
Compute and display the relative norm of the solution change given the two successive solutions.
Definition gsFlowSolverBase.hpp:182
virtual gsFlowSolverParams< T > getParams()
Retrurns the solver parameters.
Definition gsFlowSolverBase.h:176
virtual T residualRelNorm() const
Compute and return the relative residual norm for the current solution.
Definition gsFlowSolverBase.h:145
virtual const T getSolveTime() const
Returns the total time spent on solving of the linear systems.
Definition gsFlowSolverBase.h:207
virtual void applySolver(gsMatrix< T > &solution)
Solve the linear system.
Definition gsFlowSolverBase.h:106
real_t stopwatchStart()
Start measuring time (decides whether to use gsStopwatch or MPI_Wtime)
Definition gsFlowSolverBase.hpp:61
gsField< T > constructSolution(int unk) const
Construct solution field for the unknown unk for the current solution vector.
Definition gsFlowSolverBase.h:155
virtual void markDofsAsEliminatedZeros(const std::vector< gsMatrix< index_t > > &boundaryDofs, const int unk)
Eliminate given DOFs as homogeneous Dirichlet boundary.
Definition gsFlowSolverBase.hpp:195
virtual const T getSolverSetupTime() const
Returns the total time spent on linear solver setup.
Definition gsFlowSolverBase.h:204
virtual std::string getName()
Retrurns the name of the class as a string.
Definition gsFlowSolverBase.h:173
virtual void updateAssembler(bool updateSol=true)
Update the assembler with current solution.
Definition gsFlowSolverBase.h:128
void solve(const int maxIterations, const T epsilon=1e-3, const int minIterations=1)
Solve the incompressible Navier-Stokes problem.
Definition gsFlowSolverBase.hpp:122
virtual T residualRelNorm(const gsMatrix< T > &solution) const
Compute and return the relative residual norm for the given solution.
Definition gsFlowSolverBase.h:149
unsigned getIterationNumber() const
Returns the current iteration number.
Definition gsFlowSolverBase.h:192
virtual const T getAssemblyTime() const
Returns the total time spent on matrix assembly.
Definition gsFlowSolverBase.h:201
virtual void initMembers()
Initialize all members.
Definition gsFlowSolverBase.hpp:20
virtual void createOutputFile()
Create output file.
Definition gsFlowSolverBase.hpp:45
real_t stopwatchStop()
Stop measuring time (decides whether to use gsStopwatch or MPI_Wtime)
Definition gsFlowSolverBase.hpp:79
virtual void updateAssembler(const gsMatrix< T > &sol, bool updateSol=true)
Update the assembler with a given solution.
Definition gsFlowSolverBase.hpp:147
virtual void updateSizes()
Update sizes of members (when DOF numbers change, e.g. after markDofsAsEliminatedZeros()).
Definition gsFlowSolverBase.hpp:38
int checkGeoJacobian(int npts=-1, T dist=-1, T tol=-1)
Check values of jacobian near boundaries of all patches.
Definition gsFlowSolverBase.hpp:203
const gsMatrix< T > & getSolution() const
Returns the current solution vector.
Definition gsFlowSolverBase.h:185
int numDofs() const
Returns the total number of DOFs (the matrix size).
Definition gsFlowSolverBase.h:195
gsFlowLinSystSolver< T, MatOrder > * getLinSolver() const
Returns a pointer to the linear system solver.
Definition gsFlowSolverBase.h:182
virtual void initIteration()
Prepare for the solution process.
Definition gsFlowSolverBase.h:96
A class that holds all parameters needed by the incompressible flow solver.
Definition gsFlowSolverParams.h:34
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
A class that holds all parameters needed by the incompressible flow solver.
The G+Smo namespace, containing all definitions for the library.