-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmpas_vector_reconstruction.F
411 lines (336 loc) · 15.8 KB
/
mpas_vector_reconstruction.F
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
! Copyright (c) 2013, Los Alamos National Security, LLC (LANS)
! and the University Corporation for Atmospheric Research (UCAR).
!
! Unless noted otherwise source code is licensed under the BSD license.
! Additional copyright and license information can be found in the LICENSE file
! distributed with this code, or at http://mpas-dev.github.com/license.html
!
!***********************************************************************
!
! mpas_vector_reconstruction
!
!> \brief MPAS Vector reconstruction module
!> \author Xylar Asay-Davis, Todd Ringler
!> \date 03/28/13
!> \details
!> This module provides routines for performing vector reconstruction from edges to cell centers.
!
!-----------------------------------------------------------------------
module mpas_vector_reconstruction
use mpas_derived_types
use mpas_pool_routines
use mpas_constants
use mpas_rbf_interpolation
use mpas_vector_operations
implicit none
public :: mpas_init_reconstruct, mpas_reconstruct
interface mpas_reconstruct
module procedure mpas_reconstruct_1d
module procedure mpas_reconstruct_2d
end interface
contains
!***********************************************************************
!
! routine mpas_init_reconstruct
!
!> \brief MPAS Vector reconstruction initialization routine
!> \author Xylar Asay-Davis, Todd Ringler
!> \date 03/28/13
!> \details
!> Purpose: pre-compute coefficients used by the reconstruct() routine
!> Input: grid meta data
!> Output: grid % coeffs_reconstruct - coefficients used to reconstruct
!> velocity vectors at cell centers
!-----------------------------------------------------------------------
subroutine mpas_init_reconstruct(meshPool, includeHalos)!{{{
implicit none
type (mpas_pool_type), intent(in) :: &
meshPool !< Input: Mesh information
logical, optional, intent(in) :: includeHalos
! temporary arrays needed in the (to be constructed) init procedure
integer, pointer :: nCells
integer, dimension(:,:), pointer :: edgesOnCell
integer, dimension(:), pointer :: nEdgesOnCell
integer :: i, iCell, iEdge, pointCount, maxEdgeCount
real (kind=RKIND), dimension(:), pointer :: xCell, yCell, zCell, xEdge, yEdge, zEdge
real (kind=RKIND) :: r, cellCenter(3), alpha, tangentPlane(2,3)
real (kind=RKIND), allocatable, dimension(:,:) :: edgeOnCellLocations, edgeOnCellNormals, coeffs, &
edgeOnCellLocationsWork, edgeOnCellNormalsWork, coeffsWork
real(kind=RKIND), dimension(:,:), pointer :: edgeNormalVectors
real(kind=RKIND), dimension(:,:,:), pointer :: cellTangentPlane
real (kind=RKIND), dimension(:,:,:), pointer :: coeffs_reconstruct
logical, pointer :: is_periodic
real(kind=RKIND), pointer :: x_period, y_period
logical :: includeHalosLocal
call mpas_pool_get_config(meshPool, 'is_periodic', is_periodic)
call mpas_pool_get_config(meshPool, 'x_period', x_period)
call mpas_pool_get_config(meshPool, 'y_period', y_period)
if ( present(includeHalos) ) then
includeHalosLocal = includeHalos
else
includeHalosLocal = .false.
end if
!========================================================
! arrays filled and saved during init procedure
!========================================================
call mpas_pool_get_array(meshPool, 'coeffs_reconstruct', coeffs_reconstruct)
!========================================================
! temporary variables needed for init procedure
!========================================================
call mpas_pool_get_array(meshPool, 'xCell', xCell)
call mpas_pool_get_array(meshPool, 'yCell', yCell)
call mpas_pool_get_array(meshPool, 'zCell', zCell)
call mpas_pool_get_array(meshPool, 'xEdge', xEdge)
call mpas_pool_get_array(meshPool, 'yEdge', yEdge)
call mpas_pool_get_array(meshPool, 'zEdge', zEdge)
call mpas_pool_get_array(meshPool, 'nEdgesOnCell', nEdgesOnCell)
call mpas_pool_get_array(meshPool, 'edgesOnCell', edgesOnCell)
call mpas_pool_get_array(meshPool, 'edgeNormalVectors', edgeNormalVectors)
call mpas_pool_get_array(meshPool, 'cellTangentPlane', cellTangentPlane)
if ( includeHalosLocal ) then
call mpas_pool_get_dimension(meshPool, 'nCells', nCells)
else
call mpas_pool_get_dimension(meshPool, 'nCellsSolve', nCells)
end if
! init arrays
coeffs_reconstruct = 0.0
maxEdgeCount = maxval(nEdgesOnCell)
allocate(edgeOnCellLocations(maxEdgeCount,3))
allocate(edgeOnCellNormals(maxEdgeCount,3))
allocate(coeffs(maxEdgeCount,3))
! loop over all cells to be solved on this block
do iCell=1,nCells
pointCount = nEdgesOnCell(iCell)
cellCenter(1) = xCell(iCell)
cellCenter(2) = yCell(iCell)
cellCenter(3) = zCell(iCell)
do i=1,pointCount
iEdge = edgesOnCell(i,iCell)
if (is_periodic) then
edgeOnCellLocations(i,1) = mpas_fix_periodicity(xEdge(iEdge), cellCenter(1), x_period)
edgeOnCellLocations(i,2) = mpas_fix_periodicity(yEdge(iEdge), cellCenter(2), y_period)
edgeOnCellLocations(i,3) = zEdge(iEdge)
else
edgeOnCellLocations(i,1) = xEdge(iEdge)
edgeOnCellLocations(i,2) = yEdge(iEdge)
edgeOnCellLocations(i,3) = zEdge(iEdge)
end if
edgeOnCellNormals(i,:) = edgeNormalVectors(:, iEdge)
end do
alpha = 0.0
do i=1,pointCount
r = sqrt(sum((cellCenter - edgeOnCellLocations(i,:))**2))
alpha = alpha + r
enddo
alpha = alpha/pointCount
tangentPlane(1,:) = cellTangentPlane(:,1,iCell)
tangentPlane(2,:) = cellTangentPlane(:,2,iCell)
allocate(edgeOnCellLocationsWork(pointCount,3))
allocate(edgeOnCellNormalsWork(pointCount,3))
allocate(coeffsWork(pointCount,3))
edgeOnCellLocationsWork = edgeOnCellLocations(1:pointCount,:)
edgeOnCellNormalsWork = edgeOnCellNormals(1:pointCount,:)
call mpas_rbf_interp_func_3D_plane_vec_const_dir_comp_coeffs(pointCount, &
edgeOnCellLocationsWork, edgeOnCellNormalsWork, &
cellCenter, alpha, tangentPlane, coeffsWork)
coeffs(1:pointCount,:) = coeffsWork
deallocate(edgeOnCellLocationsWork)
deallocate(edgeOnCellNormalsWork)
deallocate(coeffsWork)
do i=1,pointCount
coeffs_reconstruct(:,i,iCell) = coeffs(i,:)
end do
enddo ! iCell
deallocate(edgeOnCellLocations)
deallocate(edgeOnCellNormals)
deallocate(coeffs)
end subroutine mpas_init_reconstruct!}}}
!***********************************************************************
!
! routine mpas_reconstruct_2d
!
!> \brief 2d MPAS Vector reconstruction routine
!> \author Xylar Asay-Davis, Todd Ringler
!> \date 03/28/13
!> \details
!> Purpose: reconstruct vector field at cell centers based on radial basis functions
!> Input: grid meta data and vector component data residing at cell edges
!> Output: reconstructed vector field (measured in X,Y,Z) located at cell centers
!-----------------------------------------------------------------------
subroutine mpas_reconstruct_2d(meshPool, u, uReconstructX, uReconstructY, uReconstructZ, uReconstructZonal, uReconstructMeridional, includeHalos)!{{{
implicit none
type (mpas_pool_type), intent(in) :: meshPool !< Input: Mesh information
real (kind=RKIND), dimension(:,:), intent(in) :: u !< Input: Velocity field on edges
real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructX !< Output: X Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructY !< Output: Y Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructZ !< Output: Z Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructZonal !< Output: Zonal Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructMeridional !< Output: Meridional Component of velocity reconstructed to cell centers
logical, optional, intent(in) :: includeHalos !< Input: Optional logical that allows reconstruction over halo regions
! temporary arrays needed in the compute procedure
logical :: includeHalosLocal
integer, pointer :: nCells
integer, dimension(:,:), pointer :: edgesOnCell
integer, dimension(:), pointer :: nEdgesOnCell
integer :: iCell,iEdge, i
real(kind=RKIND), dimension(:), pointer :: latCell, lonCell
real (kind=RKIND), dimension(:,:,:), pointer :: coeffs_reconstruct
logical, pointer :: on_a_sphere
real (kind=RKIND) :: clat, slat, clon, slon
if ( present(includeHalos) ) then
includeHalosLocal = includeHalos
else
includeHalosLocal = .false.
end if
! stored arrays used during compute procedure
call mpas_pool_get_array(meshPool, 'coeffs_reconstruct', coeffs_reconstruct)
! temporary variables
call mpas_pool_get_array(meshPool, 'nEdgesOnCell', nEdgesOnCell)
call mpas_pool_get_array(meshPool, 'edgesOnCell', edgesOnCell)
if ( includeHalosLocal ) then
call mpas_pool_get_dimension(meshPool, 'nCells', nCells)
else
call mpas_pool_get_dimension(meshPool, 'nCellsSolve', nCells)
end if
call mpas_pool_get_array(meshPool, 'latCell', latCell)
call mpas_pool_get_array(meshPool, 'lonCell', lonCell)
call mpas_pool_get_config(meshPool, 'on_a_sphere', on_a_sphere)
! loop over cell centers
!$omp do schedule(runtime)
do iCell = 1, nCells
! initialize the reconstructed vectors
uReconstructX(:,iCell) = 0.0
uReconstructY(:,iCell) = 0.0
uReconstructZ(:,iCell) = 0.0
! a more efficient reconstruction where rbf_values*matrix_reconstruct
! has been precomputed in coeffs_reconstruct
do i=1,nEdgesOnCell(iCell)
iEdge = edgesOnCell(i,iCell)
uReconstructX(:,iCell) = uReconstructX(:,iCell) &
+ coeffs_reconstruct(1,i,iCell) * u(:,iEdge)
uReconstructY(:,iCell) = uReconstructY(:,iCell) &
+ coeffs_reconstruct(2,i,iCell) * u(:,iEdge)
uReconstructZ(:,iCell) = uReconstructZ(:,iCell) &
+ coeffs_reconstruct(3,i,iCell) * u(:,iEdge)
enddo
enddo ! iCell
!$omp end do
call mpas_threading_barrier()
if (on_a_sphere) then
!$omp do schedule(runtime)
do iCell = 1, nCells
clat = cos(latCell(iCell))
slat = sin(latCell(iCell))
clon = cos(lonCell(iCell))
slon = sin(lonCell(iCell))
uReconstructZonal(:,iCell) = -uReconstructX(:,iCell)*slon + &
uReconstructY(:,iCell)*clon
uReconstructMeridional(:,iCell) = -(uReconstructX(:,iCell)*clon &
+ uReconstructY(:,iCell)*slon)*slat &
+ uReconstructZ(:,iCell)*clat
end do
!$omp end do
else
!$omp do schedule(runtime)
do iCell = 1, nCells
uReconstructZonal (:,iCell) = uReconstructX(:,iCell)
uReconstructMeridional(:,iCell) = uReconstructY(:,iCell)
end do
!$omp end do
end if
end subroutine mpas_reconstruct_2d!}}}
!***********************************************************************
!
! routine mpas_reconstruct_1d
!
!> \brief 1d MPAS Vector reconstruction routine
!> \author Xylar Asay-Davis, Todd Ringler, Matt Hoffman
!> \date 03/28/13
!> \details
!> Purpose: reconstruct vector field at cell centers based on radial basis functions
!> Input: grid meta data and vector component data residing at cell edges
!> Output: reconstructed vector field (measured in X,Y,Z) located at cell centers
!-----------------------------------------------------------------------
subroutine mpas_reconstruct_1d(meshPool, u, uReconstructX, uReconstructY, uReconstructZ, uReconstructZonal, uReconstructMeridional, includeHalos)!{{{
implicit none
type (mpas_pool_type), intent(in) :: meshPool !< Input: Mesh information
real (kind=RKIND), dimension(:), intent(in) :: u !< Input: Velocity field on edges
real (kind=RKIND), dimension(:), intent(out) :: uReconstructX !< Output: X Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:), intent(out) :: uReconstructY !< Output: Y Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:), intent(out) :: uReconstructZ !< Output: Z Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:), intent(out) :: uReconstructZonal !< Output: Zonal Component of velocity reconstructed to cell centers
real (kind=RKIND), dimension(:), intent(out) :: uReconstructMeridional !< Output: Meridional Component of velocity reconstructed to cell centers
logical, optional, intent(in) :: includeHalos !< Input: Logical flag that allows reconstructing over halo regions
! temporary arrays needed in the compute procedure
integer, pointer :: nCells
integer, dimension(:,:), pointer :: edgesOnCell
integer, dimension(:), pointer :: nEdgesOnCell
integer :: iCell,iEdge, i
real(kind=RKIND), dimension(:), pointer :: latCell, lonCell
real (kind=RKIND), dimension(:,:,:), pointer :: coeffs_reconstruct
logical, pointer :: on_a_sphere
logical :: includeHalosLocal
real (kind=RKIND) :: clat, slat, clon, slon
if ( present(includeHalos) ) then
includeHalosLocal = includeHalos
else
includeHalosLocal = .false.
end if
! stored arrays used during compute procedure
call mpas_pool_get_array(meshPool, 'coeffs_reconstruct', coeffs_reconstruct)
! temporary variables
call mpas_pool_get_array(meshPool, 'nEdgesOnCell', nEdgesOnCell)
call mpas_pool_get_array(meshPool, 'edgesOnCell', edgesOnCell)
if ( includeHalosLocal ) then
call mpas_pool_get_dimension(meshPool, 'nCells', nCells)
else
call mpas_pool_get_dimension(meshPool, 'nCellsSolve', nCells)
end if
call mpas_pool_get_array(meshPool, 'latCell', latCell)
call mpas_pool_get_array(meshPool, 'lonCell', lonCell)
call mpas_pool_get_config(meshPool, 'on_a_sphere', on_a_sphere)
! loop over cell centers
!$omp do schedule(runtime)
do iCell = 1, nCells
! initialize the reconstructed vectors
uReconstructX(iCell) = 0.0
uReconstructY(iCell) = 0.0
uReconstructZ(iCell) = 0.0
! a more efficient reconstruction where rbf_values*matrix_reconstruct
! has been precomputed in coeffs_reconstruct
do i=1,nEdgesOnCell(iCell)
iEdge = edgesOnCell(i,iCell)
uReconstructX(iCell) = uReconstructX(iCell) &
+ coeffs_reconstruct(1,i,iCell) * u(iEdge)
uReconstructY(iCell) = uReconstructY(iCell) &
+ coeffs_reconstruct(2,i,iCell) * u(iEdge)
uReconstructZ(iCell) = uReconstructZ(iCell) &
+ coeffs_reconstruct(3,i,iCell) * u(iEdge)
enddo
enddo ! iCell
!$omp end do
call mpas_threading_barrier()
if (on_a_sphere) then
!$omp do schedule(runtime)
do iCell = 1, nCells
clat = cos(latCell(iCell))
slat = sin(latCell(iCell))
clon = cos(lonCell(iCell))
slon = sin(lonCell(iCell))
uReconstructZonal(iCell) = -uReconstructX(iCell)*slon + &
uReconstructY(iCell)*clon
uReconstructMeridional(iCell) = -(uReconstructX(iCell)*clon &
+ uReconstructY(iCell)*slon)*slat &
+ uReconstructZ(iCell)*clat
end do
!$omp end do
else
!$omp do schedule(runtime)
do iCell = 1, nCells
uReconstructZonal (iCell) = uReconstructX(iCell)
uReconstructMeridional(iCell) = uReconstructY(iCell)
end do
!$omp end do
end if
end subroutine mpas_reconstruct_1d!}}}
end module mpas_vector_reconstruction