Skip to content

Commit 241dad8

Browse files
committed
updating nonlinear printing info (#385)
1 parent c64dc9c commit 241dad8

4 files changed

Lines changed: 60 additions & 57 deletions

File tree

applications/allenCahn_conserved/customPDE.h

Lines changed: 14 additions & 13 deletions
Original file line numberDiff line numberDiff line change
@@ -177,13 +177,13 @@ customPDE<dim, degree>::solveIncrement(bool skip_time_dependent)
177177
// parabolic and auxilary equations should also be here
178178
if (this->hasNonExplicitEquation)
179179
{
180-
bool nonlinear_it_converged = false;
181-
unsigned int nonlinear_it_index = 0;
180+
bool nonlinear_iteration_converged = false;
181+
unsigned int nonlinear_iteration_index = 0;
182182

183-
while (!nonlinear_it_converged)
183+
while (!nonlinear_iteration_converged)
184184
{
185-
nonlinear_it_converged = true; // Set to true here and will be set to false if
186-
// any variable isn't converged
185+
nonlinear_iteration_converged = true; // Set to true here and will be set to
186+
// false if any variable isn't converged
187187

188188
// Update residualSet for the non-explicitly updated variables
189189
// compute_nonexplicit_RHS()
@@ -214,16 +214,16 @@ customPDE<dim, degree>::solveIncrement(bool skip_time_dependent)
214214
this->pcout << buffer;
215215
}
216216

217-
nonlinear_it_converged =
218-
this->updateImplicitSolution(fieldIndex, nonlinear_it_index);
217+
nonlinear_iteration_converged =
218+
this->updateImplicitSolution(fieldIndex, nonlinear_iteration_index);
219219

220220
// Apply Boundary conditions
221221
this->applyBCs(fieldIndex);
222222
}
223223
else if (this->fields[fieldIndex].pdetype == AUXILIARY)
224224
{
225225
if (this->var_attributes.at(fieldIndex).is_nonlinear ||
226-
nonlinear_it_index == 0)
226+
nonlinear_iteration_index == 0)
227227
{
228228
// If the equation for this field is nonlinear, save the
229229
// old solution
@@ -281,18 +281,19 @@ customPDE<dim, degree>::solveIncrement(bool skip_time_dependent)
281281
{
282282
this->pcout << "Relative difference between nonlinear "
283283
"iterations: "
284-
<< diff << " " << nonlinear_it_index << " "
285-
<< this->currentIncrement << std::endl;
284+
<< diff << " " << nonlinear_iteration_index
285+
<< " " << this->currentIncrement
286+
<< std::endl;
286287
}
287288

288289
if (diff > MatrixFreePDE<dim, degree>::userInputs
289290
.nonlinear_solver_parameters.getToleranceValue(
290291
fieldIndex) &&
291-
nonlinear_it_index <
292+
nonlinear_iteration_index <
292293
MatrixFreePDE<dim, degree>::userInputs
293294
.nonlinear_solver_parameters.getMaxIterations())
294295
{
295-
nonlinear_it_converged = false;
296+
nonlinear_iteration_converged = false;
296297
}
297298
}
298299
else
@@ -318,7 +319,7 @@ customPDE<dim, degree>::solveIncrement(bool skip_time_dependent)
318319
}
319320
}
320321

321-
nonlinear_it_index++;
322+
nonlinear_iteration_index++;
322323
}
323324
}
324325

include/core/matrixFreePDE.h

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -288,7 +288,7 @@ class MatrixFreePDE : public Subscriptor
288288

289289
/*Method to compute an implicit timestep*/
290290
bool
291-
updateImplicitSolution(unsigned int fieldIndex, unsigned int nonlinear_it_index);
291+
updateImplicitSolution(unsigned int fieldIndex, unsigned int nonlinear_iteration_index);
292292

293293
/*Method to apply boundary conditions*/
294294
void

src/core/invM.cc

Lines changed: 0 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -42,20 +42,6 @@ MatrixFreePDE<dim, degree>::computeInvM()
4242
}
4343
}
4444

45-
// Check if invM has been found. If not, print an "error" message
46-
if (!invMscalarFound)
47-
{
48-
pcout << "matrixFreePDE.h: no PARABOLIC scalar field... hence setting "
49-
"parabolicScalarFieldIndex to 0 and marching ahead withn invM "
50-
"computation\n";
51-
}
52-
else if (!invMvectorFound)
53-
{
54-
pcout << "matrixFreePDE.h: no PARABOLIC vector field... hence setting "
55-
"parabolicVectorFieldIndex to 0 and marching ahead withn invM "
56-
"computation\n";
57-
}
58-
5945
// Initialize invM and clear its values
6046
matrixFreeObject.initialize_dof_vector(invMscalar, parabolicScalarFieldIndex);
6147
invMscalar = 0.0;

src/core/solvers/solveIncrement.cc

