timedisc_ssprk.f90
Go to the documentation of this file.
1!#############################################################################
2!# #
3!# fosite - 3D hydrodynamical simulation program #
4!# module: timedisc_ssprk.f90 #
5!# #
6!# Copyright (C) 2013 #
7!# Manuel Jung <mjung@astrophysik.uni-kiel.de> #
8!# #
9!# This program is free software; you can redistribute it and/or modify #
10!# it under the terms of the GNU General Public License as published by #
11!# the Free Software Foundation; either version 2 of the License, or (at #
12!# your option) any later version. #
13!# #
14!# This program is distributed in the hope that it will be useful, but #
15!# WITHOUT ANY WARRANTY; without even the implied warranty of #
16!# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
17!# NON INFRINGEMENT. See the GNU General Public License for more #
18!# details. #
19!# #
20!# You should have received a copy of the GNU General Public License #
21!# along with this program; if not, write to the Free Software #
22!# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
23!# #
24!#############################################################################
25
26!----------------------------------------------------------------------------!
39!----------------------------------------------------------------------------!
48 USE common_dict
49 IMPLICIT NONE
50 !--------------------------------------------------------------------------!
51 PRIVATE
53 CONTAINS
54 PROCEDURE :: inittimedisc_ssprk
55 PROCEDURE :: setbutchertableau
56 PROCEDURE :: finalize
57 END TYPE
58 CHARACTER(LEN=32), PARAMETER :: odesolver_name = "SSP Runge-Kutta method"
59
60 !--------------------------------------------------------------------------!
61 PUBLIC :: &
62 ! types
64 !--------------------------------------------------------------------------!
65
66CONTAINS
67
68 SUBROUTINE inittimedisc_ssprk(this,Mesh,Physics,config,IO)
69 IMPLICIT NONE
70 !------------------------------------------------------------------------!
71 CLASS(timedisc_ssprk), INTENT(INOUT) :: this
72 CLASS(mesh_base), INTENT(INOUT) :: Mesh
73 CLASS(physics_base), INTENT(IN) :: Physics
74 TYPE(dict_typ), POINTER :: config,IO
75 !------------------------------------------------------------------------!
76 ! set default order
77 CALL getattr(config, "order", this%order, 5)
78
79 ! set number of coefficients
80 SELECT CASE(this%GetOrder())
81 CASE(3)
82 this%m = 3
83 CASE(5)
84 this%m = 5
85 CASE DEFAULT
86 CALL this%Error("timedisc_ssprk::InitTimedisc_ssprk","time order must be 3 or 5")
87 END SELECT
88
89 CALL this%InitTimedisc(mesh,physics,config,io,ssprk,odesolver_name)
90 END SUBROUTINE inittimedisc_ssprk
91
93 SUBROUTINE setbutchertableau(this)
94 IMPLICIT NONE
95 !------------------------------------------------------------------------!
96 CLASS(timedisc_ssprk) :: this
97 !------------------------------------------------------------------------!
98 SELECT CASE(this%GetOrder())
99 CASE(3) ! SSPRK order 3
100 CALL this%Warning("timedisc_ssprk", &
101 "This 3rd order embedded SSP RK scheme has been constructed from a second " // &
102 "third order SSP RK scheme by hand! It seems to work, but I am not sure, that one " // &
103 "is allowed to do so. An optimal embedded third order SSP RK scheme is described " // &
104 "in chapter 6.3 of the main reference (see above), but still needs to be translated into " // &
105 "a butchers tableau. => Better use the 5th order scheme!")
106 this%b_high = (/ 1./6., 1./6., 2./3. /)
107 this%b_low = (/ 0.5, 0.5, 0. /)
108 this%c = (/ 0., 1., 0.5 /)
109 this%a = transpose(reshape((/ &
110 0.0, 0.0, 0.0, &
111 1.0, 0.0, 0.0, &
112 0.25, 0.25, 0.0 /),(/this%m,this%m/)))
113 CASE(5) ! SSPRK order 5
114 this%b_high = (/ 0.17279, 0.094505, 0.12947, 0.29899, 0.30424 /)
115 this%b_low = (/ 0.12293, 0.31981, -0.15316, 0.31887, 0.39155 /)
116 this%c = (/ 0., 0.36717, 0.58522, 0.44156, 0.8464 /)
117 this%a = transpose(reshape((/ &
118 0., 0., 0., 0., 0., &
119 0.36717, 0., 0., 0., 0., &
120 0.26802, 0.3172, 0., 0., 0., &
121 0.11606, 0.13735, 0.18816, 0., 0., &
122 0.11212, 0.13269, 0.18178, 0.4198, 0. /),(/this%m,this%m/)))
123 CASE DEFAULT
124 CALL this%Error("timedisc_ssprk::SetButcherTableau","only order 3 or 5 supported")
125 END SELECT
126 END SUBROUTINE setbutchertableau
127
128 SUBROUTINE finalize(this)
129 IMPLICIT NONE
130 !------------------------------------------------------------------------!
131 CLASS(timedisc_ssprk) :: this
132 !------------------------------------------------------------------------!
133 CALL this%timedisc_rkfehlberg%Finalize()
134 END SUBROUTINE
135
136END MODULE timedisc_ssprk_mod
Dictionary for generic data types.
Definition: common_dict.f90:61
base module for numerical flux functions
Definition: fluxes_base.f90:39
basic mesh module
Definition: mesh_base.f90:72
Basic physics module.
generic source terms module providing functionaly common to all source terms
integer, parameter, public ssprk
subroutines for Runge-Kutta Fehlberg method
subroutine setbutchertableau(this)
set coefficients for RK-Fehlberg schemes
character(len=32), parameter odesolver_name
subroutines for strong stability preserving (SSP) Runge Kutta methods
subroutine inittimedisc_ssprk(this, Mesh, Physics, config, IO)
mesh data structure
Definition: mesh_base.f90:122