File: test_fsunmatrix_dense_mod.f90

package info (click to toggle)
sundials 6.4.1%2Bdfsg1-3
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 79,368 kB
  • sloc: ansic: 218,700; f90: 62,503; cpp: 61,511; fortran: 5,166; python: 4,642; sh: 4,114; makefile: 562; perl: 123
file content (274 lines) | stat: -rw-r--r-- 6,961 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
! -----------------------------------------------------------------
! Programmer(s): Cody J. Balos @ LLNL
! -----------------------------------------------------------------
! SUNDIALS Copyright Start
! Copyright (c) 2002-2022, Lawrence Livermore National Security
! and Southern Methodist University.
! All rights reserved.
!
! See the top-level LICENSE and NOTICE files for details.
!
! SPDX-License-Identifier: BSD-3-Clause
! SUNDIALS Copyright End
! -----------------------------------------------------------------
! This file tests the Fortran 2003 interface to the SUNDIALS
! dense SUNMatrix implementation.
! -----------------------------------------------------------------

module test_fsunmatrix_dense
  use, intrinsic :: iso_c_binding
  use test_utilities
  implicit none

  integer(C_LONG_LONG), parameter :: N = 4

contains

  integer(C_INT) function smoke_tests() result(fails)

    !======== Inclusions ==========
    use, intrinsic :: iso_c_binding
    use fsundials_nvector_mod
    use fsundials_matrix_mod
    use fsunmatrix_dense_mod
    use fnvector_serial_mod

    !======== Declarations ========
    implicit none

    ! local variables
    type(SUNMatrix), pointer :: A, B               ! SUNMatrix
    type(N_Vector),  pointer :: x, y               ! NVectors
    real(C_DOUBLE),  pointer :: matdat(:)          ! matrix data pointer
    integer(C_LONG)          :: lenrw(1), leniw(1) ! matrix real and int work space size
    integer(C_LONG)          :: val

    fails = 0

    x => FN_VNew_Serial(N, sunctx)
    y => FN_VNew_Serial(N, sunctx)

    !===== Calls to interface =====

    ! constructor
    A => FSUNDenseMatrix(N, N, sunctx)
    if (.not. associated(A)) then
      print *,'>>> FAILED - ERROR in FSUNDenseMatrix; halting'
      stop 1
    end if

    ! misc. matrix functions
    val = FSUNMatGetID_Dense(A)
    val = FSUNDenseMatrix_Rows(A)
    val = FSUNDenseMatrix_Columns(A)
    val = FSUNDenseMatrix_LData(A)
    matdat => FSUNDenseMatrix_Data(A)
    matdat => FSUNDenseMatrix_Column(A,N)

    ! matrix operations
    B => FSUNMatClone_Dense(A)
    if (.not. associated(B)) then
      print *,'>>> FAILED - ERROR in FSUNMatClone_Dense; halting'
      stop 1
    end if
    fails = fails + FSUNMatZero_Dense(A)
    fails = fails + FSUNMatCopy_Dense(A,B)
    fails = fails + FSUNMatScaleAdd_Dense(ONE, A, B)
    fails = fails + FSUNMatScaleAddI_Dense(ONE, A)
    fails = fails + FSUNMatMatvec_Dense(A, x, y)
    fails = fails + FSUNMatSpace_Dense(A, lenrw, leniw)

    !======= Cleanup ===========
    call FSUNMatDestroy_Dense(A)
    call FSUNMatDestroy_Dense(B)
    call FN_VDestroy_Serial(x)
    call FN_VDestroy_Serial(y)

  end function smoke_tests

  integer(C_INT) function unit_tests() result(fails)
    use, intrinsic :: iso_c_binding
    use fsundials_nvector_mod
    use fsundials_matrix_mod
    use fnvector_serial_mod
    use fsunmatrix_dense_mod
    use test_sunmatrix

    implicit none

    type(SUNMatrix), pointer :: A, I
    type(N_Vector),  pointer :: x, y
    real(C_DOUBLE),  pointer :: Adata(:), Idata(:), xdata(:), ydata(:)
    integer(C_LONG)          :: ii, jj, tmp1, tmp2

    fails = 0

    A => FSUNDenseMatrix(N, N, sunctx)
    I => FSUNDenseMatrix(N, N, sunctx)
    x => FN_VNew_Serial(N, sunctx)
    y => FN_VNew_Serial(N, sunctx)

    ! fill matrix A
    Adata => FSUNDenseMatrix_Data(A)
    do jj=1, N
      do ii=1, N
        Adata((jj-1)*N + ii) = jj*(ii+jj-2)
      end do
    end do

    ! fill matrix I (identity)
    Idata => FSUNDenseMatrix_Data(I)
    do jj=1, N
        Idata((jj-1)*N + jj) = ONE
    end do

    ! fill vector x
    xdata => FN_VGetArrayPointer(x)
    do ii=1, N
      xdata(ii) = ONE / ii
    end do

    ! fill vector y
    ydata => FN_VGetArrayPointer(y)
    do ii=1, N
      tmp1 = ii-1
      tmp2 = tmp1 + N - 1
      ydata(ii) = HALF*(tmp2+1-tmp1)*(tmp1+tmp2)
    end do

    fails = fails + Test_FSUNMatGetID(A, SUNMATRIX_DENSE, 0)
    fails = fails + Test_FSUNMatClone(A, 0)
    fails = fails + Test_FSUNMatCopy(A, 0)
    fails = fails + Test_FSUNMatZero(A, 0)
    fails = fails + Test_FSUNMatScaleAdd(A, I, 0)
    fails = fails + Test_FSUNMatScaleAddI(A, I, 0)
    fails = fails + Test_FSUNMatMatvec(A, x, y, 0)
    fails = fails + Test_FSUNMatSpace(A, 0)

    ! cleanup
    call FSUNMatDestroy(A)
    call FSUNMatDestroy(I)
    call FN_VDestroy(x)
    call FN_VDestroy(y)

  end function unit_tests

