MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_lfff.t
Go to the documentation of this file.
1!> Program to extrapolate linear force-free fields in 3D Cartesian coordinates,
2!> based on exact Green function method (Chiu & Hilton 1977 ApJ 212,873).
3!>
4!> Usage:
5!> 1 In the subroutine usr_set_parameters of mod_usr.t:
6!> To extrapolate a linear force free field from a observed magnetogram
7!> prepared in a data file, e.g., 'hmiM720sxxxx.dat' replace
8!> call init_bc_fff_data('hmiM720sxxxx.dat',unit_length,unit_magneticfield)
9!> 'hmiM720sxxxx.dat' must be a binary file containing nx1,nx2,xc1,xc2,dxm1,
10!> dxm2, Bz0(nx1,nx2). Integers nx1 and nx2 give the resolution of the
11!> uniform-grid magentogram. Others are double-precision floats. xc1 and xc2
12!> are coordinates of the central point of the magnetogram. dxm1 and dxm2
13!> are the cell sizes for each direction, Bz0 is the vertical conponent
14!> of magetic field on the solar surface from observations.
15!>2 In the subroutine usr_init_one_grid of mod_usr.t,
16!> add lines like:
17!>
18!> double precision :: Bf(ixG^S,1:ndir), alpha, zshift
19!>
20!> alpha=0.d0 ! potential field
21!> !alpha=0.08d0 ! non-potential linear force-free field
22!> zshift=0.05d0 ! lift your box zshift heigher to the bottom magnetogram
23!> call calc_lin_fff(ixG^L,ix^L,Bf,x,alpha,zshift)
24!>
25!>3 Notice that the resolution of input magnetogram must be better than the best
26!> resolution of your AMR grid to have a good behavior close to the bottom layer
28 implicit none
29
30 double precision, save :: bzmax,darea
31 double precision, allocatable, save :: bz0(:,:)
32 double precision, allocatable, save :: xa1(:),xa2(:)
33 integer, save :: nx1,nx2
34
35contains
36
37 subroutine init_b_fff_data(magnetogramname,qLunit,qBunit)
39 use mod_comm_lib, only: mpistop
40
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
45 logical :: aexist
46 character(len=*), intent(in) :: magnetogramname
47 ! nx1,nx2 are numbers of cells for each direction
48 ! xc1,xc2 are coordinates of the central point of the magnetogram
49 ! dxm1,dxm2 are cell sizes for each direction
50 ! Bz0 is the 2D Bz magnetogram
51 inquire(file=magnetogramname,exist=aexist)
52 if(.not. aexist) then
53 if(mype==0) write(*,'(2a)') "can not find file:",magnetogramname
54 call mpistop("no input magnetogram----init_b_fff_data")
55 end if
56 call mpi_file_open(icomm,magnetogramname,mpi_mode_rdonly,mpi_info_null,&
57 file_handle,ierrmpi)
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)
60 allocate(bz0(nx1,nx2))
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,&
66 statuss,ierrmpi)
67 call mpi_file_close(file_handle,ierrmpi)
68 allocate(xa1(nx1))
69 allocate(xa2(nx2))
70 do i=1,nx1
71 xa1(i) = xc1 + (dble(i) - dble(nx1)/2.d0 - 0.5d0)*dxm1
72 enddo
73 do i=1,nx2
74 xa2(i) = xc2 + (dble(i) - dble(nx2)/2.d0 - 0.5d0)*dxm2
75 enddo
76 ! declare and define global variables Lunit and Bunit to be your length unit in
77 ! cm and magnetic strength unit in Gauss first
78 dxm1=dxm1/qlunit
79 dxm2=dxm2/qlunit
80 xa1=xa1/qlunit
81 xa2=xa2/qlunit
82 darea=dxm1*dxm2
83 bz0=bz0/qbunit
84 bzmax=maxval(dabs(bz0(:,:)))
85
86 ! normalize b
88 if(mype==0) then
89 print*,'magnetogram xrange:',minval(xa1),maxval(xa1)
90 print*,'magnetogram yrange:',minval(xa2),maxval(xa2)
91 end if
92
93 if(mype==0) then
94 print*,'extrapolating 3D force-free field from an observed Bz '
95 print*,'magnetogram of',nx1,'by',nx2,'pixels. Bzmax=',bzmax
96 endif
97
98 end subroutine init_b_fff_data
99
100{^ifthreed
101 subroutine calc_lin_fff(ixI^L,ixO^L,Bf,x,alpha,zshift,idir)
102 ! PURPOSE:
103 ! Calculation to determine linear FFF from the field on
104 ! the lower boundary (Chiu and Hilton 1977 ApJ 212,873).
105 ! NOTE: Only works for Cartesian coordinates
106 ! INPUT: Bf,x
107 ! OUTPUT: updated b in w
109
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)
114
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
118
119 bf=0.d0
120 twopiinv = 0.5d0/dpi*bzmax*darea
121 ! get cos and sin arrays
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))
125 ! looping Bz0 pixels
126 do ixp2=1,nx2
127 do ixp1=1,nx1
128 bigr(ixo^s)=dsqrt((x(ixo^s,1)-xa1(ixp1))**2+&
129 (x(ixo^s,2)-xa2(ixp2))**2)
130 r2=bigr**2+zk**2
131 r=dsqrt(r2)
132 r3=r**3
133 cos_ar=dcos(alpha*r)
134 sin_ar=dsin(alpha*r)
135 where(bigr/=0.d0)
136 bigr=1.d0/bigr
137 end where
138 where(r/=0.d0)
139 r=1.d0/r
140 end where
141 where(r2/=0.d0)
142 r2=1.d0/r2
143 end where
144 where(r3/=0.d0)
145 r3=1.d0/r3
146 end where
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
149 do idim=1,ndim
150 if(present(idir).and.idim/=idir) cycle
151 select case(idim)
152 case(1)
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)
155 case(2)
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)
158 case(3)
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))
161 end select
162 end do
163 end do
164 end do
165 bf(ixo^s,:)=bf(ixo^s,:)*twopiinv
166
167 end subroutine calc_lin_fff
168
169 subroutine get_potential_field_potential(ixI^L,ixO^L,potential,x,zshift)
170 ! PURPOSE:
171 ! Calculation scalar potential of potential field given
172 ! Bz at photosphere (Schmidt 1964 NASSP).
173 ! NOTE: Only works for Cartesian coordinates
174 ! INPUT: x,zshift
175 ! OUTPUT: potential
177
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)
181
182 double precision :: zk
183 integer :: ixp1,ixp2,ix^D
184
185 potential=0.d0
186 ! looping Bz0 pixels see equation (2)
187 !$OMP PARALLEL DO PRIVATE(bigr) REDUCTION(+:potential)
188 do ix3=ixomin3,ixomax3
189 zk=x(ixomin1,ixomin2,ix3,3)-xprobmin3+zshift
190 do ix2=ixomin2,ixomax2
191 do ix1=ixomin1,ixomax1
192 do ixp2=1,nx2
193 do ixp1=1,nx1
194 potential(ix^d)=potential(ix^d)+0.5d0*bz0(ixp1,ixp2)*darea/&
195 (dpi*dsqrt((x(ix^d,1)-xa1(ixp1))**2+(x(ix^d,2)-xa2(ixp2))**2+zk**2))
196 end do
197 end do
198 end do
199 end do
200 end do
201 !$OMP END PARALLEL DO
202 end subroutine get_potential_field_potential
203
204 subroutine get_potential_field_potential_sphere(ixI^L,x,potential,nth,nph,magnetogram,theta,phi,r_sphere)
205 ! PURPOSE:
206 ! Calculation scalar potential of potential field given
207 ! Bz at photosphere (Schmidt 1964 NASSP).
208 ! NOTE: Only works for spherical coordinates
209 ! OUTPUT: potential
211 integer, intent(in) :: ixI^L,nth,nph
212 real*8, intent(in) :: x(ixi^s,1:ndim)
213 ! magnetogram Br on photosphere
214 real*8, intent(in) :: magnetogram(nth,nph)
215 ! theta and phi grid of the photospheric magnetogram
216 real*8, intent(in) :: theta(nth),phi(nph)
217 ! radius of photosphere
218 real*8, intent(in) :: r_sphere
219 real*8, intent(out) :: potential(ixi^s)
220
221 real*8 :: area(nth),distance(ixi^s),dtheta_half,dphi,inv2pi
222 integer :: ix1,ix2
223
224 potential=0.d0
225 ! assume uniformly discretized theta and phi
226 dtheta_half=0.5d0*(theta(2)-theta(1))
227 dphi=phi(2)-phi(1)
228 area(1:nth)=2.d0*r_sphere**2*sin(theta(1:nth))*sin(dtheta_half)*sin(dphi)
229 inv2pi=-1.d0/(2.d0*dpi)
230
231 !$OMP PARALLEL DO PRIVATE(distance) REDUCTION(+:potential)
232 do ix2=1,nph,2
233 do ix1=1,nth,2
234 distance(ixi^s)=sqrt(x(ixi^s,1)**2+r_sphere**2-2.d0*x(ixi^s,1)*r_sphere*&
235 (sin(x(ixi^s,2))*sin(theta(ix1))*cos(phi(ix2)-x(ixi^s,3))+cos(x(ixi^s,2))*&
236 cos(theta(ix1))))
237 where(distance(ixi^s)/=0.d0)
238 distance(ixi^s)=1.d0/distance(ixi^s)
239 end where
240 potential(ixi^s)=potential(ixi^s)+inv2pi*magnetogram(ix1,ix2)*distance(ixi^s)*area(ix1)
241 end do
242 end do
243 !$OMP END PARALLEL DO
244
246
247 !> get potential magnetic field energy given normal B on all boundaries
248 subroutine potential_field_energy_mg(benergy)
249 ! NOTE: Solve Poisson equation of scalar potential using multigrid solver
250 ! Only works for 3D Cartesian coordinates
253 use mod_forest
254 use mod_geometry
255 real*8, intent(out) :: benergy
256 real*8 :: block_benergy(ixm^t), mype_benergy
257 real*8, allocatable :: tmp(:,:,:)
258 integer :: iigrid, igrid, idir
259 integer :: id,nc, lvl, ixI^L
260
262 iximin^d=ixmlo^d-1;
263 iximax^d=ixmhi^d+1;
264 allocate(tmp(ixi^s))
265 mype_benergy=0.d0
266 do iigrid=1,igridstail; igrid=igrids(iigrid);
267 block=>ps(igrid)
268 ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
269 id = igrid_to_node(igrid, mype)%node%id
270 lvl= mg%boxes(id)%lvl
271 nc = mg%box_size_lvl(lvl)
272 block_benergy=0.d0
273 do idir=1,ndim
274 call gradient(mg%boxes(id)%cc({0:nc+1},mg_iphi),ixi^l,ixm^ll,idir,tmp)
275 block_benergy(ixm^t)=block_benergy(ixm^t)+tmp(ixm^t)**2
276 end do
277 mype_benergy=mype_benergy+sum(0.5d0*block_benergy(ixm^t)*block%dvolume(ixm^t))
278 end do
279
280 call mpi_allreduce(mype_benergy,benergy,1,mpi_double_precision,mpi_sum,icomm,ierrmpi)
281
282 end subroutine potential_field_energy_mg
283
284 !> Solve Poisson equation of scalar potential using multigrid solver
287 use mod_comm_lib, only: mpistop
289 double precision :: max_residual, res
290 integer :: i, max_its
291
292 mg%operator_type = mg_laplacian
293 max_its=50
294 max_residual=1.d-8
295
296 ! Set boundary conditions
297 mg%bc(:, mg_iphi)%bc_type = mg_bc_neumann
298
299 do i=1,2*ndim
300 mg%bc(i, mg_iphi)%boundary_cond => multigrid_bc
301 end do
302
303 ! Solve laplacian(phi) = 0
304 do i=1,max_its
305 !call mg_fas_vcycle(mg, max_res=res)
306 call mg_fas_fmg(mg,.true.,max_res=res)
307 if(mype==0) write(*,*) 'mg iteration',i,'residual:', res
308 if(res < max_residual) exit
309 end do
310 if(res > max_residual) call mpistop("get potential multigrid: no convergence")
311
313
314 !> To set boundary condition on physical boundaries for mg Poisson solver
315 subroutine multigrid_bc(box, nc, iv, nb, bc_type, bc)
318 type(mg_box_t), intent(in) :: box
319 integer, intent(in) :: nc
320 integer, intent(in) :: iv !< Index of variable
321 integer, intent(in) :: nb !< number of boundary from 1 to 6 for 3D
322 integer, intent(out) :: bc_type !< Type of b.c.
323 ! mg boundary values
324 double precision, intent(out) :: bc(nc, nc)
325 ! mg coordinates of boundary layer
326 double precision :: rr(nc, nc, 3)
327 double precision :: rmina,rminb,rmaxa,rmaxb,xmina,xminb,xmaxa,xmaxb
328 ! store normal magnetic field on the boundary
329 double precision :: wbn(ixG^T)
330 double precision, allocatable :: xcoarse(:,:,:)
331 integer :: iigrid,igrid,ix^D,idir,ixbca,ixbcb,ixbcn,dlvl,wnc
332
333 bc_type = mg_bc_neumann
334
335 call mg_get_face_coords(box, nb, nc, rr)
336
337 idir=(nb+1)/2
338 bc=0.d0
339 ! send normal B on boundaries from AMRVAC to mg Poisson solver
340 ! for Neumann boundary condition
341 select case(idir)
342 case(1)
343 if(mod(nb,2)==0) then
344 ixbcn=ixmhi1
345 else
346 ixbcn=ixmlo1-1
347 end if
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);
353 block=>ps(igrid)
354 if(.not.block%is_physical_boundary(nb)) cycle
355 if(stagger_grid) then
356 !wbn(ixbcn^%1ixG^T)=block%ws(ixbcn^%1ixG^T,idir)
357 wbn(ixbcn^%1ixg^t)=block%ws(ixbcn^%1ixg^t,idir)
358 else
359 wbn(ixbcn^%1ixg^t)=half*(block%w(ixbcn^%1ixg^t,iw_mag(idir))+block%w(ixbcn+1^%1ixg^t,iw_mag(idir)))
360 end if
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
367 do ix2=1,nc
368 do ix1=1,nc
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)
372 end do
373 end do
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
376 dlvl=node(plevel_,igrid)-box%lvl
377 wnc=nc/2**dlvl
378 allocate(xcoarse(wnc,wnc,2))
379 do ix2=1,wnc
380 do ix1=1,wnc
381 xcoarse(ix1,ix2,1)=sum(block%x(1,(ix1-1)*2**dlvl+1+nghostcells:ix1*2**dlvl+nghostcells,1,2))/dble(2**dlvl)
382 xcoarse(ix1,ix2,2)=sum(block%x(1,1,(ix2-1)*2**dlvl+1+nghostcells:ix2*2**dlvl+nghostcells,3))/dble(2**dlvl)
383 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(2))
384 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(3))
385 bc(ixbca,ixbcb)=sum(wbn(ixbcn,(ix1-1)*2**dlvl+1+nghostcells:ix1*2**dlvl+nghostcells,&
386 (ix2-1)*2**dlvl+1+nghostcells:ix2*2**dlvl+nghostcells))/dble(2**(2*dlvl))
387 end do
388 end do
389 deallocate(xcoarse)
390 end if
391 end do
392 case(2)
393 if(mod(nb,2)==0) then
394 ixbcn=ixmhi2
395 else
396 ixbcn=ixmlo2-1
397 end if
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);
403 block=>ps(igrid)
404 if(.not.block%is_physical_boundary(nb)) cycle
405 if(stagger_grid) then
406 wbn(ixbcn^%2ixg^t)=block%ws(ixbcn^%2ixg^t,idir)
407 else
408 wbn(ixbcn^%2ixg^t)=half*(block%w(ixbcn^%2ixg^t,iw_mag(idir))+block%w(ixbcn+1^%2ixg^t,iw_mag(idir)))
409 end if
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
416 do ix2=1,nc
417 do ix1=1,nc
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)
421 end do
422 end do
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
425 dlvl=node(plevel_,igrid)-box%lvl
426 wnc=nc/2**dlvl
427 allocate(xcoarse(wnc,wnc,2))
428 do ix2=1,wnc
429 do ix1=1,wnc
430 xcoarse(ix1,ix2,1)=sum(block%x((ix1-1)*2**dlvl+1+nghostcells:ix1*2**dlvl+nghostcells,1,1,1))/dble(2**dlvl)
431 xcoarse(ix1,ix2,2)=sum(block%x(1,1,(ix2-1)*2**dlvl+1+nghostcells:ix2*2**dlvl+nghostcells,3))/dble(2**dlvl)
432 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(1))
433 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(3))
434 bc(ixbca,ixbcb)=sum(wbn((ix1-1)*2**dlvl+1+nghostcells:ix1*2**dlvl+nghostcells,ixbcn,&
435 (ix2-1)*2**dlvl+1+nghostcells:ix2*2**dlvl+nghostcells))/dble(2**(2*dlvl))
436 end do
437 end do
438 deallocate(xcoarse)
439 end if
440 end do
441 case(3)
442 if(mod(nb,2)==0) then
443 ixbcn=ixmhi3
444 else
445 ixbcn=ixmlo3-1
446 end if
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);
452 block=>ps(igrid)
453 if(.not.block%is_physical_boundary(nb)) cycle
454 if(stagger_grid) then
455 wbn(ixbcn^%3ixg^t)=block%ws(ixbcn^%3ixg^t,idir)
456 else
457 wbn(ixbcn^%3ixg^t)=half*(block%w(ixbcn^%3ixg^t,iw_mag(idir))+block%w(ixbcn+1^%3ixg^t,iw_mag(idir)))
458 end if
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
465 do ix2=1,nc
466 do ix1=1,nc
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)
470 end do
471 end do
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
474 dlvl=node(plevel_,igrid)-box%lvl
475 wnc=nc/2**dlvl
476 allocate(xcoarse(wnc,wnc,2))
477 do ix2=1,wnc
478 do ix1=1,wnc
479 xcoarse(ix1,ix2,1)=sum(block%x((ix1-1)*2**dlvl+1+nghostcells:ix1*2**dlvl+nghostcells,1,1,1))/dble(2**dlvl)
480 xcoarse(ix1,ix2,2)=sum(block%x(1,(ix2-1)*2**dlvl+1+nghostcells:ix2*2**dlvl+nghostcells,1,2))/dble(2**dlvl)
481 ixbca=ceiling((xcoarse(ix1,ix2,1)-rmina)/box%dr(1))
482 ixbcb=ceiling((xcoarse(ix1,ix2,2)-rminb)/box%dr(2))
483 bc(ixbca,ixbcb)=sum(wbn((ix1-1)*2**dlvl+1+nghostcells:ix1*2**dlvl+nghostcells,&
484 (ix2-1)*2**dlvl+1+nghostcells:ix2*2**dlvl+nghostcells,ixbcn))/dble(2**(2*dlvl))
485 end do
486 end do
487 deallocate(xcoarse)
488 end if
489 end do
490 end select
491
492 end subroutine multigrid_bc
493}
494end module mod_lfff
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
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
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
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.
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...
Definition mod_lfff.t:27
double precision, save darea
Definition mod_lfff.t:30
double precision, dimension(:,:), allocatable, save bz0
Definition mod_lfff.t:31
subroutine multigrid_bc(box, nc, iv, nb, bc_type, bc)
To set boundary condition on physical boundaries for mg Poisson solver.
Definition mod_lfff.t:316
integer, save nx1
Definition mod_lfff.t:33
subroutine get_potential_field_potential_sphere(ixil, x, potential, nth, nph, magnetogram, theta, phi, r_sphere)
Definition mod_lfff.t:205
subroutine init_b_fff_data(magnetogramname, qlunit, qbunit)
Definition mod_lfff.t:38
subroutine get_potential_field_potential_mg()
Solve Poisson equation of scalar potential using multigrid solver.
Definition mod_lfff.t:286
subroutine calc_lin_fff(ixil, ixol, bf, x, alpha, zshift, idir)
Definition mod_lfff.t:102
subroutine potential_field_energy_mg(benergy)
get potential magnetic field energy given normal B on all boundaries
Definition mod_lfff.t:249
integer, save nx2
Definition mod_lfff.t:33
double precision, dimension(:), allocatable, save xa2
Definition mod_lfff.t:32
double precision, save bzmax
Definition mod_lfff.t:30
subroutine get_potential_field_potential(ixil, ixol, potential, x, zshift)
Definition mod_lfff.t:170
double precision, dimension(:), allocatable, save xa1
Definition mod_lfff.t:32
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.