reconstruction_base.f90
Go to the documentation of this file.
1!############################################################################
2!# #
3!# fosite - 3D hydrodynamical simulation program #
4!# module: reconstruction_base.f90 #
5!# #
6!# Copyright (C) 2007-2012 #
7!# Tobias Illenseer <tillense@astrophysik.uni-kiel.de> #
8!# Jannes Klee <jklee@astrophysik.uni-kiel.de> #
9!# #
10!# This program is free software; you can redistribute it and/or modify #
11!# it under the terms of the GNU General Public License as published by #
12!# the Free Software Foundation; either version 2 of the License, or (at #
13!# your option) any later version. #
14!# #
15!# This program is distributed in the hope that it will be useful, but #
16!# WITHOUT ANY WARRANTY; without even the implied warranty of #
17!# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
18!# NON INFRINGEMENT. See the GNU General Public License for more #
19!# details. #
20!# #
21!# You should have received a copy of the GNU General Public License #
22!# along with this program; if not, write to the Free Software #
23!# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
24!# #
25!#############################################################################
26
27!----------------------------------------------------------------------------!
34!----------------------------------------------------------------------------!
40 USE common_dict
41 IMPLICIT NONE
42 !--------------------------------------------------------------------------!
43 PRIVATE
44 TYPE, ABSTRACT, EXTENDS (logging_base) :: reconstruction_base
46 CLASS(logging_base), ALLOCATABLE :: limiter
47 LOGICAL :: primcons
48 REAL :: limiter_param
49 CONTAINS
50 PROCEDURE :: initreconstruction
51 procedure(calculatestates), DEFERRED :: calculatestates
52 PROCEDURE :: primrecon
53 PROCEDURE :: finalize_base
54 procedure(finalize), DEFERRED :: finalize
55 END TYPE reconstruction_base
56
57 abstract INTERFACE
58 SUBROUTINE calculatestates(this,Mesh,Physics,rvar,rstates)
60 IMPLICIT NONE
61 CLASS(reconstruction_base), INTENT(INOUT) :: this
62 CLASS(mesh_base), INTENT(IN) :: Mesh
63 CLASS(physics_base), INTENT(IN) :: Physics
64 CLASS(marray_compound), INTENT(INOUT) :: rvar
65 CLASS(marray_compound), INTENT(INOUT) :: rstates
66 END SUBROUTINE
67 SUBROUTINE finalize(this)
69 IMPLICIT NONE
70 CLASS(reconstruction_base), INTENT(INOUT) :: this
71 END SUBROUTINE
72 END INTERFACE
73
74 !--------------------------------------------------------------------------!
75
78 INTEGER, PARAMETER :: constant = 1
79 INTEGER, PARAMETER :: linear = 2
80 !--------------------------------------------------------------------------!
81 PUBLIC :: &
82 ! types
84 ! constants
86 !--------------------------------------------------------------------------!
87
88CONTAINS
89
90
92 SUBROUTINE initreconstruction(this,Mesh,Physics,config,IO,rtype,rname)
93 IMPLICIT NONE
94 !------------------------------------------------------------------------!
95 CLASS(reconstruction_base), INTENT(INOUT) :: this
96 CLASS(mesh_base), INTENT(IN) :: Mesh
97 CLASS(physics_base), INTENT(IN) :: Physics
98 TYPE(dict_typ), POINTER :: config,IO
99 INTEGER :: rtype
100 CHARACTER(LEN=32) :: rname
101 !------------------------------------------------------------------------!
102 INTEGER :: variables
103 CHARACTER(LEN=32) :: infostr
104 INTEGER :: order
105 !------------------------------------------------------------------------!
106 INTENT(IN) :: config,rtype,rname
107 !------------------------------------------------------------------------!
108 CALL this%InitLogging(rtype,rname)
109
110 ! check initialization of Mesh and Physics
111 IF (.NOT.mesh%Initialized().OR..NOT.physics%Initialized()) & ! TODO: Why should this not be seperated?
112 CALL this%Error("InitFluxes","mesh and/or physics module uninitialized")
113
114 ! set general reconstruction defaults
115 CALL getattr(config, "order", order, linear)
116
117 CALL getattr(config, "variables", variables, conservative)
118 SELECT CASE(variables)
119 CASE(primitive)
120 this%primcons = .true.
121 CASE DEFAULT
122 this%primcons = .false.
123 END SELECT
124
125 ! print some information
126 CALL this%Info(" RECONSTR-> order: " // trim(this%GetName()))
127 IF (primrecon(this)) THEN
128 WRITE (infostr,'(A)') "primitive"
129 ELSE
130 WRITE (infostr,'(A)') "conservative"
131 END IF
132 CALL this%Info(" variables: " // trim(infostr))
133
134 END SUBROUTINE initreconstruction
135
136
137 PURE FUNCTION primrecon(this) RESULT(pc)
138 IMPLICIT NONE
139 !------------------------------------------------------------------------!
140 CLASS(reconstruction_base), INTENT(IN) :: this
141 LOGICAL :: pc
142 !------------------------------------------------------------------------!
143 pc = this%primcons
144 END FUNCTION primrecon
145
146
147 SUBROUTINE finalize_base(this)
148 IMPLICIT NONE
149 !------------------------------------------------------------------------!
150 CLASS(reconstruction_base), INTENT(INOUT) :: this
151 !------------------------------------------------------------------------!
152 IF (.NOT.this%Initialized()) &
153 CALL this%Error("CloseReconstruction","not initialized")
154 END SUBROUTINE finalize_base
155
elemental real function pc(tau, rho_s0, cs_inf, gamma)
Definition: agndisk.f90:515
Dictionary for generic data types.
Definition: common_dict.f90:61
type(logging_base), save this
Basic fosite module.
subroutine finalize(this)
Destructor of logging_base class.
derived class for compound of mesh arrays
basic mesh module
Definition: mesh_base.f90:72
Basic physics module.
@, public primitive
@, public conservative
base module for reconstruction process
pure logical function primrecon(this)
subroutine initreconstruction(this, Mesh, Physics, config, IO, rtype, rname)
Constructor of base reconstruction module.
integer, parameter, public constant
integer, parameter, public linear
common data structure
mesh data structure
Definition: mesh_base.f90:122