Octopus
cartesian.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
21module cartesian_oct_m
24 use debug_oct_m
25 use global_oct_m
26 use math_oct_m
30
31 implicit none
32
33 private
34 public :: &
37
38 type, extends(affine_coordinates_t) :: cartesian_t
39 private
40 contains
41 procedure :: to_cartesian => cartesian_to_cartesian
42 procedure :: from_cartesian => cartesian_from_cartesian
43 procedure :: dvector_from_cartesian => dcartesian_vector_from_cartesian
44 procedure :: zvector_from_cartesian => zcartesian_vector_from_cartesian
45 procedure :: dcovector_to_cartesian => dcartesian_covector_to_cartesian
46 procedure :: zcovector_to_cartesian => zcartesian_covector_to_cartesian
47 procedure :: det_jac => cartesian_det_jac
48 procedure :: write_info => cartesian_write_info
49 end type cartesian_t
50
51 interface cartesian_t
52 procedure cartesian_constructor
53 end interface cartesian_t
54
55contains
56
57 ! ---------------------------------------------------------
58 function cartesian_constructor(namespace, dim) result(cart)
59 type(namespace_t), intent(in) :: namespace
60 integer, intent(in) :: dim
61 class(cartesian_t), pointer :: cart
62
63 push_sub(cartesian_constructor)
64
65 safe_allocate(cart)
66
67 cart%dim = dim
68 cart%local_basis = .false.
69 cart%orthogonal = .true.
70 cart%basis = basis_vectors_t(namespace, dim, diagonal_matrix(dim, m_one))
71
73 end function cartesian_constructor
74
75 ! --------------------------------------------------------------
76 subroutine cartesian_copy(this_out, this_in)
77 type(cartesian_t), intent(inout) :: this_out
78 type(cartesian_t), intent(in) :: this_in
79
80 push_sub(cartesian_copy)
81
82 this_out%dim = this_in%dim
83 this_out%local_basis = this_in%local_basis
84
85 pop_sub(cartesian_copy)
86 end subroutine cartesian_copy
87
88 ! ---------------------------------------------------------
89 pure function cartesian_to_cartesian(this, chi) result(xx)
90 class(cartesian_t), target, intent(in) :: this
91 real(real64), intent(in) :: chi(:)
92 real(real64) :: xx(1:this%dim)
93
94 ! no PUSH_SUB, called too often
95
96 xx = chi
97
98 end function cartesian_to_cartesian
99
100 ! ---------------------------------------------------------
101 pure function cartesian_from_cartesian(this, xx) result(chi)
102 class(cartesian_t), target, intent(in) :: this
103 real(real64), intent(in) :: xx(:)
104 real(real64) :: chi(1:this%dim)
105
106 ! no PUSH_SUB, called too often
107
108 chi = xx
109
110 end function cartesian_from_cartesian
111
112 ! ---------------------------------------------------------
113 pure subroutine dcartesian_vector_from_cartesian(this, xx, vv, src)
114 class(cartesian_t), intent(in) :: this
115 real(real64), intent(in) :: xx(:)
116 real(real64), intent(inout) :: vv(:)
117 real(real64), optional, intent(in) :: src(:)
118
119 ! no PUSH_SUB, called too often
120
121 ! We are already in Cartesian coordinates, so nothing to do
122 if (present(src)) then
123 vv = src
124 end if
125
127
128 ! ---------------------------------------------------------
129 pure subroutine zcartesian_vector_from_cartesian(this, xx, vv, src)
130 class(cartesian_t), intent(in) :: this
131 real(real64), intent(in) :: xx(:)
132 complex(real64), intent(inout) :: vv(:)
133 complex(real64), optional, intent(in) :: src(:)
135 ! no PUSH_SUB, called too often
137 ! We are already in Cartesian coordinates, so nothing to do
138 if (present(src)) then
139 vv = src
140 end if
143
144 ! ---------------------------------------------------------
145 pure subroutine dcartesian_covector_to_cartesian(this, xx, cv, src)
146 class(cartesian_t), intent(in) :: this
147 real(real64), intent(in) :: xx(:)
148 real(real64), intent(inout) :: cv(:)
149 real(real64), optional, intent(in) :: src(:)
150
151 ! no PUSH_SUB, called too often
152
153 ! We are already in Cartesian coordinates, so nothing to do
154 if (present(src)) then
155 cv = src
156 end if
157
159
160 ! ---------------------------------------------------------
161 pure subroutine zcartesian_covector_to_cartesian(this, xx, cv, src)
162 class(cartesian_t), intent(in) :: this
163 real(real64), intent(in) :: xx(:)
164 complex(real64), intent(inout) :: cv(:)
165 complex(real64), optional, intent(in) :: src(:)
166
167 ! no PUSH_SUB, called too often
168
169 ! We are already in Cartesian coordinates, so nothing to do
170 if (present(src)) then
171 cv = src
172 end if
173
175
176 ! ---------------------------------------------------------
177 pure real(real64) function cartesian_det_jac(this, xx, chi) result(jdet)
178 class(cartesian_t), intent(in) :: this
179 real(real64), intent(in) :: xx(:)
180 real(real64), intent(in) :: chi(:)
181
182 ! No PUSH_SUB, called too often
183
184 jdet = m_one
185
186 end function cartesian_det_jac
187
188 ! ---------------------------------------------------------
189 subroutine cartesian_write_info(this, iunit, namespace)
190 class(cartesian_t), intent(in) :: this
191 integer, optional, intent(in) :: iunit
192 type(namespace_t), optional, intent(in) :: namespace
193
195
196 write(message(1), '(a)') ' Using Cartesian coordinates'
197 call messages_info(1, iunit=iunit, namespace=namespace)
198
199 pop_sub(cartesian_write_info)
200 end subroutine cartesian_write_info
201
202end module cartesian_oct_m
203
204!! Local Variables:
205!! mode: f90
206!! coding: utf-8
207!! End:
pure real(real64) function cartesian_det_jac(this, xx, chi)
Definition: cartesian.F90:271
pure real(real64) function, dimension(1:this%dim) cartesian_from_cartesian(this, xx)
Definition: cartesian.F90:195
class(cartesian_t) function, pointer cartesian_constructor(namespace, dim)
Definition: cartesian.F90:152
subroutine, public cartesian_copy(this_out, this_in)
Definition: cartesian.F90:170
pure subroutine dcartesian_covector_to_cartesian(this, xx, cv, src)
Definition: cartesian.F90:239
pure subroutine zcartesian_vector_from_cartesian(this, xx, vv, src)
Definition: cartesian.F90:223
pure subroutine dcartesian_vector_from_cartesian(this, xx, vv, src)
Definition: cartesian.F90:207
pure real(real64) function, dimension(1:this%dim) cartesian_to_cartesian(this, chi)
Definition: cartesian.F90:183
subroutine cartesian_write_info(this, iunit, namespace)
Definition: cartesian.F90:283
pure subroutine zcartesian_covector_to_cartesian(this, xx, cv, src)
Definition: cartesian.F90:255
real(real64), parameter, public m_one
Definition: global.F90:188
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
Vectors defining a basis in a vector space. This class provides methods to convert vector coordinates...
int true(void)