Skip to content

Commit

Permalink
Finishing touches
Browse files Browse the repository at this point in the history
added ci

added ci

things2

things

added single precision to test suite and generated golden files

added working precision to case_dicts, testing, and documentation

fixed code not passing using double working precision
  • Loading branch information
wilfonba committed May 16, 2023
1 parent ebadc94 commit 843a50b
Show file tree
Hide file tree
Showing 32 changed files with 191 additions and 184 deletions.
16 changes: 14 additions & 2 deletions .github/workflows/ci.yml
Original file line number Diff line number Diff line change
Expand Up @@ -86,17 +86,29 @@ jobs:
${{ matrix.intel-command }}
/bin/bash mfc.sh build -j $(nproc) ${{ matrix.debug }}
- name: Test Suite (Debug)
- name: Test Suite (Debug, Double Precision)
if: matrix.debug == '--debug'
run: |
${{ matrix.intel-command }}
/bin/bash mfc.sh test -j $(nproc) --debug
- name: Test Suite (No Debug)
- name: Test Suite (Debug, Single Precision)
if: matrix.debug == '--debug'
run: |
${{ matrix.intel-command }}
/bin/bash mfc.sh test -j $(nproc) --debug --single
- name: Test Suite (No Debug, Double Precision)
if: matrix.debug == '--no-debug'
run: |
${{ matrix.intel-command }}
/bin/bash mfc.sh test -j $(nproc) -a
- name: Test Suite (No Debug, Single Precision)
if: matrix.debug == '--no-debug'
run: |
${{ matrix.intel-command }}
/bin/bash mfc.sh test -j $(nproc) -a --single
self-cpu-release:
Expand Down
7 changes: 7 additions & 0 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,7 @@ option(MFC_SIMULATION "Build simulation" OFF)
option(MFC_POST_PROCESS "Build post_process" OFF)
option(MFC_DOCUMENTATION "Build documentation" OFF)
option(MFC_ALL "Build everything" OFF)
option(MFC_SINGLE_PRECISION "Build double precision" OFF)

if (MFC_BUILD_ALL)
set(MFC_PRE_PROCESS ON FORCE)
Expand Down Expand Up @@ -156,6 +157,12 @@ if (CMAKE_BUILD_TYPE STREQUAL "Debug")
add_compile_definitions(MFC_DEBUG)
endif()

if (MFC_SINGLE_PRECISION)
add_compile_definitions(MFC_SINGLE_PRECISION)
else()
add_compile_definitions(MFC_DOUBLE_PRECISION)
endif()

## === HANDLE_SOURCES
# Gather F90 source files for a given target, including common/, and preprocessing .fpp -> .f90.
# Outputs:
Expand Down
3 changes: 3 additions & 0 deletions docs/documentation/case.md
Original file line number Diff line number Diff line change
Expand Up @@ -76,8 +76,11 @@ Definition of the parameters is described in the following subsections.
| Parameter | Type | Description |
| ---: | :----: | :--- |
| `run_time_info` | Logical | Output run-time information |
| `working_precision` | Integer | Simulation working precision |

- `run_time_info` generates a text file that includes run-time information including the CFL number(s) at each time-step.
- `working_precision` sets the working precision to either (1) - single precision, or
(2) - double precision.

### 2. Computational Domain

Expand Down
3 changes: 2 additions & 1 deletion misc/run-phoenix-release-cpu.sh
Original file line number Diff line number Diff line change
Expand Up @@ -11,4 +11,5 @@
cd $SLURM_SUBMIT_DIR # Change to working directory
echo $(pwd)
. ./mfc.sh load -c p -m g
./mfc.sh test -j 12 -b mpirun -a
./mfc.sh test -j 12 -b mpirun -a --no-single
./mfc.sh test -j 12 -b mpirun -a --single
3 changes: 2 additions & 1 deletion misc/run-phoenix-release-gpu.sh
Original file line number Diff line number Diff line change
Expand Up @@ -11,4 +11,5 @@
cd $SLURM_SUBMIT_DIR # Change to working directory
echo $(pwd)
. ./mfc.sh load -c p -m g
./mfc.sh test -j 1 -b mpirun -a --gpu
./mfc.sh test -j 1 -b mpirun -a --gpu --no-single
./mfc.sh test -j 1 -b mpirun -a --gpu --single
4 changes: 2 additions & 2 deletions src/common/include/inline_conversions.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -32,11 +32,11 @@

