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 newton_msr_fd(int Nx, int Ny, msrProblem& problem) :
43 newton_argument<double>(Nx, Ny),
45 Nphi_u(problem.NUmodes + problem.liftfield.size()),
46 Nphi_p(problem.NPmodes),
47 N_BC(problem.inletIndex.rows())
48 {}
49 int operator()(const Eigen::VectorXd& x, Eigen::VectorXd& fvec) const;
50 int df(const Eigen::VectorXd& x, Eigen::MatrixXd& fjac) const;
51
52 int Nphi_u;
53 int Nphi_p;
54 int N_BC;
55
56 scalar nu;
57
58 Eigen::VectorXd BC;
60};
61
62struct newton_msr_n: public newton_argument<double>
63{
64 public:
66 newton_msr_n(int Nx, int Ny, msrProblem& problem) :
67 newton_argument<double>(Nx, Ny),
69 Nphi_flux(problem.NFluxmodes),
70 Nphi_prec1(problem.NPrecmodes(0)),
71 Nphi_prec2(problem.NPrecmodes(1)),
72 Nphi_prec3(problem.NPrecmodes(2)),
73 Nphi_prec4(problem.NPrecmodes(3)),
74 Nphi_prec5(problem.NPrecmodes(4)),
75 Nphi_prec6(problem.NPrecmodes(5)),
76 Nphi_prec7(problem.NPrecmodes(6)),
77 Nphi_prec8(problem.NPrecmodes(7)),
78 Nphi_const(problem.NCmodes),
79 d_c(problem.NCmodes),
80 nsf_c(problem.NCmodes),
81 a_c(problem.NCmodes)
82 {}
83 int operator()(const Eigen::VectorXd& n, Eigen::VectorXd& fvecn) const;
84 int df(const Eigen::VectorXd& n, Eigen::MatrixXd& fjacn) const;
85
86 int Nphi_u;
97 //fluid-dynamics constants
98 scalar nu;
99 //neutronic constants
100 scalar iv;
101 scalar l1; //lambda-ith
102 scalar l2;
103 scalar l3;
104 scalar l4;
105 scalar l5;
106 scalar l6;
107 scalar l7;
108 scalar l8;
109 scalar b1; //beta-ith
110 scalar b2;
111 scalar b3;
112 scalar b4;
113 scalar b5;
114 scalar b6;
115 scalar b7;
116 scalar b8;
117 scalar btot;
118 scalar Sc;
119
120 Eigen::VectorXd d_c;
121 Eigen::VectorXd nsf_c;
122 Eigen::VectorXd a_c;
123 Eigen::VectorXd a_tmp; //to store u time coeff
124 std::vector<SPLINTER::RBFSpline*> SPLINES_d;
125 std::vector<SPLINTER::RBFSpline*> SPLINES_nsf;
126 std::vector<SPLINTER::RBFSpline*> SPLINES_a;
128};
129
130struct newton_msr_t: public newton_argument<double>
131{
132 public:
134 newton_msr_t(int Nx, int Ny, msrProblem& problem) :
135 newton_argument<double>(Nx, Ny),
137 Nphi_T(problem.NTmodes),
138 Nphi_dec1(problem.NDecmodes(0)),
139 Nphi_dec2(problem.NDecmodes(1)),
140 Nphi_dec3(problem.NDecmodes(2)),
141 N_BCt(problem.inletIndexT.rows()),
142 Nphi_const(problem.NCmodes),
143 v_c(problem.NCmodes),
144 sp_c(problem.NCmodes),
145 txs_c(problem.NCmodes)
146 {}
147 int operator()(const Eigen::VectorXd& n, Eigen::VectorXd& fvecn) const;
148 int df(const Eigen::VectorXd& n, Eigen::MatrixXd& fjacn) const;
149
156 int N_BCt;
158 //fluid-dynamics constants
159 scalar nu;
160 //thermal constants
161 scalar cp;
162 scalar dl1;
163 scalar dl2;
164 scalar dl3;
165 scalar db1;
166 scalar db2;
167 scalar db3;
168 scalar dbtot;
169 scalar Sc;
170 scalar Pr;
171
172 Eigen::VectorXd a_tmp;
173 Eigen::VectorXd c_tmp;
174 Eigen::VectorXd BCt;
175 Eigen::VectorXd v_c;
176 Eigen::VectorXd sp_c;
177 Eigen::VectorXd txs_c;
178 std::vector<SPLINTER::RBFSpline*> SPLINES_v;
179 std::vector<SPLINTER::RBFSpline*> SPLINES_sp;
180 std::vector<SPLINTER::RBFSpline*> SPLINES_TXS;
182};
183
184
185/*---------------------------------------------------------------------------*\
186 Class reducedProblem Declaration
187\*---------------------------------------------------------------------------*/
188
190{
191
192 public:
193
194 reducedMSR();
195 explicit reducedMSR(msrProblem& problem);
197
199 Eigen::VectorXd y_old;
200 Eigen::VectorXd w_old;
201 Eigen::VectorXd z_old;
202
204 Eigen::VectorXd y;
205 Eigen::VectorXd w;
206 Eigen::VectorXd z;
207
208
210 //fluid-dynamics constants
211 scalar nu;
212 scalar rho;
213 scalar Pr;
214 scalar Sc;
215 //neutronic constants
216 scalar d;
217 scalar m;
218 scalar nsf;
219 scalar a;
220 scalar Keff;
221 scalar iv;
222 scalar l1; //lambda-ith
223 scalar l2;
224 scalar l3;
225 scalar l4;
226 scalar l5;
227 scalar l6;
228 scalar l7;
229 scalar l8;
230 scalar b1; //beta-ith
231 scalar b2;
232 scalar b3;
233 scalar b4;
234 scalar b5;
235 scalar b6;
236 scalar b7;
237 scalar b8;
238 scalar btot;
239 //thermal constants
240 scalar cp;
241 scalar dl1;
242 scalar dl2;
243 scalar dl3;
244 scalar db1;
245 scalar db2;
246 scalar db3;
247 scalar dbtot;
248 scalar sp;
249
250
251
253 List < Eigen::MatrixXd> online_solution_fd;
254 List < Eigen::MatrixXd> online_solution_n;
255 List < Eigen::MatrixXd> online_solution_t;
256 List < Eigen::MatrixXd> online_solution_C;
257
259 PtrList<volVectorField> Umodes;
260 PtrList<volScalarField> Pmodes;
261 PtrList<volScalarField> Fluxmodes;
262 PtrList<volScalarField> Prec1modes;
263 PtrList<volScalarField> Prec2modes;
264 PtrList<volScalarField> Prec3modes;
265 PtrList<volScalarField> Prec4modes;
266 PtrList<volScalarField> Prec5modes;
267 PtrList<volScalarField> Prec6modes;
268 PtrList<volScalarField> Prec7modes;
269 PtrList<volScalarField> Prec8modes;
270 PtrList<volScalarField> Tmodes;
271 PtrList<volScalarField> Dec1modes;
272 PtrList<volScalarField> Dec2modes;
273 PtrList<volScalarField> Dec3modes;
274 PtrList<volScalarField> vmodes;
275 PtrList<volScalarField> Dmodes;
276 PtrList<volScalarField> NSFmodes;
277 PtrList<volScalarField> Amodes;
278 PtrList<volScalarField> SPmodes;
279 PtrList<volScalarField> TXSmodes;
280
281
283 PtrList<volVectorField> Usnapshots;
284 PtrList<volScalarField> Psnapshots;
285 PtrList<volScalarField> Fluxsnapshots;
286 PtrList<volScalarField> Prec1snapshots;
287 PtrList<volScalarField> Prec2snapshots;
288 PtrList<volScalarField> Prec3snapshots;
289 PtrList<volScalarField> Prec4snapshots;
290 PtrList<volScalarField> Prec5snapshots;
291 PtrList<volScalarField> Prec6snapshots;
292 PtrList<volScalarField> Prec7snapshots;
293 PtrList<volScalarField> Prec8snapshots;
294 PtrList<volScalarField> Tsnapshots;
295 PtrList<volScalarField> Dec1snapshots;
296 PtrList<volScalarField> Dec2snapshots;
297 PtrList<volScalarField> Dec3snapshots;
298 PtrList<volScalarField> vsnapshots;
299 PtrList<volScalarField> Dsnapshots;
300 PtrList<volScalarField> NSFsnapshots;
301 PtrList<volScalarField> Asnapshots;
302 PtrList<volScalarField> SPsnapshots;
303 PtrList<volScalarField> TXSsnapshots;
304
306 PtrList<volVectorField> UREC;
307 PtrList<volScalarField> PREC;
308 PtrList<volScalarField> FLUXREC;
309 PtrList<volScalarField> PREC1REC;
310 PtrList<volScalarField> PREC2REC;
311 PtrList<volScalarField> PREC3REC;
312 PtrList<volScalarField> PREC4REC;
313 PtrList<volScalarField> PREC5REC;
314 PtrList<volScalarField> PREC6REC;
315 PtrList<volScalarField> PREC7REC;
316 PtrList<volScalarField> PREC8REC;
317 PtrList<volScalarField> TREC;
318 PtrList<volScalarField> DEC1REC;
319 PtrList<volScalarField> DEC2REC;
320 PtrList<volScalarField> DEC3REC;
321 PtrList<volScalarField> POWERDENSREC;
322 PtrList<volScalarField> vREC;
323 PtrList<volScalarField> DREC;
324 PtrList<volScalarField> NSFREC;
325 PtrList<volScalarField> AREC;
326 PtrList<volScalarField> SPREC;
327 PtrList<volScalarField> TXSREC;
328
333
336
354
357 bool recall = false;
358
360 int N_BC;
361 int N_BCt;
362
365
366 //Methods:
367
368
369
371 void solveOnline(Eigen::MatrixXd vel_now, Eigen::MatrixXd temp_now,
372 Eigen::VectorXd mu_online);
373
375 void clearFields();
376
382 void reconstructAP(fileName folder = "./ITHACAOutput/online_rec",
383 int printevery = 1);
385 void reconstruct_fd(fileName folder = "./ITHACAOutput/online_rec",
386 int printevery = 1);
388 void reconstruct_n(fileName folder = "./ITHACAOutput/online_rec",
389 int printevery = 1);
391 void reconstruct_t(fileName folder = "./ITHACAOutput/online_rec",
392 int printevery = 1);
394 void reconstruct_C(fileName folder = "./ITHACAOutput/online_rec",
395 int printevery = 1);
396
397
398
399 protected:
400
401 //Methods:
402
405
406
407
408};
409
410
411#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
Template object created to solve non linear problems using the Eigen library.
PtrList< volScalarField > NSFREC
Definition ReducedMSR.H:324
scalar cp
Definition ReducedMSR.H:240
scalar Keff
Definition ReducedMSR.H:220
PtrList< volScalarField > Tsnapshots
Definition ReducedMSR.H:294
PtrList< volScalarField > DEC3REC
Definition ReducedMSR.H:320
msrProblem * problem
Pointer to the FOM problem.
Definition ReducedMSR.H:335
scalar Pr
Definition ReducedMSR.H:213
PtrList< volScalarField > vREC
Definition ReducedMSR.H:322
scalar b4
Definition ReducedMSR.H:233
scalar l6
Definition ReducedMSR.H:227
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:331
PtrList< volScalarField > DEC2REC
Definition ReducedMSR.H:319
PtrList< volScalarField > Dec1modes
Definition ReducedMSR.H:271
PtrList< volScalarField > Pmodes
Definition ReducedMSR.H:260
PtrList< volVectorField > Usnapshots
List of pointers to store the snapshots for each field.
Definition ReducedMSR.H:283
PtrList< volScalarField > Psnapshots
Definition ReducedMSR.H:284
PtrList< volScalarField > FLUXREC
Definition ReducedMSR.H:308
PtrList< volVectorField > Umodes
List of pointers to store the modes for each field.
Definition ReducedMSR.H:259
PtrList< volScalarField > Prec8snapshots
Definition ReducedMSR.H:293
PtrList< volScalarField > Dec2modes
Definition ReducedMSR.H:272
scalar nu
characteristic constants of the problem
Definition ReducedMSR.H:211
PtrList< volScalarField > TXSREC
Definition ReducedMSR.H:327
bool recall
boolean variable to check if the user wants to reconstruct all the three physics of the problem
Definition ReducedMSR.H:357
PtrList< volScalarField > PREC4REC
Definition ReducedMSR.H:312
PtrList< volScalarField > Prec2modes
Definition ReducedMSR.H:263
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:291
PtrList< volScalarField > PREC3REC
Definition ReducedMSR.H:311
PtrList< volScalarField > vsnapshots
Definition ReducedMSR.H:298
scalar l7
Definition ReducedMSR.H:228
PtrList< volVectorField > UREC
Recontructed fields.
Definition ReducedMSR.H:306
List< Eigen::MatrixXd > online_solution_n
Definition ReducedMSR.H:254
scalar b2
Definition ReducedMSR.H:231
PtrList< volScalarField > Prec7modes
Definition ReducedMSR.H:268
PtrList< volScalarField > Prec6modes
Definition ReducedMSR.H:267
scalar b1
Definition ReducedMSR.H:230
PtrList< volScalarField > SPmodes
Definition ReducedMSR.H:278
void reconstruct_t(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct thermal fields
Definition ReducedMSR.C:867
PtrList< volScalarField > Prec5snapshots
Definition ReducedMSR.H:290
PtrList< volScalarField > TXSsnapshots
Definition ReducedMSR.H:303
scalar b8
Definition ReducedMSR.H:237
scalar dl2
Definition ReducedMSR.H:242
Eigen::VectorXd z
Definition ReducedMSR.H:206
PtrList< volScalarField > Dec3snapshots
Definition ReducedMSR.H:297
PtrList< volScalarField > Prec3modes
Definition ReducedMSR.H:264
scalar Sc
Definition ReducedMSR.H:214
PtrList< volScalarField > Prec3snapshots
Definition ReducedMSR.H:288
PtrList< volScalarField > Prec7snapshots
Definition ReducedMSR.H:292
scalar db1
Definition ReducedMSR.H:244
scalar nsf
Definition ReducedMSR.H:218
scalar b3
Definition ReducedMSR.H:232
scalar db2
Definition ReducedMSR.H:245
newton_msr_fd newton_object_fd
Newton object used to solve the non linear problem.
Definition ReducedMSR.H:330
PtrList< volScalarField > Amodes
Definition ReducedMSR.H:277
void loadConstants(msrProblem *problem)
Method to load all the constants needed in the ROM from ///the FOM.
PtrList< volScalarField > DEC1REC
Definition ReducedMSR.H:318
PtrList< volScalarField > NSFmodes
Definition ReducedMSR.H:276
scalar iv
Definition ReducedMSR.H:221
PtrList< volScalarField > Fluxsnapshots
Definition ReducedMSR.H:285
void reconstruct_fd(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
reconstruct fluid-dynamics
Definition ReducedMSR.C:706
scalar dbtot
Definition ReducedMSR.H:247
int N_BC
Number of parametrized boundary conditions.
Definition ReducedMSR.H:360
scalar db3
Definition ReducedMSR.H:246
PtrList< volScalarField > Dec1snapshots
Definition ReducedMSR.H:295
PtrList< volScalarField > SPsnapshots
Definition ReducedMSR.H:302
scalar l8
Definition ReducedMSR.H:229
PtrList< volScalarField > Asnapshots
Definition ReducedMSR.H:301
PtrList< volScalarField > Prec4modes
Definition ReducedMSR.H:265
PtrList< volScalarField > PREC1REC
Definition ReducedMSR.H:309
scalar btot
Definition ReducedMSR.H:238
PtrList< volScalarField > Prec1snapshots
Definition ReducedMSR.H:286
PtrList< volScalarField > Dsnapshots
Definition ReducedMSR.H:299
PtrList< volScalarField > Prec8modes
Definition ReducedMSR.H:269
List< Eigen::MatrixXd > online_solution_fd
List of Eigen matrices to store the online solution.
Definition ReducedMSR.H:253
PtrList< volScalarField > POWERDENSREC
Definition ReducedMSR.H:321
PtrList< volScalarField > PREC6REC
Definition ReducedMSR.H:314
PtrList< volScalarField > TREC
Definition ReducedMSR.H:317
Eigen::VectorXd y_old
Vector to store the previous solution during the Newton procedure.
Definition ReducedMSR.H:199
PtrList< volScalarField > NSFsnapshots
Definition ReducedMSR.H:300
PtrList< volScalarField > TXSmodes
Definition ReducedMSR.H:279
PtrList< volScalarField > Dec3modes
Definition ReducedMSR.H:273
scalar dl1
Definition ReducedMSR.H:241
int count_online_solve
Counter to count the online solutions.
Definition ReducedMSR.H:364
PtrList< volScalarField > AREC
Definition ReducedMSR.H:325
PtrList< volScalarField > vmodes
Definition ReducedMSR.H:274
PtrList< volScalarField > Prec4snapshots
Definition ReducedMSR.H:289
PtrList< volScalarField > PREC2REC
Definition ReducedMSR.H:310
Eigen::VectorXd w
Definition ReducedMSR.H:205
PtrList< volScalarField > Tmodes
Definition ReducedMSR.H:270
PtrList< volScalarField > Dmodes
Definition ReducedMSR.H:275
int Nphi_u
Number of modes for each field.
Definition ReducedMSR.H:338
List< Eigen::MatrixXd > online_solution_t
Definition ReducedMSR.H:255
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:266
scalar rho
Definition ReducedMSR.H:212
PtrList< volScalarField > Prec2snapshots
Definition ReducedMSR.H:287
newton_msr_t newton_object_t
Definition ReducedMSR.H:332
PtrList< volScalarField > PREC
Definition ReducedMSR.H:307
PtrList< volScalarField > SPREC
Definition ReducedMSR.H:326
Eigen::VectorXd z_old
Definition ReducedMSR.H:201
scalar l4
Definition ReducedMSR.H:225
scalar sp
Definition ReducedMSR.H:248
void clearFields()
Method to clear all the fields of MSR (sets the size to zero)
scalar dl3
Definition ReducedMSR.H:243
PtrList< volScalarField > Dec2snapshots
Definition ReducedMSR.H:296
scalar b6
Definition ReducedMSR.H:235
PtrList< volScalarField > Prec1modes
Definition ReducedMSR.H:262
scalar l5
Definition ReducedMSR.H:226
scalar l2
Definition ReducedMSR.H:223
scalar b5
Definition ReducedMSR.H:234
Eigen::VectorXd y
Vector to store the solution during the Newton procedure.
Definition ReducedMSR.H:204
Eigen::VectorXd w_old
Definition ReducedMSR.H:200
scalar b7
Definition ReducedMSR.H:236
scalar l1
Definition ReducedMSR.H:222
PtrList< volScalarField > PREC5REC
Definition ReducedMSR.H:313
PtrList< volScalarField > Fluxmodes
Definition ReducedMSR.H:261
PtrList< volScalarField > PREC8REC
Definition ReducedMSR.H:316
PtrList< volScalarField > DREC
Definition ReducedMSR.H:323
List< Eigen::MatrixXd > online_solution_C
Definition ReducedMSR.H:256
PtrList< volScalarField > PREC7REC
Definition ReducedMSR.H:315
scalar l3
Definition ReducedMSR.H:224
Base class for the implementation of a reduced problem.
virtual void solveOnline()
Virtual Method to perform and online Solve.
Structure to implement a newton object for a stationary MSR problem.
Definition ReducedMSR.H:39
Eigen::VectorXd BC
Definition ReducedMSR.H:58
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:42
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
Definition ReducedMSR.C:185
msrProblem * problem
Definition ReducedMSR.H:59
Eigen::VectorXd a_tmp
Definition ReducedMSR.H:123
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
Definition ReducedMSR.C:414
std::vector< SPLINTER::RBFSpline * > SPLINES_d
Definition ReducedMSR.H:124
Eigen::VectorXd d_c
Definition ReducedMSR.H:120
newton_msr_n(int Nx, int Ny, msrProblem &problem)
Definition ReducedMSR.H:66
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const
Definition ReducedMSR.C:239
std::vector< SPLINTER::RBFSpline * > SPLINES_a
Definition ReducedMSR.H:126
Eigen::VectorXd nsf_c
Definition ReducedMSR.H:121
std::vector< SPLINTER::RBFSpline * > SPLINES_nsf
Definition ReducedMSR.H:125
Eigen::VectorXd a_c
Definition ReducedMSR.H:122
msrProblem * problem
Definition ReducedMSR.H:127
Eigen::VectorXd a_tmp
Definition ReducedMSR.H:172
newton_msr_t(int Nx, int Ny, msrProblem &problem)
Definition ReducedMSR.H:134
Eigen::VectorXd sp_c
Definition ReducedMSR.H:176
std::vector< SPLINTER::RBFSpline * > SPLINES_TXS
Definition ReducedMSR.H:180
Eigen::VectorXd txs_c
Definition ReducedMSR.H:177
Eigen::VectorXd c_tmp
Definition ReducedMSR.H:173
msrProblem * problem
Definition ReducedMSR.H:181
int df(const Eigen::VectorXd &n, Eigen::MatrixXd &fjacn) const
Definition ReducedMSR.C:512
std::vector< SPLINTER::RBFSpline * > SPLINES_v
Definition ReducedMSR.H:178
Eigen::VectorXd BCt
Definition ReducedMSR.H:174
int operator()(const Eigen::VectorXd &n, Eigen::VectorXd &fvecn) const
Definition ReducedMSR.C:422
Eigen::VectorXd v_c
Definition ReducedMSR.H:175
std::vector< SPLINTER::RBFSpline * > SPLINES_sp
Definition ReducedMSR.H:179