Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fixing misc bugs #88

Merged
merged 7 commits into from
Nov 16, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
6 changes: 4 additions & 2 deletions sources/datalloc.f90
Original file line number Diff line number Diff line change
Expand Up @@ -195,7 +195,7 @@ subroutine AllocData(type)
else if (coil(icoil)%type==coil_type_spline) then
Tdof = Tdof + 1 + 3*(Splines(icoil)%NCP)
else
Tdof = Tdof + coil(icoil)%Ic + DoF(icoil)%ND
Tdof = Tdof + 1 + DoF(icoil)%ND
end if
if (DoF(icoil)%ND >= Cdof) Cdof = DoF(icoil)%ND ! find the largest ND for single coil;

Expand All @@ -209,7 +209,9 @@ subroutine AllocData(type)
SALLOCATE( xdof, (1:Ndof), zero ) ! dof vector;
SALLOCATE( dofnorm, (1:Ndof), one ) ! dof normalized value vector;
SALLOCATE( evolution, (1:Nouts+1, 0:13), zero ) !evolution array;
SALLOCATE( coilspace, (1:Nouts+1, 1:Tdof), zero ) ! all the coil parameters;
if (Tdof >= 1) then
SALLOCATE( coilspace, (1:Nouts+1, 1:Tdof), zero ) ! all the coil parameters;
endif

! determine dofnorm
if ( IsNormalize > 0 ) then
Expand Down
7 changes: 4 additions & 3 deletions sources/diagnos.f90
Original file line number Diff line number Diff line change
Expand Up @@ -63,14 +63,14 @@ SUBROUTINE diagnos
coilspace(iout, idof+1:idof+NF ) = FouCoil(icoil)%ys(1:NF) ; idof = idof + NF
coilspace(iout, idof+1:idof+NF+1) = FouCoil(icoil)%zc(0:NF) ; idof = idof + NF +1
coilspace(iout, idof+1:idof+NF ) = FouCoil(icoil)%zs(1:NF) ; idof = idof + NF
case (coil_type_spline)
case (coil_type_spline)
NCP = Splines(icoil)%NCP
coilspace(iout, idof+1:idof+NCP*3) = Splines(icoil)%Cpoints(0:3*NCP-1) ; idof = idof + 3*NCP
!!$ case default
!!$ FATAL(descent, .true., not supported coil types)
end select
enddo
!!$ FATAL( output , idof .ne. Tdof, counting error in restart )
FATAL( output , idof .ne. Tdof, counting error in restart )
endif
!-------------------------------average coil length-------------------------------------------------------
AvgLength = zero
Expand Down Expand Up @@ -317,7 +317,8 @@ SUBROUTINE diagnos

!--------------------------------calculate the stochastic Bn error----------------------------
if ( Npert .ge. 1 .and. allocated(surf(isurf)%bn) ) then


if (case_optimize .eq. 0) call perturbation(0)
call sbnormal( 0 )

if(myid .eq. 0) write(ounit, '(8X": Maximum field error after perturbations: "ES23.15)') bnormmax
Expand Down
4 changes: 2 additions & 2 deletions sources/globals.f90
Original file line number Diff line number Diff line change
Expand Up @@ -16,7 +16,7 @@ module globals
!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!


CHARACTER(10), parameter :: version='v0.17.00' ! version number
CHARACTER(10), parameter :: version='v0.17.01' ! version number


!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!
Expand Down Expand Up @@ -183,7 +183,7 @@ module globals
REAL :: weight_specw = 0.000D+00
! bnrom stochastic
REAL :: weight_sbnorm = 0.000D+00 ! Stochastic bnormal weight
INTEGER :: Npert = 10 ! Number of coil perturbations
INTEGER :: Npert = 0 ! Number of coil perturbations
INTEGER :: Nmax = 3 ! Max frequency of perturbations
REAL :: sdelta = 1.000D-02 ! Perturbation magnitude
! optimize controls
Expand Down
7 changes: 6 additions & 1 deletion sources/nissin.f90
Original file line number Diff line number Diff line change
Expand Up @@ -71,7 +71,12 @@ subroutine nisscom(ideriv)
endif

if ( coil(icoil)%Lc /= 0 ) then !if geometry is free;
call nissinDeriv1( icoil, t1N(idof+1:idof+ND), ND, NF )
if (abs(nissin) .le. machprec) then
t1N(idof+1:idof+ND) = 0 ! avoid divided-by-zero-error
else
call nissinDeriv1( icoil, t1N(idof+1:idof+ND), ND, NF )
endif

if (mnissin > 0) then ! L-M format of targets
LM_fjac(inissin+ivec, idof+1:idof+ND) = weight_nissin * t1N(idof+1:idof+ND)
ivec = ivec + 1
Expand Down
84 changes: 44 additions & 40 deletions sources/saving.f90
Original file line number Diff line number Diff line change
Expand Up @@ -174,43 +174,45 @@ subroutine saving
HWRITERA( Nteta,Nzeta , Bz , surf(plasma)%Bz(0:Nteta-1,0:Nzeta-1) )
endif

