linsolve Subroutine

public subroutine linsolve(a, b, x)

Solve linear equations A * X = B

Arguments

Type IntentOptional AttributesName
real(kind=double), intent(in), dimension(:,:):: a
real(kind=double), intent(in), dimension(:):: b
real(kind=double), intent(out), dimension(:):: x

Calls

proc~~linsolve~~CallsGraph proc~linsolve linsolve proc~matinv matinv proc~linsolve->proc~matinv dgetrs dgetrs proc~linsolve->dgetrs dgetrf dgetrf proc~linsolve->dgetrf proc~ludcmp ludcmp proc~matinv->proc~ludcmp proc~outerprod outerprod proc~ludcmp->proc~outerprod

Called by

proc~~linsolve~~CalledByGraph proc~linsolve linsolve proc~colrad colrad proc~colrad->proc~linsolve proc~fida_weights_los fida_weights_los proc~fida_weights_los->proc~colrad proc~read_equilibrium read_equilibrium proc~read_equilibrium->proc~colrad proc~attenuate attenuate proc~attenuate->proc~colrad proc~pnpa_f pnpa_f proc~pnpa_f->proc~attenuate proc~pnpa_mc pnpa_mc proc~pnpa_mc->proc~attenuate proc~npa_f npa_f proc~npa_f->proc~attenuate program~fidasim fidasim program~fidasim->proc~fida_weights_los program~fidasim->proc~read_equilibrium program~fidasim->proc~pnpa_f program~fidasim->proc~pnpa_mc program~fidasim->proc~npa_f proc~npa_mc npa_mc program~fidasim->proc~npa_mc proc~npa_weights npa_weights program~fidasim->proc~npa_weights proc~npa_mc->proc~attenuate proc~npa_weights->proc~attenuate

Contents

Source Code


Source Code

  subroutine linsolve(a,b,x)
   !+ Solve linear equations A * X = B
   real(double), dimension(:,:), intent(IN)  :: a ! assume a square
   real(double), dimension(:), intent(IN)    :: b ! assume size(b) == size(a,1)
   real(double), dimension(:), intent(OUT)   :: x ! assume size(x) == size(b)

#ifdef _USE_BLAS
   real(double), dimension(size(b),size(b)) :: lu
   integer, dimension(size(b)) :: ipiv
   integer :: n,info

   n = size(b)

   ! first factorize a
   lu(:,:) = a(:,:)
   call DGETRF(n,n,lu,n,ipiv,info)
   if (info /= 0) stop 'sub linsolve: DGETRF failed!'

   x(:) = b(:)
   call DGETRS('N',n,1,lu,n,ipiv,x,n,info)
   if (info /= 0) stop 'sub linsolve: DGETRS failed!'

#else
   real(double), dimension(size(b),size(b)) :: a_inv
   call matinv(a, a_inv)
   x = matmul(a_inv, b)!coeffs determined from states at t=0
#endif

  end subroutine linsolve