c*********************** problem name: ident ************************* c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine a1xy(x,y,u,ux,uy,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values character(len=80) :: su common /val0/k0,ku,kx,ky,kl,kuu,kux,kuy,kul, + kxu,kxx,kxy,kxl,kyu,kyx,kyy,kyl,klu,klx,kly,kll common /atest2/irl,iu(99),a,c0,c1,c2,d,ru(95),su(100) cy call setrl(rl) values(k0)=(a**2+1.0e0_rknd)*ux values(kx)=a**2+1.0e0_rknd if(irl==1) then values(kl)=ux*a*2.0e0_rknd values(klx)=a*2.0e0_rknd values(kll)=ux*2.0e0_rknd endif return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine a2xy(x,y,u,ux,uy,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values character(len=80) :: su common /val0/k0,ku,kx,ky,kl,kuu,kux,kuy,kul, + kxu,kxx,kxy,kxl,kyu,kyx,kyy,kyl,klu,klx,kly,kll common /atest2/irl,iu(99),a,c0,c1,c2,d,ru(95),su(100) cy call setrl(rl) values(k0)=(a**2+1.0e0_rknd)*uy values(ky)=a**2+1.0e0_rknd if(irl==1) then values(kl)=uy*a*2.0e0_rknd values(kly)=a*2.0e0_rknd values(kll)=uy*2.0e0_rknd endif return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine fxy(x,y,u,ux,uy,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values character(len=80) :: su common /val0/k0,ku,kx,ky,kl,kuu,kux,kuy,kul, + kxu,kxx,kxy,kxl,kyu,kyx,kyy,kyl,klu,klx,kly,kll common /atest2/irl,iu(99),a,c0,c1,c2,d,ru(95),su(100) cy call setrl(rl) values(k0)=c2*u**2+c1*u-c0 values(ku)=c2*u*2.0e0_rknd+c1 values(kuu)=c2*2.0e0_rknd if(irl==2) then values(kl)=-1.0e0_rknd elseif(irl==3) then values(kl)=u values(kul)=1.0e0_rknd elseif(irl==4) then values(kl)=u**2 values(kul)=2.0e0_rknd*u endif return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine gnxy(x,y,u,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values character(len=80) :: su common /val1/k0,ku,kl,kuu,kul,klu,kll common /atest2/irl,iu(99),a,c0,c1,c2,d,ru(95),su(100) cy call setrl(rl) c return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine gdxy(x,y,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values character(len=80) :: su common /val2/k0,kl,kll,klb,kub,kic,kim,kil common /atest2/irl,iu(99),a,c0,c1,c2,d,ru(95),su(100) cy call setrl(rl) if(itag==2) then values(k0)=d if(irl==5) values(kl)=1.0e0_rknd endif c return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine p1xy(x,y,u,ux,uy,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values common /val0/k0,ku,kx,ky,kl,kuu,kux,kuy,kul, + kxu,kxx,kxy,kxl,kyu,kyx,kyy,kyl,klu,klx,kly,kll cy uu=1.0e0_rknd e=exp(-20.0e0_rknd*(x**2+y**2)) values(k0)=e*(u-uu)**2 values(ku)=e*2.0e0_rknd*(u-uu) values(kuu)=e*2.0e0_rknd return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine p2xy(x,y,dx,dy,u,ux,uy,rl,itag,jtag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values common /val0/k0,ku,kx,ky,kl,kuu,kux,kuy,kul, + kxu,kxx,kxy,kxl,kyu,kyx,kyy,kyl,klu,klx,kly,kll cy return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine qxy(x,y,u,ux,uy,rl,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values common /val3/kf,kf1,kf2,kad cy return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine sxy(rl,s,itag,values) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) real(kind=rknd), dimension(*) :: values common /val4/jx,jy,jxs,jys,jxl,jyl,jxss,jyss,jxll,jyll, + jxsl,jysl,jxls,jyls cy return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine setrl(rl) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) character(len=80) :: su common /atest2/irl,iu(99),a,c0,c1,c2,d,ru(95),su(100) cy if(irl>0) ru(irl)=rl return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine usrcmd(vx,vy,sf,itnode,ibndry,ip,rp,sp,iu,ru,su) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) integer(kind=iknd), dimension(5,*) :: itnode integer(kind=iknd), dimension(7,*) :: ibndry integer(kind=iknd), dimension(100) :: ip,iu integer(kind=iknd), save :: len real(kind=rknd), dimension(*) :: vx,vy real(kind=rknd), dimension(2,*) :: sf real(kind=rknd), dimension(100) :: rp,ru character(len=80), dimension(100) :: sp,su character(len=80), save, dimension(20) :: file cy data len/11/ data (file(i),i= 1, 10)/ + 'n i= 1,n=irl, a= i,t=i', 1 'n i= 1,n=a, a= a,t=r', 2 'n i= 2,n=c0, a=c0,t=r', 3 'n i= 3,n=c1, a=c1,t=r', 4 'n i= 4,n=c2, a=c2,t=r', 5 'n i= 5,n=d, a= d,t=r', 6 's n=irl ,v= 1,l="a"', 7 's n=irl ,v= 2,l="c0"', 8 's n=irl ,v= 3,l="c1"', 9 's n=irl ,v= 4,l="c2"'/ data (file(i),i= 11, 11)/ + 's n=irl ,v= 5,l="d"'/ c c c enter input mode c irlsv=iu(1) call usrset(file,len,iu,ru,su) c if(irlsv/=iu(1)) then ip(7)=8 rp(1)=ru(iu(1)) endif return end c----------------------------------------------------------------------- c c piecewise lagrange triangle multi grid package c c edition 11.0 - - - june, 2012 c c----------------------------------------------------------------------- subroutine gdata(vx,vy,sf,itnode,ibndry,ip,rp,sp,iu,ru,su,sxy) cx use mthdef implicit real(kind=rknd) (a-h,o-z) implicit integer(kind=iknd) (i-n) integer(kind=iknd), dimension(5,*) :: itnode integer(kind=iknd), dimension(7,*) :: ibndry integer(kind=iknd), dimension(100) :: ip,iu integer(kind=iknd), save, dimension(13) :: ix,iy integer(kind=iknd), save :: ntf,nvf,nbf,iprob,ispd integer(kind=iknd), save :: iadapt,ifirst real(kind=rknd), dimension(*) :: vx,vy real(kind=rknd), dimension(2,*) :: sf real(kind=rknd), dimension(100) :: rp,ru real(kind=rknd), save :: hmax,grade character(len=80), dimension(100) :: sp,su cy external sxy data ix/1,3,3,1,1,-1,-1,-3,-3,-1,-1,1,0/ data iy/-1,-1,1,1,3,3,1,1,-1,-1,-3,-3,0/ data nvf,ntf,nbf/13,4,16/ data iprob,ispd,iadapt,ifirst/4,1,5,1/ data hmax,grade/0.1e0_rknd,1.5e0_rknd/ c if(ip(41)==1) then sp(2)='ident' sp(1)='ident' sp(3)='ident' sp(6)='ident_mpixxx.rw' sp(7)='ident.jnl' sp(9)='ident_mpixxx.out' ip(7)=8 iu(1)=1 ru(1)=1.0e0_rknd ru(2)=1.0e0_rknd ru(3)=0.0e0_rknd ru(4)=1.0e0_rknd ru(5)=0.0e0_rknd rp(1)=ru(iu(1)) endif c c do i=1,nvf vx(i)=real(ix(i),rknd) vy(i)=real(iy(i),rknd) enddo do i=1,nbf ibndry(1,i)=i ibndry(2,i)=i+1 if(i>=nvf) then ibndry(1,i)=(i-nvf)*3+1 ibndry(2,i)=nvf endif ibndry(3,i)=0 ibndry(4,i)=1 if(i>=nvf) ibndry(4,i)=0 ibndry(5,i)=0 ibndry(6,i)=0 ibndry(7,i)=1 ii=i-(i/3)*3 if(ii==2.and.i<=nvf) then ibndry(4,i)=2 ibndry(7,i)=i endif sf(1,i)=0.0e0_rknd sf(2,i)=0.0e0_rknd enddo ibndry(2,nvf-1)=1 c ip(1)=ntf ip(2)=nvf ip(3)=nbf ip(5)=max(ip(5),ifirst) ip(6)=iprob ip(8)=ispd ip(20)=iadapt rp(15)=hmax rp(16)=grade c c make itnode, find symmetries c call sklutl(0_iknd,vx,vy,sf,itnode,ibndry,ip,rp,iflag,sxy) call sklutl(2_iknd,vx,vy,sf,itnode,ibndry,ip,rp,iflag,sxy) return end .