Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
Gyacomo
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package Registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Antoine Cyril David Hoffmann
Gyacomo
Commits
2489ac29
Commit
2489ac29
authored
2 years ago
by
Antoine Cyril David Hoffmann
Browse files
Options
Downloads
Patches
Plain Diff
release almost ready
parent
f69bcc08
No related branches found
No related tags found
No related merge requests found
Changes
24
Hide whitespace changes
Inline
Side-by-side
Showing
4 changed files
src/stepon.F90
+151
-128
151 additions, 128 deletions
src/stepon.F90
src/tesend.F90
+1
-1
1 addition, 1 deletion
src/tesend.F90
testcases/cyclone_example/fort.90
+0
-92
0 additions, 92 deletions
testcases/cyclone_example/fort.90
wk/analysis_gyacomo.m
+9
-8
9 additions, 8 deletions
wk/analysis_gyacomo.m
with
161 additions
and
229 deletions
src/stepon.F90
+
151
−
128
View file @
2489ac29
SUBROUTINE
stepon
! Advance one time step, (num_step=4 for Runge Kutta 4 scheme)
USE
advance_field_routine
,
ONLY
:
advance_time_level
,
advance_moments
USE
basic
,
ONLY
:
nlend
USE
closure
,
ONLY
:
apply_closure_model
USE
ghosts
,
ONLY
:
update_ghosts_moments
,
update_ghosts_EM
use
mpi
,
ONLY
:
MPI_COMM_WORLD
USE
time_integration
,
ONLY
:
ntimelevel
USE
prec_const
,
ONLY
:
dp
IMPLICIT
NONE
INTEGER
::
num_step
,
ierr
LOGICAL
::
mlend
! Advance one time step, (num_step=4 for Runge Kutta 4 scheme)
USE
advance_field_routine
,
ONLY
:
advance_time_level
,
advance_moments
USE
basic
,
ONLY
:
nlend
,
chrono_advf
,
chrono_pois
,&
chrono_chck
,
chrono_clos
,
chrono_ghst
,
start_chrono
,
stop_chrono
USE
closure
,
ONLY
:
apply_closure_model
USE
ghosts
,
ONLY
:
update_ghosts_moments
,
update_ghosts_EM
use
mpi
,
ONLY
:
MPI_COMM_WORLD
USE
time_integration
,
ONLY
:
ntimelevel
USE
prec_const
,
ONLY
:
dp
IMPLICIT
NONE
INTEGER
::
num_step
,
ierr
LOGICAL
::
mlend
DO
num_step
=
1
,
ntimelevel
! eg RK4 compute successively k1, k2, k3, k4
!----- BEFORE: All fields+ghosts are updated for step = n
!----- BEFORE: All fields+ghosts are updated for step = n
! Compute right hand side from current fields
! N_rhs(N_n, nadia_n, phi_n, S_n, Tcoll_n)
CALL
assemble_RHS
...
...
@@ -25,134 +26,156 @@ SUBROUTINE stepon
! Update moments with the hierarchy RHS (step by step)
! N_n+1 = N_n + N_rhs(n)
CALL
advance_moments
CALL
start_chrono
(
chrono_advf
)
CALL
advance_moments
CALL
stop_chrono
(
chrono_advf
)
! Closure enforcement of moments
CALL
apply_closure_model
CALL
start_chrono
(
chrono_clos
)
CALL
apply_closure_model
CALL
stop_chrono
(
chrono_clos
)
! Exchanges the ghosts values of N_n+1
CALL
update_ghosts_moments
CALL
start_chrono
(
chrono_ghst
)
CALL
update_ghosts_moments
CALL
stop_chrono
(
chrono_ghst
)
! Update electrostatic potential phi_n = phi(N_n+1) and potential vect psi
CALL
solve_EM_fields
CALL
update_ghosts_EM
CALL
start_chrono
(
chrono_pois
)
CALL
solve_EM_fields
CALL
stop_chrono
(
chrono_pois
)
! Update EM ghosts
CALL
start_chrono
(
chrono_ghst
)
CALL
update_ghosts_EM
CALL
stop_chrono
(
chrono_ghst
)
!- Check before next step
CALL
checkfield_all
()
CALL
start_chrono
(
chrono_chck
)
CALL
checkfield_all
()
CALL
stop_chrono
(
chrono_chck
)
IF
(
nlend
)
EXIT
! exit do loop
CALL
MPI_BARRIER
(
MPI_COMM_WORLD
,
ierr
)
!----- AFTER: All fields are updated for step = n+1
!----- AFTER: All fields are updated for step = n+1
END
DO
CONTAINS
!!!! Basic structure to simplify stepon
SUBROUTINE
assemble_RHS
USE
moments_eq_rhs
,
ONLY
:
compute_moments_eq_rhs
USE
collision
,
ONLY
:
compute_Capj
USE
nonlinear
,
ONLY
:
compute_Sapj
USE
processing
,
ONLY
:
compute_nadiab_moments_z_gradients_and_interp
IMPLICIT
NONE
! compute auxiliary non adiabatic moments array, gradients and interp
CALL
compute_nadiab_moments_z_gradients_and_interp
! compute nonlinear term ("if linear" is included inside)
CALL
compute_Sapj
! compute collision term ("if coll, if nu >0" is included inside)
CALL
compute_Capj
! compute the moments equation rhs
CALL
compute_moments_eq_rhs
END
SUBROUTINE
assemble_RHS
SUBROUTINE
checkfield_all
! Check all the fields for inf or nan
USE
utility
,
ONLY
:
is_nan
,
is_inf
USE
basic
,
ONLY
:
t0_checkfield
,
t1_checkfield
,
tc_checkfield
USE
fields
,
ONLY
:
phi
USE
grid
,
ONLY
:
local_na
,
local_np
,
local_nj
,
local_nky
,
local_nkx
,
local_nz
,&
ngp
,
ngj
,
ngz
USE
MPI
USE
time_integration
,
ONLY
:
updatetlevel
USE
model
,
ONLY
:
LINEARITY
IMPLICIT
NONE
LOGICAL
::
checkf_
INTEGER
::
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
REAL
::
sum_
! Execution time start
CALL
cpu_time
(
t0_checkfield
)
IF
(
LINEARITY
.NE.
'linear'
)
CALL
anti_aliasing
! ensure 0 mode for 2/3 rule
IF
(
LINEARITY
.NE.
'linear'
)
CALL
enforce_symmetry
! Enforcing symmetry on kx = 0
mlend
=
.FALSE.
IF
(
.NOT.
nlend
)
THEN
sum_
=
SUM
(
ABS
(
phi
))
checkf_
=
(
is_nan
(
sum_
,
'phi'
)
.OR.
is_inf
(
sum_
,
'phi'
))
mlend
=
(
mlend
.or.
checkf_
)
CALL
MPI_ALLREDUCE
(
mlend
,
nlend
,
1
,
MPI_LOGICAL
,
MPI_LOR
,
MPI_COMM_WORLD
,
ierr
)
ENDIF
! Execution time end
CALL
cpu_time
(
t1_checkfield
)
tc_checkfield
=
tc_checkfield
+
(
t1_checkfield
-
t0_checkfield
)
END
SUBROUTINE
checkfield_all
SUBROUTINE
anti_aliasing
USE
fields
,
ONLY
:
moments
USE
time_integration
,
ONLY
:
updatetlevel
USE
grid
,
ONLY
:
local_na
,
local_np
,
local_nj
,
local_nky
,
local_nkx
,
local_nz
,&
ngp
,
ngj
,
ngz
,
AA_x
,
AA_y
IMPLICIT
NONE
INTEGER
::
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
z
:
DO
iz
=
1
,
local_nz
+
ngz
kx
:
DO
ikx
=
1
,
local_nkx
ky
:
DO
iky
=
1
,
local_nky
j
:
DO
ij
=
1
,
local_nj
+
ngj
p
:
DO
ip
=
1
,
local_np
+
ngp
a
:
DO
ia
=
1
,
local_na
moments
(
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
,
updatetlevel
)
=
&
AA_x
(
ikx
)
*
AA_y
(
iky
)
*
moments
(
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
,
updatetlevel
)
ENDDO
a
ENDDO
p
ENDDO
j
ENDDO
ky
ENDDO
kx
ENDDO
z
END
SUBROUTINE
anti_aliasing
SUBROUTINE
enforce_symmetry
! Force X(k) = X(N-k)* complex conjugate symmetry
USE
fields
,
ONLY
:
phi
,
psi
,
moments
USE
time_integration
,
ONLY
:
updatetlevel
USE
grid
,
ONLY
:
local_na
,
local_np
,
local_nj
,
total_nkx
,
local_nz
,&
ngp
,
ngj
,
ngz
,
ikx0
,
iky0
,
contains_ky0
IMPLICIT
NONE
INTEGER
::
ia
,
ip
,
ij
,
ikx
,
iz
IF
(
contains_ky0
)
THEN
! moments
! z: DO iz = 1,local_nz+ngz
! j: DO ij=1,local_nj+ngj
! p: DO ip=1,local_np+ngp
! a: DO ia=1,local_na
DO
ikx
=
2
,
total_nkx
/
2
!symmetry at ky = 0
moments
(:,:,:,
iky0
,
ikx
,:,
updatetlevel
)
=
&
CONJG
(
moments
(:,:,:,
iky0
,
total_nkx
+2
-
ikx
,:,
updatetlevel
))
END
DO
! must be real at origin and top right
moments
(:,:,:,
iky0
,
ikx0
,:,
updatetlevel
)
=
&
REAL
(
moments
(:,:,:,
iky0
,
ikx0
,:,
updatetlevel
),
dp
)
! ENDDO a
! ENDDO p
! ENDDO j
! ENDDO z
! Phi
DO
ikx
=
2
,
total_nkx
/
2
!symmetry at ky = 0
CONTAINS
!!!! Basic structure to simplify stepon
SUBROUTINE
assemble_RHS
USE
basic
,
ONLY
:
chrono_mrhs
,
chrono_sapj
,
chrono_coll
,
chrono_grad
,
chrono_napj
,
start_chrono
,
stop_chrono
USE
moments_eq_rhs
,
ONLY
:
compute_moments_eq_rhs
USE
collision
,
ONLY
:
compute_Capj
USE
nonlinear
,
ONLY
:
compute_Sapj
USE
processing
,
ONLY
:
compute_nadiab_moments
,
compute_interp_z
,
compute_gradients_z
IMPLICIT
NONE
! compute auxiliary non adiabatic moments array
CALL
start_chrono
(
chrono_napj
)
CALL
compute_nadiab_moments
CALL
stop_chrono
(
chrono_napj
)
! compute gradients and interp
CALL
start_chrono
(
chrono_grad
)
CALL
compute_gradients_z
CALL
compute_interp_z
CALL
stop_chrono
(
chrono_grad
)
! compute nonlinear term ("if linear" is included inside)
CALL
start_chrono
(
chrono_sapj
)
CALL
compute_Sapj
CALL
stop_chrono
(
chrono_sapj
)
! compute collision term ("if coll, if nu >0" is included inside)
CALL
start_chrono
(
chrono_coll
)
CALL
compute_Capj
CALL
stop_chrono
(
chrono_coll
)
! compute the moments equation rhs
CALL
start_chrono
(
chrono_mrhs
)
CALL
compute_moments_eq_rhs
CALL
stop_chrono
(
chrono_mrhs
)
END
SUBROUTINE
assemble_RHS
SUBROUTINE
checkfield_all
! Check all the fields for inf or nan
USE
utility
,
ONLY
:
is_nan
,
is_inf
USE
fields
,
ONLY
:
phi
USE
MPI
USE
model
,
ONLY
:
LINEARITY
IMPLICIT
NONE
LOGICAL
::
checkf_
REAL
::
sum_
IF
(
LINEARITY
.NE.
'linear'
)
CALL
anti_aliasing
! ensure 0 mode for 2/3 rule
IF
(
LINEARITY
.NE.
'linear'
)
CALL
enforce_symmetry
! Enforcing symmetry on kx = 0
mlend
=
.FALSE.
IF
(
.NOT.
nlend
)
THEN
sum_
=
SUM
(
ABS
(
phi
))
checkf_
=
(
is_nan
(
sum_
,
'phi'
)
.OR.
is_inf
(
sum_
,
'phi'
))
mlend
=
(
mlend
.or.
checkf_
)
CALL
MPI_ALLREDUCE
(
mlend
,
nlend
,
1
,
MPI_LOGICAL
,
MPI_LOR
,
MPI_COMM_WORLD
,
ierr
)
ENDIF
END
SUBROUTINE
checkfield_all
SUBROUTINE
anti_aliasing
USE
fields
,
ONLY
:
moments
USE
time_integration
,
ONLY
:
updatetlevel
USE
grid
,
ONLY
:
local_na
,
local_np
,
local_nj
,
local_nky
,
local_nkx
,
local_nz
,&
ngp
,
ngj
,
ngz
,
AA_x
,
AA_y
IMPLICIT
NONE
INTEGER
::
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
z
:
DO
iz
=
1
,
local_nz
+
ngz
kx
:
DO
ikx
=
1
,
local_nkx
ky
:
DO
iky
=
1
,
local_nky
j
:
DO
ij
=
1
,
local_nj
+
ngj
p
:
DO
ip
=
1
,
local_np
+
ngp
a
:
DO
ia
=
1
,
local_na
moments
(
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
,
updatetlevel
)
=
&
AA_x
(
ikx
)
*
AA_y
(
iky
)
*
moments
(
ia
,
ip
,
ij
,
iky
,
ikx
,
iz
,
updatetlevel
)
ENDDO
a
ENDDO
p
ENDDO
j
ENDDO
ky
ENDDO
kx
ENDDO
z
END
SUBROUTINE
anti_aliasing
SUBROUTINE
enforce_symmetry
! Force X(k) = X(N-k)* complex conjugate symmetry
USE
fields
,
ONLY
:
phi
,
psi
,
moments
USE
time_integration
,
ONLY
:
updatetlevel
USE
grid
,
ONLY
:
total_nkx
,
ikx0
,
iky0
,
contains_ky0
IMPLICIT
NONE
INTEGER
::
ikx
IF
(
contains_ky0
)
THEN
! moments
! z: DO iz = 1,local_nz+ngz
! j: DO ij=1,local_nj+ngj
! p: DO ip=1,local_np+ngp
! a: DO ia=1,local_na
DO
ikx
=
2
,
total_nkx
/
2
!symmetry at ky = 0
moments
(:,:,:,
iky0
,
ikx
,:,
updatetlevel
)
=
&
CONJG
(
moments
(:,:,:,
iky0
,
total_nkx
+2
-
ikx
,:,
updatetlevel
))
END
DO
! must be real at origin and top right
moments
(:,:,:,
iky0
,
ikx0
,:,
updatetlevel
)
=
&
REAL
(
moments
(:,:,:,
iky0
,
ikx0
,:,
updatetlevel
),
dp
)
! ENDDO a
! ENDDO p
! ENDDO j
! ENDDO z
! Phi
DO
ikx
=
2
,
total_nkx
/
2
!symmetry at ky = 0
phi
(
iky0
,
ikx
,:)
=
phi
(
iky0
,
total_nkx
+2
-
ikx
,:)
END
DO
! must be real at origin
phi
(
iky0
,
ikx0
,:)
=
REAL
(
phi
(
iky0
,
ikx0
,:),
dp
)
! Psi
DO
ikx
=
2
,
total_nkx
/
2
!symmetry at ky = 0
END
DO
! must be real at origin
phi
(
iky0
,
ikx0
,:)
=
REAL
(
phi
(
iky0
,
ikx0
,:),
dp
)
! Psi
DO
ikx
=
2
,
total_nkx
/
2
!symmetry at ky = 0
psi
(
iky0
,
ikx
,:)
=
psi
(
iky0
,
total_nkx
+2
-
ikx
,:)
END
DO
! must be real at origin
psi
(
iky0
,
ikx0
,:)
=
REAL
(
psi
(
iky0
,
ikx0
,:),
dp
)
ENDIF
END
SUBROUTINE
enforce_symmetry
END
DO
! must be real at origin
psi
(
iky0
,
ikx0
,:)
=
REAL
(
psi
(
iky0
,
ikx0
,:),
dp
)
ENDIF
END
SUBROUTINE
enforce_symmetry
END
SUBROUTINE
stepon
This diff is collapsed.
Click to expand it.
src/tesend.F90
+
1
−
1
View file @
2489ac29
...
...
@@ -44,7 +44,7 @@ SUBROUTINE tesend
!________________________________________________________________________________
! 4. Test on run time
CALL
cpu_time
(
tnow
)
mlend
=
(
1.1
*
(
tnow
-
start
))
.GT.
maxruntime
mlend
=
(
1.1
*
(
tnow
-
chrono_runt
%
t
start
))
.GT.
maxruntime
CALL
mpi_allreduce
(
mlend
,
nlend
,
1
,
MPI_LOGICAL
,
MPI_LOR
,
MPI_COMM_WORLD
,
ierr
)
...
...
This diff is collapsed.
Click to expand it.
testcases/cyclone_example/fort.90
deleted
100644 → 0
+
0
−
92
View file @
f69bcc08
&BASIC
nrun = 99999999
dt = 0.01
tmax = 25
maxruntime = 356400
job2load = -1
/
&GRID
pmax = 4
jmax = 2
Nx = 128
Lx = 120
Ny = 48
Ly = 120
Nz = 16
SG = .f.
Nexc = 0
/
&GEOMETRY
geom = 's-alpha'
q0 = 1.4
shear = 0.8
eps = 0.18
kappa = 1.0
s_kappa= 0.0
delta = 0.0
s_delta= 0.0
zeta = 0.0
s_zeta = 0.0
parallel_bc = 'dirichlet'
shift_y= 0.0
/
&OUTPUT_PAR
dtsave_0d = 0.1
dtsave_1d = -1
dtsave_2d = -1
dtsave_3d = 1
dtsave_5d = 10
write_doubleprecision = .f.
write_gamma = .t.
write_hf = .t.
write_phi = .t.
write_Na00 = .t.
write_Napj = .t.
write_dens = .t.
write_fvel = .t.
write_temp = .t.
/
&MODEL_PAR
! Collisionality
CLOS = 0
NL_CLOS = 0
LINEARITY = 'nonlinear'
Na = 1 ! number of species
mu_x = 1.0
mu_y = 1.0
N_HD = 4
mu_z = 2.0
HYP_V = 'hypcoll'
mu_p = 0.0
mu_j = 0.0
nu = 0.001
beta = 0.0
ADIAB_E = .t.
tau_e = 1.0
/
&SPECIES
! ions
name_ = 'ions'
tau_ = 1.0
sigma_= 1.0
q_ = 1.0
k_N_ = 2.22
k_T_ = 6.96
/
&COLLISION_PAR
collision_model = 'DG' !DG/SG/PA/LD (dougherty, sugama, pitch angle, landau)
GK_CO = .f.
INTERSPECIES = .true.
mat_file = 'gk_sugama_P_20_J_10_N_150_kpm_8.0.h5'
/
&INITIAL_CON
INIT_OPT = 'blob'
ACT_ON_MODES = 'donothing'
init_background = 0.0
init_noiselvl = 0.005
iseed = 42
/
&TIME_INTEGRATION_PAR
numerical_scheme = 'RK4'
/
This diff is collapsed.
Click to expand it.
wk/analysis_gyacomo.m
+
9
−
8
View file @
2489ac29
%% UNCOMMENT FOR TUTORIAL
gyacomodir
=
pwd
;
gyacomodir
=
gyacomodir
(
1
:
end
-
2
);
% get code directory
resdir
=
'../testcases/
zpinch
_example/'
;
%Name of the directory where the results are located
resdir
=
'../testcases/
cyclone
_example/'
;
%Name of the directory where the results are located
PARTITION
=
''
;
JOBNUMMIN
=
0
0
;
JOBNUMMAX
=
1
0
;
JOBNUMMIN
=
0
3
;
JOBNUMMAX
=
0
3
;
%%
addpath
(
genpath
([
gyacomodir
,
'matlab'
]))
% ... add
addpath
(
genpath
([
gyacomodir
,
'matlab/plot'
]))
% ... add
...
...
@@ -80,7 +80,7 @@ options.NAME = '\phi';
% options.NAME = 'Q_x';
% options.NAME = 'n_i';
% options.NAME = 'n_i-n_e';
options
.
PLAN
=
'
x
y'
;
options
.
PLAN
=
'
kxk
y'
;
% options.NAME = 'f_i';
% options.PLAN = 'sx';
options
.
COMP
=
'avg'
;
...
...
@@ -93,27 +93,28 @@ options.RESOLUTION = 256;
create_film
(
data
,
options
,
'.gif'
)
end
if
0
if
1
%% fields snapshots
% Options
options
.
INTERP
=
0
;
options
.
POLARPLOT
=
0
;
options
.
AXISEQUAL
=
0
;
options
.
NORMALIZE
=
0
;
options
.
NAME
=
'\phi'
;
%
options.NAME = '\phi';
% options.NAME = '\psi';
% options.NAME = '\omega_z';
% options.NAME = 'T_i';
% options.NAME = 'n_i';
% options.NAME = '\phi^{NZ}';
%
options.NAME = 'N_i^{00}';
options
.
NAME
=
'N_i^{00}'
;
% options.NAME = 'N_i^{00}-N_e^{00}';
% options.NAME = 's_{Ex}';
% options.NAME = 'Q_x';
% options.NAME = 'k^2n_e';
options
.
PLAN
=
'
x
y'
;
options
.
PLAN
=
'
kxk
y'
;
options
.
COMP
=
'avg'
;
options
.
TIME
=
[
15
30
50
];
% options.TIME = [49 50 51];
options
.
TIME
=
data
.
Ts3D
(
51
:
55
);
options
.
RESOLUTION
=
256
;
data
.
a
=
data
.
EPS
*
2e3
;
...
...
This diff is collapsed.
Click to expand it.
Prev
1
2
Next
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment