Skip to content

Commit

Permalink
Gfortran13 compatibility fixes
Browse files Browse the repository at this point in the history
  • Loading branch information
Esteban Ferrer committed Jan 10, 2025
1 parent c57446f commit 6e6aa11
Show file tree
Hide file tree
Showing 7 changed files with 62 additions and 32 deletions.
41 changes: 36 additions & 5 deletions Solver/src/libs/ftobject/FTDictionaryClass.f90
Original file line number Diff line number Diff line change
Expand Up @@ -441,15 +441,46 @@ FUNCTION objectForKey(self,key)
CHARACTER(LEN=*) :: key
CLASS(FTObject) , POINTER :: objectForKey
INTEGER :: h

INTEGER, EXTERNAL :: b3hs_hash_key_jenkins

CLASS(FTLinkedListRecord) , POINTER :: listRecordPtr => NULL()
CHARACTER(LEN=FTDICT_KWD_STRING_LENGTH) :: keyString
INTEGER, EXTERNAL :: b3hs_hash_key_jenkins

objectForKey => NULL()
IF(self % COUNT() == 0) RETURN

!

! -------------
! Find the hash
! -------------
!
h = b3hs_hash_key_jenkins(key,SIZE(self % entries))

IF ( self % entries(h) % COUNT() > 0 ) THEN
objectForKey => objectForKeyInList(key,self % entries(h))
!
! -----------------------
! Search through the list
! -----------------------
!
listRecordPtr => self % entries(h) % head
DO WHILE(ASSOCIATED(listRecordPtr))
!
! --------------------------------------------
! The list's recordObject is a FTKeyObjectPair
! --------------------------------------------
!
SELECT TYPE (pair => listRecordPtr % recordObject)
TYPE is (FTKeyObjectPair)
keyString = pair % keyString
IF ( TRIM(keyString) == TRIM(key) .AND. &
ASSOCIATED( pair % valueObject) ) THEN
objectForKey => pair % valueObject
EXIT
END IF
CLASS DEFAULT
END SELECT
listRecordPtr => listRecordPtr % next
END DO

END IF

