44    real(real64), 
allocatable, 
public :: half_length(:)
 
   46    integer, 
public :: n_periodic_boundaries = 0
 
   58    procedure box_parallelepiped_constructor
 
   65    integer,            
intent(in) :: dim
 
   66    real(real64),       
intent(in) :: center(dim)
 
   67    real(real64),       
intent(in) :: axes(dim, dim)
 
   68    real(real64),       
intent(in) :: length(dim)
 
   69    type(namespace_t),  
intent(in) :: namespace
 
   70    integer, 
optional,  
intent(in) :: n_periodic_boundaries
 
   71    class(box_parallelepiped_t), 
pointer :: box
 
   79    safe_allocate(box%half_length(1:dim))
 
   80    box%half_length = 
m_half * length
 
   82    if (
present(n_periodic_boundaries)) 
then 
   83      box%n_periodic_boundaries = n_periodic_boundaries
 
   85    call box_shape_init(box, namespace, dim, center, bounding_box_min=-box%half_length, bounding_box_max=box%half_length, axes=axes)
 
   87    box%bounding_box_l = box%half_length + abs(center)
 
   94    class(box_parallelepiped_t), 
intent(inout) :: box
 
   95    integer,                     
intent(in) :: dim
 
   96    real(real64),                
intent(in) :: axes(dim, dim)
 
   97    real(real64),                
intent(in) :: length(dim)
 
   98    type(namespace_t),           
intent(in) :: namespace
 
  100    real(real64) :: center(dim)
 
  105    safe_deallocate_a(box%center)
 
  106    safe_deallocate_a(box%bounding_box_l)
 
  109    box%half_length = 
m_half * length
 
  110    call box_shape_init(box, namespace, dim, center, bounding_box_min=-box%half_length, bounding_box_max=box%half_length, axes=axes)
 
  112    box%bounding_box_l = box%half_length + abs(center)
 
  120    type(box_parallelepiped_t), 
intent(inout) :: this
 
  125    safe_deallocate_a(this%half_length)
 
  132    class(box_parallelepiped_t), 
intent(in)  :: this
 
  133    integer,                     
intent(in)  :: nn
 
  134    real(real64), 
contiguous,           
intent(in)  :: xx(:,:)
 
  135    logical :: contained(1:nn)
 
  138    real(real64) :: ulimit(this%dim), llimit(this%dim)
 
  141    do idir = 1, this%dim
 
  142      if (idir <= this%n_periodic_boundaries) 
then 
  151      contained(ip) = .
true.
 
  152      do idir = 1, this%dim
 
  153        contained(ip) = contained(ip) .and. xx(ip, idir) >= llimit(idir) .and. xx(ip, idir) <= ulimit(idir)
 
  154        if (idir > this%n_periodic_boundaries) 
then 
  156          contained(ip) = contained(ip) .neqv. this%is_inside_out()
 
  170    real(real64),       
intent(in)  :: mesh_spacing(:)
 
  171    integer,            
intent(in)  :: nn
 
  172    real(real64),       
intent(in)  :: xx(:,:)
 
  173    integer, 
optional,  
intent(in)  :: number_of_layers
 
  175    logical :: surface_points(1:nn)
 
  176    integer :: idir, number_of_layers_
 
  177    real(real64)   :: shrink(this%dim), test_axis(this%dim)
 
  181    do idir = 1, this%dim
 
  183      test_axis(idir) = 
m_one 
  184      assert(all(abs(this%axes%vectors(:, idir) - test_axis(:)) < 
m_epsilon))
 
  187    number_of_layers_ = 1
 
  188    if (
present(number_of_layers)) number_of_layers_ = number_of_layers
 
  191    shrink(:) = 1 - number_of_layers_ * (1 - 
box_boundary_delta) * (mesh_spacing(:) / this%half_length(:))
 
  195      m_two*this%half_length(:)*shrink(:), namespace)
 
  200    if (
SIZE(xx, 1) == 3) 
then 
  206    safe_deallocate_p(shrinked_box)
 
  212    real(real64),                
