From b42b46dd749f22a7317bba70c0013d437b619143 Mon Sep 17 00:00:00 2001 From: Semih Akkurt Date: Fri, 20 Dec 2024 19:49:57 +0000 Subject: [PATCH] Add channel case. --- examples/channel/input.x3d | 31 +++++++++++++++ src/CMakeLists.txt | 1 + src/case/channel.f90 | 77 ++++++++++++++++++++++++++++++++++++++ src/xcompact.f90 | 4 ++ 4 files changed, 113 insertions(+) create mode 100644 examples/channel/input.x3d create mode 100644 src/case/channel.f90 diff --git a/examples/channel/input.x3d b/examples/channel/input.x3d new file mode 100644 index 00000000..322c6a64 --- /dev/null +++ b/examples/channel/input.x3d @@ -0,0 +1,31 @@ +&domain_settings +! Flow case +flow_case_name = 'channel' + +! Global domain length +L_global = 8d0, 2d0, 4d0 + +! Global number of cells in each direction +dims_global = 256, 257, 128 + +! Domain decomposition in each direction +nproc_dir = 1, 1, 2 + +! BC options are 'periodic' | 'neumann' | 'dirichlet' +BC_x = 'periodic', 'periodic' +BC_y = 'dirichlet', 'dirichlet' +BC_z = 'periodic', 'periodic' +/End + +&solver_params +Re = 4200d0 +time_intg = 'AB3' ! 'AB[1-4]' | 'RK[1-4]' +dt = 0.005d0 +n_iters = 100000 +n_output = 1000 +poisson_solver_type = 'FFT' ! 'FFT' | 'CG' +der1st_scheme = 'compact6' +der2nd_scheme = 'compact6' ! 'compact6' | 'compact6-hyperviscous' +interpl_scheme = 'optimised' ! 'classic' | 'optimised' | 'aggressive' +stagder_scheme = 'compact6' +/End diff --git a/src/CMakeLists.txt b/src/CMakeLists.txt index 8db0c706..d1b44201 100644 --- a/src/CMakeLists.txt +++ b/src/CMakeLists.txt @@ -11,6 +11,7 @@ set(SRC time_integrator.f90 vector_calculus.f90 case/base_case.f90 + case/channel.f90 case/generic.f90 case/tgv.f90 omp/backend.f90 diff --git a/src/case/channel.f90 b/src/case/channel.f90 new file mode 100644 index 00000000..3e140797 --- /dev/null +++ b/src/case/channel.f90 @@ -0,0 +1,77 @@ +module m_case_channel + use iso_fortran_env, only: stderr => error_unit + + use m_allocator, only: allocator_t, field_t + use m_base_backend, only: base_backend_t + use m_base_case, only: base_case_t + use m_common, only: dp + use m_mesh, only: mesh_t + use m_solver, only: init + + implicit none + + type, extends(base_case_t) :: case_channel_t + contains + procedure :: boundary_conditions => boundary_conditions_channel + procedure :: initial_conditions => initial_conditions_channel + procedure :: post_transeq => post_transeq_channel + procedure :: postprocess => postprocess_channel + end type case_channel_t + + interface case_channel_t + module procedure case_channel_init + end interface case_channel_t + +contains + + function case_channel_init(backend, mesh, host_allocator) result(flow_case) + implicit none + + class(base_backend_t), target, intent(inout) :: backend + type(mesh_t), target, intent(inout) :: mesh + type(allocator_t), target, intent(inout) :: host_allocator + type(case_channel_t) :: flow_case + + call flow_case%case_init(backend, mesh, host_allocator) + + end function case_channel_init + + subroutine boundary_conditions_channel(self) + implicit none + + class(case_channel_t) :: self + + end subroutine boundary_conditions_channel + + subroutine initial_conditions_channel(self) + implicit none + + class(case_channel_t) :: self + + call self%solver%u%fill(1._dp) + call self%solver%v%fill(0._dp) + call self%solver%w%fill(0._dp) + + end subroutine initial_conditions_channel + + subroutine postprocess_channel(self, t) + implicit none + + class(case_channel_t) :: self + real(dp), intent(in) :: t + + if (self%solver%mesh%par%is_root()) print *, 'time =', t + call self%print_enstrophy(self%solver%u, self%solver%v, self%solver%w) + call self%print_div_max_mean(self%solver%u, self%solver%v, self%solver%w) + + end subroutine postprocess_channel + + subroutine post_transeq_channel(self, du, dv, dw) + implicit none + + class(case_channel_t) :: self + class(field_t), intent(inout) :: du, dv, dw + + end subroutine post_transeq_channel + +end module m_case_channel diff --git a/src/xcompact.f90 b/src/xcompact.f90 index da9ac4aa..30138772 100644 --- a/src/xcompact.f90 +++ b/src/xcompact.f90 @@ -8,6 +8,7 @@ program xcompact use m_solver, only: solver_t use m_tdsops, only: tdsops_t use m_mesh + use m_case_channel, only: case_channel_t use m_case_generic, only: case_generic_t use m_case_tgv, only: case_tgv_t @@ -108,6 +109,9 @@ program xcompact if (nrank == 0) print *, 'Flow case: ', flow_case_name select case (trim(flow_case_name)) + case ('channel') + allocate (case_channel_t :: flow_case) + flow_case = case_channel_t(backend, mesh, host_allocator) case ('generic') allocate (case_generic_t :: flow_case) flow_case = case_generic_t(backend, mesh, host_allocator)