FLASH-X
Doxygen Generated Documentation From Interface Source Code
Hydro_interface.F90
Go to the documentation of this file.
1!! NOTICE
2!! Copyright 2022 UChicago Argonne, LLC and contributors
3!!
4!! Licensed under the Apache License, Version 2.0 (the "License");
5!! you may not use this file except in compliance with the License.
6!!
7!! Unless required by applicable law or agreed to in writing, software
8!! distributed under the License is distributed on an "AS IS" BASIS,
9!! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
10!! See the License for the specific language governing permissions and
11!! limitations under the License.
12!!
13!! This is the header file for the hydro module that defines its
14!! public interfaces.
16#include "constants.h"
17#include "Simulation.h"
18
19 implicit none
20
21 interface
22 subroutine Hydro_computeDt (tileDesc, &
23 x, dx, uxgrid, &
24 y, dy, uygrid, &
25 z, dz, uzgrid, &
26 blkLimits,blkLimitsGC, &
27 solnData, &
28 dt_check, dt_minloc, extraInfo )
29 use Grid_tile, ONLY : Grid_tile_t
30 implicit none
31 type(Grid_tile_t), intent(IN) :: tileDesc
32 integer,dimension(LOW:HIGH,MDIM), intent(IN) :: blkLimits,blkLimitsGC
33 real, dimension(blkLimitsGC(LOW,IAXIS):blkLimitsGC(HIGH,IAXIS)), intent(IN) :: x, dx, uxgrid
34 real, dimension(blkLimitsGC(LOW,JAXIS):blkLimitsGC(HIGH,JAXIS)), intent(IN) :: y, dy, uygrid
35 real, dimension(blkLimitsGC(LOW,KAXIS):blkLimitsGC(HIGH,KAXIS)), intent(IN) :: z, dz, uzgrid
36 real,INTENT(INOUT) :: dt_check
37 integer,INTENT(INOUT) :: dt_minloc(5)
38 real, pointer :: solnData(:,:,:,:)
39 real, OPTIONAL,intent(INOUT) :: extraInfo
40 end subroutine Hydro_computeDt
41 end interface
42
43 interface
45 end subroutine Hydro_consolidateCFL
46 end interface
47
48 interface
49
51 implicit none
52 end subroutine Hydro_prepareBuffers
53 subroutine Hydro_freeBuffers()
54 implicit none
55 end subroutine Hydro_freeBuffers
56
57
58 end interface
59
60
61 interface
62 subroutine Hydro( timeEndAdv, dt, dtOld,sweepOrder )
63 real, INTENT(IN) :: timeEndAdv, dt, dtOld
64 integer, optional, INTENT(IN) :: sweepOrder
65 end subroutine Hydro
66 end interface
67
68 interface
69 subroutine Hydro_init()
70 end subroutine Hydro_init
71 end interface
72
73 interface
74 subroutine Hydro_finalize()
75 end subroutine Hydro_finalize
76 end interface
77
78 interface
79 subroutine Hydro_detectShock(solnData, shock, blkLimits, blkLimitsGC, &
80 guardCells, &
81 primaryCoord,secondCoord,thirdCoord)
82
83 integer, intent(IN), dimension(LOW:HIGH,MDIM) :: blkLimits, blkLimitsGC
84 integer, intent(IN) :: guardCells(MDIM)
85 real, pointer,dimension(:,:,:,:) :: solnData
86#ifdef FIXEDBLOCKSIZE
87 real, intent(out),dimension(GRID_ILO_GC:GRID_IHI_GC,&
88 GRID_JLO_GC:GRID_JHI_GC,&
89 GRID_KLO_GC:GRID_KHI_GC):: shock
90 real,intent(IN),dimension(GRID_ILO_GC:GRID_IHI_GC) :: primaryCoord
91 real,intent(IN),dimension(GRID_JLO_GC:GRID_JHI_GC) :: secondCoord
92 real,intent(IN),dimension(GRID_KLO_GC:GRID_KHI_GC) :: thirdCoord
93#else
94 real,intent(out), dimension(blkLimitsGC(LOW,IAXIS):blkLimitsGC(HIGH,IAXIS),&
95 blkLimitsGC(LOW,JAXIS):blkLimitsGC(HIGH,JAXIS),&
96 blkLimitsGC(LOW,KAXIS):blkLimitsGC(HIGH,KAXIS)) :: shock
97 real,intent(IN),dimension(blkLimitsGC(LOW,IAXIS):blkLimitsGC(HIGH,IAXIS)) :: primaryCoord
98 real,intent(IN),dimension(blkLimitsGC(LOW,JAXIS):blkLimitsGC(HIGH,JAXIS)) :: secondCoord
99 real,intent(IN),dimension(blkLimitsGC(LOW,KAXIS):blkLimitsGC(HIGH,KAXIS)) :: thirdCoord
100#endif
101
102 end subroutine Hydro_detectShock
103 end interface
104
105 interface
106 subroutine Hydro_shockStrength(solnData, shock, lo,hi,loHalo,hiHalo,&
107 primaryCoord,secondCoord,thirdCoord, &
108 threshold, mode)
109 implicit none
110 integer, intent(IN), dimension(1:MDIM) :: lo,hi,loHalo,hiHalo
111 real, pointer :: solnData(:,:,:,:)
112 real,intent(inout),dimension(loHalo(IAXIS):hiHalo(IAXIS),&
113 loHalo(JAXIS):hiHalo(JAXIS),&
114 loHalo(KAXIS):hiHalo(KAXIS)) :: shock
115 real,intent(IN),dimension(loHalo(IAXIS):hiHalo(IAXIS)) :: primaryCoord
116 real,intent(IN),dimension(loHalo(JAXIS):hiHalo(JAXIS)) :: secondCoord
117 real,intent(IN),dimension(loHalo(KAXIS):hiHalo(KAXIS)) :: thirdCoord
118 real, intent(IN) :: threshold
119 integer, intent(IN) :: mode
120 end subroutine Hydro_shockStrength
121 end interface
122
123
124 interface
126
127 end subroutine Hydro_sendOutputData
128 end interface
129
130 interface
132 implicit none
134 end interface
135
136 interface
137 subroutine Hydro_mapBcType(bcTypeToApply,bcTypeFromGrid,varIndex,gridDataStruct, &
138 axis,face,idest)
139 implicit none
140 integer, intent(OUT) :: bcTypeToApply
141 integer, intent(in) :: bcTypeFromGrid,varIndex,gridDataStruct,axis,face
142 integer,intent(IN),OPTIONAL:: idest
143 end subroutine Hydro_mapBcType
144 end interface
145
146end Module Hydro_interface