(git:374b731)
Loading...
Searching...
No Matches
splines_types.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief routines for handling splines_types
10!> \par History
11!> 2001-09-21-HAF added this doc entry and changed formatting
12!> \author various
13! **************************************************************************************************
15
16 USE kinds, ONLY: dp
17#include "./base/base_uses.f90"
18
19 IMPLICIT NONE
20
21 PRIVATE
22 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'splines_types'
23
30 PUBLIC :: spline_data_type ! the data structure for spline table
31 PUBLIC :: spline_factor_type ! the multiplicative factors for splines
32
33! **************************************************************************************************
34!> \brief Data-structure that holds all needed information about
35!> a specific spline interpolation.
36!> \par History
37!> 2001-09-19-HAF added this doc entry and changed formatting
38!> \author unknown
39! **************************************************************************************************
41 INTEGER :: ref_count
42 REAL(kind=dp), POINTER :: y(:) ! the function values y(x)
43 REAL(kind=dp), POINTER :: y2(:) ! the 2nd derivative via interpolation
44 INTEGER :: n ! dimension of above arrays
45 ! not used if uniform increments
46 REAL(kind=dp) :: h ! uniform increment of x if applicable
47 REAL(kind=dp) :: invh ! inverse of h
48 REAL(kind=dp) :: h26 ! 1/6 * h**2 if uniform increments
49 ! 1/6 otherwise
50 REAL(kind=dp) :: x1 ! starting x value if uniform incr.
51 REAL(kind=dp) :: xn ! end x value if uniform incr.
52 END TYPE spline_data_type
53
54! **************************************************************************************************
56 TYPE(spline_data_type), POINTER :: spline_data
57 END TYPE spline_data_p_type
58
59! **************************************************************************************************
60 TYPE spline_data_pp_type
61 TYPE(spline_data_p_type), POINTER, DIMENSION(:) :: spl_p
62 END TYPE spline_data_pp_type
63
64! **************************************************************************************************
66 TYPE(spline_data_pp_type), POINTER, DIMENSION(:) :: spl_pp
67 INTEGER, POINTER, DIMENSION(:, :) :: spltab
69
70! **************************************************************************************************
72 REAL(kind=dp) :: rcutsq_f, cutoff
73 REAL(kind=dp), DIMENSION(:), POINTER :: rscale
74 REAL(kind=dp), DIMENSION(:), POINTER :: fscale
75 REAL(kind=dp), DIMENSION(:), POINTER :: dscale
76 END TYPE spline_factor_type
77
78CONTAINS
79
80! **************************************************************************************************
81!> \brief releases spline_env
82!> \param spline_env ...
83!> \author unknown
84! **************************************************************************************************
85 SUBROUTINE spline_env_release(spline_env)
86 TYPE(spline_environment_type), INTENT(INOUT) :: spline_env
87
88 INTEGER :: i
89 TYPE(spline_data_p_type), DIMENSION(:), POINTER :: spl_p
90
91 DEALLOCATE (spline_env%spltab)
92 DO i = 1, SIZE(spline_env%spl_pp)
93 spl_p => spline_env%spl_pp(i)%spl_p
94 CALL spline_data_p_release(spl_p)
95 END DO
96 DEALLOCATE (spline_env%spl_pp)
97
98 END SUBROUTINE spline_env_release
99
100! **************************************************************************************************
101!> \brief releases spline_data
102!> \param spline_data ...
103!> \author CJM
104! **************************************************************************************************
105 SUBROUTINE spline_data_release(spline_data)
106 TYPE(spline_data_type), POINTER :: spline_data
107
108 IF (ASSOCIATED(spline_data)) THEN
109 cpassert(spline_data%ref_count > 0)
110 spline_data%ref_count = spline_data%ref_count - 1
111 IF (spline_data%ref_count < 1) THEN
112 IF (ASSOCIATED(spline_data%y)) THEN
113 DEALLOCATE (spline_data%y)
114 END IF
115 IF (ASSOCIATED(spline_data%y2)) THEN
116 DEALLOCATE (spline_data%y2)
117 END IF
118 DEALLOCATE (spline_data)
119 END IF
120 END IF
121 END SUBROUTINE spline_data_release
122
123! **************************************************************************************************
124!> \brief releases spline_data_p
125!> \param spl_p ...
126!> \author CJM
127! **************************************************************************************************
128 SUBROUTINE spline_data_p_release(spl_p)
129 TYPE(spline_data_p_type), DIMENSION(:), POINTER :: spl_p
130
131 INTEGER :: i
132 LOGICAL :: release_kind
133
134 IF (ASSOCIATED(spl_p)) THEN
135 release_kind = .true.
136 DO i = 1, SIZE(spl_p)
137 CALL spline_data_release(spl_p(i)%spline_data)
138 release_kind = release_kind .AND. (.NOT. ASSOCIATED(spl_p(i)%spline_data))
139 END DO
140 IF (release_kind) THEN
141 DEALLOCATE (spl_p)
142 END IF
143 END IF
144
145 END SUBROUTINE spline_data_p_release
146
147! **************************************************************************************************
148!> \brief retains spline_env
149!> \param spline_data ...
150!> \author CJM
151! **************************************************************************************************
152 SUBROUTINE spline_data_retain(spline_data)
153 TYPE(spline_data_type), POINTER :: spline_data
154
155 cpassert(ASSOCIATED(spline_data))
156 cpassert(spline_data%ref_count > 0)
157 spline_data%ref_count = spline_data%ref_count + 1
158 END SUBROUTINE spline_data_retain
159
160! **************************************************************************************************
161!> \brief retains spline_data_p_type
162!> \param spl_p ...
163!> \author CJM
164! **************************************************************************************************
165 SUBROUTINE spline_data_p_retain(spl_p)
166 TYPE(spline_data_p_type), DIMENSION(:), POINTER :: spl_p
167
168 INTEGER :: i
169
170 cpassert(ASSOCIATED(spl_p))
171 DO i = 1, SIZE(spl_p)
172 CALL spline_data_retain(spl_p(i)%spline_data)
173 END DO
174 END SUBROUTINE spline_data_p_retain
175
176! **************************************************************************************************
177!> \brief Data-structure that holds all needed information about
178!> a specific spline interpolation.
179!> \param spline_env ...
180!> \param ntype ...
181!> \param ntab_in ...
182!> \par History
183!> 2001-09-19-HAF added this doc entry and changed formatting
184!> \author unknown
185! **************************************************************************************************
186 SUBROUTINE spline_env_create(spline_env, ntype, ntab_in)
187 TYPE(spline_environment_type), INTENT(OUT) :: spline_env
188 INTEGER, INTENT(IN) :: ntype
189 INTEGER, INTENT(IN), OPTIONAL :: ntab_in
190
191 CHARACTER(len=*), PARAMETER :: routinen = 'spline_env_create'
192
193 INTEGER :: handle, i, isize, j, ntab
194
195 CALL timeset(routinen, handle)
196
197 NULLIFY (spline_env%spl_pp)
198 NULLIFY (spline_env%spltab)
199 ! Allocate the number of spline data tables (upper triangular)
200 IF (PRESENT(ntab_in)) THEN
201 ntab = ntab_in
202 ELSE
203 ntab = (ntype*ntype + ntype)/2
204 END IF
205 ALLOCATE (spline_env%spl_pp(ntab))
206
207 ALLOCATE (spline_env%spltab(ntype, ntype))
208
209 DO i = 1, ntab
210 NULLIFY (spline_env%spl_pp(i)%spl_p)
211 isize = 1
212 ALLOCATE (spline_env%spl_pp(i)%spl_p(isize))
213 DO j = 1, SIZE(spline_env%spl_pp(i)%spl_p)
214 CALL spline_data_create(spline_env%spl_pp(i)%spl_p(j)%spline_data)
215 END DO
216 END DO
217
218 CALL timestop(handle)
219
220 END SUBROUTINE spline_env_create
221
222! **************************************************************************************************
223!> \brief Copy Data-structure of spline_data_p_type
224!> \param spl_p_source ...
225!> \param spl_p_dest ...
226!> \author teo 06.2007
227! **************************************************************************************************
228 SUBROUTINE spline_data_p_copy(spl_p_source, spl_p_dest)
229 TYPE(spline_data_p_type), DIMENSION(:), POINTER :: spl_p_source, spl_p_dest
230
231 INTEGER :: i, nsized, nsizes
232
233 cpassert(ASSOCIATED(spl_p_source))
234 nsizes = SIZE(spl_p_source)
235 IF (.NOT. ASSOCIATED(spl_p_dest)) THEN
236 ALLOCATE (spl_p_dest(nsizes))
237 DO i = 1, nsizes
238 NULLIFY (spl_p_dest(i)%spline_data)
239 END DO
240 ELSE
241 nsized = SIZE(spl_p_dest)
242 cpassert(nsizes == nsized)
243 DO i = 1, nsizes
244 CALL spline_data_release(spl_p_dest(i)%spline_data)
245 END DO
246 END IF
247 DO i = 1, nsizes
248 CALL spline_data_copy(spl_p_source(i)%spline_data, spl_p_dest(i)%spline_data)
249 END DO
250 END SUBROUTINE spline_data_p_copy
251
252! **************************************************************************************************
253!> \brief Copy Data-structure that constains spline table
254!> \param spline_data_source ...
255!> \param spline_data_dest ...
256!> \author teo 11.2005
257! **************************************************************************************************
258 SUBROUTINE spline_data_copy(spline_data_source, spline_data_dest)
259 TYPE(spline_data_type), POINTER :: spline_data_source, spline_data_dest
260
261 cpassert(ASSOCIATED(spline_data_source))
262 IF (.NOT. ASSOCIATED(spline_data_dest)) CALL spline_data_create(spline_data_dest)
263
264 spline_data_dest%ref_count = spline_data_source%ref_count
265 spline_data_dest%n = spline_data_source%n
266 spline_data_dest%h = spline_data_source%h
267 spline_data_dest%invh = spline_data_source%invh
268 spline_data_dest%h26 = spline_data_source%h26
269 spline_data_dest%x1 = spline_data_source%x1
270 spline_data_dest%xn = spline_data_source%xn
271 IF (ASSOCIATED(spline_data_source%y)) THEN
272 ALLOCATE (spline_data_dest%y(SIZE(spline_data_source%y)))
273 spline_data_dest%y = spline_data_source%y
274 END IF
275 IF (ASSOCIATED(spline_data_source%y2)) THEN
276 ALLOCATE (spline_data_dest%y2(SIZE(spline_data_source%y2)))
277 spline_data_dest%y2 = spline_data_source%y2
278 END IF
279 END SUBROUTINE spline_data_copy
280
281! **************************************************************************************************
282!> \brief Data-structure that constains spline table
283!> \param spline_data ...
284!> \author unknown
285! **************************************************************************************************
286 SUBROUTINE spline_data_create(spline_data)
287 TYPE(spline_data_type), POINTER :: spline_data
288
289 ALLOCATE (spline_data)
290 spline_data%ref_count = 1
291 NULLIFY (spline_data%y)
292 NULLIFY (spline_data%y2)
293 END SUBROUTINE spline_data_create
294
295! **************************************************************************************************
296!> \brief releases spline_factor
297!> \param spline_factor ...
298!> \author teo
299! **************************************************************************************************
300 SUBROUTINE spline_factor_release(spline_factor)
301 TYPE(spline_factor_type), POINTER :: spline_factor
302
303 IF (ASSOCIATED(spline_factor)) THEN
304 IF (ASSOCIATED(spline_factor%rscale)) THEN
305 DEALLOCATE (spline_factor%rscale)
306 END IF
307 IF (ASSOCIATED(spline_factor%fscale)) THEN
308 DEALLOCATE (spline_factor%fscale)
309 END IF
310 IF (ASSOCIATED(spline_factor%dscale)) THEN
311 DEALLOCATE (spline_factor%dscale)
312 END IF
313 DEALLOCATE (spline_factor)
314 END IF
315 END SUBROUTINE spline_factor_release
316
317! **************************************************************************************************
318!> \brief releases spline_factor
319!> \param spline_factor ...
320!> \author teo
321! **************************************************************************************************
322 SUBROUTINE spline_factor_create(spline_factor)
323 TYPE(spline_factor_type), POINTER :: spline_factor
324
325 cpassert(.NOT. ASSOCIATED(spline_factor))
326 ALLOCATE (spline_factor)
327 ALLOCATE (spline_factor%rscale(1))
328 ALLOCATE (spline_factor%fscale(1))
329 ALLOCATE (spline_factor%dscale(1))
330 spline_factor%rscale = 1.0_dp
331 spline_factor%fscale = 1.0_dp
332 spline_factor%dscale = 1.0_dp
333 spline_factor%rcutsq_f = 1.0_dp
334 spline_factor%cutoff = 0.0_dp
335 END SUBROUTINE spline_factor_create
336
337! **************************************************************************************************
338!> \brief releases spline_factor
339!> \param spline_factor_source ...
340!> \param spline_factor_dest ...
341!> \author teo
342! **************************************************************************************************
343 SUBROUTINE spline_factor_copy(spline_factor_source, spline_factor_dest)
344 TYPE(spline_factor_type), POINTER :: spline_factor_source, spline_factor_dest
345
346 INTEGER :: isize, jsize, ksize
347
348 IF (ASSOCIATED(spline_factor_dest)) CALL spline_factor_release(spline_factor_dest)
349 IF (ASSOCIATED(spline_factor_source)) THEN
350 isize = SIZE(spline_factor_source%rscale)
351 jsize = SIZE(spline_factor_source%fscale)
352 ksize = SIZE(spline_factor_source%dscale)
353 cpassert(isize == jsize)
354 cpassert(isize == ksize)
355 CALL spline_factor_create(spline_factor_dest)
356 spline_factor_dest%rscale = spline_factor_source%rscale
357 spline_factor_dest%fscale = spline_factor_source%fscale
358 spline_factor_dest%dscale = spline_factor_source%dscale
359 spline_factor_dest%rcutsq_f = spline_factor_source%rcutsq_f
360 spline_factor_dest%cutoff = spline_factor_source%cutoff
361 END IF
362 END SUBROUTINE spline_factor_copy
363
364END MODULE splines_types
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
routines for handling splines_types
subroutine, public spline_data_p_copy(spl_p_source, spl_p_dest)
Copy Data-structure of spline_data_p_type.
subroutine, public spline_data_retain(spline_data)
retains spline_env
subroutine, public spline_factor_release(spline_factor)
releases spline_factor
subroutine, public spline_data_p_release(spl_p)
releases spline_data_p
subroutine, public spline_factor_create(spline_factor)
releases spline_factor
subroutine, public spline_data_create(spline_data)
Data-structure that constains spline table.
subroutine, public spline_data_release(spline_data)
releases spline_data
subroutine, public spline_env_release(spline_env)
releases spline_env
subroutine, public spline_data_p_retain(spl_p)
retains spline_data_p_type
subroutine, public spline_env_create(spline_env, ntype, ntab_in)
Data-structure that holds all needed information about a specific spline interpolation.
subroutine, public spline_factor_copy(spline_factor_source, spline_factor_dest)
releases spline_factor
Data-structure that holds all needed information about a specific spline interpolation.