forked from ESCOMP/CMEPS
-
Notifications
You must be signed in to change notification settings - Fork 21
/
Copy pathmed_map_mod.F90
1455 lines (1260 loc) · 64.7 KB
/
med_map_mod.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
module med_map_mod
use med_kind_mod , only : CX=>SHR_KIND_CX, CS=>SHR_KIND_CS, CL=>SHR_KIND_CL, R8=>SHR_KIND_R8
use med_kind_mod , only : I4=>SHR_KIND_I4
use ESMF , only : ESMF_SUCCESS, ESMF_FAILURE
use ESMF , only : ESMF_LOGMSG_ERROR, ESMF_LOGMSG_INFO, ESMF_LogWrite
use ESMF , only : ESMF_Field
use med_internalstate_mod , only : InternalState, logunit, mastertask
use med_constants_mod , only : dbug_flag => med_constants_dbug_flag
use med_utils_mod , only : chkerr => med_utils_ChkErr
use perf_mod , only : t_startf, t_stopf
implicit none
private
! public routines
public :: med_map_routehandles_init
public :: med_map_rh_is_created
public :: med_map_mapnorm_init
public :: med_map_packed_field_create
public :: med_map_field_packed
public :: med_map_field_normalized
public :: med_map_field
interface med_map_routehandles_init
module procedure med_map_routehandles_initfrom_esmflds
module procedure med_map_routehandles_initfrom_fieldbundle
module procedure med_map_routehandles_initfrom_field
end interface
interface med_map_RH_is_created
module procedure med_map_RH_is_created_RH3d
module procedure med_map_RH_is_created_RH1d
end interface
type(ESMF_Field) :: uv3d_src, uv3d_dst ! needed for 3d mapping of u,v vector pairs
! private module variables
character(*), parameter :: u_FILE_u = &
__FILE__
!================================================================================
contains
!================================================================================
subroutine med_map_RouteHandles_initfrom_esmflds(gcomp, llogunit, rc)
!---------------------------------------------
! Initialize route handles in the mediator
! Assumptions:
! - Route handles are created per target field bundles NOT
! per individual fields in the bundle
! - ALL fields in the bundle are on identical grids
! - MULTIPLE route handles are going to be generated for
! given field bundle source and destination grids
! - Route handles will ONLY be created if coupling is active
! between n1 and n2
! Algorithm
! n1=source component index
! n2=destination component index
! nf=field index
! fldListFr(n)%flds(nf) is queried to determine the mapindex and mapfile
! and the appropriate route handle is created
!
! Regridding is done on a per-field basis AND only for those fields that have a
! valid mapping index for the destination component
! n = source field index field index
! destcomp = destination component index
! The fldsSrc input argument is queried for the mapping type of the field
! for the desination component
! mapindex = fldsSrc(n)%mapindex(destcomp)
! If the mapindex is 0 (there is no valid mapping) then NO mapping is done
! for the field
!---------------------------------------------
use ESMF , only : ESMF_LogWrite, ESMF_LOGMSG_INFO, ESMF_SUCCESS, ESMF_LogFlush
use ESMF , only : ESMF_GridComp, ESMF_GridCompGet, ESMF_Field
use esmFlds , only : fldListFr, ncomps, mapunset, compname
use med_methods_mod , only : med_methods_FB_getFieldN
! input/output variables
type(ESMF_GridComp) :: gcomp
integer, intent(in) :: llogunit
integer, intent(out) :: rc
! local variables
type(InternalState) :: is_local
type(ESMF_Field) :: fldsrc
type(ESMF_Field) :: flddst
integer :: n,n1,n2,m,nf
character(len=CX) :: mapfile
integer :: mapindex
logical :: mapexists = .false.
character(len=*), parameter :: subname=' (module_med_map: RouteHandles_init) '
!-----------------------------------------------------------
call t_startf('MED:'//subname)
rc = ESMF_SUCCESS
if (dbug_flag > 1) then
call ESMF_LogWrite(trim(subname)//": start", ESMF_LOGMSG_INFO)
endif
! Get the internal state from Component.
nullify(is_local%wrap)
call ESMF_GridCompGetInternalState(gcomp, is_local, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Create the necessary route handles
! First loop over source and destination components components
if (mastertask) write(logunit,*) ' '
do n1 = 1, ncomps
do n2 = 1, ncomps
if (n1 /= n2) then
if (is_local%wrap%med_coupling_active(n1,n2)) then ! If coupling is active between n1 and n2
! Get source and destination fields
call med_methods_FB_getFieldN(is_local%wrap%FBImp(n1,n1), 1, fldsrc, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call med_methods_FB_getFieldN(is_local%wrap%FBImp(n1,n2), 1, flddst, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Loop over fields
do nf = 1,size(fldListFr(n1)%flds)
! Determine the mapping type for mapping field nf from n1 to n2
mapindex = fldListFr(n1)%flds(nf)%mapindex(n2)
if (mapindex /= mapunset) then
! determine if route handle has already been created
mapexists = med_map_RH_is_created(is_local%wrap%RH,n1,n2,mapindex,rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Create route handle for target mapindex if route handle is required
! (i.e. mapindex /= mapunset) and route handle has not already been created
if (.not. mapexists) then
mapfile = trim(fldListFr(n1)%flds(nf)%mapfile(n2))
call med_map_routehandles_initfrom_field(n1, n2, fldsrc, flddst, &
mapindex, is_local%wrap%rh(n1,n2,:), mapfile=trim(mapfile), rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end if
end if ! end if mapindex is mapunset
end do ! loop over fields
end if ! if coupling is active between n1 and n2
end if ! if n1 not equal to n2
end do ! loop over n2
end do ! loop over n1
if (dbug_flag > 1) then
call ESMF_LogWrite(trim(subname)//": done", ESMF_LOGMSG_INFO)
endif
call t_stopf('MED:'//subname)
end subroutine med_map_RouteHandles_initfrom_esmflds
!================================================================================
subroutine med_map_routehandles_initfrom_fieldbundle(n1, n2, FBsrc, FBdst, mapindex, RouteHandle, rc)
use ESMF , only : ESMF_LogWrite, ESMF_LOGMSG_INFO, ESMF_SUCCESS, ESMF_LogFlush
use ESMF , only : ESMf_Field, ESMF_FieldBundle, ESMF_RouteHandle
use med_methods_mod , only : med_methods_FB_getFieldN
!---------------------------------------------
! Initialize initialize additional route handles for mapping fractions
!---------------------------------------------
! input/output variables
integer , intent(in) :: n1
integer , intent(in) :: n2
type(ESMF_FieldBundle) , intent(in) :: FBsrc
type(ESMF_FieldBundle) , intent(in) :: fBdst
integer , intent(in) :: mapindex
type(ESMF_RouteHandle) , intent(inout) :: RouteHandle(:,:,:)
integer , intent(out) :: rc
! local variables
type(ESMF_Field) :: fldsrc
type(ESMF_Field) :: flddst
character(len=*), parameter :: subname=' (module_MED_map:med_map_routehandles_initfrom_fieldbundle) '
!---------------------------------------------
rc = ESMF_SUCCESS
call t_startf('MED:'//subname)
if (dbug_flag > 1) then
call ESMF_LogWrite(trim(subname)//": start", ESMF_LOGMSG_INFO)
endif
call med_methods_FB_getFieldN(FBsrc, 1, fldsrc, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call med_methods_FB_getFieldN(FBDst, 1, flddst, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call med_map_routehandles_initfrom_field(n1, n2, fldsrc, flddst, mapindex, routehandle(n1,n2,:), rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
if (dbug_flag > 1) then
call ESMF_LogWrite(trim(subname)//": done", ESMF_LOGMSG_INFO)
endif
call t_stopf('MED:'//subname)
end subroutine med_map_routehandles_initfrom_fieldbundle
!================================================================================
subroutine med_map_routehandles_initfrom_field(n1, n2, fldsrc, flddst, mapindex, routehandles, mapfile, rc)
use ESMF , only : ESMF_RouteHandle, ESMF_RouteHandlePrint, ESMF_Field, ESMF_MAXSTR
use ESMF , only : ESMF_PoleMethod_Flag, ESMF_POLEMETHOD_ALLAVG
use ESMF , only : ESMF_FieldSMMStore, ESMF_FieldRedistStore, ESMF_FieldRegridStore
use ESMF , only : ESMF_RouteHandleIsCreated
use ESMF , only : ESMF_REGRIDMETHOD_BILINEAR, ESMF_REGRIDMETHOD_PATCH
use ESMF , only : ESMF_REGRIDMETHOD_CONSERVE, ESMF_NORMTYPE_DSTAREA, ESMF_NORMTYPE_FRACAREA
use ESMF , only : ESMF_UNMAPPEDACTION_IGNORE, ESMF_REGRIDMETHOD_NEAREST_STOD
use ESMF , only : ESMF_EXTRAPMETHOD_NEAREST_STOD
use ESMF , only : ESMF_Mesh, ESMF_MeshLoc, ESMF_MESHLOC_ELEMENT, ESMF_TYPEKIND_I4
use ESMF , only : ESMF_MeshGet, ESMF_DistGridGet, ESMF_DistGrid, ESMF_TYPEKIND_R8
use ESMF , only : ESMF_FieldGet, ESMF_FieldCreate, ESMF_FieldWrite, ESMF_FieldDestroy
use esmFlds , only : mapbilnr, mapconsf, mapconsd, mappatch, mappatch_uv3d, mapbilnr_uv3d, mapfcopy
use esmFlds , only : mapunset, mapnames, nmappers
use esmFlds , only : mapnstod, mapnstod_consd, mapnstod_consf, mapnstod_consd
use esmFlds , only : mapfillv_bilnr, mapbilnr_nstod
use esmFlds , only : ncomps, compatm, compice, compocn, compname
use esmFlds , only : coupling_mode, dststatus_print
use esmFlds , only : atm_name
use med_constants_mod , only : ispval_mask => med_constants_ispval_mask
! input/output variables
integer , intent(in) :: n1
integer , intent(in) :: n2
type(ESMF_Field) , intent(inout) :: fldsrc
type(ESMF_Field) , intent(inout) :: flddst
integer , intent(in) :: mapindex
type(ESMF_RouteHandle) , intent(inout) :: routehandles(:)
character(len=*), optional , intent(in) :: mapfile
integer , intent(out) :: rc
! local variables
type(ESMF_Mesh) :: dstmesh
type(ESMF_Field) :: dststatusfield, doffield
type(ESMF_DistGrid) :: distgrid
character(len=CS) :: string
character(len=CS) :: mapname
character(len=CL) :: fname
integer :: srcMaskValue
integer :: dstMaskValue
character(len=ESMF_MAXSTR) :: lmapfile
logical :: rhprint = .false.
integer :: ns
integer(I4), pointer :: dof(:) => null()
integer :: srcTermProcessing_Value = 0
type(ESMF_PoleMethod_Flag), parameter :: polemethod=ESMF_POLEMETHOD_ALLAVG
character(len=*), parameter :: subname=' (module_med_map: med_map_routehandles_initfrom_field) '
!---------------------------------------------
lmapfile = 'unset'
if (present(mapfile)) then
lmapfile = trim(mapfile)
end if
mapname = trim(mapnames(mapindex))
call ESMF_LogWrite(trim(subname)//": mapname "//trim(mapname), ESMF_LOGMSG_INFO)
! create a field to retrieve the dststatus field
call ESMF_FieldGet(flddst, mesh=dstmesh, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
dststatusfield = ESMF_FieldCreate(dstmesh, ESMF_TYPEKIND_I4, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
if (trim(coupling_mode) == 'cesm') then
dstMaskValue = ispval_mask
srcMaskValue = ispval_mask
if (n1 == compocn .or. n1 == compice) srcMaskValue = 0
if (n2 == compocn .or. n2 == compice) dstMaskValue = 0
else if (coupling_mode(1:4) == 'nems') then
if (n1 == compatm .and. (n2 == compocn .or. n2 == compice)) then
srcMaskValue = 1
dstMaskValue = 0
if (atm_name(1:4).eq.'datm') then
srcMaskValue = 0
endif
else if (n2 == compatm .and. (n1 == compocn .or. n1 == compice)) then
srcMaskValue = 0
dstMaskValue = 1
else if ((n1 == compocn .and. n2 == compice) .or. (n1 == compice .and. n2 == compocn)) then
srcMaskValue = 0
dstMaskValue = 0
else
! TODO: what should the condition be here?
dstMaskValue = ispval_mask
srcMaskValue = ispval_mask
end if
else if (trim(coupling_mode) == 'hafs') then
dstMaskValue = ispval_mask
srcMaskValue = ispval_mask
if (n1 == compocn .or. n1 == compice) srcMaskValue = 0
if (n2 == compocn .or. n2 == compice) dstMaskValue = 0
if (n1 == compatm .and. n2 == compocn) then
if (trim(atm_name).ne.'datm') then
srcMaskValue = 1
endif
dstMaskValue = 0
elseif (n1 == compocn .and. n2 == compatm) then
srcMaskValue = 0
dstMaskValue = ispval_mask
endif
end if
write(string,'(a)') trim(compname(n1))//' to '//trim(compname(n2))
! Create route handle
if (mapindex == mapfcopy) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH redist for '//trim(string)
end if
call ESMF_FieldRedistStore(fldsrc, flddst, routehandle=routehandles(mapfcopy), &
ignoreUnmatchedIndices = .true., rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (lmapfile /= 'unset') then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//&
' via input file '//trim(mapfile)//' for '//trim(string)
end if
call ESMF_FieldSMMStore(fldsrc, flddst, mapfile, routehandle=routehandles(mapindex), &
ignoreUnmatchedIndices=.true., &
srcTermProcessing=srcTermProcessing_Value, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (mapindex == mapbilnr .or. mapindex == mapbilnr_uv3d) then
if (.not. ESMF_RouteHandleIsCreated(routehandles(mapbilnr))) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//' for '//trim(string)
end if
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mapbilnr), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_BILINEAR, &
polemethod=polemethod, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end if
else if (mapindex == mapfillv_bilnr) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//' for '//trim(string)
end if
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mapfillv_bilnr), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_BILINEAR, &
polemethod=polemethod, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (mapindex == mapbilnr_nstod) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//' for '//trim(string)
end if
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mapbilnr_nstod), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_BILINEAR, &
extrapMethod=ESMF_EXTRAPMETHOD_NEAREST_STOD, &
polemethod=polemethod, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (mapindex == mapconsf .or. mapindex == mapnstod_consf) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//' for '//trim(string)
end if
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mapconsf), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_CONSERVE, &
normType=ESMF_NORMTYPE_FRACAREA, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (mapindex == mapconsd .or. mapindex == mapnstod_consd) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//' for '//trim(string)
end if
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mapconsd), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_CONSERVE, &
normType=ESMF_NORMTYPE_DSTAREA, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (mapindex == mappatch .or. mapindex == mappatch_uv3d) then
if (.not. ESMF_RouteHandleIsCreated(routehandles(mappatch))) then
if (mastertask) then
write(logunit,'(A)') trim(subname)//' creating RH '//trim(mapname)//' for '//trim(string)
end if
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mappatch), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_PATCH, &
polemethod=polemethod, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end if
else
if (mastertask) then
write(logunit,'(A)') trim(subname)//' mapindex '//trim(mapname)//' not supported for '//trim(string)
end if
call ESMF_LogWrite(trim(subname)//' mapindex '//trim(mapname)//' not supported ', &
ESMF_LOGMSG_ERROR, line=__LINE__, file=u_FILE_u)
rc = ESMF_FAILURE
return
end if
! Output destination status field to file if requested
if (dststatus_print) then
if (mapindex /= mapfcopy .or. lmapfile /= 'unset') then
fname = 'dststatus.'//trim(compname(n1))//'.'//trim(compname(n2))//'.'//trim(mapname)//'.nc'
call ESMF_LogWrite(trim(subname)//": writing dstStatusField to "//trim(fname), ESMF_LOGMSG_INFO)
call ESMF_FieldWrite(dststatusfield, filename=trim(fname), variableName='dststatus', &
overwrite=.true., rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! the sequence index in order to sort the dststatus field
call ESMF_MeshGet(dstmesh, elementDistgrid=distgrid, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_DistGridGet(distgrid, localDE=0, elementCount=ns, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
allocate(dof(ns))
call ESMF_DistGridGet(distgrid, localDE=0, seqIndexList=dof, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
doffield = ESMF_FieldCreate(dstmesh, dof, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_FieldWrite(doffield, fileName='dof.'//trim(compname(n2))//'.nc', variableName='dof', &
overwrite=.true., rc=rc)
deallocate(dof)
call ESMF_FieldDestroy(doffield, rc=rc, noGarbage=.true.)
end if
end if
! consd_nstod method requires a second routehandle
if (mapindex == mapnstod .or. mapindex == mapnstod_consd .or. mapindex == mapnstod_consf) then
call ESMF_FieldRegridStore(fldsrc, flddst, routehandle=routehandles(mapnstod), &
srcMaskValues=(/srcMaskValue/), &
dstMaskValues=(/dstMaskValue/), &
regridmethod=ESMF_REGRIDMETHOD_NEAREST_STOD, &
srcTermProcessing=srcTermProcessing_Value, &
ignoreDegenerate=.true., &
dstStatusField=dststatusfield, &
unmappedaction=ESMF_UNMAPPEDACTION_IGNORE, &
rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Output destination status field to file if requested
if (dststatus_print) then
fname = 'dststatus.'//trim(compname(n1))//'.'//trim(compname(n2))//'.'//trim(mapname)//'_2.nc'
call ESMF_LogWrite(trim(subname)//": writing dstStatusField to "//trim(fname), ESMF_LOGMSG_INFO)
call ESMF_FieldWrite(dststatusfield, filename=trim(fname), variableName='dststatus', overwrite=.true., rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end if
end if
! Check that a valid route handle has been created
! TODO: should this be implemented as an error check or ignored?
! if (.not. med_map_RH_is_created(routehandle ,rc=rc)) then
! string = trim(compname(n1))//"2"//trim(compname(n2))//'_weights'
! call ESMF_LogWrite(trim(subname)//trim(string)//": failed RH "//trim(mapnames(mapindex)), &
! ESMF_LOGMSG_INFO)
! endif
! Output route handle to file if requested
if (rhprint) then
if (mastertask) then
write(logunit,'(a)') trim(subname)//trim(string)//": printing RH for "//trim(mapname)
end if
call ESMF_RouteHandlePrint(routehandles(mapindex), rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
endif
call ESMF_FieldDestroy(dststatusfield, rc=rc, noGarbage=.true.)
end subroutine med_map_routehandles_initfrom_field
!================================================================================
logical function med_map_RH_is_created_RH3d(RHs,n1,n2,mapindex,rc)
use ESMF, only : ESMF_RouteHandle
! input/output variables
type(ESMF_RouteHandle) , intent(in) :: RHs(:,:,:)
integer , intent(in) :: n1
integer , intent(in) :: n2
integer , intent(in) :: mapindex
integer , intent(out) :: rc
! local variables
integer :: rc1, rc2
character(len=*), parameter :: subname=' (module_MED_map:med_map_RH_is_created_RH3d) '
!-----------------------------------------------------------
rc = ESMF_SUCCESS
med_map_RH_is_created_RH3d = med_map_RH_is_created_RH1d(RHs(n1,n2,:),mapindex,rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end function med_map_RH_is_created_RH3d
!================================================================================
logical function med_map_RH_is_created_RH1d(RHs,mapindex,rc)
use ESMF , only : ESMF_RouteHandle, ESMF_RouteHandleIsCreated
use esmFlds , only : mapconsd, mapconsf, mapnstod
use esmFlds , only : mapnstod_consd, mapnstod_consf
! input/output varaibes
type(ESMF_RouteHandle) , intent(in) :: RHs(:)
integer , intent(in) :: mapindex
integer , intent(out) :: rc
! local variables
integer :: rc1, rc2
logical :: mapexists
character(len=*), parameter :: subname=' (module_MED_map:med_map_RH_is_created_RH1d) '
!-----------------------------------------------------------
rc = ESMF_SUCCESS
rc1 = ESMF_SUCCESS
rc2 = ESMF_SUCCESS
mapexists = .false.
if (mapindex == mapnstod_consd .and. &
ESMF_RouteHandleIsCreated(RHs(mapnstod), rc=rc1) .and. &
ESMF_RouteHandleIsCreated(RHs(mapconsd), rc=rc2)) then
rc = rc1
if (chkerr(rc,__LINE__,u_FILE_u)) return
rc = rc2
if (chkerr(rc,__LINE__,u_FILE_u)) return
mapexists = .true.
else if (mapindex == mapnstod_consf .and. &
ESMF_RouteHandleIsCreated(RHs(mapnstod), rc=rc1) .and. &
ESMF_RouteHandleIsCreated(RHs(mapconsf), rc=rc2)) then
rc = rc1
if (chkerr(rc,__LINE__,u_FILE_u)) return
rc = rc2
if (chkerr(rc,__LINE__,u_FILE_u)) return
mapexists = .true.
else if (ESMF_RouteHandleIsCreated(RHs(mapindex), rc=rc1)) then
rc = rc1
if (chkerr(rc,__LINE__,u_FILE_u)) return
mapexists = .true.
end if
med_map_RH_is_created_RH1d = mapexists
end function med_map_RH_is_created_RH1d
!================================================================================
subroutine med_map_mapnorm_init(gcomp, rc)
!---------------------------------------
! Initialize unity normalization fields and do the mapping for unity normalization up front
!---------------------------------------
use ESMF , only: ESMF_LogWrite, ESMF_LOGMSG_INFO, ESMF_SUCCESS, ESMF_LogFlush
use ESMF , only: ESMF_GridComp
use ESMF , only: ESMF_Mesh, ESMF_TYPEKIND_R8, ESMF_MESHLOC_ELEMENT
use ESMF , only: ESMF_FieldBundle, ESMF_FieldBundleGet, ESMF_FieldBundleCreate
use ESMF , only: ESMF_FieldBundleIsCreated
use ESMF , only: ESMF_Field, ESMF_FieldGet, ESMF_FieldCreate, ESMF_FieldDestroy
use esmFlds , only: ncomps, nmappers, compname, mapnames
use med_constants_mod , only: czero => med_constants_czero
! input/output variables
type(ESMF_GridComp) :: gcomp
integer, intent(out) :: rc
! local variables
type(InternalState) :: is_local
integer :: n1, n2, m
real(R8), pointer :: dataptr(:) => null()
integer :: fieldCount
type(ESMF_Field), pointer :: fieldlist(:) => null()
type(ESMF_Field) :: field_src
type(ESMF_Mesh) :: mesh_src
type(ESMF_Mesh) :: mesh_dst
character(len=*),parameter :: subname=' (module_MED_MAP:MapNorm_init)'
!-----------------------------------------------------------
call t_startf('MED:'//subname)
rc = ESMF_SUCCESS
if (dbug_flag > 1) then
call ESMF_LogWrite(trim(subname)//": start", ESMF_LOGMSG_INFO)
endif
if (mastertask) then
write(logunit,*)
write(logunit,'(a)') trim(subname)//"Initializing unity map normalizations"
endif
! Get the internal state from Component.
nullify(is_local%wrap)
call ESMF_GridCompGetInternalState(gcomp, is_local, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Create the destination normalization field
do n1 = 1,ncomps
if (ESMF_FieldBundleIsCreated(is_local%wrap%FBImp(n1,n1))) then
! Get source mesh
call ESMF_FieldBundleGet(is_local%wrap%FBImp(n1,n1), fieldCount=fieldCount, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
allocate(fieldlist(fieldcount))
call ESMF_FieldBundleGet(is_local%wrap%FBImp(n1,n1), fieldlist=fieldlist, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call ESMF_FieldGet(fieldlist(1), mesh=mesh_src, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
field_src = ESMF_FieldCreate(mesh_src, ESMF_TYPEKIND_R8, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_FieldGet(field_src, farrayptr=dataPtr, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
dataptr(:) = 1.0_R8
do n2 = 1,ncomps
if ( n1 /= n2 .and. &
ESMF_FieldBundleIsCreated(is_local%wrap%FBImp(n1,n2)) .and. &
is_local%wrap%med_coupling_active(n1,n2) ) then
! Get destination mesh
call ESMF_FieldBundleGet(is_local%wrap%FBImp(n1,n2), fieldlist=fieldlist, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
call ESMF_FieldGet(fieldlist(1), mesh=mesh_dst, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Create is_local%wrap%field_NormOne(n1,n2,m)
do m = 1,nmappers
if (med_map_RH_is_created(is_local%wrap%RH,n1,n2,m,rc=rc)) then
is_local%wrap%field_NormOne(n1,n2,m) = ESMF_FieldCreate(mesh_dst, &
ESMF_TYPEKIND_R8, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_FieldGet(is_local%wrap%field_NormOne(n1,n2,m), farrayptr=dataptr, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
dataptr(:) = czero
call med_map_field( &
field_src=field_src, &
field_dst=is_local%wrap%field_NormOne(n1,n2,m), &
routehandles=is_local%wrap%RH(n1,n2,:), &
maptype=m, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
if (mastertask) then
write(logunit,'(a)') trim(subname)//' created field_NormOne for '&
//compname(n1)//'->'//compname(n2)//' with mapping '//mapnames(m)
endif
end if
end do ! end of loop over m mappers
end if ! end of if block for creating destination field
end do ! end of loop over n2
! Deallocate memory
deallocate(fieldlist)
call ESMF_FieldDestroy(field_src, rc=rc, noGarbage=.true.)
if (chkerr(rc,__LINE__,u_FILE_u)) return
end if ! end of if-block for existence of field bundle
end do ! end of loop over n1
if (dbug_flag > 1) then
call ESMF_LogWrite(trim(subname)//": done", ESMF_LOGMSG_INFO)
endif
call t_stopf('MED:'//subname)
end subroutine med_map_mapnorm_init
!================================================================================
subroutine med_map_packed_field_create(destcomp, flds_scalar_name, &
fldsSrc, FBSrc, FBDst, packed_data, rc)
use ESMF
use esmFlds , only : med_fldList_entry_type, nmappers
use esmFlds , only : ncomps, compatm, compice, compocn, compname, mapnames
use med_internalstate_mod , only : packed_data_type
! input/output variables
integer , intent(in) :: destcomp
character(len=*) , intent(in) :: flds_scalar_name
type(med_fldList_entry_type) , pointer :: fldsSrc(:) ! array over mapping types
type(ESMF_FieldBundle) , intent(in) :: FBSrc
type(ESMF_FieldBundle) , intent(inout) :: FBDst
type(packed_data_type) , intent(inout) :: packed_data(:) ! array over mapping types
integer , intent(out) :: rc
! local variables
integer :: nf, nu, ns
integer, allocatable :: npacked(:)
integer :: fieldcount
type(ESMF_Field) :: lfield
integer :: ungriddedUBound(1) ! currently the size must equal 1 for rank 2 fields
real(r8), pointer :: ptrsrc_packed(:,:) => null()
real(r8), pointer :: ptrdst_packed(:,:) => null()
integer :: lsize_src
integer :: lsize_dst
type(ESMF_Mesh) :: lmesh_src
type(ESMF_Mesh) :: lmesh_dst
integer :: mapindex
type(ESMF_Field), pointer :: fieldlist_src(:) => null()
type(ESMF_Field), pointer :: fieldlist_dst(:) => null()
character(CL), allocatable :: fieldNameList(:)
character(CS) :: mapnorm_mapindex
character(len=CX) :: tmpstr
character(len=*), parameter :: subname=' (module_MED_map:med_packed_field_create) '
!-----------------------------------------------------------
rc = ESMF_SUCCESS
! Get field count for both FBsrc and FBdst
call ESMF_FieldBundleGet(FBsrc, fieldCount=fieldCount, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
! get fields in source and destination field bundles
allocate(fieldlist_src(fieldcount))
call ESMF_FieldBundleGet(FBsrc, fieldlist=fieldlist_src, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
allocate(fieldlist_dst(fieldcount))
call ESMF_FieldBundleGet(FBdst, fieldlist=fieldlist_dst, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
! field names are the same for the source and destination field bundles
allocate(fieldnamelist(fieldcount))
call ESMF_FieldBundleGet(FBsrc, fieldnamelist=fieldnamelist, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
! Determine local size and mesh of source fields
! Allocate a source fortran pointer for the new packed field bundle
call ESMF_FieldGet(fieldlist_src(1), mesh=lmesh_src, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_MeshGet(lmesh_src, numOwnedElements=lsize_src, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Determine local size of destination fields
! Allocate a destination fortran pointer for the new packed field bundle
call ESMF_FieldGet(fieldlist_dst(1), mesh=lmesh_dst, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call ESMF_MeshGet(lmesh_dst, numOwnedElements=lsize_dst, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Gather all fields that will be mapped with a target map index into a packed field
! Calculated size of packed field based on the fact that some fields have
! ungridded dimensions and need to unwrap them into separate fields for the
! purposes of packing
if (mastertask) write(logunit,*)
! Determine the normalization type for each packed_data mapping element
! Loop over mapping types
do mapindex = 1,nmappers
mapnorm_mapindex = 'not_set'
! Loop over source field bundle
do nf = 1, fieldCount
! Loop over the fldsSrc types
do ns = 1,size(fldsSrc)
! Note that fieldnamelist is an array of names for the source fields
! The assumption is that there is only one mapping normalization
! for any given mapping type
if ( fldsSrc(ns)%mapindex(destcomp) == mapindex .and. &
trim(fldsSrc(ns)%shortname) == trim(fieldnamelist(nf))) then
! Set the normalization to the input
packed_data(mapindex)%mapnorm = fldsSrc(ns)%mapnorm(destcomp)
if (mapnorm_mapindex == 'not_set') then
mapnorm_mapindex = packed_data(mapindex)%mapnorm
write(tmpstr,*)'Map type '//trim(mapnames(mapindex)) &
//', destcomp '//trim(compname(destcomp)) &
//', mapnorm '//trim(mapnorm_mapindex) &
//' '//trim(fieldnamelist(nf))
call ESMF_LogWrite(trim(tmpstr), ESMF_LOGMSG_INFO)
else
if (mapnorm_mapindex /= packed_data(mapindex)%mapnorm) then
write(tmpstr,*)'Map type '//trim(mapnames(mapindex)) &
//', destcomp '//trim(compname(destcomp)) &
//', mapnorm '//trim(mapnorm_mapindex) &
//' set; cannot set mapnorm to '//trim(packed_data(mapindex)%mapnorm) &
//' '//trim(fieldnamelist(nf))
call ESMF_LogWrite(trim(tmpstr), ESMF_LOGMSG_INFO)
call ESMF_Finalize(endflag=ESMF_END_ABORT)
end if
end if
end if
end do
end do
end do
! Allocate memory to keep tracked of packing index for each mapping type
allocate(npacked(nmappers))
npacked(:) = 0
! Loop over mapping types
do mapindex = 1,nmappers
! Allocate the fldindex attribute of packed_indices if needed
if (.not. allocated(packed_data(mapindex)%fldindex)) then
allocate(packed_data(mapindex)%fldindex(fieldcount))
packed_data(mapindex)%fldindex(:) = -999
end if
! Loop over the fields in FBSrc
do nf = 1, fieldCount
! Loop over the fldsSrc types
do ns = 1,size(fldsSrc)
if ( fldsSrc(ns)%mapindex(destcomp) == mapindex .and. &
trim(fldsSrc(ns)%shortname) == trim(fieldnamelist(nf))) then
! Determine mapping of indices into packed field bundle
! Get source field
call ESMF_FieldGet(fieldlist_src(nf), ungriddedUBound=ungriddedUBound, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
if (ungriddedUBound(1) > 0) then
do nu = 1,ungriddedUBound(1)
npacked(mapindex) = npacked(mapindex) + 1
if (nu == 1) then
packed_data(mapindex)%fldindex(nf) = npacked(mapindex)
end if
end do
else
npacked(mapindex) = npacked(mapindex) + 1
packed_data(mapindex)%fldindex(nf) = npacked(mapindex)
end if
if (mastertask) then
write(logunit,'(5(a,2x),2x,i4)') trim(subname)//&
'Packed field: destcomp,mapping,mapnorm,fldname,index: ', &
trim(compname(destcomp)), &
trim(mapnames(mapindex)), &
trim(packed_data(mapindex)%mapnorm), &
trim(fieldnamelist(nf)), &
packed_data(mapindex)%fldindex(nf)
end if
end if! end if source field is mapped to destination field with mapindex
end do ! end loop over FBSrc fields
end do ! end loop over fldsSrc elements
if (npacked(mapindex) > 0) then
! Create the packed source field bundle for mapindex
allocate(ptrsrc_packed(npacked(mapindex), lsize_src))
packed_data(mapindex)%field_src = ESMF_FieldCreate(lmesh_src, &
ptrsrc_packed, gridToFieldMap=(/2/), meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Create the packed destination field bundle for mapindex
allocate(ptrdst_packed(npacked(mapindex), lsize_dst))
packed_data(mapindex)%field_dst = ESMF_FieldCreate(lmesh_dst, &
ptrdst_packed, gridToFieldMap=(/2/), meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
packed_data(mapindex)%field_fracsrc = ESMF_FieldCreate(lmesh_src, ESMF_TYPEKIND_R8, &
meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
packed_data(mapindex)%field_fracdst = ESMF_FieldCreate(lmesh_dst, ESMF_TYPEKIND_R8, &
meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
end if
end do ! end loop over mapindex
deallocate(npacked)
deallocate(fieldlist_src)
deallocate(fieldlist_dst)
end subroutine med_map_packed_field_create
!================================================================================
subroutine med_map_field_packed(FBSrc, FBDst, FBFracSrc, field_normOne, packed_data, routehandles, rc)
! -----------------------------------------------
! Do regridding via packed field bundles
! -----------------------------------------------
use ESMF , only : ESMF_Field, ESMF_FieldGet, ESMF_FieldIsCreated
use ESMF , only : ESMF_FieldBundle, ESMF_FieldBundleGet
use ESMF , only : ESMF_FieldBundleIsCreated
use ESMF , only : ESMF_FieldRedist, ESMF_RouteHandle
use esmFlds , only : nmappers, mapfcopy
use esmFlds , only : mappatch_uv3d, mappatch, mapbilnr_uv3d, mapbilnr
use med_internalstate_mod , only : packed_data_type
! input/output variables
type(ESMF_FieldBundle) , intent(in) :: FBSrc
type(ESMF_FieldBundle) , intent(inout) :: FBDst
type(ESMF_Field) , intent(in) :: field_normOne(:) ! array over mapping types
type(ESMF_FieldBundle) , intent(in) :: FBFracSrc ! fraction field bundle for source
type(packed_data_type) , intent(inout) :: packed_data(:) ! array over mapping types
type(ESMF_RouteHandle) , intent(inout) :: routehandles(:)
integer , intent(out) :: rc
! local variables
integer :: nf, nu, np, n
integer :: fieldcount
integer :: mapindex
integer :: ungriddedUBound(1)
real(r8), pointer :: dataptr1d(:) => null()
real(r8), pointer :: dataptr2d(:,:) => null()
real(r8), pointer :: dataptr2d_packed(:,:) => null()
type(ESMF_Field) :: lfield
type(ESMF_Field) :: field_fracsrc
type(ESMF_Field), pointer :: fieldlist_src(:) => null()
type(ESMF_Field), pointer :: fieldlist_dst(:) => null()
type(ESMF_Field) :: usrc, vsrc ! only used for 3d mapping of u,v
type(ESMF_Field) :: udst, vdst ! only used for 3d mapping of u,v
real(r8), pointer :: data_norm(:) => null()
real(r8), pointer :: data_dst(:,:) => null()
character(len=*), parameter :: subname=' (module_MED_map:med_map_field_packed) '
!-----------------------------------------------------------
call t_startf('MED:'//subname)
rc = ESMF_SUCCESS
! Get field count for both FBsrc and FBdst
if (ESMF_FieldBundleIsCreated(FBsrc)) then
call ESMF_FieldBundleGet(FBsrc, fieldCount=fieldCount, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
allocate(fieldlist_src(fieldcount))
call ESMF_FieldBundleGet(FBsrc, fieldlist=fieldlist_src, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
allocate(fieldlist_dst(fieldcount))
call ESMF_FieldBundleGet(FBdst, fieldlist=fieldlist_dst, rc=rc)
if (ChkErr(rc,__LINE__,u_FILE_u)) return
else
fieldcount=0
endif
! Loop over mapping types
do mapindex = 1,nmappers
if (ESMF_FieldIsCreated(packed_data(mapindex)%field_src)) then
if (mapindex == mappatch_uv3d) then
! For mappatch_uv3d do not use packed field bundles
call med_map_uv_cart3d(FBsrc, FBdst, routehandles, mappatch, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else if (mapindex == mapbilnr_uv3d) then
! For mapbilnr_uv3d do not use packed field bundles
call med_map_uv_cart3d(FBsrc, FBdst, routehandles, mapbilnr, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
else
! -----------------------------------
! Copy the src fields into the packed field bundle
! -----------------------------------
call t_startf('MED:'//trim(subname)//' copy from src')
! First get the pointer for the packed source data
call ESMF_FieldGet(packed_data(mapindex)%field_src, farrayptr=dataptr2d_packed, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
! Now do the copy
do nf = 1,fieldcount
! Get the indices into the packed data structure
np = packed_data(mapindex)%fldindex(nf)
if (np > 0) then
call ESMF_FieldGet(fieldlist_src(nf), ungriddedUBound=ungriddedUBound, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
if (ungriddedUBound(1) > 0) then
call ESMF_FieldGet(fieldlist_src(nf), farrayptr=dataptr2d, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
do nu = 1,ungriddedUBound(1)
dataptr2d_packed(np+nu-1,:) = dataptr2d(nu,:)
end do
else
call ESMF_FieldGet(fieldlist_src(nf), farrayptr=dataptr1d, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
dataptr2d_packed(np,:) = dataptr1d(:)
end if
end if
end do
call t_stopf('MED:'//trim(subname)//' copy from src')
! -----------------------------------
! Do the mapping
! -----------------------------------