Skip to content

Commit aa8ce58

Browse files
Merge pull request #85 from FluidNumerics/patches
Patches
2 parents 950ea58 + 3a9a053 commit aa8ce58

12 files changed

+257
-45
lines changed

CMakeLists.txt

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -115,7 +115,7 @@ endif()
115115
# Check for openmp support if offload target is provided
116116
if( SELF_ENABLE_MULTITHREADING )
117117
if( "${CMAKE_Fortran_COMPILER_ID}" STREQUAL "GNU" )
118-
set( OFFLOAD_FLAGS "-ftree-parallelize-loops=${SELF_MULTITHREADING_NTHREADS} -fopt-info-loop" )
118+
set( OFFLOAD_FLAGS "-ftree-parallelize-loops=${SELF_MULTITHREADING_NTHREADS}" ) #-fopt-info-loop"
119119
elseif( "${CMAKE_Fortran_COMPILER_ID}" STREQUAL "Intel" )
120120
set( OFFLOAD_FLAGS "-parallel -qopt-report -qopt-report-phase=par" )
121121
elseif( "${CMAKE_Fortran_COMPILER_ID}" STREQUAL "IntelLLVM" )

examples/linear_euler2d_spherical_soundwave_closeddomain.f90

Lines changed: 13 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -30,6 +30,13 @@ program LinearEuler_Example
3030
use self_LinearEuler2D
3131

3232
implicit none
33+
real(prec),parameter :: rho0 = 1.225_prec
34+
real(prec),parameter :: rhoprime = 0.01_prec
35+
real(prec),parameter :: c = 1.0_prec ! Speed of sound
36+
real(prec),parameter :: Lr = 0.06_prec
37+
real(prec),parameter :: x0 = 0.5_prec
38+
real(prec),parameter :: y0 = 0.5_prec
39+
3340
character(SELF_INTEGRATOR_LENGTH),parameter :: integrator = 'rk3'
3441
integer,parameter :: controlDegree = 7
3542
integer,parameter :: targetDegree = 15
@@ -64,16 +71,12 @@ program LinearEuler_Example
6471
! Initialize the model
6572
call modelobj%Init(mesh,geometry)
6673
modelobj%prescribed_bcs_enabled = .false. ! Disables prescribed boundary condition block for gpu accelerated implementations
67-
! modelobj%rho0 = ! optional, set the reference density
68-
! modelobj%c = ! optional set the reference sound wave speed
69-
! modelobj%g = ! optional set the gravitational acceleration (y-direction)
74+
modelobj%tecplot_enabled = .false. ! Disables tecplot output
75+
modelobj%rho0 = rho0 ! optional, set the reference density
76+
modelobj%c = c ! optional set the reference sound wave speed
7077

7178
! Set the initial condition
72-
call modelobj%solution%SetEquation(1,'d = 0.0001*exp( -ln(2)*( (x-0.5)^2 + (y-0.5)^2 )/0.0036 )') ! density
73-
call modelobj%solution%SetEquation(2,'u = 0') ! u
74-
call modelobj%solution%SetEquation(3,'v = 0') ! v
75-
call modelobj%solution%SetEquation(4,'p = 0.0001*exp( -ln(2)*( (x-0.5)^2 + (y-0.5)^2 )/0.0036 )') ! pressure
76-
call modelobj%solution%SetInteriorFromEquation(geometry,0.0_prec)
79+
call modelobj%SphericalSoundWave(rhoprime,Lr,x0,y0)
7780

7881
call modelobj%WriteModel()
7982
call modelobj%IncrementIOCounter()
@@ -90,8 +93,8 @@ program LinearEuler_Example
9093

9194
ef = modelobj%entropy
9295

93-
if(ef > e0) then
94-
print*,"Error: Final absmax greater than initial absmax! ",e0,ef
96+
if(ef /= ef) then
97+
print*,"Error: Final entropy is inf or nan",ef
9598
stop 1
9699
endif
97100
! Clean up

pyself/model2d.py

Lines changed: 58 additions & 24 deletions
Original file line numberDiff line numberDiff line change
@@ -3,15 +3,14 @@
33

44

55
# Other SELF modules
6-
import self.geometry as geometry
6+
import pyself.geometry as geometry
77

88
class model:
99
def __init__(self):
1010
self.solution = None
1111
self.pvdata = None # Pyvista data
1212
self.varnames = None
1313
self.varunits = None
14-
self.nvar = 0
1514
self.geom = geometry.semquad()
1615

1716
def load(self, hdf5File):
@@ -26,26 +25,56 @@ def load(self, hdf5File):
2625

2726
if 'controlgrid' in list(f.keys()):
2827

