Loading...
Searching...
No Matches
ReducedMSR.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 ReducedMSR_H
26#define ReducedMSR_H
27
28#include "fvCFD.H"
29#include "IOmanip.H"
30#include "ReducedProblem.H"
31#include "msrProblem.H"
32#include "ITHACAutilities.H"
33#include <Eigen/Dense>
34#include <unsupported/Eigen/NonLinearOptimization>
35#include <unsupported/Eigen/NumericalDiff>
36
38struct newton_msr_fd: public newton_argument<double>
39{
40 public:
42
43 newton_msr_fd(int Nx, int Ny, msrProblem& problem) :
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
51 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
52 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
53
54 int Nphi_u;
55 int Nphi_p;
56 int N_BC;
57
58 scalar nu;
59
60 Eigen::VectorXd BC;
62};
63
64struct newton_msr_n: public newton_argument<double>
65{
66 public:
68
69 newton_msr_n(int Nx, int Ny, msrProblem& problem) :
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
87 int operator()(const Eigen::VectorXd& n, Eigen::VectorXd& fvecn) const;
88 int df(const Eigen::VectorXd& n, Eigen::MatrixXd& fjacn) const;
89
90 int Nphi_u;
101 //fluid-dynamics constants
102 scalar nu;
103 //neutronic constants
104 scalar iv;
105 scalar l1; //lambda-ith
106 scalar l2;
107 scalar l3;
108 scalar l4;
109 scalar l5;
110 scalar l6;
111 scalar l7;
112 scalar l8;
113 scalar b1; //beta-ith
114 scalar b2;
115 scalar b3;
116 scalar b4;
117 scalar b5;
118 scalar b6;
119 scalar b7;
120 scalar b8;
121 scalar btot;
122 scalar Sc;
123
124 Eigen::VectorXd d_c;
125 Eigen::VectorXd nsf_c;
126 Eigen::VectorXd a_c;
127 Eigen::VectorXd a_tmp; //to store u time coeff
128 std::vector<SPLINTER::RBFSpline*> SPLINES_d;
129 std::vector<SPLINTER::RBFSpline*> SPLINES_nsf;
130 std::vector<SPLINTER::RBFSpline*> SPLINES_a;
132};
133
134struct newton_msr_t: public newton_argument<double>
135{
136 public:
138
139 newton_msr_t(int Nx, int Ny, msrProblem& problem) :
140 newton_argument<double>(Nx, Ny),
141 problem(& problem),
142 Nphi_T(problem.NTmodes),
143 Nphi_dec1(problem.NDecmodes(0)),
144 Nphi_dec2(problem.NDecmodes(1)),
145 Nphi_dec3(problem.NDecmodes(2)),
146 N_BCt(problem.inletIndexT.rows()),
147 Nphi_const(problem.NCmodes),
148 v_c(problem.NCmodes),
149 sp_c(problem.NCmodes),
150 txs_c(problem.NCmodes)
151 {}
152
153 int operator()(const Eigen::VectorXd& n, Eigen::VectorXd& fvecn) const;
154 int df(const Eigen::VectorXd& n, Eigen::MatrixXd& fjacn) const;
155
162 int N_BCt;
164 //fluid-dynamics constants
165 scalar nu;
166 //thermal constants
167 scalar cp;
168 scalar dl1;
169 scalar dl2;
170 scalar dl3;
171 scalar db1;
172 scalar db2;
173 scalar db3;
174 scalar dbtot;
175 scalar Sc;
176 scalar Pr;
177
178 Eigen::VectorXd a_tmp;
179 Eigen::VectorXd c_tmp;
180 Eigen::VectorXd BCt;
181 Eigen::VectorXd v_c;
182 Eigen::VectorXd sp_c;
183 Eigen::VectorXd txs_c;
184 std::vector<SPLINTER::RBFSpline*> SPLINES_v;
185 std::vector<SPLINTER::RBFSpline*> SPLINES_sp;
186 std::vector<SPLINTER::RBFSpline*> SPLINES_TXS;
188};
189
190
191/*---------------------------------------------------------------------------*\
192 Class reducedProblem Declaration
193\*---------------------------------------------------------------------------*/
194
196{
197
198 public:
199
200 reducedMSR();
201 explicit reducedMSR(msrProblem& problem);
203
205 Eigen::VectorXd y_old;
206 Eigen::VectorXd w_old;
207 Eigen::VectorXd z_old;
208
210 Eigen::VectorXd y;
211 Eigen::VectorXd w;
212 Eigen::VectorXd z;
213
214
216 //fluid-dynamics constants
217 scalar nu;
218 scalar rho;
219 scalar Pr;
220 scalar Sc;
221 //neutronic constants
222 scalar d;
223 scalar m;
224 scalar nsf;
225 scalar a;
226 scalar Keff;
227 scalar iv;
228 scalar l1; //lambda-ith
229 scalar l2;
230 scalar l3;
231 scalar l4;
232 scalar l5;
233 scalar l6;
234 scalar l7;
235 scalar l8;
236 scalar b1; //beta-ith
237 scalar b2;
238 scalar b3;
239 scalar b4;
240 scalar b5;
241 scalar b6;
242 scalar b7;
243 scalar b8;
244 scalar btot;
245 //thermal constants
246 scalar cp;
247 scalar dl1;
248 scalar dl2;
249 scalar dl3;
250 scalar db1;
251 scalar db2;
252 scalar db3;
253 scalar dbtot;
254 scalar sp;
255
256
257
259 List < Eigen::MatrixXd> online_solution_fd;
260 List < Eigen::MatrixXd> online_solution_n;
261 List < Eigen::MatrixXd> online_solution_t;
262 List < Eigen::MatrixXd> online_solution_C;
263
265 PtrList<volVectorField> Umodes;
266 PtrList<volScalarField> Pmodes;
267 PtrList<volScalarField> Fluxmodes;
268 PtrList<volScalarField> Prec1modes;
269 PtrList<volScalarField> Prec2modes;
270 PtrList<volScalarField> Prec3modes;
271 PtrList<volScalarField> Prec4modes;
272 PtrList<volScalarField> Prec5modes;
273 PtrList<volScalarField> Prec6modes;
274 PtrList<volScalarField> Prec7modes;
275 PtrList<volScalarField> Prec8modes;
276 PtrList<volScalarField> Tmodes;
277 PtrList<volScalarField> Dec1modes;
278 PtrList<volScalarField> Dec2modes;
279 PtrList<volScalarField> Dec3modes;
280 PtrList<volScalarField> vmodes;
281 PtrList<volScalarField> Dmodes;
282 PtrList<volScalarField> NSFmodes;
283 PtrList<volScalarField> Amodes;
284 PtrList<volScalarField> SPmodes;
285 PtrList<volScalarField> TXSmodes;
286
287
289 PtrList<volVectorField> Usnapshots;
290 PtrList<volScalarField> Psnapshots;
291 PtrList<volScalarField> Fluxsnapshots;
292 PtrList<volScalarField> Prec1snapshots;
293 PtrList<volScalarField> Prec2snapshots;
294 PtrList<volScalarField> Prec3snapshots;
295 PtrList<volScalarField> Prec4snapshots;
296 PtrList<volScalarField> Prec5snapshots;
297 PtrList<volScalarField> Prec6snapshots;
298 PtrList<volScalarField> Prec7snapshots;
299 PtrList<volScalarField> Prec8snapshots;
300 PtrList<volScalarField> Tsnapshots;
301 PtrList<volScalarField> Dec1snapshots;
302 PtrList<volScalarField> Dec2snapshots;
303 PtrList<volScalarField> Dec3snapshots;
304 PtrList<volScalarField> vsnapshots;
305 PtrList<volScalarField> Dsnapshots;
306 PtrList<volScalarField> NSFsnapshots;
307 PtrList<volScalarField> Asnapshots;
308 PtrList<volScalarField> SPsnapshots;
309 PtrList<volScalarField> TXSsnapshots;
310
312 PtrList<volVectorField> UREC;
313 PtrList<volScalarField> PREC;
314 PtrList<volScalarField> FLUXREC;
315 PtrList<volScalarField> PREC1REC;
316 PtrList<volScalarField> PREC2REC;
317 PtrList<volScalarField> PREC3REC;
318 PtrList<volScalarField> PREC4REC;
319 PtrList<volScalarField> PREC5REC;
320 PtrList<volScalarField> PREC6REC;
321 PtrList<volScalarField> PREC7REC;
322 PtrList<volScalarField> PREC8REC;
323 PtrList<volScalarField> TREC;
324 PtrList<volScalarField> DEC1REC;
325 PtrList<volScalarField> DEC2REC;
326 PtrList<volScalarField> DEC3REC;
327 PtrList<volScalarField> POWERDENSREC;
328 PtrList<volScalarField> vREC;
329 PtrList<volScalarField> DREC;
330 PtrList<volScalarField> NSFREC;
331 PtrList<volScalarField> AREC;
332 PtrList<volScalarField> SPREC;
333 PtrList<volScalarField> TXSREC;
334
339
342
360
363 bool recall = false;
364
366 int N_BC;
367 int N_BCt;
368
371
372 //Methods:
373
374
375
377 void solveOnline(Eigen::MatrixXd vel_now, Eigen::MatrixXd temp_now,
378 Eigen::VectorXd mu_online);
379
381 void clearFields();
382
388 void reconstructAP(fileName folder = "./ITHACAOutput/online_rec",
389 int printevery = 1);
391 void reconstruct_fd(fileName folder = "./ITHACAOutput/online_rec",
392 int printevery = 1);
394 void reconstruct_n(fileName folder = "./ITHACAOutput/online_rec",
395 int printevery = 1);
397 void reconstruct_t(fileName folder = "./ITHACAOutput/online_rec",
398 int printevery = 1);
400 void reconstruct_C(fileName folder = "./ITHACAOutput/online_rec",
401 int printevery = 1);
402
403
404
405 protected:
406
407 //Methods:
408
411
412
413
414};
415
416
417#endif
Header file of the ITHACAutilities namespace.
Header file of the reducedProblem class.
Class to implement Molten Salt Reactor multiphysics problem.
Definition msrProblem.H:36
PtrList< volScalarField > NSFREC
Definition ReducedMSR.H:330
scalar cp
Definition ReducedMSR.H:246
scalar Keff
Definition ReducedMSR.H:226
PtrList< volScalarField > Tsnapshots
Definition ReducedMSR.H:300
PtrList< volScalarField > DEC3REC
Definition ReducedMSR.H:326
msrProblem * problem
Pointer to the FOM problem.
Definition ReducedMSR.H:341
scalar Pr
Definition ReducedMSR.H:219
PtrList< volScalarField > vREC
Definition ReducedMSR.H:328
scalar b4
Definition ReducedMSR.H:239
scalar l6
Definition ReducedMSR.H:233
void reconstruct_n(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct neutronics
Definition ReducedMSR.C:751
newton_msr_n newton_object_n
Definition ReducedMSR.H:337
PtrList< volScalarField > DEC2REC
Definition ReducedMSR.H:325
PtrList< volScalarField > Dec1modes
Definition ReducedMSR.H:277
PtrList< volScalarField > Pmodes
Definition ReducedMSR.H:266
PtrList< volVectorField > Usnapshots
List of pointers to store the snapshots for each field.
Definition ReducedMSR.H:289
PtrList< volScalarField > Psnapshots
Definition ReducedMSR.H:290
PtrList< volScalarField > FLUXREC
Definition ReducedMSR.H:314
PtrList< volVectorField > Umodes
List of pointers to store the modes for each field.
Definition ReducedMSR.H:265
PtrList< volScalarField > Prec8snapshots
Definition ReducedMSR.H:299
PtrList< volScalarField > Dec2modes
Definition ReducedMSR.H:278
scalar nu
characteristic constants of the problem
Definition ReducedMSR.H:217
PtrList< volScalarField > TXSREC
Definition ReducedMSR.H:333
bool recall
boolean variable to check if the user wants to reconstruct all the three physics of the problem
Definition ReducedMSR.H:363
PtrList< volScalarField > PREC4REC
Definition ReducedMSR.H:318
PtrList< volScalarField > Prec2modes
Definition ReducedMSR.H:269
void reconstruct_C(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct temperature dependent constants
Definition ReducedMSR.C:943
PtrList< volScalarField > Prec6snapshots
Definition ReducedMSR.H:297
PtrList< volScalarField > PREC3REC
Definition ReducedMSR.H:317
PtrList< volScalarField > vsnapshots
Definition ReducedMSR.H:304
scalar l7
Definition ReducedMSR.H:234
PtrList< volVectorField > UREC
Recontructed fields.
Definition ReducedMSR.H:312
List< Eigen::MatrixXd > online_solution_n
Definition ReducedMSR.H:260
scalar b2
Definition ReducedMSR.H:237
PtrList< volScalarField > Prec7modes
Definition ReducedMSR.H:274
PtrList< volScalarField > Prec6modes
Definition ReducedMSR.H:273
scalar b1
Definition ReducedMSR.H:236
PtrList< volScalarField > SPmodes
Definition ReducedMSR.H:284
void reconstruct_t(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct thermal fields
Definition ReducedMSR.C:867
PtrList< volScalarField > Prec5snapshots
Definition ReducedMSR.H:296
PtrList< volScalarField > TXSsnapshots
Definition ReducedMSR.H:309
scalar b8
Definition ReducedMSR.H:243
scalar dl2
Definition ReducedMSR.H:248
Eigen::VectorXd z
Definition ReducedMSR.H:212
PtrList< volScalarField > Dec3snapshots
Definition ReducedMSR.H:303
PtrList< volScalarField > Prec3modes
Definition ReducedMSR.H:270
scalar Sc
Definition ReducedMSR.H:220
PtrList< volScalarField > Prec3snapshots
Definition ReducedMSR.H:294
PtrList< volScalarField > Prec7snapshots
Definition ReducedMSR.H:298
scalar db1
Definition ReducedMSR.H:250
scalar nsf
Definition ReducedMSR.H:224
scalar b3
Definition ReducedMSR.H:238
scalar db2
Definition ReducedMSR.H:251
newton_msr_fd newton_object_fd
Newton object used to solve the non linear problem.
Definition ReducedMSR.H:336
PtrList< volScalarField > Amodes
Definition ReducedMSR.H:283
void loadConstants(msrProblem *problem)
Method to load all the constants needed in the ROM from ///the FOM.
PtrList< volScalarField > DEC1REC
Definition ReducedMSR.H:324
PtrList< volScalarField > NSFmodes
Definition ReducedMSR.H:282
scalar iv
Definition ReducedMSR.H:227
PtrList< volScalarField > Fluxsnapshots
Definition ReducedMSR.H:291
void reconstruct_fd(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct fluid-dynamics
Definition ReducedMSR.C:706
scalar dbtot
Definition ReducedMSR.H:253
int N_BC
Number of parametrized boundary conditions.
Definition ReducedMSR.H:366
scalar db3
Definition ReducedMSR.H:252
PtrList< volScalarField > Dec1snapshots
Definition ReducedMSR.H:301
PtrList< volScalarField > SPsnapshots
Definition ReducedMSR.H:308
scalar l8
Definition ReducedMSR.H:235
PtrList< volScalarField > Asnapshots
Definition ReducedMSR.H:307
PtrList< volScalarField > Prec4modes
Definition ReducedMSR.H:271
PtrList< volScalarField > PREC1REC
Definition ReducedMSR.H:315
scalar btot
Definition ReducedMSR.H:244
PtrList< volScalarField > Prec1snapshots
Definition ReducedMSR.H:292
PtrList< volScalarField > Dsnapshots
Definition ReducedMSR.H:305
PtrList< volScalarField > Prec8modes
Definition ReducedMSR.H:275
List< Eigen::MatrixXd > online_solution_fd
List of Eigen matrices to store the online solution.
Definition ReducedMSR.H:259
PtrList< volScalarField > POWERDENSREC
Definition ReducedMSR.H:327
PtrList< volScalarField > PREC6REC
Definition ReducedMSR.H:320
PtrList< volScalarField > TREC
Definition ReducedMSR.H:323
Eigen::VectorXd y_old
Vector to store the previous solution during the Newton procedure.
Definition ReducedMSR.H:205
PtrList< volScalarField > NSFsnapshots
Definition ReducedMSR.H:306
PtrList< volScalarField > TXSmodes
Definition ReducedMSR.H:285
PtrList< volScalarField > Dec3modes
Definition ReducedMSR.H:279
scalar dl1
Definition ReducedMSR.H:247
int count_online_solve
Counter to count the online solutions.
Definition ReducedMSR.H:370
PtrList< volScalarField > AREC
Definition ReducedMSR.H:331
PtrList< volScalarField > vmodes
Definition ReducedMSR.H:280
PtrList< volScalarField > Prec4snapshots
Definition ReducedMSR.H:295
PtrList< volScalarField > PREC2REC
Definition ReducedMSR.H:316
Eigen::VectorXd w
Definition ReducedMSR.H:211
PtrList< volScalarField > Tmodes
Definition ReducedMSR.H:276
PtrList< volScalarField > Dmodes
Definition ReducedMSR.H:281
int Nphi_u
Number of modes for each field.
Definition ReducedMSR.H:344
List< Eigen::MatrixXd > online_solution_t
Definition ReducedMSR.H:261
void reconstructAP(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Methods to reconstruct a solution from an online solve with a PPE stabilisation technique.
Definition ReducedMSR.C:694
PtrList< volScalarField > Prec5modes
Definition ReducedMSR.H:272
scalar rho
Definition ReducedMSR.H:218
PtrList< volScalarField > Prec2snapshots
Definition ReducedMSR.H:293
newton_msr_t newton_object_t
Definition ReducedMSR.H:338
PtrList< volScalarField > PREC
Definition ReducedMSR.H:313
PtrList< volScalarField > SPREC
Definition ReducedMSR.H:332
Eigen::VectorXd z_old
Definition ReducedMSR.H:207
scalar l4
Definition ReducedMSR.H:231
scalar sp
Definition ReducedMSR.H:254
void clearFields()
Method to clear all the fields of MSR (sets the size to zero)
scalar dl3
Definition ReducedMSR.H:249
PtrList< volScalarField > Dec2snapshots
Definition ReducedMSR.H:302
scalar b6
Definition ReducedMSR.H:241
PtrList< volScalarField > Prec1modes
Definition ReducedMSR.H:268
scalar l5
Definition ReducedMSR.H:232
scalar l2
Definition ReducedMSR.H:229
scalar b5
Definition ReducedMSR.H:240
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
Definition ReducedMSR.H:210
Eigen::VectorXd w_old
Definition ReducedMSR.H:206
scalar b7
Definition ReducedMSR.H:242
scalar l1
Definition ReducedMSR.H:228
PtrList< volScalarField > PREC5REC
Definition ReducedMSR.H:319
PtrList< volScalarField > Fluxmodes
Definition ReducedMSR.H:267
PtrList< volScalarField > PREC8REC
Definition ReducedMSR.H:322
PtrList< volScalarField > DREC
Definition ReducedMSR.H:329
List< Eigen::MatrixXd > online_solution_C
Definition ReducedMSR.H:262
PtrList< volScalarField > PREC7REC
Definition ReducedMSR.H:321
scalar l3
Definition ReducedMSR.H:230
virtual void solveOnline()
Virtual Method to perform and online Solve.
reducedProblem()
Construct Null.
Structure to implement a newton object for a stationary MSR problem.
Definition ReducedMSR.H:39
Eigen::VectorXd BC
Definition ReducedMSR.H:60
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
Definition ReducedMSR.C:232
newton_msr_fd(int Nx, int Ny, msrProblem &problem)
Definition ReducedMSR.H:43
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
Definition ReducedMSR.C:185
msrProblem * problem
Definition ReducedMSR.H:61
Eigen::VectorXd a_tmp
Definition ReducedMSR.H:127
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
Definition ReducedMSR.C:414
std::vector< SPLINTER::RBFSpline * > SPLINES_d
Definition ReducedMSR.H:128
Eigen::VectorXd d_c
Definition ReducedMSR.H:124
newton_msr_n(int Nx, int Ny, msrProblem &problem)
Definition ReducedMSR.H:69
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const
Definition ReducedMSR.C:239
std::vector< SPLINTER::RBFSpline * > SPLINES_a
Definition ReducedMSR.H:130
Eigen::VectorXd nsf_c
Definition ReducedMSR.H:125
std::vector< SPLINTER::RBFSpline * > SPLINES_nsf
Definition ReducedMSR.H:129
Eigen::VectorXd a_c
Definition ReducedMSR.H:126
msrProblem * problem
Definition ReducedMSR.H:131
Eigen::VectorXd a_tmp
Definition ReducedMSR.H:178
newton_msr_t(int Nx, int Ny, msrProblem &problem)
Definition ReducedMSR.H:139
Eigen::VectorXd sp_c
Definition ReducedMSR.H:182
std::vector< SPLINTER::RBFSpline * > SPLINES_TXS
Definition ReducedMSR.H:186
Eigen::VectorXd txs_c
Definition ReducedMSR.H:183
Eigen::VectorXd c_tmp
Definition ReducedMSR.H:179
msrProblem * problem
Definition ReducedMSR.H:187
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
Definition ReducedMSR.C:512
std::vector< SPLINTER::RBFSpline * > SPLINES_v
Definition ReducedMSR.H:184
Eigen::VectorXd BCt
Definition ReducedMSR.H:180
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const
Definition ReducedMSR.C:422
Eigen::VectorXd v_c
Definition ReducedMSR.H:181
std::vector< SPLINTER::RBFSpline * > SPLINES_sp
Definition ReducedMSR.H:185