tem_createPlane Subroutine

public subroutine tem_createPlane(me, origin, vecA, vecB)

This routine creates plane definition from given origin and two vectors \verbatim vecB______ /\ | | | | | | | | | |-------------->| origin vecA \endverbatim

Arguments

Type IntentOptional Attributes Name
type(tem_plane_type), intent(out) :: me
real(kind=rk), intent(in) :: origin(3)
real(kind=rk), intent(in) :: vecA(3)
real(kind=rk), intent(in) :: vecB(3)

Calls

proc~~tem_createplane~~CallsGraph proc~tem_createplane tem_createPlane proc~cross_product3d cross_product3D proc~tem_createplane->proc~cross_product3d

Called by

proc~~tem_createplane~~CalledByGraph proc~tem_createplane tem_createPlane proc~tem_createbox tem_createBox proc~tem_createbox->proc~tem_createplane proc~tem_createbox~2 tem_createBox proc~tem_createbox~2->proc~tem_createplane proc~tem_load_onecanonicalnd tem_load_oneCanonicalND proc~tem_load_onecanonicalnd->proc~tem_createplane proc~tem_load_onecanonicalnd->proc~tem_createbox proc~tem_load_canonicalnd_vec tem_load_canonicalND_vec proc~tem_load_canonicalnd_vec->proc~tem_load_onecanonicalnd interface~tem_load_canonicalnd tem_load_canonicalND interface~tem_load_canonicalnd->proc~tem_load_onecanonicalnd interface~tem_load_canonicalnd->proc~tem_load_canonicalnd_vec proc~tem_load_shape_single tem_load_shape_single proc~tem_load_shape_single->interface~tem_load_canonicalnd proc~tem_load_shapes tem_load_shapes proc~tem_load_shapes->proc~tem_load_shape_single interface~tem_load_shape tem_load_shape interface~tem_load_shape->proc~tem_load_shape_single

Contents

Source Code


Source Code

  subroutine tem_createPlane(me, origin, vecA, vecB)
    !--------------------------------------------------------------------------!
    type(tem_plane_type), intent(out) :: me !< plane to return
    real(kind=rk), intent(in) :: origin(3) !< plane origin
    real(kind=rk), intent(in) :: vecA(3) !< 1st vector
    real(kind=rk), intent(in) :: vecB(3) !< 2nd vector
    !--------------------------------------------------------------------------!
    integer :: iTri
    !--------------------------------------------------------------------------!
    me%origin = origin
    me%vec(:,1) = vecA
    me%vec(:,2) = vecB
    ! convert plane into triangle
    do iTri = 1,2
      me%triangle(iTri)%nodes(1:3,1) = me%origin
      me%triangle(iTri)%nodes(1:3,2) = me%origin + me%vec(:,iTri)
      me%triangle(iTri)%nodes(1:3,3) = me%origin + me%vec(:,1) + me%vec(:,2)
    end do

    ! normal direction = crossproduct(vecA, vecB)
    me%unitNormal = cross_product3D(vecA, vecB)
    me%unitNormal = me%unitNormal                       &
      & / sqrt(dot_product(me%unitNormal, me%unitNormal))

  end subroutine tem_createPlane