<p><b>mpetersen@lanl.gov</b> 2012-05-29 14:18:12 -0600 (Tue, 29 May 2012)</p><p>BRANCH COMMIT: Added capability to record instantaneous transport through sections to a text file. This occurs at the frequency config_hf_transport_interval, and is independant of the normal output stream interval.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/transport_sections/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/transport_sections/src/core_ocean/Registry        2012-05-29 16:00:15 UTC (rev 1943)
+++ branches/ocean_projects/transport_sections/src/core_ocean/Registry        2012-05-29 20:18:12 UTC (rev 1944)
@@ -20,7 +20,9 @@
namelist integer io config_frames_per_outfile 0
namelist integer io config_pio_num_iotasks 0
namelist integer io config_pio_stride 1
-namelist logical io config_record_transport false
+namelist logical io config_record_transport false
+namelist logical io config_record_hf_transport false
+namelist character io config_hf_transport_interval 24:00:00
namelist character decomposition config_block_decomp_file_prefix graph.info.part.
namelist integer decomposition config_number_of_blocks 0
namelist logical decomposition config_explicit_proc_decomp .false.
@@ -313,15 +315,17 @@
var persistent real acc_u ( nVertLevels nEdges Time ) 2 o acc_u state - -
var persistent real acc_uVar ( nVertLevels nEdges Time ) 2 o acc_uVar state - -
var persistent real acc_uTransportSection ( nSections Time ) 2 o acc_uTransportSection state - -
+var persistent real hf_uTransportSection ( nSections Time ) 2 - hf_uTransportSection state - -
% Variables that contain sections for recording transport, read from grid or restart file
% To create these variables, see mpas_repo/branches/tools/transport_sections
var persistent integer sectionEdgeIndex ( maxNEdgesInSection nSections ) 0 - sectionEdgeIndex mesh - -
var persistent integer sectionEdgeSign ( maxNEdgesInSection nSections ) 0 - sectionEdgeSign mesh - -
var persistent integer nEdgesInSection ( nSections ) 0 iro nEdgesInSection mesh - -
-var persistent char sectionText ( CharLength120 nSections ) 0 - sectionText mesh - -
-var persistent char sectionAbbreviation ( CharLength8 nSections ) 0 - sectionAbbreviation mesh - -
var persistent real sectionCoord ( latLonPairs nSections ) 0 - sectionCoord mesh - -
+% 2D string variables are not yet allowed in mpas. Hopefully we can include these soon. mrp.
+%var persistent text sectionText ( CharLength120 nSections ) 0 - sectionText mesh - -
+%var persistent text sectionAbbreviation ( CharLength8 nSections ) 0 - sectionAbbreviation mesh - -
% Variables that contain sections for recording transport, not read from grid or restart file
% These variables are initialized upon start-up from the global variables.
Modified: branches/ocean_projects/transport_sections/src/core_ocean/mpas_ocn_mpas_core.F
===================================================================
--- branches/ocean_projects/transport_sections/src/core_ocean/mpas_ocn_mpas_core.F        2012-05-29 16:00:15 UTC (rev 1943)
+++ branches/ocean_projects/transport_sections/src/core_ocean/mpas_ocn_mpas_core.F        2012-05-29 20:18:12 UTC (rev 1944)
@@ -39,6 +39,7 @@
integer, parameter :: outputAlarmID = 1
integer, parameter :: restartAlarmID = 2
integer, parameter :: statsAlarmID = 3
+ integer, parameter :: hfTransportAlarmID = 4
type (timer_node), pointer :: globalDiagTimer, timeIntTimer
type (timer_node), pointer :: initDiagSolveTimer
@@ -228,6 +229,11 @@
! call mpas_add_clock_alarm(clock, statsAlarmID, alarmStartTime, alarmTimeStep, ierr=ierr)
!end if
+ ! set high frequency transport output alarm
+ call mpas_set_timeInterval(alarmTimeStep, timeString=config_hf_transport_interval, ierr=ierr)
+ alarmStartTime = startTime + alarmTimeStep
+ call mpas_add_clock_alarm(clock, hfTransportAlarmID, alarmStartTime, alarmTimeStep, ierr=ierr)
+
call mpas_get_time(curr_time=startTime, dateTimeString=startTimeStamp, ierr=ierr)
end subroutine ocn_simulation_clock_init!}}}
@@ -267,27 +273,35 @@
read(15,*) mesh % sectionCoord % array
close(15)
+!mrp 120529
+! 2D string variables are not yet allowed in mpas. Hopefully we can include these soon. mrp.
+! open (unit=15, file = "sectionAbbreviation.txt")
+! do j=1,mesh % nSections
+! read(15,'(A8)') mesh % sectionAbbreviation % array(:,j)
+! end do
+! close(15)
+
! mrp debugging output
- print *, 'mesh % nSections',mesh % nSections
- print *, 'mesh % maxNEdgesInSection',mesh % maxNEdgesInSection
- print *, 'mesh % latLonPairs',mesh % latLonPairs
- print *, 'mesh % CharLength8',mesh % CharLength8
- print *, 'mesh % CharLength120',mesh % CharLength120
+ !print *, 'mesh % nSections',mesh % nSections
+ !print *, 'mesh % maxNEdgesInSection',mesh % maxNEdgesInSection
+ !print *, 'mesh % latLonPairs',mesh % latLonPairs
+ !print *, 'mesh % CharLength8',mesh % CharLength8
+ !print *, 'mesh % CharLength120',mesh % CharLength120
- print *, 'size(mesh % nEdgesInSection',size(mesh % nEdgesInSection % array)
- print *, 'size(mesh % sectionEdgeIndex',size(mesh % sectionEdgeIndex % array,1)
- print *, 'size(mesh % sectionEdgeIndex',size(mesh % sectionEdgeIndex % array,2)
- print *, 'size(mesh % sectionEdgeSign',size(mesh % sectionEdgeSign % array,1)
- print *, 'size(mesh % sectionEdgeSign',size(mesh % sectionEdgeSign % array,2)
- print *, 'size(mesh % sectionCoord',size(mesh % sectionCoord % array,1)
- print *, 'size(mesh % sectionCoord',size(mesh % sectionCoord % array,2)
- do j=1,mesh % nSections
- print *, 'nSection = ',j
- print *, 'mesh % nEdgesInSection',mesh % nEdgesInSection % array(j)
- print *, 'mesh % sectionCoord',mesh % sectionCoord % array(:,j)
- print *, 'mesh % sectionEdgeIndex',mesh % sectionEdgeIndex % array(:,j)
- print *, 'mesh % sectionEdgeSign',mesh % sectionEdgeSign % array(:,j)
- end do
+ !print *, 'size(mesh % nEdgesInSection',size(mesh % nEdgesInSection % array)
+ !print *, 'size(mesh % sectionEdgeIndex',size(mesh % sectionEdgeIndex % array,1)
+ !print *, 'size(mesh % sectionEdgeIndex',size(mesh % sectionEdgeIndex % array,2)
+ !print *, 'size(mesh % sectionEdgeSign',size(mesh % sectionEdgeSign % array,1)
+ !print *, 'size(mesh % sectionEdgeSign',size(mesh % sectionEdgeSign % array,2)
+ !print *, 'size(mesh % sectionCoord',size(mesh % sectionCoord % array,1)
+ !print *, 'size(mesh % sectionCoord',size(mesh % sectionCoord % array,2)
+! do j=1,mesh % nSections
+ !print *, 'iSection = ',j
+ !print *, 'mesh % nEdgesInSection',mesh % nEdgesInSection % array(j)
+ !print *, 'mesh % sectionCoord',mesh % sectionCoord % array(:,j)
+ !print *, 'mesh % sectionEdgeIndex',mesh % sectionEdgeIndex % array(:,j)
+ !print *, 'mesh % sectionEdgeSign',mesh % sectionEdgeSign % array(:,j)
+! end do
! mrp debugging output end
end subroutine ocn_initialize_transport_sections!}}}
@@ -552,6 +566,11 @@
call mpas_output_state_finalize(restart_obj, domain % dminfo)
end if
+ if (mpas_is_alarm_ringing(clock, hfTransportAlarmID, ierr=ierr)) then
+ call mpas_reset_clock_alarm(clock, hfTransportAlarmID, ierr=ierr)
+ call ocn_write_hf_transport(domain,timeStamp)
+ end if
+
end do
end subroutine mpas_core_run!}}}
@@ -596,7 +615,42 @@
end if
end subroutine ocn_write_output_frame!}}}
+
+ subroutine ocn_write_hf_transport(domain,timeStamp)!{{{
+ !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
+ ! Compute diagnostic fields for a domain and write model state to output file
+ !
+ ! Input/Output: domain - contains model state; diagnostic field are computed
+ ! before returning
+ !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
+ use mpas_grid_types
+ use mpas_io_output
+
+ implicit none
+
+ type (domain_type), intent(inout) :: domain
+ character(len=StrKIND), intent(in) :: timeStamp
+
+ type (block_type), pointer :: block_ptr
+
+ if (.not.config_record_hf_transport) return
+
+ block_ptr => domain % blocklist
+ do while (associated(block_ptr))
+ call ocn_hf_transport_compute(block_ptr % mesh, block_ptr % state % time_levs(1) % state)
+ block_ptr => block_ptr % next
+ end do
+
+ block_ptr => domain % blocklist
+ do while (associated(block_ptr))
+ call ocn_hf_transport_write_to_file(block_ptr % mesh, block_ptr % state % time_levs(1) % state, domain % dminfo,timeStamp)
+ block_ptr => block_ptr % next
+ end do
+
+ end subroutine ocn_write_hf_transport!}}}
+
+
subroutine ocn_compute_output_diagnostics(state, grid)!{{{
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Compute diagnostic fields for a domain
Modified: branches/ocean_projects/transport_sections/src/core_ocean/mpas_ocn_time_average.F
===================================================================
--- branches/ocean_projects/transport_sections/src/core_ocean/mpas_ocn_time_average.F        2012-05-29 16:00:15 UTC (rev 1943)
+++ branches/ocean_projects/transport_sections/src/core_ocean/mpas_ocn_time_average.F        2012-05-29 20:18:12 UTC (rev 1944)
@@ -159,8 +159,7 @@
type (state_type), intent(inout) :: state
type (dm_info) :: dminfo
- integer :: &
- nSections
+ integer :: nSections
real, pointer :: nAccumulate
@@ -202,4 +201,83 @@
end subroutine ocn_time_average_normalize!}}}
+ subroutine ocn_hf_transport_compute(mesh, state)!{{{
+ ! Compute high freqency instantaneous transports through sections.
+
+ type (mesh_type), intent(in) :: mesh
+ type (state_type), intent(inout) :: state
+
+ integer :: iSection, iLocalEdge, iEdge, k
+
+ integer :: nSections
+ integer, dimension(:), pointer :: nLocalEdgesInSection, maxLevelEdgeTop
+ integer, dimension(:,:), pointer :: localSectionEdgeIndex, localSectionEdgeSign
+
+ real (kind=RKIND) :: tempSum, m3ps_to_Sv
+
+ real (kind=RKIND), dimension(:), pointer :: hf_uTransportSection, dvEdge, hZLevel
+
+ real (kind=RKIND), dimension(:,:), pointer :: u, h_edge
+
+ m3ps_to_Sv = 1e-6 ! m^3/sec flux to Sverdrups
+
+ nSections = mesh % nSections
+ nLocalEdgesInSection => mesh % nLocalEdgesInSection % array
+ maxLevelEdgeTop => mesh % maxLevelEdgeTop % array
+ dvEdge => mesh % dvEdge % array
+ localSectionEdgeIndex => mesh % localSectionEdgeIndex % array
+ localSectionEdgeSign => mesh % localSectionEdgeSign % array
+ hZLevel => mesh % hZLevel % array
+
+ u => state % u % array
+ h_edge => state % h_edge % array
+
+ hf_uTransportSection => state % hf_uTransportSection % array
+
+ ! Sum up transport over all edges in this processors partition,
+ ! for each section. Here we use h_edge for best precision.
+ ! To match with the matlab postprocessing code transport_sections.m,
+ ! replace h_edge(k,iEdge) with hZLevel(k).
+ hf_uTransportSection = 0.0
+ do iSection = 1,nSections
+ do iLocalEdge = 1,nLocalEdgesInSection(iSection)
+ iEdge = localSectionEdgeIndex(iLocalEdge,iSection)
+ do k=1,maxLevelEdgeTop(iEdge)
+ hf_uTransportSection(iSection) = hf_uTransportSection(iSection) &
+ + localSectionEdgeSign(iLocalEdge,iSection) * u(k,iEdge)*dvEdge(iEdge) * h_edge(k,iEdge)*m3ps_to_Sv;
+ end do
+ end do
+ end do
+
+ end subroutine ocn_hf_transport_compute!}}}
+
+ subroutine ocn_hf_transport_write_to_file(mesh,state,dminfo,timeStamp)!{{{
+ ! Sum high freqency instantaneous transports through sections across processors.
+ ! Write the resulting value to a file.
+
+ type (mesh_type), intent(in) :: mesh
+ type (state_type), intent(inout) :: state
+ type (dm_info), intent(in) :: dminfo
+ character(len=StrKIND), intent(in) :: timeStamp
+
+ integer :: nSections
+
+ real (kind=RKIND), dimension(:), pointer :: hf_uTransportSection
+ real (kind=RKIND), dimension(:), allocatable :: temp_hf_uTransportSection
+
+ nSections = mesh % nSections
+ hf_uTransportSection => state % hf_uTransportSection % array
+ allocate(temp_hf_uTransportSection(nSections))
+ temp_hf_uTransportSection = hf_uTransportSection
+ call mpas_dmpar_sum_real_array(dminfo, nSections, temp_hf_uTransportSection, hf_uTransportSection)
+ deallocate(temp_hf_uTransportSection)
+
+ if (dminfo % my_proc_id.eq.0) then
+ open(15,file='uTransportSection.txt',access='append')
+ write (15,'(a,100es14.6)') trim(timeStamp),hf_uTransportSection
+ close(15)
+ end if
+
+ end subroutine ocn_hf_transport_write_to_file!}}}
+
end module ocn_time_average
</font>
</pre>