Loading...
Searching...
No Matches
ReducedUnsteadyMSR.H
Go to the documentation of this file.
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/>.
24\*---------------------------------------------------------------------------*/
25#ifndef ReducedusMSR_H
26#define ReducedusMSR_H
27
28#include "fvCFD.H"
29#include "IOmanip.H"
30#include "ReducedMSR.H"
31#include "usmsrProblem.H"
32#include <Eigen/Dense>
33#include <unsupported/Eigen/NonLinearOptimization>
34#include <unsupported/Eigen/NumericalDiff>
35#include <fstream>
36#include <vector>
37#include <stdlib.h>
38
39struct newton_usmsr_fd: public newton_argument<double>
40{
41 public:
44 newton_argument<double>(Nx, Ny),
46 Nphi_u(problem.NUmodes + problem.liftfield.size()),
47 Nphi_p(problem.NPmodes),
48 N_BC(problem.inletIndex.rows())
49 {}
50 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
51 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
52
53 int Nphi_u;
54 int Nphi_p;
55 int N_BC;
56
57 scalar nu;
58 scalar dt;
59 Eigen::VectorXd y_old;
60 Eigen::VectorXd BC;
61
63};
64
65struct newton_usmsr_n: public newton_argument<double>
66{
67 public:
70 newton_argument<double>(Nx, Ny),
72 Nphi_flux(problem.NFluxmodes),
73 Nphi_prec1(problem.NPrecmodes(0)),
74 Nphi_prec2(problem.NPrecmodes(1)),
75 Nphi_prec3(problem.NPrecmodes(2)),
76 Nphi_prec4(problem.NPrecmodes(3)),
77 Nphi_prec5(problem.NPrecmodes(4)),
78 Nphi_prec6(problem.NPrecmodes(5)),
79 Nphi_prec7(problem.NPrecmodes(6)),
80 Nphi_prec8(problem.NPrecmodes(7)),
81 Nphi_const(problem.NCmodes),
82 d_c(problem.NCmodes),
83 nsf_c(problem.NCmodes),
84 a_c(problem.NCmodes)
85 {}
86 int operator()(const Eigen::VectorXd& n, Eigen::VectorXd& fvecn) const;
87 int df(const Eigen::VectorXd& n, Eigen::MatrixXd& fjacn) const;
97 int Nphi_u;
99 scalar nu;
100 //neutronic constants
101 scalar d;
102 scalar m;
103 scalar nsf;
104 scalar Keff;
105 scalar iv;
106 scalar l1; //lambda-ith
107 scalar l2;
108 scalar l3;
109 scalar l4;
110 scalar l5;
111 scalar l6;
112 scalar l7;
113 scalar l8;
114 scalar b1; //beta-ith
115 scalar b2;
116 scalar b3;
117 scalar b4;
118 scalar b5;
119 scalar b6;
120 scalar b7;
121 scalar b8;
122 scalar btot;
123 scalar Sc;
124
125 scalar dt;
126 Eigen::VectorXd d_c;
127 Eigen::VectorXd nsf_c;
128 Eigen::VectorXd a_c;
129 Eigen::VectorXd a_tmp; //to store u time coeff
130 Eigen::VectorXd w_old;
131 std::vector<SPLINTER::RBFSpline*> SPLINES_d;
132 std::vector<SPLINTER::RBFSpline*> SPLINES_nsf;
133 std::vector<SPLINTER::RBFSpline*> SPLINES_a;
135};
136
137struct newton_usmsr_t: public newton_argument<double>
138{
139 public:
142 newton_argument<double>(Nx, Ny),
144 Nphi_T(problem.NTmodes + problem.liftfieldT.size()),
145 Nphi_dec1(problem.NDecmodes(0)),
146 Nphi_dec2(problem.NDecmodes(1)),
147 Nphi_dec3(problem.NDecmodes(2)),
148 N_BCt(problem.inletIndexT.rows()),
149 Nphi_const(problem.NCmodes),
150 v_c(problem.NCmodes),
151 sp_c(problem.NCmodes),
152 txs_c(problem.NCmodes)
153 {}
154 int operator()(const Eigen::VectorXd& t, Eigen::VectorXd& fvect) const;
155 int df(const Eigen::VectorXd& t, Eigen::MatrixXd& fjact) const;
162 int N_BCt;
164 scalar nu;
165 //neutronic constants
166 scalar Keff;
167 //thermal constants
168 scalar sp;
169 scalar cp;
170 scalar dl1;
171 scalar dl2;
172 scalar dl3;
173 scalar db1;
174 scalar db2;
175 scalar db3;
176 scalar dbtot;
177 scalar Sc;
178 scalar Pr;
179
180 scalar dt;
181 Eigen::VectorXd a_tmp; //to store u time coeff.
182 Eigen::VectorXd c_tmp; //to store flux time coeff.
183 Eigen::VectorXd z_old;
184 Eigen::VectorXd BCt;
185 Eigen::VectorXd v_c;
186 Eigen::VectorXd sp_c;
187 Eigen::VectorXd txs_c;
188 std::vector<SPLINTER::RBFSpline*> SPLINES_v;
189 std::vector<SPLINTER::RBFSpline*> SPLINES_sp;
190 std::vector<SPLINTER::RBFSpline*> SPLINES_TXS;
191
193};
194
195/*---------------------------------------------------------------------------*\
196 Class reducedProblem Declaration
197\*---------------------------------------------------------------------------*/
199{
200 public:
201 reducedusMSR();
204
208 scalar time;
209 scalar dt;
210 scalar finalTime;
211 scalar tstart;
213
214 bool recall = false;
215 void solveOnline(Eigen::MatrixXd vel_now, Eigen::MatrixXd temp_now,
216 Eigen::VectorXd mu_online, int startSnap = 0);
217 void reconstructAP(fileName folder = "./ITHACAOutput/online_rec",
218 int printevery = 1);
219 void reconstruct_fd(fileName folder = "./ITHACAOutput/online_rec",
220 int printevery = 1);
221 void reconstruct_n(fileName folder = "./ITHACAOutput/online_rec",
222 int printevery = 1);
223 void reconstruct_t(fileName folder = "./ITHACAOutput/online_rec",
224 int printevery = 1);
225 void reconstruct_C(fileName folder = "./ITHACAOutput/online_rec",
226 int printevery = 1);
227
228
229};
230
231
232
233
234#endif
235
Template object created to solve non linear problems using the Eigen library.
virtual void solveOnline()
Virtual Method to perform and online Solve.
void reconstructAP(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
void reconstruct_fd(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
newton_usmsr_n newton_object_n
void reconstruct_t(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
usmsrProblem * problem
void reconstruct_C(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
void reconstruct_n(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
newton_usmsr_fd newton_object_fd
newton_usmsr_t newton_object_t
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
Eigen::VectorXd y_old
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
Eigen::VectorXd BC
usmsrProblem * problem
newton_usmsr_fd(int Nx, int Ny, usmsrProblem &problem)
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const
Eigen::VectorXd a_c
Eigen::VectorXd a_tmp
std::vector< SPLINTER::RBFSpline * > SPLINES_d
Eigen::VectorXd d_c
std::vector< SPLINTER::RBFSpline * > SPLINES_a
Eigen::VectorXd w_old
Eigen::VectorXd nsf_c
usmsrProblem * problem
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
newton_usmsr_n(int Nx, int Ny, usmsrProblem &problem)
std::vector< SPLINTER::RBFSpline * > SPLINES_nsf
Eigen::VectorXd txs_c
Eigen::VectorXd sp_c
Eigen::VectorXd c_tmp
Eigen::VectorXd z_old
int operator()(const Eigen::VectorXd &t, Eigen::VectorXd &fvect) const
std::vector< SPLINTER::RBFSpline * > SPLINES_TXS
Eigen::VectorXd v_c
std::vector< SPLINTER::RBFSpline * > SPLINES_v
Eigen::VectorXd a_tmp
newton_usmsr_t(int Nx, int Ny, usmsrProblem &problem)
std::vector< SPLINTER::RBFSpline * > SPLINES_sp
usmsrProblem * problem
int df(const Eigen::VectorXd &t, Eigen::MatrixXd &fjact) const
Eigen::VectorXd BCt