xref: /libCEED/tests/t302-basis-f.f90 (revision e508ec0b3d519a4b165a0015a5ee248d90be1508)
18980d4a7Sjeremylt!-----------------------------------------------------------------------
28980d4a7Sjeremylt      program test
31f9a83abSJed Brown      implicit none
4ec3da8bcSJed Brown      include 'ceed/fortran.h'
58980d4a7Sjeremylt
61f9a83abSJed Brown      integer ceed,err,i,j
7*9ac7b42eSJeremy L Thompson      integer b,p
8*9ac7b42eSJeremy L Thompson      parameter(p=4)
9*9ac7b42eSJeremy L Thompson      real*8 collograd1d(36),x2(6)
10*9ac7b42eSJeremy L Thompson      real*8 grad1d(16),qref(6)
11*9ac7b42eSJeremy L Thompson      integer*8 gradoffset,qoffset
12*9ac7b42eSJeremy L Thompson      real*8 sum
138980d4a7Sjeremylt
148980d4a7Sjeremylt      character arg*32
158980d4a7Sjeremylt
168980d4a7Sjeremylt      call getarg(1,arg)
1752bfb9bbSJeremy L Thompson
188980d4a7Sjeremylt      call ceedinit(trim(arg)//char(0),ceed,err)
198980d4a7Sjeremylt
2052bfb9bbSJeremy L Thompson!     Already collocated, GetCollocatedGrad will return grad1d
21*9ac7b42eSJeremy L Thompson      call ceedbasiscreatetensorh1lagrange(ceed,1,1,p,p,ceed_gauss_lobatto,b,&
22c8b9fe72Sjeremylt     & err)
23b300f92cSjeremylt      call ceedbasisgetcollocatedgrad(b,collograd1d,err)
24*9ac7b42eSJeremy L Thompson      call ceedbasisgetgrad1d(b,grad1d,gradoffset,err)
25*9ac7b42eSJeremy L Thompson      do i=0,p-1
26*9ac7b42eSJeremy L Thompson        do j=1,p
27*9ac7b42eSJeremy L Thompson          if (abs(collograd1d(j+p*i)-grad1d(j+p*i+gradoffset))>1.0D-13) then
28*9ac7b42eSJeremy L Thompson! LCOV_EXCL_START
29*9ac7b42eSJeremy L Thompson            write(*,*) 'Error in collocated gradient ',collograd1d(j+p*i),' != ',&
30*9ac7b42eSJeremy L Thompson     &       grad1d(j+p*i+gradoffset)
31*9ac7b42eSJeremy L Thompson! LCOV_EXCL_STOP
3252bfb9bbSJeremy L Thompson          endif
3352bfb9bbSJeremy L Thompson        enddo
3452bfb9bbSJeremy L Thompson      enddo
3552bfb9bbSJeremy L Thompson      call ceedbasisdestroy(b,err)
368980d4a7Sjeremylt
3752bfb9bbSJeremy L Thompson!     Q = P, not already collocated
38*9ac7b42eSJeremy L Thompson      call ceedbasiscreatetensorh1lagrange(ceed,1,1,p,p,ceed_gauss,b,err)
39b300f92cSjeremylt      call ceedbasisgetcollocatedgrad(b,collograd1d,err)
40*9ac7b42eSJeremy L Thompson
41*9ac7b42eSJeremy L Thompson      call ceedbasisgetqref(b,qref,qoffset,err)
42*9ac7b42eSJeremy L Thompson      do i=1,p
43*9ac7b42eSJeremy L Thompson        x2(i)=qref(i+qoffset)*qref(i+qoffset)
44*9ac7b42eSJeremy L Thompson      enddo
45*9ac7b42eSJeremy L Thompson
46*9ac7b42eSJeremy L Thompson      do i=0,p-1
47*9ac7b42eSJeremy L Thompson        sum=0
48*9ac7b42eSJeremy L Thompson        do j=1,p
49*9ac7b42eSJeremy L Thompson            sum=sum+collograd1d(j+p*i)*x2(j)
50*9ac7b42eSJeremy L Thompson        enddo
51*9ac7b42eSJeremy L Thompson        if (abs(sum-2*qref(i+1+qoffset))>1.0D-13) then
52a2546046Sjeremylt! LCOV_EXCL_START
53*9ac7b42eSJeremy L Thompson            write(*,*) 'Error in collocated gradient ',sum,' != ',&
54*9ac7b42eSJeremy L Thompson            &       2*qref(i+1+qoffset)
55de996c55Sjeremylt! LCOV_EXCL_STOP
568980d4a7Sjeremylt        endif
578980d4a7Sjeremylt      enddo
5852bfb9bbSJeremy L Thompson      call ceedbasisdestroy(b,err)
598980d4a7Sjeremylt
6052bfb9bbSJeremy L Thompson!     Q = P + 2, not already collocated
61*9ac7b42eSJeremy L Thompson      call ceedbasiscreatetensorh1lagrange(ceed,1,1,p,p+2,ceed_gauss,b,err)
62*9ac7b42eSJeremy L Thompson      call ceedbasisgetcollocatedgrad(b,collograd1d,err)
63*9ac7b42eSJeremy L Thompson
64*9ac7b42eSJeremy L Thompson      call ceedbasisgetqref(b,qref,qoffset,err)
65*9ac7b42eSJeremy L Thompson      do i=1,p+2
66*9ac7b42eSJeremy L Thompson        x2(i)=qref(i+qoffset)*qref(i+qoffset)
67*9ac7b42eSJeremy L Thompson      enddo
68*9ac7b42eSJeremy L Thompson
69*9ac7b42eSJeremy L Thompson      do i=0,p+1
70*9ac7b42eSJeremy L Thompson        sum=0
71*9ac7b42eSJeremy L Thompson        do j=1,p+2
72*9ac7b42eSJeremy L Thompson            sum=sum+collograd1d(j+(p+2)*i)*x2(j)
73*9ac7b42eSJeremy L Thompson        enddo
74*9ac7b42eSJeremy L Thompson        if (abs(sum-2*qref(i+1+qoffset))>1.0D-13) then
7552bfb9bbSJeremy L Thompson! LCOV_EXCL_START
76*9ac7b42eSJeremy L Thompson            write(*,*) 'Error in collocated gradient ',sum,' != ',&
77*9ac7b42eSJeremy L Thompson     &       2*qref(i+1+qoffset)
7852bfb9bbSJeremy L Thompson! LCOV_EXCL_STOP
7952bfb9bbSJeremy L Thompson        endif
8052bfb9bbSJeremy L Thompson      enddo
8152bfb9bbSJeremy L Thompson      call ceedbasisdestroy(b,err)
8252bfb9bbSJeremy L Thompson
838980d4a7Sjeremylt      call ceeddestroy(ceed,err)
8452bfb9bbSJeremy L Thompson
858980d4a7Sjeremylt      end
868980d4a7Sjeremylt!-----------------------------------------------------------------------