Line data Source code
1 : module edyn_esmf
2 : use shr_kind_mod, only: r8 => shr_kind_r8
3 : use cam_logfile, only: iulog
4 : use cam_abortutils, only: endrun
5 : use infnan, only: nan, assignment(=)
6 :
7 : use ESMF, only: ESMF_Grid, ESMF_Mesh, ESMF_Field, ESMF_RouteHandle
8 : use ESMF, only: ESMF_SUCCESS
9 : use ESMF, only: ESMF_KIND_R8, ESMF_KIND_I4
10 : use ESMF, only: ESMF_FieldGet
11 : use ESMF, only: ESMF_STAGGERLOC_CENTER, ESMF_FieldRegridStore, ESMF_FieldRegrid
12 : use ESMF, only: ESMF_StaggerLoc
13 : use ESMF, only: ESMF_REGRIDMETHOD_BILINEAR, ESMF_POLEMETHOD_ALLAVG
14 : use ESMF, only: ESMF_GridCreate1PeriDim, ESMF_INDEX_GLOBAL
15 : use ESMF, only: ESMF_GridAddCoord, ESMF_GridGetCoord
16 : use ESMF, only: ESMF_TYPEKIND_R8, ESMF_FieldCreate, ESMF_Array
17 : use ESMF, only: ESMF_ArraySpec, ESMF_DistGrid, ESMF_DELayout
18 : use ESMF, only: ESMF_GridGet, ESMF_ArraySpecSet
19 : use ESMF, only: ESMF_ArrayCreate
20 : use ESMF, only: ESMF_GridComp, ESMF_TERMORDER_SRCSEQ
21 : use ESMF, only: ESMF_EXTRAPMETHOD_NEAREST_IDAVG
22 : use ESMF, only: ESMF_UNMAPPEDACTION_IGNORE
23 : use ESMF, only: ESMF_GridDestroy, ESMF_FieldDestroy, ESMF_RouteHandleDestroy
24 : use ESMF, only: ESMF_Mesh, ESMF_MeshIsCreated, ESMF_MeshDestroy
25 : use ESMF, only: ESMF_MESHLOC_ELEMENT
26 : use edyn_mpi, only: mytid, ntask, ntaski, ntaskj, tasks, lon0, lon1, lat0
27 : use edyn_mpi, only: lat1, nmagtaski, nmagtaskj, mlon0, mlon1
28 : use edyn_mpi, only: mlat0,mlat1
29 : use getapex, only: gdlatdeg, gdlondeg
30 : ! dynamically allocated geo grid for Oplus transport model
31 : use edyn_geogrid, only: nlon, nlev, glon, glat
32 : use edyn_maggrid, only: gmlat, gmlon
33 : use spmd_utils, only: masterproc
34 :
35 : implicit none
36 : save
37 : private
38 :
39 : public :: edyn_esmf_update
40 : public :: edyn_esmf_final ! Clean up any edyn usage of ESMF
41 :
42 : public :: edyn_esmf_regrid_phys2geo
43 : public :: edyn_esmf_regrid_geo2phys
44 : public :: edyn_esmf_regrid_phys2mag
45 : public :: edyn_esmf_regrid_mag2phys
46 : public :: edyn_esmf_regrid_geo2mag
47 : public :: edyn_esmf_regrid_mag2geo
48 :
49 : public :: edyn_esmf_get_1dfield
50 : public :: edyn_esmf_get_2dfield ! Retrieve a pointer to 2D ESMF field data
51 : public :: edyn_esmf_get_3dfield ! Retrieve a pointer to 3D ESMF field data
52 : public :: edyn_esmf_get_2dphysfield
53 : public :: edyn_esmf_set3d_geo ! Set ESMF field with 3D geo data
54 : public :: edyn_esmf_set2d_geo ! Set ESMF field with 2D geo data
55 : public :: edyn_esmf_set3d_mag ! Set ESMF field with 3D mag field data
56 : public :: edyn_esmf_set2d_mag ! Set ESMF field with 2D mag field data
57 : public :: edyn_esmf_set3d_phys ! Set ESMF field with 3D physics field data
58 : public :: edyn_esmf_set2d_phys ! Set ESMF field with 2D physics field data
59 : public :: edyn_esmf_update_phys_mesh
60 :
61 : public :: phys_3dfld, phys_2dfld
62 : public :: geo_3dfld, geo_2dfld
63 : public :: mag_des_3dfld, mag_des_2dfld
64 : public :: mag_src_3dfld, mag_src_2dfld
65 :
66 : public :: edyn_esmf_chkerr
67 :
68 : type(ESMF_Grid) :: &
69 : mag_src_grid, & ! source grid (will not have periodic pts)
70 : mag_des_grid, & ! destination grid (will have periodic pts)
71 : geo_grid ! geographic grid for Oplus transport
72 :
73 : ! phys_mesh: Mesh representation of physics decomposition
74 : type(ESMF_Mesh), public, protected :: phys_mesh
75 :
76 : ! ESMF fields used for mapping between physics, oplus geographic, and geomagnetic grids
77 : type(ESMF_Field) :: phys_3dfld, phys_2dfld
78 : type(ESMF_Field) :: geo_3dfld, geo_2dfld
79 : type(ESMF_Field) :: mag_des_3dfld, mag_des_2dfld
80 : type(ESMF_Field) :: mag_src_3dfld, mag_src_2dfld
81 :
82 :
83 : type(ESMF_RouteHandle) :: & ! ESMF route handles for regridding
84 : routehandle_phys2geo, & ! for physics to geo 3-D regrid
85 : routehandle_geo2phys, & ! for geo to physics 3-D regrid
86 : routehandle_phys2mag, & ! for physics to mag 3-D regrid
87 : routehandle_geo2mag, & ! for geo to mag 3-D regrid
88 : routehandle_mag2geo, & ! for geo to mag 3-D regrid
89 : routehandle_phys2mag_2d, & ! for 2d geo to phys
90 : routehandle_mag2phys_2d, & ! for 2d phys to geo for AMIE fields
91 : routehandle_geo2mag_2d ! for 2d geo to mag
92 :
93 : logical, parameter :: debug = .false.
94 :
95 : integer, allocatable :: petmap(:,:,:)
96 : logical :: initialized=.false.
97 :
98 : contains
99 :
100 : !-----------------------------------------------------------------------
101 : !-----------------------------------------------------------------------
102 861312 : subroutine edyn_esmf_chkerr(subname, routine, rc)
103 : use shr_kind_mod, only: shr_kind_cl
104 :
105 : character(len=*), intent(in) :: subname
106 : character(len=*), intent(in) :: routine
107 : integer, intent(in) :: rc
108 :
109 : character(len=shr_kind_cl) :: errmsg
110 :
111 861312 : if (rc /= ESMF_SUCCESS) then
112 0 : write(errmsg, '(4a,i0)') trim(subname), ': Error return from ', trim(routine), ', rc = ', rc
113 0 : if (masterproc) then
114 0 : write(iulog, '(2a)') 'ERROR: ', trim(errmsg)
115 : end if
116 0 : call endrun(trim(errmsg))
117 : end if
118 861312 : end subroutine edyn_esmf_chkerr
119 :
120 : !-----------------------------------------------------------------------
121 : !-----------------------------------------------------------------------
122 768 : subroutine edyn_esmf_final
123 :
124 768 : call edyn_esmf_destroy_mag_objs()
125 768 : call edyn_esmf_destroy_nonmag_objs()
126 :
127 768 : end subroutine edyn_esmf_final
128 :
129 : !-----------------------------------------------------------------------
130 : !-----------------------------------------------------------------------
131 9216 : subroutine edyn_esmf_destroy_mag_objs
132 :
133 : integer :: rc ! return code for ESMF calls
134 : character(len=*), parameter :: subname = 'edyn_esmf_destroy_mag_objs'
135 :
136 768 : call ESMF_RouteHandleDestroy(routehandle_phys2mag, rc=rc )
137 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_phys2mag', rc)
138 768 : call ESMF_RouteHandleDestroy(routehandle_geo2mag, rc=rc )
139 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_geo2mag', rc)
140 768 : call ESMF_RouteHandleDestroy(routehandle_mag2geo, rc=rc )
141 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_mag2geo', rc)
142 768 : call ESMF_RouteHandleDestroy(routehandle_phys2mag_2d, rc=rc )
143 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_phys2mag_2d', rc)
144 768 : call ESMF_RouteHandleDestroy(routehandle_mag2phys_2d, rc=rc )
145 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_mag2phys_2d', rc)
146 768 : call ESMF_RouteHandleDestroy(routehandle_geo2mag_2d, rc=rc )
147 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_geo2mag_2d', rc)
148 :
149 768 : call ESMF_FieldDestroy(mag_des_3dfld, rc=rc )
150 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy mag_des_3dfld', rc)
151 768 : call ESMF_FieldDestroy(mag_des_2dfld, rc=rc )
152 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy mag_des_2dfld', rc)
153 768 : call ESMF_FieldDestroy(mag_src_3dfld, rc=rc )
154 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy mag_src_3dfld', rc)
155 768 : call ESMF_FieldDestroy(mag_src_2dfld, rc=rc )
156 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy mag_src_2dfld', rc)
157 :
158 768 : call ESMF_GridDestroy(mag_src_grid, rc=rc)
159 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridDestroy mag_src_grid', rc)
160 768 : call ESMF_GridDestroy(mag_des_grid, rc=rc)
161 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridDestroy mag_des_grid', rc)
162 :
163 768 : end subroutine edyn_esmf_destroy_mag_objs
164 :
165 : !-----------------------------------------------------------------------
166 : !-----------------------------------------------------------------------
167 6144 : subroutine edyn_esmf_destroy_nonmag_objs
168 :
169 : integer :: rc ! return code for ESMF calls
170 : character(len=*), parameter :: subname = 'edyn_esmf_destroy_nonmag_objs'
171 :
172 768 : call ESMF_RouteHandleDestroy(routehandle_phys2geo, rc=rc )
173 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_phys2geo', rc)
174 768 : call ESMF_RouteHandleDestroy(routehandle_geo2phys, rc=rc )
175 768 : call edyn_esmf_chkerr(subname, 'ESMF_RouteHandleDestroy routehandle_geo2phys', rc)
176 :
177 768 : call ESMF_FieldDestroy(phys_3dfld, rc=rc )
178 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy phys_3dfld', rc)
179 768 : call ESMF_FieldDestroy(phys_2dfld, rc=rc )
180 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy phys_2dfld', rc)
181 768 : call ESMF_FieldDestroy(geo_3dfld, rc=rc )
182 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy geo_3dfld', rc)
183 768 : call ESMF_FieldDestroy(geo_2dfld, rc=rc )
184 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldDestroy geo_2dfld', rc)
185 :
186 768 : call ESMF_GridDestroy(geo_grid, rc=rc)
187 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridDestroy geo_grid', rc)
188 768 : call ESMF_MeshDestroy(phys_mesh, rc=rc)
189 768 : call edyn_esmf_chkerr(subname, 'ESMF_MeshDestroy phys_mesh', rc)
190 :
191 768 : end subroutine edyn_esmf_destroy_nonmag_objs
192 :
193 : !-----------------------------------------------------------------------
194 : !-----------------------------------------------------------------------
195 8064 : subroutine edyn_esmf_update
196 : use getapex, only: get_apex, magfield, alonm
197 : use mo_apex, only: geomag_year_updated
198 :
199 : ! Create ESMF grids for physics, geographic (ion transport), and
200 : ! magnetic grids, and create ESMF fields as necessary on each grid.
201 : ! Define the 2d coordinates for each grid, and save an ESMF
202 : ! routehandles for regridding.
203 : !
204 : ! Local:
205 : integer :: rc ! return code for ESMF calls
206 : integer :: lbnd_destgeo(3), ubnd_destgeo(3) ! 3d bounds of dest geo grid
207 : integer :: lbnd_destmag(3), ubnd_destmag(3) ! 3d bounds of dest mag grid
208 : integer :: lbnd_srcgeo(3), ubnd_srcgeo(3) ! 3d bounds of src geo grid
209 : integer :: lbnd_srcmag(3), ubnd_srcmag(3) ! 3d bounds of src mag grid
210 8064 : real(ESMF_KIND_R8), pointer :: fptr(:,:,:)
211 8064 : integer(ESMF_KIND_I4), pointer :: factorIndexList(:,:)
212 8064 : real(ESMF_KIND_R8), pointer :: factorList(:)
213 : integer :: smm_srctermproc, smm_pipelinedep
214 :
215 : character(len=*), parameter :: subname = 'edyn_esmf_update'
216 :
217 8064 : if (.not.geomag_year_updated .and. initialized) then
218 7296 : return
219 : end if
220 :
221 768 : if (mytid<ntask) then
222 : !
223 : ! Get apex coordinates.
224 : !
225 768 : call get_apex() ! get apex coordinates (allocates alonm)
226 768 : call magfield() ! calculate magnetic field parameters
227 :
228 : endif
229 :
230 768 : smm_srctermproc = 0
231 768 : smm_pipelinedep = 16
232 :
233 768 : if (initialized) then
234 0 : call edyn_esmf_destroy_mag_objs()
235 : endif
236 768 : if (.not.initialized) then
237 : !
238 : ! Make geographic grid for phys2geo and geo2phys regridding:
239 : !
240 768 : call create_geo_grid(geo_grid) ! geo (Oplus) grid
241 : endif
242 : !
243 : ! Make magnetic grid for phys2mag regridding:
244 : !
245 768 : call create_mag_grid(mag_des_grid, 'des') ! mag destination grid
246 : !
247 : ! Make grid for mag2phys regridding:
248 : !
249 768 : call create_mag_grid(mag_src_grid, 'src')
250 : !
251 : ! Create empty fields on geographic grid or phyiscs mesh that
252 : ! will be transformed to the magnetic grid and passed as input
253 : ! to the dynamo. This does not assign any values.
254 : !
255 : ! 3d fields (inputs to edynamo) on physics mesh for phys2mag:
256 : !
257 768 : if (.not.initialized) then
258 768 : call edyn_esmf_create_physfield(phys_2dfld, phys_mesh, 'PHYS_2DFLD', 0)
259 768 : call edyn_esmf_create_physfield(phys_3dfld, phys_mesh, 'PHYS_3DFLD', nlev)
260 :
261 768 : call edyn_esmf_create_geofield(geo_2dfld, geo_grid, 'GEO_2DFLD', 0)
262 768 : call edyn_esmf_create_geofield(geo_3dfld, geo_grid, 'GEO_3DFLD', nlev)
263 : endif
264 :
265 768 : call edyn_esmf_create_magfield(mag_des_2dfld, mag_des_grid, 'MAG_DES_2DFLD', 0)
266 768 : call edyn_esmf_create_magfield(mag_des_3dfld, mag_des_grid, 'MAG_DES_3DFLD', nlev)
267 :
268 768 : call edyn_esmf_create_magfield(mag_src_2dfld, mag_src_grid, 'MAG_SRC_2DFLD', 0)
269 768 : call edyn_esmf_create_magfield(mag_src_3dfld, mag_src_grid, 'MAG_SRC_3DFLD', nlev)
270 :
271 : if (debug .and. masterproc) then
272 : !
273 : ! Get 3d bounds of source geo field:
274 : !
275 : call ESMF_FieldGet(geo_3dfld, localDe=0, farrayPtr=fptr, &
276 : computationalLBound=lbnd_srcgeo, &
277 : computationalUBound=ubnd_srcgeo, rc=rc)
278 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet, geo_3dfld', rc)
279 :
280 : write(iulog,"(2a,i4,2(', ',i4),a,i4,2(', ',i4),a)") subname, &
281 : ': Bounds of source geo field: lbnd_srcgeo = (', &
282 : lbnd_srcgeo, '), ubnd_srcgeo = (', ubnd_srcgeo,')'
283 : end if
284 :
285 : if (debug .and. masterproc) then
286 : !
287 : ! Get 3d bounds of destination mag field:
288 : !
289 : call ESMF_FieldGet(mag_des_3dfld, localDe=0, farrayPtr=fptr, &
290 : computationalLBound=lbnd_destmag, &
291 : computationalUBound=ubnd_destmag, rc=rc)
292 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet, mag_des_3dfld', rc)
293 :
294 : write(iulog,"(2a,3i4,a,3i4,' gmlon=',2f9.3)") subname, &
295 : ': Bounds of destination mag field: lbnd_destmag = ', &
296 : lbnd_destmag, ' ubnd_destmag = ', ubnd_destmag
297 : write(iulog,"(a,': lon bnd_destmag =',2i4,' gmlon = ',2f9.3)") &
298 : subname, lbnd_destmag(1), ubnd_destmag(1), &
299 : gmlon(lbnd_destmag(1)), gmlon(ubnd_destmag(1))
300 : write(iulog,"(a,': lat bnd_destmag = ',2i4,' gmlat = ',2f9.3)") &
301 : subname, lbnd_destmag(2), ubnd_destmag(2), &
302 : gmlat(lbnd_destmag(2)), gmlat(ubnd_destmag(2))
303 : !
304 : ! Get 3d bounds of source mag field:
305 : !
306 : call ESMF_FieldGet(mag_src_3dfld, localDe=0, farrayPtr=fptr, &
307 : computationalLBound=lbnd_srcmag, &
308 : computationalUBound=ubnd_srcmag, rc=rc)
309 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet, mag_src_3dfld', rc)
310 :
311 : write(iulog,"(a,2(a,i4),a)") subname, ': lon srcmag bounds = (', &
312 : lbnd_srcmag(1), ', ', ubnd_srcmag(1), ')'
313 : write(iulog,"(a,2(a,i4),a)") subname, ': lat srcmag bounds = (', &
314 : lbnd_srcmag(2), ', ', ubnd_srcmag(2), ')'
315 : !
316 : ! Get 3d bounds of destination geo field:
317 : !
318 : call ESMF_FieldGet(geo_3dfld, localDe=0, farrayPtr=fptr, &
319 : computationalLBound=lbnd_destgeo, &
320 : computationalUBound=ubnd_destgeo, rc=rc)
321 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet, geo_3dfld', rc)
322 :
323 : write(iulog,"(a,': lon bnd_destgeo=',2i4)") subname, &
324 : lbnd_destgeo(1),ubnd_destgeo(1)
325 : write(iulog,"(a,': lat bnd_destgeo=',2i4)") subname, &
326 : lbnd_destgeo(2),ubnd_destgeo(2)
327 : end if
328 :
329 : !
330 : ! Save route handles for grid transformations in both directions
331 : ! phys2mag and mag2phys. FieldRegridStore needs to be called only
332 : ! once for each transformation before the timestep loop (src and
333 : ! dest fields are still required, so just use ped here). Once inside
334 : ! the timestep loop, the same routehandle can be used for all fields
335 : ! that are regridded in the given direction.
336 : !
337 :
338 : !
339 : ! Compute and store route handle for phys2mag 2d fields:
340 : !
341 : call ESMF_FieldRegridStore(srcField=phys_2dfld, dstField=mag_des_2dfld, &
342 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
343 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
344 : extrapMethod=ESMF_EXTRAPMETHOD_NEAREST_IDAVG, &
345 : routeHandle=routehandle_phys2mag_2d, &
346 : factorIndexList=factorIndexList, &
347 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
348 768 : pipelineDepth=smm_pipelinedep, rc=rc)
349 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 2D phys2mag', rc)
350 :
351 : !
352 : ! Compute and store route handle for phys2mag 3d fields:
353 : !
354 : call ESMF_FieldRegridStore(srcField=phys_3dfld, dstField=mag_des_3dfld, &
355 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
356 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
357 : extrapMethod=ESMF_EXTRAPMETHOD_NEAREST_IDAVG, &
358 : routeHandle=routehandle_phys2mag, factorIndexList=factorIndexList, &
359 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
360 768 : pipelineDepth=smm_pipelinedep, rc=rc)
361 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 3D phys2mag', rc)
362 : !
363 : ! Compute and store route handle for mag2phys 2d (amie) fields:
364 : !
365 : call ESMF_FieldRegridStore(srcField=mag_src_2dfld, dstField=phys_2dfld,&
366 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
367 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
368 : extrapMethod=ESMF_EXTRAPMETHOD_NEAREST_IDAVG, &
369 : routeHandle=routehandle_mag2phys_2d, &
370 : factorIndexList=factorIndexList, &
371 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
372 768 : pipelineDepth=smm_pipelinedep, rc=rc)
373 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 2D mag2phys', rc)
374 768 : if (.not.initialized) then
375 : !
376 : ! Compute and store route handle for phys2geo 3d fields:
377 : !
378 : call ESMF_FieldRegridStore(srcField=phys_3dfld, dstField=geo_3dfld, &
379 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
380 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
381 : extrapMethod=ESMF_EXTRAPMETHOD_NEAREST_IDAVG, &
382 : routeHandle=routehandle_phys2geo, factorIndexList=factorIndexList, &
383 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
384 768 : pipelineDepth=smm_pipelinedep, rc=rc)
385 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 3D phys2geo', rc)
386 :
387 : !
388 : ! Compute and store route handle for geo2phys 3d fields:
389 : !
390 : call ESMF_FieldRegridStore(srcField=geo_3dfld, dstField=phys_3dfld,&
391 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
392 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
393 : extrapMethod=ESMF_EXTRAPMETHOD_NEAREST_IDAVG, &
394 : routeHandle=routehandle_geo2phys, factorIndexList=factorIndexList, &
395 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
396 768 : pipelineDepth=smm_pipelinedep, rc=rc)
397 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 3D geo2phys', rc)
398 : endif
399 : !
400 : ! Compute and store route handle for geo2mag 3d fields:
401 : !
402 : call ESMF_FieldRegridStore(srcField=geo_3dfld, dstField=mag_des_3dfld, &
403 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
404 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
405 : routeHandle=routehandle_geo2mag, factorIndexList=factorIndexList, &
406 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
407 768 : pipelineDepth=smm_pipelinedep, rc=rc)
408 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 3D geo2mag', rc)
409 : !
410 : ! Compute and store route handle for geo2mag 2d fields:
411 : !
412 : call ESMF_FieldRegridStore(srcField=geo_2dfld, dstField=mag_des_2dfld, &
413 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
414 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
415 : routeHandle=routehandle_geo2mag_2d, &
416 : factorIndexList=factorIndexList, &
417 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
418 768 : pipelineDepth=smm_pipelinedep, rc=rc)
419 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 2D geo2mag', rc)
420 :
421 : !
422 : ! Compute and store route handle for mag2geo 3d fields:
423 : !
424 : call ESMF_FieldRegridStore(srcField=mag_src_3dfld, dstField=geo_3dfld, &
425 : regridMethod=ESMF_REGRIDMETHOD_BILINEAR, &
426 : polemethod=ESMF_POLEMETHOD_ALLAVG, &
427 : routeHandle=routehandle_mag2geo, factorIndexList=factorIndexList, &
428 : factorList=factorList, srcTermProcessing=smm_srctermproc, &
429 768 : pipelineDepth=smm_pipelinedep, rc=rc)
430 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegridStore for 3D mag2geo', rc)
431 :
432 768 : initialized=.true.
433 :
434 8064 : end subroutine edyn_esmf_update
435 :
436 : !-----------------------------------------------------------------------
437 1536 : subroutine create_mag_grid(grid_out, srcdes)
438 : !
439 : ! Create ESMF geomagnetic grid, w/ lon,lat coordinates.
440 : !
441 : ! Args:
442 : type(ESMF_Grid), intent(out) :: grid_out
443 : character(len=*), intent(in) :: srcdes
444 : !
445 : ! Local:
446 : integer :: i,j,n,rc
447 1536 : real(ESMF_KIND_R8), pointer :: coordX(:,:),coordY(:,:)
448 : integer :: lbnd(2),ubnd(2)
449 3072 : integer :: nmlons_task(ntaski) ! number of lons per task
450 3072 : integer :: nmlats_task(ntaskj) ! number of lats per task
451 : character(len=*), parameter :: subname = 'create_mag_grid'
452 :
453 : !
454 : ! We are creating either a source grid or a destination grid:
455 : !
456 1536 : if (srcdes /= 'src' .and. srcdes /= 'des') then
457 0 : write(iulog,"(a)") '>>> create_mag_grid: srcdes = ''',srcdes, &
458 0 : ''' but must be either ''src'' or ''des'''
459 0 : call endrun('create_mag_grid: srcdes')
460 : end if
461 : !
462 : ! nmlons_task(nmagtaski) = number of mag lons per task in lon dim
463 : !
464 19968 : do i = 1, nmagtaski
465 121344 : loop: do n = 0, ntask - 1
466 119808 : if (tasks(n)%magtidi == i-1) then
467 18432 : nmlons_task(i) = tasks(n)%nmaglons
468 18432 : exit loop
469 : end if
470 : end do loop
471 : end do
472 : !
473 : ! Exclude periodic points (1 point fewer for mpi tasks at east end)
474 : ! for source grids (this overwrites above for eastern-most tasks):
475 : !
476 1536 : if (srcdes == 'src') then
477 295680 : do n = 0, ntask-1
478 295680 : if (tasks(n)%magtidi == nmagtaski-1) then ! east edge of proc matrix
479 24576 : nmlons_task(tasks(n)%magtidi+1) = tasks(n)%nmaglons-1
480 : end if
481 : end do
482 : end if
483 : !
484 : ! nmlats_task(nmagtaskj) = number of mag lats per task in lat dim
485 : !
486 50688 : do j = 1, nmagtaskj
487 9192960 : loop1: do n = 0, ntask-1
488 9191424 : if (tasks(n)%magtidj == j-1) then
489 49152 : nmlats_task(j) = tasks(n)%nmaglats
490 49152 : exit loop1
491 : end if
492 : end do loop1
493 : end do
494 : !
495 : ! Create curvilinear magnetic grid (both coords depend
496 : ! on both dimensions, i.e., lon(i,j),lat(i,j)):
497 : !
498 : grid_out = ESMF_GridCreate1PeriDim( &
499 : countsPerDEDim1=nmlons_task, coordDep1=(/1,2/), &
500 : countsPerDEDim2=nmlats_task, coordDep2=(/1,2/), petmap=petmap, &
501 1536 : indexflag=ESMF_INDEX_GLOBAL,rc=rc)
502 1536 : call edyn_esmf_chkerr(subname, 'ESMF_GridCreate1PeriDim', rc)
503 : !
504 : ! Allocate coordinates:
505 : !
506 1536 : call ESMF_GridAddCoord(grid_out,staggerloc=ESMF_STAGGERLOC_CENTER,rc=rc)
507 1536 : call edyn_esmf_chkerr(subname, 'ESMF_GridAddCoord', rc)
508 1536 : if (mytid<ntask) then
509 : !
510 : ! Get pointer and set mag grid longitude coordinates:
511 : !
512 : call ESMF_GridGetCoord(grid_out, coordDim=1, localDE=0, &
513 : computationalLBound=lbnd, computationalUBound=ubnd, &
514 1536 : farrayPtr=coordX, staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
515 1536 : call edyn_esmf_chkerr(subname, 'ESMF_GridGetCoord for longitude coords', rc)
516 :
517 6192 : do j = lbnd(2), ubnd(2)
518 37426 : do i = lbnd(1), ubnd(1)
519 35890 : coordX(i,j) = gdlondeg(i,j)
520 : end do
521 : end do
522 : !
523 : ! Get pointer and set mag grid latitude coordinates:
524 : !
525 : call ESMF_GridGetCoord(grid_out, coordDim=2, localDE=0, &
526 : computationalLBound=lbnd, computationalUBound=ubnd, &
527 1536 : farrayPtr=coordY, staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
528 1536 : call edyn_esmf_chkerr(subname, 'ESMF_GridGetCoord for latitude coords', rc)
529 :
530 6192 : do j = lbnd(2), ubnd(2)
531 37426 : do i = lbnd(1), ubnd(1)
532 35890 : coordY(i,j) = gdlatdeg(i,j)
533 : end do
534 : end do
535 :
536 : if (debug .and. masterproc) then
537 : write(iulog,"(3a,2i4,a,2i4,a,2i4,a,2i4)") 'Created ESMF ',srcdes, &
538 : ' mag grid: lbnd, ubnd_lon=', lbnd(1), ubnd(1), &
539 : ' mlon0,1=', mlon0, mlon1, ' lbnd, ubnd_lat=', &
540 : lbnd(2), ubnd(2), ' mlat0,1=', mlat0, mlat1
541 : end if
542 : end if
543 :
544 1536 : end subroutine create_mag_grid
545 :
546 : !-----------------------------------------------------------------------
547 768 : subroutine create_geo_grid(grid_out)
548 : !
549 : ! Args:
550 : type(ESMF_Grid),intent(out) :: grid_out
551 : !
552 : ! Local:
553 : integer :: i, j, n, rc
554 : integer :: lbnd_lat, ubnd_lat, lbnd_lon, ubnd_lon
555 : integer :: lbnd(1), ubnd(1)
556 1536 : integer :: nlons_task(ntaski) ! # lons per task
557 1536 : integer :: nlats_task(ntaskj) ! # lats per task
558 768 : real(ESMF_KIND_R8), pointer :: coordX(:), coordY(:)
559 : character(len=*), parameter :: subname = 'create_geo_grid'
560 :
561 : !
562 : ! nlons_task(ntaski) = number of geo lons per task.
563 : !
564 9984 : do i = 1, ntaski
565 60672 : loop: do n = 0, ntask-1
566 59904 : if (tasks(n)%mytidi == i-1) then
567 9216 : nlons_task(i) = tasks(n)%nlons
568 9216 : exit loop
569 : end if
570 : end do loop
571 : end do
572 : !
573 25344 : do j = 1, ntaskj
574 4596480 : loop1: do n = 0, ntask-1
575 4595712 : if (tasks(n)%mytidj == j-1) then
576 24576 : nlats_task(j) = tasks(n)%nlats
577 24576 : exit loop1
578 : end if
579 : end do loop1
580 : end do
581 : !
582 : ! Create 2d geographic source grid (with poles)
583 : grid_out = ESMF_GridCreate1PeriDim( &
584 : countsPerDEDim1=nlons_task, coordDep1=(/1/), &
585 : countsPerDEDim2=nlats_task, coordDep2=(/2/), petmap=petmap, &
586 768 : indexflag=ESMF_INDEX_GLOBAL,minIndex=(/1,1/), rc=rc)
587 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridCreate1PeriDim', rc)
588 : !
589 : ! Allocate coordinates:
590 : !
591 768 : call ESMF_GridAddCoord(grid_out, staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
592 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridAddCoord', rc)
593 : !
594 : ! Get pointer and set geo grid longitude coordinates:
595 : !
596 768 : if (mytid<ntask) then
597 : call ESMF_GridGetCoord(grid_out, coordDim=1, localDE=0, &
598 : computationalLBound=lbnd, computationalUBound=ubnd, &
599 768 : farrayPtr=coordX, staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
600 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridGetCoord for longitude coords', rc)
601 : !
602 : ! Note glon was shifted to +/-180 by sub set_geogrid (edyn_init.F90)
603 : !
604 768 : lbnd_lon = lbnd(1)
605 768 : ubnd_lon = ubnd(1)
606 9984 : do i = lbnd_lon, ubnd_lon
607 9984 : coordX(i) = glon(i) ! 1 -> 72
608 : end do
609 : !
610 : ! Get pointer and set geo grid latitude coordinates, including poles:
611 : !
612 : call ESMF_GridGetCoord(grid_out, coordDim=2, localDE=0, &
613 : computationalLBound=lbnd, computationalUBound=ubnd, &
614 768 : farrayPtr=coordY, staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
615 768 : call edyn_esmf_chkerr(subname, 'ESMF_GridGetCoord for latitude coords', rc)
616 :
617 768 : lbnd_lat = lbnd(1)
618 768 : ubnd_lat = ubnd(1)
619 3072 : do i = lbnd_lat, ubnd_lat
620 3072 : coordY(i) = glat(i)
621 : end do
622 :
623 : if (debug .and. masterproc) then
624 : write(iulog,"(2a,2i4,a,2i4,a,2i4,a,2i4)") 'Created ESMF geo_grid:', &
625 : ' lbnd,ubnd_lon=', lbnd_lon, ubnd_lon, ' lon0,1=', lon0, lon1, &
626 : ' lbnd,ubnd_lat=', lbnd_lat, ubnd_lat, ' lat0,1=', lat0, lat1
627 : write(iulog,"('coordX for geo grid = ',/,(8f10.4))") coordX
628 : write(iulog,"('coordY for geo grid = ',/,(8f10.4))") coordY
629 : end if
630 : endif
631 :
632 768 : end subroutine create_geo_grid
633 :
634 : !-----------------------------------------------------------------------
635 1536 : subroutine edyn_esmf_create_physfield(field, mesh, name, nlev)
636 : !
637 : ! Create ESMF field (2d or 3d) on physics mesh
638 : ! If nlev == 0, field is 2d (i,j), otherwise field is 3d,
639 : ! and 3rd dimension is ungridded
640 : !
641 : ! Args:
642 : integer, intent(in) :: nlev ! if nlev == 0, field is 2d (i,j)
643 : type(ESMF_Mesh), intent(in) :: mesh
644 : character(len=*), intent(in) :: name
645 : type(ESMF_Field), intent(out) :: field
646 : !
647 : ! Local:
648 : integer :: rc
649 : type(ESMF_ArraySpec) :: arrayspec
650 : character(len=*), parameter :: subname = 'edyn_esmf_create_physfield'
651 :
652 :
653 : ! Create 3d field (i,j,k), with non-distributed vertical dimension:
654 1536 : if (nlev > 0) then
655 768 : call ESMF_ArraySpecSet(arrayspec, 2, ESMF_TYPEKIND_R8, rc=rc)
656 768 : call edyn_esmf_chkerr(subname, 'ESMF_ArraySpecSet 2D', rc)
657 : field = ESMF_FieldCreate(mesh, arrayspec, &
658 : gridToFieldMap=(/2/), meshloc=ESMF_MESHLOC_ELEMENT, &
659 1536 : ungriddedLBound=(/1/), ungriddedUBound=(/nlev/), rc=rc)
660 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldCreate 2D field', rc)
661 : !
662 : ! Create 2d field (i,j):
663 : else ! create 2d field
664 768 : call ESMF_ArraySpecSet(arrayspec, 1, ESMF_TYPEKIND_R8, rc=rc)
665 768 : call edyn_esmf_chkerr(subname, 'ESMF_ArraySpecSet 2D', rc)
666 768 : field = ESMF_FieldCreate(mesh, arrayspec, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
667 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldCreate 2D field', rc)
668 : end if
669 :
670 1536 : end subroutine edyn_esmf_create_physfield
671 : !-----------------------------------------------------------------------
672 1536 : subroutine edyn_esmf_create_geofield(field, grid, name, nlev)
673 : !
674 : ! Create ESMF field (2d or 3d) on geo grid (will exclude periodic points)
675 : ! If nlev == 0, field is 2d (i,j), otherwise field is 3d,
676 : ! and 3rd dimension is ungridded
677 : !
678 : ! Args:
679 : integer, intent(in) :: nlev ! if nlev == 0, field is 2d (i,j)
680 : type(ESMF_Grid), intent(in) :: grid
681 : character(len=*), intent(in) :: name
682 : type(ESMF_Field), intent(out) :: field
683 : !
684 : ! Local:
685 : integer :: rc
686 : type(ESMF_ArraySpec) :: arrayspec
687 : character(len=*), parameter :: subname = 'edyn_esmf_create_geofield'
688 : !
689 : ! Create 3d field (i,j,k), with non-distributed vertical dimension:
690 1536 : if (nlev > 0) then
691 768 : call ESMF_ArraySpecSet(arrayspec,3,ESMF_TYPEKIND_R8,rc=rc)
692 768 : call edyn_esmf_chkerr(subname, 'ESMF_ArraySpecSet 3D', rc)
693 : field = ESMF_FieldCreate(grid, arrayspec,ungriddedLBound=(/1/), &
694 1536 : ungriddedUBound=(/nlev/),staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
695 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldCreate 3D field', rc)
696 : !
697 : ! Create 2d field (i,j):
698 : else ! create 2d field
699 768 : call ESMF_ArraySpecSet(arrayspec, 2, ESMF_TYPEKIND_R8, rc=rc)
700 768 : call edyn_esmf_chkerr(subname, 'ESMF_ArraySpecSet 2D', rc)
701 : field = ESMF_FieldCreate(grid, arrayspec,&
702 768 : staggerloc=ESMF_STAGGERLOC_CENTER,rc=rc)
703 768 : call edyn_esmf_chkerr(subname, 'ESMF_FieldCreate 2D field', rc)
704 : end if
705 1536 : end subroutine edyn_esmf_create_geofield
706 : !-----------------------------------------------------------------------
707 3072 : subroutine edyn_esmf_create_magfield(field, grid, name, nlev)
708 : !
709 : ! Create ESMF field (2d or 3d) on mag grid. This will include the
710 : ! mag periodic point, which will be zero after regridding.
711 : ! If nlev == 0, field is 2d (i,j), otherwise field is 3d,
712 : ! and 3rd dimension is ungridded
713 : !
714 : ! Args:
715 : integer, intent(in) :: nlev ! if nlev == 0, field is 2d (i,j)
716 : type(ESMF_Grid), intent(in) :: grid
717 : character(len=*), intent(in) :: name
718 : type(ESMF_Field), intent(out) :: field
719 : !
720 : ! Local:
721 : integer :: rc
722 : type(ESMF_ArraySpec) :: arrayspec
723 : type(ESMF_Array) :: array3d,array2d
724 : type(ESMF_DistGrid) :: distgrid
725 : character(len=*), parameter :: subname = 'edyn_esmf_create_magfield'
726 : !
727 : ! Get necessary information from the mag grid:
728 : call ESMF_GridGet(grid, staggerloc=ESMF_STAGGERLOC_CENTER, &
729 3072 : distgrid=distgrid,rc=rc)
730 3072 : call edyn_esmf_chkerr(subname, 'ESMF_GridGet', rc)
731 : !
732 : ! Create 3d mag field (i,j,k), with non-distributed vertical dimension:
733 : ! (add periodic point in longitude with computationalEdgeUWidth)
734 : !
735 3072 : if (nlev > 0) then
736 1536 : call ESMF_ArraySpecSet(arrayspec,3,ESMF_TYPEKIND_R8,rc=rc)
737 1536 : call edyn_esmf_chkerr(subname, 'ESMF_ArraySpecSet 3D field', rc)
738 :
739 : array3d = ESMF_ArrayCreate(arrayspec=arrayspec, &
740 : distgrid=distgrid,computationalEdgeUWidth=(/1,0/), &
741 : undistLBound=(/1/),undistUBound=(/nlev/), &
742 3072 : indexflag=ESMF_INDEX_GLOBAL,rc=rc)
743 1536 : call edyn_esmf_chkerr(subname, 'ESMF_ArrayCreate 3D field', rc)
744 :
745 : field = ESMF_FieldCreate(grid, array3d, &
746 : ungriddedLBound=(/1/), ungriddedUBound=(/nlev/), &
747 3072 : staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
748 1536 : call edyn_esmf_chkerr(subname, 'ESMF_FieldCreate 3D field', rc)
749 : !
750 : ! Create 2d mag field (i,j):
751 : ! (add periodic point in longitude with computationalEdgeUWidth)
752 : !
753 : else ! create 2d field
754 1536 : call ESMF_ArraySpecSet(arrayspec,2,ESMF_TYPEKIND_R8,rc=rc)
755 1536 : call edyn_esmf_chkerr(subname, 'ESMF_ArraySpecSet 2D field', rc)
756 :
757 : array2d = ESMF_ArrayCreate(arrayspec=arrayspec, &
758 : distgrid=distgrid,computationalEdgeUWidth=(/1,0/), &
759 1536 : indexflag=ESMF_INDEX_GLOBAL,rc=rc)
760 1536 : call edyn_esmf_chkerr(subname, 'ESMF_ArrayCreate 2D field', rc)
761 : field = ESMF_FieldCreate(grid, array2d, &
762 1536 : staggerloc=ESMF_STAGGERLOC_CENTER, rc=rc)
763 1536 : call edyn_esmf_chkerr(subname, 'ESMF_FieldCreate 2D field', rc)
764 : end if
765 3072 : end subroutine edyn_esmf_create_magfield
766 :
767 : !-----------------------------------------------------------------------
768 51072 : subroutine edyn_esmf_set3d_geo(field, fdata, ilon0, ilon1, ilat0, ilat1, ilev0, ilev1 )
769 : !
770 : ! Set values of a 3d ESMF field on geographic source grid, prior to
771 : ! geographic to physics grid transformation.
772 : ! Periodic points are excluded, geographic poles are at
773 : ! j==jspole and jnpole
774 : ! Note dimension order changes from input (k,i,j) to output (i,j,k).
775 : !
776 : ! Args:
777 : type(ESMF_Field), intent(in) :: field ! esmf fields on geo grid
778 : !
779 : ! field is input data on model subdomains (including periodic points)
780 : ! (note esmf source field excludes periodic points)
781 : !
782 : integer, intent(in) :: ilev0, ilev1, ilon0, ilon1, ilat0, ilat1
783 : real(r8), intent(in) :: fdata(ilon0:ilon1,ilat0:ilat1,ilev0:ilev1)
784 : !
785 : ! Local:
786 : integer :: i, j, k, rc
787 : integer :: lbnd(3), ubnd(3) ! 3d field bounds
788 : !
789 : ! fptr is esmf pointer (i,j,k) to 3d field, set by this subroutine
790 51072 : real(ESMF_KIND_R8), pointer :: fptr(:,:,:)
791 : character(len=*), parameter :: subname = 'edyn_esmf_set3d_geo'
792 51072 : if (mytid<ntask) then
793 : !
794 : ! Get and set pointer to the field:
795 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
796 51072 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
797 51072 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet field', rc)
798 :
799 51072 : fptr = nan
800 6690432 : do k = lbnd(3),ubnd(3)
801 26608512 : do j = lbnd(2),ubnd(2)
802 265574400 : do i = lbnd(1),ubnd(1)
803 258935040 : fptr(i,j,k) = fdata(i,j,k)
804 : end do
805 : end do
806 : end do
807 : endif
808 :
809 51072 : end subroutine edyn_esmf_set3d_geo
810 : !-----------------------------------------------------------------------
811 36480 : subroutine edyn_esmf_set2d_geo(field, fdata, ilon0, ilon1, ilat0, ilat1)
812 : !
813 : ! Set values of a 2d ESMF field on geographic source grid, prior to
814 : ! geographic to physics grid transformation. (Essentially the same
815 : ! as esmf_set3d_geo, except for 2d fields instead of 3d)
816 : ! Periodic points are excluded, geographic poles are at j==jspole and jnpole
817 : !
818 : ! Args:
819 : type(ESMF_Field), intent(in) :: field
820 : integer, intent(in) :: ilon0,ilon1,ilat0,ilat1
821 : real(r8), intent(in) :: fdata(ilon0:ilon1,ilat0:ilat1)
822 : !
823 : ! Local:
824 : integer :: i, j, rc
825 36480 : real(ESMF_KIND_R8), pointer :: fptr(:,:)
826 : integer :: lbnd(2), ubnd(2)
827 : character(len=*), parameter :: subname = 'edyn_esmf_set2d_geo'
828 36480 : if (mytid<ntask) then
829 : !
830 : ! Get pointer to the field:
831 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
832 36480 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
833 36480 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
834 : !
835 36480 : fptr(:,:) = nan
836 : !
837 : ! Set interior latitudes (excluding poles):
838 145920 : do j = lbnd(2), ubnd(2)
839 1459200 : do i = lbnd(1), ubnd(1)
840 1422720 : fptr(i,j) = fdata(i,j)
841 : end do
842 : end do
843 : endif
844 36480 : end subroutine edyn_esmf_set2d_geo
845 :
846 : !-----------------------------------------------------------------------
847 29184 : subroutine edyn_esmf_set3d_mag(field, fdata, ilon0, ilon1, ilat0, ilat1, ilev0, ilev1 )
848 : !
849 : ! Set values of a 3d ESMF field on magnetic grid, prior to magnetic to
850 : ! physics grid transformation.
851 : !
852 : ! Args:
853 : type(ESMF_Field), intent(in) :: field ! esmf fields on mag grid
854 : !
855 : ! fdata is input data on model subdomains:
856 : !
857 : integer, intent(in) :: ilev0, ilev1, ilon0, ilon1, ilat0, ilat1
858 : real(r8), intent(in) :: fdata(ilon0:ilon1,ilat0:ilat1,ilev0:ilev1)
859 : !
860 : ! Local:
861 : integer :: i, j, k, rc
862 : integer :: lbnd(3), ubnd(3) ! 3d field bounds
863 : !
864 : ! fptr is esmf pointer (i,j,k) to 3d field, set by this subroutine
865 29184 : real(ESMF_KIND_R8), pointer :: fptr(:,:,:)
866 : character(len=*), parameter :: subname = 'edyn_esmf_set3d_mag'
867 29184 : if (mytid<ntask) then
868 :
869 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
870 29184 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
871 29184 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
872 : !
873 29184 : fptr(:,:,:) = nan
874 : !
875 : ! Set ESMF pointer:
876 : !
877 3823104 : do k = lbnd(3), ubnd(3)
878 15323424 : do j = lbnd(2), ubnd(2)
879 91963040 : do i = lbnd(1), ubnd(1)
880 88169120 : fptr(i,j,k) = fdata(i,j,k)
881 : end do
882 : end do
883 : end do
884 : endif
885 :
886 29184 : end subroutine edyn_esmf_set3d_mag
887 : !-----------------------------------------------------------------------
888 :
889 0 : subroutine edyn_esmf_set2d_mag(field, fdata, ilon0, ilon1, ilat0, ilat1)
890 : !
891 : ! Set values of a 2d ESMF field on magnetic grid, prior to magnetic to
892 : ! physics grid transformation.
893 : !
894 : ! Args:
895 : type(ESMF_Field), intent(in) :: field ! esmf field on mag grid
896 : !
897 : ! fdata is input data on model subdomains:
898 : !
899 : integer, intent(in) :: ilon0, ilon1, ilat0, ilat1
900 : real(r8), intent(in) :: fdata(ilon0:ilon1,ilat0:ilat1)
901 : !
902 : ! Local:
903 : integer :: i, j, rc
904 : integer :: lbnd(2), ubnd(2) ! 2d field bounds
905 : !
906 : ! fptr is esmf pointer (i,j,k) to 2d field, set by this subroutine
907 0 : real(ESMF_KIND_R8), pointer :: fptr(:,:)
908 : character(len=*), parameter :: subname = 'edyn_esmf_set2d_mag'
909 0 : if (mytid<ntask) then
910 : !
911 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
912 0 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
913 0 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
914 : !
915 0 : fptr(:,:) = nan
916 : !
917 : ! Set ESMF pointer:
918 : !
919 0 : do j = lbnd(2), ubnd(2)
920 0 : do i = lbnd(1), ubnd(1)
921 0 : fptr(i,j) = fdata(i,j)
922 : end do
923 : end do
924 : !
925 : endif
926 :
927 0 : end subroutine edyn_esmf_set2d_mag
928 : !-----------------------------------------------------------------------
929 153216 : subroutine edyn_esmf_set3d_phys(field, fdata, ilev0, ilev1, icol0, icol1)
930 : !
931 : ! Set values ESMF field on physics mesh
932 : !
933 : ! Args:
934 : type(ESMF_Field), intent(in) :: field ! esmf field on phys grid
935 : !
936 : ! fdata is input data on model subdomains:
937 : !
938 : integer, intent(in) :: ilev0, ilev1, icol0, icol1
939 : real(r8), intent(in) :: fdata(ilev0:ilev1,icol0:icol1)
940 : !
941 : ! Local:
942 : integer :: i, k, rc
943 : integer :: lbnd(2), ubnd(2) ! 2d field bounds
944 : !
945 : ! fptr is esmf pointer to 2d field, set by this subroutine
946 153216 : real(ESMF_KIND_R8), pointer :: fptr(:,:)
947 : character(len=*), parameter :: subname = 'edyn_esmf_set2d_phys'
948 : !
949 : ! Fields loop:
950 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
951 : computationalLBound=lbnd, computationalUBound=ubnd, &
952 153216 : rc=rc)
953 153216 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
954 : !
955 153216 : fptr(:,:) = nan
956 : !
957 : ! Set ESMF pointer:
958 : !
959 5668992 : do i = lbnd(2), ubnd(2)
960 722719872 : do k = lbnd(1), ubnd(1)
961 722566656 : fptr(k,i) = fdata(k,i)
962 : end do
963 : end do
964 153216 : end subroutine edyn_esmf_set3d_phys
965 : !-----------------------------------------------------------------------
966 : !
967 0 : subroutine edyn_esmf_set2d_phys(field, fdata, icol0, icol1)
968 : !
969 : ! Set values of a ESMF field on phys mesh
970 : !
971 : ! Args:
972 : type(ESMF_Field), intent(in) :: field ! esmf field on mag grid
973 : !
974 : ! fdata is input data on model subdomains:
975 : !
976 : integer, intent(in) :: icol0, icol1
977 : real(r8), intent(in) :: fdata(icol0:icol1)
978 : !
979 : ! Local:
980 : integer :: i, rc
981 : integer :: lbnd(1), ubnd(1) ! 2d field bounds
982 : !
983 : ! fptr is esmf 1d field, set by this subroutine
984 0 : real(ESMF_KIND_R8), pointer :: fptr(:)
985 : character(len=*), parameter :: subname = 'edyn_esmf_set1d_phys'
986 : !
987 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
988 0 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
989 0 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
990 : !
991 0 : fptr(:) = nan
992 : !
993 : ! Set ESMF pointer:
994 : !
995 0 : do i = lbnd(1), ubnd(1)
996 0 : fptr(i) = fdata(i)
997 : end do
998 : !
999 0 : end subroutine edyn_esmf_set2d_phys
1000 : !-----------------------------------------------------------------------
1001 196992 : subroutine edyn_esmf_get_3dfield(field, data, i0,i1,j0,j1,k0,k1 )
1002 : !
1003 : ! Get pointer to 3d esmf field (i,j,k):
1004 : !
1005 : ! Args:
1006 : integer, intent(in) :: i0,i1,j0,j1,k0,k1
1007 : type(ESMF_field), intent(in) :: field
1008 : real(r8), intent(out) :: data(i0:i1,j0:j1,k0:k1)
1009 : !
1010 : ! Local:
1011 196992 : real(r8), pointer :: fptr(:,:,:)
1012 : integer :: rc, lbnd(3), ubnd(3)
1013 : character(len=*), parameter :: subname = 'edyn_esmf_get_3dfield'
1014 :
1015 196992 : if (mytid<ntask) then
1016 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
1017 196992 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
1018 196992 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
1019 :
1020 196992 : if (i0/=lbnd(1).or.i1/=ubnd(1).or.j0/=lbnd(2).or.j1/=ubnd(2).or.k0/=lbnd(3).or.k1/=ubnd(3)) then
1021 0 : call endrun(subname//' array bnds do not match')
1022 : endif
1023 :
1024 951011142 : data(:,:,:) = fptr(:,:,:)
1025 : endif
1026 :
1027 196992 : end subroutine edyn_esmf_get_3dfield
1028 : !-----------------------------------------------------------------------
1029 36480 : subroutine edyn_esmf_get_2dfield(field, data, i0,i1,j0,j1 )
1030 : !
1031 : ! Get pointer to 2d esmf field (i,j):
1032 : !
1033 : ! Args:
1034 : integer, intent(in) :: i0,i1,j0,j1
1035 : type(ESMF_field), intent(in) :: field
1036 : real(r8), intent(out) :: data(i0:i1,j0:j1)
1037 : !
1038 : ! Local:
1039 36480 : real(r8), pointer :: fptr(:,:)
1040 : integer :: rc, lbnd(2), ubnd(2)
1041 : character(len=*), parameter :: subname = 'edyn_esmf_get_2dfield'
1042 36480 : if (mytid<ntask) then
1043 :
1044 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
1045 36480 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
1046 36480 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
1047 :
1048 36480 : if (i0/=lbnd(1).or.i1/=ubnd(1).or.j0/=lbnd(2).or.j1/=ubnd(2)) then
1049 0 : call endrun(subname//' array bnds do not match')
1050 : endif
1051 :
1052 893475 : data(:,:) = fptr(:,:)
1053 : endif
1054 36480 : end subroutine edyn_esmf_get_2dfield
1055 : !-----------------------------------------------------------------------
1056 36480 : subroutine edyn_esmf_get_2dphysfield(field, data, i0,i1,j0,j1 )
1057 : !
1058 : ! Get pointer to 2d esmf field (i,j):
1059 : !
1060 : ! Args:
1061 : integer, intent(in) :: i0,i1,j0,j1
1062 : type(ESMF_field), intent(in) :: field
1063 : real(r8), intent(out) :: data(i0:i1,j0:j1)
1064 : !
1065 : ! Local:
1066 36480 : real(r8), pointer :: fptr(:,:)
1067 : integer :: rc, lbnd(2), ubnd(2)
1068 : character(len=*), parameter :: subname = 'edyn_esmf_get_2dphysfield'
1069 :
1070 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
1071 36480 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
1072 36480 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
1073 :
1074 36480 : if (i0/=lbnd(1).or.i1/=ubnd(1).or.j0/=lbnd(2).or.j1/=ubnd(2)) then
1075 0 : call endrun(subname//' array bnds do not match')
1076 : endif
1077 :
1078 172076160 : data(:,:) = fptr(:,:)
1079 :
1080 36480 : end subroutine edyn_esmf_get_2dphysfield
1081 : !-----------------------------------------------------------------------
1082 0 : subroutine edyn_esmf_get_1dfield(field, data, i0,i1 )
1083 : !
1084 : ! Get pointer to 2d esmf field (i,j):
1085 : !
1086 : ! Args:
1087 : integer, intent(in) :: i0,i1
1088 : type(ESMF_field), intent(in) :: field
1089 : real(r8), intent(out) :: data(i0:i1)
1090 : !
1091 : ! Local:
1092 0 : real(r8), pointer :: fptr(:)
1093 : integer :: rc, lbnd(1), ubnd(1)
1094 : character(len=*), parameter :: subname = 'edyn_esmf_get_1dfield'
1095 :
1096 : call ESMF_FieldGet(field, localDe=0, farrayPtr=fptr, &
1097 0 : computationalLBound=lbnd, computationalUBound=ubnd, rc=rc)
1098 0 : call edyn_esmf_chkerr(subname, 'ESMF_FieldGet', rc)
1099 :
1100 0 : if (i0/=lbnd(1).or.i1/=ubnd(1)) then
1101 0 : call endrun(subname//' array bnds do not match')
1102 : endif
1103 :
1104 0 : data(:) = fptr(:)
1105 :
1106 0 : end subroutine edyn_esmf_get_1dfield
1107 : !-----------------------------------------------------------------------
1108 21888 : subroutine edyn_esmf_regrid_phys2mag(srcfield, dstfield, ndim)
1109 : !
1110 : ! Args:
1111 : integer :: ndim
1112 : type(ESMF_Field), intent(inout) :: srcfield, dstfield
1113 : !
1114 : ! Local:
1115 : integer :: rc
1116 : character(len=*), parameter :: subname = 'edyn_esmf_regrid_phys2mag'
1117 : !
1118 21888 : if (ndim == 2) then
1119 : !
1120 : ! Do sparse matrix multiply for 2d phys2mag.
1121 : !
1122 : call ESMF_FieldRegrid(srcfield, dstfield, routehandle_phys2mag_2d, &
1123 0 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1124 0 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid phys2mag 2D', rc)
1125 : else ! 3d geo2mag
1126 : !
1127 : ! Do sparse matrix multiply for 3d geo2mag.
1128 : !
1129 : call ESMF_FieldRegrid(srcfield, dstfield, routehandle_phys2mag, &
1130 21888 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1131 21888 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid phys2mag 3D', rc)
1132 : end if
1133 21888 : end subroutine edyn_esmf_regrid_phys2mag
1134 : !-----------------------------------------------------------------------
1135 0 : subroutine edyn_esmf_regrid_mag2phys(srcfield, dstfield, ndim)
1136 : !
1137 : ! Args:
1138 : type(ESMF_Field), intent(inout) :: srcfield, dstfield
1139 : integer :: ndim
1140 : !
1141 : ! Local:
1142 : integer :: rc
1143 : character(len=*), parameter :: subname = 'edyn_esmf_regrid_mag2phys'
1144 : !
1145 0 : if (ndim == 2) then
1146 : call ESMF_FieldRegrid(srcfield, dstfield, routehandle_mag2phys_2d, &
1147 0 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1148 0 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid mag2phys 2D', rc)
1149 : else
1150 : ! call ESMF_FieldRegrid(srcfield, dstfield, routehandle_mag2phys, &
1151 : ! termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1152 : ! call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid mag2phys 3D', rc)
1153 0 : call endrun(subname//': routehandle_mag2phys not implemented')
1154 : end if
1155 0 : end subroutine edyn_esmf_regrid_mag2phys
1156 : !-----------------------------------------------------------------------
1157 131328 : subroutine edyn_esmf_regrid_phys2geo(srcfield, dstfield, ndim)
1158 : !
1159 : ! Args:
1160 : integer :: ndim
1161 : type(ESMF_Field), intent(inout) :: srcfield, dstfield
1162 : !
1163 : ! Local:
1164 : integer :: rc
1165 : character(len=*), parameter :: subname = 'edyn_esmf_regrid_phys2geo'
1166 : !
1167 131328 : if (ndim == 2) then
1168 : ! call ESMF_FieldRegrid( srcfield, dstfield, routehandle_phys2geo_2d, &
1169 : ! termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1170 : ! call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid phys2geo 2D', rc)
1171 0 : call endrun(subname//': routehandle_phys2geo_2d not implemented')
1172 : else ! 3d phys2geo
1173 : !
1174 : ! Do sparse matrix multiply for 3d phys2geo.
1175 : !
1176 : call ESMF_FieldRegrid( srcfield, dstfield, routehandle_phys2geo, &
1177 131328 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1178 131328 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid phys2geo 3D', rc)
1179 : end if
1180 131328 : end subroutine edyn_esmf_regrid_phys2geo
1181 : !-----------------------------------------------------------------------
1182 36480 : subroutine edyn_esmf_regrid_geo2phys(srcfield, dstfield, ndim)
1183 : !
1184 : ! Args:
1185 : integer :: ndim
1186 : type(ESMF_Field), intent(inout) :: srcfield, dstfield
1187 : !
1188 : ! Local:
1189 : integer :: rc
1190 : character(len=*), parameter :: subname = 'edyn_esmf_regrid_geo2phys'
1191 : !
1192 36480 : if (ndim == 2) then
1193 : ! call ESMF_FieldRegrid(srcfield, dstfield, routehandle_geo2phys_2d, &
1194 : ! termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1195 : ! call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid geo2phys 2D', rc)
1196 0 : call endrun(subname//': routehandle_geo2phys_2d not implemented')
1197 : else
1198 : call ESMF_FieldRegrid( srcfield, dstfield, routehandle_geo2phys, &
1199 36480 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1200 36480 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid geo2phys 3D', rc)
1201 : end if
1202 36480 : end subroutine edyn_esmf_regrid_geo2phys
1203 : !-----------------------------------------------------------------------
1204 51072 : subroutine edyn_esmf_regrid_geo2mag(srcfield, dstfield, ndim)
1205 : !
1206 : ! Args:
1207 : integer :: ndim
1208 : type(ESMF_Field), intent(inout) :: srcfield, dstfield
1209 : !
1210 : ! Local:
1211 : integer :: rc
1212 : character(len=*), parameter :: subname = 'edyn_esmf_regrid_geo2mag'
1213 : !
1214 51072 : if (ndim == 2) then
1215 : !
1216 : ! Do sparse matrix multiply for 2d geo2mag.
1217 : !
1218 : call ESMF_FieldRegrid(srcfield, dstfield, routehandle_geo2mag_2d, &
1219 36480 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1220 36480 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid geo2mag 2D', rc)
1221 : else ! 3d geo2mag
1222 : !
1223 : ! Do sparse matrix multiply for 3d geo2mag.
1224 : !
1225 : call ESMF_FieldRegrid(srcfield, dstfield, routehandle_geo2mag, &
1226 14592 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1227 14592 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid geo2mag 3D', rc)
1228 : end if
1229 51072 : end subroutine edyn_esmf_regrid_geo2mag
1230 : !-----------------------------------------------------------------------
1231 :
1232 29184 : subroutine edyn_esmf_regrid_mag2geo(srcfield, dstfield, ndim)
1233 : !
1234 : ! Args:
1235 : integer :: ndim
1236 : type(ESMF_Field), intent(inout) :: srcfield, dstfield
1237 : !
1238 : ! Local:
1239 : integer :: rc
1240 : character(len=*), parameter :: subname = 'edyn_esmf_regrid_mag2geo'
1241 : !
1242 29184 : if (ndim == 2) then
1243 : ! call ESMF_FieldRegrid(srcfield, dstfield, routehandle_mag2geo_2d, &
1244 : ! termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1245 : ! call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid geo2mag 2D', rc)
1246 0 : call endrun(subname//': routehandle_mag2geo_2d not implemented')
1247 : else ! 3d geo2mag
1248 : !
1249 : ! Do sparse matrix multiply for 3d geo2mag.
1250 : !
1251 : call ESMF_FieldRegrid(srcfield, dstfield, routehandle_mag2geo, &
1252 29184 : termorderflag=ESMF_TERMORDER_SRCSEQ, rc=rc)
1253 29184 : call edyn_esmf_chkerr(subname, 'ESMF_FieldRegrid geo2mag 3D', rc)
1254 : end if
1255 29184 : end subroutine edyn_esmf_regrid_mag2geo
1256 :
1257 :
1258 : !-----------------------------------------------------------------------
1259 :
1260 768 : subroutine edyn_esmf_update_phys_mesh(new_phys_mesh)
1261 :
1262 : ! Dummy argument
1263 : type(ESMF_Mesh), intent(in) :: new_phys_mesh
1264 :
1265 : integer :: petcnt, i,j
1266 :
1267 : ! Ignore return code here as all we need is an attempt to reclaim memory
1268 768 : if (ESMF_MeshIsCreated(phys_mesh)) then
1269 0 : call ESMF_MeshDestroy(phys_mesh)
1270 : end if
1271 :
1272 768 : phys_mesh = new_phys_mesh
1273 :
1274 768 : if (.not. allocated(petmap)) then
1275 3072 : allocate(petmap(ntaski,ntaskj,1))
1276 : endif
1277 :
1278 768 : petcnt = 0
1279 25344 : do j = 1,ntaskj
1280 320256 : do i = 1,ntaski
1281 294912 : petmap(i,j,1) = petcnt
1282 319488 : petcnt = petcnt+1
1283 : end do
1284 : end do
1285 :
1286 768 : end subroutine edyn_esmf_update_phys_mesh
1287 :
1288 : end module edyn_esmf
|