Loading...
Searching...
No Matches
ITHACAparameters.H
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser 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 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Class
25 ITHACAparameters
26Description
27 Some parameters for input output informations
28SourceFiles
29 ITHACAparameters.H
30\*---------------------------------------------------------------------------*/
31#ifndef Parameters_H
32#define Parameters_H
33
34#include <iostream>
35#include "fvCFD.H"
36#include "ITHACAassert.H"
37
40class ITHACAparameters
41{
42 public:
43
44 ITHACAparameters(const fvMesh& mesh, Time& localTime);
45
48
50 label precision;
51
53 bool exportPython;
54 bool exportMatlab;
55 bool exportTxt;
56 bool exportNpy;
57 bool debug;
58 bool warnings;
59 bool correctBC;
60
62 std::_Ios_Fmtflags outytpe;
63
65 IOdictionary* ITHACAdict;
66
68 Time& runTime;
69
71 const fvMesh& mesh;
72
74 static ITHACAparameters* instance;
75
76 public:
85 static ITHACAparameters* getInstance(const fvMesh& mesh, Time& localTime);
86
92 static ITHACAparameters* getInstance();
93
94
97
98
99 private:
100
101};
102
103#endif
Implementation of the assert function for ITHACA-FV.
Time & runTime
runTime defined locally
static ITHACAparameters * instance
Pointer to the ITHACAparameters class.
word eigensolver
type of eigensolver used in the eigenvalue decomposition can be either be eigen or spectra
label precision
precision of the output Market Matrix objects (i.e. reduced matrices, eigenvalues,...
IOdictionary * ITHACAdict
Dictionary for input objects from file.
const fvMesh & mesh
Mesh object defined locally.
~ITHACAparameters()=delete
Delete empty constructor.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
std::_Ios_Fmtflags outytpe
type of output format can be fixed or scientific