-
Notifications
You must be signed in to change notification settings - Fork 122
/
Copy pathatmosphere.F90
2543 lines (2229 loc) · 103 KB
/
atmosphere.F90
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
!***********************************************************************
!* GNU Lesser General Public License
!*
!* This file is part of the FV3 dynamical core.
!*
!* The FV3 dynamical core is free software: you can redistribute it
!* and/or modify it under the terms of the
!* GNU Lesser General Public License as published by the
!* Free Software Foundation, either version 3 of the License, or
!* (at your option) any later version.
!*
!* The FV3 dynamical core is distributed in the hope that it will be
!* useful, but WITHOUT ANYWARRANTY; without even the implied warranty
!* of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
!* See the GNU General Public License for more details.
!*
!* You should have received a copy of the GNU Lesser General Public
!* License along with the FV3 dynamical core.
!* If not, see <http://www.gnu.org/licenses/>.
!***********************************************************************
!>@brief The module 'atmosphere' provides the interface for the
!! Cubed-Sphere FV dynamical core
module atmosphere_mod
! <table>
! <tr>
! <th>Module Name</th>
! <th>Functions Included</th>
! </tr>
! <tr>
! <td>block_control_mod</td>
! <td>block_control_type</td>
! </tr>
! <tr>
! <td>constants_mod</td>
! <td>cp_air, rdgas, grav, rvgas, kappa, pstd_mks</td>
! </tr>
! <tr>
! <td>field_manager_mod</td>
! <td>MODEL_ATMOS</td>
! </tr>
! <tr>
! <td>fms_mod</td>
! <td>file_exist, error_mesg, FATAL, check_nml_error,
! stdlog,write_version_number, mpp_clock_id,
! mpp_clock_begin, mpp_clock_end, CLOCK_SUBCOMPONENT,
! clock_flag_default</td>
! </tr>
! <tr>
! <td>fms2_io_mod</td>
! <td>file_exists</td>
! <tr>
! <td>fv_arrays_mod</td>
! <td>fv_atmos_type, R_GRID</td>
! </tr>
! <tr>
! <td>fv_control_mod</td>
! <td>fv_init, fv_end, ngrids</td>
! </tr>
! <tr>
! <td>fv_diagnostics_mod</td>
! <td>fv_diag_init, fv_diag, fv_time, prt_maxmin, prt_height</td>
! </tr>
! <tr>
! <td>fv_dynamics_mod</td>
! <td>fv_dynamics</td>
! </tr>
! <tr>
! <td>fv_eta_mod</td>
! <td>get_eta_level</td>
! </tr>
! <tr>
! <td>fv_fill_mod</td>
! <td>fill_gfs</td>
! </tr>
! <tr>
! <td>fv_mp_mod</td>
! <td>switch_current_Atm</td>
! </tr>
! <tr>
! <td>fv_nesting_mod</td>
! <td>twoway_nesting</td>
! </tr>
! <tr>
! <td>fv_nggps_diags_mod</td>
! <td>fv_nggps_diag_init, fv_nggps_diag</td>
! </tr>
! <tr>
! <td>fv_nwp_nudge_mod</td>
! <td>fv_nwp_nudge_init, fv_nwp_nudge_end, do_adiabatic_init</td>
! </tr>
! <tr>
! <td>fv_restart_mod</td>
! <td>fv_restart, fv_write_restart</td>
! </tr>
! <tr>
! <td>fv_sg_mod</td>
! <td>fv_subgrid_z</td>
! </tr>
! <tr>
! <td>fv_timing_mod</td>
! <td>timing_on, timing_off</td>
! </tr>
! <tr>
! <td>fv_update_phys_mod</td>
! <td>fv_update_phys</td>
! </tr>
! <tr>
! <td>IPD_typedefs_mod</td>
! <td>IPD_data_type, kind_phys => IPD_kind_phys</td>
! </tr>
! <tr>
! <td>mpp_mod</td>
! <td>mpp_error, stdout, FATAL, NOTE, input_nml_file, mpp_root_pe,
! mpp_npes, mpp_pe, mpp_chksum,mpp_get_current_pelist,
! mpp_set_current_pelist</td>
! </tr>
! <tr>
! <td>mpp_domains_mod</td>
! <td>mpp_get_data_domain, mpp_get_compute_domain, domain2d, mpp_update_domains</td>
! </tr>
! <tr>
! <td>mpp_parameter_mod</td>
! <td>EUPDATE, WUPDATE, SUPDATE, NUPDATE</td>
! </tr>
! <tr>
! <td>time_manager_mod</td>
! <td>time_type, get_time, set_time, operator(+), operator(-)</td>
! </tr>
! <tr>
! <td>tracer_manager_mod</td>
! <td>get_tracer_index, get_number_tracers, NO_TRACER</td>
! </tr>
! <tr>
! <td>xgrid_mod</td>
! <td>grid_box_type</td>
! </tr>
! </table>
#include <fms_platform.h>
!-----------------
! FMS modules:
!-----------------
use block_control_mod, only: block_control_type
#ifdef OVERLOAD_R4
use constantsR4_mod, only: cp_air, rdgas, grav, rvgas, kappa, pstd_mks
#else
use constants_mod, only: cp_air, rdgas, grav, rvgas, kappa, pstd_mks
#endif
use time_manager_mod, only: time_type, get_time, set_time, operator(+), &
operator(-), operator(/), time_type_to_real
use fms_mod, only: error_mesg, FATAL, &
check_nml_error, stdlog, &
write_version_number, &
mpp_clock_id, mpp_clock_begin, &
mpp_clock_end, CLOCK_SUBCOMPONENT, &
clock_flag_default
use fms2_io_mod, only: file_exists
use mpp_mod, only: mpp_error, stdout, FATAL, WARNING, NOTE, &
input_nml_file, mpp_root_pe, &
mpp_npes, mpp_pe, mpp_chksum, &
mpp_get_current_pelist, &
mpp_set_current_pelist, &
mpp_sync, mpp_sync_self, mpp_send, mpp_recv, mpp_broadcast
use mpp_parameter_mod, only: EUPDATE, WUPDATE, SUPDATE, NUPDATE
use mpp_domains_mod, only: CENTER, CORNER, NORTH, EAST, WEST, SOUTH
use mpp_domains_mod, only: domain2d, mpp_update_domains, mpp_global_field
use mpp_domains_mod, only: mpp_get_data_domain, mpp_get_compute_domain, mpp_get_global_domain
use xgrid_mod, only: grid_box_type
use field_manager_mod, only: MODEL_ATMOS
use tracer_manager_mod, only: get_tracer_index, get_number_tracers, &
NO_TRACER, get_tracer_names
use DYCORE_typedefs, only: DYCORE_data_type
#ifdef GFS_TYPES
use GFS_typedefs, only: IPD_data_type => GFS_data_type, IPD_control_type => GFS_control_type, kind_phys
#else
use IPD_typedefs, only: IPD_data_type, IPD_control_type, kind_phys => IPD_kind_phys
#endif
use fv_iau_mod, only: IAU_external_data_type
#ifdef MULTI_GASES
use multi_gases_mod, only: virq, virq_max, num_gas, ri, cpi
#endif
!-----------------
! FV core modules:
!-----------------
use fv_arrays_mod, only: fv_atmos_type, R_GRID, fv_grid_bounds_type, phys_diag_type
use fv_control_mod, only: fv_control_init, fv_end, ngrids
use fv_eta_mod, only: get_eta_level
use fv_fill_mod, only: fill_gfs
use fv_dynamics_mod, only: fv_dynamics
use fv_nesting_mod, only: twoway_nesting
use boundary_mod, only: fill_nested_grid
use fv_diagnostics_mod, only: fv_diag_init, fv_diag, fv_time, prt_maxmin, prt_height
use fv_nggps_diags_mod, only: fv_nggps_diag_init, fv_nggps_diag, fv_nggps_tavg
use fv_restart_mod, only: fv_restart, fv_write_restart
use fv_timing_mod, only: timing_on, timing_off
use fv_mp_mod, only: is_master, tile_fine
use fv_sg_mod, only: fv_subgrid_z
use fv_update_phys_mod, only: fv_update_phys
use fv_io_mod, only: fv_io_register_nudge_restart
use fv_nwp_nudge_mod, only: fv_nwp_nudge_init, fv_nwp_nudge_end, do_adiabatic_init
use fv_regional_mod, only: start_regional_restart, read_new_bc_data, &
a_step, p_step, current_time_in_seconds
use fv_grid_utils_mod, only: g_sum
use coarse_graining_mod, only: coarse_graining_init
use coarse_grained_diagnostics_mod, only: fv_coarse_diag_init, fv_coarse_diag
use coarse_grained_restart_files_mod, only: fv_coarse_restart_init
use diag_manager_mod, only: send_data
implicit none
private
!--- driver routines
public :: atmosphere_init, atmosphere_end, atmosphere_restart, &
atmosphere_dynamics, atmosphere_state_update
!--- utility routines
public :: atmosphere_resolution, atmosphere_grid_bdry, &
atmosphere_grid_ctr, atmosphere_domain, &
atmosphere_control_data, atmosphere_pref, &
atmosphere_diag_axes, atmosphere_etalvls, &
atmosphere_hgt, atmosphere_scalar_field_halo, &
! experimental APIs not ready for use
! atmosphere_tracer_postinit, &
atmosphere_diss_est, & ! dissipation estimate for SKEB
atmosphere_get_bottom_layer, &
atmosphere_nggps_diag, &
get_bottom_mass, get_bottom_wind, &
get_stock_pe, set_atmosphere_pelist, &
get_nth_domain_info
!--- physics/radiation data exchange routines
public :: atmos_phys_driver_statein
!--- coupling data exchange routines
public :: atmosphere_fill_nest_cpl
!-----------------------------------------------------------------------
! version number of this module
! Include variable "version" to be written to log file.
#include<file_version.h>
character(len=20) :: mod_name = 'fvGFS/atmosphere_mod'
!---- private data ----
type (time_type) :: Time_step_atmos
public Atm, mygrid, p_split, dt_atmos ! Share over to moving nest functions.
!These are convenience variables for local use only, and are set to values in Atm%
real :: dt_atmos
real :: zvir
integer :: npx, npy, npz, ncnst, pnats
integer :: isc, iec, jsc, jec
integer :: isd, ied, jsd, jed
integer :: nq ! number of transported tracers
integer :: sec, seconds, days
integer :: id_dynam, id_fv_diag, id_subgridz
logical :: cold_start = .false. ! used in initial condition
integer, dimension(:), allocatable :: id_tracerdt_dyn
integer :: sphum, liq_wat, rainwat, ice_wat, snowwat, graupel, hailwat, cld_amt ! condensate species tracer indices
integer :: mygrid = 1
integer :: p_split = 1
integer, allocatable :: pelist(:)
logical, allocatable :: grids_on_this_pe(:)
type(fv_atmos_type), allocatable, target :: Atm(:)
integer :: id_udt_dyn, id_vdt_dyn
real, parameter:: w0_big = 60. ! to prevent negative w-tracer diffusion
!---dynamics tendencies for use in fv_subgrid_z and during fv_update_phys
real, allocatable, dimension(:,:,:) :: u_dt, v_dt, t_dt, qv_dt
real, allocatable :: pref(:,:), dum1d(:)
logical :: first_diag = .true.
contains
!>@brief The subroutine 'atmosphere_init' is an API to initialize the FV3 dynamical core,
!! including the grid structures, memory, initial state (self-initialization or restart),
!! and diagnostics.
subroutine atmosphere_init (Time_init, Time, Time_step, Grid_box, area)
use ccpp_static_api, only: ccpp_physics_init
use CCPP_data, only: ccpp_suite, &
cdata => cdata_tile, &
GFDL_interstitial
#ifdef OPENMP
use omp_lib
#endif
type (time_type), intent(in) :: Time_init, Time, Time_step
type(grid_box_type), intent(inout) :: Grid_box
real(kind=kind_phys), pointer, dimension(:,:), intent(inout) :: area
!--- local variables ---
integer :: i, n
! integer :: itrac
logical :: do_atmos_nudge
character(len=32) :: tracer_name, tracer_units
real :: ps1, ps2
integer :: nthreads, ierr
integer :: nlunit = 9999
character (len = 64) :: fn_nml = 'input.nml'
! DH* 20210326
! This is a temporary workaround until the implementation
! of the generic interface to call GFDL or CCPP physics is
! completed. We need the name of the CCPP suite here in order
! to run the adiabatic init with fast physics turned on. All
! other vaiables are ignored (set to the same default value
! as in fv3atm's atmos_model.F90.
integer :: io
integer :: blocksize = 1
logical :: chksum_debug = .false.
logical :: dycore_only = .false.
logical :: debug = .false.
logical :: sync = .false.
logical :: ignore_rst_cksum = .false.
real :: avg_max_length = 3600.
namelist /atmos_model_nml/ blocksize, chksum_debug, dycore_only, debug, sync, ccpp_suite, avg_max_length, &
ignore_rst_cksum
! *DH 20210326
!For regional
a_step = 0
current_time_in_seconds = time_type_to_real( Time - Time_init )
if (mpp_pe() == 0) write(*,"('atmosphere_init: current_time_seconds = ',f9.1)")current_time_in_seconds
call timing_on('ATMOS_INIT')
allocate(pelist(mpp_npes()))
call mpp_get_current_pelist(pelist)
zvir = rvgas/rdgas - 1.
!---- compute physics/atmos time step in seconds ----
Time_step_atmos = Time_step
call get_time (Time_step_atmos, sec)
dt_atmos = real(sec)
!----- initialize FV dynamical core -----
!NOTE do we still need the second file_exist call?
cold_start = (.not.file_exists('INPUT/fv_core.res.nc') .and. .not.file_exists('INPUT/fv_core.res.tile1.nc'))
call fv_control_init( Atm, dt_atmos, mygrid, grids_on_this_pe, p_split ) ! allocates Atm components; sets mygrid
Atm(mygrid)%Time_init = Time_init
if(Atm(mygrid)%flagstruct%warm_start) then
a_step = nint(current_time_in_seconds/dt_atmos)
if (Atm(mygrid)%coarse_graining%write_coarse_restart_files .or. &
Atm(mygrid)%coarse_graining%write_coarse_diagnostics) then
call coarse_graining_init(Atm(mygrid)%flagstruct%npx, Atm(mygrid)%npz, &
Atm(mygrid)%layout, Atm(mygrid)%bd%is, Atm(mygrid)%bd%ie, &
Atm(mygrid)%bd%js, Atm(mygrid)%bd%je, Atm(mygrid)%coarse_graining%factor, &
Atm(mygrid)%coarse_graining%nx_coarse, &
Atm(mygrid)%coarse_graining%strategy, &
Atm(mygrid)%coarse_graining%domain)
endif
endif
!----- write version and namelist to log file -----
call write_version_number ( mod_name, version )
!-----------------------------------
npx = Atm(mygrid)%npx
npy = Atm(mygrid)%npy
npz = Atm(mygrid)%npz
ncnst = Atm(mygrid)%ncnst
pnats = Atm(mygrid)%flagstruct%pnats
isc = Atm(mygrid)%bd%isc
iec = Atm(mygrid)%bd%iec
jsc = Atm(mygrid)%bd%jsc
jec = Atm(mygrid)%bd%jec
isd = isc - Atm(mygrid)%bd%ng
ied = iec + Atm(mygrid)%bd%ng
jsd = jsc - Atm(mygrid)%bd%ng
jed = jec + Atm(mygrid)%bd%ng
nq = ncnst-pnats
sphum = get_tracer_index (MODEL_ATMOS, 'sphum' )
liq_wat = get_tracer_index (MODEL_ATMOS, 'liq_wat' )
ice_wat = get_tracer_index (MODEL_ATMOS, 'ice_wat' )
rainwat = get_tracer_index (MODEL_ATMOS, 'rainwat' )
snowwat = get_tracer_index (MODEL_ATMOS, 'snowwat' )
graupel = get_tracer_index (MODEL_ATMOS, 'graupel' )
hailwat = get_tracer_index (MODEL_ATMOS, 'hailwat' )
cld_amt = get_tracer_index (MODEL_ATMOS, 'cld_amt')
if (max(sphum,liq_wat,ice_wat,rainwat,snowwat,graupel,hailwat) > Atm(mygrid)%flagstruct%nwat) then
call mpp_error (FATAL,' atmosphere_init: condensate species are not first in the list of &
&tracers defined in the field_table')
endif
! Allocate grid variables to be used to calculate gradient in 2nd order flux exchange
! This data is only needed for the COARSEST grid.
!call switch_current_Atm(Atm(mygrid))
allocate(Grid_box%dx ( isc:iec , jsc:jec+1))
allocate(Grid_box%dy ( isc:iec+1, jsc:jec ))
allocate(Grid_box%area ( isc:iec , jsc:jec ))
allocate(Grid_box%edge_w( jsc:jec+1))
allocate(Grid_box%edge_e( jsc:jec+1))
allocate(Grid_box%edge_s( isc:iec+1 ))
allocate(Grid_box%edge_n( isc:iec+1 ))
allocate(Grid_box%en1 (3, isc:iec , jsc:jec+1))
allocate(Grid_box%en2 (3, isc:iec+1, jsc:jec ))
allocate(Grid_box%vlon (3, isc:iec , jsc:jec ))
allocate(Grid_box%vlat (3, isc:iec , jsc:jec ))
Grid_box%dx ( isc:iec , jsc:jec+1) = Atm(mygrid)%gridstruct%dx ( isc:iec, jsc:jec+1)
Grid_box%dy ( isc:iec+1, jsc:jec ) = Atm(mygrid)%gridstruct%dy ( isc:iec+1, jsc:jec )
Grid_box%area ( isc:iec , jsc:jec ) = Atm(mygrid)%gridstruct%area ( isc:iec , jsc:jec )
Grid_box%edge_w( jsc:jec+1) = Atm(mygrid)%gridstruct%edge_w( jsc:jec+1)
Grid_box%edge_e( jsc:jec+1) = Atm(mygrid)%gridstruct%edge_e( jsc:jec+1)
Grid_box%edge_s( isc:iec+1 ) = Atm(mygrid)%gridstruct%edge_s( isc:iec+1)
Grid_box%edge_n( isc:iec+1 ) = Atm(mygrid)%gridstruct%edge_n( isc:iec+1)
Grid_box%en1 (:, isc:iec , jsc:jec+1) = Atm(mygrid)%gridstruct%en1 (:, isc:iec , jsc:jec+1)
Grid_box%en2 (:, isc:iec+1, jsc:jec ) = Atm(mygrid)%gridstruct%en2 (:, isc:iec+1, jsc:jec )
do i = 1,3
Grid_box%vlon (i, isc:iec , jsc:jec ) = Atm(mygrid)%gridstruct%vlon (isc:iec , jsc:jec, i )
Grid_box%vlat (i, isc:iec , jsc:jec ) = Atm(mygrid)%gridstruct%vlat (isc:iec , jsc:jec, i )
enddo
allocate (area(isc:iec , jsc:jec ))
area(isc:iec,jsc:jec) = Atm(mygrid)%gridstruct%area_64(isc:iec,jsc:jec)
!----- allocate and zero out the dynamics (and accumulated) tendencies
allocate( u_dt(isd:ied,jsd:jed,npz), &
v_dt(isd:ied,jsd:jed,npz), &
t_dt(isc:iec,jsc:jec,npz), &
qv_dt(isc:iec,jsc:jec,npz) )
!--- allocate pref
allocate(pref(npz+1,2), dum1d(npz+1))
! DH* 20210326
! First, read atmos_model_nml namelist section - this is a workaround to avoid
! unnecessary additional changes to the input namelists, in anticipation of the
! implementation of a generic interface for GFDL and CCPP fast physics soon
read(input_nml_file, nml=atmos_model_nml, iostat=io)
ierr = check_nml_error(io, 'atmos_model_nml')
!write(0,'(a)') "It's me, and my physics suite is '" // trim(ccpp_suite) // "'"
! *DH 20210326
call fv_restart(Atm(mygrid)%domain, Atm, seconds, days, cold_start, Atm(mygrid)%gridstruct%grid_type, mygrid)
fv_time = Time
!----- initialize atmos_axes and fv_dynamics diagnostics
!I've had trouble getting this to work with multiple grids at a time; worth revisiting?
call fv_diag_init(Atm(mygrid:mygrid), Atm(mygrid)%atmos_axes, Time, npx, npy, npz, Atm(mygrid)%flagstruct%p_ref)
if (Atm(mygrid)%coarse_graining%write_coarse_diagnostics) then
call fv_coarse_diag_init(Atm, Time, Atm(mygrid)%atmos_axes(3), &
Atm(mygrid)%atmos_axes(4), Atm(mygrid)%coarse_graining)
endif
if (Atm(mygrid)%coarse_graining%write_coarse_restart_files) then
call fv_coarse_restart_init(Atm(mygrid)%npz, Atm(mygrid)%flagstruct%nt_prog, &
Atm(mygrid)%flagstruct%nt_phys, Atm(mygrid)%flagstruct%hydrostatic, &
Atm(mygrid)%flagstruct%hybrid_z, Atm(mygrid)%flagstruct%fv_land, &
Atm(mygrid)%coarse_graining%write_coarse_dgrid_vel_rst, &
Atm(mygrid)%coarse_graining%write_coarse_agrid_vel_rst, &
Atm(mygrid)%coarse_graining%restart)
endif
!---------- reference profile -----------
ps1 = 101325.
ps2 = 81060.
pref(npz+1,1) = ps1
pref(npz+1,2) = ps2
call get_eta_level ( npz, ps1, pref(1,1), dum1d, Atm(mygrid)%ak, Atm(mygrid)%bk )
call get_eta_level ( npz, ps2, pref(1,2), dum1d, Atm(mygrid)%ak, Atm(mygrid)%bk )
! --- initialize clocks for dynamics, physics_down and physics_up
id_dynam = mpp_clock_id ('FV dy-core', flags = clock_flag_default, grain=CLOCK_SUBCOMPONENT )
id_subgridz = mpp_clock_id ('FV subgrid_z',flags = clock_flag_default, grain=CLOCK_SUBCOMPONENT )
id_fv_diag = mpp_clock_id ('FV Diag', flags = clock_flag_default, grain=CLOCK_SUBCOMPONENT )
call timing_off('ATMOS_INIT')
! Do CCPP fast physics initialization before call to adiabatic_init (since this calls fv_dynamics)
! For fast physics running over the entire domain, block
! and thread number are not used; set to safe values
cdata%blk_no = 1
cdata%thrd_no = 1
cdata%thrd_cnt = 1
! Create shared data type for fast and slow physics, one for each thread
#ifdef OPENMP
nthreads = omp_get_max_threads()
#else
nthreads = 1
#endif
! Create interstitial data type for fast physics; for multi-gases physics,
! pass q(:,:,:,1:num_gas) as qvi, otherwise pass q(:,:,:,1:1) as 4D array
call GFDL_interstitial%create(Atm(mygrid)%bd%is, Atm(mygrid)%bd%ie, Atm(mygrid)%bd%isd, Atm(mygrid)%bd%ied, &
Atm(mygrid)%bd%js, Atm(mygrid)%bd%je, Atm(mygrid)%bd%jsd, Atm(mygrid)%bd%jed, &
Atm(mygrid)%npz, Atm(mygrid)%ng, &
dt_atmos, p_split, Atm(mygrid)%flagstruct%k_split, &
zvir, Atm(mygrid)%flagstruct%p_ref, Atm(mygrid)%ak, Atm(mygrid)%bk, &
liq_wat>0, ice_wat>0, rainwat>0, snowwat>0, graupel>0, &
cld_amt>0, kappa, Atm(mygrid)%flagstruct%hydrostatic, &
Atm(mygrid)%flagstruct%do_sat_adj, &
Atm(mygrid)%delp, Atm(mygrid)%delz, Atm(mygrid)%gridstruct%area_64, &
Atm(mygrid)%peln, Atm(mygrid)%phis, Atm(mygrid)%pkz, Atm(mygrid)%pt, &
#ifdef MULTI_GASES
Atm(mygrid)%q(:,:,:,1:max(1,num_gas)), &
#else
Atm(mygrid)%q(:,:,:,1:1), &
#endif
Atm(mygrid)%q(:,:,:,sphum), Atm(mygrid)%q(:,:,:,liq_wat), &
Atm(mygrid)%q(:,:,:,ice_wat), Atm(mygrid)%q(:,:,:,rainwat), &
Atm(mygrid)%q(:,:,:,snowwat), Atm(mygrid)%q(:,:,:,graupel), &
Atm(mygrid)%q(:,:,:,cld_amt), Atm(mygrid)%q_con, nthreads, &
Atm(mygrid)%flagstruct%nwat, &
#ifdef MULTI_GASES
ngas=num_gas, rilist=ri, cpilist=cpi, &
#endif
mpirank=mpp_pe(), mpiroot=mpp_root_pe())
if (Atm(mygrid)%flagstruct%do_sat_adj) then
! Initialize fast physics
call ccpp_physics_init(cdata, suite_name=trim(ccpp_suite), group_name="fast_physics", ierr=ierr)
if (ierr/=0) then
cdata%errmsg = ' atmosphere_dynamics: error in ccpp_physics_init for group fast_physics: ' // trim(cdata%errmsg)
call mpp_error (FATAL, cdata%errmsg)
endif
endif
! --- initiate the start for a restarted regional forecast
if ( Atm(mygrid)%gridstruct%regional .and. Atm(mygrid)%flagstruct%warm_start ) then
call start_regional_restart(Atm(1), &
isc, iec, jsc, jec, &
isd, ied, jsd, jed )
endif
if ( Atm(mygrid)%flagstruct%nudge ) then
call fv_nwp_nudge_init( Time, Atm(mygrid)%atmos_axes, npz, zvir, Atm(mygrid)%ak, Atm(mygrid)%bk, Atm(mygrid)%ts, &
Atm(mygrid)%phis, Atm(mygrid)%gridstruct, Atm(mygrid)%ks, Atm(mygrid)%npx, Atm(mygrid)%neststruct, Atm(mygrid)%bd)
call mpp_error(NOTE, 'NWP nudging is active')
endif
call fv_io_register_nudge_restart ( Atm )
if ( Atm(mygrid)%flagstruct%na_init>0 ) then
if ( .not. Atm(mygrid)%flagstruct%hydrostatic ) then
call prt_maxmin('Before adi: W', Atm(mygrid)%w, isc, iec, jsc, jec, Atm(mygrid)%ng, npz, 1.)
endif
call adiabatic_init(zvir,Atm(mygrid)%flagstruct%nudge_dz, time)
if ( .not. Atm(mygrid)%flagstruct%hydrostatic ) then
call prt_maxmin('After adi: W', Atm(mygrid)%w, isc, iec, jsc, jec, Atm(mygrid)%ng, npz, 1.)
! Not nested?
call prt_height('na_ini Z500', isc,iec, jsc,jec, 3, npz, 500.E2, Atm(mygrid)%phis, Atm(mygrid)%delz, &
Atm(mygrid)%peln, Atm(mygrid)%gridstruct%area_64(isc:iec,jsc:jec), Atm(mygrid)%gridstruct%agrid_64(isc:iec,jsc:jec,2))
endif
else
call mpp_error(NOTE,'No adiabatic initialization correction in use')
endif
#ifdef DEBUG
call fv_diag(Atm(mygrid:mygrid), zvir, Time, -1)
if (Atm(mygrid)%coarse_graining%write_coarse_diagnostics) then
call fv_coarse_diag(Atm(mygrid:mygrid), fv_time)
endif
#endif
end subroutine atmosphere_init
!>@brief The subroutine 'p_adi' computes (ps, pk, pe, peln, pkz)
!! given (ptop, delp).
subroutine p_adi(km, ng, ifirst, ilast, jfirst, jlast, ptop, &
delp, pt, ps, pe, peln, pk, pkz, hydrostatic)
! Given (ptop, delp) computes (ps, pk, pe, peln, pkz)
! Input:
integer, intent(in):: km, ng
integer, intent(in):: ifirst, ilast !< Longitude strip
integer, intent(in):: jfirst, jlast !< Latitude strip
logical, intent(in):: hydrostatic
real, intent(in):: ptop
real, intent(in):: pt(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
real, intent(in):: delp(ifirst-ng:ilast+ng,jfirst-ng:jlast+ng, km)
! Output:
real, intent(out) :: ps(ifirst-ng:ilast+ng, jfirst-ng:jlast+ng)
real, intent(out) :: pk(ifirst:ilast, jfirst:jlast, km+1)
real, intent(out) :: pe(ifirst-1:ilast+1,km+1,jfirst-1:jlast+1) !< Ghosted Edge pressure
real, intent(out) :: peln(ifirst:ilast, km+1, jfirst:jlast) !< Edge pressure
real, intent(out) :: pkz(ifirst:ilast, jfirst:jlast, km)
! Local
real pek
integer i, j, k
pek = ptop ** kappa
!$OMP parallel do default (none) &
!$OMP shared (ifirst,ilast,jfirst,jlast,km,ptop,pek,pe,pk, &
!$OMP ps,delp,peln,hydrostatic,pkz) &
!$OMP private (j, i, k)
do j=jfirst,jlast
do i=ifirst,ilast
pe(i,1,j) = ptop
pk(i,j,1) = pek
enddo
do k=2,km+1
do i=ifirst,ilast
pe(i,k,j) = pe(i,k-1,j) + delp(i,j,k-1)
peln(i,k,j) = log(pe(i,k,j))
pk(i,j,k) = exp( kappa*peln(i,k,j) )
enddo
enddo
do i=ifirst,ilast
ps(i,j) = pe(i,km+1,j)
enddo
if ( hydrostatic ) then
do k=1,km
do i=ifirst,ilast
pkz(i,j,k) = (pk(i,j,k+1)-pk(i,j,k))/(kappa*(peln(i,k+1,j)-peln(i,k,j)))
enddo
enddo
endif
enddo
end subroutine p_adi
!>@brief The subroutine 'atmosphere_dynamics' is an API for the main driver
!! of the FV3 dynamical core responsible for executing a "dynamics" step.
subroutine atmosphere_dynamics ( Time )
type(time_type),intent(in) :: Time
integer :: n, psc, atmos_time_step
integer :: k, w_diff, nt_dyn, n_split_loc, seconds, days
logical :: used
real :: rdt
type(time_type) :: atmos_time
!---- Call FV dynamics -----
call mpp_clock_begin (id_dynam)
n = mygrid
call get_time (time, seconds, days)
! if (seconds < 10800 .and. days == 0) then
! n_split_loc = (Atm(n)%flagstruct%n_split * 3) / 2
if (seconds < nint(3600*Atm(n)%flagstruct%fhouri) .and. Atm(n)%flagstruct%fac_n_spl > 1.0) then
n_split_loc = nint(Atm(n)%flagstruct%n_split * Atm(n)%flagstruct%fac_n_spl)
else
n_split_loc = Atm(n)%flagstruct%n_split
endif
! write(0,*)' before calling init n_split_loc=',n_split_loc,' seconds=',seconds,' days=',days,&
! ' n_split=',Atm(mytile)%flagstruct%n_split,' mytile=',mytile
a_step = a_step + 1
!
!*** If this is a regional run then read in the next boundary data when it is time.
!
if(Atm(n)%flagstruct%regional)then
call read_new_bc_data(Atm(n), Time, Time_step_atmos, p_split, &
isd, ied, jsd, jed )
endif
do psc=1,abs(p_split)
p_step = psc
call timing_on('fv_dynamics')
!uc/vc only need be same on coarse grid? However BCs do need to be the same
call fv_dynamics(npx, npy, npz, nq, Atm(n)%ng, dt_atmos/real(abs(p_split)),&
Atm(n)%flagstruct%consv_te, Atm(n)%flagstruct%fill, &
Atm(n)%flagstruct%reproduce_sum, kappa, cp_air, zvir, &
Atm(n)%ptop, Atm(n)%ks, nq, &
n_split_loc, Atm(n)%flagstruct%q_split, &
! Atm(n)%flagstruct%n_split, Atm(n)%flagstruct%q_split, &
Atm(n)%u, Atm(n)%v, Atm(n)%w, Atm(n)%delz, &
Atm(n)%flagstruct%hydrostatic, &
Atm(n)%pt , Atm(n)%delp, Atm(n)%q, Atm(n)%ps, &
Atm(n)%pe, Atm(n)%pk, Atm(n)%peln, &
Atm(n)%pkz, Atm(n)%phis, Atm(n)%q_con, &
Atm(n)%omga, Atm(n)%ua, Atm(n)%va, Atm(n)%uc, &
Atm(n)%vc, Atm(n)%ak, Atm(n)%bk, Atm(n)%mfx, &
Atm(n)%mfy , Atm(n)%cx, Atm(n)%cy, Atm(n)%ze0, &
Atm(n)%flagstruct%hybrid_z, &
Atm(n)%gridstruct, Atm(n)%flagstruct, &
Atm(n)%neststruct, Atm(n)%idiag, Atm(n)%bd, &
Atm(n)%parent_grid, Atm(n)%domain,Atm(n)%diss_est, &
Atm(n)%inline_mp)
call timing_off('fv_dynamics')
if (ngrids > 1 .and. (psc < p_split .or. p_split < 0)) then
call mpp_sync()
call timing_on('TWOWAY_UPDATE')
call twoway_nesting(Atm, ngrids, grids_on_this_pe, zvir, fv_time, mygrid)
call timing_off('TWOWAY_UPDATE')
endif
enddo !p_split
call mpp_clock_end (id_dynam)
!-----------------------------------------------------
!--- COMPUTE SUBGRID Z
!-----------------------------------------------------
!--- zero out tendencies
call mpp_clock_begin (id_subgridz)
u_dt(:,:,:) = 0. ! These are updated by fv_subgrid_z
v_dt(:,:,:) = 0.
! t_dt is used for two different purposes:
! 1 - to calculate the diagnostic temperature tendency from fv_subgrid_z
! 2 - as an accumulator for the IAU increment and physics tendency
! because of this, it will need to be zeroed out after the diagnostic is calculated
t_dt(:,:,:) = Atm(n)%pt(isc:iec,jsc:jec,:)
qv_dt(:,:,:) = Atm(n)%q (isc:iec,jsc:jec,:,sphum)
rdt = 1./dt_atmos
w_diff = get_tracer_index (MODEL_ATMOS, 'w_diff' )
if ( Atm(n)%flagstruct%fv_sg_adj > 0 ) then
nt_dyn = nq
if ( w_diff /= NO_TRACER ) then
nt_dyn = nq - 1
endif
call fv_subgrid_z(isd, ied, jsd, jed, isc, iec, jsc, jec, Atm(n)%npz, &
nt_dyn, dt_atmos, Atm(n)%flagstruct%fv_sg_adj, &
Atm(n)%flagstruct%nwat, Atm(n)%delp, Atm(n)%pe, &
Atm(n)%peln, Atm(n)%pkz, Atm(n)%pt, Atm(n)%q, &
Atm(n)%ua, Atm(n)%va, Atm(n)%flagstruct%hydrostatic,&
Atm(n)%w, Atm(n)%delz, u_dt, v_dt, t_dt, &
Atm(n)%flagstruct%n_sponge)
endif
#ifdef USE_Q_DT
if ( .not. Atm(n)%flagstruct%hydrostatic .and. w_diff /= NO_TRACER ) then
!$OMP parallel do default (none) &
!$OMP shared (isc, iec, jsc, jec, w_diff, n, Atm, q_dt) &
!$OMP private (k)
do k=1, Atm(n)%npz
Atm(n)%q(isc:iec,jsc:jec,k,w_diff) = Atm(n)%w(isc:iec,jsc:jec,k) + w0_big
q_dt(:,:,k,w_diff) = 0.
enddo
endif
#endif
if (Atm(1)%idiag%id_u_dt_sg > 0) then
used = send_data(Atm(1)%idiag%id_u_dt_sg, u_dt(isc:iec,jsc:jec,:), fv_time)
endif
if (Atm(1)%idiag%id_v_dt_sg > 0) then
used = send_data(Atm(1)%idiag%id_v_dt_sg, v_dt(isc:iec,jsc:jec,:), fv_time)
endif
if (Atm(1)%idiag%id_t_dt_sg > 0) then
t_dt(:,:,:) = rdt*(Atm(1)%pt(isc:iec,jsc:jec,:) - t_dt(:,:,:))
used = send_data(Atm(1)%idiag%id_t_dt_sg, t_dt, fv_time)
endif
if (Atm(1)%idiag%id_qv_dt_sg > 0) then
qv_dt(:,:,:) = rdt*(Atm(1)%q(isc:iec,jsc:jec,:,sphum) - qv_dt(:,:,:))
used = send_data(Atm(1)%idiag%id_qv_dt_sg, qv_dt, fv_time)
endif
! zero out t_dt for use as an accumulator
t_dt = 0.
call mpp_clock_end (id_subgridz)
end subroutine atmosphere_dynamics
!>@brief The subroutine 'atmosphere_end' is an API for the termination of the
!! FV3 dynamical core responsible for writing out a restart and final diagnostic state.
subroutine atmosphere_end (Time, Grid_box, restart_endfcst)
use ccpp_static_api, only: ccpp_physics_finalize
use CCPP_data, only: ccpp_suite
use CCPP_data, only: cdata => cdata_tile
type (time_type), intent(in) :: Time
type(grid_box_type), intent(inout) :: Grid_box
logical, intent(in) :: restart_endfcst
integer :: ierr
if (Atm(mygrid)%flagstruct%do_sat_adj) then
! Finalize fast physics
call ccpp_physics_finalize(cdata, suite_name=trim(ccpp_suite), group_name="fast_physics", ierr=ierr)
if (ierr/=0) then
cdata%errmsg = ' atmosphere_dynamics: error in ccpp_physics_finalize for group fast_physics: ' // trim(cdata%errmsg)
call mpp_error (FATAL, cdata%errmsg)
endif
endif
! initialize domains for writing global physics data
if ( Atm(mygrid)%flagstruct%nudge ) call fv_nwp_nudge_end
if (first_diag) then
call timing_on('FV_DIAG')
call fv_diag(Atm(mygrid:mygrid), zvir, fv_time, Atm(mygrid)%flagstruct%print_freq)
call fv_nggps_diag_init(Atm(mygrid:mygrid), Atm(mygrid)%atmos_axes, fv_time)
call fv_nggps_diag(Atm(mygrid:mygrid), zvir, fv_time)
if (Atm(mygrid)%coarse_graining%write_coarse_diagnostics) then
call fv_coarse_diag(Atm(mygrid:mygrid), fv_time)
endif
first_diag = .false.
call timing_off('FV_DIAG')
endif
call fv_end(Atm, mygrid, restart_endfcst)
deallocate (Atm)
deallocate( u_dt, v_dt, t_dt, qv_dt, pref, dum1d )
end subroutine atmosphere_end
!>@brief The subroutine 'atmosphere_restart' is an API to save restart information
!! at a given timestamp.
!>@detail This API is used to provide intermediate restart capability to the atmospheric
!! driver.
subroutine atmosphere_restart(timestamp)
character(len=*), intent(in) :: timestamp
call fv_write_restart(Atm(mygrid), timestamp)
end subroutine atmosphere_restart
!>@brief The subroutine 'atmospehre_resolution' is an API to return the local
!! extents of the current MPI-rank or the global extents of the current
!! cubed-sphere tile.
subroutine atmosphere_resolution (i_size, j_size, global)
integer, intent(out) :: i_size, j_size
logical, intent(in), optional :: global
logical :: local
local = .true.
if( PRESENT(global) ) local = .NOT.global
if( local ) then
i_size = iec - isc + 1
j_size = jec - jsc + 1
else
i_size = npx - 1
j_size = npy - 1
endif
end subroutine atmosphere_resolution
!>@brief The subroutine 'atmosphere_pref' is an API to return the reference pressure.
subroutine atmosphere_pref (p_ref)
real, dimension(:,:), intent(inout) :: p_ref
p_ref = pref
end subroutine atmosphere_pref
subroutine atmosphere_control_data (i1, i2, j1, j2, kt, p_hydro, hydro, tile_of_mosaic, global_tile_num)
integer, intent(out) :: i1, i2, j1, j2, kt
logical, intent(out), optional :: p_hydro, hydro
integer, intent(out), optional :: tile_of_mosaic, global_tile_num
i1 = Atm(mygrid)%bd%isc
i2 = Atm(mygrid)%bd%iec
j1 = Atm(mygrid)%bd%jsc
j2 = Atm(mygrid)%bd%jec
kt = Atm(mygrid)%npz
if (present(global_tile_num)) global_tile_num = Atm(mygrid)%global_tile
if (present(tile_of_mosaic)) tile_of_mosaic = Atm(mygrid)%tile_of_mosaic
if (present(p_hydro)) p_hydro = Atm(mygrid)%flagstruct%phys_hydrostatic
if (present( hydro)) hydro = Atm(mygrid)%flagstruct%hydrostatic
end subroutine atmosphere_control_data
!>@brief The subroutine 'atmosphere_grid_ctr' is an API that returns the longitude and
!! latitude cell centers of the current MPI-rank.
subroutine atmosphere_grid_ctr (lon, lat)
real(kind=kind_phys), intent(out) :: lon(:,:), lat(:,:) !< Unit: radian
! Local data:
integer i,j
do j=jsc,jec
do i=isc,iec
lon(i-isc+1,j-jsc+1) = Atm(mygrid)%gridstruct%agrid_64(i,j,1)
lat(i-isc+1,j-jsc+1) = Atm(mygrid)%gridstruct%agrid_64(i,j,2)
enddo
enddo
end subroutine atmosphere_grid_ctr
!>@brief The subroutine 'atmosphere_grid_bdry' is an API to returns the
!! longitude and latitude finite volume edges (grid box) for the current MPI-rank.
subroutine atmosphere_grid_bdry (blon, blat, global)
real(kind=kind_phys), intent(out) :: blon(:,:), blat(:,:) !< Unit: radian
logical, intent(in), optional :: global
! Local data:
integer i,j
if( PRESENT(global) ) then
if (global) call mpp_error(FATAL, '==> global grid is no longer available &
& in the Cubed Sphere')
endif
do j=jsc,jec+1
do i=isc,iec+1
blon(i-isc+1,j-jsc+1) = Atm(mygrid)%gridstruct%grid(i,j,1)
blat(i-isc+1,j-jsc+1) = Atm(mygrid)%gridstruct%grid(i,j,2)
enddo
enddo
end subroutine atmosphere_grid_bdry
subroutine set_atmosphere_pelist ()
call mpp_set_current_pelist(Atm(mygrid)%pelist, no_sync=.TRUE.)
end subroutine set_atmosphere_pelist
subroutine get_nth_domain_info(n, layout, nx, ny, pelist)
integer, intent(in) :: n
integer, intent(out) :: layout(2)
integer, intent(out) :: nx, ny
integer, pointer, intent(out) :: pelist(:)
layout(1:2) = Atm(n)%layout(1:2)
nx = Atm(n)%npx -1
ny = Atm(n)%npy -1
pelist => Atm(n)%pelist
end subroutine get_nth_domain_info
!>@brief The subroutine 'atmosphere_domain' is an API to return
!! the "domain2d" variable associated with the coupling grid and the
!! decomposition for the current cubed-sphere tile.
!>@detail Coupling is done using the mass/temperature grid with no halos.
subroutine atmosphere_domain ( fv_domain, rd_domain, layout, regional, nested, &
ngrids_atmos, mygrid_atmos, pelist )
type(domain2d), intent(out) :: fv_domain, rd_domain
integer, intent(out) :: layout(2)
logical, intent(out) :: regional
logical, intent(out) :: nested
integer, intent(out) :: ngrids_atmos
integer, intent(out) :: mygrid_atmos
integer, pointer, intent(out) :: pelist(:)
integer :: n
fv_domain = Atm(mygrid)%domain_for_coupler
rd_domain = Atm(mygrid)%domain_for_read
layout(1:2) = Atm(mygrid)%layout(1:2)
regional = Atm(mygrid)%flagstruct%regional
nested = ngrids > 1
ngrids_atmos = ngrids
mygrid_atmos = mygrid
call set_atmosphere_pelist()
pelist => Atm(mygrid)%pelist
end subroutine atmosphere_domain
!>@brief The subroutine 'atmosphere_diag_axes' is an API to return the axis indices
!! for the atmospheric (mass) grid.
subroutine atmosphere_diag_axes ( axes )
integer, intent(out) :: axes (:)
!----- returns the axis indices for the atmospheric (mass) grid -----
if ( size(axes(:)) < 0 .or. size(axes(:)) > 4 ) call error_mesg ( &
'get_atmosphere_axes in atmosphere_mod', &
'size of argument is incorrect', FATAL )
axes (1:size(axes(:))) = Atm(mygrid)%atmos_axes (1:size(axes(:)))
end subroutine atmosphere_diag_axes
!>@brief The subroutine 'atmosphere_etalvls' is an API to return the ak/bk
!! pairs used to compute the eta or pressure levels.
!>@detail By default, the vertical dimension assumes the standard FV3 convention
!! of TOA (k=1) to Surface (k=npz).
subroutine atmosphere_etalvls (ak, bk, flip)
real(kind=kind_phys), pointer, dimension(:), intent(inout) :: ak, bk
logical, intent(in) :: flip !< control vertical index flipping
allocate(ak(npz+1))
allocate(bk(npz+1))
if (flip) then
ak(1:npz+1) = Atm(mygrid)%ak(npz+1:1:-1)