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 :: det_jac => curv_briggs_det_jac
50 procedure :: write_info => curv_briggs_write_info
51 procedure :: surface_element => curv_briggs_surface_element
53 end type curv_briggs_t
54
55 interface curv_briggs_t
56 procedure curv_briggs_constructor
57 end interface curv_briggs_t
58
59contains
60
61 ! ---------------------------------------------------------
62 function curv_briggs_constructor(namespace, dim, lsize, spacing) result(briggs)
63 type(namespace_t), intent(in) :: namespace
64 integer, intent(in) :: dim
65 real(real64), intent(in) :: lsize(1:dim)
66 real(real64), intent(in) :: spacing(1:dim)
67 class(curv_briggs_t), pointer :: briggs
68
69 integer :: idim
70
72
73 safe_allocate(briggs)
74
75 briggs%local_basis = .true.
76 briggs%orthogonal = .true. ! This needs to be checked
77 briggs%dim = dim
78 safe_allocate(briggs%lsize(1:dim))
79 briggs%lsize(1:dim) = lsize(1:dim)
80
81 call parse_variable(namespace, 'CurvBriggsBeta', m_half, briggs%beta)
82
83 if (briggs%beta < m_zero .or. briggs%beta > m_one) then
84 message(1) = 'The parameter "CurvBriggsBeta" must lie between 0 and 1.'
85 call messages_fatal(1, namespace=namespace)
86 end if
87
88 briggs%min_mesh_scaling_product = m_one
89 do idim = 1, briggs%dim
90 ! corresponds to the distance of grid points at [+spacing/2,-spacing/2]
91 briggs%min_mesh_scaling_product = briggs%min_mesh_scaling_product * (m_one / &
92 (m_one - briggs%lsize(idim) * briggs%beta / (m_pi * spacing(idim)) * sin(m_pi * spacing(idim) / briggs%lsize(idim))))
93 end do
94
96 end function curv_briggs_constructor
97
98 ! ---------------------------------------------------------
99 subroutine curv_briggs_finalize(this)
100 type(curv_briggs_t), intent(inout) :: this
101
102 push_sub(curv_briggs_finalize)
103
104 safe_deallocate_a(this%lsize)
105
106 pop_sub(curv_briggs_finalize)
107 end subroutine curv_briggs_finalize
108
109 ! ---------------------------------------------------------
110 subroutine curv_briggs_copy(this_out, this_in)
111 type(curv_briggs_t), intent(inout) :: this_out
112 type(curv_briggs_t), intent(in) :: this_in
113
114 push_sub(curv_briggs_copy)
115
116 safe_allocate_source_a(this_out%lsize, this_in%lsize)
117 this_out%beta = this_in%beta
118
120 end subroutine curv_briggs_copy
121
122 ! ---------------------------------------------------------
123 pure function curv_briggs_to_cartesian(this, chi) result(xx)
124 class(curv_briggs_t), target, intent(in) :: this
125 real(real64), intent(in) :: chi(:)
126 real(real64) :: xx(1:this%dim)
127
128 ! no PUSH_SUB, called too often
129
130 xx = chi - this%lsize*this%beta/(m_two*m_pi)*sin(m_two*m_pi*chi/this%lsize)
131
132 end function curv_briggs_to_cartesian
133
134 ! ---------------------------------------------------------
135 function curv_briggs_from_cartesian(this, xx) result(chi)
136 class(curv_briggs_t), target, intent(in) :: this
137 real(real64), intent(in) :: xx(:)
138 real(real64) :: chi(1:this%dim)
139
140 ! no PUSH_SUB, called too often
142 chi = m_zero
143 message(1) = "Internal error in curv_briggs_from_cartesian"
146 end function curv_briggs_from_cartesian
147
148 ! ---------------------------------------------------------
149 pure real(real64) function curv_briggs_det_jac(this, xx, chi) result(jdet)
150 class(curv_briggs_t), intent(in) :: this
151 real(real64), intent(in) :: xx(:)
152 real(real64), intent(in) :: chi(:)
153
154 integer :: idim
155 real(real64) :: jac(1:this%dim)
156
157 ! no PUSH_SUB, called too often
158
159 ! Jacobian is diagonal in this method
160 do idim = 1, this%dim
161 jac(idim) = m_one - this%beta*cos(m_two*m_pi*chi(idim)/this%lsize(idim))
162 end do
163 jdet = product(jac)
164
165 end function curv_briggs_det_jac
166
167 ! ---------------------------------------------------------
168 subroutine curv_briggs_write_info(this, iunit, namespace)
169 class(curv_briggs_t), intent(in) :: this
170 integer, optional, intent(in) :: iunit
171 type(namespace_t), optional, intent(in) :: namespace
172
173 push_sub(curv_briggs_write_info)
174
175 write(message(1), '(a)') ' Curvilinear Method = briggs'
176 call messages_info(1, iunit=iunit, namespace=namespace)
177
179 end subroutine curv_briggs_write_info
180
181 ! ---------------------------------------------------------
182 real(real64) function curv_briggs_surface_element(this, idir) result(ds)
183 class(curv_briggs_t), intent(in) :: this
184 integer, intent(in) :: idir
185
186 ds = m_zero
187 message(1) = 'Surface element with briggs curvilinear coordinates not implemented'
188 call messages_fatal(1)
189
190 end function curv_briggs_surface_element
191
193
194!! Local Variables:
195!! mode: f90
196!! coding: utf-8
197!! 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,...
class(curv_briggs_t) function, pointer curv_briggs_constructor(namespace, dim, lsize, spacing)
pure real(real64) function curv_briggs_det_jac(this, xx, chi)
subroutine, public curv_briggs_copy(this_out, this_in)
real(real64) function curv_briggs_surface_element(this, idir)
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:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
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
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
abstract class to describe coordinate systems
int true(void)