@@ -118,21 +118,32 @@ SNESSolver::SNESSolver(Options* opts)
118118 maxits((*options)["max_nonlinear_iterations"]
119119 .doc(" Maximum number of nonlinear iterations per SNES solve" )
120120 .withDefault(50 )),
121+ maxf((*options)["maxf"]
122+ .doc(" Maximum number of function evaluations per SNES solve" )
123+ .withDefault(10000 )),
121124 lower_its((*options)["lower_its"]
122125 .doc(" Iterations below which the next timestep is increased" )
123126 .withDefault(static_cast <int >(maxits * 0.5 ))),
124127 upper_its((*options)["upper_its"]
125128 .doc(" Iterations above which the next timestep is reduced" )
126129 .withDefault(static_cast <int >(maxits * 0.8 ))),
127130 timestep_factor_on_failure((*options)["timestep_factor_on_failure"]
128- .doc(" Multiply timestep on convergence failure" )
129- .withDefault(0.5 )),
130- timestep_factor_on_upper_its((*options)["timestep_factor_on_upper_its"]
131- .doc(" Multiply timestep if iterations exceed upper_its" )
132- .withDefault(0.9 )),
133- timestep_factor_on_lower_its((*options)["timestep_factor_on_lower_its"]
134- .doc(" Multiply timestep if iterations are below lower_its" )
135- .withDefault(1.4 )),
131+ .doc(" Multiply timestep on convergence failure" )
132+ .withDefault(0.5 )),
133+ timestep_factor_on_upper_its(
134+ (*options)["timestep_factor_on_upper_its"]
135+ .doc(" Multiply timestep if iterations exceed upper_its" )
136+ .withDefault(0.9 )),
137+ timestep_factor_on_lower_its(
138+ (*options)["timestep_factor_on_lower_its"]
139+ .doc(" Multiply timestep if iterations are below lower_its" )
140+ .withDefault(1.4 )),
141+ pid_controller(
142+ (*options)["pid_controller"].doc(" Use PID controller?" ).withDefault(false )),
143+ target_its((*options)["target_its"].doc(" Target snes iterations" ).withDefault(7 )),
144+ kP((*options)["kP"].doc(" Proportional PID parameter" ).withDefault(0.7 )),
145+ kI((*options)["kI"].doc(" Integral PID parameter" ).withDefault(0.3 )),
146+ kD((*options)["kD"].doc(" Derivative PID parameter" ).withDefault(0.2 )),
136147 diagnose(
137148 (*options)["diagnose"].doc(" Print additional diagnostics" ).withDefault(false )),
138149 diagnose_failures((*options)["diagnose_failures"]
@@ -340,6 +351,7 @@ int SNESSolver::init() {
340351
341352 // Set size of Matrix on each processor to nlocal x nlocal
342353 MatCreate (BoutComm::get (), &Jfd);
354+ MatSetOption (Jfd, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
343355 MatSetSizes (Jfd, nlocal, nlocal, PETSC_DETERMINE, PETSC_DETERMINE);
344356 MatSetFromOptions (Jfd);
345357 // Determine which row/columns of the matrix are locally owned
@@ -632,14 +644,21 @@ int SNESSolver::init() {
632644 // Note: If the 'Amat' Jacobian is matrix free, SNESComputeJacobian
633645 // always updates its reference 'u' vector every nonlinear iteration
634646 SNESSetLagJacobian (snes, lag_jacobian);
635- // Set Jacobian and preconditioner to persist across time steps
636- SNESSetLagJacobianPersists (snes, PETSC_TRUE);
637- SNESSetLagPreconditionerPersists (snes, PETSC_TRUE);
647+ if (pid_controller) {
648+ nl_its_prev = target_its;
649+ nl_its_prev2 = target_its;
650+ SNESSetLagJacobianPersists (snes, PETSC_FALSE);
651+ SNESSetLagPreconditionerPersists (snes, PETSC_FALSE);
652+ } else {
653+ // Set Jacobian and preconditioner to persist across time steps
654+ SNESSetLagJacobianPersists (snes, PETSC_TRUE);
655+ SNESSetLagPreconditionerPersists (snes, PETSC_TRUE);
656+ }
638657 SNESSetLagPreconditioner (snes, 1 ); // Rebuild when Jacobian is rebuilt
639658 }
640659
641660 // Set tolerances
642- SNESSetTolerances (snes, atol, rtol, stol, maxits, PETSC_DEFAULT );
661+ SNESSetTolerances (snes, atol, rtol, stol, maxits, maxf );
643662
644663 // Force SNES to take at least one nonlinear iteration.
645664 // This may prevent the solver from getting stuck in false steady state conditions
@@ -843,6 +862,10 @@ int SNESSolver::run() {
843862 VecAXPBY (snes_x, -beta, (1 . + beta), x1);
844863 }
845864
865+ if (pid_controller) {
866+ SNESSetLagJacobian (snes, lag_jacobian);
867+ }
868+
846869 // Run the solver
847870 PetscErrorCode ierr = SNESSolve (snes, nullptr , snes_x);
848871
@@ -900,6 +923,7 @@ int SNESSolver::run() {
900923 updateColoring ();
901924 jacobian_pruned = false ; // Reset flag. Will be set after pruning.
902925 }
926+
903927 if (saved_jacobian_lag == 0 ) {
904928 // This triggers a Jacobian recalculation
905929 SNESGetLagJacobian (snes, &saved_jacobian_lag);
@@ -1033,44 +1057,61 @@ int SNESSolver::run() {
10331057#endif // PETSC_VERSION_GE(3,20,0)
10341058
10351059 if (looping) {
1036- // Consider changing the timestep.
1037- // Note: The preconditioner depends on the timestep,
1038- // so if it is not recalculated the it will be less
1039- // effective.
1040- if ((nl_its <= lower_its) && (timestep < max_timestep)
1041- && (steps_since_snes_failure > 2 )) {
1042- // Increase timestep slightly
1043- timestep *= timestep_factor_on_lower_its;
1044-
1045- if (timestep > max_timestep) {
1046- timestep = max_timestep;
1060+
1061+ if (pid_controller) {
1062+ // Changing the timestep.
1063+ // Note: The preconditioner depends on the timestep,
1064+ // so we recalculate the jacobian and the preconditioner
1065+ // every time the timestep changes
1066+
1067+ timestep = pid (timestep, nl_its);
1068+
1069+ // NOTE(malamast): Do we really need this?
1070+ // Recompute Jacobian (for now)
1071+ if (saved_jacobian_lag == 0 ) {
1072+ SNESGetLagJacobian (snes, &saved_jacobian_lag);
1073+ SNESSetLagJacobian (snes, 1 );
10471074 }
10481075
1049- // Note: Setting the SNESJacobianFn to NULL retains
1050- // previously set evaluation function.
1051- //
1052- // The SNES Jacobian is a combination of the RHS Jacobian
1053- // and a factor involving the timestep.
1054- // Depends on equation_form
1055- // -> Probably call SNESSetJacobian(snes, Jfd, Jfd, NULL, fdcoloring);
1076+ } else {
1077+
1078+ // Consider changing the timestep.
1079+ // Note: The preconditioner depends on the timestep,
1080+ // so if it is not recalculated the it will be less
1081+ // effective.
1082+ if ((nl_its <= lower_its) && (timestep < max_timestep)
1083+ && (steps_since_snes_failure > 2 )) {
1084+ // Increase timestep slightly
1085+ timestep *= timestep_factor_on_lower_its;
1086+
1087+ timestep = std::min (timestep, max_timestep);
1088+
1089+ // Note: Setting the SNESJacobianFn to NULL retains
1090+ // previously set evaluation function.
1091+ //
1092+ // The SNES Jacobian is a combination of the RHS Jacobian
1093+ // and a factor involving the timestep.
1094+ // Depends on equation_form
1095+ // -> Probably call SNESSetJacobian(snes, Jfd, Jfd, NULL, fdcoloring);
1096+
1097+ if (static_cast <BoutReal>(lin_its) / nl_its > 4 ) {
1098+ // Recompute Jacobian (for now)
1099+ if (saved_jacobian_lag == 0 ) {
1100+ SNESGetLagJacobian (snes, &saved_jacobian_lag);
1101+ SNESSetLagJacobian (snes, 1 );
1102+ }
1103+ }
10561104
1057- if (static_cast <BoutReal>(lin_its) / nl_its > 4 ) {
1058- // Recompute Jacobian (for now)
1105+ } else if (nl_its >= upper_its) {
1106+ // Reduce timestep slightly
1107+ timestep *= timestep_factor_on_upper_its;
1108+
1109+ // Recompute Jacobian
10591110 if (saved_jacobian_lag == 0 ) {
10601111 SNESGetLagJacobian (snes, &saved_jacobian_lag);
10611112 SNESSetLagJacobian (snes, 1 );
10621113 }
10631114 }
1064-
1065- } else if (nl_its >= upper_its) {
1066- // Reduce timestep slightly
1067- timestep *= timestep_factor_on_upper_its;
1068-
1069- // Recompute Jacobian
1070- if (saved_jacobian_lag == 0 ) {
1071- SNESGetLagJacobian (snes, &saved_jacobian_lag);
1072- SNESSetLagJacobian (snes, 1 );
1073- }
10741115 }
10751116 }
10761117 snes_failures = 0 ;
@@ -1360,7 +1401,10 @@ void SNESSolver::updateColoring() {
13601401 // Re-calculate the coloring
13611402 MatColoring coloring = NULL ;
13621403 MatColoringCreate (Jfd, &coloring);
1363- MatColoringSetType (coloring, MATCOLORINGSL);
1404+ // MatColoringSetType(coloring, MATCOLORINGSL); // Serial algorithm. Better for smale-to-medium size problems.
1405+ MatColoringSetType (
1406+ coloring, MATCOLORINGGREEDY); // Parallel algorith. Better for large parallel runs
1407+ // MatColoringSetType(coloring, MATCOLORINGJP); // This didn't work
13641408 MatColoringSetFromOptions (coloring);
13651409
13661410 // Calculate new index sets
@@ -1386,4 +1430,25 @@ void SNESSolver::updateColoring() {
13861430 }
13871431}
13881432
1433+ BoutReal SNESSolver::pid (BoutReal timestep, int nl_its) {
1434+
1435+ /* ---------- multiplicative PID factors ---------- */
1436+ const BoutReal facP = std::pow (double (target_its) / double (nl_its), kP );
1437+ const BoutReal facI = std::pow (double (nl_its_prev) / double (nl_its), kI );
1438+ const BoutReal facD = std::pow (double (nl_its_prev) * double (nl_its_prev)
1439+ / double (nl_its) / double (nl_its_prev2),
1440+ kD );
1441+
1442+ // clamp groth factor to avoid huge changes
1443+ const BoutReal fac = std::clamp (facP * facI * facD, 0.2 , 5.0 );
1444+
1445+ /* ---------- update timestep and history ---------- */
1446+ const BoutReal dt_new = std::min (timestep * fac, max_timestep);
1447+
1448+ nl_its_prev2 = nl_its_prev;
1449+ nl_its_prev = nl_its;
1450+
1451+ return dt_new;
1452+ }
1453+
13891454#endif // BOUT_HAS_PETSC
0 commit comments