44 real(real64),
allocatable :: lsize(:)
56 procedure curv_briggs_constructor
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
75 briggs%local_basis = .
true.
76 briggs%orthogonal = .
true.
78 safe_allocate(briggs%lsize(1:dim))
79 briggs%lsize(1:dim) = lsize(1:dim)
83 if (briggs%beta <
m_zero .or. briggs%beta >
m_one)
then
84 message(1) =
'The parameter "CurvBriggsBeta" must lie between 0 and 1.'
88 briggs%min_mesh_scaling_product =
m_one
89 do idim = 1, briggs%dim
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))))
100 type(curv_briggs_t),
intent(inout) :: this
104 safe_deallocate_a(this%lsize)
111 type(curv_briggs_t),
intent(inout) :: this_out
112 type(curv_briggs_t),
intent(in) :: this_in
116 safe_allocate_source_a(this_out%lsize, this_in%lsize)
117 this_out%beta = this_in%beta
124 class(curv_briggs_t),
target,
intent(in) :: this
125 real(real64),
intent(in) :: chi(:)
126 real(real64) :: xx(1:this%dim)
136 class(curv_briggs_t),
target,
intent(in) :: this
137 real(real64),
intent(in) :: xx(:)
138 real(real64) :: chi(1:this%dim)
143 message(1) =
"Internal error in curv_briggs_from_cartesian"
151 real(real64),
intent(in) :: xx(:)
152 real(real64),
intent(in) :: chi(:)
155 real(real64) :: jac(1:this%dim)
160 do idim = 1, this%dim
170 integer,
optional,
intent(in) :: iunit
171 type(
namespace_t),
optional,
intent(in) :: namespace
175 write(
message(1),
'(a)')
' Curvilinear Method = briggs'
184 integer,
intent(in) :: idir
187 message(1) =
'Surface element with briggs curvilinear coordinates not implemented'
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
real(real64), parameter, public m_zero
real(real64), parameter, public m_pi
some mathematical constants
real(real64), parameter, public m_half
real(real64), parameter, public m_one
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
abstract class to describe coordinate systems