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:
43
45 newton_argument<double>(Nx, Ny),
47 Nphi_u(problem.NUmodes + problem.liftfield.size()),
48 Nphi_p(problem.NPmodes),
49 N_BC(problem.inletIndex.rows())
50 {}
51
52 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
53 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
54
55 int Nphi_u;
56 int Nphi_p;
57 int N_BC;
58
59 scalar nu;
60 scalar dt;
61 Eigen::VectorXd y_old;
62 Eigen::VectorXd BC;
63
65};
66
67struct newton_usmsr_n: public newton_argument<double>
68{
69 public:
71
73 newton_argument<double>(Nx, Ny),
75 Nphi_flux(problem.NFluxmodes),
76 Nphi_prec1(problem.NPrecmodes(0)),
77 Nphi_prec2(problem.NPrecmodes(1)),
78 Nphi_prec3(problem.NPrecmodes(2)),
79 Nphi_prec4(problem.NPrecmodes(3)),
80 Nphi_prec5(problem.NPrecmodes(4)),
81 Nphi_prec6(problem.NPrecmodes(5)),
82 Nphi_prec7(problem.NPrecmodes(6)),
83 Nphi_prec8(problem.NPrecmodes(7)),
84 Nphi_const(problem.NCmodes),
85 d_c(problem.NCmodes),
86 nsf_c(problem.NCmodes),
87 a_c(problem.NCmodes)
88 {}
89
90 int operator()(const Eigen::VectorXd& n, Eigen::VectorXd& fvecn) const;
91 int df(const Eigen::VectorXd& n, Eigen::MatrixXd& fjacn) const;
103 scalar nu;
104 //neutronic constants
105 scalar d;
106 scalar m;
107 scalar nsf;
108 scalar Keff;
109 scalar iv;
110 scalar l1; //lambda-ith
111 scalar l2;
112 scalar l3;
113 scalar l4;
114 scalar l5;
115 scalar l6;
116 scalar l7;
117 scalar l8;
118 scalar b1; //beta-ith
119 scalar b2;
120 scalar b3;
121 scalar b4;
122 scalar b5;
123 scalar b6;
124 scalar b7;
125 scalar b8;
126 scalar btot;
127 scalar Sc;
128
129 scalar dt;
130 Eigen::VectorXd d_c;
131 Eigen::VectorXd nsf_c;
132 Eigen::VectorXd a_c;
133 Eigen::VectorXd a_tmp; //to store u time coeff
134 Eigen::VectorXd w_old;
135 std::vector<SPLINTER::RBFSpline*> SPLINES_d;
136 std::vector<SPLINTER::RBFSpline*> SPLINES_nsf;
137 std::vector<SPLINTER::RBFSpline*> SPLINES_a;
139};
140
141struct newton_usmsr_t: public newton_argument<double>
142{
143 public:
145
147 newton_argument<double>(Nx, Ny),
148 problem(& problem),
149 Nphi_T(problem.NTmodes + problem.liftfieldT.size()),
150 Nphi_dec1(problem.NDecmodes(0)),
151 Nphi_dec2(problem.NDecmodes(1)),
152 Nphi_dec3(problem.NDecmodes(2)),
153 N_BCt(problem.inletIndexT.rows()),
154 Nphi_const(problem.NCmodes),
155 v_c(problem.NCmodes),
156 sp_c(problem.NCmodes),
157 txs_c(problem.NCmodes)
158 {}
159
160 int operator()(const Eigen::VectorXd& t, Eigen::VectorXd& fvect) const;
161 int df(const Eigen::VectorXd& t, Eigen::MatrixXd& fjact) const;
168 int N_BCt;
170 scalar nu;
171 //neutronic constants
172 scalar Keff;
173 //thermal constants
174 scalar sp;
175 scalar cp;
176 scalar dl1;
177 scalar dl2;
178 scalar dl3;
179 scalar db1;
180 scalar db2;
181 scalar db3;
182 scalar dbtot;
183 scalar Sc;
184 scalar Pr;
185
186 scalar dt;
187 Eigen::VectorXd a_tmp; //to store u time coeff.
188 Eigen::VectorXd c_tmp; //to store flux time coeff.
189 Eigen::VectorXd z_old;
190 Eigen::VectorXd BCt;
191 Eigen::VectorXd v_c;
192 Eigen::VectorXd sp_c;
193 Eigen::VectorXd txs_c;
194 std::vector<SPLINTER::RBFSpline*> SPLINES_v;
195 std::vector<SPLINTER::RBFSpline*> SPLINES_sp;
196 std::vector<SPLINTER::RBFSpline*> SPLINES_TXS;
197
199};
200
201/*---------------------------------------------------------------------------*\
202 Class reducedProblem Declaration
203\*---------------------------------------------------------------------------*/
205{
206 public:
207 reducedusMSR();
210
214 scalar time;
215 scalar dt;
216 scalar finalTime;
217 scalar tstart;
219
220 bool recall = false;
221 void solveOnline(Eigen::MatrixXd vel_now, Eigen::MatrixXd temp_now,
222 Eigen::VectorXd mu_online, int startSnap = 0);
223 void reconstructAP(fileName folder = "./ITHACAOutput/online_rec",
224 int printevery = 1);
225 void reconstruct_fd(fileName folder = "./ITHACAOutput/online_rec",
226 int printevery = 1);
227 void reconstruct_n(fileName folder = "./ITHACAOutput/online_rec",
228 int printevery = 1);
229 void reconstruct_t(fileName folder = "./ITHACAOutput/online_rec",
230 int printevery = 1);
231 void reconstruct_C(fileName folder = "./ITHACAOutput/online_rec",
232 int printevery = 1);
233
234
235};
236
237
238
239
240#endif
241
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