timestep_module.f90
Go to the documentation of this file.
1 !> @file timestep.f90
2 !!
3 !! The error file code for this file is ***W43***.
4 !!
5 !! Contains the module \ref timestep_module
6 !!
7 !! @brief Timestepping subroutines: \ref timestep, \ref exp_timestep etc.
8 
9 
10 !>
11 !! @brief The timestep_module contains all timestep subroutines and newton-raphsons
12 !! to solve the ODE system.
13 !!
14 !! This module contains several different treatments for the timestep and for
15 !! solving the ODE with different solvers (so far implicit euler and gear).
16 !! The gear solver develops its own timestep in \ref gear_module. It is only
17 !! affected by the subroutine \ref restrict_timestep.
18 !!
19 #include "macros.h"
23  implicit none
24 
25  !
26  ! Public and private fields and methods of the module.
27  !
28  public:: &
30  private:: &
33 
34 contains
35 
36 !>
37 !! @brief Network timestepping along the thermodynamic trajectory,
38 !! specified by {ztime, ztemp, zdens}[zsteps]
39 !!
40 !! This subroutine calls more specific subroutines to calculate the
41 !! timestep. This also depends on the conditions. For example, in NSE it
42 !! call \ref ye_timestep, while it calls \ref abchange_timestep otherwise.
43 !!
44 !! @note Gear solver will only be affected by the last call of this subroutine
45 !! of \ref restrict_timestep
46 !!
47 !! \b Edited:
48 !! - OK: 13.01.14
49 !! - OK: 7.11.16
50 !! - OK: 27.08.17 - merged in the rest of timestepping subroutines
51 !! - MR: 27.06.18 - introduced additional subroutine "restrict_timestep"
52 !! .
53 subroutine timestep (ctime,temp,dens,entropy,Y,dYdt,Ye,delYe,h,evmode)
54  use global_class, only: net_size
56  use gear_module, only: get_timestep
57  implicit none
58 
59  real(r_kind), intent(in) :: ctime !< current time
60  real(r_kind), intent(in) :: temp !< temperature [GK]
61  real(r_kind), intent(in) :: dens !< density [g/cm^3]
62  real(r_kind), intent(in) :: entropy !< [kB/baryon]
63  real(r_kind),dimension(net_size),intent(in):: Y !< abundances
64  real(r_kind),dimension(net_size),intent(in):: dYdt !< change of abundances
65  real(r_kind), intent(in) :: Ye !< electron fraction
66  real(r_kind), intent(in) :: delYe !< change of elec.fraction (Ye-Ye_p)
67  real(r_kind), intent(inout) :: h !< resulting timestep
68  integer, intent(in) :: evmode!< zone evolution mode
69  !
70  real(r_kind) :: h_old
71  real(r_kind) :: rad
72  real(r_kind) :: tdens =1e24 ! interpolated density: initialize to huge value
73  real(r_kind) :: ttemp =1e24 ! interpolated temperature: the same
74 
75  info_entry("timestep")
76 
77  !-- save old timestep
78  h_old = h
79 
80  ! MR: even without the if statement, the gear solver does not care
81  ! about the timestep, so I implemented it here
82  if (solver .ne. 1) then
83  !-- attempt to increase the timestep
84  h = h*timestep_max
85 
86  if (evmode.eq.em_nse) then
87  !-- timescale of Ye changes
88  h= min(h, ye_timestep(h, ye, delye))
89  if (h .ne. h) call raise_exception("Timestep got NaN after ye_timestep.",&
90  "timestep",430003)
91 
92  else
93  !-- timescale of abundance change (h_abchange)
94  h= min(h, abchange_timestep(h,y,dydt))
95  if (h .ne. h) call raise_exception("Timestep got NaN after abchange_timestep.",&
96  "timestep",430004)
97  endif
98 
99  !-- hydrodynamic timescale (h_hydro)
100  !-- Only do this in case temperature is not at lower boundary
101  if (temp .gt. freeze_rate_temp*1.000001) then
102  h= min(h, hydro_timestep(h,ctime,temp,dens,entropy,ttemp,tdens,ye))
103  if (h .ne. h) call raise_exception("Timestep got NaN after hydro_timestep.",&
104  "timestep",430005)
105  end if
106  elseif (solver .eq. 1) then
107  ! Also gear should be regulated by the restriction
108  h = get_timestep()
109  ! Get hydro quantities
110  call update_hydro(ctime+h, ttemp, tdens, rad, entropy, ye)
111  end if
112 
113 
114  ! Restrict the timestep to match specific values
115  call restrict_timestep(ctime,h,temp,dens,ttemp,tdens)
116  if (h .ne. h) call raise_exception("Timestep got NaN after restrict_timestep.",&
117  "timestep",430006)
118 
119  info_exit("timestep")
120 
121 end subroutine timestep
122 
123 
124 
125 !>
126 !! @brief Restricting the previously selected timestep
127 !!
128 !! Possible restrictions are given by the termination criterion and
129 !! the trajectory datapoints as well as custom selected snapshots.
130 !! Furthermore, it is restricted to hit the NSE transition temperature
131 !! (\ref parameter_class::nsetemp_cold).
132 !!
133 !! @note The timestep of the gear solver, which is evolved independently in
134 !! \ref gear_module is also effected by this subroutine.
135 !!
136 !! \b Edited:
137 !! - MR: 27.06.18
138 !! - MR: 22.12.20 - Removed "const" parameter
139 !! .
140 subroutine restrict_timestep(ctime,h,temp,dens,ttemp,tdens)
141  use parameter_class, only: &
150  use global_class, only: heating_switch
152  use parser_module
153  implicit none
154  real(r_kind), intent(in) :: ctime !< current time
155  real(r_kind), intent(in) :: temp !< temperature [GK]
156  real(r_kind), intent(in) :: dens !< density [g/cm^3]
157  real(r_kind), intent(in) :: ttemp !< temperature [GK] of the next step
158  real(r_kind), intent(in) :: tdens !< density [g/cm^3] of the next step
159  real(r_kind), intent(inout) :: h !< resulting timestep
160  real(r_kind) :: h_in !< ingoing timestep
161  real(r_kind) :: tnext !< Next time of the trajectory
162  real(r_kind) :: temp_tmp !< Temperature on hypothetical timestep
163  real(r_kind) :: dens_tmp !< Density on hypothetical timestep
164  integer :: i !< Loop variable
165  logical :: converged !< Convergence indicator for analytic mode
166 
167  ! Store the ingoing timestep for later comparison
168  h_in = h
169 
170  !-- cut by termination criterion
171  select case (termination_criterion)
172  case (0) ! by the end of trajectory file (regardless of timestep_traj_limit)
173  if (ctime+h.gt.ztime(zsteps)) h = ztime(zsteps) - ctime
174  case (1) ! by final_time
175  if (ctime+h.gt.final_time) h = final_time - ctime
176  case (2) ! by final_temp
177  if (ttemp .lt. final_temp*0.999d0) then
178  if ((trim(adjustl(trajectory_mode)) .eq. 'from_file')) then
179  ! Try to interpolate if still within the trajectory
180  if ((ctime.lt.ztime(zsteps)) .and. (.not. (heating_switch))) then
181  ! Interpolate the NSE temperature
182  tnext = ctime+h
183  call inverse_interp1d(zsteps,ztemp,final_temp*0.999d0,ztime,tnext)
184  ! Check that it worked, sometimes for very flat trajectories it has
185  ! rounding errors
186  call interp1d (zsteps,ztime,tnext,ztemp,temp_tmp)
187  if (temp_tmp .lt. final_temp) h = max(tnext-ctime,1e-15) ! it worked
188  if (temp_tmp .ge. final_temp) h = min(1.1*max(tnext-ctime,1e-15),h) ! make the timestep larger
189  else
190  ! otherwise, basic interpolation
191  h = h * (temp-final_temp)/(temp-ttemp)
192  end if
193  elseif ((trim(adjustl(trajectory_mode)) .eq. 'analytic')) then
194  ! Find the time of T=final_temp
195  call find_start_value(t9_analytic,final_temp*0.999d0,ctime+h,converged,tnext)
196  if (converged .and. (tnext .gt. ctime)) h = max(tnext-ctime,1e-15)
197  end if
198  end if
199  case (3) ! by final_dens
200  if (tdens.le.final_dens*0.999d0) then
201  if ((trim(adjustl(trajectory_mode)) .eq. 'from_file')) then
202  ! Try to interpolate if still within the trajectory
203  if ((ctime.lt.ztime(zsteps))) then
204  ! Interpolate the NSE temperature
205  tnext = ctime+h
206  call inverse_interp1d(zsteps,zdens,final_dens*0.999d0,ztime,tnext)
207  ! Check that it worked, sometimes for very flat trajectories it has
208  ! rounding errors
209  call interp1d (zsteps,ztime,tnext,zdens,temp_tmp)
210  if (temp_tmp .lt. final_dens) h = max(tnext-ctime,1e-15) ! it worked
211  if (temp_tmp .ge. final_dens) h = min(1.1*max(tnext-ctime,1e-15),h) ! make the timestep larger
212  else
213  ! otherwise, basic interpolation
214  h = h * (dens-final_dens)/(dens-tdens)
215  end if
216  elseif ((trim(adjustl(trajectory_mode)) .eq. 'analytic')) then
217  ! Find the time of rho=final_dens
218  call find_start_value(rho_analytic,final_dens*0.999d0,ctime+h,converged,tnext)
219  if (converged .and. (tnext .gt. ctime)) h = max(tnext-ctime,1e-15)
220  end if
221  end if
222  endselect
223 
224 
225  !-- Prevent NSE overshooting
226  if ((evolution_mode .eq. em_nse) .and. (ttemp .lt. nsetemp_cold)) then
227  if (trajectory_mode .eq. "from_file") then
228  ! No heating
229  if ((.not. (heating_switch)) .and. (.not. (ctime+h .gt.ztime(zsteps)))) then
230  if (ttemp .lt. nsetemp_cold*(1d0-1d-9)) then
231  ! Interpolate the NSE temperature
232  tnext = ctime+h
233  call inverse_interp1d(zsteps,ztemp,nsetemp_cold*(1d0-1d-9),ztime,tnext)
234  ! Check that it worked, sometimes for very flat trajectories it has
235  ! rounding errors
236  call interp1d (zsteps,ztime,tnext,ztemp,temp_tmp)
237  if (temp_tmp .lt. nsetemp_cold) h = min(max(tnext-ctime,1e-15),h) ! it worked
238  if (temp_tmp .ge. nsetemp_cold) h = min(1.1*max(tnext-ctime,1e-15),h) ! make the timestep larger
239  end if
240  end if
241  else if (trajectory_mode .eq. "analytic") then
242  ! No heating
243  if (.not. (heating_switch)) then
244  if (ttemp .lt. nsetemp_cold*(1d0-1d-9)) then
245  ! Find the time of T=T_nsecold
246  call find_start_value(t9_analytic,nsetemp_cold*(1d0-1d-9),ctime+h,converged,tnext)
247  if (converged .and. (tnext .gt. ctime)) h = max(tnext-ctime,1e-15)
248  end if
249  end if
250  end if
251  end if
252 
253 
254  !-- Prevent density overshooting for heating switch
255  if ((heating_mode .gt. 0) .and. (.not. heating_switch) .and. (tdens .lt. heating_density)) then
256  if (trajectory_mode .eq. "from_file") then
257  if (tdens .lt. heating_density*(1d0-1d-8)) then
258  ! Interpolate the NSE temperature
259  tnext = ctime+h
260  call inverse_interp1d(zsteps,zdens,heating_density*(1d0-1d-8),ztime,tnext)
261  ! Check that it worked, sometimes for very flat trajectories it has
262  ! rounding errors
263  call interp1d (zsteps,ztime,tnext,zdens,dens_tmp)
264  if (dens_tmp .lt. heating_density) h = min(max(tnext-ctime,1e-15),h) ! it worked
265  if (dens_tmp .ge. heating_density) h = min(1.1*max(tnext-ctime,1e-15),h) ! make the timestep larger
266  end if
267  else if (trajectory_mode .eq. "analytic") then
268  if (tdens .lt. heating_density*(1d0-1d-8)) then
269  ! Find the time of T=T_nsecold
270  call find_start_value(rho_analytic,heating_density*(1d0-1d-8),ctime+h,converged,tnext)
271  if (converged .and. (tnext .gt. ctime)) h = max(tnext-ctime,1e-15)
272  end if
273  end if
274  end if
275 
276  !-- if timestep is limited by trajectory steps, cut to the next grid point
277  if (timestep_traj_limit .and. &
278  & (trim(adjustl(trajectory_mode)) .eq. 'from_file')) then
279  if (ctime.lt.ztime(zsteps)) then
280  find_next_time: do i=1,zsteps-1
281  if (ztime(i).gt.ctime) exit find_next_time
282  enddo find_next_time
283  tnext= ztime(i)
284  if ((ctime.lt.tnext).and.(ctime+h.gt.tnext)) then
285  ! To ensure that the timestep is not to small
286  h = max(tnext-ctime,1e-15)
287  end if
288  endif
289  endif
290 
291  !-- Jump at least at the last trajectory point
292  if ((trim(adjustl(trajectory_mode)) .eq. 'from_file')) then
293  if ((ctime.lt.ztime(zsteps)) .and. (ctime+h .gt.ztime(zsteps))) then
294  ! To ensure that the timestep is not to small
295  h = max(ztime(zsteps)-ctime,1e-15)
296  endif
297  endif
298 
299  !-- Restrict timestep to match custom snapshots
300  if ((custom_snapshots) .or. (h_custom_snapshots)) then
301  ! Find the next time, should be cheap for short arrays
302  next_time: do i=1,snapshot_amount
303  if (snapshot_time(i) .gt. ctime) exit next_time
304  enddo next_time
305  ! Be careful with the array bounds
306  if (i .gt. snapshot_amount) i = snapshot_amount
307  ! Restrict the timestep if necessary
308  tnext=snapshot_time(i)
309  if ((ctime.lt.tnext).and.(ctime+h.gt.tnext)) then
310  ! To ensure that the timestep is not to small
311  h = max(tnext-ctime,1e-15)
312  end if
313  end if
314 
315  ! Also gear should be regulated by the restriction so change the timestep
316  if ((solver == 1) .and. (h_in .ne. h)) then
317  call set_timestep(h)
318  end if
319 end subroutine restrict_timestep
320 
321 
322 
323 !>
324 !! @brief Estimate the timestep from the change of electron fraction
325 !!
326 !! In NSE, only (slow) weak reactions are taken into account when calculating the
327 !! timestep. Therefore, only the change in \f$ y_e \f$ is considered here and
328 !! the timestep is usually larger than in the full network mode.
329 !! The timestep is calculated by:
330 !! \f[ h = h_\mathrm{old} \cdot \eta \frac{y_e}{|\Delta y_e|} \f]
331 !! with \f$ \eta \f$ being the parameter \ref parameter_class::timestep_factor.
332 !!
333 !! \b Created: OK 28.08.2017
334 !!
335 function ye_timestep(h,Ye,delYe)
337 implicit none
338 real(r_kind):: ye_timestep
339 real(r_kind),intent(in):: ye !< electron fraction
340 real(r_kind),intent(in):: delye !< change of the electron fraction
341 real(r_kind),intent(in):: h !< timescale on which this change happens
342 !
343  if (abs(delye) .gt. 1d-15) then
344  ye_timestep = h * timestep_factor * ye/abs(delye)
345  else
346  ye_timestep = h
347  end if
348 
349 end function ye_timestep
350 
351 
352 !>
353 !! @brief Estimate the timestep from the abundances change
354 !!
355 !! In the network evolution mode, the timestep is calculated based on the
356 !! most changing isotope. The timestep is the maximum of :
357 !! \f[ h = \eta \cdot \left|\frac{Y}{\Delta Y} \right| \f]
358 !! for all nuclei above a certain threshold abundance
359 !! \ref parameter_class::timestep_Ymin. The factor \f$ \eta \f$ is determined
360 !! by the parameter \ref parameter_class::timestep_factor.
361 !!
362 !! \b Created: OK 28.08.2017
363 !!
364 function abchange_timestep(h,Y,dYdt)
366 use global_class, only: net_size
367 implicit none
369 real(r_kind),dimension(net_size),intent(in):: y !< abundances
370 real(r_kind),dimension(net_size),intent(in):: dydt !< change of abundances
371 real(r_kind),intent(in):: h !< timescale on which this change happens
372 !
373 integer:: i
374 real(r_kind):: h1
375 
376  h1= h
377  do i = 1, net_size
378  if (y(i) .le. timestep_ymin) cycle
379  if (dydt(i) .eq. 0.d0) cycle
380  if (dabs(dydt(i)*h1) .gt. dabs(y(i))*timestep_factor) then
381  h1= timestep_factor * dabs(y(i)/dydt(i))
382  endif
383  end do
385 
386 end function abchange_timestep
387 
388 
389 !>
390 !! @brief Estimate the hydro timestep
391 !!
392 !! The timestep is restricted by the change of the hydrodynamic conditions
393 !! (i.e., temperature and density). The maximum change of them is determined
394 !! by the parameter \ref parameter_class::timestep_hydro_factor.
395 !! If the temperature or density change to rapidly the timestep will be decreased
396 !! and ultimately, if the timestep is lower than \f$ 10^{-24}s \f$ an exception
397 !! will be raised and the code terminates. The subroutine is also responsible
398 !! to update the next temperature (ttemp) and density (tdens) that are used
399 !! in \ref restrict_timestep.
400 !!
401 !! @note This subroutine restricts only the timestep for the implicit euler solver,
402 !! not for the gear solver.
403 !!
404 !! \b Created: OK 28.08.2017
405 !!
406 function hydro_timestep(h1,ctime,temp,dens,entropy,ttemp,tdens,Ye)
407 use parameter_class, only: heating_mode
409 use gear_module, only: get_timestep
410 use global_class, only: heating_switch
411 use expansion_module,only: rad
412 implicit none
413 real(r_kind):: hydro_timestep
414 real(r_kind),intent(in):: h1 !< timescale on which this change happens
415 real(r_kind),intent(in):: ctime !< current time
416 real(r_kind),intent(in):: temp !< temperature at ctime
417 real(r_kind),intent(in):: dens !< density at ctime
418 real(r_kind),intent(in):: entropy !< [kB/baryon]
419 real(r_kind),intent(inout):: ttemp !< computed temperature at ctime + h
420 real(r_kind),intent(inout):: tdens !< computed density at ctime + h
421 real(r_kind),intent(in):: ye !< electron fraction
422 !
423 real(r_kind) :: ttime,h
424 character(200) :: e_message2 !< Error message
425 
426 
427 h= h1
428 hydro_loop: do
429  if (h .ne. h) call raise_exception("Timestep got NaN.","hydro_timestep",430007)
430 
431  ! Suggested next time
432  ttime = ctime + h
433  ! Initialize old values
434  ttemp = temp
435  tdens = dens
436  ! Get temperature and density in the next step
437  call update_hydro(ttime, ttemp, tdens, rad, entropy, ye)
438 
439  ! check density (and temperature if no heating) increment
440  ! is small enough
441  if (heating_switch) then
442  if (dabs(tdens-dens)/dens.lt.timestep_hydro_factor) exit hydro_loop
443  h = h * timestep_hydro_factor * dens / abs(tdens-dens) &
444  * .999 ! black magic
445  else
446  if ((dabs(tdens-dens)/dens.lt.timestep_hydro_factor) .and. &
447  (dabs(ttemp-temp)/temp.lt.timestep_hydro_factor)) exit hydro_loop
448  h = h * timestep_hydro_factor &
449  / max(abs(tdens-dens)/dens,abs(ttemp-temp)/temp) &
450  * .999 ! black magic
451  endif
452 
453  ! if timestep too small, complain and exit
454  if (h.lt.1e-24) then
455  ! Prepare a meaningful message
456  ! Say something about the time where it occured
457  e_message2 = 'time: t ='//trim(adjustl(num_to_str(ctime)))//', h ='//trim(adjustl(num_to_str(h)))//new_line('A')
458  e_message2 = trim(adjustl(e_message2))//'dens: d(t)='//trim(adjustl(num_to_str(dens)))//', d(t+h)='//trim(adjustl(num_to_str(tdens)))//new_line('A')
459  e_message2 = trim(adjustl(e_message2))//'temp: d(t)='//trim(adjustl(num_to_str(temp)))//', d(t+h)='//trim(adjustl(num_to_str(ttemp)))//new_line('A')
460 
461  ! Raise the exception
462  call raise_exception(&
463  "Timestep too small (<1e-24). This happened because the hydrodynamic"//new_line('A')//&
464  'conditions were changing to fast. Check the hydrodynamic conditions'//new_line('A')//&
465  'or choose the "timestep_hydro_factor" parameter larger ('//trim(adjustl(num_to_str(timestep_hydro_factor)))//').'//&
466  new_line('A')//&
467  trim(adjustl(e_message2))&
468  ,'hydro_timestep',430008)
469  endif
470 enddo hydro_loop
471 
472 
474 
475 end function hydro_timestep
476 
477 
478 !>
479 !! Selects the timestep
480 !!
481 !! This subroutine is the main interface to the \ref driver.f90 .
482 !! It calls \ref timestep and initializes the timestep in case
483 !! of the first call.
484 !!
485 !! \b Edited:
486 !! - 17.12.2016
487 !! .
492 use nucstuff_class, only: el_ab
493 use jacobian_class, only: abchange
494 use gear_module, only: init_dydt
495 implicit none
496 logical, intent(inout):: init
497 
498  if (init) then
499 
500  ! calculate dYdt from reaction rates, needed for timestep
501  call el_ab (y, ye)
502  call abchange (time, t9, rhob, ye, rkm, y, dydt, evolution_mode)
503 
504  if(solver==1) call init_dydt(dydt)
505 
506  ! adjust the stepsize
508  endif
509 
510  call el_ab (y, ye)
511  call timestep (time,t9,rhob,ent,y,dydt,ye,ye-ye_p,stepsize,evolution_mode)
512  init = .false.
513 
514 end subroutine select_optimistic_timestep
515 
516 
517 
518 !>
519 !! @brief Returns temperature, density, and radius for a given time "time_i".
520 !!
521 !! The time, entropy, and electron fraction are only inputs and are not changed.
522 !! The update depends on the trajectory mode of the run and can be either
523 !! interpolated from a trajectory, got from the expansion, or from an analytic
524 !! expression.
525 !!
526 !! \b Edited:
527 !! - MR: 15.01.2020 - Moved it from advance_implicit_euler/gear to this new subroutine
528 !! - MR: 10.02.2023 - Gave an optional output for the temperature
529 !! .
530 subroutine update_hydro(time_i, temperature, density, radius, entropy, efraction, T_tr)
531  use single_zone_vars, only: y, t9_p, time_p
535  use inter_module, only: interp1d
537  use expansion_module, only: expansion,vel
539  use parser_module, only: parse_string
540  implicit none
541  real(r_kind),intent(in) :: time_i !< Time [s]
542  real(r_kind),intent(in) :: entropy !< Entropy [kB/nuc]
543  real(r_kind),intent(in) :: efraction !< Electron fraction
544  real(r_kind),intent(inout) :: temperature !< temperature at time_i
545  real(r_kind),intent(inout) :: density !< density at time_i
546  real(r_kind),intent(inout) :: radius !< radius at time_i
547  real(r_kind),intent(out),optional :: t_tr !< temperature at time_i, ignoring heating
548  real(r_kind) :: ysum, yasum
549  real(r_kind) :: htmp !< Temporary stepsize
550 
551 
552  ! Calculate necessary sums for expansion
553  ysum = sum(y(1:net_size))
554  yasum = sum(y(1:net_size)*isotope(1:net_size)%mass)
555  ! yzsum = sum(Y(1:net_size)*isotope(1:net_size)%p_nr)
556 
557  ! interpolate temperature, density and radius
558  if (trim(adjustl(trajectory_mode)) .eq. "analytic") then
559  if (.not. (heating_switch)) then
560  ! Update temperature from analytic expression
561  temperature = parse_string(t9_analytic,time_i)
562  else
563  temperature = t9_p
564  if (present(t_tr)) then
565  t_tr = parse_string(t9_analytic,time_i)
566  end if
567  end if
568  ! Update density and radius from analytic expression
569  density = parse_string(rho_analytic,time_i)
570  radius = parse_string(rkm_analytic,time_i)
571  ! Avoid underflows/overflows
572  if (temperature .le. freeze_rate_temp) temperature = freeze_rate_temp
573  if (density .le. 1d-99) density = 1d-99
574  if (radius .ge. 1d99) radius = 1d99
575 
576  else if (trim(adjustl(trajectory_mode)) .eq. "from_file") then
577  if (time_i.gt.ztime(zsteps)) then
578  htmp = time_i - time_p
579  call expansion(time_p,htmp,efraction,density,radius,temperature,entropy,vel,yasum/ysum)
580  if (present(t_tr)) then
581  t_tr = temperature
582  end if
583  if (heating_switch) then
584  temperature = t9_p
585  end if
586  else
587  if (.not. (heating_switch)) then
588  call interp1d (zsteps,ztime,time_i,ztemp,temperature)
589  else
590  if (present(t_tr)) then
591  call interp1d (zsteps,ztime,time_i,ztemp,temperature)
592  t_tr = temperature
593  end if
594  temperature = t9_p
595  end if
596 
597  call interp1d (zsteps,ztime,time_i,zdens,density)
598  call interp1d (zsteps,ztime,time_i,zrad,radius,itype=itype_linear)
599  end if
600  end if
601 
602  ! Generally only allow temperatures above 1e-6 GK
603  temperature = max(temperature,freeze_rate_temp)
604 
605 end subroutine update_hydro
606 
607 
608 !>
609 !! Advance system to the next step for the implicit Euler method
610 !!
611 !! This subroutine solves one timestep for the implicit Euler method
612 !! via a newton-raphson (NR) scheme. The maximum amount of iterations for
613 !! the NR is controlled by the parameter \ref parameter_class::nr_maxcount.
614 !! The NR is successful if the mass conservation deviates less than
615 !! parameter_class::nr_tol. If the maximum iterations have been reached without
616 !! convergence, it repeats the calculation with half of the timestep.
617 !! This is done parameter_class::adapt_stepsize_maxcount many times. If after
618 !! this amount of iterations still no convergence could be achieved it will
619 !! raise an error.
620 !!
621 !! @see [Winteler 2013](https://edoc.unibas.ch/29895/)
622 !!
623 !! \b Edited:
624 !! - OK: 15.06.2017 - Made into a separate subroutine
625 !! - MR: 06.07.2022 - Let it crash if the result does not converge properly
626 !! .
627 subroutine advance_implicit_euler(cnt)
632 use global_class, only: net_size, isotope,ipro,ineu,&
634 use hydro_trajectory,only: zsteps,ztime
636 use winnse_module, only: pf,winnse_guess
637 use nucstuff_class, only: el_ab,masscalc
638 use jacobian_class, only: jacobi_init
639 use pardiso_class, only: netsolve
642 use error_msg_class, only: num_to_str
644 implicit none
645 integer, intent(in) :: cnt !< global iteration counter
646 !
647 integer :: k !< Current iteration of adapt stepsize
648 integer :: nr_count !< Current iteration of the newton-raphson
649 real(r_kind) :: m_tot !< Mass conservation
650 type(timmes_eos_state) :: state !< State variable for Helmholtz EOS
651 integer :: eos_status !< Status of the eos, 0=working, 0!=failing
652 logical :: exit_condition!< Helper variable that keeps track of the exit condition of the NR
653 real(r_kind) :: t9_conv !< Convergence of the temperature
654  ! Check the mass conservation and kill the calculation if something goes odd.
655  ! Mass conservation should not be within 0.001% to nr_tol
656  ! simulataneously with a small timestep.
657  ! Note: Other networks such as SkyNet rescale the abundances at this point.
658  call masscalc(y_p, m_tot)
659  if ((1.00001*dabs(m_tot-1d0) .ge. nr_tol) .and. (stepsize .le. 1e-15)) then
660  call raise_exception("The result did not converge properly. "//new_line("A")//&
661  "Mass conservation was "//num_to_str(dabs(m_tot-1d0))//"."//new_line("A")//&
662  "Stepsize was "//num_to_str(stepsize)//"s."//new_line("A")//&
663  "This error could be a result of inconsistent reaction rates, "//new_line("A")//&
664  "if this is not the case,"//&
665  "try to use different values of nr_tol, nr_maxcount, or timestep_factor.",&
666  'advance_implicit_euler',430011)
667  end if
668 
669  adapt_stepsize: do k=0,adapt_stepsize_maxcount
670 
671  ! compute some handy reductions and auxiliary variables
672  y = y_p
673  rkm = rkm_p
674  time = time_p + stepsize
675  ent = ent_p
676  t9 = t9_p
677 
678  call el_ab (y, ye)
679  ! Get temperature, density, and radius
680  call update_hydro(time, t9, rhob, rkm, ent, ye, t9h)
681 
682  !-- Newton-Rhapson iterative solve ---------
683  newton_raphson: do nr_count=1,nr_maxcount
684  exit_condition = .true.
685 
686  ! calculate jacobian and rhs of system of DE (f())
688 
689  ! solve the linearized system
690  call netsolve(f, y)
691 
692  ! thermodynamics update
693  call el_ab(y, ye)
694 
695  ! Check if thinks are reasonable or
696  ! went terribly wrong
697  if ((ye .lt. 0) .or. ( ye .gt. 1)) then
698  if (verbose_level .gt. 1) then
699  write(*,*) "Ye went wrong, obtained Ye = "//num_to_str(ye)
700  end if
701  exit_condition = .false.
702  exit newton_raphson
703  end if
704 
705  if (heating_switch) then
706  ! Update NSE if heating is turned on
707  if (evolution_mode .eq. em_nse) then
708  call winnse_guess( t9, rhob, ye, y(ineu), y(ipro), y)
709  end if
710  ! update entropy from nuclear heating, then update temperature
712  end if
713 
714  ! mass conservation check
715  call masscalc(y, m_tot)
716 
717  ! periodic output
718  call output_nr_diagnostic(cnt,k,nr_count,time,t9,stepsize,m_tot,y_p,y)
719 
720  ! exit if mass is conserved sufficiently well
721  ! if (((nr_count.ge.nr_mincount) .and. (dabs(m_tot-1d0).lt.nr_tol)) ) &
722  if (((.not. heating_switch) .and. (nr_count.ge.nr_mincount) .and. (dabs(m_tot-1d0).lt.nr_tol)) &
723  .or. &
724  ((heating_switch) .and. (nr_count.ge.nr_mincount) .and. (dabs(m_tot-1d0).lt.nr_tol) .and. &
725  ((t9_conv .lt. heating_t9_tol))) .or. &
726  ((heating_switch) .and. (nr_count.ge.nr_mincount) .and. (dabs(m_tot-1d0).lt.nr_tol) .and. &
727  k .eq. adapt_stepsize_maxcount)) &
728  exit newton_raphson
729 
730 
731  end do newton_raphson
732 
733 
734  ! See if one has to repeat the timestep
735  exit_condition = exit_condition .and. (nr_count.le.nr_maxcount)
736  if (heating_switch) then
737  ! Dont be too strict in the last adapt stepsize iteration
738  if (k .ne. adapt_stepsize_maxcount) then
739  exit_condition = exit_condition .and. (max(t9,t9_p)/min(t9,t9_p)-1 .lt. timestep_hydro_factor)
740  end if
741  end if
742  if (exit_condition) exit adapt_stepsize
743 
744  ! Check if the maxcount is reached and try one last time
745  ! to get it to convergence
746  if (k .eq. adapt_stepsize_maxcount) then
748  else
749  stepsize = 0.5*stepsize
750  end if
751  end do adapt_stepsize
752 
753  ! Update entropy
754  if (.not. (heating_switch)) then
755  if (trim(adjustl(trajectory_mode)) .eq. "from_file") then
756  if (time.le.ztime(zsteps) .and. (t9 .gt. freeze_rate_temp*1.0000001)) then
757  ! update entropy using the temperature
758  state%abar = sum(y(1:net_size)*isotope(1:net_size)%mass) &
759  / sum(y(1:net_size))
760  call timmes_eos(ink,t9*1.d9,rhob,ye,state,eos_status)
761  if(eos_status.ne.0) call raise_exception("An error occured in the EOS.",&
762  'advance_implicit_euler',430009)
763  ent = state%s
764  else
765  ! adiabatic expansion
766  ent = ent_p
767  end if
768  else if (trim(adjustl(trajectory_mode)) .eq. "analytic") then
769  if (t9 .gt. freeze_rate_temp*1.0000001) then
770  ! update entropy using the temperature
771  state%abar = sum(y(1:net_size)*isotope(1:net_size)%mass) &
772  / sum(y(1:net_size))
773  call timmes_eos(ink,t9*1.d9,rhob,ye,state,eos_status)
774  if(eos_status.ne.0) call raise_exception("An error occured in the EOS.",&
775  'advance_implicit_euler',430009)
776  ent = state%s
777  end if
778  end if
779  end if
780 
781 
782  ! if still not converged, complain and exit
783  if (k.gt.adapt_stepsize_maxcount) then
784  call raise_exception("Max. number of steps reached. Probably try to change "//new_line('A')//&
785  'the "adapt_stepsize_maxcount", "nr_maxcount", or "nr_tol" parameter.',&
786  'advance_implicit_euler',430010)
787  end if
788 
789  ! Update NSE abundances, only in NSE evolution mode
790  if (evolution_mode .eq. em_nse) then
791  ! Update at all?
792  if (nse_calc_every .ne. 0) then
793  ! Update every nse_calc_step
794  if (modulo(cnt,nse_calc_every) .eq. 0) then
795  call winnse_guess( t9, rhob, ye, y(ineu), y(ipro), y)
796  end if
797  end if
798  end if
799 
800 end subroutine advance_implicit_euler
801 
802 
803 !>
804 !! Advance system to the next step using the Gear method
805 !! (using \ref gear_module)
806 !!
807 !! This subroutine is similar to \ref advance_implicit_euler. However, no
808 !! additional check for convergence is done after the newton-raphson. This is
809 !! not necessary as the Gear solver has its estimate of a possible error
810 !! in the next step, given by evolving higher orders and comparing to lower ones.
811 !!
812 !! @see \ref gear_module, [Martin, D. 2017](https://tuprints.ulb.tu-darmstadt.de/6301/)
813 !!
814 !! \b Edited:
815 !! - DM 15.06.2017
816 !! - OK 27.08.2017
817 !! - MR 22.12.2020
818 !! .
819 subroutine advance_gear(cnt)
827 use hydro_trajectory,only: zsteps,ztime
829 use winnse_module, only: pf,winnse_guess
830 use nucstuff_class, only: el_ab,masscalc
831 use jacobian_class, only: jacobi_init
832 use pardiso_class, only: netsolve
836 use gear_module, only: get_time, get_timestep, get_l1, &
841 implicit none
842 integer, intent(in) :: cnt !< global iteration counter
843 !
844 integer :: k !< Dummy for nr_diagnostic
845 real(r_kind) :: m_tot !< Mass conservation
846 !
847 real(r_kind),dimension(net_size) :: predictor
848 real(r_kind),dimension(net_size) :: delta
849 real(r_kind),dimension(net_size) :: ydiff
850 real(r_kind),dimension(net_size) :: olol
851 
852 type(timmes_eos_state) :: state
853 integer :: eos_status
854 logical :: exit_condition!< Helper variable that keeps
855  ! track of the exit condition of the NR
856 real(r_kind) :: t9_conv !< Convergence of the temperature
857 
858  adapt_stepsize: do k=0,adapt_stepsize_maxcount
859  ! compute some handy reductions and auxiliary variables
861 
862  olol = get_solution()
863  y = y_p
864  rkm = rkm_p
865  time = time_p + stepsize
866  ent = ent_p
867 
868  call el_ab (y, ye)
869 
870  ! Get temperature, density, and radius
871  call update_hydro(time, t9, rhob, rkm, ent, ye, t9h)
872 
873  ! shift tau list and insert h, set time
874  call shift_tau
875  time = get_time()
876 
877  ! predictor
878  call nordsieck_product
879 
880  ! calculate xi vector
881  call set_xi
882 
883  ! calculate l vector
884  call set_l
885 
886  predictor = get_predictor_y()
887  y = predictor
888 
889  !-- Newton-Rhapson iterative solve ---------
890  newton_raphson: do gear_nr_count=1,gear_nr_maxcount
891  ! Set a variable that takes track of problems
892  exit_condition = .true.
893 
894  ! calculate jacobian and rhs of system of DE (f())
896 
897  ! solve the system
898  call netsolve(f, delta)
899 
900  ! add correction to the abundance array
901  y= y + delta
902 
903  ! thermodynamics update
904  call el_ab(y, ye)
905 
906  ! Check if thinks are reasonable or
907  ! went terribly wrong
908  if (ye .lt. 0 .or. ye .gt. 1) then
909  exit_condition = .false.
910  exit newton_raphson
911  end if
912 
913 
914  if (heating_switch) then
915  ! Update NSE if heating is turned on
916  if (evolution_mode .eq. em_nse) then
917  call winnse_guess( t9, rhob, ye, y(ineu), y(ipro), y)
918  end if
920  end if
921 
922  ! mass conservation check
923  call masscalc(y, m_tot)
924 
925  ! periodic output
927 
928  ! exit if mass is conserved sufficiently well
929  ! Also check the temperature for the heating case
930  if (((.not. heating_switch) .and. maxval(dabs(delta)/max(y,gear_escale)).lt.gear_nr_eps .and. &
931  dabs(m_tot-1.0d0).lt.gear_nr_eps .and. gear_nr_count .ge. gear_nr_mincount) .or. &
932  (heating_switch .and. maxval(dabs(delta)/max(y,gear_escale)).lt.gear_nr_eps .and. &
933  dabs(m_tot-1.0d0).lt.gear_nr_eps .and. gear_nr_count .ge. gear_nr_mincount .and.&
934  ((t9_conv .lt. heating_t9_tol))) .or. &
935  ((heating_switch) .and. maxval(dabs(delta)/max(y,gear_escale)).lt.gear_nr_eps .and. &
936  dabs(m_tot-1.0d0).lt.gear_nr_eps .and. &
937  (k .eq. adapt_stepsize_maxcount))) then
938  exit newton_raphson
939  end if
940 
941  end do newton_raphson
942 
943  ! Gear usually converges quite well without adapting the time step
944  if (.not. gear_ignore_adapt_stepsize) then
945  ! if previous loop didn't have problems converging ...
946  ! This basically never happens, but it's good to have it
947  ! Dont do it in NSE, because it will always fail
948  exit_condition = exit_condition .and. ((gear_nr_count.le.gear_nr_maxcount) .or. (evolution_mode .eq. em_nse))
949  if (heating_switch) then
950  ! Dont be too strict in the last adapt stepsize iteration
951  if (k .ne. adapt_stepsize_maxcount) then
952  exit_condition = exit_condition .and. (max(t9,t9_p)/min(t9,t9_p)-1 .lt. timestep_hydro_factor)
953  end if
954  end if
955  if (exit_condition) exit adapt_stepsize
956  else
957  exit adapt_stepsize
958  end if
959 
960  ! Check if the maxcount is reached and try one last time
961  ! to get it to convergence
962  if (k .eq. adapt_stepsize_maxcount) then
964  else
965  stepsize = 0.5*stepsize
966  end if
968 
969  end do adapt_stepsize
970 
971  ! if still not converged, complain and exit
972  if (k.gt.adapt_stepsize_maxcount) then
973  call raise_exception("Max. number of steps reached. Probably try to change "//new_line('A')//&
974  'the "adapt_stepsize_maxcount", "gear_nr_maxcount", or "gear_nr_eps" parameter.',&
975  'advance_gear',430010)
976  end if
977 
978 
979  ! Update entropy
980  if (.not. (heating_switch)) then
981  if (trim(adjustl(trajectory_mode)) .eq. "from_file") then
982  if (time.le.ztime(zsteps) .and. (t9 .gt. freeze_rate_temp*1.0000001)) then
983  ! update entropy using the temperature
984  state%abar = sum(y(1:net_size)*isotope(1:net_size)%mass) &
985  / sum(y(1:net_size))
986  call timmes_eos(ink,t9*1.d9,rhob,ye,state,eos_status)
987  if(eos_status.ne.0) call raise_exception("An error occured in the EOS.",&
988  'advance_implicit_euler',430009)
989  ent = state%s
990  else
991  ! adiabatic expansion
992  ent = ent_p
993  end if
994  else if (trim(adjustl(trajectory_mode)) .eq. "analytic") then
995  if (t9 .gt. freeze_rate_temp*1.0000001) then
996  ! update entropy using the temperature
997  state%abar = sum(y(1:net_size)*isotope(1:net_size)%mass) &
998  / sum(y(1:net_size))
999  call timmes_eos(ink,t9*1.d9,rhob,ye,state,eos_status)
1000  if(eos_status.ne.0) call raise_exception("An error occured in the EOS.",&
1001  'advance_implicit_euler',430009)
1002  ent = state%s
1003  end if
1004  end if
1005  end if
1006 
1007  ! Update NSE abundances, only in NSE evolution mode
1008  if (evolution_mode .eq. em_nse) then
1009  ! Update at all?
1010  if (nse_calc_every .ne. 0) then
1011  ! Update every nse_calc_step
1012  if (modulo(cnt,nse_calc_every) .eq. 0) then
1013  call winnse_guess( t9, rhob, ye, y(ineu), y(ipro), y)
1014  end if
1015  end if
1016  end if
1017 
1018  ydiff = y-predictor
1019 
1020  ! set the new result based on ydiff obtained in NR scheme
1021  call set_new_result(ydiff)
1022 
1023  ! calculate errors and possible stepsizes and take q(+-1) with largest stepsize
1024  call prepare_next_step
1025 
1026 end subroutine advance_gear
1027 
1028 end module timestep_module
timestep_module
The timestep_module contains all timestep subroutines and newton-raphsons to solve the ODE system.
Definition: timestep_module.f90:20
single_zone_vars::y
real(r_kind), dimension(:), allocatable y
Definition: single_zone_vars.f90:21
ls_timmes_eos_module::ink
integer, parameter ink
Definition: ls_timmes_eos_module.f:14
winnse_module::pf
real(r_kind), dimension(:), allocatable, public pf
Partition function for a given temperature.
Definition: winnse_module.f90:38
gear_module::get_predictor_dydt
real(r_kind) function, dimension(net_size), public get_predictor_dydt()
Determines the predicted dydt_next.
Definition: gear_module.f90:332
timestep_module::select_optimistic_timestep
subroutine, public select_optimistic_timestep(init)
Selects the timestep.
Definition: timestep_module.f90:489
gear_module::get_time
real(r_kind) function, public get_time()
The current time.
Definition: gear_module.f90:196
winnse_module
Module to calculate NSE.
Definition: winnse_module.f90:34
global_class::heating_switch
logical, public heating_switch
Variables related to nuclear heating.
Definition: global_class.f90:41
parameter_class::freeze_rate_temp
real(r_kind) freeze_rate_temp
Tmperature at which rates get frozen (for reacl. rates this should be 1d-2GK)
Definition: parameter_class.f90:84
parameter_class::timestep_ymin
real(r_kind) timestep_ymin
Lower limit of the abundance to contribute to the timestep calculation, default value is 1....
Definition: parameter_class.f90:135
parameter_class::nse_calc_every
integer nse_calc_every
Compute NSE abundances every x step.
Definition: parameter_class.f90:89
winnse_module::winnse_guess
subroutine, public winnse_guess(t9, rho, ye, yn_guess, yp_guess, ysol)
Calculate NSE composition with an initial guess.
Definition: winnse_module.f90:116
global_class::ipro
integer, public ipro
index of alphas, neutrons and protons
Definition: global_class.f90:94
expansion_module::expansion
subroutine expansion(t, dt, ye, dens, rad_in, temp, entropy, vel_in, abar)
Returns temperature, radius, density, and entropy after the trajectory has ended.
Definition: expansion_module.f90:240
analysis::snapshot_amount
integer snapshot_amount
Amount of custom snapshots.
Definition: analysis.f90:31
ls_timmes_eos_module
Definition: ls_timmes_eos_module.f:2
inter_module
Module inter_module with interpolation routines.
Definition: inter_module.f90:13
inter_module::inverse_interp1d
subroutine, public inverse_interp1d(n, xp, xb, yp, res, flin, itype)
The inverse of the 1D interpolation function.
Definition: inter_module.f90:931
single_zone_vars::rhob
real(r_kind) rhob
Definition: single_zone_vars.f90:17
parameter_class::timestep_factor
real(r_kind) timestep_factor
Factor for the change of the timestep (see nu in Winteler 2012 Eq. 2.49). Default value is 1....
Definition: parameter_class.f90:133
parameter_class::gear_nr_maxcount
integer gear_nr_maxcount
Maximum newton-raphson iterations for gear solver.
Definition: parameter_class.f90:141
analysis::snapshot_time
real, dimension(:), allocatable snapshot_time
point in time for custom snapshot
Definition: analysis.f90:30
timestep_module::advance_implicit_euler
subroutine, public advance_implicit_euler(cnt)
Advance system to the next step for the implicit Euler method.
Definition: timestep_module.f90:628
single_zone_vars::t9
real(r_kind) t9
Definition: single_zone_vars.f90:15
nuclear_heating
Takes care of the temperature and entropy updates due to the nuclear energy release.
Definition: nuclear_heating.f90:11
single_zone_vars::rkm
real(r_kind) rkm
Definition: single_zone_vars.f90:18
parameter_class::timestep_max
real(r_kind) timestep_max
Maximum factor for the change of the timestep. The new timestep is only allowed to be timestep_max * ...
Definition: parameter_class.f90:132
single_zone_vars::ent
real(r_kind) ent
Definition: single_zone_vars.f90:20
gear_module::get_solution_at
subroutine, public get_solution_at(time_inter, Y_inter)
Determines the solution at a given point in time for [ti-h, ti].
Definition: gear_module.f90:346
parameter_class::gear_nr_eps
real(r_kind) gear_nr_eps
Convergence criterion for the newton-raphson of the gear solver.
Definition: parameter_class.f90:139
gear_module::get_timestep
real(r_kind) function, public get_timestep()
The current timestep.
Definition: gear_module.f90:210
error_msg_class
Error handling routines.
Definition: error_msg_class.f90:16
expansion_module
Subroutines to handle parametric evolution of hydrodynamic quantities after the final timestep of the...
Definition: expansion_module.f90:16
timestep_module::update_hydro
subroutine, private update_hydro(time_i, temperature, density, radius, entropy, efraction, T_tr)
Returns temperature, density, and radius for a given time "time_i".
Definition: timestep_module.f90:531
parameter_class::termination_criterion
integer termination_criterion
condition to terminate the simulation ([0]=trajectory_file, 1=final_time, 2=final_temp,...
Definition: parameter_class.f90:127
expansion_module::rad
real(r_kind) rad
radial distance [km]
Definition: expansion_module.f90:33
gear_module::set_xi
subroutine, public set_xi
Function to calculate xi (needed for getting l in the corrector step)
Definition: gear_module.f90:421
parameter_class::gear_escale
real(r_kind) gear_escale
Normalization cutoff for gear solver, similar to timestep_Ymin for Euler.
Definition: parameter_class.f90:137
single_zone_vars::t9h_p
real(r_kind) t9h_p
Temperature [GK] storage for heating mode 2.
Definition: single_zone_vars.f90:16
gear_module::gear_nr_count
integer, public gear_nr_count
counter for NR steps
Definition: gear_module.f90:68
parameter_class::nsetemp_cold
real(r_kind) nsetemp_cold
T [GK] for the nse->network switch.
Definition: parameter_class.f90:119
nucstuff_class::masscalc
subroutine, public masscalc(Y, m_tot)
Total mass used to check the mass conservation.
Definition: nucstuff_class.f90:384
single_zone_vars::evolution_mode
integer evolution_mode
NSE, network hot/cold, etc.
Definition: single_zone_vars.f90:23
timestep_module::timestep
subroutine timestep(ctime, temp, dens, entropy, Y, dYdt, Ye, delYe, h, evmode)
Network timestepping along the thermodynamic trajectory, specified by {ztime, ztemp,...
Definition: timestep_module.f90:54
pardiso_class::netsolve
subroutine, public netsolve(rhs, res)
The solver core.
Definition: pardiso_class.f90:247
ls_timmes_eos_module::timmes_eos
subroutine timmes_eos(var, vin, d, Ye, state, status)
Definition: ls_timmes_eos_module.f:49
parameter_class::final_dens
real(r_kind) final_dens
termination density [g/cm3]
Definition: parameter_class.f90:131
parameter_class::nr_maxcount
integer nr_maxcount
– Newton-Raphson iterative loop parameters
Definition: parameter_class.f90:201
error_msg_class::raise_exception
subroutine, public raise_exception(msg, sub, error_code)
Raise a exception with a given error message.
Definition: error_msg_class.f90:245
parameter_class::initial_stepsize
real(r_kind) initial_stepsize
this value is used as a stepsize at initial step
Definition: parameter_class.f90:128
parameter_class::nr_mincount
integer nr_mincount
Minimum iterations in NR.
Definition: parameter_class.f90:202
parameter_class::timestep_traj_limit
logical timestep_traj_limit
Should the timestep be limited by the timestep of the trajectory.
Definition: parameter_class.f90:144
parameter_class::timestep_hydro_factor
real(r_kind) timestep_hydro_factor
Factor for the maximum change of the hydrodynamic quantities (density and temperature)
Definition: parameter_class.f90:134
single_zone_vars::rkm_p
real(r_kind) rkm_p
Radius of the outflow [km].
Definition: single_zone_vars.f90:18
timestep_module::hydro_timestep
real(r_kind) function, private hydro_timestep(h1, ctime, temp, dens, entropy, ttemp, tdens, Ye)
Estimate the hydro timestep.
Definition: timestep_module.f90:407
timestep_module::ye_timestep
real(r_kind) function, private ye_timestep(h, Ye, delYe)
Estimate the timestep from the change of electron fraction.
Definition: timestep_module.f90:336
jacobian_class::abchange
subroutine, public abchange(time, itemp, rho, Ye, rkm, Y, dYdt, evolution_mode)
This subroutine calculates the change in abundances (dYdt) due to reaction rates at given temperature...
Definition: jacobian_class.f90:289
gear_module::get_solution
real(r_kind) function, dimension(net_size), public get_solution()
Gives a copy of the solution at the current time.
Definition: gear_module.f90:297
parameter_class::gear_nr_mincount
integer gear_nr_mincount
Minimum newton-raphson iterations for gear solver.
Definition: parameter_class.f90:142
parameter_class::gear_ignore_adapt_stepsize
logical gear_ignore_adapt_stepsize
Flag whether gear should ignore the adapt stepsize loop.
Definition: parameter_class.f90:143
gear_module::set_new_result
subroutine, public set_new_result(ydiff)
Updates the solution (corrector step)
Definition: gear_module.f90:476
global_class::ineu
integer, public ineu
Definition: global_class.f90:94
parameter_class::rho_analytic
character(max_fname_len) rho_analytic
analytic density [g/cm3]
Definition: parameter_class.f90:187
global_class::isotope
type(isotope_type), dimension(:), allocatable, public isotope
all nuclides used in the network
Definition: global_class.f90:34
single_zone_vars::time_p
real(r_kind) time_p
Time.
Definition: single_zone_vars.f90:13
hydro_trajectory::ztime
real(r_kind), dimension(:), allocatable, public ztime
time information from trajectory
Definition: hydro_trajectory.f90:19
nucstuff_class
Module with helper functions such as rate tabulation, partition functions, and theoretical weak rates...
Definition: nucstuff_class.f90:11
parameter_class::rkm_analytic
character(max_fname_len) rkm_analytic
analytic radial scale [km]
Definition: parameter_class.f90:188
inter_module::interp1d
subroutine, public interp1d(n, xp, xb, yp, res, flin, itype)
Interface for 1D interpolation routines.
Definition: inter_module.f90:99
parameter_class::heating_density
real(r_kind) heating_density
Density at which nuclear heating will be switched on (-1) to always include heating.
Definition: parameter_class.f90:124
gear_module::set_timestep
subroutine, public set_timestep(timestep)
Set the current timestep.
Definition: gear_module.f90:233
parameter_class::heating_t9_tol
real(r_kind) heating_t9_tol
Convergence parameter for the temperature in the heating mode.
Definition: parameter_class.f90:149
single_zone_vars::ent_p
real(r_kind) ent_p
Entropy [kB/baryon].
Definition: single_zone_vars.f90:20
single_zone_vars::f
real(r_kind), dimension(:), allocatable f
Definition: single_zone_vars.f90:22
parser_module
Subroutines for equation parsing in case of an analytic trajectory or luminosity mode.
Definition: parser_module.f90:22
hydro_trajectory::zrad
real(r_kind), dimension(:), allocatable, public zrad
radii from trajectory
Definition: hydro_trajectory.f90:23
parser_module::find_start_value
subroutine, public find_start_value(input_string, eq_value, initial_guess, converged, result)
Finds a value of the variable for a given y-value.
Definition: parser_module.f90:856
expansion_module::vel
real(r_kind) vel
velocity [km/s]
Definition: expansion_module.f90:32
gear_module::get_l1
real(r_kind) function, public get_l1()
The current l_1 (in Fortran: l(2))
Definition: gear_module.f90:281
parameter_class::solver
integer solver
solver flag (0 - implicit Euler, 1 - Gear's method, ...), is integer as it is faster than comparing s...
Definition: parameter_class.f90:147
nucstuff_class::el_ab
subroutine, public el_ab(Y, Ye)
Computes the electron fraction.
Definition: nucstuff_class.f90:219
parser_module::parse_string
real(r_kind) function, public parse_string(input_string, var_value)
Takes a string and evaluates the expression.
Definition: parser_module.f90:807
global_class::net_size
integer, public net_size
total number of isotopes (network size)
Definition: global_class.f90:93
gear_module::get_predictor_y
real(r_kind) function, dimension(net_size), public get_predictor_y()
Determines the predicted y_next.
Definition: gear_module.f90:314
nuclear_heating::nuclear_heating_update
subroutine, public nuclear_heating_update(nr_count, rho, Ye, pf, Y_p, Y, entropy_p, entropy, T9_p, T9, T9_tr_p, T9_tr, timestep, T9conv)
Modifies the temperature and entropy to account for nuclear heating.
Definition: nuclear_heating.f90:320
parameter_class::adapt_stepsize_maxcount
integer adapt_stepsize_maxcount
max. iterations in adapting the stepsize
Definition: parameter_class.f90:204
parameter_class::trajectory_mode
character *20 trajectory_mode
determines how trajectory is calculated (from_file|analytic|expand)
Definition: parameter_class.f90:90
error_msg_class::num_to_str
character(:) function, allocatable, public num_to_str(num)
Converts a given real to a string with format "(1pE10.2)".
Definition: error_msg_class.f90:205
parameter_class::heating_mode
integer heating_mode
Mode for heating: 0 - no heating, 1 - heating using an entropy equation, 2 - heating from the energy ...
Definition: parameter_class.f90:148
analysis
Contains various routines for analysis and diagnostic output.
Definition: analysis.f90:9
single_zone_vars::time
real(r_kind) time
Definition: single_zone_vars.f90:13
pardiso_class
Contains subroutines for sparse matrix assembly and the solver core.
Definition: pardiso_class.f90:24
single_zone_vars::stepsize
real(r_kind) stepsize
Stepsize.
Definition: single_zone_vars.f90:14
r_kind
#define r_kind
Definition: macros.h:46
parameter_class::h_custom_snapshots
logical h_custom_snapshots
Same, but in hdf5 format.
Definition: parameter_class.f90:108
timestep_module::restrict_timestep
subroutine, private restrict_timestep(ctime, h, temp, dens, ttemp, tdens)
Restricting the previously selected timestep.
Definition: timestep_module.f90:141
gear_module::set_l
subroutine, public set_l
Function to calculate (needed for the corrector step)
Definition: gear_module.f90:448
hydro_trajectory::zsteps
integer zsteps
number of timesteps in the hydro trajectory
Definition: hydro_trajectory.f90:17
hydro_trajectory
Contains arrays representing thermodynamic conditions from hydro trajectory file.
Definition: hydro_trajectory.f90:14
hydro_trajectory::zdens
real(r_kind), dimension(:), allocatable, public zdens
density information from trajectory
Definition: hydro_trajectory.f90:21
parameter_class::nr_tol
real(r_kind) nr_tol
exit NR if tolerance less than this value
Definition: parameter_class.f90:203
single_zone_vars::dydt
real(r_kind), dimension(:), allocatable dydt
Time derivative of the abundances.
Definition: single_zone_vars.f90:22
gear_module::prepare_next_step
subroutine, public prepare_next_step
Stepsize and order control to pepare the next step.
Definition: gear_module.f90:503
jacobian_class::jacobi_init
subroutine, public jacobi_init(time, itemp, rho, rkm, Y, Y_p, dYdt, rhs, h, evolution_mode)
This subroutine calculates the entries in the jacobian and the right-hand side of the numerical integ...
Definition: jacobian_class.f90:523
single_zone_vars::y_p
real(r_kind), dimension(:), allocatable y_p
Abundances.
Definition: single_zone_vars.f90:21
parameter_class::final_time
real(r_kind) final_time
termination time in seconds
Definition: parameter_class.f90:129
jacobian_class
Contains subroutines to assemble the system Jacobian and changes in the abundances.
Definition: jacobian_class.f90:12
parameter_class::final_temp
real(r_kind) final_temp
termination temperature [GK]
Definition: parameter_class.f90:130
global_class
Contains types and objects shared between multiple modules.
Definition: global_class.f90:19
parameter_class::t9_analytic
character(max_fname_len) t9_analytic
analytic temperature [T9]
Definition: parameter_class.f90:186
single_zone_vars::ye
real(r_kind) ye
Definition: single_zone_vars.f90:19
analysis::output_nr_diagnostic
subroutine, public output_nr_diagnostic(cnt, k, nr_c, ctime, T9, stepsize, mtot, Y_p, Y)
Output of the Newton-Raphson loop diagnostics.
Definition: analysis.f90:754
ls_timmes_eos_module::timmes_eos_state
Definition: ls_timmes_eos_module.f:23
timestep_module::abchange_timestep
real(r_kind) function, private abchange_timestep(h, Y, dYdt)
Estimate the timestep from the abundances change.
Definition: timestep_module.f90:365
hydro_trajectory::ztemp
real(r_kind), dimension(:), allocatable, public ztemp
temperature information from trajectory
Definition: hydro_trajectory.f90:20
gear_module::nordsieck_product
subroutine, public nordsieck_product
Nordsieck vector-matrix product.
Definition: gear_module.f90:393
parameter_class::custom_snapshots
logical custom_snapshots
If true, a file must be provided with numbers in days. Snapshots will be created for these points in ...
Definition: parameter_class.f90:107
single_zone_vars
Simulation variables for a single zone model.
Definition: single_zone_vars.f90:11
gear_module::revert_timestep
subroutine, public revert_timestep(timestep)
Revert a timestep.
Definition: gear_module.f90:262
gear_module::shift_tau
subroutine, public shift_tau
Shifts tau vector and controls the current time.
Definition: gear_module.f90:371
timestep_module::advance_gear
subroutine, public advance_gear(cnt)
Advance system to the next step using the Gear method (using gear_module)
Definition: timestep_module.f90:820
gear_module
gear_module contains adaptive high-order Gear solver
Definition: gear_module.f90:37
single_zone_vars::t9_p
real(r_kind) t9_p
Temperature [GK].
Definition: single_zone_vars.f90:15
gear_module::init_dydt
subroutine, public init_dydt(dYdt)
Initialize dYdt component of the Nordsieck vector.
Definition: gear_module.f90:181
parameter_class
Contains all runtime parameters as well as phys and math constants.
Definition: parameter_class.f90:24
single_zone_vars::ye_p
real(r_kind) ye_p
Electron fraction [mol/g].
Definition: single_zone_vars.f90:19
single_zone_vars::t9h
real(r_kind) t9h
Definition: single_zone_vars.f90:16