Lines changed: 45 additions & 29 deletions
Original file line numberDiff line numberDiff line change
@@ -66,19 +66,14 @@ MatrixFreePDE<dim, degree>::solveIncrement(bool skip_time_dependent)
6666
// parabolic and auxilary equations should also be here
6767
if (hasNonExplicitEquation)
6868
{
69-
bool nonlinear_it_converged = false;
70-
unsigned int nonlinear_it_index = 0;
69+
bool nonlinear_iteration_converged = false;
70+
unsigned int nonlinear_iteration_index = 0;
7171

72-
while (!nonlinear_it_converged)
72+
while (!nonlinear_iteration_converged)
7373
{
74-
nonlinear_it_converged = true; // Set to true here and will be set to false if
75-
// any variable isn't converged
74+
nonlinear_iteration_converged = true;
7675

7776
// Update residualSet for the non-explicitly updated variables
78-
// compute_nonexplicit_RHS()
79-
// Ideally, I'd just do this for the non-explicit variables, but for
80-
// now I'll do all of them this is a little redundant, but hopefully
81-
// not too terrible
8277
computeNonexplicitRHS();
8378

8479
for (const auto &[fieldIndex, variable] : var_attributes)
@@ -102,15 +97,15 @@ MatrixFreePDE<dim, degree>::solveIncrement(bool skip_time_dependent)
10297
pcout << buffer;
10398
}
10499

105-
nonlinear_it_converged =
106-
updateImplicitSolution(fieldIndex, nonlinear_it_index);
100+
nonlinear_iteration_converged =
101+
updateImplicitSolution(fieldIndex, nonlinear_iteration_index);
107102

108103
// Apply Boundary conditions
109104
applyBCs(fieldIndex);
110105
}
111106
else if (fields[fieldIndex].pdetype == AUXILIARY)
112107
{
113-
if (variable.is_nonlinear || nonlinear_it_index == 0)
108+
if (variable.is_nonlinear || nonlinear_iteration_index == 0)
114109
{
115110
// If the equation for this field is nonlinear, save the old
116111
// solution
@@ -164,19 +159,26 @@ MatrixFreePDE<dim, degree>::solveIncrement(bool skip_time_dependent)
164159
}
165160
if (currentIncrement % userInputs.skip_print_steps == 0)
166161
{
167-
pcout << "Relative difference between nonlinear "
168-
"iterations: "
169-
<< diff << " " << nonlinear_it_index << " "
170-
<< currentIncrement << "\n";
162+
snprintf(buffer,
163+
sizeof(buffer),
164+
" field '%2s' [nonlinear solve] current "
165+
"increment: %u, nonlinear "
166+
"iteration: "
167+
"%u, dU: %12.6e\n",
168+
fields[fieldIndex].name.c_str(),
169+
currentIncrement,
170+
nonlinear_iteration_index,
171+
diff);
172+
pcout << buffer;
171173
}
172174

173175
if (diff > userInputs.nonlinear_solver_parameters
174176
.getToleranceValue(fieldIndex) &&
175-
nonlinear_it_index <
177+
nonlinear_iteration_index <
176178
userInputs.nonlinear_solver_parameters
177179
.getMaxIterations())
178180
{
179-
nonlinear_it_converged = false;
181+
nonlinear_iteration_converged = false;
180182
}
181183
}
182184
else
@@ -201,7 +203,7 @@ MatrixFreePDE<dim, degree>::solveIncrement(bool skip_time_dependent)
201203
}
202204
}
203205

204-
nonlinear_it_index++;
206+
nonlinear_iteration_index++;
205207
}
206208
}
207209

@@ -286,12 +288,12 @@ MatrixFreePDE<dim, degree>::updateExplicitSolution(unsigned int fieldIndex)
286288
template <int dim, int degree>
287289
bool
288290
MatrixFreePDE<dim, degree>::updateImplicitSolution(unsigned int fieldIndex,
289-
unsigned int nonlinear_it_index)
291+
unsigned int nonlinear_iteration_index)
290292
{
291293
char buffer[200];
292294

293295
// Assume convergence criterion is met, unless otherwise proven later on.
294-
bool nonlinear_it_converged = true;
296+
bool nonlinear_iteration_converged = true;
295297

296298
// Apply Dirichlet BC's. This clears the residual where we want to apply Dirichlet BCs,
297299
// otherwise the solver sees a positive residual
@@ -469,18 +471,32 @@ MatrixFreePDE<dim, degree>::updateImplicitSolution(unsigned int fieldIndex,
469471
}
470472
if (currentIncrement % userInputs.skip_print_steps == 0)
471473
{
472-
pcout << "Relative difference between nonlinear "
473-
"iterations: "
474-
<< diff << " " << nonlinear_it_index << " " << currentIncrement
475-
<< "\n";
474+
snprintf(buffer,
475+
sizeof(buffer),
476+
" field '%2s' [nonlinear solve] current increment: %u, nonlinear "
477+
"iteration: "
478+
"%u, dU: %12.6e\n",
479+
fields[fieldIndex].name.c_str(),
480+
currentIncrement,
481+
nonlinear_iteration_index,
482+
diff);
483+
pcout << buffer;
476484
}
477485

478486
if (diff >
479487
userInputs.nonlinear_solver_parameters.getToleranceValue(fieldIndex) &&
480-
nonlinear_it_index <
488+
nonlinear_iteration_index <
481489
userInputs.nonlinear_solver_parameters.getMaxIterations())
482490
{
483-
nonlinear_it_converged = false;
491+
nonlinear_iteration_converged = false;
492+
}
493+
else if (diff >
494+
userInputs.nonlinear_solver_parameters.getToleranceValue(fieldIndex))
495+
{
496+
pcout << "\nWarning: nonlinear solver did not converge as "
497+
"per set tolerances. consider increasing the "
498+
"maximum number of iterations or decreasing the "
499+
"solver tolerance.\n";
484500
}
485501
}
486502
else
@@ -492,7 +508,7 @@ MatrixFreePDE<dim, degree>::updateImplicitSolution(unsigned int fieldIndex,
492508
}
493509
else
494510
{
495-
if (nonlinear_it_index == 0)
511+
if (nonlinear_iteration_index == 0)
496512
{
497513
if (fields[fieldIndex].type == SCALAR)
498514
{
@@ -532,5 +548,5 @@ MatrixFreePDE<dim, degree>::updateImplicitSolution(unsigned int fieldIndex,
532548
}
533549
}
534550

535-
return nonlinear_it_converged;
551+
return nonlinear_iteration_converged;
536552
}

0 commit comments

Comments
 (0)