if (mpp_lim .and. (num_fluids > 1)) then
c = (1._wp/gamma + 1._wp)* &
(pres + pi_inf)/rho
(pres + pi_inf/(gamma+1._wp))/rho
else
c = &
(1._wp/gamma + 1._wp)* &
(pres + pi_inf)/ &
(pres + pi_inf/(gamma + 1._wp))/ &
(rho*(1._wp - adv(num_fluids)))
end if
else
Expand Down
15 changes: 14 additions & 1 deletion src/common/m_constants.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,20 @@

module m_constants

use m_precision_select
use mpi

integer, parameter :: sp = selected_real_kind(6, 37)
integer, parameter :: dp = selected_real_kind(15, 307)

#ifdef MFC_DOUBLE_PRECISION
integer, parameter :: wp = dp
integer, parameter :: mpi_p = MPI_DOUBLE_PRECISION
#endif

#ifdef MFC_SINGLE_PRECISION
integer, parameter :: wp = sp
integer, parameter :: mpi_p = MPI_REAL
#endif

character, parameter :: dflt_char = ' ' !< Default string value

Expand Down
2 changes: 1 addition & 1 deletion src/common/m_derived_types.f90
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,7 @@
module m_derived_types

use m_constants !< Constants
use m_precision_select


implicit none

Expand Down
2 changes: 1 addition & 1 deletion src/common/m_eigen_solver.f90
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,7 @@
!! modifications for compatibility.
module m_eigen_solver

use m_precision_select
use m_constants

implicit none

Expand Down
20 changes: 10 additions & 10 deletions src/common/m_helper.f90
Original file line number Diff line number Diff line change
Expand Up @@ -92,11 +92,11 @@ end subroutine s_compute_finite_difference_coefficients ! --------------
!! @param ntmp is the output number bubble density
subroutine s_comp_n_from_prim(vftmp, Rtmp, ntmp, weights)
!$acc routine seq
real(kind(0._wp)), intent(IN) :: vftmp
real(kind(0._wp)), dimension(nb), intent(IN) :: Rtmp
real(kind(0._wp)), intent(OUT) :: ntmp
real(kind(0._wp)) :: R3
real(kind(0._wp)), dimension(nb) :: weights
real(wp), intent(IN) :: vftmp
real(wp), dimension(nb), intent(IN) :: Rtmp
real(wp), intent(OUT) :: ntmp
real(wp) :: R3
real(wp), dimension(nb) :: weights

R3 = dot_product(weights, Rtmp**3._wp)
ntmp = (3._wp/(4._wp*pi))*vftmp/R3
Expand All @@ -105,11 +105,11 @@ end subroutine s_comp_n_from_prim

subroutine s_comp_n_from_cons(vftmp, nRtmp, ntmp, weights)
!$acc routine seq
real(kind(0._wp)), intent(IN) :: vftmp
real(kind(0._wp)), dimension(nb), intent(IN) :: nRtmp
real(kind(0._wp)), intent(OUT) :: ntmp
real(kind(0._wp)) :: nR3
real(kind(0._wp)), dimension(nb) :: weights
real(wp), intent(IN) :: vftmp
real(wp), dimension(nb), intent(IN) :: nRtmp
real(wp), intent(OUT) :: ntmp
real(wp) :: nR3
real(wp), dimension(nb) :: weights

nR3 = dot_product(weights, nRtmp**3._wp)
ntmp = sqrt((4._wp*pi/3._wp)*nR3/vftmp)
Expand Down
23 changes: 0 additions & 23 deletions src/common/m_precision_select.f90

This file was deleted.

7 changes: 3 additions & 4 deletions src/common/m_variables_conversion.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -158,7 +158,7 @@ contains

pres = ( &
energy - &
0.5_wp*(mom**2._wp)/rho - &
(5._wp * (10._wp ** (-1)))*(mom**2._wp)/rho - &
pi_inf - E_e &
)/gamma
end if
Expand Down Expand Up @@ -235,7 +235,6 @@ contains

real(wp), optional, dimension(2), intent(OUT) :: Re_K
real(wp), dimension(num_fluids) :: alpha_rho_K, alpha_K
real(wp), optional, dimension(2), intent(OUT) :: Re_K

