xref: /libCEED/tests/t302-basis-f.f90 (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1!-----------------------------------------------------------------------
2      program test
3      implicit none
4      include 'ceed/fortran.h'
5
6      integer ceed,err,i,j
7      integer b,p
8      parameter(p=4)
9      real*8 collograd1d(36),x2(6)
10      real*8 grad1d(16),qref(6)
11      integer*8 gradoffset,qoffset
12      real*8 sum
13
14      character arg*32
15
16      call getarg(1,arg)
17
18      call ceedinit(trim(arg)//char(0),ceed,err)
19
20!     Already collocated, GetCollocatedGrad will return grad1d
21      call ceedbasiscreatetensorh1lagrange(ceed,1,1,p,p,ceed_gauss_lobatto,b,&
22     & err)
23      call ceedbasisgetcollocatedgrad(b,collograd1d,err)
24      call ceedbasisgetgrad1d(b,grad1d,gradoffset,err)
25      do i=0,p-1
26        do j=1,p
27          if (abs(collograd1d(j+p*i)-grad1d(j+p*i+gradoffset))>1.0D-13) then
28! LCOV_EXCL_START
29            write(*,*) 'Error in collocated gradient ',collograd1d(j+p*i),' != ',&
30     &       grad1d(j+p*i+gradoffset)
31! LCOV_EXCL_STOP
32          endif
33        enddo
34      enddo
35      call ceedbasisdestroy(b,err)
36
37!     Q = P, not already collocated
38      call ceedbasiscreatetensorh1lagrange(ceed,1,1,p,p,ceed_gauss,b,err)
39      call ceedbasisgetcollocatedgrad(b,collograd1d,err)
40
41      call ceedbasisgetqref(b,qref,qoffset,err)
42      do i=1,p
43        x2(i)=qref(i+qoffset)*qref(i+qoffset)
44      enddo
45
46      do i=0,p-1
47        sum=0
48        do j=1,p
49            sum=sum+collograd1d(j+p*i)*x2(j)
50        enddo
51        if (abs(sum-2*qref(i+1+qoffset))>1.0D-13) then
52! LCOV_EXCL_START
53            write(*,*) 'Error in collocated gradient ',sum,' != ',&
54            &       2*qref(i+1+qoffset)
55! LCOV_EXCL_STOP
56        endif
57      enddo
58      call ceedbasisdestroy(b,err)
59
60!     Q = P + 2, not already collocated
61      call ceedbasiscreatetensorh1lagrange(ceed,1,1,p,p+2,ceed_gauss,b,err)
62      call ceedbasisgetcollocatedgrad(b,collograd1d,err)
63
64      call ceedbasisgetqref(b,qref,qoffset,err)
65      do i=1,p+2
66        x2(i)=qref(i+qoffset)*qref(i+qoffset)
67      enddo
68
69      do i=0,p+1
70        sum=0
71        do j=1,p+2
72            sum=sum+collograd1d(j+(p+2)*i)*x2(j)
73        enddo
74        if (abs(sum-2*qref(i+1+qoffset))>1.0D-13) then
75! LCOV_EXCL_START
76            write(*,*) 'Error in collocated gradient ',sum,' != ',&
77     &       2*qref(i+1+qoffset)
78! LCOV_EXCL_STOP
79        endif
80      enddo
81      call ceedbasisdestroy(b,err)
82
83      call ceeddestroy(ceed,err)
84
85      end
86!-----------------------------------------------------------------------