Skip to content

Commit

Permalink
Avoid stress calc at t=0 and use correct arrays for stress objs
Browse files Browse the repository at this point in the history
  • Loading branch information
bpatel2107 committed Oct 24, 2024
1 parent a368356 commit 404bd32
Show file tree
Hide file tree
Showing 3 changed files with 14 additions and 23 deletions.
31 changes: 11 additions & 20 deletions cgyro/src/cgyro_flux.f90
Original file line number Diff line number Diff line change
Expand Up @@ -220,27 +220,18 @@ subroutine cgyro_flux

if (stress_flag .eq. 1) then
stress(:, :, :, :) = 0.0
if (i_time .ne. 0) then

! Set up stress for each field
do i_field=1, n_field
call cgyro_nl_fftw_comm1_async
call cgyro_nl_fftw_comm2_async_stress(i_field)
call cgyro_nl_fftw
call cgyro_nl_fftw_comm1_r_stress(i_field)
end do

! Set up phi stress
call cgyro_nl_fftw_comm1_async
call cgyro_nl_fftw_comm2_async_stress(1)
call cgyro_nl_fftw
call cgyro_nl_fftw_comm1_r_stress(1)
! velocity + select ky=0

! Set up A|| stress
call cgyro_nl_fftw_comm1_async
call cgyro_nl_fftw_comm2_async_stress(2)
call cgyro_nl_fftw
call cgyro_nl_fftw_comm1_r_stress(2)

! Set up B|| stress
call cgyro_nl_fftw_comm1_async
call cgyro_nl_fftw_comm2_async_stress(3)
call cgyro_nl_fftw
call cgyro_nl_fftw_comm1_r_stress(3)

! stress(kx-theta, vel_coord, ky, phi)
end if
! stress(kx-theta, vel_coord, ky, phi)
stress_integrated_loc = 0.0
iv_loc = 0
do iv=nv1,nv2
Expand Down
2 changes: 1 addition & 1 deletion cgyro/src/cgyro_globals.F90
Original file line number Diff line number Diff line change
Expand Up @@ -250,7 +250,7 @@ module cgyro_globals
(/'bin.cgyro.kxky_n','bin.cgyro.kxky_e','bin.cgyro.kxky_v'/)
character(len=19), dimension(3) :: binfile_kxky_field = &
(/'bin.cgyro.kxky_phi ','bin.cgyro.kxky_apar','bin.cgyro.kxky_bpar'/)
character(len=19), dimension(3) :: binfile_stress = &
character(len=21), dimension(3) :: binfile_stress = &
(/'bin.cgyro.stress_phi ','bin.cgyro.stress_apar','bin.cgyro.stress_bpar'/)
character(len=20), dimension(3) :: binfile_lky_flux = &
(/'bin.cgyro.lky_flux_n','bin.cgyro.lky_flux_e','bin.cgyro.lky_flux_v'/)
Expand Down
4 changes: 2 additions & 2 deletions cgyro/src/cgyro_init_manager.F90
Original file line number Diff line number Diff line change
Expand Up @@ -309,8 +309,8 @@ subroutine cgyro_init_manager

if (stress_flag == 1) then
allocate(stress(nc,nv_loc,nt1:nt2,n_field))
allocate(stress_integrated_loc(n_radial,nt_loc,nt1:nt2,n_field))
allocate(stress_integrated(n_radial,nt_loc,nt1:nt2,n_field))
allocate(stress_integrated_loc(n_radial,n_theta,nt1:nt2,n_field))
allocate(stress_integrated(n_radial,n_theta,nt1:nt2,n_field))
#if defined(OMPGPU)
!$omp target enter data map(alloc:stress,stress_integrated,stress_integrated_loc)
#elif defined(_OPENACC)
Expand Down

0 comments on commit 404bd32

Please sign in to comment.