real(wp), optional, intent(OUT) :: G_K
real(wp), optional, dimension(num_fluids), intent(IN) :: G
Expand Down Expand Up @@ -510,7 +509,7 @@ contains
real(wp), dimension(num_fluids), intent(IN) :: alpha_rho_K, alpha_K !<
!! Partial densities and volume fractions

real(kind(0d0)), dimension(2), intent(OUT) :: Re_K
real(wp), dimension(2), intent(OUT) :: Re_K

integer, intent(IN) :: k, l, r
integer :: i, j !< Generic loop iterators
Expand Down Expand Up @@ -731,7 +730,7 @@ contains

integer :: i, j, k, l !< Generic loop iterators

real(kind(0._wp)) :: ntmp
real(wp) :: ntmp

if (bubbles) then
allocate(nRtmp(nb))
Expand Down
3 changes: 0 additions & 3 deletions src/post_process/m_checker.f90
Original file line number Diff line number Diff line change
Expand Up @@ -65,9 +65,6 @@ subroutine s_check_inputs()
! Constraints on model equations and number of fluids in the flow
elseif (all(model_eqns /= (/1, 2, 3, 4/))) then
call s_mpi_abort('Unsupported value of model_eqns. Exiting ...')
elseif (wp == single_precision .and. precision == 2) then
call s_mpi_abort('Unsupported combination of working precision'// &
'and silo precision')
elseif (num_fluids /= dflt_int &
.and. &
(num_fluids < 1 .or. num_fluids > num_fluids)) then
Expand Down
10 changes: 5 additions & 5 deletions src/post_process/m_data_output.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -140,7 +140,7 @@ contains
if (n == 0) then
allocate (q_root_sf(0:m_root, 0:0, 0:0))
if (precision == 1 .and. wp == double_precision) then
if (precision == 1) then
allocate (q_root_sf_s(0:m_root, 0:0, 0:0))
end if
end if
Expand Down Expand Up @@ -769,10 +769,10 @@ contains
! and write it to the formatted database master file.
if (n == 0) then
if (precision == 1 .and. wp == double_precision) then
if (precision == 1 .and. wp == dp) then
x_cc_s = real(x_cc, sp)
q_sf_s = real(q_sf, sp)
elseif (precision == 1 .and. wp == single_precision) then
elseif (precision == 1 .and. wp == sp) then
x_cc_s = x_cc
q_sf_s = q_sf
end if
Expand Down Expand Up @@ -866,7 +866,7 @@ contains
! Finally, each of the local processor(s) proceeds to write
! the flow variable data that it is responsible for to the
! formatted database slave file.
if (wp == double_precision) then
if (wp == dp) then
if (precision == 1) then
do i = -offset_x%beg, m + offset_x%end
do j = -offset_y%beg, n + offset_y%end
Expand Down Expand Up @@ -895,7 +895,7 @@ contains
end do
end if
end if
elseif (wp == single_precision) then
elseif (wp == dp) then
do i = -offset_x%beg, m + offset_x%end
do j = -offset_y%beg, n + offset_y%end
do k = -offset_z%beg, p + offset_z%end
Expand Down
52 changes: 26 additions & 26 deletions src/post_process/m_global_parameters.f90
Original file line number Diff line number Diff line change
Expand Up @@ -629,28 +629,28 @@ end subroutine s_initialize_global_parameters_module ! --------------------
subroutine s_initialize_nonpoly

integer :: ir
real(kind(0._wp)) :: rhol0
real(kind(0._wp)) :: pl0
real(kind(0._wp)) :: uu
real(kind(0._wp)) :: D_m
real(kind(0._wp)) :: temp
real(kind(0._wp)) :: omega_ref
real(kind(0._wp)), dimension(Nb) :: chi_vw0
real(kind(0._wp)), dimension(Nb) :: cp_m0
real(kind(0._wp)), dimension(Nb) :: k_m0
real(kind(0._wp)), dimension(Nb) :: rho_m0
real(kind(0._wp)), dimension(Nb) :: x_vw
real(wp) :: rhol0
real(wp) :: pl0
real(wp) :: uu
real(wp) :: D_m
real(wp) :: temp
real(wp) :: omega_ref
real(wp), dimension(Nb) :: chi_vw0
real(wp), dimension(Nb) :: cp_m0
real(wp), dimension(Nb) :: k_m0
real(wp), dimension(Nb) :: rho_m0
real(wp), dimension(Nb) :: x_vw

