38 double complex,
allocatable :: flm(:,:),Alm(:,:),Blm(:,:)
39 double precision,
allocatable :: Rlm(:,:), xrg(:)
40 double precision,
public ::
r_s=2.5d0,
r_0=1.d0
41 integer,
allocatable :: lmaxarray(:)
43 logical,
public ::
trunc=.false.
56 double precision,
allocatable :: b_r0(:,:)
57 double precision,
allocatable :: theta(:),phi(:),cfwm(:)
58 double precision :: rsl,xrl,dxr
59 integer :: xm,ym,
l,m,amode,file_handle,il,ir,nlarr,nsh
60 integer,
dimension(MPI_STATUS_SIZE) :: statuss
62 character(len=*) :: mapname
63 character(len=80) :: fharmcoef
65 fharmcoef=mapname//
'.coef'
66 inquire(file=fharmcoef, exist=aexist)
69 call mpi_file_open(mpi_comm_self,fharmcoef,mpi_mode_rdonly, &
70 mpi_info_null,file_handle,
ierrmpi)
71 call mpi_file_read(file_handle,
lmax,1,mpi_integer,statuss,
ierrmpi)
73 call mpi_file_read(file_handle,flm,(
lmax+1)*(
lmax+1),&
74 mpi_double_complex,statuss,
ierrmpi)
75 call mpi_file_close(file_handle,
ierrmpi)
85 inquire(file=mapname,exist=aexist)
87 if(
mype==0)
write(*,
'(2a)')
"can not find file:",mapname
88 call mpistop(
"no input magnetogram found")
90 call mpi_file_open(mpi_comm_self,mapname,mpi_mode_rdonly,mpi_info_null,&
92 call mpi_file_read(file_handle,xm,1,mpi_integer,statuss,
ierrmpi)
93 call mpi_file_read(file_handle,ym,1,mpi_integer,statuss,
ierrmpi)
99 call mpi_file_read(file_handle,theta,ym,mpi_double_precision,&
101 call mpi_file_read(file_handle,phi,xm,mpi_double_precision,&
103 call mpi_file_read(file_handle,b_r0,xm*ym,mpi_double_precision,&
105 call mpi_file_close(file_handle,
ierrmpi)
106 print*,
'nphi,ntheta',xm,ym
107 print*,
'theta range:',minval(theta),maxval(theta)
108 print*,
'phi range:',minval(phi),maxval(phi)
109 print*,
'Brmax,Brmin',maxval(b_r0),minval(b_r0)
113 call coef(b_r0,xm,ym,dcos(theta),dsin(theta),cfwm)
117 amode=ior(mpi_mode_create,mpi_mode_wronly)
118 call mpi_file_open(mpi_comm_self,fharmcoef,amode, &
119 mpi_info_null,file_handle,
ierrmpi)
120 call mpi_file_write(file_handle,
lmax,1,mpi_integer,statuss,
ierrmpi)
121 call mpi_file_write(file_handle,flm,(
lmax+1)*(
lmax+1),&
122 mpi_double_complex,statuss,
ierrmpi)
123 call mpi_file_close(file_handle,
ierrmpi)
129 if(
npe>1)
call mpi_bcast(flm,(
lmax+1)*(
lmax+1),mpi_double_complex,0,&
134 allocate(lmaxarray(nlarr))
138 dxr=(
r_s-
r_0)/dble(nlarr-1)
140 xrg(ir)=dxr*dble(ir-1)+
r_0
159 rlm(
l,m)=dsqrt(dble(
l**2-m**2)/dble(4*
l**2-1))
160 blm(
l,m)=-flm(
l,m)/(1.d0+dble(
l)+dble(
l)*rsl)
161 alm(
l,m)=-rsl*blm(
l,m)
170 integer,
intent(in) :: ym
171 double precision,
intent(in) :: miu(ym)
172 double precision,
intent(out) :: cfwm(ym)
174 double precision,
dimension(ym) :: Pl,Pm2,Pm1,Pprime,sintheta
175 double precision :: lr
178 sintheta=dsqrt(1.d0-miu**2)
185 pl=(2.d0-lr)*pm1*miu-(1.d0-lr)*pm2
190 pprime=(dble(ym)*pl)/sintheta**2
191 cfwm=2.d0/(sintheta*pprime)**2
196 subroutine coef(b_r0,xm,ym,miu,mius,cfwm)
199 integer,
intent(in) :: xm,ym
200 double precision,
intent(in) :: b_r0(xm,ym),cfwm(ym),miu(ym),mius(ym)
202 double complex :: Bm(0:xm-1,0:ym-1)
203 double precision,
dimension(xm) :: fftmr,fftmi
204 double precision,
dimension(0:lmax) :: N_mm
205 double precision,
dimension(ym) :: P_lm1,P_lm2,old_Pmm,P_l
206 double precision :: mr,lr,c1,c2
207 integer :: l,m,i,j,stat
211 fftmr=b_r0(:,i)/dble(xm)
213 call fft(fftmr,fftmi,xm,xm,xm,-1)
214 bm(:,i-1)=(fftmr+(0.d0,1.d0)*fftmi)
216 n_mm(0)=1.d0/dsqrt(4.d0*dpi)
218 n_mm(m)=-n_mm(m-1)*dsqrt(1.d0+1.d0/dble(2*m))
222 p_lm1=p_lm2*miu*dsqrt(3.d0)
224 flm(0,0)=sum(bm(0,:)*p_lm2*cfwm)
226 flm(1,0)=sum(bm(0,:)*p_lm1*cfwm)
229 c1=dsqrt(4.d0-1.d0/lr**2)
230 c2=-(1.d0-1.d0/lr)*dsqrt((2.d0*lr+1.d0)/(2.d0*lr-3.d0))
231 p_l=c1*miu*p_lm1+c2*p_lm2
233 flm(l,0)=sum(bm(0,:)*p_l*cfwm)
244 p_lm2=old_pmm*mius*n_mm(m)/n_mm(m-1)
245 p_lm1=p_lm2*miu*dsqrt(dble(2*m+3))
249 flm(m,m)=sum(bm(m,:)*p_lm2*cfwm)
251 if(m<
lmax) flm(m+1,m)=sum(bm(m,:)*p_lm1*cfwm)
255 c1=dsqrt((4.d0*lr**2-1.d0)/(lr**2-mr**2))
256 c2=-dsqrt(((2.d0*lr+1.d0)*((lr-1.d0)**2-mr**2))/((2.d0*lr-3.d0)*(lr**2-&
258 p_l=c1*miu*p_lm1+c2*p_lm2
259 flm(l,m)=sum(bm(m,:)*p_l*cfwm)
267 subroutine pfss(ixI^L,ixO^L,Bpf,x)
270 integer,
intent(in) :: ixi^
l,ixo^
l
271 double precision,
intent(in) :: x(ixi^s,1:
ndim)
272 double precision,
intent(out) :: bpf(ixi^s,1:
ndir)
274 double complex :: bt(0:
lmax,0:
lmax,ixomin1:ixomax1)
275 double precision :: phase(ixi^s,1:
ndir),bpfiv(ixomin3:ixomax3,ixomin2:ixomax2)
276 double precision :: miu(ixomin2:ixomax2),mius(ixomin2:ixomax2),xr
277 integer ::
l,m,ix^
d,j,l1,l2,ntheta,nphi,ir,qlmax
280 nphi=ixomax3-ixomin3+1
281 ntheta=ixomax2-ixomin2+1
282 miu(ixomin2:ixomax2)=dcos(x(ixomin1,ixomax2:ixomin2:-1,ixomin3,2))
283 mius(ixomin2:ixomax2)=dsin(x(ixomin1,ixomax2:ixomin2:-1,ixomin3,2))
284 do ix1=ixomin1,ixomax1
285 xr=x(ix1,ixomin2,ixomin3,1)
287 do ir=1,
size(lmaxarray)
290 if(ir>
size(lmaxarray)) ir=
size(lmaxarray)
298 bt(
l,m,ix1)=alm(
l,m)*dble(
l)*xr**(
l-1)-blm(
l,m)*dble(
l+1)*xr**(-
l-2)
302 ixomin3:ixomax3,3),miu,mius,nphi,ntheta,bpfiv,qlmax)
303 do ix3=ixomin3,ixomax3
304 do ix2=ixomin2,ixomax2
305 bpf(ix1,ix2,ix3,1)=bpfiv(ix3,ixomax2-ix2+ixomin2)
312 bt(
l,m,ix1)=-rlm(
l+1,m)*dble(
l+2)*&
313 (alm(
l+1,m)*xr**
l+blm(
l+1,m)*xr**(-
l-3))
314 else if (
l>=1 .and.
l<=
lmax-1)
then
315 bt(
l,m,ix1)=rlm(
l,m)*&
316 dble(
l-1)*(alm(
l-1,m)*xr**(
l-2)+blm(
l-1,m)*&
317 xr**(-
l-1))-rlm(
l+1,m)*dble(
l+2)*&
318 (alm(
l+1,m)*xr**
l+blm(
l+1,m)*xr**(-
l-3))
320 bt(
l,m,ix1)=rlm(
l,m)*&
321 dble(
l-1)*(alm(
l-1,m)*xr**(
l-2)+blm(
l-1,m)*xr**(-
l-1))
326 ixomin3:ixomax3,3),miu,mius,nphi,ntheta,bpfiv,qlmax)
327 do ix3=ixomin3,ixomax3
328 do ix2=ixomin2,ixomax2
329 bpf(ix1,ix2,ix3,2)=bpfiv(ix3,ixomax2-ix2+ixomin2)/mius(&
337 bt(
l,m,ix1)=(0.d0,1.d0)*m*(alm(
l,m)*xr**(
l-1)+blm(
l,m)*xr**(-
l-2))
341 ixomin3:ixomax3,3),miu,mius,nphi,ntheta,bpfiv,qlmax)
342 do ix3=ixomin3,ixomax3
343 do ix2=ixomin2,ixomax2
344 bpf(ix1,ix2,ix3,3)=bpfiv(ix3,ixomax2-ix2+ixomin2)/mius(&
364 integer,
intent(in) :: nphi,ntheta,qlmax
365 double complex,
intent(in) :: Bt(0:lmax,0:lmax)
366 double precision,
intent(in) :: phi(nphi),miu(ntheta),mius(ntheta)
367 double precision,
intent(out) :: Bpf(nphi,ntheta)
369 double precision,
dimension(0:lmax,0:lmax) :: cp,phase,Bamp
370 double precision,
dimension(ntheta) :: cp_1_0,cp_l_0,cp_lm1_0,cp_lm2_0,cp_m_m
371 double precision,
dimension(ntheta) :: cp_1_m,cp_l_m,cp_lm1_m,cp_lm2_m,cp_mp1_m
372 double precision :: angpart(nphi)
373 double precision :: ld,md,c1,c2,cp_0_0
374 integer :: l,m,iph,ith
378 phase=atan2(dimag(bt),dble(bt))
382 cp_0_0=dsqrt(1.d0/(4.d0*dpi))
384 bpf=bpf+bamp(0,0)*dcos(phase(0,0))*cp_0_0
387 cp_1_0=dsqrt(3.d0)*miu*cp_0_0
389 bpf(iph,:)=bpf(iph,:)+bamp(1,0)*dcos(phase(1,0))*cp_1_0
399 c1=dsqrt(4.d0*ld**2-1.d0)/ld
400 c2=dsqrt((2.d0*ld+1.d0)/(2.d0*ld-3.d0))*(ld-1.d0)/ld
401 cp_l_0=c1*miu*cp_lm1_0-c2*cp_lm2_0
403 bpf(iph,:)=bpf(iph,:)+bamp(l,0)*dcos(phase(l,0))*cp_l_0
412 cp_m_m=-dsqrt(1.d0+1.d0/(2.d0*md))*mius*cp_m_m
414 angpart(iph)=dcos(md*phi(iph)+phase(m,m))
418 bpf(iph,ith)=bpf(iph,ith)+bamp(m,m)*angpart(iph)*cp_m_m(ith)
424 cp_mp1_m=dsqrt(2.d0*md+3.d0)*miu*cp_m_m
425 angpart=dcos(md*phi+phase(m+1,m))
428 bpf(iph,ith)=bpf(iph,ith)+bamp(m+1,m)*angpart(iph)*cp_mp1_m(ith)
441 c1=dsqrt((4.d0*ld**2-1.d0)/(ld**2-md**2))
442 c2=dsqrt((2.d0*ld+1.d0)*((ld-1.d0)**2-md**2)/(2.d0*ld-3.d0)/(ld**2-md**2))
443 cp_l_m=c1*miu*cp_lm1_m-c2*cp_lm2_m
444 angpart=dcos(md*phi+phase(l,m))
447 bpf(iph,ith)=bpf(iph,ith)+bamp(l,m)*angpart(iph)*cp_l_m(ith)
456 subroutine fft(a,b,ntot,n,nspan,isn)
500 double precision :: a(:),b(:)
504 dimension nfac(11),np(209)
506 dimension at(23),ck(23),bt(23),sk(23)
507 double precision :: c72,s72,s120,rad,radf,sd,cd,ak,bk,c1
508 double precision :: s1,aj,bj,akp,ajp,ajm,akm,bkp,bkm,bjp,bjm,aa
509 double precision :: bb,sk,ck,at,bt,s3,c3,s2,c2
510 integer :: i,ii,maxp,maxf,n,inc,isn,nt,ntot,ks,nspan,kspan,nn,jc,jf,m
511 integer :: k,j,jj,nfac,kt,np,kk,k1,k2,k3,k4,kspnn
525 c72=0.30901699437494742d0
526 s72=0.95105651629515357d0
527 s120=0.86602540378443865d0
528 rad=6.2831853071796d0
529 if(isn .ge. 0)
go to 10
539 radf=rad*dble(jc)*0.5d0
549 20
if(k-(k/16)*16 .eq. 0)
go to 15
556 30
if(mod(k,jj) .eq. 0)
go to 25
559 if(jj .le. k)
go to 30
560 if(k .gt. 4)
go to 40
565 40
if(k-(k/4)*4 .ne. 0)
go to 50
571 60
if(mod(k,j) .ne. 0)
go to 70
576 if(j .le. k)
go to 60
577 80
if(kt .eq. 0)
go to 100
582 if(j .ne. 0)
go to 90
584 100 sd=radf/dble(kspan)
589 if(nfac(i) .ne. 2)
go to 400
601 if(kk .le. nn)
go to 210
603 if(kk .le. jc)
go to 210
604 if(kk .gt. kspan)
go to 800
615 if(kk .lt. nt)
go to 230
619 if(kk .gt. k2)
go to 230
622 c1=2.d0-(ak**2+s1**2)
626 if(kk .lt. k2)
go to 230
629 if(kk .le. jc+jc)
go to 220
642 aj=(a(k1)-a(k2))*s120
643 bj=(b(k1)-b(k2))*s120
649 if(kk .lt. nn)
go to 320
651 if(kk .le. kspan)
go to 320
654 400
if(nfac(i) .ne. 4)
go to 600
674 if(isn .lt. 0)
go to 450
679 if(s1 .eq. 0.d0)
go to 460
680 430 a(k1)=akp*c1-bkp*s1
687 if(kk .le. nt)
go to 420
688 440 c2=c1-(cd*c1+sd*s1)
690 c1=2.d0-(c2**2+s1**2)
698 if(kk .le. kspan)
go to 420
700 if(kk .le. jc)
go to 410
701 if(kspan .eq. jc)
go to 800
707 if(s1 .ne. 0)
go to 430
715 if(kk .le. nt)
go to 420
753 if(kk .lt. nn)
go to 520
755 if(kk .le. kspan)
go to 520
761 if(k .eq. 3)
go to 320
762 if(k .eq. 5)
go to 510
763 if(k .eq. jf)
go to 640
768 if(jf .gt. maxf)
go to 998
772 630 ck(j)=ck(k)*c1+sk(k)*s1
773 sk(j)=ck(k)*s1-sk(k)*c1
778 if(j .lt. k)
go to 630
797 if(k1 .lt. k2)
go to 650
818 if(jj .gt. jf) jj=jj-jf
819 if(k .lt. jf)
go to 670
826 if(j .lt. k)
go to 660
828 if(kk .le. nn)
go to 640
830 if(kk .le. kspan)
go to 640
832 700
if(i .eq. m)
go to 800
843 if(kk .le. nt)
go to 730
848 if(kk .le. kspnn)
go to 730
851 c1=2.d0-(c2**2+s1**2)
855 if(kk .le. kspan)
go to 720
857 if(kk .le. jc+jc)
go to 710
862 if(kt .eq. 0)
go to 890
867 810 np(j+1)=np(j)/nfac(j)
868 np(k)=np(k+1)*nfac(j)
871 if(j .lt. k)
go to 810
877 if(n .ne. ntot)
go to 850
887 if(k2 .lt. ks)
go to 820
891 if(k2 .gt. np(j))
go to 830
893 840
if(kk .lt. k2)
go to 820
896 if(k2 .lt. ks)
go to 840
897 if(kk .lt. ks)
go to 830
910 if(kk .lt. k)
go to 860
913 if(kk .lt. nt)
go to 850
916 if(k2 .lt. ks)
go to 850
920 if(k2 .gt. np(j))
go to 870
922 880
if(kk .lt. k2)
go to 850
925 if(k2 .lt. ks)
go to 880
926 if(kk .lt. ks)
go to 870
928 890
if(2*kt+1 .ge. m)
return
933 900 nfac(j)=nfac(j)*nfac(j+1)
935 if(j .ne. kt)
go to 900
938 if(nn .gt. maxp)
go to 998
947 if(jj .ge. k2)
go to 902
953 if(j .le. nn)
go to 904
960 if(kk .ne. j)
go to 910
964 if(kk .lt. 0)
go to 914
965 if(kk .ne. j)
go to 910
967 if(j .ne. nn)
go to 914
972 if(np(j) .lt. 0)
go to 924
975 if(jj .gt. maxf) kspan=maxf
985 if(k1 .ne. kk)
go to 928
993 if(k1 .ne. kk)
go to 936
995 if(k .ne. j)
go to 932
1002 if(k1 .ne. kk)
go to 940
1003 if(jj .ne. 0)
go to 926
1004 if(j .ne. 1)
go to 924
1008 if(nt .ge. 0)
go to 924
1014 999
format(44h0array bounds exceeded within
subroutine fft)
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer icomm
The MPI communicator.
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 ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
module mod_pfss.t – potential field source surface model PURPOSE : to extrapolate global potential ma...
subroutine fft(a, b, ntot, n, nspan, isn)
subroutine coef(b_r0, xm, ym, miu, mius, cfwm)
subroutine, public pfss(ixIL, ixOL, Bpf, x)
subroutine inv_sph_transform(Bt, phi, miu, mius, nphi, ntheta, Bpf, qlmax)
double precision, public r_0
subroutine cfweights(ym, miu, cfwm)
subroutine, public harm_coef(mapname)
double precision, public r_s