Skip to content
GitLab
Next
Projects
Groups
Snippets
Help
Loading...
Help
See what's new at GitLab
4
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
A
afivo
Project overview
Project overview
Details
Activity
Releases
Repository
Repository
Files
Commits
Branches
Tags
Contributors
Graph
Compare
Locked Files
Issues
1
Issues
1
List
Boards
Labels
Service Desk
Milestones
Iterations
Merge Requests
0
Merge Requests
0
Requirements
Requirements
List
CI / CD
CI / CD
Pipelines
Jobs
Schedules
Test Cases
Security & Compliance
Security & Compliance
Dependency List
License Compliance
Operations
Operations
Incidents
Environments
Packages & Registries
Packages & Registries
Container Registry
Analytics
Analytics
CI / CD
Code Review
Insights
Issue
Repository
Value Stream
Wiki
Wiki
Snippets
Snippets
Members
Members
Collapse sidebar
Close sidebar
Activity
Graph
Create a new issue
Jobs
Commits
Issue Boards
Open sidebar
MD-CWI-NL
afivo
Commits
e9430ad0
Commit
e9430ad0
authored
Apr 07, 2020
by
Jannis Teunissen
Browse files
Options
Browse Files
Download
Email Patches
Plain Diff
Pass both i_step and n_steps to forward_euler method
parent
47899d9a
Changes
3
Hide whitespace changes
Inline
Side-by-side
Showing
3 changed files
with
18 additions
and
11 deletions
+18
-11
examples/advection_v2.f90
examples/advection_v2.f90
+3
-2
examples/euler_gas_dynamics.f90
examples/euler_gas_dynamics.f90
+3
-2
src/m_af_advance.f90
src/m_af_advance.f90
+12
-7
No files found.
examples/advection_v2.f90
View file @
e9430ad0
...
...
@@ -207,7 +207,8 @@ contains
end
select
end
function
solution
subroutine
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
s_deriv
,
s_prev
,
s_out
,
istep
)
subroutine
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
s_deriv
,
s_prev
,
s_out
,
&
i_step
,
n_steps
)
type
(
af_t
),
intent
(
inout
)
::
tree
real
(
dp
),
intent
(
in
)
::
dt
real
(
dp
),
intent
(
in
)
::
time
...
...
@@ -215,7 +216,7 @@ contains
integer
,
intent
(
in
)
::
s_deriv
integer
,
intent
(
in
)
::
s_prev
integer
,
intent
(
in
)
::
s_out
integer
,
intent
(
in
)
::
i
step
integer
,
intent
(
in
)
::
i
_step
,
n_steps
real
(
dp
)
::
wmax
(
NDIM
)
call
flux_generic_tree
(
tree
,
1
,
[
i_phi
+
s_deriv
],
[
i_flux
],
wmax
,
&
...
...
examples/euler_gas_dynamics.f90
View file @
e9430ad0
...
...
@@ -162,7 +162,8 @@ program KT_euler
contains
subroutine
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
s_deriv
,
s_prev
,
s_out
,
istep
)
subroutine
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
s_deriv
,
s_prev
,
s_out
,
&
i_step
,
n_steps
)
type
(
af_t
),
intent
(
inout
)
::
tree
real
(
dp
),
intent
(
in
)
::
dt
real
(
dp
),
intent
(
out
)
::
dt_lim
...
...
@@ -170,7 +171,7 @@ contains
integer
,
intent
(
in
)
::
s_deriv
integer
,
intent
(
in
)
::
s_prev
integer
,
intent
(
in
)
::
s_out
integer
,
intent
(
in
)
::
i
step
integer
,
intent
(
in
)
::
i
_step
,
n_steps
real
(
dp
)
::
wmax
(
NDIM
)
call
flux_generic_tree
(
tree
,
n_vars
,
variables
+
s_deriv
,
fluxes
,
wmax
,
&
...
...
src/m_af_advance.f90
View file @
e9430ad0
...
...
@@ -19,7 +19,8 @@ module m_af_advance
integer
,
parameter
::
req_copies
(
n_integrators
)
=
af_advance_num_steps
interface
subroutine
subr_feuler
(
tree
,
dt
,
dt_lim
,
time
,
s_deriv
,
s_prev
,
s_out
,
istep
)
subroutine
subr_feuler
(
tree
,
dt
,
dt_lim
,
time
,
s_deriv
,
s_prev
,
s_out
,
&
i_step
,
n_steps
)
import
type
(
af_t
),
intent
(
inout
)
::
tree
real
(
dp
),
intent
(
in
)
::
dt
!< Time step
...
...
@@ -28,7 +29,8 @@ module m_af_advance
integer
,
intent
(
in
)
::
s_deriv
!< State to compute derivatives from
integer
,
intent
(
in
)
::
s_prev
!< Previous state
integer
,
intent
(
in
)
::
s_out
!< Output state
integer
,
intent
(
in
)
::
istep
!< Step of the integrator
integer
,
intent
(
in
)
::
i_step
!< Step of the integrator
integer
,
intent
(
in
)
::
n_steps
!< Total number of steps
end
subroutine
subr_feuler
end
interface
...
...
@@ -45,6 +47,7 @@ contains
integer
,
intent
(
in
)
::
i_cc
(:)
!< Index of cell-centered variables
integer
,
intent
(
in
)
::
time_integrator
procedure
(
subr_feuler
)
::
forward_euler
integer
::
n_steps
if
(
time_integrator
<
1
.or.
time_integrator
>
n_integrators
)
&
error stop
"Invalid time integrator"
...
...
@@ -52,18 +55,20 @@ contains
if
(
any
(
tree
%
cc_num_copies
(
i_cc
)
<
req_copies
(
time_integrator
)))
&
error stop
"Not enough copies available"
n_steps
=
af_advance_num_steps
(
time_integrator
)
select
case
(
time_integrator
)
case
(
af_forward_euler
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
0
,
0
,
0
,
1
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
0
,
0
,
0
,
1
,
n_steps
)
time
=
time
+
dt
case
(
af_midpoint_method
)
call
forward_euler
(
tree
,
0.5_dp
*
dt
,
dt_lim
,
time
,
0
,
0
,
1
,
1
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
1
,
0
,
0
,
2
)
call
forward_euler
(
tree
,
0.5_dp
*
dt
,
dt_lim
,
time
,
0
,
0
,
1
,
1
,
n_steps
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
1
,
0
,
0
,
2
,
n_steps
)
time
=
time
+
dt
case
(
af_heuns_method
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
0
,
0
,
1
,
1
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
0
,
0
,
1
,
1
,
n_steps
)
time
=
time
+
dt
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
1
,
1
,
1
,
2
)
call
forward_euler
(
tree
,
dt
,
dt_lim
,
time
,
1
,
1
,
1
,
2
,
n_steps
)
call
combine_substeps
(
tree
,
i_cc
,
2
,
[
0
,
1
],
[
0.5_dp
,
0.5_dp
],
0
)
end
select
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment