-
Notifications
You must be signed in to change notification settings - Fork 378
/
Copy pathmpas_ocn_diagnostics.F
4726 lines (4052 loc) · 184 KB
/
mpas_ocn_diagnostics.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! Copyright (c) 2013, Los Alamos National Security, LLC (LANS)
! and the University Corporation for Atmospheric Research (UCAR).
!
! Unless noted otherwise source code is licensed under the BSD license.
! Additional copyright and license information can be found in the LICENSE file
! distributed with this code, or at http://mpas-dev.github.io/license.html
!
!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
!
! ocn_diagnostics
!
!> \brief MPAS ocean diagnostics driver
!> \author Mark Petersen
!> \date 23 September 2011
!> \details
!> This module contains the routines for computing
!> diagnostic variables, and other quantities such as vertAleTransportTop.
!
!-----------------------------------------------------------------------
module ocn_diagnostics
use mpas_timer
use mpas_derived_types
use mpas_pool_routines
use mpas_constants
use mpas_threading
use mpas_vector_reconstruction
use mpas_stream_manager
use mpas_io_units
use ocn_constants
use ocn_config
use ocn_gm
use ocn_equation_of_state
use ocn_thick_ale
use ocn_diagnostics_variables
use ocn_mesh
use ocn_surface_land_ice_fluxes
use ocn_vertical_advection
use ocn_subgrid
implicit none
private
save
!--------------------------------------------------------------------
!
! Public parameters
!
!--------------------------------------------------------------------
!--------------------------------------------------------------------
!
! Public member functions
!
!--------------------------------------------------------------------
public :: ocn_diagnostic_solve, &
ocn_diagnostic_solve_layerThicknessEdge, &
ocn_diagnostic_solve_wctEdge, &
ocn_relativeVorticity_circulation, &
ocn_vert_transport_velocity_top, &
ocn_fuperp, &
ocn_filter_btr_mode_tend_vel, &
ocn_reconstruct_eddy_vectors, &
ocn_compute_mixing_input_fields, &
ocn_validate_state, &
ocn_build_log_filename, &
ocn_diagnostics_init
!--------------------------------------------------------------------
!
! Private module variables
!
!--------------------------------------------------------------------
integer :: ke_cell_flag, ke_vertex_flag
real (kind=RKIND) :: fCoef
real (kind=RKIND) :: landIceTopDragCoeff
real (kind=RKIND), pointer :: coef_3rd_order
! Methods for computing thickness at edges for flux calculations
integer :: &
thickEdgeFluxChoice, &! choice of thickness flux type
thickEdgeDragChoice ! choice of thickness drag type
integer, parameter :: &
thickEdgeFluxCenter = 1, &! use mean thickness of cell neighbors
thickEdgeFluxUpwind = 2, &! use upwind cell thickness at edge
thickEdgeFluxConstant = 3 ! use constant thickness in time from
! RestingThickness, for linear test problems
integer, parameter :: &
thickEdgeDragCenter = 1, &! use mean thickness of cell neighbors
thickEdgeDragHarMean = 2 ! use harmonic mean cell thickness at edge
!***********************************************************************
contains
!***********************************************************************
!
! routine ocn_diagnostic_solve
!
!> \brief Computes diagnostic variables
!> \author Mark Petersen
!> \date 23 September 2011
!> \details
!> This routine computes the diagnostic variables for the ocean
!
!-----------------------------------------------------------------------
subroutine ocn_diagnostic_solve(dt, statePool, forcingPool, meshPool, verticalMeshPool, scratchPool, tracersPool, &
timeLevelIn, full)!{{{
real (kind=RKIND), intent(in) :: dt !< Input: Time step
type (mpas_pool_type), intent(in) :: statePool !< Input: State information
type (mpas_pool_type), intent(in) :: forcingPool !< Input: Forcing information
type (mpas_pool_type), intent(in) :: meshPool !< Input: mesh information
type (mpas_pool_type), intent(in) :: verticalMeshPool !< Input: vertical mesh information
type (mpas_pool_type), intent(in) :: scratchPool !< Input: scratch variables
type (mpas_pool_type), intent(in) :: tracersPool !< Input: tracer fields
integer, intent(in), optional :: timeLevelIn !< Input: Time level in state
logical, intent(in), optional :: full !< Input: Run all computations?
integer :: iEdge, iCell, iTracer, k, kmin, kmax
integer :: err, nCells, nEdges, nVertices, numTracers
integer :: timeLevel
real (kind=RKIND) :: coef_3rd_order
real (kind=RKIND), dimension(:), pointer :: &
ssh, seaIcePressure, atmosphericPressure
real (kind=RKIND), dimension(:,:), pointer :: &
layerThickness, layerThicknessLag, normalVelocity, restingThickness
real (kind=RKIND), dimension(:,:,:), pointer :: activeTracers
integer, pointer :: indexTemperaturePtr, indexSalinityPtr
integer :: indexTemperature, indexSalinity
logical :: full_compute = .true.
real (kind=RKIND), dimension(:), pointer :: &
frazilSurfacePressure, landIcePressure, landIceDraft, landIceFraction
integer, dimension(:), pointer :: &
landIceFloatingMask
call mpas_timer_start('diagnostic solve')
if ( present(full) ) then
full_compute = full
else
full_compute = .true.
end if
if (present(timeLevelIn)) then
timeLevel = timeLevelIn
else
timeLevel = 1
end if
call mpas_pool_get_dimension(tracersPool, 'index_temperature', &
indexTemperaturePtr)
call mpas_pool_get_dimension(tracersPool, 'index_salinity', &
indexSalinityPtr)
indexTemperature = indexTemperaturePtr
indexSalinity = indexSalinityPtr
call mpas_pool_get_array(tracersPool, 'activeTracers', activeTracers, timeLevel)
call mpas_pool_get_array(statePool, 'layerThickness', layerThickness, timeLevel)
call mpas_pool_get_array(statePool, 'normalVelocity', normalVelocity, timeLevel)
call mpas_pool_get_array(statePool, 'ssh', ssh, timeLevel)
call mpas_pool_get_array(forcingPool, 'seaIcePressure', seaIcePressure)
call mpas_pool_get_array(forcingPool, 'atmosphericPressure', atmosphericPressure)
call mpas_pool_get_array(forcingPool, 'frazilSurfacePressure', frazilSurfacePressure)
call mpas_pool_get_array(verticalMeshPool, 'restingThickness', restingThickness)
if (landIcePressureOn) then
call mpas_pool_get_array(forcingPool, 'landIcePressure', landIcePressure)
call mpas_pool_get_array(forcingPool, 'landIceDraft', landIceDraft)
call mpas_pool_get_array(forcingPool, 'landIceFloatingMask', landIceFloatingMask)
call mpas_pool_get_array(forcingPool, 'landIceFraction', landIceFraction)
if (nCellsAll.gt.0 .AND. landIceFloatingMask(1) == -1) then
call mpas_log_write('landIceFloatingMask contains the default value which likely indicates that this field is missing in the initial condition file (e.g. because it is meant for an older E3SM version).', &
MPAS_LOG_CRIT)
endif
end if
nEdges = nEdgesAll
nCells = nCellsAll
nVertices = nVerticesAll
coef_3rd_order = config_coef_3rd_order
#ifdef MPAS_OPENACC
!! Listing new inputs required by the first routine where they are used
!! ocn_diagnostic_solve_layerThicknessEdge :: layerThickness, (not on device)
!! normalVelocity (not on device)
!! ocn_relativeVorticity_circulation ::
!! ocn_diagnostic_solve_vortVel :: normalTransportVelocity, (diagnostics array)
!! normalGMBolusVelocity (diagnostics array)
!! ocn_diagnostic_solve_kineticEnergy ::
!! ocn_diagnostic_solve_vorticity ::
!! ocn_diagnostic_solve_surface_pressure :: atmosphericPressure, (not on device)
!! seaIcePressure, (not on device)
!! frazilSurfacePressure, (not on device)
!! landIcePressure (not on device)
!! ocn_diagnostic_solve_z_coordinates ::
!! ocn_equation_of_state_density :: activeTracers, (not on device)
!! tracersSurfaceValue, (diagnostics array)
!! maxLevelCell, (mesh array)
!! refBottomDepth (mesh array)
!! ocn_diagnostic_solve_pressure ::
!! ocn_diagnostic_solve_richardson :: edgeAreaFractionOfCell (mesh array)
!! tracersSurfaceValue calc :: minLevelCell (mesh array)
!! normalVelocitySurfaceLayer calc :: minLevelEdgeBot (mesh array)
!! surfaceFluxAttenuationCoefficientRunoff :: layerThickEdgeFlux (diagnostics array)
!! surfaceFluxAttenuationCoefficientSubglacialRunoff :: layerThickEdgeFlux (diagnostics array)
!! ocn_diagnostic_solve_ssh :: ssh, (not on device)
!! landIceDraft (not on device)
! !$acc enter data copyin(layerThickness, normalVelocity)
! !$acc update device (normalTransportVelocity, &
! !$acc normalGMBolusVelocity)
! !$acc enter data copyin(atmosphericPressure, seaIcePressure)
! !$acc enter data copyin(activeTracers, ssh)
! !$acc enter data copyin(activeTracers)
! !$acc update device(tracersSurfaceValue)
! if ( associated(frazilSurfacePressure) ) then
! !$acc enter data copyin(frazilSurfacePressure)
! endif
! if (landIcePressureOn) then
! !$acc enter data copyin(landIcePressure)
! !$acc enter data copyin(landIceDraft)
! endif
#endif
!
! Z-coordinates
! This section must be placed in the code before computing the density.
!
! inputs : layerThickness
! outputs : zMid, zTop, ssh)
call ocn_diagnostic_solve_z_coordinates(layerThickness, zMid, zTop, ssh)
! inputs: layerThickness, normalVelocity
! output: layerThickEdgeMean, layerThickEdgeDrag, layerThickEdgeFlux
call ocn_diagnostic_solve_layerThicknessEdge(normalVelocity, &
layerThickness, restingThickness, &
ssh)
! inputs: normalVelocity
! outputs: relativeVorticity, circulation
call ocn_relativeVorticity_circulation(relativeVorticity, circulation, normalVelocity, err)
! inputs: relativeVorticity, layerThickEdgeFlux, normalVelocity, normalTransportVelocity,
! normalGMBolusVelocity
! outputs: vertTransportVelocityTop, vertGMBolusVelocityTop, relativeVorticityCell,
! divergence, kineticEnergyCell, tangentialVelocity, vertMLEBolusVelocityTop
! in and out: vertVelocityTop
call ocn_diagnostic_solve_vortVel(relativeVorticity, &
layerThickEdgeFlux, normalVelocity, normalTransportVelocity, &
vertVelocityTop, vertTransportVelocityTop, relativeVorticityCell, &
divergence, kineticEnergyCell, tangentialVelocity)
!compute vertical velocity due to gm bolus velocity if requested
if (config_use_gm) then
call ocn_diagnostic_solve_GMvel(layerThickEdgeFlux, normalGMBolusVelocity, vertGMBolusVelocityTop)
end if
!compute vertical velocity due to submesoscale eddies (MLE) velocity if requested
if (config_submesoscale_enable) then
call ocn_diagnostic_solve_MLEvel(layerThickEdgeFlux, normalMLEVelocity, vertMLEBolusVelocityTop)
end if
if ( configVertAdvMethod == vertAdvRemap ) then
! Overwrite vertVelocityTop with exact Eulerian solution from remapping
call mpas_pool_get_array(statePool, 'layerThicknessLag', &
layerThicknessLag, timeLevel)
call ocn_diagnostic_solve_vertVel_remap(layerThickness, &
layerThicknessLag, dt, vertVelocityTop)
! Overwrite vertical tracer budget terms with new vertVelocityTop
if ( config_compute_active_tracer_budgets ) then
numTracers = size(activeTracers, dim=1)
do iTracer = 1, numTracers
!$omp parallel
!$omp do schedule(runtime) private(k,kmin,kmax)
do iCell = 1, nCells
kmin = minLevelCell(iCell)
kmax = maxLevelCell(iCell)
activeTracerVerticalAdvectionTopFlux(iTracer, :, iCell) = 0.0_RKIND
do k = kmin+1, kmax
activeTracerVerticalAdvectionTopFlux(iTracer, k, iCell) = &
min(0.0_RKIND,vertVelocityTop(k, iCell))* &
activeTracers(iTracer, k-1, iCell) &
+ max(0.0_RKIND,vertVelocityTop(k,iCell))* &
activeTracers(iTracer, k, iCell)
enddo
do k = kmin, kmax-1
activeTracerVerticalAdvectionTendency(iTracer, k, iCell) = ( &
activeTracerVerticalAdvectionTopFlux(iTracer, k+1, iCell) &
- activeTracerVerticalAdvectionTopFlux(iTracer, k, iCell) &
) / layerThickness(k, iCell)
enddo
activeTracerVerticalAdvectionTendency(iTracer, kmax, iCell) = &
- activeTracerVerticalAdvectionTopFlux(iTracer, kmax, iCell) &
/ layerThickness(kmax, iCell)
enddo
!$omp end do
!$omp end parallel
enddo
endif
endif
!
! Compute kinetic energy
!
! inputs :: normalVelocity
! outputs :: kineticEnergyCell
if ( ke_vertex_flag == 1 ) then
call ocn_diagnostic_solve_kineticEnergy(normalVelocity, kineticEnergyCell)
end if
! inputs :: normalVelocity, tangentialVelocity, layerThickness, relativeVorticity
! outputs :: normRelVortEdge, normPlanetVortEdge,
! normalizedRelativeVorticityCell
call ocn_diagnostic_solve_vorticity(dt, normalVelocity, tangentialVelocity, &
layerThickness, relativeVorticity, &
normRelVortEdge, normPlanetVortEdge, &
normalizedRelativeVorticityCell)
!
! Surface pressure
! This section must be placed in the code before computing the density.
!
! inputs: atmosphericPressure, seaIcePressure, frazilSurfacePressure, landIcePressure
! outputs: surfacePressure
call ocn_diagnostic_solve_surface_pressure(forcingPool, atmosphericPressure, &
seaIcePressure, surfacePressure)
!
! equation of state
!
! Only need densities on 0, 1, and 2 halo cells
nCells = nCellsAll
! inputs: activeTracers, tracersSurfaceValue,maxLevelCell,refBottomDepth
! outputs: density, potentialDensity, displacedDensity, thermalExpansionCoeff, salineContractionCoeff
! compute in-place density
call ocn_equation_of_state_density(statePool, meshPool, activeTracers, tracersSurfaceValue, &
nCells, 0, 'relative', density, err, &
thermExpCoeff, salineContractCoeff, &
timeLevelIn=timeLevel)
! compute potentialDensity, the density displaced adiabatically to the mid-depth of top layer.
call ocn_equation_of_state_density(statePool, meshPool, activeTracers, tracersSurfaceValue, &
nCells, 1, 'absolute', potentialDensity, &
err, timeLevelIn=timeLevel)
! compute displacedDensity, density displaced adiabatically to the mid-depth one layer deeper.
! That is, layer k has been displaced to the depth of layer k+1.
call ocn_equation_of_state_density(statePool, meshPool, activeTracers, tracersSurfaceValue, &
nCells, 1, 'relative', displacedDensity, &
err, timeLevelIn=timeLevel)
!
! Pressure
! This section must be placed in the code after computing the density.
!
! inputs: layerThickness, density, surfacePressure
! outputs: montgomeryPotential, pressure
call ocn_diagnostic_solve_pressure(layerThickness, density, &
surfacePressure, montgomeryPotential, pressure)
nCells = nCellsAll
if ( full_compute ) then
!
! Brunt-Vaisala frequency and gradient Richardson number
!
! inputs: displacedDensity, density, zMid, normalVelocity, edgeAreaFractionOfCell
! outputs: BruntVaisalaFreqTop, RiTopOfCell
call ocn_diagnostic_solve_richardson(displacedDensity, density, zMid, &
normalVelocity, edgeAreaFractionOfCell, BruntVaisalaFreqTop, RiTopOfCell)
end if ! full_compute
!
! extrapolate tracer values to ocean surface
! this eventually be a modelled process
! at present, just copy k=1 tracer values onto surface values
! field will be updated below is better approximations are available
!TDR need to consider how to handle tracersSurfaceValues
#ifdef MPAS_OPENACC
! tracersSurfaceValue calc
! inputs: activeTracers, minLevelCell
! outputs: tracersSurfaceValue
!$acc kernels present(tracersSurfaceValue, activeTracers, minLevelCell)
#else
!$omp parallel
!$omp do schedule(runtime)
#endif
#ifdef MPAS_OPENACC
!$acc loop independent
#endif
do iCell = 1, nCells
do iTracer = 1,size(tracersSurfaceValue,1)
tracersSurfaceValue(iTracer, iCell) = activeTracers(iTracer, minLevelCell(iCell), iCell)
enddo
end do
#ifdef MPAS_OPENACC
!$acc end kernels
#else
!$omp end do
#endif
#ifdef MPAS_OPENACC
! normalVelocitySurfaceLayer calc
! inputs: normalVelocity, minLevelEdgeBot
! outputs: normalVelocitySurfaceLayer
!$acc parallel loop gang vector &
!$acc present(normalVelocity, normalVelocitySurfaceLayer, minLevelEdgeBot)
#else
!$omp do schedule(runtime)
#endif
do iEdge = 1, nEdges
normalVelocitySurfaceLayer(iEdge) = normalVelocity(minLevelEdgeBot(iEdge), iEdge)
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
if ( full_compute ) then
!
! average tracer values over the ocean surface layer
! the ocean surface layer is generally assumed to be about 0.1 of the boundary layer depth
!
! inputs: layerThickness, activeTracers, layerThicknessEdgeFlux, normalVelocity
! outputs: tracersSurfaceLayerValue, indexSurfaceLayerDepth, normalVelocitySurfaceLayer,
! surfaceFluxAttenuationCoefficient, surfaceFluxAttenuationCoefficientRunoff
call ocn_diagnostic_solve_surfaceLayer(layerThickness, activeTracers, &
layerThickEdgeFlux, normalVelocity, tracersSurfaceLayerValue, &
indexSurfaceLayerDepth, normalVelocitySurfaceLayer, &
sfcFlxAttCoeff, surfaceFluxAttenuationCoefficientRunoff)
if (trim(config_subglacial_runoff_mode) == 'data') then
! outputs: surfaceFluxAttenuationCoefficientSubglacialRunoff
call ocn_diagnostic_solve_surfaceLayer_subglacialRunoff(surfaceFluxAttenuationCoefficientSubglacialRunoff)
end if
end if ! full_compute
!
! compute fields needed to compute land-ice fluxes, either in the ocean model or in the coupler
!
! inputs: layerThickness, normalVelocity, landIceFraction, landIceFloatingMask
! outputs:
if (landIcePressureOn) then
call ocn_compute_land_ice_flux_input_fields(layerThickness, normalVelocity, activeTracers, &
landIceFraction, landIceFloatingMask, timeLevel, indexTemperature, indexSalinity)
end if
if ( full_compute ) then
!
! inputs: ssh, seaIcePressure, landIceDraft
! outputs: pressureAdjustedSSH, gradSSH
if (landIcePressureOn) then
call ocn_diagnostic_solve_ssh(forcingPool, ssh, seaIcePressure, &
pressureAdjustedSSH, gradSSH, landIceDraft=landIceDraft )
else
call ocn_diagnostic_solve_ssh(forcingPool, ssh, seaIcePressure, &
pressureAdjustedSSH, gradSSH)
end if
end if ! full_compute
#ifdef MPAS_OPENACC
!! Listing outputs by first routine that modifies variable
!! ocn_diagnostic_solve_layerThicknessEdge :: layerThickEdgeDrag, (diagnostics array)
!! layerThickEdgeFlux, (diagnostics array)
!! layerThickEdgeMean (diagnostics array)
!! ocn_relativeVorticity_circulation :: relativeVorticity, (diagnostics array)
!! circulation (diagnostics array)
!! ocn_diagnostic_solve_vortVel :: vertTransportVelocityTop, (diagnostics array)
!! vertGMBolusVelocityTop, (diagnostics array)
!! relativeVorticityCell, (diagnostics array)
!! divergence, (diagnostics array)
!! kineticEnergyCell, (diagnostics array)
!! tangentialVelocity (diagnostics array)
!! vertVelocityTop (diagnostics array)
!! ocn_diagnostic_solve_kineticEnergy ::
!! ocn_diagnostic_solve_vorticity :: normRelVortEdge, (diagnostics array)
!! normPlanetVortEdge, (diagnostics array)
!! normalizedRelativeVorticityCell (diagnostics array)
!! ocn_diagnostic_solve_surface_pressure :: surfacePressure (diagnostics array)
!! ocn_diagnostic_solve_z_coordinates :: zMid, (diagnostics array)
!! zTop,(diagnostics array)
!! ssh (not on device)
!! ocn_equation_of_state_density :: density, (diagnostics array)
!! potentialDensity, (diagnostics array)
!! displacedDensity, (diagnostics array)
!! thermalExpansionCoeff, (diagnostics array)
!! salineContractionCoeff (diagnostics array)
!! ocn_diagnostic_solve_pressure :: montgomeryPotential, (diagnostics array)
!! pressure (diagnostics array)
!! ocn_diagnostic_solve_richardson :: BruntVaisalaFreqTop, (diagnostics array)
!! RiTopOfCell (diagnostics array)
!! tracersSurfaceValue calc :: tracersSurfaceValue (diagnostics array)
!! normalVelocitySurfaceLayer calc :: normalVelocitySurfaceLayer (diagnostics array)
!! surfaceFluxAttenuationCoefficientRunoff :: tracersSurfaceLayerValue, (diagnostics array)
!! indexSurfaceLayerDepth, (diagnostics array)
!! normalVelocitySurfaceLayer, (diagnostics array)
!! surfaceFluxAttenuationCoefficient, ????
!! surfaceFluxAttenuationCoefficientRunoff (diagnostics array)
!! surfaceFluxAttenuationCoefficientSubglacialRunoff :: tracersSurfaceLayerValue, (diagnostics array)
!! indexSurfaceLayerDepth, (diagnostics array)
!! normalVelocitySurfaceLayer, (diagnostics array)
!! surfaceFluxAttenuationCoefficient, ????
!! surfaceFluxAttenuationCoefficientSubglacialRunoff (diagnostics array)
!! ocn_diagnostic_solve_ssh :: pressureAdjustedSSH, (diagnostics array)
!! gradSSH (diagnostics array)
! !$acc update host(layerThickEdgeDrag, &
! !$acc layerThickEdgeFlux, &
! !$acc layerThickEdgeMean)
! !$acc update host(relativeVorticity, circulation)
! !$acc update host(vertTransportVelocityTop, &
! !$acc vertGMBolusVelocityTop, &
! !$acc relativeVorticityCell, &
! !$acc divergence, &
! !$acc kineticEnergyCell, &
! !$acc tangentialVelocity, &
! !$acc vertVelocityTop)
! !$acc update host(normRelVortEdge, normPlanetVortEdge, &
! !$acc normalizedRelativeVorticityCell)
! !$acc update host (surfacePressure)
! !$acc update host(zMid, zTop)
! !$acc exit data copyout(ssh)
! !$acc update host(density, potentialDensity, displacedDensity)
! !$acc update host(thermExpCoeff, &
! !$acc& salineContractCoeff)
! !$acc update host(tracersSurfaceValue)
! !$acc update host(normalVelocitySurfaceLayer)
! !$acc update host(montgomeryPotential, pressure)
! if ( full_compute ) then
! !$acc update host(RiTopOfCell, &
! !$acc BruntVaisalaFreqTop)
! !$acc update host(tracersSurfaceLayerValue, &
! !$acc indexSurfaceLayerDepth, &
! !$acc normalVelocitySurfaceLayer, &
! !$acc sfcFlxAttCoeff, &
! !$acc surfaceFluxAttenuationCoefficientRunoff)
! if (trim(config_subglacial_runoff_mode) == 'data') then
! !$acc update host(surfaceFluxAttenuationCoefficientSubglacialRunoff)
! end if
! end if ! full_compute
! !$acc exit data delete (atmosphericPressure, seaIcePressure)
! if ( associated(frazilSurfacePressure) ) then
! !$acc exit data delete(frazilSurfacePressure)
! endif
! if (landIcePressureOn) then
! !$acc exit data delete(landIcePressure)
! !$acc exit data delete(landIceDraft)
! endif
! !$acc exit data delete (layerThickness, normalVelocity, &
! !$acc activeTracers)
#endif
call mpas_timer_stop('diagnostic solve')
end subroutine ocn_diagnostic_solve!}}}
!***********************************************************************
!
! routine ocn_relativeVorticity_circulation
!
!> \brief Computes relative vorticity and circulation
!> \author Mark Petersen, Doug Jacobsen, Todd Ringler
!> \date November 2013
!> \details
!> Computes relative vorticity and circulation
!
!-----------------------------------------------------------------------
subroutine ocn_relativeVorticity_circulation(relativeVorticity, circulation, normalVelocity, err)!{{{
!-----------------------------------------------------------------
!
! input variables
!
!-----------------------------------------------------------------
real (kind=RKIND), dimension(:,:), intent(in) :: &
normalVelocity
!-----------------------------------------------------------------
!
! output variables
!
!-----------------------------------------------------------------
real (kind=RKIND), dimension(:,:), intent(out) :: &
relativeVorticity
real (kind=RKIND), dimension(:,:), intent(out) :: &
circulation
integer, intent(out) :: err !< Output: error flag
!-----------------------------------------------------------------
!
! local variables
!
!-----------------------------------------------------------------
integer :: iVertex, iEdge, i, k
real (kind=RKIND) :: invAreaTri1, r_tmp
err = 0
#ifdef MPAS_OPENACC
!$acc parallel loop collapse(2) present(circulation, relativeVorticity)
#else
!$omp parallel
!$omp do schedule(runtime)
#endif
do iVertex = 1, nVerticesAll
do k = 1, nVertLevels
circulation(k, iVertex) = 0.0_RKIND
relativeVorticity(k, iVertex) = 0.0_RKIND
enddo
enddo
#ifndef MPAS_OPENACC
!$omp end do
#endif
#ifdef MPAS_OPENACC
!$acc parallel loop &
!$acc present(circulation, relativeVorticity, areaTriangle, edgesOnVertex, &
!$acc maxLevelVertexBot, dcEdge, normalVelocity, edgeSignOnVertex, &
!$acc minLevelVertexTop) &
!$acc private(invAreaTri1, iVertex, iEdge, k, r_tmp)
#else
!$omp do schedule(runtime) private(invAreaTri1, i, iEdge, k, r_tmp)
#endif
do iVertex = 1, nVerticesAll
invAreaTri1 = 1.0_RKIND / areaTriangle(iVertex)
do i = 1, vertexDegree
iEdge = edgesOnVertex(i, iVertex)
do k = minLevelVertexTop(iVertex), maxLevelVertexBot(iVertex)
r_tmp = dcEdge(iEdge) * normalVelocity(k, iEdge)
circulation(k, iVertex) = circulation(k, iVertex) + edgeSignOnVertex(i, iVertex) * r_tmp
relativeVorticity(k, iVertex) = relativeVorticity(k, iVertex) + edgeSignOnVertex(i, iVertex) * r_tmp * invAreaTri1
end do
end do
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
!--------------------------------------------------------------------
end subroutine ocn_relativeVorticity_circulation!}}}
!***********************************************************************
!
! routine ocn_diagnostic_solve_wctEdge
!
!> \brief Computes layer thickness at edge
!> \author Carolyn Begeman
!> \date July 2023
!> \details
!> This routine computes the water column thickness on edges, wctEdge, which
!> can either be the mean or an upwind value, depending on user input.
!
!-----------------------------------------------------------------------
subroutine ocn_diagnostic_solve_wctEdge(btrVelocity, deltaSSH, &
baroclinicThickness, wctEdge)!{{{
!-----------------------------------------------------------------
! input variables
!-----------------------------------------------------------------
real (kind=RKIND), dimension(:), intent(in) :: &
btrVelocity, &!< barotropic velocity
deltaSSH, &!< ssh change from baroclinic step to barotropic
!< substep
baroclinicThickness !< water column thickness on edges at the last
!< baroclinic step (centered or upwinded)
!-----------------------------------------------------------------
! output variables
!-----------------------------------------------------------------
! Outputs are the shared variables from diagnostic_variables:
real (kind=RKIND), dimension(:), intent(inout) :: &
wctEdge !< water column thickness on edge
!-----------------------------------------------------------------
! local variables
!-----------------------------------------------------------------
integer :: &
iCell, iEdge, &! edge and cell loop indices
cell1, cell2 ! neighbor cell indices across an edge
! End preamble
!-----------------------------------------------------------------
! Begin code
! Compute edge flux based on option set on init
select case (thickEdgeFluxChoice)
case (thickEdgeFluxCenter)
! Use centered (mean) thickness as flux value
! Compute mean layer thickness at edge
#ifdef MPAS_OPENACC
!$acc enter data copyin(baroclinicThickness, deltaSSH, wctEdge)
!$acc parallel loop &
!$acc present(baroclinicThickness, deltaSSH, wctEdge, cellsOnEdge) &
!$acc private(iEdge, cell1, cell2)
#else
!$omp parallel
!$omp do schedule(runtime) private (cell1, cell2)
#endif
do iEdge = 1, nEdgesAll
if ( maxLevelEdgeTop(iEdge) == 0 ) cycle
cell1 = cellsOnEdge(1,iEdge)
cell2 = cellsOnEdge(2,iEdge)
wctEdge(iEdge) = 0.5_RKIND*(deltaSSH(cell1) + deltaSSH(cell2)) &
+ baroclinicThickness(iEdge)
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
#ifdef MPAS_OPENACC
!$acc exit data copyout(wctEdge) delete(baroclinicThickness, deltaSSH)
#endif
case (thickEdgeFluxUpwind)
! Use upwind thickness as the edge flux value
#ifdef MPAS_OPENACC
!$acc enter data copyin(baroclinicThickness, deltaSSH, btrVelocity, wctEdge)
!$acc parallel loop &
!$acc present(baroclinicThickness, btrVelocity, deltaSSH, cellsOnEdge, maxLevelEdgeTop, wctEdge) &
!$acc private(cell1, cell2, iEdge)
#else
!$omp parallel
!$omp do schedule(runtime) private(cell1, cell2)
#endif
do iEdge = 1, nEdgesAll
if ( maxLevelEdgeTop(iEdge) == 0 ) cycle
cell1 = cellsOnEdge(1, iEdge)
cell2 = cellsOnEdge(2, iEdge)
if (btrVelocity(iEdge) > 0.0_RKIND) then
wctEdge(iEdge) = deltaSSH(cell1) + baroclinicThickness(iEdge)
elseif (btrVelocity(iEdge) < 0.0_RKIND) then
wctEdge(iEdge) = deltaSSH(cell2) + baroclinicThickness(iEdge)
else
wctEdge(iEdge) = max(deltaSSH(cell1), deltaSSH(cell2)) &
+ baroclinicThickness(iEdge)
end if
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
#ifdef MPAS_OPENACC
!$acc exit data copyout(wctEdge) delete(baroclinicThickness, deltaSSH, btrVelocity)
#endif
case default
! Should have been caught on init
call mpas_log_write('Thickness flux option unknown', &
MPAS_LOG_CRIT)
end select
end subroutine ocn_diagnostic_solve_wctEdge!}}}
!***********************************************************************
!
! routine ocn_diagnostic_solve_layerThicknessEdge
!
!> \brief Computes layer thickness at edge
!> \author Matt Turner
!> \date October 2020
!> \details
!> This routine computes the diagnostic variables layerThickEdgeDrag,
!> layerThickEdgeFlux, and layerThickEdgeMean. The Mean is a simple
!> mean across the edge but the Flux value can either be the mean or
!> an upwind value, depending on user input.
!
!-----------------------------------------------------------------------
subroutine ocn_diagnostic_solve_layerThicknessEdge(normalVelocity, &
layerThickness, restingThickness, &
ssh)!{{{
!-----------------------------------------------------------------
! input variables
!-----------------------------------------------------------------
real (kind=RKIND), dimension(:,:), intent(in) :: &
normalVelocity, &!< [in] transport
layerThickness, &!< [in] layer thickness at cell center
restingThickness !< [in] initial layer thickness at cell center
real (kind=RKIND), dimension(:), intent(in) :: &
ssh !< [in] sea surface height at cell center
!-----------------------------------------------------------------
! output variables
!-----------------------------------------------------------------
! Outputs are the shared variables from diagnostic_variables:
!real (kind=RKIND), dimension(:,:), intent(out) :: &
! layerThickEdgeDrag, & !< [out] layer thickness on edge for drag
! layerThickEdgeMean, & !< [out] centered layer thickness on edge
! layerThickEdgeFlux !< [out] flux-related thickness on edge
!-----------------------------------------------------------------
! local variables
!-----------------------------------------------------------------
integer :: &
iEdge, k, &! edge and vertical loop indices
cell1, cell2, &! neighbor cell indices across an edge
kmin, kmax ! min,max active vertical levels
! End preamble
!-----------------------------------------------------------------
! Begin code
! Compute mean layer thickness at edge
#ifdef MPAS_OPENACC
!$acc parallel loop &
!$acc present(layerThickness, layerThickEdgeMean, &
!$acc minLevelEdgeBot, maxLevelEdgeTop, cellsOnEdge) &
!$acc private(k, kmin, kmax, cell1, cell2)
#else
!$omp parallel
!$omp do schedule(runtime) private(k, kmin, kmax, cell1, cell2)
#endif
do iEdge = 1, nEdgesAll
kmin = minLevelEdgeBot(iEdge)
kmax = maxLevelEdgeTop(iEdge)
cell1 = cellsOnEdge(1,iEdge)
cell2 = cellsOnEdge(2,iEdge)
do k = 1,nVertLevels
! initialize layerThicknessEdgeMean to avoid divide by
! zero and NaN problems.
layerThickEdgeMean(k,iEdge) = -1.0e34_RKIND
end do
do k = kmin,kmax
! central differenced
layerThickEdgeMean(k,iEdge) = 0.5_RKIND * &
(layerThickness(k,cell1) + &
layerThickness(k,cell2))
end do
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
! Compute edge flux based on option set on init
select case (thickEdgeFluxChoice)
case (thickEdgeFluxCenter)
! Use centered (mean) thickness as flux value
if ( config_use_subgrid_wetting_drying ) then
#ifdef MPAS_OPENACC
!$acc update host(ssh)
#endif
call ocn_subgrid_layerThickEdgeFlux_center(ssh, layerThickEdgeFlux)
#ifdef MPAS_OPENACC
!$acc update device(layerThickEdgeFlux)
#endif
else
#ifdef MPAS_OPENACC
!$acc parallel loop collapse(2) &
!$acc present(layerThickEdgeFlux, layerThickEdgeMean)
#else
!$omp parallel
!$omp do schedule(runtime) private(k)
#endif
do iEdge = 1, nEdgesAll
do k = 1,nVertLevels
layerThickEdgeFlux(k,iEdge) = &
layerThickEdgeMean(k,iEdge)
end do
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
end if
case (thickEdgeFluxUpwind)
! Use upwind thickness as the edge flux value
if (config_use_subgrid_wetting_drying) then
#ifdef MPAS_OPENACC
!$acc update host(ssh, normalVelocity, layerThickness)
#endif
call ocn_subgrid_layerThickEdgeFlux_upwind(ssh, normalVelocity, layerThickness, layerThickEdgeFlux)
#ifdef MPAS_OPENACC
!$acc update device(layerThickEdgeFlux)
#endif
else
#ifdef MPAS_OPENACC
!$acc parallel loop &
!$acc present(normalVelocity, layerThickness, &
!$acc minLevelEdgeBot, maxLevelEdgeTop, &
!$acc layerThickEdgeFlux, cellsOnEdge) &
!$acc private(k, kmin, kmax, cell1, cell2)
#else
!$omp parallel
!$omp do schedule(runtime) private(k, kmin, kmax, cell1, cell2)
#endif
do iEdge = 1, nEdgesAll
kmin = minLevelEdgeBot(iEdge)
kmax = maxLevelEdgeTop(iEdge)
cell1 = cellsOnEdge(1,iEdge)
cell2 = cellsOnEdge(2,iEdge)
do k=1,nVertLevels
! initialize layerThicknessEdgeFlux to avoid divide by
! zero and NaN problems.
layerThickEdgeFlux(k,iEdge) = -1.0e34_RKIND
end do
do k = kmin,kmax
if (normalVelocity(k,iEdge) > 0.0_RKIND) then
layerThickEdgeFlux(k,iEdge) = layerThickness(k,cell1)
elseif (normalVelocity(k,iEdge) < 0.0_RKIND) then
layerThickEdgeFlux(k,iEdge) = layerThickness(k,cell2)
else
layerThickEdgeFlux(k,iEdge) = &
max(layerThickness(k,cell1), &
layerThickness(k,cell2))
end if
end do
end do
#ifndef MPAS_OPENACC
!$omp end do
!$omp end parallel
#endif
end if
case (thickEdgeFluxConstant)
! Use linearized version H*u where H is constant in time
#ifdef MPAS_OPENACC
!$acc parallel loop &
!$acc present(restingThickness, &
!$acc minLevelEdgeBot, maxLevelEdgeTop, &
!$acc layerThickEdgeFlux, cellsOnEdge) &
!$acc private(k, kmin, kmax, cell1, cell2)
#else
!$omp parallel
!$omp do schedule(runtime) private(k, kmin, kmax, cell1, cell2)