Octopus
curv_briggs.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
25
28 use debug_oct_m
29 use global_oct_m
32 use parser_oct_m
34
35 implicit none
36
37 private
38 public :: &
41
42 type, extends(coordinate_system_t) :: curv_briggs_t
43 private
44 real(real64), allocatable :: lsize(:)
45 real(real64) :: beta
46 contains
47 procedure :: to_cartesian => curv_briggs_to_cartesian
48 procedure :: from_cartesian => curv_briggs_from_cartesian
49 procedure :: write_info => curv_briggs_write_info
50 procedure :: surface_element => curv_briggs_surface_element
51 procedure :: jacobian => curv_briggs_jacobian
52 procedure :: jacobian_inverse => curv_briggs_jacobian_inverse
54 end type curv_briggs_t
55
56 interface curv_briggs_t
57 procedure curv_briggs_constructor
58 end interface curv_briggs_t
59
60contains
61
62 ! ---------------------------------------------------------
63 function curv_briggs_constructor(namespace, dim, lsize, spacing) result(briggs)
64 type(namespace_t), intent(in) :: namespace
65 integer, intent(in) :: dim
66 real(real64), intent(in) :: lsize(1:dim)
67 real(real64), intent(in) :: spacing(1:dim)
68 class(curv_briggs_t), pointer :: briggs
69
70 integer :: idim
71
73
74 safe_allocate(briggs)
75
76 briggs%local_basis = .true.
77 briggs%orthogonal = .true. ! This needs to be checked
78 briggs%dim = dim
79 safe_allocate(briggs%lsize(1:dim))
80 briggs%lsize(1:dim) = lsize(1:dim)
81
82 call parse_variable(namespace, 'CurvBriggsBeta', m_half, briggs%beta)
83
84 if (briggs%beta < m_zero .or. briggs%beta > m_one) then
85 message(1) = 'The parameter "CurvBriggsBeta" must lie between 0 and 1.'
86 call messages_fatal(1, namespace=namespace)
87 end if
88
89 briggs%min_mesh_scaling_product = m_one
90 do idim = 1, briggs%dim
91 ! corresponds to the distance of grid points at [+spacing/2,-spacing/2]
92 briggs%min_mesh_scaling_product = briggs%min_mesh_scaling_product * (m_one / &
93 (m_one - briggs%lsize(idim) * briggs%beta / (m_pi * spacing(idim)) * sin(m_pi * spacing(idim) / briggs%lsize(idim))))
94 end do
95
97 end function curv_briggs_constructor
98
99 ! ---------------------------------------------------------
100 subroutine curv_briggs_finalize(this)
101 type(curv_briggs_t), intent(inout) :: this
102
103 push_sub(curv_briggs_finalize)
104
105 safe_deallocate_a(this%lsize)
106
107 pop_sub(curv_briggs_finalize)
108 end subroutine curv_briggs_finalize
109
110 ! ---------------------------------------------------------
111 subroutine curv_briggs_copy(this_out, this_in)
112 type(curv_briggs_t), intent(inout) :: this_out
113 type(curv_briggs_t), intent(in) :: this_in
114
115 push_sub(curv_briggs_copy)
116
117 safe_allocate_source_a(this_out%lsize, this_in%lsize)
118 this_out%beta = this_in%beta
120 pop_sub(curv_briggs_copy)
121 end subroutine curv_briggs_copy
122
123 ! ---------------------------------------------------------
124 pure function curv_briggs_to_cartesian(this, chi) result(xx)
125 class(curv_briggs_t), target, intent(in) :: this
126 real(real64), intent(in) :: chi(:)
127 real(real64) :: xx(1:this%dim)
128
129 ! no PUSH_SUB, called too often
130
131 xx = chi - this%lsize*this%beta/(m_two*m_pi)*sin(m_two*m_pi*chi/this%lsize)
132
133 end function curv_briggs_to_cartesian
134
135 ! ---------------------------------------------------------
136 function curv_briggs_from_cartesian(this, xx) result(chi)
137 class(curv_briggs_t), target, intent(in) :: this
138 real(real64), intent(in) :: xx(:)
139 real(real64) :: chi(1:this%dim)
141 ! no PUSH_SUB, called too often
143 chi = m_zero
144 message(1) = "Internal error in curv_briggs_from_cartesian"
147 end function curv_briggs_from_cartesian
148
149 ! ---------------------------------------------------------
150 subroutine curv_briggs_write_info(this, iunit, namespace)
151 class(curv_briggs_t), intent(in) :: this
152 integer, optional, intent(in) :: iunit
153 type(namespace_t), optional, intent(in) :: namespace
154
155 push_sub(curv_briggs_write_info)
157 write(message(1), '(a)') ' Curvilinear Method = briggs'
158 call messages_info(1, iunit=iunit, namespace=namespace)
159
161 end subroutine curv_briggs_write_info
162
163 ! ---------------------------------------------------------
164 real(real64) function curv_briggs_surface_element(this, idir) result(ds)
165 class(curv_briggs_t), intent(in) :: this
166 integer, intent(in) :: idir
167
168 ds = m_zero
169 message(1) = 'Surface element with briggs curvilinear coordinates not implemented'
170 call messages_fatal(1)
171
172 end function curv_briggs_surface_element
173
174 ! ---------------------------------------------------------
175 function curv_briggs_jacobian(this, chi) result(jacobian)
176 class(curv_briggs_t), intent(in) :: this
177 real(real64), intent(in) :: chi(:)
178 real(real64) :: jacobian(1:this%dim, 1:this%dim)
179
180 integer :: idim
181
182 ! Jacobian is diagonal in this method
183 jacobian(:, :) = m_zero
184 do idim = 1, this%dim
185 jacobian(idim, idim) = m_one - this%beta*cos(m_two*m_pi*chi(idim)/this%lsize(idim))
186 end do
187
188 end function curv_briggs_jacobian
189
190 ! ---------------------------------------------------------
191 function curv_briggs_jacobian_inverse(this, chi) result(jacobian_inverse)
192 class(curv_briggs_t), intent(in) :: this
193 real(real64), intent(in) :: chi(:)
194 real(real64) :: jacobian_inverse(1:this%dim, 1:this%dim)
195
196 integer :: idim
197
198 jacobian_inverse = this%jacobian(chi)
199 do idim = 1, this%dim
200 jacobian_inverse(idim, idim) = m_one/jacobian_inverse(idim, idim)
201 end do
202
205end module curv_briggs_oct_m
206
207!! Local Variables:
208!! mode: f90
209!! coding: utf-8
210!! End:
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
This module implements the curvilinear coordinates given in E.L. Briggs, D.J. Sullivan,...
real(real64) function, dimension(1:this%dim, 1:this%dim) curv_briggs_jacobian(this, chi)
class(curv_briggs_t) function, pointer curv_briggs_constructor(namespace, dim, lsize, spacing)
subroutine, public curv_briggs_copy(this_out, this_in)
real(real64) function curv_briggs_surface_element(this, idir)
real(real64) function, dimension(1:this%dim, 1:this%dim) curv_briggs_jacobian_inverse(this, chi)
subroutine curv_briggs_finalize(this)
pure real(real64) function, dimension(1:this%dim) curv_briggs_to_cartesian(this, chi)
subroutine curv_briggs_write_info(this, iunit, namespace)
real(real64) function, dimension(1:this%dim) curv_briggs_from_cartesian(this, xx)
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
abstract class to describe coordinate systems
int true(void)