Initialization of the mean sprectrum with N Gaussian
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n_gauss | number of Gaussian |
||
real(kind=xp), | intent(inout), | dimension(3*n_gauss) | :: | params | params to optimize |
|
integer, | intent(in) | :: | dim_v | dimension along v axis |
||
real(kind=xp), | intent(in), | dimension(dim_v) | :: | line | spectrum |
|
real(kind=xp), | intent(in) | :: | amp_fact_init | times max amplitude of additional Gaussian |
||
real(kind=xp), | intent(in) | :: | sig_init | dispersion of additional Gaussian |
||
integer, | intent(in) | :: | maxiter | Max number of iteration |
||
integer, | intent(in) | :: | m | number of corrections used in the limited memory matrix by LBFGS-B |
||
integer, | intent(in) | :: | iprint | print option |
subroutine init_spectrum(n_gauss, params, dim_v, line, amp_fact_init, sig_init, maxiter, m, iprint)
!! Initialization of the mean sprectrum with N Gaussian
implicit none
integer, intent(in) :: n_gauss !! number of Gaussian
integer, intent(in) :: dim_v !! dimension along v axis
integer, intent(in) :: maxiter !! Max number of iteration
integer, intent(in) :: m !! number of corrections used in the limited memory matrix by LBFGS-B
integer, intent(in) :: iprint !! print option
real(xp), intent(in), dimension(dim_v) :: line !! spectrum
real(xp), intent(in) :: amp_fact_init !! times max amplitude of additional Gaussian
real(xp), intent(in) :: sig_init !! dispersion of additional Gaussian
real(xp), intent(inout), dimension(3*n_gauss) :: params !! params to optimize
integer :: i, j, k, p
real(xp), dimension(:), allocatable :: lb, ub
real(xp), dimension(dim_v) :: model, residual
real(xp), dimension(:), allocatable :: x
do i=1, n_gauss
allocate(lb(3*i), ub(3*i))
model = 0._xp
residual = 0._xp
lb = 0._xp; ub=0._xp
call init_bounds(line, i, dim_v, lb, ub)
do j=1, i
do k=1, dim_v
model(k) = model(k) + gaussian(k, params(1+(3*(j-1))), params(2+(3*(j-1))), params(3+(3*(j-1))))
end do
enddo
residual = model - line
allocate(x(3*i))
x = 0._xp
do p=1, 3*(i-1)
x(p) = params(p);
end do
x(2+(3*(i-1))) = minloc(residual, dim_v)
x(1+(3*(i-1))) = line(int(x(2+(3*(i-1))))) * amp_fact_init
x(3+(3*(i-1))) = sig_init;
call minimize_spec(3*i, m, x, lb, ub, line, dim_v, i, maxiter, iprint)
do p=1, 3*i
params(p) = x(p);
end do
deallocate(x)
deallocate(lb, ub)
enddo
end subroutine init_spectrum