forked from idaholab/moose
/
TransientMultiApp.C
546 lines (431 loc) · 17.5 KB
/
TransientMultiApp.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
/****************************************************************/
/* DO NOT MODIFY THIS HEADER */
/* MOOSE - Multiphysics Object Oriented Simulation Environment */
/* */
/* (c) 2010 Battelle Energy Alliance, LLC */
/* ALL RIGHTS RESERVED */
/* */
/* Prepared by Battelle Energy Alliance, LLC */
/* Under Contract No. DE-AC07-05ID14517 */
/* With the U. S. Department of Energy */
/* */
/* See COPYRIGHT for full restrictions */
/****************************************************************/
// MOOSE includes
#include "TransientMultiApp.h"
#include "AllLocalDofIndicesThread.h"
#include "AuxiliarySystem.h"
#include "Console.h"
#include "LayeredSideFluxAverage.h"
#include "MooseMesh.h"
#include "Output.h"
#include "TimeStepper.h"
#include "Transient.h"
#include "libmesh/mesh_tools.h"
template <>
InputParameters
validParams<TransientMultiApp>()
{
InputParameters params = validParams<MultiApp>();
params += validParams<TransientInterface>();
params.addParam<bool>("sub_cycling",
false,
"Set to true to allow this MultiApp to take smaller "
"timesteps than the rest of the simulation. More "
"than one timestep will be performed for each "
"'master' timestep");
params.addParam<bool>("interpolate_transfers",
false,
"Only valid when sub_cycling. This allows "
"transferred values to be interpolated "
"over the time frame the MultiApp is "
"executing over when sub_cycling");
params.addParam<bool>("detect_steady_state",
false,
"If true then while sub_cycling a steady state check will be "
"done. In this mode output will only be done once the "
"MultiApp reaches the target time or steady state is reached");
params.addParam<Real>("steady_state_tol",
1e-8,
"The relative difference between the new "
"solution and the old solution that will be "
"considered to be at steady state");
params.addParam<bool>("output_sub_cycles", false, "If true then every sub-cycle will be output.");
params.addParam<bool>(
"print_sub_cycles", true, "Toggle the display of sub-cycles on the screen.");
params.addParam<unsigned int>(
"max_failures", 0, "Maximum number of solve failures tolerated while sub_cycling.");
params.addParam<bool>("tolerate_failure",
false,
"If true this MultiApp won't participate in dt "
"decisions and will always be fast-forwarded to "
"the current time.");
params.addParam<bool>(
"catch_up",
false,
"If true this will allow failed solves to attempt to 'catch up' using smaller timesteps.");
params.addParam<Real>("max_catch_up_steps",
2,
"Maximum number of steps to allow an app to take "
"when trying to catch back up after a failed "
"solve.");
return params;
}
TransientMultiApp::TransientMultiApp(const InputParameters & parameters)
: MultiApp(parameters),
_sub_cycling(getParam<bool>("sub_cycling")),
_interpolate_transfers(getParam<bool>("interpolate_transfers")),
_detect_steady_state(getParam<bool>("detect_steady_state")),
_steady_state_tol(getParam<Real>("steady_state_tol")),
_output_sub_cycles(getParam<bool>("output_sub_cycles")),
_max_failures(getParam<unsigned int>("max_failures")),
_tolerate_failure(getParam<bool>("tolerate_failure")),
_failures(0),
_catch_up(getParam<bool>("catch_up")),
_max_catch_up_steps(getParam<Real>("max_catch_up_steps")),
_first(declareRecoverableData<bool>("first", true)),
_auto_advance(false),
_print_sub_cycles(getParam<bool>("print_sub_cycles"))
{
// Transfer interpolation only makes sense for sub-cycling solves
if (_interpolate_transfers && !_sub_cycling)
mooseError("MultiApp ",
name(),
" is set to interpolate_transfers but is not sub_cycling! That is not valid!");
// Subcycling overrides catch up, we don't want to confuse users by allowing them to set both.
if (_sub_cycling && _catch_up)
mooseError("MultiApp ",
name(),
" sub_cycling and catch_up cannot both be set to true simultaneously.");
}
NumericVector<Number> &
TransientMultiApp::appTransferVector(unsigned int app, std::string var_name)
{
if (std::find(_transferred_vars.begin(), _transferred_vars.end(), var_name) ==
_transferred_vars.end())
_transferred_vars.push_back(var_name);
if (_interpolate_transfers)
return appProblemBase(app).getAuxiliarySystem().system().get_vector("transfer");
return appProblemBase(app).getAuxiliarySystem().solution();
}
void
TransientMultiApp::initialSetup()
{
MultiApp::initialSetup();
if (!_has_an_app)
return;
Moose::ScopedCommSwapper swapper(_my_comm);
if (_has_an_app)
{
_transient_executioners.resize(_my_num_apps);
// Grab Transient Executioners from each app
for (unsigned int i = 0; i < _my_num_apps; i++)
setupApp(i);
}
}
bool
TransientMultiApp::solveStep(Real dt, Real target_time, bool auto_advance)
{
if (!_has_an_app)
return true;
_auto_advance = auto_advance;
_console << "Solving MultiApp " << name() << std::endl;
// "target_time" must always be in global time
target_time += _app.getGlobalTimeOffset();
Moose::ScopedCommSwapper swapper(_my_comm);
bool return_value = true;
// Make sure we swap back the communicator regardless of how this routine is exited
try
{
int rank;
int ierr;
ierr = MPI_Comm_rank(_orig_comm, &rank);
mooseCheckMPIErr(ierr);
for (unsigned int i = 0; i < _my_num_apps; i++)
{
FEProblemBase & problem = appProblemBase(_first_local_app + i);
Transient * ex = _transient_executioners[i];
// The App might have a different local time from the rest of the problem
Real app_time_offset = _apps[i]->getGlobalTimeOffset();
if ((ex->getTime() + app_time_offset) + 2e-14 >=
target_time) // Maybe this MultiApp was already solved
continue;
if (_sub_cycling)
{
Real time_old = ex->getTime() + app_time_offset;
if (_interpolate_transfers)
{
AuxiliarySystem & aux_system = problem.getAuxiliarySystem();
System & libmesh_aux_system = aux_system.system();
NumericVector<Number> & solution = *libmesh_aux_system.solution;
NumericVector<Number> & transfer_old = libmesh_aux_system.get_vector("transfer_old");
solution.close();
// Save off the current auxiliary solution
transfer_old = solution;
transfer_old.close();
// Snag all of the local dof indices for all of these variables
AllLocalDofIndicesThread aldit(libmesh_aux_system, _transferred_vars);
ConstElemRange & elem_range = *problem.mesh().getActiveLocalElementRange();
Threads::parallel_reduce(elem_range, aldit);
_transferred_dofs = aldit._all_dof_indices;
}
// Disable/enable output for sub cycling
problem.allowOutput(_output_sub_cycles); // disables all outputs, including console
problem.allowOutput<Console>(_print_sub_cycles); // re-enables Console to print, if desired
ex->setTargetTime(target_time - app_time_offset);
// unsigned int failures = 0;
bool at_steady = false;
if (_first && !_app.isRecovering())
problem.advanceState();
bool local_first = _first;
// Now do all of the solves we need
while ((!at_steady && ex->getTime() + app_time_offset + 2e-14 < target_time) ||
!ex->lastSolveConverged())
{
if (local_first != true)
ex->incrementStepOrReject();
local_first = false;
ex->preStep();
ex->computeDT();
if (_interpolate_transfers)
{
// See what time this executioner is going to go to.
Real future_time = ex->getTime() + app_time_offset + ex->getDT();
// How far along we are towards the target time:
Real step_percent = (future_time - time_old) / (target_time - time_old);
Real one_minus_step_percent = 1.0 - step_percent;
// Do the interpolation for each variable that was transferred to
FEProblemBase & problem = appProblemBase(_first_local_app + i);
AuxiliarySystem & aux_system = problem.getAuxiliarySystem();
System & libmesh_aux_system = aux_system.system();
NumericVector<Number> & solution = *libmesh_aux_system.solution;
NumericVector<Number> & transfer = libmesh_aux_system.get_vector("transfer");
NumericVector<Number> & transfer_old = libmesh_aux_system.get_vector("transfer_old");
solution.close(); // Just to be sure
transfer.close();
transfer_old.close();
for (const auto & dof : _transferred_dofs)
{
solution.set(dof,
(transfer_old(dof) * one_minus_step_percent) +
(transfer(dof) * step_percent));
// solution.set(dof, transfer_old(dof));
// solution.set(dof, transfer(dof));
// solution.set(dof, 1);
}
solution.close();
}
ex->takeStep();
bool converged = ex->lastSolveConverged();
if (!converged)
{
mooseWarning(
"While sub_cycling ", name(), _first_local_app + i, " failed to converge!\n");
_failures++;
if (_failures > _max_failures)
{
std::stringstream oss;
oss << "While sub_cycling " << name() << _first_local_app << i << " REALLY failed!";
throw MultiAppSolveFailure(oss.str());
}
}
Real solution_change_norm = ex->getSolutionChangeNorm();
if (_detect_steady_state)
_console << "Solution change norm: " << solution_change_norm << std::endl;
if (converged && _detect_steady_state && solution_change_norm < _steady_state_tol)
{
_console << "Detected Steady State! Fast-forwarding to " << target_time << std::endl;
at_steady = true;
// Indicate that the next output call (occurs in ex->endStep()) should output,
// regardless of intervals etc...
problem.forceOutput();
// Clean up the end
ex->endStep(target_time - app_time_offset);
ex->postStep();
}
else
{
ex->endStep();
ex->postStep();
}
}
// If we were looking for a steady state, but didn't reach one, we still need to output one
// more time, regardless of interval
if (!at_steady)
problem.outputStep(EXEC_FORCED);
} // sub_cycling
else if (_tolerate_failure)
{
ex->takeStep(dt);
ex->endStep(target_time - app_time_offset);
ex->postStep();
}
else
{
_console << "Solving Normal Step!" << std::endl;
if (_first && !_app.isRecovering())
problem.advanceState();
if (auto_advance)
problem.allowOutput(true);
ex->takeStep(dt);
if (auto_advance)
{
ex->endStep();
ex->postStep();
if (!ex->lastSolveConverged())
{
mooseWarning(name(), _first_local_app + i, " failed to converge!\n");
if (_catch_up)
{
_console << "Starting Catch Up!" << std::endl;
bool caught_up = false;
unsigned int catch_up_step = 0;
Real catch_up_dt = dt / 2;
while (!caught_up && catch_up_step < _max_catch_up_steps)
{
Moose::err << "Solving " << name() << "catch up step " << catch_up_step
<< std::endl;
ex->incrementStepOrReject();
ex->computeDT();
ex->takeStep(catch_up_dt); // Cut the timestep in half to try two half-step solves
if (ex->lastSolveConverged())
{
if (ex->getTime() + app_time_offset +
ex->timestepTol() * std::abs(ex->getTime()) >=
target_time)
{
problem.outputStep(EXEC_FORCED);
caught_up = true;
}
}
else
catch_up_dt /= 2.0;
ex->endStep();
ex->postStep();
catch_up_step++;
}
if (!caught_up)
throw MultiAppSolveFailure(name() + " Failed to catch up!\n");
}
} // if (!ex->lastSolveConverged())
} // if (auto_advance)
else if (!ex->lastSolveConverged())
throw MultiAppSolveFailure(name() + " failed to converge");
}
// Re-enable all output (it may of been disabled by sub-cycling)
problem.allowOutput(true);
}
_first = false;
_console << "Successfully Solved MultiApp " << name() << "." << std::endl;
}
catch (MultiAppSolveFailure & e)
{
mooseWarning(e.what());
_console << "Failed to Solve MultiApp " << name() << ", attempting to recover." << std::endl;
return_value = false;
}
_transferred_vars.clear();
return return_value;
}
void
TransientMultiApp::incrementTStep()
{
if (!_sub_cycling)
{
for (unsigned int i = 0; i < _my_num_apps; i++)
{
Transient * ex = _transient_executioners[i];
ex->incrementStepOrReject();
}
}
}
void
TransientMultiApp::finishStep()
{
if (!_sub_cycling)
{
for (unsigned int i = 0; i < _my_num_apps; i++)
{
Transient * ex = _transient_executioners[i];
ex->endStep();
ex->postStep();
}
}
}
bool
TransientMultiApp::needsRestoration()
{
return _sub_cycling || _catch_up || _auto_advance || _tolerate_failure || _detect_steady_state;
}
Real
TransientMultiApp::computeDT()
{
if (_sub_cycling) // Bow out of the timestep selection dance
return std::numeric_limits<Real>::max();
Real smallest_dt = std::numeric_limits<Real>::max();
if (_has_an_app)
{
Moose::ScopedCommSwapper swapper(_my_comm);
for (unsigned int i = 0; i < _my_num_apps; i++)
{
Transient * ex = _transient_executioners[i];
ex->computeDT();
Real dt = ex->getDT();
smallest_dt = std::min(dt, smallest_dt);
}
}
if (_tolerate_failure) // Bow out of the timestep selection dance, we do this down here because we
// need to call computeConstrainedDT at least once for these
// executioners...
return std::numeric_limits<Real>::max();
_communicator.min(smallest_dt);
return smallest_dt;
}
void TransientMultiApp::resetApp(
unsigned int global_app,
Real /*time*/) // FIXME: Note that we are passing in time but also grabbing it below
{
if (hasLocalApp(global_app))
{
unsigned int local_app = globalAppToLocal(global_app);
// Grab the current time the App is at so we can start the new one at the same place
Real time =
_transient_executioners[local_app]->getTime() + _apps[local_app]->getGlobalTimeOffset();
// Reset the Multiapp
MultiApp::resetApp(global_app, time);
Moose::ScopedCommSwapper swapper(_my_comm);
// Setup the app, disable the output so that the initial condition does not output
// When an app is reset the initial condition was effectively already output before reset
FEProblemBase & problem = appProblemBase(local_app);
problem.allowOutput(false);
setupApp(local_app, time);
problem.allowOutput(true);
}
}
void TransientMultiApp::setupApp(unsigned int i, Real /*time*/) // FIXME: Should we be passing time?
{
auto & app = _apps[i];
Transient * ex = dynamic_cast<Transient *>(app->getExecutioner());
if (!ex)
mooseError("MultiApp ", name(), " is not using a Transient Executioner!");
// Get the FEProblemBase for the current MultiApp
FEProblemBase & problem = appProblemBase(_first_local_app + i);
// Update the file numbers for the outputs from the parent application
app->getOutputWarehouse().setFileNumbers(_app.getOutputFileNumbers());
// Call initialization method of Executioner (Note, this preforms the output of the initial time
// step, if desired)
ex->init();
if (_interpolate_transfers)
{
AuxiliarySystem & aux_system = problem.getAuxiliarySystem();
System & libmesh_aux_system = aux_system.system();
// We'll store a copy of the auxiliary system's solution at the old time in here
libmesh_aux_system.add_vector("transfer_old", false);
// This will be where we'll transfer the value to for the "target" time
libmesh_aux_system.add_vector("transfer", false);
}
ex->preExecute();
if (!_app.isRecovering())
problem.advanceState();
_transient_executioners[i] = ex;
}