forked from geoschem/FVdycoreCubed_GridComp
-
Notifications
You must be signed in to change notification settings - Fork 0
/
AppGridCreate.F90
610 lines (507 loc) · 22.9 KB
/
AppGridCreate.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
! routine to create the global edges and centers of the cubed-sphere grid
! but not the ESMF grid
subroutine AppCSEdgeCreateF(IM_WORLD, LonEdge,LatEdge, LonCenter, LatCenter, rc)
#include "MAPL_Generic.h"
use ESMF
use MAPL
use MAPL_ConstantsMod, only : pi=> MAPL_PI_R8
use fv_arrays_mod, only: REAL4, REAL8, R_GRID
use fv_grid_utils_mod, only: gnomonic_grids, cell_center2, direct_transform
use fv_grid_tools_mod, only: mirror_grid
use FV_StateMod, only: FV_Atm
implicit none
! !ARGUMENTS:
integer, intent(IN) :: IM_WORLD
integer, optional, intent(OUT) :: rc
real(ESMF_KIND_R8), intent(inout) :: LonEdge(IM_World+1,IM_World+1,6)
real(ESMF_KIND_R8), intent(inout) :: LatEdge(IM_World+1,IM_World+1,6)
real(ESMF_KIND_R8), optional, intent(inout) :: LonCenter(IM_World,IM_World)
real(ESMF_KIND_R8), optional, intent(inout) :: LatCenter(IM_World,IM_World)
! ErrLog variables
!-----------------
integer :: STATUS
character(len=ESMF_MAXSTR), parameter :: Iam="AppCSEdgeCreateF"
! Local variables
!-----------------
#ifdef EIGHT_BYTE
integer, parameter:: f_p = selected_real_kind(15) ! same as 12 on Altix
#else
! Higher precisions for grid geometrical factors:
integer, parameter:: f_p = selected_real_kind(20)
#endif
integer, parameter :: grid_type = 0
integer :: npts
integer :: ntiles=6
integer :: ndims=2
integer :: I, J, N
integer :: IG, JG
real(ESMF_KIND_R8) :: alocs(2)
real(R_GRID), allocatable :: grid_global(:,:,:,:)
integer :: L
npts = IM_World
allocate( grid_global(npts+1,npts+1,ndims,ntiles) )
call gnomonic_grids(grid_type, npts, grid_global(:,:,1,1), grid_global(:,:,2,1))
! mirror_grid assumes that the tile=1 is centered on equator and greenwich meridian Lon[-pi,pi]
call mirror_grid(grid_global, 0, npts+1, npts+1, 2, 6)
if (allocated(FV_Atm)) then
if (.not. FV_Atm(1)%flagstruct%do_schmidt) &
grid_global(:,:,1,:) = grid_global(:,:,1,:) - PI/FV_Atm(1)%flagstruct%shift_fac
else
grid_global(:,:,1,:) = grid_global(:,:,1,:) - PI/18.0
end if
where(grid_global(:,:,1,:) < 0.) &
grid_global(:,:,1,:) = grid_global(:,:,1,:) + 2.* PI
! Keep Equator and Greenwich exact
!---------------------------------
where(abs(grid_global(:,:,:,1)) < 1.d-10) grid_global(:,:,:,1) = 0.0
!---------------------------------
! Clean Up Corners
!---------------------------------
grid_global( 1,1:npts+1,:,2)=grid_global(npts+1,1:npts+1,:,1)
grid_global( 1,1:npts+1,:,3)=grid_global(npts+1:1:-1,npts+1,:,1)
grid_global(1:npts+1,npts+1,:,5)=grid_global(1,npts+1:1:-1,:,1)
grid_global(1:npts+1,npts+1,:,6)=grid_global(1:npts+1,1,:,1)
grid_global(1:npts+1, 1,:,3)=grid_global(1:npts+1,npts+1,:,2)
grid_global(1:npts+1, 1,:,4)=grid_global(npts+1,npts+1:1:-1,:,2)
grid_global(npts+1,1:npts+1,:,6)=grid_global(npts+1:1:-1,1,:,2)
grid_global( 1,1:npts+1,:,4)=grid_global(npts+1,1:npts+1,:,3)
grid_global( 1,1:npts+1,:,5)=grid_global(npts+1:1:-1,npts+1,:,3)
grid_global(npts+1,1:npts+1,:,3)=grid_global(1,1:npts+1,:,4)
grid_global(1:npts+1, 1,:,5)=grid_global(1:npts+1,npts+1,:,4)
grid_global(1:npts+1, 1,:,6)=grid_global(npts+1,npts+1:1:-1,:,4)
grid_global( 1,1:npts+1,:,6)=grid_global(npts+1,1:npts+1,:,5)
!------------------------
! Schmidt transformation:
!------------------------
if (allocated(FV_Atm)) then
if ( FV_Atm(1)%flagstruct%do_schmidt ) then
do n=1,ntiles
call direct_transform(FV_Atm(1)%flagstruct%stretch_fac, 1, npts+1, 1, npts+1, FV_Atm(1)%flagstruct%target_lon, FV_Atm(1)%flagstruct%target_lat, &
n, grid_global(1:npts+1,1:npts+1,1,n), grid_global(1:npts+1,1:npts+1,2,n))
enddo
endif
end if
if (present(LonCenter) .and. present(LatCenter)) then
do n=1,ntiles
do j=1,npts
do i=1,npts
call cell_center2(grid_global(i,j, 1:2,n), grid_global(i+1,j, 1:2,n), &
grid_global(i,j+1,1:2,n), grid_global(i+1,j+1,1:2,n), &
alocs)
jg = (n-1)*npts + j
LonCenter(i,jg) = alocs(1)
LatCenter(i,jg) = alocs(2)
enddo
enddo
enddo
end if
LonEdge = grid_global(:,:,1,:)
LatEdge = grid_global(:,:,2,:)
deallocate( grid_global )
RETURN_(ESMF_SUCCESS)
end subroutine AppCSEdgeCreateF
!!!!!!!!!!!!!!!%%%%%%%%%%%%%%%%%%%%%%%!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
function AppGridCreateF(IM_WORLD, JM_WORLD, LM, NX, NY, rc) result(esmfgrid)
#include "MAPL_Generic.h"
#define DEALLOCGLOB_(A) if(associated(A))then;A=0;if(MAPL_ShmInitialized)then; call MAPL_DeAllocNodeArray(A,rc=STATUS);else; deallocate(A,stat=STATUS);endif;VERIFY_(STATUS);NULLIFY(A);endif
use ESMF
use MAPL, pi=> MAPL_PI_R8
use fv_arrays_mod, only: REAL4, REAL8, R_GRID
use fv_grid_utils_mod, only: gnomonic_grids, cell_center2, direct_transform
use fv_grid_tools_mod, only: mirror_grid
use FV_StateMod, only: FV_Atm
implicit none
! !ARGUMENTS:
integer, intent(IN) :: IM_WORLD, JM_WORLD, LM
integer, intent(IN) :: NX, NY
integer, optional, intent(OUT) :: rc
type (ESMF_Grid) :: esmfgrid
! ErrLog variables
!-----------------
integer :: STATUS
character(len=ESMF_MAXSTR), parameter :: Iam="AppGridCreateF"
! Local variables
!-----------------
type(ESMF_DistGrid) :: distGrid
integer :: DECOUNT(2)
integer :: GLOBAL_DE_START(2)
integer :: NX2, NY2
integer :: NPES, NPES_X, NPES_Y
integer :: NPX, NPY
integer :: isg, ieg
integer :: jsg, jeg
integer :: is, ie
integer :: js, je
integer :: myTile
integer :: npts
integer :: ntiles, grid_type
integer :: ndims=2
integer :: I, J, N
integer :: IG, JG
real(REAL8) :: deglon=0.0
type (ESMF_Array), pointer :: coords(:)
real(REAL8), pointer :: lons(:,:)
real(REAL8), pointer :: lats(:,:)
type (ESMF_Array), target :: tarray(2)
real(REAL8) :: alocs(2), dx, dy
integer :: L
integer, allocatable :: IMS(:), JMS(:)
character(len=ESMF_MAXSTR) :: gridname, FMT, FMTIM, FMTJM
type (ESMF_VM) :: VM
type (ESMF_Config) :: cf
real(REAL8), allocatable :: gridCornerLats(:)
real(REAL8), allocatable :: gridCornerLons(:)
integer :: idx
real(R_GRID), pointer :: grid_global(:,:,:,:) => null()
logical :: AppGridStandAlone
! We need the VM to create the grid ??
!------------------------------------
AppGridStandAlone = .false.
call ESMF_VMGetCurrent(vm, rc=STATUS)
VERIFY_(STATUS)
! Check FV3 grid_type
! -------------------
cf = ESMF_ConfigCreate(rc=rc)
call ESMF_ConfigLoadFile( cf, 'fvcore_layout.rc', rc = rc )
call ESMF_ConfigGetAttribute( cf, ntiles, label = 'ntiles:', default=6, rc = rc )
call ESMF_ConfigGetAttribute( cf, grid_type, label = 'grid_type:', default=0, rc = rc )
call ESMF_ConfigDestroy(cf, rc=rc)
! Give the IMS and JMS the MAPL default distribution
! --------------------------------------------------
allocate( IMS(0:NX-1) )
allocate( JMS(0:NY-1) )
call GET_INT_FORMAT(IM_WORLD, FMTIM)
call GET_INT_FORMAT(JM_WORLD, FMTJM)
FMT = '(A,' // trim(FMTIM) //',A,' // trim(FMTJM) // ',A)'
if (grid_type <= 3) then
! Cubed-Sphere
call DecomposeDim ( IM_WORLD , IMS , NX )
call DecomposeDim ( JM_WORLD/6, JMS(0:NY/6 -1) , NY/6 )
do n=2,6
JMS((n-1)*NY/6 : n*NY/6 -1) = JMS(0:NY/6 -1)
enddo
write(gridname,trim(FMT)) 'PE',IM_WORLD,'x',JM_WORLD,'-CF'
else
! Doubly Periodic Cartesian
call DecomposeDim ( IM_WORLD , IMS , NX )
call DecomposeDim ( JM_WORLD , JMS , NY )
write(gridname,trim(FMT)) 'PE',IM_WORLD,'x',JM_WORLD,'-DP'
endif
! We should have a much simpler create with the next ESMF grid design
!--------------------------------------------------------------------
esmfgrid = ESMF_GridCreate( &
name=gridname, &
countsPerDEDim1=IMS, &
countsPerDEDim2=JMS, &
indexFlag = ESMF_INDEX_DELOCAL,&
gridEdgeLWidth = (/0,0/), &
gridEdgeUWidth = (/0,0/), &
coordDep1 = (/1,2/), &
coordDep2 = (/1,2/), &
rc=status)
VERIFY_(STATUS)
! Allocate coords at default stagger location
call ESMF_GridAddCoord(esmfgrid, rc=status)
VERIFY_(STATUS)
call ESMF_AttributeSet(esmfgrid, name='GRID_LM', value=LM, rc=status)
VERIFY_(STATUS)
if (grid_type <= 3) then
call ESMF_AttributeSet(esmfgrid, 'GridType', 'Cubed-Sphere', rc=STATUS)
VERIFY_(STATUS)
else
call ESMF_AttributeSet(esmfgrid, 'GridType', 'Doubly-Periodic', rc=STATUS)
VERIFY_(STATUS)
endif
! -------------
NPES_X = size(IMS)
NPES_Y = size(JMS)
NPES = NPES_X+NPES_Y
call MAPL_GRID_INTERIOR(esmfgrid,isg,ieg,jsg,jeg)
npx = IM_WORLD
npy = JM_WORLD
myTile = jsg/(npy/ntiles)
is = isg
ie = ieg
js = jsg - myTile*(npy/ntiles)
je = jeg - myTile*(npy/ntiles)
npts = (npy/ntiles)
if (npts /= npx) then
print*, 'Error npts /= npx', npts, npx
STATUS=1
endif
VERIFY_(STATUS)
if (.not.allocated(FV_Atm)) then
AppGridStandAlone = .true.
if(MAPL_ShmInitialized) then
call MAPL_AllocNodeArray(grid_global,Shp=(/npts+1,npts+1,ndims,ntiles/),rc=status)
VERIFY_(STATUS)
else
allocate( grid_global(npts+1,npts+1,ndims,ntiles) )
end if
if (grid_type <= 3) then
! Cubed-Sphere
call gnomonic_grids(grid_type, npts, grid_global(:,:,1,1), grid_global(:,:,2,1))
! mirror_grid assumes that the tile=1 is centered on equator and greenwich meridian Lon[-pi,pi]
call mirror_grid(grid_global, 0, npts+1, npts+1, 2, 6)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! ! MAT BMA You cannot refer to FV_Atm here because line 288 says this !
! ! is only for FV_Atm not allocated !
! !
! if (.not. FV_Atm(1)%flagstruct%do_schmidt) & !
! grid_global(:,:,1,:) = grid_global(:,:,1,:) - PI/FV_Atm(1)%flagstruct%shift_fac !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
grid_global(:,:,1,:) = grid_global(:,:,1,:) - PI/18.0
where(grid_global(:,:,1,:) < 0.) &
grid_global(:,:,1,:) = grid_global(:,:,1,:) + 2.* PI
! Keep Equator and Greenwich exact
!---------------------------------
where(abs(grid_global(:,:,:,1)) < 1.d-10) grid_global(:,:,:,1) = 0.0
!---------------------------------
! Clean Up Corners
!---------------------------------
grid_global( 1,1:npts+1,:,2)=grid_global(npts+1,1:npts+1,:,1)
grid_global( 1,1:npts+1,:,3)=grid_global(npts+1:1:-1,npts+1,:,1)
grid_global(1:npts+1,npts+1,:,5)=grid_global(1,npts+1:1:-1,:,1)
grid_global(1:npts+1,npts+1,:,6)=grid_global(1:npts+1,1,:,1)
grid_global(1:npts+1, 1,:,3)=grid_global(1:npts+1,npts+1,:,2)
grid_global(1:npts+1, 1,:,4)=grid_global(npts+1,npts+1:1:-1,:,2)
grid_global(npts+1,1:npts+1,:,6)=grid_global(npts+1:1:-1,1,:,2)
grid_global( 1,1:npts+1,:,4)=grid_global(npts+1,1:npts+1,:,3)
grid_global( 1,1:npts+1,:,5)=grid_global(npts+1:1:-1,npts+1,:,3)
grid_global(npts+1,1:npts+1,:,3)=grid_global(1,1:npts+1,:,4)
grid_global(1:npts+1, 1,:,5)=grid_global(1:npts+1,npts+1,:,4)
grid_global(1:npts+1, 1,:,6)=grid_global(npts+1,npts+1:1:-1,:,4)
grid_global( 1,1:npts+1,:,6)=grid_global(npts+1,1:npts+1,:,5)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! ! MAT BMA You cannot refer to FV_Atm here because line 288 says this !
! ! is only for FV_Atm not allocated !
! !
! !------------------------ !
! ! Schmidt transformation: !
! !------------------------ !
! if ( FV_Atm(1)%flagstruct%do_schmidt ) then !
! do n=1,ntiles !
! call direct_transform(FV_Atm(1)%flagstruct%stretch_fac, 1, npts+1, 1, npts+1, FV_Atm(1)%flagstruct%target_lon, FV_Atm(1)%flagstruct%target_lat, & !
! n, grid_global(1:npts+1,1:npts+1,1,n), grid_global(1:npts+1,1:npts+1,2,n)) !
! enddo !
! endif !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
else
if (MAPL_AM_I_ROOT()) write(*,*)'AppGridCreate can not make DP grid by itself'
_ASSERT(.false.,'needs informative message')
end if
end if
! Fill lat/lons at cell center locations
! Retrieve the coordinates so we can set them
call ESMF_GridGetCoord(esmfgrid, coordDim=1, localDE=0, &
staggerloc=ESMF_STAGGERLOC_CENTER, &
farrayPtr=lons, rc=status)
VERIFY_(STATUS)
call ESMF_GridGetCoord(esmfgrid, coordDim=2, localDE=0, &
staggerloc=ESMF_STAGGERLOC_CENTER, &
farrayPtr=lats, rc=status)
VERIFY_(STATUS)
! Save corners into the esmfgrid
!---------------------------
allocate(gridCornerLons((size(lons,1)+1)*(size(lons,2)+1)), stat=status)
VERIFY_(STATUS)
allocate(gridCornerLats((size(lats,1)+1)*(size(lats,2)+1)), stat=status)
VERIFY_(STATUS)
idx = 0
do jg=jsg,jeg+1
do ig=isg,ieg+1
i=ig
j=jg-myTile*npts
idx = idx + 1
if (AppGridStandAlone) then
gridCornerLons(idx) = grid_global(i, j, 1, myTile+1)
gridCornerLats(idx) = grid_global(i, j, 2, myTile+1)
else
gridCornerLons(idx) = FV_Atm(1)%gridstruct%grid(i,j,1)
gridCornerLats(idx) = FV_Atm(1)%gridstruct%grid(i,j,2)
end if
end do
end do
call ESMF_AttributeSet(esmfgrid, name='GridCornerLons:', &
itemCount = idx, valueList=gridCornerLons, rc=status)
VERIFY_(STATUS)
call ESMF_AttributeSet(esmfgrid, name='GridCornerLats:', &
itemCount = idx, valueList=gridCornerLats, rc=status)
VERIFY_(STATUS)
deallocate(gridCornerLats, stat=status)
VERIFY_(STATUS)
deallocate(gridCornerLons, stat=status)
VERIFY_(STATUS)
do jg=jsg,jeg
do ig=isg,ieg
i=ig
j=jg-myTile*npts
if (AppGridStandAlone) then
call cell_center2(grid_global(i,j, 1:2,myTile+1), grid_global(i+1,j, 1:2,myTile+1), &
grid_global(i,j+1,1:2,myTile+1), grid_global(i+1,j+1,1:2,myTile+1), &
alocs)
else
call cell_center2(dble(FV_Atm(1)%gridstruct%grid(i,j, 1:2)), &
dble(FV_Atm(1)%gridstruct%grid(i+1,j, 1:2)), &
dble(FV_Atm(1)%gridstruct%grid(i,j+1,1:2)), &
dble(FV_Atm(1)%gridstruct%grid(i+1,j+1,1:2)), &
alocs)
end if
i=ig-isg+1
j=jg-jsg+1
lons(i,j) = alocs(1)
lats(i,j) = alocs(2)
enddo
enddo
if (AppGridStandAlone) then
DEALLOCGLOB_(grid_global)
end if
call MAPL_MemUtilsWrite(VM, trim(Iam), RC=STATUS )
VERIFY_(STATUS)
RETURN_(ESMF_SUCCESS)
end function AppGridCreateF
subroutine AppGridCreate (META, esmfgrid, RC)
#include "MAPL_Generic.h"
use ESMF
use MAPL
use fv_arrays_mod, only: REAL4, REAL8, R_GRID
implicit none
! !ARGUMENTS:
type(MAPL_MetaComp), intent(INOUT) :: META
type (ESMF_Grid), intent( OUT) :: esmfgrid
integer, optional, intent( OUT) :: rc
! ErrLog variables
!-----------------
integer :: STATUS
character(len=ESMF_MAXSTR), parameter :: Iam="AppGridCreate"
! Local variables
!-----------------
integer :: IM_WORLD
integer :: JM_WORLD
integer :: LM
integer :: NX, NY
character(len=ESMF_MAXSTR) :: tmpname, gridname
character(len=ESMF_MAXSTR) :: FMT, FMTIM, FMTJM
real(REAL8) :: deglon=0.0
interface
function AppGridCreateF(IM_WORLD, JM_WORLD, LM, NX, NY, rc)
use ESMF
implicit none
type(ESMF_Grid) :: AppGridCreateF
integer, intent(IN) :: IM_WORLD, JM_WORLD, LM
integer, intent(IN) :: NX, NY
integer, optional, intent(OUT) :: rc
end function AppGridCreateF
end interface
! Get Decomposition from CF
!--------------------------
! !RESOURCE_ITEM: none :: Processing elements in 1st dimension
call MAPL_GetResource( META, NX, label ='NX:', default=1, rc = status )
VERIFY_(STATUS)
! !RESOURCE_ITEM: none :: Processing elements in 2nd dimension
call MAPL_GetResource( META, NY, label ='NY:', default=1, rc = status )
VERIFY_(STATUS)
! Get World problem size from CF
!-------------------------------
! !RESOURCE_ITEM: none :: Grid size in 1st dimension
call MAPL_GetResource( META, IM_WORLD, 'AGCM_IM:', rc = status )
VERIFY_(STATUS)
! !RESOURCE_ITEM: none :: Grid size in 2nd dimension
call MAPL_GetResource( META, JM_WORLD, 'AGCM_JM:', rc = status )
VERIFY_(STATUS)
! !RESOURCE_ITEM: none :: Grid size in 3rd dimension
call MAPL_GetResource( META, LM, 'AGCM_LM:', default=1, rc = status )
VERIFY_(STATUS)
! The grid's name is optional
!----------------------------
call GET_INT_FORMAT(IM_WORLD, FMTIM)
call GET_INT_FORMAT(JM_WORLD, FMTJM)
FMT = '(A,' // trim(FMTIM) //',A,' // trim(FMTJM) // ',A)'
write(tmpname,trim(FMT)) 'PE',IM_WORLD,'x',JM_WORLD,'-CF'
! !RESOURCE_ITEM: none :: Optional grid name
call MAPL_GetResource( META, GRIDNAME, 'AGCM_GRIDNAME:', default=trim(tmpname), rc = status )
VERIFY_(STATUS)
esmfgrid = AppGridCreateF(IM_WORLD, JM_WORLD, LM, NX, NY, STATUS)
VERIFY_(STATUS)
RETURN_(ESMF_SUCCESS)
end subroutine AppGridCreate
subroutine GET_INT_FORMAT(N, FMT)
integer :: N
character(len=*) :: FMT
IF(N < 10) THEN
FMT = 'I1'
ELSE IF (N< 100) THEN
FMT = 'I2'
ELSE IF (N< 1000) THEN
FMT = 'I3'
ELSE IF (N< 10000) THEN
FMT = 'I4'
else
FMT = 'I5'
end IF
end subroutine GET_INT_FORMAT
subroutine DecomposeDim ( dim_world,dim,NDEs )
!
! From FMS/MPP
!
implicit none
integer dim_world, NDEs
integer dim(0:NDEs-1)
integer :: is,ie,isg,ieg
integer :: ndiv,ndivs,imax,ndmax,ndmirror,n
integer :: ibegin(0:NDEs-1)
integer :: iend(0:NDEs-1)
logical :: symmetrize
logical :: even, odd
even(n) = (mod(n,2).EQ.0)
odd (n) = (mod(n,2).EQ.1)
isg = 1
ieg = dim_world
ndivs = NDEs
is = isg
n = 0
do ndiv=0,ndivs-1
!modified for mirror-symmetry
!original line
! ie = is + CEILING( float(ieg-is+1)/(ndivs-ndiv) ) - 1
!problem of dividing nx points into n domains maintaining symmetry
!i.e nx=18 n=4 4554 and 5445 are solutions but 4455 is not.
!this will always work for nx even n even or odd
!this will always work for nx odd, n odd
!this will never work for nx odd, n even: for this case we supersede the mirror calculation
! symmetrize = .NOT. ( mod(ndivs,2).EQ.0 .AND. mod(ieg-isg+1,2).EQ.1 )
!nx even n odd fails if n>nx/2
symmetrize = ( even(ndivs) .AND. even(ieg-isg+1) ) .OR. &
( odd(ndivs) .AND. odd(ieg-isg+1) ) .OR. &
( odd(ndivs) .AND. even(ieg-isg+1) .AND. ndivs.LT.(ieg-isg+1)/2 )
!mirror domains are stored in the list and retrieved if required.
if( ndiv.EQ.0 )then
!initialize max points and max domains
imax = ieg
ndmax = ndivs
end if
!do bottom half of decomposition, going over the midpoint for odd ndivs
if( ndiv.LT.(ndivs-1)/2+1 )then
!domain is sized by dividing remaining points by remaining domains
ie = is + CEILING( REAL(imax-is+1)/(ndmax-ndiv) ) - 1
ndmirror = (ndivs-1) - ndiv !mirror domain
if( ndmirror.GT.ndiv .AND. symmetrize )then !only for domains over the midpoint
!mirror extents, the max(,) is to eliminate overlaps
ibegin(ndmirror) = max( isg+ieg-ie, ie+1 )
iend(ndmirror) = max( isg+ieg-is, ie+1 )
imax = ibegin(ndmirror) - 1
ndmax = ndmax - 1
end if
else
if( symmetrize )then
!do top half of decomposition by retrieving saved values
is = ibegin(ndiv)
ie = iend(ndiv)
else
ie = is + CEILING( REAL(imax-is+1)/(ndmax-ndiv) ) - 1
end if
end if
dim(ndiv) = ie-is+1
is = ie + 1
end do
end subroutine DecomposeDim