Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix OMP race conditions and some chkcsum calls #142

Merged
merged 7 commits into from
Mar 9, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion src/core/MOM.F90
Original file line number Diff line number Diff line change
Expand Up @@ -2897,7 +2897,7 @@ subroutine extract_surface_state(CS, sfc_state)


if (allocated(sfc_state%melt_potential)) then
!$OMP parallel do default(shared)
!$OMP parallel do default(shared) private(depth_ml, dh, T_freeze, depth, delT)
do j=js,je
do i=is,ie
depth(i) = 0.0
Expand Down
2 changes: 1 addition & 1 deletion src/diagnostics/MOM_diagnostics.F90
Original file line number Diff line number Diff line change
Expand Up @@ -600,7 +600,7 @@ subroutine calculate_diagnostic_fields(u, v, h, uh, vh, tv, ADp, CDp, p_surf, &
if (CS%id_rhopot2 > 0) call post_data(CS%id_rhopot2, Rcv, CS%diag)
endif
if (CS%id_rhoinsitu > 0) then
!$OMP parallel do default(none) shared(tv,Rcv,is,ie,js,je,nz,pressure_1d,h,GV)
!$OMP parallel do default(none) shared(tv,Rcv,is,ie,js,je,nz,h,GV) private(pressure_1d)
do j=js,je
pressure_1d(:) = 0. ! Start at p=0 Pa at surface
do k=1,nz
Expand Down
3 changes: 2 additions & 1 deletion src/parameterizations/lateral/MOM_MEKE.F90
Original file line number Diff line number Diff line change
Expand Up @@ -286,7 +286,8 @@ subroutine step_forward_MEKE(MEKE, h, SN_u, SN_v, visc, dt, G, GV, US, CS, hu, h
! Calculates bottomFac2, barotrFac2 and LmixScale
call MEKE_lengthScales(CS, MEKE, G, GV, US, SN_u, SN_v, MEKE%MEKE, bottomFac2, barotrFac2, LmixScale)
if (CS%debug) then
call uvchksum("MEKE drag_vel_[uv]", drag_vel_u, drag_vel_v, G%HI, scale=US%Z_to_m*US%s_to_T)
if (CS%visc_drag) &
call uvchksum("MEKE drag_vel_[uv]", drag_vel_u, drag_vel_v, G%HI, scale=US%Z_to_m*US%s_to_T)
call hchksum(mass, 'MEKE mass',G%HI,haloshift=1, scale=US%R_to_kg_m3*US%Z_to_m)
call hchksum(drag_rate_visc, 'MEKE drag_rate_visc', G%HI, scale=US%L_T_to_m_s)
call hchksum(bottomFac2, 'MEKE bottomFac2', G%HI)
Expand Down
20 changes: 10 additions & 10 deletions src/parameterizations/vertical/MOM_CVMix_KPP.F90
Original file line number Diff line number Diff line change
Expand Up @@ -962,16 +962,16 @@ subroutine KPP_compute_BLD(CS, G, GV, US, h, Temp, Salt, u, v, EOS, uStar, buoyF
buoy_scale = US%L_to_m**2*US%s_to_T**3

! loop over horizontal points on processor
!$OMP parallel do default(none) private(surfFricVel, iFaceHeight, hcorr, dh, cellHeight, &
!$OMP surfBuoyFlux, U_H, V_H, u, v, Coriolis, pRef, SLdepth_0d, &
!$OMP ksfc, surfHtemp, surfHsalt, surfHu, surfHv, surfHuS, &
!$OMP surfHvS, hTot, delH, surftemp, surfsalt, surfu, surfv, &
!$OMP surfUs, surfVs, Uk, Vk, deltaU2, km1, kk, pres_1D, &
!$OMP Temp_1D, salt_1D, surfBuoyFlux2, MLD_GUESS, LA, rho_1D, &
!$OMP deltarho, N2_1d, ws_1d, LangEnhVT2, enhvt2, wst, &
!$OMP BulkRi_1d, zBottomMinusOffset) &
!$OMP shared(G, GV, CS, US, uStar, h, buoy_scale, buoyFlux, &
!$OMP Temp, Salt, waves, EOS, GoRho)
!GOMP parallel do default(none) private(surfFricVel, iFaceHeight, hcorr, dh, cellHeight, &
!GOMP surfBuoyFlux, U_H, V_H, u, v, Coriolis, pRef, SLdepth_0d, &
!GOMP ksfc, surfHtemp, surfHsalt, surfHu, surfHv, surfHuS, &
!GOMP surfHvS, hTot, delH, surftemp, surfsalt, surfu, surfv, &
!GOMP surfUs, surfVs, Uk, Vk, deltaU2, km1, kk, pres_1D, &
!GOMP Temp_1D, salt_1D, surfBuoyFlux2, MLD_GUESS, LA, rho_1D, &
!GOMP deltarho, N2_1d, ws_1d, LangEnhVT2, enhvt2, wst, &
!GOMP BulkRi_1d, zBottomMinusOffset) &
!GOMP shared(G, GV, CS, US, uStar, h, buoy_scale, buoyFlux, &
!GOMP Temp, Salt, waves, EOS, GoRho)
do j = G%jsc, G%jec
do i = G%isc, G%iec

Expand Down
1 change: 1 addition & 0 deletions src/parameterizations/vertical/MOM_diabatic_driver.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1332,6 +1332,7 @@ subroutine diabatic_ALE(u, v, h, tv, Hml, fluxes, visc, ADp, CDp, dt, Time_end,
nkmb = GV%nk_rho_varies
h_neglect = GV%H_subroundoff ; h_neglect2 = h_neglect*h_neglect
Kd_heat(:,:,:) = 0.0 ; Kd_salt(:,:,:) = 0.0
ea_s(:,:,:) = 0.0; eb_s(:,:,:) = 0.0; ea_t(:,:,:) = 0.0; eb_t(:,:,:) = 0.0

showCallTree = callTree_showQuery()
if (showCallTree) call callTree_enter("diabatic_ALE(), MOM_diabatic_driver.F90")
Expand Down
5 changes: 3 additions & 2 deletions src/parameterizations/vertical/MOM_set_viscosity.F90
Original file line number Diff line number Diff line change
Expand Up @@ -365,10 +365,11 @@ subroutine set_viscous_BBL(u, v, h, tv, visc, G, GV, US, CS, symmetrize)

if (.not.use_BBL_EOS) Rml_vel(:,:) = 0.0

!$OMP parallel do default(private) shared(u,v,h,tv,visc,G,GV,US,CS,Rml,is,ie,js,je,nz,nkmb, &
!$OMP parallel do default(private) shared(u,v,h,tv,visc,G,GV,US,CS,Rml,nz,nkmb, &
!$OMP nkml,Isq,Ieq,Jsq,Jeq,h_neglect,Rho0x400_G,C2pi_3, &
!$OMP U_bg_sq,cdrag_sqrt_Z,cdrag_sqrt,K2,use_BBL_EOS, &
!$OMP OBC,maxitt,Vol_quit,D_u,D_v,mask_u,mask_v)
!$OMP OBC,maxitt,D_u,D_v,mask_u,mask_v) &
!$OMP firstprivate(Vol_quit)
do j=Jsq,Jeq ; do m=1,2

if (m==1) then
Expand Down
6 changes: 3 additions & 3 deletions src/tracer/MOM_tracer_advect.F90
Original file line number Diff line number Diff line change
Expand Up @@ -250,9 +250,9 @@ subroutine advect_tracer(h_end, uhtr, vhtr, OBC, dt, G, GV, US, CS, Reg, &
isv = isv + stencil ; iev = iev - stencil
jsv = jsv + stencil ; jev = jev - stencil

!$OMP parallel do default(none) shared(nz,domore_k,x_first,Tr,hprev,uhr,uh_neglect, &
!$OMP OBC,domore_u,ntr,Idt,isv,iev,jsv,jev,stencil, &
!$OMP G,GV,CS,vhr,vh_neglect,domore_v,US)
!GOMP parallel do default(none) shared(nz,domore_k,x_first,Tr,hprev,uhr,uh_neglect, &
!GOMP OBC,domore_u,ntr,Idt,isv,iev,jsv,jev,stencil, &
!GOMP G,GV,CS,vhr,vh_neglect,domore_v,US)

! To ensure positive definiteness of the thickness at each iteration, the
! mass fluxes out of each layer are checked each step, and limited to keep
Expand Down