-
Notifications
You must be signed in to change notification settings - Fork 137
/
xgrid.F90
5455 lines (4815 loc) · 213 KB
/
xgrid.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 GFDL Flexible Modeling System (FMS).
!*
!* FMS 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.
!*
!* FMS is distributed in the hope that it will be useful, but WITHOUT
!* ANY WARRANTY; 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 FMS. If not, see <http://www.gnu.org/licenses/>.
!***********************************************************************
!> @defgroup xgrid_mod xgrid_mod
!> @ingroup exchange
!> @brief Implements exchange grids for coupled models running on multiple processors
!> @author Michael Winton, Zhi Liang
!!
!! An exchange grid is formed from the union of
!! the bounding lines of the two (logically rectangular) participating
!! grids. The exchange grid is therefore the coarsest grid that is a
!! refinement of both participating grids. Exchange grids are used for
!! two purposes by coupled models:
!! 1. conservative interpolation of fields
!! between models uses the exchange grid cell areas as weights and
!! 2. the surface flux calculation takes place on the exchange grid thereby
!! using the finest scale data available.
!! <TT>xgrid_mod</TT> uses a NetCDF grid
!! specification file containing the grid cell overlaps in combination with
!! the @link ftp://ftp.gfdl.gov/pub/vb/mpp/mpp_domains.F90 @endlink domain
!! decomposition information to determine
!! the grid and processor connectivities.
!!
!!
!! xgrid_mod - implements exchange grids. An exchange grid is the grid whose
!! boundary set is the union of the boundaries of the participating
!! grids. The exchange grid is the coarsest grid that is a
!! refinement of each of the participating grids. Every exchange
!! grid cell is a subarea of one and only one cell in each of the
!! participating grids. The exchange grid has two purposes:
!!
!! (1) The exchange cell areas are used as weights for
!! conservative interpolation between model grids.
!!
!! (2) Computation of surface fluxes takes place on it,
!! thereby using the finest scale data obtainable.
!!
!! The exchange cells are the 2D intersections between cells of the
!! participating grids. They are computed elsewhere and are
!! read here from a NetCDF grid file as a sequence of quintuples
!! (i and j on each of two grids and the cell area).
!!
!! Each processing element (PE) computes a subdomain of each of the
!! participating grids as well as a subset of the exchange cells.
!! The geographic regions corresponding to these subdomains will,
!! in general, not be the same so communication must occur between
!! the PEs. The scheme for doing this is as follows. A distinction
!! is drawn between the participating grids. There is a single
!! "side 1" grid and it does not have partitions (sub-grid surface
!! types). There are one or more "side 2" grids and they may have
!! more than 1 partition. In standard usage, the atmosphere grid is
!! on side 1 and the land and sea ice grids are on side 2. The set
!! of exchange cells computed on a PE corresponds to its side 2
!! geographic region(s). Communication between the PEs takes place
!! on the side 1 grid. Note: this scheme does not generally allow
!! reproduction of answers across varying PE counts. This is
!! because, in the side 1 "get", exchange cells are first summed
!! locally onto a side 1 grid, then these side 1 contributions are
!! further summed after they have been communicated to their target
!! PE. For the make_exchange_reproduce option, a special side 1 get
!! is used. This get communicates individual exchange cells. The
!! cells are summed in the order they appear in the grid spec. file.
!!
!! <TT>xgrid_mod</TT> reads a NetCDF grid specification file to determine the
!! grid and processor connectivities. The exchange grids are defined
!! by a sequence of quintuples: the <TT>i/j</TT> indices of the intersecting
!! cells of the two participating grids and their areal overlap.
!! The names of the five fields are generated automatically from the
!! three character ids of the participating grids. For example, if
!! the side one grid id is "ATM" and the side two grid id is "OCN",
!! <TT>xgrid_mod</TT> expects to find the following five fields in the grid
!! specification file: <TT>I_ATM_ATMxOCN, J_ATM_ATMxOCN, I_OCN_ATMxOCN,
!! J_OCN_ATMxOCN, and AREA_ATMxOCN</TT>. These fields may be generated
!! by the <TT>make_xgrids</TT> utility.
!> @addtogroup xgrid_mod
!> @{
module xgrid_mod
use fms_mod, only: check_nml_error, &
error_mesg, FATAL, NOTE, stdlog, &
write_version_number, lowercase, string
use mpp_mod, only: mpp_npes, mpp_pe, mpp_root_pe, mpp_send, mpp_recv, &
mpp_sync_self, stdout, mpp_max, EVENT_RECV, &
mpp_get_current_pelist, mpp_clock_id, mpp_min, &
mpp_alltoall, &
mpp_clock_begin, mpp_clock_end, MPP_CLOCK_SYNC, &
COMM_TAG_1, COMM_TAG_2, COMM_TAG_3, COMM_TAG_4, &
COMM_TAG_5, COMM_TAG_6, COMM_TAG_7, COMM_TAG_8, &
COMM_TAG_9, COMM_TAG_10
use mpp_mod, only: input_nml_file, mpp_set_current_pelist, mpp_sum, mpp_sync
use mpp_domains_mod, only: mpp_get_compute_domain, mpp_get_compute_domains, &
Domain2d, mpp_global_sum, mpp_update_domains, &
mpp_modify_domain, mpp_get_data_domain, XUPDATE, &
YUPDATE, mpp_get_current_ntile, mpp_get_tile_id, &
mpp_get_ntile_count, mpp_get_tile_list, &
mpp_get_global_domain, Domain1d, &
mpp_deallocate_domain, mpp_define_domains, &
mpp_get_domain_npes, mpp_get_domain_root_pe, &
mpp_domain_is_initialized, mpp_broadcast_domain, &
mpp_get_domain_pelist, mpp_compute_extent, &
domainUG, mpp_get_ug_compute_domains, &
mpp_get_ug_domains_index, mpp_get_ug_domain_grid_index, &
mpp_get_ug_domain_tile_list, mpp_pass_sg_to_ug
use constants_mod, only: PI, RADIUS
use mosaic2_mod, only: get_mosaic_xgrid, get_mosaic_xgrid_size, &
get_mosaic_ntiles, get_mosaic_ncontacts, &
get_mosaic_contact, get_mosaic_grid_sizes, &
get_mosaic_tile_grid
use stock_constants_mod, only: ISTOCK_TOP, ISTOCK_BOTTOM, ISTOCK_SIDE, STOCK_NAMES, &
STOCK_UNITS, NELEMS, stocks_file, stock_type
use gradient_mod, only: gradient_cubic
use fms2_io_mod, only: FmsNetcdfFile_t, open_file, variable_exists, close_file
use fms2_io_mod, only: FmsNetcdfDomainFile_t, read_data, get_dimension_size
use fms2_io_mod, only: get_variable_units, dimension_exists
use platform_mod, only: r8_kind, i8_kind, FMS_FILE_LEN
implicit none
private
public xmap_type, setup_xmap, set_frac_area, put_to_xgrid, get_from_xgrid, &
xgrid_count, some, conservation_check, xgrid_init, &
AREA_ATM_SPHERE, AREA_OCN_SPHERE, &
AREA_ATM_MODEL, AREA_OCN_MODEL, &
get_ocean_model_area_elements, grid_box_type, &
get_xmap_grid_area, put_to_xgrid_ug, get_from_xgrid_ug, &
set_frac_area_ug
!--- paramters that determine the remapping method
integer, parameter :: FIRST_ORDER = 1
integer, parameter :: SECOND_ORDER = 2
integer, parameter :: VERSION1 = 1 !< grid spec file
integer, parameter :: VERSION2 = 2 !< mosaic grid file
integer, parameter :: MAX_FIELDS = 80
logical :: make_exchange_reproduce = .false. !< Set to .true. to make <TT>xgrid_mod</TT> reproduce answers on different
!! numbers of PEs. This option has a considerable performance impact.
!< exactly same on different # PEs
character(len=64) :: interp_method = 'first_order' !< Exchange grid interpolation method.
!! It has two options: "first_order", "second_order".
logical :: debug_stocks = .false.
logical :: xgrid_clocks_on = .false.
logical :: monotonic_exchange = .false.
integer :: nsubset = 0 !< Number of processors to read exchange grid information. Those processors
!! that read the exchange grid information will send data to other processors
!! to prepare for flux exchange. Default value is 0. When nsubset is 0, each
!! processor will read part of the exchange grid information. The purpose of
!! this namelist is to improve performance of setup_xmap when running on
!! higher processor count and solve receiving size mismatch issue on high
!! processor count. Try to set nsubset = mpp_npes/MPI_rank_per_node.
logical :: do_alltoall = .true.
logical :: do_alltoallv = .false.
logical :: use_mpp_io = .false.!< use_mpp_io Default = .false. When true, uses mpp_io for IO.
!! When false, uses fms2_io for IO.
!> @brief xgrid nml
namelist /xgrid_nml/ make_exchange_reproduce, interp_method, debug_stocks, xgrid_clocks_on, &
monotonic_exchange, nsubset, do_alltoall, do_alltoallv, &
use_mpp_io
integer :: remapping_method
!> Area elements used inside each model
real(r8_kind), allocatable, dimension(:,:) :: AREA_ATM_MODEL, AREA_LND_MODEL, AREA_OCN_MODEL
!> Area elements based on a the spherical model used by the ICE layer
real(r8_kind), allocatable, dimension(:,:) :: AREA_ATM_SPHERE, AREA_LND_SPHERE, AREA_OCN_SPHERE
!> @}
!> @brief Scatters data from model grid onto exchange grid.
!!
!> Example usage:
!! @code{.F90}
!! call put_to_xgrid(d, grid_id, x, xmap, remap_order)
!! @endcode
!!
!> @ingroup xgrid_mod
interface put_to_xgrid
module procedure put_side1_to_xgrid
module procedure put_side2_to_xgrid
end interface
!> @brief Sums data from exchange grid to model grid.
!!
!> <br>Example usage:
!! @code{.F90}
!! call get_from_xgrid(d, grid_id, x, xmap)
!! @endcode
!> @ingroup xgrid_mod
interface get_from_xgrid
module procedure get_side1_from_xgrid
module procedure get_side2_from_xgrid
end interface
!> @brief @ref put_to_xgrid for unstructured grids.
!!
!> Scatters data from unstructured grid onto exchange grid.
!> @ingroup xgrid_mod
interface put_to_xgrid_ug
module procedure put_side1_to_xgrid_ug
module procedure put_side2_to_xgrid_ug
end interface
!> @brief @ref get_from_xgrid for unstructured grids.
!!
!> Sums data from exchange grid to model grid.
!> @ingroup xgrid_mod
interface get_from_xgrid_ug
module procedure get_side2_from_xgrid_ug
module procedure get_side1_from_xgrid_ug
end interface
!> @brief Sets sub-grid area and numbering in the given exchange grid.
!> @ingroup xgrid_mod
interface set_frac_area
module procedure set_frac_area_sg
module procedure set_frac_area_ug
end interface
!> @brief Returns three numbers which are the global sum of a variable.
!! @details Returns three numbers which are the global sum of a
!! variable (1) on its home model grid, (2) after interpolation to the other
!! side grid(s), and (3) after re_interpolation back onto its home side grid(s).
!! Conservation_check must be called by all PEs to work properly.
!!
!! @param d real(r8_kind) data
!! @param grid_id 3 character grid ID
!! @param[inout] xmap exchange grid
!! @param[out] global sum of a variable on home model grid, after side grid interpolation and after
!! reinterpolation
!!
!! <br>Example usage:
!! @code{.F90}
!! call conservation_check(d, grid_id, xmap,remap_order)
!! @endcode
!> @ingroup xgrid_mod
interface conservation_check
module procedure conservation_check_side1
module procedure conservation_check_side2
end interface
!> For an unstructured grid, returns three numbers which are the global sum of a
!! variable (1) on its home model grid, (2) after interpolation to the other
!! side grid(s), and (3) after re_interpolation back onto its home side grid(s).
!> @ingroup xgrid_mod
interface conservation_check_ug
module procedure conservation_check_ug_side1
module procedure conservation_check_ug_side2
end interface
!> Private type for cell indices and data in the exchange grid
!> @ingroup xgrid_mod
type xcell_type
integer :: i1 !< indices of cell in model arrays on both sides
integer :: j1 !< indices of cell in model arrays on both sides
integer :: i2 !< indices of cell in model arrays on both sides
integer :: j2 !< indices of cell in model arrays on both sides
integer :: l1, l2
integer :: recv_pos !< position in the receive buffer.
integer :: pe !< other side pe that has this cell
integer :: tile !< tile index of side 1 mosaic.
real(r8_kind) :: area !< geographic area of exchange cell
! real(r8_kind) :: area1_ratio !(= x_area/grid1_area), will be added in the future to improve efficiency
! real(r8_kind) :: area2_ratio !(= x_area/grid2_area), will be added in the future to improve efficiency
real(r8_kind) :: di !< Weight for the gradient of flux
real(r8_kind) :: dj !< Weight for the gradient of flux
real(r8_kind) :: scale
end type xcell_type
!> Type to hold pointers for grid boxes
!> @ingroup xgrid_mod
type grid_box_type
real(r8_kind), dimension(:,:), pointer :: dx => NULL()
real(r8_kind), dimension(:,:), pointer :: dy => NULL()
real(r8_kind), dimension(:,:), pointer :: area => NULL()
real(r8_kind), dimension(:), pointer :: edge_w => NULL()
real(r8_kind), dimension(:), pointer :: edge_e => NULL()
real(r8_kind), dimension(:), pointer :: edge_s => NULL()
real(r8_kind), dimension(:), pointer :: edge_n => NULL()
real(r8_kind), dimension(:,:,:), pointer :: en1 => NULL()
real(r8_kind), dimension(:,:,:), pointer :: en2 => NULL()
real(r8_kind), dimension(:,:,:), pointer :: vlon => NULL()
real(r8_kind), dimension(:,:,:), pointer :: vlat => NULL()
end type grid_box_type
!> Private type to hold all data needed from given grid for an exchange grid
!> @ingroup xgrid_mod
type grid_type
character(len=3) :: id !< grid identifier
integer :: npes !< number of processor on this grid.
logical :: on_this_pe !< indicate the domain is defined on this pe
integer :: root_pe !< indicate the root pe of the domain
integer, pointer, dimension(:) :: pelist !< pelist of the domain
integer :: ntile !< number of tiles in mosaic
integer :: ni !< max of global size of all the tiles
integer :: nj !< max of global size of all the tiles
integer, pointer, dimension(:) :: tile =>NULL() !< tile id ( pe index )
integer, pointer, dimension(:) :: is =>NULL() !< domain - i-range (pe index)
integer, pointer, dimension(:) :: ie =>NULL() !< domain - i-range (pe index)
integer, pointer, dimension(:) :: js =>NULL() !< domain - j-range (pe index)
integer, pointer, dimension(:) :: je =>NULL() !< domain - j-range (pe index)
integer, pointer :: is_me =>NULL() !< my domain - i-range
integer, pointer :: ie_me =>NULL() !< my domain - i-range
integer, pointer :: js_me =>NULL() !< my domain - j-range
integer, pointer :: je_me =>NULL() !< my domain - j-range
integer :: isd_me !< my data domain - i-range
integer :: ied_me !< my data domain - i-range
integer :: jsd_me !< my data domain - j-range
integer :: jed_me !< my data domain - j-range
integer :: nxd_me !< data domain size
integer :: nyd_me !< data domain size
integer :: nxc_me !< compute domain size
integer :: nyc_me !< compute domain size
integer, pointer :: tile_me !< my tile id
integer :: im !< global domain range
integer :: jm !< global domain range
integer :: km !< global domain range
real(r8_kind), pointer, dimension(:) :: lon =>NULL() !< center of global grids
real(r8_kind), pointer, dimension(:) :: lat =>NULL() !< center of global grids
real(r8_kind), pointer, dimension(:,:) :: geolon=>NULL() !< geographical grid center
real(r8_kind), pointer, dimension(:,:) :: geolat=>NULL() !< geographical grid center
real(r8_kind), pointer, dimension(:,:,:) :: frac_area =>NULL() !< partition fractions
real(r8_kind), pointer, dimension(:,:) :: area =>NULL() !< cell area
real(r8_kind), pointer, dimension(:,:) :: area_inv =>NULL() !< 1 / area for normalization
integer :: first !< xgrid index range
integer :: last !< xgrid index range
integer :: first_get !< xgrid index range for get_2_from_xgrid
integer :: last_get !< xgrid index range for get_2_from_xgrid
integer :: size !< # xcell patterns
type(xcell_type), pointer :: x(:) =>NULL() !< xcell patterns
integer :: size_repro !< # side 1 patterns for repro
type(xcell_type), pointer :: x_repro(:) =>NULL() !< side 1 patterns for repro
type(Domain2d) :: domain !< used for conservation checks
type(Domain2d) :: domain_with_halo !< used for second order remapping
logical :: is_latlon !< indicate if the grid is lat-lon grid or not.
type(grid_box_type) :: box !< used for second order remapping.
!--- The following is for land unstruct domain
logical :: is_ug
integer :: nxl_me
integer, pointer :: ls_me =>NULL() !< unstruct domain
integer, pointer :: le_me =>NULL() !< unstruct domain
integer, pointer, dimension(:) :: ls =>NULL(), le =>NULL()
integer, pointer :: gs_me =>NULL(), ge_me =>NULL()
integer, pointer, dimension(:) :: gs =>NULL(), ge =>NULL()
integer, pointer, dimension(:) :: l_index =>NULL()
type(DomainUG) :: ug_domain
end type grid_type
!> Private type for exchange grid data
!> @ingroup xgrid_mod
type x1_type
integer :: i, j
real(r8_kind) :: area !< (= geographic area * frac_area)
! real(r8_kind) :: area_ratio !(= x1_area/grid1_area) ! will be added in the future to improve efficiency
real(r8_kind) :: di !< weight for the gradient of flux
real(r8_kind) :: dj !< weight for the gradient of flux
integer :: tile !< tile index of side 1 mosaic.
integer :: pos
end type x1_type
!> Private type for exchange grid data
!> @ingroup xgrid_mod
type x2_type
integer :: i, j, l, k, pos
real(r8_kind) :: area !< geographic area of exchange cell
! real(r8_kind) :: area_ratio !(=x2_area/grid2_area ) ! will be added in the future to improve efficiency
end type x2_type
!> Private type for overlap exchange grid data
!> @ingroup xgrid_mod
type overlap_type
integer :: count
integer :: pe
integer :: buffer_pos
integer, allocatable :: i(:)
integer, allocatable :: j(:)
integer, allocatable :: g(:)
integer, allocatable :: xLoc(:)
integer, allocatable :: tile(:)
real(r8_kind), allocatable :: di(:)
real(r8_kind), allocatable :: dj(:)
end type overlap_type
!> Private type used for exchange grid communication
!> @ingroup xgrid_mod
type comm_type
integer :: nsend, nrecv
integer :: sendsize, recvsize
integer, pointer, dimension(:) :: unpack_ind=>NULL()
type(overlap_type), pointer, dimension(:) :: send=>NULL()
type(overlap_type), pointer, dimension(:) :: recv=>NULL()
end type comm_type
!> @brief Type for an exchange grid, holds pointers to included grids and any necessary data.
!> @ingroup xgrid_mod
type xmap_type
private
integer :: size !< # of exchange grid cells with area > 0 on this pe
integer :: size_put1 !< # of exchange grid cells for put_1_to_xgrid
integer :: size_get2 !< # of exchange grid cells for get_2_to_xgrid
integer :: me, npes, root_pe
logical, pointer, dimension(:) :: your1my2 =>NULL()!< true if side 1 domain on
!! indexed pe overlaps side 2
!! domain on this pe
logical, pointer, dimension(:) :: your2my1 =>NULL() !< true if a side 2 domain on
!! indexed pe overlaps side 1
!! domain on this pe
integer, pointer, dimension(:) :: your2my1_size=>NULL() !< number of exchange grid of
!! a side 2 domain on
!! indexed pe overlaps side 1
!! domain on this pe
type (grid_type), pointer, dimension(:) :: grids =>NULL() !< 1st grid is side 1;
!! rest on side 2
!
! Description of the individual exchange grid cells (index is cell #)
!
type(x1_type), pointer, dimension(:) :: x1 =>NULL() !< side 1 info
type(x1_type), pointer, dimension(:) :: x1_put =>NULL() !< side 1 info
type(x2_type), pointer, dimension(:) :: x2 =>NULL() !< side 2 info
type(x2_type), pointer, dimension(:) :: x2_get =>NULL() !< side 2 info
integer, pointer, dimension(:) :: send_count_repro =>NULL()
integer, pointer, dimension(:) :: recv_count_repro =>NULL()
integer :: send_count_repro_tot !< sum(send_count_repro)
integer :: recv_count_repro_tot !< sum(recv_count_repro)
integer :: version !< version of xgrids. version=VERSION! is for grid_spec file
!! and version=VERSION2 is for mosaic grid.
integer, pointer, dimension(:) :: ind_get1 =>NULL() !< indx for side1 get and side2 put.
integer, pointer, dimension(:) :: ind_put1 =>NULL() !< indx for side1 put and side 2get.
type(comm_type), pointer :: put1 =>NULL() !< for put_1_to_xgrid
type(comm_type), pointer :: get1 =>NULL() !< for get_1_from_xgrid
type(comm_type), pointer :: get1_repro =>NULL()!< for get_1_from_xgrid_repro
end type xmap_type
!> @addtogroup stock_constants_mod
!> @{
!-----------------------------------------------------------------------
! Include variable "version" to be written to log file.
#include<file_version.h>
real(r8_kind), parameter :: EPS = 1.0e-10_r8_kind
real(r8_kind), parameter :: LARGE_NUMBER = 1.e20_r8_kind
logical :: module_is_initialized = .FALSE.
integer :: id_put_1_to_xgrid_order_1 = 0
integer :: id_put_1_to_xgrid_order_2 = 0
integer :: id_get_1_from_xgrid = 0
integer :: id_get_1_from_xgrid_repro = 0
integer :: id_get_2_from_xgrid = 0
integer :: id_put_2_to_xgrid = 0
integer :: id_setup_xmap = 0
integer :: id_load_xgrid1, id_load_xgrid2, id_load_xgrid3
integer :: id_load_xgrid4, id_load_xgrid5
integer :: id_load_xgrid, id_set_comm, id_regen, id_conservation_check
! The following is for nested model
integer :: nnest=0, tile_nest, tile_parent
integer :: is_nest=0, ie_nest=0, js_nest=0, je_nest=0
integer :: is_parent=0, ie_parent=0, js_parent=0, je_parent=0
!> @}
! The following is required to compute stocks of water, heat, ...
!> @ingroup xgrid_mod
interface stock_move
module procedure stock_move_3d, stock_move_2d
end interface
!> @ingroup xgrid_mod
interface stock_move_ug
module procedure stock_move_ug_3d
end interface
public stock_move, stock_type, stock_print, get_index_range, stock_integrate_2d
public FIRST_ORDER, SECOND_ORDER, stock_move_ug
!> @ingroup xgrid_mod
interface get_area_elements
module procedure get_area_elements_fms2_io
end interface
!> @ingroup xgrid_mod
interface get_nest_contact
module procedure get_nest_contact_fms2_io
end interface
contains
!> @addtogroup xgrid_mod
!> @{
!#######################################################################
!> @return logical in_box
logical function in_box(i, j, is, ie, js, je)
integer, intent(in) :: i, j, is, ie, js, je
in_box = (i>=is) .and. (i<=ie) .and. (j>=js) .and. (j<=je)
end function in_box
!#######################################################################
!> @brief Initialize the xgrid_mod.
!! @details Initialization routine for the xgrid module. It reads the xgrid_nml,
!! writes the version information and xgrid_nml to the log file.
subroutine xgrid_init(remap_method)
integer, intent(out) :: remap_method !< exchange grid interpolation method. It has four possible values:
!! FIRST_ORDER (=1), SECOND_ORDER(=2).
integer :: iunit, ierr, io, out_unit
if (module_is_initialized) return
module_is_initialized = .TRUE.
read (input_nml_file, xgrid_nml, iostat=io)
ierr = check_nml_error ( io, 'xgrid_nml' )
!--------- write version number and namelist ------------------
call write_version_number("XGRID_MOD", version)
iunit = stdlog ( )
out_unit = stdout()
if ( mpp_pe() == mpp_root_pe() ) write (iunit,nml=xgrid_nml)
if (use_mpp_io) then
! FATAL error if trying to use mpp_io
call error_mesg('xgrid_init', &
'MPP_IO is no longer supported. Please remove use_mpp_io from namelists',&
FATAL)
endif
!--------- check interp_method has suitable value
!--- when monotonic_exchange is true, interp_method must be second order.
select case(trim(interp_method))
case('first_order')
remap_method = FIRST_ORDER
if( monotonic_exchange ) call error_mesg('xgrid_mod', &
'xgrid_nml monotonic_exchange must be .false. when interp_method = first_order', FATAL)
write(out_unit,*)"NOTE from xgrid_mod: use first_order conservative exchange"
case('second_order')
if(monotonic_exchange) then
write(out_unit,*)"NOTE from xgrid_mod: use monotonic second_order conservative exchange"
else
write(out_unit,*)"NOTE from xgrid_mod: use second_order conservative exchange"
endif
remap_method = SECOND_ORDER
case default
call error_mesg('xgrid_mod', ' nml interp_method = ' //trim(interp_method)// &
' is not a valid namelist option', FATAL)
end select
if(xgrid_clocks_on) then
id_put_1_to_xgrid_order_1 = mpp_clock_id("put_1_to_xgrid_order_1", flags=MPP_CLOCK_SYNC)
id_put_1_to_xgrid_order_2 = mpp_clock_id("put_1_to_xgrid_order_2", flags=MPP_CLOCK_SYNC)
id_get_1_from_xgrid = mpp_clock_id("get_1_from_xgrid", flags=MPP_CLOCK_SYNC)
id_get_1_from_xgrid_repro = mpp_clock_id("get_1_from_xgrid_repro", flags=MPP_CLOCK_SYNC)
id_get_2_from_xgrid = mpp_clock_id("get_2_from_xgrid", flags=MPP_CLOCK_SYNC)
id_put_2_to_xgrid = mpp_clock_id("put_2_to_xgrid", flags=MPP_CLOCK_SYNC)
id_setup_xmap = mpp_clock_id("setup_xmap", flags=MPP_CLOCK_SYNC)
id_set_comm = mpp_clock_id("set_comm")
id_regen = mpp_clock_id("regen")
id_conservation_check = mpp_clock_id("conservation_check")
id_load_xgrid = mpp_clock_id("load_xgrid")
id_load_xgrid1 = mpp_clock_id("load_xgrid1")
id_load_xgrid2 = mpp_clock_id("load_xgrid2")
id_load_xgrid3 = mpp_clock_id("load_xgrid3")
id_load_xgrid4 = mpp_clock_id("load_xgrid4")
id_load_xgrid5 = mpp_clock_id("load_xgrid5")
endif
remapping_method = remap_method
end subroutine xgrid_init
!#######################################################################
subroutine load_xgrid (xmap, grid, grid_file, grid1_id, grid_id, tile1, tile2, use_higher_order)
type(xmap_type), intent(inout) :: xmap
type(grid_type), intent(inout) :: grid
character(len=*), intent(in) :: grid_file
character(len=3), intent(in) :: grid1_id, grid_id
integer, intent(in) :: tile1, tile2
logical, intent(in) :: use_higher_order
integer, pointer, dimension(:) :: i1=>NULL(), j1=>NULL()
integer, pointer, dimension(:) :: i2=>NULL(), j2=>NULL()
real(r8_kind), pointer, dimension(:) :: di=>NULL(), dj=>NULL()
real(r8_kind), pointer, dimension(:) :: area =>NULL()
integer, pointer, dimension(:) :: i1_tmp=>NULL(), j1_tmp=>NULL()
integer, pointer, dimension(:) :: i2_tmp=>NULL(), j2_tmp=>NULL()
real(r8_kind), pointer, dimension(:) :: di_tmp=>NULL(), dj_tmp=>NULL()
real(r8_kind), pointer, dimension(:) :: area_tmp =>NULL()
integer, pointer, dimension(:) :: i1_side1=>NULL(), j1_side1=>NULL()
integer, pointer, dimension(:) :: i2_side1=>NULL(), j2_side1=>NULL()
real(r8_kind), pointer, dimension(:) :: di_side1=>NULL(), dj_side1=>NULL()
real(r8_kind), pointer, dimension(:) :: area_side1 =>NULL()
real(r8_kind), allocatable, dimension(:,:) :: tmp
real(r8_kind), allocatable, dimension(:) :: send_buffer, recv_buffer
type (grid_type), pointer, save :: grid1 =>NULL()
integer :: l, ll, ll_repro, p, nxgrid, size_prev
type(xcell_type), allocatable :: x_local(:)
integer :: size_repro, out_unit
logical :: scale_exist = .false.
logical :: is_distribute = .false.
real(r8_kind), allocatable, dimension(:) :: scale
real(r8_kind) :: garea
integer :: npes, isc, iec, nxgrid_local, pe, nxgrid_local_orig
integer :: nxgrid1, nxgrid2, nset1, nset2, ndivs, cur_ind
integer :: pos, nsend, nrecv, l1, l2, n, mypos
integer :: start(4), nread(4)
logical :: found
character(len=128) :: attvalue
integer, dimension(0:xmap%npes-1) :: pelist
logical, dimension(0:xmap%npes-1) :: subset_rootpe
integer, dimension(0:xmap%npes-1) :: nsend1, nsend2, nrecv1, nrecv2
integer, dimension(0:xmap%npes-1) :: send_cnt, recv_cnt
integer, dimension(0:xmap%npes-1) :: send_buffer_pos, recv_buffer_pos
integer, dimension(0:xmap%npes-1) :: ibegin, iend, pebegin, peend
integer, dimension(2*xmap%npes) :: ibuf1, ibuf2
integer, dimension(0:xmap%npes-1) :: pos_x, y2m1_size
integer, allocatable, dimension(:) :: y2m1_pe
integer, pointer, save :: iarray(:), jarray(:)
integer, allocatable, save :: pos_s(:)
integer, pointer, dimension(:) :: iarray2(:)=>NULL(), jarray2(:)=>NULL()
logical :: last_grid
integer :: nxgrid1_old
integer :: lll
type(FmsNetcdfFile_t) :: fileobj
if(.not. open_file(fileobj, grid_file, 'read' )) then
call error_mesg('xgrid_mod(load_xgrid)', 'Error in opening file '//trim(grid_file), FATAL)
endif
scale_exist = .false.
grid1 => xmap%grids(1)
out_unit = stdout()
npes = xmap%npes
pe = mpp_pe()
mypos = mpp_pe()-mpp_root_pe()
call mpp_get_current_pelist(pelist)
!--- make sure npes = pelist(npes-1) - pelist(0) + 1
if( npes .NE. pelist(npes-1) - pelist(0) + 1 ) then
print*, "npes =", npes, ", pelist(npes-1)=", pelist(npes-1), ", pelist(0)=", pelist(0)
call error_mesg('xgrid_mod', 'npes .NE. pelist(npes-1) - pelist(0)', FATAL)
endif
select case(xmap%version)
case(VERSION1)
nxgrid = 0
if (dimension_exists(fileobj, 'i_'//lowercase(grid1_id)//'X'//lowercase(grid_id))) then
call get_dimension_size(fileobj, 'i_'//lowercase(grid1_id)//'X'//lowercase(grid_id), nxgrid)
endif
if(nxgrid .LE. 0) return
case(VERSION2)
!--- max_size is the exchange grid size between super grid.
nxgrid = get_mosaic_xgrid_size(fileobj)
if(nxgrid .LE. 0) return
end select
!--- define a domain to read exchange grid.
if(nxgrid > npes) then
ndivs = npes
if(nsubset >0 .AND. nsubset < npes) ndivs = nsubset
call mpp_compute_extent( 1, nxgrid, ndivs, ibegin, iend)
if(npes == ndivs) then
p = mpp_pe()-mpp_root_pe()
isc = ibegin(p)
iec = iend(p)
subset_rootpe(:) = .true.
else
isc = 0; iec = -1
call mpp_compute_extent(pelist(0), pelist(npes-1), ndivs, pebegin, peend)
do n = 0, ndivs-1
if(pe == pebegin(n)) then
isc = ibegin(n)
iec = iend(n)
exit
endif
enddo
cur_ind = 0
subset_rootpe(:) = .false.
do n = 0, npes-1
if(pelist(n) == pebegin(cur_ind)) then
subset_rootpe(n) = .true.
cur_ind = cur_ind+1
if(cur_ind == ndivs) exit
endif
enddo
endif
is_distribute = .true.
else
is_distribute = .false.
isc = 1; iec = nxgrid
endif
nset1 = 5
nset2 = 5
if(use_higher_order) then
nset1 = nset1 + 2
nset2 = nset2 + 2
end if
if(scale_exist) nset2 = nset1 + 1
call mpp_clock_begin(id_load_xgrid1)
if(iec .GE. isc) then
nxgrid_local = iec - isc + 1
allocate(i1_tmp(isc:iec), j1_tmp(isc:iec), i2_tmp(isc:iec), j2_tmp(isc:iec), area_tmp(isc:iec) )
if(use_higher_order) allocate(di_tmp(isc:iec), dj_tmp(isc:iec))
start = 1; nread = 1
select case(xmap%version)
case(VERSION1)
start(1) = isc; nread(1) = nxgrid_local
allocate(tmp(nxgrid_local,1))
call read_data(fileobj, 'I_'//grid1_id//'_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
i1_tmp = int(tmp(:,1))
call read_data(fileobj, 'J_'//grid1_id//'_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
j1_tmp = int(tmp(:,1))
call read_data(fileobj, 'I_'//grid_id//'_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
i2_tmp = int(tmp(:,1))
call read_data(fileobj, 'J_'//grid_id//'_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
j2_tmp = int(tmp(:,1))
call read_data(fileobj, 'AREA_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
area_tmp = tmp(:,1)
if(use_higher_order) then
call read_data(fileobj, 'DI_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
di_tmp = tmp(:,1)
call read_data(fileobj, 'DJ_'//grid1_id//'x'//grid_id, tmp, corner=start, edge_lengths=nread)
dj_tmp = tmp(:,1)
end if
deallocate(tmp)
case(VERSION2)
nread(1) = 2; start(2) = isc; nread(2) = nxgrid_local
allocate(tmp(2, isc:iec))
call read_data(fileobj, "tile1_cell", tmp, corner=start, edge_lengths=nread)
i1_tmp(isc:iec) = int(tmp(1, isc:iec))
j1_tmp(isc:iec) = int(tmp(2, isc:iec))
call read_data(fileobj, "tile2_cell", tmp, corner=start, edge_lengths=nread)
i2_tmp(isc:iec) = int(tmp(1, isc:iec))
j2_tmp(isc:iec) = int(tmp(2, isc:iec))
if(use_higher_order) then
call read_data(fileobj, "tile1_distance", tmp, corner=start, edge_lengths=nread)
di_tmp(isc:iec) = tmp(1, isc:iec)
dj_tmp(isc:iec) = tmp(2, isc:iec)
end if
start = 1; nread = 1
start(1) = isc; nread(1) = nxgrid_local
deallocate(tmp)
allocate(tmp(isc:iec,1) )
call read_data(fileobj, "xgrid_area", tmp(:,1:1), corner=start, edge_lengths=nread)
! check the units of "xgrid_area
call get_variable_units(fileobj, "xgrid_area", attvalue)
if( trim(attvalue) == 'm2' ) then
garea = 4.0_r8_kind * PI * RADIUS * RADIUS;
area_tmp = tmp(:,1)/garea
else if( trim(attvalue) == 'none' ) then
area_tmp = tmp(:,1)
else
call error_mesg('xgrid_mod', 'In file '//trim(grid_file)//', xgrid_area units = '// &
trim(attvalue)//' should be "m2" or "none"', FATAL)
endif
!--- if field "scale" exist, read this field. Normally this
!--- field only exist in landXocean exchange grid cell.
if(grid1_id == 'LND' .AND. grid_id == 'OCN') then
if(variable_exists(fileobj, "scale")) then
allocate(scale(isc:iec))
write(out_unit, *)"NOTE from load_xgrid(xgrid_mod): field 'scale' exist in the file "// &
& trim(grid_file)//", this field will be read and the exchange grid cell area will be"// &
& " multiplied by scale"
call read_data(fileobj, "scale", tmp, corner=start, edge_lengths=nread)
scale = tmp(:,1)
scale_exist = .true.
endif
endif
deallocate(tmp)
end select
!---z1l: The following change is for the situation that some processor is masked out.
!---loop through all the pe to see if side 1 and side of each exchange grid is on some processor
nxgrid_local_orig = nxgrid_local
allocate(i1(isc:iec), j1(isc:iec), i2(isc:iec), j2(isc:iec), area(isc:iec) )
if(use_higher_order) allocate(di(isc:iec), dj(isc:iec))
pos = isc-1
do l = isc, iec
found = .false.
!--- first check if the exchange grid is on one of side 1 processor
do p = 0, npes - 1
if(grid1%tile(p) == tile1) then
if(in_box_nbr(i1_tmp(l), j1_tmp(l), grid1, p)) then
found = .true.
exit
endif
endif
enddo
!--- Then check if the exchange grid is on one of side 2 processor
if( found ) then
do p = 0, npes - 1
if(grid%tile(p) == tile2) then
if (in_box_nbr(i2_tmp(l), j2_tmp(l), grid, p)) then
pos = pos+1
i1(pos) = i1_tmp(l)
j1(pos) = j1_tmp(l)
i2(pos) = i2_tmp(l)
j2(pos) = j2_tmp(l)
area(pos) = area_tmp(l)
if(use_higher_order) then
di(pos) = di_tmp(l)
dj(pos) = dj_tmp(l)
endif
exit
endif
endif
enddo
endif
enddo
deallocate(i1_tmp, i2_tmp, j1_tmp, j2_tmp, area_tmp)
if(use_higher_order) deallocate( di_tmp, dj_tmp)
iec = pos
if(iec .GE. isc) then
nxgrid_local = iec - isc + 1
else
nxgrid_local = 0
endif
else
nxgrid_local = 0
nxgrid_local_orig = 0
endif
call close_file(fileobj)
call mpp_clock_end(id_load_xgrid1)
if(is_distribute) then
!--- Since the xgrid is distributed according to side 2 grid. Send all the xgrid to its own side 2.
!--- Also need to send the xgrid to its own side 1 for the reproducing ability between processor count.
!--- first find out number of points need to send to other pe and fill the send buffer.
nsend1(:) = 0; nrecv1(:) = 0
nsend2(:) = 0; nrecv2(:) = 0
ibuf1(:)= 0; ibuf2(:)= 0
call mpp_clock_begin(id_load_xgrid2)
if(nxgrid_local>0) then
allocate( send_buffer(nxgrid_local * (nset1+nset2)) )
pos = 0
do p = 0, npes - 1
send_buffer_pos(p) = pos
if(grid%tile(p) == tile2) then
do l = isc, iec
if(in_box_nbr(i2(l), j2(l), grid, p) ) then
nsend2(p) = nsend2(p) + 1
send_buffer(pos+1) = real(i1(l), r8_kind)
send_buffer(pos+2) = real(j1(l), r8_kind)
send_buffer(pos+3) = real(i2(l), r8_kind)
send_buffer(pos+4) = real(j2(l), r8_kind)
send_buffer(pos+5) = area(l)
if(use_higher_order) then
send_buffer(pos+6) = di(l)
send_buffer(pos+7) = dj(l)
endif
if(scale_exist) send_buffer(pos+nset2) = scale(l)
pos = pos + nset2
endif
enddo
endif
if(grid1%tile(p) == tile1) then
do l = isc, iec
if(in_box_nbr(i1(l), j1(l), grid1, p)) then
nsend1(p) = nsend1(p) + 1
send_buffer(pos+1) = real(i1(l), r8_kind)
send_buffer(pos+2) = real(j1(l), r8_kind)
send_buffer(pos+3) = real(i2(l), r8_kind)
send_buffer(pos+4) = real(j2(l), r8_kind)
send_buffer(pos+5) = area(l)
if(use_higher_order) then
send_buffer(pos+6) = di(l)
send_buffer(pos+7) = dj(l)
endif
pos = pos + nset1
endif
enddo
endif
enddo
endif
call mpp_clock_end(id_load_xgrid2)
!--- send the size of the data on side 1 to be sent over.
call mpp_clock_begin(id_load_xgrid3)
if (do_alltoall) then
do p = 0, npes-1
ibuf1(2*p+1) = nsend1(p)
ibuf1(2*p+2) = nsend2(p)
enddo
call mpp_alltoall(ibuf1, 2, ibuf2, 2)
else
do n = 0, npes-1
p = mod(mypos+npes-n, npes)
if(.not. subset_rootpe(p)) cycle
call mpp_recv( ibuf2(2*p+1), glen=2, from_pe=pelist(p), block=.FALSE., tag=COMM_TAG_1)
enddo
if(nxgrid_local_orig>0) then
do n = 0, npes-1
p = mod(mypos+n, npes)
ibuf1(2*p+1) = nsend1(p)
ibuf1(2*p+2) = nsend2(p)
call mpp_send( ibuf1(2*p+1), plen=2, to_pe=pelist(p), tag=COMM_TAG_1)
enddo
endif
call mpp_sync_self(check=EVENT_RECV)
endif
do p = 0, npes-1
nrecv1(p) = ibuf2(2*p+1)
nrecv2(p) = ibuf2(2*p+2)
enddo
if(.not. do_alltoall) call mpp_sync_self()
call mpp_clock_end(id_load_xgrid3)
call mpp_clock_begin(id_load_xgrid4)
pos = 0
do p = 0, npes - 1
recv_buffer_pos(p) = pos
pos = pos + nrecv1(p) * nset1 + nrecv2(p) * nset2
end do
!--- now get the data
nxgrid1 = sum(nrecv1)
nxgrid2 = sum(nrecv2)
if(nxgrid1>0 .OR. nxgrid2>0) allocate(recv_buffer(nxgrid1*nset1+nxgrid2*nset2))
if (do_alltoallv) then
! Construct the send and receive counters
send_cnt(:) = nset1 * nsend1(:) + nset2 * nsend2(:)
recv_cnt(:) = nset1 * nrecv1(:) + nset2 * nrecv2(:)
call mpp_alltoall(send_buffer, send_cnt, send_buffer_pos, &
recv_buffer, recv_cnt, recv_buffer_pos)
else
do n = 0, npes-1
p = mod(mypos+npes-n, npes)
nrecv = nrecv1(p)*nset1+nrecv2(p)*nset2
if(nrecv==0) cycle
pos = recv_buffer_pos(p)
call mpp_recv(recv_buffer(pos+1), glen=nrecv, from_pe=pelist(p), &
block=.FALSE., tag=COMM_TAG_2)
end do
do n = 0, npes-1
p = mod(mypos+n, npes)
nsend = nsend1(p)*nset1 + nsend2(p)*nset2
if(nsend==0) cycle
pos = send_buffer_pos(p)
call mpp_send(send_buffer(pos+1), plen=nsend, to_pe=pelist(p), &
tag=COMM_TAG_2)
end do
call mpp_sync_self(check=EVENT_RECV)
end if
call mpp_clock_end(id_load_xgrid4)
!--- unpack buffer.
if( nxgrid_local>0) then
deallocate(i1,j1,i2,j2,area)
endif
allocate(i1(nxgrid2), j1(nxgrid2))
allocate(i2(nxgrid2), j2(nxgrid2))
allocate(area(nxgrid2))
allocate(i1_side1(nxgrid1), j1_side1(nxgrid1))
allocate(i2_side1(nxgrid1), j2_side1(nxgrid1))
allocate(area_side1(nxgrid1))
if(use_higher_order) then
if(nxgrid_local>0) deallocate(di,dj)
allocate(di (nxgrid2), dj (nxgrid2))
allocate(di_side1(nxgrid1), dj_side1(nxgrid1))
endif
if(scale_exist) then
if(nxgrid_local>0)deallocate(scale)