28 use,
intrinsic :: iso_fortran_env
55 integer,
public,
parameter :: &
56 POISSON_FFT_KERNEL_NONE = -1, &
66 type(fourier_space_op_t) :: coulb
68 real(real64) :: soft_coulb_param
71 real(real64),
parameter :: TOL_VANISHING_Q = 1e-6_real64
74 subroutine poisson_fft_init(this, namespace, space, cube, kernel, soft_coulb_param, fullcube)
75 type(poisson_fft_t),
intent(out) :: this
76 type(namespace_t),
intent(in) :: namespace
77 class(space_t),
intent(in) :: space
78 type(cube_t),
intent(inout) :: cube
79 integer,
intent(in) :: kernel
80 real(real64),
optional,
intent(in) :: soft_coulb_param
81 type(cube_t),
optional,
intent(in) :: fullcube
88 safe_allocate(this%coulb%qq(1:space%dim))
89 this%coulb%qq(1:space%periodic_dim) =
m_zero
90 this%coulb%qq(space%periodic_dim+1:space%dim) = 1e-5_real64
91 this%coulb%singularity =
m_zero
102 type(namespace_t),
intent(in) :: namespace
103 class(space_t),
intent(in) :: space
104 type(cube_t),
intent(in) :: cube
105 type(fourier_space_op_t),
intent(inout) :: coulb
106 integer,
intent(in) :: kernel
107 real(real64),
optional,
intent(in) :: soft_coulb_param
108 type(cube_t),
optional,
intent(in) :: fullcube
114 message(1) =
"The screened Coulomb potential is only implemented in 3D for PoissonFFTKernel=fft_nocut."
121 if (.not.
present(fullcube))
then
122 message(1) =
"Hockney's FFT-kernel needs cube of full unit cell "
125 if (.not.
allocated(fullcube%fft))
then
126 message(1) =
"Hockney's FFT-kernel needs PoissonSolver=fft"
133 select case (space%dim)
135 assert(
present(soft_coulb_param))
142 message(1) =
"Invalid Poisson FFT kernel for 1D."
155 message(1) =
"Invalid Poisson FFT kernel for 2D."
177 message(1) =
"Invalid Poisson FFT kernel for 3D."
187 subroutine get_cutoff(namespace, default_r_c, r_c)
189 real(real64),
intent(in) :: default_r_c
190 real(real64),
intent(out) :: r_c
201 call messages_write(
'Poisson cutoff radius is larger than cell size.', new_line = .
true.)
202 call messages_write(
'You can see electrons in neighboring cell(s).')
212 type(
cube_t),
intent(in) :: cube
215 integer :: ix, iy, iz, ixx(3), db(3), n1, n2, n3, lx, ly, lz
216 real(real64) :: temp(3), modg2, inv_four_mu2, ecut
217 real(real64) :: gg(3)
218 real(real64) :: singularity_term
219 real(real64),
allocatable :: fft_coulb_fs(:,:,:)
223 db(1:3) = cube%rs_n_global(1:3)
231 n1 = max(1, cube%fs_n(1))
232 n2 = max(1, cube%fs_n(2))
233 n3 = max(1, cube%fs_n(3))
239 singularity_term =
m_four *
m_pi * inv_four_mu2 * coulb%beta
243 singularity_term = coulb%singularity*coulb%alpha +
m_four*
m_pi*inv_four_mu2 * coulb%beta
245 singularity_term = coulb%singularity
251 safe_allocate(fft_coulb_fs(1:n1, 1:n2, 1:n3))
254 temp(1:3) =
m_two*
m_pi/(db(1:3)*cube%spacing(1:3))
263 iz = cube%fs_istart(3) + lz - 1
265 iy = cube%fs_istart(2) + ly - 1
268 ix = cube%fs_istart(1) + lx - 1
274 if(modg2 >
m_two * ecut * 1.001_real64) cycle
277 if (cube%fft%library.eq.
fftlib_nfft) modg2=cube%Lfs(ix,1)**2 + cube%Lfs(iy,2)**2 + cube%Lfs(iz,3)**2
279 if (modg2 > tol_vanishing_q)
then
283 fft_coulb_fs(lx, ly, lz) =
m_four *
m_pi / modg2 * (coulb%alpha + coulb%beta *
exp(-modg2*inv_four_mu2))
285 fft_coulb_fs(lx, ly, lz) =
m_four *
m_pi / modg2 * coulb%beta * (
m_one-
exp(-modg2 * inv_four_mu2))
289 fft_coulb_fs(lx, ly, lz) =
m_four *
m_pi / modg2 * coulb%alpha
291 fft_coulb_fs(lx, ly, lz) =
m_four *
m_pi / modg2
295 fft_coulb_fs(lx, ly, lz) = singularity_term
304 safe_deallocate_a(fft_coulb_fs)
315 type(
cube_t),
intent(in) :: cube
317 type(
cube_t),
intent(in) :: fullcube
319 integer :: ix, iy, iz, ixx(3), db(3), nfs(3), nrs(3), nfs_s(3), nrs_s(3), dnrs(3)
320 real(real64) :: temp(3), modg2, weight
321 real(real64) :: gg(3)
322 real(real64),
allocatable :: fft_Coulb_small_RS(:,:,:)
323 real(real64),
allocatable :: fft_Coulb_RS(:,:,:)
324 complex(real64),
allocatable :: fft_Coulb_small_FS(:,:,:), fft_Coulb_FS(:,:,:)
331 nfs(1:3) = fullcube%fs_n_global(1:3)
332 nrs(1:3) = fullcube%rs_n_global(1:3)
334 safe_allocate(fft_coulb_fs(1:nfs(1),1:nfs(2),1:nfs(3)))
335 safe_allocate(fft_coulb_rs(1:nrs(1),1:nrs(2),1:nrs(3)))
338 nfs_s(1:3) = cube%fs_n_global(1:3)
339 nrs_s(1:3) = cube%rs_n_global(1:3)
341 safe_allocate(fft_coulb_small_fs(1:nfs_s(1),1:nfs_s(2),1:nfs_s(3)))
342 safe_allocate(fft_coulb_small_rs(1:nrs_s(1),1:nrs_s(2),1:nrs_s(3)))
347 db(1:3) = fullcube%rs_n_global(1:3)
348 temp(1:3) =
m_two*
m_pi/(db(1:3)*cube%spacing(1:3))
359 if (abs(modg2) > tol_vanishing_q)
then
360 fft_coulb_fs(ix, iy, iz) =
m_one/modg2
362 fft_coulb_fs(ix, iy, iz) =
m_zero
370 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
375 fft_coulb_fs(ix, iy, iz) = weight*fft_coulb_fs(ix, iy, iz)
390 if (iz > nrs_s(3)/2+1) ixx(3) = ixx(3) + dnrs(3)
393 if (iy > nrs_s(2)/2+1) ixx(2) = ixx(2) + dnrs(2)
396 if (ix > nrs_s(1)/2+1) ixx(1) = ixx(1) + dnrs(1)
397 fft_coulb_small_rs(ix, iy, iz) = fft_coulb_rs(ixx(1),ixx(2),ixx(3))
402 call dfft_forward(cube%fft, fft_coulb_small_rs, fft_coulb_small_fs)
404 fft_coulb_small_rs(1:nfs_s(1),1:nfs_s(2),1:nfs_s(3)) = &
405 real( fft_Coulb_small_FS(1:nfs_s(1),1:nfs_s(2),1:nfs_s(3)), real64)
411 fft_coulb_small_rs(cube%fs_istart(1):cube%fs_istart(1)+cube%fs_n(1), &
412 cube%fs_istart(2):cube%fs_istart(2)+cube%fs_n(2), &
413 cube%fs_istart(3):cube%fs_istart(3)+cube%fs_n(3)))
415 safe_deallocate_a(fft_coulb_fs)
416 safe_deallocate_a(fft_coulb_rs)
417 safe_deallocate_a(fft_coulb_small_fs)
418 safe_deallocate_a(fft_coulb_small_rs)
428 type(
cube_t),
intent(in) :: cube
431 integer :: ix, iy, iz, ixx(3), db(3)
432 integer :: lx, ly, lz, n1, n2, n3
433 real(real64) :: temp(3), modg2, ecut, weight
434 real(real64) :: gpar, gz, r_c, gg(3), default_r_c
435 real(real64),
allocatable :: fft_coulb_FS(:,:,:)
439 db(1:3) = cube%rs_n_global(1:3)
444 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
458 default_r_c = db(3)*cube%spacing(3)/
m_two
461 n1 = max(1, cube%fs_n(1))
462 n2 = max(1, cube%fs_n(2))
463 n3 = max(1, cube%fs_n(3))
465 safe_allocate(fft_coulb_fs(1:n1, 1:n2, 1:n3))
468 temp(1:3) =
m_two*
m_pi/(db(1:3)*cube%spacing(1:3))
473 iz = cube%fs_istart(3) + lz - 1
476 iy = cube%fs_istart(2) + ly - 1
479 ix = cube%fs_istart(1) + lx - 1
484 if(sum(gg(1:2)**2) >
m_two*ecut*1.001_real64) cycle
486 if (abs(modg2) > tol_vanishing_q)
then
488 gpar =
hypot(gg(1), gg(2))
492 fft_coulb_fs(lx, ly, lz) = -
m_half*r_c**2
494 fft_coulb_fs(lx, ly, lz) = weight*fft_coulb_fs(lx, ly, lz)
502 safe_deallocate_a(fft_coulb_fs)
512 class(
space_t),
intent(in) :: space
513 type(
cube_t),
intent(in) :: cube
517 real(real64),
allocatable :: x(:), y(:)
518 integer :: ix, iy, iz, ixx(3), db(3), k, ngp
519 integer :: lx, ly, lz, n1, n2, n3, lxx(3)
520 real(real64) :: temp(3), modg2, xmax, weight
521 real(real64) :: gperp, gx, gy, gz, r_c, gg(3), default_r_c
522 real(real64),
allocatable :: fft_coulb_FS(:,:,:)
529 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
532 db(1:3) = cube%rs_n_global(1:3)
534 default_r_c = maxval(db(2:3)*cube%spacing(2:3)/
m_two)
537 n1 = max(1, cube%fs_n(1))
538 n2 = max(1, cube%fs_n(2))
539 n3 = max(1, cube%fs_n(3))
541 safe_allocate(fft_coulb_fs(1:n1, 1:n2, 1:n3))
544 temp(1:3) =
m_two*
m_pi/(db(1:3)*cube%spacing(1:3))
546 if (.not. space%is_periodic())
then
548 safe_allocate(x(1:ngp))
549 safe_allocate(y(1:ngp))
555 ix = cube%fs_istart(1) + lx - 1
557 lxx(1) = ixx(1) - cube%fs_istart(1) + 1
560 if (.not. space%is_periodic())
then
562 xmax = norm2(temp(2:3)*db(2:3))/2
564 x(k) = (k-1)*(xmax/(ngp-1))
566 maxval(norm2(cube%latt%rlattice_primitive(:, 2:3), dim=1)))
572 iy = cube%fs_istart(2) + ly - 1
574 lxx(2) = ixx(2) - cube%fs_istart(2) + 1
576 iz = cube%fs_istart(3) + lz - 1
578 lxx(3) = ixx(3) - cube%fs_istart(3) + 1
582 if (abs(modg2) > tol_vanishing_q)
then
583 gperp =
hypot(gg(2), gg(3))
584 if (space%periodic_dim == 1)
then
585 if (gperp > r_c)
then
586 fft_coulb_fs(lx, ly, lz) =
m_zero
590 else if (.not. space%is_periodic())
then
594 fft_coulb_fs(lx, ly, lz) =
spline_eval(cylinder_cutoff_f, gperp)
597 fft_coulb_fs(lx, ly, lz) = fft_coulb_fs(lx, -lxx(2) + 1, lz)
600 fft_coulb_fs(lx, ly, lz) = fft_coulb_fs(lx, ly, -lxx(3) + 1)
603 fft_coulb_fs(lx, ly, lz) = fft_coulb_fs(lx, -lxx(2) + 1, -lxx(3) + 1)
608 if (space%periodic_dim == 1)
then
610 else if (.not. space%is_periodic())
then
612 norm2(cube%latt%rlattice_primitive(:, 1)), maxval(norm2(cube%latt%rlattice_primitive(:, 2:3), dim=1)))
616 fft_coulb_fs(lx, ly, lz) = weight*fft_coulb_fs(lx, ly, lz)
620 if (.not. space%is_periodic())
then
627 safe_deallocate_a(fft_coulb_fs)
639 type(
cube_t),
intent(in) :: cube
640 integer,
intent(in) :: kernel
642 logical,
intent(in) :: is_periodic
644 integer :: ix, iy, iz, ixx(3), db(3), lx, ly, lz, n1, n2, n3
645 real(real64) :: temp(3), modg2, ecut, weight
646 real(real64) :: r_c, gg(3), default_r_c
647 real(real64),
allocatable :: fft_coulb_FS(:,:,:)
648 real(real64) :: axis(3,3)
655 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
658 db(1:3) = cube%rs_n_global(1:3)
664 axis(:,ix) = cube%latt%rlattice_primitive(:, ix) * cube%spacing(ix) * db(ix) /
m_two
673 if (.not. cube%latt%nonorthogonal)
then
674 temp(1:3) = axis(:, ix)
675 default_r_c = min(default_r_c, norm2(temp(1:3)))
684 temp = temp / norm2(temp)
685 default_r_c = min(default_r_c, dot_product(temp, axis(:, ix)-axis(:, iy)))
691 n1 = max(1, cube%fs_n(1))
692 n2 = max(1, cube%fs_n(2))
693 n3 = max(1, cube%fs_n(3))
697 safe_allocate(fft_coulb_fs(1:n1,1:n2,1:n3))
700 temp(1:3) =
m_two*
m_pi/(db(1:3)*cube%spacing(1:3))
705 iz = cube%fs_istart(3) + lz - 1
708 iy = cube%fs_istart(2) + ly - 1
711 ix = cube%fs_istart(1) + lx - 1
718 modg2=cube%Lfs(ix,1)**2+cube%Lfs(iy,2)**2+cube%Lfs(iz,3)**2
719 r_c = default_r_c*
m_two
723 if(modg2 >
m_two*ecut*1.001_real64 .and. is_periodic) cycle
725 if (abs(modg2) > tol_vanishing_q)
then
730 fft_coulb_fs(lx, ly, lz) = weight/modg2
735 fft_coulb_fs(lx, ly, lz) = weight*r_c**2/
m_two
737 fft_coulb_fs(lx, ly, lz) =
m_zero
746 safe_deallocate_a(fft_coulb_fs)
756 type(
cube_t),
intent(in) :: cube
760 integer :: i, ix, iy, ixx(2), db(2), npoints
761 real(real64) :: temp(2), vec, r_c, maxf, dk, default_r_c, weight
762 real(real64),
allocatable :: x(:), y(:)
763 real(real64),
allocatable :: fft_coulb_FS(:,:,:)
770 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
772 db(1:2) = cube%rs_n_global(1:2)
774 default_r_c = maxval(db(1:2)*cube%spacing(1:2)/
m_two)
780 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
782 temp(1:2) =
m_two*
m_pi/(db(1:2)*cube%spacing(1:2))
784 maxf = r_c * norm2(temp(1:2)*db(1:2))/2
786 npoints = nint(maxf/dk)
787 safe_allocate(x(1:npoints))
788 safe_allocate(y(1:npoints))
792 x(i) = (i-1) * maxf / (npoints-1)
797 do iy = 1, cube%fs_n_global(2)
799 do ix = 1, cube%fs_n_global(1)
801 vec = norm2(temp(1:2)*ixx(1:2))
803 if (vec*r_c >= x(npoints))
then
804 fft_coulb_fs(ix, iy, 1) = weight * y(npoints)
805 else if (vec >
m_zero)
then
808 fft_coulb_fs(ix, iy, 1) = weight *
m_two *
m_pi * r_c
815 safe_deallocate_a(fft_coulb_fs)
828 type(
cube_t),
intent(in) :: cube
831 integer :: ix, iy, ixx(2), db(2)
832 real(real64) :: temp(2), r_c, gx, gy, default_r_c, weight
833 real(real64),
allocatable :: fft_coulb_FS(:,:,:)
840 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
843 db(1:2) = cube%rs_n_global(1:2)
845 default_r_c = db(2)*cube%spacing(2)/
m_two
849 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
851 temp(1:2) =
m_two*
m_pi/(db(1:2)*cube%spacing(1:2))
855 do iy = 2, cube%fs_n_global(2)
864 do ix = 2, cube%fs_n_global(1)
873 safe_deallocate_a(fft_coulb_fs)
883 type(
cube_t),
intent(in) :: cube
886 integer :: ix, iy, ixx(2), db(2)
887 real(real64) :: temp(2), vec, weight
888 real(real64),
allocatable :: fft_coulb_FS(:,:,:)
895 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
897 db(1:2) = cube%rs_n_global(1:2)
900 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
902 temp(1:2) =
m_two*
m_pi/(db(1:2)*cube%spacing(1:2))
904 do iy = 1, cube%fs_n_global(2)
906 do ix = 1, cube%fs_n_global(1)
908 vec =
sqrt((temp(1) * ixx(1))**2 + (temp(2) * ixx(2))**2)
909 if (vec >
m_zero) fft_coulb_fs(ix, iy, 1) =
m_two *
m_pi / vec * weight
915 safe_deallocate_a(fft_coulb_fs)
923 type(
cube_t),
intent(in) :: cube
925 real(real64),
intent(in) :: poisson_soft_coulomb_param
929 real(real64),
allocatable :: fft_coulb_fs(:, :, :), weight
936 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
938 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
942 do ix = 1, cube%fs_n_global(1)
944 g = (ixx + coulb%qq(1))*
m_two*
m_pi/abs(cube%latt%rlattice(1,1))
945 if (abs(g) > tol_vanishing_q)
then
946 fft_coulb_fs(ix, 1, 1) =
m_two *
loct_bessel_k0(poisson_soft_coulomb_param*abs(g)) * weight
948 fft_coulb_fs(ix, 1, 1) = coulb%singularity *
m_two * weight
953 safe_deallocate_a(fft_coulb_fs)
963 type(
cube_t),
intent(in) :: cube
965 real(real64),
intent(in) :: poisson_soft_coulomb_param
967 integer :: box(1), ixx(1), ix
968 real(real64) :: temp(1), g, r_c, default_r_c, weight
969 real(real64),
allocatable :: fft_coulb_fs(:, :, :)
976 if(coulb%alpha >
m_epsilon) weight = weight * coulb%alpha
978 box(1:1) = cube%rs_n_global(1:1)
980 default_r_c = box(1)*cube%spacing(1)/
m_two
983 safe_allocate(fft_coulb_fs(1:cube%fs_n_global(1), 1:cube%fs_n_global(2), 1:cube%fs_n_global(3)))
985 temp(1:1) =
m_two*
m_pi/(box(1:1)*cube%spacing(1:1))
988 do ix = 1, cube%fs_n_global(1)
995 safe_deallocate_a(fft_coulb_fs)
1015#include "poisson_fft_inc.F90"
1017#include "complex.F90"
1018#include "poisson_fft_inc.F90"
Some operations may be done for one spline-function, or for an array of them.
double hypot(double __x, double __y) __attribute__((__nothrow__
double log(double __x) __attribute__((__nothrow__
double exp(double __x) __attribute__((__nothrow__
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
real(real64) function, public fft_get_ecut_from_box(box_dim, fs_istart, latt, gspacing, periodic_dim, qq)
Given an fft box (fixed by the real-space grid), it returns the cutoff energy of the sphere that fits...
pure integer function, public pad_feq(ii, nn, mode)
convert between array index and G-vector
integer, parameter, public fftlib_nfft
pure subroutine, public fft_gg_transform(gg_in, temp, periodic_dim, latt, qq, gg, modg2)
Convert FFT grid index into the Cartesian reciprocal-space vector .
subroutine, public fourier_space_op_end(this)
subroutine, public dfourier_space_op_init(this, cube, op, in_device)
real(real64), parameter, public m_two
real(real64), parameter, public m_huge
real(real64), parameter, public m_zero
real(real64), parameter, public m_four
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_fourth
real(real64), parameter, public m_epsilon
real(real64), parameter, public m_half
real(real64), parameter, public m_one
This module is intended to contain "only mathematical" functions and procedures.
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
This module defines the meshes, which are used in Octopus.
subroutine, public messages_warning(no_lines, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
real(real64) function, public poisson_cutoff_3d_2d(p, z, r)
real(real64) function, public poisson_cutoff_3d_1d(x, p, rmax)
real(real64) function, public poisson_cutoff_3d_0d(x, r)
integer, parameter, public poisson_fft_kernel_hockney
subroutine poisson_fft_build_2d_0d(namespace, cube, coulb)
A. Castro et al., Phys. Rev. B 80, 033102 (2009)
subroutine poisson_fft_build_3d_1d(namespace, space, cube, coulb)
C. A. Rozzi et al., Phys. Rev. B 73, 205119 (2006), Table I.
subroutine poisson_fft_build_3d_2d(namespace, cube, coulb)
C. A. Rozzi et al., Phys. Rev. B 73, 205119 (2006), Table I.
integer, parameter, public poisson_fft_kernel_nocut
integer, parameter, public poisson_fft_kernel_cyl
subroutine poisson_fft_build_2d_1d(namespace, cube, coulb)
A. Castro et al., Phys. Rev. B 80, 033102 (2009)
subroutine poisson_fft_build_1d_0d(namespace, cube, coulb, poisson_soft_coulomb_param)
subroutine, public zpoisson_fft_solve(this, mesh, cube, pot, rho, mesh_cube_map, average_to_zero, kernel, sm)
subroutine poisson_fft_build_3d_0d(namespace, cube, kernel, coulb, is_periodic)
C. A. Rozzi et al., Phys. Rev. B 73, 205119 (2006), Table I.
subroutine get_cutoff(namespace, default_r_c, r_c)
subroutine poisson_fft_build_3d_3d(cube, coulb)
Compute the Coulomb kernel in reciprocal space, for a 3D FFT grid.
subroutine poisson_fft_build_1d_1d(cube, coulb, poisson_soft_coulomb_param)
subroutine, public poisson_fft_get_kernel(namespace, space, cube, coulb, kernel, soft_coulb_param, fullcube)
subroutine, public poisson_fft_end(this)
subroutine poisson_fft_build_2d_2d(cube, coulb)
A. Castro et al., Phys. Rev. B 80, 033102 (2009)
subroutine, public poisson_fft_init(this, namespace, space, cube, kernel, soft_coulb_param, fullcube)
integer, parameter, public poisson_fft_kernel_pla
integer, parameter, public poisson_fft_kernel_corrected
integer, parameter, public poisson_fft_kernel_sph
subroutine, public dpoisson_fft_solve(this, mesh, cube, pot, rho, mesh_cube_map, average_to_zero, kernel, sm)
subroutine poisson_fft_build_3d_3d_hockney(cube, coulb, fullcube)
Kernel for Hockneys algorithm that solves the poisson equation in a small box while respecting the pe...
subroutine, public spline_fit(nrc, rofi, ffit, spl, threshold)
real(real64) function, public spline_eval(spl, x)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
the basic spline datatype