Skip to content

Commit bf3d7a6

Browse files
mark-petersendouglasjacobsen
authored andcommitted
Fix unitialized pointer bug.
Before, highFreqThicknessNew would be passed into ocn_vert_transport_velocity_top even when it was an unitialized pointer. Change highFreqThicknessNew to be an optional variable. Otherwise, the pgi compiler always errors out on this line in debug mode.
1 parent 73c5c4d commit bf3d7a6

File tree

4 files changed

+50
-20
lines changed

4 files changed

+50
-20
lines changed

src/core_ocean/mode_forward/mpas_ocn_time_integration_rk4.F

Lines changed: 23 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -387,18 +387,32 @@ subroutine ocn_time_integrator_rk4(domain, dt)!{{{
387387
call mpas_pool_get_array(provisStatePool, 'normalVelocity', normalVelocityProvis, 1)
388388
call mpas_pool_get_array(provisStatePool, 'highFreqThickness', highFreqThicknessProvis, 1)
389389

390-
! advection of u uses u, while advection of layerThickness and tracers use uTransport.
391-
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
392-
layerThicknessCur,layerThicknessEdge, normalVelocityProvis, &
393-
sshCur, highFreqThicknessProvis, rk_substep_weights(rk_step), &
394-
vertTransportVelocityTop, err)
390+
! advection of u uses u, while advection of layerThickness and tracers use normalTransportVelocity.
391+
if (associated(highFreqThicknessProvis)) then
392+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
393+
layerThicknessCur,layerThicknessEdge, normalVelocityProvis, &
394+
sshCur, rk_substep_weights(rk_step), &
395+
vertTransportVelocityTop, err, highFreqThicknessProvis)
396+
else
397+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
398+
layerThicknessCur,layerThicknessEdge, normalVelocityProvis, &
399+
sshCur, rk_substep_weights(rk_step), &
400+
vertTransportVelocityTop, err)
401+
endif
395402

396403
call ocn_tend_vel(tendPool, provisStatePool, forcingPool, diagnosticsPool, meshPool, scratchPool, 1)
397404

398-
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
399-
layerThicknessCur, layerThicknessEdge, uTransport, &
400-
sshCur, highFreqThicknessProvis, rk_substep_weights(rk_step), &
401-
vertTransportVelocityTop, err)
405+
if (associated(highFreqThicknessProvis)) then
406+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
407+
layerThicknessCur, layerThicknessEdge, uTransport, &
408+
sshCur, rk_substep_weights(rk_step), &
409+
vertTransportVelocityTop, err, highFreqThicknessProvis)
410+
else
411+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
412+
layerThicknessCur, layerThicknessEdge, uTransport, &
413+
sshCur, rk_substep_weights(rk_step), &
414+
vertTransportVelocityTop, err)
415+
endif
402416

403417
call ocn_tend_thick(tendPool, forcingPool, diagnosticsPool, meshPool)
404418

src/core_ocean/mode_forward/mpas_ocn_time_integration_split.F

Lines changed: 18 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -406,9 +406,15 @@ subroutine ocn_time_integrator_split(domain, dt)!{{{
406406
407407
! compute vertTransportVelocityTop. Use u (rather than uTransport) for momentum advection.
408408
! Use the most recent time level available.
409-
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
410-
layerThicknessCur, layerThicknessEdge, normalVelocityCur, &
411-
sshCur, highFreqThicknessNew, dt, vertTransportVelocityTop, err)
409+
if (associated(highFreqThicknessNew)) then
410+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
411+
layerThicknessCur, layerThicknessEdge, normalVelocityCur, &
412+
sshCur, dt, vertTransportVelocityTop, err, highFreqThicknessNew)
413+
else
414+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
415+
layerThicknessCur, layerThicknessEdge, normalVelocityCur, &
416+
sshCur, dt, vertTransportVelocityTop, err)
417+
endif
412418
413419
call ocn_tend_vel(tendPool, statePool, forcingPool, diagnosticsPool, meshPool, scratchPool, stage1_tend_time)
414420
@@ -1158,9 +1164,15 @@ subroutine ocn_time_integrator_split(domain, dt)!{{{
11581164
! compute vertTransportVelocityTop. Use uTransport for advection of layerThickness and tracers.
11591165
! Use time level 1 values of layerThickness and layerThicknessEdge because
11601166
! layerThickness has not yet been computed for time level 2.
1161-
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
1162-
layerThicknessCur, layerThicknessEdge, uTransport, &
1163-
sshCur, highFreqThicknessNew, dt, vertTransportVelocityTop, err)
1167+
if (associated(highFreqThicknessNew)) then
1168+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
1169+
layerThicknessCur, layerThicknessEdge, uTransport, &
1170+
sshCur, dt, vertTransportVelocityTop, err, highFreqThicknessNew)
1171+
else
1172+
call ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, &
1173+
layerThicknessCur, layerThicknessEdge, uTransport, &
1174+
sshCur, dt, vertTransportVelocityTop, err)
1175+
endif
11641176

