append_intpMatrixLSF Subroutine

private subroutine append_intpMatrixLSF(me, order, QQ, nDims, nSources, cxDirRK, neighDir, pos, success)

This routine builds up the matrix for least square fit used in linear and quadratic interpolation.

Compute interpolation matrix for least square fit using stencil direction of available sources The parent of target childs coord is 0,0,0 so we could just use of stencil%cxDir to build up this matrix entries Every row in matrix is evaluated with coord of source element

Arguments

Type IntentOptional Attributes Name
type(tem_intpMatrixLSF_type), intent(inout) :: me

intpMatrix for LSF fill

integer, intent(inout) :: order

interpolation order calculated for current element depending on nSources if quadratic LSF matrix is singular fall back to linear

integer, intent(in) :: QQ

Number of stencil directions

integer, intent(in) :: nDims

Number of dimensions

integer, intent(in) :: nSources

Number of sources from coarser found

real(kind=rk), intent(in) :: cxDirRK(3,QQ)

Stencil directions

integer, intent(in) :: neighDir(nSources)

direction in which sources are found

integer, intent(out) :: pos

Pointer to position of interpolation matrix in growing array of matrix

logical, intent(out) :: success

success if false if matrix is singular reduce interpolation order


Calls

proc~~append_intpmatrixlsf~~CallsGraph proc~append_intpmatrixlsf append_intpMatrixLSF interface~append~24 append proc~append_intpmatrixlsf->interface~append~24 proc~build_matrixlsf_linearintp build_matrixLSF_linearIntp proc~append_intpmatrixlsf->proc~build_matrixlsf_linearintp proc~build_matrixlsf_quadintp build_matrixLSF_quadIntp proc~append_intpmatrixlsf->proc~build_matrixlsf_quadintp proc~append_arrayga2d_real append_arrayga2d_real interface~append~24->proc~append_arrayga2d_real proc~append_singlega2d_real append_singlega2d_real interface~append~24->proc~append_singlega2d_real proc~alloc_matrix alloc_matrix proc~build_matrixlsf_linearintp->proc~alloc_matrix proc~invert_matrix invert_matrix proc~build_matrixlsf_linearintp->proc~invert_matrix proc~polylinear_1d polyLinear_1D proc~build_matrixlsf_linearintp->proc~polylinear_1d proc~polylinear_2d polyLinear_2D proc~build_matrixlsf_linearintp->proc~polylinear_2d proc~polylinear_3d polyLinear_3D proc~build_matrixlsf_linearintp->proc~polylinear_3d proc~tem_abort tem_abort proc~build_matrixlsf_linearintp->proc~tem_abort proc~tem_matrix_dump tem_matrix_dump proc~build_matrixlsf_linearintp->proc~tem_matrix_dump proc~build_matrixlsf_quadintp->proc~alloc_matrix proc~build_matrixlsf_quadintp->proc~invert_matrix proc~polyquadratic_1d polyQuadratic_1D proc~build_matrixlsf_quadintp->proc~polyquadratic_1d proc~polyquadratic_2d polyQuadratic_2D proc~build_matrixlsf_quadintp->proc~polyquadratic_2d proc~polyquadratic_3d polyQuadratic_3D proc~build_matrixlsf_quadintp->proc~polyquadratic_3d proc~build_matrixlsf_quadintp->proc~tem_abort proc~build_matrixlsf_quadintp->proc~tem_matrix_dump proc~alloc_matrix->proc~tem_abort interface~expand~22 expand proc~append_arrayga2d_real->interface~expand~22 proc~append_singlega2d_real->interface~expand~22 proc~invert_matrix->proc~tem_abort dgetrf dgetrf proc~invert_matrix->dgetrf dgetri dgetri proc~invert_matrix->dgetri mpi_abort mpi_abort proc~tem_abort->mpi_abort proc~expand_ga2d_real expand_ga2d_real interface~expand~22->proc~expand_ga2d_real

Called by

proc~~append_intpmatrixlsf~~CalledByGraph proc~append_intpmatrixlsf append_intpMatrixLSF interface~append~9 append interface~append~9->proc~append_intpmatrixlsf

Source Code

  subroutine append_intpMatrixLSF(me, order, QQ, nDims, nSources, cxDirRK, &
    &                             neighDir, pos, success)
    ! --------------------------------------------------------------------------
    !> intpMatrix for LSF fill
    type(tem_intpMatrixLSF_type), intent(inout) :: me
    !> interpolation order calculated for current element depending on nSources
    !! if quadratic LSF matrix is singular fall back to linear
    integer, intent(inout) :: order
    !> Number of stencil directions
    integer, intent(in) :: QQ
    !> Number of dimensions
    integer, intent(in) :: nDims
    !> Number of sources from coarser found
    integer, intent(in) :: nSources
    !> Stencil directions
    real(kind=rk), intent(in) :: cxDirRK(3,QQ)
    !> direction in which sources are found
    integer, intent(in) :: neighDir(nSources)
    !> Pointer to position of interpolation matrix in growing array of matrix
    integer, intent(out) :: pos
    !> success if false if matrix is singular reduce interpolation order
    logical, intent(out) :: success
    ! --------------------------------------------------------------------------
    integer :: hashID
    type(tem_matrix_type) :: matLSF_tmp
    logical :: wasAdded
    integer :: iNeigh, iSrc
    ! --------------------------------------------------------------------------
write(dbgUnit(4),"(A,I0)") 'Inside append least square fit matrix for intp ' &
  &                     //'order: ', order
    ! hashID to identify unique combination of available nSources
    hashID = 0
    do iSrc = 1, nSources
      iNeigh = neighDir(iSrc)
      hashID = ibset(hashID, iNeigh)
    end do
write(dbgUnit(4),"(A,I0)") '   hashID: ', hashID
    call append( me       = me%ID,    &
      &          val      = hashID,   &
      &          pos      = pos,      &
      &          wasAdded = wasAdded  )
write(dbgUnit(4),"(A,L5)") '   wasAdded: ', wasAdded
write(dbgUnit(4),"(A,I5)") '        pos: ', pos

    if (wasAdded) then

      select case(order)
      case(1) ! linear
        call build_matrixLSF_linearIntp( me       = matLSF_tmp, &
          &                              QQ       = QQ,         &
          &                              nDims    = nDims,      &
          &                              nSources = nSources,   &
          &                              cxDirRK  = cxDirRK,    &
          &                              neighDir = neighDir,   &
          &                              nCoeffs  = me%nCoeffs, &
          &                              success  = success     )
      case(2) ! quadratic
        call build_matrixLSF_quadIntp( me       = matLSF_tmp, &
          &                            QQ       = QQ,         &
          &                            nDims    = nDims,      &
          &                            nSources = nSources,   &
          &                            cxDirRK  = cxDirRK,    &
          &                            neighDir = neighDir,   &
          &                            nCoeffs  = me%nCoeffs, &
          &                            success  = success     )
      case default
        write(logUnit(1),*) 'Unsupported interpolation order to build matrix ' &
          &                 //'for LSF: ', order
      end select

      call append( me = me%isInvertible, val = success )
      call append( me = me%matArray, val = matLSF_tmp )
    else
      ! hashID already exist then return success depends on matrix isInvertible
      if (me%isInvertible%val(pos)) then
        success = .true.
      else
        success = .false.
      end if
    end if
write(dbgUnit(4),"(A,L5)") '     success  : ', success

  end subroutine append_intpMatrixLSF