MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_twofl_phys.t
Go to the documentation of this file.
1!> Magneto-hydrodynamics module
3
4#include "amrvac.h"
5
7 use mod_global_parameters, only: std_len
12 use mod_comm_lib, only: mpistop
13
14 implicit none
15 private
16 !! E_c = E_kin + E_mag + E_int
17 !! E_n = E_kin + E_int
18 integer, public, parameter :: eq_energy_tot=2
19 !! E_c = E_int
20 !! E_n = E_int
21 integer, public, parameter :: eq_energy_int=1
22 !! E_n, E_c are calculated from density as c_adiab rho^gamma
23 !! No energy equation => no variable assigned for it
24 integer, public, parameter :: eq_energy_none=0
25 !! E_c = E_kin + E_int
26 !! E_n = E_kin + E_int
27 integer, public, parameter :: eq_energy_ki=3
28
29 integer, public, protected :: twofl_eq_energy = eq_energy_tot
30
31 !> Whether hyperdiffusivity is used
32 logical, public, protected :: twofl_hyperdiffusivity = .false.
33 logical, public, protected :: twofl_dump_hyperdiffusivity_coef = .false.
34 double precision, public, protected, allocatable :: c_shk(:)
35 double precision, public, protected, allocatable :: c_hyp(:)
36
37 !> Whether thermal conduction is used
38 logical, public, protected :: twofl_thermal_conduction_c = .false.
39 !> type of TC used: 1: adapted module (mhd implementation), 2: adapted module (hd implementation)
40 integer, parameter, private :: mhd_tc =1
41 integer, parameter, private :: hd_tc =2
42 integer, protected :: use_twofl_tc_c = mhd_tc
43
44 !> Whether radiative cooling is added
45 logical, public, protected :: twofl_radiative_cooling_c = .false.
46 type(rc_fluid), public, allocatable :: rc_fl_c
47
48 !> Whether viscosity is added
49 logical, public, protected :: twofl_viscosity = .false.
50
51 !> Whether gravity is added: common flag for charges and neutrals
52 logical, public, protected :: twofl_gravity = .false.
53
54 !> whether dump full variables (when splitting is used) in a separate dat file
55 logical, public, protected :: twofl_dump_full_vars = .false.
56
57 !> Whether Hall-MHD is used
58 logical, public, protected :: twofl_hall = .false.
59
60 type(tc_fluid), public, allocatable :: tc_fl_c
61 type(te_fluid), public, allocatable :: te_fl_c
62
63 type(tc_fluid), allocatable :: tc_fl_n
64 logical, public, protected :: twofl_thermal_conduction_n = .false.
65 logical, public, protected :: twofl_radiative_cooling_n = .false.
66 type(rc_fluid), allocatable :: rc_fl_n
67
68 !> Whether TRAC method is used
69 logical, public, protected :: twofl_trac = .false.
70
71 !> Whether GLM-MHD is used
72 logical, public, protected :: twofl_glm = .false.
73
74 !> Which TRAC method is used
75 integer, public, protected :: twofl_trac_type=1
76
77 !> Height of the mask used in the TRAC method
78 double precision, public, protected :: twofl_trac_mask = 0.d0
79
80 !> Whether divB cleaning sources are added splitting from fluid solver
81 logical, public, protected :: source_split_divb = .false.
82
83 !> GLM-MHD parameter: ratio of the diffusive and advective time scales for div b
84 !> taking values within [0, 1]
85 double precision, public :: twofl_glm_alpha = 0.5d0
86
87 !> MHD fourth order
88 logical, public, protected :: twofl_4th_order = .false.
89
90 !> Index of the density (in the w array)
91 integer, public :: rho_c_
92
93 !> Indices of the momentum density
94 integer, allocatable, public :: mom_c(:)
95 !> Indices of the momentum density for the form of better vectorization
96 integer, public, protected :: ^c&m^C_
97
98 !> Index of the energy density (-1 if not present)
99 integer, public :: e_c_=-1
100 !> Indices of the magnetic field for the form of better vectorization
101 integer, public, protected :: ^c&b^C_
102
103 !> Index of the cutoff temperature for the TRAC method
104 integer, public :: tcoff_c_
105 integer, public :: tweight_c_
106
107 !> Indices of the GLM psi
108 integer, public, protected :: psi_
109
110 !> equi vars flags
111 logical, public :: has_equi_rho_c0 = .false.
112 logical, public :: has_equi_pe_c0 = .false.
113
114 !> equi vars indices in the state%equi_vars array
115 integer, public :: equi_rho_c0_ = -1
116 integer, public :: equi_pe_c0_ = -1
117 logical, public :: twofl_equi_thermal_c = .false.
118
119 logical, public :: twofl_equi_thermal = .false.
120 !neutrals:
121
122 integer, public :: rho_n_
123 integer, allocatable, public :: mom_n(:)
124 integer, public :: e_n_
125 integer, public :: tcoff_n_
126 integer, public :: tweight_n_
127 logical, public :: has_equi_rho_n0 = .false.
128 logical, public :: has_equi_pe_n0 = .false.
129 integer, public :: equi_rho_n0_ = -1
130 integer, public :: equi_pe_n0_ = -1
131
132 ! related to collisions:
133 !> collisional alpha
134 double precision, public :: twofl_alpha_coll = 0d0
135 logical, public :: twofl_alpha_coll_constant = .true.
136 !> whether include thermal exchange collisional terms
137 logical, public :: twofl_coll_inc_te = .true.
138 !> whether include ionization/recombination inelastic collisional terms
139 logical, public :: twofl_coll_inc_ionrec = .false.
140 logical, public :: twofl_equi_thermal_n = .false.
141 double precision, public :: dtcollpar = -1d0 !negative value does not impose restriction on the timestep
142 !> whether dump collisional terms in a separte dat file
143 logical, public, protected :: twofl_dump_coll_terms = .false.
144
145 ! TODO Helium abundance not used, radiative cooling init uses it
146 ! not in parameters list anymore
147 double precision, public, protected :: he_abundance = 0d0
148 ! two fluid is only H plasma
149 double precision, public, protected :: rc = 2d0
150 double precision, public, protected :: rn = 1d0
151
152 !> The adiabatic index
153 double precision, public :: twofl_gamma = 5.d0/3.0d0
154
155 !> The adiabatic constant
156 double precision, public :: twofl_adiab = 1.0d0
157
158 !> The MHD resistivity
159 double precision, public :: twofl_eta = 0.0d0
160
161 !> The MHD hyper-resistivity
162 double precision, public :: twofl_eta_hyper = 0.0d0
163
164 !> The MHD Hall coefficient
165 double precision, public :: twofl_etah = 0.0d0
166
167 !> The small_est allowed energy
168 double precision, protected :: small_e
169
170 !> Method type to clean divergence of B
171 character(len=std_len), public, protected :: typedivbfix = 'linde'
172
173 !> Method type of constrained transport
174 character(len=std_len), public, protected :: type_ct = 'uct_contact'
175
176 !> Whether divB is computed with a fourth order approximation
177 integer, public, protected :: twofl_divb_nth = 1
178
179 !> Method type in a integer for good performance
180 integer :: type_divb
181
182 !> Coefficient of diffusive divB cleaning
183 double precision :: divbdiff = 0.8d0
184
185 !> Update all equations due to divB cleaning
186 character(len=std_len) :: typedivbdiff = 'all'
187
188 !> clean initial divB
189 logical, public :: clean_initial_divb = .false.
190
191 !> Add divB wave in Roe solver
192 logical, public :: divbwave = .true.
193
194 !> To control divB=0 fix for boundary
195 logical, public, protected :: boundary_divbfix(2*^nd)=.true.
196
197 !> To skip * layer of ghost cells during divB=0 fix for boundary
198 integer, public, protected :: boundary_divbfix_skip(2*^nd)=0
199
200 !> B0 field is force-free
201 logical, public, protected :: b0field_forcefree=.true.
202
203 logical :: twofl_cbounds_species = .true.
204
205 !> added from modules: gravity
206 !> source split or not
207 logical :: grav_split= .false.
208
209 !> gamma minus one and its inverse
210 double precision :: gamma_1, inv_gamma_1
211
212 ! DivB cleaning methods
213 integer, parameter :: divb_none = 0
214 integer, parameter :: divb_multigrid = -1
215 integer, parameter :: divb_glm = 1
216 integer, parameter :: divb_powel = 2
217 integer, parameter :: divb_janhunen = 3
218 integer, parameter :: divb_linde = 4
219 integer, parameter :: divb_lindejanhunen = 5
220 integer, parameter :: divb_lindepowel = 6
221 integer, parameter :: divb_lindeglm = 7
222 integer, parameter :: divb_ct = 8
223
224 ! Public methods
225 public :: twofl_phys_init
226 public :: twofl_to_conserved
227 public :: twofl_to_primitive
228 public :: get_divb
229 public :: get_rhoc_tot
230 public :: twofl_get_v_c_idim
231 ! TODO needed for the roe, see if can be used for n
233 public :: get_rhon_tot
234 public :: get_alpha_coll
235 public :: get_gamma_ion_rec
236 public :: twofl_get_v_n_idim
237 public :: get_current
238 public :: twofl_get_pthermal_c
239 public :: twofl_get_pthermal_n
240 public :: twofl_face_to_center
241 public :: get_normalized_divb
243 public :: usr_mask_gamma_ion_rec
244 public :: usr_mask_alpha
245
246 {^nooned
248 }
249
250 abstract interface
251
252 subroutine implicit_mult_factor_subroutine(ixI^L, ixO^L, step_dt, JJ, res)
253 integer, intent(in) :: ixi^l, ixo^l
254 double precision, intent(in) :: step_dt
255 double precision, intent(in) :: jj(ixi^s)
256 double precision, intent(out) :: res(ixi^s)
257
258 end subroutine implicit_mult_factor_subroutine
259
260 subroutine mask_subroutine(ixI^L,ixO^L,w,x,res)
262 integer, intent(in) :: ixi^l, ixo^l
263 double precision, intent(in) :: x(ixi^s,1:ndim)
264 double precision, intent(in) :: w(ixi^s,1:nw)
265 double precision, intent(inout) :: res(ixi^s)
266 end subroutine mask_subroutine
267
268 subroutine mask_subroutine2(ixI^L,ixO^L,w,x,res1, res2)
270 integer, intent(in) :: ixi^l, ixo^l
271 double precision, intent(in) :: x(ixi^s,1:ndim)
272 double precision, intent(in) :: w(ixi^s,1:nw)
273 double precision, intent(inout) :: res1(ixi^s),res2(ixi^s)
274 end subroutine mask_subroutine2
275
276 end interface
277
278 procedure(implicit_mult_factor_subroutine), pointer :: calc_mult_factor => null()
279 integer, protected :: twofl_implicit_calc_mult_method = 1
280 procedure(mask_subroutine), pointer :: usr_mask_alpha => null()
281 procedure(mask_subroutine2), pointer :: usr_mask_gamma_ion_rec => null()
282
283contains
284
285 !> Read this module"s parameters from a file
286 subroutine twofl_read_params(files)
288 character(len=*), intent(in) :: files(:)
289 integer :: n
290
291 namelist /twofl_list/ twofl_eq_energy, twofl_gamma, twofl_adiab,&
295 typedivbdiff, type_ct, divbwave, si_unit, b0field,&
302 twofl_dump_coll_terms,twofl_implicit_calc_mult_method,&
305 twofl_trac, twofl_trac_type, twofl_trac_mask,twofl_cbounds_species
306
307 do n = 1, size(files)
308 open(unitpar, file=trim(files(n)), status="old")
309 read(unitpar, twofl_list, end=111)
310111 close(unitpar)
311 end do
312
313 end subroutine twofl_read_params
314
315 subroutine twofl_init_hyper(files)
318 character(len=*), intent(in) :: files(:)
319 integer :: n
320
321 namelist /hyperdiffusivity_list/ c_shk, c_hyp
322
323 do n = 1, size(files)
324 open(unitpar, file=trim(files(n)), status="old")
325 read(unitpar, hyperdiffusivity_list, end=113)
326113 close(unitpar)
327 end do
328
329 call hyperdiffusivity_init()
330
331 !!DEBUG
332 if(mype .eq. 0) then
333 print*, "Using Hyperdiffusivity"
334 print*, "C_SHK ", c_shk(:)
335 print*, "C_HYP ", c_hyp(:)
336 endif
337
338 end subroutine twofl_init_hyper
339
340 !> Write this module's parameters to a snapsoht
341 subroutine twofl_write_info(fh)
343 integer, intent(in) :: fh
344 integer, parameter :: n_par = 1
345 double precision :: values(n_par)
346 character(len=name_len) :: names(n_par)
347 integer, dimension(MPI_STATUS_SIZE) :: st
348 integer :: er
349
350 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
351
352 names(1) = "gamma"
353 values(1) = twofl_gamma
354 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
355 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
356 end subroutine twofl_write_info
357
358 subroutine twofl_phys_init()
363 !use mod_gravity, only: gravity_init
366 {^nooned
368 }
369 integer :: itr, idir
370
371 call twofl_read_params(par_files)
372 physics_type = "twofl"
373 if (twofl_cbounds_species) then
374 number_species = 2
375 endif
376 phys_energy=.true.
377 !> Solve total energy equation or not
378 ! for the two fluid the true value means
379 ! E_charges = E_mag + E_kin_charges + E_int_charges
380 ! E_neutrals = E_kin_neutrals + E_int_neutrals
381 phys_total_energy=.false.
382
383 !> Solve internal energy instead of total energy
384 ! for the two fluid the true value means
385 ! E_charges = E_int_charges
386 ! E_neutrals = E_int_neutrals
387 phys_internal_e=.false.
388
389 ! For the two fluid phys_energy=.true. and phys_internal_e=.false. and phys_total_energy = .false. means
390 ! E_charges = E_kin_charges + E_int_charges
391 ! E_neutrals = E_kin_neutrals + E_int_neutrals
392 phys_gamma = twofl_gamma
393
395 phys_internal_e = .true.
396 elseif(twofl_eq_energy == eq_energy_tot) then
397 phys_total_energy = .true.
398 elseif(twofl_eq_energy == eq_energy_none) then
399 phys_energy = .false.
400 endif
401
404
405 if(.not. phys_energy) then
408 if(mype==0) write(*,*) 'WARNING: set twofl_thermal_conduction_n=F when twofl_energy=F'
409 end if
412 if(mype==0) write(*,*) 'WARNING: set twofl_radiative_cooling_n=F when twofl_energy=F'
413 end if
416 if(mype==0) write(*,*) 'WARNING: set twofl_thermal_conduction_c=F when twofl_energy=F'
417 end if
420 if(mype==0) write(*,*) 'WARNING: set twofl_radiative_cooling_c=F when twofl_energy=F'
421 end if
422 if(twofl_trac) then
423 twofl_trac=.false.
424 if(mype==0) write(*,*) 'WARNING: set twofl_trac=F when twofl_energy=F'
425 end if
426 end if
427 {^ifoned
428 if(twofl_trac .and. twofl_trac_type .gt. 1) then
430 if(mype==0) write(*,*) 'WARNING: set twofl_trac_type=1 for 1D simulation'
431 end if
432 }
433 if(twofl_trac .and. twofl_trac_type .le. 3) then
434 twofl_trac_mask=bigdouble
435 if(mype==0) write(*,*) 'WARNING: set twofl_trac_mask==bigdouble for global TRAC method'
436 end if
438
439 ! set default gamma for polytropic/isothermal process
440 if(ndim==1) typedivbfix='none'
441 select case (typedivbfix)
442 case ('none')
443 type_divb = divb_none
444 {^nooned
445 case ('multigrid')
446 type_divb = divb_multigrid
447 use_multigrid = .true.
448 mg%operator_type = mg_laplacian
449 phys_global_source_after => twofl_clean_divb_multigrid
450 }
451 case ('glm')
452 twofl_glm = .true.
453 need_global_cmax = .true.
454 type_divb = divb_glm
455 case ('powel', 'powell')
456 type_divb = divb_powel
457 case ('janhunen')
458 type_divb = divb_janhunen
459 case ('linde')
460 type_divb = divb_linde
461 case ('lindejanhunen')
462 type_divb = divb_lindejanhunen
463 case ('lindepowel')
464 type_divb = divb_lindepowel
465 case ('lindeglm')
466 twofl_glm = .true.
467 need_global_cmax = .true.
468 type_divb = divb_lindeglm
469 case ('ct')
470 type_divb = divb_ct
471 stagger_grid = .true.
472 case default
473 call mpistop('Unknown divB fix')
474 end select
475
476 select case(type_ct)
477 case('average')
478 transverse_ghost_cells = 1
479 case('uct_contact')
480 transverse_ghost_cells = 1
481 case('uct_hll')
482 transverse_ghost_cells = 2
483 case default
484 call mpistop('choose average, uct_contact,or uct_hll for type_ct!')
485 end select
486
487 allocate(start_indices(number_species))
488 allocate(stop_indices(number_species))
489 start_indices(1)=1
490 !allocate charges first and the same order as in mhd module
491 rho_c_ = var_set_fluxvar("rho_c", "rho_c")
492 !set variables from mod_variables to point to charges vars
493 iw_rho = rho_c_
494
495 allocate(mom_c(ndir))
496 do idir=1,ndir
497 mom_c(idir) = var_set_fluxvar("m_c","v_c",idir)
498 enddo
499 m^c_=mom_c(^c);
500
501 allocate(iw_mom(ndir))
502 iw_mom(1:ndir) = mom_c(1:ndir)
503
504 ! Set index of energy variable
505 if (phys_energy) then
506 e_c_ = var_set_fluxvar("e_c", "p_c")
507 iw_e = e_c_
508 else
509 e_c_ = -1
510 end if
511
512 ! ambipolar sts assumes mag and energy charges are continuous
513 allocate(mag(ndir))
514 mag(:) = var_set_bfield(ndir)
515 b^c_=mag(^c);
516
517 if (twofl_glm) then
518 psi_ = var_set_fluxvar('psi', 'psi', need_bc=.false.)
519 else
520 psi_ = -1
521 end if
522
523 ! set cutoff temperature when using the TRAC method, as well as an auxiliary weight
524 tweight_c_ = -1
525 if(twofl_trac) then
526 tcoff_c_ = var_set_wextra()
527 iw_tcoff = tcoff_c_
528 if(twofl_trac_type > 2) then
529 tweight_c_ = var_set_wextra()
530 endif
531 else
532 tcoff_c_ = -1
533 end if
534
535 !now allocate neutrals
536
537 ! TODO so far number_species is only used to treat them differently
538 ! in the solvers (different cbounds)
539 if (twofl_cbounds_species) then
540 stop_indices(1)=nwflux
541 start_indices(2)=nwflux+1
542 endif
543
544 ! Determine flux variables
545 rho_n_ = var_set_fluxvar("rho_n", "rho_n")
546 allocate(mom_n(ndir))
547 do idir=1,ndir
548 mom_n(idir) = var_set_fluxvar("m_n","v_n",idir)
549 enddo
550 if (phys_energy) then
551 e_n_ = var_set_fluxvar("e_n", "p_n")
552 else
553 e_n_ = -1
554 end if
555
556 tweight_n_ = -1
557 if(twofl_trac) then
558 tcoff_n_ = var_set_wextra()
559 if(twofl_trac_type > 2) then
560 tweight_n_ = var_set_wextra()
561 endif
562 else
563 tcoff_n_ = -1
564 end if
565
566 stop_indices(number_species)=nwflux
567
568 ! set indices of equi vars and update number_equi_vars
570 if(has_equi_rho_n0) then
573 endif
574 if(has_equi_pe_n0) then
577 phys_equi_pe=.true.
578 endif
579 if(has_equi_rho_c0) then
582 iw_equi_rho = equi_rho_c0_
583 endif
584 if(has_equi_pe_c0) then
587 iw_equi_p = equi_pe_c0_
588 phys_equi_pe=.true.
589 endif
590
591 ! set number of variables which need update ghostcells
592 nwgc=nwflux+nwaux
593
594 ! determine number of stagger variables
595 nws=ndim
596
597 ! Check whether custom flux types have been defined
598 if (.not. allocated(flux_type)) then
599 allocate(flux_type(ndir, nw))
600 flux_type = flux_default
601 else if (any(shape(flux_type) /= [ndir, nw])) then
602 call mpistop("phys_check error: flux_type has wrong shape")
603 end if
604
605 if(ndim>1) then
606 if(twofl_glm) then
607 flux_type(:,psi_)=flux_special
608 do idir=1,ndir
609 flux_type(idir,mag(idir))=flux_special
610 end do
611 else
612 do idir=1,ndir
613 flux_type(idir,mag(idir))=flux_tvdlf
614 end do
615 end if
616 end if
617
618 phys_get_dt => twofl_get_dt
619 phys_get_cmax => twofl_get_cmax
620 phys_get_a2max => twofl_get_a2max
621 !phys_get_tcutoff => twofl_get_tcutoff_c
622 if(twofl_cbounds_species) then
623 if (mype .eq. 0) print*, "Using different cbounds for each species nspecies = ", number_species
624 phys_get_cbounds => twofl_get_cbounds_species
625 phys_get_h_speed => twofl_get_h_speed_species
626 else
627 if (mype .eq. 0) print*, "Using same cbounds for all species"
628 phys_get_cbounds => twofl_get_cbounds_one
629 phys_get_h_speed => twofl_get_h_speed_one
630 endif
631 phys_get_flux => twofl_get_flux
632 phys_add_source_geom => twofl_add_source_geom
633 phys_add_source => twofl_add_source
634 phys_to_conserved => twofl_to_conserved
635 phys_to_primitive => twofl_to_primitive
636 phys_check_params => twofl_check_params
637 phys_check_w => twofl_check_w
638 phys_write_info => twofl_write_info
639 phys_handle_small_values => twofl_handle_small_values
640 !set equilibrium variables for the new grid
641 if(number_equi_vars>0) then
642 phys_set_equi_vars => set_equi_vars_grid
643 endif
644 ! convert_type is not known here, so associate the corresp. subroutine in check_params
645 if(type_divb==divb_glm) then
646 phys_modify_wlr => twofl_modify_wlr
647 end if
648
649 ! if using ct stagger grid, boundary divb=0 is not done here
650 if(stagger_grid) then
651 select case(type_ct)
652 case('average')
653 transverse_ghost_cells = 1
654 phys_get_ct_velocity => twofl_get_ct_velocity_average
655 phys_update_faces => twofl_update_faces_average
656 case('uct_contact')
657 transverse_ghost_cells = 1
658 phys_get_ct_velocity => twofl_get_ct_velocity_contact
659 phys_update_faces => twofl_update_faces_contact
660 case('uct_hll')
661 transverse_ghost_cells = 2
662 phys_get_ct_velocity => twofl_get_ct_velocity_hll
663 phys_update_faces => twofl_update_faces_hll
664 case default
665 call mpistop('choose average, uct_contact,or uct_hll for type_ct!')
666 end select
667 phys_face_to_center => twofl_face_to_center
668 phys_modify_wlr => twofl_modify_wlr
669 else if(ndim>1) then
670 phys_boundary_adjust => twofl_boundary_adjust
671 end if
672
673 {^nooned
674 ! clean initial divb
676 }
677
678 ! derive units from basic units
679 call twofl_physical_units()
680
681 if(.not. phys_energy .and. (twofl_thermal_conduction_c&
683 call mpistop("thermal conduction needs twofl_energy=T")
684 end if
685
686 ! initialize thermal conduction module
689 call sts_init()
691 endif
693 allocate(tc_fl_c)
694 if(has_equi_pe_c0 .and. has_equi_rho_c0) then
695 tc_fl_c%get_temperature_from_eint => twofl_get_temperature_from_eint_c_with_equi
696 if(phys_internal_e) then
697 tc_fl_c%get_temperature_from_conserved => twofl_get_temperature_from_eint_c_with_equi
698 else
699 if(twofl_eq_energy == eq_energy_ki) then
700 tc_fl_c%get_temperature_from_conserved => twofl_get_temperature_from_eki_c_with_equi
701 else
702 tc_fl_c%get_temperature_from_conserved => twofl_get_temperature_from_etot_c_with_equi
703 endif
704 endif
705 if(twofl_equi_thermal_c) then
706 tc_fl_c%has_equi = .true.
707 tc_fl_c%get_temperature_equi => twofl_get_temperature_c_equi
708 tc_fl_c%get_rho_equi => twofl_get_rho_c_equi
709 else
710 tc_fl_c%has_equi = .false.
711 endif
712 else
713 if(phys_internal_e) then
714 tc_fl_c%get_temperature_from_conserved => twofl_get_temperature_from_eint_c
715 else
716 if(twofl_eq_energy == eq_energy_ki) then
717 tc_fl_c%get_temperature_from_conserved => twofl_get_temperature_from_eki_c
718 else
719 tc_fl_c%get_temperature_from_conserved => twofl_get_temperature_from_etot_c
720 endif
721 endif
722 tc_fl_c%get_temperature_from_eint => twofl_get_temperature_from_eint_c
723 endif
724 if(use_twofl_tc_c .eq. mhd_tc) then
725 call tc_get_mhd_params(tc_fl_c,tc_c_params_read_mhd)
726 call add_sts_method(twofl_get_tc_dt_mhd_c,twofl_sts_set_source_tc_c_mhd,e_c_,1,e_c_,1,.false.)
727 else if(use_twofl_tc_c .eq. hd_tc) then
728 call tc_get_hd_params(tc_fl_c,tc_c_params_read_hd)
729 call add_sts_method(twofl_get_tc_dt_hd_c,twofl_sts_set_source_tc_c_hd,e_c_,1,e_c_,1,.false.)
730 endif
731 if(.not. phys_internal_e) then
732 call set_conversion_methods_to_head(twofl_e_to_ei_c, twofl_ei_to_e_c)
733 endif
734 call set_error_handling_to_head(twofl_tc_handle_small_e_c)
735 tc_fl_c%get_rho => get_rhoc_tot
736 tc_fl_c%e_ = e_c_
737 tc_fl_c%Tcoff_ = tcoff_c_
738 end if
740 allocate(tc_fl_n)
741 call tc_get_hd_params(tc_fl_n,tc_n_params_read_hd)
742 if(has_equi_pe_n0 .and. has_equi_rho_n0) then
743 tc_fl_n%get_temperature_from_eint => twofl_get_temperature_from_eint_n_with_equi
744 if(twofl_equi_thermal_n) then
745 tc_fl_n%has_equi = .true.
746 tc_fl_n%get_temperature_equi => twofl_get_temperature_n_equi
747 tc_fl_n%get_rho_equi => twofl_get_rho_n_equi
748 else
749 tc_fl_n%has_equi = .false.
750 endif
751 else
752 tc_fl_n%get_temperature_from_eint => twofl_get_temperature_from_eint_n
753 endif
754 if(phys_internal_e) then
755 if(has_equi_pe_n0 .and. has_equi_rho_n0) then
756 tc_fl_n%get_temperature_from_conserved => twofl_get_temperature_from_eint_n_with_equi
757 else
758 tc_fl_n%get_temperature_from_conserved => twofl_get_temperature_from_eint_n
759 endif
760 call add_sts_method(twofl_get_tc_dt_hd_n,twofl_sts_set_source_tc_n_hd,e_n_,1,e_n_,1,.false.)
761 else
762 if(has_equi_pe_n0 .and. has_equi_rho_n0) then
763 tc_fl_n%get_temperature_from_conserved => twofl_get_temperature_from_etot_n_with_equi
764 else
765 tc_fl_n%get_temperature_from_conserved => twofl_get_temperature_from_etot_n
766 endif
767 call add_sts_method(twofl_get_tc_dt_hd_n,twofl_sts_set_source_tc_n_hd,e_n_,1,e_n_,1,.false.)
768 call set_conversion_methods_to_head(twofl_e_to_ei_n, twofl_ei_to_e_n)
769 endif
770 call set_error_handling_to_head(twofl_tc_handle_small_e_n)
771 tc_fl_n%get_rho => get_rhon_tot
772 tc_fl_n%e_ = e_n_
773 tc_fl_n%Tcoff_ = tcoff_n_
774 end if
775
776
777 if(.not. phys_energy .and. (twofl_radiative_cooling_c&
778 .or. twofl_radiative_cooling_n)) then
779 call mpistop("radiative cooling needs twofl_energy=T")
780 end if
781
782 if(twofl_equi_thermal .and. (.not. has_equi_pe_c0 .or. .not. has_equi_pe_n0)) then
783 call mpistop("twofl_equi_thermal=T has_equi_pe_n0 and has _equi_pe_c0=T")
784 endif
785
786 ! initialize thermal conduction module
789 ! Initialize radiative cooling module
792 allocate(rc_fl_c)
793 call radiative_cooling_init(rc_fl_c,rc_params_read_c)
794 rc_fl_c%get_rho => get_rhoc_tot
795 rc_fl_c%get_pthermal => twofl_get_pthermal_c
796 rc_fl_c%get_var_Rfactor => rfactor_c
797 rc_fl_c%e_ = e_c_
798 rc_fl_c%Tcoff_ = tcoff_c_
800 rc_fl_c%has_equi = .true.
801 rc_fl_c%get_rho_equi => twofl_get_rho_c_equi
802 rc_fl_c%get_pthermal_equi => twofl_get_pe_c_equi
803 else
804 rc_fl_c%has_equi = .false.
805 end if
806 end if
807 end if
808 allocate(te_fl_c)
809 te_fl_c%get_rho=> get_rhoc_tot
810 te_fl_c%get_pthermal=> twofl_get_pthermal_c
811 te_fl_c%get_var_Rfactor => rfactor_c
812{^ifthreed
813 phys_te_images => twofl_te_images
814}
815
816 ! Initialize viscosity module
817 !if (twofl_viscosity) call viscosity_init(phys_wider_stencil)
818
819 ! Initialize gravity module
820 if(twofl_gravity) then
821 ! call gravity_init()
822 call grav_params_read(par_files)
823 end if
824
825 ! Initialize particles module
826 ! For Hall, we need one more reconstructed layer since currents are computed
827 ! in getflux: assuming one additional ghost layer (two for FOURTHORDER) was
828 ! added in nghostcells.
829 if (twofl_hall) then
830 if (twofl_4th_order) then
831 phys_wider_stencil = 2
832 else
833 phys_wider_stencil = 1
834 end if
835 end if
836
838 allocate(c_shk(1:nwflux))
839 allocate(c_hyp(1:nwflux))
840 call twofl_init_hyper(par_files)
841 end if
842
843 end subroutine twofl_phys_init
844
845{^ifthreed
846 subroutine twofl_te_images
849
850 select case(convert_type)
851 case('EIvtiCCmpi','EIvtuCCmpi')
853 case('ESvtiCCmpi','ESvtuCCmpi')
855 case('SIvtiCCmpi','SIvtuCCmpi')
857 case('WIvtiCCmpi','WIvtuCCmpi')
859 case default
860 call mpistop("Error in synthesize emission: Unknown convert_type")
861 end select
862 end subroutine twofl_te_images
863}
864
865 ! wrappers for STS functions in thermal_conductivity module
866 ! which take as argument the tc_fluid (defined in the physics module)
867 subroutine twofl_sts_set_source_tc_c_mhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
871 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
872 double precision, intent(in) :: x(ixi^s,1:ndim)
873 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
874 double precision, intent(in) :: my_dt
875 logical, intent(in) :: fix_conserve_at_step
876 call sts_set_source_tc_mhd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl_c)
877 end subroutine twofl_sts_set_source_tc_c_mhd
878
879 subroutine twofl_sts_set_source_tc_c_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
883 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
884 double precision, intent(in) :: x(ixi^s,1:ndim)
885 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
886 double precision, intent(in) :: my_dt
887 logical, intent(in) :: fix_conserve_at_step
888 call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl_c)
889 end subroutine twofl_sts_set_source_tc_c_hd
890
891 function twofl_get_tc_dt_mhd_c(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
892 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
893 !where tc_k_para_i=tc_k_para*B_i**2/B**2
894 !and T=p/rho
897
898 integer, intent(in) :: ixi^l, ixo^l
899 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
900 double precision, intent(in) :: w(ixi^s,1:nw)
901 double precision :: dtnew
902
903 dtnew=get_tc_dt_mhd(w,ixi^l,ixo^l,dx^d,x,tc_fl_c)
904 end function twofl_get_tc_dt_mhd_c
905
906 function twofl_get_tc_dt_hd_c(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
907 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
908 !where tc_k_para_i=tc_k_para*B_i**2/B**2
909 !and T=p/rho
912
913 integer, intent(in) :: ixi^l, ixo^l
914 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
915 double precision, intent(in) :: w(ixi^s,1:nw)
916 double precision :: dtnew
917
918 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl_c)
919 end function twofl_get_tc_dt_hd_c
920
921 subroutine twofl_tc_handle_small_e_c(w, x, ixI^L, ixO^L, step)
924
925 integer, intent(in) :: ixi^l,ixo^l
926 double precision, intent(inout) :: w(ixi^s,1:nw)
927 double precision, intent(in) :: x(ixi^s,1:ndim)
928 integer, intent(in) :: step
929
930 character(len=140) :: error_msg
931
932 write(error_msg,"(a,i3)") "Charges thermal conduction step ", step
933 call twofl_handle_small_ei_c(w,x,ixi^l,ixo^l,e_c_,error_msg)
934 end subroutine twofl_tc_handle_small_e_c
935
936 subroutine twofl_sts_set_source_tc_n_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
940 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
941 double precision, intent(in) :: x(ixi^s,1:ndim)
942 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
943 double precision, intent(in) :: my_dt
944 logical, intent(in) :: fix_conserve_at_step
945 call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl_n)
946 end subroutine twofl_sts_set_source_tc_n_hd
947
948 subroutine twofl_tc_handle_small_e_n(w, x, ixI^L, ixO^L, step)
950
951 integer, intent(in) :: ixi^l,ixo^l
952 double precision, intent(inout) :: w(ixi^s,1:nw)
953 double precision, intent(in) :: x(ixi^s,1:ndim)
954 integer, intent(in) :: step
955
956 character(len=140) :: error_msg
957
958 write(error_msg,"(a,i3)") "Neutral thermal conduction step ", step
959 call twofl_handle_small_ei_n(w,x,ixi^l,ixo^l,e_n_,error_msg)
960 end subroutine twofl_tc_handle_small_e_n
961
962 function twofl_get_tc_dt_hd_n(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
963 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
964 !where tc_k_para_i=tc_k_para*B_i**2/B**2
965 !and T=p/rho
968
969 integer, intent(in) :: ixi^l, ixo^l
970 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
971 double precision, intent(in) :: w(ixi^s,1:nw)
972 double precision :: dtnew
973
974 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl_n)
975 end function twofl_get_tc_dt_hd_n
976
977 subroutine tc_n_params_read_hd(fl)
980 type(tc_fluid), intent(inout) :: fl
981 integer :: n
982 logical :: tc_saturate=.false.
983 double precision :: tc_k_para=0d0
984
985 namelist /tc_n_list/ tc_saturate, tc_k_para
986
987 do n = 1, size(par_files)
988 open(unitpar, file=trim(par_files(n)), status="old")
989 read(unitpar, tc_n_list, end=111)
990111 close(unitpar)
991 end do
992 fl%tc_saturate = tc_saturate
993 fl%tc_k_para = tc_k_para
994
995 end subroutine tc_n_params_read_hd
996
997 subroutine rc_params_read_n(fl)
999 use mod_constants, only: bigdouble
1000 type(rc_fluid), intent(inout) :: fl
1001 integer :: n
1002 ! list parameters
1003 integer :: ncool = 4000
1004 double precision :: cfrac=0.1d0
1005
1006 !> Name of cooling curve
1007 character(len=std_len) :: coolcurve='JCorona'
1008
1009 !> Name of cooling method
1010 character(len=std_len) :: coolmethod='exact'
1011
1012 !> Fixed temperature not lower than tlow
1013 logical :: tfix=.false.
1014
1015 !> Lower limit of temperature
1016 double precision :: tlow=bigdouble
1017
1018 !> Add cooling source in a split way (.true.) or un-split way (.false.)
1019 logical :: rc_split=.false.
1020
1021 namelist /rc_list_n/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
1022
1023 do n = 1, size(par_files)
1024 open(unitpar, file=trim(par_files(n)), status="old")
1025 read(unitpar, rc_list_n, end=111)
1026111 close(unitpar)
1027 end do
1028
1029 fl%ncool=ncool
1030 fl%coolcurve=coolcurve
1031 fl%coolmethod=coolmethod
1032 fl%tlow=tlow
1033 fl%Tfix=tfix
1034 fl%rc_split=rc_split
1035 fl%cfrac=cfrac
1036 end subroutine rc_params_read_n
1037
1038 !end wrappers
1039
1040 ! fill in tc_fluid fields from namelist
1041 subroutine tc_c_params_read_mhd(fl)
1043 type(tc_fluid), intent(inout) :: fl
1044
1045 integer :: n
1046
1047 ! list parameters
1048 logical :: tc_perpendicular=.false.
1049 logical :: tc_saturate=.false.
1050 double precision :: tc_k_para=0d0
1051 double precision :: tc_k_perp=0d0
1052 character(len=std_len) :: tc_slope_limiter="MC"
1053
1054 namelist /tc_c_list/ tc_perpendicular, tc_saturate, tc_slope_limiter, tc_k_para, tc_k_perp
1055 do n = 1, size(par_files)
1056 open(unitpar, file=trim(par_files(n)), status="old")
1057 read(unitpar, tc_c_list, end=111)
1058111 close(unitpar)
1059 end do
1060
1061 fl%tc_perpendicular = tc_perpendicular
1062 fl%tc_saturate = tc_saturate
1063 fl%tc_k_para = tc_k_para
1064 fl%tc_k_perp = tc_k_perp
1065 select case(tc_slope_limiter)
1066 case ('no','none')
1067 fl%tc_slope_limiter = 0
1068 case ('MC')
1069 ! montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
1070 fl%tc_slope_limiter = 1
1071 case('minmod')
1072 ! minmod limiter
1073 fl%tc_slope_limiter = 2
1074 case ('superbee')
1075 ! Roes superbee limiter (eq.3.51i)
1076 fl%tc_slope_limiter = 3
1077 case ('koren')
1078 ! Barry Koren Right variant
1079 fl%tc_slope_limiter = 4
1080 case default
1081 call mpistop("Unknown tc_slope_limiter, choose MC, minmod")
1082 end select
1083 end subroutine tc_c_params_read_mhd
1084
1085 subroutine tc_c_params_read_hd(fl)
1088 type(tc_fluid), intent(inout) :: fl
1089 integer :: n
1090 logical :: tc_saturate=.false.
1091 double precision :: tc_k_para=0d0
1092
1093 namelist /tc_c_list/ tc_saturate, tc_k_para
1094
1095 do n = 1, size(par_files)
1096 open(unitpar, file=trim(par_files(n)), status="old")
1097 read(unitpar, tc_c_list, end=111)
1098111 close(unitpar)
1099 end do
1100 fl%tc_saturate = tc_saturate
1101 fl%tc_k_para = tc_k_para
1102
1103 end subroutine tc_c_params_read_hd
1104
1105!! end th cond
1106
1107!!rad cool
1108 subroutine rc_params_read_c(fl)
1110 use mod_constants, only: bigdouble
1111 type(rc_fluid), intent(inout) :: fl
1112 integer :: n
1113 ! list parameters
1114 integer :: ncool = 4000
1115 double precision :: cfrac=0.1d0
1116
1117 !> Name of cooling curve
1118 character(len=std_len) :: coolcurve='JCcorona'
1119
1120 !> Name of cooling method
1121 character(len=std_len) :: coolmethod='exact'
1122
1123 !> Fixed temperature not lower than tlow
1124 logical :: tfix=.false.
1125
1126 !> Lower limit of temperature
1127 double precision :: tlow=bigdouble
1128
1129 !> Add cooling source in a split way (.true.) or un-split way (.false.)
1130 logical :: rc_split=.false.
1131
1132
1133 namelist /rc_list_c/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
1134
1135 do n = 1, size(par_files)
1136 open(unitpar, file=trim(par_files(n)), status="old")
1137 read(unitpar, rc_list_c, end=111)
1138111 close(unitpar)
1139 end do
1140
1141 fl%ncool=ncool
1142 fl%coolcurve=coolcurve
1143 fl%coolmethod=coolmethod
1144 fl%tlow=tlow
1145 fl%Tfix=tfix
1146 fl%rc_split=rc_split
1147 fl%cfrac=cfrac
1148 end subroutine rc_params_read_c
1149
1150!! end rad cool
1151
1152 !> sets the equilibrium variables
1153 subroutine set_equi_vars_grid_faces(igrid,x,ixI^L,ixO^L)
1156 use mod_usr_methods
1157 integer, intent(in) :: igrid, ixi^l, ixo^l
1158 double precision, intent(in) :: x(ixi^s,1:ndim)
1159
1160 double precision :: delx(ixi^s,1:ndim)
1161 double precision :: xc(ixi^s,1:ndim),xshift^d
1162 integer :: idims, ixc^l, hxo^l, ix, idims2
1163
1164 if(slab_uniform)then
1165 ^d&delx(ixi^s,^d)=rnode(rpdx^d_,igrid)\
1166 else
1167 ! for all non-cartesian and stretched cartesian coordinates
1168 delx(ixi^s,1:ndim)=ps(igrid)%dx(ixi^s,1:ndim)
1169 endif
1170
1171
1172 do idims=1,ndim
1173 hxo^l=ixo^l-kr(idims,^d);
1174 if(stagger_grid) then
1175 ! ct needs all transverse cells
1176 ixcmax^d=ixomax^d+nghostcells-nghostcells*kr(idims,^d); ixcmin^d=hxomin^d-nghostcells+nghostcells*kr(idims,^d);
1177 else
1178 ! ixC is centered index in the idims direction from ixOmin-1/2 to ixOmax+1/2
1179 ixcmax^d=ixomax^d; ixcmin^d=hxomin^d;
1180 end if
1181 ! always xshift=0 or 1/2
1182 xshift^d=half*(one-kr(^d,idims));
1183 do idims2=1,ndim
1184 select case(idims2)
1185 {case(^d)
1186 do ix = ixc^lim^d
1187 ! xshift=half: this is the cell center coordinate
1188 ! xshift=0: this is the cell edge i+1/2 coordinate
1189 xc(ix^d%ixC^s,^d)=x(ix^d%ixC^s,^d)+(half-xshift^d)*delx(ix^d%ixC^s,^d)
1190 end do\}
1191 end select
1192 end do
1193 call usr_set_equi_vars(ixi^l,ixc^l,xc,ps(igrid)%equi_vars(ixi^s,1:number_equi_vars,idims))
1194 end do
1195
1196 end subroutine set_equi_vars_grid_faces
1197
1198 !> sets the equilibrium variables
1199 subroutine set_equi_vars_grid(igrid)
1201 use mod_usr_methods
1202
1203 integer, intent(in) :: igrid
1204
1205 !values at the center
1206 call usr_set_equi_vars(ixg^ll,ixg^ll,ps(igrid)%x,ps(igrid)%equi_vars(ixg^t,1:number_equi_vars,0))
1207
1208 !values at the interfaces
1209 call set_equi_vars_grid_faces(igrid,ps(igrid)%x,ixg^ll,ixm^ll)
1210
1211 end subroutine set_equi_vars_grid
1212
1213 ! w, wnew conserved
1214 function convert_vars_splitting(ixI^L,ixO^L, w, x, nwc) result(wnew)
1216 integer, intent(in) :: ixi^l,ixo^l, nwc
1217 double precision, intent(in) :: w(ixi^s, 1:nw)
1218 double precision, intent(in) :: x(ixi^s,1:ndim)
1219 double precision :: wnew(ixo^s, 1:nwc)
1220 double precision :: rho(ixi^s)
1221
1222 call get_rhon_tot(w,x,ixi^l,ixo^l,rho(ixi^s))
1223 wnew(ixo^s,rho_n_) = rho(ixo^s)
1224 wnew(ixo^s,mom_n(:)) = w(ixo^s,mom_n(:))
1225 call get_rhoc_tot(w,x,ixi^l,ixo^l,rho(ixi^s))
1226 wnew(ixo^s,rho_c_) = rho(ixo^s)
1227 wnew(ixo^s,mom_c(:)) = w(ixo^s,mom_c(:))
1228
1229 if (b0field) then
1230 ! add background magnetic field B0 to B
1231 wnew(ixo^s,mag(:))=w(ixo^s,mag(:))+block%B0(ixo^s,:,0)
1232 else
1233 wnew(ixo^s,mag(:))=w(ixo^s,mag(:))
1234 end if
1235
1236 if(phys_energy) then
1237 wnew(ixo^s,e_n_) = w(ixo^s,e_n_)
1238 if(has_equi_pe_n0) then
1239 wnew(ixo^s,e_n_) = wnew(ixo^s,e_n_) + block%equi_vars(ixo^s,equi_pe_n0_,0)* inv_gamma_1
1240 endif
1241 wnew(ixo^s,e_c_) = w(ixo^s,e_c_)
1242 if(has_equi_pe_c0) then
1243 wnew(ixo^s,e_c_) = wnew(ixo^s,e_c_) + block%equi_vars(ixo^s,equi_pe_c0_,0)* inv_gamma_1
1244 endif
1245 if(b0field .and. phys_total_energy) then
1246 wnew(ixo^s,e_c_)=wnew(ixo^s,e_c_)+0.5d0*sum(block%B0(ixo^s,:,0)**2,dim=ndim+1) &
1247 + sum(w(ixo^s,mag(:))*block%B0(ixo^s,:,0),dim=ndim+1)
1248 endif
1249 endif
1250
1251 end function convert_vars_splitting
1252
1253 !> copied from mod_gravity
1254 subroutine grav_params_read(files)
1256 character(len=*), intent(in) :: files(:)
1257 integer :: n
1258
1259 namelist /grav_list/ grav_split
1260
1261 do n = 1, size(files)
1262 open(unitpar, file=trim(files(n)), status="old")
1263 read(unitpar, grav_list, end=111)
1264111 close(unitpar)
1265 end do
1266
1267 end subroutine grav_params_read
1268
1269 subroutine associate_dump_hyper()
1272 integer :: ii
1273 do ii = 1,ndim
1274 if(ii==1) then
1275 call add_convert_method(dump_hyperdiffusivity_coef_x, nw, cons_wnames(1:nw), "hyper_x")
1276 elseif(ii==2) then
1277 call add_convert_method(dump_hyperdiffusivity_coef_y, nw, cons_wnames(1:nw), "hyper_y")
1278 else
1279 call add_convert_method(dump_hyperdiffusivity_coef_z, nw, cons_wnames(1:nw), "hyper_z")
1280 endif
1281 enddo
1282 end subroutine associate_dump_hyper
1283
1284 subroutine twofl_check_params
1286 use mod_usr_methods
1288
1289 ! after user parameter setting
1290 gamma_1=twofl_gamma-1.d0
1291 if (.not. phys_energy) then
1292 if (twofl_gamma <= 0.0d0) call mpistop ("Error: twofl_gamma <= 0")
1293 if (twofl_adiab < 0.0d0) call mpistop ("Error: twofl_adiab < 0")
1295 else
1296 if (twofl_gamma <= 0.0d0 .or. twofl_gamma == 1.0d0) &
1297 call mpistop ("Error: twofl_gamma <= 0 or twofl_gamma == 1")
1298 inv_gamma_1=1.d0/gamma_1
1299 small_e = small_pressure * inv_gamma_1
1300 end if
1301
1302 ! this has to be done here as use_imex_scheme is not set in init subroutine,
1303 ! but here it is
1304 if(use_imex_scheme) then
1305 if(has_collisions()) then
1306 ! implicit collisional terms update
1307 phys_implicit_update => twofl_implicit_coll_terms_update
1308 phys_evaluate_implicit => twofl_evaluate_implicit
1309 if(mype .eq. 1) then
1310 print*, "IMPLICIT UPDATE with calc_mult_factor", twofl_implicit_calc_mult_method
1311 endif
1312 if(twofl_implicit_calc_mult_method == 1) then
1313 calc_mult_factor => calc_mult_factor1
1314 else
1315 calc_mult_factor => calc_mult_factor2
1316 endif
1317 endif
1318 else
1319 ! check dtcoll par for explicit implementation of the coll. terms
1320 if(dtcollpar .le. 0d0 .or. dtcollpar .ge. 1d0) then
1321 if (mype .eq. 0) print*, "Explicit update of coll terms requires 0<dtcollpar<1, dtcollpar set to 0.8."
1322 dtcollpar = 0.8
1323 endif
1324
1325 endif
1326! if(H_ion_fr == 0d0 .and. He_ion_fr == 0d0) then
1327! call mpistop("H_ion_fr or He_ion_fr must be > 0 or use hd module")
1328! endif
1329! if(H_ion_fr == 1d0 .and. He_ion_fr == 1d0) then
1330! call mpistop("H_ion_fr or He_ion_fr must be < 1 or use mhd module")
1331! endif
1332 if (number_equi_vars > 0 .and. .not. associated(usr_set_equi_vars)) then
1333 call mpistop("usr_set_equi_vars has to be implemented in the user file")
1334 endif
1335 if(convert .or. autoconvert) then
1336 if(convert_type .eq. 'dat_generic_mpi') then
1337 if(twofl_dump_full_vars) then
1338 if(mype .eq. 0) print*, " add conversion method: split -> full "
1339 call add_convert_method(convert_vars_splitting, nw, cons_wnames, "new")
1340 endif
1341 if(twofl_dump_coll_terms) then
1342 if(mype .eq. 0) print*, " add conversion method: dump coll terms "
1343 call add_convert_method(dump_coll_terms, 3, (/"alpha ", "gamma_rec", "gamma_ion"/), "_coll")
1344 endif
1346 if(mype .eq. 0) print*, " add conversion method: dump hyperdiffusivity coeff. "
1347 call associate_dump_hyper()
1348 endif
1349 endif
1350 endif
1351 end subroutine twofl_check_params
1352
1353 subroutine twofl_physical_units()
1355 double precision :: mp,kb,miu0,c_lightspeed
1356 double precision :: a,b
1357 ! Derive scaling units
1358 if(si_unit) then
1359 mp=mp_si
1360 kb=kb_si
1361 miu0=miu0_si
1362 c_lightspeed=c_si
1363 else
1364 mp=mp_cgs
1365 kb=kb_cgs
1366 miu0=4.d0*dpi
1367 c_lightspeed=const_c
1368 end if
1369
1370 a=1d0
1371 b=1d0
1372 rc=2d0
1373 rn=1d0
1374 ! assume unit_length is alway specified
1375 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
1376 if(unit_density/=1.d0) then
1378 else if(unit_numberdensity/=1.d0) then
1380 end if
1381 if(unit_temperature/=1.d0) then
1386 else if(unit_magneticfield/=1.d0) then
1391 else if(unit_pressure/=1.d0) then
1396 else if(unit_velocity/=1.d0) then
1401 else if(unit_time/=1.d0) then
1406 end if
1407 else if(unit_temperature/=1.d0) then
1408 ! units of temperature and velocity are dependent
1409 if(unit_magneticfield/=1.d0) then
1415 else if(unit_pressure/=1.d0) then
1421 end if
1422 else if(unit_magneticfield/=1.d0) then
1423 ! units of magnetic field and pressure are dependent
1424 if(unit_velocity/=1.d0) then
1430 else if(unit_time/=0.d0) then
1436 end if
1437 else if(unit_pressure/=1.d0) then
1438 if(unit_velocity/=1.d0) then
1444 else if(unit_time/=0.d0) then
1450 end if
1451 end if
1452 ! Additional units needed for the particles
1453 c_norm=c_lightspeed/unit_velocity
1455 if (.not. si_unit) unit_charge = unit_charge*const_c
1457 end subroutine twofl_physical_units
1458
1459 subroutine twofl_check_w(primitive,ixI^L,ixO^L,w,flag)
1461
1462 logical, intent(in) :: primitive
1463 integer, intent(in) :: ixi^l, ixo^l
1464 double precision, intent(in) :: w(ixi^s,nw)
1465 double precision :: tmp(ixi^s)
1466 logical, intent(inout) :: flag(ixi^s,1:nw)
1467
1468 flag=.false.
1469
1470 if(has_equi_rho_n0) then
1471 tmp(ixo^s) = w(ixo^s,rho_n_) + block%equi_vars(ixo^s,equi_rho_n0_,0)
1472 else
1473 tmp(ixo^s) = w(ixo^s,rho_n_)
1474 endif
1475 where(tmp(ixo^s) < small_density) flag(ixo^s,rho_n_) = .true.
1476 if(has_equi_rho_c0) then
1477 tmp(ixo^s) = w(ixo^s,rho_c_) + block%equi_vars(ixo^s,equi_rho_c0_,0)
1478 else
1479 tmp(ixo^s) = w(ixo^s,rho_c_)
1480 endif
1481 where(tmp(ixo^s) < small_density) flag(ixo^s,rho_c_) = .true.
1482 if(phys_energy) then
1483 if(primitive) then
1484 tmp(ixo^s) = w(ixo^s,e_n_)
1485 if(has_equi_pe_n0) then
1486 tmp(ixo^s) = tmp(ixo^s)+block%equi_vars(ixo^s,equi_pe_n0_,0)
1487 endif
1488 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_n_) = .true.
1489 tmp(ixo^s) = w(ixo^s,e_c_)
1490 if(has_equi_pe_c0) then
1491 tmp(ixo^s) = tmp(ixo^s)+block%equi_vars(ixo^s,equi_pe_c0_,0)
1492 endif
1493 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_c_) = .true.
1494 else
1495 if(phys_internal_e) then
1496 tmp(ixo^s)=w(ixo^s,e_n_)
1497 if(has_equi_pe_n0) then
1498 tmp(ixo^s) = tmp(ixo^s)+block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
1499 endif
1500 where(tmp(ixo^s) < small_e) flag(ixo^s,e_n_) = .true.
1501 tmp(ixo^s)=w(ixo^s,e_c_)
1502 if(has_equi_pe_c0) then
1503 tmp(ixo^s) = tmp(ixo^s)+block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
1504 endif
1505 where(tmp(ixo^s) < small_e) flag(ixo^s,e_c_) = .true.
1506 else
1507 !neutrals
1508 tmp(ixo^s)=w(ixo^s,e_n_)-&
1509 twofl_kin_en_n(w,ixi^l,ixo^l)
1510 if(has_equi_pe_n0) then
1511 tmp(ixo^s) = tmp(ixo^s)+block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
1512 endif
1513 where(tmp(ixo^s) < small_e) flag(ixo^s,e_n_) = .true.
1514 if(phys_total_energy) then
1515 tmp(ixo^s)=w(ixo^s,e_c_)-&
1516 twofl_kin_en_c(w,ixi^l,ixo^l)-twofl_mag_en(w,ixi^l,ixo^l)
1517 else
1518 tmp(ixo^s)=w(ixo^s,e_c_)-&
1519 twofl_kin_en_c(w,ixi^l,ixo^l)
1520 end if
1521 if(has_equi_pe_c0) then
1522 tmp(ixo^s) = tmp(ixo^s)+block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
1523 endif
1524 where(tmp(ixo^s) < small_e) flag(ixo^s,e_c_) = .true.
1525 end if
1526 endif
1527 end if
1528
1529 end subroutine twofl_check_w
1530
1531 !> Transform primitive variables into conservative ones
1532 subroutine twofl_to_conserved(ixI^L,ixO^L,w,x)
1534 integer, intent(in) :: ixi^l, ixo^l
1535 double precision, intent(inout) :: w(ixi^s, nw)
1536 double precision, intent(in) :: x(ixi^s, 1:ndim)
1537 integer :: idir
1538 double precision :: rhoc(ixi^s)
1539 double precision :: rhon(ixi^s)
1540
1541 !if (fix_small_values) then
1542 ! call twofl_handle_small_values(.true., w, x, ixI^L, ixO^L, 'twofl_to_conserved')
1543 !end if
1544
1545 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
1546 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
1547
1548 ! Calculate total energy from pressure, kinetic and magnetic energy
1549 if(phys_energy) then
1550 if(phys_internal_e) then
1551 w(ixo^s,e_n_)=w(ixo^s,e_n_)*inv_gamma_1
1552 w(ixo^s,e_c_)=w(ixo^s,e_c_)*inv_gamma_1
1553 else
1554 w(ixo^s,e_n_)=w(ixo^s,e_n_)*inv_gamma_1&
1555 +half*sum(w(ixo^s,mom_n(:))**2,dim=ndim+1)*rhon(ixo^s)
1556 if(phys_total_energy) then
1557 w(ixo^s,e_c_)=w(ixo^s,e_c_)*inv_gamma_1&
1558 +half*sum(w(ixo^s,mom_c(:))**2,dim=ndim+1)*rhoc(ixo^s)&
1559 +twofl_mag_en(w, ixi^l, ixo^l)
1560 else
1561 ! kinetic energy + internal energy is evolved
1562 w(ixo^s,e_c_)=w(ixo^s,e_c_)*inv_gamma_1&
1563 +half*sum(w(ixo^s,mom_c(:))**2,dim=ndim+1)*rhoc(ixo^s)
1564 end if
1565 end if
1566 end if
1567
1568 ! Convert velocity to momentum
1569 do idir = 1, ndir
1570 w(ixo^s, mom_n(idir)) = rhon(ixo^s) * w(ixo^s, mom_n(idir))
1571 w(ixo^s, mom_c(idir)) = rhoc(ixo^s) * w(ixo^s, mom_c(idir))
1572 end do
1573 end subroutine twofl_to_conserved
1574
1575 !> Transform conservative variables into primitive ones
1576 subroutine twofl_to_primitive(ixI^L,ixO^L,w,x)
1578 integer, intent(in) :: ixi^l, ixo^l
1579 double precision, intent(inout) :: w(ixi^s, nw)
1580 double precision, intent(in) :: x(ixi^s, 1:ndim)
1581 integer :: idir
1582 double precision :: rhoc(ixi^s)
1583 double precision :: rhon(ixi^s)
1584
1585 if (fix_small_values) then
1586 call twofl_handle_small_values(.false., w, x, ixi^l, ixo^l, 'twofl_to_primitive')
1587 end if
1588
1589 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
1590 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
1591
1592 if(phys_energy) then
1593 if(phys_internal_e) then
1594 w(ixo^s,e_n_)=w(ixo^s,e_n_)*gamma_1
1595 w(ixo^s,e_c_)=w(ixo^s,e_c_)*gamma_1
1596 else
1597 ! neutrals evolved energy = ke + e_int
1598 w(ixo^s,e_n_)=gamma_1*(w(ixo^s,e_n_)&
1599 -twofl_kin_en_n(w,ixi^l,ixo^l))
1600 ! charges
1601 if(phys_total_energy) then
1602 ! evolved energy = ke + e_int + e_mag
1603 w(ixo^s,e_c_)=gamma_1*(w(ixo^s,e_c_)&
1604 -twofl_kin_en_c(w,ixi^l,ixo^l)&
1605 -twofl_mag_en(w,ixi^l,ixo^l))
1606 else
1607 ! evolved energy = ke + e_int
1608 w(ixo^s,e_c_)=gamma_1*(w(ixo^s,e_c_)&
1609 -twofl_kin_en_c(w,ixi^l,ixo^l))
1610 end if
1611 end if
1612 end if
1613
1614 ! Convert momentum to velocity
1615 do idir = 1, ndir
1616 w(ixo^s, mom_c(idir)) = w(ixo^s, mom_c(idir))/rhoc(ixo^s)
1617 w(ixo^s, mom_n(idir)) = w(ixo^s, mom_n(idir))/rhon(ixo^s)
1618 end do
1619
1620 end subroutine twofl_to_primitive
1621
1622!!USED IN TC
1623 !> Transform internal energy to total energy
1624 subroutine twofl_ei_to_e_c(ixI^L,ixO^L,w,x)
1626 integer, intent(in) :: ixi^l, ixo^l
1627 double precision, intent(inout) :: w(ixi^s, nw)
1628 double precision, intent(in) :: x(ixi^s, 1:ndim)
1629
1630 ! Calculate total energy from internal, kinetic and magnetic energy
1631 if(twofl_eq_energy == eq_energy_ki) then
1632 w(ixo^s,e_c_)=w(ixo^s,e_c_)&
1633 +twofl_kin_en_c(w,ixi^l,ixo^l)
1634 else
1635 w(ixo^s,e_c_)=w(ixo^s,e_c_)&
1636 +twofl_kin_en_c(w,ixi^l,ixo^l)&
1637 +twofl_mag_en(w,ixi^l,ixo^l)
1638 endif
1639 end subroutine twofl_ei_to_e_c
1640
1641 !> Transform total energy to internal energy
1642 subroutine twofl_e_to_ei_c(ixI^L,ixO^L,w,x)
1644 integer, intent(in) :: ixi^l, ixo^l
1645 double precision, intent(inout) :: w(ixi^s, nw)
1646 double precision, intent(in) :: x(ixi^s, 1:ndim)
1647
1648 if(twofl_eq_energy == eq_energy_ki) then
1649 w(ixo^s,e_c_)=w(ixo^s,e_c_)&
1650 -twofl_kin_en_c(w,ixi^l,ixo^l)
1651 else
1652 ! Calculate ei = e - ek - eb
1653 w(ixo^s,e_c_)=w(ixo^s,e_c_)&
1654 -twofl_kin_en_c(w,ixi^l,ixo^l)&
1655 -twofl_mag_en(w,ixi^l,ixo^l)
1656 endif
1657 end subroutine twofl_e_to_ei_c
1658
1659 !Neutrals
1660 subroutine twofl_ei_to_e_n(ixI^L,ixO^L,w,x)
1662 integer, intent(in) :: ixi^l, ixo^l
1663 double precision, intent(inout) :: w(ixi^s, nw)
1664 double precision, intent(in) :: x(ixi^s, 1:ndim)
1665
1666 ! Calculate total energy from internal and kinetic energy
1667
1668 w(ixo^s,e_n_)=w(ixo^s,e_n_)+twofl_kin_en_n(w,ixi^l,ixo^l)
1669
1670 end subroutine twofl_ei_to_e_n
1671
1672 !> Transform total energy to internal energy
1673 subroutine twofl_e_to_ei_n(ixI^L,ixO^L,w,x)
1675 integer, intent(in) :: ixi^l, ixo^l
1676 double precision, intent(inout) :: w(ixi^s, nw)
1677 double precision, intent(in) :: x(ixi^s, 1:ndim)
1678
1679 ! Calculate ei = e - ek
1680 w(ixo^s,e_n_)=w(ixo^s,e_n_)-twofl_kin_en_n(w,ixi^l,ixo^l)
1681
1682 call twofl_handle_small_ei_n(w,x,ixi^l,ixo^l,e_n_,"e_to_ei_n")
1683 end subroutine twofl_e_to_ei_n
1684
1685 subroutine twofl_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1688 logical, intent(in) :: primitive
1689 integer, intent(in) :: ixi^l,ixo^l
1690 double precision, intent(inout) :: w(ixi^s,1:nw)
1691 double precision, intent(in) :: x(ixi^s,1:ndim)
1692 character(len=*), intent(in) :: subname
1693
1694 integer :: idir
1695 logical :: flag(ixi^s,1:nw)
1696 double precision :: tmp2(ixi^s)
1697 double precision :: tmp1(ixi^s)
1698
1699 call twofl_check_w(primitive, ixi^l, ixo^l, w, flag)
1700
1701 if(any(flag)) then
1702 select case (small_values_method)
1703 case ("replace")
1704 if(has_equi_rho_c0) then
1705 where(flag(ixo^s,rho_c_)) w(ixo^s,rho_c_) = &
1706 small_density-block%equi_vars(ixo^s,equi_rho_c0_,0)
1707 else
1708 where(flag(ixo^s,rho_c_)) w(ixo^s,rho_c_) = small_density
1709 end if
1710 if(has_equi_rho_n0) then
1711 where(flag(ixo^s,rho_n_)) w(ixo^s,rho_n_) = &
1712 small_density-block%equi_vars(ixo^s,equi_rho_n0_,0)
1713 else
1714 where(flag(ixo^s,rho_n_)) w(ixo^s,rho_n_) = small_density
1715 end if
1716 do idir = 1, ndir
1717 if(small_values_fix_iw(mom_n(idir))) then
1718 where(flag(ixo^s,rho_n_)) w(ixo^s, mom_n(idir)) = 0.0d0
1719 end if
1720 if(small_values_fix_iw(mom_c(idir))) then
1721 where(flag(ixo^s,rho_c_)) w(ixo^s, mom_c(idir)) = 0.0d0
1722 end if
1723 end do
1724
1725 if(phys_energy) then
1726 if(primitive) then
1727 if(has_equi_pe_n0) then
1728 tmp1(ixo^s) = small_pressure - &
1729 block%equi_vars(ixo^s,equi_pe_n0_,0)
1730 else
1731 tmp1(ixo^s) = small_pressure
1732 end if
1733 if(has_equi_pe_c0) then
1734 tmp2(ixo^s) = small_e - &
1735 block%equi_vars(ixo^s,equi_pe_c0_,0)
1736 else
1737 tmp2(ixo^s) = small_pressure
1738 end if
1739 else
1740 ! conserved
1741 if(has_equi_pe_n0) then
1742 tmp1(ixo^s) = small_e - &
1743 block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
1744 else
1745 tmp1(ixo^s) = small_e
1746 end if
1747 if(has_equi_pe_c0) then
1748 tmp2(ixo^s) = small_e - &
1749 block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
1750 else
1751 tmp2(ixo^s) = small_e
1752 end if
1753 if(phys_internal_e) then
1754 where(flag(ixo^s,e_n_))
1755 w(ixo^s,e_n_)=tmp1(ixo^s)
1756 end where
1757 where(flag(ixo^s,e_c_))
1758 w(ixo^s,e_c_)=tmp2(ixo^s)
1759 end where
1760 else
1761 where(flag(ixo^s,e_n_))
1762 w(ixo^s,e_n_) = tmp1(ixo^s)+&
1763 twofl_kin_en_n(w,ixi^l,ixo^l)
1764 end where
1765 if(phys_total_energy) then
1766 where(flag(ixo^s,e_c_))
1767 w(ixo^s,e_c_) = tmp2(ixo^s)+&
1768 twofl_kin_en_c(w,ixi^l,ixo^l)+&
1769 twofl_mag_en(w,ixi^l,ixo^l)
1770 end where
1771 else
1772 where(flag(ixo^s,e_c_))
1773 w(ixo^s,e_c_) = tmp2(ixo^s)+&
1774 twofl_kin_en_c(w,ixi^l,ixo^l)
1775 end where
1776 end if
1777 end if
1778 end if
1779 end if
1780 case ("average")
1781 call small_values_average(ixi^l, ixo^l, w, x, flag)
1782 case default
1783 if(.not.primitive) then
1784 !convert w to primitive
1785 ! Calculate pressure = (gamma-1) * (e-ek-eb)
1786 if(phys_energy) then
1787 if(phys_internal_e) then
1788 w(ixo^s,e_c_)=w(ixo^s,e_c_)*gamma_1
1789 w(ixo^s,e_n_)=w(ixo^s,e_n_)*gamma_1
1790 else
1791 w(ixo^s,e_n_)=gamma_1*(w(ixo^s,e_n_)&
1792 -twofl_kin_en_n(w,ixi^l,ixo^l))
1793 if(phys_total_energy) then
1794 w(ixo^s,e_c_)=gamma_1*(w(ixo^s,e_c_)&
1795 -twofl_kin_en_c(w,ixi^l,ixo^l)&
1796 -twofl_mag_en(w,ixi^l,ixo^l))
1797 else
1798 w(ixo^s,e_c_)=gamma_1*(w(ixo^s,e_c_)&
1799 -twofl_kin_en_c(w,ixi^l,ixo^l))
1800
1801 end if
1802 end if
1803 end if
1804 ! Convert momentum to velocity
1805 if(has_equi_rho_n0) then
1806 tmp1(ixo^s) = w(ixo^s,rho_n_) + block%equi_vars(ixo^s,equi_rho_n0_,0)
1807 else
1808 tmp1(ixo^s) = w(ixo^s,rho_n_)
1809 end if
1810
1811 if(has_equi_rho_c0) then
1812 tmp2(ixo^s) = w(ixo^s,rho_c_) + block%equi_vars(ixo^s,equi_rho_c0_,0)
1813 else
1814 tmp2(ixo^s) = w(ixo^s,rho_c_)
1815 end if
1816 do idir = 1, ndir
1817 w(ixo^s, mom_n(idir)) = w(ixo^s, mom_n(idir))/tmp1(ixo^s)
1818 w(ixo^s, mom_c(idir)) = w(ixo^s, mom_c(idir))/tmp2(ixo^s)
1819 end do
1820 end if
1821 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1822 end select
1823 end if
1824 end subroutine twofl_handle_small_values
1825
1826 !> Calculate cmax_idim=csound+abs(v_idim) within ixO^L
1827 subroutine twofl_get_cmax(w,x,ixI^L,ixO^L,idim,cmax)
1829
1830 integer, intent(in) :: ixi^l, ixo^l, idim
1831 ! w in primitive form
1832 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
1833 double precision, intent(inout) :: cmax(ixi^s)
1834 double precision :: cmax2(ixi^s),rhon(ixi^s)
1835
1836 call twofl_get_csound_c_idim(w,x,ixi^l,ixo^l,idim,cmax)
1837 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
1838 if(phys_energy) then
1839 if(has_equi_pe_n0) then
1840 cmax2(ixo^s)=sqrt(twofl_gamma*(w(ixo^s,e_n_)+&
1841 block%equi_vars(ixo^s,equi_pe_n0_,b0i))/rhon(ixo^s))
1842 else
1843 cmax2(ixo^s)=sqrt(twofl_gamma*w(ixo^s,e_n_)/rhon(ixo^s))
1844 end if
1845 else
1846 cmax2(ixo^s)=sqrt(twofl_gamma*twofl_adiab*rhon(ixo^s)**gamma_1)
1847 end if
1848 cmax(ixo^s)=max(abs(w(ixo^s,mom_n(idim)))+cmax2(ixo^s),&
1849 abs(w(ixo^s,mom_c(idim)))+cmax(ixo^s))
1850
1851 end subroutine twofl_get_cmax
1852
1853 subroutine twofl_get_a2max(w,x,ixI^L,ixO^L,a2max)
1855
1856 integer, intent(in) :: ixi^l, ixo^l
1857 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
1858 double precision, intent(inout) :: a2max(ndim)
1859 double precision :: a2(ixi^s,ndim,nw)
1860 integer :: gxo^l,hxo^l,jxo^l,kxo^l,i,j
1861
1862 a2=zero
1863 do i = 1,ndim
1864 !> 4th order
1865 hxo^l=ixo^l-kr(i,^d);
1866 gxo^l=hxo^l-kr(i,^d);
1867 jxo^l=ixo^l+kr(i,^d);
1868 kxo^l=jxo^l+kr(i,^d);
1869 a2(ixo^s,i,1:nw)=abs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
1870 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
1871 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
1872 end do
1873 end subroutine twofl_get_a2max
1874
1875 ! COPIED from hd/moh_hd_phys
1876 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
1877 subroutine twofl_get_tcutoff_n(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
1879 integer, intent(in) :: ixi^l,ixo^l
1880 double precision, intent(in) :: x(ixi^s,1:ndim),w(ixi^s,1:nw)
1881 double precision, intent(out) :: tco_local, tmax_local
1882
1883 double precision, parameter :: delta=0.25d0
1884 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
1885 integer :: jxo^l,hxo^l
1886 logical :: lrlt(ixi^s)
1887
1888 {^ifoned
1889 ! reuse lts as rhon
1890 call get_rhon_tot(w,x,ixi^l,ixi^l,lts)
1891 tmp1(ixi^s)=w(ixi^s,e_n_)-0.5d0*sum(w(ixi^s,mom_n(:))**2,dim=ndim+1)/lts(ixi^s)
1892 te(ixi^s)=tmp1(ixi^s)/lts(ixi^s)*(twofl_gamma-1.d0)
1893
1894 tmax_local=maxval(te(ixo^s))
1895
1896 hxo^l=ixo^l-1;
1897 jxo^l=ixo^l+1;
1898 lts(ixo^s)=0.5d0*abs(te(jxo^s)-te(hxo^s))/te(ixo^s)
1899 lrlt=.false.
1900 where(lts(ixo^s) > delta)
1901 lrlt(ixo^s)=.true.
1902 end where
1903 tco_local=zero
1904 if(any(lrlt(ixo^s))) then
1905 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
1906 end if
1907 }
1908 end subroutine twofl_get_tcutoff_n
1909
1910 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
1911 subroutine twofl_get_tcutoff_c(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
1913 use mod_geometry
1914 integer, intent(in) :: ixi^l,ixo^l
1915 double precision, intent(in) :: x(ixi^s,1:ndim)
1916 double precision, intent(inout) :: w(ixi^s,1:nw)
1917 double precision, intent(out) :: tco_local,tmax_local
1918
1919 double precision, parameter :: trac_delta=0.25d0
1920 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
1921 double precision, dimension(ixI^S,1:ndir) :: bunitvec
1922 double precision, dimension(ixI^S,1:ndim) :: gradt
1923 double precision :: bdir(ndim)
1924 double precision :: ltr(ixi^s),ltrc,ltrp,altr(ixi^s)
1925 integer :: idims,jxo^l,hxo^l,ixa^d,ixb^d
1926 integer :: jxp^l,hxp^l,ixp^l
1927 logical :: lrlt(ixi^s)
1928
1929 ! reuse lts as rhoc
1930 call get_rhoc_tot(w,x,ixi^l,ixi^l,lts)
1931 if(phys_internal_e) then
1932 tmp1(ixi^s)=w(ixi^s,e_c_)
1933 else
1934 tmp1(ixi^s)=w(ixi^s,e_c_)-0.5d0*(sum(w(ixi^s,mom_c(:))**2,dim=ndim+1)/&
1935 lts(ixi^s)+sum(w(ixi^s,mag(:))**2,dim=ndim+1))
1936 end if
1937 te(ixi^s)=tmp1(ixi^s)/lts(ixi^s)*(twofl_gamma-1.d0)
1938 tmax_local=maxval(te(ixo^s))
1939
1940 {^ifoned
1941 select case(twofl_trac_type)
1942 case(0)
1943 !> test case, fixed cutoff temperature
1944 w(ixi^s,tcoff_c_)=2.5d5/unit_temperature
1945 case(1)
1946 hxo^l=ixo^l-1;
1947 jxo^l=ixo^l+1;
1948 lts(ixo^s)=0.5d0*abs(te(jxo^s)-te(hxo^s))/te(ixo^s)
1949 lrlt=.false.
1950 where(lts(ixo^s) > trac_delta)
1951 lrlt(ixo^s)=.true.
1952 end where
1953 if(any(lrlt(ixo^s))) then
1954 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
1955 end if
1956 case(2)
1957 !> iijima et al. 2021, LTRAC method
1958 ltrc=1.5d0
1959 ltrp=2.5d0
1960 ixp^l=ixo^l^ladd1;
1961 hxo^l=ixo^l-1;
1962 jxo^l=ixo^l+1;
1963 hxp^l=ixp^l-1;
1964 jxp^l=ixp^l+1;
1965 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
1966 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
1967 w(ixo^s,tcoff_c_)=te(ixo^s)*&
1968 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
1969 case default
1970 call mpistop("twofl_trac_type not allowed for 1D simulation")
1971 end select
1972 }
1973 {^nooned
1974 select case(twofl_trac_type)
1975 case(0)
1976 !> test case, fixed cutoff temperature
1977 w(ixi^s,tcoff_c_)=2.5d5/unit_temperature
1978 case(1,4,6)
1979 ! temperature gradient at cell centers
1980 do idims=1,ndim
1981 call gradient(te,ixi^l,ixo^l,idims,tmp1)
1982 gradt(ixo^s,idims)=tmp1(ixo^s)
1983 end do
1984 ! B vector
1985 if(b0field) then
1986 bunitvec(ixo^s,:)=w(ixo^s,iw_mag(:))+block%B0(ixo^s,:,0)
1987 else
1988 bunitvec(ixo^s,:)=w(ixo^s,iw_mag(:))
1989 end if
1990 if(twofl_trac_type .gt. 1) then
1991 ! B direction at cell center
1992 bdir=zero
1993 {do ixa^d=0,1\}
1994 ixb^d=(ixomin^d+ixomax^d-1)/2+ixa^d;
1995 bdir(1:ndim)=bdir(1:ndim)+bunitvec(ixb^d,1:ndim)
1996 {end do\}
1997 if(sum(bdir(:)**2) .gt. zero) then
1998 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
1999 end if
2000 block%special_values(3:ndim+2)=bdir(1:ndim)
2001 end if
2002 tmp1(ixo^s)=dsqrt(sum(bunitvec(ixo^s,:)**2,dim=ndim+1))
2003 where(tmp1(ixo^s)/=0.d0)
2004 tmp1(ixo^s)=1.d0/tmp1(ixo^s)
2005 elsewhere
2006 tmp1(ixo^s)=bigdouble
2007 end where
2008 ! b unit vector: magnetic field direction vector
2009 do idims=1,ndim
2010 bunitvec(ixo^s,idims)=bunitvec(ixo^s,idims)*tmp1(ixo^s)
2011 end do
2012 ! temperature length scale inversed
2013 lts(ixo^s)=abs(sum(gradt(ixo^s,1:ndim)*bunitvec(ixo^s,1:ndim),dim=ndim+1))/te(ixo^s)
2014 ! fraction of cells size to temperature length scale
2015 if(slab_uniform) then
2016 lts(ixo^s)=minval(dxlevel)*lts(ixo^s)
2017 else
2018 lts(ixo^s)=minval(block%ds(ixo^s,:),dim=ndim+1)*lts(ixo^s)
2019 end if
2020 lrlt=.false.
2021 where(lts(ixo^s) > trac_delta)
2022 lrlt(ixo^s)=.true.
2023 end where
2024 if(any(lrlt(ixo^s))) then
2025 block%special_values(1)=maxval(te(ixo^s), mask=lrlt(ixo^s))
2026 else
2027 block%special_values(1)=zero
2028 end if
2029 block%special_values(2)=tmax_local
2030 case(2)
2031 !> iijima et al. 2021, LTRAC method
2032 ltrc=1.5d0
2033 ltrp=4.d0
2034 ixp^l=ixo^l^ladd1;
2035 ! temperature gradient at cell centers
2036 do idims=1,ndim
2037 call gradient(te,ixi^l,ixp^l,idims,tmp1)
2038 gradt(ixp^s,idims)=tmp1(ixp^s)
2039 end do
2040 ! B vector
2041 if(b0field) then
2042 bunitvec(ixp^s,:)=w(ixp^s,iw_mag(:))+block%B0(ixp^s,:,0)
2043 else
2044 bunitvec(ixp^s,:)=w(ixp^s,iw_mag(:))
2045 end if
2046 tmp1(ixp^s)=dsqrt(sum(bunitvec(ixp^s,:)**2,dim=ndim+1))
2047 where(tmp1(ixp^s)/=0.d0)
2048 tmp1(ixp^s)=1.d0/tmp1(ixp^s)
2049 elsewhere
2050 tmp1(ixp^s)=bigdouble
2051 end where
2052 ! b unit vector: magnetic field direction vector
2053 do idims=1,ndim
2054 bunitvec(ixp^s,idims)=bunitvec(ixp^s,idims)*tmp1(ixp^s)
2055 end do
2056 ! temperature length scale inversed
2057 lts(ixp^s)=abs(sum(gradt(ixp^s,1:ndim)*bunitvec(ixp^s,1:ndim),dim=ndim+1))/te(ixp^s)
2058 ! fraction of cells size to temperature length scale
2059 if(slab_uniform) then
2060 lts(ixp^s)=minval(dxlevel)*lts(ixp^s)
2061 else
2062 lts(ixp^s)=minval(block%ds(ixp^s,:),dim=ndim+1)*lts(ixp^s)
2063 end if
2064 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
2065
2066 altr(ixi^s)=zero
2067 do idims=1,ndim
2068 hxo^l=ixo^l-kr(idims,^d);
2069 jxo^l=ixo^l+kr(idims,^d);
2070 altr(ixo^s)=altr(ixo^s) &
2071 +0.25*(ltr(hxo^s)+two*ltr(ixo^s)+ltr(jxo^s))*bunitvec(ixo^s,idims)**2
2072 w(ixo^s,tcoff_c_)=te(ixo^s)*altr(ixo^s)**(0.4*ltrp)
2073 end do
2074 case(3,5)
2075 !> do nothing here
2076 case default
2077 call mpistop("unknown twofl_trac_type")
2078 end select
2079 }
2080 end subroutine twofl_get_tcutoff_c
2081
2082 !> get H speed for H-correction to fix the carbuncle problem at grid-aligned shock front
2083 subroutine twofl_get_h_speed_one(wprim,x,ixI^L,ixO^L,idim,Hspeed)
2085
2086 integer, intent(in) :: ixi^l, ixo^l, idim
2087 double precision, intent(in) :: wprim(ixi^s, nw)
2088 double precision, intent(in) :: x(ixi^s,1:ndim)
2089 double precision, intent(out) :: hspeed(ixi^s,1:number_species)
2090
2091 double precision :: csound(ixi^s,ndim),tmp(ixi^s)
2092 integer :: jxc^l, ixc^l, ixa^l, id, ix^d
2093
2094 hspeed=0.d0
2095 ixa^l=ixo^l^ladd1;
2096 do id=1,ndim
2097 call twofl_get_csound_prim(wprim,x,ixi^l,ixa^l,id,tmp)
2098 csound(ixa^s,id)=tmp(ixa^s)
2099 end do
2100 ixcmax^d=ixomax^d;
2101 ixcmin^d=ixomin^d+kr(idim,^d)-1;
2102 jxcmax^d=ixcmax^d+kr(idim,^d);
2103 jxcmin^d=ixcmin^d+kr(idim,^d);
2104 hspeed(ixc^s,1)=0.5d0*abs(&
2105 0.5d0 * (wprim(jxc^s,mom_c(idim))+ wprim(jxc^s,mom_n(idim))) &
2106 +csound(jxc^s,idim)- &
2107 0.5d0 * (wprim(ixc^s,mom_c(idim)) + wprim(ixc^s,mom_n(idim)))&
2108 +csound(ixc^s,idim))
2109
2110 do id=1,ndim
2111 if(id==idim) cycle
2112 ixamax^d=ixcmax^d+kr(id,^d);
2113 ixamin^d=ixcmin^d+kr(id,^d);
2114 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(&
2115 0.5d0 * (wprim(ixa^s,mom_c(id)) + wprim(ixa^s,mom_n(id)))&
2116 +csound(ixa^s,id)-&
2117 0.5d0 * (wprim(ixc^s,mom_c(id)) + wprim(ixc^s,mom_n(id)))&
2118 +csound(ixc^s,id)))
2119
2120
2121 ixamax^d=ixcmax^d-kr(id,^d);
2122 ixamin^d=ixcmin^d-kr(id,^d);
2123 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(&
2124 0.5d0 * (wprim(ixc^s,mom_c(id)) + wprim(ixc^s,mom_n(id)))&
2125 +csound(ixc^s,id)-&
2126 0.5d0 * (wprim(ixa^s,mom_c(id)) + wprim(ixa^s,mom_n(id)))&
2127 +csound(ixa^s,id)))
2128
2129 end do
2130
2131 do id=1,ndim
2132 if(id==idim) cycle
2133 ixamax^d=jxcmax^d+kr(id,^d);
2134 ixamin^d=jxcmin^d+kr(id,^d);
2135 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(&
2136 0.5d0 * (wprim(ixa^s,mom_c(id)) + wprim(ixa^s,mom_n(id)))&
2137 +csound(ixa^s,id)-&
2138 0.5d0 * (wprim(jxc^s,mom_c(id)) + wprim(jxc^s,mom_n(id)))&
2139 +csound(jxc^s,id)))
2140 ixamax^d=jxcmax^d-kr(id,^d);
2141 ixamin^d=jxcmin^d-kr(id,^d);
2142 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(&
2143 0.5d0 * (wprim(jxc^s,mom_c(id)) + wprim(jxc^s,mom_n(id)))&
2144 +csound(jxc^s,id)-&
2145 0.5d0 * (wprim(ixa^s,mom_c(id)) + wprim(ixa^s,mom_n(id)))&
2146 +csound(ixa^s,id)))
2147 end do
2148
2149 end subroutine twofl_get_h_speed_one
2150
2151 !> get H speed for H-correction to fix the carbuncle problem at grid-aligned shock front
2152 subroutine twofl_get_h_speed_species(wprim,x,ixI^L,ixO^L,idim,Hspeed)
2154
2155 integer, intent(in) :: ixi^l, ixo^l, idim
2156 double precision, intent(in) :: wprim(ixi^s, nw)
2157 double precision, intent(in) :: x(ixi^s,1:ndim)
2158 double precision, intent(out) :: hspeed(ixi^s,1:number_species)
2159
2160 double precision :: csound(ixi^s,ndim),tmp(ixi^s)
2161 integer :: jxc^l, ixc^l, ixa^l, id, ix^d
2162
2163 hspeed=0.d0
2164 ! charges
2165 ixa^l=ixo^l^ladd1;
2166 do id=1,ndim
2167 call twofl_get_csound_prim_c(wprim,x,ixi^l,ixa^l,id,tmp)
2168 csound(ixa^s,id)=tmp(ixa^s)
2169 end do
2170 ixcmax^d=ixomax^d;
2171 ixcmin^d=ixomin^d+kr(idim,^d)-1;
2172 jxcmax^d=ixcmax^d+kr(idim,^d);
2173 jxcmin^d=ixcmin^d+kr(idim,^d);
2174 hspeed(ixc^s,1)=0.5d0*abs(wprim(jxc^s,mom_c(idim))+csound(jxc^s,idim)-wprim(ixc^s,mom_c(idim))+csound(ixc^s,idim))
2175
2176 do id=1,ndim
2177 if(id==idim) cycle
2178 ixamax^d=ixcmax^d+kr(id,^d);
2179 ixamin^d=ixcmin^d+kr(id,^d);
2180 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(wprim(ixa^s,mom_c(id))+csound(ixa^s,id)-wprim(ixc^s,mom_c(id))+csound(ixc^s,id)))
2181 ixamax^d=ixcmax^d-kr(id,^d);
2182 ixamin^d=ixcmin^d-kr(id,^d);
2183 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(wprim(ixc^s,mom_c(id))+csound(ixc^s,id)-wprim(ixa^s,mom_c(id))+csound(ixa^s,id)))
2184 end do
2185
2186 do id=1,ndim
2187 if(id==idim) cycle
2188 ixamax^d=jxcmax^d+kr(id,^d);
2189 ixamin^d=jxcmin^d+kr(id,^d);
2190 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(wprim(ixa^s,mom_c(id))+csound(ixa^s,id)-wprim(jxc^s,mom_c(id))+csound(jxc^s,id)))
2191 ixamax^d=jxcmax^d-kr(id,^d);
2192 ixamin^d=jxcmin^d-kr(id,^d);
2193 hspeed(ixc^s,1)=max(hspeed(ixc^s,1),0.5d0*abs(wprim(jxc^s,mom_c(id))+csound(jxc^s,id)-wprim(ixa^s,mom_c(id))+csound(ixa^s,id)))
2194 end do
2195
2196 ! neutrals
2197 ixa^l=ixo^l^ladd1;
2198 do id=1,ndim
2199 call twofl_get_csound_prim_n(wprim,x,ixi^l,ixa^l,id,tmp)
2200 csound(ixa^s,id)=tmp(ixa^s)
2201 end do
2202 ixcmax^d=ixomax^d;
2203 ixcmin^d=ixomin^d+kr(idim,^d)-1;
2204 jxcmax^d=ixcmax^d+kr(idim,^d);
2205 jxcmin^d=ixcmin^d+kr(idim,^d);
2206 hspeed(ixc^s,2)=0.5d0*abs(wprim(jxc^s,mom_n(idim))+csound(jxc^s,idim)-wprim(ixc^s,mom_n(idim))+csound(ixc^s,idim))
2207
2208 do id=1,ndim
2209 if(id==idim) cycle
2210 ixamax^d=ixcmax^d+kr(id,^d);
2211 ixamin^d=ixcmin^d+kr(id,^d);
2212 hspeed(ixc^s,2)=max(hspeed(ixc^s,2),0.5d0*abs(wprim(ixa^s,mom_n(id))+csound(ixa^s,id)-wprim(ixc^s,mom_n(id))+csound(ixc^s,id)))
2213 ixamax^d=ixcmax^d-kr(id,^d);
2214 ixamin^d=ixcmin^d-kr(id,^d);
2215 hspeed(ixc^s,2)=max(hspeed(ixc^s,2),0.5d0*abs(wprim(ixc^s,mom_n(id))+csound(ixc^s,id)-wprim(ixa^s,mom_n(id))+csound(ixa^s,id)))
2216 end do
2217
2218 do id=1,ndim
2219 if(id==idim) cycle
2220 ixamax^d=jxcmax^d+kr(id,^d);
2221 ixamin^d=jxcmin^d+kr(id,^d);
2222 hspeed(ixc^s,2)=max(hspeed(ixc^s,2),0.5d0*abs(wprim(ixa^s,mom_n(id))+csound(ixa^s,id)-wprim(jxc^s,mom_n(id))+csound(jxc^s,id)))
2223 ixamax^d=jxcmax^d-kr(id,^d);
2224 ixamin^d=jxcmin^d-kr(id,^d);
2225 hspeed(ixc^s,2)=max(hspeed(ixc^s,2),0.5d0*abs(wprim(jxc^s,mom_n(id))+csound(jxc^s,id)-wprim(ixa^s,mom_n(id))+csound(ixa^s,id)))
2226 end do
2227
2228 end subroutine twofl_get_h_speed_species
2229
2230 !> Estimating bounds for the minimum and maximum signal velocities
2231 subroutine twofl_get_cbounds_one(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
2234
2235 integer, intent(in) :: ixi^l, ixo^l, idim
2236 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
2237 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
2238 double precision, intent(in) :: x(ixi^s,1:ndim)
2239 double precision, intent(inout) :: cmax(ixi^s,number_species)
2240 double precision, intent(inout), optional :: cmin(ixi^s,number_species)
2241 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
2242
2243 double precision :: wmean(ixi^s,nw)
2244 double precision :: rhon(ixi^s)
2245 double precision :: rhoc(ixi^s)
2246 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
2247 integer :: ix^d
2248
2249 select case (boundspeed)
2250 case (1)
2251 ! This implements formula (10.52) from "Riemann Solvers and Numerical
2252 ! Methods for Fluid Dynamics" by Toro.
2253 call get_rhoc_tot(wlp,x,ixi^l,ixo^l,rhoc)
2254 call get_rhon_tot(wlp,x,ixi^l,ixo^l,rhon)
2255 tmp1(ixo^s)=sqrt(abs(rhoc(ixo^s) +rhon(ixo^s)))
2256
2257 call get_rhoc_tot(wrp,x,ixi^l,ixo^l,rhoc)
2258 call get_rhon_tot(wrp,x,ixi^l,ixo^l,rhon)
2259 tmp2(ixo^s)=sqrt(abs(rhoc(ixo^s) +rhon(ixo^s)))
2260
2261 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
2262 umean(ixo^s)=(0.5*(wlp(ixo^s,mom_n(idim))+wlp(ixo^s,mom_c(idim)))*tmp1(ixo^s) + &
2263 0.5*(wrp(ixo^s,mom_n(idim))+wrp(ixo^s,mom_c(idim)))*tmp2(ixo^s))*tmp3(ixo^s)
2264 call twofl_get_csound_prim(wlp,x,ixi^l,ixo^l,idim,csoundl)
2265 call twofl_get_csound_prim(wrp,x,ixi^l,ixo^l,idim,csoundr)
2266
2267 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)**2+tmp2(ixo^s)*csoundr(ixo^s)**2)*tmp3(ixo^s)+&
2268 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2*(&
2269 0.5*(wrp(ixo^s,mom_n(idim))+wrp(ixo^s,mom_c(idim)))- &
2270 0.5*(wlp(ixo^s,mom_n(idim))+wlp(ixo^s,mom_c(idim))))**2
2271 dmean(ixo^s)=sqrt(dmean(ixo^s))
2272 if(present(cmin)) then
2273 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
2274 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
2275 if(h_correction) then
2276 {do ix^db=ixomin^db,ixomax^db\}
2277 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
2278 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
2279 {end do\}
2280 end if
2281 else
2282 cmax(ixo^s,1)=abs(umean(ixo^s))+dmean(ixo^s)
2283 end if
2284 case (2)
2285 ! typeboundspeed=='cmaxmean'
2286 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
2287 call get_rhon_tot(wmean,x,ixi^l,ixo^l,rhon)
2288 tmp2(ixo^s)=wmean(ixo^s,mom_n(idim))/rhon(ixo^s)
2289 call get_rhoc_tot(wmean,x,ixi^l,ixo^l,rhoc)
2290 tmp1(ixo^s)=wmean(ixo^s,mom_c(idim))/rhoc(ixo^s)
2291 call twofl_get_csound(wmean,x,ixi^l,ixo^l,idim,csoundr)
2292 if(present(cmin)) then
2293 cmax(ixo^s,1)=max(max(abs(tmp2(ixo^s)), abs(tmp1(ixo^s)) ) +csoundr(ixo^s),zero)
2294 cmin(ixo^s,1)=min(min(abs(tmp2(ixo^s)), abs(tmp1(ixo^s)) ) -csoundr(ixo^s),zero)
2295 if(h_correction) then
2296 {do ix^db=ixomin^db,ixomax^db\}
2297 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
2298 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
2299 {end do\}
2300 end if
2301 else
2302 cmax(ixo^s,1)= max(abs(tmp2(ixo^s)),abs(tmp1(ixo^s)))+csoundr(ixo^s)
2303 end if
2304 case (3)
2305 ! Miyoshi 2005 JCP 208, 315 equation (67)
2306 call twofl_get_csound(wlp,x,ixi^l,ixo^l,idim,csoundl)
2307 call twofl_get_csound(wrp,x,ixi^l,ixo^l,idim,csoundr)
2308 csoundl(ixo^s)=max(csoundl(ixo^s),csoundr(ixo^s))
2309 if(present(cmin)) then
2310 cmin(ixo^s,1)=min(0.5*(wlp(ixo^s,mom_c(idim))+ wrp(ixo^s,mom_n(idim))),&
2311 0.5*(wrp(ixo^s,mom_c(idim))+ wrp(ixo^s,mom_n(idim))))-csoundl(ixo^s)
2312 cmax(ixo^s,1)=max(0.5*(wlp(ixo^s,mom_c(idim))+ wrp(ixo^s,mom_n(idim))),&
2313 0.5*(wrp(ixo^s,mom_c(idim))+ wrp(ixo^s,mom_n(idim))))+csoundl(ixo^s)
2314 if(h_correction) then
2315 {do ix^db=ixomin^db,ixomax^db\}
2316 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
2317 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
2318 {end do\}
2319 end if
2320 else
2321 cmax(ixo^s,1)=max(0.5*(wlp(ixo^s,mom_c(idim))+ wrp(ixo^s,mom_n(idim))),&
2322 0.5*(wrp(ixo^s,mom_c(idim))+ wrp(ixo^s,mom_n(idim))))+csoundl(ixo^s)
2323 end if
2324 end select
2325
2326 end subroutine twofl_get_cbounds_one
2327
2328 !> Calculate fast magnetosonic wave speed
2329 subroutine twofl_get_csound_prim_c(w,x,ixI^L,ixO^L,idim,csound)
2331
2332 integer, intent(in) :: ixi^l, ixo^l, idim
2333 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
2334 double precision, intent(out):: csound(ixi^s)
2335 double precision :: cfast2(ixi^s), avmincs2(ixi^s), b2(ixi^s), kmax
2336 double precision :: inv_rho(ixo^s)
2337 double precision :: rhoc(ixi^s)
2338
2339 integer :: ix1,ix2
2340
2341
2342 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
2343 inv_rho(ixo^s)=1.d0/rhoc(ixo^s)
2344
2345 if(phys_energy) then
2346 call twofl_get_pthermal_c_primitive(w,x,ixi^l,ixo^l,csound)
2347 csound(ixo^s)=twofl_gamma*csound(ixo^s)/rhoc(ixo^s)
2348 else
2349 call twofl_get_csound2_adiab_c(w,x,ixi^l,ixo^l,csound)
2350 endif
2351
2352 ! store |B|^2 in v
2353 b2(ixo^s) = twofl_mag_en_all(w,ixi^l,ixo^l)
2354 cfast2(ixo^s) = b2(ixo^s) * inv_rho(ixo^s)+csound(ixo^s)
2355 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
2356 * twofl_mag_i_all(w,ixi^l,ixo^l,idim)**2 &
2357 * inv_rho(ixo^s)
2358
2359 where(avmincs2(ixo^s)<zero)
2360 avmincs2(ixo^s)=zero
2361 end where
2362
2363 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
2364
2365 if (.not. twofl_hall) then
2366 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
2367 else
2368 ! take the Hall velocity into account:
2369 ! most simple estimate, high k limit:
2370 ! largest wavenumber supported by grid: Nyquist (in practise can reduce by some factor)
2371 kmax = dpi/min({dxlevel(^d)},bigdouble)*half
2372 csound(ixo^s) = max(sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s))), &
2373 twofl_etah * sqrt(b2(ixo^s))*inv_rho(ixo^s)*kmax)
2374 end if
2375
2376 end subroutine twofl_get_csound_prim_c
2377
2378 !> Calculate fast magnetosonic wave speed
2379 subroutine twofl_get_csound_prim_n(w,x,ixI^L,ixO^L,idim,csound)
2381
2382 integer, intent(in) :: ixi^l, ixo^l, idim
2383 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
2384 double precision, intent(out):: csound(ixi^s)
2385 double precision :: rhon(ixi^s)
2386
2387 if(phys_energy) then
2388 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
2389 call twofl_get_pthermal_n_primitive(w,x,ixi^l,ixo^l,csound)
2390 csound(ixo^s)=twofl_gamma*csound(ixo^s)/rhon(ixo^s)
2391 else
2392 call twofl_get_csound2_adiab_n(w,x,ixi^l,ixo^l,csound)
2393 endif
2394 csound(ixo^s) = sqrt(csound(ixo^s))
2395
2396 end subroutine twofl_get_csound_prim_n
2397
2398 !> Estimating bounds for the minimum and maximum signal velocities
2399 subroutine twofl_get_cbounds_species(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
2402 use mod_variables
2403
2404 integer, intent(in) :: ixi^l, ixo^l, idim
2405 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
2406 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
2407 double precision, intent(in) :: x(ixi^s,1:ndim)
2408 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
2409 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
2410 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
2411
2412 double precision :: wmean(ixi^s,nw)
2413 double precision :: rho(ixi^s)
2414 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
2415 integer :: ix^d
2416
2417 select case (boundspeed)
2418 case (1)
2419 ! This implements formula (10.52) from "Riemann Solvers and Numerical
2420 ! Methods for Fluid Dynamics" by Toro.
2421 ! charges
2422 call get_rhoc_tot(wlp,x,ixi^l,ixo^l,rho)
2423 tmp1(ixo^s)=sqrt(abs(rho(ixo^s)))
2424
2425 call get_rhoc_tot(wrp,x,ixi^l,ixo^l,rho)
2426 tmp2(ixo^s)=sqrt(abs(rho(ixo^s)))
2427
2428 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
2429 umean(ixo^s)=(wlp(ixo^s,mom_c(idim))*tmp1(ixo^s)+wrp(ixo^s,mom_c(idim))*tmp2(ixo^s))*tmp3(ixo^s)
2430 call twofl_get_csound_prim_c(wlp,x,ixi^l,ixo^l,idim,csoundl)
2431 call twofl_get_csound_prim_c(wrp,x,ixi^l,ixo^l,idim,csoundr)
2432
2433
2434 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)**2+tmp2(ixo^s)*csoundr(ixo^s)**2)*tmp3(ixo^s)+&
2435 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2*&
2436 (wrp(ixo^s,mom_c(idim)) - wlp(ixo^s,mom_c(idim)))**2
2437 dmean(ixo^s)=sqrt(dmean(ixo^s))
2438 if(present(cmin)) then
2439 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
2440 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
2441 if(h_correction) then
2442 {do ix^db=ixomin^db,ixomax^db\}
2443 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
2444 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
2445 {end do\}
2446 end if
2447 else
2448 cmax(ixo^s,1)=abs(umean(ixo^s))+dmean(ixo^s)
2449 end if
2450
2451 ! neutrals
2452
2453 call get_rhon_tot(wlp,x,ixi^l,ixo^l,rho)
2454 tmp1(ixo^s)=sqrt(abs(rho(ixo^s)))
2455
2456 call get_rhon_tot(wrp,x,ixi^l,ixo^l,rho)
2457 tmp2(ixo^s)=sqrt(abs(rho(ixo^s)))
2458
2459 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
2460 umean(ixo^s)=(wlp(ixo^s,mom_n(idim))*tmp1(ixo^s)+wrp(ixo^s,mom_n(idim))*tmp2(ixo^s))*tmp3(ixo^s)
2461 call twofl_get_csound_prim_n(wlp,x,ixi^l,ixo^l,idim,csoundl)
2462 call twofl_get_csound_prim_n(wrp,x,ixi^l,ixo^l,idim,csoundr)
2463
2464
2465 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)**2+tmp2(ixo^s)*csoundr(ixo^s)**2)*tmp3(ixo^s)+&
2466 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2*&
2467 (wrp(ixo^s,mom_n(idim)) - wlp(ixo^s,mom_n(idim)))**2
2468 dmean(ixo^s)=sqrt(dmean(ixo^s))
2469 if(present(cmin)) then
2470 cmin(ixo^s,2)=umean(ixo^s)-dmean(ixo^s)
2471 cmax(ixo^s,2)=umean(ixo^s)+dmean(ixo^s)
2472 if(h_correction) then
2473 {do ix^db=ixomin^db,ixomax^db\}
2474 cmin(ix^d,2)=sign(one,cmin(ix^d,2))*max(abs(cmin(ix^d,2)),hspeed(ix^d,2))
2475 cmax(ix^d,2)=sign(one,cmax(ix^d,2))*max(abs(cmax(ix^d,2)),hspeed(ix^d,2))
2476 {end do\}
2477 end if
2478 else
2479 cmax(ixo^s,2)=abs(umean(ixo^s))+dmean(ixo^s)
2480 end if
2481
2482 case (2)
2483 ! typeboundspeed=='cmaxmean'
2484 wmean(ixo^s,1:nwflux)=0.5d0*(wlp(ixo^s,1:nwflux)+wrp(ixo^s,1:nwflux))
2485 ! charges
2486 tmp1(ixo^s)=wmean(ixo^s,mom_c(idim))
2487 call twofl_get_csound_c_idim(wmean,x,ixi^l,ixo^l,idim,csoundr)
2488 if(present(cmin)) then
2489 cmax(ixo^s,1)=max(abs(tmp1(ixo^s))+csoundr(ixo^s),zero)
2490 cmin(ixo^s,1)=min(abs(tmp1(ixo^s))-csoundr(ixo^s),zero)
2491 if(h_correction) then
2492 {do ix^db=ixomin^db,ixomax^db\}
2493 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
2494 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
2495 {end do\}
2496 end if
2497 else
2498 cmax(ixo^s,1)=abs(tmp1(ixo^s))+csoundr(ixo^s)
2499 end if
2500 !neutrals
2501
2502 tmp1(ixo^s)=wmean(ixo^s,mom_n(idim))
2503 call twofl_get_csound_n(wmean,x,ixi^l,ixo^l,csoundr)
2504 if(present(cmin)) then
2505 cmax(ixo^s,2)=max(abs(tmp1(ixo^s))+csoundr(ixo^s),zero)
2506 cmin(ixo^s,2)=min(abs(tmp1(ixo^s))-csoundr(ixo^s),zero)
2507 if(h_correction) then
2508 {do ix^db=ixomin^db,ixomax^db\}
2509 cmin(ix^d,2)=sign(one,cmin(ix^d,2))*max(abs(cmin(ix^d,2)),hspeed(ix^d,2))
2510 cmax(ix^d,2)=sign(one,cmax(ix^d,2))*max(abs(cmax(ix^d,2)),hspeed(ix^d,2))
2511 {end do\}
2512 end if
2513 else
2514 cmax(ixo^s,2)= abs(tmp1(ixo^s))+csoundr(ixo^s)
2515 end if
2516 case (3)
2517 ! Miyoshi 2005 JCP 208, 315 equation (67)
2518 call twofl_get_csound_c_idim(wlp,x,ixi^l,ixo^l,idim,csoundl)
2519 call twofl_get_csound_c_idim(wrp,x,ixi^l,ixo^l,idim,csoundr)
2520 csoundl(ixo^s)=max(csoundl(ixo^s),csoundr(ixo^s))
2521 if(present(cmin)) then
2522 cmin(ixo^s,1)=min(wlp(ixo^s,mom_c(idim)),wrp(ixo^s,mom_c(idim)))-csoundl(ixo^s)
2523 cmax(ixo^s,1)=max(wlp(ixo^s,mom_c(idim)),wrp(ixo^s,mom_c(idim)))+csoundl(ixo^s)
2524 if(h_correction) then
2525 {do ix^db=ixomin^db,ixomax^db\}
2526 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
2527 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
2528 {end do\}
2529 end if
2530 else
2531 cmax(ixo^s,1)=max(wlp(ixo^s,mom_c(idim)),wrp(ixo^s,mom_c(idim)))+csoundl(ixo^s)
2532 end if
2533 call twofl_get_csound_n(wlp,x,ixi^l,ixo^l,csoundl)
2534 call twofl_get_csound_n(wrp,x,ixi^l,ixo^l,csoundr)
2535 csoundl(ixo^s)=max(csoundl(ixo^s),csoundr(ixo^s))
2536 if(present(cmin)) then
2537 cmin(ixo^s,2)=min(wlp(ixo^s,mom_n(idim)),wrp(ixo^s,mom_n(idim)))-csoundl(ixo^s)
2538 cmax(ixo^s,2)=max(wlp(ixo^s,mom_n(idim)),wrp(ixo^s,mom_n(idim)))+csoundl(ixo^s)
2539 if(h_correction) then
2540 {do ix^db=ixomin^db,ixomax^db\}
2541 cmin(ix^d,2)=sign(one,cmin(ix^d,2))*max(abs(cmin(ix^d,1)),hspeed(ix^d,2))
2542 cmax(ix^d,2)=sign(one,cmax(ix^d,2))*max(abs(cmax(ix^d,1)),hspeed(ix^d,2))
2543 {end do\}
2544 end if
2545 else
2546 cmax(ixo^s,2)=max(wlp(ixo^s,mom_n(idim)),wrp(ixo^s,mom_n(idim)))+csoundl(ixo^s)
2547 end if
2548
2549 end select
2550
2551 end subroutine twofl_get_cbounds_species
2552
2553 !> prepare velocities for ct methods
2554 subroutine twofl_get_ct_velocity_average(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
2556
2557 integer, intent(in) :: ixi^l, ixo^l, idim
2558 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
2559 double precision, intent(in) :: cmax(ixi^s)
2560 double precision, intent(in), optional :: cmin(ixi^s)
2561 type(ct_velocity), intent(inout):: vcts
2562
2563 end subroutine twofl_get_ct_velocity_average
2564
2565 subroutine twofl_get_ct_velocity_contact(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
2567
2568 integer, intent(in) :: ixi^l, ixo^l, idim
2569 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
2570 double precision, intent(in) :: cmax(ixi^s)
2571 double precision, intent(in), optional :: cmin(ixi^s)
2572 type(ct_velocity), intent(inout):: vcts
2573
2574 if(.not.allocated(vcts%vnorm)) allocate(vcts%vnorm(ixi^s,1:ndim))
2575 ! get average normal velocity at cell faces
2576 vcts%vnorm(ixo^s,idim)=0.5d0*(wlp(ixo^s,mom_c(idim))+wrp(ixo^s,mom_c(idim)))
2577
2578 end subroutine twofl_get_ct_velocity_contact
2579
2580 subroutine twofl_get_ct_velocity_hll(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
2582
2583 integer, intent(in) :: ixi^l, ixo^l, idim
2584 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
2585 double precision, intent(in) :: cmax(ixi^s)
2586 double precision, intent(in), optional :: cmin(ixi^s)
2587 type(ct_velocity), intent(inout):: vcts
2588
2589 integer :: idime,idimn
2590
2591 if(.not.allocated(vcts%vbarC)) then
2592 allocate(vcts%vbarC(ixi^s,1:ndir,2),vcts%vbarLC(ixi^s,1:ndir,2),vcts%vbarRC(ixi^s,1:ndir,2))
2593 allocate(vcts%cbarmin(ixi^s,1:ndim),vcts%cbarmax(ixi^s,1:ndim))
2594 end if
2595 ! Store magnitude of characteristics
2596 if(present(cmin)) then
2597 vcts%cbarmin(ixo^s,idim)=max(-cmin(ixo^s),zero)
2598 vcts%cbarmax(ixo^s,idim)=max( cmax(ixo^s),zero)
2599 else
2600 vcts%cbarmax(ixo^s,idim)=max( cmax(ixo^s),zero)
2601 vcts%cbarmin(ixo^s,idim)=vcts%cbarmax(ixo^s,idim)
2602 end if
2603
2604 idimn=mod(idim,ndir)+1 ! 'Next' direction
2605 idime=mod(idim+1,ndir)+1 ! Electric field direction
2606 ! Store velocities
2607 vcts%vbarLC(ixo^s,idim,1)=wlp(ixo^s,mom_c(idimn))
2608 vcts%vbarRC(ixo^s,idim,1)=wrp(ixo^s,mom_c(idimn))
2609 vcts%vbarC(ixo^s,idim,1)=(vcts%cbarmax(ixo^s,idim)*vcts%vbarLC(ixo^s,idim,1) &
2610 +vcts%cbarmin(ixo^s,idim)*vcts%vbarRC(ixo^s,idim,1))&
2611 /(vcts%cbarmax(ixo^s,idim)+vcts%cbarmin(ixo^s,idim))
2612
2613 vcts%vbarLC(ixo^s,idim,2)=wlp(ixo^s,mom_c(idime))
2614 vcts%vbarRC(ixo^s,idim,2)=wrp(ixo^s,mom_c(idime))
2615 vcts%vbarC(ixo^s,idim,2)=(vcts%cbarmax(ixo^s,idim)*vcts%vbarLC(ixo^s,idim,2) &
2616 +vcts%cbarmin(ixo^s,idim)*vcts%vbarRC(ixo^s,idim,1))&
2617 /(vcts%cbarmax(ixo^s,idim)+vcts%cbarmin(ixo^s,idim))
2618
2619 end subroutine twofl_get_ct_velocity_hll
2620
2621 subroutine twofl_get_csound_c_idim(w,x,ixI^L,ixO^L,idim,csound)
2623
2624 integer, intent(in) :: ixi^l, ixo^l, idim
2625 ! w in primitive form
2626 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
2627 double precision, intent(out):: csound(ixi^s)
2628 double precision :: cfast2(ixi^s), avmincs2(ixi^s), b2(ixi^s), kmax
2629 double precision :: inv_rho(ixo^s)
2630 double precision :: tmp(ixi^s)
2631#if (!defined(ONE_FLUID) || ONE_FLUID==0) && (defined(A_TOT) && A_TOT == 1)
2632 double precision :: rhon(ixi^s)
2633#endif
2634 call get_rhoc_tot(w,x,ixi^l,ixo^l,tmp)
2635#if (!defined(ONE_FLUID) || ONE_FLUID==0) && (defined(A_TOT) && A_TOT == 1)
2636 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
2637 inv_rho(ixo^s) = 1d0/(rhon(ixo^s)+tmp(ixo^s))
2638#else
2639 inv_rho(ixo^s)=1.d0/tmp(ixo^s)
2640#endif
2641
2642 if(phys_energy) then
2643 csound(ixo^s)=twofl_gamma*w(ixo^s,e_c_)*inv_rho(ixo^s)
2644 else
2645 csound(ixo^s)=twofl_gamma*twofl_adiab*tmp(ixo^s)**gamma_1
2646 end if
2647
2648 ! store |B|^2 in v
2649 b2(ixo^s) = twofl_mag_en_all(w,ixi^l,ixo^l)
2650
2651 cfast2(ixo^s) = b2(ixo^s) * inv_rho(ixo^s)+csound(ixo^s)
2652 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
2653 * twofl_mag_i_all(w,ixi^l,ixo^l,idim)**2 &
2654 * inv_rho(ixo^s)
2655
2656 where(avmincs2(ixo^s)<zero)
2657 avmincs2(ixo^s)=zero
2658 end where
2659
2660 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
2661
2662 if (.not. twofl_hall) then
2663 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
2664 else
2665 ! take the Hall velocity into account:
2666 ! most simple estimate, high k limit:
2667 ! largest wavenumber supported by grid: Nyquist (in practise can reduce by some factor)
2668 kmax = dpi/min({dxlevel(^d)},bigdouble)*half
2669 csound(ixo^s) = max(sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s))), &
2670 twofl_etah * sqrt(b2(ixo^s))*inv_rho(ixo^s)*kmax)
2671 end if
2672
2673 end subroutine twofl_get_csound_c_idim
2674
2675 !> Calculate fast magnetosonic wave speed when cbounds_species=false
2676 subroutine twofl_get_csound_prim(w,x,ixI^L,ixO^L,idim,csound)
2678
2679 integer, intent(in) :: ixi^l, ixo^l, idim
2680 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
2681 double precision, intent(out):: csound(ixi^s)
2682 double precision :: cfast2(ixi^s), avmincs2(ixi^s), b2(ixi^s), kmax
2683 double precision :: inv_rho(ixo^s)
2684 double precision :: rhoc(ixi^s)
2685#if (defined(A_TOT) && A_TOT == 1)
2686 double precision :: rhon(ixi^s)
2687#endif
2688 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
2689#if (defined(A_TOT) && A_TOT == 1)
2690 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
2691 inv_rho(ixo^s) = 1d0/(rhon(ixo^s)+rhoc(ixo^s))
2692#else
2693 inv_rho(ixo^s)=1.d0/rhoc(ixo^s)
2694#endif
2695
2696 call twofl_get_csound2_primitive(w,x,ixi^l,ixo^l,csound)
2697
2698 ! store |B|^2 in v
2699 b2(ixo^s) = twofl_mag_en_all(w,ixi^l,ixo^l)
2700 cfast2(ixo^s) = b2(ixo^s) * inv_rho(ixo^s)+csound(ixo^s)
2701 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
2702 * twofl_mag_i_all(w,ixi^l,ixo^l,idim)**2 &
2703 * inv_rho(ixo^s)
2704
2705 where(avmincs2(ixo^s)<zero)
2706 avmincs2(ixo^s)=zero
2707 end where
2708
2709 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
2710
2711 if (.not. twofl_hall) then
2712 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
2713 else
2714 ! take the Hall velocity into account:
2715 ! most simple estimate, high k limit:
2716 ! largest wavenumber supported by grid: Nyquist (in practise can reduce by some factor)
2717 kmax = dpi/min({dxlevel(^d)},bigdouble)*half
2718 csound(ixo^s) = max(sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s))), &
2719 twofl_etah * sqrt(b2(ixo^s))*inv_rho(ixo^s)*kmax)
2720 end if
2721
2722 contains
2723 !TODO copy it inside
2724 subroutine twofl_get_csound2_primitive(w,x,ixI^L,ixO^L,csound2)
2726 integer, intent(in) :: ixI^L, ixO^L
2727 double precision, intent(in) :: w(ixI^S,nw)
2728 double precision, intent(in) :: x(ixI^S,1:ndim)
2729 double precision, intent(out) :: csound2(ixI^S)
2730 double precision :: pth_c(ixI^S)
2731 double precision :: pth_n(ixI^S)
2732
2733 if(phys_energy) then
2734 call twofl_get_pthermal_c_primitive(w,x,ixi^l,ixo^l,pth_c)
2735 call twofl_get_pthermal_n_primitive(w,x,ixi^l,ixo^l,pth_n)
2736 call twofl_get_csound2_from_pthermal(w,x,ixi^l,ixo^l,pth_c,pth_n,csound2)
2737 else
2738 call twofl_get_csound2_adiab(w,x,ixi^l,ixo^l,csound2)
2739 endif
2740 end subroutine twofl_get_csound2_primitive
2741
2742 end subroutine twofl_get_csound_prim
2743
2744 subroutine twofl_get_csound2(w,x,ixI^L,ixO^L,csound2)
2746 integer, intent(in) :: ixI^L, ixO^L
2747 double precision, intent(in) :: w(ixI^S,nw)
2748 double precision, intent(in) :: x(ixI^S,1:ndim)
2749 double precision, intent(out) :: csound2(ixI^S)
2750 double precision :: pth_c(ixI^S)
2751 double precision :: pth_n(ixI^S)
2752
2753 if(phys_energy) then
2754 call twofl_get_pthermal_c(w,x,ixi^l,ixo^l,pth_c)
2755 call twofl_get_pthermal_n(w,x,ixi^l,ixo^l,pth_n)
2756 call twofl_get_csound2_from_pthermal(w,x,ixi^l,ixo^l,pth_c,pth_n,csound2)
2757 else
2758 call twofl_get_csound2_adiab(w,x,ixi^l,ixo^l,csound2)
2759 endif
2760 end subroutine twofl_get_csound2
2761
2762 subroutine twofl_get_csound2_adiab(w,x,ixI^L,ixO^L,csound2)
2764 integer, intent(in) :: ixI^L, ixO^L
2765 double precision, intent(in) :: w(ixI^S,nw)
2766 double precision, intent(in) :: x(ixI^S,1:ndim)
2767 double precision, intent(out) :: csound2(ixI^S)
2768 double precision :: rhoc(ixI^S)
2769 double precision :: rhon(ixI^S)
2770
2771 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
2772 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
2773 csound2(ixo^s)=twofl_gamma*twofl_adiab*&
2774 max((rhoc(ixo^s)**twofl_gamma + rhon(ixo^s)**twofl_gamma)/(rhoc(ixo^s)+ rhon(ixo^s)),&
2775 rhon(ixo^s)**gamma_1,rhoc(ixo^s)**gamma_1)
2776 end subroutine twofl_get_csound2_adiab
2777
2778 subroutine twofl_get_csound(w,x,ixI^L,ixO^L,idim,csound)
2780
2781 integer, intent(in) :: ixI^L, ixO^L, idim
2782 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
2783 double precision, intent(out):: csound(ixI^S)
2784 double precision :: cfast2(ixI^S), AvMinCs2(ixI^S), b2(ixI^S), kmax
2785 double precision :: inv_rho(ixO^S)
2786 double precision :: rhoc(ixI^S)
2787#if (defined(A_TOT) && A_TOT == 1)
2788 double precision :: rhon(ixI^S)
2789#endif
2790 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
2791#if (defined(A_TOT) && A_TOT == 1)
2792 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
2793 inv_rho(ixo^s) = 1d0/(rhon(ixo^s)+rhoc(ixo^s))
2794#else
2795 inv_rho(ixo^s)=1.d0/rhoc(ixo^s)
2796#endif
2797
2798 call twofl_get_csound2(w,x,ixi^l,ixo^l,csound)
2799
2800 ! store |B|^2 in v
2801 b2(ixo^s) = twofl_mag_en_all(w,ixi^l,ixo^l)
2802
2803 cfast2(ixo^s) = b2(ixo^s) * inv_rho(ixo^s)+csound(ixo^s)
2804 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
2805 * twofl_mag_i_all(w,ixi^l,ixo^l,idim)**2 &
2806 * inv_rho(ixo^s)
2807
2808 where(avmincs2(ixo^s)<zero)
2809 avmincs2(ixo^s)=zero
2810 end where
2811
2812 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
2813
2814 if (.not. twofl_hall) then
2815 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
2816 else
2817 ! take the Hall velocity into account:
2818 ! most simple estimate, high k limit:
2819 ! largest wavenumber supported by grid: Nyquist (in practise can reduce by some factor)
2820 kmax = dpi/min({dxlevel(^d)},bigdouble)*half
2821 csound(ixo^s) = max(sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s))), &
2822 twofl_etah * sqrt(b2(ixo^s))*inv_rho(ixo^s)*kmax)
2823 end if
2824
2825 end subroutine twofl_get_csound
2826
2827 subroutine twofl_get_csound2_from_pthermal(w,x,ixI^L,ixO^L,pth_c,pth_n,csound2)
2829 integer, intent(in) :: ixI^L, ixO^L
2830 double precision, intent(in) :: w(ixI^S,nw)
2831 double precision, intent(in) :: x(ixI^S,1:ndim)
2832 double precision, intent(in) :: pth_c(ixI^S)
2833 double precision, intent(in) :: pth_n(ixI^S)
2834 double precision, intent(out) :: csound2(ixI^S)
2835 double precision :: csound1(ixI^S),rhon(ixI^S),rhoc(ixI^S)
2836
2837 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
2838 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
2839#if !defined(C_TOT) || C_TOT == 0
2840 csound2(ixo^s)=twofl_gamma*max((pth_c(ixo^s) + pth_n(ixo^s))/(rhoc(ixo^s) + rhon(ixo^s)),&
2841 pth_n(ixo^s)/rhon(ixo^s), pth_c(ixo^s)/rhoc(ixo^s))
2842#else
2843 csound2(ixo^s)=twofl_gamma*(csound2(ixo^s) + csound1(ixo^s))/(rhoc(ixo^s) + rhon(ixo^s))
2844
2845#endif
2846 end subroutine twofl_get_csound2_from_pthermal
2847
2848! end cbounds_species=false
2849
2850 subroutine twofl_get_csound_n(w,x,ixI^L,ixO^L,csound)
2852
2853 integer, intent(in) :: ixI^L, ixO^L
2854 double precision, intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
2855 double precision, intent(out):: csound(ixI^S)
2856 double precision :: pe_n1(ixI^S)
2857 call twofl_get_csound2_n_from_conserved(w,x,ixi^l,ixo^l,csound)
2858 csound(ixo^s) = sqrt(csound(ixo^s))
2859 end subroutine twofl_get_csound_n
2860
2861 !> separate routines so that it is faster
2862 !> Calculate temperature=p/rho when in e_ the internal energy is stored
2863 subroutine twofl_get_temperature_from_eint_n(w, x, ixI^L, ixO^L, res)
2865 integer, intent(in) :: ixI^L, ixO^L
2866 double precision, intent(in) :: w(ixI^S, 1:nw)
2867 double precision, intent(in) :: x(ixI^S, 1:ndim)
2868 double precision, intent(out):: res(ixI^S)
2869
2870 res(ixo^s) = 1d0/rn * gamma_1 * w(ixo^s, e_n_) /w(ixo^s,rho_n_)
2871
2872 end subroutine twofl_get_temperature_from_eint_n
2873
2874 subroutine twofl_get_temperature_from_eint_n_with_equi(w, x, ixI^L, ixO^L, res)
2876 integer, intent(in) :: ixI^L, ixO^L
2877 double precision, intent(in) :: w(ixI^S, 1:nw)
2878 double precision, intent(in) :: x(ixI^S, 1:ndim)
2879 double precision, intent(out):: res(ixI^S)
2880
2881 res(ixo^s) = 1d0/rn * (gamma_1 * w(ixo^s, e_n_) + block%equi_vars(ixo^s,equi_pe_n0_,b0i)) /&
2882 (w(ixo^s,rho_n_) +block%equi_vars(ixo^s,equi_rho_n0_,b0i))
2883 end subroutine twofl_get_temperature_from_eint_n_with_equi
2884
2885! subroutine twofl_get_temperature_n_pert_from_tot(Te, ixI^L, ixO^L, res)
2886! use mod_global_parameters
2887! integer, intent(in) :: ixI^L, ixO^L
2888! double precision, intent(in) :: Te(ixI^S)
2889! double precision, intent(out):: res(ixI^S)
2890! res(ixO^S) = Te(ixO^S) -1d0/Rn * &
2891! block%equi_vars(ixO^S,equi_pe_n0_,0)/block%equi_vars(ixO^S,equi_rho_n0_,0)
2892! end subroutine twofl_get_temperature_n_pert_from_tot
2893
2894 subroutine twofl_get_temperature_n_equi(w,x, ixI^L, ixO^L, res)
2896 integer, intent(in) :: ixI^L, ixO^L
2897 double precision, intent(in) :: w(ixI^S, 1:nw)
2898 double precision, intent(in) :: x(ixI^S, 1:ndim)
2899 double precision, intent(out):: res(ixI^S)
2900 res(ixo^s) = 1d0/rn * &
2901 block%equi_vars(ixo^s,equi_pe_n0_,b0i)/block%equi_vars(ixo^s,equi_rho_n0_,b0i)
2902 end subroutine twofl_get_temperature_n_equi
2903
2904 subroutine twofl_get_rho_n_equi(w, x,ixI^L, ixO^L, res)
2906 integer, intent(in) :: ixI^L, ixO^L
2907 double precision, intent(in) :: w(ixI^S, 1:nw)
2908 double precision, intent(in) :: x(ixI^S, 1:ndim)
2909 double precision, intent(out):: res(ixI^S)
2910 res(ixo^s) = block%equi_vars(ixo^s,equi_rho_n0_,b0i)
2911 end subroutine twofl_get_rho_n_equi
2912
2913 subroutine twofl_get_pe_n_equi(w, x, ixI^L, ixO^L, res)
2915 integer, intent(in) :: ixI^L, ixO^L
2916 double precision, intent(in) :: w(ixI^S, 1:nw)
2917 double precision, intent(in) :: x(ixI^S, 1:ndim)
2918 double precision, intent(out):: res(ixI^S)
2919 res(ixo^s) = block%equi_vars(ixo^s,equi_pe_n0_,b0i)
2920 end subroutine twofl_get_pe_n_equi
2921
2922 !> Calculate temperature=p/rho when in e_ the total energy is stored
2923 !> this does not check the values of twofl_energy and twofl_internal_e,
2924 !> twofl_energy = .true. and twofl_internal_e = .false.
2925 !> also check small_values is avoided
2926 subroutine twofl_get_temperature_from_etot_n(w, x, ixI^L, ixO^L, res)
2928 integer, intent(in) :: ixI^L, ixO^L
2929 double precision, intent(in) :: w(ixI^S, 1:nw)
2930 double precision, intent(in) :: x(ixI^S, 1:ndim)
2931 double precision, intent(out):: res(ixI^S)
2932 res(ixo^s)=1d0/rn * (gamma_1*(w(ixo^s,e_n_)&
2933 - twofl_kin_en_n(w,ixi^l,ixo^l)))/w(ixo^s,rho_n_)
2934 end subroutine twofl_get_temperature_from_etot_n
2935
2936 subroutine twofl_get_temperature_from_etot_n_with_equi(w, x, ixI^L, ixO^L, res)
2938 integer, intent(in) :: ixI^L, ixO^L
2939 double precision, intent(in) :: w(ixI^S, 1:nw)
2940 double precision, intent(in) :: x(ixI^S, 1:ndim)
2941 double precision, intent(out):: res(ixI^S)
2942 res(ixo^s)=1d0/rn * (gamma_1*(w(ixo^s,e_n_)&
2943 - twofl_kin_en_n(w,ixi^l,ixo^l)) + block%equi_vars(ixo^s,equi_pe_n0_,b0i))&
2944 /(w(ixo^s,rho_n_) +block%equi_vars(ixo^s,equi_rho_n0_,b0i))
2945
2946 end subroutine twofl_get_temperature_from_etot_n_with_equi
2947
2948 !> separate routines so that it is faster
2949 !> Calculate temperature=p/rho when in e_ the internal energy is stored
2950 subroutine twofl_get_temperature_from_eint_c(w, x, ixI^L, ixO^L, res)
2952 integer, intent(in) :: ixI^L, ixO^L
2953 double precision, intent(in) :: w(ixI^S, 1:nw)
2954 double precision, intent(in) :: x(ixI^S, 1:ndim)
2955 double precision, intent(out):: res(ixI^S)
2956
2957 res(ixo^s) = 1d0/rc * gamma_1 * w(ixo^s, e_c_) /w(ixo^s,rho_c_)
2958
2959 end subroutine twofl_get_temperature_from_eint_c
2960
2961 subroutine twofl_get_temperature_from_eint_c_with_equi(w, x, ixI^L, ixO^L, res)
2963 integer, intent(in) :: ixI^L, ixO^L
2964 double precision, intent(in) :: w(ixI^S, 1:nw)
2965 double precision, intent(in) :: x(ixI^S, 1:ndim)
2966 double precision, intent(out):: res(ixI^S)
2967 res(ixo^s) = 1d0/rc * (gamma_1 * w(ixo^s, e_c_) + block%equi_vars(ixo^s,equi_pe_c0_,b0i)) /&
2968 (w(ixo^s,rho_c_) +block%equi_vars(ixo^s,equi_rho_c0_,b0i))
2969 end subroutine twofl_get_temperature_from_eint_c_with_equi
2970
2971! subroutine twofl_get_temperature_c_pert_from_tot(Te, ixI^L, ixO^L, res)
2972! use mod_global_parameters
2973! integer, intent(in) :: ixI^L, ixO^L
2974! double precision, intent(in) :: Te(ixI^S)
2975! double precision, intent(out):: res(ixI^S)
2976! res(ixO^S) = Te(ixO^S) -1d0/Rc * &
2977! block%equi_vars(ixO^S,equi_pe_c0_,0)/block%equi_vars(ixO^S,equi_rho_c0_,0)
2978! end subroutine twofl_get_temperature_c_pert_from_tot
2979
2980 subroutine twofl_get_temperature_c_equi(w,x, ixI^L, ixO^L, res)
2982 integer, intent(in) :: ixI^L, ixO^L
2983 double precision, intent(in) :: w(ixI^S, 1:nw)
2984 double precision, intent(in) :: x(ixI^S, 1:ndim)
2985 double precision, intent(out):: res(ixI^S)
2986 res(ixo^s) = 1d0/rc * &
2987 block%equi_vars(ixo^s,equi_pe_c0_,b0i)/block%equi_vars(ixo^s,equi_rho_c0_,b0i)
2988 end subroutine twofl_get_temperature_c_equi
2989
2990 subroutine twofl_get_rho_c_equi(w, x, ixI^L, ixO^L, res)
2992 integer, intent(in) :: ixI^L, ixO^L
2993 double precision, intent(in) :: w(ixI^S, 1:nw)
2994 double precision, intent(in) :: x(ixI^S, 1:ndim)
2995 double precision, intent(out):: res(ixI^S)
2996 res(ixo^s) = block%equi_vars(ixo^s,equi_rho_c0_,b0i)
2997 end subroutine twofl_get_rho_c_equi
2998
2999 subroutine twofl_get_pe_c_equi(w,x, ixI^L, ixO^L, res)
3001 integer, intent(in) :: ixI^L, ixO^L
3002 double precision, intent(in) :: w(ixI^S, 1:nw)
3003 double precision, intent(in) :: x(ixI^S, 1:ndim)
3004 double precision, intent(out):: res(ixI^S)
3005 res(ixo^s) = block%equi_vars(ixo^s,equi_pe_c0_,b0i)
3006 end subroutine twofl_get_pe_c_equi
3007
3008 !> Calculate temperature=p/rho when in e_ the total energy is stored
3009 !> this does not check the values of twofl_energy and twofl_internal_e,
3010 !> twofl_energy = .true. and twofl_internal_e = .false.
3011 !> also check small_values is avoided
3012 subroutine twofl_get_temperature_from_etot_c(w, x, ixI^L, ixO^L, res)
3014 integer, intent(in) :: ixI^L, ixO^L
3015 double precision, intent(in) :: w(ixI^S, 1:nw)
3016 double precision, intent(in) :: x(ixI^S, 1:ndim)
3017 double precision, intent(out):: res(ixI^S)
3018 res(ixo^s)=1d0/rc * (gamma_1*(w(ixo^s,e_c_)&
3019 - twofl_kin_en_c(w,ixi^l,ixo^l)&
3020 - twofl_mag_en(w,ixi^l,ixo^l)))/w(ixo^s,rho_c_)
3021 end subroutine twofl_get_temperature_from_etot_c
3022 subroutine twofl_get_temperature_from_eki_c(w, x, ixI^L, ixO^L, res)
3024 integer, intent(in) :: ixI^L, ixO^L
3025 double precision, intent(in) :: w(ixI^S, 1:nw)
3026 double precision, intent(in) :: x(ixI^S, 1:ndim)
3027 double precision, intent(out):: res(ixI^S)
3028 res(ixo^s)=1d0/rc * (gamma_1*(w(ixo^s,e_c_)&
3029 - twofl_kin_en_c(w,ixi^l,ixo^l)))/w(ixo^s,rho_c_)
3030 end subroutine twofl_get_temperature_from_eki_c
3031
3032 subroutine twofl_get_temperature_from_etot_c_with_equi(w, x, ixI^L, ixO^L, res)
3034 integer, intent(in) :: ixI^L, ixO^L
3035 double precision, intent(in) :: w(ixI^S, 1:nw)
3036 double precision, intent(in) :: x(ixI^S, 1:ndim)
3037 double precision, intent(out):: res(ixI^S)
3038 res(ixo^s)=1d0/rc * (gamma_1*(w(ixo^s,e_c_)&
3039 - twofl_kin_en_c(w,ixi^l,ixo^l)&
3040 - twofl_mag_en(w,ixi^l,ixo^l)) + block%equi_vars(ixo^s,equi_pe_c0_,b0i))&
3041 /(w(ixo^s,rho_c_) +block%equi_vars(ixo^s,equi_rho_c0_,b0i))
3042
3043 end subroutine twofl_get_temperature_from_etot_c_with_equi
3044
3045 subroutine twofl_get_temperature_from_eki_c_with_equi(w, x, ixI^L, ixO^L, res)
3047 integer, intent(in) :: ixI^L, ixO^L
3048 double precision, intent(in) :: w(ixI^S, 1:nw)
3049 double precision, intent(in) :: x(ixI^S, 1:ndim)
3050 double precision, intent(out):: res(ixI^S)
3051 res(ixo^s)=1d0/rc * (gamma_1*(w(ixo^s,e_c_)&
3052 - twofl_kin_en_c(w,ixi^l,ixo^l)) + block%equi_vars(ixo^s,equi_pe_c0_,b0i))&
3053 /(w(ixo^s,rho_c_) +block%equi_vars(ixo^s,equi_rho_c0_,b0i))
3054
3055 end subroutine twofl_get_temperature_from_eki_c_with_equi
3056
3057 subroutine twofl_get_csound2_adiab_n(w,x,ixI^L,ixO^L,csound2)
3059 integer, intent(in) :: ixI^L, ixO^L
3060 double precision, intent(in) :: w(ixI^S,nw)
3061 double precision, intent(in) :: x(ixI^S,1:ndim)
3062 double precision, intent(out) :: csound2(ixI^S)
3063 double precision :: rhon(ixI^S)
3064
3065 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3066 csound2(ixo^s)=twofl_gamma*twofl_adiab*rhon(ixo^s)**gamma_1
3067
3068 end subroutine twofl_get_csound2_adiab_n
3069
3070 subroutine twofl_get_csound2_n_from_conserved(w,x,ixI^L,ixO^L,csound2)
3072 integer, intent(in) :: ixI^L, ixO^L
3073 double precision, intent(in) :: w(ixI^S,nw)
3074 double precision, intent(in) :: x(ixI^S,1:ndim)
3075 double precision, intent(out) :: csound2(ixI^S)
3076 double precision :: rhon(ixI^S)
3077
3078 if(phys_energy) then
3079 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3080 call twofl_get_pthermal_n(w,x,ixi^l,ixo^l,csound2)
3081 csound2(ixo^s)=twofl_gamma*csound2(ixo^s)/rhon(ixo^s)
3082 else
3083 call twofl_get_csound2_adiab_n(w,x,ixi^l,ixo^l,csound2)
3084 endif
3085 end subroutine twofl_get_csound2_n_from_conserved
3086
3087 !! TO DELETE
3088 subroutine twofl_get_csound2_n_from_primitive(w,x,ixI^L,ixO^L,csound2)
3090 integer, intent(in) :: ixI^L, ixO^L
3091 double precision, intent(in) :: w(ixI^S,nw)
3092 double precision, intent(in) :: x(ixI^S,1:ndim)
3093 double precision, intent(out) :: csound2(ixI^S)
3094 double precision :: rhon(ixI^S)
3095
3096 if(phys_energy) then
3097 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3098 call twofl_get_pthermal_n_primitive(w,x,ixi^l,ixo^l,csound2)
3099 csound2(ixo^s)=twofl_gamma*csound2(ixo^s)/rhon(ixo^s)
3100 else
3101 call twofl_get_csound2_adiab_n(w,x,ixi^l,ixo^l,csound2)
3102 endif
3103 end subroutine twofl_get_csound2_n_from_primitive
3104
3105 subroutine twofl_get_csound2_adiab_c(w,x,ixI^L,ixO^L,csound2)
3107 integer, intent(in) :: ixI^L, ixO^L
3108 double precision, intent(in) :: w(ixI^S,nw)
3109 double precision, intent(in) :: x(ixI^S,1:ndim)
3110 double precision, intent(out) :: csound2(ixI^S)
3111 double precision :: rhoc(ixI^S)
3112
3113 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
3114 csound2(ixo^s)=twofl_gamma*twofl_adiab* rhoc(ixo^s)**gamma_1
3115
3116 end subroutine twofl_get_csound2_adiab_c
3117
3118 subroutine twofl_get_csound2_c_from_conserved(w,x,ixI^L,ixO^L,csound2)
3120 integer, intent(in) :: ixi^l, ixo^l
3121 double precision, intent(in) :: w(ixi^s,nw)
3122 double precision, intent(in) :: x(ixi^s,1:ndim)
3123 double precision, intent(out) :: csound2(ixi^s)
3124 double precision :: rhoc(ixi^s)
3125
3126 if(phys_energy) then
3127 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
3128 call twofl_get_pthermal_c(w,x,ixi^l,ixo^l,csound2)
3129 csound2(ixo^s)=twofl_gamma*csound2(ixo^s)/rhoc(ixo^s)
3130 else
3131 call twofl_get_csound2_adiab_c(w,x,ixi^l,ixo^l,csound2)
3132 endif
3134
3135 !> Calculate fluxes within ixO^L.
3136 subroutine twofl_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
3138 use mod_geometry
3139
3140 integer, intent(in) :: ixi^l, ixo^l, idim
3141 ! conservative w
3142 double precision, intent(in) :: wc(ixi^s,nw)
3143 ! primitive w
3144 double precision, intent(in) :: w(ixi^s,nw)
3145 double precision, intent(in) :: x(ixi^s,1:ndim)
3146 double precision,intent(out) :: f(ixi^s,nwflux)
3147
3148 double precision :: pgas(ixo^s), ptotal(ixo^s),tmp(ixi^s)
3149 double precision, allocatable:: vhall(:^d&,:)
3150 integer :: idirmin, iw, idir, jdir, kdir
3151
3152 ! value at the interfaces, idim = block%iw0 --> b0i
3153 ! reuse tmp, used afterwards
3154 ! value at the interface so we can't put momentum
3155 call get_rhoc_tot(w,x,ixi^l,ixo^l,tmp)
3156 ! Get flux of density
3157 f(ixo^s,rho_c_)=w(ixo^s,mom_c(idim))*tmp(ixo^s)
3158 ! pgas is time dependent only
3159 if(phys_energy) then
3160 pgas(ixo^s)=w(ixo^s,e_c_)
3161 else
3162 pgas(ixo^s)=twofl_adiab*tmp(ixo^s)**twofl_gamma
3163 if(has_equi_pe_c0) then
3164 pgas(ixo^s)=pgas(ixo^s)-block%equi_vars(ixo^s,equi_pe_c0_,b0i)
3165 end if
3166 end if
3167
3168 if (twofl_hall) then
3169 allocate(vhall(ixi^s,1:ndir))
3170 call twofl_getv_hall(w,x,ixi^l,ixo^l,vhall)
3171 end if
3172
3173 if(b0field) tmp(ixo^s)=sum(block%B0(ixo^s,:,idim)*w(ixo^s,mag(:)),dim=ndim+1)
3174
3175 ptotal(ixo^s) = pgas(ixo^s) + 0.5d0*sum(w(ixo^s, mag(:))**2, dim=ndim+1)
3176
3177 ! Get flux of momentum
3178 ! f_i[m_k]=v_i*m_k-b_k*b_i [+ptotal if i==k]
3179 do idir=1,ndir
3180 if(idim==idir) then
3181 f(ixo^s,mom_c(idir))=ptotal(ixo^s)-w(ixo^s,mag(idim))*w(ixo^s,mag(idir))
3182 if(b0field) f(ixo^s,mom_c(idir))=f(ixo^s,mom_c(idir))+tmp(ixo^s)
3183 else
3184 f(ixo^s,mom_c(idir))= -w(ixo^s,mag(idir))*w(ixo^s,mag(idim))
3185 end if
3186 if (b0field) then
3187 f(ixo^s,mom_c(idir))=f(ixo^s,mom_c(idir))&
3188 -w(ixo^s,mag(idir))*block%B0(ixo^s,idim,idim)&
3189 -w(ixo^s,mag(idim))*block%B0(ixo^s,idir,idim)
3190 end if
3191 f(ixo^s,mom_c(idir))=f(ixo^s,mom_c(idir))+w(ixo^s,mom_c(idim))*wc(ixo^s,mom_c(idir))
3192 end do
3193
3194 ! Get flux of energy
3195 ! f_i[e]=v_i*e+v_i*ptotal-b_i*(b_k*v_k)
3196 if(phys_energy) then
3197 if (phys_internal_e) then
3198 f(ixo^s,e_c_)=w(ixo^s,mom_c(idim))*wc(ixo^s,e_c_)
3199 else if(twofl_eq_energy == eq_energy_ki) then
3200
3201 f(ixo^s,e_c_)=w(ixo^s,mom_c(idim))*(wc(ixo^s,e_c_)+pgas(ixo^s))
3202 else
3203 f(ixo^s,e_c_)=w(ixo^s,mom_c(idim))*(wc(ixo^s,e_c_)+ptotal(ixo^s))&
3204 -w(ixo^s,mag(idim))*sum(w(ixo^s,mag(:))*w(ixo^s,mom_c(:)),dim=ndim+1)
3205
3206 if (b0field) then
3207 f(ixo^s,e_c_) = f(ixo^s,e_c_) &
3208 + w(ixo^s,mom_c(idim)) * tmp(ixo^s) &
3209 - sum(w(ixo^s,mom_c(:))*w(ixo^s,mag(:)),dim=ndim+1) * block%B0(ixo^s,idim,idim)
3210 end if
3211
3212 if (twofl_hall) then
3213 ! f_i[e]= f_i[e] + vHall_i*(b_k*b_k) - b_i*(vHall_k*b_k)
3214 if (twofl_etah>zero) then
3215 f(ixo^s,e_c_) = f(ixo^s,e_c_) + vhall(ixo^s,idim) * &
3216 sum(w(ixo^s, mag(:))**2,dim=ndim+1) &
3217 - w(ixo^s,mag(idim)) * sum(vhall(ixo^s,:)*w(ixo^s,mag(:)),dim=ndim+1)
3218 if (b0field) then
3219 f(ixo^s,e_c_) = f(ixo^s,e_c_) &
3220 + vhall(ixo^s,idim) * tmp(ixo^s) &
3221 - sum(vhall(ixo^s,:)*w(ixo^s,mag(:)),dim=ndim+1) * block%B0(ixo^s,idim,idim)
3222 end if
3223 end if
3224 end if
3225 end if !total_energy
3226 ! add flux of equilibrium internal energy corresponding to pe_c0
3227 if(has_equi_pe_c0) then
3228#if !defined(E_RM_W0) || E_RM_W0 == 1
3229 f(ixo^s,e_c_)= f(ixo^s,e_c_) &
3230 + w(ixo^s,mom_c(idim)) * block%equi_vars(ixo^s,equi_pe_c0_,idim) * inv_gamma_1
3231#else
3232 if(phys_internal_e) then
3233 f(ixo^s,e_c_)= f(ixo^s,e_c_) &
3234 + w(ixo^s,mom_c(idim)) * block%equi_vars(ixo^s,equi_pe_c0_,idim) * inv_gamma_1
3235 else
3236 f(ixo^s,e_c_)= f(ixo^s,e_c_) &
3237 + w(ixo^s,mom_c(idim)) * block%equi_vars(ixo^s,equi_pe_c0_,idim) * twofl_gamma * inv_gamma_1
3238 end if
3239#endif
3240 end if
3241 end if !phys_energy
3242
3243 ! compute flux of magnetic field
3244 ! f_i[b_k]=v_i*b_k-v_k*b_i
3245 do idir=1,ndir
3246 if (idim==idir) then
3247 ! f_i[b_i] should be exactly 0, so we do not use the transport flux
3248 if (twofl_glm) then
3249 f(ixo^s,mag(idir))=w(ixo^s,psi_)
3250 else
3251 f(ixo^s,mag(idir))=zero
3252 end if
3253 else
3254 f(ixo^s,mag(idir))=w(ixo^s,mom_c(idim))*w(ixo^s,mag(idir))-w(ixo^s,mag(idim))*w(ixo^s,mom_c(idir))
3255
3256 if (b0field) then
3257 f(ixo^s,mag(idir))=f(ixo^s,mag(idir))&
3258 +w(ixo^s,mom_c(idim))*block%B0(ixo^s,idir,idim)&
3259 -w(ixo^s,mom_c(idir))*block%B0(ixo^s,idim,idim)
3260 end if
3261
3262 if (twofl_hall) then
3263 ! f_i[b_k] = f_i[b_k] + vHall_i*b_k - vHall_k*b_i
3264 if (twofl_etah>zero) then
3265 if (b0field) then
3266 f(ixo^s,mag(idir)) = f(ixo^s,mag(idir)) &
3267 - vhall(ixo^s,idir)*(w(ixo^s,mag(idim))+block%B0(ixo^s,idim,idim)) &
3268 + vhall(ixo^s,idim)*(w(ixo^s,mag(idir))+block%B0(ixo^s,idir,idim))
3269 else
3270 f(ixo^s,mag(idir)) = f(ixo^s,mag(idir)) &
3271 - vhall(ixo^s,idir)*w(ixo^s,mag(idim)) &
3272 + vhall(ixo^s,idim)*w(ixo^s,mag(idir))
3273 end if
3274 end if
3275 end if
3276
3277 end if
3278 end do
3279
3280 if (twofl_glm) then
3281 !f_i[psi]=Ch^2*b_{i} Eq. 24e and Eq. 38c Dedner et al 2002 JCP, 175, 645
3282 f(ixo^s,psi_) = cmax_global**2*w(ixo^s,mag(idim))
3283 end if
3284
3285 if (twofl_hall) then
3286 deallocate(vhall)
3287 end if
3288
3289 !!neutrals
3290 call get_rhon_tot(w,x,ixi^l,ixo^l,tmp)
3291 f(ixo^s,rho_n_)=w(ixo^s,mom_n(idim))*tmp(ixo^s)
3292 if(phys_energy) then
3293 pgas(ixo^s) = w(ixo^s, e_n_)
3294 else
3295 pgas(ixo^s)=twofl_adiab*tmp(ixo^s)**twofl_gamma
3296 if(has_equi_pe_n0) then
3297 pgas(ixo^s)=pgas(ixo^s)-block%equi_vars(ixo^s,equi_pe_n0_,b0i)
3298 end if
3299 end if
3300 ! Momentum flux is v_i*m_i, +p in direction idim
3301 do idir = 1, ndir
3302 !if(idim==idir) then
3303 ! f(ixO^S,mom_c(idir)) = pgas(ixO^S)
3304 !else
3305 ! f(ixO^S,mom_c(idir)) = 0.0d0
3306 !end if
3307 !f(ixO^S,mom_c(idir))=f(ixO^S,mom_c(idir))+w(ixO^S,mom_c(idim))*wC(ixO^S,mom_c(idir))
3308 f(ixo^s, mom_n(idir)) = w(ixo^s,mom_n(idim)) * wc(ixo^s, mom_n(idir))
3309 end do
3310
3311 f(ixo^s, mom_n(idim)) = f(ixo^s, mom_n(idim)) + pgas(ixo^s)
3312
3313 if(phys_energy) then
3314 !reuse pgas for storing a in the term: div (u_n * a) and make multiplication at the end
3315 pgas(ixo^s) = wc(ixo^s,e_n_)
3316 if(.not. phys_internal_e) then
3317 ! add pressure perturbation
3318 pgas(ixo^s) = pgas(ixo^s) + w(ixo^s,e_n_)
3319 end if
3320 ! add flux of equilibrium internal energy corresponding to pe_n0
3321 if(has_equi_pe_n0) then
3322#if !defined(E_RM_W0) || E_RM_W0 == 1
3323 pgas(ixo^s) = pgas(ixo^s) + block%equi_vars(ixo^s,equi_pe_n0_,idim) * inv_gamma_1
3324#else
3325 pgas(ixo^s) = pgas(ixo^s) + block%equi_vars(ixo^s,equi_pe_n0_,idim) * twofl_gamma * inv_gamma_1
3326#endif
3327 end if
3328 ! add u_n * a in the flux
3329 f(ixo^s, e_n_) = w(ixo^s,mom_n(idim)) * pgas(ixo^s)
3330 end if
3331
3332 end subroutine twofl_get_flux
3333
3334 !> w[iws]=w[iws]+qdt*S[iws,wCT] where S is the source based on wCT within ixO
3335 subroutine twofl_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
3339 !use mod_gravity, only: gravity_add_source
3340
3341 integer, intent(in) :: ixi^l, ixo^l
3342 double precision, intent(in) :: qdt,dtfactor
3343 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw),x(ixi^s,1:ndim)
3344 double precision, intent(inout) :: w(ixi^s,1:nw)
3345 logical, intent(in) :: qsourcesplit
3346 logical, intent(inout) :: active
3347
3348 if (.not. qsourcesplit) then
3349 ! Source for solving internal energy
3350 if(phys_internal_e) then
3351 active = .true.
3352 call internal_energy_add_source_n(qdt,ixi^l,ixo^l,wct,w,x)
3353 call internal_energy_add_source_c(qdt,ixi^l,ixo^l,wct,w,x,e_c_)
3354 else
3355#if !defined(E_RM_W0) || E_RM_W0==1
3356 ! add -p0 div v source terms when equi are present
3357 if(has_equi_pe_n0) then
3358 active = .true.
3359 call add_pe_n0_divv(qdt,ixi^l,ixo^l,wct,w,x)
3360 endif
3361 if(has_equi_pe_c0) then
3362 active = .true.
3363 call add_pe_c0_divv(qdt,ixi^l,ixo^l,wct,w,x)
3364 endif
3365#endif
3366 if(twofl_eq_energy == eq_energy_ki) then
3367 active = .true.
3368 call add_source_lorentz_work(qdt,ixi^l,ixo^l,w,wct,x)
3369 endif
3370 endif
3371
3372 ! Source for B0 splitting
3373 if (b0field) then
3374 active = .true.
3375 call add_source_b0split(qdt,ixi^l,ixo^l,wct,w,x)
3376 end if
3377
3378 ! Sources for resistivity in eqs. for e, B1, B2 and B3
3379 if (abs(twofl_eta)>smalldouble)then
3380 active = .true.
3381 call add_source_res2(qdt,ixi^l,ixo^l,wct,w,x)
3382 end if
3383
3384 if (twofl_eta_hyper>0.d0)then
3385 active = .true.
3386 call add_source_hyperres(qdt,ixi^l,ixo^l,wct,w,x)
3387 end if
3388 !it is not added in a split manner
3389 if(.not. use_imex_scheme .and. has_collisions()) then
3390 active = .true.
3391 call twofl_explicit_coll_terms_update(qdt,ixi^l,ixo^l,w,wct,x)
3392 endif
3393
3394 if(twofl_hyperdiffusivity) then
3395 active = .true.
3396 call add_source_hyperdiffusive(qdt,ixi^l,ixo^l,w,wct,x)
3397 endif
3398
3399 end if
3400
3401 {^nooned
3402 if(source_split_divb .eqv. qsourcesplit) then
3403 ! Sources related to div B
3404 select case (type_divb)
3405 case (divb_none)
3406 ! Do nothing
3407 case (divb_glm)
3408 active = .true.
3409 call add_source_glm(qdt,ixi^l,ixo^l,wct,w,x)
3410 case (divb_powel)
3411 active = .true.
3412 call add_source_powel(qdt,ixi^l,ixo^l,wct,w,x)
3413 case (divb_janhunen)
3414 active = .true.
3415 call add_source_janhunen(qdt,ixi^l,ixo^l,wct,w,x)
3416 case (divb_linde)
3417 active = .true.
3418 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
3419 case (divb_lindejanhunen)
3420 active = .true.
3421 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
3422 call add_source_janhunen(qdt,ixi^l,ixo^l,wct,w,x)
3423 case (divb_lindepowel)
3424 active = .true.
3425 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
3426 call add_source_powel(qdt,ixi^l,ixo^l,wct,w,x)
3427 case (divb_lindeglm)
3428 active = .true.
3429 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
3430 call add_source_glm(qdt,ixi^l,ixo^l,wct,w,x)
3431 case (divb_ct)
3432 continue ! Do nothing
3433 case (divb_multigrid)
3434 continue ! Do nothing
3435 case default
3436 call mpistop('Unknown divB fix')
3437 end select
3438 end if
3439 }
3440
3442 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
3443 w,x,qsourcesplit,active,rc_fl_c)
3444 end if
3446 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
3447 w,x,qsourcesplit,active,rc_fl_n)
3448 end if
3449!
3450! if(twofl_viscosity) then
3451! call viscosity_add_source(qdt,ixI^L,ixO^L,wCT,wCTprim,&
3452! w,x,phys_energy,qsourcesplit,active)
3453! end if
3454!
3455 if(twofl_gravity) then
3456 call gravity_add_source(qdt,ixi^l,ixo^l,wct,&
3457 w,x,twofl_eq_energy .eq. eq_energy_ki .or. phys_total_energy,qsourcesplit,active)
3458 end if
3459
3460 end subroutine twofl_add_source
3461
3462 subroutine add_pe_n0_divv(qdt,ixI^L,ixO^L,wCT,w,x)
3464 use mod_geometry
3465
3466 integer, intent(in) :: ixi^l, ixo^l
3467 double precision, intent(in) :: qdt
3468 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3469 double precision, intent(inout) :: w(ixi^s,1:nw)
3470 double precision :: v(ixi^s,1:ndir)
3471
3472 call twofl_get_v_n(wct,x,ixi^l,ixi^l,v)
3473 call add_geom_pdivv(qdt,ixi^l,ixo^l,v,-block%equi_vars(ixi^s,equi_pe_n0_,0),w,x,e_n_)
3474
3475 end subroutine add_pe_n0_divv
3476
3477 subroutine add_pe_c0_divv(qdt,ixI^L,ixO^L,wCT,w,x)
3479 use mod_geometry
3480
3481 integer, intent(in) :: ixi^l, ixo^l
3482 double precision, intent(in) :: qdt
3483 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3484 double precision, intent(inout) :: w(ixi^s,1:nw)
3485 double precision :: v(ixi^s,1:ndir)
3486
3487 call twofl_get_v_c(wct,x,ixi^l,ixi^l,v)
3488 call add_geom_pdivv(qdt,ixi^l,ixo^l,v,-block%equi_vars(ixi^s,equi_pe_c0_,0),w,x,e_c_)
3489
3490 end subroutine add_pe_c0_divv
3491
3492 subroutine add_geom_pdivv(qdt,ixI^L,ixO^L,v,p,w,x,ind)
3494 use mod_geometry
3495 integer, intent(in) :: ixi^l, ixo^l,ind
3496 double precision, intent(in) :: qdt
3497 double precision, intent(in) :: p(ixi^s), v(ixi^s,1:ndir), x(ixi^s,1:ndim)
3498 double precision, intent(inout) :: w(ixi^s,1:nw)
3499 double precision :: divv(ixi^s)
3500
3501 if(slab_uniform) then
3502 if(nghostcells .gt. 2) then
3503 call divvector(v,ixi^l,ixo^l,divv,3)
3504 else
3505 call divvector(v,ixi^l,ixo^l,divv,2)
3506 end if
3507 else
3508 call divvector(v,ixi^l,ixo^l,divv)
3509 end if
3510 w(ixo^s,ind)=w(ixo^s,ind)+qdt*p(ixo^s)*divv(ixo^s)
3511 end subroutine add_geom_pdivv
3512
3513 !> Compute the Lorentz force (JxB)
3514 subroutine get_lorentz(ixI^L,ixO^L,w,JxB)
3516 integer, intent(in) :: ixi^l, ixo^l
3517 double precision, intent(in) :: w(ixi^s,1:nw)
3518 double precision, intent(inout) :: jxb(ixi^s,3)
3519 double precision :: a(ixi^s,3), b(ixi^s,3), tmp(ixi^s,3)
3520 integer :: idir, idirmin
3521 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
3522 double precision :: current(ixi^s,7-2*ndir:3)
3523
3524 b=0.0d0
3525 do idir = 1, ndir
3526 b(ixo^s, idir) = twofl_mag_i_all(w, ixi^l, ixo^l,idir)
3527 end do
3528
3529 ! store J current in a
3530 call get_current(w,ixi^l,ixo^l,idirmin,current)
3531
3532 a=0.0d0
3533 do idir=7-2*ndir,3
3534 a(ixo^s,idir)=current(ixo^s,idir)
3535 end do
3536
3537 call cross_product(ixi^l,ixo^l,a,b,jxb)
3538 end subroutine get_lorentz
3539
3540 subroutine add_source_lorentz_work(qdt,ixI^L,ixO^L,w,wCT,x)
3542 integer, intent(in) :: ixi^l, ixo^l
3543 double precision, intent(in) :: qdt
3544 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3545 double precision, intent(inout) :: w(ixi^s,1:nw)
3546 double precision :: a(ixi^s,3), b(ixi^s,1:ndir)
3547
3548 call get_lorentz(ixi^l, ixo^l,wct,a)
3549 call twofl_get_v_c(wct,x,ixi^l,ixo^l,b)
3550 w(ixo^s,e_c_)=w(ixo^s,e_c_)+qdt*sum(a(ixo^s,1:ndir)*b(ixo^s,1:ndir),dim=ndim+1)
3551
3552 end subroutine add_source_lorentz_work
3553
3554 !> Calculate v_n vector
3555 subroutine twofl_get_v_n(w,x,ixI^L,ixO^L,v)
3557
3558 integer, intent(in) :: ixi^l, ixo^l
3559 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
3560 double precision, intent(out) :: v(ixi^s,ndir)
3561 double precision :: rhon(ixi^s)
3562 integer :: idir
3563
3564 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3565
3566 do idir=1,ndir
3567 v(ixo^s,idir) = w(ixo^s, mom_n(idir)) / rhon(ixo^s)
3568 end do
3569
3570 end subroutine twofl_get_v_n
3571
3572 subroutine get_rhon_tot(w,x,ixI^L,ixO^L,rhon)
3574 integer, intent(in) :: ixi^l, ixo^l
3575 double precision, intent(in) :: w(ixi^s,1:nw), x(ixi^s,1:ndim)
3576 double precision, intent(out) :: rhon(ixi^s)
3577 if(has_equi_rho_n0) then
3578 rhon(ixo^s) = w(ixo^s,rho_n_) + block%equi_vars(ixo^s,equi_rho_n0_,b0i)
3579 else
3580 rhon(ixo^s) = w(ixo^s,rho_n_)
3581 endif
3582
3583 end subroutine get_rhon_tot
3584
3585 subroutine twofl_get_pthermal_n(w,x,ixI^L,ixO^L,pth)
3588 integer, intent(in) :: ixi^l, ixo^l
3589 double precision, intent(in) :: w(ixi^s,1:nw)
3590 double precision, intent(in) :: x(ixi^s,1:ndim)
3591 double precision, intent(out) :: pth(ixi^s)
3592
3593 integer :: ix^d, iw
3594
3595 if(phys_energy) then
3596 if(phys_internal_e) then
3597 pth(ixo^s)=gamma_1*w(ixo^s,e_n_)
3598 else
3599 pth(ixo^s)=gamma_1*(w(ixo^s,e_n_)&
3600 - twofl_kin_en_n(w,ixi^l,ixo^l))
3601 end if
3602 if(has_equi_pe_n0) then
3603 pth(ixo^s) = pth(ixo^s) + block%equi_vars(ixo^s,equi_pe_n0_,b0i)
3604 endif
3605 else
3606 call get_rhon_tot(w,x,ixi^l,ixo^l,pth)
3607 pth(ixo^s)=twofl_adiab*pth(ixo^s)**twofl_gamma
3608 end if
3609
3610 if (fix_small_values) then
3611 {do ix^db= ixo^lim^db\}
3612 if(pth(ix^d)<small_pressure) then
3613 pth(ix^d)=small_pressure
3614 end if
3615 {enddo^d&\}
3616 else if (check_small_values) then
3617 {do ix^db= ixo^lim^db\}
3618 if(pth(ix^d)<small_pressure) then
3619 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
3620 " encountered when call twofl_get_pthermal_n"
3621 write(*,*) "Iteration: ", it, " Time: ", global_time
3622 write(*,*) "Location: ", x(ix^d,:)
3623 write(*,*) "Cell number: ", ix^d
3624 do iw=1,nw
3625 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
3626 end do
3627 ! use erroneous arithmetic operation to crash the run
3628 if(trace_small_values) write(*,*) sqrt(pth(ix^d)-bigdouble)
3629 write(*,*) "Saving status at the previous time step"
3630 crash=.true.
3631 end if
3632 {enddo^d&\}
3633 end if
3634
3635 end subroutine twofl_get_pthermal_n
3636
3637 subroutine twofl_get_pthermal_n_primitive(w,x,ixI^L,ixO^L,pth)
3639 integer, intent(in) :: ixi^l, ixo^l
3640 double precision, intent(in) :: w(ixi^s,1:nw)
3641 double precision, intent(in) :: x(ixi^s,1:ndim)
3642 double precision, intent(out) :: pth(ixi^s)
3643
3644 if(phys_energy) then
3645 if(has_equi_pe_n0) then
3646 pth(ixo^s) = w(ixo^s,e_n_) + block%equi_vars(ixo^s,equi_pe_n0_,b0i)
3647 else
3648 pth(ixo^s) = w(ixo^s,e_n_)
3649 endif
3650 else
3651 call get_rhon_tot(w,x,ixi^l,ixo^l,pth)
3652 pth(ixo^s)=twofl_adiab*pth(ixo^s)**twofl_gamma
3653 end if
3654 end subroutine twofl_get_pthermal_n_primitive
3655
3656 !> Calculate v component
3657 subroutine twofl_get_v_n_idim(w,x,ixI^L,ixO^L,idim,v)
3659
3660 integer, intent(in) :: ixi^l, ixo^l, idim
3661 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
3662 double precision, intent(out) :: v(ixi^s)
3663 double precision :: rhon(ixi^s)
3664
3665 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3666 v(ixo^s) = w(ixo^s, mom_n(idim)) / rhon(ixo^s)
3667
3668 end subroutine twofl_get_v_n_idim
3669
3670 subroutine internal_energy_add_source_n(qdt,ixI^L,ixO^L,wCT,w,x)
3672 use mod_geometry
3673
3674 integer, intent(in) :: ixi^l, ixo^l
3675 double precision, intent(in) :: qdt
3676 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3677 double precision, intent(inout) :: w(ixi^s,1:nw)
3678 double precision :: pth(ixi^s),v(ixi^s,1:ndir),divv(ixi^s)
3679
3680 call twofl_get_pthermal_n(wct,x,ixi^l,ixo^l,pth)
3681 call twofl_get_v_n(wct,x,ixi^l,ixi^l,v)
3682 call add_geom_pdivv(qdt,ixi^l,ixo^l,v,-pth,w,x,e_n_)
3683
3684 if(fix_small_values .and. .not. has_equi_pe_n0) then
3685 call twofl_handle_small_ei_n(w,x,ixi^l,ixo^l,e_n_,'internal_energy_add_source')
3686 end if
3687 end subroutine internal_energy_add_source_n
3688
3689 !> Calculate v_c vector
3690 subroutine twofl_get_v_c(w,x,ixI^L,ixO^L,v)
3692
3693 integer, intent(in) :: ixi^l, ixo^l
3694 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
3695 double precision, intent(out) :: v(ixi^s,ndir)
3696 double precision :: rhoc(ixi^s)
3697 integer :: idir
3698
3699 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
3700 do idir=1,ndir
3701 v(ixo^s,idir) = w(ixo^s, mom_c(idir)) / rhoc(ixo^s)
3702 end do
3703
3704 end subroutine twofl_get_v_c
3705
3706 subroutine get_rhoc_tot(w,x,ixI^L,ixO^L,rhoc)
3708 integer, intent(in) :: ixi^l, ixo^l
3709 double precision, intent(in) :: w(ixi^s,1:nw), x(ixi^s,1:ndim)
3710 double precision, intent(out) :: rhoc(ixi^s)
3711 if(has_equi_rho_c0) then
3712 rhoc(ixo^s) = w(ixo^s,rho_c_) + block%equi_vars(ixo^s,equi_rho_c0_,b0i)
3713 else
3714 rhoc(ixo^s) = w(ixo^s,rho_c_)
3715 endif
3716
3717 end subroutine get_rhoc_tot
3718
3719 subroutine twofl_get_pthermal_c(w,x,ixI^L,ixO^L,pth)
3722 integer, intent(in) :: ixi^l, ixo^l
3723 double precision, intent(in) :: w(ixi^s,1:nw)
3724 double precision, intent(in) :: x(ixi^s,1:ndim)
3725 double precision, intent(out) :: pth(ixi^s)
3726 integer :: ix^d, iw
3727
3728 if(phys_energy) then
3729 if(phys_internal_e) then
3730 pth(ixo^s)=gamma_1*w(ixo^s,e_c_)
3731 elseif(phys_total_energy) then
3732 pth(ixo^s)=gamma_1*(w(ixo^s,e_c_)&
3733 - twofl_kin_en_c(w,ixi^l,ixo^l)&
3734 - twofl_mag_en(w,ixi^l,ixo^l))
3735 else
3736 pth(ixo^s)=gamma_1*(w(ixo^s,e_c_)&
3737 - twofl_kin_en_c(w,ixi^l,ixo^l))
3738 end if
3739 if(has_equi_pe_c0) then
3740 pth(ixo^s) = pth(ixo^s) + block%equi_vars(ixo^s,equi_pe_c0_,b0i)
3741 endif
3742 else
3743 call get_rhoc_tot(w,x,ixi^l,ixo^l,pth)
3744 pth(ixo^s)=twofl_adiab*pth(ixo^s)**twofl_gamma
3745 end if
3746
3747 if (fix_small_values) then
3748 {do ix^db= ixo^lim^db\}
3749 if(pth(ix^d)<small_pressure) then
3750 pth(ix^d)=small_pressure
3751 end if
3752 {enddo^d&\}
3753 else if (check_small_values) then
3754 {do ix^db= ixo^lim^db\}
3755 if(pth(ix^d)<small_pressure) then
3756 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
3757 " encountered when call twofl_get_pe_c1"
3758 write(*,*) "Iteration: ", it, " Time: ", global_time
3759 write(*,*) "Location: ", x(ix^d,:)
3760 write(*,*) "Cell number: ", ix^d
3761 do iw=1,nw
3762 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
3763 end do
3764 ! use erroneous arithmetic operation to crash the run
3765 if(trace_small_values) write(*,*) sqrt(pth(ix^d)-bigdouble)
3766 write(*,*) "Saving status at the previous time step"
3767 crash=.true.
3768 end if
3769 {enddo^d&\}
3770 end if
3771
3772 end subroutine twofl_get_pthermal_c
3773
3774 subroutine twofl_get_pthermal_c_primitive(w,x,ixI^L,ixO^L,pth)
3776 integer, intent(in) :: ixi^l, ixo^l
3777 double precision, intent(in) :: w(ixi^s,1:nw)
3778 double precision, intent(in) :: x(ixi^s,1:ndim)
3779 double precision, intent(out) :: pth(ixi^s)
3780
3781 if(phys_energy) then
3782 if(has_equi_pe_c0) then
3783 pth(ixo^s) = w(ixo^s,e_c_) + block%equi_vars(ixo^s,equi_pe_c0_,b0i)
3784 else
3785 pth(ixo^s) = w(ixo^s,e_c_)
3786 endif
3787 else
3788 call get_rhoc_tot(w,x,ixi^l,ixo^l,pth)
3789 pth(ixo^s)=twofl_adiab*pth(ixo^s)**twofl_gamma
3790 end if
3791 end subroutine twofl_get_pthermal_c_primitive
3792
3793 !> Calculate v_c component
3794 subroutine twofl_get_v_c_idim(w,x,ixI^L,ixO^L,idim,v)
3796
3797 integer, intent(in) :: ixi^l, ixo^l, idim
3798 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
3799 double precision, intent(out) :: v(ixi^s)
3800 double precision :: rhoc(ixi^s)
3801
3802 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
3803 v(ixo^s) = w(ixo^s, mom_c(idim)) / rhoc(ixo^s)
3804
3805 end subroutine twofl_get_v_c_idim
3806
3807 subroutine internal_energy_add_source_c(qdt,ixI^L,ixO^L,wCT,w,x,ie)
3809 use mod_geometry
3810
3811 integer, intent(in) :: ixi^l, ixo^l,ie
3812 double precision, intent(in) :: qdt
3813 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3814 double precision, intent(inout) :: w(ixi^s,1:nw)
3815 double precision :: pth(ixi^s),v(ixi^s,1:ndir),divv(ixi^s)
3816
3817 call twofl_get_pthermal_c(wct,x,ixi^l,ixo^l,pth)
3818 call twofl_get_v_c(wct,x,ixi^l,ixi^l,v)
3819 call add_geom_pdivv(qdt,ixi^l,ixo^l,v,-pth,w,x,ie)
3820 if(fix_small_values .and. .not. has_equi_pe_c0) then
3821 call twofl_handle_small_ei_c(w,x,ixi^l,ixo^l,ie,'internal_energy_add_source')
3822 end if
3823 end subroutine internal_energy_add_source_c
3824
3825 !> handle small or negative internal energy
3826 subroutine twofl_handle_small_ei_c(w, x, ixI^L, ixO^L, ie, subname)
3829 integer, intent(in) :: ixi^l,ixo^l, ie
3830 double precision, intent(inout) :: w(ixi^s,1:nw)
3831 double precision, intent(in) :: x(ixi^s,1:ndim)
3832 character(len=*), intent(in) :: subname
3833
3834 integer :: idir
3835 logical :: flag(ixi^s,1:nw)
3836 double precision :: rhoc(ixi^s)
3837 double precision :: rhon(ixi^s)
3838
3839 flag=.false.
3840 if(has_equi_pe_c0) then
3841 where(w(ixo^s,ie)+block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1<small_e)&
3842 flag(ixo^s,ie)=.true.
3843 else
3844 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
3845 endif
3846 if(any(flag(ixo^s,ie))) then
3847 select case (small_values_method)
3848 case ("replace")
3849 if(has_equi_pe_c0) then
3850 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e - &
3851 block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
3852 else
3853 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
3854 endif
3855 case ("average")
3856 call small_values_average(ixi^l, ixo^l, w, x, flag, ie)
3857 case default
3858 ! small values error shows primitive variables
3859 ! to_primitive subroutine cannot be used as this error handling
3860 ! is also used in TC where e_to_ei is explicitly called
3861 w(ixo^s,e_n_)=w(ixo^s,e_n_)*gamma_1
3862 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3863 w(ixo^s,e_c_)=w(ixo^s,e_c_)*gamma_1
3864 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
3865 do idir = 1, ndir
3866 w(ixo^s, mom_n(idir)) = w(ixo^s, mom_n(idir))/rhon(ixo^s)
3867 w(ixo^s, mom_c(idir)) = w(ixo^s, mom_c(idir))/rhoc(ixo^s)
3868 end do
3869 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
3870 end select
3871 end if
3872
3873 end subroutine twofl_handle_small_ei_c
3874
3875 !> handle small or negative internal energy
3876 subroutine twofl_handle_small_ei_n(w, x, ixI^L, ixO^L, ie, subname)
3879 integer, intent(in) :: ixi^l,ixo^l, ie
3880 double precision, intent(inout) :: w(ixi^s,1:nw)
3881 double precision, intent(in) :: x(ixi^s,1:ndim)
3882 character(len=*), intent(in) :: subname
3883
3884 integer :: idir
3885 logical :: flag(ixi^s,1:nw)
3886 double precision :: rhoc(ixi^s)
3887 double precision :: rhon(ixi^s)
3888
3889 flag=.false.
3890 if(has_equi_pe_n0) then
3891 where(w(ixo^s,ie)+block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1<small_e)&
3892 flag(ixo^s,ie)=.true.
3893 else
3894 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
3895 endif
3896 if(any(flag(ixo^s,ie))) then
3897 select case (small_values_method)
3898 case ("replace")
3899 if(has_equi_pe_n0) then
3900 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e - &
3901 block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
3902 else
3903 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
3904 endif
3905 case ("average")
3906 call small_values_average(ixi^l, ixo^l, w, x, flag, ie)
3907 case default
3908 ! small values error shows primitive variables
3909 w(ixo^s,e_n_)=w(ixo^s,e_n_)*gamma_1
3910 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
3911 w(ixo^s,e_c_)=w(ixo^s,e_c_)*gamma_1
3912 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
3913 do idir = 1, ndir
3914 w(ixo^s, mom_n(idir)) = w(ixo^s, mom_n(idir))/rhon(ixo^s)
3915 w(ixo^s, mom_c(idir)) = w(ixo^s, mom_c(idir))/rhoc(ixo^s)
3916 end do
3917 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
3918 end select
3919 end if
3920
3921 end subroutine twofl_handle_small_ei_n
3922
3923 !> Source terms after split off time-independent magnetic field
3924 subroutine add_source_b0split(qdt,ixI^L,ixO^L,wCT,w,x)
3926
3927 integer, intent(in) :: ixi^l, ixo^l
3928 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3929 double precision, intent(inout) :: w(ixi^s,1:nw)
3930
3931 double precision :: a(ixi^s,3), b(ixi^s,3), axb(ixi^s,3)
3932 integer :: idir
3933
3934 a=0.d0
3935 b=0.d0
3936 ! for force-free field J0xB0 =0
3937 if(.not.b0field_forcefree) then
3938 ! store B0 magnetic field in b
3939 b(ixo^s,1:ndir)=block%B0(ixo^s,1:ndir,0)
3940
3941 ! store J0 current in a
3942 do idir=7-2*ndir,3
3943 a(ixo^s,idir)=block%J0(ixo^s,idir)
3944 end do
3945 call cross_product(ixi^l,ixo^l,a,b,axb)
3946 axb(ixo^s,:)=axb(ixo^s,:)*qdt
3947 ! add J0xB0 source term in momentum equations
3948 w(ixo^s,mom_c(1:ndir))=w(ixo^s,mom_c(1:ndir))+axb(ixo^s,1:ndir)
3949 end if
3950
3951 if(phys_total_energy) then
3952 a=0.d0
3953 ! for free-free field -(vxB0) dot J0 =0
3954 b(ixo^s,:)=wct(ixo^s,mag(:))
3955 ! store full magnetic field B0+B1 in b
3956 if(.not.b0field_forcefree) b(ixo^s,:)=b(ixo^s,:)+block%B0(ixo^s,:,0)
3957 ! store velocity in a
3958 do idir=1,ndir
3959 call twofl_get_v_c_idim(wct,x,ixi^l,ixo^l,idir,a(ixi^s,idir))
3960 end do
3961 call cross_product(ixi^l,ixo^l,a,b,axb)
3962 axb(ixo^s,:)=axb(ixo^s,:)*qdt
3963 ! add -(vxB) dot J0 source term in energy equation
3964 do idir=7-2*ndir,3
3965 w(ixo^s,e_c_)=w(ixo^s,e_c_)-axb(ixo^s,idir)*block%J0(ixo^s,idir)
3966 end do
3967 end if
3968
3969 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_B0')
3970
3971 end subroutine add_source_b0split
3972
3973 !> Add resistive source to w within ixO Uses 3 point stencil (1 neighbour) in
3974 !> each direction, non-conservative. If the fourthorder precompiler flag is
3975 !> set, uses fourth order central difference for the laplacian. Then the
3976 !> stencil is 5 (2 neighbours).
3977 subroutine add_source_res1(qdt,ixI^L,ixO^L,wCT,w,x)
3979 use mod_usr_methods
3980 use mod_geometry
3981
3982 integer, intent(in) :: ixi^l, ixo^l
3983 double precision, intent(in) :: qdt
3984 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
3985 double precision, intent(inout) :: w(ixi^s,1:nw)
3986 integer :: ixa^l,idir,jdir,kdir,idirmin,idim,jxo^l,hxo^l,ix
3987 integer :: lxo^l, kxo^l
3988
3989 double precision :: tmp(ixi^s),tmp2(ixi^s)
3990
3991 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
3992 double precision :: current(ixi^s,7-2*ndir:3),eta(ixi^s)
3993 double precision :: gradeta(ixi^s,1:ndim), bf(ixi^s,1:ndir)
3994
3995 ! Calculating resistive sources involve one extra layer
3996 if (twofl_4th_order) then
3997 ixa^l=ixo^l^ladd2;
3998 else
3999 ixa^l=ixo^l^ladd1;
4000 end if
4001
4002 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
4003 call mpistop("Error in add_source_res1: Non-conforming input limits")
4004
4005 ! Calculate current density and idirmin
4006 call get_current(wct,ixi^l,ixo^l,idirmin,current)
4007
4008 if (twofl_eta>zero)then
4009 eta(ixa^s)=twofl_eta
4010 gradeta(ixo^s,1:ndim)=zero
4011 else
4012 call usr_special_resistivity(wct,ixi^l,ixa^l,idirmin,x,current,eta)
4013 ! assumes that eta is not function of current?
4014 do idim=1,ndim
4015 call gradient(eta,ixi^l,ixo^l,idim,tmp)
4016 gradeta(ixo^s,idim)=tmp(ixo^s)
4017 end do
4018 end if
4019
4020 if(b0field) then
4021 bf(ixi^s,1:ndir)=wct(ixi^s,mag(1:ndir))+block%B0(ixi^s,1:ndir,0)
4022 else
4023 bf(ixi^s,1:ndir)=wct(ixi^s,mag(1:ndir))
4024 end if
4025
4026 do idir=1,ndir
4027 ! Put B_idir into tmp2 and eta*Laplace B_idir into tmp
4028 if (twofl_4th_order) then
4029 tmp(ixo^s)=zero
4030 tmp2(ixi^s)=bf(ixi^s,idir)
4031 do idim=1,ndim
4032 lxo^l=ixo^l+2*kr(idim,^d);
4033 jxo^l=ixo^l+kr(idim,^d);
4034 hxo^l=ixo^l-kr(idim,^d);
4035 kxo^l=ixo^l-2*kr(idim,^d);
4036 tmp(ixo^s)=tmp(ixo^s)+&
4037 (-tmp2(lxo^s)+16.0d0*tmp2(jxo^s)-30.0d0*tmp2(ixo^s)+16.0d0*tmp2(hxo^s)-tmp2(kxo^s)) &
4038 /(12.0d0 * dxlevel(idim)**2)
4039 end do
4040 else
4041 tmp(ixo^s)=zero
4042 tmp2(ixi^s)=bf(ixi^s,idir)
4043 do idim=1,ndim
4044 jxo^l=ixo^l+kr(idim,^d);
4045 hxo^l=ixo^l-kr(idim,^d);
4046 tmp(ixo^s)=tmp(ixo^s)+&
4047 (tmp2(jxo^s)-2.0d0*tmp2(ixo^s)+tmp2(hxo^s))/dxlevel(idim)**2
4048 end do
4049 end if
4050
4051 ! Multiply by eta
4052 tmp(ixo^s)=tmp(ixo^s)*eta(ixo^s)
4053
4054 ! Subtract grad(eta) x J = eps_ijk d_j eta J_k if eta is non-constant
4055 if (twofl_eta<zero)then
4056 do jdir=1,ndim; do kdir=idirmin,3
4057 if (lvc(idir,jdir,kdir)/=0)then
4058 if (lvc(idir,jdir,kdir)==1)then
4059 tmp(ixo^s)=tmp(ixo^s)-gradeta(ixo^s,jdir)*current(ixo^s,kdir)
4060 else
4061 tmp(ixo^s)=tmp(ixo^s)+gradeta(ixo^s,jdir)*current(ixo^s,kdir)
4062 end if
4063 end if
4064 end do; end do
4065 end if
4066
4067 ! Add sources related to eta*laplB-grad(eta) x J to B and e
4068 w(ixo^s,mag(idir))=w(ixo^s,mag(idir))+qdt*tmp(ixo^s)
4069 if (phys_total_energy) then
4070 w(ixo^s,e_c_)=w(ixo^s,e_c_)+qdt*tmp(ixo^s)*bf(ixo^s,idir)
4071 end if
4072 end do ! idir
4073
4074 if (phys_energy) then
4075 ! de/dt+=eta*J**2
4076 w(ixo^s,e_c_)=w(ixo^s,e_c_)+qdt*eta(ixo^s)*sum(current(ixo^s,:)**2,dim=ndim+1)
4077 end if
4078
4079 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_res1')
4080
4081 end subroutine add_source_res1
4082
4083 !> Add resistive source to w within ixO
4084 !> Uses 5 point stencil (2 neighbours) in each direction, conservative
4085 subroutine add_source_res2(qdt,ixI^L,ixO^L,wCT,w,x)
4087 use mod_usr_methods
4088 use mod_geometry
4089
4090 integer, intent(in) :: ixi^l, ixo^l
4091 double precision, intent(in) :: qdt
4092 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
4093 double precision, intent(inout) :: w(ixi^s,1:nw)
4094
4095 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
4096 double precision :: current(ixi^s,7-2*ndir:3),eta(ixi^s),curlj(ixi^s,1:3)
4097 double precision :: tmpvec(ixi^s,1:3),tmp(ixo^s)
4098 integer :: ixa^l,idir,idirmin,idirmin1
4099
4100 ixa^l=ixo^l^ladd2;
4101
4102 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
4103 call mpistop("Error in add_source_res2: Non-conforming input limits")
4104
4105 ixa^l=ixo^l^ladd1;
4106 ! Calculate current density within ixL: J=curl B, thus J_i=eps_ijk*d_j B_k
4107 ! Determine exact value of idirmin while doing the loop.
4108 call get_current(wct,ixi^l,ixa^l,idirmin,current)
4109
4110 if (twofl_eta>zero)then
4111 eta(ixa^s)=twofl_eta
4112 else
4113 call usr_special_resistivity(wct,ixi^l,ixa^l,idirmin,x,current,eta)
4114 end if
4115
4116 ! dB/dt= -curl(J*eta), thus B_i=B_i-eps_ijk d_j Jeta_k
4117 tmpvec(ixa^s,1:ndir)=zero
4118 do idir=idirmin,3
4119 tmpvec(ixa^s,idir)=current(ixa^s,idir)*eta(ixa^s)
4120 end do
4121 curlj=0.d0
4122 call curlvector(tmpvec,ixi^l,ixo^l,curlj,idirmin1,1,3)
4123 if(stagger_grid.and.ndim==2.and.ndir==3) then
4124 ! if 2.5D
4125 w(ixo^s,mag(ndir)) = w(ixo^s,mag(ndir))-qdt*curlj(ixo^s,ndir)
4126 else
4127 w(ixo^s,mag(1:ndir)) = w(ixo^s,mag(1:ndir))-qdt*curlj(ixo^s,1:ndir)
4128 end if
4129
4130 if(phys_energy) then
4131 if(phys_total_energy) then
4132 ! de/dt= +div(B x Jeta) = eta J^2 - B dot curl(eta J)
4133 ! de1/dt= eta J^2 - B1 dot curl(eta J)
4134 w(ixo^s,e_c_)=w(ixo^s,e_c_)+qdt*(eta(ixo^s)*sum(current(ixo^s,:)**2,dim=ndim+1)-&
4135 sum(wct(ixo^s,mag(1:ndir))*curlj(ixo^s,1:ndir),dim=ndim+1))
4136 else
4137 ! add eta*J**2 source term in the internal energy equation
4138 w(ixo^s,e_c_)=w(ixo^s,e_c_)+qdt*eta(ixo^s)*sum(current(ixo^s,:)**2,dim=ndim+1)
4139 end if
4140
4141 end if
4142
4143 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_res2')
4144 end subroutine add_source_res2
4145
4146 !> Add Hyper-resistive source to w within ixO
4147 !> Uses 9 point stencil (4 neighbours) in each direction.
4148 subroutine add_source_hyperres(qdt,ixI^L,ixO^L,wCT,w,x)
4150 use mod_geometry
4151
4152 integer, intent(in) :: ixi^l, ixo^l
4153 double precision, intent(in) :: qdt
4154 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
4155 double precision, intent(inout) :: w(ixi^s,1:nw)
4156 !.. local ..
4157 double precision :: current(ixi^s,7-2*ndir:3)
4158 double precision :: tmpvec(ixi^s,1:3),tmpvec2(ixi^s,1:3),tmp(ixi^s),ehyper(ixi^s,1:3)
4159 integer :: ixa^l,idir,jdir,kdir,idirmin,idirmin1
4160
4161 ixa^l=ixo^l^ladd3;
4162 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
4163 call mpistop("Error in add_source_hyperres: Non-conforming input limits")
4164
4165 call get_current(wct,ixi^l,ixa^l,idirmin,current)
4166 tmpvec(ixa^s,1:ndir)=zero
4167 do jdir=idirmin,3
4168 tmpvec(ixa^s,jdir)=current(ixa^s,jdir)
4169 end do
4170
4171 ixa^l=ixo^l^ladd2;
4172 call curlvector(tmpvec,ixi^l,ixa^l,tmpvec2,idirmin1,1,3)
4173
4174 ixa^l=ixo^l^ladd1;
4175 tmpvec(ixa^s,1:ndir)=zero
4176 call curlvector(tmpvec2,ixi^l,ixa^l,tmpvec,idirmin1,1,3)
4177 ehyper(ixa^s,1:ndir) = - tmpvec(ixa^s,1:ndir)*twofl_eta_hyper
4178
4179 ixa^l=ixo^l;
4180 tmpvec2(ixa^s,1:ndir)=zero
4181 call curlvector(ehyper,ixi^l,ixa^l,tmpvec2,idirmin1,1,3)
4182
4183 do idir=1,ndir
4184 w(ixo^s,mag(idir)) = w(ixo^s,mag(idir))-tmpvec2(ixo^s,idir)*qdt
4185 end do
4186
4187 if (phys_energy) then
4188 ! de/dt= +div(B x Ehyper)
4189 ixa^l=ixo^l^ladd1;
4190 tmpvec2(ixa^s,1:ndir)=zero
4191 do idir=1,ndir; do jdir=1,ndir; do kdir=idirmin,3
4192 tmpvec2(ixa^s,idir) = tmpvec(ixa^s,idir)&
4193 + lvc(idir,jdir,kdir)*wct(ixa^s,mag(jdir))*ehyper(ixa^s,kdir)
4194 end do; end do; end do
4195 tmp(ixo^s)=zero
4196 call divvector(tmpvec2,ixi^l,ixo^l,tmp)
4197 w(ixo^s,e_c_)=w(ixo^s,e_c_)+tmp(ixo^s)*qdt
4198 end if
4199
4200 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_hyperres')
4201
4202 end subroutine add_source_hyperres
4203
4204 subroutine add_source_glm(qdt,ixI^L,ixO^L,wCT,w,x)
4205 ! Add divB related sources to w within ixO
4206 ! corresponding to Dedner JCP 2002, 175, 645 _equation 24_
4207 ! giving the EGLM-MHD scheme
4209 use mod_geometry
4210
4211 integer, intent(in) :: ixi^l, ixo^l
4212 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
4213 double precision, intent(inout) :: w(ixi^s,1:nw)
4214 double precision:: divb(ixi^s)
4215 integer :: idim,idir
4216 double precision :: gradpsi(ixi^s)
4217
4218 ! We calculate now div B
4219 call get_divb(wct,ixi^l,ixo^l,divb, twofl_divb_nth)
4220
4221 ! dPsi/dt = - Ch^2/Cp^2 Psi
4222 if (twofl_glm_alpha < zero) then
4223 w(ixo^s,psi_) = abs(twofl_glm_alpha)*wct(ixo^s,psi_)
4224 else
4225 ! implicit update of Psi variable
4226 ! equation (27) in Mignone 2010 J. Com. Phys. 229, 2117
4227 if(slab_uniform) then
4228 w(ixo^s,psi_) = dexp(-qdt*cmax_global*twofl_glm_alpha/minval(dxlevel(:)))*w(ixo^s,psi_)
4229 else
4230 w(ixo^s,psi_) = dexp(-qdt*cmax_global*twofl_glm_alpha/minval(block%ds(ixo^s,:),dim=ndim+1))*w(ixo^s,psi_)
4231 end if
4232 end if
4233
4234 ! gradient of Psi
4235 do idim=1,ndim
4236 select case(typegrad)
4237 case("central")
4238 call gradient(wct(ixi^s,psi_),ixi^l,ixo^l,idim,gradpsi)
4239 case("limited")
4240 call gradientl(wct(ixi^s,psi_),ixi^l,ixo^l,idim,gradpsi)
4241 end select
4242 if (phys_total_energy) then
4243 ! e = e -qdt (b . grad(Psi))
4244 w(ixo^s,e_c_) = w(ixo^s,e_c_)-qdt*wct(ixo^s,mag(idim))*gradpsi(ixo^s)
4245 end if
4246 end do
4247
4248 ! m = m - qdt b div b
4249 do idir=1,ndir
4250 w(ixo^s,mom_c(idir))=w(ixo^s,mom_c(idir))-qdt*twofl_mag_i_all(w,ixi^l,ixo^l,idir)*divb(ixo^s)
4251 end do
4252
4253 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_glm')
4254
4255 end subroutine add_source_glm
4256
4257 !> Add divB related sources to w within ixO corresponding to Powel
4258 subroutine add_source_powel(qdt,ixI^L,ixO^L,wCT,w,x)
4260
4261 integer, intent(in) :: ixi^l, ixo^l
4262 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
4263 double precision, intent(inout) :: w(ixi^s,1:nw)
4264 double precision :: divb(ixi^s),v(ixi^s,1:ndir)
4265 integer :: idir
4266
4267 ! We calculate now div B
4268 call get_divb(wct,ixi^l,ixo^l,divb, twofl_divb_nth)
4269
4270 ! calculate velocity
4271 call twofl_get_v_c(wct,x,ixi^l,ixo^l,v)
4272
4273 if (phys_total_energy) then
4274 ! e = e - qdt (v . b) * div b
4275 w(ixo^s,e_c_)=w(ixo^s,e_c_)-&
4276 qdt*sum(v(ixo^s,:)*wct(ixo^s,mag(:)),dim=ndim+1)*divb(ixo^s)
4277 end if
4278
4279 ! b = b - qdt v * div b
4280 do idir=1,ndir
4281 w(ixo^s,mag(idir))=w(ixo^s,mag(idir))-qdt*v(ixo^s,idir)*divb(ixo^s)
4282 end do
4283
4284 ! m = m - qdt b div b
4285 do idir=1,ndir
4286 w(ixo^s,mom_c(idir))=w(ixo^s,mom_c(idir))-qdt*twofl_mag_i_all(w,ixi^l,ixo^l,idir)*divb(ixo^s)
4287 end do
4288
4289 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_powel')
4290
4291 end subroutine add_source_powel
4292
4293 subroutine add_source_janhunen(qdt,ixI^L,ixO^L,wCT,w,x)
4294 ! Add divB related sources to w within ixO
4295 ! corresponding to Janhunen, just the term in the induction equation.
4297
4298 integer, intent(in) :: ixi^l, ixo^l
4299 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
4300 double precision, intent(inout) :: w(ixi^s,1:nw)
4301 double precision :: divb(ixi^s),vel(ixi^s)
4302 integer :: idir
4303
4304 ! We calculate now div B
4305 call get_divb(wct,ixi^l,ixo^l,divb, twofl_divb_nth)
4306
4307 ! b = b - qdt v * div b
4308 do idir=1,ndir
4309 call twofl_get_v_c_idim(wct,x,ixi^l,ixo^l,idir,vel)
4310 w(ixo^s,mag(idir))=w(ixo^s,mag(idir))-qdt*vel(ixo^s)*divb(ixo^s)
4311 end do
4312
4313 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_janhunen')
4314
4315 end subroutine add_source_janhunen
4316
4317 subroutine add_source_linde(qdt,ixI^L,ixO^L,wCT,w,x)
4318 ! Add Linde's divB related sources to wnew within ixO
4320 use mod_geometry
4321
4322 integer, intent(in) :: ixi^l, ixo^l
4323 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
4324 double precision, intent(inout) :: w(ixi^s,1:nw)
4325 integer :: idim, idir, ixp^l, i^d, iside
4326 double precision :: divb(ixi^s),graddivb(ixi^s)
4327 logical, dimension(-1:1^D&) :: leveljump
4328
4329 ! Calculate div B
4330 ixp^l=ixo^l^ladd1;
4331 call get_divb(wct,ixi^l,ixp^l,divb, twofl_divb_nth)
4332
4333 ! for AMR stability, retreat one cell layer from the boarders of level jump
4334 {do i^db=-1,1\}
4335 if(i^d==0|.and.) cycle
4336 if(neighbor_type(i^d,block%igrid)==2 .or. neighbor_type(i^d,block%igrid)==4) then
4337 leveljump(i^d)=.true.
4338 else
4339 leveljump(i^d)=.false.
4340 end if
4341 {end do\}
4342
4343 ixp^l=ixo^l;
4344 do idim=1,ndim
4345 select case(idim)
4346 {case(^d)
4347 do iside=1,2
4348 i^dd=kr(^dd,^d)*(2*iside-3);
4349 if (leveljump(i^dd)) then
4350 if (iside==1) then
4351 ixpmin^d=ixomin^d-i^d
4352 else
4353 ixpmax^d=ixomax^d-i^d
4354 end if
4355 end if
4356 end do
4357 \}
4358 end select
4359 end do
4360
4361 ! Add Linde's diffusive terms
4362 do idim=1,ndim
4363 ! Calculate grad_idim(divb)
4364 select case(typegrad)
4365 case("central")
4366 call gradient(divb,ixi^l,ixp^l,idim,graddivb)
4367 case("limited")
4368 call gradientl(divb,ixi^l,ixp^l,idim,graddivb)
4369 end select
4370
4371 ! Multiply by Linde's eta*dt = divbdiff*(c_max*dx)*dt = divbdiff*dx**2
4372 if (slab_uniform) then
4373 graddivb(ixp^s)=graddivb(ixp^s)*divbdiff/(^d&1.0d0/dxlevel(^d)**2+)
4374 else
4375 graddivb(ixp^s)=graddivb(ixp^s)*divbdiff &
4376 /(^d&1.0d0/block%ds(ixp^s,^d)**2+)
4377 end if
4378
4379 w(ixp^s,mag(idim))=w(ixp^s,mag(idim))+graddivb(ixp^s)
4380
4381 if (typedivbdiff=='all' .and. phys_total_energy) then
4382 ! e += B_idim*eta*grad_idim(divb)
4383 w(ixp^s,e_c_)=w(ixp^s,e_c_)+wct(ixp^s,mag(idim))*graddivb(ixp^s)
4384 end if
4385 end do
4386
4387 if (fix_small_values) call twofl_handle_small_values(.false.,w,x,ixi^l,ixo^l,'add_source_linde')
4388
4389 end subroutine add_source_linde
4390
4391
4392 !> get dimensionless div B = |divB| * volume / area / |B|
4393 subroutine get_normalized_divb(w,ixI^L,ixO^L,divb)
4394
4396
4397 integer, intent(in) :: ixi^l, ixo^l
4398 double precision, intent(in) :: w(ixi^s,1:nw)
4399 double precision :: divb(ixi^s), dsurface(ixi^s)
4400
4401 double precision :: invb(ixo^s)
4402 integer :: ixa^l,idims
4403
4404 call get_divb(w,ixi^l,ixo^l,divb)
4405 invb(ixo^s)=sqrt(twofl_mag_en_all(w,ixi^l,ixo^l))
4406 where(invb(ixo^s)/=0.d0)
4407 invb(ixo^s)=1.d0/invb(ixo^s)
4408 end where
4409 if(slab_uniform) then
4410 divb(ixo^s)=0.5d0*abs(divb(ixo^s))*invb(ixo^s)/sum(1.d0/dxlevel(:))
4411 else
4412 ixamin^d=ixomin^d-1;
4413 ixamax^d=ixomax^d-1;
4414 dsurface(ixo^s)= sum(block%surfaceC(ixo^s,:),dim=ndim+1)
4415 do idims=1,ndim
4416 ixa^l=ixo^l-kr(idims,^d);
4417 dsurface(ixo^s)=dsurface(ixo^s)+block%surfaceC(ixa^s,idims)
4418 end do
4419 divb(ixo^s)=abs(divb(ixo^s))*invb(ixo^s)*&
4420 block%dvolume(ixo^s)/dsurface(ixo^s)
4421 end if
4422
4423 end subroutine get_normalized_divb
4424
4425 !> Calculate idirmin and the idirmin:3 components of the common current array
4426 !> make sure that dxlevel(^D) is set correctly.
4427 subroutine get_current(w,ixI^L,ixO^L,idirmin,current)
4429 use mod_geometry
4430
4431 integer, intent(in) :: ixo^l, ixi^l
4432 double precision, intent(in) :: w(ixi^s,1:nw)
4433 integer, intent(out) :: idirmin
4434 integer :: idir, idirmin0
4435
4436 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
4437 double precision :: current(ixi^s,7-2*ndir:3),bvec(ixi^s,1:ndir)
4438
4439 idirmin0 = 7-2*ndir
4440
4441 bvec(ixi^s,1:ndir)=w(ixi^s,mag(1:ndir))
4442
4443 call curlvector(bvec,ixi^l,ixo^l,current,idirmin,idirmin0,ndir)
4444
4445 if(b0field) current(ixo^s,idirmin0:3)=current(ixo^s,idirmin0:3)+&
4446 block%J0(ixo^s,idirmin0:3)
4447
4448 end subroutine get_current
4449
4450 ! copied from gravity
4451 !> w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
4452 subroutine gravity_add_source(qdt,ixI^L,ixO^L,wCT,w,x,&
4453 energy,qsourcesplit,active)
4455 use mod_usr_methods
4456
4457 integer, intent(in) :: ixi^l, ixo^l
4458 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
4459 double precision, intent(in) :: wct(ixi^s,1:nw)
4460 double precision, intent(inout) :: w(ixi^s,1:nw)
4461 logical, intent(in) :: energy,qsourcesplit
4462 logical, intent(inout) :: active
4463 double precision :: vel(ixi^s)
4464 integer :: idim
4465
4466 double precision :: gravity_field(ixi^s,ndim)
4467
4468 if(qsourcesplit .eqv. grav_split) then
4469 active = .true.
4470
4471 if (.not. associated(usr_gravity)) then
4472 write(*,*) "mod_usr.t: please point usr_gravity to a subroutine"
4473 write(*,*) "like the phys_gravity in mod_usr_methods.t"
4474 call mpistop("gravity_add_source: usr_gravity not defined")
4475 else
4476 call usr_gravity(ixi^l,ixo^l,wct,x,gravity_field)
4477 end if
4478
4479 do idim = 1, ndim
4480 w(ixo^s,mom_n(idim)) = w(ixo^s,mom_n(idim)) &
4481 + qdt * gravity_field(ixo^s,idim) * wct(ixo^s,rho_n_)
4482 w(ixo^s,mom_c(idim)) = w(ixo^s,mom_c(idim)) &
4483 + qdt * gravity_field(ixo^s,idim) * wct(ixo^s,rho_c_)
4484 if(energy) then
4485#if !defined(E_RM_W0) || E_RM_W0 == 1
4486 call twofl_get_v_n_idim(wct,x,ixi^l,ixo^l,idim,vel)
4487 w(ixo^s,e_n_)=w(ixo^s,e_n_) &
4488 + qdt * gravity_field(ixo^s,idim) * vel(ixo^s) * wct(ixo^s,rho_n_)
4489 call twofl_get_v_c_idim(wct,x,ixi^l,ixo^l,idim,vel)
4490 w(ixo^s,e_c_)=w(ixo^s,e_c_) &
4491 + qdt * gravity_field(ixo^s,idim) * vel(ixo^s) * wct(ixo^s,rho_c_)
4492#else
4493 w(ixo^s,e_n_)=w(ixo^s,e_n_) &
4494 + qdt * gravity_field(ixo^s,idim) * wct(ixo^s,mom_n(idim))
4495 w(ixo^s,e_c_)=w(ixo^s,e_c_) &
4496 + qdt * gravity_field(ixo^s,idim) * wct(ixo^s,mom_c(idim))
4497#endif
4498
4499
4500 end if
4501 end do
4502 end if
4503
4504 end subroutine gravity_add_source
4505
4506 subroutine gravity_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
4508 use mod_usr_methods
4509
4510 integer, intent(in) :: ixi^l, ixo^l
4511 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim), w(ixi^s,1:nw)
4512 double precision, intent(inout) :: dtnew
4513
4514 double precision :: dxinv(1:ndim), max_grav
4515 integer :: idim
4516
4517 double precision :: gravity_field(ixi^s,ndim)
4518
4519 ^d&dxinv(^d)=one/dx^d;
4520
4521 if(.not. associated(usr_gravity)) then
4522 write(*,*) "mod_usr.t: please point usr_gravity to a subroutine"
4523 write(*,*) "like the phys_gravity in mod_usr_methods.t"
4524 call mpistop("gravity_get_dt: usr_gravity not defined")
4525 else
4526 call usr_gravity(ixi^l,ixo^l,w,x,gravity_field)
4527 end if
4528
4529 do idim = 1, ndim
4530 max_grav = maxval(abs(gravity_field(ixo^s,idim)))
4531 max_grav = max(max_grav, epsilon(1.0d0))
4532 dtnew = min(dtnew, 1.0d0 / sqrt(max_grav * dxinv(idim)))
4533 end do
4534
4535 end subroutine gravity_get_dt
4536
4537 !> If resistivity is not zero, check diffusion time limit for dt
4538 subroutine twofl_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
4540 use mod_usr_methods
4542 !use mod_viscosity, only: viscosity_get_dt
4543 !use mod_gravity, only: gravity_get_dt
4544
4545 integer, intent(in) :: ixi^l, ixo^l
4546 double precision, intent(inout) :: dtnew
4547 double precision, intent(in) :: dx^d
4548 double precision, intent(in) :: w(ixi^s,1:nw)
4549 double precision, intent(in) :: x(ixi^s,1:ndim)
4550
4551 integer :: idirmin,idim
4552 double precision :: dxarr(ndim)
4553 double precision :: current(ixi^s,7-2*ndir:3),eta(ixi^s)
4554
4555 dtnew = bigdouble
4556
4557 ^d&dxarr(^d)=dx^d;
4558 if (twofl_eta>zero)then
4559 dtnew=dtdiffpar*minval(dxarr(1:ndim))**2/twofl_eta
4560 else if (twofl_eta<zero)then
4561 call get_current(w,ixi^l,ixo^l,idirmin,current)
4562 call usr_special_resistivity(w,ixi^l,ixo^l,idirmin,x,current,eta)
4563 dtnew=bigdouble
4564 do idim=1,ndim
4565 if(slab_uniform) then
4566 dtnew=min(dtnew,&
4567 dtdiffpar/(smalldouble+maxval(eta(ixo^s)/dxarr(idim)**2)))
4568 else
4569 dtnew=min(dtnew,&
4570 dtdiffpar/(smalldouble+maxval(eta(ixo^s)/block%ds(ixo^s,idim)**2)))
4571 end if
4572 end do
4573 end if
4574
4575 if(twofl_eta_hyper>zero) then
4576 if(slab_uniform) then
4577 dtnew=min(dtdiffpar*minval(dxarr(1:ndim))**4/twofl_eta_hyper,dtnew)
4578 else
4579 dtnew=min(dtdiffpar*minval(block%ds(ixo^s,1:ndim))**4/twofl_eta_hyper,dtnew)
4580 end if
4581 end if
4582
4583 ! the timestep related to coll terms: 1/(rho_n rho_c alpha)
4584 if(dtcollpar>0d0 .and. has_collisions()) then
4585 call coll_get_dt(w,x,ixi^l,ixo^l,dtnew)
4586 endif
4587
4589 call cooling_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x,rc_fl_c)
4590 end if
4592 call cooling_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x,rc_fl_n)
4593 end if
4594!
4595! if(twofl_viscosity) then
4596! call viscosity_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
4597! end if
4598!
4599 if(twofl_gravity) then
4600 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
4601 end if
4602 if(twofl_hyperdiffusivity) then
4603 call hyperdiffusivity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
4604 end if
4605
4606
4607 end subroutine twofl_get_dt
4608
4609 pure function has_collisions() result(res)
4610 logical :: res
4611 res = .not. twofl_alpha_coll_constant .or. twofl_alpha_coll >0d0
4612 end function has_collisions
4613
4614 subroutine coll_get_dt(w,x,ixI^L,ixO^L,dtnew)
4616 integer, intent(in) :: ixi^l, ixo^l
4617 double precision, intent(in) :: w(ixi^s,1:nw)
4618 double precision, intent(in) :: x(ixi^s,1:ndim)
4619 double precision, intent(inout) :: dtnew
4620
4621 double precision :: rhon(ixi^s), rhoc(ixi^s), alpha(ixi^s)
4622 double precision, allocatable :: gamma_rec(:^d&), gamma_ion(:^D&)
4623 double precision :: max_coll_rate
4624
4625 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
4626 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
4627
4628 call get_alpha_coll(ixi^l, ixo^l, w, x, alpha)
4629 max_coll_rate = maxval(alpha(ixo^s) * max(rhon(ixo^s), rhoc(ixo^s)))
4630
4631 if(twofl_coll_inc_ionrec) then
4632 allocate(gamma_ion(ixi^s), gamma_rec(ixi^s))
4633 call get_gamma_ion_rec(ixi^l, ixo^l, w, x, gamma_rec, gamma_ion)
4634 max_coll_rate=max(max_coll_rate, maxval(gamma_ion(ixo^s)), maxval(gamma_rec(ixo^s)))
4635 deallocate(gamma_ion, gamma_rec)
4636 endif
4637 dtnew = min(dtcollpar/max_coll_rate, dtnew)
4638
4639 end subroutine coll_get_dt
4640
4641 ! Add geometrical source terms to w
4642 subroutine twofl_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
4644 use mod_geometry
4645
4646 integer, intent(in) :: ixi^l, ixo^l
4647 double precision, intent(in) :: qdt, dtfactor,x(ixi^s,1:ndim)
4648 double precision, intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw), w(ixi^s,1:nw)
4649
4650 integer :: iw,idir, h1x^l{^nooned, h2x^l}
4651 double precision :: tmp(ixi^s),tmp1(ixi^s),tmp2(ixi^s),rho(ixi^s)
4652
4653 integer :: mr_,mphi_ ! Polar var. names
4654 integer :: br_,bphi_
4655
4656 ! charges
4657
4658 mr_=mom_c(1); mphi_=mom_c(1)-1+phi_ ! Polar var. names
4659 br_=mag(1); bphi_=mag(1)-1+phi_
4660 call get_rhoc_tot(wct,x,ixi^l,ixo^l,rho)
4661
4662 select case (coordinate)
4663 case (cylindrical)
4664 call twofl_get_p_c_total(wct,x,ixi^l,ixo^l,tmp)
4665
4666 if(phi_>0) then
4667 w(ixo^s,mr_)=w(ixo^s,mr_)+qdt/x(ixo^s,1)*(tmp(ixo^s)-&
4668 wct(ixo^s,bphi_)**2+wct(ixo^s,mphi_)**2/rho(ixo^s))
4669 w(ixo^s,mphi_)=w(ixo^s,mphi_)+qdt/x(ixo^s,1)*(&
4670 -wct(ixo^s,mphi_)*wct(ixo^s,mr_)/rho(ixo^s) &
4671 +wct(ixo^s,bphi_)*wct(ixo^s,br_))
4672 if(.not.stagger_grid) then
4673 w(ixo^s,bphi_)=w(ixo^s,bphi_)+qdt/x(ixo^s,1)*&
4674 (wct(ixo^s,bphi_)*wct(ixo^s,mr_) &
4675 -wct(ixo^s,br_)*wct(ixo^s,mphi_)) &
4676 /rho(ixo^s)
4677 end if
4678 else
4679 w(ixo^s,mr_)=w(ixo^s,mr_)+qdt/x(ixo^s,1)*tmp(ixo^s)
4680 end if
4681 if(twofl_glm) w(ixo^s,br_)=w(ixo^s,br_)+qdt*wct(ixo^s,psi_)/x(ixo^s,1)
4682 case (spherical)
4683 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
4684 call twofl_get_p_c_total(wct,x,ixi^l,ixo^l,tmp1)
4685 tmp(ixo^s)=tmp1(ixo^s)
4686 if(b0field) then
4687 tmp2(ixo^s)=sum(block%B0(ixo^s,:,0)*wct(ixo^s,mag(:)),dim=ndim+1)
4688 tmp(ixo^s)=tmp(ixo^s)+tmp2(ixo^s)
4689 end if
4690 ! m1
4691 tmp(ixo^s)=tmp(ixo^s)*x(ixo^s,1) &
4692 *(block%surfaceC(ixo^s,1)-block%surfaceC(h1x^s,1))/block%dvolume(ixo^s)
4693 if(ndir>1) then
4694 do idir=2,ndir
4695 tmp(ixo^s)=tmp(ixo^s)+wct(ixo^s,mom_c(idir))**2/rho(ixo^s)-wct(ixo^s,mag(idir))**2
4696 if(b0field) tmp(ixo^s)=tmp(ixo^s)-2.0d0*block%B0(ixo^s,idir,0)*wct(ixo^s,mag(idir))
4697 end do
4698 end if
4699 w(ixo^s,mom_c(1))=w(ixo^s,mom_c(1))+qdt*tmp(ixo^s)/x(ixo^s,1)
4700 ! b1
4701 if(twofl_glm) then
4702 w(ixo^s,mag(1))=w(ixo^s,mag(1))+qdt/x(ixo^s,1)*2.0d0*wct(ixo^s,psi_)
4703 end if
4704
4705 {^nooned
4706 ! m2
4707 tmp(ixo^s)=tmp1(ixo^s)
4708 if(b0field) then
4709 tmp(ixo^s)=tmp(ixo^s)+tmp2(ixo^s)
4710 end if
4711 ! This will make hydrostatic p=const an exact solution
4712 w(ixo^s,mom_c(2))=w(ixo^s,mom_c(2))+qdt*tmp(ixo^s) &
4713 *(block%surfaceC(ixo^s,2)-block%surfaceC(h2x^s,2)) &
4714 /block%dvolume(ixo^s)
4715 tmp(ixo^s)=-(wct(ixo^s,mom_c(1))*wct(ixo^s,mom_c(2))/rho(ixo^s) &
4716 -wct(ixo^s,mag(1))*wct(ixo^s,mag(2)))
4717 if (b0field) then
4718 tmp(ixo^s)=tmp(ixo^s)+block%B0(ixo^s,1,0)*wct(ixo^s,mag(2)) &
4719 +wct(ixo^s,mag(1))*block%B0(ixo^s,2,0)
4720 end if
4721 if(ndir==3) then
4722 tmp(ixo^s)=tmp(ixo^s)+(wct(ixo^s,mom_c(3))**2/rho(ixo^s) &
4723 -wct(ixo^s,mag(3))**2)*dcos(x(ixo^s,2))/dsin(x(ixo^s,2))
4724 if (b0field) then
4725 tmp(ixo^s)=tmp(ixo^s)-2.0d0*block%B0(ixo^s,3,0)*wct(ixo^s,mag(3))&
4726 *dcos(x(ixo^s,2))/dsin(x(ixo^s,2))
4727 end if
4728 end if
4729 w(ixo^s,mom_c(2))=w(ixo^s,mom_c(2))+qdt*tmp(ixo^s)/x(ixo^s,1)
4730 ! b2
4731 if(.not.stagger_grid) then
4732 tmp(ixo^s)=(wct(ixo^s,mom_c(1))*wct(ixo^s,mag(2)) &
4733 -wct(ixo^s,mom_c(2))*wct(ixo^s,mag(1)))/rho(ixo^s)
4734 if(b0field) then
4735 tmp(ixo^s)=tmp(ixo^s)+(wct(ixo^s,mom_c(1))*block%B0(ixo^s,2,0) &
4736 -wct(ixo^s,mom_c(2))*block%B0(ixo^s,1,0))/rho(ixo^s)
4737 end if
4738 if(twofl_glm) then
4739 tmp(ixo^s)=tmp(ixo^s) &
4740 + dcos(x(ixo^s,2))/dsin(x(ixo^s,2))*wct(ixo^s,psi_)
4741 end if
4742 w(ixo^s,mag(2))=w(ixo^s,mag(2))+qdt*tmp(ixo^s)/x(ixo^s,1)
4743 end if
4744 }
4745
4746 if(ndir==3) then
4747 ! m3
4748 tmp(ixo^s)=-(wct(ixo^s,mom_c(3))*wct(ixo^s,mom_c(1))/rho(ixo^s) &
4749 -wct(ixo^s,mag(3))*wct(ixo^s,mag(1))) {^nooned &
4750 -(wct(ixo^s,mom_c(2))*wct(ixo^s,mom_c(3))/rho(ixo^s) &
4751 -wct(ixo^s,mag(2))*wct(ixo^s,mag(3))) &
4752 *dcos(x(ixo^s,2))/dsin(x(ixo^s,2)) }
4753 if (b0field) then
4754 tmp(ixo^s)=tmp(ixo^s)+block%B0(ixo^s,1,0)*wct(ixo^s,mag(3)) &
4755 +wct(ixo^s,mag(1))*block%B0(ixo^s,3,0) {^nooned &
4756 +(block%B0(ixo^s,2,0)*wct(ixo^s,mag(3)) &
4757 +wct(ixo^s,mag(2))*block%B0(ixo^s,3,0)) &
4758 *dcos(x(ixo^s,2))/dsin(x(ixo^s,2)) }
4759 end if
4760 w(ixo^s,mom_c(3))=w(ixo^s,mom_c(3))+qdt*tmp(ixo^s)/x(ixo^s,1)
4761 ! b3
4762 if(.not.stagger_grid) then
4763 tmp(ixo^s)=(wct(ixo^s,mom_c(1))*wct(ixo^s,mag(3)) &
4764 -wct(ixo^s,mom_c(3))*wct(ixo^s,mag(1)))/rho(ixo^s) {^nooned &
4765 -(wct(ixo^s,mom_c(3))*wct(ixo^s,mag(2)) &
4766 -wct(ixo^s,mom_c(2))*wct(ixo^s,mag(3)))*dcos(x(ixo^s,2)) &
4767 /(rho(ixo^s)*dsin(x(ixo^s,2))) }
4768 if (b0field) then
4769 tmp(ixo^s)=tmp(ixo^s)+(wct(ixo^s,mom_c(1))*block%B0(ixo^s,3,0) &
4770 -wct(ixo^s,mom_c(3))*block%B0(ixo^s,1,0))/rho(ixo^s){^nooned &
4771 -(wct(ixo^s,mom_c(3))*block%B0(ixo^s,2,0) &
4772 -wct(ixo^s,mom_c(2))*block%B0(ixo^s,3,0))*dcos(x(ixo^s,2)) &
4773 /(rho(ixo^s)*dsin(x(ixo^s,2))) }
4774 end if
4775 w(ixo^s,mag(3))=w(ixo^s,mag(3))+qdt*tmp(ixo^s)/x(ixo^s,1)
4776 end if
4777 end if
4778 end select
4779
4780 ! neutrals
4781 !TODO no dust: see and implement them from hd/mod_hd_phys !
4782 !uncomment cartesian expansion
4783 call get_rhon_tot(wct,x,ixi^l,ixo^l,rho)
4784 call twofl_get_pthermal_n(wct, x, ixi^l, ixo^l, tmp1)
4785
4786 select case (coordinate)
4787! case(Cartesian_expansion)
4788! !the user provides the functions of exp_factor and del_exp_factor
4789! if(associated(usr_set_surface)) call usr_set_surface(ixI^L,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
4790! tmp(ixO^S) = tmp1(ixO^S)*del_exp_factor(ixO^S)/exp_factor(ixO^S)
4791! w(ixO^S,mom(1)) = w(ixO^S,mom(1)) + qdt*tmp(ixO^S)
4792
4793 case (cylindrical)
4794 mr_ = mom_n(r_)
4795 if (phi_ > 0) then
4796 where (rho(ixo^s) > 0d0)
4797 tmp(ixo^s) = tmp1(ixo^s) + wct(ixo^s, mphi_)**2 / rho(ixo^s)
4798 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * tmp(ixo^s) / x(ixo^s, r_)
4799 end where
4800 ! s[mphi]=(-mphi*mr/rho)/radius
4801 where (rho(ixo^s) > 0d0)
4802 tmp(ixo^s) = -wct(ixo^s, mphi_) * wct(ixo^s, mr_) / rho(ixo^s)
4803 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * tmp(ixo^s) / x(ixo^s, r_)
4804 end where
4805 else
4806 ! s[mr]=2pthermal/radius
4807 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * tmp1(ixo^s) / x(ixo^s, r_)
4808 end if
4809 case (spherical)
4810 if(phi_>0) mphi_ = mom_n(phi_)
4811 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
4812 ! s[mr]=((mtheta**2+mphi**2)/rho+2*p)/r
4813 tmp(ixo^s) = tmp1(ixo^s) * x(ixo^s, 1) &
4814 *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
4815 /block%dvolume(ixo^s)
4816 if (ndir > 1) then
4817 do idir = 2, ndir
4818 tmp(ixo^s) = tmp(ixo^s) + wct(ixo^s, mom_n(idir))**2 / rho(ixo^s)
4819 end do
4820 end if
4821 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * tmp(ixo^s) / x(ixo^s, 1)
4822
4823 {^nooned
4824 ! s[mtheta]=-(mr*mtheta/rho)/r+cot(theta)*(mphi**2/rho+p)/r
4825 tmp(ixo^s) = tmp1(ixo^s) * x(ixo^s, 1) &
4826 * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
4827 / block%dvolume(ixo^s)
4828 if (ndir == 3) then
4829 tmp(ixo^s) = tmp(ixo^s) + (wct(ixo^s, mom_n(3))**2 / rho(ixo^s)) / tan(x(ixo^s, 2))
4830 end if
4831 tmp(ixo^s) = tmp(ixo^s) - (wct(ixo^s, mom_n(2)) * wct(ixo^s, mr_)) / rho(ixo^s)
4832 w(ixo^s, mom_n(2)) = w(ixo^s, mom_n(2)) + qdt * tmp(ixo^s) / x(ixo^s, 1)
4833
4834 if (ndir == 3) then
4835 ! s[mphi]=-(mphi*mr/rho)/r-cot(theta)*(mtheta*mphi/rho)/r
4836 tmp(ixo^s) = -(wct(ixo^s, mom_n(3)) * wct(ixo^s, mr_)) / rho(ixo^s)&
4837 - (wct(ixo^s, mom_n(2)) * wct(ixo^s, mom_n(3))) / rho(ixo^s) / tan(x(ixo^s, 2))
4838 w(ixo^s, mom_n(3)) = w(ixo^s, mom_n(3)) + qdt * tmp(ixo^s) / x(ixo^s, 1)
4839 end if
4840 }
4841 end select
4842
4843 contains
4844 subroutine twofl_get_p_c_total(w,x,ixI^L,ixO^L,p)
4846
4847 integer, intent(in) :: ixI^L, ixO^L
4848 double precision, intent(in) :: w(ixI^S,nw)
4849 double precision, intent(in) :: x(ixI^S,1:ndim)
4850 double precision, intent(out) :: p(ixI^S)
4851
4852 call twofl_get_pthermal_c(w,x,ixi^l,ixo^l,p)
4853
4854 p(ixo^s) = p(ixo^s) + 0.5d0 * sum(w(ixo^s, mag(:))**2, dim=ndim+1)
4855
4856 end subroutine twofl_get_p_c_total
4857
4858 end subroutine twofl_add_source_geom
4859
4860 subroutine twofl_get_temp_c_pert_from_etot(w, x, ixI^L, ixO^L, res)
4862 integer, intent(in) :: ixI^L, ixO^L
4863 double precision, intent(in) :: w(ixI^S, 1:nw)
4864 double precision, intent(in) :: x(ixI^S, 1:ndim)
4865 double precision, intent(out):: res(ixI^S)
4866
4867 ! store pe1 in res
4868 res(ixo^s)=(gamma_1*(w(ixo^s,e_c_)&
4869 - twofl_kin_en_c(w,ixi^l,ixo^l)&
4870 - twofl_mag_en(w,ixi^l,ixo^l)))
4871 if(has_equi_pe_c0) then
4872 res(ixo^s) = res(ixo^s) + block%equi_vars(ixo^s,equi_pe_c0_,b0i)
4873 if(has_equi_rho_c0) then
4874 res(ixo^s) = res(ixo^s)/(rc * (w(ixo^s,rho_c_)+ block%equi_vars(ixo^s,equi_rho_c0_,b0i))) - &
4875 block%equi_vars(ixo^s,equi_pe_c0_,b0i)/(rc * block%equi_vars(ixo^s,equi_rho_c0_,b0i))
4876 else
4877 ! infinite equi temperature with p0 and 0 density
4878 res(ixo^s) = 0d0
4879 endif
4880 else
4881 res(ixo^s) = res(ixo^s)/(rc * w(ixo^s,rho_c_))
4882 endif
4883
4884 end subroutine twofl_get_temp_c_pert_from_etot
4885
4886 !> Compute 2 times total magnetic energy
4887 function twofl_mag_en_all(w, ixI^L, ixO^L) result(mge)
4889 integer, intent(in) :: ixI^L, ixO^L
4890 double precision, intent(in) :: w(ixI^S, nw)
4891 double precision :: mge(ixO^S)
4892
4893 if (b0field) then
4894 mge(ixo^s) = sum((w(ixo^s, mag(:))+block%B0(ixo^s,:,b0i))**2, dim=ndim+1)
4895 else
4896 mge(ixo^s) = sum(w(ixo^s, mag(:))**2, dim=ndim+1)
4897 end if
4898 end function twofl_mag_en_all
4899
4900 !> Compute full magnetic field by direction
4901 function twofl_mag_i_all(w, ixI^L, ixO^L,idir) result(mgf)
4903 integer, intent(in) :: ixI^L, ixO^L, idir
4904 double precision, intent(in) :: w(ixI^S, nw)
4905 double precision :: mgf(ixO^S)
4906
4907 if (b0field) then
4908 mgf(ixo^s) = w(ixo^s, mag(idir))+block%B0(ixo^s,idir,b0i)
4909 else
4910 mgf(ixo^s) = w(ixo^s, mag(idir))
4911 end if
4912 end function twofl_mag_i_all
4913
4914 !> Compute evolving magnetic energy
4915 function twofl_mag_en(w, ixI^L, ixO^L) result(mge)
4916 use mod_global_parameters, only: nw, ndim
4917 integer, intent(in) :: ixI^L, ixO^L
4918 double precision, intent(in) :: w(ixI^S, nw)
4919 double precision :: mge(ixO^S)
4920
4921 mge(ixo^s) = 0.5d0 * sum(w(ixo^s, mag(:))**2, dim=ndim+1)
4922 end function twofl_mag_en
4923
4924 !> compute kinetic energy of neutrals
4925 function twofl_kin_en_n(w, ixI^L, ixO^L) result(ke)
4926 use mod_global_parameters, only: nw, ndim,block
4927 integer, intent(in) :: ixI^L, ixO^L
4928 double precision, intent(in) :: w(ixI^S, nw)
4929 double precision :: ke(ixO^S)
4930
4931 if(has_equi_rho_n0) then
4932 ke(ixo^s) = 0.5d0 * sum(w(ixo^s, mom_n(:))**2, dim=ndim+1) / (w(ixo^s, rho_n_) + block%equi_vars(ixo^s,equi_rho_n0_,0))
4933 else
4934 ke(ixo^s) = 0.5d0 * sum(w(ixo^s, mom_n(:))**2, dim=ndim+1) / w(ixo^s, rho_n_)
4935 endif
4936
4937 end function twofl_kin_en_n
4938
4939 subroutine twofl_get_temp_n_pert_from_etot(w, x, ixI^L, ixO^L, res)
4941 integer, intent(in) :: ixI^L, ixO^L
4942 double precision, intent(in) :: w(ixI^S, 1:nw)
4943 double precision, intent(in) :: x(ixI^S, 1:ndim)
4944 double precision, intent(out):: res(ixI^S)
4945
4946 ! store pe1 in res
4947 res(ixo^s)=(gamma_1*(w(ixo^s,e_c_)- twofl_kin_en_c(w,ixi^l,ixo^l)))
4948 if(has_equi_pe_n0) then
4949 res(ixo^s) = res(ixo^s) + block%equi_vars(ixo^s,equi_pe_n0_,b0i)
4950 if(has_equi_rho_n0) then
4951 res(ixo^s) = res(ixo^s)/(rn * (w(ixo^s,rho_n_)+ block%equi_vars(ixo^s,equi_rho_n0_,b0i))) - &
4952 block%equi_vars(ixo^s,equi_pe_n0_,b0i)/(rn * block%equi_vars(ixo^s,equi_rho_n0_,b0i))
4953 else
4954 ! infinite equi temperature with p0 and 0 density
4955 res(ixo^s) = 0d0
4956 endif
4957 else
4958 res(ixo^s) = res(ixo^s)/(rn * w(ixo^s,rho_n_))
4959 endif
4960
4961 end subroutine twofl_get_temp_n_pert_from_etot
4962
4963 !> compute kinetic energy of charges
4964 !> w are conserved variables
4965 function twofl_kin_en_c(w, ixI^L, ixO^L) result(ke)
4966 use mod_global_parameters, only: nw, ndim,block
4967 integer, intent(in) :: ixI^L, ixO^L
4968 double precision, intent(in) :: w(ixI^S, nw)
4969 double precision :: ke(ixO^S)
4970
4971 if(has_equi_rho_c0) then
4972 ke(ixo^s) = 0.5d0 * sum(w(ixo^s, mom_c(:))**2, dim=ndim+1) / (w(ixo^s, rho_c_) + block%equi_vars(ixo^s,equi_rho_c0_,0))
4973 else
4974 ke(ixo^s) = 0.5d0 * sum(w(ixo^s, mom_c(:))**2, dim=ndim+1) / w(ixo^s, rho_c_)
4975 endif
4976 end function twofl_kin_en_c
4977
4978 subroutine twofl_getv_hall(w,x,ixI^L,ixO^L,vHall)
4980
4981 integer, intent(in) :: ixI^L, ixO^L
4982 double precision, intent(in) :: w(ixI^S,nw)
4983 double precision, intent(in) :: x(ixI^S,1:ndim)
4984 double precision, intent(inout) :: vHall(ixI^S,1:3)
4985
4986 integer :: idir, idirmin
4987 double precision :: current(ixI^S,7-2*ndir:3)
4988 double precision :: rho(ixI^S)
4989
4990 call get_rhoc_tot(w,x,ixi^l,ixo^l,rho)
4991 ! Calculate current density and idirmin
4992 call get_current(w,ixi^l,ixo^l,idirmin,current)
4993 vhall(ixo^s,1:3) = zero
4994 vhall(ixo^s,idirmin:3) = - twofl_etah*current(ixo^s,idirmin:3)
4995 do idir = idirmin, 3
4996 vhall(ixo^s,idir) = vhall(ixo^s,idir)/rho(ixo^s)
4997 end do
4998
4999 end subroutine twofl_getv_hall
5000
5001! the following not used
5002! subroutine twofl_getdt_Hall(w,x,ixI^L,ixO^L,dx^D,dthall)
5003! use mod_global_parameters
5004!
5005! integer, intent(in) :: ixI^L, ixO^L
5006! double precision, intent(in) :: dx^D
5007! double precision, intent(in) :: w(ixI^S,1:nw)
5008! double precision, intent(in) :: x(ixI^S,1:ndim)
5009! double precision, intent(out) :: dthall
5010! !.. local ..
5011! double precision :: dxarr(ndim)
5012! double precision :: bmag(ixI^S)
5013!
5014! dthall=bigdouble
5015!
5016! ! because we have that in cmax now:
5017! return
5018!
5019! ^D&dxarr(^D)=dx^D;
5020!
5021! if (.not. B0field) then
5022! bmag(ixO^S)=sqrt(sum(w(ixO^S,mag(:))**2, dim=ndim+1))
5023! bmag(ixO^S)=sqrt(sum((w(ixO^S,mag(:)) + block%B0(ixO^S,1:ndir,b0i))**2))
5024! end if
5025!
5026! if(slab_uniform) then
5027! dthall=dtdiffpar*minval(dxarr(1:ndim))**2.0d0/(twofl_etah*maxval(bmag(ixO^S)/w(ixO^S,rho_c_)))
5028! else
5029! dthall=dtdiffpar*minval(block%ds(ixO^S,1:ndim))**2.0d0/(twofl_etah*maxval(bmag(ixO^S)/w(ixO^S,rho_c_)))
5030! end if
5031!
5032! end subroutine twofl_getdt_Hall
5033
5034 subroutine twofl_modify_wlr(ixI^L,ixO^L,qt,wLC,wRC,wLp,wRp,s,idir)
5036 use mod_usr_methods
5037 integer, intent(in) :: ixI^L, ixO^L, idir
5038 double precision, intent(in) :: qt
5039 double precision, intent(inout) :: wLC(ixI^S,1:nw), wRC(ixI^S,1:nw)
5040 double precision, intent(inout) :: wLp(ixI^S,1:nw), wRp(ixI^S,1:nw)
5041 type(state) :: s
5042 double precision :: dB(ixI^S), dPsi(ixI^S)
5043
5044 if(stagger_grid) then
5045 wlc(ixo^s,mag(idir))=s%ws(ixo^s,idir)
5046 wrc(ixo^s,mag(idir))=s%ws(ixo^s,idir)
5047 wlp(ixo^s,mag(idir))=s%ws(ixo^s,idir)
5048 wrp(ixo^s,mag(idir))=s%ws(ixo^s,idir)
5049 else
5050 ! Solve the Riemann problem for the linear 2x2 system for normal
5051 ! B-field and GLM_Psi according to Dedner 2002:
5052 ! This implements eq. (42) in Dedner et al. 2002 JcP 175
5053 ! Gives the Riemann solution on the interface
5054 ! for the normal B component and Psi in the GLM-MHD system.
5055 ! 23/04/2013 Oliver Porth
5056 db(ixo^s) = wrp(ixo^s,mag(idir)) - wlp(ixo^s,mag(idir))
5057 dpsi(ixo^s) = wrp(ixo^s,psi_) - wlp(ixo^s,psi_)
5058
5059 wlp(ixo^s,mag(idir)) = 0.5d0 * (wrp(ixo^s,mag(idir)) + wlp(ixo^s,mag(idir))) &
5060 - 0.5d0/cmax_global * dpsi(ixo^s)
5061 wlp(ixo^s,psi_) = 0.5d0 * (wrp(ixo^s,psi_) + wlp(ixo^s,psi_)) &
5062 - 0.5d0*cmax_global * db(ixo^s)
5063
5064 wrp(ixo^s,mag(idir)) = wlp(ixo^s,mag(idir))
5065 wrp(ixo^s,psi_) = wlp(ixo^s,psi_)
5066
5067 if(phys_total_energy) then
5068 wrc(ixo^s,e_c_)=wrc(ixo^s,e_c_)-half*wrc(ixo^s,mag(idir))**2
5069 wlc(ixo^s,e_c_)=wlc(ixo^s,e_c_)-half*wlc(ixo^s,mag(idir))**2
5070 end if
5071 wrc(ixo^s,mag(idir)) = wlp(ixo^s,mag(idir))
5072 wrc(ixo^s,psi_) = wlp(ixo^s,psi_)
5073 wlc(ixo^s,mag(idir)) = wlp(ixo^s,mag(idir))
5074 wlc(ixo^s,psi_) = wlp(ixo^s,psi_)
5075 ! modify total energy according to the change of magnetic field
5076 if(phys_total_energy) then
5077 wrc(ixo^s,e_c_)=wrc(ixo^s,e_c_)+half*wrc(ixo^s,mag(idir))**2
5078 wlc(ixo^s,e_c_)=wlc(ixo^s,e_c_)+half*wlc(ixo^s,mag(idir))**2
5079 end if
5080 end if
5081
5082 if(associated(usr_set_wlr)) call usr_set_wlr(ixi^l,ixo^l,qt,wlc,wrc,wlp,wrp,s,idir)
5083
5084 end subroutine twofl_modify_wlr
5085
5086 subroutine twofl_boundary_adjust(igrid,psb)
5088 integer, intent(in) :: igrid
5089 type(state), target :: psb(max_blocks)
5090
5091 integer :: iB, idims, iside, ixO^L, i^D
5092
5093 block=>ps(igrid)
5094 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
5095 do idims=1,ndim
5096 ! to avoid using as yet unknown corner info in more than 1D, we
5097 ! fill only interior mesh ranges of the ghost cell ranges at first,
5098 ! and progressively enlarge the ranges to include corners later
5099 do iside=1,2
5100 i^d=kr(^d,idims)*(2*iside-3);
5101 if (neighbor_type(i^d,igrid)/=1) cycle
5102 ib=(idims-1)*2+iside
5103 if(.not.boundary_divbfix(ib)) cycle
5104 if(any(typeboundary(:,ib)==bc_special)) then
5105 ! MF nonlinear force-free B field extrapolation and data driven
5106 ! require normal B of the first ghost cell layer to be untouched by
5107 ! fixdivB=0 process, set boundary_divbfix_skip(iB)=1 in par file
5108 select case (idims)
5109 {case (^d)
5110 if (iside==2) then
5111 ! maximal boundary
5112 ixomin^dd=ixghi^d+1-nghostcells+boundary_divbfix_skip(2*^d)^d%ixOmin^dd=ixglo^dd;
5113 ixomax^dd=ixghi^dd;
5114 else
5115 ! minimal boundary
5116 ixomin^dd=ixglo^dd;
5117 ixomax^dd=ixglo^d-1+nghostcells-boundary_divbfix_skip(2*^d-1)^d%ixOmax^dd=ixghi^dd;
5118 end if \}
5119 end select
5120 call fixdivb_boundary(ixg^ll,ixo^l,psb(igrid)%w,psb(igrid)%x,ib)
5121 end if
5122 end do
5123 end do
5124
5125 end subroutine twofl_boundary_adjust
5126
5127 subroutine fixdivb_boundary(ixG^L,ixO^L,w,x,iB)
5129
5130 integer, intent(in) :: ixG^L,ixO^L,iB
5131 double precision, intent(inout) :: w(ixG^S,1:nw)
5132 double precision, intent(in) :: x(ixG^S,1:ndim)
5133
5134 double precision :: dx1x2,dx1x3,dx2x1,dx2x3,dx3x1,dx3x2
5135 integer :: ix^D,ixF^L
5136
5137 select case(ib)
5138 case(1)
5139 ! 2nd order CD for divB=0 to set normal B component better
5140 {^iftwod
5141 ixfmin1=ixomin1+1
5142 ixfmax1=ixomax1+1
5143 ixfmin2=ixomin2+1
5144 ixfmax2=ixomax2-1
5145 if(slab_uniform) then
5146 dx1x2=dxlevel(1)/dxlevel(2)
5147 do ix1=ixfmax1,ixfmin1,-1
5148 w(ix1-1,ixfmin2:ixfmax2,mag(1))=w(ix1+1,ixfmin2:ixfmax2,mag(1)) &
5149 +dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))-&
5150 w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))
5151 enddo
5152 else
5153 do ix1=ixfmax1,ixfmin1,-1
5154 w(ix1-1,ixfmin2:ixfmax2,mag(1))=( (w(ix1+1,ixfmin2:ixfmax2,mag(1))+&
5155 w(ix1,ixfmin2:ixfmax2,mag(1)))*block%surfaceC(ix1,ixfmin2:ixfmax2,1)&
5156 +(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))+w(ix1,ixfmin2:ixfmax2,mag(2)))*&
5157 block%surfaceC(ix1,ixfmin2:ixfmax2,2)&
5158 -(w(ix1,ixfmin2:ixfmax2,mag(2))+w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))*&
5159 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,2) )&
5160 /block%surfaceC(ix1-1,ixfmin2:ixfmax2,1)-w(ix1,ixfmin2:ixfmax2,mag(1))
5161 end do
5162 end if
5163 }
5164 {^ifthreed
5165 ixfmin1=ixomin1+1
5166 ixfmax1=ixomax1+1
5167 ixfmin2=ixomin2+1
5168 ixfmax2=ixomax2-1
5169 ixfmin3=ixomin3+1
5170 ixfmax3=ixomax3-1
5171 if(slab_uniform) then
5172 dx1x2=dxlevel(1)/dxlevel(2)
5173 dx1x3=dxlevel(1)/dxlevel(3)
5174 do ix1=ixfmax1,ixfmin1,-1
5175 w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
5176 w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)) &
5177 +dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))-&
5178 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2))) &
5179 +dx1x3*(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))-&
5180 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))
5181 end do
5182 else
5183 do ix1=ixfmax1,ixfmin1,-1
5184 w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
5185 ( (w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))+&
5186 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)))*&
5187 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)&
5188 +(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))+&
5189 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2)))*&
5190 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,2)&
5191 -(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2))+&
5192 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2)))*&
5193 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,2)&
5194 +(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))+&
5195 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3)))*&
5196 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,3)&
5197 -(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3))+&
5198 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))*&
5199 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,3) )&
5200 /block%surfaceC(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)-&
5201 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))
5202 end do
5203 end if
5204 }
5205 case(2)
5206 {^iftwod
5207 ixfmin1=ixomin1-1
5208 ixfmax1=ixomax1-1
5209 ixfmin2=ixomin2+1
5210 ixfmax2=ixomax2-1
5211 if(slab_uniform) then
5212 dx1x2=dxlevel(1)/dxlevel(2)
5213 do ix1=ixfmin1,ixfmax1
5214 w(ix1+1,ixfmin2:ixfmax2,mag(1))=w(ix1-1,ixfmin2:ixfmax2,mag(1)) &
5215 -dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))-&
5216 w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))
5217 enddo
5218 else
5219 do ix1=ixfmin1,ixfmax1
5220 w(ix1+1,ixfmin2:ixfmax2,mag(1))=( (w(ix1-1,ixfmin2:ixfmax2,mag(1))+&
5221 w(ix1,ixfmin2:ixfmax2,mag(1)))*block%surfaceC(ix1-1,ixfmin2:ixfmax2,1)&
5222 -(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))+w(ix1,ixfmin2:ixfmax2,mag(2)))*&
5223 block%surfaceC(ix1,ixfmin2:ixfmax2,2)&
5224 +(w(ix1,ixfmin2:ixfmax2,mag(2))+w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))*&
5225 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,2) )&
5226 /block%surfaceC(ix1,ixfmin2:ixfmax2,1)-w(ix1,ixfmin2:ixfmax2,mag(1))
5227 end do
5228 end if
5229 }
5230 {^ifthreed
5231 ixfmin1=ixomin1-1
5232 ixfmax1=ixomax1-1
5233 ixfmin2=ixomin2+1
5234 ixfmax2=ixomax2-1
5235 ixfmin3=ixomin3+1
5236 ixfmax3=ixomax3-1
5237 if(slab_uniform) then
5238 dx1x2=dxlevel(1)/dxlevel(2)
5239 dx1x3=dxlevel(1)/dxlevel(3)
5240 do ix1=ixfmin1,ixfmax1
5241 w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
5242 w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)) &
5243 -dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))-&
5244 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2))) &
5245 -dx1x3*(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))-&
5246 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))
5247 end do
5248 else
5249 do ix1=ixfmin1,ixfmax1
5250 w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
5251 ( (w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))+&
5252 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)))*&
5253 block%surfaceC(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)&
5254 -(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))+&
5255 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2)))*&
5256 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,2)&
5257 +(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2))+&
5258 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2)))*&
5259 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,2)&
5260 -(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))+&
5261 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3)))*&
5262 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,3)&
5263 +(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3))+&
5264 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))*&
5265 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,3) )&
5266 /block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)-&
5267 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))
5268 end do
5269 end if
5270 }
5271 case(3)
5272 {^iftwod
5273 ixfmin1=ixomin1+1
5274 ixfmax1=ixomax1-1
5275 ixfmin2=ixomin2+1
5276 ixfmax2=ixomax2+1
5277 if(slab_uniform) then
5278 dx2x1=dxlevel(2)/dxlevel(1)
5279 do ix2=ixfmax2,ixfmin2,-1
5280 w(ixfmin1:ixfmax1,ix2-1,mag(2))=w(ixfmin1:ixfmax1,ix2+1,mag(2)) &
5281 +dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))-&
5282 w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))
5283 enddo
5284 else
5285 do ix2=ixfmax2,ixfmin2,-1
5286 w(ixfmin1:ixfmax1,ix2-1,mag(2))=( (w(ixfmin1:ixfmax1,ix2+1,mag(2))+&
5287 w(ixfmin1:ixfmax1,ix2,mag(2)))*block%surfaceC(ixfmin1:ixfmax1,ix2,2)&
5288 +(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))+w(ixfmin1:ixfmax1,ix2,mag(1)))*&
5289 block%surfaceC(ixfmin1:ixfmax1,ix2,1)&
5290 -(w(ixfmin1:ixfmax1,ix2,mag(1))+w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))*&
5291 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,1) )&
5292 /block%surfaceC(ixfmin1:ixfmax1,ix2-1,2)-w(ixfmin1:ixfmax1,ix2,mag(2))
5293 end do
5294 end if
5295 }
5296 {^ifthreed
5297 ixfmin1=ixomin1+1
5298 ixfmax1=ixomax1-1
5299 ixfmin3=ixomin3+1
5300 ixfmax3=ixomax3-1
5301 ixfmin2=ixomin2+1
5302 ixfmax2=ixomax2+1
5303 if(slab_uniform) then
5304 dx2x1=dxlevel(2)/dxlevel(1)
5305 dx2x3=dxlevel(2)/dxlevel(3)
5306 do ix2=ixfmax2,ixfmin2,-1
5307 w(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,mag(2))=w(ixfmin1:ixfmax1,&
5308 ix2+1,ixfmin3:ixfmax3,mag(2)) &
5309 +dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))-&
5310 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1))) &
5311 +dx2x3*(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))-&
5312 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))
5313 end do
5314 else
5315 do ix2=ixfmax2,ixfmin2,-1
5316 w(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,mag(2))=&
5317 ( (w(ixfmin1:ixfmax1,ix2+1,ixfmin3:ixfmax3,mag(2))+&
5318 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2)))*&
5319 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,2)&
5320 +(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))+&
5321 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1)))*&
5322 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,1)&
5323 -(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1))+&
5324 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1)))*&
5325 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,1)&
5326 +(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))+&
5327 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3)))*&
5328 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,3)&
5329 -(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3))+&
5330 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))*&
5331 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,3) )&
5332 /block%surfaceC(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,2)-&
5333 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2))
5334 end do
5335 end if
5336 }
5337 case(4)
5338 {^iftwod
5339 ixfmin1=ixomin1+1
5340 ixfmax1=ixomax1-1
5341 ixfmin2=ixomin2-1
5342 ixfmax2=ixomax2-1
5343 if(slab_uniform) then
5344 dx2x1=dxlevel(2)/dxlevel(1)
5345 do ix2=ixfmin2,ixfmax2
5346 w(ixfmin1:ixfmax1,ix2+1,mag(2))=w(ixfmin1:ixfmax1,ix2-1,mag(2)) &
5347 -dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))-&
5348 w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))
5349 end do
5350 else
5351 do ix2=ixfmin2,ixfmax2
5352 w(ixfmin1:ixfmax1,ix2+1,mag(2))=( (w(ixfmin1:ixfmax1,ix2-1,mag(2))+&
5353 w(ixfmin1:ixfmax1,ix2,mag(2)))*block%surfaceC(ixfmin1:ixfmax1,ix2-1,2)&
5354 -(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))+w(ixfmin1:ixfmax1,ix2,mag(1)))*&
5355 block%surfaceC(ixfmin1:ixfmax1,ix2,1)&
5356 +(w(ixfmin1:ixfmax1,ix2,mag(1))+w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))*&
5357 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,1) )&
5358 /block%surfaceC(ixfmin1:ixfmax1,ix2,2)-w(ixfmin1:ixfmax1,ix2,mag(2))
5359 end do
5360 end if
5361 }
5362 {^ifthreed
5363 ixfmin1=ixomin1+1
5364 ixfmax1=ixomax1-1
5365 ixfmin3=ixomin3+1
5366 ixfmax3=ixomax3-1
5367 ixfmin2=ixomin2-1
5368 ixfmax2=ixomax2-1
5369 if(slab_uniform) then
5370 dx2x1=dxlevel(2)/dxlevel(1)
5371 dx2x3=dxlevel(2)/dxlevel(3)
5372 do ix2=ixfmin2,ixfmax2
5373 w(ixfmin1:ixfmax1,ix2+1,ixfmin3:ixfmax3,mag(2))=w(ixfmin1:ixfmax1,&
5374 ix2-1,ixfmin3:ixfmax3,mag(2)) &
5375 -dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))-&
5376 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1))) &
5377 -dx2x3*(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))-&
5378 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))
5379 end do
5380 else
5381 do ix2=ixfmin2,ixfmax2
5382 w(ixfmin1:ixfmax1,ix2+1,ixfmin3:ixfmax3,mag(2))=&
5383 ( (w(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,mag(2))+&
5384 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2)))*&
5385 block%surfaceC(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,2)&
5386 -(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))+&
5387 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1)))*&
5388 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,1)&
5389 +(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1))+&
5390 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1)))*&
5391 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,1)&
5392 -(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))+&
5393 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3)))*&
5394 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,3)&
5395 +(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3))+&
5396 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))*&
5397 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,3) )&
5398 /block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,2)-&
5399 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2))
5400 end do
5401 end if
5402 }
5403 {^ifthreed
5404 case(5)
5405 ixfmin1=ixomin1+1
5406 ixfmax1=ixomax1-1
5407 ixfmin2=ixomin2+1
5408 ixfmax2=ixomax2-1
5409 ixfmin3=ixomin3+1
5410 ixfmax3=ixomax3+1
5411 if(slab_uniform) then
5412 dx3x1=dxlevel(3)/dxlevel(1)
5413 dx3x2=dxlevel(3)/dxlevel(2)
5414 do ix3=ixfmax3,ixfmin3,-1
5415 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,mag(3))=w(ixfmin1:ixfmax1,&
5416 ixfmin2:ixfmax2,ix3+1,mag(3)) &
5417 +dx3x1*(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))-&
5418 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1))) &
5419 +dx3x2*(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))-&
5420 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))
5421 end do
5422 else
5423 do ix3=ixfmax3,ixfmin3,-1
5424 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,mag(3))=&
5425 ( (w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3+1,mag(3))+&
5426 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3)))*&
5427 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,3)&
5428 +(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))+&
5429 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1)))*&
5430 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,1)&
5431 -(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1))+&
5432 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1)))*&
5433 block%surfaceC(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,1)&
5434 +(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))+&
5435 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2)))*&
5436 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,2)&
5437 -(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2))+&
5438 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))*&
5439 block%surfaceC(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,2) )&
5440 /block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,3)-&
5441 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3))
5442 end do
5443 end if
5444 case(6)
5445 ixfmin1=ixomin1+1
5446 ixfmax1=ixomax1-1
5447 ixfmin2=ixomin2+1
5448 ixfmax2=ixomax2-1
5449 ixfmin3=ixomin3-1
5450 ixfmax3=ixomax3-1
5451 if(slab_uniform) then
5452 dx3x1=dxlevel(3)/dxlevel(1)
5453 dx3x2=dxlevel(3)/dxlevel(2)
5454 do ix3=ixfmin3,ixfmax3
5455 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3+1,mag(3))=w(ixfmin1:ixfmax1,&
5456 ixfmin2:ixfmax2,ix3-1,mag(3)) &
5457 -dx3x1*(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))-&
5458 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1))) &
5459 -dx3x2*(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))-&
5460 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))
5461 end do
5462 else
5463 do ix3=ixfmin3,ixfmax3
5464 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3+1,mag(3))=&
5465 ( (w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,mag(3))+&
5466 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3)))*&
5467 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,3)&
5468 -(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))+&
5469 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1)))*&
5470 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,1)&
5471 +(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1))+&
5472 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1)))*&
5473 block%surfaceC(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,1)&
5474 -(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))+&
5475 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2)))*&
5476 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,2)&
5477 +(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2))+&
5478 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))*&
5479 block%surfaceC(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,2) )&
5480 /block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,3)-&
5481 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3))
5482 end do
5483 end if
5484 }
5485 case default
5486 call mpistop("Special boundary is not defined for this region")
5487 end select
5488
5489 end subroutine fixdivb_boundary
5490
5491 {^nooned
5492 subroutine twofl_clean_divb_multigrid(qdt, qt, active)
5493 use mod_forest
5496 use mod_geometry
5497
5498 double precision, intent(in) :: qdt !< Current time step
5499 double precision, intent(in) :: qt !< Current time
5500 logical, intent(inout) :: active !< Output if the source is active
5501 integer :: iigrid, igrid, id
5502 integer :: n, nc, lvl, ix^l, ixc^l, idim
5503 type(tree_node), pointer :: pnode
5504 double precision :: tmp(ixg^t), grad(ixg^t, ndim)
5505 double precision :: res
5506 double precision, parameter :: max_residual = 1d-3
5507 double precision, parameter :: residual_reduction = 1d-10
5508 integer, parameter :: max_its = 50
5509 double precision :: residual_it(max_its), max_divb
5510
5511 mg%operator_type = mg_laplacian
5512
5513 ! Set boundary conditions
5514 do n = 1, 2*ndim
5515 idim = (n+1)/2
5516 select case (typeboundary(mag(idim), n))
5517 case (bc_symm)
5518 ! d/dx B = 0, take phi = 0
5519 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
5520 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
5521 case (bc_asymm)
5522 ! B = 0, so grad(phi) = 0
5523 mg%bc(n, mg_iphi)%bc_type = mg_bc_neumann
5524 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
5525 case (bc_cont)
5526 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
5527 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
5528 case (bc_special)
5529 ! Assume Dirichlet boundary conditions, derivative zero
5530 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
5531 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
5532 case (bc_periodic)
5533 ! Nothing to do here
5534 case default
5535 print *, "divb_multigrid warning: unknown b.c.: ", &
5536 typeboundary(mag(idim), n)
5537 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
5538 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
5539 end select
5540 end do
5541
5542 ix^l=ixm^ll^ladd1;
5543 max_divb = 0.0d0
5544
5545 ! Store divergence of B as right-hand side
5546 do iigrid = 1, igridstail
5547 igrid = igrids(iigrid);
5548 pnode => igrid_to_node(igrid, mype)%node
5549 id = pnode%id
5550 lvl = mg%boxes(id)%lvl
5551 nc = mg%box_size_lvl(lvl)
5552
5553 ! Geometry subroutines expect this to be set
5554 block => ps(igrid)
5555 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
5556
5557 call get_divb(ps(igrid)%w(ixg^t, 1:nw), ixg^ll, ixm^ll, tmp, &
5559 mg%boxes(id)%cc({1:nc}, mg_irhs) = tmp(ixm^t)
5560 max_divb = max(max_divb, maxval(abs(tmp(ixm^t))))
5561 end do
5562
5563 ! Solve laplacian(phi) = divB
5564 if(stagger_grid) then
5565 call mpi_allreduce(mpi_in_place, max_divb, 1, mpi_double_precision, &
5566 mpi_max, icomm, ierrmpi)
5567
5568 if (mype == 0) print *, "Performing multigrid divB cleaning"
5569 if (mype == 0) print *, "iteration vs residual"
5570 ! Solve laplacian(phi) = divB
5571 do n = 1, max_its
5572 call mg_fas_fmg(mg, n>1, max_res=residual_it(n))
5573 if (mype == 0) write(*, "(I4,E11.3)") n, residual_it(n)
5574 if (residual_it(n) < residual_reduction * max_divb) exit
5575 end do
5576 if (mype == 0 .and. n > max_its) then
5577 print *, "divb_multigrid warning: not fully converged"
5578 print *, "current amplitude of divb: ", residual_it(max_its)
5579 print *, "multigrid smallest grid: ", &
5580 mg%domain_size_lvl(:, mg%lowest_lvl)
5581 print *, "note: smallest grid ideally has <= 8 cells"
5582 print *, "multigrid dx/dy/dz ratio: ", mg%dr(:, 1)/mg%dr(1, 1)
5583 print *, "note: dx/dy/dz should be similar"
5584 end if
5585 else
5586 do n = 1, max_its
5587 call mg_fas_vcycle(mg, max_res=res)
5588 if (res < max_residual) exit
5589 end do
5590 if (res > max_residual) call mpistop("divb_multigrid: no convergence")
5591 end if
5592
5593
5594 ! Correct the magnetic field
5595 do iigrid = 1, igridstail
5596 igrid = igrids(iigrid);
5597 pnode => igrid_to_node(igrid, mype)%node
5598 id = pnode%id
5599
5600 ! Geometry subroutines expect this to be set
5601 block => ps(igrid)
5602 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
5603
5604 ! Compute the gradient of phi
5605 tmp(ix^s) = mg%boxes(id)%cc({:,}, mg_iphi)
5606
5607 if(stagger_grid) then
5608 do idim =1, ndim
5609 ixcmin^d=ixmlo^d-kr(idim,^d);
5610 ixcmax^d=ixmhi^d;
5611 call gradientf(tmp,ps(igrid)%x,ixg^ll,ixc^l,idim,grad(ixg^t,idim))
5612 ! Apply the correction B* = B - gradient(phi)
5613 ps(igrid)%ws(ixc^s,idim)=ps(igrid)%ws(ixc^s,idim)-grad(ixc^s,idim)
5614 end do
5615 ! store cell-center magnetic energy
5616 tmp(ixm^t) = sum(ps(igrid)%w(ixm^t, mag(1:ndim))**2, dim=ndim+1)
5617 ! change cell-center magnetic field
5618 call twofl_face_to_center(ixm^ll,ps(igrid))
5619 else
5620 do idim = 1, ndim
5621 call gradient(tmp,ixg^ll,ixm^ll,idim,grad(ixg^t, idim))
5622 end do
5623 ! store cell-center magnetic energy
5624 tmp(ixm^t) = sum(ps(igrid)%w(ixm^t, mag(1:ndim))**2, dim=ndim+1)
5625 ! Apply the correction B* = B - gradient(phi)
5626 ps(igrid)%w(ixm^t, mag(1:ndim)) = &
5627 ps(igrid)%w(ixm^t, mag(1:ndim)) - grad(ixm^t, :)
5628 end if
5629
5630 if(phys_total_energy) then
5631 ! Determine magnetic energy difference
5632 tmp(ixm^t) = 0.5_dp * (sum(ps(igrid)%w(ixm^t, &
5633 mag(1:ndim))**2, dim=ndim+1) - tmp(ixm^t))
5634 ! Keep thermal pressure the same
5635 ps(igrid)%w(ixm^t, e_c_) = ps(igrid)%w(ixm^t, e_c_) + tmp(ixm^t)
5636 end if
5637 end do
5638
5639 active = .true.
5640
5641 end subroutine twofl_clean_divb_multigrid
5642 }
5643 !> get electric field though averaging neighors to update faces in CT
5644 subroutine twofl_update_faces_average(ixI^L,ixO^L,qt,qdt,wp,fC,fE,sCT,s,vcts)
5646 use mod_usr_methods
5647
5648 integer, intent(in) :: ixi^l, ixo^l
5649 double precision, intent(in) :: qt,qdt
5650 ! cell-center primitive variables
5651 double precision, intent(in) :: wp(ixi^s,1:nw)
5652 type(state) :: sct, s
5653 type(ct_velocity) :: vcts
5654 double precision, intent(in) :: fc(ixi^s,1:nwflux,1:ndim)
5655 double precision, intent(inout) :: fe(ixi^s,sdim:3)
5656
5657 double precision :: circ(ixi^s,1:ndim)
5658 ! non-ideal electric field on cell edges
5659 double precision, dimension(ixI^S,sdim:3) :: e_resi, e_ambi
5660 integer :: ix^d,ixc^l,ixa^l,i1kr^d,i2kr^d
5661 integer :: idim1,idim2,idir,iwdim1,iwdim2
5662
5663 associate(bfaces=>s%ws,x=>s%x)
5664
5665 ! Calculate contribution to FEM of each edge,
5666 ! that is, estimate value of line integral of
5667 ! electric field in the positive idir direction.
5668
5669 ! if there is resistivity, get eta J
5670 if(twofl_eta/=zero) call get_resistive_electric_field(ixi^l,ixo^l,wp,sct,s,e_resi)
5671
5672 do idim1=1,ndim
5673 iwdim1 = mag(idim1)
5674 i1kr^d=kr(idim1,^d);
5675 do idim2=1,ndim
5676 iwdim2 = mag(idim2)
5677 i2kr^d=kr(idim2,^d);
5678 do idir=sdim,3! Direction of line integral
5679 ! Allow only even permutations
5680 if (lvc(idim1,idim2,idir)==1) then
5681 ixcmax^d=ixomax^d;
5682 ixcmin^d=ixomin^d+kr(idir,^d)-1;
5683 ! average cell-face electric field to cell edges
5684 {do ix^db=ixcmin^db,ixcmax^db\}
5685 fe(ix^d,idir)=quarter*&
5686 (fc(ix^d,iwdim1,idim2)+fc({ix^d+i1kr^d},iwdim1,idim2)&
5687 -fc(ix^d,iwdim2,idim1)-fc({ix^d+i2kr^d},iwdim2,idim1))
5688 ! add resistive electric field at cell edges E=-vxB+eta J
5689 if(twofl_eta/=zero) fe(ix^d,idir)=fe(ix^d,idir)+e_resi(ix^d,idir)
5690 ! times time step and edge length
5691 fe(ix^d,idir)=fe(ix^d,idir)*qdt*s%dsC(ix^d,idir)
5692 {end do\}
5693 end if
5694 end do
5695 end do
5696 end do
5697
5698 ! allow user to change inductive electric field, especially for boundary driven applications
5699 if(associated(usr_set_electric_field)) &
5700 call usr_set_electric_field(ixi^l,ixo^l,qt,qdt,fe,sct)
5701
5702 circ(ixi^s,1:ndim)=zero
5703
5704 ! Calculate circulation on each face
5705 do idim1=1,ndim ! Coordinate perpendicular to face
5706 ixcmax^d=ixomax^d;
5707 ixcmin^d=ixomin^d-kr(idim1,^d);
5708 do idim2=1,ndim
5709 ixa^l=ixc^l-kr(idim2,^d);
5710 do idir=sdim,3 ! Direction of line integral
5711 ! Assemble indices
5712 if(lvc(idim1,idim2,idir)==1) then
5713 ! Add line integrals in direction idir
5714 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
5715 +(fe(ixc^s,idir)&
5716 -fe(ixa^s,idir))
5717 else if(lvc(idim1,idim2,idir)==-1) then
5718 ! Add line integrals in direction idir
5719 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
5720 -(fe(ixc^s,idir)&
5721 -fe(ixa^s,idir))
5722 end if
5723 end do
5724 end do
5725 {do ix^db=ixcmin^db,ixcmax^db\}
5726 ! Divide by the area of the face to get dB/dt
5727 if(s%surfaceC(ix^d,idim1) > smalldouble) then
5728 ! Time update cell-face magnetic field component
5729 bfaces(ix^d,idim1)=bfaces(ix^d,idim1)-circ(ix^d,idim1)/s%surfaceC(ix^d,idim1)
5730 end if
5731 {end do\}
5732 end do
5733
5734 end associate
5735
5736 end subroutine twofl_update_faces_average
5737
5738 !> update faces using UCT contact mode by Gardiner and Stone 2005 JCP 205, 509
5739 subroutine twofl_update_faces_contact(ixI^L,ixO^L,qt,qdt,wp,fC,fE,sCT,s,vcts)
5741 use mod_usr_methods
5742 use mod_geometry
5743
5744 integer, intent(in) :: ixi^l, ixo^l
5745 double precision, intent(in) :: qt, qdt
5746 ! cell-center primitive variables
5747 double precision, intent(in) :: wp(ixi^s,1:nw)
5748 type(state) :: sct, s
5749 type(ct_velocity) :: vcts
5750 double precision, intent(in) :: fc(ixi^s,1:nwflux,1:ndim)
5751 double precision, intent(inout) :: fe(ixi^s,sdim:3)
5752
5753 double precision :: circ(ixi^s,1:ndim)
5754 ! electric field at cell centers
5755 double precision :: ecc(ixi^s,sdim:3)
5756 double precision :: ein(ixi^s,sdim:3)
5757 ! gradient of E at left and right side of a cell face
5758 double precision :: el(ixi^s),er(ixi^s)
5759 ! gradient of E at left and right side of a cell corner
5760 double precision :: elc,erc
5761 ! non-ideal electric field on cell edges
5762 double precision, dimension(ixI^S,sdim:3) :: e_resi, e_ambi
5763 ! current on cell edges
5764 double precision :: jce(ixi^s,sdim:3)
5765 ! location at cell faces
5766 double precision :: xs(ixgs^t,1:ndim)
5767 double precision :: gradi(ixgs^t)
5768 integer :: ixc^l,ixa^l
5769 integer :: idim1,idim2,idir,iwdim1,iwdim2,ix^d,i1kr^d,i2kr^d
5770
5771 associate(bfaces=>s%ws,x=>s%x,w=>s%w,vnorm=>vcts%vnorm,wcts=>sct%ws)
5772
5773 ! if there is resistivity, get eta J
5774 if(twofl_eta/=zero) call get_resistive_electric_field(ixi^l,ixo^l,wp,sct,s,e_resi)
5775
5776 if(b0field) then
5777 {do ix^db=iximin^db,iximax^db\}
5778 ! Calculate electric field at cell centers
5779 {^ifthreed
5780 ecc(ix^d,1)=(wp(ix^d,b2_)+block%B0(ix^d,2,0))*wp(ix^d,m3_)-(wp(ix^d,b3_)+block%B0(ix^d,3,0))*wp(ix^d,m2_)
5781 ecc(ix^d,2)=(wp(ix^d,b3_)+block%B0(ix^d,3,0))*wp(ix^d,m1_)-(wp(ix^d,b1_)+block%B0(ix^d,1,0))*wp(ix^d,m3_)
5782 ecc(ix^d,3)=(wp(ix^d,b1_)+block%B0(ix^d,1,0))*wp(ix^d,m2_)-(wp(ix^d,b2_)+block%B0(ix^d,2,0))*wp(ix^d,m1_)
5783 }
5784 {^iftwod
5785 ecc(ix^d,3)=wp(ix^d,b1_)*wp(ix^d,m2_)-wp(ix^d,b2_)*wp(ix^d,m1_)
5786 }
5787 {^ifoned
5788 ecc(ix^d,3)=0.d0
5789 }
5790 {end do\}
5791 else
5792 {do ix^db=iximin^db,iximax^db\}
5793 ! Calculate electric field at cell centers
5794 {^ifthreed
5795 ecc(ix^d,1)=wp(ix^d,b2_)*wp(ix^d,m3_)-wp(ix^d,b3_)*wp(ix^d,m2_)
5796 ecc(ix^d,2)=wp(ix^d,b3_)*wp(ix^d,m1_)-wp(ix^d,b1_)*wp(ix^d,m3_)
5797 ecc(ix^d,3)=wp(ix^d,b1_)*wp(ix^d,m2_)-wp(ix^d,b2_)*wp(ix^d,m1_)
5798 }
5799 {^iftwod
5800 ecc(ix^d,3)=wp(ix^d,b1_)*wp(ix^d,m2_)-wp(ix^d,b2_)*wp(ix^d,m1_)
5801 }
5802 {^ifoned
5803 ecc(ix^d,3)=0.d0
5804 }
5805 {end do\}
5806 end if
5807
5808 ! Calculate contribution to FEM of each edge,
5809 ! that is, estimate value of line integral of
5810 ! electric field in the positive idir direction.
5811 ! evaluate electric field along cell edges according to equation (41)
5812 do idim1=1,ndim
5813 iwdim1 = mag(idim1)
5814 i1kr^d=kr(idim1,^d);
5815 do idim2=1,ndim
5816 iwdim2 = mag(idim2)
5817 i2kr^d=kr(idim2,^d);
5818 do idir=sdim,3 ! Direction of line integral
5819 ! Allow only even permutations
5820 if (lvc(idim1,idim2,idir)==1) then
5821 ixcmax^d=ixomax^d;
5822 ixcmin^d=ixomin^d+kr(idir,^d)-1;
5823 ! Assemble indices
5824 ! average cell-face electric field to cell edges
5825 {do ix^db=ixcmin^db,ixcmax^db\}
5826 fe(ix^d,idir)=quarter*&
5827 (fc(ix^d,iwdim1,idim2)+fc({ix^d+i1kr^d},iwdim1,idim2)&
5828 -fc(ix^d,iwdim2,idim1)-fc({ix^d+i2kr^d},iwdim2,idim1))
5829 {end do\}
5830 ! add slope in idim2 direction from equation (50)
5831 ixamin^d=ixcmin^d;
5832 ixamax^d=ixcmax^d+i1kr^d;
5833 {do ix^db=ixamin^db,ixamax^db\}
5834 el(ix^d)=fc(ix^d,iwdim1,idim2)-ecc(ix^d,idir)
5835 er(ix^d)=fc(ix^d,iwdim1,idim2)-ecc({ix^d+i2kr^d},idir)
5836 {end do\}
5837 {!dir$ ivdep
5838 do ix^db=ixcmin^db,ixcmax^db\}
5839 if(vnorm(ix^d,idim1)>0.d0) then
5840 elc=el(ix^d)
5841 else if(vnorm(ix^d,idim1)<0.d0) then
5842 elc=el({ix^d+i1kr^d})
5843 else
5844 elc=0.5d0*(el(ix^d)+el({ix^d+i1kr^d}))
5845 end if
5846 if(vnorm({ix^d+i2kr^d},idim1)>0.d0) then
5847 erc=er(ix^d)
5848 else if(vnorm({ix^d+i2kr^d},idim1)<0.d0) then
5849 erc=er({ix^d+i1kr^d})
5850 else
5851 erc=0.5d0*(er(ix^d)+er({ix^d+i1kr^d}))
5852 end if
5853 fe(ix^d,idir)=fe(ix^d,idir)+0.25d0*(elc+erc)
5854 {end do\}
5855
5856 ! add slope in idim1 direction from equation (50)
5857 ixamin^d=ixcmin^d;
5858 ixamax^d=ixcmax^d+i2kr^d;
5859 {do ix^db=ixamin^db,ixamax^db\}
5860 el(ix^d)=-fc(ix^d,iwdim2,idim1)-ecc(ix^d,idir)
5861 er(ix^d)=-fc(ix^d,iwdim2,idim1)-ecc({ix^d+i1kr^d},idir)
5862 {end do\}
5863 {!dir$ ivdep
5864 do ix^db=ixcmin^db,ixcmax^db\}
5865 if(vnorm(ix^d,idim2)>0.d0) then
5866 elc=el(ix^d)
5867 else if(vnorm(ix^d,idim2)<0.d0) then
5868 elc=el({ix^d+i2kr^d})
5869 else
5870 elc=0.5d0*(el(ix^d)+el({ix^d+i2kr^d}))
5871 end if
5872 if(vnorm({ix^d+i1kr^d},idim2)>0.d0) then
5873 erc=er(ix^d)
5874 else if(vnorm({ix^d+i1kr^d},idim2)<0.d0) then
5875 erc=er({ix^d+i2kr^d})
5876 else
5877 erc=0.5d0*(er(ix^d)+er({ix^d+i2kr^d}))
5878 end if
5879 fe(ix^d,idir)=fe(ix^d,idir)+0.25d0*(elc+erc)
5880 ! add resistive electric field at cell edges E=-vxB+eta J
5881 if(twofl_eta/=zero) fe(ix^d,idir)=fe(ix^d,idir)+e_resi(ix^d,idir)
5882 ! times time step and edge length
5883 fe(ix^d,idir)=fe(ix^d,idir)*qdt*s%dsC(ix^d,idir)
5884 {end do\}
5885 end if
5886 end do
5887 end do
5888 end do
5889
5890 ! allow user to change inductive electric field, especially for boundary driven applications
5891 if(associated(usr_set_electric_field)) &
5892 call usr_set_electric_field(ixi^l,ixo^l,qt,qdt,fe,sct)
5893
5894 circ(ixi^s,1:ndim)=zero
5895
5896 ! Calculate circulation on each face
5897 do idim1=1,ndim ! Coordinate perpendicular to face
5898 ixcmax^d=ixomax^d;
5899 ixcmin^d=ixomin^d-kr(idim1,^d);
5900 do idim2=1,ndim
5901 ixa^l=ixc^l-kr(idim2,^d);
5902 do idir=sdim,3 ! Direction of line integral
5903 ! Assemble indices
5904 if(lvc(idim1,idim2,idir)==1) then
5905 ! Add line integrals in direction idir
5906 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
5907 +(fe(ixc^s,idir)&
5908 -fe(ixa^s,idir))
5909 else if(lvc(idim1,idim2,idir)==-1) then
5910 ! Add line integrals in direction idir
5911 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
5912 -(fe(ixc^s,idir)&
5913 -fe(ixa^s,idir))
5914 end if
5915 end do
5916 end do
5917 {do ix^db=ixcmin^db,ixcmax^db\}
5918 ! Divide by the area of the face to get dB/dt
5919 if(s%surfaceC(ix^d,idim1) > smalldouble) then
5920 ! Time update cell-face magnetic field component
5921 bfaces(ix^d,idim1)=bfaces(ix^d,idim1)-circ(ix^d,idim1)/s%surfaceC(ix^d,idim1)
5922 end if
5923 {end do\}
5924 end do
5925
5926 end associate
5927
5928 end subroutine twofl_update_faces_contact
5929
5930 !> update faces
5931 subroutine twofl_update_faces_hll(ixI^L,ixO^L,qt,qdt,wp,fC,fE,sCT,s,vcts)
5933 use mod_usr_methods
5935
5936 integer, intent(in) :: ixi^l, ixo^l
5937 double precision, intent(in) :: qt, qdt
5938 ! cell-center primitive variables
5939 double precision, intent(in) :: wp(ixi^s,1:nw)
5940 type(state) :: sct, s
5941 type(ct_velocity) :: vcts
5942 double precision, intent(in) :: fc(ixi^s,1:nwflux,1:ndim)
5943 double precision, intent(inout) :: fe(ixi^s,sdim:3)
5944
5945 double precision :: vtill(ixi^s,2)
5946 double precision :: vtilr(ixi^s,2)
5947 double precision :: bfacetot(ixi^s,ndim)
5948 double precision :: btill(ixi^s,ndim)
5949 double precision :: btilr(ixi^s,ndim)
5950 double precision :: cp(ixi^s,2)
5951 double precision :: cm(ixi^s,2)
5952 double precision :: circ(ixi^s,1:ndim)
5953 ! non-ideal electric field on cell edges
5954 double precision, dimension(ixI^S,sdim:3) :: e_resi, e_ambi
5955 integer :: hxc^l,ixc^l,ixcp^l,jxc^l,ixcm^l
5956 integer :: idim1,idim2,idir,ix^d
5957
5958 associate(bfaces=>s%ws,bfacesct=>sct%ws,x=>s%x,vbarc=>vcts%vbarC,cbarmin=>vcts%cbarmin,&
5959 cbarmax=>vcts%cbarmax)
5960
5961 ! Calculate contribution to FEM of each edge,
5962 ! that is, estimate value of line integral of
5963 ! electric field in the positive idir direction.
5964
5965 ! Loop over components of electric field
5966
5967 ! idir: electric field component we need to calculate
5968 ! idim1: directions in which we already performed the reconstruction
5969 ! idim2: directions in which we perform the reconstruction
5970
5971 ! if there is resistivity, get eta J
5972 if(twofl_eta/=zero) call get_resistive_electric_field(ixi^l,ixo^l,wp,sct,s,e_resi)
5973
5974 do idir=sdim,3
5975 ! Indices
5976 ! idir: electric field component
5977 ! idim1: one surface
5978 ! idim2: the other surface
5979 ! cyclic permutation: idim1,idim2,idir=1,2,3
5980 ! Velocity components on the surface
5981 ! follow cyclic premutations:
5982 ! Sx(1),Sx(2)=y,z ; Sy(1),Sy(2)=z,x ; Sz(1),Sz(2)=x,y
5983
5984 ixcmax^d=ixomax^d;
5985 ixcmin^d=ixomin^d-1+kr(idir,^d);
5986
5987 ! Set indices and directions
5988 idim1=mod(idir,3)+1
5989 idim2=mod(idir+1,3)+1
5990
5991 jxc^l=ixc^l+kr(idim1,^d);
5992 ixcp^l=ixc^l+kr(idim2,^d);
5993
5994 ! Reconstruct transverse transport velocities
5995 call reconstruct(ixi^l,ixc^l,idim2,vbarc(ixi^s,idim1,1),&
5996 vtill(ixi^s,2),vtilr(ixi^s,2))
5997
5998 call reconstruct(ixi^l,ixc^l,idim1,vbarc(ixi^s,idim2,2),&
5999 vtill(ixi^s,1),vtilr(ixi^s,1))
6000
6001 ! Reconstruct magnetic fields
6002 ! Eventhough the arrays are larger, reconstruct works with
6003 ! the limits ixG.
6004 if(b0field) then
6005 bfacetot(ixi^s,idim1)=bfacesct(ixi^s,idim1)+block%B0(ixi^s,idim1,idim1)
6006 bfacetot(ixi^s,idim2)=bfacesct(ixi^s,idim2)+block%B0(ixi^s,idim2,idim2)
6007 else
6008 bfacetot(ixi^s,idim1)=bfacesct(ixi^s,idim1)
6009 bfacetot(ixi^s,idim2)=bfacesct(ixi^s,idim2)
6010 end if
6011 call reconstruct(ixi^l,ixc^l,idim2,bfacetot(ixi^s,idim1),&
6012 btill(ixi^s,idim1),btilr(ixi^s,idim1))
6013
6014 call reconstruct(ixi^l,ixc^l,idim1,bfacetot(ixi^s,idim2),&
6015 btill(ixi^s,idim2),btilr(ixi^s,idim2))
6016
6017 ! Take the maximum characteristic
6018
6019 cm(ixc^s,1)=max(cbarmin(ixcp^s,idim1),cbarmin(ixc^s,idim1))
6020 cp(ixc^s,1)=max(cbarmax(ixcp^s,idim1),cbarmax(ixc^s,idim1))
6021
6022 cm(ixc^s,2)=max(cbarmin(jxc^s,idim2),cbarmin(ixc^s,idim2))
6023 cp(ixc^s,2)=max(cbarmax(jxc^s,idim2),cbarmax(ixc^s,idim2))
6024
6025
6026 ! Calculate eletric field
6027 fe(ixc^s,idir)=-(cp(ixc^s,1)*vtill(ixc^s,1)*btill(ixc^s,idim2) &
6028 + cm(ixc^s,1)*vtilr(ixc^s,1)*btilr(ixc^s,idim2) &
6029 - cp(ixc^s,1)*cm(ixc^s,1)*(btilr(ixc^s,idim2)-btill(ixc^s,idim2)))&
6030 /(cp(ixc^s,1)+cm(ixc^s,1)) &
6031 +(cp(ixc^s,2)*vtill(ixc^s,2)*btill(ixc^s,idim1) &
6032 + cm(ixc^s,2)*vtilr(ixc^s,2)*btilr(ixc^s,idim1) &
6033 - cp(ixc^s,2)*cm(ixc^s,2)*(btilr(ixc^s,idim1)-btill(ixc^s,idim1)))&
6034 /(cp(ixc^s,2)+cm(ixc^s,2))
6035
6036 ! add resistive electric field at cell edges E=-vxB+eta J
6037 if(twofl_eta/=zero) fe(ixc^s,idir)=fe(ixc^s,idir)+e_resi(ixc^s,idir)
6038 fe(ixc^s,idir)=qdt*s%dsC(ixc^s,idir)*fe(ixc^s,idir)
6039
6040 if (.not.slab) then
6041 where(abs(x(ixc^s,r_)+half*dxlevel(r_)).lt.1.0d-9)
6042 fe(ixc^s,idir)=zero
6043 end where
6044 end if
6045
6046 end do
6047
6048 ! allow user to change inductive electric field, especially for boundary driven applications
6049 if(associated(usr_set_electric_field)) &
6050 call usr_set_electric_field(ixi^l,ixo^l,qt,qdt,fe,sct)
6051
6052 circ(ixi^s,1:ndim)=zero
6053
6054 ! Calculate circulation on each face: interal(fE dot dl)
6055 do idim1=1,ndim ! Coordinate perpendicular to face
6056 ixcmax^d=ixomax^d;
6057 ixcmin^d=ixomin^d-kr(idim1,^d);
6058 do idim2=1,ndim
6059 do idir=sdim,3 ! Direction of line integral
6060 ! Assemble indices
6061 if(lvc(idim1,idim2,idir)/=0) then
6062 hxc^l=ixc^l-kr(idim2,^d);
6063 ! Add line integrals in direction idir
6064 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
6065 +lvc(idim1,idim2,idir)&
6066 *(fe(ixc^s,idir)&
6067 -fe(hxc^s,idir))
6068 end if
6069 end do
6070 end do
6071 {do ix^db=ixcmin^db,ixcmax^db\}
6072 ! Divide by the area of the face to get dB/dt
6073 if(s%surfaceC(ix^d,idim1) > smalldouble) then
6074 ! Time update cell-face magnetic field component
6075 bfaces(ix^d,idim1)=bfaces(ix^d,idim1)-circ(ix^d,idim1)/s%surfaceC(ix^d,idim1)
6076 end if
6077 {end do\}
6078 end do
6079
6080 end associate
6081 end subroutine twofl_update_faces_hll
6082
6083 !> calculate eta J at cell edges
6084 subroutine get_resistive_electric_field(ixI^L,ixO^L,wp,sCT,s,jce)
6086 use mod_usr_methods
6087 use mod_geometry
6088
6089 integer, intent(in) :: ixi^l, ixo^l
6090 ! cell-center primitive variables
6091 double precision, intent(in) :: wp(ixi^s,1:nw)
6092 type(state), intent(in) :: sct, s
6093 ! current on cell edges
6094 double precision :: jce(ixi^s,sdim:3)
6095
6096 ! current on cell centers
6097 double precision :: jcc(ixi^s,7-2*ndir:3)
6098 ! location at cell faces
6099 double precision :: xs(ixgs^t,1:ndim)
6100 ! resistivity
6101 double precision :: eta(ixi^s)
6102 double precision :: gradi(ixgs^t)
6103 integer :: ix^d,ixc^l,ixa^l,ixb^l,idir,idirmin,idim1,idim2
6104
6105 associate(x=>s%x,dx=>s%dx,w=>s%w,wct=>sct%w,wcts=>sct%ws)
6106 ! calculate current density at cell edges
6107 jce=0.d0
6108 do idim1=1,ndim
6109 do idim2=1,ndim
6110 do idir=sdim,3
6111 if (lvc(idim1,idim2,idir)==0) cycle
6112 ixcmax^d=ixomax^d;
6113 ixcmin^d=ixomin^d+kr(idir,^d)-1;
6114 ixbmax^d=ixcmax^d-kr(idir,^d)+1;
6115 ixbmin^d=ixcmin^d;
6116 ! current at transverse faces
6117 xs(ixb^s,:)=x(ixb^s,:)
6118 xs(ixb^s,idim2)=x(ixb^s,idim2)+half*dx(ixb^s,idim2)
6119 call gradientf(wcts(ixgs^t,idim2),xs,ixgs^ll,ixc^l,idim1,gradi,2)
6120 if (lvc(idim1,idim2,idir)==1) then
6121 jce(ixc^s,idir)=jce(ixc^s,idir)+gradi(ixc^s)
6122 else
6123 jce(ixc^s,idir)=jce(ixc^s,idir)-gradi(ixc^s)
6124 end if
6125 end do
6126 end do
6127 end do
6128 ! get resistivity
6129 if(twofl_eta>zero)then
6130 jce(ixi^s,:)=jce(ixi^s,:)*twofl_eta
6131 else
6132 ixa^l=ixo^l^ladd1;
6133 call get_current(wct,ixi^l,ixa^l,idirmin,jcc)
6134 call usr_special_resistivity(wp,ixi^l,ixa^l,idirmin,x,jcc,eta)
6135 ! calcuate eta on cell edges
6136 do idir=sdim,3
6137 ixcmax^d=ixomax^d;
6138 ixcmin^d=ixomin^d+kr(idir,^d)-1;
6139 jcc(ixc^s,idir)=0.d0
6140 {do ix^db=0,1\}
6141 if({ ix^d==1 .and. ^d==idir | .or.}) cycle
6142 ixamin^d=ixcmin^d+ix^d;
6143 ixamax^d=ixcmax^d+ix^d;
6144 jcc(ixc^s,idir)=jcc(ixc^s,idir)+eta(ixa^s)
6145 {end do\}
6146 jcc(ixc^s,idir)=jcc(ixc^s,idir)*0.25d0
6147 jce(ixc^s,idir)=jce(ixc^s,idir)*jcc(ixc^s,idir)
6148 end do
6149 end if
6150
6151 end associate
6152 end subroutine get_resistive_electric_field
6153
6154 !> calculate cell-center values from face-center values
6155 subroutine twofl_face_to_center(ixO^L,s)
6157 ! Non-staggered interpolation range
6158 integer, intent(in) :: ixo^l
6159 type(state) :: s
6160
6161 integer :: fxo^l, gxo^l, hxo^l, jxo^l, kxo^l, idim
6162
6163 associate(w=>s%w, ws=>s%ws)
6164
6165 ! calculate cell-center values from face-center values in 2nd order
6166 do idim=1,ndim
6167 ! Displace index to the left
6168 ! Even if ixI^L is the full size of the w arrays, this is ok
6169 ! because the staggered arrays have an additional place to the left.
6170 hxo^l=ixo^l-kr(idim,^d);
6171 ! Interpolate to cell barycentre using arithmetic average
6172 ! This might be done better later, to make the method less diffusive.
6173 w(ixo^s,mag(idim))=half/s%surface(ixo^s,idim)*&
6174 (ws(ixo^s,idim)*s%surfaceC(ixo^s,idim)&
6175 +ws(hxo^s,idim)*s%surfaceC(hxo^s,idim))
6176 end do
6177
6178 ! calculate cell-center values from face-center values in 4th order
6179 !do idim=1,ndim
6180 ! gxO^L=ixO^L-2*kr(idim,^D);
6181 ! hxO^L=ixO^L-kr(idim,^D);
6182 ! jxO^L=ixO^L+kr(idim,^D);
6183
6184 ! ! Interpolate to cell barycentre using fourth order central formula
6185 ! w(ixO^S,mag(idim))=(0.0625d0/s%surface(ixO^S,idim))*&
6186 ! ( -ws(gxO^S,idim)*s%surfaceC(gxO^S,idim) &
6187 ! +9.0d0*ws(hxO^S,idim)*s%surfaceC(hxO^S,idim) &
6188 ! +9.0d0*ws(ixO^S,idim)*s%surfaceC(ixO^S,idim) &
6189 ! -ws(jxO^S,idim)*s%surfaceC(jxO^S,idim) )
6190 !end do
6191
6192 ! calculate cell-center values from face-center values in 6th order
6193 !do idim=1,ndim
6194 ! fxO^L=ixO^L-3*kr(idim,^D);
6195 ! gxO^L=ixO^L-2*kr(idim,^D);
6196 ! hxO^L=ixO^L-kr(idim,^D);
6197 ! jxO^L=ixO^L+kr(idim,^D);
6198 ! kxO^L=ixO^L+2*kr(idim,^D);
6199
6200 ! ! Interpolate to cell barycentre using sixth order central formula
6201 ! w(ixO^S,mag(idim))=(0.00390625d0/s%surface(ixO^S,idim))* &
6202 ! ( +3.0d0*ws(fxO^S,idim)*s%surfaceC(fxO^S,idim) &
6203 ! -25.0d0*ws(gxO^S,idim)*s%surfaceC(gxO^S,idim) &
6204 ! +150.0d0*ws(hxO^S,idim)*s%surfaceC(hxO^S,idim) &
6205 ! +150.0d0*ws(ixO^S,idim)*s%surfaceC(ixO^S,idim) &
6206 ! -25.0d0*ws(jxO^S,idim)*s%surfaceC(jxO^S,idim) &
6207 ! +3.0d0*ws(kxO^S,idim)*s%surfaceC(kxO^S,idim) )
6208 !end do
6209
6210 end associate
6211
6212 end subroutine twofl_face_to_center
6213
6214 !> calculate magnetic field from vector potential
6215 subroutine b_from_vector_potential(ixIs^L, ixI^L, ixO^L, ws, x)
6218
6219 integer, intent(in) :: ixis^l, ixi^l, ixo^l
6220 double precision, intent(inout) :: ws(ixis^s,1:nws)
6221 double precision, intent(in) :: x(ixi^s,1:ndim)
6222
6223 double precision :: adummy(ixis^s,1:3)
6224
6225 call b_from_vector_potentiala(ixis^l, ixi^l, ixo^l, ws, x, adummy)
6226
6227 end subroutine b_from_vector_potential
6228
6229 subroutine hyperdiffusivity_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
6232 integer, intent(in) :: ixi^l, ixo^l
6233 double precision, intent(in) :: w(ixi^s,1:nw)
6234 double precision, intent(in) :: x(ixi^s,1:ndim)
6235 double precision, intent(in) :: dx^d
6236 double precision, intent(inout) :: dtnew
6237
6238 double precision :: nu(ixi^s),tmp(ixi^s),rho(ixi^s),temp(ixi^s)
6239 double precision :: divv(ixi^s,1:ndim)
6240 double precision :: vel(ixi^s,1:ndir)
6241 double precision :: csound(ixi^s),csound_dim(ixi^s,1:ndim)
6242 double precision :: dxarr(ndim)
6243 double precision :: maxcoef
6244 integer :: ixoo^l, hxb^l, hx^l, ii, jj
6245
6246
6247 ^d&dxarr(^d)=dx^d;
6248 maxcoef = smalldouble
6249
6250 ! charges
6251 call twofl_get_v_c(w,x,ixi^l,ixi^l,vel)
6252 call get_rhoc_tot(w,x,ixi^l,ixi^l,rho)
6253 call twofl_get_csound2_c_from_conserved(w,x,ixi^l,ixi^l,csound)
6254 csound(ixi^s) = sqrt(csound(ixi^s)) + sqrt(twofl_mag_en_all(w,ixi^l,ixi^l) /rho(ixi^s))
6255 csound(ixi^s) = csound(ixi^s) + sqrt(sum(vel(ixi^s,1:ndir)**2 ,dim=ndim+1))
6256 do ii=1,ndim
6257 call div_vel_coeff(ixi^l, ixoo^l, vel, ii, divv(ixi^s,ii))
6258 hxmin^d=iximin^d+1;
6259 hxmax^d=iximax^d-1;
6260 hxb^l=hx^l-kr(ii,^d);
6261 csound_dim(hx^s,ii) = (csound(hxb^s)+csound(hx^s))/2d0
6262 enddo
6263 call twofl_get_temp_c_pert_from_etot(w, x, ixi^l, ixi^l, temp)
6264 do ii=1,ndim
6265 !TODO the following is copied
6266 !rho_c
6267 call hyp_coeff(ixi^l, ixoo^l, w(ixi^s,rho_c_), ii, tmp(ixi^s))
6268 nu(ixo^s) = c_hyp(rho_c_) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6269 c_shk(rho_c_) * (dxlevel(ii)**2) *divv(ixo^s,ii)
6270 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6271
6272 !TH c
6273 call hyp_coeff(ixi^l, ixoo^l, temp(ixi^s), ii, tmp(ixi^s))
6274 nu(ixo^s) = c_hyp(e_c_) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6275 c_shk(e_c_) * (dxlevel(ii)**2) *divv(ixo^s,ii)
6276 nu(ixo^s) = nu(ixo^s) * rho(ixo^s) * rc/(twofl_gamma-1d0)
6277 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6278
6279 !visc c
6280 do jj=1,ndir
6281 call hyp_coeff(ixi^l, ixoo^l, vel(ixi^s,jj), ii, tmp(ixi^s))
6282 nu(ixo^s) = c_hyp(mom_c(jj)) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6283 c_shk(mom_c(jj)) * (dxlevel(ii)**2) *divv(ixo^s,ii)
6284 nu(ixo^s) = nu(ixo^s) * rho(ixo^s)
6285 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6286 enddo
6287
6288 ! Ohmic
6289 do jj=1,ndir
6290 if(ii .ne. jj) then
6291 call hyp_coeff(ixi^l, ixoo^l, w(ixi^s,mag(jj)), ii, tmp(ixi^s))
6292 nu(ixo^s) = c_hyp(mag(jj)) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6293 c_shk(mag(jj)) * (dxlevel(ii)**2) *divv(ixo^s,ii)
6294 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6295 endif
6296 enddo
6297
6298 enddo
6299
6300 !TODO the following is copied, as charges, and as in add_source!
6301 ! neutrals
6302 call twofl_get_v_n(w,x,ixi^l,ixi^l,vel)
6303 call twofl_get_csound_n(w,x,ixi^l,ixi^l,csound)
6304 csound(ixi^s) = csound(ixi^s) + sqrt(sum(vel(ixi^s,1:ndir)**2 ,dim=ndim+1))
6305 do ii=1,ndim
6306 call div_vel_coeff(ixi^l, ixoo^l, vel, ii, divv(ixi^s,ii))
6307 hxmin^d=iximin^d+1;
6308 hxmax^d=iximax^d-1;
6309 hxb^l=hx^l-kr(ii,^d);
6310 csound_dim(hx^s,ii) = (csound(hxb^s)+csound(hx^s))/2d0
6311 enddo
6312 call get_rhon_tot(w,x,ixi^l,ixo^l,rho)
6313 call twofl_get_temp_n_pert_from_etot(w, x, ixi^l, ixi^l, temp)
6314 do ii=1,ndim
6315 !rho_n
6316 call hyp_coeff(ixi^l, ixoo^l, w(ixi^s,rho_n_), ii, tmp(ixi^s))
6317 nu(ixo^s) = c_hyp(rho_n_) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6318 c_shk(rho_n_) * (dxlevel(ii)**2) *divv(ixoo^s,ii)
6319 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6320
6321 !TH n
6322 call hyp_coeff(ixi^l, ixoo^l, temp(ixi^s), ii, tmp(ixi^s))
6323 nu(ixo^s) = c_hyp(e_n_) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6324 c_shk(e_n_) * (dxlevel(ii)**2) *divv(ixo^s,ii)
6325 nu(ixo^s) = nu(ixo^s) * rho(ixo^s) * rn/(twofl_gamma-1d0)
6326 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6327
6328 !visc n
6329 do jj=1,ndir
6330 call hyp_coeff(ixi^l, ixoo^l, vel(ixi^s,jj), ii, tmp(ixi^s))
6331 nu(ixo^s) = c_hyp(mom_n(jj)) * csound_dim(ixo^s,ii) * dxlevel(ii) * tmp(ixo^s) + &
6332 c_shk(mom_n(jj)) * (dxlevel(ii)**2) *divv(ixo^s,ii)
6333 nu(ixo^s) = nu(ixo^s) * rho(ixo^s)
6334 maxcoef = max(maxcoef,maxval(nu(ixo^s)))
6335 enddo
6336 enddo
6337
6338 dtnew=min(dtdiffpar*minval(dxarr(1:ndim))**2/maxcoef,dtnew)
6339 end subroutine hyperdiffusivity_get_dt
6340
6341 subroutine add_source_hyperdiffusive(qdt,ixI^L,ixO^L,w,wCT,x)
6344
6345 integer, intent(in) :: ixi^l, ixo^l
6346 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
6347 double precision, intent(inout) :: w(ixi^s,1:nw)
6348 double precision, intent(in) :: wct(ixi^s,1:nw)
6349
6350 double precision :: divv(ixi^s,1:ndim)
6351 double precision :: vel(ixi^s,1:ndir)
6352 double precision :: csound(ixi^s),csound_dim(ixi^s,1:ndim)
6353 integer :: ii,ixoo^l,hxb^l,hx^l
6354 double precision :: rho(ixi^s)
6355
6356 call twofl_get_v_c(wct,x,ixi^l,ixi^l,vel)
6357 call get_rhoc_tot(wct,x,ixi^l,ixi^l,rho)
6358 call twofl_get_csound2_c_from_conserved(wct,x,ixi^l,ixi^l,csound)
6359 csound(ixi^s) = sqrt(csound(ixi^s)) + sqrt(twofl_mag_en_all(wct,ixi^l,ixi^l) /rho(ixi^s))
6360 csound(ixi^s) = csound(ixi^s) + sqrt(sum(vel(ixi^s,1:ndir)**2 ,dim=ndim+1))
6361 do ii=1,ndim
6362 call div_vel_coeff(ixi^l, ixoo^l, vel, ii, divv(ixi^s,ii))
6363 hxmin^d=iximin^d+1;
6364 hxmax^d=iximax^d-1;
6365 hxb^l=hx^l-kr(ii,^d);
6366 csound_dim(hx^s,ii) = (csound(hxb^s)+csound(hx^s))/2d0
6367 enddo
6369 call add_viscosity_hyper_source(rho,mom_c(1), e_c_)
6370 call add_th_cond_c_hyper_source(rho)
6371 call add_ohmic_hyper_source()
6372
6373 call twofl_get_v_n(wct,x,ixi^l,ixi^l,vel)
6374 call twofl_get_csound_n(wct,x,ixi^l,ixi^l,csound)
6375 csound(ixi^s) = csound(ixi^s) + sqrt(sum(vel(ixi^s,1:ndir)**2 ,dim=ndim+1))
6376 do ii=1,ndim
6377 call div_vel_coeff(ixi^l, ixoo^l, vel, ii, divv(ixi^s,ii))
6378 hxmin^d=iximin^d+1;
6379 hxmax^d=iximax^d-1;
6380 hxb^l=hx^l-kr(ii,^d);
6381 csound_dim(hx^s,ii) = (csound(hxb^s)+csound(hx^s))/2d0
6382 enddo
6384 call get_rhon_tot(wct,x,ixi^l,ixi^l,rho)
6385 call add_viscosity_hyper_source(rho,mom_n(1), e_n_)
6386 call add_th_cond_n_hyper_source(rho)
6387
6388 contains
6389
6390 subroutine add_density_hyper_source(index_rho)
6391 integer, intent(in) :: index_rho
6392
6393 double precision :: nu(ixI^S), tmp(ixI^S)
6394
6395 do ii=1,ndim
6396 call hyp_coeff(ixi^l, ixoo^l, wct(ixi^s,index_rho), ii, tmp(ixi^s))
6397 nu(ixoo^s) = c_hyp(index_rho) * csound_dim(ixoo^s,ii) * dxlevel(ii) * tmp(ixoo^s) + &
6398 c_shk(index_rho) * (dxlevel(ii)**2) *divv(ixoo^s,ii)
6399 !print*, "IXOO HYP ", ixOO^L, " IDIMM ", ii
6400 call second_same_deriv(ixi^l, ixoo^l, nu(ixi^s), wct(ixi^s,index_rho), ii, tmp)
6401
6402 w(ixo^s,index_rho) = w(ixo^s,index_rho) + qdt * tmp(ixo^s)
6403 !print*, "RHO ", index_rho, maxval(abs(tmp(ixO^S)))
6404 enddo
6405 end subroutine add_density_hyper_source
6406
6407 subroutine add_th_cond_c_hyper_source(var2)
6408 double precision, intent(in) :: var2(ixI^S)
6409 double precision :: nu(ixI^S), tmp(ixI^S), var(ixI^S)
6410 call twofl_get_temp_c_pert_from_etot(wct, x, ixi^l, ixi^l, var)
6411 do ii=1,ndim
6412 call hyp_coeff(ixi^l, ixoo^l, var(ixi^s), ii, tmp(ixi^s))
6413 nu(ixoo^s) = c_hyp(e_c_) * csound_dim(ixoo^s,ii) * dxlevel(ii) * tmp(ixoo^s) + &
6414 c_shk(e_c_) * (dxlevel(ii)**2) *divv(ixoo^s,ii)
6415 call second_same_deriv2(ixi^l, ixoo^l, nu(ixi^s), var2(ixi^s) ,var(ixi^s), ii, tmp)
6416 w(ixo^s,e_c_) = w(ixo^s,e_c_) + qdt * tmp(ixo^s) * rc/(twofl_gamma-1d0)
6417 !print*, "TH C ", maxval(abs(tmp(ixO^S)))
6418 enddo
6419 end subroutine add_th_cond_c_hyper_source
6420
6421 subroutine add_th_cond_n_hyper_source(var2)
6422 double precision, intent(in) :: var2(ixI^S)
6423 double precision :: nu(ixI^S), tmp(ixI^S), var(ixI^S)
6424 call twofl_get_temp_n_pert_from_etot(wct, x, ixi^l, ixi^l, var)
6425 do ii=1,ndim
6426 call hyp_coeff(ixi^l, ixoo^l, var(ixi^s), ii, tmp(ixi^s))
6427 nu(ixoo^s) = c_hyp(e_n_) * csound_dim(ixoo^s,ii) * dxlevel(ii) * tmp(ixoo^s) + &
6428 c_shk(e_n_) * (dxlevel(ii)**2) *divv(ixoo^s,ii)
6429 call second_same_deriv2(ixi^l, ixoo^l, nu(ixi^s), var2(ixi^s) ,var(ixi^s), ii, tmp)
6430 w(ixo^s,e_n_) = w(ixo^s,e_n_) + qdt * tmp(ixo^s) * rn/(twofl_gamma-1d0)
6431 !print*, "TH N ", maxval(abs(tmp(ixO^S)))
6432 enddo
6433 end subroutine add_th_cond_n_hyper_source
6434
6435 subroutine add_viscosity_hyper_source(rho,index_mom1, index_e)
6436 double precision, intent(in) :: rho(ixI^S)
6437 integer, intent(in) :: index_mom1, index_e
6438
6439 double precision :: nu(ixI^S,1:ndir,1:ndim), tmp(ixI^S),tmp2(ixI^S)
6440 integer :: jj
6441
6442 do jj=1,ndir
6443 do ii=1,ndim
6444 call hyp_coeff(ixi^l, ixoo^l, vel(ixi^s,jj), ii, tmp(ixi^s))
6445 nu(ixoo^s,jj,ii) = c_hyp(index_mom1-1+jj) * csound_dim(ixoo^s,ii) * dxlevel(ii) * tmp(ixoo^s) + &
6446 c_shk(index_mom1-1+jj) * (dxlevel(ii)**2) *divv(ixoo^s,ii)
6447 enddo
6448 enddo
6449
6450 do jj=1,ndir
6451 do ii=1,ndim
6452 call second_same_deriv2(ixi^l, ixoo^l, nu(ixi^s,jj,ii), rho(ixi^s), vel(ixi^s,jj), ii, tmp)
6453 call second_same_deriv2(ixi^l, ixoo^l, nu(ixi^s,jj,ii), wct(ixi^s,index_mom1-1+jj), vel(ixi^s,jj), ii, tmp2)
6454 if(ii .eq. jj) then
6455 w(ixo^s,index_mom1-1+jj) = w(ixo^s,index_mom1-1+jj) + qdt * tmp(ixo^s)
6456 w(ixo^s,index_e) = w(ixo^s,index_e) + qdt * tmp2(ixo^s)
6457
6458 else
6459 w(ixo^s,index_mom1-1+jj) = w(ixo^s,index_mom1-1+jj) + 0.5*qdt * tmp(ixo^s)
6460 w(ixo^s,index_e) = w(ixo^s,index_e) + 0.5*qdt * tmp2(ixo^s)
6461 call second_cross_deriv2(ixi^l, ixoo^l, nu(ixi^s,ii,jj), rho(ixi^s), vel(ixi^s,ii), jj, ii, tmp)
6462 w(ixo^s,index_mom1-1+jj) = w(ixo^s,index_mom1-1+jj) + 0.5*qdt * tmp(ixo^s)
6463 call second_cross_deriv2(ixi^l, ixoo^l, nu(ixi^s,jj,ii), wct(ixi^s,index_mom1-1+jj), vel(ixi^s,jj), ii, jj, tmp2)
6464 w(ixo^s,index_e) = w(ixo^s,index_e) + 0.5*qdt * tmp2(ixo^s)
6465 endif
6466
6467 enddo
6468 enddo
6469
6470 end subroutine add_viscosity_hyper_source
6471
6472 subroutine add_ohmic_hyper_source()
6473 double precision :: nu(ixI^S,1:ndir,1:ndim), tmp(ixI^S)
6474 integer :: jj
6475
6476 do jj=1,ndir
6477 do ii=1,ndim
6478 if(ii .ne. jj) then
6479 call hyp_coeff(ixi^l, ixoo^l, wct(ixi^s,mag(jj)), ii, tmp(ixi^s))
6480 nu(ixoo^s,jj,ii) = c_hyp(mag(jj)) * csound_dim(ixoo^s,ii) * dxlevel(ii) * tmp(ixoo^s) + &
6481 c_shk(mag(jj)) * (dxlevel(ii)**2) *divv(ixoo^s,ii)
6482 endif
6483 enddo
6484 enddo
6485
6486 do jj=1,ndir
6487 do ii=1,ndim
6488 if(ii .ne. jj) then
6489 !mag field
6490 call second_same_deriv(ixi^l, ixoo^l, nu(ixi^s,jj,ii), wct(ixi^s,mag(jj)), ii, tmp)
6491 w(ixo^s,mag(jj)) = w(ixo^s,mag(jj)) + qdt * tmp(ixo^s)
6492 call second_cross_deriv(ixi^l, ixoo^l, nu(ixi^s,ii,jj), wct(ixi^s,mag(ii)), jj, ii, tmp)
6493 w(ixo^s,mag(jj)) = w(ixo^s,mag(jj)) + qdt * tmp(ixo^s)
6494 !in the total energy
6495 call second_same_deriv(ixi^l, ixoo^l, nu(ixi^s,jj,ii), wct(ixi^s,mag(jj)), ii, tmp)
6496 w(ixo^s,e_c_) = w(ixo^s,e_c_) + qdt * tmp(ixo^s)
6497 call second_cross_deriv2(ixi^l, ixoo^l, nu(ixi^s,ii,jj), wct(ixi^s,mag(jj)), wct(ixi^s,mag(ii)), jj, ii, tmp)
6498 w(ixo^s,e_c_) = w(ixo^s,e_c_) + qdt * tmp(ixo^s)
6499 endif
6500
6501 enddo
6502 enddo
6503
6504 end subroutine add_ohmic_hyper_source
6505
6506 end subroutine add_source_hyperdiffusive
6507
6508 function dump_hyperdiffusivity_coef_x(ixI^L,ixO^L, w, x, nwc) result(wnew)
6511 integer, intent(in) :: ixI^L, ixO^L, nwc
6512 double precision, intent(in) :: w(ixI^S, 1:nw)
6513 double precision, intent(in) :: x(ixI^S,1:ndim)
6514 double precision :: wnew(ixO^S, 1:nwc)
6515
6516 if(nw .ne. nwc) call mpistop("nw != nwc")
6517 wnew(ixo^s,1:nw) = dump_hyperdiffusivity_coef_dim(ixi^l,ixo^l, w, x, 1)
6518
6519 end function dump_hyperdiffusivity_coef_x
6520
6521 function dump_hyperdiffusivity_coef_y(ixI^L,ixO^L, w, x, nwc) result(wnew)
6524 integer, intent(in) :: ixI^L, ixO^L, nwc
6525 double precision, intent(in) :: w(ixI^S, 1:nw)
6526 double precision, intent(in) :: x(ixI^S,1:ndim)
6527 double precision :: wnew(ixO^S, 1:nwc)
6528
6529 if(nw .ne. nwc) call mpistop("nw != nwc")
6530 wnew(ixo^s,1:nw) = dump_hyperdiffusivity_coef_dim(ixi^l,ixo^l, w, x, 2)
6531
6532 end function dump_hyperdiffusivity_coef_y
6533
6534 function dump_hyperdiffusivity_coef_z(ixI^L,ixO^L, w, x, nwc) result(wnew)
6537 integer, intent(in) :: ixI^L, ixO^L, nwc
6538 double precision, intent(in) :: w(ixI^S, 1:nw)
6539 double precision, intent(in) :: x(ixI^S,1:ndim)
6540 double precision :: wnew(ixO^S, 1:nwc)
6541
6542 if(nw .ne. nwc) call mpistop("nw != nwc")
6543 wnew(ixo^s,1:nw) = dump_hyperdiffusivity_coef_dim(ixi^l,ixo^l, w, x, 3)
6544
6545 end function dump_hyperdiffusivity_coef_z
6546
6547 function dump_hyperdiffusivity_coef_dim(ixI^L,ixOP^L, w, x, ii) result(wnew)
6550 integer, intent(in) :: ixI^L, ixOP^L, ii
6551 double precision, intent(in) :: w(ixI^S, 1:nw)
6552 double precision, intent(in) :: x(ixI^S,1:ndim)
6553 double precision :: wnew(ixOP^S, 1:nw)
6554
6555 double precision :: nu(ixI^S),tmp(ixI^S),rho(ixI^S),temp(ixI^S)
6556 double precision :: divv(ixI^S)
6557 double precision :: vel(ixI^S,1:ndir)
6558 double precision :: csound(ixI^S),csound_dim(ixI^S)
6559 double precision :: dxarr(ndim)
6560 integer :: ixOO^L, hxb^L, hx^L, jj, ixO^L
6561
6562 ! this is done because of save_physical_boundary = true
6563 ixomin^d=max(ixopmin^d,iximin^d+3);
6564 ixomax^d=min(ixopmax^d,iximax^d-3);
6565
6566 wnew(ixop^s,1:nw) = 0d0
6567
6568 ! charges
6569 call twofl_get_temp_c_pert_from_etot(w, x, ixi^l, ixi^l, temp)
6570 call twofl_get_v_c(w,x,ixi^l,ixi^l,vel)
6571 call get_rhoc_tot(w,x,ixi^l,ixi^l,rho)
6572 call twofl_get_csound2_c_from_conserved(w,x,ixi^l,ixi^l,csound)
6573 csound(ixi^s) = sqrt(csound(ixi^s)) + sqrt(twofl_mag_en_all(w,ixi^l,ixi^l) /rho(ixi^s))
6574 csound(ixi^s) = csound(ixi^s) + sqrt(sum(vel(ixi^s,1:ndir)**2 ,dim=ndim+1))
6575 !for dim
6576 call div_vel_coeff(ixi^l, ixoo^l, vel, ii, divv(ixi^s))
6577 hxmin^d=iximin^d+1;
6578 hxmax^d=iximax^d-1;
6579 hxb^l=hx^l-kr(ii,^d);
6580 csound_dim(hx^s) = (csound(hxb^s)+csound(hx^s))/2d0
6581
6582 !TODO the following is copied
6583 !rho_c
6584 call hyp_coeff(ixi^l, ixoo^l, w(ixi^s,rho_c_), ii, tmp(ixi^s))
6585 nu(ixo^s) = c_hyp(rho_c_) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6586 c_shk(rho_c_) * (dxlevel(ii)**2) *divv(ixo^s)
6587
6588 wnew(ixo^s,rho_c_) = nu(ixo^s)
6589
6590 !TH c
6591 call hyp_coeff(ixi^l, ixoo^l, temp(ixi^s), ii, tmp(ixi^s))
6592 nu(ixo^s) = c_hyp(e_c_) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6593 c_shk(e_c_) * (dxlevel(ii)**2) *divv(ixo^s)
6594 nu(ixo^s) = nu(ixo^s) * rho(ixo^s) * rc/(twofl_gamma-1d0)
6595 wnew(ixo^s,e_c_) = nu(ixo^s)
6596
6597 !visc c
6598 do jj=1,ndir
6599 call hyp_coeff(ixi^l, ixoo^l, vel(ixi^s,jj), ii, tmp(ixi^s))
6600 nu(ixo^s) = c_hyp(mom_c(jj)) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6601 c_shk(mom_c(jj)) * (dxlevel(ii)**2) *divv(ixo^s)
6602 nu(ixo^s) = nu(ixo^s) * rho(ixo^s)
6603 wnew(ixo^s,mom_c(jj)) = nu(ixo^s)
6604 enddo
6605
6606 ! Ohmic
6607 do jj=1,ndir
6608 if(ii .ne. jj) then
6609 call hyp_coeff(ixi^l, ixoo^l, w(ixi^s,mag(jj)), ii, tmp(ixi^s))
6610 nu(ixo^s) = c_hyp(mag(jj)) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6611 c_shk(mag(jj)) * (dxlevel(ii)**2) *divv(ixo^s)
6612 wnew(ixo^s,mag(jj)) = nu(ixo^s)
6613 endif
6614 enddo
6615
6616 !end for dim
6617
6618 ! neutrals
6619 call get_rhon_tot(w,x,ixi^l,ixo^l,rho)
6620 call twofl_get_temp_n_pert_from_etot(w, x, ixi^l, ixi^l, temp)
6621 call twofl_get_v_n(w,x,ixi^l,ixi^l,vel)
6622 call twofl_get_csound_n(w,x,ixi^l,ixi^l,csound)
6623 csound(ixi^s) = csound(ixi^s) + sqrt(sum(vel(ixi^s,1:ndir)**2 ,dim=ndim+1))
6624 !for dim
6625 call div_vel_coeff(ixi^l, ixoo^l, vel, ii, divv(ixi^s))
6626 hxb^l=ixoo^l-kr(ii,^d);
6627 csound_dim(ixoo^s) = (csound(hxb^s)+csound(ixoo^s))/2d0
6628 !rho_n
6629 call hyp_coeff(ixi^l, ixoo^l, w(ixi^s,rho_n_), ii, tmp(ixi^s))
6630 nu(ixo^s) = c_hyp(rho_n_) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6631 c_shk(rho_n_) * (dxlevel(ii)**2) *divv(ixoo^s)
6632 wnew(ixo^s,rho_n_) = nu(ixo^s)
6633
6634 !TH n
6635 call hyp_coeff(ixi^l, ixoo^l, temp(ixi^s), ii, tmp(ixi^s))
6636 nu(ixo^s) = c_hyp(e_n_) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6637 c_shk(e_n_) * (dxlevel(ii)**2) *divv(ixo^s)
6638 nu(ixo^s) = nu(ixo^s) * rho(ixo^s) * rn/(twofl_gamma-1d0)
6639 wnew(ixo^s,e_n_) = nu(ixo^s)
6640
6641 !visc n
6642 do jj=1,ndir
6643 call hyp_coeff(ixi^l, ixoo^l, vel(ixi^s,jj), ii, tmp(ixi^s))
6644 nu(ixo^s) = c_hyp(mom_n(jj)) * csound_dim(ixo^s) * dxlevel(ii) * tmp(ixo^s) + &
6645 c_shk(mom_n(jj)) * (dxlevel(ii)**2) *divv(ixo^s)
6646 nu(ixo^s) = nu(ixo^s) * rho(ixo^s)
6647 wnew(ixo^s,mom_n(jj)) = nu(ixo^s)
6648 enddo
6649 !end for dim
6650
6651 end function dump_hyperdiffusivity_coef_dim
6652
6653 function dump_coll_terms(ixI^L,ixO^L, w, x, nwc) result(wnew)
6655 integer, intent(in) :: ixI^L,ixO^L, nwc
6656 double precision, intent(in) :: w(ixI^S, 1:nw)
6657 double precision, intent(in) :: x(ixI^S,1:ndim)
6658 double precision :: wnew(ixO^S, 1:nwc)
6659 double precision :: tmp(ixI^S),tmp2(ixI^S)
6660
6661 call get_alpha_coll(ixi^l, ixo^l, w, x, tmp(ixi^s))
6662 wnew(ixo^s,1)= tmp(ixo^s)
6663 call get_gamma_ion_rec(ixi^l, ixo^l, w, x, tmp(ixi^s), tmp2(ixi^s))
6664 wnew(ixo^s,2)= tmp(ixo^s)
6665 wnew(ixo^s,3)= tmp2(ixo^s)
6666
6667 end function dump_coll_terms
6668
6669 subroutine get_gamma_ion_rec(ixI^L, ixO^L, w, x, gamma_rec, gamma_ion)
6671
6672 integer, intent(in) :: ixi^l, ixo^l
6673 double precision, intent(in) :: w(ixi^s,1:nw)
6674 double precision, intent(in) :: x(ixi^s,1:ndim)
6675 double precision, intent(out) :: gamma_rec(ixi^s),gamma_ion(ixi^s)
6676 ! calculations are done in S.I. units
6677 double precision, parameter :: a = 2.91e-14, & !m3/s
6678 k = 0.39, &
6679 xx = 0.232, &
6680 eion = 13.6 ! eV
6681 double precision, parameter :: echarge=1.6022d-19 !C
6682 double precision :: rho(ixi^s), tmp(ixi^s)
6683
6684 call twofl_get_pthermal_c(w,x,ixi^l,ixo^l,tmp)
6685 call get_rhoc_tot(w,x,ixi^l,ixo^l,rho)
6686 tmp(ixo^s) = tmp(ixo^s)/(rc * rho(ixo^s))
6687
6688 !transform to SI units
6689 tmp(ixo^s) = tmp(ixo^s) * unit_temperature * kb_si/echarge !* BK/ECHARGE means K to eV
6690 !number electrons rho_c = n_e * MH, in normalized units MH=1 and n = rho
6691 rho(ixo^s) = rho(ixo^s) * unit_numberdensity
6692 if(.not. si_unit) then
6693 !1/cm^3 = 1e6/m^3
6694 rho(ixo^s) = rho(ixo^s) * 1d6
6695 endif
6696 gamma_rec(ixo^s) = rho(ixo^s) /sqrt(tmp(ixo^s)) * 2.6e-19
6697 gamma_ion(ixo^s) = ((rho(ixo^s) * a) /(xx + eion/tmp(ixo^s))) * ((eion/tmp(ixo^s))**k) * exp(-eion/tmp(ixo^s))
6698 ! see Voronov table: valid for temp min = 1eV(approx 11605 K), Temp max = 20KeV
6699 !to normalized
6700 gamma_rec(ixo^s) = gamma_rec(ixo^s) * unit_time
6701 gamma_ion(ixo^s) = gamma_ion(ixo^s) * unit_time
6702
6703 if (associated(usr_mask_gamma_ion_rec)) then
6704 call usr_mask_gamma_ion_rec(ixi^l,ixo^l,w,x,gamma_ion, gamma_rec)
6705 end if
6706 end subroutine get_gamma_ion_rec
6707
6708 subroutine get_alpha_coll(ixI^L, ixO^L, w, x, alpha)
6710 integer, intent(in) :: ixi^l, ixo^l
6711 double precision, intent(in) :: w(ixi^s,1:nw)
6712 double precision, intent(in) :: x(ixi^s,1:ndim)
6713 double precision, intent(out) :: alpha(ixi^s)
6715 alpha(ixo^s) = twofl_alpha_coll
6716 else
6717 call get_alpha_coll_plasma(ixi^l, ixo^l, w, x, alpha)
6718 endif
6719 if (associated(usr_mask_alpha)) then
6720 call usr_mask_alpha(ixi^l,ixo^l,w,x,alpha)
6721 end if
6722 end subroutine get_alpha_coll
6723
6724 subroutine get_alpha_coll_plasma(ixI^L, ixO^L, w, x, alpha)
6726 integer, intent(in) :: ixi^l, ixo^l
6727 double precision, intent(in) :: w(ixi^s,1:nw)
6728 double precision, intent(in) :: x(ixi^s,1:ndim)
6729 double precision, intent(out) :: alpha(ixi^s)
6730 double precision :: pe(ixi^s),rho(ixi^s), tmp(ixi^s), tmp2(ixi^s)
6731
6732 double precision :: sigma_in = 1e-19 ! m^2
6733 ! make calculation in SI physical units
6734
6735 call twofl_get_pthermal_c(w,x,ixi^l,ixo^l,pe)
6736 call get_rhoc_tot(w,x,ixi^l,ixo^l,rho)
6737 tmp(ixo^s) = pe(ixo^s)/(rc * rho(ixo^s))
6738 call twofl_get_pthermal_n(w,x,ixi^l,ixo^l,pe)
6739 call get_rhon_tot(w,x,ixi^l,ixo^l,rho)
6740 tmp2(ixo^s) = pe(ixo^s)/(rn * rho(ixo^s))
6741 alpha(ixo^s) = (2d0/(mp_si**(3d0/2) * sqrt(dpi))*sqrt(0.5*(tmp(ixo^s)+tmp2(ixo^s))*unit_temperature*kb_si) * sigma_in)*unit_time * unit_density
6742 if(.not. si_unit) then
6743 alpha(ixo^s) = alpha(ixo^s) * 1d3 ! this comes from unit_density: g/cm^3 = 1e-3 kg/m^3
6744 endif
6745
6746 end subroutine get_alpha_coll_plasma
6747
6748 subroutine calc_mult_factor1(ixI^L, ixO^L, step_dt, JJ, res)
6749 integer, intent(in) :: ixi^l, ixo^l
6750 double precision, intent(in) :: step_dt
6751 double precision, intent(in) :: jj(ixi^s)
6752 double precision, intent(out) :: res(ixi^s)
6753
6754 res(ixo^s) = step_dt/(1d0 + step_dt * jj(ixo^s))
6755
6756 end subroutine calc_mult_factor1
6757
6758 subroutine calc_mult_factor2(ixI^L, ixO^L, step_dt, JJ, res)
6759 integer, intent(in) :: ixi^l, ixo^l
6760 double precision, intent(in) :: step_dt
6761 double precision, intent(in) :: jj(ixi^s)
6762 double precision, intent(out) :: res(ixi^s)
6763
6764 res(ixo^s) = (1d0 - exp(-step_dt * jj(ixo^s)))/jj(ixo^s)
6765
6766 end subroutine calc_mult_factor2
6767
6768 subroutine advance_implicit_grid(ixI^L, ixO^L, w, wout, x, dtfactor,qdt)
6770 integer, intent(in) :: ixi^l, ixo^l
6771 double precision, intent(in) :: qdt
6772 double precision, intent(in) :: dtfactor
6773 double precision, intent(in) :: w(ixi^s,1:nw)
6774 double precision, intent(in) :: x(ixi^s,1:ndim)
6775 double precision, intent(out) :: wout(ixi^s,1:nw)
6776
6777 integer :: idir
6778 double precision :: tmp(ixi^s),tmp1(ixi^s),tmp2(ixi^s),tmp3(ixi^s),tmp4(ixi^s),tmp5(ixi^s)
6779 double precision :: v_c(ixi^s,ndir), v_n(ixi^s,ndir)
6780 double precision :: rhon(ixi^s), rhoc(ixi^s), alpha(ixi^s)
6781 double precision, allocatable :: gamma_rec(:^d&), gamma_ion(:^D&)
6782
6783 !TODO latest changes sets already wout to w in implicit update (see where psb=psa)
6784 ! commment out setting mag and density when they are not modified here
6785
6786 ! copy vars at the indices which are not updated here: mag. field
6787 wout(ixo^s,mag(:)) = w(ixo^s,mag(:))
6788
6789 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
6790 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
6791 !update density
6792 if(twofl_coll_inc_ionrec) then
6793 allocate(gamma_ion(ixi^s), gamma_rec(ixi^s))
6794 call get_gamma_ion_rec(ixi^l, ixo^l, w, x, gamma_rec, gamma_ion)
6795 tmp2(ixo^s) = gamma_rec(ixo^s) + gamma_ion(ixo^s)
6796 call calc_mult_factor(ixi^l, ixo^l, dtfactor * qdt, tmp2, tmp3)
6797 tmp(ixo^s) = (-gamma_ion(ixo^s) * rhon(ixo^s) + &
6798 gamma_rec(ixo^s) * rhoc(ixo^s))
6799 wout(ixo^s,rho_n_) = w(ixo^s,rho_n_) + tmp(ixo^s) * tmp3(ixo^s)
6800 wout(ixo^s,rho_c_) = w(ixo^s,rho_c_) - tmp(ixo^s) * tmp3(ixo^s)
6801 else
6802 wout(ixo^s,rho_n_) = w(ixo^s,rho_n_)
6803 wout(ixo^s,rho_c_) = w(ixo^s,rho_c_)
6804 endif
6805
6806 call get_alpha_coll(ixi^l, ixo^l, w, x, alpha)
6807
6808 !-J11 + J12 for momentum and kinetic energy
6809 tmp2(ixo^s) = alpha(ixo^s) * (rhon(ixo^s) + rhoc(ixo^s))
6810 if(twofl_coll_inc_ionrec) then
6811 tmp2(ixo^s) = tmp2(ixo^s) + gamma_ion(ixo^s) + gamma_rec(ixo^s)
6812 endif
6813 call calc_mult_factor(ixi^l, ixo^l, dtfactor * qdt, tmp2, tmp3)
6814
6815 ! momentum update
6816 do idir=1,ndir
6817
6818 tmp(ixo^s) = alpha(ixo^s)* (-rhoc(ixo^s) * w(ixo^s,mom_n(idir)) + rhon(ixo^s) * w(ixo^s,mom_c(idir)))
6819 if(twofl_coll_inc_ionrec) then
6820 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s) * w(ixo^s,mom_n(idir)) + gamma_rec(ixo^s) * w(ixo^s,mom_c(idir))
6821 endif
6822
6823 wout(ixo^s,mom_n(idir)) = w(ixo^s,mom_n(idir)) + tmp(ixo^s) * tmp3(ixo^s)
6824 wout(ixo^s,mom_c(idir)) = w(ixo^s,mom_c(idir)) - tmp(ixo^s) * tmp3(ixo^s)
6825 enddo
6826
6827 ! energy update
6828
6829 ! kinetic energy update
6830 if(.not. phys_internal_e) then
6831 ! E_tot includes kinetic energy
6832 tmp1(ixo^s) = twofl_kin_en_n(w,ixi^l,ixo^l)
6833 tmp2(ixo^s) = twofl_kin_en_c(w,ixi^l,ixo^l)
6834 tmp4(ixo^s) = w(ixo^s,e_n_) - tmp1(ixo^s) !E_tot - E_kin
6835 tmp5(ixo^s) = w(ixo^s,e_c_) - tmp2(ixo^s)
6836 if(phys_total_energy) then
6837 tmp5(ixo^s) = tmp5(ixo^s) - twofl_mag_en(w,ixi^l,ixo^l)
6838 endif
6839
6840 !!implicit update
6841 tmp(ixo^s) = alpha(ixo^s)*(-rhoc(ixo^s) * tmp1(ixo^s) + rhon(ixo^s) * tmp2(ixo^s))
6842 if(twofl_coll_inc_ionrec) then
6843 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s) * tmp1(ixo^s) + gamma_rec(ixo^s) * tmp2(ixo^s)
6844 endif
6845
6846 wout(ixo^s,e_n_) = w(ixo^s,e_n_) + tmp(ixo^s) * tmp3(ixo^s)
6847 wout(ixo^s,e_c_) = w(ixo^s,e_c_) - tmp(ixo^s) * tmp3(ixo^s)
6848
6849 else
6850 tmp4(ixo^s) = w(ixo^s,e_n_)
6851 tmp5(ixo^s) = w(ixo^s,e_c_)
6852 ! calculate velocities, using the already updated variables
6853 call twofl_get_v_n(wout,x,ixi^l,ixo^l,v_n)
6854 call twofl_get_v_c(wout,x,ixi^l,ixo^l,v_c)
6855 tmp1(ixo^s) = alpha(ixo^s) * rhoc(ixo^s) * rhon(ixo^s)
6856 tmp2(ixo^s) = tmp1(ixo^s)
6857 if(twofl_coll_inc_ionrec) then
6858 tmp1(ixo^s) = tmp1(ixo^s) + rhoc(ixo^s) * gamma_rec(ixo^s)
6859 tmp2(ixo^s) = tmp2(ixo^s) + rhon(ixo^s) * gamma_ion(ixo^s)
6860 endif
6861
6862 tmp(ixo^s) = 0.5d0 * sum((v_c(ixo^s,1:ndir) - v_n(ixo^s,1:ndir))**2, dim=ndim+1) &
6863 * dtfactor * qdt
6864 wout(ixo^s,e_n_) = w(ixo^s,e_n_) + tmp(ixo^s)*tmp1(ixo^s)
6865 wout(ixo^s,e_c_) = w(ixo^s,e_c_) + tmp(ixo^s)*tmp2(ixo^s)
6866 endif
6867
6868 !update internal energy
6869 if(twofl_coll_inc_te) then
6870 if(has_equi_pe_n0) then
6871 tmp2(ixo^s)= block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
6872 endif
6873 if(has_equi_pe_c0) then
6874 tmp3(ixo^s)=block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
6875 endif
6876 if (twofl_equi_thermal) then
6877 tmp(ixo^s) = alpha(ixo^s) *(-1d0/rn*(rhoc(ixo^s) * tmp4(ixo^s) + &
6878 tmp2(ixo^s)*w(ixo^s,rho_c_)) + 1d0/rc*(rhon(ixo^s) * tmp5(ixo^s) +&
6879 tmp3(ixo^s)*w(ixo^s,rho_n_)))
6880 endif
6881 if(has_equi_pe_n0) then
6882 tmp4(ixo^s) = tmp2(ixo^s) + tmp4(ixo^s)
6883 endif
6884 if(has_equi_pe_c0) then
6885 tmp5(ixo^s) = tmp3(ixo^s) + tmp5(ixo^s)
6886 endif
6887 if (.not. twofl_equi_thermal) then
6888 tmp(ixo^s) = alpha(ixo^s) *(-rhoc(ixo^s)/rn * tmp4(ixo^s) + rhon(ixo^s)/rc * tmp5(ixo^s))
6889 endif
6890 tmp2(ixo^s) = alpha(ixo^s) * (rhon(ixo^s)/rc + rhoc(ixo^s)/rn)
6891 if(twofl_coll_inc_ionrec) then
6892 tmp2(ixo^s) = tmp2(ixo^s) + gamma_rec(ixo^s)/rc + gamma_ion(ixo^s)/rn
6893 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s)/rn * tmp4(ixo^s) + gamma_rec(ixo^s)/rc * tmp5(ixo^s)
6894 endif
6895 call calc_mult_factor(ixi^l, ixo^l, dtfactor * qdt, tmp2, tmp3)
6896 wout(ixo^s,e_n_) = wout(ixo^s,e_n_)+tmp(ixo^s)*tmp3(ixo^s)
6897 wout(ixo^s,e_c_) = wout(ixo^s,e_c_)-tmp(ixo^s)*tmp3(ixo^s)
6898 endif
6899 if(twofl_coll_inc_ionrec) then
6900 deallocate(gamma_ion, gamma_rec)
6901 endif
6902 end subroutine advance_implicit_grid
6903
6904 !> Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
6905 subroutine twofl_implicit_coll_terms_update(dtfactor,qdt,qtC,psb,psa)
6908
6909 type(state), target :: psa(max_blocks)
6910 type(state), target :: psb(max_blocks)
6911 double precision, intent(in) :: qdt
6912 double precision, intent(in) :: qtc
6913 double precision, intent(in) :: dtfactor
6914
6915 integer :: iigrid, igrid
6916 !print*, "IMPL call ", it
6917
6918 call getbc(global_time,0.d0,psa,1,nw)
6919 !$OMP PARALLEL DO PRIVATE(igrid)
6920 do iigrid=1,igridstail; igrid=igrids(iigrid);
6921 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
6922 block=>psa(igrid)
6923 call advance_implicit_grid(ixg^ll, ixg^ll, psa(igrid)%w, psb(igrid)%w, psa(igrid)%x, dtfactor,qdt)
6924 end do
6925 !$OMP END PARALLEL DO
6926
6927 end subroutine twofl_implicit_coll_terms_update
6928
6929 !> inplace update of psa==>F_im(psa)
6930 subroutine twofl_evaluate_implicit(qtC,psa)
6932 type(state), target :: psa(max_blocks)
6933 double precision, intent(in) :: qtc
6934
6935 integer :: iigrid, igrid, level
6936
6937 !$OMP PARALLEL DO PRIVATE(igrid)
6938 do iigrid=1,igridstail; igrid=igrids(iigrid);
6939 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
6940 block=>psa(igrid)
6941 call coll_terms(ixg^ll,ixm^ll,psa(igrid)%w,psa(igrid)%x)
6942 end do
6943 !$OMP END PARALLEL DO
6944
6945 end subroutine twofl_evaluate_implicit
6946
6947 subroutine coll_terms(ixI^L,ixO^L,w,x)
6949 integer, intent(in) :: ixi^l, ixo^l
6950 double precision, intent(inout) :: w(ixi^s, 1:nw)
6951 double precision, intent(in) :: x(ixi^s,1:ndim)
6952
6953 integer :: idir
6954 double precision :: tmp(ixi^s),tmp1(ixi^s),tmp2(ixi^s),tmp3(ixi^s),tmp4(ixi^s),tmp5(ixi^s)
6955 !double precision :: v_c(ixI^S,ndir), v_n(ixI^S,ndir)
6956 double precision, allocatable :: v_c(:^d&,:), v_n(:^D&,:)
6957 double precision, allocatable :: rho_c1(:^d&), rho_n1(:^D&)
6958 double precision :: rhon(ixi^s), rhoc(ixi^s), alpha(ixi^s)
6959 double precision, allocatable :: gamma_rec(:^d&), gamma_ion(:^D&)
6960
6961 ! copy density before overwrite
6962 if(twofl_equi_thermal) then
6963 allocate(rho_n1(ixi^s), rho_c1(ixi^s))
6964 rho_n1(ixo^s) = w(ixo^s,rho_n_)
6965 rho_c1(ixo^s) = w(ixo^s,rho_c_)
6966 endif
6967
6968 ! get total density before overwrite density
6969 call get_rhon_tot(w,x,ixi^l,ixo^l,rhon)
6970 call get_rhoc_tot(w,x,ixi^l,ixo^l,rhoc)
6971 if(phys_internal_e) then
6972 ! get velocity before overwrite momentum
6973 allocate(v_n(ixi^s,ndir), v_c(ixi^s,ndir))
6974 call twofl_get_v_n(w,x,ixi^l,ixo^l,v_n)
6975 call twofl_get_v_c(w,x,ixi^l,ixo^l,v_c)
6976 else
6977 ! get ke before overwrite density and momentum
6978 tmp1(ixo^s) = twofl_kin_en_n(w,ixi^l,ixo^l)
6979 tmp2(ixo^s) = twofl_kin_en_c(w,ixi^l,ixo^l)
6980 endif
6981
6982 !update density
6983 if(twofl_coll_inc_ionrec) then
6984 allocate(gamma_ion(ixi^s), gamma_rec(ixi^s))
6985 call get_gamma_ion_rec(ixi^l, ixo^l, w, x, gamma_rec, gamma_ion)
6986 tmp(ixo^s) = -gamma_ion(ixo^s) * rhon(ixo^s) + &
6987 gamma_rec(ixo^s) * rhoc(ixo^s)
6988 w(ixo^s,rho_n_) = tmp(ixo^s)
6989 w(ixo^s,rho_c_) = -tmp(ixo^s)
6990 else
6991 w(ixo^s,rho_n_) = 0d0
6992 w(ixo^s,rho_c_) = 0d0
6993
6994 endif
6995
6996 call get_alpha_coll(ixi^l, ixo^l, w, x, alpha)
6997
6998 ! momentum update
6999 do idir=1,ndir
7000
7001 tmp(ixo^s) = alpha(ixo^s)* (-rhoc(ixo^s) * w(ixo^s,mom_n(idir)) + rhon(ixo^s) * w(ixo^s,mom_c(idir)))
7002 if(twofl_coll_inc_ionrec) then
7003 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s) * w(ixo^s,mom_n(idir)) + gamma_rec(ixo^s) * w(ixo^s,mom_c(idir))
7004 endif
7005
7006 w(ixo^s,mom_n(idir)) = tmp(ixo^s)
7007 w(ixo^s,mom_c(idir)) = -tmp(ixo^s)
7008 enddo
7009
7010 ! energy update
7011
7012 ! kinetic energy update
7013 if(.not. phys_internal_e) then
7014 ! E_tot includes kinetic energy
7015 tmp4(ixo^s) = w(ixo^s,e_n_) - tmp1(ixo^s) !E_tot - E_kin
7016 tmp5(ixo^s) = w(ixo^s,e_c_) - tmp2(ixo^s)
7017 if(phys_total_energy) then
7018 tmp5(ixo^s) = tmp5(ixo^s) - twofl_mag_en(w,ixi^l,ixo^l)
7019 endif
7020 ! tmp4 = eint_n, tmp5 = eint_c
7021 ! tmp1 = ke_n, tmp2 = ke_c
7022 tmp(ixo^s) = alpha(ixo^s)*(-rhoc(ixo^s) * tmp1(ixo^s) + rhon(ixo^s) * tmp2(ixo^s))
7023 if(twofl_coll_inc_ionrec) then
7024 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s) * tmp1(ixo^s) + gamma_rec(ixo^s) * tmp2(ixo^s)
7025 endif
7026
7027 w(ixo^s,e_n_) = tmp(ixo^s)
7028 w(ixo^s,e_c_) = -tmp(ixo^s)
7029
7030 else
7031 tmp4(ixo^s) = w(ixo^s,e_n_)
7032 tmp5(ixo^s) = w(ixo^s,e_c_)
7033 tmp1(ixo^s) = alpha(ixo^s) * rhoc(ixo^s) * rhon(ixo^s)
7034 tmp2(ixo^s) = tmp1(ixo^s)
7035 if(twofl_coll_inc_ionrec) then
7036 tmp1(ixo^s) = tmp1(ixo^s) + rhoc(ixo^s) * gamma_rec(ixo^s)
7037 tmp2(ixo^s) = tmp2(ixo^s) + rhon(ixo^s) * gamma_ion(ixo^s)
7038 endif
7039
7040 tmp(ixo^s) = 0.5d0 * sum((v_c(ixo^s,1:ndir) - v_n(ixo^s,1:ndir))**2, dim=ndim+1)
7041 w(ixo^s,e_n_) = tmp(ixo^s)*tmp1(ixo^s)
7042 w(ixo^s,e_c_) = tmp(ixo^s)*tmp2(ixo^s)
7043 endif
7044
7045 !update internal energy
7046 if(twofl_coll_inc_te) then
7047
7048 if(has_equi_pe_n0) then
7049 tmp2(ixo^s)= block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
7050 endif
7051 if(has_equi_pe_c0) then
7052 tmp3(ixo^s)=block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
7053 endif
7054 if (twofl_equi_thermal) then
7055 tmp(ixo^s) = alpha(ixo^s) *(-1d0/rn*(rhoc(ixo^s) * tmp4(ixo^s) + &
7056 tmp2(ixo^s)*rho_c1(ixo^s)) + 1d0/rc*(rhon(ixo^s) * tmp5(ixo^s) +&
7057 tmp3(ixo^s)*rho_n1(ixo^s)))
7058 endif
7059 if(has_equi_pe_n0) then
7060 tmp4(ixo^s) = tmp2(ixo^s) + tmp4(ixo^s)
7061 endif
7062 if(has_equi_pe_c0) then
7063 tmp5(ixo^s) = tmp3(ixo^s) + tmp5(ixo^s)
7064 endif
7065 if (.not. twofl_equi_thermal) then
7066 tmp(ixo^s) = alpha(ixo^s) *(-rhoc(ixo^s)/rn * tmp4(ixo^s) + rhon(ixo^s)/rc * tmp5(ixo^s))
7067 endif
7068
7069 if(twofl_coll_inc_ionrec) then
7070 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s)/rn * tmp4(ixo^s) + gamma_rec(ixo^s)/rc * tmp5(ixo^s)
7071 endif
7072
7073 w(ixo^s,e_n_) = w(ixo^s,e_n_)+tmp(ixo^s)
7074 w(ixo^s,e_c_) = w(ixo^s,e_c_)-tmp(ixo^s)
7075 endif
7076 if(twofl_coll_inc_ionrec) then
7077 deallocate(gamma_ion, gamma_rec)
7078 endif
7079 if(phys_internal_e) then
7080 deallocate(v_n, v_c)
7081 endif
7082 if(twofl_equi_thermal) then
7083 deallocate(rho_n1, rho_c1)
7084 endif
7085 !set contribution to mag field
7086 w(ixo^s,mag(1:ndir)) = 0d0
7087
7088 end subroutine coll_terms
7089
7090 subroutine twofl_explicit_coll_terms_update(qdt,ixI^L,ixO^L,w,wCT,x)
7092
7093 integer, intent(in) :: ixi^l, ixo^l
7094 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
7095 double precision, intent(inout) :: w(ixi^s,1:nw)
7096 double precision, intent(in) :: wct(ixi^s,1:nw)
7097
7098 integer :: idir
7099 double precision :: tmp(ixi^s),tmp1(ixi^s),tmp2(ixi^s),tmp3(ixi^s),tmp4(ixi^s),tmp5(ixi^s)
7100 double precision :: v_c(ixi^s,ndir), v_n(ixi^s,ndir)
7101 double precision :: rhon(ixi^s), rhoc(ixi^s), alpha(ixi^s)
7102 double precision, allocatable :: gamma_rec(:^d&), gamma_ion(:^D&)
7103
7104 call get_rhon_tot(wct,x,ixi^l,ixo^l,rhon)
7105 call get_rhoc_tot(wct,x,ixi^l,ixo^l,rhoc)
7106 !update density
7107 if(twofl_coll_inc_ionrec) then
7108 allocate(gamma_ion(ixi^s), gamma_rec(ixi^s))
7109 call get_gamma_ion_rec(ixi^l, ixo^l, wct, x, gamma_rec, gamma_ion)
7110 tmp(ixo^s) = qdt *(-gamma_ion(ixo^s) * rhon(ixo^s) + &
7111 gamma_rec(ixo^s) * rhoc(ixo^s))
7112 w(ixo^s,rho_n_) = w(ixo^s,rho_n_) + tmp(ixo^s)
7113 w(ixo^s,rho_c_) = w(ixo^s,rho_c_) - tmp(ixo^s)
7114 endif
7115
7116 call get_alpha_coll(ixi^l, ixo^l, wct, x, alpha)
7117
7118 ! momentum update
7119 do idir=1,ndir
7120
7121 tmp(ixo^s) = alpha(ixo^s)* (-rhoc(ixo^s) * wct(ixo^s,mom_n(idir)) + rhon(ixo^s) * wct(ixo^s,mom_c(idir)))
7122 if(twofl_coll_inc_ionrec) then
7123 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s) * wct(ixo^s,mom_n(idir)) + gamma_rec(ixo^s) * wct(ixo^s,mom_c(idir))
7124 endif
7125 tmp(ixo^s) =tmp(ixo^s) * qdt
7126
7127 w(ixo^s,mom_n(idir)) = w(ixo^s,mom_n(idir)) + tmp(ixo^s)
7128 w(ixo^s,mom_c(idir)) = w(ixo^s,mom_c(idir)) - tmp(ixo^s)
7129 enddo
7130
7131 ! energy update
7132
7133 ! kinetic energy update
7134 if(.not. phys_internal_e) then
7135 ! E_tot includes kinetic energy
7136 tmp1(ixo^s) = twofl_kin_en_n(wct,ixi^l,ixo^l)
7137 tmp2(ixo^s) = twofl_kin_en_c(wct,ixi^l,ixo^l)
7138 tmp4(ixo^s) = wct(ixo^s,e_n_) - tmp1(ixo^s) !E_tot - E_kin
7139 tmp5(ixo^s) = wct(ixo^s,e_c_) - tmp2(ixo^s)
7140 if(phys_total_energy) then
7141 tmp5(ixo^s) = tmp5(ixo^s) - twofl_mag_en(wct,ixi^l,ixo^l)
7142 endif
7143
7144 tmp(ixo^s) = alpha(ixo^s)*(-rhoc(ixo^s) * tmp1(ixo^s) + rhon(ixo^s) * tmp2(ixo^s))
7145 if(twofl_coll_inc_ionrec) then
7146 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s) * tmp1(ixo^s) + gamma_rec(ixo^s) * tmp2(ixo^s)
7147 endif
7148 tmp(ixo^s) =tmp(ixo^s) * qdt
7149
7150 w(ixo^s,e_n_) = w(ixo^s,e_n_) + tmp(ixo^s)
7151 w(ixo^s,e_c_) = w(ixo^s,e_c_) - tmp(ixo^s)
7152
7153 else
7154 tmp4(ixo^s) = w(ixo^s,e_n_)
7155 tmp5(ixo^s) = w(ixo^s,e_c_)
7156 call twofl_get_v_n(wct,x,ixi^l,ixo^l,v_n)
7157 call twofl_get_v_c(wct,x,ixi^l,ixo^l,v_c)
7158 tmp1(ixo^s) = alpha(ixo^s) * rhoc(ixo^s) * rhon(ixo^s)
7159 tmp2(ixo^s) = tmp1(ixo^s)
7160 if(twofl_coll_inc_ionrec) then
7161 tmp1(ixo^s) = tmp1(ixo^s) + rhoc(ixo^s) * gamma_rec(ixo^s)
7162 tmp2(ixo^s) = tmp2(ixo^s) + rhon(ixo^s) * gamma_ion(ixo^s)
7163 endif
7164
7165 tmp(ixo^s) = 0.5d0 * sum((v_c(ixo^s,1:ndir) - v_n(ixo^s,1:ndir))**2, dim=ndim+1) * qdt
7166 w(ixo^s,e_n_) = w(ixo^s,e_n_) + tmp(ixo^s)*tmp1(ixo^s)
7167 w(ixo^s,e_c_) = w(ixo^s,e_c_) + tmp(ixo^s)*tmp2(ixo^s)
7168 endif
7169
7170 !update internal energy
7171 if(twofl_coll_inc_te) then
7172 if(has_equi_pe_n0) then
7173 tmp2(ixo^s)= block%equi_vars(ixo^s,equi_pe_n0_,0)*inv_gamma_1
7174 endif
7175 if(has_equi_pe_c0) then
7176 tmp3(ixo^s)=block%equi_vars(ixo^s,equi_pe_c0_,0)*inv_gamma_1
7177 endif
7178 if (twofl_equi_thermal) then
7179 tmp(ixo^s) = alpha(ixo^s) *(-1d0/rn*(rhoc(ixo^s) * tmp4(ixo^s) + &
7180 tmp2(ixo^s)*wct(ixo^s,rho_c_)) + 1d0/rc*(rhon(ixo^s) * tmp5(ixo^s) +&
7181 tmp3(ixo^s)*wct(ixo^s,rho_n_)))
7182 endif
7183 if(has_equi_pe_n0) then
7184 tmp4(ixo^s) = tmp2(ixo^s) + tmp4(ixo^s)
7185 endif
7186 if(has_equi_pe_c0) then
7187 tmp5(ixo^s) = tmp3(ixo^s) + tmp5(ixo^s)
7188 endif
7189 if (.not. twofl_equi_thermal) then
7190 tmp(ixo^s) = alpha(ixo^s) *(-rhoc(ixo^s)/rn * tmp4(ixo^s) + rhon(ixo^s)/rc * tmp5(ixo^s))
7191 endif
7192
7193 if(twofl_coll_inc_ionrec) then
7194 tmp(ixo^s) = tmp(ixo^s) - gamma_ion(ixo^s)/rn * tmp4(ixo^s) + gamma_rec(ixo^s)/rc * tmp5(ixo^s)
7195 endif
7196
7197 tmp(ixo^s) =tmp(ixo^s) * qdt
7198
7199 w(ixo^s,e_n_) = w(ixo^s,e_n_)+tmp(ixo^s)
7200 w(ixo^s,e_c_) = w(ixo^s,e_c_)-tmp(ixo^s)
7201 endif
7202 if(twofl_coll_inc_ionrec) then
7203 deallocate(gamma_ion, gamma_rec)
7204 endif
7205 end subroutine twofl_explicit_coll_terms_update
7206
7207 subroutine rfactor_c(w,x,ixI^L,ixO^L,Rfactor)
7209 integer, intent(in) :: ixi^l, ixo^l
7210 double precision, intent(in) :: w(ixi^s,1:nw)
7211 double precision, intent(in) :: x(ixi^s,1:ndim)
7212 double precision, intent(out):: rfactor(ixi^s)
7213
7214 rfactor(ixo^s)=rc
7215
7216 end subroutine rfactor_c
7217
7218end module mod_twofl_phys
subroutine twofl_get_p_c_total(w, x, ixil, ixol, p)
subroutine twofl_get_csound2_primitive(w, x, ixil, ixol, csound2)
subroutine add_density_hyper_source(index_rho)
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter bigdouble
A very large real number.
subroutine reconstruct(ixil, ixcl, idir, q, ql, qr)
Reconstruct scalar q within ixO^L to 1/2 dx in direction idir Return both left and right reconstructe...
subroutine b_from_vector_potentiala(ixisl, ixil, ixol, ws, x, a)
calculate magnetic field from vector potential A at cell edges
subroutine add_convert_method(phys_convert_vars, nwc, dataset_names, file_suffix)
Definition mod_convert.t:59
Module for flux conservation near refinement boundaries.
Module with basic grid data structures.
Definition mod_forest.t:2
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Definition mod_forest.t:32
subroutine, public get_divb(w, ixil, ixol, divb, nth_in)
Calculate div B within ixO.
integer, dimension(:), allocatable, public mag
Indices of the magnetic field.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
subroutine divvector(qvec, ixil, ixol, divq, nth_in)
integer coordinate
Definition mod_geometry.t:7
integer, parameter spherical
integer, parameter cylindrical
subroutine curlvector(qvec, ixil, ixol, curlvec, idirmin, idirmin0, ndir0, fourthorder)
Calculate curl of a vector qvec within ixL Options to employ standard second order CD evaluations use...
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
subroutine gradientf(q, x, ixil, ixol, idir, gradq, nth_in, pm_in)
subroutine gradientl(q, ixil, ixol, idir, gradq)
update ghost cells of all blocks including physical boundaries
subroutine getbc(time, qdt, psb, nwstart, nwbc)
do update ghost cells of all blocks including physical boundaries
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.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
double precision dtdiffpar
For resistive MHD, the time step is also limited by the diffusion time: .
character(len=std_len) typegrad
double precision unit_charge
Physical scaling factor for charge.
integer ixghi
Upper index of grid block arrays.
pure subroutine cross_product(ixil, ixol, a, b, axb)
Cross product of two vectors.
integer, dimension(3, 3, 3) lvc
Levi-Civita tensor.
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
double precision global_time
The global simulation time.
double precision unit_mass
Physical scaling factor for mass.
logical use_imex_scheme
whether IMEX in use or not
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer it
Number of time steps taken.
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
double precision unit_numberdensity
Physical scaling factor for number density.
character(len=std_len) convert_type
Which format to use when converting.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
logical stagger_grid
True for using stagger grid.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer icomm
The MPI communicator.
double precision bdip
amplitude of background dipolar, quadrupolar, octupolar, user's field
integer b0i
background magnetic field location indicator
integer mype
The rank of the current MPI task.
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 autoconvert
If true, already convert to output format during the run.
logical slab
Cartesian geometry or not.
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision unit_magneticfield
Physical scaling factor for magnetic field.
double precision unit_velocity
Physical scaling factor for velocity.
double precision c_norm
Normalised speed of light.
logical b0field
split magnetic field as background B0 field
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision unit_temperature
Physical scaling factor for temperature.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer nghostcells
Number of ghost cells surrounding a grid.
integer, parameter sdim
starting dimension for electric field
logical phys_trac
Use TRAC for MHD or 1D HD.
logical need_global_cmax
need global maximal wave speed
logical convert
If true and restart_from_file is given, convert snapshots to other file formats.
logical fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical use_multigrid
Use multigrid (only available in 2D and 3D)
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer max_blocks
The maximum number of grid blocks in a processor.
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
integer number_equi_vars
number of equilibrium set variables, besides the mag field
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
integer, parameter ixglo
Lower index of grid block arrays (always 1)
Subroutines for Roe-type Riemann solver for HD.
subroutine second_same_deriv(ixil, ixol, nu_hyper, var, idimm, res)
subroutine div_vel_coeff(ixil, ixol, vel, idimm, nu_vel)
subroutine second_cross_deriv2(ixil, ixol, nu_hyper, var2, var, idimm, idimm2, res)
subroutine hyp_coeff(ixil, ixol, var, idimm, nu_hyp)
subroutine second_cross_deriv(ixil, ixol, nu_hyper, var, idimm, idimm2, res)
subroutine second_same_deriv2(ixil, ixol, nu_hyper, var2, var, idimm, res)
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
type(mg_t) mg
Data structure containing the multigrid tree.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
module radiative cooling – add optically thin radiative cooling for HD and MHD
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine cooling_get_dt(w, ixil, ixol, dtnew, dxd, x, fl)
subroutine radiative_cooling_init(fl, read_params)
subroutine radiative_cooling_add_source(qdt, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active, fl)
Module for handling problematic values in simulations, such as negative pressures.
subroutine, public small_values_average(ixil, ixol, w, x, w_flag, windex)
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_error(wprim, x, ixil, ixol, w_flag, subname)
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
character(len=20), public small_values_method
How to handle small values.
Generic supertimestepping method 1) in amrvac.par in sts_list set the following parameters which have...
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...
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...
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.
Thermal conduction for HD and MHD or RHD and RMHD or twofl (plasma-neutral) module Adaptation of mod_...
subroutine, public tc_get_hd_params(fl, read_hd_params)
Init TC coefficients: HD case.
double precision function, public get_tc_dt_mhd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (mhd implementation) Note: for multi-D MHD (1D MHD will use HD f...
double precision function, public get_tc_dt_hd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (hd implementation) Note: also used in 1D MHD (or for neutrals i...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_hd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
subroutine, public sts_set_source_tc_mhd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
anisotropic thermal conduction with slope limited symmetric scheme Sharma 2007 Journal of Computation...
subroutine, public tc_get_mhd_params(fl, read_mhd_params)
Init TC coefficients: MHD case.
subroutine get_euv_image(qunit, fl)
subroutine get_sxr_image(qunit, fl)
subroutine get_euv_spectrum(qunit, fl)
subroutine get_whitelight_image(qunit, fl)
Magneto-hydrodynamics module.
logical, public twofl_coll_inc_ionrec
whether include ionization/recombination inelastic collisional terms
subroutine, public b_from_vector_potential(ixisl, ixil, ixol, ws, x)
calculate magnetic field from vector potential
logical, public, protected twofl_dump_full_vars
whether dump full variables (when splitting is used) in a separate dat file
subroutine, public get_gamma_ion_rec(ixil, ixol, w, x, gamma_rec, gamma_ion)
subroutine, public twofl_get_v_c_idim(w, x, ixil, ixol, idim, v)
Calculate v_c component.
double precision, public, protected rn
logical, public clean_initial_divb
clean initial divB
double precision, public twofl_eta_hyper
The MHD hyper-resistivity.
subroutine, public twofl_get_csound2_c_from_conserved(w, x, ixil, ixol, csound2)
double precision, public twofl_eta
The MHD resistivity.
integer, public, protected twofl_trac_type
Which TRAC method is used
logical, public has_equi_pe_c0
type(tc_fluid), allocatable, public tc_fl_c
logical, public twofl_alpha_coll_constant
double precision, dimension(:), allocatable, public, protected c_shk
subroutine, public twofl_face_to_center(ixol, s)
calculate cell-center values from face-center values
logical, public, protected twofl_dump_hyperdiffusivity_coef
double precision, public twofl_glm_alpha
GLM-MHD parameter: ratio of the diffusive and advective time scales for div b taking values within [0...
integer, parameter, public eq_energy_ki
logical, public, protected twofl_thermal_conduction_n
subroutine, public twofl_phys_init()
subroutine, public get_rhon_tot(w, x, ixil, ixol, rhon)
logical, public, protected twofl_thermal_conduction_c
Whether thermal conduction is used.
double precision, public twofl_adiab
The adiabatic constant.
logical, public twofl_equi_thermal_c
character(len=std_len), public, protected type_ct
Method type of constrained transport.
integer, public tweight_c_
subroutine, public twofl_get_pthermal_c(w, x, ixil, ixol, pth)
logical, public, protected twofl_radiative_cooling_n
integer, parameter, public eq_energy_none
integer, public e_n_
type(te_fluid), allocatable, public te_fl_c
procedure(mask_subroutine2), pointer, public usr_mask_gamma_ion_rec
double precision, public, protected rc
logical, public, protected twofl_dump_coll_terms
whether dump collisional terms in a separte dat file
logical, public twofl_equi_thermal_n
logical, public, protected twofl_radiative_cooling_c
Whether radiative cooling is added.
logical, public, protected b0field_forcefree
B0 field is force-free.
integer, public e_c_
Index of the energy density (-1 if not present)
integer, public equi_rho_n0_
integer, public, protected c
Indices of the momentum density for the form of better vectorization.
subroutine, public twofl_get_v_n_idim(w, x, ixil, ixol, idim, v)
Calculate v component.
integer, parameter, public eq_energy_int
integer, dimension(:), allocatable, public mom_c
Indices of the momentum density.
logical, public twofl_coll_inc_te
whether include thermal exchange collisional terms
logical, public has_equi_rho_c0
equi vars flags
logical, public, protected twofl_viscosity
Whether viscosity is added.
double precision, public dtcollpar
logical, public divbwave
Add divB wave in Roe solver.
subroutine, public twofl_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
logical, public, protected twofl_4th_order
MHD fourth order.
integer, public tcoff_n_
double precision, dimension(:), allocatable, public, protected c_hyp
integer, public equi_rho_c0_
equi vars indices in the stateequi_vars array
logical, public, protected twofl_glm
Whether GLM-MHD is used.
double precision, public twofl_alpha_coll
collisional alpha
logical, public, protected twofl_trac
Whether TRAC method is used.
double precision, public twofl_etah
The MHD Hall coefficient.
subroutine, public get_alpha_coll(ixil, ixol, w, x, alpha)
logical, dimension(2 *^nd), public, protected boundary_divbfix
To control divB=0 fix for boundary.
integer, public, protected b
integer, public equi_pe_c0_
integer, parameter, public eq_energy_tot
integer, dimension(:), allocatable, public mom_n
logical, public, protected twofl_gravity
Whether gravity is added: common flag for charges and neutrals.
integer, public tcoff_c_
Index of the cutoff temperature for the TRAC method.
subroutine, public twofl_clean_divb_multigrid(qdt, qt, active)
double precision, public, protected twofl_trac_mask
Height of the mask used in the TRAC method.
logical, public has_equi_pe_n0
procedure(mask_subroutine), pointer, public usr_mask_alpha
integer, public, protected twofl_divb_nth
Whether divB is computed with a fourth order approximation.
integer, public rho_c_
Index of the density (in the w array)
subroutine, public get_normalized_divb(w, ixil, ixol, divb)
get dimensionless div B = |divB| * volume / area / |B|
type(rc_fluid), allocatable, public rc_fl_c
subroutine, public get_current(w, ixil, ixol, idirmin, current)
Calculate idirmin and the idirmin:3 components of the common current array make sure that dxlevel(^D)...
subroutine, public twofl_to_conserved(ixil, ixol, w, x)
Transform primitive variables into conservative ones.
logical, public twofl_equi_thermal
integer, public, protected c_
logical, public has_equi_rho_n0
subroutine, public get_rhoc_tot(w, x, ixil, ixol, rhoc)
integer, public rho_n_
subroutine, public twofl_get_pthermal_n(w, x, ixil, ixol, pth)
logical, public, protected twofl_hyperdiffusivity
Whether hyperdiffusivity is used.
integer, public, protected twofl_eq_energy
integer, dimension(2 *^nd), public, protected boundary_divbfix_skip
To skip * layer of ghost cells during divB=0 fix for boundary.
integer, public, protected m
character(len=std_len), public, protected typedivbfix
Method type to clean divergence of B.
double precision, public twofl_gamma
The adiabatic index.
integer, public equi_pe_n0_
logical, public, protected twofl_hall
Whether Hall-MHD is used.
integer, public tweight_n_
integer, public, protected psi_
Indices of the GLM psi.
logical, public, protected source_split_divb
Whether divB cleaning sources are added splitting from fluid solver.
Module with all the methods that users can customize in AMRVAC.
procedure(special_resistivity), pointer usr_special_resistivity
procedure(phys_gravity), pointer usr_gravity
procedure(set_equi_vars), pointer usr_set_equi_vars
procedure(set_electric_field), pointer usr_set_electric_field
procedure(set_wlr), pointer usr_set_wlr
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
The module add viscous source terms and check time step.
subroutine viscosity_init(phys_wider_stencil)
Initialize the module.
procedure(sub_add_source), pointer, public viscosity_add_source
The data structure that contains information about a tree node/grid block.
Definition mod_forest.t:11