41 double precision,
intent(in) :: qLunit,qBunit
42 double precision :: xc1,xc2,dxm1,dxm2
43 integer,
dimension(MPI_STATUS_SIZE) :: statuss
44 integer :: file_handle,i
46 character(len=*),
intent(in) :: magnetogramname
51 inquire(file=magnetogramname,exist=aexist)
53 if(
mype==0)
write(*,
'(2a)')
"can not find file:",magnetogramname
54 call mpistop(
"no input magnetogram----init_b_fff_data")
56 call mpi_file_open(
icomm,magnetogramname,mpi_mode_rdonly,mpi_info_null,&
58 call mpi_file_read_all(file_handle,
nx1,1,mpi_integer,statuss,
ierrmpi)
59 call mpi_file_read_all(file_handle,
nx2,1,mpi_integer,statuss,
ierrmpi)
61 call mpi_file_read_all(file_handle,xc1,1,mpi_double_precision,statuss,
ierrmpi)
62 call mpi_file_read_all(file_handle,xc2,1,mpi_double_precision,statuss,
ierrmpi)
63 call mpi_file_read_all(file_handle,dxm1,1,mpi_double_precision,statuss,
ierrmpi)
64 call mpi_file_read_all(file_handle,dxm2,1,mpi_double_precision,statuss,
ierrmpi)
65 call mpi_file_read_all(file_handle,
bz0,
nx1*
nx2,mpi_double_precision,&
67 call mpi_file_close(file_handle,
ierrmpi)
71 xa1(i) = xc1 + (dble(i) - dble(
nx1)/2.d0 - 0.5d0)*dxm1
74 xa2(i) = xc2 + (dble(i) - dble(
nx2)/2.d0 - 0.5d0)*dxm2
89 print*,
'magnetogram xrange:',minval(
xa1),maxval(
xa1)
90 print*,
'magnetogram yrange:',minval(
xa2),maxval(
xa2)
94 print*,
'extrapolating 3D force-free field from an observed Bz '
95 print*,
'magnetogram of',
nx1,
'by',
nx2,
'pixels. Bzmax=',
bzmax
110 integer,
intent(in) :: ixI^L, ixO^L
111 integer,
optional,
intent(in) :: idir
112 double precision,
intent(in) :: x(ixI^S,1:ndim),alpha,zshift
113 double precision,
intent(inout) :: Bf(ixI^S,1:ndir)
115 double precision,
dimension(ixO^S) :: cos_az,sin_az,zk,bigr,r,r2,r3,cos_ar,sin_ar,g,dgdz
116 double precision :: gx(ixO^S,1:ndim),twopiinv
117 integer :: idim,ixp1,ixp2
122 zk(ixo^s)=x(ixo^s,3)-xprobmin3+zshift
123 cos_az(ixo^s)=dcos(alpha*zk(ixo^s))
124 sin_az(ixo^s)=dsin(alpha*zk(ixo^s))
128 bigr(ixo^s)=dsqrt((x(ixo^s,1)-
xa1(ixp1))**2+&
129 (x(ixo^s,2)-
xa2(ixp2))**2)
147 g=(zk*cos_ar*r-cos_az)*bigr
148 dgdz=(cos_ar*(r-zk**2*r3)-alpha*zk**2*sin_ar*r2+alpha*sin_az)*bigr
150 if(
present(idir).and.idim/=idir) cycle
153 bf(ixo^s,1)=bf(ixo^s,1)+
bz0(ixp1,ixp2)*((x(ixo^s,1)-
xa1(ixp1))*dgdz(ixo^s)&
154 +alpha*g(ixo^s)*(x(ixo^s,2)-
xa2(ixp2)))*bigr(ixo^s)
156 bf(ixo^s,2)=bf(ixo^s,2)+
bz0(ixp1,ixp2)*((x(ixo^s,2)-
xa2(ixp2))*dgdz(ixo^s)&
157 -alpha*g(ixo^s)*(x(ixo^s,1)-
xa1(ixp1)))*bigr(ixo^s)
159 bf(ixo^s,3)=bf(ixo^s,3)+
bz0(ixp1,ixp2)*(zk(ixo^s)*cos_ar(ixo^s)*r3(ixo^s)+alpha*&
160 zk(ixo^s)*sin_ar(ixo^s)*r2(ixo^s))
165 bf(ixo^s,:)=bf(ixo^s,:)*twopiinv
318 type(mg_box_t),
intent(in) :: box
319 integer,
intent(in) :: nc
320 integer,
intent(in) :: iv
321 integer,
intent(in) :: nb
322 integer,
intent(out) :: bc_type
324 double precision,
intent(out) :: bc(nc, nc)
326 double precision :: rr(nc, nc, 3)
327 double precision :: rmina,rminb,rmaxa,rmaxb,xmina,xminb,xmaxa,xmaxb
329 double precision :: wbn(ixG^T)
330 double precision,
allocatable :: xcoarse(:,:,:)
331 integer :: iigrid,igrid,ix^D,idir,ixbca,ixbcb,ixbcn,dlvl,wnc
333 bc_type = mg_bc_neumann
335 call mg_get_face_coords(box, nb, nc, rr)
343 if(mod(nb,2)==0)
then
348 rmina=rr(1,1,2)-0.5d0*box%dr(2)
349 rmaxa=rr(nc,1,2)+0.5d0*box%dr(2)
350 rminb=rr(1,1,3)-0.5d0*box%dr(3)
351 rmaxb=rr(1,nc,3)+0.5d0*box%dr(3)
352 do iigrid=1,igridstail; igrid=igrids(iigrid);
354 if(.not.
block%is_physical_boundary(nb)) cycle
357 wbn(ixbcn^%1ixg^t)=
block%ws(ixbcn^%1ixg^t,idir)
359 wbn(ixbcn^%1ixg^t)=half*(
block%w(ixbcn^%1ixg^t,iw_mag(idir))+
block%w(ixbcn+1^%1ixg^t,iw_mag(idir)))
361 xmina=
block%x(1,1,1,2)-0.5d0*
rnode(rpdx2_,igrid)
362 xmaxa=
block%x(1,ixghi2,1,2)+0.5d0*
rnode(rpdx2_,igrid)
363 xminb=
block%x(1,1,1,3)-0.5d0*
rnode(rpdx3_,igrid)
364 xmaxb=
block%x(1,1,ixghi3,3)+0.5d0*
rnode(rpdx3_,igrid)
365 if(xmina<rr(1,1,2) .and. xmaxa>rr(nc,1,2) .and. &
366 xminb<rr(1,1,3) .and. xmaxb>rr(1,nc,3))
then
369 ixbca=ceiling((rr(ix1,ix2,2)-xmina)/
rnode(rpdx2_,igrid))
370 ixbcb=ceiling((rr(ix1,ix2,3)-xminb)/
rnode(rpdx3_,igrid))
371 bc(ix1,ix2)=wbn(ixbcn,ixbca,ixbcb)
374 else if(
block%x(1,ixmlo2,1,2)>rmina .and.
block%x(1,ixmhi2,1,2)<rmaxa .and. &
375 block%x(1,1,ixmlo3,3)>rminb .and.
block%x(1,1,ixmhi3,3)<rmaxb)
then
378 allocate(xcoarse(wnc,wnc,2))
383 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(2))
384 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(3))
393 if(mod(nb,2)==0)
then
398 rmina=rr(1,1,1)-0.5d0*box%dr(1)
399 rmaxa=rr(nc,1,1)+0.5d0*box%dr(1)
400 rminb=rr(1,1,3)-0.5d0*box%dr(3)
401 rmaxb=rr(1,nc,3)+0.5d0*box%dr(3)
402 do iigrid=1,igridstail; igrid=igrids(iigrid);
404 if(.not.
block%is_physical_boundary(nb)) cycle
406 wbn(ixbcn^%2ixg^t)=
block%ws(ixbcn^%2ixg^t,idir)
408 wbn(ixbcn^%2ixg^t)=half*(
block%w(ixbcn^%2ixg^t,iw_mag(idir))+
block%w(ixbcn+1^%2ixg^t,iw_mag(idir)))
410 xmina=
block%x(1,1,1,1)-0.5d0*
rnode(rpdx1_,igrid)
411 xmaxa=
block%x(ixghi1,1,1,1)+0.5d0*
rnode(rpdx1_,igrid)
412 xminb=
block%x(1,1,1,3)-0.5d0*
rnode(rpdx3_,igrid)
413 xmaxb=
block%x(1,1,ixghi3,3)+0.5d0*
rnode(rpdx3_,igrid)
414 if(xmina<rr(1,1,1) .and. xmaxa>rr(nc,1,1) .and. &
415 xminb<rr(1,1,3) .and. xmaxb>rr(1,nc,3))
then
418 ixbca=ceiling((rr(ix1,ix2,1)-xmina)/
rnode(rpdx1_,igrid))
419 ixbcb=ceiling((rr(ix1,ix2,3)-xminb)/
rnode(rpdx3_,igrid))
420 bc(ix1,ix2)=wbn(ixbca,ixbcn,ixbcb)
423 else if(
block%x(ixmlo1,1,1,1)>rmina .and.
block%x(ixmhi1,1,1,1)<rmaxa .and. &
424 block%x(1,1,ixmlo3,3)>rminb .and.
block%x(1,1,ixmhi3,3)<rmaxb)
then
427 allocate(xcoarse(wnc,wnc,2))
432 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(1))
433 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(3))
442 if(mod(nb,2)==0)
then
447 rmina=rr(1,1,1)-0.5d0*box%dr(1)
448 rmaxa=rr(nc,1,1)+0.5d0*box%dr(1)
449 rminb=rr(1,1,2)-0.5d0*box%dr(2)
450 rmaxb=rr(1,nc,2)+0.5d0*box%dr(2)
451 do iigrid=1,igridstail; igrid=igrids(iigrid);
453 if(.not.
block%is_physical_boundary(nb)) cycle
455 wbn(ixbcn^%3ixg^t)=
block%ws(ixbcn^%3ixg^t,idir)
457 wbn(ixbcn^%3ixg^t)=half*(
block%w(ixbcn^%3ixg^t,iw_mag(idir))+
block%w(ixbcn+1^%3ixg^t,iw_mag(idir)))
459 xmina=
block%x(1,1,1,1)-0.5d0*
rnode(rpdx1_,igrid)
460 xmaxa=
block%x(ixghi1,1,1,1)+0.5d0*
rnode(rpdx1_,igrid)
461 xminb=
block%x(1,1,1,2)-0.5d0*
rnode(rpdx2_,igrid)
462 xmaxb=
block%x(1,ixghi2,1,2)+0.5d0*
rnode(rpdx2_,igrid)
463 if(xmina<rr(1,1,1) .and. xmaxa>rr(nc,1,1) .and. &
464 xminb<rr(1,1,2) .and. xmaxb>rr(1,nc,2))
then
467 ixbca=ceiling((rr(ix1,ix2,1)-xmina)/
rnode(rpdx1_,igrid))
468 ixbcb=ceiling((rr(ix1,ix2,2)-xminb)/
rnode(rpdx2_,igrid))
469 bc(ix1,ix2)=wbn(ixbca,ixbcb,ixbcn)
472 else if(
block%x(ixmlo1,1,1,1)>rmina .and.
block%x(ixmhi1,1,1,1)<rmaxa .and. &
473 block%x(1,ixmlo2,1,2)>rminb .and.
block%x(1,ixmhi2,1,2)<rmaxb)
then
476 allocate(xcoarse(wnc,wnc,2))
481 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(1))
482 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(2))