MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
Loading...
Searching...
No Matches
mod_mf_phys.t
Go to the documentation of this file.
1!> Magnetofriction module
3 use mod_global_parameters, only: std_len
5 use mod_comm_lib, only: mpistop
6
7
8 implicit none
9 private
10
11 !> viscosity coefficient in s cm^-2 for solar corona (Cheung 2012 ApJ)
12 double precision, public :: mf_nu = 1.d-15
13
14 !> maximal limit of magnetofrictional velocity in cm s^-1 (Pomoell 2019 A&A)
15 double precision, public :: mf_vmax = 3.d6
16
17 !> decay scale of frictional velocity near boundaries
18 double precision, public :: mf_decay_scale(2*^nd)=0.d0
19
20 !> GLM-MHD parameter: ratio of the diffusive and advective time scales for div b
21 !> taking values within [0, 1]
22 double precision, public :: mf_glm_alpha = 0.5d0
23
24 !> The resistivity
25 double precision, public :: mf_eta = 0.0d0
26
27 !> The hyper-resistivity
28 double precision, public :: mf_eta_hyper = 0.0d0
29
30 !> Coefficient of diffusive divB cleaning
31 double precision :: divbdiff = 0.8d0
32
33 !> Helium abundance over Hydrogen
34 double precision, public, protected :: he_abundance=0.1d0
35
36 !> Method type in a integer for good performance
37 integer :: type_divb
38
39 !> Indices of the momentum density
40 integer, allocatable, public, protected :: mom(:)
41
42 !> Indices of the momentum density for the form of better vectorization
43 integer, public, protected :: ^c&m^C_
44
45 !> Indices of the magnetic field for the form of better vectorization
46 integer, public, protected :: ^c&b^C_
47
48 !> Indices of the GLM psi
49 integer, public, protected :: psi_
50
51 !> To skip * layer of ghost cells during divB=0 fix for boundary
52 integer, public, protected :: boundary_divbfix_skip(2*^nd)=0
53
54 ! DivB cleaning methods
55 integer, parameter :: divb_none = 0
56 integer, parameter :: divb_multigrid = -1
57 integer, parameter :: divb_glm = 1
58 integer, parameter :: divb_powel = 2
59 integer, parameter :: divb_janhunen = 3
60 integer, parameter :: divb_linde = 4
61 integer, parameter :: divb_lindejanhunen = 5
62 integer, parameter :: divb_lindepowel = 6
63 integer, parameter :: divb_lindeglm = 7
64 integer, parameter :: divb_ct = 8
65
66 !> Whether particles module is added
67 logical, public, protected :: mf_particles = .false.
68
69 !> Whether GLM-MHD is used
70 logical, public, protected :: mf_glm = .false.
71
72 !> Whether divB cleaning sources are added splitting from fluid solver
73 logical, public, protected :: source_split_divb = .false.
74
75 !> MHD fourth order
76 logical, public, protected :: mf_4th_order = .false.
77
78 !> set to true if need to record electric field on cell edges
79 logical, public, protected :: mf_record_electric_field = .false.
80
81 !> Whether divB is computed with a fourth order approximation
82 integer, public, protected :: mf_divb_nth = 1
83
84 !> To control divB=0 fix for boundary
85 logical, public, protected :: boundary_divbfix(2*^nd)=.true.
86
87 !> Use a compact way to add resistivity
88 logical :: compactres = .false.
89
90 !> Add divB wave in Roe solver
91 logical, public :: divbwave = .true.
92
93 !> clean divb in the initial condition
94 logical, public, protected :: clean_initial_divb=.false.
95
96 !> Method type to clean divergence of B
97 character(len=std_len), public, protected :: typedivbfix = 'ct'
98
99 !> Method type of constrained transport
100 character(len=std_len), public, protected :: type_ct = 'average'
101
102 !> Update all equations due to divB cleaning
103 character(len=std_len) :: typedivbdiff = 'all'
104
105
106 ! Public methods
107 public :: mf_phys_init
108 public :: mf_to_conserved
109 public :: mf_to_primitive
110 public :: mf_face_to_center
111 public :: get_divb
112 public :: get_current
113 public :: get_normalized_divb
116 {^nooned
118 }
119
120contains
121
122 !> Read this module's parameters from a file
123 subroutine mf_read_params(files)
125 use mod_particles, only: particles_eta
126 character(len=*), intent(in) :: files(:)
127 integer :: n
128
129 namelist /mf_list/ mf_nu, mf_vmax, mf_decay_scale, &
131 particles_eta, mf_record_electric_field,&
133 typedivbdiff, type_ct, compactres, divbwave, he_abundance, si_unit, &
136
137 do n = 1, size(files)
138 open(unitpar, file=trim(files(n)), status="old")
139 read(unitpar, mf_list, end=111)
140111 close(unitpar)
141 end do
142
143 end subroutine mf_read_params
144
145 !> Write this module's parameters to a snapsoht
146 subroutine mf_write_info(fh)
148 integer, intent(in) :: fh
149 integer, parameter :: n_par = 1
150 double precision :: values(n_par)
151 character(len=name_len) :: names(n_par)
152 integer, dimension(MPI_STATUS_SIZE) :: st
153 integer :: er
154
155 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
156
157 names(1) = "nu"
158 values(1) = mf_nu
159 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
160 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
161 end subroutine mf_write_info
162
163 subroutine mf_phys_init()
165 use mod_physics
166 use mod_particles, only: particles_init, particles_eta, particles_etah
167 {^nooned
169 }
170
171 integer :: itr, idir
172
173 call mf_read_params(par_files)
174
175 physics_type = "mf"
176 phys_energy = .false.
178
179 if(ndim==1) typedivbfix='none'
180 select case (typedivbfix)
181 case ('none')
182 type_divb = divb_none
183 {^nooned
184 case ('multigrid')
185 type_divb = divb_multigrid
186 use_multigrid = .true.
187 mg%operator_type = mg_laplacian
189 }
190 case ('glm')
191 mf_glm = .true.
192 need_global_cmax = .true.
193 type_divb = divb_glm
194 case ('powel', 'powell')
195 type_divb = divb_powel
196 case ('janhunen')
197 type_divb = divb_janhunen
198 case ('linde')
199 type_divb = divb_linde
200 case ('lindejanhunen')
201 type_divb = divb_lindejanhunen
202 case ('lindepowel')
203 type_divb = divb_lindepowel
204 case ('lindeglm')
205 mf_glm = .true.
206 need_global_cmax = .true.
207 type_divb = divb_lindeglm
208 case ('ct')
209 type_divb = divb_ct
210 stagger_grid = .true.
211 case default
212 call mpistop('Unknown divB fix')
213 end select
214
215 ! set velocity field as flux variables
216 allocate(mom(ndir))
217 mom(:) = var_set_momentum(ndir)
218 m^c_=mom(^c);
219
220 ! set magnetic field as flux variables
221 allocate(mag(ndir))
222 mag(:) = var_set_bfield(ndir)
223 b^c_=mag(^c);
224
225 ! start with magnetic field and skip velocity when update ghostcells
226 iwstart=mag(1)
227 allocate(start_indices(number_species),stop_indices(number_species))
228 ! set the index of the first flux variable for species 1
229 start_indices(1)=mag(1)
230
231 if (mf_glm) then
232 psi_ = var_set_fluxvar('psi', 'psi', need_bc=.false.)
233 else
234 psi_ = -1
235 end if
236
237 ! set number of variables which need update ghostcells
238 nwgc=nwflux-ndir
239
240 ! set the index of the last flux variable for species 1
241 stop_indices(1)=nwflux
242
243 ! determine number of stagger variables
244 nws=ndim
245
246 nvector = 2 ! No. vector vars
247 allocate(iw_vector(nvector))
248 iw_vector(1) = mom(1) - 1 ! TODO: why like this?
249 iw_vector(2) = mag(1) - 1 ! TODO: why like this?
250
251 ! Check whether custom flux types have been defined
252 if (.not. allocated(flux_type)) then
253 allocate(flux_type(ndir, nw))
255 else if (any(shape(flux_type) /= [ndir, nw])) then
256 call mpistop("phys_check error: flux_type has wrong shape")
257 end if
258 do idir=1,ndir
259 if(ndim>1) flux_type(idir,mag(idir))=flux_tvdlf
260 end do
261 if(mf_glm .and. ndim>1) flux_type(:,psi_)=flux_tvdlf
262
263 phys_get_dt => mf_get_dt
264 phys_get_cmax => mf_get_cmax
265 phys_get_cbounds => mf_get_cbounds
266 phys_get_flux => mf_get_flux
267 phys_add_source_geom => mf_add_source_geom
268 phys_add_source => mf_add_source
271 phys_check_params => mf_check_params
272 phys_write_info => mf_write_info
273 phys_special_advance => mf_velocity_update
274
275 if(type_divb==divb_glm) then
276 phys_modify_wlr => mf_modify_wlr
277 end if
278
279 ! pass to global variable to record electric field
281
282 ! Initialize particles module
283 if(mf_particles) then
284 call particles_init()
285 ! never allow Hall effects in particles when doing magnetofrictional
286 particles_etah=0.0d0
287 if(mype==0) then
288 write(*,*) '*****Using particles: with mf_eta, mf_eta_hyper :', mf_eta, mf_eta_hyper
289 write(*,*) '*****Using particles: particles_eta, particles_etah :', particles_eta, particles_etah
290 end if
291 end if
292
293 ! if using ct stagger grid, boundary divb=0 is not done here
294 if(stagger_grid) then
296 select case(type_ct)
297 case('average')
299 phys_get_ct_velocity => mf_get_ct_velocity_average
300 phys_update_faces => mf_update_faces_average
301 case('uct_contact')
303 phys_get_ct_velocity => mf_get_ct_velocity_contact
304 phys_update_faces => mf_update_faces_contact
305 case('uct_hll')
307 phys_get_ct_velocity => mf_get_ct_velocity_hll
308 phys_update_faces => mf_update_faces_hll
309 case default
310 call mpistop('choose average, uct_contact,or uct_hll for type_ct!')
311 end select
312 else if(ndim>1) then
313 phys_boundary_adjust => mf_boundary_adjust
314 end if
315
316 {^nooned
317 ! clean initial divb
319 }
320
321 ! derive units from basic units
322 call mf_physical_units()
323
324 end subroutine mf_phys_init
325
326 subroutine mf_check_params
328
329 end subroutine mf_check_params
330
331 subroutine mf_physical_units()
333 double precision :: mp,kb,miu0,c_lightspeed
334 double precision :: a,b
335 ! Derive scaling units
336 if(si_unit) then
337 mp=mp_si
338 kb=kb_si
339 miu0=miu0_si
340 c_lightspeed=c_si
341 else
342 mp=mp_cgs
343 kb=kb_cgs
344 miu0=4.d0*dpi
345 c_lightspeed=const_c
346 end if
347 a=1.d0+4.d0*he_abundance
348 b=2.d0+3.d0*he_abundance
349 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
350 if(unit_density/=1.d0) then
352 else if(unit_numberdensity/=1.d0) then
354 end if
355 if(unit_temperature/=1.d0) then
359 if(unit_length/=1.d0) then
361 else if(unit_time/=1.d0) then
363 end if
364 else if(unit_magneticfield/=1.d0) then
368 if(unit_length/=1.d0) then
370 else if(unit_time/=1.d0) then
372 end if
373 else if(unit_pressure/=1.d0) then
377 if(unit_length/=1.d0) then
379 else if(unit_time/=1.d0) then
381 end if
382 else if(unit_velocity/=1.d0) then
386 if(unit_length/=1.d0) then
388 else if(unit_time/=1.d0) then
390 end if
391 else if(unit_time/=1.d0) then
396 end if
397 else if(unit_temperature/=1.d0) then
398 ! units of temperature and velocity are dependent
399 if(unit_magneticfield/=1.d0) then
404 if(unit_length/=1.d0) then
406 else if(unit_time/=1.d0) then
408 end if
409 else if(unit_pressure/=1.d0) then
414 if(unit_length/=1.d0) then
416 else if(unit_time/=1.d0) then
418 end if
419 end if
420 else if(unit_magneticfield/=1.d0) then
421 ! units of magnetic field and pressure are dependent
422 if(unit_velocity/=1.d0) then
427 if(unit_length/=1.d0) then
429 else if(unit_time/=1.d0) then
431 end if
432 else if(unit_time/=0.d0) then
438 end if
439 else if(unit_pressure/=1.d0) then
440 if(unit_velocity/=1.d0) then
445 if(unit_length/=1.d0) then
447 else if(unit_time/=1.d0) then
449 end if
450 else if(unit_time/=0.d0) then
456 end if
457 end if
458 ! Additional units needed for the particles
459 c_norm=c_lightspeed/unit_velocity
461 if (.not. si_unit) unit_charge = unit_charge*const_c
463
464 ! get dimensionless mf nu
466
467 ! get dimensionless maximal mf velocity limit
469
470 end subroutine mf_physical_units
471
472 !> Transform primitive variables into conservative ones
473 subroutine mf_to_conserved(ixI^L,ixO^L,w,x)
475 integer, intent(in) :: ixi^l, ixo^l
476 double precision, intent(inout) :: w(ixi^s, nw)
477 double precision, intent(in) :: x(ixi^s, 1:ndim)
478
479 ! nothing to do for mf
480 end subroutine mf_to_conserved
481
482 !> Transform conservative variables into primitive ones
483 subroutine mf_to_primitive(ixI^L,ixO^L,w,x)
485 integer, intent(in) :: ixi^l, ixo^l
486 double precision, intent(inout) :: w(ixi^s, nw)
487 double precision, intent(in) :: x(ixi^s, 1:ndim)
488
489 ! nothing to do for mf
490 end subroutine mf_to_primitive
491
492 !> Calculate cmax_idim=csound+abs(v_idim) within ixO^L
493 subroutine mf_get_cmax(w,x,ixI^L,ixO^L,idim,cmax)
495
496 integer, intent(in) :: ixi^l, ixo^l, idim
497 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
498 double precision, intent(inout) :: cmax(ixi^s)
499
500 cmax(ixo^s)=abs(w(ixo^s,mom(idim)))+one
501
502 end subroutine mf_get_cmax
503
504 !> Estimating bounds for the minimum and maximum signal velocities
505 subroutine mf_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
507 use mod_variables
508
509 integer, intent(in) :: ixi^l, ixo^l, idim
510 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
511 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
512 double precision, intent(in) :: x(ixi^s,1:ndim)
513 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
514 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
515 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
516
517 double precision, dimension(ixI^S) :: tmp1
518
519 tmp1(ixo^s)=0.5d0*(wlc(ixo^s,mom(idim))+wrc(ixo^s,mom(idim)))
520 if(present(cmin)) then
521 cmax(ixo^s,1)=max(tmp1(ixo^s)+one,zero)
522 cmin(ixo^s,1)=min(tmp1(ixo^s)-one,zero)
523 else
524 cmax(ixo^s,1)=abs(tmp1(ixo^s))+one
525 end if
526
527 end subroutine mf_get_cbounds
528
529 !> prepare velocities for ct methods
530 subroutine mf_get_ct_velocity_average(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
532
533 integer, intent(in) :: ixi^l, ixo^l, idim
534 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
535 double precision, intent(in) :: cmax(ixi^s)
536 double precision, intent(in), optional :: cmin(ixi^s)
537 type(ct_velocity), intent(inout):: vcts
538
539 end subroutine mf_get_ct_velocity_average
540
541 subroutine mf_get_ct_velocity_contact(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
543
544 integer, intent(in) :: ixi^l, ixo^l, idim
545 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
546 double precision, intent(in) :: cmax(ixi^s)
547 double precision, intent(in), optional :: cmin(ixi^s)
548 type(ct_velocity), intent(inout):: vcts
549
550 ! calculate velocities related to different UCT schemes
551 if(.not.allocated(vcts%vnorm)) allocate(vcts%vnorm(ixi^s,1:ndim))
552 ! get average normal velocity at cell faces
553 vcts%vnorm(ixo^s,idim)=0.5d0*(wlp(ixo^s,mom(idim))+wrp(ixo^s,mom(idim)))
554
555 end subroutine mf_get_ct_velocity_contact
556
557 subroutine mf_get_ct_velocity_hll(vcts,wLp,wRp,ixI^L,ixO^L,idim,cmax,cmin)
559
560 integer, intent(in) :: ixi^l, ixo^l, idim
561 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
562 double precision, intent(in) :: cmax(ixi^s)
563 double precision, intent(in), optional :: cmin(ixi^s)
564 type(ct_velocity), intent(inout):: vcts
565
566 integer :: idime,idimn
567
568 if(.not.allocated(vcts%vbarC)) then
569 allocate(vcts%vbarC(ixi^s,1:ndir,2),vcts%vbarLC(ixi^s,1:ndir,2),vcts%vbarRC(ixi^s,1:ndir,2))
570 allocate(vcts%cbarmin(ixi^s,1:ndim),vcts%cbarmax(ixi^s,1:ndim))
571 end if
572 ! Store magnitude of characteristics
573 if(present(cmin)) then
574 vcts%cbarmin(ixo^s,idim)=max(-cmin(ixo^s),zero)
575 vcts%cbarmax(ixo^s,idim)=max( cmax(ixo^s),zero)
576 else
577 vcts%cbarmax(ixo^s,idim)=max( cmax(ixo^s),zero)
578 vcts%cbarmin(ixo^s,idim)=vcts%cbarmax(ixo^s,idim)
579 end if
580
581 idimn=mod(idim,ndir)+1 ! 'Next' direction
582 idime=mod(idim+1,ndir)+1 ! Electric field direction
583 ! Store velocities
584 vcts%vbarLC(ixo^s,idim,1)=wlp(ixo^s,mom(idimn))
585 vcts%vbarRC(ixo^s,idim,1)=wrp(ixo^s,mom(idimn))
586 vcts%vbarC(ixo^s,idim,1)=(vcts%cbarmax(ixo^s,idim)*vcts%vbarLC(ixo^s,idim,1) &
587 +vcts%cbarmin(ixo^s,idim)*vcts%vbarRC(ixo^s,idim,1))&
588 /(vcts%cbarmax(ixo^s,idim)+vcts%cbarmin(ixo^s,idim))
589
590 vcts%vbarLC(ixo^s,idim,2)=wlp(ixo^s,mom(idime))
591 vcts%vbarRC(ixo^s,idim,2)=wrp(ixo^s,mom(idime))
592 vcts%vbarC(ixo^s,idim,2)=(vcts%cbarmax(ixo^s,idim)*vcts%vbarLC(ixo^s,idim,2) &
593 +vcts%cbarmin(ixo^s,idim)*vcts%vbarRC(ixo^s,idim,1))&
594 /(vcts%cbarmax(ixo^s,idim)+vcts%cbarmin(ixo^s,idim))
595
596 end subroutine mf_get_ct_velocity_hll
597
598 !> Calculate fluxes within ixO^L.
599 subroutine mf_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
602 use mod_geometry
603
604 integer, intent(in) :: ixi^l, ixo^l, idim
605 ! conservative w
606 double precision, intent(in) :: wc(ixi^s,nw)
607 ! primitive w
608 double precision, intent(in) :: w(ixi^s,nw)
609 double precision, intent(in) :: x(ixi^s,1:ndim)
610 double precision,intent(out) :: f(ixi^s,nwflux)
611
612 integer :: ix^d
613
614 {do ix^db=ixomin^db,ixomax^db\}
615 ! compute flux of magnetic field
616 ! f_i[b_k]=v_i*b_k-v_k*b_i
617 ^c&f(ix^d,b^c_)=w(ix^d,mom(idim))*w(ix^d,b^c_)-w(ix^d,mag(idim))*w(ix^d,m^c_)\
618 {end do\}
619 if(mf_glm) then
620 {do ix^db=ixomin^db,ixomax^db\}
621 f(ix^d,mag(idim))=w(ix^d,psi_)
622 !f_i[psi]=Ch^2*b_{i} Eq. 24e and Eq. 38c Dedner et al 2002 JCP, 175, 645
623 f(ix^d,psi_)=cmax_global**2*w(ix^d,mag(idim))
624 {end do\}
625 end if
626
627 end subroutine mf_get_flux
628
629 !> Add global source terms to update frictional velocity on complete domain
630 subroutine mf_velocity_update(qt,psa)
633 double precision, intent(in) :: qt !< Current time
634 type(state), target :: psa(max_blocks) !< Compute based on this state
635
636 integer :: iigrid,igrid
637 logical :: stagger_flag
638 logical :: firstmf=.true.
639
640 if(firstmf) then
641 ! point bc mpi datatype to partial type for velocity field
649 firstmf=.false.
650 end if
651
652 !$OMP PARALLEL DO PRIVATE(igrid)
653 do iigrid=1,igridstail_active; igrid=igrids_active(iigrid);
654 block=>psa(igrid)
655 call frictional_velocity(psa(igrid)%w,psa(igrid)%x,ixg^ll,ixm^ll)
656 end do
657 !$OMP END PARALLEL DO
658
659 ! only update mf velocity in ghost cells
660 stagger_flag=stagger_grid
667 stagger_grid=.false.
668 bcphys=.false.
669 call getbc(qt,0.d0,psa,mom(1),ndir)
670 bcphys=.true.
677 stagger_grid=stagger_flag
678
679 end subroutine mf_velocity_update
680
681 !> w[iws]=w[iws]+qdt*S[iws,wCT] where S is the source based on wCT within ixO
682 subroutine mf_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
684
685 integer, intent(in) :: ixi^l, ixo^l
686 double precision, intent(in) :: qdt, dtfactor
687 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:ndim)
688 double precision, intent(inout) :: w(ixi^s,1:nw)
689 logical, intent(in) :: qsourcesplit
690 logical, intent(inout) :: active
691
692 if (.not. qsourcesplit) then
693
694 ! Sources for resistivity in eqs. for e, B1, B2 and B3
695 if (abs(mf_eta)>smalldouble)then
696 active = .true.
697 call add_source_res2(qdt,ixi^l,ixo^l,wct,w,x)
698 end if
699
700 if (mf_eta_hyper>0.d0)then
701 active = .true.
702 call add_source_hyperres(qdt,ixi^l,ixo^l,wct,w,x)
703 end if
704
705 end if
706
707 {^nooned
708 if(source_split_divb .eqv. qsourcesplit) then
709 ! Sources related to div B
710 select case (type_divb)
711 case (divb_none)
712 ! Do nothing
713 case (divb_glm)
714 active = .true.
715 call add_source_glm(qdt,ixi^l,ixo^l,wct,w,x)
716 case (divb_powel)
717 active = .true.
718 call add_source_powel(qdt,ixi^l,ixo^l,wct,w,x)
719 case (divb_janhunen)
720 active = .true.
721 call add_source_janhunen(qdt,ixi^l,ixo^l,wct,w,x)
722 case (divb_linde)
723 active = .true.
724 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
725 case (divb_lindejanhunen)
726 active = .true.
727 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
728 call add_source_janhunen(qdt,ixi^l,ixo^l,wct,w,x)
729 case (divb_lindepowel)
730 active = .true.
731 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
732 call add_source_powel(qdt,ixi^l,ixo^l,wct,w,x)
733 case (divb_lindeglm)
734 active = .true.
735 call add_source_linde(qdt,ixi^l,ixo^l,wct,w,x)
736 call add_source_glm(qdt,ixi^l,ixo^l,wct,w,x)
737 case (divb_ct)
738 continue ! Do nothing
739 case (divb_multigrid)
740 continue ! Do nothing
741 case default
742 call mpistop('Unknown divB fix')
743 end select
744 end if
745 }
746
747 end subroutine mf_add_source
748
749 subroutine frictional_velocity(w,x,ixI^L,ixO^L)
751
752 integer, intent(in) :: ixi^l, ixo^l
753 double precision, intent(in) :: x(ixi^s,1:ndim)
754 double precision, intent(inout) :: w(ixi^s,1:nw)
755
756 double precision :: xmin(ndim),xmax(ndim)
757 double precision :: decay(ixo^s)
758 double precision :: current(ixi^s,7-2*ndir:3),tmp(ixo^s)
759 integer :: ix^d, idirmin,idir,jdir,kdir
760
761 call get_current(w,ixi^l,ixo^l,idirmin,current)
762 ! extrapolate current for the outmost layer,
763 ! because extrapolation of B at boundaries introduces artificial current
764{
765 if(block%is_physical_boundary(2*^d-1)) then
766 if(slab) then
767 ! exclude boundary 5 in 3D Cartesian
768 if(2*^d-1==5.and.ndim==3) then
769 else
770 current(ixomin^d^d%ixO^s,:)=current(ixomin^d+1^d%ixO^s,:)
771 end if
772 else
773 ! exclude boundary 1 spherical/cylindrical
774 if(2*^d-1>1) then
775 current(ixomin^d^d%ixO^s,:)=current(ixomin^d+1^d%ixO^s,:)
776 end if
777 end if
778 end if
779 if(block%is_physical_boundary(2*^d)) then
780 current(ixomax^d^d%ixO^s,:)=current(ixomax^d-1^d%ixO^s,:)
781 end if
782\}
783 w(ixo^s,mom(:))=0.d0
784 ! calculate Lorentz force
785 do idir=1,ndir; do jdir=idirmin,3; do kdir=1,ndir
786 if(lvc(idir,jdir,kdir)/=0) then
787 if(lvc(idir,jdir,kdir)==1) then
788 w(ixo^s,mom(idir))=w(ixo^s,mom(idir))+current(ixo^s,jdir)*w(ixo^s,mag(kdir))
789 else
790 w(ixo^s,mom(idir))=w(ixo^s,mom(idir))-current(ixo^s,jdir)*w(ixo^s,mag(kdir))
791 end if
792 end if
793 end do; end do; end do
794
795 tmp(ixo^s)=sum(w(ixo^s,mag(:))**2,dim=ndim+1)
796 ! frictional coefficient
797 where(tmp(ixo^s)/=0.d0)
798 tmp(ixo^s)=1.d0/(tmp(ixo^s)*mf_nu)
799 endwhere
800
801 do idir=1,ndir
802 w(ixo^s,mom(idir))=w(ixo^s,mom(idir))*tmp(ixo^s)
803 end do
804
805 ! decay frictional velocity near selected boundaries
806 ^d&xmin(^d)=xprobmin^d\
807 ^d&xmax(^d)=xprobmax^d\
808 decay(ixo^s)=1.d0
809 do idir=1,ndim
810 if(mf_decay_scale(2*idir-1)>0.d0) decay(ixo^s)=decay(ixo^s)*&
811 (1.d0-exp(-(x(ixo^s,idir)-xmin(idir))/mf_decay_scale(2*idir-1)))
812 if(mf_decay_scale(2*idir)>0.d0) decay(ixo^s)=decay(ixo^s)*&
813 (1.d0-exp((x(ixo^s,idir)-xmax(idir))/mf_decay_scale(2*idir)))
814 end do
815
816 ! saturate mf velocity at mf_vmax
817 tmp(ixo^s)=sqrt(sum(w(ixo^s,mom(:))**2,dim=ndim+1))/mf_vmax+1.d-12
818 tmp(ixo^s)=dtanh(tmp(ixo^s))/tmp(ixo^s)
819 do idir=1,ndir
820 w(ixo^s,mom(idir))=w(ixo^s,mom(idir))*tmp(ixo^s)*decay(ixo^s)
821 end do
822
823 end subroutine frictional_velocity
824
825 !> Add resistive source to w within ixO Uses 3 point stencil (1 neighbour) in
826 !> each direction, non-conservative. If the fourthorder precompiler flag is
827 !> set, uses fourth order central difference for the laplacian. Then the
828 !> stencil is 5 (2 neighbours).
829 subroutine add_source_res1(qdt,ixI^L,ixO^L,wCT,w,x)
832 use mod_geometry
833
834 integer, intent(in) :: ixi^l, ixo^l
835 double precision, intent(in) :: qdt
836 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
837 double precision, intent(inout) :: w(ixi^s,1:nw)
838
839 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
840 double precision :: current(ixi^s,7-2*ndir:3),eta(ixi^s)
841 double precision :: gradeta(ixi^s,1:ndim), bf(ixi^s,1:ndir)
842 double precision :: tmp(ixi^s),tmp2(ixi^s)
843 integer :: ixa^l,idir,jdir,kdir,idirmin,idim,jxo^l,hxo^l,ix
844 integer :: lxo^l, kxo^l
845
846 ! Calculating resistive sources involve one extra layer
847 if (mf_4th_order) then
848 ixa^l=ixo^l^ladd2;
849 else
850 ixa^l=ixo^l^ladd1;
851 end if
852
853 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
854 call mpistop("Error in add_source_res1: Non-conforming input limits")
855
856 ! Calculate current density and idirmin
857 call get_current(wct,ixi^l,ixo^l,idirmin,current)
858
859 if (mf_eta>zero)then
860 eta(ixa^s)=mf_eta
861 gradeta(ixo^s,1:ndim)=zero
862 else
863 call usr_special_resistivity(wct,ixi^l,ixa^l,idirmin,x,current,eta)
864 ! assumes that eta is not function of current?
865 do idim=1,ndim
866 call gradient(eta,ixi^l,ixo^l,idim,tmp)
867 gradeta(ixo^s,idim)=tmp(ixo^s)
868 end do
869 end if
870
871 bf(ixi^s,1:ndir)=wct(ixi^s,mag(1:ndir))
872
873 do idir=1,ndir
874 ! Put B_idir into tmp2 and eta*Laplace B_idir into tmp
875 if (mf_4th_order) then
876 tmp(ixo^s)=zero
877 tmp2(ixi^s)=bf(ixi^s,idir)
878 do idim=1,ndim
879 lxo^l=ixo^l+2*kr(idim,^d);
880 jxo^l=ixo^l+kr(idim,^d);
881 hxo^l=ixo^l-kr(idim,^d);
882 kxo^l=ixo^l-2*kr(idim,^d);
883 tmp(ixo^s)=tmp(ixo^s)+&
884 (-tmp2(lxo^s)+16.0d0*tmp2(jxo^s)-30.0d0*tmp2(ixo^s)+16.0d0*tmp2(hxo^s)-tmp2(kxo^s)) &
885 /(12.0d0 * dxlevel(idim)**2)
886 end do
887 else
888 tmp(ixo^s)=zero
889 tmp2(ixi^s)=bf(ixi^s,idir)
890 do idim=1,ndim
891 jxo^l=ixo^l+kr(idim,^d);
892 hxo^l=ixo^l-kr(idim,^d);
893 tmp(ixo^s)=tmp(ixo^s)+&
894 (tmp2(jxo^s)-2.0d0*tmp2(ixo^s)+tmp2(hxo^s))/dxlevel(idim)**2
895 end do
896 end if
897
898 ! Multiply by eta
899 tmp(ixo^s)=tmp(ixo^s)*eta(ixo^s)
900
901 ! Subtract grad(eta) x J = eps_ijk d_j eta J_k if eta is non-constant
902 if (mf_eta<zero)then
903 do jdir=1,ndim; do kdir=idirmin,3
904 if (lvc(idir,jdir,kdir)/=0)then
905 if (lvc(idir,jdir,kdir)==1)then
906 tmp(ixo^s)=tmp(ixo^s)-gradeta(ixo^s,jdir)*current(ixo^s,kdir)
907 else
908 tmp(ixo^s)=tmp(ixo^s)+gradeta(ixo^s,jdir)*current(ixo^s,kdir)
909 end if
910 end if
911 end do; end do
912 end if
913
914 ! Add sources related to eta*laplB-grad(eta) x J to B and e
915 w(ixo^s,mag(idir))=w(ixo^s,mag(idir))+qdt*tmp(ixo^s)
916
917 end do ! idir
918
919 end subroutine add_source_res1
920
921 !> Add resistive source to w within ixO
922 !> Uses 5 point stencil (2 neighbours) in each direction, conservative
923 subroutine add_source_res2(qdt,ixI^L,ixO^L,wCT,w,x)
926 use mod_geometry
927
928 integer, intent(in) :: ixi^l, ixo^l
929 double precision, intent(in) :: qdt
930 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
931 double precision, intent(inout) :: w(ixi^s,1:nw)
932
933 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
934 double precision :: current(ixi^s,7-2*ndir:3),eta(ixi^s),curlj(ixi^s,1:3)
935 double precision :: tmpvec(ixi^s,1:3)
936 integer :: ixa^l,idir,idirmin,idirmin1
937
938 ! resistive term already added as an electric field component
939 if(stagger_grid .and. ndim==ndir) return
940
941 ixa^l=ixo^l^ladd2;
942
943 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
944 call mpistop("Error in add_source_res2: Non-conforming input limits")
945
946 ixa^l=ixo^l^ladd1;
947 ! Calculate current density within ixL: J=curl B, thus J_i=eps_ijk*d_j B_k
948 ! Determine exact value of idirmin while doing the loop.
949 call get_current(wct,ixi^l,ixa^l,idirmin,current)
950
951 if (mf_eta>zero)then
952 eta(ixa^s)=mf_eta
953 else
954 call usr_special_resistivity(wct,ixi^l,ixa^l,idirmin,x,current,eta)
955 end if
956
957 ! dB/dt= -curl(J*eta), thus B_i=B_i-eps_ijk d_j Jeta_k
958 tmpvec(ixa^s,1:ndir)=zero
959 do idir=idirmin,3
960 tmpvec(ixa^s,idir)=current(ixa^s,idir)*eta(ixa^s)
961 end do
962 curlj=0.d0
963 call curlvector(tmpvec,ixi^l,ixo^l,curlj,idirmin1,1,3)
964 if(stagger_grid) then
965 if(ndim==2.and.ndir==3) then
966 ! if 2.5D
967 w(ixo^s,mag(ndir)) = w(ixo^s,mag(ndir))-qdt*curlj(ixo^s,ndir)
968 end if
969 else
970 w(ixo^s,mag(1:ndir)) = w(ixo^s,mag(1:ndir))-qdt*curlj(ixo^s,1:ndir)
971 end if
972
973 end subroutine add_source_res2
974
975 !> Add Hyper-resistive source to w within ixO
976 !> Uses 9 point stencil (4 neighbours) in each direction.
977 subroutine add_source_hyperres(qdt,ixI^L,ixO^L,wCT,w,x)
979 use mod_geometry
980
981 integer, intent(in) :: ixi^l, ixo^l
982 double precision, intent(in) :: qdt
983 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
984 double precision, intent(inout) :: w(ixi^s,1:nw)
985 !.. local ..
986 double precision :: current(ixi^s,7-2*ndir:3)
987 double precision :: tmpvec(ixi^s,1:3),tmpvec2(ixi^s,1:3),tmp(ixi^s),ehyper(ixi^s,1:3)
988 integer :: ixa^l,idir,jdir,kdir,idirmin,idirmin1
989
990 ixa^l=ixo^l^ladd3;
991 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
992 call mpistop("Error in add_source_hyperres: Non-conforming input limits")
993
994 call get_current(wct,ixi^l,ixa^l,idirmin,current)
995 tmpvec(ixa^s,1:ndir)=zero
996 do jdir=idirmin,3
997 tmpvec(ixa^s,jdir)=current(ixa^s,jdir)
998 end do
999
1000 ixa^l=ixo^l^ladd2;
1001 call curlvector(tmpvec,ixi^l,ixa^l,tmpvec2,idirmin1,1,3)
1002
1003 ixa^l=ixo^l^ladd1;
1004 tmpvec(ixa^s,1:ndir)=zero
1005 call curlvector(tmpvec2,ixi^l,ixa^l,tmpvec,idirmin1,1,3)
1006 ehyper(ixa^s,1:ndir) = - tmpvec(ixa^s,1:ndir)*mf_eta_hyper
1007
1008 ixa^l=ixo^l;
1009 tmpvec2(ixa^s,1:ndir)=zero
1010 call curlvector(ehyper,ixi^l,ixa^l,tmpvec2,idirmin1,1,3)
1011
1012 do idir=1,ndir
1013 w(ixo^s,mag(idir)) = w(ixo^s,mag(idir))-tmpvec2(ixo^s,idir)*qdt
1014 end do
1015
1016 end subroutine add_source_hyperres
1017
1018 subroutine add_source_glm(qdt,ixI^L,ixO^L,wCT,w,x)
1019 ! Add divB related sources to w within ixO
1020 ! corresponding to Dedner JCP 2002, 175, 645 _equation 24_
1021 ! giving the EGLM-MHD scheme
1023 use mod_geometry
1024
1025 integer, intent(in) :: ixi^l, ixo^l
1026 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1027 double precision, intent(inout) :: w(ixi^s,1:nw)
1028 double precision:: divb(ixi^s)
1029 double precision :: gradpsi(ixi^s)
1030 integer :: idim,idir
1031
1032 ! We calculate now div B
1033 call get_divb(wct,ixi^l,ixo^l,divb,mf_divb_nth)
1034
1035 ! dPsi/dt = - Ch^2/Cp^2 Psi
1036 if (mf_glm_alpha < zero) then
1037 w(ixo^s,psi_) = abs(mf_glm_alpha)*wct(ixo^s,psi_)
1038 else
1039 ! implicit update of Psi variable
1040 ! equation (27) in Mignone 2010 J. Com. Phys. 229, 2117
1041 if(slab_uniform) then
1042 w(ixo^s,psi_) = dexp(-qdt*cmax_global*mf_glm_alpha/minval(dxlevel(:)))*w(ixo^s,psi_)
1043 else
1044 w(ixo^s,psi_) = dexp(-qdt*cmax_global*mf_glm_alpha/minval(block%ds(ixo^s,:),dim=ndim+1))*w(ixo^s,psi_)
1045 end if
1046 end if
1047
1048 ! gradient of Psi
1049 do idim=1,ndim
1050 select case(typegrad)
1051 case("central")
1052 call gradient(wct(ixi^s,psi_),ixi^l,ixo^l,idim,gradpsi)
1053 case("limited")
1054 call gradientl(wct(ixi^s,psi_),ixi^l,ixo^l,idim,gradpsi)
1055 end select
1056 end do
1057 end subroutine add_source_glm
1058
1059 !> Add divB related sources to w within ixO corresponding to Powel
1060 subroutine add_source_powel(qdt,ixI^L,ixO^L,wCT,w,x)
1062
1063 integer, intent(in) :: ixi^l, ixo^l
1064 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1065 double precision, intent(inout) :: w(ixi^s,1:nw)
1066
1067 double precision :: divb(ixi^s)
1068 integer :: idir
1069
1070 ! We calculate now div B
1071 call get_divb(wct,ixi^l,ixo^l,divb,mf_divb_nth)
1072
1073 ! b = b - qdt v * div b
1074 do idir=1,ndir
1075 w(ixo^s,mag(idir))=w(ixo^s,mag(idir))-qdt*wct(ixo^s,mom(idir))*divb(ixo^s)
1076 end do
1077
1078 end subroutine add_source_powel
1079
1080 subroutine add_source_janhunen(qdt,ixI^L,ixO^L,wCT,w,x)
1081 ! Add divB related sources to w within ixO
1082 ! corresponding to Janhunen, just the term in the induction equation.
1084
1085 integer, intent(in) :: ixi^l, ixo^l
1086 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1087 double precision, intent(inout) :: w(ixi^s,1:nw)
1088 double precision :: divb(ixi^s)
1089 integer :: idir
1090
1091 ! We calculate now div B
1092 call get_divb(wct,ixi^l,ixo^l,divb,mf_divb_nth)
1093
1094 ! b = b - qdt v * div b
1095 do idir=1,ndir
1096 w(ixo^s,mag(idir))=w(ixo^s,mag(idir))-qdt*wct(ixo^s,mom(idir))*divb(ixo^s)
1097 end do
1098
1099 end subroutine add_source_janhunen
1100
1101 subroutine add_source_linde(qdt,ixI^L,ixO^L,wCT,w,x)
1102 ! Add Linde's divB related sources to wnew within ixO
1104 use mod_geometry
1105
1106 integer, intent(in) :: ixi^l, ixo^l
1107 double precision, intent(in) :: qdt, wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1108 double precision, intent(inout) :: w(ixi^s,1:nw)
1109 double precision :: divb(ixi^s),graddivb(ixi^s)
1110 integer :: idim, idir, ixp^l, i^d, iside
1111 logical, dimension(-1:1^D&) :: leveljump
1112
1113 ! Calculate div B
1114 ixp^l=ixo^l^ladd1;
1115 call get_divb(wct,ixi^l,ixp^l,divb,mf_divb_nth)
1116
1117 ! for AMR stability, retreat one cell layer from the boarders of level jump
1118 {do i^db=-1,1\}
1119 if(i^d==0|.and.) cycle
1120 if(neighbor_type(i^d,block%igrid)==2 .or. neighbor_type(i^d,block%igrid)==4) then
1121 leveljump(i^d)=.true.
1122 else
1123 leveljump(i^d)=.false.
1124 end if
1125 {end do\}
1126
1127 ixp^l=ixo^l;
1128 do idim=1,ndim
1129 select case(idim)
1130 {case(^d)
1131 do iside=1,2
1132 i^dd=kr(^dd,^d)*(2*iside-3);
1133 if (leveljump(i^dd)) then
1134 if (iside==1) then
1135 ixpmin^d=ixomin^d-i^d
1136 else
1137 ixpmax^d=ixomax^d-i^d
1138 end if
1139 end if
1140 end do
1141 \}
1142 end select
1143 end do
1144
1145 ! Add Linde's diffusive terms
1146 do idim=1,ndim
1147 ! Calculate grad_idim(divb)
1148 select case(typegrad)
1149 case("central")
1150 call gradient(divb,ixi^l,ixp^l,idim,graddivb)
1151 case("limited")
1152 call gradientl(divb,ixi^l,ixp^l,idim,graddivb)
1153 end select
1154
1155 ! Multiply by Linde's eta*dt = divbdiff*(c_max*dx)*dt = divbdiff*dx**2
1156 if (slab_uniform) then
1157 graddivb(ixp^s)=graddivb(ixp^s)*divbdiff/(^d&1.0d0/dxlevel(^d)**2+)
1158 else
1159 graddivb(ixp^s)=graddivb(ixp^s)*divbdiff &
1160 /(^d&1.0d0/block%ds(ixp^s,^d)**2+)
1161 end if
1162
1163 w(ixp^s,mag(idim))=w(ixp^s,mag(idim))+graddivb(ixp^s)
1164 end do
1165
1166 end subroutine add_source_linde
1167
1168
1169 !> get dimensionless div B = |divB| * volume / area / |B|
1170 subroutine get_normalized_divb(w,ixI^L,ixO^L,divb)
1171
1173
1174 integer, intent(in) :: ixi^l, ixo^l
1175 double precision, intent(in) :: w(ixi^s,1:nw)
1176 double precision :: divb(ixi^s), dsurface(ixi^s)
1177
1178 double precision :: invb(ixo^s)
1179 integer :: ixa^l,idims
1180
1181 call get_divb(w,ixi^l,ixo^l,divb)
1182 invb(ixo^s)=sqrt(sum(w(ixo^s, mag(:))**2, dim=ndim+1))
1183 where(invb(ixo^s)/=0.d0)
1184 invb(ixo^s)=1.d0/invb(ixo^s)
1185 end where
1186 if(slab_uniform) then
1187 divb(ixo^s)=0.5d0*abs(divb(ixo^s))*invb(ixo^s)/sum(1.d0/dxlevel(:))
1188 else
1189 ixamin^d=ixomin^d-1;
1190 ixamax^d=ixomax^d-1;
1191 dsurface(ixo^s)= sum(block%surfaceC(ixo^s,:),dim=ndim+1)
1192 do idims=1,ndim
1193 ixa^l=ixo^l-kr(idims,^d);
1194 dsurface(ixo^s)=dsurface(ixo^s)+block%surfaceC(ixa^s,idims)
1195 end do
1196 divb(ixo^s)=abs(divb(ixo^s))*invb(ixo^s)*&
1197 block%dvolume(ixo^s)/dsurface(ixo^s)
1198 end if
1199
1200 end subroutine get_normalized_divb
1201
1202 !> Calculate idirmin and the idirmin:3 components of the common current array
1203 !> make sure that dxlevel(^D) is set correctly.
1204 subroutine get_current(w,ixI^L,ixO^L,idirmin,current,fourthorder)
1206 use mod_geometry
1207
1208 integer, intent(in) :: ixo^l, ixi^l
1209 double precision, intent(in) :: w(ixi^s,1:nw)
1210 integer, intent(out) :: idirmin
1211 logical, intent(in), optional :: fourthorder
1212
1213 ! For ndir=2 only 3rd component of J can exist, ndir=1 is impossible for MHD
1214 double precision :: current(ixi^s,7-2*ndir:3)
1215 integer :: idir, idirmin0
1216
1217 idirmin0 = 7-2*ndir
1218
1219 if(present(fourthorder)) then
1220 call curlvector(w(ixi^s,mag(1:ndir)),ixi^l,ixo^l,current,idirmin,idirmin0,ndir,fourthorder)
1221 else
1222 call curlvector(w(ixi^s,mag(1:ndir)),ixi^l,ixo^l,current,idirmin,idirmin0,ndir)
1223 end if
1224
1225 end subroutine get_current
1226
1227 !> If resistivity is not zero, check diffusion time limit for dt
1228 subroutine mf_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
1230 use mod_usr_methods
1231
1232 integer, intent(in) :: ixi^l, ixo^l
1233 double precision, intent(inout) :: dtnew
1234 double precision, intent(in) :: dx^d
1235 double precision, intent(in) :: w(ixi^s,1:nw)
1236 double precision, intent(in) :: x(ixi^s,1:ndim)
1237
1238 double precision :: dxarr(ndim)
1239 double precision :: current(ixi^s,7-2*ndir:3),eta(ixi^s)
1240 integer :: idirmin,idim
1241
1242 dtnew = bigdouble
1243
1244 ^d&dxarr(^d)=dx^d;
1245 if (mf_eta>zero)then
1246 dtnew=dtdiffpar*minval(dxarr(1:ndim))**2/mf_eta
1247 else if (mf_eta<zero)then
1248 call get_current(w,ixi^l,ixo^l,idirmin,current)
1249 call usr_special_resistivity(w,ixi^l,ixo^l,idirmin,x,current,eta)
1250 dtnew=bigdouble
1251 do idim=1,ndim
1252 if(slab_uniform) then
1253 dtnew=min(dtnew,&
1254 dtdiffpar/(smalldouble+maxval(eta(ixo^s)/dxarr(idim)**2)))
1255 else
1256 dtnew=min(dtnew,&
1257 dtdiffpar/(smalldouble+maxval(eta(ixo^s)/block%ds(ixo^s,idim)**2)))
1258 end if
1259 end do
1260 end if
1261
1262 if(mf_eta_hyper>zero) then
1263 if(slab_uniform) then
1264 dtnew=min(dtdiffpar*minval(dxarr(1:ndim))**4/mf_eta_hyper,dtnew)
1265 else
1266 dtnew=min(dtdiffpar*minval(block%ds(ixo^s,1:ndim))**4/mf_eta_hyper,dtnew)
1267 end if
1268 end if
1269 end subroutine mf_get_dt
1270
1271 ! Add geometrical source terms to w
1272 subroutine mf_add_source_geom(qdt,dtfactor, ixI^L,ixO^L,wCT,wprim,w,x)
1274 use mod_geometry
1275
1276 integer, intent(in) :: ixi^l, ixo^l
1277 double precision, intent(in) :: qdt, dtfactor, x(ixi^s,1:ndim)
1278 double precision, intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw), w(ixi^s,1:nw)
1279
1280 double precision :: tmp,invr,cs2
1281 integer :: iw,idir,ix^d
1282
1283 integer :: mr_,mphi_ ! Polar var. names
1284 integer :: br_,bphi_
1285
1286 mr_=mom(1); mphi_=mom(1)-1+phi_ ! Polar var. names
1287 br_=mag(1); bphi_=mag(1)-1+phi_
1288
1289 select case (coordinate)
1290 case (cylindrical)
1291 if(phi_>0) then
1292 if(.not.stagger_grid) then
1293 {do ix^db=ixomin^db,ixomax^db\}
1294 w(ix^d,bphi_)=w(ix^d,bphi_)+qdt/x(ix^d,1)*&
1295 (wct(ix^d,bphi_)*wct(ix^d,mr_) &
1296 -wct(ix^d,br_)*wct(ix^d,mphi_))
1297 {end do\}
1298 end if
1299 end if
1300 if(mf_glm) w(ixo^s,br_)=w(ixo^s,br_)+qdt*wct(ixo^s,psi_)/x(ixo^s,1)
1301 case (spherical)
1302 {do ix^db=ixomin^db,ixomax^db\}
1303 invr=qdt/x(ix^d,1)
1304 ! b1
1305 if(mf_glm) then
1306 w(ix^d,mag(1))=w(ix^d,mag(1))+invr*2.0d0*wct(ix^d,psi_)
1307 end if
1308
1309 {^nooned
1310 ! b2
1311 if(.not.stagger_grid) then
1312 tmp=wct(ix^d,mom(1))*wct(ix^d,mag(2))-wct(ix^d,mom(2))*wct(ix^d,mag(1))
1313 cs2=1.d0/tan(x(ix^d,2))
1314 if(mf_glm) then
1315 tmp=tmp+cs2*wct(ix^d,psi_)
1316 end if
1317 w(ix^d,mag(2))=w(ix^d,mag(2))+tmp*invr
1318 end if
1319 }
1320
1321 if(ndir==3) then
1322 {^ifoned
1323 ! b3
1324 if(.not.stagger_grid) then
1325 w(ix^d,mag(3))=w(ix^d,mag(3))+invr*&
1326 (wct(ix^d,mom(1))*wct(ix^d,mag(3)) &
1327 -wct(ix^d,mom(3))*wct(ix^d,mag(1)))
1328 end if
1329 }
1330 {^nooned
1331 ! b3
1332 if(.not.stagger_grid) then
1333 w(ix^d,mag(3))=w(ix^d,mag(3))+invr*&
1334 (wct(ix^d,mom(1))*wct(ix^d,mag(3)) &
1335 -wct(ix^d,mom(3))*wct(ix^d,mag(1)) &
1336 -(wct(ix^d,mom(3))*wct(ix^d,mag(2)) &
1337 -wct(ix^d,mom(2))*wct(ix^d,mag(3)))*cs2)
1338 end if
1339 }
1340 end if
1341 {end do\}
1342 end select
1343
1344 end subroutine mf_add_source_geom
1345
1346 subroutine mf_modify_wlr(ixI^L,ixO^L,qt,wLC,wRC,wLp,wRp,s,idir)
1348 use mod_usr_methods
1349 integer, intent(in) :: ixi^l, ixo^l, idir
1350 double precision, intent(in) :: qt
1351 double precision, intent(inout) :: wlc(ixi^s,1:nw), wrc(ixi^s,1:nw)
1352 double precision, intent(inout) :: wlp(ixi^s,1:nw), wrp(ixi^s,1:nw)
1353 type(state) :: s
1354 double precision :: db(ixi^s), dpsi(ixi^s)
1355
1356 ! Solve the Riemann problem for the linear 2x2 system for normal
1357 ! B-field and GLM_Psi according to Dedner 2002:
1358 ! This implements eq. (42) in Dedner et al. 2002 JcP 175
1359 ! Gives the Riemann solution on the interface
1360 ! for the normal B component and Psi in the GLM-MHD system.
1361 ! 23/04/2013 Oliver Porth
1362 db(ixo^s) = wrp(ixo^s,mag(idir)) - wlp(ixo^s,mag(idir))
1363 dpsi(ixo^s) = wrp(ixo^s,psi_) - wlp(ixo^s,psi_)
1364
1365 wlp(ixo^s,mag(idir)) = 0.5d0 * (wrp(ixo^s,mag(idir)) + wlp(ixo^s,mag(idir))) &
1366 - 0.5d0/cmax_global * dpsi(ixo^s)
1367 wlp(ixo^s,psi_) = 0.5d0 * (wrp(ixo^s,psi_) + wlp(ixo^s,psi_)) &
1368 - 0.5d0*cmax_global * db(ixo^s)
1369
1370 wrp(ixo^s,mag(idir)) = wlp(ixo^s,mag(idir))
1371 wrp(ixo^s,psi_) = wlp(ixo^s,psi_)
1372
1373 end subroutine mf_modify_wlr
1374
1375 subroutine mf_boundary_adjust(igrid,psb)
1377 integer, intent(in) :: igrid
1378 type(state), target :: psb(max_blocks)
1379
1380 integer :: ib, idims, iside, ixo^l, i^d
1381
1382 block=>ps(igrid)
1383 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
1384 do idims=1,ndim
1385 ! to avoid using as yet unknown corner info in more than 1D, we
1386 ! fill only interior mesh ranges of the ghost cell ranges at first,
1387 ! and progressively enlarge the ranges to include corners later
1388 do iside=1,2
1389 i^d=kr(^d,idims)*(2*iside-3);
1390 if (neighbor_type(i^d,igrid)/=1) cycle
1391 ib=(idims-1)*2+iside
1392 if(.not.boundary_divbfix(ib)) cycle
1393 if(any(typeboundary(:,ib)==bc_special)) then
1394 ! MF nonlinear force-free B field extrapolation and data driven
1395 ! require normal B of the first ghost cell layer to be untouched by
1396 ! fixdivB=0 process, set boundary_divbfix_skip(iB)=1 in par file
1397 select case (idims)
1398 {case (^d)
1399 if (iside==2) then
1400 ! maximal boundary
1401 ixomin^dd=ixghi^d+1-nghostcells+boundary_divbfix_skip(2*^d)^d%ixOmin^dd=ixglo^dd;
1402 ixomax^dd=ixghi^dd;
1403 else
1404 ! minimal boundary
1405 ixomin^dd=ixglo^dd;
1406 ixomax^dd=ixglo^d-1+nghostcells-boundary_divbfix_skip(2*^d-1)^d%ixOmax^dd=ixghi^dd;
1407 end if \}
1408 end select
1409 call fixdivb_boundary(ixg^ll,ixo^l,psb(igrid)%w,psb(igrid)%x,ib)
1410 end if
1411 end do
1412 end do
1413
1414 end subroutine mf_boundary_adjust
1415
1416 subroutine fixdivb_boundary(ixG^L,ixO^L,w,x,iB)
1418
1419 integer, intent(in) :: ixg^l,ixo^l,ib
1420 double precision, intent(inout) :: w(ixg^s,1:nw)
1421 double precision, intent(in) :: x(ixg^s,1:ndim)
1422
1423 double precision :: dx1x2,dx1x3,dx2x1,dx2x3,dx3x1,dx3x2
1424 integer :: ix^d,ixf^l
1425
1426 select case(ib)
1427 case(1)
1428 ! 2nd order CD for divB=0 to set normal B component better
1429 {^iftwod
1430 ixfmin1=ixomin1+1
1431 ixfmax1=ixomax1+1
1432 ixfmin2=ixomin2+1
1433 ixfmax2=ixomax2-1
1434 if(slab_uniform) then
1435 dx1x2=dxlevel(1)/dxlevel(2)
1436 do ix1=ixfmax1,ixfmin1,-1
1437 w(ix1-1,ixfmin2:ixfmax2,mag(1))=w(ix1+1,ixfmin2:ixfmax2,mag(1)) &
1438 +dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))-&
1439 w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))
1440 enddo
1441 else
1442 do ix1=ixfmax1,ixfmin1,-1
1443 w(ix1-1,ixfmin2:ixfmax2,mag(1))=( (w(ix1+1,ixfmin2:ixfmax2,mag(1))+&
1444 w(ix1,ixfmin2:ixfmax2,mag(1)))*block%surfaceC(ix1,ixfmin2:ixfmax2,1)&
1445 +(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))+w(ix1,ixfmin2:ixfmax2,mag(2)))*&
1446 block%surfaceC(ix1,ixfmin2:ixfmax2,2)&
1447 -(w(ix1,ixfmin2:ixfmax2,mag(2))+w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))*&
1448 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,2) )&
1449 /block%surfaceC(ix1-1,ixfmin2:ixfmax2,1)-w(ix1,ixfmin2:ixfmax2,mag(1))
1450 end do
1451 end if
1452 }
1453 {^ifthreed
1454 ixfmin1=ixomin1+1
1455 ixfmax1=ixomax1+1
1456 ixfmin2=ixomin2+1
1457 ixfmax2=ixomax2-1
1458 ixfmin3=ixomin3+1
1459 ixfmax3=ixomax3-1
1460 if(slab_uniform) then
1461 dx1x2=dxlevel(1)/dxlevel(2)
1462 dx1x3=dxlevel(1)/dxlevel(3)
1463 do ix1=ixfmax1,ixfmin1,-1
1464 w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
1465 w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)) &
1466 +dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))-&
1467 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2))) &
1468 +dx1x3*(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))-&
1469 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))
1470 end do
1471 else
1472 do ix1=ixfmax1,ixfmin1,-1
1473 w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
1474 ( (w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))+&
1475 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)))*&
1476 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)&
1477 +(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))+&
1478 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2)))*&
1479 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,2)&
1480 -(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2))+&
1481 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2)))*&
1482 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,2)&
1483 +(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))+&
1484 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3)))*&
1485 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,3)&
1486 -(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3))+&
1487 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))*&
1488 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,3) )&
1489 /block%surfaceC(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)-&
1490 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))
1491 end do
1492 end if
1493 }
1494 case(2)
1495 {^iftwod
1496 ixfmin1=ixomin1-1
1497 ixfmax1=ixomax1-1
1498 ixfmin2=ixomin2+1
1499 ixfmax2=ixomax2-1
1500 if(slab_uniform) then
1501 dx1x2=dxlevel(1)/dxlevel(2)
1502 do ix1=ixfmin1,ixfmax1
1503 w(ix1+1,ixfmin2:ixfmax2,mag(1))=w(ix1-1,ixfmin2:ixfmax2,mag(1)) &
1504 -dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))-&
1505 w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))
1506 enddo
1507 else
1508 do ix1=ixfmin1,ixfmax1
1509 w(ix1+1,ixfmin2:ixfmax2,mag(1))=( (w(ix1-1,ixfmin2:ixfmax2,mag(1))+&
1510 w(ix1,ixfmin2:ixfmax2,mag(1)))*block%surfaceC(ix1-1,ixfmin2:ixfmax2,1)&
1511 -(w(ix1,ixfmin2+1:ixfmax2+1,mag(2))+w(ix1,ixfmin2:ixfmax2,mag(2)))*&
1512 block%surfaceC(ix1,ixfmin2:ixfmax2,2)&
1513 +(w(ix1,ixfmin2:ixfmax2,mag(2))+w(ix1,ixfmin2-1:ixfmax2-1,mag(2)))*&
1514 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,2) )&
1515 /block%surfaceC(ix1,ixfmin2:ixfmax2,1)-w(ix1,ixfmin2:ixfmax2,mag(1))
1516 end do
1517 end if
1518 }
1519 {^ifthreed
1520 ixfmin1=ixomin1-1
1521 ixfmax1=ixomax1-1
1522 ixfmin2=ixomin2+1
1523 ixfmax2=ixomax2-1
1524 ixfmin3=ixomin3+1
1525 ixfmax3=ixomax3-1
1526 if(slab_uniform) then
1527 dx1x2=dxlevel(1)/dxlevel(2)
1528 dx1x3=dxlevel(1)/dxlevel(3)
1529 do ix1=ixfmin1,ixfmax1
1530 w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
1531 w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)) &
1532 -dx1x2*(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))-&
1533 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2))) &
1534 -dx1x3*(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))-&
1535 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))
1536 end do
1537 else
1538 do ix1=ixfmin1,ixfmax1
1539 w(ix1+1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))=&
1540 ( (w(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))+&
1541 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1)))*&
1542 block%surfaceC(ix1-1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)&
1543 -(w(ix1,ixfmin2+1:ixfmax2+1,ixfmin3:ixfmax3,mag(2))+&
1544 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2)))*&
1545 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,2)&
1546 +(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(2))+&
1547 w(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,mag(2)))*&
1548 block%surfaceC(ix1,ixfmin2-1:ixfmax2-1,ixfmin3:ixfmax3,2)&
1549 -(w(ix1,ixfmin2:ixfmax2,ixfmin3+1:ixfmax3+1,mag(3))+&
1550 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3)))*&
1551 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,3)&
1552 +(w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(3))+&
1553 w(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,mag(3)))*&
1554 block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3-1:ixfmax3-1,3) )&
1555 /block%surfaceC(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,1)-&
1556 w(ix1,ixfmin2:ixfmax2,ixfmin3:ixfmax3,mag(1))
1557 end do
1558 end if
1559 }
1560 case(3)
1561 {^iftwod
1562 ixfmin1=ixomin1+1
1563 ixfmax1=ixomax1-1
1564 ixfmin2=ixomin2+1
1565 ixfmax2=ixomax2+1
1566 if(slab_uniform) then
1567 dx2x1=dxlevel(2)/dxlevel(1)
1568 do ix2=ixfmax2,ixfmin2,-1
1569 w(ixfmin1:ixfmax1,ix2-1,mag(2))=w(ixfmin1:ixfmax1,ix2+1,mag(2)) &
1570 +dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))-&
1571 w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))
1572 enddo
1573 else
1574 do ix2=ixfmax2,ixfmin2,-1
1575 w(ixfmin1:ixfmax1,ix2-1,mag(2))=( (w(ixfmin1:ixfmax1,ix2+1,mag(2))+&
1576 w(ixfmin1:ixfmax1,ix2,mag(2)))*block%surfaceC(ixfmin1:ixfmax1,ix2,2)&
1577 +(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))+w(ixfmin1:ixfmax1,ix2,mag(1)))*&
1578 block%surfaceC(ixfmin1:ixfmax1,ix2,1)&
1579 -(w(ixfmin1:ixfmax1,ix2,mag(1))+w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))*&
1580 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,1) )&
1581 /block%surfaceC(ixfmin1:ixfmax1,ix2-1,2)-w(ixfmin1:ixfmax1,ix2,mag(2))
1582 end do
1583 end if
1584 }
1585 {^ifthreed
1586 ixfmin1=ixomin1+1
1587 ixfmax1=ixomax1-1
1588 ixfmin3=ixomin3+1
1589 ixfmax3=ixomax3-1
1590 ixfmin2=ixomin2+1
1591 ixfmax2=ixomax2+1
1592 if(slab_uniform) then
1593 dx2x1=dxlevel(2)/dxlevel(1)
1594 dx2x3=dxlevel(2)/dxlevel(3)
1595 do ix2=ixfmax2,ixfmin2,-1
1596 w(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,mag(2))=w(ixfmin1:ixfmax1,&
1597 ix2+1,ixfmin3:ixfmax3,mag(2)) &
1598 +dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))-&
1599 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1))) &
1600 +dx2x3*(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))-&
1601 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))
1602 end do
1603 else
1604 do ix2=ixfmax2,ixfmin2,-1
1605 w(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,mag(2))=&
1606 ( (w(ixfmin1:ixfmax1,ix2+1,ixfmin3:ixfmax3,mag(2))+&
1607 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2)))*&
1608 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,2)&
1609 +(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))+&
1610 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1)))*&
1611 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,1)&
1612 -(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1))+&
1613 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1)))*&
1614 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,1)&
1615 +(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))+&
1616 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3)))*&
1617 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,3)&
1618 -(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3))+&
1619 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))*&
1620 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,3) )&
1621 /block%surfaceC(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,2)-&
1622 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2))
1623 end do
1624 end if
1625 }
1626 case(4)
1627 {^iftwod
1628 ixfmin1=ixomin1+1
1629 ixfmax1=ixomax1-1
1630 ixfmin2=ixomin2-1
1631 ixfmax2=ixomax2-1
1632 if(slab_uniform) then
1633 dx2x1=dxlevel(2)/dxlevel(1)
1634 do ix2=ixfmin2,ixfmax2
1635 w(ixfmin1:ixfmax1,ix2+1,mag(2))=w(ixfmin1:ixfmax1,ix2-1,mag(2)) &
1636 -dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))-&
1637 w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))
1638 end do
1639 else
1640 do ix2=ixfmin2,ixfmax2
1641 w(ixfmin1:ixfmax1,ix2+1,mag(2))=( (w(ixfmin1:ixfmax1,ix2-1,mag(2))+&
1642 w(ixfmin1:ixfmax1,ix2,mag(2)))*block%surfaceC(ixfmin1:ixfmax1,ix2-1,2)&
1643 -(w(ixfmin1+1:ixfmax1+1,ix2,mag(1))+w(ixfmin1:ixfmax1,ix2,mag(1)))*&
1644 block%surfaceC(ixfmin1:ixfmax1,ix2,1)&
1645 +(w(ixfmin1:ixfmax1,ix2,mag(1))+w(ixfmin1-1:ixfmax1-1,ix2,mag(1)))*&
1646 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,1) )&
1647 /block%surfaceC(ixfmin1:ixfmax1,ix2,2)-w(ixfmin1:ixfmax1,ix2,mag(2))
1648 end do
1649 end if
1650 }
1651 {^ifthreed
1652 ixfmin1=ixomin1+1
1653 ixfmax1=ixomax1-1
1654 ixfmin3=ixomin3+1
1655 ixfmax3=ixomax3-1
1656 ixfmin2=ixomin2-1
1657 ixfmax2=ixomax2-1
1658 if(slab_uniform) then
1659 dx2x1=dxlevel(2)/dxlevel(1)
1660 dx2x3=dxlevel(2)/dxlevel(3)
1661 do ix2=ixfmin2,ixfmax2
1662 w(ixfmin1:ixfmax1,ix2+1,ixfmin3:ixfmax3,mag(2))=w(ixfmin1:ixfmax1,&
1663 ix2-1,ixfmin3:ixfmax3,mag(2)) &
1664 -dx2x1*(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))-&
1665 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1))) &
1666 -dx2x3*(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))-&
1667 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))
1668 end do
1669 else
1670 do ix2=ixfmin2,ixfmax2
1671 w(ixfmin1:ixfmax1,ix2+1,ixfmin3:ixfmax3,mag(2))=&
1672 ( (w(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,mag(2))+&
1673 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2)))*&
1674 block%surfaceC(ixfmin1:ixfmax1,ix2-1,ixfmin3:ixfmax3,2)&
1675 -(w(ixfmin1+1:ixfmax1+1,ix2,ixfmin3:ixfmax3,mag(1))+&
1676 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1)))*&
1677 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,1)&
1678 +(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(1))+&
1679 w(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,mag(1)))*&
1680 block%surfaceC(ixfmin1-1:ixfmax1-1,ix2,ixfmin3:ixfmax3,1)&
1681 -(w(ixfmin1:ixfmax1,ix2,ixfmin3+1:ixfmax3+1,mag(3))+&
1682 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3)))*&
1683 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,3)&
1684 +(w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(3))+&
1685 w(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,mag(3)))*&
1686 block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3-1:ixfmax3-1,3) )&
1687 /block%surfaceC(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,2)-&
1688 w(ixfmin1:ixfmax1,ix2,ixfmin3:ixfmax3,mag(2))
1689 end do
1690 end if
1691 }
1692 {^ifthreed
1693 case(5)
1694 ixfmin1=ixomin1+1
1695 ixfmax1=ixomax1-1
1696 ixfmin2=ixomin2+1
1697 ixfmax2=ixomax2-1
1698 ixfmin3=ixomin3+1
1699 ixfmax3=ixomax3+1
1700 if(slab_uniform) then
1701 dx3x1=dxlevel(3)/dxlevel(1)
1702 dx3x2=dxlevel(3)/dxlevel(2)
1703 do ix3=ixfmax3,ixfmin3,-1
1704 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,mag(3))=w(ixfmin1:ixfmax1,&
1705 ixfmin2:ixfmax2,ix3+1,mag(3)) &
1706 +dx3x1*(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))-&
1707 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1))) &
1708 +dx3x2*(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))-&
1709 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))
1710 end do
1711 else
1712 do ix3=ixfmax3,ixfmin3,-1
1713 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,mag(3))=&
1714 ( (w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3+1,mag(3))+&
1715 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3)))*&
1716 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,3)&
1717 +(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))+&
1718 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1)))*&
1719 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,1)&
1720 -(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1))+&
1721 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1)))*&
1722 block%surfaceC(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,1)&
1723 +(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))+&
1724 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2)))*&
1725 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,2)&
1726 -(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2))+&
1727 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))*&
1728 block%surfaceC(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,2) )&
1729 /block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,3)-&
1730 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3))
1731 end do
1732 end if
1733 case(6)
1734 ixfmin1=ixomin1+1
1735 ixfmax1=ixomax1-1
1736 ixfmin2=ixomin2+1
1737 ixfmax2=ixomax2-1
1738 ixfmin3=ixomin3-1
1739 ixfmax3=ixomax3-1
1740 if(slab_uniform) then
1741 dx3x1=dxlevel(3)/dxlevel(1)
1742 dx3x2=dxlevel(3)/dxlevel(2)
1743 do ix3=ixfmin3,ixfmax3
1744 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3+1,mag(3))=w(ixfmin1:ixfmax1,&
1745 ixfmin2:ixfmax2,ix3-1,mag(3)) &
1746 -dx3x1*(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))-&
1747 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1))) &
1748 -dx3x2*(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))-&
1749 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))
1750 end do
1751 else
1752 do ix3=ixfmin3,ixfmax3
1753 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3+1,mag(3))=&
1754 ( (w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,mag(3))+&
1755 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3)))*&
1756 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3-1,3)&
1757 -(w(ixfmin1+1:ixfmax1+1,ixfmin2:ixfmax2,ix3,mag(1))+&
1758 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1)))*&
1759 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,1)&
1760 +(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(1))+&
1761 w(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,mag(1)))*&
1762 block%surfaceC(ixfmin1-1:ixfmax1-1,ixfmin2:ixfmax2,ix3,1)&
1763 -(w(ixfmin1:ixfmax1,ixfmin2+1:ixfmax2+1,ix3,mag(2))+&
1764 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2)))*&
1765 block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,2)&
1766 +(w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(2))+&
1767 w(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,mag(2)))*&
1768 block%surfaceC(ixfmin1:ixfmax1,ixfmin2-1:ixfmax2-1,ix3,2) )&
1769 /block%surfaceC(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,3)-&
1770 w(ixfmin1:ixfmax1,ixfmin2:ixfmax2,ix3,mag(3))
1771 end do
1772 end if
1773 }
1774 case default
1775 call mpistop("Special boundary is not defined for this region")
1776 end select
1777
1778 end subroutine fixdivb_boundary
1779
1780 {^nooned
1781 subroutine mf_clean_divb_multigrid(qdt, qt, active)
1782 use mod_forest
1785 use mod_geometry
1786
1787 double precision, intent(in) :: qdt !< Current time step
1788 double precision, intent(in) :: qt !< Current time
1789 logical, intent(inout) :: active !< Output if the source is active
1790
1791 double precision :: tmp(ixg^t), grad(ixg^t, ndim)
1792 double precision :: res
1793 double precision, parameter :: max_residual = 1d-3
1794 double precision, parameter :: residual_reduction = 1d-10
1795 integer :: id
1796 integer, parameter :: max_its = 50
1797 double precision :: residual_it(max_its), max_divb
1798 integer :: iigrid, igrid
1799 integer :: n, nc, lvl, ix^l, ixc^l, idim
1800 type(tree_node), pointer :: pnode
1801
1802 mg%operator_type = mg_laplacian
1803
1804 ! Set boundary conditions
1805 do n = 1, 2*ndim
1806 idim = (n+1)/2
1807 select case (typeboundary(mag(idim), n))
1808 case (bc_symm)
1809 ! d/dx B = 0, take phi = 0
1810 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
1811 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
1812 case (bc_asymm)
1813 ! B = 0, so grad(phi) = 0
1814 mg%bc(n, mg_iphi)%bc_type = mg_bc_neumann
1815 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
1816 case (bc_cont)
1817 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
1818 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
1819 case (bc_special)
1820 ! Assume Dirichlet boundary conditions, derivative zero
1821 mg%bc(n, mg_iphi)%bc_type = mg_bc_neumann
1822 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
1823 case (bc_periodic)
1824 ! Nothing to do here
1825 case default
1826 write(*,*) "mf_clean_divb_multigrid warning: unknown boundary type"
1827 mg%bc(n, mg_iphi)%bc_type = mg_bc_dirichlet
1828 mg%bc(n, mg_iphi)%bc_value = 0.0_dp
1829 end select
1830 end do
1831
1832 ix^l=ixm^ll^ladd1;
1833 max_divb = 0.0d0
1834
1835 ! Store divergence of B as right-hand side
1836 do iigrid = 1, igridstail
1837 igrid = igrids(iigrid);
1838 pnode => igrid_to_node(igrid, mype)%node
1839 id = pnode%id
1840 lvl = mg%boxes(id)%lvl
1841 nc = mg%box_size_lvl(lvl)
1842
1843 ! Geometry subroutines expect this to be set
1844 block => ps(igrid)
1845 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
1846
1847 call get_divb(ps(igrid)%w(ixg^t, 1:nw), ixg^ll, ixm^ll, tmp, &
1849 mg%boxes(id)%cc({1:nc}, mg_irhs) = tmp(ixm^t)
1850 max_divb = max(max_divb, maxval(abs(tmp(ixm^t))))
1851 end do
1852
1853 ! Solve laplacian(phi) = divB
1854 if(stagger_grid) then
1855 call mpi_allreduce(mpi_in_place, max_divb, 1, mpi_double_precision, &
1856 mpi_max, icomm, ierrmpi)
1857
1858 if (mype == 0) print *, "Performing multigrid divB cleaning"
1859 if (mype == 0) print *, "iteration vs residual"
1860 ! Solve laplacian(phi) = divB
1861 do n = 1, max_its
1862 call mg_fas_fmg(mg, n>1, max_res=residual_it(n))
1863 if (mype == 0) write(*, "(I4,E11.3)") n, residual_it(n)
1864 if (residual_it(n) < residual_reduction * max_divb) exit
1865 end do
1866 if (mype == 0 .and. n > max_its) then
1867 print *, "divb_multigrid warning: not fully converged"
1868 print *, "current amplitude of divb: ", residual_it(max_its)
1869 print *, "multigrid smallest grid: ", &
1870 mg%domain_size_lvl(:, mg%lowest_lvl)
1871 print *, "note: smallest grid ideally has <= 8 cells"
1872 print *, "multigrid dx/dy/dz ratio: ", mg%dr(:, 1)/mg%dr(1, 1)
1873 print *, "note: dx/dy/dz should be similar"
1874 end if
1875 else
1876 do n = 1, max_its
1877 call mg_fas_vcycle(mg, max_res=res)
1878 if (res < max_residual) exit
1879 end do
1880 if (res > max_residual) call mpistop("divb_multigrid: no convergence")
1881 end if
1882
1883
1884 ! Correct the magnetic field
1885 do iigrid = 1, igridstail
1886 igrid = igrids(iigrid);
1887 pnode => igrid_to_node(igrid, mype)%node
1888 id = pnode%id
1889
1890 ! Geometry subroutines expect this to be set
1891 block => ps(igrid)
1892 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
1893
1894 ! Compute the gradient of phi
1895 tmp(ix^s) = mg%boxes(id)%cc({:,}, mg_iphi)
1896
1897 if(stagger_grid) then
1898 do idim =1, ndim
1899 ixcmin^d=ixmlo^d-kr(idim,^d);
1900 ixcmax^d=ixmhi^d;
1901 call gradientf(tmp,ps(igrid)%x,ixg^ll,ixc^l,idim,grad(ixg^t,idim))
1902 ps(igrid)%ws(ixc^s,idim)=ps(igrid)%ws(ixc^s,idim)-grad(ixc^s,idim)
1903 end do
1904 call mf_face_to_center(ixm^ll,ps(igrid))
1905 else
1906 do idim = 1, ndim
1907 call gradient(tmp,ixg^ll,ixm^ll,idim,grad(ixg^t, idim))
1908 end do
1909 ! Apply the correction B* = B - gradient(phi)
1910 ps(igrid)%w(ixm^t, mag(1:ndim)) = &
1911 ps(igrid)%w(ixm^t, mag(1:ndim)) - grad(ixm^t, :)
1912 end if
1913
1914 end do
1915
1916 active = .true.
1917
1918 end subroutine mf_clean_divb_multigrid
1919 }
1920
1921 !> get electric field though averaging neighors to update faces in CT
1922 subroutine mf_update_faces_average(ixI^L,ixO^L,qt,qdt,wp,fC,fE,sCT,s,vcts)
1924 use mod_usr_methods
1925
1926 integer, intent(in) :: ixi^l, ixo^l
1927 double precision, intent(in) :: qt, qdt
1928 ! cell-center primitive variables
1929 double precision, intent(in) :: wp(ixi^s,1:nw)
1930 type(state) :: sct, s
1931 type(ct_velocity) :: vcts
1932 double precision, intent(in) :: fc(ixi^s,1:nwflux,1:ndim)
1933 double precision, intent(inout) :: fe(ixi^s,sdim:3)
1934
1935 double precision :: circ(ixi^s,1:ndim)
1936 double precision :: e_resi(ixi^s,sdim:3)
1937 integer :: ix^d,ixc^l,ixa^l,i1kr^d,i2kr^d
1938 integer :: idim1,idim2,idir,iwdim1,iwdim2
1939
1940 associate(bfaces=>s%ws,x=>s%x)
1941
1942 ! Calculate contribution to FEM of each edge,
1943 ! that is, estimate value of line integral of
1944 ! electric field in the positive idir direction.
1945
1946 ! if there is resistivity, get eta J
1947 if(mf_eta/=zero) call get_resistive_electric_field(ixi^l,ixo^l,sct,s,e_resi)
1948
1949 do idim1=1,ndim
1950 iwdim1 = mag(idim1)
1951 i1kr^d=kr(idim1,^d);
1952 do idim2=1,ndim
1953 iwdim2 = mag(idim2)
1954 i2kr^d=kr(idim2,^d);
1955 do idir=sdim,3! Direction of line integral
1956 ! Allow only even permutations
1957 if (lvc(idim1,idim2,idir)==1) then
1958 ixcmax^d=ixomax^d;
1959 ixcmin^d=ixomin^d+kr(idir,^d)-1;
1960 ! average cell-face electric field to cell edges
1961 {do ix^db=ixcmin^db,ixcmax^db\}
1962 fe(ix^d,idir)=quarter*&
1963 (fc(ix^d,iwdim1,idim2)+fc({ix^d+i1kr^d},iwdim1,idim2)&
1964 -fc(ix^d,iwdim2,idim1)-fc({ix^d+i2kr^d},iwdim2,idim1))
1965 ! add resistive electric field at cell edges E=-vxB+eta J
1966 if(mf_eta/=zero) fe(ix^d,idir)=fe(ix^d,idir)+e_resi(ix^d,idir)
1967 if(record_electric_field) s%we(ix^d,idir)=fe(ix^d,idir)
1968 ! times time step and edge length
1969 fe(ix^d,idir)=fe(ix^d,idir)*qdt*s%dsC(ix^d,idir)
1970 {end do\}
1971 end if
1972 end do
1973 end do
1974 end do
1975
1976 ! allow user to change inductive electric field, especially for boundary driven applications
1977 if(associated(usr_set_electric_field)) &
1978 call usr_set_electric_field(ixi^l,ixo^l,qt,qdt,fe,sct)
1979
1980 circ(ixi^s,1:ndim)=zero
1981
1982 ! Calculate circulation on each face
1983
1984 do idim1=1,ndim ! Coordinate perpendicular to face
1985 ixcmax^d=ixomax^d;
1986 ixcmin^d=ixomin^d-kr(idim1,^d);
1987 do idim2=1,ndim
1988 ixa^l=ixc^l-kr(idim2,^d);
1989 do idir=sdim,3 ! Direction of line integral
1990 ! Assemble indices
1991 if(lvc(idim1,idim2,idir)==1) then
1992 ! Add line integrals in direction idir
1993 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
1994 +(fe(ixc^s,idir)&
1995 -fe(ixa^s,idir))
1996 else if(lvc(idim1,idim2,idir)==-1) then
1997 ! Add line integrals in direction idir
1998 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
1999 -(fe(ixc^s,idir)&
2000 -fe(ixa^s,idir))
2001 end if
2002 end do
2003 end do
2004 ! Divide by the area of the face to get dB/dt
2005 circ(ixc^s,idim1)=circ(ixc^s,idim1)/s%surfaceC(ixc^s,idim1)
2006 ! Time update
2007 bfaces(ixc^s,idim1)=bfaces(ixc^s,idim1)-circ(ixc^s,idim1)
2008 end do
2009
2010 end associate
2011
2012 end subroutine mf_update_faces_average
2013
2014 !> update faces using UCT contact mode by Gardiner and Stone 2005 JCP 205, 509
2015 subroutine mf_update_faces_contact(ixI^L,ixO^L,qt,qdt,wp,fC,fE,sCT,s,vcts)
2017 use mod_usr_methods
2018
2019 integer, intent(in) :: ixi^l, ixo^l
2020 double precision, intent(in) :: qt, qdt
2021 ! cell-center primitive variables
2022 double precision, intent(in) :: wp(ixi^s,1:nw)
2023 type(state) :: sct, s
2024 type(ct_velocity) :: vcts
2025 double precision, intent(in) :: fc(ixi^s,1:nwflux,1:ndim)
2026 double precision, intent(inout) :: fe(ixi^s,sdim:3)
2027
2028 double precision :: circ(ixi^s,1:ndim)
2029 ! electric field at cell centers
2030 double precision :: ecc(ixi^s,sdim:3)
2031 ! gradient of E at left and right side of a cell face
2032 double precision :: el(ixi^s),er(ixi^s)
2033 ! gradient of E at left and right side of a cell corner
2034 double precision :: elc(ixi^s),erc(ixi^s)
2035 ! current on cell edges
2036 double precision :: jce(ixi^s,sdim:3)
2037 integer :: hxc^l,ixc^l,jxc^l,ixa^l,ixb^l
2038 integer :: idim1,idim2,idir,iwdim1,iwdim2
2039
2040 associate(bfaces=>s%ws,x=>s%x,w=>s%w,vnorm=>vcts%vnorm)
2041
2042 ecc=0.d0
2043 ! Calculate electric field at cell centers
2044 do idim1=1,ndim; do idim2=1,ndim; do idir=sdim,3
2045 if(lvc(idim1,idim2,idir)==1)then
2046 ecc(ixi^s,idir)=ecc(ixi^s,idir)+wp(ixi^s,mag(idim1))*wp(ixi^s,mom(idim2))
2047 else if(lvc(idim1,idim2,idir)==-1) then
2048 ecc(ixi^s,idir)=ecc(ixi^s,idir)-wp(ixi^s,mag(idim1))*wp(ixi^s,mom(idim2))
2049 endif
2050 enddo; enddo; enddo
2051
2052 ! if there is resistivity, get eta J
2053 if(mf_eta/=zero) call get_resistive_electric_field(ixi^l,ixo^l,sct,s,jce)
2054
2055 ! Calculate contribution to FEM of each edge,
2056 ! that is, estimate value of line integral of
2057 ! electric field in the positive idir direction.
2058 ! evaluate electric field along cell edges according to equation (41)
2059 do idim1=1,ndim
2060 iwdim1 = mag(idim1)
2061 do idim2=1,ndim
2062 iwdim2 = mag(idim2)
2063 do idir=sdim,3 ! Direction of line integral
2064 ! Allow only even permutations
2065 if (lvc(idim1,idim2,idir)==1) then
2066 ixcmax^d=ixomax^d;
2067 ixcmin^d=ixomin^d+kr(idir,^d)-1;
2068 ! Assemble indices
2069 jxc^l=ixc^l+kr(idim1,^d);
2070 hxc^l=ixc^l+kr(idim2,^d);
2071 ! average cell-face electric field to cell edges
2072 fe(ixc^s,idir)=quarter*&
2073 (fc(ixc^s,iwdim1,idim2)+fc(jxc^s,iwdim1,idim2)&
2074 -fc(ixc^s,iwdim2,idim1)-fc(hxc^s,iwdim2,idim1))
2075
2076 ! add slope in idim2 direction from equation (50)
2077 ixamin^d=ixcmin^d;
2078 ixamax^d=ixcmax^d+kr(idim1,^d);
2079 el(ixa^s)=fc(ixa^s,iwdim1,idim2)-ecc(ixa^s,idir)
2080 hxc^l=ixa^l+kr(idim2,^d);
2081 er(ixa^s)=fc(ixa^s,iwdim1,idim2)-ecc(hxc^s,idir)
2082 where(vnorm(ixc^s,idim1)>0.d0)
2083 elc(ixc^s)=el(ixc^s)
2084 else where(vnorm(ixc^s,idim1)<0.d0)
2085 elc(ixc^s)=el(jxc^s)
2086 else where
2087 elc(ixc^s)=0.5d0*(el(ixc^s)+el(jxc^s))
2088 end where
2089 hxc^l=ixc^l+kr(idim2,^d);
2090 where(vnorm(hxc^s,idim1)>0.d0)
2091 erc(ixc^s)=er(ixc^s)
2092 else where(vnorm(hxc^s,idim1)<0.d0)
2093 erc(ixc^s)=er(jxc^s)
2094 else where
2095 erc(ixc^s)=0.5d0*(er(ixc^s)+er(jxc^s))
2096 end where
2097 fe(ixc^s,idir)=fe(ixc^s,idir)+0.25d0*(elc(ixc^s)+erc(ixc^s))
2098
2099 ! add slope in idim1 direction from equation (50)
2100 jxc^l=ixc^l+kr(idim2,^d);
2101 ixamin^d=ixcmin^d;
2102 ixamax^d=ixcmax^d+kr(idim2,^d);
2103 el(ixa^s)=-fc(ixa^s,iwdim2,idim1)-ecc(ixa^s,idir)
2104 hxc^l=ixa^l+kr(idim1,^d);
2105 er(ixa^s)=-fc(ixa^s,iwdim2,idim1)-ecc(hxc^s,idir)
2106 where(vnorm(ixc^s,idim2)>0.d0)
2107 elc(ixc^s)=el(ixc^s)
2108 else where(vnorm(ixc^s,idim2)<0.d0)
2109 elc(ixc^s)=el(jxc^s)
2110 else where
2111 elc(ixc^s)=0.5d0*(el(ixc^s)+el(jxc^s))
2112 end where
2113 hxc^l=ixc^l+kr(idim1,^d);
2114 where(vnorm(hxc^s,idim2)>0.d0)
2115 erc(ixc^s)=er(ixc^s)
2116 else where(vnorm(hxc^s,idim2)<0.d0)
2117 erc(ixc^s)=er(jxc^s)
2118 else where
2119 erc(ixc^s)=0.5d0*(er(ixc^s)+er(jxc^s))
2120 end where
2121 fe(ixc^s,idir)=fe(ixc^s,idir)+0.25d0*(elc(ixc^s)+erc(ixc^s))
2122
2123 ! add current component of electric field at cell edges E=-vxB+eta J
2124 if(mf_eta/=zero) fe(ixc^s,idir)=fe(ixc^s,idir)+jce(ixc^s,idir)
2125
2126 if(record_electric_field) s%we(ixc^s,idir)=fe(ixc^s,idir)
2127 ! times time step and edge length
2128 fe(ixc^s,idir)=fe(ixc^s,idir)*qdt*s%dsC(ixc^s,idir)
2129 if (.not.slab) then
2130 where(abs(x(ixc^s,r_)+half*dxlevel(r_))<1.0d-9)
2131 fe(ixc^s,idir)=zero
2132 end where
2133 end if
2134 end if
2135 end do
2136 end do
2137 end do
2138
2139 ! allow user to change inductive electric field, especially for boundary driven applications
2140 if(associated(usr_set_electric_field)) &
2141 call usr_set_electric_field(ixi^l,ixo^l,qt,qdt,fe,sct)
2142
2143 circ(ixi^s,1:ndim)=zero
2144
2145 ! Calculate circulation on each face
2146 do idim1=1,ndim ! Coordinate perpendicular to face
2147 ixcmax^d=ixomax^d;
2148 ixcmin^d=ixomin^d-kr(idim1,^d);
2149 do idim2=1,ndim
2150 do idir=sdim,3 ! Direction of line integral
2151 ! Assemble indices
2152 if(lvc(idim1,idim2,idir)/=0) then
2153 hxc^l=ixc^l-kr(idim2,^d);
2154 ! Add line integrals in direction idir
2155 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
2156 +lvc(idim1,idim2,idir)&
2157 *(fe(ixc^s,idir)&
2158 -fe(hxc^s,idir))
2159 end if
2160 end do
2161 end do
2162 ! Divide by the area of the face to get dB/dt
2163 where(s%surfaceC(ixc^s,idim1) > 1.0d-9*s%dvolume(ixc^s))
2164 circ(ixc^s,idim1)=circ(ixc^s,idim1)/s%surfaceC(ixc^s,idim1)
2165 elsewhere
2166 circ(ixc^s,idim1)=zero
2167 end where
2168 ! Time update cell-face magnetic field component
2169 bfaces(ixc^s,idim1)=bfaces(ixc^s,idim1)-circ(ixc^s,idim1)
2170 end do
2171
2172 end associate
2173
2174 end subroutine mf_update_faces_contact
2175
2176 !> update faces
2177 subroutine mf_update_faces_hll(ixI^L,ixO^L,qt,qdt,wp,fC,fE,sCT,s,vcts)
2179 use mod_usr_methods
2181
2182 integer, intent(in) :: ixi^l, ixo^l
2183 double precision, intent(in) :: qt, qdt
2184 ! cell-center primitive variables
2185 double precision, intent(in) :: wp(ixi^s,1:nw)
2186 type(state) :: sct, s
2187 type(ct_velocity) :: vcts
2188 double precision, intent(in) :: fc(ixi^s,1:nwflux,1:ndim)
2189 double precision, intent(inout) :: fe(ixi^s,sdim:3)
2190
2191 double precision :: vtill(ixi^s,2)
2192 double precision :: vtilr(ixi^s,2)
2193 double precision :: btill(ixi^s,ndim)
2194 double precision :: btilr(ixi^s,ndim)
2195 double precision :: cp(ixi^s,2)
2196 double precision :: cm(ixi^s,2)
2197 double precision :: circ(ixi^s,1:ndim)
2198 ! current on cell edges
2199 double precision :: jce(ixi^s,sdim:3)
2200 integer :: hxc^l,ixc^l,ixcp^l,jxc^l,ixcm^l
2201 integer :: idim1,idim2,idir
2202
2203 associate(bfaces=>s%ws,bfacesct=>sct%ws,x=>s%x,vbarc=>vcts%vbarC,cbarmin=>vcts%cbarmin,&
2204 cbarmax=>vcts%cbarmax)
2205
2206 ! Calculate contribution to FEM of each edge,
2207 ! that is, estimate value of line integral of
2208 ! electric field in the positive idir direction.
2209
2210 ! Loop over components of electric field
2211
2212 ! idir: electric field component we need to calculate
2213 ! idim1: directions in which we already performed the reconstruction
2214 ! idim2: directions in which we perform the reconstruction
2215
2216 ! if there is resistivity, get eta J
2217 if(mf_eta/=zero) call get_resistive_electric_field(ixi^l,ixo^l,sct,s,jce)
2218
2219 do idir=sdim,3
2220 ! Indices
2221 ! idir: electric field component
2222 ! idim1: one surface
2223 ! idim2: the other surface
2224 ! cyclic permutation: idim1,idim2,idir=1,2,3
2225 ! Velocity components on the surface
2226 ! follow cyclic premutations:
2227 ! Sx(1),Sx(2)=y,z ; Sy(1),Sy(2)=z,x ; Sz(1),Sz(2)=x,y
2228
2229 ixcmax^d=ixomax^d;
2230 ixcmin^d=ixomin^d-1+kr(idir,^d);
2231
2232 ! Set indices and directions
2233 idim1=mod(idir,3)+1
2234 idim2=mod(idir+1,3)+1
2235
2236 jxc^l=ixc^l+kr(idim1,^d);
2237 ixcp^l=ixc^l+kr(idim2,^d);
2238
2239 ! Reconstruct transverse transport velocities
2240 call reconstruct(ixi^l,ixc^l,idim2,vbarc(ixi^s,idim1,1),&
2241 vtill(ixi^s,2),vtilr(ixi^s,2))
2242
2243 call reconstruct(ixi^l,ixc^l,idim1,vbarc(ixi^s,idim2,2),&
2244 vtill(ixi^s,1),vtilr(ixi^s,1))
2245
2246 ! Reconstruct magnetic fields
2247 ! Eventhough the arrays are larger, reconstruct works with
2248 ! the limits ixG.
2249 call reconstruct(ixi^l,ixc^l,idim2,bfacesct(ixi^s,idim1),&
2250 btill(ixi^s,idim1),btilr(ixi^s,idim1))
2251
2252 call reconstruct(ixi^l,ixc^l,idim1,bfacesct(ixi^s,idim2),&
2253 btill(ixi^s,idim2),btilr(ixi^s,idim2))
2254
2255 ! Take the maximum characteristic
2256
2257 cm(ixc^s,1)=max(cbarmin(ixcp^s,idim1),cbarmin(ixc^s,idim1))
2258 cp(ixc^s,1)=max(cbarmax(ixcp^s,idim1),cbarmax(ixc^s,idim1))
2259
2260 cm(ixc^s,2)=max(cbarmin(jxc^s,idim2),cbarmin(ixc^s,idim2))
2261 cp(ixc^s,2)=max(cbarmax(jxc^s,idim2),cbarmax(ixc^s,idim2))
2262
2263
2264 ! Calculate eletric field
2265 fe(ixc^s,idir)=-(cp(ixc^s,1)*vtill(ixc^s,1)*btill(ixc^s,idim2) &
2266 + cm(ixc^s,1)*vtilr(ixc^s,1)*btilr(ixc^s,idim2) &
2267 - cp(ixc^s,1)*cm(ixc^s,1)*(btilr(ixc^s,idim2)-btill(ixc^s,idim2)))&
2268 /(cp(ixc^s,1)+cm(ixc^s,1)) &
2269 +(cp(ixc^s,2)*vtill(ixc^s,2)*btill(ixc^s,idim1) &
2270 + cm(ixc^s,2)*vtilr(ixc^s,2)*btilr(ixc^s,idim1) &
2271 - cp(ixc^s,2)*cm(ixc^s,2)*(btilr(ixc^s,idim1)-btill(ixc^s,idim1)))&
2272 /(cp(ixc^s,2)+cm(ixc^s,2))
2273
2274 ! add current component of electric field at cell edges E=-vxB+eta J
2275 if(mf_eta/=zero) fe(ixc^s,idir)=fe(ixc^s,idir)+jce(ixc^s,idir)
2276
2277 if(record_electric_field) s%we(ixc^s,idir)=fe(ixc^s,idir)
2278 ! times time step and edge length
2279 fe(ixc^s,idir)=qdt*s%dsC(ixc^s,idir)*fe(ixc^s,idir)
2280
2281 if (.not.slab) then
2282 where(abs(x(ixc^s,r_)+half*dxlevel(r_)).lt.1.0d-9)
2283 fe(ixc^s,idir)=zero
2284 end where
2285 end if
2286
2287 end do
2288
2289 ! allow user to change inductive electric field, especially for boundary driven applications
2290 if(associated(usr_set_electric_field)) &
2291 call usr_set_electric_field(ixi^l,ixo^l,qt,qdt,fe,sct)
2292
2293 circ(ixi^s,1:ndim)=zero
2294
2295 ! Calculate circulation on each face: interal(fE dot dl)
2296
2297 do idim1=1,ndim ! Coordinate perpendicular to face
2298 ixcmax^d=ixomax^d;
2299 ixcmin^d=ixomin^d-kr(idim1,^d);
2300 do idim2=1,ndim
2301 do idir=sdim,3 ! Direction of line integral
2302 ! Assemble indices
2303 if(lvc(idim1,idim2,idir)/=0) then
2304 hxc^l=ixc^l-kr(idim2,^d);
2305 ! Add line integrals in direction idir
2306 circ(ixc^s,idim1)=circ(ixc^s,idim1)&
2307 +lvc(idim1,idim2,idir)&
2308 *(fe(ixc^s,idir)&
2309 -fe(hxc^s,idir))
2310 end if
2311 end do
2312 end do
2313 ! Divide by the area of the face to get dB/dt
2314 where(s%surfaceC(ixc^s,idim1) > 1.0d-9*s%dvolume(ixc^s))
2315 circ(ixc^s,idim1)=circ(ixc^s,idim1)/s%surfaceC(ixc^s,idim1)
2316 elsewhere
2317 circ(ixc^s,idim1)=zero
2318 end where
2319 ! Time update
2320 bfaces(ixc^s,idim1)=bfaces(ixc^s,idim1)-circ(ixc^s,idim1)
2321 end do
2322
2323 end associate
2324 end subroutine mf_update_faces_hll
2325
2326 !> calculate eta J at cell edges
2327 subroutine get_resistive_electric_field(ixI^L,ixO^L,sCT,s,jce)
2329 use mod_usr_methods
2330 use mod_geometry
2331
2332 integer, intent(in) :: ixi^l, ixo^l
2333 type(state), intent(in) :: sct, s
2334 ! current on cell edges
2335 double precision :: jce(ixi^s,sdim:3)
2336
2337 ! current on cell centers
2338 double precision :: jcc(ixi^s,7-2*ndir:3)
2339 ! location at cell faces
2340 double precision :: xs(ixgs^t,1:ndim)
2341 ! resistivity
2342 double precision :: eta(ixi^s)
2343 double precision :: gradi(ixgs^t)
2344 integer :: ix^d,ixc^l,ixa^l,ixb^l,idir,idirmin,idim1,idim2
2345
2346 associate(x=>s%x,dx=>s%dx,w=>s%w,wct=>sct%w,wcts=>sct%ws)
2347 ! calculate current density at cell edges
2348 jce=0.d0
2349 do idim1=1,ndim
2350 do idim2=1,ndim
2351 do idir=sdim,3
2352 if (lvc(idim1,idim2,idir)==0) cycle
2353 ixcmax^d=ixomax^d+1;
2354 ixcmin^d=ixomin^d+kr(idir,^d)-2;
2355 ixbmax^d=ixcmax^d-kr(idir,^d)+1;
2356 ixbmin^d=ixcmin^d;
2357 ! current at transverse faces
2358 xs(ixb^s,:)=x(ixb^s,:)
2359 xs(ixb^s,idim2)=x(ixb^s,idim2)+half*dx(ixb^s,idim2)
2360 call gradientf(wcts(ixgs^t,idim2),xs,ixgs^ll,ixc^l,idim1,gradi)
2361 if (lvc(idim1,idim2,idir)==1) then
2362 jce(ixc^s,idir)=jce(ixc^s,idir)+gradi(ixc^s)
2363 else
2364 jce(ixc^s,idir)=jce(ixc^s,idir)-gradi(ixc^s)
2365 end if
2366 end do
2367 end do
2368 end do
2369 ! get resistivity
2370 if(mf_eta>zero)then
2371 jce(ixi^s,:)=jce(ixi^s,:)*mf_eta
2372 else
2373 ixa^l=ixo^l^ladd1;
2374 call get_current(wct,ixi^l,ixa^l,idirmin,jcc)
2375 call usr_special_resistivity(wct,ixi^l,ixa^l,idirmin,x,jcc,eta)
2376 ! calcuate eta on cell edges
2377 do idir=sdim,3
2378 ixcmax^d=ixomax^d;
2379 ixcmin^d=ixomin^d+kr(idir,^d)-1;
2380 jcc(ixc^s,idir)=0.d0
2381 {do ix^db=0,1\}
2382 if({ ix^d==1 .and. ^d==idir | .or.}) cycle
2383 ixamin^d=ixcmin^d+ix^d;
2384 ixamax^d=ixcmax^d+ix^d;
2385 jcc(ixc^s,idir)=jcc(ixc^s,idir)+eta(ixa^s)
2386 {end do\}
2387 jcc(ixc^s,idir)=jcc(ixc^s,idir)*0.25d0
2388 jce(ixc^s,idir)=jce(ixc^s,idir)*jcc(ixc^s,idir)
2389 enddo
2390 end if
2391
2392 end associate
2393 end subroutine get_resistive_electric_field
2394
2395 !> calculate cell-center values from face-center values
2396 subroutine mf_face_to_center(ixO^L,s)
2398 ! Non-staggered interpolation range
2399 integer, intent(in) :: ixo^l
2400 type(state) :: s
2401
2402 integer :: ix^d
2403
2404 ! calculate cell-center values from face-center values in 2nd order
2405 {!dir$ ivdep
2406 do ix^db=ixomin^db,ixomax^db\}
2407 {^ifthreed
2408 s%w(ix^d,mag(1))=half/s%surface(ix^d,1)*(s%ws(ix^d,1)*s%surfaceC(ix^d,1)&
2409 +s%ws(ix1-1,ix2,ix3,1)*s%surfaceC(ix1-1,ix2,ix3,1))
2410 s%w(ix^d,mag(2))=half/s%surface(ix^d,2)*(s%ws(ix^d,2)*s%surfaceC(ix^d,2)&
2411 +s%ws(ix1,ix2-1,ix3,2)*s%surfaceC(ix1,ix2-1,ix3,2))
2412 s%w(ix^d,mag(3))=half/s%surface(ix^d,3)*(s%ws(ix^d,3)*s%surfaceC(ix^d,3)&
2413 +s%ws(ix1,ix2,ix3-1,3)*s%surfaceC(ix1,ix2,ix3-1,3))
2414 }
2415 {^iftwod
2416 s%w(ix^d,mag(1))=half/s%surface(ix^d,1)*(s%ws(ix^d,1)*s%surfaceC(ix^d,1)&
2417 +s%ws(ix1-1,ix2,1)*s%surfaceC(ix1-1,ix2,1))
2418 s%w(ix^d,mag(2))=half/s%surface(ix^d,2)*(s%ws(ix^d,2)*s%surfaceC(ix^d,2)&
2419 +s%ws(ix1,ix2-1,2)*s%surfaceC(ix1,ix2-1,2))
2420 }
2421 {end do\}
2422
2423 ! calculate cell-center values from face-center values in 4th order
2424 !do idim=1,ndim
2425 ! gxO^L=ixO^L-2*kr(idim,^D);
2426 ! hxO^L=ixO^L-kr(idim,^D);
2427 ! jxO^L=ixO^L+kr(idim,^D);
2428
2429 ! ! Interpolate to cell barycentre using fourth order central formula
2430 ! w(ixO^S,mag(idim))=(0.0625d0/s%surface(ixO^S,idim))*&
2431 ! ( -ws(gxO^S,idim)*s%surfaceC(gxO^S,idim) &
2432 ! +9.0d0*ws(hxO^S,idim)*s%surfaceC(hxO^S,idim) &
2433 ! +9.0d0*ws(ixO^S,idim)*s%surfaceC(ixO^S,idim) &
2434 ! -ws(jxO^S,idim)*s%surfaceC(jxO^S,idim) )
2435 !end do
2436
2437 ! calculate cell-center values from face-center values in 6th order
2438 !do idim=1,ndim
2439 ! fxO^L=ixO^L-3*kr(idim,^D);
2440 ! gxO^L=ixO^L-2*kr(idim,^D);
2441 ! hxO^L=ixO^L-kr(idim,^D);
2442 ! jxO^L=ixO^L+kr(idim,^D);
2443 ! kxO^L=ixO^L+2*kr(idim,^D);
2444
2445 ! ! Interpolate to cell barycentre using sixth order central formula
2446 ! w(ixO^S,mag(idim))=(0.00390625d0/s%surface(ixO^S,idim))* &
2447 ! ( +3.0d0*ws(fxO^S,idim)*s%surfaceC(fxO^S,idim) &
2448 ! -25.0d0*ws(gxO^S,idim)*s%surfaceC(gxO^S,idim) &
2449 ! +150.0d0*ws(hxO^S,idim)*s%surfaceC(hxO^S,idim) &
2450 ! +150.0d0*ws(ixO^S,idim)*s%surfaceC(ixO^S,idim) &
2451 ! -25.0d0*ws(jxO^S,idim)*s%surfaceC(jxO^S,idim) &
2452 ! +3.0d0*ws(kxO^S,idim)*s%surfaceC(kxO^S,idim) )
2453 !end do
2454
2455 end subroutine mf_face_to_center
2456
2457 !> calculate magnetic field from vector potential
2458 subroutine b_from_vector_potential(ixIs^L, ixI^L, ixO^L, ws, x)
2461
2462 integer, intent(in) :: ixis^l, ixi^l, ixo^l
2463 double precision, intent(inout) :: ws(ixis^s,1:nws)
2464 double precision, intent(in) :: x(ixi^s,1:ndim)
2465
2466 double precision :: adummy(ixis^s,1:3)
2467
2468 call b_from_vector_potentiala(ixis^l, ixi^l, ixo^l, ws, x, adummy)
2469
2470 end subroutine b_from_vector_potential
2471
2474
2475 double precision :: sum_jbb_ipe, sum_j_ipe, sum_l_ipe, f_i_ipe, volume_pe
2476 double precision :: sum_jbb, sum_j, sum_l, f_i, volume, cw_sin_theta
2477 integer :: iigrid, igrid, ix^d
2478 integer :: amode, istatus(mpi_status_size)
2479 integer, save :: fhmf
2480 logical :: patchwi(ixg^t)
2481 logical, save :: logmfopened=.false.
2482 character(len=800) :: filename,filehead
2483 character(len=800) :: line,datastr
2484
2485 sum_jbb_ipe = 0.d0
2486 sum_j_ipe = 0.d0
2487 sum_l_ipe = 0.d0
2488 f_i_ipe = 0.d0
2489 volume_pe=0.d0
2490 do iigrid=1,igridstail; igrid=igrids(iigrid);
2491 block=>ps(igrid)
2492 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
2493 call mask_inner(ixg^ll,ixm^ll,ps(igrid)%w,ps(igrid)%x,patchwi,volume_pe)
2494 sum_jbb_ipe = sum_jbb_ipe+integral_grid_mf(ixg^ll,ixm^ll,ps(igrid)%w,&
2495 ps(igrid)%x,1,patchwi)
2496 sum_j_ipe = sum_j_ipe+integral_grid_mf(ixg^ll,ixm^ll,ps(igrid)%w,&
2497 ps(igrid)%x,2,patchwi)
2498 f_i_ipe=f_i_ipe+integral_grid_mf(ixg^ll,ixm^ll,ps(igrid)%w,&
2499 ps(igrid)%x,3,patchwi)
2500 sum_l_ipe = sum_l_ipe+integral_grid_mf(ixg^ll,ixm^ll,ps(igrid)%w,&
2501 ps(igrid)%x,4,patchwi)
2502 end do
2503 call mpi_reduce(sum_jbb_ipe,sum_jbb,1,mpi_double_precision,&
2504 mpi_sum,0,icomm,ierrmpi)
2505 call mpi_reduce(sum_j_ipe,sum_j,1,mpi_double_precision,mpi_sum,0,&
2506 icomm,ierrmpi)
2507 call mpi_reduce(f_i_ipe,f_i,1,mpi_double_precision,mpi_sum,0,&
2508 icomm,ierrmpi)
2509 call mpi_reduce(sum_l_ipe,sum_l,1,mpi_double_precision,mpi_sum,0,&
2510 icomm,ierrmpi)
2511 call mpi_reduce(volume_pe,volume,1,mpi_double_precision,mpi_sum,0,&
2512 icomm,ierrmpi)
2513
2514 if(mype==0) then
2515 ! current- and volume-weighted average of the sine of the angle
2516 ! between the magnetic field B and the current density J
2517 cw_sin_theta = sum_jbb/sum_j
2518 ! volume-weighted average of the absolute value of the fractional
2519 ! magnetic flux change
2520 f_i = f_i/volume
2521 sum_j=sum_j/volume
2522 sum_l=sum_l/volume
2523 if(.not.logmfopened) then
2524 ! generate filename
2525 write(filename,"(a,a)") trim(base_filename), "_mf_metrics.csv"
2526
2527 ! Delete the log when not doing a restart run
2528 if (restart_from_file == undefined) then
2529 open(unit=20,file=trim(filename),status='replace')
2530 close(20, status='delete')
2531 end if
2532
2533 amode=ior(mpi_mode_create,mpi_mode_wronly)
2534 amode=ior(amode,mpi_mode_append)
2535 call mpi_file_open(mpi_comm_self,filename,amode,mpi_info_null,fhmf,ierrmpi)
2536 logmfopened=.true.
2537 filehead=" it, time, CW_sin_theta, current, Lorenz_force, f_i"
2538 if (it == 0) then
2539 call mpi_file_write(fhmf,filehead,len_trim(filehead), &
2540 mpi_character,istatus,ierrmpi)
2541 call mpi_file_write(fhmf,achar(10),1,mpi_character,istatus,ierrmpi)
2542 end if
2543 end if
2544 line=''
2545 write(datastr,'(i8,a)') it,','
2546 line=trim(line)//trim(datastr)
2547 write(datastr,'(es13.6,a)') global_time,','
2548 line=trim(line)//trim(datastr)
2549 write(datastr,'(es13.6,a)') cw_sin_theta,','
2550 line=trim(line)//trim(datastr)
2551 write(datastr,'(es13.6,a)') sum_j,','
2552 line=trim(line)//trim(datastr)
2553 write(datastr,'(es13.6,a)') sum_l,','
2554 line=trim(line)//trim(datastr)
2555 write(datastr,'(es13.6)') f_i
2556 line=trim(line)//trim(datastr)//new_line('A')
2557 call mpi_file_write(fhmf,line,len_trim(line),mpi_character,istatus,ierrmpi)
2558 if(.not.time_advance) then
2559 call mpi_file_close(fhmf,ierrmpi)
2560 logmfopened=.false.
2561 end if
2562 end if
2563
2564 end subroutine record_force_free_metrics
2565
2566 subroutine mask_inner(ixI^L,ixO^L,w,x,patchwi,volume_pe)
2568
2569 integer, intent(in) :: ixi^l,ixo^l
2570 double precision, intent(in):: w(ixi^s,nw),x(ixi^s,1:ndim)
2571 double precision, intent(inout) :: volume_pe
2572 logical, intent(out) :: patchwi(ixi^s)
2573
2574 double precision :: xo^l
2575 integer :: ix^d
2576
2577 {xomin^d = xprobmin^d + 0.05d0*(xprobmax^d-xprobmin^d)\}
2578 {xomax^d = xprobmax^d - 0.05d0*(xprobmax^d-xprobmin^d)\}
2579 if(slab) then
2580 xomin^nd = xprobmin^nd
2581 else
2582 xomin1 = xprobmin1
2583 end if
2584
2585 {do ix^db=ixomin^db,ixomax^db\}
2586 if({ x(ix^dd,^d) > xomin^d .and. x(ix^dd,^d) < xomax^d | .and. }) then
2587 patchwi(ix^d)=.true.
2588 volume_pe=volume_pe+block%dvolume(ix^d)
2589 else
2590 patchwi(ix^d)=.false.
2591 endif
2592 {end do\}
2593
2594 end subroutine mask_inner
2595
2596 function integral_grid_mf(ixI^L,ixO^L,w,x,iw,patchwi)
2598
2599 integer, intent(in) :: ixi^l,ixo^l,iw
2600 double precision, intent(in) :: x(ixi^s,1:ndim)
2601 double precision :: w(ixi^s,nw+nwauxio)
2602 logical, intent(in) :: patchwi(ixi^s)
2603
2604 double precision, dimension(ixI^S,7-2*ndir:3) :: current
2605 double precision, dimension(ixI^S,1:ndir) :: bvec,qvec
2606 double precision :: integral_grid_mf,tmp(ixi^s),bm2
2607 integer :: ix^d,idirmin0,idirmin,idir,jdir,kdir
2608
2609 integral_grid_mf=0.d0
2610 select case(iw)
2611 case(1)
2612 ! Sum(|JxB|/|B|*dvolume)
2613 bvec(ixo^s,:)=w(ixo^s,mag(:))
2614 call get_current(w,ixi^l,ixo^l,idirmin,current)
2615 ! calculate Lorentz force
2616 qvec(ixo^s,1:ndir)=zero
2617 do idir=1,ndir; do jdir=idirmin,3; do kdir=1,ndir
2618 if(lvc(idir,jdir,kdir)/=0)then
2619 if(lvc(idir,jdir,kdir)==1)then
2620 qvec(ixo^s,idir)=qvec(ixo^s,idir)+current(ixo^s,jdir)*bvec(ixo^s,kdir)
2621 else
2622 qvec(ixo^s,idir)=qvec(ixo^s,idir)-current(ixo^s,jdir)*bvec(ixo^s,kdir)
2623 end if
2624 end if
2625 end do; end do; end do
2626
2627 {do ix^db=ixomin^db,ixomax^db\}
2628 if(patchwi(ix^d)) then
2629 bm2=sum(bvec(ix^d,:)**2)
2630 if(bm2/=0.d0) bm2=1.d0/bm2
2631 integral_grid_mf=integral_grid_mf+sqrt(sum(qvec(ix^d,:)**2)*&
2632 bm2)*block%dvolume(ix^d)
2633 end if
2634 {end do\}
2635 case(2)
2636 ! Sum(|J|*dvolume)
2637 call get_current(w,ixi^l,ixo^l,idirmin,current)
2638 {do ix^db=ixomin^db,ixomax^db\}
2639 if(patchwi(ix^d)) integral_grid_mf=integral_grid_mf+sqrt(sum(current(ix^d,:)**2))*&
2640 block%dvolume(ix^d)
2641 {end do\}
2642 case(3)
2643 ! f_i solenoidal property of B: (dvolume |div B|)/(dsurface |B|)
2644 ! Sum(f_i*dvolume)
2645 call get_normalized_divb(w,ixi^l,ixo^l,tmp)
2646 {do ix^db=ixomin^db,ixomax^db\}
2647 if(patchwi(ix^d)) integral_grid_mf=integral_grid_mf+tmp(ix^d)*block%dvolume(ix^d)
2648 {end do\}
2649 case(4)
2650 ! Sum(|JxB|*dvolume)
2651 bvec(ixo^s,:)=w(ixo^s,mag(:))
2652 call get_current(w,ixi^l,ixo^l,idirmin,current)
2653 ! calculate Lorentz force
2654 qvec(ixo^s,1:ndir)=zero
2655 do idir=1,ndir; do jdir=idirmin,3; do kdir=1,ndir
2656 if(lvc(idir,jdir,kdir)/=0)then
2657 if(lvc(idir,jdir,kdir)==1)then
2658 qvec(ixo^s,idir)=qvec(ixo^s,idir)+current(ixo^s,jdir)*bvec(ixo^s,kdir)
2659 else
2660 qvec(ixo^s,idir)=qvec(ixo^s,idir)-current(ixo^s,jdir)*bvec(ixo^s,kdir)
2661 end if
2662 end if
2663 end do; end do; end do
2664
2665 {do ix^db=ixomin^db,ixomax^db\}
2666 if(patchwi(ix^d)) integral_grid_mf=integral_grid_mf+&
2667 sqrt(sum(qvec(ix^d,:)**2))*block%dvolume(ix^d)
2668 {end do\}
2669 end select
2670 return
2671 end function integral_grid_mf
2672
2673end module mod_mf_phys
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
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
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
integer coordinate
Definition mod_geometry.t:7
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
integer, dimension(0:3^d &), target type_recv_p_p1
integer, dimension( :^d &), pointer type_recv_r
integer, dimension(-1:1^d &), target type_send_srl_p1
subroutine getbc(time, qdt, psb, nwstart, nwbc)
do update ghost cells of all blocks including physical boundaries
integer, dimension(0:3^d &), target type_send_p_f
integer, dimension( :^d &), pointer type_send_p
integer, dimension( :^d &), pointer type_send_srl
integer, dimension(-1:1^d &), target type_send_r_p1
subroutine create_bc_mpi_datatype(nwstart, nwbc)
integer, dimension(0:3^d &), target type_recv_r_f
integer, dimension(-1:1^d &), target type_recv_srl_f
integer, dimension(-1:1^d &), target type_send_r_f
integer, dimension(-1:1^d &), target type_send_srl_f
integer, dimension(0:3^d &), target type_send_p_p1
integer, dimension( :^d &), pointer type_send_r
integer, dimension(0:3^d &), target type_recv_r_p1
integer, dimension( :^d &), pointer type_recv_p
integer, dimension(-1:1^d &), target type_recv_srl_p1
integer, dimension(0:3^d &), target type_recv_p_f
integer, dimension( :^d &), pointer type_recv_srl
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
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.
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.
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.
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
logical use_particles
Use particles module or not.
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 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 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.
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision unit_velocity
Physical scaling factor for velocity.
logical time_advance
do time evolving
character(len=std_len) restart_from_file
If not 'unavailable', resume from snapshot with this base file name.
double precision c_norm
Normalised speed of light.
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
character(len= *), parameter undefined
logical need_global_cmax
need global maximal wave speed
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)
character(len=std_len) base_filename
Base file name for simulation output, which will be followed by a number.
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:
logical record_electric_field
True for record electric field.
integer, parameter ixglo
Lower index of grid block arrays (always 1)
Magnetofriction module.
Definition mod_mf_phys.t:2
integer, public, protected b
Definition mod_mf_phys.t:46
double precision, public mf_nu
viscosity coefficient in s cm^-2 for solar corona (Cheung 2012 ApJ)
Definition mod_mf_phys.t:12
subroutine, public mf_phys_init()
integer, public, protected mf_divb_nth
Whether divB is computed with a fourth order approximation.
Definition mod_mf_phys.t:82
subroutine, public b_from_vector_potential(ixisl, ixil, ixol, ws, x)
calculate magnetic field from vector potential
logical, public, protected mf_record_electric_field
set to true if need to record electric field on cell edges
Definition mod_mf_phys.t:79
logical, public, protected clean_initial_divb
clean divb in the initial condition
Definition mod_mf_phys.t:94
logical, public divbwave
Add divB wave in Roe solver.
Definition mod_mf_phys.t:91
subroutine, public mf_face_to_center(ixol, s)
calculate cell-center values from face-center values
double precision, public mf_vmax
maximal limit of magnetofrictional velocity in cm s^-1 (Pomoell 2019 A&A)
Definition mod_mf_phys.t:15
logical, public, protected mf_glm
Whether GLM-MHD is used.
Definition mod_mf_phys.t:70
double precision, public mf_eta_hyper
The hyper-resistivity.
Definition mod_mf_phys.t:28
double precision, public mf_glm_alpha
GLM-MHD parameter: ratio of the diffusive and advective time scales for div b taking values within [0...
Definition mod_mf_phys.t:22
integer, public, protected c
Indices of the momentum density for the form of better vectorization.
Definition mod_mf_phys.t:43
integer, public, protected c_
Definition mod_mf_phys.t:43
logical, dimension(2 *^nd), public, protected boundary_divbfix
To control divB=0 fix for boundary.
Definition mod_mf_phys.t:85
character(len=std_len), public, protected type_ct
Method type of constrained transport.
character(len=std_len), public, protected typedivbfix
Method type to clean divergence of B.
Definition mod_mf_phys.t:97
subroutine, public get_normalized_divb(w, ixil, ixol, divb)
get dimensionless div B = |divB| * volume / area / |B|
integer, dimension(2 *^nd), public, protected boundary_divbfix_skip
To skip * layer of ghost cells during divB=0 fix for boundary.
Definition mod_mf_phys.t:52
double precision, dimension(2 *^nd), public mf_decay_scale
decay scale of frictional velocity near boundaries
Definition mod_mf_phys.t:18
subroutine, public record_force_free_metrics()
subroutine, public mf_to_conserved(ixil, ixol, w, x)
Transform primitive variables into conservative ones.
logical, public, protected mf_4th_order
MHD fourth order.
Definition mod_mf_phys.t:76
subroutine, public get_current(w, ixil, ixol, idirmin, current, fourthorder)
Calculate idirmin and the idirmin:3 components of the common current array make sure that dxlevel(^D)...
integer, public, protected psi_
Indices of the GLM psi.
Definition mod_mf_phys.t:49
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
Definition mod_mf_phys.t:40
subroutine, public mf_clean_divb_multigrid(qdt, qt, active)
logical, public, protected mf_particles
Whether particles module is added.
Definition mod_mf_phys.t:67
integer, public, protected m
Definition mod_mf_phys.t:43
subroutine, public mf_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
double precision, public mf_eta
The resistivity.
Definition mod_mf_phys.t:25
logical, public, protected source_split_divb
Whether divB cleaning sources are added splitting from fluid solver.
Definition mod_mf_phys.t:73
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
Definition mod_mf_phys.t:34
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.
Module containing all the particle routines.
subroutine particles_init()
Initialize particle data and parameters.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_get_ct_velocity), pointer phys_get_ct_velocity
Definition mod_physics.t:80
procedure(sub_convert), pointer phys_to_primitive
Definition mod_physics.t:57
integer, parameter flux_tvdlf
Indicates the flux should be treated with tvdlf.
Definition mod_physics.t:26
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:78
procedure(sub_get_flux), pointer phys_get_flux
Definition mod_physics.t:65
procedure(sub_get_cbounds), pointer phys_get_cbounds
Definition mod_physics.t:64
procedure(sub_get_dt), pointer phys_get_dt
Definition mod_physics.t:68
procedure(sub_add_source_geom), pointer phys_add_source_geom
Definition mod_physics.t:69
procedure(sub_check_params), pointer phys_check_params
Definition mod_physics.t:54
integer, parameter flux_default
Indicates a normal flux.
Definition mod_physics.t:24
integer transverse_ghost_cells
extra ghost cells in the transverse dimensions for fluxes of CT
Definition mod_physics.t:34
integer, dimension(:, :), allocatable flux_type
Array per direction per variable, which can be used to specify that certain fluxes have to be treated...
Definition mod_physics.t:21
procedure(sub_update_faces), pointer phys_update_faces
Definition mod_physics.t:81
procedure(sub_convert), pointer phys_to_conserved
Definition mod_physics.t:56
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:52
procedure(sub_clean_divb), pointer phys_clean_divb
Definition mod_physics.t:85
procedure(sub_boundary_adjust), pointer phys_boundary_adjust
Definition mod_physics.t:77
procedure(sub_global_source), pointer phys_global_source_after
Definition mod_physics.t:71
procedure(sub_add_source), pointer phys_add_source
Definition mod_physics.t:70
procedure(sub_face_to_center), pointer phys_face_to_center
Definition mod_physics.t:82
procedure(sub_modify_wlr), pointer phys_modify_wlr
Definition mod_physics.t:58
procedure(sub_get_cmax), pointer phys_get_cmax
Definition mod_physics.t:59
logical phys_energy
Solve energy equation or not.
Definition mod_physics.t:37
procedure(sub_special_advance), pointer phys_special_advance
Definition mod_physics.t:72
Module with all the methods that users can customize in AMRVAC.
procedure(special_resistivity), pointer usr_special_resistivity
procedure(set_electric_field), pointer usr_set_electric_field
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 data structure that contains information about a tree node/grid block.
Definition mod_forest.t:11