forked from idaholab/moose
/
PetscOutput.C
260 lines (221 loc) · 8.74 KB
/
PetscOutput.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
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html
// MOOSE includes
#include "PetscOutput.h"
#include "FEProblem.h"
#include "NonlinearSystem.h"
#include "libmesh/libmesh_common.h"
#include "libmesh/petsc_nonlinear_solver.h"
defineLegacyParams(PetscOutput);
InputParameters
PetscOutput::validParams()
{
InputParameters params = Output::validParams();
// Toggled for outputting nonlinear and linear residuals, only if we have PETSc
#ifdef LIBMESH_HAVE_PETSC
params.addParam<bool>(
"output_linear", false, "Specifies whether output occurs on each linear residual evaluation");
params.addParam<bool>("output_nonlinear",
false,
"Specifies whether output occurs on each nonlinear residual evaluation");
// **** DEPRECATED PARAMETERS ****
params.addDeprecatedParam<bool>(
"linear_residuals",
false,
"Specifies whether output occurs on each linear residual evaluation",
"Please use 'output_linear' to get this behavior.");
params.addDeprecatedParam<bool>(
"nonlinear_residuals",
false,
"Specifies whether output occurs on each nonlinear residual evaluation",
"Please use 'output_nonlinear' to get this behavior.");
// Psuedo time step divisors
params.addParam<Real>(
"nonlinear_residual_dt_divisor",
1000,
"Number of divisions applied to time step when outputting non-linear residuals");
params.addParam<Real>(
"linear_residual_dt_divisor",
1000,
"Number of divisions applied to time step when outputting linear residuals");
// Start times for residual output
params.addParam<Real>(
"linear_residual_start_time",
"Specifies a start time to begin output on each linear residual evaluation");
params.addParam<Real>(
"nonlinear_residual_start_time",
"Specifies a start time to begin output on each nonlinear residual evaluation");
// End time for residual output
/* Note, No default is given here so that in Peacock giant numbers do not show up by default, the
* defaults are set in the initialization list */
params.addParam<Real>("linear_residual_end_time",
"Specifies an end time to begin output on each linear residual evaluation");
params.addParam<Real>(
"nonlinear_residual_end_time",
"Specifies an end time to begin output on each nonlinear residual evaluation");
params.addParamNamesToGroup("linear_residuals nonlinear_residuals linear_residual_start_time "
"nonlinear_residual_start_time linear_residual_end_time "
"nonlinear_residual_end_time nonlinear_residual_dt_divisor "
"linear_residual_dt_divisor",
"PETSc");
#endif
return params;
}
PetscOutput::PetscOutput(const InputParameters & parameters)
: Output(parameters),
_nonlinear_iter(0),
_linear_iter(0),
_on_linear_residual(false),
_on_nonlinear_residual(false),
_nonlinear_dt_divisor(getParam<Real>("nonlinear_residual_dt_divisor")),
_linear_dt_divisor(getParam<Real>("linear_residual_dt_divisor")),
_nonlinear_start_time(-std::numeric_limits<Real>::max()),
_linear_start_time(-std::numeric_limits<Real>::max()),
_nonlinear_end_time(std::numeric_limits<Real>::max()),
_linear_end_time(std::numeric_limits<Real>::max())
{
// Output toggle support
if (getParam<bool>("output_linear"))
_execute_on.push_back("linear");
if (getParam<bool>("output_nonlinear"))
_execute_on.push_back("nonlinear");
// **** DEPRECATED PARAMETER SUPPORT ****
if (getParam<bool>("linear_residuals"))
_execute_on.push_back("linear");
if (getParam<bool>("nonlinear_residuals"))
_execute_on.push_back("nonlinear");
// Nonlinear residual start-time supplied by user
if (isParamValid("nonlinear_residual_start_time"))
{
_nonlinear_start_time = getParam<Real>("nonlinear_residual_start_time");
_execute_on.push_back("nonlinear");
}
// Nonlinear residual end-time supplied by user
if (isParamValid("nonlinear_residual_end_time"))
_nonlinear_end_time = getParam<Real>("nonlinear_residual_end_time");
// Linear residual start-time supplied by user
if (isParamValid("linear_residual_start_time"))
{
_linear_start_time = getParam<Real>("linear_residual_start_time");
_execute_on.push_back("linear");
}
// Linear residual end-time supplied by user
if (isParamValid("linear_residual_end_time"))
_linear_end_time = getParam<Real>("linear_residual_end_time");
}
void
PetscOutput::solveSetup()
{
// Only execute if PETSc exists
#ifdef LIBMESH_HAVE_PETSC
// Extract the non-linear and linear solvers from PETSc
NonlinearSystemBase & nl = _problem_ptr->getNonlinearSystemBase();
SNES snes = nl.getSNES();
KSP ksp;
SNESGetKSP(snes, &ksp);
// Update the pseudo times
_nonlinear_time = _time_old; // non-linear time starts with the previous time step
if (_dt != 0)
_nonlinear_dt = _dt / _nonlinear_dt_divisor; // set the pseudo non-linear timestep as fraction
// of real timestep for transient executioners
else
_nonlinear_dt = 1. / _nonlinear_dt_divisor; // set the pseudo non-linear timestep for steady
// executioners (here _dt==0)
_linear_dt = _nonlinear_dt / _linear_dt_divisor; // set the pseudo linear timestep
// Set the PETSc monitor functions
if (_execute_on.contains(EXEC_NONLINEAR) &&
(_time >= _nonlinear_start_time - _t_tol && _time <= _nonlinear_end_time + _t_tol))
{
PetscErrorCode ierr = SNESMonitorSet(snes, petscNonlinearOutput, this, PETSC_NULL);
CHKERRABORT(_communicator.get(), ierr);
}
if (_execute_on.contains(EXEC_LINEAR) &&
(_time >= _linear_start_time - _t_tol && _time <= _linear_end_time + _t_tol))
{
PetscErrorCode ierr = KSPMonitorSet(ksp, petscLinearOutput, this, PETSC_NULL);
CHKERRABORT(_communicator.get(), ierr);
}
#endif
}
// Only define the monitor functions if PETSc exists
#ifdef LIBMESH_HAVE_PETSC
PetscErrorCode
PetscOutput::petscNonlinearOutput(SNES, PetscInt its, PetscReal norm, void * void_ptr)
{
// Get the outputter object
PetscOutput * ptr = static_cast<PetscOutput *>(void_ptr);
// Update the pseudo times
ptr->_nonlinear_time += ptr->_nonlinear_dt;
ptr->_linear_time = ptr->_nonlinear_time;
// Set the current norm and iteration number
ptr->_norm = norm;
ptr->_nonlinear_iter = its;
// Set the flag indicating that output is occurring on the non-linear residual
ptr->_on_nonlinear_residual = true;
// Perform the output
ptr->outputStep(EXEC_NONLINEAR);
/**
* This is one of three locations where we explicitly flush the output buffers during a
* simulation:
* PetscOutput::petscNonlinearOutput()
* PetscOutput::petscLinearOutput()
* OutputWarehouse::outputStep()
*
* All other Console output _should_ be using newlines to avoid covering buffer errors
* and to avoid excessive I/O. This call is necessary. In the PETSc callback the
* context bypasses the OutputWarehouse.
*/
ptr->_app.getOutputWarehouse().flushConsoleBuffer();
// Reset the non-linear output flag and the simulation time
ptr->_on_nonlinear_residual = false;
// Done
return 0;
}
PetscErrorCode
PetscOutput::petscLinearOutput(KSP, PetscInt its, PetscReal norm, void * void_ptr)
{
// Get the Outputter object
PetscOutput * ptr = static_cast<PetscOutput *>(void_ptr);
// Update the pseudo time
ptr->_linear_time += ptr->_linear_dt;
// Set the current norm and iteration number
ptr->_norm = norm;
ptr->_linear_iter = its;
// Set the flag indicating that output is occurring on the non-linear residual
ptr->_on_linear_residual = true;
// Perform the output
ptr->outputStep(EXEC_LINEAR);
/**
* This is one of three locations where we explicitly flush the output buffers during a
* simulation:
* PetscOutput::petscNonlinearOutput()
* PetscOutput::petscLinearOutput()
* OutputWarehouse::outputStep()
*
* All other Console output _should_ be using newlines to avoid covering buffer errors
* and to avoid excessive I/O. This call is necessary. In the PETSc callback the
* context bypasses the OutputWarehouse.
*/
ptr->_app.getOutputWarehouse().flushConsoleBuffer();
// Reset the linear output flag and the simulation time
ptr->_on_linear_residual = false;
// Done
return 0;
}
#endif
Real
PetscOutput::time()
{
if (_on_nonlinear_residual)
return _nonlinear_time;
else if (_on_linear_residual)
return _linear_time;
else
return Output::time();
}