• Home
  • Features
  • Pricing
  • Docs
  • Announcements
  • Sign In

gyrokinetics / gs2 / 1991098022

19 Aug 2025 09:13AM UTC coverage: 10.671% (+0.007%) from 10.664%
1991098022

push

gitlab-ci

David Dickinson
Merged in experimental/add_override_flags_to_configs (pull request #1138)

4766 of 44662 relevant lines covered (10.67%)

124981.59 hits per line

Source File
Press 'n' to go to next uncovered line, 'b' for previous

0.0
/src/split_nonlinear_terms.fpp
1
!> A module to deal with advancing the nonlinear term
2
!> separately from the linear terms.
3
module split_nonlinear_terms
4
  use abstract_config, only: abstract_config_type, CONFIG_MAX_NAME_LEN
5
  use rk_schemes, only: rk_scheme_type
6
  implicit none
7
  private
8

9
  public :: init_split_nonlinear_terms, finish_split_nonlinear_terms, set_split_nonlinear_terms_config
10
  public :: advance_nonlinear_term, advance_nonadiabatic_dfn, strang_split
11
  public :: split_nonlinear_terms_config_type, get_split_nonlinear_terms_config
12

13
  logical :: initialized = .false.
14

15
  !> Statistics for the number of steps taken in the integrators
16
  integer :: nfailed_steps = 0, nrhs_calls = 0, nsuccessful_steps = 0
17

18
  !> Interface to describe the call signatures for routines to calculate the
19
  !> field from gnew (invert_field_func) and to calculate the source term used
20
  !> in the advance_nonlinear_term methods. Used to allow us to pass in these
21
  !> methods, giving us a way to test the integrators on different problems.
22
  interface
23
     subroutine invert_field_func (g_in, phi, apar, bpar, gf_lo, local_only)
24
       implicit none
25
       complex, dimension (:,:,:), intent(in) :: g_in
26
       complex, dimension (:,:,:), intent (out) :: phi, apar, bpar
27
       logical, intent(in), optional :: gf_lo, local_only
28
     end subroutine invert_field_func
29

30
     subroutine source_term_func(g_in, g1, phi, apar, bpar, max_vel_local, &
31
          need_to_adjust, calculate_cfl_limit)
32
       implicit none
33
       complex, dimension (:,:,:), intent (in) :: g_in
34
       complex, dimension (:,:,:), intent (out) :: g1
35
       complex, dimension (:,:,:), intent (in) :: phi, apar, bpar
36
       real, intent(out) :: max_vel_local
37
       logical, intent(in) :: need_to_adjust
38
       logical, intent(in), optional :: calculate_cfl_limit
39
     end subroutine source_term_func
40
  end interface
41

42
  ! Related to inputs
43
  integer :: split_method_switch
44
  integer, parameter :: split_method_ab3 = 1, split_method_backwards_euler = 2, &
45
       split_method_rk = 3, split_method_picard = 4
46
  type(rk_scheme_type) :: the_rk_scheme
47
  logical :: show_statistics, advance_nonadiabatic_dfn, strang_split
48
  real :: convergence_tolerance, time_step_safety_factor
49
  real :: relative_tolerance, absolute_tolerance
50

51
  !> Used to represent the input configuration of split_nonlinear_terms
52
  type, extends(abstract_config_type) :: split_nonlinear_terms_config_type
53
     ! namelist : split_nonlinear_terms_knobs
54
     ! indexed : false
55
     !> The absolute tolerance used in error controlled schemes. Attempt to adjust
56
     !> time step to keep error below this tolerance.
57
     real :: absolute_tolerance = 1.0e-3
58
     !> If true then nonlinear integrator advances the nonadiabatic dfn, h,
59
     !> rather than the modified distribution function, g.
60
     logical :: advance_nonadiabatic_dfn = .false.
61
     !> The tolerance used in convergence checks. Currently only used
62
     !> for halting the fixed-point iteration in the beuler method
63
     real :: convergence_tolerance = 1.0e-1
64
     !> The relative tolerance used in error controlled schemes. Attempt to adjust
65
     !> time step to keep error below this tolerance.
66
     real :: relative_tolerance = 1.0e-2
67
     !> Choose which specific RK scheme to use in the RK method.
68
     !> Generally if the tolerances are small a higher order scheme
69
     !> will be more effective than low order. At loose tolerances
70
     !> low order schemes are generally more efficient
71
     !> Valid options include:
72
     !>
73
     !>  - 'cashkarp' -- Cash-Karp scheme of order 4(5)
74
     !>  - 'default' -- Same as heun
75
     !>  - 'heun' -- Heun scheme of order 1(2)
76
     !>
77
     !> See [[rk_schemes]] for other options.
78
     character(len = 20) :: rk_method = 'default'
79
     !> If true then reports the number of right hand side (NL term) calculations,
80
     !> the number of failed internal steps and the number of successful steps at the
81
     !> end of each linear size step.
82
     logical :: show_statistics = .false.
83
     !> What algorithm do we use to advance the nonlinear term if
84
     !> split_nonlinear is true.
85
     !> Valid options include:
86
     !>
87
     !>  - 'AB3' -- Adams-Bashforth (up to) 3rd order. CLF controlled time step.
88
     !>  - 'beuler' -- Backwards Euler. Relative and absolute error controlled time step.
89
     !>            Newton iteration controlled by convergence_tolerance.
90
     !>  - 'default' -- The same as 'RK'.
91
     !>  - 'RK' -- RK scheme with embedded error estimate. Exact scheme can be switched out
92
     !>            using rk_method switch. Relative and absolute error controlled time step.
93
     !>  - 'picard' -- Use simple Picard iteration. Not generally recommended but available
94
     !>            for experimentation.
95
     !>
96
     character(len = 20) :: split_method = 'default'
97
     !> If true then indicates that we want to use a strang split approach where
98
     !> we advance the NL operator by a half step, linear by a full step and then
99
     !> NL by another half step. This _should_ be second order accurate whilst the
100
     !> regular split is only first order accurate.
101
     logical :: strang_split = .false.
102
     !> Multiplies the new time step calculated from either the error or
103
     !> cfl limits. Smaller values are more conservative but may help avoid
104
     !> repeated violation of the error/cfl limits and hence reduce the
105
     !> number of failed steps.
106
     real :: time_step_safety_factor = 0.9
107
#include "split_nonlinear_terms_overrides_and_bound_auto_gen.inc"
108
  end type split_nonlinear_terms_config_type
109

110
  type(split_nonlinear_terms_config_type) :: split_nonlinear_terms_config
111

112
contains
113

114
  !> Initialises the split_nonlinear_terms module. Primarily just reading input
115
  subroutine init_split_nonlinear_terms(split_nonlinear_terms_config_in)
×
116
    implicit none
117
    type(split_nonlinear_terms_config_type), intent(in), optional :: split_nonlinear_terms_config_in
118

119
    if (initialized) return
×
120
    initialized = .true.
×
121

122
    call read_parameters(split_nonlinear_terms_config_in)
×
123

124
  end subroutine init_split_nonlinear_terms
125

126
  !> Read the split_nonlinear_terms namelist
127
  subroutine read_parameters(split_nonlinear_terms_config_in)
×
128
    use file_utils, only: error_unit
129
    use text_options, only: text_option, get_option_value
130
    use rk_schemes, only: get_rk_schemes_as_text_options, get_rk_scheme_by_id
131
    implicit none
132
    type(split_nonlinear_terms_config_type), intent(in), optional :: split_nonlinear_terms_config_in
133
    type(text_option), dimension (*), parameter :: split_method_opts = &
134
         [ &
135
         text_option('default', split_method_rk), &
136
         text_option('AB3', split_method_ab3), &
137
         text_option('RK', split_method_rk), &
138
         text_option('beuler', split_method_backwards_euler),  &
139
         text_option('picard', split_method_picard)  &
140
         ]
141
    character(len = 20) :: rk_method, split_method
142
    integer :: ierr, rk_method_switch
143

144
    if (present(split_nonlinear_terms_config_in)) split_nonlinear_terms_config = split_nonlinear_terms_config_in
×
145

146
    call split_nonlinear_terms_config%init(name = 'split_nonlinear_terms_knobs', requires_index = .false.)
×
147

148
    ! Copy out internal values into module level parameters
149
    associate(self => split_nonlinear_terms_config)
150
#include "split_nonlinear_terms_copy_out_auto_gen.inc"
151
    end associate
152

153
    ierr = error_unit()
×
154

155
    call get_option_value &
156
         (split_method, split_method_opts, split_method_switch, &
157
         ierr, "split_method in split_nonlinear_terms_knobs",.true.)
×
158

159
    call get_option_value &
160
         (rk_method, get_rk_schemes_as_text_options(), rk_method_switch, &
161
         ierr, "rk_method in split_nonlinear_terms_knobs",.true.)
×
162
    the_rk_scheme = get_rk_scheme_by_id(rk_method_switch)
×
163
  end subroutine read_parameters
×
164

165
  !> Reset the module, freeing memory etc.
166
  subroutine finish_split_nonlinear_terms
×
167
    implicit none
168
    initialized = .false.
×
169
    call split_nonlinear_terms_config%reset()
×
170
  end subroutine finish_split_nonlinear_terms
×
171

172
  !> Reset the step statistics
173
  subroutine reset_step_statistics
×
174
    nfailed_steps = 0 ; nrhs_calls = 0 ; nsuccessful_steps = 0
×
175
  end subroutine reset_step_statistics
×
176

177
  subroutine report_step_statistics(istep, unit_in)
×
178
    use iso_fortran_env, only: output_unit
179
    use optionals, only: get_option_with_default
180
    implicit none
181
    integer, intent(in) :: istep
182
    integer, intent(in), optional :: unit_in
183
    integer :: unit
184
    unit = get_option_with_default(unit_in, output_unit)
×
185
    write(unit,'("Iteration : ",I0," Number of internal steps : ",I0," Number of failed steps : ",I0," Number of RHS calls : ",I0)') &
186
         istep, nsuccessful_steps, nfailed_steps, nrhs_calls
×
187
  end subroutine report_step_statistics
×
188

189
  !> Advances dg/dt = NL(g,chi) from g_state, chi from t -> t+dt by calling
190
  !> specific requested method.
191
  !> On output g_state contains new g at next step, chi is left unchanged
192
  subroutine advance_nonlinear_term(g_state, istep, phi, apar, bpar, &
×
193
       dt, fields_func, source_func)
194
    use dist_fn_arrays, only: g_adjust, to_g_gs2, from_g_gs2
195
    use mp, only: mp_abort, get_mp_times, proc0
196
    use job_manage, only: time_message
197
    use nonlinear_terms, only: time_add_explicit_terms, time_add_explicit_terms_mpi
198
    use run_parameters, only: has_phi, has_apar, has_bpar
199
    use array_utils, only: copy
200
    implicit none
201
    integer, intent(in) :: istep
202
    complex, dimension (:,:,:), intent(in out) :: g_state
203
    complex, dimension (:,:,:), intent(in) :: phi, apar, bpar
204
    real, intent(in) :: dt
205
    complex, dimension (:,:,:), allocatable :: phinew, aparnew, bparnew
×
206
    real :: mp_total_after, mp_total
207
    procedure(invert_field_func) :: fields_func
208
    procedure(source_term_func) :: source_func
209

210
    call time_message(.false., time_add_explicit_terms, 'Explicit terms')
×
211
    call get_mp_times(total_time = mp_total)
×
212
    call reset_step_statistics
×
213

214
    ! Copy current fields to local copies
215
    allocate(phinew, mold = phi)
×
216
    allocate(aparnew, mold = apar)
×
217
    allocate(bparnew, mold = bpar)
×
218
    if (has_phi) call copy(phi, phinew)
×
219
    if (has_apar) call copy(apar, aparnew)
×
220
    if (has_bpar) call copy(bpar, bparnew)
×
221

222
    ! If we want to advance h rather than g then adjust
223
    ! here to get h. Note we rely on the call site also
224
    ! providing the correct fields_func to calculate fields
225
    ! from h rather than g. A different option would be
226
    ! to require the call site to always pass a way to calculate
227
    ! from g and a way from h and we then decide here which to use.
228
    if (advance_nonadiabatic_dfn) then
×
229
       call g_adjust(g_state, phinew, bparnew, direction = from_g_gs2)
×
230
    end if
231

232
    ! Call requested method for advancing nonlinear term
233
    select case(split_method_switch)
×
234
    case(split_method_ab3)
235
       call advance_nonlinear_term_ab3(g_state, phinew, aparnew, bparnew, dt, &
×
236
            fields_func, source_func)
×
237
    case(split_method_backwards_euler)
238
       call advance_nonlinear_term_beuler(g_state, phinew, aparnew, bparnew, dt, &
×
239
            fields_func, source_func)
×
240
    case(split_method_rk)
241
       call advance_nonlinear_term_rk(g_state, phinew, aparnew, bparnew, dt, &
×
242
            fields_func, source_func)
×
243
    case(split_method_picard)
244
       call advance_nonlinear_term_picard(g_state, phinew, aparnew, bparnew, dt, &
×
245
            fields_func, source_func)
×
246
    case default
247
       call mp_abort("Invalid split_method_switch", .true.)
×
248
    end select
249

250
    if (proc0 .and. show_statistics) call report_step_statistics(istep)
×
251

252
    ! If we advanced h then we need to reconstruct g from this
253
    if (advance_nonadiabatic_dfn) then
×
254
       call g_adjust(g_state, phi, bpar, direction = to_g_gs2)
×
255
    end if
256

257
    ! Now that we've advanced the explicit terms by dt and set the
258
    ! current g_state to this new result. Note that we do not update the
259
    ! fields, instead resetting them to their original values. The
260
    ! update to the fields is calculated as a part of the implicit
261
    ! field solve and including the NL contribution here leads to
262
    ! "double counting".
263

264
    call time_message(.false., time_add_explicit_terms, 'Explicit terms')
×
265
    call get_mp_times(total_time = mp_total_after)
×
266
    time_add_explicit_terms_mpi = time_add_explicit_terms_mpi + (mp_total_after - mp_total)
×
267
  end subroutine advance_nonlinear_term
×
268

269
  !> Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using
270
  !> AB scheme of orders up to 3rd.
271
  subroutine advance_nonlinear_term_ab3(g_state, phinew, aparnew, bparnew, dt, &
×
272
       fields_func, source_func)
273
    use theta_grid, only: ntgrid
274
    use gs2_time, only: get_adams_bashforth_coefficients
275
    use gs2_time, only: code_dt_min, dt_not_set, code_dt_prev1, code_dt_prev2
276
    use mp, only: mp_abort, max_allreduce
277
    use gs2_layouts, only: g_lo
278
    use array_utils, only: copy, zero_array
279
    implicit none
280
    complex, dimension (-ntgrid:, :, g_lo%llim_proc:), intent(in out) :: g_state
281
    complex, dimension (-ntgrid:,:,:), intent (in out) :: phinew, aparnew, bparnew
282
    !> Arrays to hold the nonlinear source history similar to gexp_1/2/3
283
    complex, dimension(:, :, :), allocatable :: gnl_1, gnl_2, gnl_3
×
284
    real, intent(in) :: dt
285
    real :: internal_time, time_step, max_vel, half_time
286
    real, dimension(:), allocatable :: ab_coeffs
×
287
    integer :: counter, iglo
288
    procedure(invert_field_func) :: fields_func
289
    procedure(source_term_func) :: source_func
290

291
    ! Initialise internal variables
292
    internal_time = 0.0
×
293
    counter = 1
×
294

295
    allocate(gnl_1(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc)) ; call zero_array(gnl_1)
×
296
    allocate(gnl_2(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc)) ; call zero_array(gnl_2)
×
297
    allocate(gnl_3(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc)) ; call zero_array(gnl_3)
×
298

299
    ! Reset time step sizes to reflect fact we're restarting the explicit scheme
300
    code_dt_prev1 = dt_not_set
×
301
    code_dt_prev2 = dt_not_set
×
302

303
    ! Keep taking explicit steps until we have advanced by dt
304
    do while (internal_time < dt)
×
305
       ! Get the current nonlinear source term and cfl limit
306
       call source_func(g_state, gnl_1, phinew, aparnew, bparnew, &
×
307
            max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
308
            calculate_cfl_limit = .true.)
×
309
       nrhs_calls = nrhs_calls + 1
×
310

311
       ! Get the maximum velocity across all processors
312
       ! and convert to time_step estimate.
313
       call max_allreduce(max_vel)
×
314
       time_step = 1.0 / max_vel
×
315
       time_step = time_step * time_step_safety_factor
×
316

317
       ! If the cfl time step is too small then abort
318
       if (time_step < code_dt_min) call mp_abort("Time step has become too small in advance_nonlinear_term.", .true.)
×
319

320
       ! Now make sure we're don't overshoot the target time by limiting time step.
321
       time_step = min(time_step, dt - internal_time)
×
322

323
       ! Get the Adams-Bashforth coefficients for the current collection of time step sizes
324
       ab_coeffs = get_adams_bashforth_coefficients(target_dt = time_step)
×
325

326
       half_time = 0.5 * time_step
×
327

328
       ! Advance solution by time_step
329
       select case(counter)
×
330
       case(1)
331
          !$OMP PARALLEL DO DEFAULT(none) &
332
          !$OMP PRIVATE(iglo) &
333
          !$OMP SHARED(g_lo, g_state, half_time, ab_coeffs, gnl_1) &
334
          !$OMP SCHEDULE(static)
335
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
336
             g_state(:, :, iglo) = g_state(:, :, iglo) + half_time *( &
×
337
                  ab_coeffs(1) * gnl_1(:, :, iglo)       &
×
338
                  )
×
339
          end do
340
          !$OMP END PARALLEL DO
341
       case(2)
342
          !$OMP PARALLEL DO DEFAULT(none) &
343
          !$OMP PRIVATE(iglo) &
344
          !$OMP SHARED(g_lo, g_state, half_time, ab_coeffs, gnl_1, gnl_2) &
345
          !$OMP SCHEDULE(static)
346
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
347
             g_state(:, :, iglo) = g_state(:, :, iglo) + half_time *( &
×
348
                  ab_coeffs(1) * gnl_1(:, :, iglo) +     &
×
349
                  ab_coeffs(2) * gnl_2(:, :, iglo)       &
×
350
                  )
×
351
          end do
352
          !$OMP END PARALLEL DO
353
       case default
354
          !$OMP PARALLEL DO DEFAULT(none) &
355
          !$OMP PRIVATE(iglo) &
356
          !$OMP SHARED(g_lo, g_state, half_time, ab_coeffs, gnl_1, gnl_2, gnl_3) &
357
          !$OMP SCHEDULE(static)
358
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
359
             g_state(:, :, iglo) = g_state(:, :, iglo) + half_time *( &
×
360
                  ab_coeffs(1) * gnl_1(:, :, iglo) +     &
×
361
                  ab_coeffs(2) * gnl_2(:, :, iglo) +     &
×
362
                  ab_coeffs(3) * gnl_3(:, :, iglo)       &
×
363
                  )
×
364
          end do
365
          !$OMP END PARALLEL DO
366
       end select
367

368
       ! Calculate the fields consistent with the new g
369
       call calculate_fields(fields_func, g_state, phinew, aparnew, bparnew)
×
370

371
       ! Increment internal state
372
       internal_time = internal_time + time_step
×
373
       counter = counter + 1
×
374
       nsuccessful_steps = nsuccessful_steps + 1
×
375

376
       ! Shift time history along one point
377
       code_dt_prev2 = code_dt_prev1
×
378
       code_dt_prev1 = time_step
×
379
       call copy(gnl_2, gnl_3)
×
380
       call copy(gnl_1, gnl_2)
×
381
    end do
382

383
  end subroutine advance_nonlinear_term_ab3
×
384

385
  !> Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using
386
  !> backwards Euler with fixed-point iteration.
387
  subroutine advance_nonlinear_term_beuler(g_state, phinew, aparnew, bparnew, dt, &
×
388
       fields_func, source_func)
389
    use theta_grid, only: ntgrid
390
    use gs2_time, only: code_dt_min
391
    use mp, only: mp_abort, max_allreduce, proc0
392
    use run_parameters, only: immediate_reset
393
    use gs2_layouts, only: g_lo
394
    use array_utils, only: copy, gs2_max_abs, zero_array
395
    use warning_helpers, only: exactly_equal
396
    implicit none
397
    complex, dimension (-ntgrid:, :, g_lo%llim_proc:), intent(in out) :: g_state
398
    complex, dimension (-ntgrid:,:,:), intent (in out) :: phinew, aparnew, bparnew
399
    complex, dimension(:, :, :), allocatable :: g_cur, g_local, gnl_1, gnl_2
×
400
    real, intent(in) :: dt
401
    real :: internal_time, time_step, max_vel
402
    integer :: counter, iglo
403
    integer, parameter :: iteration_limit = 20
404
    real, parameter :: time_adjust_factor = 2
405
    logical, parameter :: debug = .false.
406
    real :: inverse_error, error, actual_dt, half_time, cfl_time_step
407
    real, dimension(3) :: errors
408
    real, save :: time_step_last = -1
409
    procedure(invert_field_func) :: fields_func
410
    procedure(source_term_func) :: source_func
411
    logical :: reject_step, cfl_limited
412
    complex, dimension(:, :, :), allocatable :: error_array
×
413

414
    ! Initialise internal variables
415
    if (time_step_last < 0) time_step_last = dt
×
416
    internal_time = 0.0
×
417
    time_step = time_step_last
×
418

419
    allocate(error_array(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
420
    allocate(g_cur(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
421
    allocate(g_local(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
422
    allocate(gnl_1(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
423
    allocate(gnl_2(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
424
    call copy(g_state, g_cur)
×
425
    call zero_array(error_array) ; call zero_array(g_local)
×
426
    call zero_array(gnl_1); call zero_array(gnl_2)
×
427

428
    ! Keep on advancing until we have advanced by dt
429
    do while (internal_time < dt)
×
430
       ! Reset counter and error
431
       counter = 0
×
432
       error = 1.0
×
433
       errors = 0.0
×
434

435
       ! Limit the actual time step that we take to avoid overshooting (i.e.
436
       ! to keep internal_time <= dt).
437
       actual_dt = min(time_step, dt - internal_time)
×
438
       half_time = 0.5 * actual_dt
×
439

440
       ! Calculate the current NL source term, for use in error
441
       ! estimates later.
442
       call source_func(g_state, gnl_2, phinew, aparnew, bparnew, &
×
443
            max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
444
            calculate_cfl_limit = .false.)
×
445
       nrhs_calls = nrhs_calls + 1
×
446

447
       ! Try to take a step of size actual_dt using fixed point iteration with
448
       ! backwards Euler.
449
       do while (error > convergence_tolerance .and. counter < iteration_limit)
×
450

451
          ! Get the current nonlinear source term and cfl limit -- note we don't use the
452
          ! cfl limit here
453
          call source_func(g_state, gnl_1, phinew, aparnew, bparnew, &
×
454
               max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
455
               calculate_cfl_limit = .false.)
×
456
          nrhs_calls = nrhs_calls + 1
×
457

458
          ! Calculate the estimate of the new g
459
          !$OMP PARALLEL DO DEFAULT(none) &
460
          !$OMP PRIVATE(iglo) &
461
          !$OMP SHARED(g_lo, g_local, g_cur, half_time, gnl_1) &
462
          !$OMP SCHEDULE(static)
463
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
464
             g_local(:, :, iglo) = g_cur(:, :, iglo) + half_time * gnl_1(:, :, iglo)
×
465
          end do
466
          !$OMP END PARALLEL DO
467

468
          ! Calculate the maximum absolute error and the maximum value of the previous solution
469

470
          ! Note we _could_ merge the following loop witth the gs2_max_abs kernel
471
          ! so that we only loop over the error_array once (and actually don't
472
          ! need to explicitly store the full array)
473
          !$OMP PARALLEL DO DEFAULT(none) &
474
          !$OMP PRIVATE(iglo) &
475
          !$OMP SHARED(g_lo, error_array, g_state, g_local) &
476
          !$OMP SCHEDULE(static)
477
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
478
             error_array(:, :, iglo) = g_state(:, :, iglo) - g_local(:, :, iglo)
×
479
          end do
480
          !$OMP END PARALLEL DO
481
          errors(1) = gs2_max_abs(error_array)
×
482
          errors(2) = gs2_max_abs(g_local)
×
483
          call max_allreduce(errors)
×
484

485
          ! Calculate an approximate relative error
486
          error = errors(1) / errors(2)
×
487

488
          ! Update internals -- increment counter and update the 'previous' solution with the new one.
489
          counter = counter + 1
×
490
          call copy(g_local, g_state)
×
491

492
          ! Update the fields consistently with the current solution
493
          call fields_func(g_state, phinew, aparnew, bparnew, local_only = .true.)
×
494
       end do
495

496
       ! Now we need to decide if the last step was successful or not. There are two exit conditions, either the
497
       ! error was small enough (good step) or the iteration limit was exceeded (bad step).
498
       if (counter >= iteration_limit) then
×
499
          ! If it was a bad step then throw away these changes and reduce the time step
500
          if (proc0 .and. debug) print*,'Too many iterations, reducing step from',time_step,'to',time_step/time_adjust_factor,'with internal_time = ',internal_time
501

502
          ! Reduce the time step
503
          time_step = time_step / time_adjust_factor
×
504
          ! If the time step is now too small then abort
505
          if (time_step < code_dt_min) call mp_abort("Time step has become too small in advance_nonlinear_term.", .true.)
×
506

507
          ! If we've taken too many iterations then we want to drop the time step and try
508
          ! again. To try again we must ensure that gnew has been reset to g, and the fields
509
          ! are consistent with this.
510
          reject_step = .true.
×
511
       else
512
          ! If we have met the convergence condition then we should estimate the error
513
          ! on the solution and see if we need to reject the step or how we should adjust
514
          ! the step size.
515
          ! To estimate the error we compare the backward Euler update with an approximate
516
          ! Crank-Nicolson update.
517
          ! g_new_beuler = g + 0.5 * actual_dt * gnl_1(g_new_beuler)
518
          ! g_new_cn = g + 0.5 * actual_dt * (gnl_1(g_new_beuler) + gnl_1(g))/2
519
          ! abs_error  = |g_new_beuler - g_new_cn|
520
          !      = 0.5 * actual_dt * 0.5 * (gnl_1(g_new_beuler) - gnl_1(g)
521
          ! Where gnl_1(h) is the source func evaluated with gnew=h.
522

523
          ! Here we calculate gnl_1(g_new_beuler) exactly. We could probably
524
          ! rely on using the existing value of gnl_1(g_local) as, whilst this is
525
          ! strictly from the previous iteration of the fixed-point scheme, we know
526
          ! it is within the convergence tolerance of the true value.
527
          call source_func(g_state, gnl_1, phinew, aparnew, bparnew, &
×
528
               max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
529
               calculate_cfl_limit = .true.)
×
530
          nrhs_calls = nrhs_calls + 1
×
531

532
          ! Note we _could_ merge the following loop witth the gs2_max_abs kernel
533
          ! so that we only loop over the error_array once (and actually don't
534
          ! need to explicitly store the full array)
535
          !$OMP PARALLEL DO DEFAULT(none) &
536
          !$OMP PRIVATE(iglo) &
537
          !$OMP SHARED(g_lo, error_array, gnl_1, gnl_2) &
538
          !$OMP SCHEDULE(static)
539
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
540
             error_array(:, :, iglo) = gnl_1(:, :, iglo) - gnl_2(:, :, iglo)
×
541
          end do
542
          !$OMP END PARALLEL DO
543
          errors(1) = gs2_max_abs(error_array) * (0.5 * half_time)
×
544
          errors(3) = max_vel
×
545
          call max_allreduce(errors)
×
546
          inverse_error = get_inverse_error(errors)
×
547
          cfl_time_step = 1.0 / errors(3)
×
548

549
          !Calculate new time step based on error
550
          time_step = min(actual_dt * sqrt(inverse_error), cfl_time_step)
×
551
          cfl_limited = time_step < actual_dt .and. exactly_equal(time_step, cfl_time_step)
×
552
          time_step = time_step * time_step_safety_factor
×
553

554
          ! If error is too large reject the step
555
          if ( inverse_error < 1.0 .or. (cfl_limited .and. immediate_reset)) then
×
556
             if (proc0 .and. debug) print*,'Attempted beuler step failed due to error or cfl, adapting time step'
557
             ! Reset as we're retrying the step.
558
             reject_step = .true.
×
559
          else
560
             reject_step = .false.
×
561
          end if
562
       end if
563

564
       if (reject_step) then
×
565
          ! If the time step is now too small then abort
566
          if (time_step < code_dt_min) call mp_abort("Time step has become too small in advance_nonlinear_term.", .true.)
×
567
          call copy(g_cur, g_state)
×
568
          ! If we store the initial field then we could perhaps just copy here instead.
569
          ! This could be faster as it avoids communications involved in the velocity integration.
570
          call calculate_fields(fields_func, g_state, phinew, aparnew, bparnew)
×
571
          nfailed_steps = nfailed_steps + 1
×
572
       else
573
          ! If we have met the error condition and not taken too many steps
574
          ! then g_state represents our new solution so make sure we update g to reflect this
575
          internal_time = internal_time + actual_dt
×
576
          call copy(g_state, g_cur)
×
577
          nsuccessful_steps = nsuccessful_steps + 1
×
578
       end if
579

580
    end do
581

582
    ! Store the final time step size so we can start from this in the next callg
583
    time_step_last = time_step
×
584

585
    if(proc0.and.debug) print*,'Explicit done in ',counter,'steps with final error',inverse_error,'and step',time_step,'final time',internal_time
586

587
  end subroutine advance_nonlinear_term_beuler
×
588

589
  !> Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using
590
  !> an Runge-Kutta scheme with embedded error estimate for error control. This method
591
  !> is a wrapper to the true generic RK implementation.
592
  subroutine advance_nonlinear_term_rk(g_state, phinew, aparnew, bparnew, dt, &
×
593
       fields_func, source_func)
594
    implicit none
595
    complex, dimension (:,:,:), intent(in out) :: g_state
596
    complex, dimension (:,:,:), intent (in out) :: phinew, aparnew, bparnew
597
    real, intent(in) :: dt
598
    procedure(invert_field_func) :: fields_func
599
    procedure(source_term_func) :: source_func
600

601
    call advance_nonlinear_term_rk_implementation( &
602
         g_state, phinew, aparnew, bparnew, dt, &
×
603
         fields_func, source_func, the_rk_scheme)
×
604
  end subroutine advance_nonlinear_term_rk
×
605

606
  !> Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using
607
  !> an Runge-Kutta scheme with embedded error estimate for error control
608
  subroutine advance_nonlinear_term_rk_implementation(g_state, phinew, &
×
609
       aparnew, bparnew, dt, fields_func, source_func, scheme)
×
610
    use theta_grid, only: ntgrid
611
    use gs2_time, only: code_dt_min
612
    use mp, only: mp_abort, proc0, max_allreduce
613
    use gs2_layouts, only: g_lo
614
    use rk_schemes, only: rk_scheme_type
615
    use run_parameters, only: immediate_reset
616
    use array_utils, only: copy, gs2_max_abs, zero_array
617
    use warning_helpers, only: exactly_equal
618
    implicit none
619
    complex, dimension (-ntgrid:, :, g_lo%llim_proc:), intent(in out) :: g_state
620
    complex, dimension (-ntgrid:,:,:), intent (in out) :: phinew, aparnew, bparnew
621
    complex, dimension(:, :, :), allocatable :: g_cur, gnl_1
×
622
    real, intent(in) :: dt
623
    type(rk_scheme_type), intent(in) :: scheme
624
    complex, dimension (:, :, :, :), allocatable :: stages
×
625
    real :: internal_time, time_step, max_vel, inverse_error
626
    real :: new_time_step, cfl_time_step
627
    integer :: counter, nstages
628
    logical, parameter :: debug = .false.
629
    integer :: istage, istage_sub, iglo
630
    real, dimension(:), allocatable :: solution_coeffs, error_coeffs
×
631
    real, dimension(3) :: errors
632
    real, save :: time_step_last = -1
633
    procedure(invert_field_func) :: fields_func
634
    procedure(source_term_func) :: source_func
635
    logical :: cfl_limited
636
    real :: half_time
637
    real, dimension(:), allocatable :: stage_coeffs
×
638

639
    ! Initialise internal variables
640
    if (time_step_last < 0) time_step_last = dt
×
641
    internal_time = 0.0
×
642
    counter = 1
×
643
    time_step = time_step_last
×
644

645
    ! Copy the relevant coefficients for forming the solution to reduce
646
    ! later duplication.
647
    if (scheme%follow_high_order) then
×
648
       allocate(solution_coeffs, source = scheme%high_order_coeffs)
×
649
    else
650
       allocate(solution_coeffs, source = scheme%lower_order_coeffs)
×
651
    end if
652

653
    ! Get the error coefficients
654
    allocate(error_coeffs, source = scheme%high_order_coeffs - scheme%lower_order_coeffs)
×
655

656
    nstages = scheme%number_of_stages
×
657
    allocate(stages(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc, nstages))
×
658
    allocate(gnl_1(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
659
    allocate(g_cur(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
660
    call copy(g_state, g_cur) ; call zero_array(stages) ; call zero_array(gnl_1)
×
661

662
    ! Keep taking explicit steps until we have advanced by dt
663
    do while (internal_time < dt)
×
664
       if(proc0.and.debug) print*,'Taking step with size',time_step,'at time',internal_time
665
       half_time = 0.5 * time_step
×
666

667
       do istage = 1, nstages
×
668
          call copy(g_cur, g_state)
×
669
          stage_coeffs = half_time * scheme%coeffs(:, istage)
×
670

671
          !$OMP PARALLEL DO DEFAULT(none) &
672
          !$OMP PRIVATE(istage_sub, iglo) &
673
          !$OMP SHARED(g_lo, g_state, half_time, scheme, stages, istage, stage_coeffs) &
674
          !$OMP COLLAPSE(2) &
675
          !$OMP SCHEDULE(static)
676
          do istage_sub = 1, istage-1
×
677
             do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
678
                g_state(:, :, iglo) = g_state(:, :, iglo) + stage_coeffs(istage_sub) * stages(:, :, iglo, istage_sub)
×
679
             end do
680
          end do
681
          !$OMP END PARALLEL DO
682
          call calculate_fields(fields_func, g_state, phinew, aparnew, bparnew)
×
683

684
          ! Get the current nonlinear source term
685
          call source_func(g_state,  stages(:, :, :, istage), &
×
686
               phinew, aparnew, bparnew, &
×
687
               max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
688
               calculate_cfl_limit = (istage == 1))
×
689

690
          nrhs_calls = nrhs_calls + 1
×
691
          if (istage == 1) errors(3) = max_vel
×
692

693
       end do
694

695
       ! One can directly form the error by simply forming error =
696
       ! time_step * 0.5 * sum(stages_i * (high_order_coeffs_i -
697
       ! low_order_coeffs_i)). Here we store this in gnl_1
698
       ! without the time_step * 0.5 factor (to avoid extra
699
       ! array operations). We multiply the final error by
700
       ! time_step/2 to ensure full consistency.
701
       ! Could we merge this OpenMP parallel region with the
702
       ! other loop directly below to avoid some overhead?
703
       !$OMP PARALLEL DO DEFAULT(none) &
704
       !$OMP PRIVATE(iglo) &
705
       !$OMP SHARED(g_lo, gnl_1, error_coeffs, stages) &
706
       !$OMP SCHEDULE(static)
707
       do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
708
          gnl_1(:, :, iglo) = error_coeffs(1) * stages(:, :, iglo, 1)
×
709
       end do
710
       !$OMP END PARALLEL DO
711

712
       !$OMP PARALLEL DO DEFAULT(none) &
713
       !$OMP PRIVATE(iglo, istage) &
714
       !$OMP SHARED(g_lo, gnl_1, error_coeffs, nstages, stages) &
715
       !$OMP COLLAPSE(2) &
716
       !$OMP SCHEDULE(static)
717
       do istage = 2, nstages
×
718
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
719
             gnl_1(:, :, iglo) = gnl_1(:, :, iglo) + error_coeffs(istage) * stages(:, :, iglo, istage)
×
720
          end do
721
       end do
722
       !$OMP END PARALLEL DO
723

724
       ! Store the absolute error
725
       errors(1) = gs2_max_abs(gnl_1) * half_time
×
726
       ! Estimate the absolution size of the solution using
727
       ! the old solution
728
       errors(2) = gs2_max_abs(g_cur)
×
729
       call max_allreduce(errors)
×
730
       inverse_error = get_inverse_error(errors)
×
731
       cfl_time_step = 1.0 / errors(3)
×
732

733
       !Calculate the new_time_step
734
       new_time_step = min(time_step * (inverse_error)**(1.0/scheme%order), cfl_time_step)
×
735
       cfl_limited = new_time_step < time_step .and. exactly_equal(new_time_step, cfl_time_step)
×
736
       new_time_step = new_time_step * time_step_safety_factor
×
737

738
       if(proc0 .and. debug) print*,'Inverse error is ',inverse_error,'at internal_time', &
739
            internal_time,'with step',time_step,'new time step is',new_time_step, &
740
            'counter=',counter
741

742
       ! If the cfl time step is too small then abort
743
       if (new_time_step < code_dt_min) &
×
744
            call mp_abort("Time step has become too small in advance_nonlinear_term.", .true.)
×
745

746
       ! If the error is too large we need to retry with a smaller step, so don't
747
       ! update anything
748
       if (inverse_error < 1.0 .or. (cfl_limited .and. immediate_reset)) then
×
749
          if(proc0.and.debug) print*,'Time step with size',time_step, &
750
               'failed at internal_time',internal_time,'and counter',counter, &
751
               'retrying with step',new_time_step
752
          nfailed_steps = nfailed_steps + 1
×
753
       else
754
          stage_coeffs = solution_coeffs * half_time
×
755
          !Update the solution. Note we could probably merge the two
756
          !OpenMP loops into a single OpenMP parallel region and then
757
          !just add OMP DO to each loop. This would save a bit of overhead.
758
          !$OMP PARALLEL DO DEFAULT(none) &
759
          !$OMP PRIVATE(iglo) &
760
          !$OMP SHARED(g_lo, g_cur, solution_coeffs, stage_coeffs, stages) &
761
          !$OMP SCHEDULE(static)
762
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
763
             g_cur(:, :, iglo) = g_cur(:, :, iglo) + stage_coeffs(1) * stages(:, :, iglo, 1)
×
764
          end do
765
          !$OMP END PARALLEL DO
766

767
          !$OMP PARALLEL DO DEFAULT(none) &
768
          !$OMP PRIVATE(istage, iglo) &
769
          !$OMP SHARED(g_lo, g_cur, solution_coeffs, nstages, half_time, &
770
          !$OMP stages, stage_coeffs) &
771
          !$OMP COLLAPSE(2) &
772
          !$OMP SCHEDULE(static)
773
          do istage = 2, nstages
×
774
             do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
775
                g_cur(:, :, iglo) = g_cur(:, :, iglo) + stage_coeffs(istage) * stages(:, :, iglo, istage)
×
776
             end do
777
          end do
778
          !$OMP END PARALLEL DO
779

780
          internal_time = internal_time + time_step
×
781
          nsuccessful_steps = nsuccessful_steps + 1
×
782
       end if
783

784
       !Now update the time step, capped to dt
785
       time_step = min(new_time_step, dt)
×
786

787
       ! Make sure we don't overshoot the target time by limiting time step.
788
       if(internal_time < dt) then
×
789
          time_step = min(time_step, dt - internal_time)
×
790
          counter = counter + 1
×
791
       end if
792

793
    end do
794

795
    if(proc0.and.debug) print*,'Completed rk step in',counter,'steps to',internal_time
796

797
    ! We may be able to avoid this copy and field solve on exit if we restructure the main
798
    ! loop slightly so g_state represents the new solution. This may require use of g_work
799
    ! or similar to hold intermediate values.
800
    call copy(g_cur, g_state)
×
801

802
    time_step_last = time_step
×
803
  end subroutine advance_nonlinear_term_rk_implementation
×
804

805
  !> Advances dg/dt = NL(g,chi) from current g_state, chi from t -> t+dt by using
806
  !> a simple Picard iteration scheme
807
  subroutine advance_nonlinear_term_picard(g_state, phinew, &
×
808
       aparnew, bparnew, dt, fields_func, source_func)
×
809
    use theta_grid, only: ntgrid
810
    use gs2_time, only: code_dt_min
811
    use mp, only: mp_abort, proc0, max_allreduce
812
    use gs2_layouts, only: g_lo
813
    use array_utils, only: copy, gs2_max_abs, zero_array
814
    implicit none
815
    complex, dimension (-ntgrid:, :, g_lo%llim_proc:), intent(in out) :: g_state
816
    complex, dimension (-ntgrid:,:,:), intent (in out) :: phinew, aparnew, bparnew
817
    complex, dimension(:, :, :), allocatable :: gnl_1
×
818
    real, intent(in) :: dt
819
    real :: internal_time, time_step, max_vel, cfl_limit
820
    integer :: counter, ipic, iglo
821
    logical, parameter :: debug = .false.
822
    real, dimension(4) :: errors
823
    real, save :: time_step_last = -1
824
    procedure(invert_field_func) :: fields_func
825
    procedure(source_term_func) :: source_func
826
    complex, dimension(:, :, :), allocatable :: g_start
×
827
    integer, parameter :: npic_limit = 15, npic_low_limit = 3
828
    real, parameter :: time_step_factor = 2.0
829
    logical, parameter :: check_cfl = .true.
830
    logical :: converged
831
    real :: half_time
832

833
    ! Initialise internal variables
834
    if (time_step_last < 0) time_step_last = dt
×
835
    internal_time = 0.0
×
836
    counter = 1
×
837
    time_step = time_step_last
×
838

839
    allocate(g_start(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
840
    allocate(gnl_1(-ntgrid:ntgrid, 2, g_lo%llim_proc:g_lo%ulim_alloc))
×
841
    call zero_array(g_start) ; call zero_array(gnl_1)
×
842
    ! Keep taking explicit steps until we have advanced by dt
843
    do while (internal_time < dt)
×
844
       call copy(g_state, g_start)
×
845
       half_time = 0.5 * time_step
×
846

847
       if(proc0.and.debug) print*,'Taking step with size',time_step,'at time',internal_time
848
       ipic = 1
×
849

850
       call calculate_fields(fields_func, g_state, phinew, aparnew, bparnew)
×
851

852
       ! Get the current nonlinear source term
853
       call source_func(g_state, gnl_1, &
×
854
            phinew, aparnew, bparnew, &
×
855
            max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
856
            calculate_cfl_limit = check_cfl)
×
857
       nrhs_calls = nrhs_calls + 1
×
858

859
       errors(1) = 0.0
×
860
       errors(2) = gs2_max_abs(gnl_1) * half_time
×
861
       errors(3) = gs2_max_abs(g_start)
×
862
       errors(4) = max_vel
×
863
       call max_allreduce(errors)
×
864

865
       ! Check for CFL limit - we can simply rescale the change in g and carry on
866
       ! so here we simply rescale the effective time step,
867
       cfl_limit = 1.0 / errors(4)
×
868
       if (time_step > cfl_limit .and. check_cfl) then
×
869
          ! Drop the effective time step a little from the limit
870
          cfl_limit = 0.9 * cfl_limit
×
871
          if(proc0.and.debug) print*,'Dropping time step due to cfl from',time_step,'to',cfl_limit
872
          errors(2) = errors(2) * cfl_limit / time_step
×
873
          time_step = cfl_limit
×
874
          half_time = time_step * 0.5
×
875
       end if
876

877
       converged = picard_converged(errors)
×
878

879
       !$OMP PARALLEL DO DEFAULT(none) &
880
       !$OMP PRIVATE(iglo) &
881
       !$OMP SHARED(g_state, g_start, gnl_1, half_time, g_lo)&
882
       !$OMP SCHEDULE(static)
883
       do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
884
          g_state(:, :, iglo) = g_start(:, :, iglo) + gnl_1(:, :, iglo) * half_time
×
885
       end do
886
       !$OMP END PARALLEL DO
887

888
       do while(.not. converged)
×
889
          ! Note the first step (i.e. call before this loop) is included in
890
          ! ipic here (i.e. it counts towards our approach towards npic_limit).
891
          ipic = ipic + 1
×
892
          if (ipic > npic_limit) exit
×
893

894
          call calculate_fields(fields_func, g_state, phinew, aparnew, bparnew)
×
895

896
          ! Get the current nonlinear source term
897
          call source_func(g_state, gnl_1, &
×
898
               phinew, aparnew, bparnew, &
×
899
               max_vel, need_to_adjust = .not. advance_nonadiabatic_dfn, &
900
               calculate_cfl_limit = check_cfl)
×
901
          nrhs_calls = nrhs_calls + 1
×
902

903
          errors(1) = errors(2)
×
904
          errors(2) = gs2_max_abs(gnl_1) * half_time
×
905
          errors(4) = max_vel
×
906
          call max_allreduce(errors)
×
907
          cfl_limit = 1.0 / errors(4)
×
908

909
          converged = picard_converged(errors)
×
910

911
          if(proc0.and.debug) print*,'ipic',ipic,converged,'|',errors
912

913
          if (time_step > cfl_limit .and. .not. converged) exit
×
914

915
          !$OMP PARALLEL DO DEFAULT(none) &
916
          !$OMP PRIVATE(iglo) &
917
          !$OMP SHARED(g_state, g_start, gnl_1, half_time, g_lo)&
918
          !$OMP SCHEDULE(static)
919
          do iglo = g_lo%llim_proc, g_lo%ulim_proc
×
920
             g_state(:, :, iglo) = g_start(:, :, iglo) + gnl_1(:, :, iglo) * half_time
×
921
          end do
922
          !$OMP END PARALLEL DO
923
       end do
924

925
       if (.not. converged) then
×
926
          ! If we didn't converge in the limit then reset state
927
          ! to the start of the step and reduce the time step.
928
          call copy(g_start, g_state)
×
929
          time_step = time_step / time_step_factor
×
930
          nfailed_steps = nfailed_steps + 1
×
931

932
          ! If the time step is too small then abort
933
          if (time_step < code_dt_min) &
×
934
               call mp_abort("Time step has become too small in advance_nonlinear_term.", .true.)
×
935
       else
936
          internal_time = internal_time + time_step
×
937
          nsuccessful_steps = nsuccessful_steps + 1
×
938
          ! If we converged quickly then increase the time step
939
          if (ipic <= npic_low_limit) time_step = time_step * time_step_factor
×
940
       end if
941

942
       !Now update the time step, capped to dt
943
       time_step = min(time_step, dt)
×
944

945
       ! Make sure we don't overshoot the target time by limiting time step.
946
       if(internal_time < dt) then
×
947
          time_step = min(time_step, dt - internal_time)
×
948
          counter = counter + 1
×
949
       end if
950

951
    end do
952

953
    if(proc0.and.debug) print*,'Completed Picard step in',counter,'steps to',internal_time,'last time step',time_step
954

955
    time_step_last = time_step
×
956
  contains
957
    pure logical function picard_converged(errors) result(converged)
×
958
      implicit none
959
      real, dimension(4), intent(in) :: errors
960
      converged = &
961
           !Change in solution less than abs tol
962
           (errors(2) < absolute_tolerance) .or.  &
963
           ! Change in the change, less than abs tol
964
           (abs(errors(2) - errors(1)) < absolute_tolerance) .or. &
965
           ! Change is less than the relative tolerance
966
           (errors(2) < relative_tolerance * errors(3)) .or. &
967
           ! Change in the change less than the relative tolerance
968
           (abs(errors(2) - errors(1)) < relative_tolerance * errors(3))
×
969
    end function picard_converged
×
970

971
  end subroutine advance_nonlinear_term_picard
972

973
  !> Returns the normalised inverse error estimate. In other words
974
  !> the error tolerance, rtol*errors(2) + atol, divided by the error
975
  !> estimate (plus a small number to avoid divide by zero).
976
  !>
977
  !> Assumes that errors contains the global magnitude of the error
978
  !> estimate in the first element and the global maximum of the solution
979
  !> in the second element.
980
  real pure function get_inverse_error(errors) result(inverse_error)
×
981
    implicit none
982
    real, dimension(:), intent(in) :: errors
983
    inverse_error = (relative_tolerance * errors(2) + absolute_tolerance) / &
×
984
         (errors(1) + epsilon(errors(1)))
×
985
  end function get_inverse_error
×
986

987
  !> Calculates the fields consistent with the passed distribution function
988
  subroutine calculate_fields(fields_func, g_in, phi_out, apar_out, bpar_out)
×
989
    use nonlinear_terms, only: time_add_explicit_terms_field
990
    use job_manage, only: time_message
991
    procedure(invert_field_func) :: fields_func
992
    complex, dimension(:, :, :), intent(in) :: g_in
993
    complex, dimension (:,:,:), intent (out) :: phi_out, apar_out, bpar_out
994
    call time_message(.false., time_add_explicit_terms_field, 'Explicit terms - field')
×
995
    call fields_func(g_in, phi_out, apar_out, bpar_out, local_only = .true.)
×
996
    call time_message(.false., time_add_explicit_terms_field, 'Explicit terms - field')
×
997
  end subroutine calculate_fields
×
998

999
#include "split_nonlinear_terms_auto_gen.inc"
1000
end module split_nonlinear_terms
×
STATUS · Troubleshooting · Open an Issue · Sales · Support · CAREERS · ENTERPRISE · START FREE · SCHEDULE DEMO
ANNOUNCEMENTS · TWITTER · TOS & SLA · Supported CI Services · What's a CI service? · Automated Testing

© 2026 Coveralls, Inc