11651177
call ocn_tend_thick(tendPool, forcingPool, diagnosticsPool, meshPool)
11661178

src/core_ocean/shared/mpas_ocn_diagnostics.F

Lines changed: 7 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -574,7 +574,7 @@ end subroutine ocn_diagnostic_solve!}}}
574574
!
575575
!-----------------------------------------------------------------------
576576
subroutine ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, oldLayerThickness, layerThicknessEdge, &
577-
normalVelocity, oldSSH, newHighFreqThickness, dt, vertTransportVelocityTop, err)!{{{
577+
normalVelocity, oldSSH, dt, vertTransportVelocityTop, err, newHighFreqThickness)!{{{
578578

579579
!-----------------------------------------------------------------
580580
!
@@ -600,7 +600,7 @@ subroutine ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, oldLayerT
600600
real (kind=RKIND), dimension(:), intent(in) :: &
601601
oldSSH !< Input: sea surface height at old time
602602

603-
real (kind=RKIND), dimension(:,:), intent(in) :: &
603+
real (kind=RKIND), dimension(:,:), intent(in), optional :: &
604604
newHighFreqThickness !< Input: high frequency thickness. Alters ALE thickness.
605605

606606
real (kind=RKIND), intent(in) :: &
@@ -684,7 +684,11 @@ subroutine ocn_vert_transport_velocity_top(meshPool, verticalMeshPool, oldLayerT
684684
!
685685
! Compute desired thickness at new time
686686
!
687-
call ocn_ALE_thickness(meshPool, verticalMeshPool, oldSSH, div_hu_btr, newHighFreqThickness, dt, ALE_thickness, err)
687+
if (present(newHighFreqThickness)) then
688+
call ocn_ALE_thickness(meshPool, verticalMeshPool, oldSSH, div_hu_btr, dt, ALE_thickness, err, newHighFreqThickness)
689+
else
690+
call ocn_ALE_thickness(meshPool, verticalMeshPool, oldSSH, div_hu_btr, dt, ALE_thickness, err)
691+
endif
688692

689693
!
690694
! Vertical transport through layer interfaces

src/core_ocean/shared/mpas_ocn_thick_ale.F

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -70,7 +70,7 @@ module ocn_thick_ale
7070
!> (z-tilde), and imposes a minimum layer thickness.
7171
!
7272
!-----------------------------------------------------------------------
73-
subroutine ocn_ALE_thickness(meshPool, verticalMeshPool, oldSSH, div_hu_btr, newHighFreqThickness, dt, ALE_thickness, err)!{{{
73+
subroutine ocn_ALE_thickness(meshPool, verticalMeshPool, oldSSH, div_hu_btr, dt, ALE_thickness, err, newHighFreqThickness)!{{{
7474

7575
!-----------------------------------------------------------------
7676
!
@@ -88,7 +88,7 @@ subroutine ocn_ALE_thickness(meshPool, verticalMeshPool, oldSSH, div_hu_btr, new
8888
oldSSH, &!< Input: sea surface height at old time
8989
div_hu_btr !< Input: thickness-weighted barotropic divergence
9090

91-
real (kind=RKIND), dimension(:,:), intent(in) :: &
91+
real (kind=RKIND), dimension(:,:), intent(in), optional :: &
9292
newHighFreqThickness !< Input: high frequency thickness. Alters ALE thickness.
9393

9494
real (kind=RKIND), intent(in) :: &

0 commit comments

Comments
 (0)