timedisc_dormand_prince.f90
Go to the documentation of this file.
1 !#############################################################################
2 !# #
3 !# fosite - 3D hydrodynamical simulation program #
4 !# module: timedisc_dormand_prince.f90 #
5 !# #
6 !# Copyright (C) 2013 #
7 !# Björn Sperling <sperling@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 !----------------------------------------------------------------------------!
42  USE mesh_base_mod
43  USE fluxes_base_mod
48  USE common_dict
49  IMPLICIT NONE
50  !--------------------------------------------------------------------------!
51  PRIVATE
53  CONTAINS
55  PROCEDURE :: setbutchertableau
56  PROCEDURE :: finalize
58  CHARACTER(LEN=32), PARAMETER :: odesolver_name = "Dormand-Prince method"
59 
60  !--------------------------------------------------------------------------!
61  PUBLIC :: &
62  ! types
64  !--------------------------------------------------------------------------!
65 
66 CONTAINS
67 
68  SUBROUTINE inittimedisc_dormand_prince(this,Mesh,Physics,config,IO)
69  IMPLICIT NONE
70  !------------------------------------------------------------------------!
71  CLASS(timedisc_dormand_prince), 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(5)
82  this%m = 7
83  CASE DEFAULT
84  CALL this%Error("InitTimedisc_dormand_prince","only order 5 supported")
85  END SELECT
86 
87  CALL this%InitTimedisc(mesh,physics,config,io,dormand_prince,odesolver_name)
88  END SUBROUTINE inittimedisc_dormand_prince
89 
91  SUBROUTINE setbutchertableau(this)
92  IMPLICIT NONE
93  !------------------------------------------------------------------------!
94  CLASS(timedisc_dormand_prince) :: this
95  !------------------------------------------------------------------------!
96  SELECT CASE(this%GetOrder())
97  CASE(5)
98  this%b_high = (/ 35./384.,0.,500./1113.,125./192.,-2187./6784.,11./84.,0. /)
99  this%b_low = (/ 5179./57600., 0., 7571./16695., 393./640., &
100  -92097./339200., 187./2100., 1./40. /)
101  this%c = (/ 0., .2, 3./10., 4./5., 8./9., 1., 1. /)
102  this%a = transpose(reshape((/ &
103  0., 0., 0., 0., 0., 0., 0., &
104  .2, 0., 0., 0., 0., 0., 0., &
105  3./40., 9./40., 0., 0., 0., 0., 0., &
106  44./45., -56./15., 32./9., 0., 0., 0., 0., &
107  19372./6561., -25360./2187., 64448./6561., -212./729., 0., 0., 0.,&
108  9017./3168.,-355./33., 46732./5247.,49./176.,-5103./18656.,0.,0., &
109  35./384.,0.,500./1113.,125./192.,-2187./6784.,11./84.,0. /),&
110  (/this%m,this%m/)))
111  CASE DEFAULT
112  CALL this%Error("timedisc_dormand_prince::SetButcherTableau","only order 5 supported")
113  END SELECT
114  END SUBROUTINE setbutchertableau
115 
116  SUBROUTINE finalize(this)
117  IMPLICIT NONE
118  !-----------------------------------------------------------------------!
119  CLASS(timedisc_dormand_prince) :: this
120  !-----------------------------------------------------------------------!
121  CALL this%timedisc_rkfehlberg%Finalize()
122  END SUBROUTINE
generic source terms module providing functionaly common to all source terms
character(len=32), parameter odesolver_name
subroutine setbutchertableau(this)
set coefficients for dormand_prince scheme
subroutines for Dormand-Prince method
subroutine inittimedisc_dormand_prince(this, Mesh, Physics, config, IO)
Basic physics module.
Dictionary for generic data types.
Definition: common_dict.f90:61
subroutines for Runge-Kutta Fehlberg method
integer, parameter, public dormand_prince
base module for numerical flux functions
Definition: fluxes_base.f90:39