boundary_absorbing.f90
Go to the documentation of this file.
1 !#############################################################################
2 !# #
3 !# fosite - 3D hydrodynamical simulation program #
4 !# module: boundary_absorbing.f90 #
5 !# #
6 !# Copyright (C) 2009-2014 #
7 !# Tobias Illenseer <tillense@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 !----------------------------------------------------------------------------!
38 !----------------------------------------------------------------------------!
41  USE mesh_base_mod
44  USE common_dict
45  IMPLICIT NONE
46  !--------------------------------------------------------------------------!
47  PRIVATE
49  REAL, DIMENSION(:,:,:), ALLOCATABLE :: xvar, & !< characteristic variables for absorbing bc
50  lambda
51  CONTAINS
53  PROCEDURE :: setboundarydata
54  PROCEDURE :: finalize
55  END TYPE boundary_absorbing
56  CHARACTER(LEN=32), PARAMETER :: boundcond_name = "absorbing"
57  !--------------------------------------------------------------------------!
58  PUBLIC :: &
60  !--------------------------------------------------------------------------!
61 
62 CONTAINS
63 
68  SUBROUTINE initboundary_absorbing(this,Mesh,Physics,dir,config)
69  IMPLICIT NONE
70  !------------------------------------------------------------------------!
71  CLASS(Boundary_absorbing), INTENT(INOUT) :: this
72  CLASS(Mesh_base), INTENT(IN) :: Mesh
73  CLASS(Physics_base), INTENT(IN) :: Physics
74  TYPE(Dict_TYP), POINTER, INTENT(IN) :: config
75  INTEGER, INTENT(IN) :: dir
76  !------------------------------------------------------------------------!
77  INTEGER :: err
78  !------------------------------------------------------------------------!
79  CALL this%InitBoundary(mesh,physics,absorbing,boundcond_name,dir,config)
80  ! check if physics supports absorbing boundary conditions
81  IF (.NOT.physics%supports_absorbing) &
82  CALL this%Error("InitBoundary_absorbing", &
83  "boundary condition not supported for this type of physics")
84  SELECT CASE(this%direction%GetType())
85  CASE(west,east)
86  ALLOCATE(this%xvar(mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,physics%VNUM), &
87  this%lambda(mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,physics%VNUM), &
88  stat=err)
89  CASE(south,north)
90  ALLOCATE(this%xvar(mesh%IMIN:mesh%IMAX,mesh%KMIN:mesh%KMAX,physics%VNUM), &
91  this%lambda(mesh%IMIN:mesh%IMAX,mesh%KMIN:mesh%KMAX,physics%VNUM), &
92  stat=err)
93  CASE(bottom,top)
94  ALLOCATE(this%xvar(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,physics%VNUM), &
95  this%lambda(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,physics%VNUM), &
96  stat=err)
97  END SELECT
98  IF (err.NE.0) &
99  CALL this%Error("InitBoundary_absorbing", "Unable to allocate memory.")
100  ! initialize the data
101  this%xvar(:,:,:) = 0.0
102  this%lambda(:,:,:) = 0.0
103  END SUBROUTINE initboundary_absorbing
104 
105 
114  PURE SUBROUTINE setboundarydata(this,Mesh,Physics,time,pvar)
115  IMPLICIT NONE
116  !------------------------------------------------------------------------!
117  CLASS(boundary_absorbing),INTENT(INOUT) :: this
118  CLASS(mesh_base), INTENT(IN) :: mesh
119  CLASS(physics_base), INTENT(IN) :: physics
120  REAL, INTENT(IN) :: time
121  CLASS(marray_compound), INTENT(INOUT) :: pvar
122  !------------------------------------------------------------------------!
123  INTEGER :: i,j,k
124  !------------------------------------------------------------------------!
125  SELECT CASE(this%direction%GetType())
126  CASE(west)
127  ! get characteristic variables
128  CALL physics%CalculateCharSystemX(mesh,mesh%IMIN,mesh%IMIN+1,pvar,this%lambda,this%xvar)
129  ! set characteristic variables to zero for all incomming waves
130  WHERE (this%lambda(:,:,:).GE.0.0)
131  this%xvar(:,:,:) = 0.0
132  END WHERE
133  ! transform back to primitive variables at the boundary
134  CALL physics%CalculateBoundaryDataX(mesh,mesh%IMIN,mesh%IMIN-1,this%xvar,pvar)
135  ! copy data to outer ghost cells
136  DO i=2,mesh%GINUM
137  pvar%data4d(mesh%IMIN-i,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,:) &
138  = pvar%data4d(mesh%IMIN-1,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,:)
139  END DO
140  CASE(east)
141  ! get characteristic variables
142  CALL physics%CalculateCharSystemX(mesh,mesh%IMAX,mesh%IMAX-1,pvar,this%lambda,this%xvar)
143  ! set characteristic variables to zero for all incomming waves
144  WHERE (this%lambda(:,:,:).LE.0.0)
145  this%xvar(:,:,:) = 0.0
146  END WHERE
147  ! transform back to primitive variables at the boundary
148  CALL physics%CalculateBoundaryDataX(mesh,mesh%IMAX,mesh%IMAX+1,this%xvar,pvar)
149  ! copy data to outer ghost cells
150  DO i=2,mesh%GINUM
151  pvar%data4d(mesh%IMAX+i,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,:) &
152  = pvar%data4d(mesh%IMAX+1,mesh%JMIN:mesh%JMAX,mesh%KMIN:mesh%KMAX,:)
153  END DO
154  CASE(south)
155  ! get characteristic variables
156  CALL physics%CalculateCharSystemY(mesh,mesh%JMIN,mesh%JMIN+1,pvar,this%lambda,this%xvar)
157  ! set characteristic variables to zero for all incomming waves
158  WHERE (this%lambda(:,:,:).GE.0.0)
159  this%xvar(:,:,:) = 0.0
160  END WHERE
161  ! transform back to primitive variables at the boundary
162  CALL physics%CalculateBoundaryDataY(mesh,mesh%JMIN,mesh%JMIN-1,this%xvar,pvar)
163  ! copy data to outer ghost cells
164  DO j=2,mesh%GJNUM
165  pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN-j,mesh%KMIN:mesh%KMAX,:) &
166  = pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN-1,mesh%KMIN:mesh%KMAX,:)
167  END DO
168  CASE(north)
169  ! get characteristic variables
170  CALL physics%CalculateCharSystemY(mesh,mesh%JMAX,mesh%JMAX-1,pvar,this%lambda,this%xvar)
171  ! set characteristic variables to zero for all incomming waves
172  WHERE (this%lambda(:,:,:).LE.0.0)
173  this%xvar(:,:,:) = 0.0
174  END WHERE
175  ! transform back to primitive variables at the boundary
176  CALL physics%CalculateBoundaryDataY(mesh,mesh%JMAX,mesh%JMAX+1,this%xvar,pvar)
177  ! copy data to outer ghost cells
178  DO j=2,mesh%GJNUM
179  pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMAX+j,mesh%KMIN:mesh%KMAX,:) &
180  = pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMAX+1,mesh%KMIN:mesh%KMAX,:)
181  END DO
182  CASE(bottom)
183  ! get characteristic variables
184  CALL physics%CalculateCharSystemZ(mesh,mesh%KMIN,mesh%KMIN+1,pvar,this%lambda,this%xvar)
185  ! set characteristic variables to zero for all incomming waves
186  WHERE (this%lambda(:,:,:).GE.0.0)
187  this%xvar(:,:,:) = 0.0
188  END WHERE
189  ! transform back to primitive variables at the boundary
190  CALL physics%CalculateBoundaryDataZ(mesh,mesh%KMIN,mesh%KMIN-1,this%xvar,pvar)
191  ! copy data to outer ghost cells
192  DO k=2,mesh%GKNUM
193  pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN-k,:) &
194  = pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMIN-1,:)
195  END DO
196  CASE(top)
197  ! get characteristic variables
198  CALL physics%CalculateCharSystemZ(mesh,mesh%KMAX,mesh%KMAX-1,pvar,this%lambda,this%xvar)
199  ! set characteristic variables to zero for all incomming waves
200  WHERE (this%lambda(:,:,:).LE.0.0)
201  this%xvar(:,:,:) = 0.0
202  END WHERE
203  ! transform back to primitive variables at the boundary
204  CALL physics%CalculateBoundaryDataZ(mesh,mesh%KMAX,mesh%KMAX+1,this%xvar,pvar)
205  ! copy data to outer ghost cells
206  DO k=2,mesh%GKNUM
207  pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMAX+k,:) &
208  = pvar%data4d(mesh%IMIN:mesh%IMAX,mesh%JMIN:mesh%JMAX,mesh%KMAX+1,:)
209  END DO
210  END SELECT
211  END SUBROUTINE setboundarydata
212 
213 
215  SUBROUTINE finalize(this)
216  IMPLICIT NONE
217  !------------------------------------------------------------------------!
218  CLASS(Boundary_absorbing), INTENT(INOUT) :: this
219  !------------------------------------------------------------------------!
220  DEALLOCATE(this%xvar,this%lambda)
221  CALL this%Finalize_base()
222  END SUBROUTINE finalize
223 
224 END MODULE boundary_absorbing_mod
subroutine finalize(this)
Destructor for absorbing boundary conditions.
type(logging_base), save this
derived class for compound of mesh arrays
pure subroutine setboundarydata(this, Mesh, Physics, time, pvar)
Applies the absorbing boundary condition.
character(len=32), parameter boundcond_name
subroutine initboundary_absorbing(this, Mesh, Physics, dir, config)
Constructor for absorbing boundary conditions.
Boundary module for absorbing (non-reflecting) conditions.
named integer constants for flavour of state vectors
Basic physics module.
Dictionary for generic data types.
Definition: common_dict.f90:61