29-
s = f['controlgrid/solution']
30-
self.solution = []
31-
i = 0
32-
for k in s.keys():
33-
if k == 'metadata':
34-
for v in s[f"{k}/name"].keys():
35-
self.solution.append( {"name":s[f"{k}/name/{v}"][()][0].decode('utf-8'), "units":s[f"{k}/units/{v}"][()][0].decode('utf-8'), 'data': None} )
36-
self.nvar = len(self.solution)
28+
controlgrid = f['controlgrid']
29+
for group_name in controlgrid.keys():
30+
31+
if( group_name == 'geometry' ):
32+
continue
33+
34+
group = controlgrid[group_name]
35+
# Create a list to hold data for this group
36+
setattr(self, group_name, [])
37+
group_data = getattr(self, group_name)
38+
print(f"Loading {group_name} group")
39+
40+
# Load metadata information
41+
if( 'metadata' in list(group.keys()) ):
42+
for v in group[f"metadata/name"].keys():
43+
44+
name = group[f"metadata/name/{v}"].asstr()[()][0]
45+
try:
46+
units = group[f"metadata/units/{v}"].asstr()[()][0]
47+
except:
48+
units = "error"
49+
50+
group_data.append({
51+
"name": name,
52+
"units": units,
53+
'data': None
54+
})
3755
else:
38-
d = s[k]
39-
N = d.shape[2]
40-
# Find index for this field
41-
i=0
42-
for sol in self.solution:
43-
if sol['name'] == k:
44-
break
45-
else:
46-
i+=1
47-
48-
self.solution[i]['data'] = da.from_array(d, chunks=(self.geom.daskChunkSize,N,N))
56+
print(f"Error: /controlgrid/{group_name}/metadata group not found in {hdf5File}.")
57+
return 1
58+
59+
for k in group.keys():
60+
k_decoded = k.encode('utf-8').decode('utf-8')
61+
if k == 'metadata':
62+
continue
63+
else:
64+
print(f"Loading {k_decoded} field")
65+
# Load the actual data
66+
d = group[k]
67+
N = d.shape[2]
68+
69+
# Find index for this field
70+
i = 0
71+
for sol in group_data:
72+
if sol['name'] == k_decoded:
73+
break
74+
else:
75+
i += 1
76+
77+
group_data[i]['data'] = da.from_array(d, chunks=(self.geom.daskChunkSize, N, N))
4978

5079
self.generate_pyvista()
5180

@@ -97,8 +126,13 @@ def generate_pyvista(self):
97126

98127
# Load fields into pvdata
99128
k = 0
100-
for var in self.solution:
101-
self.pvdata.point_data.set_array(var['data'].flatten(),var['name'])
102-
k+=1
129+
for attr in self.__dict__:
130+
if not attr in ['pvdata','varnames','varunits','geom']:
131+
controlgroup = getattr(self, attr)
132+
#print(f"Loading {attr} into pvdata")
133+
for var in controlgroup:
134+
# print(f"Loading {var['name']} into pvdata")
135+
self.pvdata.point_data.set_array(var['data'].flatten(),var['name'])
136+
k+=1
103137

104138
print(self.pvdata)

