xref: /libCEED/tests/t323-basis-f.f90 (revision e735508cdfa08f2ead7f8e6e8f825d0fa851c858)
1!-----------------------------------------------------------------------
2!
3! Header with common subroutine
4!
5      include 't320-basis-f.h'
6!-----------------------------------------------------------------------
7      subroutine feval(x1,x2,val)
8      real*8 x1,x2,val
9
10      val=x1*x1+x2*x2+x1*x2+1
11
12      end
13!-----------------------------------------------------------------------
14      subroutine dfeval(x1,x2,val)
15      real*8 x1,x2,val
16
17      val=2*x1+x2
18
19      end
20!-----------------------------------------------------------------------
21      program test
22      implicit none
23      include 'ceed/fortran.h'
24
25      integer ceed,err
26      integer input,output
27      integer p,q,d
28      parameter(p=6)
29      parameter(q=4)
30      parameter(d=2)
31
32      real*8 qref(d*q)
33      real*8 qweight(q)
34      real*8 interp(p*q)
35      real*8 grad(d*p*q)
36      real*8 xq(d*q)
37      real*8 xr(d*p)
38      real*8 iinput(p)
39      real*8 ooutput(d*q)
40      real*8 val,diff
41      real*8 x1,x2
42      integer*8 ioffset,ooffset
43
44      integer b,i
45
46      character arg*32
47
48      xq=(/2.d-1,6.d-1,1.d0/3.d0,2.d-1,2.d-1,2.d-1,   1.d0/3.d0,6.d-1/)
49      xr=(/0.d0,5.d-1,1.d0,0.d0,5.d-1,0.d0,0.d0,0.d0,   0.d0,5.d-1,5.d-1,1.d0/)
50
51      call getarg(1,arg)
52
53      call buildmats(qref,qweight,interp,grad)
54
55      call ceedinit(trim(arg)//char(0),ceed,err)
56
57      call ceedbasiscreateh1(ceed,ceed_triangle,1,p,q,interp,grad,qref,qweight,&
58     & b,err)
59
60      do i=1,p
61        x1=xr(0*p+i)
62        x2=xr(1*p+i)
63        call feval(x1,x2,val)
64        iinput(i)=val
65      enddo
66
67      call ceedvectorcreate(ceed,p,input,err)
68      ioffset=0
69      call ceedvectorsetarray(input,ceed_mem_host,ceed_use_pointer,iinput,&
70     & ioffset,err)
71      call ceedvectorcreate(ceed,q*d,output,err)
72      call ceedvectorsetvalue(output,0.d0,err)
73
74      call ceedbasisapply(b,1,ceed_notranspose,ceed_eval_grad,input,output,err)
75
76      call ceedvectorgetarrayread(output,ceed_mem_host,ooutput,ooffset,err)
77      do i=1,q
78        x1=xq(0*q+i)
79        x2=xq(1*q+i)
80        call dfeval(x1,x2,val)
81        diff=val-ooutput(0*q+i+ooffset)
82        if (abs(diff)>1.0d-10) then
83! LCOV_EXCL_START
84          write(*,'(A,I1,A,F12.8,A,F12.8)')  '[',i,'] ',ooutput(i+ooffset),&
85     &     ' != ',val
86! LCOV_EXCL_STOP
87        endif
88        call dfeval(x2,x1,val)
89        diff=val-ooutput(1*q+i+ooffset)
90        if (abs(diff)>1.0d-10) then
91! LCOV_EXCL_START
92          write(*,'(A,I1,A,F12.8,A,F12.8)')  '[',i,'] ',ooutput(i+ooffset),&
93     &     ' != ',val
94! LCOV_EXCL_STOP
95        endif
96      enddo
97      call ceedvectorrestorearrayread(output,ooutput,ooffset,err)
98
99      call ceedvectordestroy(input,err)
100      call ceedvectordestroy(output,err)
101      call ceedbasisdestroy(b,err)
102      call ceeddestroy(ceed,err)
103
104      end
105!-----------------------------------------------------------------------
106