31 double precision,
allocatable,
save ::
bz0(:,:)
32 double precision,
allocatable,
save ::
xa1(:),
xa2(:)
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
178 integer,
intent(in) :: ixI^L, ixO^L
179 double precision,
intent(in) :: x(ixI^S,1:ndim),zshift
180 double precision,
intent(inout) :: potential(ixI^S)
182 double precision,
dimension(ixO^S) :: zk,bigr
185 zk(ixo^s)=x(ixo^s,3)-xprobmin3+zshift
191 bigr(ixo^s)=dsqrt((x(ixo^s,1)-
xa1(ixp1))**2+&
192 (x(ixo^s,2)-
xa2(ixp2))**2+&
194 potential(ixo^s)=potential(ixo^s)+0.5d0*
bz0(ixp1,ixp2)/bigr*
darea/dpi
207 integer,
intent(in) :: ixI^L,nth,nph
208 real*8,
intent(in) :: x(ixi^s,1:
ndim)
210 real*8,
intent(in) :: magnetogram(nth,nph)
212 real*8,
intent(in) :: theta(nth),phi(nph)
214 real*8,
intent(in) :: r_sphere
215 real*8,
intent(out) :: potential(ixi^s)
217 real*8 :: area(nth),distance(ixi^s),dtheta_half,dphi,inv2pi
222 dtheta_half=0.5d0*(theta(2)-theta(1))
224 area(1:nth)=2.d0*r_sphere**2*sin(theta(1:nth))*sin(dtheta_half)*sin(dphi)
225 inv2pi=-1.d0/(2.d0*dpi)
230 distance(ixi^s)=sqrt(x(ixi^s,1)**2+r_sphere**2-2.d0*x(ixi^s,1)*r_sphere*&
231 (sin(x(ixi^s,2))*sin(theta(ix1))*cos(phi(ix2)-x(ixi^s,3))+cos(x(ixi^s,2))*&
233 where(distance(ixi^s)/=0.d0)
234 distance(ixi^s)=1.d0/distance(ixi^s)
236 potential(ixi^s)=potential(ixi^s)+inv2pi*magnetogram(ix1,ix2)*distance(ixi^s)*area(ix1)
251 real*8,
intent(out) :: benergy
252 real*8 :: block_benergy(
ixm^t), mype_benergy
253 real*8,
allocatable :: tmp(:,:,:)
254 integer :: iigrid, igrid, idir
255 integer :: id,nc, lvl, ixI^L
262 do iigrid=1,igridstail; igrid=igrids(iigrid);
266 lvl=
mg%boxes(id)%lvl
267 nc =
mg%box_size_lvl(lvl)
271 block_benergy(
ixm^t)=block_benergy(
ixm^t)+tmp(
ixm^t)**2
273 mype_benergy=mype_benergy+sum(0.5d0*block_benergy(
ixm^t)*
block%dvolume(
ixm^t))
276 call mpi_allreduce(mype_benergy,benergy,1,mpi_double_precision,mpi_sum,
icomm,
ierrmpi)
285 double precision :: max_residual, res
286 integer :: i, max_its
288 mg%operator_type = mg_laplacian
293 mg%bc(:, mg_iphi)%bc_type = mg_bc_neumann
302 call mg_fas_fmg(
mg,.true.,max_res=res)
303 if(
mype==0)
write(*,*)
'mg iteration',i,
'residual:', res
304 if(res < max_residual)
exit
306 if(res > max_residual)
call mpistop(
"get potential multigrid: no convergence")
314 type(mg_box_t),
intent(in) :: box
315 integer,
intent(in) :: nc
316 integer,
intent(in) :: iv
317 integer,
intent(in) :: nb
318 integer,
intent(out) :: bc_type
320 double precision,
intent(out) :: bc(nc, nc)
322 double precision :: rr(nc, nc, 3)
323 double precision :: rmina,rminb,rmaxa,rmaxb,xmina,xminb,xmaxa,xmaxb
325 double precision :: wbn(ixG^T)
326 double precision,
allocatable :: xcoarse(:,:,:)
327 integer :: iigrid,igrid,ix^D,idir,ixbca,ixbcb,ixbcn,dlvl,wnc
329 bc_type = mg_bc_neumann
331 call mg_get_face_coords(box, nb, nc, rr)
339 if(mod(nb,2)==0)
then
344 rmina=rr(1,1,2)-0.5d0*box%dr(2)
345 rmaxa=rr(nc,1,2)+0.5d0*box%dr(2)
346 rminb=rr(1,1,3)-0.5d0*box%dr(3)
347 rmaxb=rr(1,nc,3)+0.5d0*box%dr(3)
348 do iigrid=1,igridstail; igrid=igrids(iigrid);
350 if(.not.
block%is_physical_boundary(nb)) cycle
353 wbn(ixbcn^%1ixg^t)=
block%ws(ixbcn^%1ixg^t,idir)
355 wbn(ixbcn^%1ixg^t)=half*(
block%w(ixbcn^%1ixg^t,iw_mag(idir))+
block%w(ixbcn+1^%1ixg^t,iw_mag(idir)))
357 xmina=
block%x(1,1,1,2)-0.5d0*
rnode(rpdx2_,igrid)
358 xmaxa=
block%x(1,ixghi2,1,2)+0.5d0*
rnode(rpdx2_,igrid)
359 xminb=
block%x(1,1,1,3)-0.5d0*
rnode(rpdx3_,igrid)
360 xmaxb=
block%x(1,1,ixghi3,3)+0.5d0*
rnode(rpdx3_,igrid)
361 if(xmina<rr(1,1,2) .and. xmaxa>rr(nc,1,2) .and. &
362 xminb<rr(1,1,3) .and. xmaxb>rr(1,nc,3))
then
365 ixbca=ceiling((rr(ix1,ix2,2)-xmina)/
rnode(rpdx2_,igrid))
366 ixbcb=ceiling((rr(ix1,ix2,3)-xminb)/
rnode(rpdx3_,igrid))
367 bc(ix1,ix2)=wbn(ixbcn,ixbca,ixbcb)
370 else if(
block%x(1,ixmlo2,1,2)>rmina .and.
block%x(1,ixmhi2,1,2)<rmaxa .and. &
371 block%x(1,1,ixmlo3,3)>rminb .and.
block%x(1,1,ixmhi3,3)<rmaxb)
then
374 allocate(xcoarse(wnc,wnc,2))
379 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(2))
380 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(3))
389 if(mod(nb,2)==0)
then
394 rmina=rr(1,1,1)-0.5d0*box%dr(1)
395 rmaxa=rr(nc,1,1)+0.5d0*box%dr(1)
396 rminb=rr(1,1,3)-0.5d0*box%dr(3)
397 rmaxb=rr(1,nc,3)+0.5d0*box%dr(3)
398 do iigrid=1,igridstail; igrid=igrids(iigrid);
400 if(.not.
block%is_physical_boundary(nb)) cycle
402 wbn(ixbcn^%2ixg^t)=
block%ws(ixbcn^%2ixg^t,idir)
404 wbn(ixbcn^%2ixg^t)=half*(
block%w(ixbcn^%2ixg^t,iw_mag(idir))+
block%w(ixbcn+1^%2ixg^t,iw_mag(idir)))
406 xmina=
block%x(1,1,1,1)-0.5d0*
rnode(rpdx1_,igrid)
407 xmaxa=
block%x(ixghi1,1,1,1)+0.5d0*
rnode(rpdx1_,igrid)
408 xminb=
block%x(1,1,1,3)-0.5d0*
rnode(rpdx3_,igrid)
409 xmaxb=
block%x(1,1,ixghi3,3)+0.5d0*
rnode(rpdx3_,igrid)
410 if(xmina<rr(1,1,1) .and. xmaxa>rr(nc,1,1) .and. &
411 xminb<rr(1,1,3) .and. xmaxb>rr(1,nc,3))
then
414 ixbca=ceiling((rr(ix1,ix2,1)-xmina)/
rnode(rpdx1_,igrid))
415 ixbcb=ceiling((rr(ix1,ix2,3)-xminb)/
rnode(rpdx3_,igrid))
416 bc(ix1,ix2)=wbn(ixbca,ixbcn,ixbcb)
419 else if(
block%x(ixmlo1,1,1,1)>rmina .and.
block%x(ixmhi1,1,1,1)<rmaxa .and. &
420 block%x(1,1,ixmlo3,3)>rminb .and.
block%x(1,1,ixmhi3,3)<rmaxb)
then
423 allocate(xcoarse(wnc,wnc,2))
428 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(1))
429 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(3))
438 if(mod(nb,2)==0)
then
443 rmina=rr(1,1,1)-0.5d0*box%dr(1)
444 rmaxa=rr(nc,1,1)+0.5d0*box%dr(1)
445 rminb=rr(1,1,2)-0.5d0*box%dr(2)
446 rmaxb=rr(1,nc,2)+0.5d0*box%dr(2)
447 do iigrid=1,igridstail; igrid=igrids(iigrid);
449 if(.not.
block%is_physical_boundary(nb)) cycle
451 wbn(ixbcn^%3ixg^t)=
block%ws(ixbcn^%3ixg^t,idir)
453 wbn(ixbcn^%3ixg^t)=half*(
block%w(ixbcn^%3ixg^t,iw_mag(idir))+
block%w(ixbcn+1^%3ixg^t,iw_mag(idir)))
455 xmina=
block%x(1,1,1,1)-0.5d0*
rnode(rpdx1_,igrid)
456 xmaxa=
block%x(ixghi1,1,1,1)+0.5d0*
rnode(rpdx1_,igrid)
457 xminb=
block%x(1,1,1,2)-0.5d0*
rnode(rpdx2_,igrid)
458 xmaxb=
block%x(1,ixghi2,1,2)+0.5d0*
rnode(rpdx2_,igrid)
459 if(xmina<rr(1,1,1) .and. xmaxa>rr(nc,1,1) .and. &
460 xminb<rr(1,1,2) .and. xmaxb>rr(1,nc,2))
then
463 ixbca=ceiling((rr(ix1,ix2,1)-xmina)/
rnode(rpdx1_,igrid))
464 ixbcb=ceiling((rr(ix1,ix2,2)-xminb)/
rnode(rpdx2_,igrid))
465 bc(ix1,ix2)=wbn(ixbca,ixbcb,ixbcn)
468 else if(
block%x(ixmlo1,1,1,1)>rmina .and.
block%x(ixmhi1,1,1,1)<rmaxa .and. &
469 block%x(1,ixmlo2,1,2)>rminb .and.
block%x(1,ixmhi2,1,2)<rmaxb)
then
472 allocate(xcoarse(wnc,wnc,2))
477 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(1))
478 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(2))
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with basic grid data structures.
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Module with geometry-related routines (e.g., divergence, curl)
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer, parameter plevel_
double precision, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer nghostcells
Number of ghost cells surrounding a grid.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
integer, dimension(:,:), allocatable node
Program to extrapolate linear force-free fields in 3D Cartesian coordinates, based on exact Green fun...
double precision, save darea
double precision, dimension(:,:), allocatable, save bz0
subroutine get_potential_field_potential(ixIL, ixOL, potential, x, zshift)
subroutine multigrid_bc(box, nc, iv, nb, bc_type, bc)
To set boundary condition on physical boundaries for mg Poisson solver.
subroutine get_potential_field_potential_mg()
Solve Poisson equation of scalar potential using multigrid solver.
subroutine potential_field_energy_mg(benergy)
get potential magnetic field energy given normal B on all boundaries
subroutine init_b_fff_data(magnetogramname, qLunit, qBunit)
double precision, dimension(:), allocatable, save xa2
subroutine calc_lin_fff(ixIL, ixOL, Bf, x, alpha, zshift, idir)
double precision, save bzmax
subroutine get_potential_field_potential_sphere(ixIL, x, potential, nth, nph, magnetogram, theta, phi, r_sphere)
double precision, dimension(:), allocatable, save xa1
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.