tgv.f90 Source File


This file depends on

sourcefile~~tgv.f90~~EfferentGraph sourcefile~tgv.f90 tgv.f90 sourcefile~allocator.f90~2 allocator.f90 sourcefile~tgv.f90->sourcefile~allocator.f90~2 sourcefile~backend.f90 backend.f90 sourcefile~tgv.f90->sourcefile~backend.f90 sourcefile~base_case.f90 base_case.f90 sourcefile~tgv.f90->sourcefile~base_case.f90 sourcefile~common.f90~3 common.f90 sourcefile~tgv.f90->sourcefile~common.f90~3 sourcefile~field.f90 field.f90 sourcefile~tgv.f90->sourcefile~field.f90 sourcefile~mesh.f90 mesh.f90 sourcefile~tgv.f90->sourcefile~mesh.f90 sourcefile~solver.f90 solver.f90 sourcefile~tgv.f90->sourcefile~solver.f90 sourcefile~allocator.f90~2->sourcefile~common.f90~3 sourcefile~allocator.f90~2->sourcefile~field.f90 sourcefile~allocator.f90~2->sourcefile~mesh.f90 sourcefile~backend.f90->sourcefile~allocator.f90~2 sourcefile~backend.f90->sourcefile~common.f90~3 sourcefile~backend.f90->sourcefile~field.f90 sourcefile~backend.f90->sourcefile~mesh.f90 sourcefile~poisson_fft.f90 poisson_fft.f90 sourcefile~backend.f90->sourcefile~poisson_fft.f90 sourcefile~tdsops.f90 tdsops.f90 sourcefile~backend.f90->sourcefile~tdsops.f90 sourcefile~base_case.f90->sourcefile~allocator.f90~2 sourcefile~base_case.f90->sourcefile~backend.f90 sourcefile~base_case.f90->sourcefile~common.f90~3 sourcefile~base_case.f90->sourcefile~field.f90 sourcefile~base_case.f90->sourcefile~mesh.f90 sourcefile~base_case.f90->sourcefile~solver.f90 sourcefile~field.f90->sourcefile~common.f90~3 sourcefile~mesh.f90->sourcefile~common.f90~3 sourcefile~mesh.f90->sourcefile~field.f90 sourcefile~decomp_2decompfft.f90 decomp_2decompfft.f90 sourcefile~mesh.f90->sourcefile~decomp_2decompfft.f90 sourcefile~mesh_content.f90 mesh_content.f90 sourcefile~mesh.f90->sourcefile~mesh_content.f90 sourcefile~solver.f90->sourcefile~allocator.f90~2 sourcefile~solver.f90->sourcefile~backend.f90 sourcefile~solver.f90->sourcefile~common.f90~3 sourcefile~solver.f90->sourcefile~field.f90 sourcefile~solver.f90->sourcefile~mesh.f90 sourcefile~config.f90 config.f90 sourcefile~solver.f90->sourcefile~config.f90 sourcefile~solver.f90->sourcefile~tdsops.f90 sourcefile~time_integrator.f90 time_integrator.f90 sourcefile~solver.f90->sourcefile~time_integrator.f90 sourcefile~vector_calculus.f90 vector_calculus.f90 sourcefile~solver.f90->sourcefile~vector_calculus.f90 sourcefile~config.f90->sourcefile~common.f90~3 sourcefile~decomp_2decompfft.f90->sourcefile~mesh_content.f90 sourcefile~mesh_content.f90->sourcefile~common.f90~3 sourcefile~poisson_fft.f90->sourcefile~common.f90~3 sourcefile~poisson_fft.f90->sourcefile~field.f90 sourcefile~poisson_fft.f90->sourcefile~mesh.f90 sourcefile~poisson_fft.f90->sourcefile~tdsops.f90 sourcefile~tdsops.f90->sourcefile~common.f90~3 sourcefile~time_integrator.f90->sourcefile~allocator.f90~2 sourcefile~time_integrator.f90->sourcefile~backend.f90 sourcefile~time_integrator.f90->sourcefile~common.f90~3 sourcefile~time_integrator.f90->sourcefile~field.f90 sourcefile~vector_calculus.f90->sourcefile~allocator.f90~2 sourcefile~vector_calculus.f90->sourcefile~backend.f90 sourcefile~vector_calculus.f90->sourcefile~common.f90~3 sourcefile~vector_calculus.f90->sourcefile~field.f90 sourcefile~vector_calculus.f90->sourcefile~tdsops.f90

Files dependent on this one

sourcefile~~tgv.f90~~AfferentGraph sourcefile~tgv.f90 tgv.f90 sourcefile~xcompact.f90 xcompact.f90 sourcefile~xcompact.f90->sourcefile~tgv.f90

Source Code

module m_case_tgv
  use iso_fortran_env, only: stderr => error_unit

  use m_allocator, only: allocator_t
  use m_base_backend, only: base_backend_t
  use m_base_case, only: base_case_t
  use m_common, only: dp, VERT, DIR_C
  use m_field, only: field_t
  use m_mesh, only: mesh_t
  use m_solver, only: init

  implicit none

  type, extends(base_case_t) :: case_tgv_t
  contains
    procedure :: boundary_conditions => boundary_conditions_tgv
    procedure :: initial_conditions => initial_conditions_tgv
    procedure :: forcings => forcings_tgv
    procedure :: pre_correction => pre_correction_tgv
    procedure :: postprocess => postprocess_tgv
  end type case_tgv_t

  interface case_tgv_t
    module procedure case_tgv_init
  end interface case_tgv_t

contains

  function case_tgv_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_tgv_t) :: flow_case

    call flow_case%case_init(backend, mesh, host_allocator)

  end function case_tgv_init

  subroutine initial_conditions_tgv(self)
    implicit none

    class(case_tgv_t) :: self

    call self%set_init(self%solver%u, u_func)
    call self%set_init(self%solver%v, v_func)
    call self%solver%w%fill(0._dp)

    call self%solver%u%set_data_loc(VERT)
    call self%solver%v%set_data_loc(VERT)
    call self%solver%w%set_data_loc(VERT)

  end subroutine initial_conditions_tgv

  pure function u_func(coords) result(r)
    implicit none

    real(dp), intent(in) :: coords(3)
    real(dp) :: r

    r = sin(coords(1))*cos(coords(2))*cos(coords(3))
  end function u_func

  pure function v_func(coords) result(r)
    implicit none

    real(dp), intent(in) :: coords(3)
    real(dp) :: r

    r = -cos(coords(1))*sin(coords(2))*cos(coords(3))
  end function v_func

  subroutine boundary_conditions_tgv(self)
    implicit none

    class(case_tgv_t) :: self

    ! do nothing for TGV case
  end subroutine boundary_conditions_tgv

  subroutine forcings_tgv(self, du, dv, dw, iter)
    implicit none

    class(case_tgv_t) :: self
    class(field_t), intent(inout) :: du, dv, dw
    integer, intent(in) :: iter

    ! do nothing for TGV case
  end subroutine forcings_tgv

  subroutine pre_correction_tgv(self, u, v, w)
    implicit none

    class(case_tgv_t) :: self
    class(field_t), intent(inout) :: u, v, w

    ! do nothing for TGV case
  end subroutine pre_correction_tgv

  subroutine postprocess_tgv(self, iter, t)
    implicit none

    class(case_tgv_t) :: self
    integer, intent(in) :: iter
    real(dp), intent(in) :: t

    if (self%solver%mesh%par%is_root()) then
      print *, 'time =', t, 'iteration =', iter
    end if

    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_tgv

end module m_case_tgv