(git:b279b6b)
cell_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 Handles all functions related to the CELL
10 !> \par History
11 !> 11.2008 Teodoro Laino [tlaino] - deeply cleaning cell_type from units
12 !> 10.2014 Moved many routines from cell_types.F here.
13 !> \author Matthias KracK (16.01.2002, based on a earlier version of CJM, JGH)
14 ! **************************************************************************************************
15 MODULE cell_types
16  USE cp_units, ONLY: cp_unit_to_cp2k
17  USE kinds, ONLY: dp
18  USE mathconstants, ONLY: degree
19  USE mathlib, ONLY: angle
20 #include "../base/base_uses.f90"
21 
22  IMPLICIT NONE
23 
24  PRIVATE
25 
26  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'cell_types'
27 
28  ! Impose cell symmetry
29  INTEGER, PARAMETER, PUBLIC :: cell_sym_none = 0, &
30  cell_sym_triclinic = 1, &
31  cell_sym_monoclinic = 2, &
40  cell_sym_cubic = 11
41 
42  INTEGER, PARAMETER, PUBLIC :: use_perd_x = 0, &
43  use_perd_y = 1, &
44  use_perd_z = 2, &
45  use_perd_xy = 3, &
46  use_perd_xz = 4, &
47  use_perd_yz = 5, &
48  use_perd_xyz = 6, &
49  use_perd_none = 7
50 
51 ! **************************************************************************************************
52 !> \brief Type defining parameters related to the simulation cell
53 !> \version 1.0
54 ! **************************************************************************************************
55  TYPE cell_type
56  CHARACTER(LEN=12) :: tag = "CELL"
57  INTEGER :: ref_count = -1, &
58  symmetry_id = use_perd_none
59  LOGICAL :: orthorhombic = .false. ! actually means a diagonal hmat
60  REAL(kind=dp) :: deth = 0.0_dp
61  INTEGER, DIMENSION(3) :: perd = -1
62  REAL(kind=dp), DIMENSION(3, 3) :: hmat = 0.0_dp, &
63  h_inv = 0.0_dp
64  END TYPE cell_type
65 
66  TYPE cell_p_type
67  TYPE(cell_type), POINTER :: cell => null()
68  END TYPE cell_p_type
69 
70  ! Public data types
71  PUBLIC :: cell_type, &
72  cell_p_type
73 
74  ! Public subroutines
75  PUBLIC :: cell_clone, &
76  cell_copy, &
77  cell_release, &
78  cell_retain, &
79  get_cell, &
81 
82 #if defined (__PLUMED2)
83  PUBLIC :: pbc_cp2k_plumed_getset_cell
84 #endif
85 
86  ! Public functions
87  PUBLIC :: plane_distance, &
88  pbc, &
91 
92  INTERFACE pbc
93  MODULE PROCEDURE pbc1, pbc2, pbc3, pbc4
94  END INTERFACE
95 
96 CONTAINS
97 
98 ! **************************************************************************************************
99 !> \brief Clone cell variable
100 !> \param cell_in Cell variable to be clone
101 !> \param cell_out Cloned cell variable
102 !> \param tag Optional new tag for cloned cell variable
103 !> \par History
104 !> - Optional tag added (17.05.2023, MK)
105 ! **************************************************************************************************
106  SUBROUTINE cell_clone(cell_in, cell_out, tag)
107 
108  TYPE(cell_type), POINTER :: cell_in, cell_out
109  CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: tag
110 
111  cell_out = cell_in
112  cell_out%ref_count = 1
113  IF (PRESENT(tag)) cell_out%tag = tag
114 
115  END SUBROUTINE cell_clone
116 
117 ! **************************************************************************************************
118 !> \brief Copy cell variable
119 !> \param cell_in Cell variable to be copied
120 !> \param cell_out Copy of cell variable
121 !> \param tag Optional new tag
122 !> \par History
123 !> - Optional tag added (17.05.2023, MK)
124 ! **************************************************************************************************
125  SUBROUTINE cell_copy(cell_in, cell_out, tag)
126 
127  TYPE(cell_type), POINTER :: cell_in, cell_out
128  CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: tag
129 
130  cell_out%deth = cell_in%deth
131  cell_out%perd = cell_in%perd
132  cell_out%hmat = cell_in%hmat
133  cell_out%h_inv = cell_in%h_inv
134  cell_out%orthorhombic = cell_in%orthorhombic
135  cell_out%symmetry_id = cell_in%symmetry_id
136  IF (PRESENT(tag)) THEN
137  cell_out%tag = tag
138  ELSE
139  cell_out%tag = cell_in%tag
140  END IF
141 
142  END SUBROUTINE cell_copy
143 
144 ! **************************************************************************************************
145 !> \brief Read cell info from a line (parsed from a file)
146 !> \param input_line ...
147 !> \param cell_itimes ...
148 !> \param cell_time ...
149 !> \param h ...
150 !> \param vol ...
151 !> \date 19.02.2008
152 !> \author Teodoro Laino [tlaino] - University of Zurich
153 !> \version 1.0
154 ! **************************************************************************************************
155  SUBROUTINE parse_cell_line(input_line, cell_itimes, cell_time, h, vol)
156 
157  CHARACTER(LEN=*), INTENT(IN) :: input_line
158  INTEGER, INTENT(OUT) :: cell_itimes
159  REAL(kind=dp), INTENT(OUT) :: cell_time
160  REAL(kind=dp), DIMENSION(3, 3), INTENT(OUT) :: h
161  REAL(kind=dp), INTENT(OUT) :: vol
162 
163  INTEGER :: i, j
164 
165  READ (input_line, *) cell_itimes, cell_time, &
166  h(1, 1), h(2, 1), h(3, 1), h(1, 2), h(2, 2), h(3, 2), h(1, 3), h(2, 3), h(3, 3), vol
167  DO i = 1, 3
168  DO j = 1, 3
169  h(j, i) = cp_unit_to_cp2k(h(j, i), "angstrom")
170  END DO
171  END DO
172 
173  END SUBROUTINE parse_cell_line
174 
175 ! **************************************************************************************************
176 !> \brief Get informations about a simulation cell.
177 !> \param cell ...
178 !> \param alpha ...
179 !> \param beta ...
180 !> \param gamma ...
181 !> \param deth ...
182 !> \param orthorhombic ...
183 !> \param abc ...
184 !> \param periodic ...
185 !> \param h ...
186 !> \param h_inv ...
187 !> \param symmetry_id ...
188 !> \param tag ...
189 !> \date 16.01.2002
190 !> \author Matthias Krack
191 !> \version 1.0
192 ! **************************************************************************************************
193  SUBROUTINE get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, &
194  h, h_inv, symmetry_id, tag)
195 
196  TYPE(cell_type), POINTER :: cell
197  REAL(kind=dp), INTENT(OUT), OPTIONAL :: alpha, beta, gamma, deth
198  LOGICAL, INTENT(OUT), OPTIONAL :: orthorhombic
199  REAL(kind=dp), DIMENSION(3), INTENT(OUT), OPTIONAL :: abc
200  INTEGER, DIMENSION(3), INTENT(OUT), OPTIONAL :: periodic
201  REAL(kind=dp), DIMENSION(3, 3), INTENT(OUT), &
202  OPTIONAL :: h, h_inv
203  INTEGER, INTENT(OUT), OPTIONAL :: symmetry_id
204  CHARACTER(LEN=*), INTENT(OUT), OPTIONAL :: tag
205 
206  cpassert(ASSOCIATED(cell))
207 
208  IF (PRESENT(deth)) deth = cell%deth ! the volume
209  IF (PRESENT(orthorhombic)) orthorhombic = cell%orthorhombic
210  IF (PRESENT(periodic)) periodic(:) = cell%perd(:)
211  IF (PRESENT(h)) h(:, :) = cell%hmat(:, :)
212  IF (PRESENT(h_inv)) h_inv(:, :) = cell%h_inv(:, :)
213 
214  ! Calculate the lengths of the cell vectors a, b, and c
215  IF (PRESENT(abc)) THEN
216  abc(1) = sqrt(cell%hmat(1, 1)*cell%hmat(1, 1) + &
217  cell%hmat(2, 1)*cell%hmat(2, 1) + &
218  cell%hmat(3, 1)*cell%hmat(3, 1))
219  abc(2) = sqrt(cell%hmat(1, 2)*cell%hmat(1, 2) + &
220  cell%hmat(2, 2)*cell%hmat(2, 2) + &
221  cell%hmat(3, 2)*cell%hmat(3, 2))
222  abc(3) = sqrt(cell%hmat(1, 3)*cell%hmat(1, 3) + &
223  cell%hmat(2, 3)*cell%hmat(2, 3) + &
224  cell%hmat(3, 3)*cell%hmat(3, 3))
225  END IF
226 
227  ! Angles between the cell vectors a, b, and c
228  ! alpha = <(b,c)
229  IF (PRESENT(alpha)) alpha = angle(cell%hmat(:, 2), cell%hmat(:, 3))*degree
230  ! beta = <(a,c)
231  IF (PRESENT(beta)) beta = angle(cell%hmat(:, 1), cell%hmat(:, 3))*degree
232  ! gamma = <(a,b)
233  IF (PRESENT(gamma)) gamma = angle(cell%hmat(:, 1), cell%hmat(:, 2))*degree
234  IF (PRESENT(symmetry_id)) symmetry_id = cell%symmetry_id
235  IF (PRESENT(tag)) tag = cell%tag
236 
237  END SUBROUTINE get_cell
238 
239 ! **************************************************************************************************
240 !> \brief Calculate the distance between two lattice planes as defined by
241 !> a triple of Miller indices (hkl).
242 !> \param h ...
243 !> \param k ...
244 !> \param l ...
245 !> \param cell ...
246 !> \return ...
247 !> \date 18.11.2004
248 !> \author Matthias Krack
249 !> \version 1.0
250 ! **************************************************************************************************
251  FUNCTION plane_distance(h, k, l, cell) RESULT(distance)
252 
253  INTEGER, INTENT(IN) :: h, k, l
254  TYPE(cell_type), POINTER :: cell
255  REAL(kind=dp) :: distance
256 
257  REAL(kind=dp) :: a, alpha, b, beta, c, cosa, cosb, cosg, &
258  d, gamma, x, y, z
259  REAL(kind=dp), DIMENSION(3) :: abc
260 
261  x = real(h, kind=dp)
262  y = real(k, kind=dp)
263  z = real(l, kind=dp)
264 
265  CALL get_cell(cell=cell, abc=abc)
266 
267  a = abc(1)
268  b = abc(2)
269  c = abc(3)
270 
271  IF (cell%orthorhombic) THEN
272 
273  d = (x/a)**2 + (y/b)**2 + (z/c)**2
274 
275  ELSE
276 
277  CALL get_cell(cell=cell, &
278  alpha=alpha, &
279  beta=beta, &
280  gamma=gamma)
281 
282  alpha = alpha/degree
283  beta = beta/degree
284  gamma = gamma/degree
285 
286  cosa = cos(alpha)
287  cosb = cos(beta)
288  cosg = cos(gamma)
289 
290  d = ((x*b*c*sin(alpha))**2 + &
291  (y*c*a*sin(beta))**2 + &
292  (z*a*b*sin(gamma))**2 + &
293  2.0_dp*a*b*c*(x*y*c*(cosa*cosb - cosg) + &
294  z*x*b*(cosg*cosa - cosb) + &
295  y*z*a*(cosb*cosg - cosa)))/ &
296  ((a*b*c)**2*(1.0_dp - cosa**2 - cosb**2 - cosg**2 + &
297  2.0_dp*cosa*cosb*cosg))
298 
299  END IF
300 
301  distance = 1.0_dp/sqrt(d)
302 
303  END FUNCTION plane_distance
304 
305 ! **************************************************************************************************
306 !> \brief Apply the periodic boundary conditions defined by a simulation
307 !> cell to a position vector r.
308 !> \param r ...
309 !> \param cell ...
310 !> \return ...
311 !> \date 16.01.2002
312 !> \author Matthias Krack
313 !> \version 1.0
314 ! **************************************************************************************************
315  FUNCTION pbc1(r, cell) RESULT(r_pbc)
316 
317  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: r
318  TYPE(cell_type), POINTER :: cell
319  REAL(kind=dp), DIMENSION(3) :: r_pbc
320 
321  REAL(kind=dp), DIMENSION(3) :: s
322 
323  cpassert(ASSOCIATED(cell))
324 
325  IF (cell%orthorhombic) THEN
326  r_pbc(1) = r(1) - cell%hmat(1, 1)*cell%perd(1)*anint(cell%h_inv(1, 1)*r(1))
327  r_pbc(2) = r(2) - cell%hmat(2, 2)*cell%perd(2)*anint(cell%h_inv(2, 2)*r(2))
328  r_pbc(3) = r(3) - cell%hmat(3, 3)*cell%perd(3)*anint(cell%h_inv(3, 3)*r(3))
329  ELSE
330  s(1) = cell%h_inv(1, 1)*r(1) + cell%h_inv(1, 2)*r(2) + cell%h_inv(1, 3)*r(3)
331  s(2) = cell%h_inv(2, 1)*r(1) + cell%h_inv(2, 2)*r(2) + cell%h_inv(2, 3)*r(3)
332  s(3) = cell%h_inv(3, 1)*r(1) + cell%h_inv(3, 2)*r(2) + cell%h_inv(3, 3)*r(3)
333  s(1) = s(1) - cell%perd(1)*anint(s(1))
334  s(2) = s(2) - cell%perd(2)*anint(s(2))
335  s(3) = s(3) - cell%perd(3)*anint(s(3))
336  r_pbc(1) = cell%hmat(1, 1)*s(1) + cell%hmat(1, 2)*s(2) + cell%hmat(1, 3)*s(3)
337  r_pbc(2) = cell%hmat(2, 1)*s(1) + cell%hmat(2, 2)*s(2) + cell%hmat(2, 3)*s(3)
338  r_pbc(3) = cell%hmat(3, 1)*s(1) + cell%hmat(3, 2)*s(2) + cell%hmat(3, 3)*s(3)
339  END IF
340 
341  END FUNCTION pbc1
342 
343 ! **************************************************************************************************
344 !> \brief Apply the periodic boundary conditions defined by a simulation
345 !> cell to a position vector r subtracting nl from the periodic images
346 !> \param r ...
347 !> \param cell ...
348 !> \param nl ...
349 !> \return ...
350 !> \date 16.01.2002
351 !> \author Matthias Krack
352 !> \version 1.0
353 ! **************************************************************************************************
354  FUNCTION pbc2(r, cell, nl) RESULT(r_pbc)
355 
356  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: r
357  TYPE(cell_type), POINTER :: cell
358  INTEGER, DIMENSION(3), INTENT(IN) :: nl
359  REAL(kind=dp), DIMENSION(3) :: r_pbc
360 
361  REAL(kind=dp), DIMENSION(3) :: s
362 
363  cpassert(ASSOCIATED(cell))
364 
365  IF (cell%orthorhombic) THEN
366  r_pbc(1) = r(1) - cell%hmat(1, 1)*cell%perd(1)* &
367  REAL(nint(cell%h_inv(1, 1)*r(1)) - nl(1), dp)
368  r_pbc(2) = r(2) - cell%hmat(2, 2)*cell%perd(2)* &
369  REAL(nint(cell%h_inv(2, 2)*r(2)) - nl(2), dp)
370  r_pbc(3) = r(3) - cell%hmat(3, 3)*cell%perd(3)* &
371  REAL(nint(cell%h_inv(3, 3)*r(3)) - nl(3), dp)
372  ELSE
373  s(1) = cell%h_inv(1, 1)*r(1) + cell%h_inv(1, 2)*r(2) + cell%h_inv(1, 3)*r(3)
374  s(2) = cell%h_inv(2, 1)*r(1) + cell%h_inv(2, 2)*r(2) + cell%h_inv(2, 3)*r(3)
375  s(3) = cell%h_inv(3, 1)*r(1) + cell%h_inv(3, 2)*r(2) + cell%h_inv(3, 3)*r(3)
376  s(1) = s(1) - cell%perd(1)*real(nint(s(1)) - nl(1), dp)
377  s(2) = s(2) - cell%perd(2)*real(nint(s(2)) - nl(2), dp)
378  s(3) = s(3) - cell%perd(3)*real(nint(s(3)) - nl(3), dp)
379  r_pbc(1) = cell%hmat(1, 1)*s(1) + cell%hmat(1, 2)*s(2) + cell%hmat(1, 3)*s(3)
380  r_pbc(2) = cell%hmat(2, 1)*s(1) + cell%hmat(2, 2)*s(2) + cell%hmat(2, 3)*s(3)
381  r_pbc(3) = cell%hmat(3, 1)*s(1) + cell%hmat(3, 2)*s(2) + cell%hmat(3, 3)*s(3)
382  END IF
383 
384  END FUNCTION pbc2
385 
386 ! **************************************************************************************************
387 !> \brief Apply the periodic boundary conditions defined by the simulation
388 !> cell cell to the vector pointing from atom a to atom b.
389 !> \param ra ...
390 !> \param rb ...
391 !> \param cell ...
392 !> \return ...
393 !> \date 11.03.2004
394 !> \author Matthias Krack
395 !> \version 1.0
396 ! **************************************************************************************************
397  FUNCTION pbc3(ra, rb, cell) RESULT(rab_pbc)
398 
399  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: ra, rb
400  TYPE(cell_type), POINTER :: cell
401  REAL(kind=dp), DIMENSION(3) :: rab_pbc
402 
403  INTEGER :: icell, jcell, kcell
404  INTEGER, DIMENSION(3) :: periodic
405  REAL(kind=dp) :: rab2, rab2_pbc
406  REAL(kind=dp), DIMENSION(3) :: r, ra_pbc, rab, rb_image, rb_pbc, s2r
407 
408  CALL get_cell(cell=cell, periodic=periodic)
409 
410  ra_pbc(:) = pbc(ra(:), cell)
411  rb_pbc(:) = pbc(rb(:), cell)
412 
413  rab2_pbc = huge(1.0_dp)
414 
415  DO icell = -periodic(1), periodic(1)
416  DO jcell = -periodic(2), periodic(2)
417  DO kcell = -periodic(3), periodic(3)
418  r = real((/icell, jcell, kcell/), dp)
419  CALL scaled_to_real(s2r, r, cell)
420  rb_image(:) = rb_pbc(:) + s2r
421  rab(:) = rb_image(:) - ra_pbc(:)
422  rab2 = rab(1)*rab(1) + rab(2)*rab(2) + rab(3)*rab(3)
423  IF (rab2 < rab2_pbc) THEN
424  rab2_pbc = rab2
425  rab_pbc(:) = rab(:)
426  END IF
427  END DO
428  END DO
429  END DO
430 
431  END FUNCTION pbc3
432 
433  !if positive_range == true, r(i) (or s(i)) in range [0, hmat(i,i)],
434  !else, r(i) (s(i)) in range [-hmat(i,i)/2, hmat(i,i)/2]
435 ! **************************************************************************************************
436 !> \brief ...
437 !> \param r ...
438 !> \param cell ...
439 !> \param positive_range ...
440 !> \return ...
441 ! **************************************************************************************************
442  FUNCTION pbc4(r, cell, positive_range) RESULT(r_pbc)
443 
444  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: r
445  TYPE(cell_type), POINTER :: cell
446  LOGICAL :: positive_range
447  REAL(kind=dp), DIMENSION(3) :: r_pbc
448 
449  REAL(kind=dp), DIMENSION(3) :: s
450 
451  cpassert(ASSOCIATED(cell))
452 
453  IF (positive_range) THEN
454  IF (cell%orthorhombic) THEN
455  r_pbc(1) = r(1) - cell%hmat(1, 1)*cell%perd(1)*floor(cell%h_inv(1, 1)*r(1))
456  r_pbc(2) = r(2) - cell%hmat(2, 2)*cell%perd(2)*floor(cell%h_inv(2, 2)*r(2))
457  r_pbc(3) = r(3) - cell%hmat(3, 3)*cell%perd(3)*floor(cell%h_inv(3, 3)*r(3))
458  ELSE
459  s(1) = cell%h_inv(1, 1)*r(1) + cell%h_inv(1, 2)*r(2) + cell%h_inv(1, 3)*r(3)
460  s(2) = cell%h_inv(2, 1)*r(1) + cell%h_inv(2, 2)*r(2) + cell%h_inv(2, 3)*r(3)
461  s(3) = cell%h_inv(3, 1)*r(1) + cell%h_inv(3, 2)*r(2) + cell%h_inv(3, 3)*r(3)
462  s(1) = s(1) - cell%perd(1)*floor(s(1))
463  s(2) = s(2) - cell%perd(2)*floor(s(2))
464  s(3) = s(3) - cell%perd(3)*floor(s(3))
465  r_pbc(1) = cell%hmat(1, 1)*s(1) + cell%hmat(1, 2)*s(2) + cell%hmat(1, 3)*s(3)
466  r_pbc(2) = cell%hmat(2, 1)*s(1) + cell%hmat(2, 2)*s(2) + cell%hmat(2, 3)*s(3)
467  r_pbc(3) = cell%hmat(3, 1)*s(1) + cell%hmat(3, 2)*s(2) + cell%hmat(3, 3)*s(3)
468  END IF
469  ELSE
470  r_pbc = pbc1(r, cell)
471  END IF
472 
473  END FUNCTION pbc4
474 
475 ! **************************************************************************************************
476 !> \brief Transform real to scaled cell coordinates.
477 !> s=h_inv*r
478 !> \param s ...
479 !> \param r ...
480 !> \param cell ...
481 !> \date 16.01.2002
482 !> \author Matthias Krack
483 !> \version 1.0
484 ! **************************************************************************************************
485  SUBROUTINE real_to_scaled(s, r, cell)
486 
487  REAL(kind=dp), DIMENSION(3), INTENT(OUT) :: s
488  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: r
489  TYPE(cell_type), POINTER :: cell
490 
491  cpassert(ASSOCIATED(cell))
492 
493  IF (cell%orthorhombic) THEN
494  s(1) = cell%h_inv(1, 1)*r(1)
495  s(2) = cell%h_inv(2, 2)*r(2)
496  s(3) = cell%h_inv(3, 3)*r(3)
497  ELSE
498  s(1) = cell%h_inv(1, 1)*r(1) + cell%h_inv(1, 2)*r(2) + cell%h_inv(1, 3)*r(3)
499  s(2) = cell%h_inv(2, 1)*r(1) + cell%h_inv(2, 2)*r(2) + cell%h_inv(2, 3)*r(3)
500  s(3) = cell%h_inv(3, 1)*r(1) + cell%h_inv(3, 2)*r(2) + cell%h_inv(3, 3)*r(3)
501  END IF
502 
503  END SUBROUTINE real_to_scaled
504 
505 ! **************************************************************************************************
506 !> \brief Transform scaled cell coordinates real coordinates.
507 !> r=h*s
508 !> \param r ...
509 !> \param s ...
510 !> \param cell ...
511 !> \date 16.01.2002
512 !> \author Matthias Krack
513 !> \version 1.0
514 ! **************************************************************************************************
515  SUBROUTINE scaled_to_real(r, s, cell)
516 
517  REAL(kind=dp), DIMENSION(3), INTENT(OUT) :: r
518  REAL(kind=dp), DIMENSION(3), INTENT(IN) :: s
519  TYPE(cell_type), POINTER :: cell
520 
521  cpassert(ASSOCIATED(cell))
522 
523  IF (cell%orthorhombic) THEN
524  r(1) = cell%hmat(1, 1)*s(1)
525  r(2) = cell%hmat(2, 2)*s(2)
526  r(3) = cell%hmat(3, 3)*s(3)
527  ELSE
528  r(1) = cell%hmat(1, 1)*s(1) + cell%hmat(1, 2)*s(2) + cell%hmat(1, 3)*s(3)
529  r(2) = cell%hmat(2, 1)*s(1) + cell%hmat(2, 2)*s(2) + cell%hmat(2, 3)*s(3)
530  r(3) = cell%hmat(3, 1)*s(1) + cell%hmat(3, 2)*s(2) + cell%hmat(3, 3)*s(3)
531  END IF
532 
533  END SUBROUTINE scaled_to_real
534 ! **************************************************************************************************
535 !> \brief retains the given cell (see doc/ReferenceCounting.html)
536 !> \param cell the cell to retain
537 !> \par History
538 !> 09.2003 created [fawzi]
539 !> \author Fawzi Mohamed
540 ! **************************************************************************************************
541  SUBROUTINE cell_retain(cell)
542 
543  TYPE(cell_type), POINTER :: cell
544 
545  cpassert(ASSOCIATED(cell))
546  cpassert(cell%ref_count > 0)
547  cell%ref_count = cell%ref_count + 1
548 
549  END SUBROUTINE cell_retain
550 
551 ! **************************************************************************************************
552 !> \brief releases the given cell (see doc/ReferenceCounting.html)
553 !> \param cell the cell to release
554 !> \par History
555 !> 09.2003 created [fawzi]
556 !> \author Fawzi Mohamed
557 ! **************************************************************************************************
558  SUBROUTINE cell_release(cell)
559 
560  TYPE(cell_type), POINTER :: cell
561 
562  IF (ASSOCIATED(cell)) THEN
563  cpassert(cell%ref_count > 0)
564  cell%ref_count = cell%ref_count - 1
565  IF (cell%ref_count == 0) THEN
566  DEALLOCATE (cell)
567  END IF
568  NULLIFY (cell)
569  END IF
570 
571  END SUBROUTINE cell_release
572 
573 #if defined (__PLUMED2)
574 ! **************************************************************************************************
575 !> \brief For the interface with plumed, pass a cell pointer and retrieve it
576 !> later. It's a hack, but avoids passing the cell back and forth
577 !> across the Fortran/C++ interface
578 !> \param cell ...
579 !> \param set ...
580 !> \date 28.02.2013
581 !> \author RK
582 !> \version 1.0
583 ! **************************************************************************************************
584  SUBROUTINE pbc_cp2k_plumed_getset_cell(cell, set)
585 
586  TYPE(cell_type), POINTER :: cell
587  LOGICAL :: set
588 
589  TYPE(cell_type), POINTER, SAVE :: stored_cell
590 
591  IF (set) THEN
592  stored_cell => cell
593  ELSE
594  cell => stored_cell
595  END IF
596 
597  END SUBROUTINE pbc_cp2k_plumed_getset_cell
598 #endif
599 
600 END MODULE cell_types
subroutine pbc(r, r_pbc, s, s_pbc, a, b, c, alpha, beta, gamma, debug, info, pbc0, h, hinv)
...
Definition: dumpdcd.F:1203
Handles all functions related to the CELL.
Definition: cell_types.F:15
subroutine, public scaled_to_real(r, s, cell)
Transform scaled cell coordinates real coordinates. r=h*s.
Definition: cell_types.F:516
integer, parameter, public use_perd_xyz
Definition: cell_types.F:42
subroutine, public parse_cell_line(input_line, cell_itimes, cell_time, h, vol)
Read cell info from a line (parsed from a file)
Definition: cell_types.F:156
integer, parameter, public cell_sym_monoclinic
Definition: cell_types.F:29
integer, parameter, public use_perd_y
Definition: cell_types.F:42
integer, parameter, public cell_sym_triclinic
Definition: cell_types.F:29
integer, parameter, public cell_sym_tetragonal_ab
Definition: cell_types.F:29
integer, parameter, public use_perd_xz
Definition: cell_types.F:42
integer, parameter, public cell_sym_rhombohedral
Definition: cell_types.F:29
subroutine, public real_to_scaled(s, r, cell)
Transform real to scaled cell coordinates. s=h_inv*r.
Definition: cell_types.F:486
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
Definition: cell_types.F:559
integer, parameter, public use_perd_x
Definition: cell_types.F:42
subroutine, public cell_clone(cell_in, cell_out, tag)
Clone cell variable.
Definition: cell_types.F:107
integer, parameter, public cell_sym_tetragonal_ac
Definition: cell_types.F:29
integer, parameter, public use_perd_z
Definition: cell_types.F:42
integer, parameter, public use_perd_yz
Definition: cell_types.F:42
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition: cell_types.F:195
integer, parameter, public use_perd_none
Definition: cell_types.F:42
subroutine, public cell_retain(cell)
retains the given cell (see doc/ReferenceCounting.html)
Definition: cell_types.F:542
integer, parameter, public cell_sym_hexagonal_gamma_60
Definition: cell_types.F:29
integer, parameter, public cell_sym_orthorhombic
Definition: cell_types.F:29
integer, parameter, public cell_sym_none
Definition: cell_types.F:29
integer, parameter, public cell_sym_hexagonal_gamma_120
Definition: cell_types.F:29
subroutine, public cell_copy(cell_in, cell_out, tag)
Copy cell variable.
Definition: cell_types.F:126
integer, parameter, public cell_sym_monoclinic_gamma_ab
Definition: cell_types.F:29
integer, parameter, public cell_sym_cubic
Definition: cell_types.F:29
integer, parameter, public use_perd_xy
Definition: cell_types.F:42
integer, parameter, public cell_sym_tetragonal_bc
Definition: cell_types.F:29
real(kind=dp) function, public plane_distance(h, k, l, cell)
Calculate the distance between two lattice planes as defined by a triple of Miller indices (hkl).
Definition: cell_types.F:252
unit conversion facility
Definition: cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition: cp_units.F:1150
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Definition: gamma.F:15
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Definition of mathematical constants and functions.
Definition: mathconstants.F:16
real(kind=dp), parameter, public degree
Collection of simple mathematical functions and subroutines.
Definition: mathlib.F:15
pure real(kind=dp) function, public angle(a, b)
Calculation of the angle between the vectors a and b. The angle is returned in radians.
Definition: mathlib.F:175