From b90afc04b61c052197cfc452c4be1914641f7628 Mon Sep 17 00:00:00 2001 From: dreamer2368 Date: Sat, 24 Dec 2022 17:45:25 -0600 Subject: [PATCH] StokesSecondWall levelset class. --- CMakeLists.txt | 2 + include/StokesSecondWallLevelset.f90 | 72 ++++++++++++++ src/RegionImpl.f90 | 4 + src/StokesSecondWallLevelsetImpl.f90 | 138 +++++++++++++++++++++++++++ 4 files changed, 216 insertions(+) create mode 100644 include/StokesSecondWallLevelset.f90 create mode 100644 src/StokesSecondWallLevelsetImpl.f90 diff --git a/CMakeLists.txt b/CMakeLists.txt index 37b2906..a1a1435 100644 --- a/CMakeLists.txt +++ b/CMakeLists.txt @@ -216,6 +216,8 @@ set(magudiObj_SOURCES ${PROJECT_SOURCE_DIR}/src/ImmersedBoundaryImpl.f90 ${PROJECT_SOURCE_DIR}/include/SinusoidalWallLevelset.f90 ${PROJECT_SOURCE_DIR}/src/SinusoidalWallLevelsetImpl.f90 + ${PROJECT_SOURCE_DIR}/include/StokesSecondWallLevelset.f90 + ${PROJECT_SOURCE_DIR}/src/StokesSecondWallLevelsetImpl.f90 # C source files. ${PROJECT_SOURCE_DIR}/src/PLOT3DFormat.c diff --git a/include/StokesSecondWallLevelset.f90 b/include/StokesSecondWallLevelset.f90 new file mode 100644 index 0000000..68f1b7b --- /dev/null +++ b/include/StokesSecondWallLevelset.f90 @@ -0,0 +1,72 @@ +#include "config.h" + +module StokesSecondWallLevelset_mod + + use LevelsetFactory_mod, only : t_LevelsetFactory + + implicit none + + type, extends(t_LevelsetFactory), public :: t_StokesSecondWallLevelset + + real(SCALAR_KIND), allocatable :: levelsetLoc(:), & ! this point lies on level set plane. + levelsetNormal(:), & ! normal direction of the plane. + levelsetDirection(:), & ! oscillation direction. + levelsetAmp, & ! amplitude of oscillation. + levelsetPeriod ! period of oscillation. + + contains + + procedure, pass :: setup => setupStokesSecondWallLevelset + procedure, pass :: cleanup => cleanupStokesSecondWallLevelset + procedure, pass :: updateLevelset => updateStokesSecondWallLevelset + + end type t_StokesSecondWallLevelset + + interface + + subroutine setupStokesSecondWallLevelset(this, grids, states) + + use Grid_mod, only : t_Grid + use State_mod, only : t_State + + import :: t_StokesSecondWallLevelset + + class(t_StokesSecondWallLevelset) :: this + class(t_Grid), intent(in) :: grids(:) + class(t_State) :: states(:) + + end subroutine setupStokesSecondWallLevelset + + end interface + + interface + + subroutine cleanupStokesSecondWallLevelset(this) + + import :: t_StokesSecondWallLevelset + + class(t_StokesSecondWallLevelset) :: this + + end subroutine cleanupStokesSecondWallLevelset + + end interface + + interface + + subroutine updateStokesSecondWallLevelset(this, mode, grids, states) + + use Grid_mod, only : t_Grid + use State_mod, only : t_State + + import :: t_StokesSecondWallLevelset + + class(t_StokesSecondWallLevelset) :: this + integer, intent(in) :: mode + class(t_Grid), intent(in) :: grids(:) + class(t_State) :: states(:) + + end subroutine updateStokesSecondWallLevelset + + end interface + +end module StokesSecondWallLevelset_mod diff --git a/src/RegionImpl.f90 b/src/RegionImpl.f90 index 7fb6952..4c645e2 100644 --- a/src/RegionImpl.f90 +++ b/src/RegionImpl.f90 @@ -2055,6 +2055,7 @@ subroutine connectLevelsetFactory(this) use ImmersedBoundaryPatch_mod, only : t_ImmersedBoundaryPatch use LevelsetFactory_mod, only : t_LevelsetFactory use SinusoidalWallLevelset_mod, only : t_SinusoidalWallLevelset + use StokesSecondWallLevelset_mod, only : t_StokesSecondWallLevelset ! <<< Internal modules >>> use InputHelper, only : getRequiredOption @@ -2100,6 +2101,9 @@ subroutine connectLevelsetFactory(this) case ('sinusoidal_wall') allocate(t_SinusoidalWallLevelset :: this%levelsetFactory) + case ('stokes_second_wall') + allocate(t_StokesSecondWallLevelset :: this%levelsetFactory) + case default this%levelsetType = "" diff --git a/src/StokesSecondWallLevelsetImpl.f90 b/src/StokesSecondWallLevelsetImpl.f90 new file mode 100644 index 0000000..8f0aed5 --- /dev/null +++ b/src/StokesSecondWallLevelsetImpl.f90 @@ -0,0 +1,138 @@ +#include "config.h" + +subroutine setupStokesSecondWallLevelset(this, grids, states) + + ! <<< Derived types >>> + use StokesSecondWallLevelset_mod, only : t_StokesSecondWallLevelset + use Grid_mod, only : t_Grid + use State_mod, only : t_State + + ! <<< Internal modules >>> + use InputHelper, only : getRequiredOption + + implicit none + + ! <<< Arguments >>> + class(t_StokesSecondWallLevelset) :: this + class(t_Grid), intent(in) :: grids(:) + class(t_State) :: states(:) + + ! <<< Local variables >>> + integer, parameter :: wp = SCALAR_KIND + real(wp), parameter :: pi = 4.0_wp * atan(1.0_wp) + integer :: i + real(wp) :: buf + character(len = STRING_LENGTH) :: key + + call this%cleanup() + + assert(size(states) == size(grids)) + allocate(this%levelsetLoc(grids(1)%nDimensions)) + allocate(this%levelsetNormal(grids(1)%nDimensions)) + allocate(this%levelsetDirection(grids(1)%nDimensions)) + + do i = 1, grids(1)%nDimensions + write(key, '(A,I1.1)') "immersed_boundary/location", i + call getRequiredOption(key, this%levelsetLoc(i), grids(1)%comm) + + write(key, '(A,I1.1)') "immersed_boundary/normal", i + call getRequiredOption(key, this%levelsetNormal(i), grids(1)%comm) + + write(key, '(A,I1.1)') "immersed_boundary/oscillation_direction", i + call getRequiredOption(key, this%levelsetDirection(i), grids(1)%comm) + end do + + ! normalize normal. + buf = 0.0_wp + do i = 1, grids(1)%nDimensions + buf = buf + this%levelsetNormal(i) * this%levelsetNormal(i) + end do + this%levelsetNormal = this%levelsetNormal / sqrt(buf) + + ! normalize direction. + buf = 0.0_wp + do i = 1, grids(1)%nDimensions + buf = buf + this%levelsetDirection(i) * this%levelsetDirection(i) + end do + this%levelsetDirection = this%levelsetDirection / sqrt(buf) + + call getRequiredOption("immersed_boundary/amplitude", this%levelsetAmp, grids(1)%comm) + call getRequiredOption("immersed_boundary/period", this%levelsetPeriod, grids(1)%comm) + +end subroutine setupStokesSecondWallLevelset + +subroutine cleanupStokesSecondWallLevelset(this) + + ! <<< Derived types >>> + use StokesSecondWallLevelset_mod, only : t_StokesSecondWallLevelset + + ! <<< Arguments >>> + class(t_StokesSecondWallLevelset) :: this + + SAFE_DEALLOCATE(this%levelsetLoc) + SAFE_DEALLOCATE(this%levelsetNormal) + SAFE_DEALLOCATE(this%levelsetDirection) + +end subroutine cleanupStokesSecondWallLevelset + +subroutine updateStokesSecondWallLevelset(this, mode, grids, states) + + ! <<< Derived types >>> + use StokesSecondWallLevelset_mod, only : t_StokesSecondWallLevelset + use Grid_mod, only : t_Grid + use State_mod, only : t_State + + ! <<< Enumerations >>> + use Region_enum, only : FORWARD + + ! <<< Internal modules >>> + use MPITimingsHelper, only : startTiming, endTiming + + implicit none + + ! <<< Arguments >>> + class(t_StokesSecondWallLevelset) :: this + integer, intent(in) :: mode + class(t_Grid), intent(in) :: grids(:) + class(t_State) :: states(:) + + ! <<< Local variables >>> + integer, parameter :: wp = SCALAR_KIND + real(wp), parameter :: pi = 4.0_wp * atan(1.0_wp) + real(wp) :: timeFactor, timeDerivativeFactor + real(wp), dimension(grids(1)%nDimensions) :: loc, vel + integer :: i, j + + call startTiming("updateStokesSecondWallLevelset") + + if (mode .ne. FORWARD) then + return + end if + + timeFactor = cos(2.0_wp * pi * states(1)%time / this%levelsetPeriod) + timeDerivativeFactor = -2.0_wp * pi / this%levelsetPeriod & + * sin(2.0_wp * pi * states(1)%time / this%levelsetPeriod) + loc = this%levelsetLoc + timeFactor * this%levelsetAmp * this%levelsetDirection + vel = timeDerivativeFactor * this%levelsetAmp * this%levelsetDirection + + do i = 1, size(states) + !NOTE: you cannot simply cycle here, since even in the same grid, some + !processors may not have ibm patch. + !if (.not. states(i)%ibmPatchExists) cycle + + assert(size(states(i)%levelset, 1) == size(states(i)%conservedVariables, 1)) + assert(size(states(i)%levelset, 2) == 1) + + states(i)%levelset = 0.0_wp + do j = 1, grids(i)%nDimensions + states(i)%levelset(:, 1) = states(i)%levelset(:, 1) & + + this%levelsetNormal(j) * (grids(i)%coordinates(:, j) - loc(j)) + + states(i)%levelsetNormal(:, j) = this%levelsetNormal(j) + states(i)%objectVelocity(:, j) = vel(j) + end do + end do + + call endTiming("updateStokesSecondWallLevelset") + +end subroutine updateStokesSecondWallLevelset