setup.py

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -11,7 +11,8 @@
1111
packages=['pyself'],
1212
install_requires=['h5py>=3.7.0',
1313
'dask',
14-
'pyvista'],
14+
'pyvista',
15+
'imageio[ffmpeg]'],
1516
classifiers=[
1617
'Development Status :: 1 - Planning',
1718
'Intended Audience :: Science/Research',

src/SELF_DGModel2D_t.f90

Lines changed: 39 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -64,6 +64,7 @@ module SELF_DGModel2D_t
6464
procedure :: SourceMethod => sourcemethod_DGModel2D_t
6565
procedure :: SetBoundaryCondition => setboundarycondition_DGModel2D_t
6666
procedure :: SetGradientBoundaryCondition => setgradientboundarycondition_DGModel2D_t
67+
procedure :: ReportMetrics => ReportMetrics_DGModel2D_t
6768

6869
procedure :: UpdateSolution => UpdateSolution_DGModel2D_t
6970

@@ -149,6 +150,43 @@ subroutine Free_DGModel2D_t(this)
149150

150151
endsubroutine Free_DGModel2D_t
151152

153+
subroutine ReportMetrics_DGModel2D_t(this)
154+
!! Base method for reporting the entropy of a model
155+
!! to stdout. Only override this procedure if additional
156+
!! reporting is needed. Alternatively, if you think
157+
!! additional reporting would be valuable for all models,
158+
!! open a pull request with modifications to this base
159+
!! method.
160+
implicit none
161+
class(DGModel2D_t),intent(inout) :: this
162+
! Local
163+
character(len=20) :: modelTime
164+
character(len=20) :: minv,maxv
165+
character(len=:),allocatable :: str
166+
integer :: ivar
167+
168+
! Copy the time and entropy to a string
169+
write(modelTime,"(ES16.7E3)") this%t
170+
171+
do ivar = 1,this%nvar
172+
write(maxv,"(ES16.7E3)") maxval(this%solution%interior(:,:,:,ivar))
173+
write(minv,"(ES16.7E3)") minval(this%solution%interior(:,:,:,ivar))
174+
175+
! Write the output to STDOUT
176+
open(output_unit,ENCODING='utf-8')
177+
write(output_unit,'(1x, A," : ")',ADVANCE='no') __FILE__
178+
str = 'tᵢ ='//trim(modelTime)
179+
write(output_unit,'(A)',ADVANCE='no') str
180+
str = ' | min('//trim(this%solution%meta(ivar)%name)// &
181+
'), max('//trim(this%solution%meta(ivar)%name)//') = '// &
182+
minv//" , "//maxv
183+
write(output_unit,'(A)',ADVANCE='yes') str
184+
enddo
185+
186+
call this%ReportUserMetrics()
187+
188+
endsubroutine ReportMetrics_DGModel2D_t
189+
152190
subroutine SetSolutionFromEqn_DGModel2D_t(this,eqn)
153191
implicit none
154192
class(DGModel2D_t),intent(inout) :: this
@@ -306,7 +344,7 @@ subroutine CalculateEntropy_DGModel2D_t(this)
306344
do iel = 1,this%geometry%nelem
307345
do j = 1,this%solution%interp%N+1
308346
do i = 1,this%solution%interp%N+1
309-
jac = this%geometry%J%interior(i,j,iel,1)
347+
jac = abs(this%geometry%J%interior(i,j,iel,1))
310348
s = this%solution%interior(i,j,iel,1:this%nvar)
311349
e = e+this%entropy_func(s)*jac
312350
enddo

src/SELF_DGModel3D_t.f90

Lines changed: 39 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -62,6 +62,7 @@ module SELF_DGModel3D_t
6262
procedure :: SourceMethod => sourcemethod_DGModel3D_t
6363
procedure :: SetBoundaryCondition => setboundarycondition_DGModel3D_t
6464
procedure :: SetGradientBoundaryCondition => setgradientboundarycondition_DGModel3D_t
65+
procedure :: ReportMetrics => ReportMetrics_DGModel3D_t
6566

6667
procedure :: UpdateSolution => UpdateSolution_DGModel3D_t
6768

@@ -147,6 +148,43 @@ subroutine Free_DGModel3D_t(this)
147148

148149
endsubroutine Free_DGModel3D_t
149150

151+
subroutine ReportMetrics_DGModel3D_t(this)
152+
!! Base method for reporting the entropy of a model
153+
!! to stdout. Only override this procedure if additional
154+
!! reporting is needed. Alternatively, if you think
155+
!! additional reporting would be valuable for all models,
156+
!! open a pull request with modifications to this base
157+
!! method.
158+
implicit none
159+
class(DGModel3D_t),intent(inout) :: this
160+
! Local
161+
character(len=20) :: modelTime
162+
character(len=20) :: minv,maxv
163+
character(len=:),allocatable :: str
164+
integer :: ivar
165+
166+
! Copy the time and entropy to a string
167+
write(modelTime,"(ES16.7E3)") this%t
168+
169+
do ivar = 1,this%nvar
170+
write(maxv,"(ES16.7E3)") maxval(this%solution%interior(:,:,:,:,ivar))
171+
write(minv,"(ES16.7E3)") minval(this%solution%interior(:,:,:,:,ivar))
172+
173+
! Write the output to STDOUT
174+
open(output_unit,ENCODING='utf-8')
175+
write(output_unit,'(1x, A," : ")',ADVANCE='no') __FILE__
176+
str = 'tᵢ ='//trim(modelTime)
177+
write(output_unit,'(A)',ADVANCE='no') str
178+
str = ' | min('//trim(this%solution%meta(ivar)%name)// &
179+
'), max('//trim(this%solution%meta(ivar)%name)//') = '// &
180+
minv//" , "//maxv
181+
write(output_unit,'(A)',ADVANCE='yes') str
182+
enddo
183+
184+
call this%ReportUserMetrics()
185+
186+
endsubroutine ReportMetrics_DGModel3D_t
187+
150188
subroutine SetSolutionFromEqn_DGModel3D_t(this,eqn)
151189
implicit none
152190
class(DGModel3D_t),intent(inout) :: this
@@ -305,7 +343,7 @@ subroutine CalculateEntropy_DGModel3D_t(this)
305343
do k = 1,this%solution%interp%N+1
306344
do j = 1,this%solution%interp%N+1
307345
do i = 1,this%solution%interp%N+1
308-
jac = this%geometry%J%interior(i,j,k,iel,1)
346+
jac = abs(this%geometry%J%interior(i,j,k,iel,1))
309347
s = this%solution%interior(i,j,k,iel,1:this%nvar)
310348
e = e+this%entropy_func(s)*jac
311349
enddo

src/SELF_LinearEuler2D_t.f90

Lines changed: 47 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -74,6 +74,7 @@ module self_LinearEuler2D_t
7474
procedure :: flux2d => flux2d_LinearEuler2D_t
7575
procedure :: riemannflux2d => riemannflux2d_LinearEuler2D_t
7676
!procedure :: source2d => source2d_LinearEuler2D_t
77+
procedure :: SphericalSoundWave => SphericalSoundWave_LinearEuler2D_t
7778

7879
endtype LinearEuler2D_t
7980

@@ -100,7 +101,7 @@ subroutine SetMetadata_LinearEuler2D_t(this)
100101
call this%solution%SetName(3,"v") ! y-velocity component
101102
call this%solution%SetUnits(3,"m⋅s⁻¹")
102103

103-
call this%solution%SetName(4,"pressure") ! Pressure
104+
call this%solution%SetName(4,"P") ! Pressure
104105
call this%solution%SetUnits(4,"kg⋅m⁻¹⋅s⁻²")
105106

106107
endsubroutine SetMetadata_LinearEuler2D_t
@@ -188,4 +189,49 @@ pure function riemannflux2d_LinearEuler2D_t(this,sL,sR,dsdx,nhat) result(flux)
188189

189190
endfunction riemannflux2d_LinearEuler2D_t
190191

192+
subroutine SphericalSoundWave_LinearEuler2D_t(this,rhoprime,Lr,x0,y0)
193+
!! This subroutine sets the initial condition for a weak blast wave
194+
!! problem. The initial condition is given by
195+
!!
196+
!! \begin{equation}
197+
!! \begin{aligned}
198+
!! \rho &= \rho_0 + \rho' \exp\left( -\ln(2) \frac{(x-x_0)^2 + (y-y_0)^2}{L_r^2} \right)
199+
!! u &= 0 \\
200+
!! v &= 0 \\
201+
!! E &= \frac{P_0}{\gamma - 1} + E \exp\left( -\ln(2) \frac{(x-x_0)^2 + (y-y_0)^2}{L_e^2} \right)
202+
!! \end{aligned}
203+
!! \end{equation}
204+
!!
205+
implicit none
206+
class(LinearEuler2D_t),intent(inout) :: this
207+
real(prec),intent(in) :: rhoprime,Lr,x0,y0
208+
! Local
209+
integer :: i,j,iEl
210+
real(prec) :: x,y,rho,r,E
211+
212+
print*,__FILE__," : Configuring weak blast wave initial condition. "
213+
print*,__FILE__," : rhoprime = ",rhoprime
214+
print*,__FILE__," : Lr = ",Lr
215+
print*,__FILE__," : x0 = ",x0
216+
print*,__FILE__," : y0 = ",y0
217+
218+
do concurrent(i=1:this%solution%N+1,j=1:this%solution%N+1, &
219+
iel=1:this%mesh%nElem)
220+
x = this%geometry%x%interior(i,j,iEl,1,1)-x0
221+
y = this%geometry%x%interior(i,j,iEl,1,2)-y0
222+
r = sqrt(x**2+y**2)
223+
224+
rho = (rhoprime)*exp(-log(2.0_prec)*r**2/Lr**2)
225+
226+
this%solution%interior(i,j,iEl,1) = rho
227+
this%solution%interior(i,j,iEl,2) = 0.0_prec
228+
this%solution%interior(i,j,iEl,3) = 0.0_prec
229+
this%solution%interior(i,j,iEl,4) = rho*this%c*this%c
230+
231+
enddo
232+
233+
call this%ReportMetrics()
234+
235+
endsubroutine SphericalSoundWave_LinearEuler2D_t
236+
191237
endmodule self_LinearEuler2D_t

src/SELF_Mesh_2D_t.f90

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -225,6 +225,8 @@ subroutine ResetBoundaryConditionType_Mesh2D_t(this,bcid)
225225
enddo
226226
enddo
227227

228+
call this%UpdateDevice()
229+
228230
endsubroutine ResetBoundaryConditionType_Mesh2D_t
229231

230232
subroutine UniformStructuredMesh_Mesh2D_t(this,nxPerTile,nyPerTile,nTileX,nTileY,dx,dy,bcids)

src/SELF_Mesh_3D_t.f90

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -301,6 +301,8 @@ subroutine ResetBoundaryConditionType_Mesh3D_t(this,bcid)
301301
enddo
302302
enddo
303303

304+
call this%UpdateDevice()
305+
304306
endsubroutine ResetBoundaryConditionType_Mesh3D_t
305307

306308
subroutine RecalculateFlip_Mesh3D_t(this)

0 commit comments

Comments
 (0)