-
Notifications
You must be signed in to change notification settings - Fork 562
/
Copy pathwav_import_export.F90
1711 lines (1510 loc) · 66.5 KB
/
wav_import_export.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
!> @file wav_import_export
!!
!> Manage the import/export state and fields
!!
!> @details Contains the public routines to advertise and realize
!! the import and export fields and the public routines to fill
!! the import and export fields within the ESMF States.
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
module wav_import_export
use ESMF
use NUOPC
use NUOPC_Model
use wav_shr_flags
use wav_kind_mod , only : r8 => shr_kind_r8, r4 => shr_kind_r4, i4 => shr_kind_i4
use wav_kind_mod , only : CL => shr_kind_cl, CS => shr_kind_cs
use wav_shr_mod , only : ymd2date
use wav_shr_mod , only : chkerr
use wav_shr_mod , only : state_diagnose, state_reset, state_getfldptr, state_fldchk
use wav_shr_mod , only : wav_coupling_to_cice, nwav_elev_spectrum, merge_import, dbug_flag, multigrid, unstr_mesh
use constants , only : grav, tpi, dwat, dair
use w3parall , only : init_get_isea
implicit none
private ! except
public :: advertise_fields !< @public create a list of fields and advertise them
public :: realize_fields !< @public realize a list of advertised fields
public :: import_fields !< @public fill WW3 fields using values in the import state
public :: export_fields !< @public fill values in the export state using WW3 fields
public :: CalcRoughl !< @public calculate the roughness length
private :: fldlist_add !< @private add a field name to a list of field names
private :: fldlist_realize !< @private realize a field in a list of field names
private :: set_importmask !< @private set the import mask when merge_import is true
private :: check_globaldata !< @private write values in a field to a netCDF file for debugging
private :: readfromfile !< @private read values from a file
interface FillGlobalInput
module procedure fillglobal_with_import
module procedure fillglobal_with_merge_import
end interface FillGlobalInput
type fld_list_type !< @private a structure for the list of fields
character(len=128) :: stdname !< a standard field name
integer :: ungridded_lbound = 0 !< the ungridded dimension lower bound
integer :: ungridded_ubound = 0 !< the ugridded dimension upper bound
end type fld_list_type
integer, parameter :: fldsMax = 100 !< the maximum allowed number of fields in a state
integer :: fldsToWav_num = 0 !< initial value of the number of fields sent to the wave model
integer :: fldsFrWav_num = 0 !< initial value of the number of fields sent from the wave model
type (fld_list_type) :: fldsToWav(fldsMax) !< a structure containing the list of fields to the wave model
type (fld_list_type) :: fldsFrWav(fldsMax) !< a structure containing the list of fields from the wave model
real(r4), allocatable :: import_mask(:) !< the mask for valid import data
real(r8), parameter :: zero = 0.0_r8 !< a named constant
#ifdef W3_CESMCOUPLED
logical :: cesmcoupled = .true. !< logical defining a CESM use case
#else
logical :: cesmcoupled = .false. !< logical defining a non-CESM use case (UWM)
#endif
integer, public :: nseal_cpl !< the number of local sea points on a processor, exclusive
!! of the ghost points. For non-PDLIB cases, this is nseal
character(*),parameter :: u_FILE_u = & !< a character string for an ESMF log message
__FILE__
!===============================================================================
contains
!===============================================================================
!> Set up the list of exchanged field to be advertised
!!
!> @details Called by InitializAdvertise, a list of standard field names to or
!! from the wave model is created and then advertised in either the import or
!! export state. A field with name set by the configuration variable ScalarFieldName
!! and size of ScalarFieldCount is added to the list of fields in the export state
!! and is used by CMEPS to write mediator history and restart fields as 2D arrays
!!
!! @param importState an ESMF_State for the import
!! @param exportState an ESMF_State for the export
!! @param[in] flds_scalar_name the name of the scalar field
!! @param[out] rc a return code
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine advertise_fields(importState, ExportState, flds_scalar_name, rc)
! input/output variables
type(ESMF_State) :: importState
type(ESMF_State) :: exportState
character(len=*) , intent(in) :: flds_scalar_name
integer , intent(out) :: rc
! local variables
integer :: n, num
character(len=2) :: fvalue
character(len=*), parameter :: subname='(wav_import_export:advertise_fields)'
!-------------------------------------------------------------------------------
rc = ESMF_SUCCESS
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' called', ESMF_LOGMSG_INFO)
!--------------------------------
! Advertise import fields
!--------------------------------
!call fldlist_add(fldsToWav_num, fldsToWav, 'So_h' )
call fldlist_add(fldsToWav_num, fldsToWav, 'Si_ifrac' )
call fldlist_add(fldsToWav_num, fldsToWav, 'So_u' )
call fldlist_add(fldsToWav_num, fldsToWav, 'So_v' )
call fldlist_add(fldsToWav_num, fldsToWav, 'So_t' )
call fldlist_add(fldsToWav_num, fldsToWav, 'Sa_tbot' )
if (cesmcoupled) then
call fldlist_add(fldsToWav_num, fldsToWav, 'Sa_u' )
call fldlist_add(fldsToWav_num, fldsToWav, 'Sa_v' )
call fldlist_add(fldsToWav_num, fldsToWav, 'So_bldepth' )
else
call fldlist_add(fldsToWav_num, fldsToWav, 'Sa_u10m' )
call fldlist_add(fldsToWav_num, fldsToWav, 'Sa_v10m' )
end if
if (wav_coupling_to_cice) then
call fldlist_add(fldsToWav_num, fldsToWav, 'Si_thick' )
call fldlist_add(fldsToWav_num, fldsToWav, 'Si_floediam')
end if
do n = 1,fldsToWav_num
call NUOPC_Advertise(importState, standardName=fldsToWav(n)%stdname, &
TransferOfferGeomObject='will provide', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end do
!--------------------------------
! Advertise export fields
!--------------------------------
if (.not. unstr_mesh) then
call fldlist_add(fldsFrWav_num, fldsFrWav, trim(flds_scalar_name))
end if
if (cesmcoupled) then
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_lamult' )
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_lasl' )
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_ustokes')
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_vstokes')
else
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_z0')
end if
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_pstokes_x', ungridded_lbound=1, ungridded_ubound=3)
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_pstokes_y', ungridded_lbound=1, ungridded_ubound=3)
! AA TODO: In the above fldlist_add calls, we are passing hardcoded ungridded_ubound values (3) because, USSPF(2)
! is not initialized yet. It is set during w3init which gets called at a later phase (realize). A permanent solution
! will be implemented soon based on receiving USSP and USSPF from the coupler instead of the mod_def file. This will
! also ensure compatibility with the ocean component since ocean will also receive these from the coupler.
if (wav_coupling_to_cice) then
call fldlist_add(fldsFrWav_num, fldsFrWav, 'Sw_elevation_spectrum', &
ungridded_lbound=1, ungridded_ubound=nwav_elev_spectrum)
end if
do n = 1,fldsFrWav_num
call NUOPC_Advertise(exportState, standardName=fldsFrWav(n)%stdname, &
TransferOfferGeomObject='will provide', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end do
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' done', ESMF_LOGMSG_INFO)
end subroutine advertise_fields
!===============================================================================
!> Realize the advertised fields
!!
!> @details Called by InitializeRealize, realize the advertised fields on the mesh
!! and set all initial values to zero
!!
!! @param gcomp an ESMF_GridComp object
!! @param mesh an ESMF_Mesh object
!! @param[in] flds_scalar_name the name of the scalar field
!! @param[in] flds_scalar_num the number of scalar fields
!! @param[out] rc a return code
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine realize_fields(gcomp, mesh, flds_scalar_name, flds_scalar_num, rc)
! input/output variables
type(ESMF_GridComp) :: gcomp
type(ESMF_Mesh) :: mesh
character(len=*) , intent(in) :: flds_scalar_name
integer , intent(in) :: flds_scalar_num
integer , intent(out) :: rc
! local variables
type(ESMF_State) :: importState
type(ESMF_State) :: exportState
character(len=*), parameter :: subname='(wav_import_export:realize_fields)'
!---------------------------------------------------------------------------
rc = ESMF_SUCCESS
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' called', ESMF_LOGMSG_INFO)
call NUOPC_ModelGet(gcomp, importState=importState, exportState=exportState, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call fldlist_realize( &
state=ExportState, &
fldList=fldsFrWav, &
numflds=fldsFrWav_num, &
flds_scalar_name=flds_scalar_name, &
flds_scalar_num=flds_scalar_num, &
tag=subname//':WW3Export',&
mesh=mesh, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call fldlist_realize( &
state=importState, &
fldList=fldsToWav, &
numflds=fldsToWav_num, &
flds_scalar_name=flds_scalar_name, &
flds_scalar_num=flds_scalar_num, &
tag=subname//':WW3Import',&
mesh=mesh, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_reset(ExportState, zero, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_reset(ImportState, zero, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (dbug_flag > 5) then
call state_diagnose(exportState, 'after state_reset', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' done', ESMF_LOGMSG_INFO)
end subroutine realize_fields
!===============================================================================
!> Fill WW3 fields with values from the import state
!!
!> @details Called by ModelAdvance, a global field for each connected field is
!! created in SetGlobalInput and used to fill the internal WW3 global variables in
!! FillGlobalInput. Optionally, the WW3 field can be created by merging with a
!! provided field in cases where the WW3 model domain extends outside the source
!! domain
!!
!! @param[inout] gcomp an ESMF_GridComp object
!! @param[in] time0 the starting time of ModelAdvance
!! @param[in] timen the ending time of ModelAdvance
!! @param[out] rc return code
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine import_fields( gcomp, time0, timen, rc )
!---------------------------------------------------------------------------
! Obtain the wave input from the mediator
!---------------------------------------------------------------------------
use w3gdatmd , only: nsea, nseal, MAPSTA, NX, NY, w3setg
use w3idatmd , only: CX0, CY0, CXN, CYN, DT0, DTN, ICEI, WLEV, INFLAGS1, ICEP1, ICEP5
use w3idatmd , only: TC0, TCN, TLN, TIN, TI1, TI5, TW0, TWN, WX0, WY0, WXN, WYN
use w3idatmd , only: UX0, UY0, UXN, UYN, TU0, TUN
use w3idatmd , only: tfn, w3seti
use w3odatmd , only: w3seto
use w3wdatmd , only: time, w3setw
#ifdef W3_CESMCOUPLED
use w3idatmd , only: HSL
#else
use wmupdtmd , only: wmupd2
use wmmdatmd , only: wmsetm
use wmmdatmd , only: mdse, mdst, nrgrd, inpmap
#ifdef W3_MPI
use wmmdatmd , only: mpi_comm_grd
#endif
#endif
! input/output variables
type(ESMF_GridComp) , intent(inout) :: gcomp
integer , intent(in) :: time0(2), timen(2)
integer , intent(out) :: rc
! Local variables
type(ESMF_State) :: importState
type(ESMF_VM) :: vm
type(ESMF_Clock) :: clock
real(r4) :: global_data(nsea)
real(r4), allocatable :: global_data2(:)
real(r4) :: def_value
character(len=10) :: uwnd
character(len=10) :: vwnd
integer :: imod, j, jmod
integer :: mpi_comm_null = -1
real(r4), allocatable :: wxdata(:) ! only needed if merge_import
real(r4), allocatable :: wydata(:) ! only needed if merge_import
character(len=CL) :: msgString
character(len=*), parameter :: subname='(wav_import_export:import_fields)'
!---------------------------------------------------------------------------
rc = ESMF_SUCCESS
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' called', ESMF_LOGMSG_INFO)
if (cesmcoupled) then
uwnd = 'Sa_u'
vwnd = 'Sa_v'
else
uwnd = 'Sa_u10m'
vwnd = 'Sa_v10m'
end if
! Get import state, clock and vm
call ESMF_GridCompGet(gcomp, clock=clock, importState=importState, vm=vm, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (dbug_flag > 5) then
call state_diagnose(importState, 'at import ', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
! input fields associated with W3FLDG calls in ww3_shel.ftn
! fill both the lower (0) and upper (N) bound data with the same values
! fill with special values as default, these should not be used in practice
! set time for input data to time0 and timen (shouldn't matter)
def_value = 0.0_r4
#ifndef W3_CESMCOUPLED
call w3setg ( 1, mdse, mdst )
call w3seti ( 1, mdse, mdst )
#endif
! ---------------
! INFLAGS1(1)
! ---------------
if (INFLAGS1(1)) then
TLN = timen
WLEV(:,:) = def_value ! water level
if (state_fldchk(importState, 'So_h')) then
call SetGlobalInput(importState, 'So_h', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, WLEV)
end if
endif
! ---------------
! INFLAGS1(2) - ocn current fields
! ---------------
if (INFLAGS1(2)) then
TC0 = time0 ! times for ocn current fields
TCN = timen
CX0(:,:) = def_value ! ocn u current
CXN(:,:) = def_value
if (state_fldchk(importState, 'So_u')) then
call SetGlobalInput(importState, 'So_u', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, CX0)
call FillGlobalInput(global_data, CXN)
end if
CY0(:,:) = def_value ! ocn v current
CYN(:,:) = def_value
if (state_fldchk(importState, 'So_v')) then
call SetGlobalInput(importState, 'So_v', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, CY0)
call FillGlobalInput(global_data, CYN)
end if
end if
! ---------------
! INFLAGS1(3) - atm wind/temp fields
! ---------------
if (INFLAGS1(3)) then
TW0 = time0 ! times for atm wind/temp fields.
TWN = timen
if (merge_import) then
! set mask using u-wind field if merge_import; assume all import fields
! will have same missing overlap region
! import_mask memory will be allocate in set_importmask
call set_importmask(importState, clock, trim(uwnd), vm, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
allocate(wxdata(nsea))
allocate(wydata(nsea))
call readfromfile('WND', wxdata, wydata, time0, timen, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (dbug_flag > 10) then
call check_globaldata(gcomp, 'wxdata', wxdata, nsea, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call check_globaldata(gcomp, 'wydata', wydata, nsea, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call check_globaldata(gcomp, 'import_mask', import_mask, nsea, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
end if
! atm u wind
WX0(:,:) = def_value
WXN(:,:) = def_value
if (state_fldchk(importState, trim(uwnd))) then
call SetGlobalInput(importState, trim(uwnd), vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (merge_import) then
call FillGlobalInput(global_data, import_mask, wxdata, WX0)
call FillGlobalInput(global_data, import_mask, wxdata, WXN)
if (dbug_flag > 10) then
call check_globaldata(gcomp, 'wx0', wx0, nx*ny, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
else
call FillGlobalInput(global_data, WX0)
call FillGlobalInput(global_data, WXN)
end if
end if
! atm v wind
WY0(:,:) = def_value
WYN(:,:) = def_value
if (state_fldchk(importState, trim(vwnd))) then
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call SetGlobalInput(importState, trim(vwnd), vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
if (merge_import) then
call FillGlobalInput(global_data, import_mask, wydata, WY0)
call FillGlobalInput(global_data, import_mask, wydata, WYN)
if (dbug_flag > 10) then
call check_globaldata(gcomp, 'wy0', wy0, nx*ny, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
else
call FillGlobalInput(global_data, WY0)
call FillGlobalInput(global_data, WYN)
end if
end if
! air temp - ocn temp
DT0(:,:) = def_value
DTN(:,:) = def_value
if ((state_fldchk(importState, 'So_t')) .and. (state_fldchk(importState, 'Sa_tbot'))) then
allocate(global_data2(nsea))
call SetGlobalInput(importState, 'Sa_tbot', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call SetGlobalInput(importState, 'So_t', vm, global_data2, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
! So_tbot - So_t
global_data = global_data - global_data2
call FillGlobalInput(global_data, DT0)
call FillGlobalInput(global_data, DTN)
deallocate(global_data2)
end if
! Deallocate memory for merge_import
if (merge_import) then
deallocate(wxdata)
deallocate(wydata)
end if
end if
! ---------------
! INFLAGS1(4) - ice fraction field
! ---------------
if (INFLAGS1(4)) then
TIN = timen ! time for ice field
ICEI(:,:) = def_value ! ice frac
if (state_fldchk(importState, 'Si_ifrac')) then
call SetGlobalInput(importState, 'Si_ifrac', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, ICEI)
end if
end if
#ifdef W3_CESMCOUPLED
! ---------------
! ocean boundary layer depth - always assume that this is being imported for CESM
! ---------------
call SetGlobalInput(importState, 'So_bldepth', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
! ocn mixing layer depth
global_data = max(global_data, 5.)*0.2
call FillGlobalInput(global_data, HSL)
#endif
! ---------------
! INFLAGS1(5) - atm momentum fields
! ---------------
if (INFLAGS1(5)) then
TU0 = time0 ! times for atm momentum fields.
TUN = timen
UX0(:,:) = def_value ! atm u momentum
UXN(:,:) = def_value
if (state_fldchk(importState, 'Faxa_taux')) then
call SetGlobalInput(importState, 'Faxa_taux', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, UX0)
call FillGlobalInput(global_data, UXN)
end if
UY0(:,:) = def_value ! atm v momentum
UYN(:,:) = def_value
if (state_fldchk(importState, 'Faxa_tauy')) then
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call SetGlobalInput(importState, 'Faxa_tauy', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, UY0)
call FillGlobalInput(global_data, UYN)
end if
end if
! ---------------
! INFLAGS1(-7)
! ---------------
if (INFLAGS1(-7)) then
TI1 = timen ! time for ice field
ICEP1(:,:) = def_value ! ice thickness
if (state_fldchk(importState, 'Si_thick')) then
call SetGlobalInput(importState, 'Si_thick', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, ICEP1)
end if
end if
! ---------------
! INFLAGS1(-3)
! ---------------
if (INFLAGS1(-3)) then
TI5 = timen ! time for ice field
ICEP5(:,:) = def_value ! ice floe size
if (state_fldchk(importState, 'Si_floediam')) then
call SetGlobalInput(importState, 'Si_floediam', vm, global_data, rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call FillGlobalInput(global_data, ICEP5)
end if
end if
#ifndef W3_CESMCOUPLED
if (multigrid) then
do j = lbound(inflags1,1),ubound(inflags1,1)
if (inflags1(j)) then
do imod = 1,nrgrd
tfn(:,j) = timen(:)
call w3setg ( imod, mdse, mdst )
call w3setw ( imod, mdse, mdst )
call w3seti ( imod, mdse, mdst )
call w3seto ( imod, mdse, mdst )
call wmsetm ( imod, mdse, mdst )
#ifdef W3_MPI
if ( mpi_comm_grd .eq. mpi_comm_null ) cycle
#endif
!TODO: when is this active? jmod = -999
jmod = inpmap(imod,j)
if ( jmod.lt.0 .and. jmod.ne.-999 ) then
call wmupd2( imod, j, jmod, rc )
if (ChkErr(rc,__LINE__,u_FILE_u)) return
endif
end do
end if
end do
end if
#endif
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' done', ESMF_LOGMSG_INFO)
end subroutine import_fields
!===============================================================================
!> Fill the export state with values from WW3 fields
!!
!> @details Called by ModelAdvance, fill or compute the values in the export state.
!!
!! @param gcomp an ESMF_GridComp object
!! @param[out] rc return code
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine export_fields (gcomp, rc)
!---------------------------------------------------------------------------
! Create the export state
!---------------------------------------------------------------------------
use wav_kind_mod, only : R8 => SHR_KIND_R8
use w3adatmd , only : USSX, USSY, USSP
use w3adatmd , only : w3seta
use w3idatmd , only : w3seti
use w3wdatmd , only : va, w3setw
use w3odatmd , only : w3seto, naproc, iaproc
use w3gdatmd , only : nseal, mapsf, MAPSTA, USSPF, NK, w3setg
use w3iogomd , only : CALC_U3STOKES
#ifdef W3_CESMCOUPLED
use w3wdatmd , only : ASF, UST
use w3adatmd , only : USSHX, USSHY, UD, HS
use w3idatmd , only : HSL
#else
use wmmdatmd , only : mdse, mdst, wmsetm
#endif
! input/output/variables
type(ESMF_GridComp) :: gcomp
integer , intent(out) :: rc
! Local variables
#ifdef W3_CESMCOUPLED
real(R8) :: fillvalue = 1.0e30_R8 ! special missing value
real :: sww, langmt, lasl, laslpj, alphal
#else
real(R8) :: fillvalue = zero ! special missing value
#endif
type(ESMF_State) :: exportState
integer :: n, jsea, isea, ix, iy, ib
real(r8), pointer :: z0rlen(:)
real(r8), pointer :: charno(:)
real(r8), pointer :: wbcuru(:)
real(r8), pointer :: wbcurv(:)
real(r8), pointer :: wbcurp(:)
real(r8), pointer :: sxxn(:)
real(r8), pointer :: sxyn(:)
real(r8), pointer :: syyn(:)
real(r8), pointer :: sw_lamult(:)
real(r8), pointer :: sw_lasl(:)
real(r8), pointer :: sw_ustokes(:)
real(r8), pointer :: sw_vstokes(:)
! d2 is location, d1 is frequency - nwav_elev_spectrum frequencies will be used
real(r8), pointer :: wave_elevation_spectrum(:,:)
! Partitioned stokes drift
real(r8), pointer :: sw_pstokes_x(:,:)
real(r8), pointer :: sw_pstokes_y(:,:)
character(len=*), parameter :: subname='(wav_import_export:export_fields)'
!---------------------------------------------------------------------------
rc = ESMF_SUCCESS
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' called', ESMF_LOGMSG_INFO)
! Get export state
call NUOPC_ModelGet(gcomp, exportState=exportState, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
#ifndef W3_CESMCOUPLED
call w3setg ( 1, mdse, mdst )
call w3setw ( 1, mdse, mdst )
call w3seta ( 1, mdse, mdst )
call w3seti ( 1, mdse, mdst )
call w3seto ( 1, mdse, mdst )
if (multigrid) then
call wmsetm ( 1, mdse, mdst )
end if
#else
if (state_fldchk(exportState, 'Sw_lamult')) then
call state_getfldptr(exportState, 'Sw_lamult', sw_lamult, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
sw_lamult(:) = fillvalue
do jsea=1, nseal_cpl
call init_get_isea(isea, jsea)
ix = mapsf(isea,1)
iy = mapsf(isea,2)
if (mapsta(iy,ix) == 1 .and. HS(jsea) > zero .and. &
sqrt(USSX(jsea)**2+USSY(jsea)**2)>zero .and. sqrt(USSHX(jsea)**2+USSHY(jsea)**2)>zero ) then
sww = atan2(USSHY(jsea),USSHX(jsea)) - UD(isea)
alphal = atan( sin(sww) / ( &
2.5 * UST(isea)*ASF(isea)*sqrt(dair/dwat) &
/ max(1.e-14_r8, sqrt(USSX(jsea)**2+USSY(jsea)**2)) &
* log(max(1.0, abs(1.25*HSL(ix,iy)/HS(jsea)))) &
+ cos(sww) ) &
)
lasl = sqrt(ust(isea) * asf(isea) * sqrt(dair/dwat) &
/ sqrt(usshx(jsea)**2 + usshy(jsea)**2 ))
laslpj = lasl * sqrt(abs(cos(alphal)) &
/ abs(cos(sww-alphal)))
sw_lamult(jsea) = min(5.0, abs(cos(alphal)) * &
sqrt(1.0+(1.5*laslpj)**(-2)+(5.4_r8*laslpj)**(-4)))
else
sw_lamult(jsea) = 1.
endif
enddo
end if
if (state_fldchk(exportState, 'Sw_lasl')) then
call state_getfldptr(exportState, 'Sw_lasl', sw_lasl, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
sw_lasl(:) = fillvalue
do jsea=1, nseal
isea = iaproc + (jsea-1)*naproc
ix = mapsf(isea,1)
iy = mapsf(isea,2)
if (mapsta(iy,ix) == 1) then
! note: an arbitrary minimum value of 0.2 is set to avoid zero
! Langmuir number which may result from zero surface friction
! velocity but may cause unphysically strong Langmuir mixing
sw_lasl(jsea) = max(0.2, sqrt(UST(isea)*ASF(isea)*sqrt(dair/dwat) &
/ max(1.e-14, sqrt(USSHX(jsea)**2+USSHY(jsea)**2))))
else
sw_lasl(jsea) = 1.e6
endif
enddo
end if
#endif
! surface stokes drift
if (state_fldchk(exportState, 'Sw_ustokes')) then
call state_getfldptr(exportState, 'Sw_ustokes', sw_ustokes, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
sw_ustokes(:) = fillvalue
do jsea=1, nseal_cpl
call init_get_isea(isea, jsea)
ix = mapsf(isea,1)
iy = mapsf(isea,2)
if (mapsta(iy,ix) == 1) then
sw_ustokes(jsea) = USSX(jsea)
else
sw_ustokes(jsea) = 0.
endif
enddo
end if
if (state_fldchk(exportState, 'Sw_vstokes')) then
call state_getfldptr(exportState, 'Sw_vstokes', sw_vstokes, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
sw_vstokes(:) = fillvalue
do jsea=1, nseal_cpl
call init_get_isea(isea, jsea)
ix = mapsf(isea,1)
iy = mapsf(isea,2)
if (mapsta(iy,ix) == 1) then
sw_vstokes(jsea) = USSY(jsea)
else
sw_vstokes(jsea) = 0.
endif
enddo
end if
if (state_fldchk(exportState, 'Sw_ch')) then
call state_getfldptr(exportState, 'charno', charno, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call CalcCharnk(charno)
endif
if (state_fldchk(exportState, 'Sw_z0')) then
call state_getfldptr(exportState, 'Sw_z0', z0rlen, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call CalcRoughl(z0rlen)
endif
if ( state_fldchk(exportState, 'wbcuru') .and. &
state_fldchk(exportState, 'wbcurv') .and. &
state_fldchk(exportState, 'wbcurp')) then
call state_getfldptr(exportState, 'wbcuru', wbcuru, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_getfldptr(exportState, 'wbcurv', wbcurv, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_getfldptr(exportState, 'wbcurp', wbcurp, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call CalcBotcur( va, wbcuru, wbcurv, wbcurp)
end if
if ( state_fldchk(exportState, 'wavsuu') .and. &
state_fldchk(exportState, 'wavsuv') .and. &
state_fldchk(exportState, 'wavsvv')) then
call state_getfldptr(exportState, 'sxxn', sxxn, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_getfldptr(exportState, 'sxyn', sxyn, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_getfldptr(exportState, 'syyn', syyn, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call CalcRadstr2D( va, sxxn, sxyn, syyn)
end if
if (wav_coupling_to_cice) then
call state_getfldptr(exportState, 'Sw_elevation_spectrum', wave_elevation_spectrum, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
! Initialize wave elevation spectrum
wave_elevation_spectrum(:,:) = fillvalue
call CalcEF(va, wave_elevation_spectrum)
end if
if ( state_fldchk(exportState, 'Sw_pstokes_x') .and. &
state_fldchk(exportState, 'Sw_pstokes_y') )then
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_getfldptr(exportState, 'Sw_pstokes_x', sw_pstokes_x, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call state_getfldptr(exportState, 'Sw_pstokes_y', sw_pstokes_y, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
sw_pstokes_x(:,:) = fillvalue
sw_pstokes_y(:,:) = fillvalue
if (USSPF(1) > 0) then ! Partitioned Stokes drift computation is turned on in mod_def file.
call CALC_U3STOKES(va, 2)
do ib = 1, USSPF(2)
do jsea = 1, nseal_cpl
call init_get_isea(isea, jsea)
ix = mapsf(isea,1)
iy = mapsf(isea,2)
sw_pstokes_x(ib,jsea) = ussp(jsea,ib)
sw_pstokes_y(ib,jsea) = ussp(jsea,nk+ib)
enddo
end do
end if
endif
if (dbug_flag > 5) then
call state_diagnose(exportState, 'at export ', rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
end subroutine export_fields
!===============================================================================
!> Add a fieldname to a list of fields in a state
!!
!! @param[inout] num a counter for added fields
!! @param[inout] fldlist a structure for the standard name and ungridded dims
!! @param[in] stdname a standard field name
!! @param[in] ungridded_lbound the lower bound of an ungridded dimension
!! @param[in] ungridded_ubound the upper bound of an ungridded dimension
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine fldlist_add(num, fldlist, stdname, ungridded_lbound, ungridded_ubound)
integer, intent(inout) :: num
type(fld_list_type), intent(inout) :: fldlist(:)
character(len=*), intent(in) :: stdname
integer, optional, intent(in) :: ungridded_lbound
integer, optional, intent(in) :: ungridded_ubound
! local variables
character(len=*), parameter :: subname='(wav_import_export:fldlist_add)'
!-------------------------------------------------------------------------------
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' called', ESMF_LOGMSG_INFO)
! Set up a list of field information
num = num + 1
if (num > fldsMax) then
call ESMF_LogWrite(trim(subname)//": ERROR num > fldsMax "//trim(stdname), &
ESMF_LOGMSG_ERROR, line=__LINE__, file=__FILE__)
return
endif
fldlist(num)%stdname = trim(stdname)
if (present(ungridded_lbound) .and. present(ungridded_ubound)) then
fldlist(num)%ungridded_lbound = ungridded_lbound
fldlist(num)%ungridded_ubound = ungridded_ubound
end if
end subroutine fldlist_add
!===============================================================================
!> Realize a list of fields in a state
!!
!> @details For a connected field in a State, create an ESMF_Field object of
!! the required dimensionality on the ESMF_Mesh. Remove any unconnected fields from
!! the State. For a scalar field, create a field of dimensionality (1:flds_scalar_num)
!!
!! @param[inout] state an ESMF_State object
!! @param[in] fldlist a list of fields in the State
!! @param[in] numflds the number of fields in the state
!! @param[in] flds_scalar_name the name of the scalar field
!! @param[in] flds_scalar_num the count of scalar fields
!! @param[in] tag a character string for logging
!! @param[in] mesh an ESMF_Mesh object
!! @param[inout] rc a return code
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine fldlist_realize(state, fldList, numflds, flds_scalar_name, flds_scalar_num, mesh, tag, rc)
! input/output variables
type(ESMF_State) , intent(inout) :: state
type(fld_list_type) , intent(in) :: fldList(:)
integer , intent(in) :: numflds
character(len=*) , intent(in) :: flds_scalar_name
integer , intent(in) :: flds_scalar_num
character(len=*) , intent(in) :: tag
type(ESMF_Mesh) , intent(in) :: mesh
integer , intent(inout) :: rc
! local variables
integer :: n
type(ESMF_Field) :: field
character(len=80) :: stdname
character(len=*),parameter :: subname='(wav_import_export:fldlist_realize)'
! ----------------------------------------------
rc = ESMF_SUCCESS
if (dbug_flag > 5) call ESMF_LogWrite(trim(subname)//' called', ESMF_LOGMSG_INFO)
do n = 1, numflds
stdname = fldList(n)%stdname
if (NUOPC_IsConnected(state, fieldName=stdname)) then
if (stdname == trim(flds_scalar_name)) then
call ESMF_LogWrite(trim(subname)//trim(tag)//" Field = "//trim(stdname)//" is connected on root pe", &
ESMF_LOGMSG_INFO)
! Create the scalar field
call SetScalarField(field, flds_scalar_name, flds_scalar_num, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
else
call ESMF_LogWrite(trim(subname)//trim(tag)//" Field = "//trim(stdname)//" is connected using mesh", &
ESMF_LOGMSG_INFO)
! Create the field
if (fldlist(n)%ungridded_lbound > 0 .and. fldlist(n)%ungridded_ubound > 0) then
call ESMF_LogWrite(trim(subname)//trim(tag)//" Field = "//trim(stdname)//" has ungridded dimension", &
ESMF_LOGMSG_INFO)
field = ESMF_FieldCreate(mesh, ESMF_TYPEKIND_R8, name=stdname, meshloc=ESMF_MESHLOC_ELEMENT, &
ungriddedLbound=(/fldlist(n)%ungridded_lbound/), &
ungriddedUbound=(/fldlist(n)%ungridded_ubound/), &
gridToFieldMap=(/2/), rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
else
field = ESMF_FieldCreate(mesh, ESMF_TYPEKIND_R8, name=stdname, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
end if ! if not scalar field
! NOW call NUOPC_Realize
call NUOPC_Realize(state, field=field, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
else
if (stdname /= trim(flds_scalar_name)) then
call ESMF_LogWrite(subname // trim(tag) // " Field = "// trim(stdname) // " is not connected.", &
ESMF_LOGMSG_INFO)
call ESMF_StateRemove(state, (/stdname/), rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
end if
end if
end do
contains !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!> Create a field with scalar data on the root pe
!!
!! @param[inout] field an ESMF_Field
!! @param[in] flds_scalar_name the scalar field name
!! @param[in[ flds_scalar_num the dimnsionality of the scalar field
!! @param[inout] rc a return code
!!
!> @author mvertens@ucar.edu, Denise.Worthen@noaa.gov
!> @date 01-05-2022
subroutine SetScalarField(field, flds_scalar_name, flds_scalar_num, rc)
! ----------------------------------------------
! create a field with scalar data on the root pe
! ----------------------------------------------
type(ESMF_Field) , intent(inout) :: field
character(len=*) , intent(in) :: flds_scalar_name
integer , intent(in) :: flds_scalar_num
integer , intent(inout) :: rc
! local variables
type(ESMF_Distgrid) :: distgrid
type(ESMF_Grid) :: grid
character(len=*), parameter :: subname='(wav_import_export:SetScalarField)'
! ----------------------------------------------
rc = ESMF_SUCCESS
! create a DistGrid with a single index space element, which gets mapped onto DE 0.
distgrid = ESMF_DistGridCreate(minIndex=(/1/), maxIndex=(/1/), rc=rc)
if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
grid = ESMF_GridCreate(distgrid, rc=rc)
if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
field = ESMF_FieldCreate(name=trim(flds_scalar_name), grid=grid, typekind=ESMF_TYPEKIND_R8, &
ungriddedLBound=(/1/), ungriddedUBound=(/flds_scalar_num/), gridToFieldMap=(/2/), rc=rc) ! num of scalar values
if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
end subroutine SetScalarField
end subroutine fldlist_realize
!===============================================================================
!> Calculate Charnok parameter for export
!!
!> @details TODO:
!!
!! @param[inout] chkn a 1-D pointer to a field on a mesh
!!
!> @author T. J. Campbell, NRL
!> @date 09-Aug-2017
subroutine CalcCharnk ( chkn )
! Calculate Charnok for export
use w3gdatmd, only : nseal, nk, nth, sig, mapsf, mapsta, nspec
use w3adatmd, only : cg, wn, charn, u10, u10d
use w3wdatmd, only : va
use w3odatmd, only : naproc, iaproc
#ifdef W3_ST3
use w3src3md, only : w3spr3
#endif
#ifdef W3_ST4
use w3src4md, only : w3spr4
#endif
! input/output variables
real(ESMF_KIND_R8), pointer :: chkn(:) ! 1D Charnock export field pointer
! local variables
integer :: isea, jsea, ix, iy
real :: emean, fmean, fmean1, wnmean, amax, ustar, ustdr
real :: tauwx, tauwy, cd, z0, fmeanws, dlwmean
logical :: llws(nspec)
logical, save :: firstCall = .true.
!----------------------------------------------------------------------