Octopus
affine_coordinates.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
25 use debug_oct_m
26 use global_oct_m
28 use math_oct_m
32 use unit_oct_m
34
35 implicit none
36
37 private
38 public :: &
41
43 private
44 type(basis_vectors_t), public :: basis
45 contains
46 procedure :: to_cartesian => affine_coordinates_to_cartesian
47 procedure :: from_cartesian => affine_coordinates_from_cartesian
48 procedure :: dvector_from_cartesian => daffine_coordinates_vector_from_cartesian
49 procedure :: zvector_from_cartesian => zaffine_coordinates_vector_from_cartesian
50 procedure :: dcovector_to_cartesian => daffine_coordinates_covector_to_cartesian
51 procedure :: zcovector_to_cartesian => zaffine_coordinates_covector_to_cartesian
52 procedure :: det_jac => affine_coordinates_det_jac
53 procedure :: write_info => affine_coordinates_write_info
54 procedure :: surface_element => affine_coordinates_surface_element
56
57 interface affine_coordinates_t
58 procedure affine_coordinates_constructor
59 end interface affine_coordinates_t
60
61contains
62
63 function affine_coordinates_constructor(namespace, dim, basis_vectors) result(affine)
64 type(namespace_t), intent(in) :: namespace
65 integer, intent(in) :: dim
66 real(real64), intent(in) :: basis_vectors(1:dim, 1:dim)
67 class(affine_coordinates_t), pointer :: affine
68
70
71 safe_allocate(affine)
72
73 affine%dim = dim
74 affine%local_basis = .false.
75 affine%basis = basis_vectors_t(namespace, dim, basis_vectors)
76 affine%orthogonal = affine%basis%orthogonal
77
80
81 ! --------------------------------------------------------------
82 subroutine affine_coordinates_copy(this_out, this_in)
83 type(affine_coordinates_t), intent(inout) :: this_out
84 type(affine_coordinates_t), intent(in) :: this_in
85
87
88 this_out%dim = this_in%dim
89 this_out%local_basis = this_in%local_basis
90 this_out%basis = this_in%basis
91 this_out%orthogonal = this_in%orthogonal
92
94 end subroutine affine_coordinates_copy
95
96 ! ---------------------------------------------------------
97 pure function affine_coordinates_to_cartesian(this, chi) result(xx)
98 class(affine_coordinates_t), target, intent(in) :: this
99 real(real64), intent(in) :: chi(:)
100 real(real64) :: xx(1:this%dim)
101
102 ! no PUSH_SUB, called too often
103
104 xx(:) = this%basis%to_cartesian(chi(:))
105
107
108 ! ---------------------------------------------------------
109 pure function affine_coordinates_from_cartesian(this, xx) result(chi)
110 class(affine_coordinates_t), target, intent(in) :: this
111 real(real64), intent(in) :: xx(:)
112 real(real64) :: chi(1:this%dim)
113
114 ! no PUSH_SUB, called too often
116 chi(:) = this%basis%from_cartesian(xx(:))
117
119
120 ! ---------------------------------------------------------
121 real(real64) function affine_coordinates_det_jac(this, xx, chi) result(jdet)
122 class(affine_coordinates_t), intent(in) :: this
123 real(real64), intent(in) :: xx(:)
124 real(real64), intent(in) :: chi(:)
125
126 real(real64) :: jac(1:this%dim, 1:this%dim)
127
128 ! No PUSH_SUB, called too often
129
130 jac(1:this%dim, 1:this%dim) = this%basis%vectors(1:this%dim,1:this%dim)
131 jdet = lalg_determinant(this%dim, jac, preserve_mat = .false.)
132
133 end function affine_coordinates_det_jac
134
135 ! ---------------------------------------------------------
136 subroutine affine_coordinates_write_info(this, iunit, namespace)
137 class(affine_coordinates_t), intent(in) :: this
138 integer, optional, intent(in) :: iunit
139 type(namespace_t), optional, intent(in) :: namespace
141 integer :: idir1, idir2
145 write(message(1), '(a)') ' Using affine coordinates'
146 write(message(2), '(a,a,a)') ' Basis vectors [', trim(units_abbrev(units_out%length)), ']:'
147 call messages_info(2, iunit=iunit, namespace=namespace)
148 do idir1 = 1, this%dim
149 write(message(2), '(4x,a,99(f8.3,1x),a)') '(', &
150 (units_from_atomic(units_out%length, this%basis%vectors(idir2, idir1)), idir2 = 1, this%dim - 1), &
151 units_from_atomic(units_out%length, this%basis%vectors(this%dim, idir1)), ')'
152 call messages_info(1, iunit=iunit, namespace=namespace)
153 end do
154
157
158 ! ---------------------------------------------------------
159 real(real64) function affine_coordinates_surface_element(this, idir) result(ds)
160 class(affine_coordinates_t), intent(in) :: this
161 integer, intent(in) :: idir
162
164
165 select case (this%dim)
166 case (3)
167 select case (idir)
168 case (1)
169 ds = norm2(dcross_product(this%basis%vectors(1:3, 2), this%basis%vectors(1:3, 3)))
170 case (2)
171 ds = norm2(dcross_product(this%basis%vectors(1:3, 3), this%basis%vectors(1:3, 1)))
172 case (3)
173 ds = norm2(dcross_product(this%basis%vectors(1:3, 1), this%basis%vectors(1:3, 2)))
174 end select
176 case (2)
177 select case (idir)
178 case (3)
179 ds = this%basis%vectors(1, 1)*this%basis%vectors(2, 2) - this%basis%vectors(2, 1)*this%basis%vectors(1, 2)
180 case default
181 ! We can only get the surface element along z, as the only existing plane is the xy plane
182 assert(.false.)
183 end select
184
185 case default
186 ! We only know how to do this for 2D and 3D
187 assert(.false.)
188 end select
189
192
193#include "undef.F90"
194#include "real.F90"
195#include "affine_coordinates_inc.F90"
196
197#include "undef.F90"
198#include "complex.F90"
199#include "affine_coordinates_inc.F90"
200
203!! Local Variables:
204!! mode: f90
205!! coding: utf-8
206!! End:
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
Definition: lalg_adv.F90:189
pure subroutine zaffine_coordinates_vector_from_cartesian(this, xx, vv, src)
subroutine affine_coordinates_write_info(this, iunit, namespace)
pure subroutine zaffine_coordinates_covector_to_cartesian(this, xx, cv, src)
pure real(real64) function, dimension(1:this%dim) affine_coordinates_to_cartesian(this, chi)
pure subroutine daffine_coordinates_vector_from_cartesian(this, xx, vv, src)
class(affine_coordinates_t) function, pointer affine_coordinates_constructor(namespace, dim, basis_vectors)
subroutine, public affine_coordinates_copy(this_out, this_in)
pure real(real64) function, dimension(1:this%dim) affine_coordinates_from_cartesian(this, xx)
real(real64) function affine_coordinates_surface_element(this, idir)
pure subroutine daffine_coordinates_covector_to_cartesian(this, xx, cv, src)
real(real64) function affine_coordinates_det_jac(this, xx, chi)
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
Vectors defining a basis in a vector space. This class provides methods to convert vector coordinates...
abstract class to describe coordinate systems