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