end module

program main
  !======== Inclusions ==========
  use, intrinsic :: iso_c_binding
  use test_fsunmatrix_dense

  !======== Declarations ========
  implicit none
  integer(C_INT) :: fails = 0

  !============== Introduction =============
  print *, 'Dense SUNMatrix Fortran 2003 interface test'

  call Test_Init(c_null_ptr)

  fails = unit_tests()
  if (fails /= 0) then
    print *, 'FAILURE: n unit tests failed'
    stop 1
  else
    print *, 'SUCCESS: all unit tests passed'
  end if

  call Test_Finalize()

end program main

! exported functions used by test_sunmatrix
integer(C_INT) function check_matrix(A, B, tol) result(fails)
  use, intrinsic :: iso_c_binding
  use fsundials_matrix_mod
  use fsunmatrix_dense_mod
  use test_utilities

  implicit none

  type(SUNMatrix) :: A, B
  real(C_DOUBLE)  :: tol
  real(C_DOUBLE), pointer :: Adata(:), Bdata(:)
  integer(C_LONG) :: Aldata, Bldata, i

  fails = 0

  ! get data pointers
  Adata => FSUNDenseMatrix_Data(A)
  Bdata => FSUNDenseMatrix_Data(B)

  ! get and check data lengths
  Aldata = FSUNDenseMatrix_LData(A)
  Bldata = FSUNDenseMatrix_LData(B)

  if (Aldata /= Bldata) then
    print *, ">>> ERROR: check_matrix: Different data array lengths"
    fails = 1
    return
  end if

  ! compare data
  do i=1, Aldata
    fails = fails + FNEQTOL(Adata(i), Bdata(i), tol)
  end do

end function check_matrix

integer(C_INT) function check_matrix_entry(A, c, tol) result(fails)
  use, intrinsic :: iso_c_binding
  use fsundials_matrix_mod
  use fsunmatrix_dense_mod
  use test_utilities

  implicit none

  type(SUNMatrix) :: A
  real(C_DOUBLE)  :: c, tol
  real(C_DOUBLE), pointer :: Adata(:)
  integer(C_LONG) :: Aldata, i

  fails = 0

  ! get data pointers
  Adata => FSUNDenseMatrix_Data(A)

  ! get and check data lengths
  Aldata = FSUNDenseMatrix_LData(A)

  ! compare data
  do i=1, Aldata
    fails = fails + FNEQTOL(Adata(i), c, tol)
  end do

  if (fails > ZERO) then
    print *, ">>> ERROR: check_matrix_entry failures: "
    do i=1, Aldata
      if (FNEQTOL(Adata(i), c, tol) /= 0) then
        write(*,'(A,I0,A,E14.7,A,E14.7)') &
          "Adata[ ", i, "] =", Adata(i) ," c = ", c
      end if
    end do
  end if

end function check_matrix_entry

logical function is_square(A) result(res)
  use, intrinsic :: iso_c_binding
  use fsundials_matrix_mod
  use fsunmatrix_dense_mod

  implicit none

  type(SUNMatrix) :: A

  if (FSUNDenseMatrix_Rows(A) == FSUNDenseMatrix_Columns(A)) then
    res = .true.
  else
    res = .false.
  end if

end function is_square