! liquid physical properties
real(kind(0._wp)) :: mul0, ss, pv, gamma_v, M_v, mu_v
real(wp) :: mul0, ss, pv, gamma_v, M_v, mu_v

! gas physical properties
real(kind(0._wp)) :: gamma_m, gamma_n, M_n, mu_n
real(wp) :: gamma_m, gamma_n, M_n, mu_n

! polytropic index used to compute isothermal natural frequency
real(kind(0._wp)), parameter :: k_poly = 1._wp
real(wp), parameter :: k_poly = 1._wp
! universal gas constant
real(kind(0._wp)), parameter :: Ru = 8314._wp
real(wp), parameter :: Ru = 8314._wp

rhol0 = rhoref
pl0 = pref
Expand Down Expand Up @@ -749,13 +749,13 @@ end subroutine s_initialize_nonpoly
!> Subroutine to compute the transfer coefficient for non-polytropic gas modeling
subroutine s_transcoeff(omega, peclet, Re_trans, Im_trans)

real(kind(0._wp)), intent(IN) :: omega
real(kind(0._wp)), intent(IN) :: peclet
real(kind(0._wp)), intent(OUT) :: Re_trans
real(kind(0._wp)), intent(OUT) :: Im_trans
real(wp), intent(IN) :: omega
real(wp), intent(IN) :: peclet
real(wp), intent(OUT) :: Re_trans
real(wp), intent(OUT) :: Im_trans
complex :: trans, c1, c2, c3
complex :: imag = (0., 1.)
real(kind(0._wp)) :: f_transcoeff
real(wp) :: f_transcoeff

c1 = imag*omega*peclet
c2 = CSQRT(c1)
Expand Down Expand Up @@ -844,12 +844,12 @@ subroutine s_simpson(Npt)

integer, intent(IN) :: Npt
integer :: ir
real(kind(0._wp)) :: R0mn
real(kind(0._wp)) :: R0mx
real(kind(0._wp)) :: dphi
real(kind(0._wp)) :: tmp
real(kind(0._wp)) :: sd
real(kind(0._wp)), dimension(Npt) :: phi
real(wp) :: R0mn
real(wp) :: R0mx
real(wp) :: dphi
real(wp) :: tmp
real(wp) :: sd
real(wp), dimension(Npt) :: phi

! nondiml. min. & max. initial radii for numerical quadrature
!sd = 0.05D0
Expand Down
10 changes: 2 additions & 8 deletions src/pre_process/m_assign_variables.f90
Original file line number Diff line number Diff line change
Expand Up @@ -238,15 +238,12 @@ subroutine s_assign_patch_species_primitive_variables_bubbles(patch_id, j, k, l,
real(wp) :: orig_pi_inf
real(wp) :: muR, muV

real(wp), dimension(sys_size) :: orig_prim_vf !<
real(wp), dimension(sys_size) :: orig_prim_vf !< Vector to hold original values of cell for smoothing purposes
real(wp), dimension(int(E_idx - mom_idx%beg)) :: vel !< velocity
real(wp) :: pres !< pressure
real(wp) :: x_centroid, y_centroid
real(wp) :: epsilon, beta

real(kind(0d0)), dimension(sys_size) :: orig_prim_vf !<
!! Vector to hold original values of cell for smoothing purposes

integer :: i !< Generic loop iterator
integer :: smooth_patch_id

Expand Down Expand Up @@ -648,15 +645,12 @@ subroutine s_assign_patch_species_primitive_variables(patch_id, j, k, l, &
!! function, respectively, obtained from the combination of primitive
!! variables of the current and smoothing patches

real(wp), dimension(sys_size) :: orig_prim_vf !<
real(wp), dimension(sys_size) :: orig_prim_vf !< Vector to hold original values of cell for smoothing purposes
real(wp), dimension(int(E_idx - mom_idx%beg)) :: vel !< velocity
real(wp) :: pres !< pressure
real(wp) :: x_centroid, y_centroid
real(wp) :: epsilon, beta

real(kind(0d0)), dimension(sys_size) :: orig_prim_vf !<
! Vector to hold original values of cell for smoothing purposes

integer :: smooth_patch_id

integer :: i !< generic loop iterator
Expand Down
Loading

0 comments on commit 843a50b

Please sign in to comment.