Skip to content

Commit f0ef368

Browse files
Merge pull request #294 from loganoz/292-acoustics-read-base-flow-from-file
292 acoustics read base flow from file
2 parents 35bd6a4 + 2552b6c commit f0ef368

File tree

4 files changed

+195
-18
lines changed

4 files changed

+195
-18
lines changed

Solver/src/libs/discretization/DGSEMClass.f90

Lines changed: 1 addition & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -144,9 +144,6 @@ SUBROUTINE ConstructDGSem( self, meshFileName_, controlVariables, &
144144
logical :: useWeightsPartition ! Partitioning mesh using DOF of elements as weights
145145
logical :: genMonitor
146146
logical :: isReconstruct=.false.
147-
#if defined(ACOUSTIC)
148-
real(kind=RP) :: QbaseUniform(1:NCONSB)
149-
#endif
150147
character(len=*), parameter :: TWOD_OFFSET_DIR_KEY = "2d mesh offset direction"
151148
procedure(UserDefinedInitialCondition_f) :: UserDefinedInitialCondition
152149
#if (!defined(NAVIERSTOKES))
@@ -401,9 +398,7 @@ SUBROUTINE ConstructDGSem( self, meshFileName_, controlVariables, &
401398
! --------------------
402399
!
403400
#if defined(ACOUSTIC)
404-
! start by default with no flow conditions
405-
QbaseUniform = [1.0_RP,0.0_RP,0.0_RP,0.0_RP,1.0_RP/(dimensionless % gammaM2),1.0_RP/POW2(dimensionless%Mach)]
406-
call self % mesh % SetUniformBaseFlow(QbaseUniform)
401+
call self % mesh % InitializeBaseFlow(controlVariables)
407402
call self % mesh % ProlongBaseSolutionToFaces(NCONSB)
408403
#ifdef _HAS_MPI_
409404
!$omp single

Solver/src/libs/io/SolutionFile.f90

Lines changed: 9 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -39,7 +39,8 @@ module SolutionFile
3939
#endif
4040

4141
private
42-
public :: MESH_FILE, SOLUTION_FILE, SOLUTION_AND_GRADIENTS_FILE, STATS_FILE, ZONE_MESH_FILE
42+
public :: MESH_FILE, SOLUTION_FILE, SOLUTION_AND_GRADIENTS_FILE, ZONE_MESH_FILE
43+
public :: STATS_FILE, STATS_AND_GRADIENTS_FILE
4344
public :: ZONE_SOLUTION_FILE, ZONE_SOLUTION_AND_DOT_FILE
4445
public :: SOLUTION_AND_SENSOR_FILE, SOLUTION_AND_GRADIENTS_AND_SENSOR_FILE
4546
public :: BEGINNING_DATA
@@ -61,11 +62,12 @@ module SolutionFile
6162
integer, parameter :: SOLUTION_FILE = 2
6263
integer, parameter :: SOLUTION_AND_GRADIENTS_FILE = 3
6364
integer, parameter :: STATS_FILE = 4
64-
integer, parameter :: ZONE_MESH_FILE = 5
65-
integer, parameter :: ZONE_SOLUTION_FILE = 6
66-
integer, parameter :: ZONE_SOLUTION_AND_DOT_FILE = 7
67-
integer, parameter :: SOLUTION_AND_SENSOR_FILE = 8
68-
integer, parameter :: SOLUTION_AND_GRADIENTS_AND_SENSOR_FILE = 9
65+
integer, parameter :: STATS_AND_GRADIENTS_FILE = 5
66+
integer, parameter :: ZONE_MESH_FILE = 6
67+
integer, parameter :: ZONE_SOLUTION_FILE = 7
68+
integer, parameter :: ZONE_SOLUTION_AND_DOT_FILE = 8
69+
integer, parameter :: SOLUTION_AND_SENSOR_FILE = 9
70+
integer, parameter :: SOLUTION_AND_GRADIENTS_AND_SENSOR_FILE = 10
6971

7072
integer, parameter :: SOLFILE_STR_LEN = 128
7173
integer, parameter :: END_OF_FILE = 99
@@ -133,6 +135,7 @@ subroutine CreateNewSolutionFile(name, type_, nodes, no_of_elements, iter, time,
133135
case(SOLUTION_FILE)
134136
case(SOLUTION_AND_GRADIENTS_FILE)
135137
case(STATS_FILE)
138+
case(STATS_AND_GRADIENTS_FILE)
136139
case(ZONE_MESH_FILE)
137140
case(ZONE_SOLUTION_FILE)
138141
case(ZONE_SOLUTION_AND_DOT_FILE)

Solver/src/libs/mesh/HexMesh.f90

Lines changed: 184 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -138,8 +138,9 @@ MODULE HexMeshClass
138138
procedure :: LoadSolutionForRestart => HexMesh_LoadSolutionForRestart
139139
procedure :: WriteCoordFile
140140
#if defined(ACOUSTIC)
141+
procedure :: InitializeBaseFlow => HexMesh_InitializeBaseFlow
141142
procedure :: SetUniformBaseFlow => HexMesh_SetUniformBaseFlow
142-
! procedure :: LoadBaseFlowSolution => HexMesh_LoadBaseFlowSolution
143+
procedure :: LoadBaseFlowSolution => HexMesh_LoadBaseFlowSolution
143144
procedure :: ProlongBaseSolutionToFaces => HexMesh_ProlongBaseSolutionToFaces
144145
procedure :: UpdateMPIFacesBaseSolution => HexMesh_UpdateMPIFacesBaseSolution
145146
procedure :: GatherMPIFacesBaseSolution => HexMesh_GatherMPIFacesBaseSolution
@@ -3554,7 +3555,7 @@ subroutine HexMesh_SaveStatistics(self, iter, time, name, saveGradients)
35543555
! Local variables
35553556
! ---------------
35563557
!
3557-
integer :: fid, eID
3558+
integer :: fid, eID, fileType
35583559
integer :: no_stat_s
35593560
integer(kind=AddrInt) :: pos
35603561
real(kind=RP) :: refs(NO_OF_SAVED_REFS)
@@ -3568,11 +3569,13 @@ subroutine HexMesh_SaveStatistics(self, iter, time, name, saveGradients)
35683569
refs(V_REF) = refValues % V
35693570
refs(T_REF) = refValues % T
35703571
refs(MACH_REF) = dimensionless % Mach
3571-
refs(RE_REF) = dimensionless % Re
3572+
! refs(RE_REF) = dimensionless % Re
35723573

35733574
! Create new file
35743575
! ---------------
3575-
call CreateNewSolutionFile(trim(name),STATS_FILE, self % nodeType, self % no_of_allElements, iter, time, refs)
3576+
fileType = STATS_FILE
3577+
if ( saveGradients .and. computeGradients ) fileType = STATS_AND_GRADIENTS_FILE
3578+
call CreateNewSolutionFile(trim(name), fileType, self % nodeType, self % no_of_allElements, iter, time, refs)
35763579
!
35773580
! Write arrays
35783581
! ------------
@@ -3583,7 +3586,6 @@ subroutine HexMesh_SaveStatistics(self, iter, time, name, saveGradients)
35833586
no_stat_s = 9
35843587
call writeArray(fid, e % storage % stats % data(1:no_stat_s,:,:,:), position=pos)
35853588
allocate(Q(NCONS, 0:e % Nxyz(1), 0:e % Nxyz(2), 0:e % Nxyz(3)))
3586-
! write(fid) e%storage%stats%data(7:,:,:,:)
35873589
Q(1:NCONS,:,:,:) = e % storage % stats % data(no_stat_s+1:no_stat_s+NCONS,:,:,:)
35883590
write(fid) Q
35893591
deallocate(Q)
@@ -4038,6 +4040,99 @@ END SUBROUTINE HexMesh_LoadSolution
40384040
!///////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
40394041
!
40404042
#if defined(ACOUSTIC)
4043+
4044+
subroutine HexMesh_InitializeBaseFlow(self, controlVariables)
4045+
use FTObjectClass
4046+
use mainKeywordsModule
4047+
use FileReadingUtilities, only: getRealArrayFromString
4048+
Implicit None
4049+
CLASS(HexMesh) :: self
4050+
class(FTValueDictionary) :: controlVariables
4051+
4052+
4053+
class(FTObject), pointer :: obj, obj2
4054+
character(len=LINE_LENGTH) :: qBaseMode
4055+
character(len=LINE_LENGTH) :: fileName
4056+
real(kind=RP) :: QbaseUniform(1:NCONSB)
4057+
4058+
CHARACTER(LEN=KEYWORD_LENGTH) :: qBaseKey = "qBase"
4059+
CHARACTER(LEN=KEYWORD_LENGTH) :: qBaseFileNameKey = "qBase file name"
4060+
CHARACTER(LEN=KEYWORD_LENGTH) :: qBaseVectorKey = "qBase vector"
4061+
character(len=LINE_LENGTH) :: qBaseByFile = 'file'
4062+
character(len=LINE_LENGTH) :: qBaseByUniformField = 'uniform'
4063+
4064+
! In the control file, the keyword 'qbase' is mandatory:
4065+
! qbase = file/uniform
4066+
! * When qbase is given by a file, the keyword 'qbase file name' is mandatory:
4067+
! qbase file name = path/to/file.stats.hsol
4068+
! * When qbase is given by a uniform field, the keyword 'qbase vector' is mandatory:
4069+
! qbase vector = [1.0_RP,0.0_RP,0.0_RP,0.0_RP,1.0_RP]
4070+
4071+
4072+
! Check which type of qBase we have: file or uniform
4073+
call toLower(qBaseKey)
4074+
obj => controlVariables % objectForKey(trim(qBaseKey))
4075+
if ( associated(obj) ) then
4076+
4077+
qBaseMode = controlVariables % stringValueForKey(trim(qBaseKey), requestedLength = LINE_LENGTH)
4078+
call ToLower(qBaseMode)
4079+
4080+
if ( trim(qBaseMode) .eq. trim(qBaseByFile) ) then
4081+
!
4082+
! Read Qbase from file
4083+
!
4084+
! Check that the user has specified the file to read from
4085+
call toLower(qBaseFileNameKey)
4086+
obj2 => controlVariables % objectForKey(trim(qBaseFileNameKey))
4087+
if ( .not. associated(obj2) ) then
4088+
print *, trim(qBaseKey), " = ", trim(qBaseMode), ", but no file specified. Use:"
4089+
print *, trim(qBaseFileNameKey), " = filename"
4090+
errorMessage(STD_OUT)
4091+
error stop
4092+
end if
4093+
! Load the base flow from the specified stats file
4094+
fileName = controlVariables % stringValueForKey(qBaseFileNameKey,requestedLength = LINE_LENGTH)
4095+
call self % LoadBaseFlowSolution(fileName)
4096+
elseif ( trim(qBaseMode) .eq. trim(qbaseByUniformField) ) then
4097+
!
4098+
! Read Qbase uniform field from control file
4099+
!
4100+
! Check that the user has specified the field
4101+
call toLower(qBaseVectorKey)
4102+
obj2 => controlVariables % objectForKey(trim(qBaseVectorKey))
4103+
if ( .not. associated(obj2) ) then
4104+
print *, trim(qBaseKey), " = ", trim(qBaseMode), ", but no vector specified. Use:"
4105+
print *, trim(qBaseVectorKey), " = [1.0_RP,0.0_RP,0.0_RP,0.0_RP,1.0_RP]"
4106+
errorMessage(STD_OUT)
4107+
error stop
4108+
end if
4109+
! Read the field
4110+
QbaseUniform = 0.0_RP
4111+
QbaseUniform = GetRealArrayFromString( controlVariables % StringValueForKey(qBaseVectorKey,requestedLength = LINE_LENGTH))
4112+
! Set uniform field
4113+
call self % SetUniformBaseFlow(QbaseUniform)
4114+
else
4115+
print*, 'Unknown qBase mode "',trim(qBaseMode),'".'
4116+
print*, "Implemented modes are:"
4117+
print*, " * ", trim(qBaseByFile)
4118+
print*, " * ", trim(qbaseByUniformField)
4119+
errorMessage(STD_OUT)
4120+
error stop
4121+
end if
4122+
4123+
else
4124+
! Keyword not present:
4125+
print*, 'Argument qBase mode is mandatory'
4126+
print*, "Implemented modes are:"
4127+
print*, " * ", trim(qBaseByFile)
4128+
print*, " * ", trim(qbaseByUniformField)
4129+
errorMessage(STD_OUT)
4130+
error stop
4131+
4132+
end if
4133+
4134+
end subroutine HexMesh_InitializeBaseFlow
4135+
40414136
Subroutine HexMesh_SetUniformBaseFlow(self,Q_in)
40424137
Implicit None
40434138
CLASS(HexMesh) :: self
@@ -4056,6 +4151,90 @@ Subroutine HexMesh_SetUniformBaseFlow(self,Q_in)
40564151
!
40574152
End Subroutine HexMesh_SetUniformBaseFlow
40584153
!
4154+
!////////////////////////////////////////////////////////////////////////
4155+
!
4156+
Subroutine HexMesh_LoadBaseFlowSolution(self, fileName)
4157+
use VariableConversion_CAA, only: PressureBaseFlow
4158+
Implicit None
4159+
CLASS(HexMesh) :: self
4160+
character(len=*) :: fileName
4161+
! ---------------
4162+
! Local variables
4163+
! ---------------
4164+
INTEGER :: fID, eID, fileType, no_of_elements, nodetype
4165+
integer :: i, j, k
4166+
integer(kind=AddrInt) :: pos
4167+
integer :: no_stat_s, no_stats_read
4168+
real(kind=RP), allocatable :: Q(:,:,:,:)
4169+
4170+
!
4171+
! Get the file type
4172+
! -----------------
4173+
fileType = getSolutionFileType(trim(fileName))
4174+
4175+
if ( (fileType .ne. STATS_FILE) .and. (fileType .ne. STATS_AND_GRADIENTS_FILE) ) then
4176+
print*, "The selected file is not a statistics file"
4177+
errorMessage(STD_OUT)
4178+
error stop
4179+
end if
4180+
!
4181+
! Get the node type
4182+
! -----------------
4183+
nodeType = getSolutionFileNodeType(trim(fileName))
4184+
4185+
if ( nodeType .ne. self % nodeType ) then
4186+
print*, "WARNING: Solution file uses a different discretization nodes than the mesh."
4187+
print*, "Add restart polorder = (Pol order in your restart file) in the control file if you want interpolation routines to be used."
4188+
print*, "If restart polorder is not specified the values in the original set of nodes are loaded into the new nodes without interpolation."
4189+
errorMessage(STD_OUT)
4190+
end if
4191+
!
4192+
! Read the number of elements
4193+
! ---------------------------
4194+
no_of_elements = getSolutionFileNoOfElements(trim(fileName))
4195+
if ( no_of_elements .ne. self % no_of_allElements ) then
4196+
write(STD_OUT,'(A,A)') "The number of elements stored in the restart file ", &
4197+
"do not match that of the mesh file"
4198+
errorMessage(STD_OUT)
4199+
error stop
4200+
end if
4201+
4202+
!
4203+
! Read elements data
4204+
! ------------------
4205+
fID = putSolutionFileInReadDataMode(trim(fileName))
4206+
no_stat_s = 9
4207+
no_stats_read = no_stat_s + NCONS
4208+
if (fileType .eq. STATS_AND_GRADIENTS_FILE) no_stats_read = no_stats_read + NGRAD*NDIM
4209+
do eID = 1, size(self % elements)
4210+
associate( e => self % elements(eID) )
4211+
pos = POS_INIT_DATA + (e % globID-1)*5_AddrInt*SIZEOF_INT + 1_AddrInt*no_stats_read*e % offsetIO*SIZEOF_RP
4212+
pos = pos + 5_AddrInt*SIZEOF_INT ! This is to skip the reading of the dimensions and shape in writeArray
4213+
4214+
! Read and initialize velocity
4215+
allocate(Q(1:no_stat_s, 0:e % Nxyz(1), 0:e % Nxyz(2), 0:e % Nxyz(3)))
4216+
read(fID, pos=pos) Q
4217+
! WARNING: 1,2,3 should match U,V,W in StatisticsMonitor
4218+
self % elements(eID) % storage % Qbase(ICAAU,:,:,:) = Q(1,:,:,:)
4219+
self % elements(eID) % storage % Qbase(ICAAV,:,:,:) = Q(2,:,:,:)
4220+
self % elements(eID) % storage % Qbase(ICAAW,:,:,:) = Q(3,:,:,:)
4221+
deallocate(Q)
4222+
! Read NCONS variables
4223+
allocate(Q(1:NCONS, 0:e % Nxyz(1), 0:e % Nxyz(2), 0:e % Nxyz(3)))
4224+
read(fID) Q
4225+
! Initialize density
4226+
self % elements(eID) % storage % Qbase(ICAARHO,:,:,:) = Q(IRHO,:,:,:)
4227+
! Read and initialize pressure
4228+
do k = 0, e % Nxyz(3) ; do j = 0, e % Nxyz(2) ; do i = 0, e % Nxyz(1)
4229+
self % elements(eID) % storage % Qbase(ICAAP,i,j,k) = PressureBaseFlow(Q(:,i,j,k))
4230+
end do ; end do ; end do
4231+
deallocate(Q)
4232+
end associate
4233+
end do
4234+
4235+
4236+
End Subroutine HexMesh_LoadBaseFlowSolution
4237+
40594238
!////////////////////////////////////////////////////////////////////////
40604239
!
40614240
subroutine HexMesh_ProlongBaseSolutionToFaces(self, nEqn)

Solver/src/libs/physics/acoustic/VariableConversion_caa.f90

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -6,7 +6,7 @@ module VariableConversion_CAA
66
implicit none
77

88
private
9-
public Pressure, PressureDot
9+
public Pressure, PressureDot, PressureBaseFlow
1010
public NSGradientVariables_STATE
1111
! public getPrimitiveVariables
1212
public getVelocityGradients, getTemperatureGradient, getConservativeGradients

0 commit comments

Comments
 (0)