File: test_nvector.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 (183 lines) | stat: -rw-r--r-- 4,649 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
! -----------------------------------------------------------------
! Programmer(s): Cody J. Balos @ LLNL
! -----------------------------------------------------------------
! Acknowledgements: These testing routines are based on
! test_nvector.c written by David Gardner and Slaven Peles @ 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
! -----------------------------------------------------------------
! These test functions are designed to check the fortran interface
! to an NVECTOR module implementation. It does not test every
! function. It tests the N_VMake constructor, one standard vector
! operation (N_VConst), N_VGetArrayPointer, and one fused operation.
! -----------------------------------------------------------------

module test_fnvector
  use, intrinsic :: iso_c_binding
  use fsundials_nvector_mod
  use fsundials_types_mod
  use test_utilities
  implicit none

  integer(C_INT), external :: check_ans
  logical, external        :: has_data

contains


integer(C_INT) function Test_FN_VMake(X, local_length, myid) &
    result(failure)
  implicit none

  type(N_Vector)  :: X
  integer(C_LONG) :: local_length
  integer(C_INT)  :: myid

  if (.not. has_data(X)) then
    print *, '(I4)', '>>> FAILED test -- FN_VMake, Proc ', myid
    print *, '    vector data is not associated'
    failure = 1
    return
  end if

  if (myid == 0) then
    print *, 'PASSED test -- FN_VMake'
  end if

  failure = 0
end function Test_FN_VMake


!! ----------------------------------------------------------------------
!! NOTE: This routine depends on FN_VConst to check vector data.
!! ----------------------------------------------------------------------
integer(C_INT) function Test_FN_VGetArrayPointer(W, local_length, myid) &
    result(failure)
  implicit none

  type(N_Vector)  :: W
  integer(C_LONG) :: local_length
  integer(C_INT)  :: myid

  ! check vector data
  if (.not. has_data(W)) then
    print *, '>>> FAILED test -- FN_VGetArrayPointer, Proc ', myid
    print *, '    Vector data == NULL \n\n'
    failure = 1
    return;
  end if

  call FN_VConst(NEG_HALF, W)
  failure = check_ans(NEG_HALF, W, local_length)

  if (failure > 0) then
    print *, '(I2)', '>>> FAILED test -- FN_VGetArrayPointer, Proc ', myid
    print *, '    Failed FN_VConst check \n\n'
    failure = 1
    return
  end if

  if (myid == 0) then
    print *, 'PASSED test -- FN_VConst'
    print *, 'PASSED test -- FN_VGetArrayPointer'
  end if

  failure = 0
end function Test_FN_VGetArrayPointer


integer(C_INT) function Test_FN_VLinearCombination(X, local_length, myid) &
    result(failure)

  type(N_Vector)          :: X
  integer(C_LONG)         :: local_length
  integer(C_INT)          :: myid, ierr
  type(N_Vector), pointer :: Y1, Y2, Y3
  type(c_ptr), target     :: V(3)
  type(c_ptr)             :: Vptr
  real(C_DOUBLE)          :: c(3)

  failure = 0

  ! create vectors for testing
  Y1 => FN_VClone(X)
  Y2 => FN_VClone(X)
  Y3 => FN_VClone(X)

  ! set vectors in vector array
  V(1) = c_loc(Y1)
  V(2) = c_loc(Y2)
  V(3) = c_loc(Y3)
  Vptr = c_loc(V)

  ! initialize c values
  c = ZERO

  !
  ! Case 1a: V[0] = a V[0], FN_VScale
  !

  ! fill vector data
  call FN_VConst(TWO, Y1)

  ! set scaling factors
  c =  HALF

  ierr = FN_VLinearCombination(1, c, Vptr, Y1)

  ! Y1 should be vector of +1
  if (ierr == 0) then
    failure = check_ans(ONE, Y1, local_length)
  else
    failure = 1
  end if

  if (failure > 0) then
    print *, '(I4)', '>>> FAILED test -- FN_VLinearCombination Case 1a, Proc ', myid
  else if (myid == 0) then
    print *, 'PASSED test -- FN_VLinearCombination Case 1a'
  end if

  !
  ! Case 3a: V[0] = V[0] + b V[1] + c V[2]
  !

  call FN_VConst(TWO, Y1)
  call FN_VConst(NEG_TWO, Y2)
  call FN_VConst(NEG_ONE, Y3)

  c(1) = ONE
  c(2) = HALF
  c(3) = NEG_TWO

  ierr = FN_VLinearCombination(3, c, Vptr, Y1)

  ! Y1 should be vector of +3
  if (ierr == 0) then
    failure = check_ans(TWO+ONE, Y1, local_length)
  else
    failure = 1
  end if

  if (failure > 0) then
    print *, '(I4)', '>>> FAILED test -- FN_VLinearCombination Case 3a, Proc ', myid
  else if (myid == 0) then
    print *, 'PASSED test -- FN_VLinearCombination Case 3a'
  end if

  ! Free vectors
  call FN_VDestroy(Y1);
  call FN_VDestroy(Y2);
  call FN_VDestroy(Y3);

end function Test_FN_VLinearCombination

end module