END FUNCTION ObjectForKey
Expand Down
24 changes: 10 additions & 14 deletions Solver/src/libs/mesh/FaceClass.f90
Original file line number Diff line number Diff line change
Expand Up @@ -875,33 +875,31 @@ subroutine Face_ProjectGradientFluxToElements(self, nEqn, Hflux, whichElements,
do side = 1, 2
select case ( whichElements(side) )
case (1) ! Prolong from left element
associate(unStar => self % storage(1) % unStar)
select case ( self % projectionType(1) )
case (0)
unStar(:,:,:,:) = Hflux
self % storage(1) % unStar(:,:,:,:) = Hflux
case (1)
unStar(:,:,:,:) = 0.0
self % storage(1) % unStar(:,:,:,:) = 0.0
do j = 0, self % NelLeft(2) ; do l = 0, self % Nf(1) ; do i = 0, self % NelLeft(1)
unStar(:,:,i,j) = unStar(:,:,i,j) + Tset(self % Nf(1), self % NfLeft(1)) % T(i,l) * Hflux(:,:,l,j)
self % storage(1) % unStar(:,:,i,j) = self % storage(1) % unStar(:,:,i,j) + Tset(self % Nf(1), self % NfLeft(1)) % T(i,l) * Hflux(:,:,l,j)
end do ; end do ; end do

case (2)
unStar(:,:,:,:) = 0.0
self % storage(1) % unStar(:,:,:,:) = 0.0
do l = 0, self % Nf(2) ; do j = 0, self % NelLeft(2) ; do i = 0, self % NelLeft(1)
unStar(:,:,i,j) = unStar(:,:,i,j) + Tset(self % Nf(2), self % NfLeft(2)) % T(j,l) * Hflux(:,:,i,l)
self % storage(1) % unStar(:,:,i,j) = self % storage(1) % unStar(:,:,i,j) + Tset(self % Nf(2), self % NfLeft(2)) % T(j,l) * Hflux(:,:,i,l)
end do ; end do ; end do

case (3)
unStar(:,:,:,:) = 0.0
self % storage(1) % unStar(:,:,:,:) = 0.0
do l = 0, self % Nf(2) ; do j = 0, self % NfLeft(2)
do m = 0, self % Nf(1) ; do i = 0, self % NfLeft(1)
unStar(:,:,i,j) = unStar(:,:,i,j) + Tset(self % Nf(1), self % NfLeft(1)) % T(i,m) &
self % storage(1) % unStar(:,:,i,j) = self % storage(1) % unStar(:,:,i,j) + Tset(self % Nf(1), self % NfLeft(1)) % T(i,m) &
* Tset(self % Nf(2), self % NfLeft(2)) % T(j,l) &
* Hflux(:,:,m,l)
end do ; end do
end do ; end do
end select
end associate

case (2) ! Prolong from right element
!
Expand Down Expand Up @@ -939,19 +937,17 @@ subroutine Face_ProjectGradientFluxToElements(self, nEqn, Hflux, whichElements,
! 2nd stage: Rotation
! *********
!
associate(unStar => self % storage(2) % unStar)
do j = 0, self % NfRight(2) ; do i = 0, self % NfRight(1)
call leftIndexes2Right(i,j,self % NfRight(1), self % NfRight(2), self % rotation, ii, jj)
unStar(:,:,ii,jj) = HstarAux(:,:,i,j)
self % storage(2) % unStar(:,:,ii,jj) = HstarAux(:,:,i,j)
end do ; end do
!
! *********
! 3rd stage: Multiplication by a factor (inversion usually)
! *********
!
unStar = factor * unStar
self % storage(2) % unStar = factor * self % storage(2) % unStar

end associate
end select
end do

Expand Down Expand Up @@ -1225,4 +1221,4 @@ elemental subroutine Face_Assign(to,from)
to % geom = from % geom
to % storage = from % storage
end subroutine Face_Assign
end Module FaceClass
end Module FaceClass
2 changes: 1 addition & 1 deletion Solver/src/libs/mesh/HexMesh.f90
Original file line number Diff line number Diff line change
Expand Up @@ -67,7 +67,7 @@ MODULE HexMeshClass
type(Element), dimension(:), allocatable :: elements
type(MPI_FacesSet_t) :: MPIfaces
type(IBM_type) :: IBM
class(Zone_t), dimension(:), allocatable :: zones
type(Zone_t), dimension(:), allocatable :: zones
logical :: child = .FALSE. ! Is this a (multigrid) child mesh? default .FALSE.
logical :: meshIs2D = .FALSE. ! Is this a 2D mesh? default .FALSE.
integer :: dir2D = 0 ! If it is in fact a 2D mesh, dir 2D stores the global direction IX, IY or IZ
Expand Down
10 changes: 5 additions & 5 deletions Solver/src/libs/mesh/ZoneClass.f90
Original file line number Diff line number Diff line change
Expand Up @@ -36,8 +36,8 @@ module ZoneClass
! ------------------------------------------
subroutine ConstructZones( faces , zones )
implicit none
class(Face), target :: faces(:)
class(Zone_t), allocatable :: zones(:)
type(Face), target :: faces(:)
type(Zone_t), allocatable :: zones(:)
!
! ---------------
! Local variables
Expand Down Expand Up @@ -82,8 +82,8 @@ end subroutine ConstructZones
! ------------------------------------------------------------------------------------------
subroutine ReassignZones( faces , zones )
implicit none
class(Face), target :: faces(:)
class(Zone_t) :: zones(:)
type(Face), target :: faces(:)
type(Zone_t) :: zones(:)
!
! ---------------
! Local variables
Expand Down Expand Up @@ -196,7 +196,7 @@ end subroutine Zone_AssignFaces
function AllZoneNames(no_of_zones, zones)
implicit none
integer, intent(in) :: no_of_zones
class(Zone_t), intent(in) :: zones(no_of_zones)
type(Zone_t), intent(in) :: zones(no_of_zones)
character(len=LINE_LENGTH) :: AllZoneNames(no_of_zones)
!
! ---------------
Expand Down
5 changes: 4 additions & 1 deletion Solver/src/libs/ml/Clustering.f90
Original file line number Diff line number Diff line change
Expand Up @@ -775,8 +775,11 @@ subroutine GMM_final(self)
nullify(self % g % logdet)

! Required for gfortran?
! The commented out line makes gfortran-13 crash. Either the code shouldnt reach this point at all
! meaning that the initialization of the object is wrong, or there is a bug in gfortran with
! checking the associated status of a pointer.
if (self % initialized) then
if (associated(self % g % storage)) deallocate(self % g % storage)
!if (associated(self % g % storage)) deallocate(self % g % storage)
if (allocated(self % g % tcov)) deallocate(self % g % tcov)
if (allocated(self % prob)) deallocate(self % prob)
end if
Expand Down
10 changes: 4 additions & 6 deletions Solver/src/libs/physics/incns/PhysicsStorage_iNS.f90
Original file line number Diff line number Diff line change
Expand Up @@ -258,6 +258,8 @@ SUBROUTINE ConstructPhysicsStorage_iNS( controlVariables, Lref, timeref, success
array = getRealArrayFromString( controlVariables % StringValueForKey(GRAVITY_DIRECTION_KEY,&
KEYWORD_LENGTH))
dimensionless_ % gravity_dir = array(1:3) / norm2(array(1:3))
else
dimensionless_ % gravity_dir = (/0.0_RP,0.0_RP,0.0_RP/)
end if

refValues_ % g0 = controlVariables % DoublePrecisionValueForKey(GRAVITY_ACCELERATION_KEY)
Expand Down Expand Up @@ -453,19 +455,15 @@ SUBROUTINE CheckPhysics_iNSInputIntegrity( controlVariables, success )

end if

else
if ( controlVariables % ContainsKey(GRAVITY_ACCELERATION_KEY) ) then
elseif ( controlVariables % ContainsKey(GRAVITY_ACCELERATION_KEY) ) then
print*, "Gravity acceleration requires gravity direction."
print*, "Specify gravity direction with:"
print*, " ", GRAVITY_DIRECTION_KEY, " = [x,y,z]"
errorMessage(STD_OUT)
error stop

else
else
call controlVariables % AddValueForKey("0.0d0", GRAVITY_ACCELERATION_KEY)

end if

end if


Expand Down
2 changes: 2 additions & 0 deletions Solver/src/libs/physics/multiphase/PhysicsStorage_MU.f90
Original file line number Diff line number Diff line change
Expand Up @@ -247,6 +247,8 @@ SUBROUTINE ConstructPhysicsStorage_MU( controlVariables, Lref, timeref, success
else
dimensionless_ % gravity_dir = 0.0_RP
end if
else
dimensionless_ % gravity_dir = 0.0_RP
end if

if ( almostEqual(abs(refValues_ % g0), 0.0_RP) ) then
Expand Down

0 comments on commit 6e6aa11

Please sign in to comment.