-
Notifications
You must be signed in to change notification settings - Fork 136
/
mpp_update_domains2D.fh
642 lines (582 loc) · 30.1 KB
/
mpp_update_domains2D.fh
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
! -*-f90-*-
!***********************************************************************
!* 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/>.
!***********************************************************************
!> @addtogroup mpp_domains_mod
!> @{
!> Updates data domain of 2D field whose computational domains have been computed
subroutine MPP_UPDATE_DOMAINS_2D_( field, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count)
MPP_TYPE_, intent(inout) :: field(:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags
logical, intent(in), optional :: complete
integer, intent(in), optional :: position
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
MPP_TYPE_ :: field3D(size(field,1),size(field,2),1)
pointer( ptr, field3D )
ptr = LOC(field)
call mpp_update_domains( field3D, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count )
return
end subroutine MPP_UPDATE_DOMAINS_2D_
!> Updates data domain of 3D field whose computational domains have been computed
subroutine MPP_UPDATE_DOMAINS_3D_( field, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count)
MPP_TYPE_, intent(inout) :: field(:,:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags
logical, intent(in), optional :: complete
integer, intent(in), optional :: position
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo ! specify halo region to be updated.
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
integer :: update_position, update_whalo, update_ehalo, update_shalo, update_nhalo, ntile
integer(i8_kind),dimension(MAX_DOMAIN_FIELDS, MAX_TILES),save :: f_addrs=-9999
integer :: tile, max_ntile
character(len=3) :: text
logical :: set_mismatch, is_complete
logical :: do_update
integer, save :: isize=0, jsize=0, ke=0, l_size=0, list=0
integer, save :: pos, whalosz, ehalosz, shalosz, nhalosz
MPP_TYPE_ :: d_type
type(overlapSpec), pointer :: update => NULL()
type(overlapSpec), pointer :: check => NULL()
if(present(whalo)) then
update_whalo = whalo
if(abs(update_whalo) > domain%whalo ) call mpp_error(FATAL, "MPP_UPDATE_3D: "// &
"optional argument whalo should not be larger than the whalo when define domain.")
else
update_whalo = domain%whalo
end if
if(present(ehalo)) then
update_ehalo = ehalo
if(abs(update_ehalo) > domain%ehalo ) call mpp_error(FATAL, "MPP_UPDATE_3D: "// &
"optional argument ehalo should not be larger than the ehalo when define domain.")
else
update_ehalo = domain%ehalo
end if
if(present(shalo)) then
update_shalo = shalo
if(abs(update_shalo) > domain%shalo ) call mpp_error(FATAL, "MPP_UPDATE_3D: "// &
"optional argument shalo should not be larger than the shalo when define domain.")
else
update_shalo = domain%shalo
end if
if(present(nhalo)) then
update_nhalo = nhalo
if(abs(update_nhalo) > domain%nhalo ) call mpp_error(FATAL, "MPP_UPDATE_3D: "// &
"optional argument nhalo should not be larger than the nhalo when define domain.")
else
update_nhalo = domain%nhalo
end if
!--- when there is NINETY or MINUS_NINETY rotation for some contact, the salar data can not be on E or N-cell,
if(present(position)) then
if(domain%rotated_ninety .AND. ( position == EAST .OR. position == NORTH ) ) &
call mpp_error(FATAL, 'MPP_UPDATE_3D: hen there is NINETY or MINUS_NINETY rotation, ' // &
'can not use scalar version update_domain for data on E or N-cell' )
end if
max_ntile = domain%max_ntile_pe
ntile = size(domain%x(:))
is_complete = .true.
if(PRESENT(complete)) then
is_complete = complete
end if
tile = 1
if(max_ntile>1) then
if(ntile>MAX_TILES) then
write( text,'(i2)' ) MAX_TILES
call mpp_error(FATAL,'MPP_UPDATE_3D: MAX_TILES='//text//' is less than number of tiles on this pe.' )
endif
if(.NOT. present(tile_count) ) call mpp_error(FATAL, "MPP_UPDATE_3D: "// &
"optional argument tile_count should be present when number of tiles on this pe is more than 1")
tile = tile_count
end if
do_update = (tile == ntile) .AND. is_complete
list = list+1
if(list > MAX_DOMAIN_FIELDS)then
write( text,'(i2)' ) MAX_DOMAIN_FIELDS
call mpp_error(FATAL,'MPP_UPDATE_3D: MAX_DOMAIN_FIELDS='//text//' exceeded for group update.' )
endif
f_addrs(list, tile) = LOC(field)
update_position = CENTER
if(present(position)) update_position = position
if(list == 1 .AND. tile == 1 )then
isize=size(field,1); jsize=size(field,2); ke = size(field,3); pos = update_position
whalosz = update_whalo; ehalosz = update_ehalo; shalosz = update_shalo; nhalosz = update_nhalo
else
set_mismatch = .false.
set_mismatch = set_mismatch .OR. (isize /= size(field,1))
set_mismatch = set_mismatch .OR. (jsize /= size(field,2))
set_mismatch = set_mismatch .OR. (ke /= size(field,3))
set_mismatch = set_mismatch .OR. (update_position /= pos)
set_mismatch = set_mismatch .OR. (update_whalo /= whalosz)
set_mismatch = set_mismatch .OR. (update_ehalo /= ehalosz)
set_mismatch = set_mismatch .OR. (update_shalo /= shalosz)
set_mismatch = set_mismatch .OR. (update_nhalo /= nhalosz)
if(set_mismatch)then
write( text,'(i2)' ) list
call mpp_error(FATAL,'MPP_UPDATE_3D: Incompatible field at count '//text//' for group update.' )
endif
endif
if(is_complete) then
l_size = list
list = 0
end if
if(do_update )then
if( domain_update_is_needed(domain, update_whalo, update_ehalo, update_shalo, update_nhalo) )then
if(debug_update_level .NE. NO_CHECK) then
check => search_check_overlap(domain, update_position)
if(ASSOCIATED(check) ) then
call mpp_do_check(f_addrs(1:l_size,1:ntile), domain, check, d_type, ke, flags, name )
endif
endif
update => search_update_overlap(domain, update_whalo, update_ehalo, update_shalo, update_nhalo, &
& update_position)
!call mpp_do_update( f_addrs(1:l_size,1:ntile), domain, update, d_type, ke, &
! b_addrs(1:l_size,1:ntile), bsize, flags)
if ( PRESENT ( flags ) ) then
call mpp_do_update( f_addrs(1:l_size,1:ntile), domain, update, d_type, ke, flags )
else
call mpp_do_update( f_addrs(1:l_size,1:ntile), domain, update, d_type, ke )
endif
end if
l_size=0; f_addrs=-9999; isize=0; jsize=0; ke=0
endif
return
end subroutine MPP_UPDATE_DOMAINS_3D_
!> Updates data domain of 4D field whose computational domains have been computed
subroutine MPP_UPDATE_DOMAINS_4D_( field, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count )
MPP_TYPE_, intent(inout) :: field(:,:,:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags
logical, intent(in), optional :: complete
integer, intent(in), optional :: position
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
MPP_TYPE_ :: field3D(size(field,1),size(field,2),size(field,3)*size(field,4))
pointer( ptr, field3D )
ptr = LOC(field)
call mpp_update_domains( field3D, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count)
return
end subroutine MPP_UPDATE_DOMAINS_4D_
!> Updates data domain of 5D field whose computational domains have been computed
subroutine MPP_UPDATE_DOMAINS_5D_( field, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count )
MPP_TYPE_, intent(inout) :: field(:,:,:,:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags
logical, intent(in), optional :: complete
integer, intent(in), optional :: position
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
MPP_TYPE_ :: field3D(size(field,1),size(field,2),size(field,3)*size(field,4)*size(field,5))
pointer( ptr, field3D )
ptr = LOC(field)
call mpp_update_domains( field3D, domain, flags, complete, position, &
whalo, ehalo, shalo, nhalo, name, tile_count )
return
end subroutine MPP_UPDATE_DOMAINS_5D_
subroutine MPP_REDISTRIBUTE_2D_( domain_in, field_in, domain_out, field_out, complete, free, list_size, &
& dc_handle, position )
type(domain2D), intent(in) :: domain_in, domain_out
MPP_TYPE_, intent(in) :: field_in (:,:)
MPP_TYPE_, intent(out) :: field_out(:,:)
logical, intent(in), optional :: complete, free
integer, intent(in), optional :: list_size
integer, intent(in), optional :: position
MPP_TYPE_ :: field3D_in (size(field_in, 1),size(field_in, 2),1)
MPP_TYPE_ :: field3D_out(size(field_out,1),size(field_out,2),1)
type(DomainCommunicator2D),pointer,optional :: dc_handle
pointer( ptr_in, field3D_in )
pointer( ptr_out, field3D_out )
field_out = 0
ptr_in = 0
ptr_out = 0
if(domain_in%initialized) ptr_in = LOC(field_in )
if(domain_out%initialized) ptr_out = LOC(field_out)
call mpp_redistribute( domain_in, field3D_in, domain_out, field3D_out, complete, free, list_size, &
& dc_handle, position )
return
end subroutine MPP_REDISTRIBUTE_2D_
subroutine MPP_REDISTRIBUTE_3D_( domain_in, field_in, domain_out, field_out, complete, free, list_size, &
& dc_handle, position )
type(domain2D), intent(in) :: domain_in, domain_out
MPP_TYPE_, intent(in) :: field_in (:,:,:)
MPP_TYPE_, intent(out) :: field_out(:,:,:)
logical, intent(in), optional :: complete, free
integer, intent(in), optional :: list_size
integer, intent(in), optional :: position
type(DomainCommunicator2D),pointer,optional :: dc_handle
type(DomainCommunicator2D),pointer,save :: d_comm =>NULL()
logical :: do_redist,free_comm
integer :: lsize
integer(i8_kind),dimension(MAX_DOMAIN_FIELDS),save :: l_addrs_in=-9999, l_addrs_out=-9999
integer, save :: isize_in=0,jsize_in=0,ke_in=0,l_size=0
integer, save :: isize_out=0,jsize_out=0,ke_out=0
logical :: set_mismatch
integer :: ke
character(len=2) :: text
MPP_TYPE_ :: d_type
integer(i8_kind) :: floc_in, floc_out
floc_in = 0
floc_out = 0
if(domain_in%initialized) floc_in = LOC(field_in)
if(domain_out%initialized) floc_out = LOC(field_out)
if(present(position)) then
if(position .NE. CENTER) call mpp_error( FATAL, &
'MPP_REDISTRIBUTE_3Dold_: only position = CENTER is implemented, contact author')
endif
do_redist=.true.; if(PRESENT(complete))do_redist=complete
free_comm=.false.; if(PRESENT(free))free_comm=free
if(free_comm)then
l_addrs_in(1) = floc_in; l_addrs_out(1) = floc_out
if(l_addrs_out(1)>0)then
ke = size(field_out,3)
else
ke = size(field_in,3)
end if
lsize=1; if(PRESENT(list_size))lsize=list_size
call mpp_redistribute_free_comm(domain_in,l_addrs_in(1),domain_out,l_addrs_out(1),ke,lsize)
else
l_size = l_size+1
if(l_size > MAX_DOMAIN_FIELDS)then
write( text,'(i2)' ) MAX_DOMAIN_FIELDS
call mpp_error(FATAL,'MPP_REDISTRIBUTE_3D: MAX_DOMAIN_FIELDS='//text//' exceeded for group redistribute.' )
end if
l_addrs_in(l_size) = floc_in; l_addrs_out(l_size) = floc_out
if(l_size == 1)then
if(l_addrs_in(l_size) > 0)then
isize_in=size(field_in,1); jsize_in=size(field_in,2); ke_in = size(field_in,3)
end if
if(l_addrs_out(l_size) > 0)then
isize_out=size(field_out,1); jsize_out=size(field_out,2); ke_out = size(field_out,3)
endif
else
set_mismatch = .false.
set_mismatch = l_addrs_in(l_size) == 0 .AND. l_addrs_in(l_size-1) /= 0
set_mismatch = set_mismatch .OR. (l_addrs_in(l_size) > 0 .AND. l_addrs_in(l_size-1) == 0)
set_mismatch = set_mismatch .OR. (l_addrs_out(l_size) == 0 .AND. l_addrs_out(l_size-1) /= 0)
set_mismatch = set_mismatch .OR. (l_addrs_out(l_size) > 0 .AND. l_addrs_out(l_size-1) == 0)
if(l_addrs_in(l_size) > 0)then
set_mismatch = set_mismatch .OR. (isize_in /= size(field_in,1))
set_mismatch = set_mismatch .OR. (jsize_in /= size(field_in,2))
set_mismatch = set_mismatch .OR. (ke_in /= size(field_in,3))
endif
if(l_addrs_out(l_size) > 0)then
set_mismatch = set_mismatch .OR. (isize_out /= size(field_out,1))
set_mismatch = set_mismatch .OR. (jsize_out /= size(field_out,2))
set_mismatch = set_mismatch .OR. (ke_out /= size(field_out,3))
endif
if(set_mismatch)then
write( text,'(i2)' ) l_size
call mpp_error(FATAL,'MPP_REDISTRIBUTE_3D: Incompatible field at count '// &
& text//' for group redistribute.' )
endif
endif
if(do_redist)then
if(PRESENT(dc_handle))d_comm =>dc_handle ! User has kept pointer to d_comm
if(.not.ASSOCIATED(d_comm))then ! d_comm needs initialization or lookup
d_comm =>mpp_redistribute_init_comm(domain_in,l_addrs_in(1:l_size),domain_out,l_addrs_out(1:l_size), &
isize_in,jsize_in,ke_in,isize_out,jsize_out,ke_out)
if(PRESENT(dc_handle))dc_handle =>d_comm ! User wants to keep pointer to d_comm
endif
call mpp_do_redistribute( l_addrs_in(1:l_size), l_addrs_out(1:l_size), d_comm, d_type )
l_size=0; l_addrs_in=-9999; l_addrs_out=-9999
isize_in=0; jsize_in=0; ke_in=0
isize_out=0; jsize_out=0; ke_out=0
d_comm =>NULL()
endif
endif
end subroutine MPP_REDISTRIBUTE_3D_
subroutine MPP_REDISTRIBUTE_4D_( domain_in, field_in, domain_out, field_out, complete, free, list_size, &
& dc_handle, position )
type(domain2D), intent(in) :: domain_in, domain_out
MPP_TYPE_, intent(in) :: field_in (:,:,:,:)
MPP_TYPE_, intent(out) :: field_out(:,:,:,:)
logical, intent(in), optional :: complete, free
integer, intent(in), optional :: list_size
integer, intent(in), optional :: position
MPP_TYPE_ :: field3D_in (size(field_in, 1),size(field_in, 2),size(field_in ,3)*size(field_in ,4))
MPP_TYPE_ :: field3D_out(size(field_out,1),size(field_out,2),size(field_out,3)*size(field_out,4))
type(DomainCommunicator2D),pointer,optional :: dc_handle
pointer( ptr_in, field3D_in )
pointer( ptr_out, field3D_out )
field_out = 0
ptr_in = 0
ptr_out = 0
if(domain_in%initialized) ptr_in = LOC(field_in )
if(domain_out%initialized) ptr_out = LOC(field_out)
call mpp_redistribute( domain_in, field3D_in, domain_out, field3D_out, complete, free, list_size, &
& dc_handle, position )
return
end subroutine MPP_REDISTRIBUTE_4D_
subroutine MPP_REDISTRIBUTE_5D_( domain_in, field_in, domain_out, field_out, complete, free, list_size, &
& dc_handle, position )
type(domain2D), intent(in) :: domain_in, domain_out
MPP_TYPE_, intent(in) :: field_in (:,:,:,:,:)
MPP_TYPE_, intent(out) :: field_out(:,:,:,:,:)
logical, intent(in), optional :: complete, free
integer, intent(in), optional :: list_size
integer, intent(in), optional :: position
MPP_TYPE_ :: field3D_in (size(field_in, 1),size(field_in, 2), &
& size(field_in ,3)*size(field_in,4)*size(field_in ,5))
MPP_TYPE_ :: field3D_out(size(field_out,1),size(field_out,2), &
& size(field_out,3)*size(field_out,4)*size(field_out,5))
type(DomainCommunicator2D),pointer,optional :: dc_handle
pointer( ptr_in, field3D_in )
pointer( ptr_out, field3D_out )
field_out = 0
ptr_in = 0
ptr_out = 0
if(domain_in%initialized) ptr_in = LOC(field_in )
if(domain_out%initialized) ptr_out = LOC(field_out)
call mpp_redistribute( domain_in, field3D_in, domain_out, field3D_out, complete, free, list_size, &
& dc_handle, position )
return
end subroutine MPP_REDISTRIBUTE_5D_
#ifdef VECTOR_FIELD_
!VECTOR_FIELD_ is set to false for MPP_TYPE_ integer.
!vector fields
subroutine MPP_UPDATE_DOMAINS_2D_V_( fieldx, fieldy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count)
!updates data domain of 2D field whose computational domains have been computed
MPP_TYPE_, intent(inout) :: fieldx(:,:), fieldy(:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags, gridtype
logical, intent(in), optional :: complete
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
MPP_TYPE_ :: field3Dx(size(fieldx,1),size(fieldx,2),1)
MPP_TYPE_ :: field3Dy(size(fieldy,1),size(fieldy,2),1)
pointer( ptrx, field3Dx )
pointer( ptry, field3Dy )
ptrx = LOC(fieldx)
ptry = LOC(fieldy)
call mpp_update_domains( field3Dx, field3Dy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count)
return
end subroutine MPP_UPDATE_DOMAINS_2D_V_
subroutine MPP_UPDATE_DOMAINS_3D_V_( fieldx, fieldy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count)
!updates data domain of 3D field whose computational domains have been computed
MPP_TYPE_, intent(inout) :: fieldx(:,:,:), fieldy(:,:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags, gridtype
logical, intent(in), optional :: complete
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
integer :: update_whalo, update_ehalo, update_shalo, update_nhalo, ntile
integer :: grid_offset_type
logical :: exchange_uv
integer(i8_kind),dimension(MAX_DOMAIN_FIELDS, MAX_TILES),save :: f_addrsx=-9999, f_addrsy=-9999
logical :: do_update, is_complete
integer, save :: isize(2)=0,jsize(2)=0,ke=0,l_size=0, offset_type=0, list=0
integer, save :: whalosz, ehalosz, shalosz, nhalosz
integer :: tile, max_ntile
integer :: position_x, position_y
logical :: set_mismatch
character(len=3) :: text
MPP_TYPE_ :: d_type
type(overlapSpec), pointer :: updatex => NULL()
type(overlapSpec), pointer :: updatey => NULL()
type(overlapSpec), pointer :: checkx => NULL()
type(overlapSpec), pointer :: checky => NULL()
if(present(whalo)) then
update_whalo = whalo
if(abs(update_whalo) > domain%whalo ) call mpp_error(FATAL, "MPP_UPDATE_3D_V: "// &
"optional argument whalo should not be larger than the whalo when define domain.")
else
update_whalo = domain%whalo
end if
if(present(ehalo)) then
update_ehalo = ehalo
if(abs(update_ehalo) > domain%ehalo ) call mpp_error(FATAL, "MPP_UPDATE_3D_V: "// &
"optional argument ehalo should not be larger than the ehalo when define domain.")
else
update_ehalo = domain%ehalo
end if
if(present(shalo)) then
update_shalo = shalo
if(abs(update_shalo) > domain%shalo ) call mpp_error(FATAL, "MPP_UPDATE_3D_V: "// &
"optional argument shalo should not be larger than the shalo when define domain.")
else
update_shalo = domain%shalo
end if
if(present(nhalo)) then
update_nhalo = nhalo
if(abs(update_nhalo) > domain%nhalo ) call mpp_error(FATAL, "MPP_UPDATE_3D_V: "// &
"optional argument nhalo should not be larger than the nhalo when define domain.")
else
update_nhalo = domain%nhalo
end if
grid_offset_type = AGRID
if( PRESENT(gridtype) ) grid_offset_type = gridtype
exchange_uv = .false.
if(grid_offset_type == DGRID_NE) then
exchange_uv = .true.
grid_offset_type = CGRID_NE
else if( grid_offset_type == DGRID_SW ) then
exchange_uv = .true.
grid_offset_type = CGRID_SW
end if
max_ntile = domain%max_ntile_pe
ntile = size(domain%x(:))
is_complete = .true.
if(PRESENT(complete)) then
is_complete = complete
end if
tile = 1
if(max_ntile>1) then
if(ntile>MAX_TILES) then
write( text,'(i2)' ) MAX_TILES
call mpp_error(FATAL,'MPP_UPDATE_3D_V: MAX_TILES='//text//' is less than number of tiles on this pe.' )
endif
if(.NOT. present(tile_count) ) call mpp_error(FATAL, "MPP_UPDATE_3D_V: "// &
"optional argument tile_count should be present when number of tiles on some pe is more than 1")
tile = tile_count
end if
do_update = (tile == ntile) .AND. is_complete
list = list+1
if(list > MAX_DOMAIN_FIELDS)then
write( text,'(i2)' ) MAX_DOMAIN_FIELDS
call mpp_error(FATAL,'MPP_UPDATE_3D_V: MAX_DOMAIN_FIELDS='//text//' exceeded for group update.' )
endif
f_addrsx(list, tile) = LOC(fieldx)
f_addrsy(list, tile) = LOC(fieldy)
if(list == 1 .AND. tile == 1)then
isize(1)=size(fieldx,1); jsize(1)=size(fieldx,2); ke = size(fieldx,3)
isize(2)=size(fieldy,1); jsize(2)=size(fieldy,2)
offset_type = grid_offset_type
whalosz = update_whalo; ehalosz = update_ehalo; shalosz = update_shalo; nhalosz = update_nhalo
else
set_mismatch = .false.
set_mismatch = set_mismatch .OR. (isize(1) /= size(fieldx,1))
set_mismatch = set_mismatch .OR. (jsize(1) /= size(fieldx,2))
set_mismatch = set_mismatch .OR. (ke /= size(fieldx,3))
set_mismatch = set_mismatch .OR. (isize(2) /= size(fieldy,1))
set_mismatch = set_mismatch .OR. (jsize(2) /= size(fieldy,2))
set_mismatch = set_mismatch .OR. (ke /= size(fieldy,3))
set_mismatch = set_mismatch .OR. (grid_offset_type /= offset_type)
set_mismatch = set_mismatch .OR. (update_whalo /= whalosz)
set_mismatch = set_mismatch .OR. (update_ehalo /= ehalosz)
set_mismatch = set_mismatch .OR. (update_shalo /= shalosz)
set_mismatch = set_mismatch .OR. (update_nhalo /= nhalosz)
if(set_mismatch)then
write( text,'(i2)' ) list
call mpp_error(FATAL,'MPP_UPDATE_3D_V: Incompatible field at count '//text//' for group vector update.' )
end if
end if
if(is_complete) then
l_size = list
list = 0
end if
if(do_update)then
if( domain_update_is_needed(domain, update_whalo, update_ehalo, update_shalo, update_nhalo) )then
select case(grid_offset_type)
case (AGRID)
position_x = CENTER
position_y = CENTER
case (BGRID_NE, BGRID_SW)
position_x = CORNER
position_y = CORNER
case (CGRID_NE, CGRID_SW)
position_x = EAST
position_y = NORTH
case default
call mpp_error(FATAL, "mpp_update_domains2D.h: invalid value of grid_offset_type")
end select
if(debug_update_level .NE. NO_CHECK) then
checkx => search_check_overlap(domain, position_x)
checky => search_check_overlap(domain, position_y)
if(ASSOCIATED(checkx)) then
if(exchange_uv) then
call mpp_do_check(f_addrsx(1:l_size,1:ntile),f_addrsy(1:l_size,1:ntile), domain, &
checky, checkx, d_type, ke, flags, name)
else
call mpp_do_check(f_addrsx(1:l_size,1:ntile),f_addrsy(1:l_size,1:ntile), domain, &
checkx, checky, d_type, ke, flags, name)
end if
endif
endif
updatex => search_update_overlap(domain, update_whalo, update_ehalo, update_shalo, update_nhalo, &
& position_x)
updatey => search_update_overlap(domain, update_whalo, update_ehalo, update_shalo, update_nhalo, &
& position_y)
if(exchange_uv) then
call mpp_do_update(f_addrsx(1:l_size,1:ntile),f_addrsy(1:l_size,1:ntile), domain, updatey, updatex, &
d_type,ke, grid_offset_type, flags)
else
call mpp_do_update(f_addrsx(1:l_size,1:ntile),f_addrsy(1:l_size,1:ntile), domain, updatex, updatey, &
d_type,ke,grid_offset_type, flags)
end if
end if
l_size=0; f_addrsx=-9999; f_addrsy=-9999; isize=0; jsize=0; ke=0
end if
return
end subroutine MPP_UPDATE_DOMAINS_3D_V_
subroutine MPP_UPDATE_DOMAINS_4D_V_( fieldx, fieldy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count )
!updates data domain of 4D field whose computational domains have been computed
MPP_TYPE_, intent(inout) :: fieldx(:,:,:,:), fieldy(:,:,:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags, gridtype
logical, intent(in), optional :: complete
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
MPP_TYPE_ :: field3Dx(size(fieldx,1),size(fieldx,2),size(fieldx,3)*size(fieldx,4))
MPP_TYPE_ :: field3Dy(size(fieldy,1),size(fieldy,2),size(fieldy,3)*size(fieldy,4))
pointer( ptrx, field3Dx )
pointer( ptry, field3Dy )
ptrx = LOC(fieldx)
ptry = LOC(fieldy)
call mpp_update_domains( field3Dx, field3Dy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count)
return
end subroutine MPP_UPDATE_DOMAINS_4D_V_
subroutine MPP_UPDATE_DOMAINS_5D_V_( fieldx, fieldy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count )
!updates data domain of 5D field whose computational domains have been computed
MPP_TYPE_, intent(inout) :: fieldx(:,:,:,:,:), fieldy(:,:,:,:,:)
type(domain2D), intent(inout) :: domain
integer, intent(in), optional :: flags, gridtype
logical, intent(in), optional :: complete
integer, intent(in), optional :: whalo, ehalo, shalo, nhalo
character(len=*), intent(in), optional :: name
integer, intent(in), optional :: tile_count
MPP_TYPE_ :: field3Dx(size(fieldx,1),size(fieldx,2),size(fieldx,3)*size(fieldx,4)*size(fieldx,5))
MPP_TYPE_ :: field3Dy(size(fieldy,1),size(fieldy,2),size(fieldy,3)*size(fieldy,4)*size(fieldy,5))
pointer( ptrx, field3Dx )
pointer( ptry, field3Dy )
ptrx = LOC(fieldx)
ptry = LOC(fieldy)
call mpp_update_domains( field3Dx, field3Dy, domain, flags, gridtype, complete, &
whalo, ehalo, shalo, nhalo, name, tile_count)
return
end subroutine MPP_UPDATE_DOMAINS_5D_V_
#endif /* VECTOR_FIELD_ */
!> @}