Octopus
coordinate_system.F90
Go to the documentation of this file.
1!! Copyright (C) 2021 M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use, intrinsic :: iso_fortran_env
24 implicit none
25
26 private
27 public :: coordinate_system_t
28
35 type, abstract :: coordinate_system_t
36 logical :: local_basis
38 logical :: orthogonal
39 integer :: dim
40 real(real64) :: min_mesh_scaling_product
41 contains
42 generic :: vector_from_cartesian => dvector_from_cartesian, zvector_from_cartesian
43
44 procedure :: dvector_from_cartesian => dcoordinate_system_vector_from_cartesian
45
46 procedure :: zvector_from_cartesian => zcoordinate_system_vector_from_cartesian
47
48 generic :: covector_to_cartesian => dcovector_to_cartesian, zcovector_to_cartesian
49
50 procedure :: dcovector_to_cartesian => dcoordinate_system_covector_to_cartesian
51
52 procedure :: zcovector_to_cartesian => zcoordinate_system_covector_to_cartesian
53
54 procedure(coordinate_system_to_cartesian), deferred :: to_cartesian
55
56 procedure(coordinate_system_from_cartesian), deferred :: from_cartesian
57
58 procedure(coordinate_system_det_jac), deferred :: det_jac
59
60 procedure(coordinate_system_write_info), deferred :: write_info
61
62 procedure(coordinates_surface_element), deferred :: surface_element
63
64 end type coordinate_system_t
65
66 abstract interface
67 ! ---------------------------------------------------------
70 function coordinate_system_to_cartesian(this, chi) result(xx)
72 import real64
73 class(coordinate_system_t), target, intent(in) :: this
74 real(real64), intent(in) :: chi(:)
75 real(real64) :: xx(1:this%dim)
77
78 ! ---------------------------------------------------------
80 function coordinate_system_from_cartesian(this, xx) result(chi)
82 import real64
83 class(coordinate_system_t), target, intent(in) :: this
84 real(real64), intent(in) :: xx(:)
85 real(real64) :: chi(1:this%dim)
87
88 ! ---------------------------------------------------------
89 real(real64) function coordinate_system_det_jac(this, xx, chi) result(jdet)
91 import real64
92 class(coordinate_system_t), intent(in) :: this
93 real(real64), intent(in) :: xx(:)
94 real(real64), intent(in) :: chi(:)
95 end function coordinate_system_det_jac
96
97 ! ---------------------------------------------------------
98 subroutine coordinate_system_write_info(this, iunit, namespace)
100 import namespace_t
101 class(coordinate_system_t), intent(in) :: this
102 integer, optional, intent(in) :: iunit
103 type(namespace_t), optional, intent(in) :: namespace
104 end subroutine coordinate_system_write_info
105
106 ! ---------------------------------------------------------
107 real(real64) function coordinates_surface_element(this, idir) result(ds)
109 import real64
110 class(coordinate_system_t), intent(in) :: this
111 integer, intent(in) :: idir
112 end function coordinates_surface_element
113 end interface
115contains
116
117 ! ---------------------------------------------------------
118 subroutine dcoordinate_system_vector_from_cartesian(this, xx, vv, src)
119 class(coordinate_system_t), intent(in) :: this
120 real(real64), intent(in) :: xx(:)
121 real(real64), intent(inout) :: vv(:)
122 real(real64), optional, intent(in) :: src(:)
123
124 ! Classes that want to provide this method must override it.
125 assert(.false.)
126
129 ! ---------------------------------------------------------
130 subroutine zcoordinate_system_vector_from_cartesian(this, xx, vv, src)
131 class(coordinate_system_t), intent(in) :: this
132 real(real64), intent(in) :: xx(:)
133 complex(real64), intent(inout) :: vv(:)
134 complex(real64), optional, intent(in) :: src(:)
136 ! Classes that want to provide this method must override it.
137 ! Note that when the basis is local, this transformantion depends on the
138 ! position
139 assert(.false.)
140
142
143 ! ---------------------------------------------------------
144 subroutine dcoordinate_system_covector_to_cartesian(this, xx, cv, src)
145 class(coordinate_system_t), intent(in) :: this
146 real(real64), intent(in) :: xx(:)
147 real(real64), intent(inout) :: cv(:)
148 real(real64), optional, intent(in) :: src(:)
150 ! Classes that want to provide this method must override it.
151 assert(.false.)
152
154
155 ! ---------------------------------------------------------
156 subroutine zcoordinate_system_covector_to_cartesian(this, xx, cv, src)
157 class(coordinate_system_t), intent(in) :: this
158 real(real64), intent(in) :: xx(:)
159 complex(real64), intent(inout) :: cv(:)
160 complex(real64), optional, intent(in) :: src(:)
161
162 ! Classes that want to provide this method must override it.
163 assert(.false.)
164
166
168
169!! Local Variables:
170!! mode: f90
171!! coding: utf-8
172!! End:
Convert Cartesian coordinates to coordinates in this coordinate system.
Convert coordinates given in this coordinate system to Cartesian coordinates.
subroutine dcoordinate_system_covector_to_cartesian(this, xx, cv, src)
subroutine zcoordinate_system_covector_to_cartesian(this, xx, cv, src)
subroutine dcoordinate_system_vector_from_cartesian(this, xx, vv, src)
subroutine zcoordinate_system_vector_from_cartesian(this, xx, vv, src)
abstract class to describe coordinate systems