MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_supertimestepping.t
Go to the documentation of this file.
1!> Generic supertimestepping method
2!> which can be used for multiple source terms in the governing equations
3!> where each term introduces a sts_term type and all sts terms are stored in a pointer list
4!>
5!> user control is in
6!> 1) amrvac.par where namelist sts_list sets the following parameters which have default values:
7!> sts_dtpar=0.5,sts_ncycles=1000,sts_method=1,sourcetype_sts=2
8!> These parametes are shared for all the terms added
9!>
10!> Any STS-handled term is done via
11!> 2) program in the code a term with the subroutine add_sts_method
12!> This method takes as parameters a function which calculates the explicit timestep
13!> associated with the term and a subroutine which sets the source term
14!> for the variables startVar:endVar=startVar+nflux-1
15!> flux conservation (fixconserve) is done for the variables as specified by ixChangeStart, ixChangeN, ixChangeFixC
16!>
17!> 3) one can have hooks before_first_cycle, after_last_cycle (e.g. conversion from e_tot to e_int before first sts cycle
18!> and back from e_int to e_tot after the last STS cycle as used for the thermal conductivity module)
19!> add those hooks just afterwards with the subroutine set_conversion_methods_to_head
20!>
21!> 4) add the hook for error handling (e.g check small values in the thermal conductivity module)
22!> call set_error_handling_to_head which takes as parameter a subroutine
23!> this error handling subroutine is called before setting BC
24!>
26 use mod_geometry
27 use mod_comm_lib, only: mpistop
28 implicit none
29 private
30
31 public :: is_sts_initialized
32 public :: sts_init
34 public :: sts_add_source
35 public :: set_dt_sts_ncycles
37
38 ! input parameters from parameter file
39 !> the coefficient that multiplies the sts dt
40 double precision :: sts_dtpar=0.5d0
41
42 !The following is used only for method 2, not input parameter
43 double precision,parameter :: nu_sts = 0.5d0
44 !> the maximum number of subcycles
45 integer :: sts_ncycles=1000
46 integer :: sts_method = 1
47 integer, parameter :: sourcetype_sts_prior =0
48 integer, parameter :: sourcetype_sts_after =1
49 integer, parameter :: sourcetype_sts_split =2
51 !> Whether to conserve fluxes at the current partial step
52 logical :: fix_conserve_at_step = .true.
53 logical :: sts_initialized = .false.
54
55 abstract interface
56
57 !>interface for setting sources in the derived type
58 subroutine subr1(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
60 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
61 double precision, intent(in) :: x(ixi^s,1:ndim)
62 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
63 double precision, intent(in) :: my_dt
64 logical, intent(in) :: fix_conserve_at_step
65 end subroutine subr1
66
67 !>interface for the function which gets the timestep in dtnew in the derived type
68 function subr2(w,ixG^L,ix^L,dx^D,x) result(dtnew)
70 integer, intent(in) :: ixg^l, ix^l
71 double precision, intent(in) :: dx^d, x(ixg^s,1:ndim)
72 double precision, intent(in) :: w(ixg^s,1:nw)
73 double precision :: dtnew
74 end function subr2
75
76 !>interface for error handling subroutine in the derived type
77 subroutine subr_e(w, x, ixI^L, ixO^L, step)
80 integer, intent(in) :: ixi^l,ixo^l
81 double precision, intent(inout) :: w(ixi^s,1:nw)
82 double precision, intent(in) :: x(ixi^s,1:ndim)
83 integer, intent(in) :: step
84 end subroutine subr_e
85
86 !>interface for the subroutines before_first_cycle and after_last_cycle in the derived type
87 subroutine subr5(ixI^L, ixO^L, w, x)
89 integer, intent(in) :: ixi^l, ixo^l
90 double precision, intent(in) :: x(ixi^s,1:ndim)
91 double precision, intent(inout) :: w(ixi^s,1:nw)
92 end subroutine subr5
93
94 !>for the subroutines in this module, which do not depend on the term, but
95 !>on the parameter sts_method = 1/2 in the parameter file
96 !>sts_add_source
97 subroutine subr3(dt)
98 double precision,intent(in) :: dt
99 end subroutine subr3
100
101 !>sts_get_ncycles
102 function subr4(dt,dtnew,dt_modified) result(s)
103 double precision,intent(in) :: dtnew
104 double precision,intent(inout) :: dt
105 logical,intent(inout) :: dt_modified
106 integer :: s
107 end function subr4
108
109 end interface
110
111 type sts_term
112
113 double precision :: dt_expl
114 integer, public :: s
115
116 !>types used for send/recv ghosts, see mod_ghostcells_update
117 integer, dimension(-1:1^D&) :: type_send_srl_sts_1, type_recv_srl_sts_1
118 integer, dimension(-1:1^D&) :: type_send_r_sts_1
119 integer, dimension( 0:3^D&) :: type_recv_r_sts_1
120 integer, dimension( 0:3^D&) :: type_recv_p_sts_1, type_send_p_sts_1
121
122 integer, dimension(-1:1^D&) :: type_send_srl_sts_2, type_recv_srl_sts_2
123 integer, dimension(-1:1^D&) :: type_send_r_sts_2
124 integer, dimension( 0:3^D&) :: type_recv_r_sts_2
125 integer, dimension( 0:3^D&) :: type_recv_p_sts_2, type_send_p_sts_2
126
127 integer :: startvar
128 integer :: endvar
129 !> number of flux involved in STS update
130 integer :: nflux
131 integer :: startwbc
132 integer :: nwbc
133 logical :: types_initialized
134 logical :: evolve_magnetic_field
135 procedure(subr1), pointer, nopass :: sts_set_sources
136 procedure(subr2), pointer, nopass :: sts_getdt
137 procedure(subr5), pointer, nopass :: sts_before_first_cycle, sts_after_last_cycle
138 procedure(subr_e), pointer, nopass :: sts_handle_errors
139 type(sts_term), pointer :: next
140
141 end type sts_term
142
143 type(sts_term), pointer :: head_sts_terms
144 !The following two subroutine/function pointers
145 !make the difference between the two STS methods implemented
146 procedure(subr3), pointer :: sts_add_source
147 procedure(subr4), pointer :: sts_get_ncycles
148
149contains
150
151 !> Initialize sts module
152 subroutine sts_init()
154 use mod_physics
155 if(.not. sts_initialized) then
156 nullify(head_sts_terms)
157 call sts_params_read(par_files)
158 sts_dtpar=sts_dtpar/dble(ndim)
159 sts_initialized = .true.
160 if(sts_method .eq. 1) then
161 sts_add_source => sts_add_source1
162 sts_get_ncycles => sts_get_ncycles1
163 else if(sts_method .eq. 2) then
164 sts_add_source => sts_add_source2
165 sts_get_ncycles => sts_get_ncycles2
166 else
167 call mpistop("Unknown sts method")
168 end if
169 endif
170
171 end subroutine sts_init
172
173 pure function is_sts_initialized() result(res)
174 logical :: res
175 if (sts_initialized) then
176 res = associated(head_sts_terms)
177 else
178 res = .false.
179 endif
180 end function is_sts_initialized
181
182 !> Read module parameters from par file
183 subroutine sts_params_read(files)
185 character(len=*), intent(in) :: files(:)
186 integer :: n
187
188 namelist /sts_list/ sts_dtpar,sts_ncycles,sts_method,sourcetype_sts
189
190 do n = 1, size(files)
191 open(unitpar, file=trim(files(n)), status="old")
192 read(unitpar, sts_list, end=111)
193111 close(unitpar)
194 end do
195
196 end subroutine sts_params_read
197
198 !> subroutine which added programatically a term to be calculated using STS
199 !> Params:
200 !> sts_getdt function calculates the explicit timestep for this term
201 !> sts_set_sources subroutine sets the source term
202 !> startVar, nflux (making endVar=startVar+nflux-1) indices of start and number of the variables that need fix conservation
203 !> startwbc, nwbc indices of start and number of the variables that need ghost cell exchange
204 !> These terms implemented by an element of the derived type sts_term are put in a linked list
205 subroutine add_sts_method(sts_getdt, sts_set_sources, startVar, nflux, startwbc, nwbc, evolve_B)
208
209 integer, intent(in) :: startvar, nflux, startwbc, nwbc
210 logical, intent(in) :: evolve_b
211
212 interface
213
214 subroutine sts_set_sources(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
217 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
218 double precision, intent(in) :: x(ixi^s,1:ndim)
219 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
220 double precision, intent(in) :: my_dt
221 logical, intent(in) :: fix_conserve_at_step
222 end subroutine sts_set_sources
223
224 function sts_getdt(w,ixG^L,ix^L,dx^D,x) result(dtnew)
226 integer, intent(in) :: ixg^l, ix^l
227 double precision, intent(in) :: dx^d, x(ixg^s,1:ndim)
228 double precision, intent(in) :: w(ixg^s,1:nw)
229 double precision :: dtnew
230 end function sts_getdt
231
232 end interface
233
234 type(sts_term), pointer :: temp
235 allocate(temp)
236
237 temp%sts_getdt => sts_getdt
238 temp%sts_set_sources => sts_set_sources
239 temp%sts_before_first_cycle => null()
240 temp%sts_after_last_cycle => null()
241 temp%sts_handle_errors => null()
242 temp%startVar = startvar
243 temp%endVar= startvar+nflux-1
244 temp%nflux = nflux
245 temp%startwbc = startwbc
246 temp%nwbc = nwbc
247 temp%types_initialized = .false.
248 temp%evolve_magnetic_field=evolve_b
249
250 temp%next => head_sts_terms
251 head_sts_terms => temp
252
253 end subroutine add_sts_method
254
255 !> Set the hooks called before the first cycle and after the last cycle in the STS update
256 !> This method should be called after add_sts_method. The hooks are added to the last term added with this subroutine
257 !> Params: sts_before_first_cycle, sts_after_last_cycle subroutines which implement the hooks called before first cycle and after last cycle
258 subroutine set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
259 interface
260 subroutine sts_before_first_cycle(ixI^L, ixO^L, w, x)
262 integer, intent(in) :: ixi^l, ixo^l
263 double precision, intent(in) :: x(ixi^s,1:ndim)
264 double precision, intent(inout) :: w(ixi^s,1:nw)
265 end subroutine sts_before_first_cycle
266
267 subroutine sts_after_last_cycle(ixI^L, ixO^L, w, x)
269 integer, intent(in) :: ixi^l, ixo^l
270 double precision, intent(in) :: x(ixi^s,1:ndim)
271 double precision, intent(inout) :: w(ixi^s,1:nw)
272 end subroutine sts_after_last_cycle
273 end interface
274
275 head_sts_terms%sts_before_first_cycle => sts_before_first_cycle
276 head_sts_terms%sts_after_last_cycle => sts_after_last_cycle
277
278 end subroutine set_conversion_methods_to_head
279
280 !> Set the hook of error handling in the STS update. This method is called before updating the BC.
281 !> This method should be called after add_sts_method. The hook is added to the last term added with this subroutine.
282 !> Param: sts_error_handing the subroutine which handles the errors
283 subroutine set_error_handling_to_head(sts_error_handling)
284 interface
285 subroutine sts_error_handling(w, x, ixI^L, ixO^L, step)
288 integer, intent(in) :: ixi^l,ixo^l
289 double precision, intent(inout) :: w(ixi^s,1:nw)
290 double precision, intent(in) :: x(ixi^s,1:ndim)
291 integer, intent(in) :: step
292 end subroutine sts_error_handling
293 end interface
294 head_sts_terms%sts_handle_errors => sts_error_handling
295
296 end subroutine set_error_handling_to_head
297
298 !> method used to set the number of cycles for the STS1 method
299 function sts_get_ncycles1(dt,dtnew,dt_modified) result(is)
300 double precision,intent(in) :: dtnew
301 double precision,intent(inout) :: dt
302 logical,intent(inout) :: dt_modified
303 integer :: is
304
305 double precision :: ss
306
307 !!ss is now limit of dt because of sts_ncycles
308 ss = dtnew*((2.d0*sts_ncycles+1)**2-9.d0)/16.d0
309 if(dt>ss) then
310 dt_modified = .true.
311 dt = ss
312 is = sts_ncycles
313 else
314 ss = dt/dtnew
315 ! get number of sub-steps of supertime stepping (Meyer 2012 MNRAS 422,2102)
316 if(ss .le. 1.d0) then
317 is=1
318 else
319 is=ceiling((dsqrt(9.d0+16.d0*ss)-1.d0)*0.5d0)
320 is=is/2*2+1
321 end if
322 end if
323
324 end function sts_get_ncycles1
325
326 !> method used to set the number of cycles for the STS2 method
327 function sts_get_ncycles2(dt,dtnew,dt_modified) result(is)
328 double precision,intent(in) :: dtnew
329 double precision,intent(inout) :: dt
330 logical,intent(inout) :: dt_modified
331 integer :: is
332
333 double precision :: ss,rr
334 integer:: ncycles
335
336 rr = dt/dtnew
337 !print*, dt, " --DTEXPL-- ", dtnew, ", rr ",rr
338 ncycles = sts_ncycles
339 !print*, "NCYCLES BEFORE ",ncycles
340 ss=sum_chev(nu_sts,ncycles,rr)
341 !print*, "NCYCLES AFTER ",ncycles
342 is = ncycles
343 !print*, "SUMCHEV ", ss, " NCYCLES ", is
344 if(ss < rr) then
345 dt_modified = .true.
346 dt = ss * dtnew
347 endif
348
349 end function sts_get_ncycles2
350
351 !> This sets the explicit dt and calculates the number of cycles for each of the terms implemented with STS.
352 function set_dt_sts_ncycles(my_dt) result(dt_modified)
354
355 double precision,intent(inout) :: my_dt
356 double precision :: my_dt1
357 logical :: dt_modified, dt_modified1, dt_modified2
358
359 double precision :: dtnew,dtmin_mype
360 double precision :: dx^d, ss
361 integer:: iigrid, igrid, ncycles
362 type(sts_term), pointer :: temp,oldtemp
363 nullify(oldtemp)
364 temp => head_sts_terms
365 dt_modified = .false.
366 do while(associated(temp))
367 dt_modified2 = .false.
368 dtmin_mype=bigdouble
369 !$OMP PARALLEL DO PRIVATE(igrid,dx^D) REDUCTION(min:dtmin_mype)
370 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
371 ! maybe the following global variables are needed in get_dt!
372 ! next few lines ensure correct usage of routines like divvector etc
373 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
374 block=>ps(igrid)
375 ! end maybe the following global variables are needed in get_dt!!!!!!!
376 dx^d=rnode(rpdx^d_,igrid);
377 dtmin_mype=min(dtmin_mype, sts_dtpar * temp%sts_getdt(ps(igrid)%w,ixg^ll,ixm^ll,dx^d,ps(igrid)%x))
378 end do
379 !$OMP END PARALLEL DO
380 call mpi_allreduce(dtmin_mype,dtnew,1,mpi_double_precision,mpi_min,icomm,ierrmpi)
381 temp%s = sts_get_ncycles(my_dt,dtnew,dt_modified2)
382
383 !print*, "NCYCLES ", temp%s, dt_modified2, my_dt, dtnew
384 temp%dt_expl = dtnew
385
386 ! Note that as for some term it may happen that the dt is modified: it may be reduced if the
387 ! number of cycles is overpassed, the list has to be reiterated to update ncycles for previous
388 ! terms which did not modify dt TODO add pointer to previous and loop backward to update
389 if(dt_modified2) then
390 dt_modified = .true.
391 !reiterate all the other sts elements and recalculate s
392 oldtemp => head_sts_terms
393 my_dt1 = my_dt
394 dt_modified1 = .false.
395 do while(.not. associated(oldtemp,temp))
396 oldtemp%s = sts_get_ncycles(my_dt1,oldtemp%dt_expl,dt_modified1)
397 !check dt is not modified again, and this should not happen, except for bug in sts_get_ncycles1,2
398 if(dt_modified1) call mpistop("sts dt modified twice")
399 oldtemp=>oldtemp%next
400 end do
401 end if
402 temp=>temp%next
403
404 end do
405
406 end function set_dt_sts_ncycles
407
408 pure FUNCTION chev(j,nu,N)
409 use mod_constants
410
411 double precision, INTENT(IN) :: nu
412 INTEGER, INTENT(IN) :: j, n
413 double precision :: chev
414
415 chev = 1d0 / ((-1d0 + nu)*cos(((2d0*j - 1d0) / n)* (dpi/2d0)) + 1d0 + nu)
416
417 END FUNCTION chev
418
419 FUNCTION sum_chev(nu,N,limMax)
420 double precision, intent(in) :: nu,limmax
421 integer, intent(inout) :: n
422 double precision :: sum_chev, tmp
423
424 integer :: j
425
426 j=1
427 sum_chev = 0d0
428 do while (j < n .and. sum_chev < limmax)
429 sum_chev = sum_chev + chev(j,nu,n)
430 j=j+1
431 enddo
432 n=j-1
433 END FUNCTION sum_chev
434
435 !> Iterates all the terms implemented with STS and adds the sources
436 !> STS method 2 implementation
437 subroutine sts_add_source2(my_dt)
438 ! Turlough Downes 2006,2007
442 use mod_physics
443
444 double precision, intent(in) :: my_dt
445 double precision, allocatable :: bj(:)
446 double precision :: sumbj,dtj
447
448 integer:: iigrid, igrid, j, ixc^l
449 logical :: stagger_flag=.false., prolong_flag=.false., coarsen_flag=.false.
450 type(sts_term), pointer :: temp
451
452 ! do not fill physical boundary conditions
453 bcphys=.false.
454
455 fix_conserve_at_step = time_advance .and. levmax>levmin
456
457 temp => head_sts_terms
458 do while(associated(temp))
459
460 if(.not.temp%evolve_magnetic_field) then
461 ! not do fix conserve and getbc for staggered values
462 stagger_flag=stagger_grid
463 stagger_grid=.false.
464 else if(stagger_grid) then
465 ixcmax^d=ixmhi^d;
466 ixcmin^d=ixmlo^d-1;
467 end if
468
469 call init_comm_fix_conserve(1,ndim,temp%nflux)
470
471 if(associated(temp%sts_before_first_cycle)) then
472 prolong_flag=prolongprimitive
473 coarsen_flag=coarsenprimitive
474 prolongprimitive=.false.
475 coarsenprimitive=.false.
476 do iigrid=1,igridstail; igrid=igrids(iigrid);
477 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
478 block=>ps(igrid)
479 call temp%sts_before_first_cycle(ixg^ll,ixg^ll,ps(igrid)%w,ps(igrid)%x)
480 end do
481 end if
482
483 allocate(bj(1:temp%s))
484 do j=1,temp%s
485 bj(j) = chev(j,nu_sts,sts_ncycles)
486 end do
487
488 type_send_srl=>temp%type_send_srl_sts_1
489 type_recv_srl=>temp%type_recv_srl_sts_1
490 type_send_r=>temp%type_send_r_sts_1
491 type_recv_r=>temp%type_recv_r_sts_1
492 type_send_p=>temp%type_send_p_sts_1
493 type_recv_p=>temp%type_recv_p_sts_1
494
495 if(.not. temp%types_initialized) then
496 call create_bc_mpi_datatype(temp%startwbc,temp%nwbc)
497 if(temp%nflux>temp%nwbc) then
498 ! prepare types for the changed no-need-ghost-update variables in the last getbc
499 type_send_srl=>temp%type_send_srl_sts_2
500 type_recv_srl=>temp%type_recv_srl_sts_2
501 type_send_r=>temp%type_send_r_sts_2
502 type_recv_r=>temp%type_recv_r_sts_2
503 type_send_p=>temp%type_send_p_sts_2
504 type_recv_p=>temp%type_recv_p_sts_2
505 call create_bc_mpi_datatype(temp%startVar,temp%nflux)
506 type_send_srl=>temp%type_send_srl_sts_1
507 type_recv_srl=>temp%type_recv_srl_sts_1
508 type_send_r=>temp%type_send_r_sts_1
509 type_recv_r=>temp%type_recv_r_sts_1
510 type_send_p=>temp%type_send_p_sts_1
511 type_recv_p=>temp%type_recv_p_sts_1
512 end if
513 temp%types_initialized = .true.
514 end if
515
516 sumbj=0.d0
517 do j=1,temp%s
518 if(j .eq. temp%s .and. (sumbj + bj(j)) * temp%dt_expl > my_dt) then
519 dtj = my_dt - sumbj * temp%dt_expl
520 else
521 dtj = bj(j)* temp%dt_expl
522 end if
523 sumbj = sumbj + bj(j)
524 if(stagger_grid) then
525 !$OMP PARALLEL DO PRIVATE(igrid)
526 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
527 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
528 block=>ps(igrid)
529 call temp%sts_set_sources(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x,ps1(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
530 if(temp%nflux>ndir) then
531 ps(igrid)%w(ixm^t,temp%startVar)=ps(igrid)%w(ixm^t,temp%startVar)+dtj*ps1(igrid)%w(ixm^t,temp%startVar)
532 end if
533 ps(igrid)%ws(ixc^s,1:nws)=ps(igrid)%ws(ixc^s,1:nws)+dtj*ps1(igrid)%w(ixc^s,iw_mag(1:nws))
534 call phys_face_to_center(ixm^ll,ps(igrid))
535 end do
536 !$OMP END PARALLEL DO
537 else
538 !$OMP PARALLEL DO PRIVATE(igrid)
539 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
540 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
541 block=>ps(igrid)
542 call temp%sts_set_sources(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x,ps1(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
543 ps(igrid)%w(ixm^t,temp%startVar:temp%endVar)=ps(igrid)%w(ixm^t,temp%startVar:temp%endVar)+&
544 dtj*ps1(igrid)%w(ixm^t,temp%startVar:temp%endVar)
545 end do
546 !$OMP END PARALLEL DO
547 end if
548 !fix conserve the fluxes set in the STS method
549 if(fix_conserve_at_step) then
550 call recvflux(1,ndim)
551 call sendflux(1,ndim)
552 call fix_conserve(ps,1,ndim,temp%startVar,temp%nflux)
553 if(stagger_grid) then
554 call fix_edges(ps,1,ndim)
555 ! fill the cell-center values from the updated staggered variables
556 !$OMP PARALLEL DO PRIVATE(igrid)
557 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
558 call phys_face_to_center(ixg^ll,ps(igrid))
559 end do
560 !$OMP END PARALLEL DO
561 end if
562 end if
563 if(associated(temp%sts_handle_errors)) then
564 !$OMP PARALLEL DO PRIVATE(igrid)
565 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
566 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
567 block=>ps(igrid)
568 call temp%sts_handle_errors(ps(igrid)%w,ps(igrid)%x,ixg^ll,ixm^ll,j)
569 end do
570 !$OMP END PARALLEL DO
571 end if
572
573 if(temp%nflux>temp%nwbc.and.temp%s==j) then
574 ! include the changed no-need-ghost-update variables in the last getbc
575 type_send_srl=>temp%type_send_srl_sts_2
576 type_recv_srl=>temp%type_recv_srl_sts_2
577 type_send_r=>temp%type_send_r_sts_2
578 type_recv_r=>temp%type_recv_r_sts_2
579 type_send_p=>temp%type_send_p_sts_2
580 type_recv_p=>temp%type_recv_p_sts_2
581 call getbc(global_time,0.d0,ps,temp%startVar,temp%nflux)
582 else
583 call getbc(global_time,0.d0,ps,temp%startwbc,temp%nwbc)
584 end if
585 end do
586
587 if(associated(temp%sts_after_last_cycle)) then
588 do iigrid=1,igridstail; igrid=igrids(iigrid);
589 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
590 block=>ps(igrid)
591 call temp%sts_after_last_cycle(ixg^ll,ixg^ll,ps(igrid)%w,ps(igrid)%x)
592 end do
593 prolongprimitive = prolong_flag
594 coarsenprimitive = coarsen_flag
595 end if
596 deallocate(bj)
597
598 if(.not.temp%evolve_magnetic_field) then
599 ! restore stagger_grid value
600 stagger_grid=stagger_flag
601 end if
602
603 temp=>temp%next
604 end do
605
606 if(associated(head_sts_terms)) then
607 ! point bc mpi data type back to full type for (M)HD
614 end if
615
616 bcphys=.true.
617 end subroutine sts_add_source2
618
619 !> Iterates all the terms implemented with STS and adds the sources
620 !> STS method 1 implementation
621 subroutine sts_add_source1(my_dt)
622 ! Meyer 2012 MNRAS 422,2102
626 use mod_physics
628
629 double precision, intent(in) :: my_dt
630 double precision :: dtj
631 double precision :: omega1,cmu,cmut,cnu,cnut,one_mu_nu
632 double precision, allocatable :: bj(:)
633 integer:: iigrid, igrid, j, ixc^l, ixgext^l
634 logical :: evenstep, stagger_flag=.false., prolong_flag=.false., coarsen_flag=.false., total_energy_flag=.true.
635 type(sts_term), pointer :: temp
636 type(state), dimension(:), pointer :: tmpps1, tmpps2
637
638 ! do not fill physical boundary conditions
639 bcphys=.false.
640
641 fix_conserve_at_step = time_advance .and. levmax>levmin
642
643 temp => head_sts_terms
644 do while(associated(temp))
645
646 if(.not.temp%evolve_magnetic_field) then
647 ! not do fix conserve and getbc for staggered values
648 stagger_flag=stagger_grid
649 stagger_grid=.false.
650 else if(stagger_grid) then
651 ixcmax^d=ixmhi^d;
652 ixcmin^d=ixmlo^d-1;
653 end if
654
655 call init_comm_fix_conserve(1,ndim,temp%nflux)
656
657 if(associated(temp%sts_before_first_cycle)) then
658 prolong_flag = prolongprimitive
659 coarsen_flag = coarsenprimitive
660 prolongprimitive = .false.
661 coarsenprimitive = .false.
662 total_energy_flag=phys_total_energy
663 phys_total_energy=.false.
664 !$OMP PARALLEL DO PRIVATE(igrid)
665 do iigrid=1,igridstail; igrid=igrids(iigrid);
666 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
667 block=>ps(igrid)
668 call temp%sts_before_first_cycle(ixg^ll,ixg^ll,ps(igrid)%w,ps(igrid)%x)
669 if(.not. allocated(ps2(igrid)%w)) allocate(ps2(igrid)%w(ixg^t,1:nw))
670 if(.not. allocated(ps3(igrid)%w)) allocate(ps3(igrid)%w(ixg^t,1:nw))
671 if(.not. allocated(ps4(igrid)%w)) allocate(ps4(igrid)%w(ixg^t,1:nw))
672 ps1(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
673 ps2(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
674 end do
675 !$OMP END PARALLEL DO
676 else
677 if(stagger_grid) then
678 ixgext^l=ixg^ll^ladd1;
679 !$OMP PARALLEL DO PRIVATE(igrid)
680 do iigrid=1,igridstail; igrid=igrids(iigrid);
681 if(.not. allocated(ps2(igrid)%w)) then
682 call alloc_state(igrid, ps2(igrid), ixg^ll, ixgext^l, .false.)
683 end if
684 if(.not. allocated(ps3(igrid)%w)) allocate(ps3(igrid)%w(ixg^t,1:nw))
685 if(.not. allocated(ps4(igrid)%w)) allocate(ps4(igrid)%w(ixg^t,1:nw))
686 ps1(igrid)%w=ps(igrid)%w
687 ps2(igrid)%w=ps(igrid)%w
688 ps1(igrid)%ws=ps(igrid)%ws
689 ps2(igrid)%ws=ps(igrid)%ws
690 end do
691 !$OMP END PARALLEL DO
692 else
693 !$OMP PARALLEL DO PRIVATE(igrid)
694 do iigrid=1,igridstail; igrid=igrids(iigrid);
695 if(.not. allocated(ps2(igrid)%w)) allocate(ps2(igrid)%w(ixg^t,1:nw))
696 if(.not. allocated(ps3(igrid)%w)) allocate(ps3(igrid)%w(ixg^t,1:nw))
697 if(.not. allocated(ps4(igrid)%w)) allocate(ps4(igrid)%w(ixg^t,1:nw))
698 ps1(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
699 ps2(igrid)%w(ixg^t,1:nw)=ps(igrid)%w(ixg^t,1:nw)
700 end do
701 !$OMP END PARALLEL DO
702 end if
703 end if
704
705 allocate(bj(0:temp%s))
706 bj(0)=1.d0/3.d0
707 bj(1)=bj(0)
708 if(temp%s>1) then
709 omega1=4.d0/dble(temp%s**2+temp%s-2)
710 cmut=omega1/3.d0
711 else
712 omega1=0.d0
713 cmut=1.d0
714 end if
715
716 type_send_srl=>temp%type_send_srl_sts_1
717 type_recv_srl=>temp%type_recv_srl_sts_1
718 type_send_r=>temp%type_send_r_sts_1
719 type_recv_r=>temp%type_recv_r_sts_1
720 type_send_p=>temp%type_send_p_sts_1
721 type_recv_p=>temp%type_recv_p_sts_1
722
723 if(.not. temp%types_initialized) then
724 call create_bc_mpi_datatype(temp%startwbc,temp%nwbc)
725 if(temp%nflux>temp%nwbc) then
726 ! prepare types for the changed no-need-ghost-update variables in the last getbc
727 type_send_srl=>temp%type_send_srl_sts_2
728 type_recv_srl=>temp%type_recv_srl_sts_2
729 type_send_r=>temp%type_send_r_sts_2
730 type_recv_r=>temp%type_recv_r_sts_2
731 type_send_p=>temp%type_send_p_sts_2
732 type_recv_p=>temp%type_recv_p_sts_2
733 call create_bc_mpi_datatype(temp%startVar,temp%nflux)
734 type_send_srl=>temp%type_send_srl_sts_1
735 type_recv_srl=>temp%type_recv_srl_sts_1
736 type_send_r=>temp%type_send_r_sts_1
737 type_recv_r=>temp%type_recv_r_sts_1
738 type_send_p=>temp%type_send_p_sts_1
739 type_recv_p=>temp%type_recv_p_sts_1
740 end if
741 temp%types_initialized = .true.
742 end if
743 dtj = cmut*my_dt
744 if(stagger_grid) then
745 !$OMP PARALLEL DO PRIVATE(igrid)
746 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
747 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
748 block=>ps(igrid)
749 ps4(igrid)%w=zero
750 call temp%sts_set_sources(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
751 !!!eq solved: dU/dt = S, ps3 is stored S^n
752 ps3(igrid)%w(ixc^s,temp%startVar:temp%endVar) = my_dt * ps4(igrid)%w(ixc^s,temp%startVar:temp%endVar)
753 if(temp%nflux>ndir) then
754 ps1(igrid)%w(ixm^t,temp%startVar) = ps1(igrid)%w(ixm^t,temp%startVar) + cmut * ps3(igrid)%w(ixm^t,temp%startVar)
755 end if
756 ps1(igrid)%ws(ixc^s,1:nws) = ps1(igrid)%ws(ixc^s,1:nws) + cmut * ps3(igrid)%w(ixc^s,iw_mag(1:nws))
757 call phys_face_to_center(ixm^ll,ps1(igrid))
758 end do
759 !$OMP END PARALLEL DO
760 else
761 !$OMP PARALLEL DO PRIVATE(igrid)
762 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
763 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
764 block=>ps(igrid)
765 call temp%sts_set_sources(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
766 !!!eq solved: dU/dt = S, ps3 is stored S^n
767 ps3(igrid)%w(ixm^t,temp%startVar:temp%endVar) = my_dt * ps4(igrid)%w(ixm^t,temp%startVar:temp%endVar)
768 ps1(igrid)%w(ixm^t,temp%startVar:temp%endVar) = ps1(igrid)%w(ixm^t,temp%startVar:temp%endVar) + &
769 cmut * ps3(igrid)%w(ixm^t,temp%startVar:temp%endVar)
770 end do
771 !$OMP END PARALLEL DO
772 end if
773 if(fix_conserve_at_step) then
774 call recvflux(1,ndim)
775 call sendflux(1,ndim)
776 call fix_conserve(ps1,1,ndim,temp%startVar,temp%nflux)
777 if(stagger_grid) then
778 call fix_edges(ps1,1,ndim)
779 ! fill the cell-center values from the updated staggered variables
780 !$OMP PARALLEL DO PRIVATE(igrid)
781 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
782 call phys_face_to_center(ixg^ll,ps1(igrid))
783 end do
784 !$OMP END PARALLEL DO
785 end if
786 end if
787 ! fix conservation of AMR grid by replacing flux from finer neighbors
788 if(associated(temp%sts_handle_errors)) then
789 !$OMP PARALLEL DO PRIVATE(igrid)
790 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
791 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
792 block=>ps(igrid)
793 call temp%sts_handle_errors(ps1(igrid)%w,ps1(igrid)%x,ixg^ll,ixm^ll,1)
794 end do
795 !$OMP END PARALLEL DO
796 end if
797 if(temp%nflux>temp%nwbc.and.temp%s==1) then
798 ! include the changed no-need-ghost-update variables in the last getbc
799 type_send_srl=>temp%type_send_srl_sts_2
800 type_recv_srl=>temp%type_recv_srl_sts_2
801 type_send_r=>temp%type_send_r_sts_2
802 type_recv_r=>temp%type_recv_r_sts_2
803 type_send_p=>temp%type_send_p_sts_2
804 type_recv_p=>temp%type_recv_p_sts_2
805 call getbc(global_time,0.d0,ps1,temp%startVar,temp%nflux)
806 else
807 call getbc(global_time,0.d0,ps1,temp%startwbc,temp%nwbc)
808 end if
809 !!first step end
810
811 evenstep=.true.
812
813 tmpps2=>ps1
814
815 do j=2,temp%s
816 bj(j)=dble(j**2+j-2)/dble(2*j*(j+1))
817 cmu=dble(2*j-1)/dble(j)*bj(j)/bj(j-1)
818 cmut=omega1*cmu
819 cnu=dble(1-j)/dble(j)*bj(j)/bj(j-2)
820 cnut=(bj(j-1)-1.d0)*cmut
821 one_mu_nu=1.d0-cmu-cnu
822 if(evenstep) then
823 tmpps1=>ps1
824 tmpps2=>ps2
825 else
826 tmpps1=>ps2
827 tmpps2=>ps1
828 end if
829
830 dtj = cmut*my_dt
831 if(stagger_grid) then
832 !$OMP PARALLEL DO PRIVATE(igrid)
833 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
834 ! maybe the following global variables are needed in set_sources
835 ! next few lines ensure correct usage of routines like divvector etc
836 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
837 block=>ps(igrid)
838 ! end maybe the following global variables are needed in set_sources
839 call temp%sts_set_sources(ixg^ll,ixm^ll,tmpps1(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
840 if(temp%nflux>ndir) then
841 tmpps2(igrid)%w(ixm^t,temp%startVar)=cmu*tmpps1(igrid)%w(ixm^t,temp%startVar)+&
842 cnu*tmpps2(igrid)%w(ixm^t,temp%startVar)+one_mu_nu*ps(igrid)%w(ixm^t,temp%startVar)+&
843 dtj*ps4(igrid)%w(ixm^t,temp%startVar)+cnut*ps3(igrid)%w(ixm^t,temp%startVar)
844 end if
845 tmpps2(igrid)%ws(ixc^s,1:nws)=cmu*tmpps1(igrid)%ws(ixc^s,1:nws)+&
846 cnu*tmpps2(igrid)%ws(ixc^s,1:nws)+one_mu_nu*ps(igrid)%ws(ixc^s,1:nws)+&
847 dtj*ps4(igrid)%w(ixc^s,iw_mag(1:nws))+cnut*ps3(igrid)%w(ixc^s,iw_mag(1:nws))
848 call phys_face_to_center(ixm^ll,tmpps2(igrid))
849 end do
850 !$OMP END PARALLEL DO
851 else
852 !$OMP PARALLEL DO PRIVATE(igrid)
853 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
854 ! maybe the following global variables are needed in set_sources
855 ! next few lines ensure correct usage of routines like divvector etc
856 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
857 block=>ps(igrid)
858 ! end maybe the following global variables are needed in set_sources
859 call temp%sts_set_sources(ixg^ll,ixm^ll,tmpps1(igrid)%w,ps(igrid)%x,ps4(igrid)%w,fix_conserve_at_step,dtj,igrid,temp%nflux)
860 tmpps2(igrid)%w(ixm^t,temp%startVar:temp%endVar)=cmu*tmpps1(igrid)%w(ixm^t,temp%startVar:temp%endVar)+&
861 cnu*tmpps2(igrid)%w(ixm^t,temp%startVar:temp%endVar)+one_mu_nu*ps(igrid)%w(ixm^t,temp%startVar:temp%endVar)+&
862 dtj*ps4(igrid)%w(ixm^t,temp%startVar:temp%endVar)+cnut*ps3(igrid)%w(ixm^t,temp%startVar:temp%endVar)
863 end do
864 !$OMP END PARALLEL DO
865 end if
866 if(fix_conserve_at_step) then
867 call recvflux(1,ndim)
868 call sendflux(1,ndim)
869 call fix_conserve(tmpps2,1,ndim,temp%startVar,temp%nflux)
870 if(stagger_grid) then
871 call fix_edges(tmpps2,1,ndim)
872 ! fill the cell-center values from the updated staggered variables
873 !$OMP PARALLEL DO PRIVATE(igrid)
874 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
875 call phys_face_to_center(ixg^ll,tmpps2(igrid))
876 end do
877 !$OMP END PARALLEL DO
878 end if
879 end if
880 if(associated(temp%sts_handle_errors)) then
881 !$OMP PARALLEL DO PRIVATE(igrid)
882 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
883 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
884 block=>ps(igrid)
885 call temp%sts_handle_errors(tmpps2(igrid)%w,ps(igrid)%x,ixg^ll,ixm^ll,j)
886 end do
887 !$OMP END PARALLEL DO
888 end if
889 if(temp%nflux>temp%nwbc.and.temp%s==j) then
890 ! include the changed no-need-ghost-update variables in the last getbc
891 type_send_srl=>temp%type_send_srl_sts_2
892 type_recv_srl=>temp%type_recv_srl_sts_2
893 type_send_r=>temp%type_send_r_sts_2
894 type_recv_r=>temp%type_recv_r_sts_2
895 type_send_p=>temp%type_send_p_sts_2
896 type_recv_p=>temp%type_recv_p_sts_2
897 call getbc(global_time,0.d0,tmpps2,temp%startVar,temp%nflux)
898 else
899 call getbc(global_time,0.d0,tmpps2,temp%startwbc,temp%nwbc)
900 end if
901 evenstep=.not.evenstep
902 end do
903
904 if(associated(temp%sts_after_last_cycle)) then
905 !$OMP PARALLEL DO PRIVATE(igrid)
906 do iigrid=1,igridstail; igrid=igrids(iigrid);
907 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
908 block=>ps(igrid)
909 ps(igrid)%w(ixg^t,temp%startVar:temp%endVar)=tmpps2(igrid)%w(ixg^t,temp%startVar:temp%endVar)
910 call temp%sts_after_last_cycle(ixg^ll,ixg^ll,ps(igrid)%w,ps(igrid)%x)
911 end do
912 !$OMP END PARALLEL DO
913 phys_total_energy=total_energy_flag
914 prolongprimitive = prolong_flag
915 coarsenprimitive = coarsen_flag
916 else
917 if(stagger_grid) then
918 !$OMP PARALLEL DO PRIVATE(igrid)
919 do iigrid=1,igridstail; igrid=igrids(iigrid);
920 ps(igrid)%w(ixg^t,temp%startVar:temp%endVar)=tmpps2(igrid)%w(ixg^t,temp%startVar:temp%endVar)
921 ps(igrid)%ws=tmpps2(igrid)%ws
922 end do
923 !$OMP END PARALLEL DO
924 else
925 !$OMP PARALLEL DO PRIVATE(igrid)
926 do iigrid=1,igridstail; igrid=igrids(iigrid);
927 ps(igrid)%w(ixg^t,temp%startVar:temp%endVar)=tmpps2(igrid)%w(ixg^t,temp%startVar:temp%endVar)
928 end do
929 !$OMP END PARALLEL DO
930 end if
931 end if
932
933 deallocate(bj)
934
935 if(.not.temp%evolve_magnetic_field) then
936 ! restore stagger_grid value
937 stagger_grid=stagger_flag
938 end if
939
940 temp=>temp%next
941 end do
942
943 if(associated(head_sts_terms)) then
944 ! point bc mpi data type back to full type for (M)HD
951 end if
952
953 bcphys=.true.
954 end subroutine sts_add_source1
955
956end module mod_supertimestepping
subroutine, public alloc_state(igrid, s, ixgl, ixgextl, alloc_once_for_ps)
allocate memory to physical state of igrid node
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter dpi
Pi.
Module for flux conservation near refinement boundaries.
subroutine, public init_comm_fix_conserve(idimlim, nwfluxin)
subroutine, public fix_edges(psuse, idimlim)
subroutine, public recvflux(idimlim)
subroutine, public sendflux(idimlim)
subroutine, public fix_conserve(psb, idimlim, nw0, nwfluxin)
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
update ghost cells of all blocks including physical boundaries
integer, dimension( :^d &), pointer type_recv_r
integer, dimension(^nd, 0:3) l
subroutine getbc(time, qdt, psb, nwstart, nwbc)
do update ghost cells of all blocks including physical boundaries
integer, dimension(0:3^d &), target type_send_p_f
integer, dimension( :^d &), pointer type_send_p
integer, dimension( :^d &), pointer type_send_srl
subroutine create_bc_mpi_datatype(nwstart, nwbc)
integer, dimension(0:3^d &), target type_recv_r_f
integer, dimension(-1:1^d &), target type_recv_srl_f
integer, dimension(-1:1^d &), target type_send_r_f
integer, dimension(-1:1^d &), target type_send_srl_f
integer, dimension( :^d &), pointer type_send_r
integer, dimension( :^d &), pointer type_recv_p
integer, dimension(0:3^d &), target type_recv_p_f
integer, dimension( :^d &), pointer type_recv_srl
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
integer, parameter unitpar
file handle for IO
double precision global_time
The global simulation time.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer icomm
The MPI communicator.
logical coarsenprimitive
coarsen primitive variables in level-jump ghost cells
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
logical time_advance
do time evolving
logical prolongprimitive
prolongate primitive variables in level-jump ghost cells
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
logical phys_total_energy
Solve total energy equation or not.
Definition mod_physics.t:40
procedure(sub_face_to_center), pointer phys_face_to_center
Definition mod_physics.t:79
Module for handling problematic values in simulations, such as negative pressures.
Generic supertimestepping method which can be used for multiple source terms in the governing equatio...
pure logical function, public is_sts_initialized()
integer, parameter, public sourcetype_sts_prior
logical function, public set_dt_sts_ncycles(my_dt)
This sets the explicit dt and calculates the number of cycles for each of the terms implemented with ...
subroutine, public add_sts_method(sts_getdt, sts_set_sources, startvar, nflux, startwbc, nwbc, evolve_b)
subroutine which added programatically a term to be calculated using STS Params: sts_getdt function c...
type(sts_term), pointer head_sts_terms
subroutine, public set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
Set the hooks called before the first cycle and after the last cycle in the STS update This method sh...
integer, parameter, public sourcetype_sts_split
subroutine, public set_error_handling_to_head(sts_error_handling)
Set the hook of error handling in the STS update. This method is called before updating the BC....
subroutine, public sts_init()
Initialize sts module.
integer, parameter, public sourcetype_sts_after
procedure(subr3), pointer, public sts_add_source