intent(in)  :: point_coordinates(:)
 
  213    real(real64),                
intent(in)  :: mesh_spacing(:)
 
  214    real(real64),                
intent(out) :: normal_vector(:)
 
  215    real(real64),                
intent(out) :: surface_element
 
  217    real(real64)   :: vector_norm
 
  218    integer :: idir, first_index, second_index, third_index, face_counter
 
  230    normal_vector(:) = (point_coordinates(:) - this%center(:)) / this%half_length(:)
 
  235    do idir = 1, this%dim
 
  237      first_index = this%dim - mod(idir, this%dim)
 
  238      second_index = this%dim - mod(idir + 1, this%dim)
 
  239      third_index = this%dim - mod(idir + 2, this%dim)
 
  240      if (abs(normal_vector(first_index)) >= abs(normal_vector(second_index)) .and. &
 
  241        abs(normal_vector(first_index)) >= abs(normal_vector(third_index))) 
then 
  243        surface_element = mesh_spacing(second_index) * mesh_spacing(third_index)
 
  244        face_counter = face_counter + 1
 
  246        if (abs(normal_vector(second_index)) / abs(normal_vector(first_index)) < &
 
  247          m_one - 
m_half * mesh_spacing(second_index) / this%half_length(first_index)) 
then 
  248          normal_vector(second_index) = 
m_zero 
  250        if (abs(normal_vector(third_index)) / abs(normal_vector(first_index)) < &
 
  251          m_one - 
m_half * mesh_spacing(third_index) / this%half_length(first_index)) 
then 
  252          normal_vector(third_index) = 
m_zero 
  257    if (face_counter == 3) surface_element = surface_element * (
m_three / 
m_four)
 
  261    vector_norm = norm2(normal_vector(:))
 
  264    normal_vector(:) = normal_vector(:) / vector_norm
 
  272    integer,           
optional, 
intent(in) :: iunit
 
  273    type(
namespace_t), 
optional, 
intent(in) :: namespace
 
  279    write(
message(1),
'(2x,a)') 
'Type = parallelepiped' 
  291    type(
unit_t),                
intent(in) :: unit_length
 
  297    write(
info, 
'(a,a,a,99(f11.6,a))') 
'BoxShape = parallelepiped; Lengths [', trim(
units_abbrev(unit_length)),
'] = [', &
 
real(real64), parameter, public box_boundary_delta
 
class(box_parallelepiped_t) function, pointer box_parallelepiped_constructor(dim, center, axes, length, namespace, n_periodic_boundaries)
 
subroutine box_parallelepiped_shape_get_surface_point_info(this, point_coordinates, mesh_spacing, normal_vector, surface_element)
 
subroutine box_parallelepiped_regenerate(box, dim, axes, length, namespace)
 
subroutine box_parallelepiped_finalize(this)
 
subroutine box_parallelepiped_write_info(this, iunit, namespace)
 
logical function, dimension(1:nn) box_parallelepiped_shape_get_surface_points(this, namespace, mesh_spacing, nn, xx, number_of_layers)
Get a mask for the grid points telling which of them are surface points.
 
character(len=box_info_len) function box_parallelepiped_short_info(this, unit_length)
 
logical function, dimension(1:nn) box_parallelepiped_shape_contains_points(this, nn, xx)
 
subroutine, public box_shape_init(this, namespace, dim, center, bounding_box_min, bounding_box_max, axes)
 
subroutine, public box_shape_end(this)
 
real(real64), parameter, public m_two
 
real(real64), parameter, public m_zero
 
real(real64), parameter, public m_four
 
real(real64), parameter, public m_epsilon
 
real(real64), parameter, public m_half
 
real(real64), parameter, public m_one
 
real(real64), parameter, public m_three
 
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
 
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
 
character(len=20) pure function, public units_abbrev(this)
 
This module defines the unit system, used for input and output.
 
type(unit_system_t), public units_out
 
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
 
Base class for more specialized boxes that are defined by a shape and have a center and basis vectors...