NSmax = 0
do icoil = 1, Ncoils
if(coil(icoil)%NS .gt. NSmax) NSmax = coil(icoil)%NS
enddo
SALLOCATE(tempvar, (1:Ncoils, 0:NSmax) , 0.0)

! Save coil data
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%xx(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, xx , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%yy(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, yy , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%zz(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, zz , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%xt(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, xt , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%yt(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, yt , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%zt(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, zt , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
if (weight_sbnorm > 0) then
NSmax = 0
do icoil = 1, Ncoils
if(coil(icoil)%NS .gt. NSmax) NSmax = coil(icoil)%NS
enddo
SALLOCATE(tempvar, (1:Ncoils, 0:NSmax) , 0.0)

! Save coil data
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%xx(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, xx , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%yy(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, yy , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%zz(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, zz , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%xt(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, xt , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%yt(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, yt , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
do icoil = 1, Ncoils
tempvar(icoil,0:coil(icoil)%NS) = coil(icoil)%zt(0:coil(icoil)%NS)
enddo
HWRITERA( Ncoils, NSmax+1, zt , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0
endif

! Save filamentary body force loads
if (filforce .eq. 1) then
Expand Down Expand Up @@ -278,17 +280,19 @@ subroutine saving
enddo
HWRITERA( Ncoils, NSmax+1, bfbz , tempvar(1:Ncoils,0:NSmax) )
tempvar(1:Ncoils,0:NSmax) = 0.0

DALLOCATE(tempvar)
endif

DALLOCATE(tempvar)

HWRITEIV( 1 , iout , iout )
HWRITERV( 1 , Inorm , Inorm )
HWRITERV( 1 , Gnorm , Gnorm )
HWRITERV( 1 , Mnorm , Mnorm )
HWRITERV( 1 , overlap , overlap )
HWRITERA( iout, 14 , evolution , evolution(1:iout, 0:13) )
HWRITERA( iout, Tdof , coilspace , coilspace(1:iout, 1:Tdof) )
if (allocated(coilspace)) then
HWRITERA( iout, Tdof , coilspace , coilspace(1:iout, 1:Tdof) )
endif

if (allocated(deriv)) then
HWRITERA( Ndof, 12 , deriv , deriv(1:Ndof, 0:11) )
Expand Down
5 changes: 2 additions & 3 deletions sources/solvers.f90
Original file line number Diff line number Diff line change
Expand Up @@ -495,7 +495,6 @@ subroutine costfun(ideriv)
!if(myid==0) write(ounit, '("-------i=9, chi = "ES23.15)') chi
! nissin
if (weight_nissin > 0.0_dp) then

call nisscom(ideriv)
chi = chi + weight_nissin * nissin
if ( ideriv == 1 ) then
Expand Down Expand Up @@ -725,7 +724,7 @@ subroutine normweight

!-!-!-!-!-!-!-!-!-!-sbnorm-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!-!

if( weight_sbnorm >= 0.0_dp ) then
if( weight_sbnorm > 0.0_dp ) then

call sbnormal(0)
if (abs(bnormavg) > machprec) weight_sbnorm = weight_sbnorm / bnormavg
Expand Down Expand Up @@ -853,7 +852,7 @@ subroutine output (mark)
coilspace(iout, idof+1:idof+3*NCP) = Splines(icoil)%Cpoints(0:3*NCP-1) ; idof = idof + 3*NCP
end select
enddo
!!$ FATAL( output , idof .ne. Tdof, counting error in restart )
FATAL( output , idof .ne. Tdof, counting error in restart )
endif

if(mod(iout,save_freq) .eq. 0) call saving
Expand Down
2 changes: 1 addition & 1 deletion sources/specinp.f90
Original file line number Diff line number Diff line change
Expand Up @@ -109,7 +109,7 @@ SUBROUTINE specinp
do imn = 1, surf(isurf)%Nfou
write(wunit,1010) surf(isurf)%bin(imn)/Nfp_raw, surf(isurf)%bim(imn), surf(isurf)%Rbc(imn), &
surf(isurf)%bin(imn)/Nfp_raw, surf(isurf)%bim(imn), surf(isurf)%Zbs(imn), &
surf(isurf)%bin(imn/Nfp_raw), surf(isurf)%bim(imn), surf(isurf)%Rbs(imn), &
surf(isurf)%bin(imn)/Nfp_raw, surf(isurf)%bim(imn), surf(isurf)%Rbs(imn), &
surf(isurf)%bin(imn)/Nfp_raw, surf(isurf)%bim(imn), surf(isurf)%Zbc(imn) ! wall is read as plasma boundary
enddo
do imn = 1, Nbf
Expand Down