moindres_carres_lineaire Subroutine

public subroutine moindres_carres_lineaire(nb_var, nb_pts, hij, beta, Jf)

Note

Function that returns the parameters of a function that approximates a data set. The parameters determination is achieved by linear least squares approximation.

Arguments

Type IntentOptional Attributes Name
integer(kind=I4), intent(in) :: nb_var

number of parameters to be determined

integer(kind=I4), intent(in) :: nb_pts

number of points for function evaluation

real(kind=R8), intent(in), dimension(1:nb_pts) :: hij

vector of evaluation points

real(kind=R8), intent(out), dimension(1:nb_var) :: beta

parameters vector

real(kind=R8), intent(in), dimension(1:nb_pts, 1:nb_var) :: Jf

Jacobian


Calls

proc~~moindres_carres_lineaire~~CallsGraph proc~moindres_carres_lineaire moindres_carres_lineaire proc~choldc choldc proc~moindres_carres_lineaire->proc~choldc proc~cholsl cholsl proc~moindres_carres_lineaire->proc~cholsl

Called by

proc~~moindres_carres_lineaire~~CalledByGraph proc~moindres_carres_lineaire moindres_carres_lineaire proc~least_squares_tcheby least_squares_tcheby proc~least_squares_tcheby->proc~moindres_carres_lineaire program~test_least test_least program~test_least->proc~moindres_carres_lineaire program~test_tchebychev test_tchebychev program~test_tchebychev->proc~least_squares_tcheby

Source Code

   subroutine moindres_carres_lineaire(nb_var, nb_pts, hij, beta, Jf)
   !================================================================================================
   !< @note Function that returns the parameters of a function that approximates a data set. The parameters
   !<        determination is achieved by linear least squares approximation.
   !<
   !<  @endnote
   !------------------------------------------------------------------------------------------------
   implicit none
   integer(kind=I4), intent(in )                                  :: nb_var   !! *number of parameters to be determined*
   integer(kind=I4), intent(in )                                  :: nb_pts   !! *number of points for function evaluation*
   real   (kind=R8), intent(in ), dimension(1:nb_pts)             :: hij      !! *vector of evaluation points*
   real   (kind=R8), intent(in ), dimension(1:nb_pts, 1:nb_var)   :: Jf       !! *Jacobian*
   real   (kind=R8), intent(out), dimension(1:nb_var)             :: beta     !! *parameters vector*

      integer(kind=I4) :: info

      real(kind=R8),              dimension(1:nb_var)           :: scd, pt
      real(kind=R8),              dimension(1:nb_var, 1:nb_var) :: mat
      real(kind=R8), allocatable, dimension(:,:)                :: JfT

      allocate( JfT(1:nb_var, 1:nb_pts) )

      JfT(1:nb_var, 1:nb_pts) = transpose(  Jf(1:nb_pts, 1:nb_var) )
      mat(1:nb_var, 1:nb_var)  =   matmul( JfT(1:nb_var, 1:nb_pts),  Jf(1:nb_pts, 1:nb_var) )
      scd(1:nb_var)            =   matmul( JfT(1:nb_var, 1:nb_pts), hij(1:nb_pts) )

      call choldc(a = mat(1:nb_var, 1:nb_var),  & !
                  n = nb_var,                   & !
                 np = nb_var,                   & !
                  p = pt(1:nb_var),             & !
               info = info)                       !

      call cholsl(a = mat(1:nb_var, 1:nb_var),  & !
                  n = nb_var,                   & !
                 np = nb_var,                   & !
                  p =  pt(1:nb_var),            & !
                  b = scd(1:nb_var),            & !
                  x = scd(1:nb_var),            & !
                  info = info)                    !

      beta(1:nb_var) = scd(1:nb_var)

      deallocate( JfT )

   return
   endsubroutine moindres_carres_lineaire