FE-Project
Loading...
Searching...
No Matches
scale_cubedsphere_coord_cnv.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> Module common / Coordinate conversion with cubed-sphere projection
3!!
4!! @par Description
5!! A module to provide coordinate conversions with an equiangular gnomonic cubed-sphere projection
6!!
7!!
8!! @par Reference
9!! - Nair et al. 2015:
10!! A Discontinuous Galerkin Transport Scheme on the Cubed Sphere.
11!! Monthly Weather Review, 133, 814–828.
12!! (Appendix A)
13!! - Yin et al. 2017:
14!! Parallel numerical simulation of the thermal convection in the Earth’s outer core on the cubed-sphere.
15!! Geophysical Journal International, 209, 1934–1954.
16!! (Appendix A)
17!!
18!! @author Yuta Kawai, Team SCALE
19!!
20#include "scaleFElib.h"
22 !-----------------------------------------------------------------------------
23 !
24 !++ used modules
25 !
26 use scale_const, only: &
27 pi => const_pi, &
28 eps => const_eps, &
29 rplanet => const_radius
30 use scale_precision
31 use scale_prc
32 use scale_io
33 use scale_prc
34
35 !-----------------------------------------------------------------------------
36 implicit none
37 private
38
39 !-----------------------------------------------------------------------------
40 !
41 !++ Public type & procedure
42 !
50
52 module procedure :: cubedspherecoordcnv_cs2localorthvec_alpha_0
53 module procedure :: cubedspherecoordcnv_cs2localorthvec_alpha_1
54 end interface
56
58 module procedure :: cubedspherecoordcnv_localorth2csvec_alpha_0
59 module procedure :: cubedspherecoordcnv_localorth2csvec_alpha_1
60 end interface
62
63 !-----------------------------------------------------------------------------
64 !
65 !++ Public parameters & variables
66 !
67 !-----------------------------------------------------------------------------
68
69 !
70 !++ Private type & procedure
71 !
72 !-----------------------------------------------------------------------------
73
74 !-----------------------------------------------------------------------------
75 !
76 !++ Private parameters & variables
77 !
78 !-----------------------------------------------------------------------------
79
80
81contains
82!> Calculate longitude and latitude coordinates from local coordinates using the central angles in an equiangular gnomonic cubed-sphere projection
83!!
84!OCL SERIAL
86 panelID, alpha, beta, gam, Np, & ! (in)
87 lon, lat ) ! (out)
90 implicit none
91
92 integer, intent(in) :: panelid !< Panel ID of cubed-sphere coordinates
93 integer, intent(in) :: np !< Array size
94 real(rp), intent(in) :: alpha(np) !< Local coordinate using the central angles [rad]
95 real(rp), intent(in) :: beta (np) !< Local coordinate using the central angles [rad]
96 real(rp), intent(in) :: gam(np) !< A factor of r/a
97 real(rp), intent(out) :: lon(np) !< Longitude coordinate [rad]
98 real(rp), intent(out) :: lat(np) !< Latitude coordinate [rad]
99
100 real(rp) :: cartpos(np,3)
101 real(rp) :: geogpos(np,3)
102
103 integer :: p
104 !-----------------------------------------------------------------------------
105
106 select case( panelid )
107 case( 1, 2, 3, 4 )
108 !$omp parallel
109 !$omp do
110 do p=1, np
111 lon(p) = alpha(p) + 0.5_rp * pi * dble(panelid - 1)
112 lat(p) = atan( tan( beta(p) ) * cos( alpha(p) ) )
113 end do
114 !$omp workshare
115 where( lon(:) < 0.0_rp )
116 lon(:) = lon(:) + 2.0_rp * pi
117 end where
118 !$omp end workshare
119 !$omp end parallel
120 case(5, 6)
121 call cubedspherecoordcnv_cs2cartpos( panelid, alpha, beta, gam, np, & ! (in)
122 cartpos(:,1), cartpos(:,2), cartpos(:,3) ) ! (out)
123
124 call geographiccoordcnv_orth_to_geo_pos( cartpos(:,:), np, & ! (in)
125 geogpos(:,:) ) ! (out)
126
127 !$omp parallel
128 !$omp do
129 do p=1, np
130 lon(p) = geogpos(p,1)
131 lat(p) = geogpos(p,2)
132 end do
133 !$omp workshare
134 where( alpha < 0.0_rp )
135 lon(:) = lon(:) + 2.0_rp * pi
136 end where
137 !$omp end workshare
138 !$omp end parallel
139 case default
140 log_error("CubedSphereCoordCnv_CS2LonLatPos",'(a,i2,a)') "panelID ", panelid, " is invalid. Check!"
141 call prc_abort
142 end select
143
144 return
146
147!> Covert the components of a vector in local coordinates with an equiangular gnomonic cubed-sphere projection to those in longitude and latitude coordinates
148!!
149!OCL SERIAL
151 panelID, alpha, beta, gam, Np, & ! (in)
152 vecalpha, vecbeta, & ! (in)
153 veclon, veclat, & ! (out)
154 lat ) ! (in, optional)
155
156 use scale_const, only: &
157 eps => const_eps
158 implicit none
159
160 integer, intent(in) :: panelid !< Panel ID of cubed-sphere coordinates
161 integer, intent(in) :: np !< Array size
162 real(rp), intent(in) :: alpha(np) !< Local coordinate using the central angles [rad]
163 real(rp), intent(in) :: beta (np) !< Local coordinate using the central angles [rad]
164 real(rp), intent(in) :: gam(np) !< A factor of RPlanet / r
165 real(dp), intent(in) :: vecalpha(np) !< A component of vector in the alpha-coordinate
166 real(dp), intent(in) :: vecbeta (np) !< A component of vector in the beta-coordinate
167 real(rp), intent(out) :: veclon(np) !< A component of vector in the longitude-coordinate
168 real(rp), intent(out) :: veclat(np) !< A component of vector in the latitude-coordinate
169 real(rp), intent(in), optional :: lat(np) !< latitude [rad]
170
171 integer :: p
172 real(rp) :: x ,y, del2
173 real(rp) :: s
174
175 real(rp) :: radius
176 real(rp) :: cos_lat(np)
177 !-----------------------------------------------------------------------------
178
179 if (present(lat)) then
180 !$omp parallel do
181 do p=1, np
182 cos_lat(p) = cos(lat(p))
183 end do
184 end if
185
186 select case( panelid )
187 case( 1, 2, 3, 4 )
188 if (.not. present(lat)) then
189 !$omp parallel do
190 do p=1, np
191 cos_lat(p) = cos( atan( tan( beta(p) ) * cos( alpha(p) ) ) )
192 end do
193 end if
194
195 !$omp parallel do private( X, Y, del2, radius )
196 do p=1, np
197 x = tan( alpha(p) )
198 y = tan( beta(p) )
199 del2 = 1.0_rp + x**2 + y**2
200 radius = rplanet * gam(p)
201
202 veclon(p) = vecalpha(p) * cos_lat(p) * radius
203 veclat(p) = ( - x * y * vecalpha(p) + ( 1.0_rp + y**2 ) * vecbeta(p) ) &
204 * radius * sqrt( 1.0_rp + x**2 ) / del2
205 end do
206 case ( 5 )
207 s = 1.0_rp
208 case ( 6 )
209 s = - 1.0_rp
210 case default
211 log_error("CubedSphereCoordCnv_CS2LonLatVec",'(a,i2,a)') "panelID ", panelid, " is invalid. Check!"
212 call prc_abort
213 end select
214
215 select case( panelid )
216 case( 5, 6 )
217 !$omp parallel do private( X, Y, del2, radius )
218 do p=1, np
219 x = tan( alpha(p) )
220 y = tan( beta(p) )
221 del2 = 1.0_rp + x**2 + y**2
222 radius = s * rplanet * gam(p)
223
224 if (.not. present(lat)) then
225 cos_lat(p) = cos( atan( sign(1.0_rp, s) / max( sqrt( x**2 + y**2 ), eps ) ) )
226 end if
227
228 veclon(p) = (- y * ( 1.0_rp + x**2 ) * vecalpha(p) + x * ( 1.0_rp + y**2 ) * vecbeta(p) ) &
229 * radius / max( x**2 + y**2, eps ) * cos_lat(p)
230 veclat(p) = (- x * ( 1.0_rp + x**2 ) * vecalpha(p) - y * ( 1.0_rp + y**2 ) * vecbeta(p) ) &
231 * radius / ( del2 * ( max( sqrt( x**2 + y**2 ), eps ) ) )
232 end do
233 end select
234
235 return
237
238!> Calculate local coordinates using the central angles in an equiangular gnomonic cubed-sphere projection from longitude and latitude coordinates
239!!
240!OCL SERIAL
242 panelID, lon, lat, Np, & ! (in)
243 alpha, beta ) ! (out)
244
245 implicit none
246
247 integer, intent(in) :: panelid !< Panel ID of cubed-sphere coordinates
248 integer, intent(in) :: np !< Array size
249 real(rp), intent(in) :: lon(np) !< Longitude coordinate [rad]
250 real(rp), intent(in) :: lat(np) !< Latitude coordinate [rad]
251 real(rp), intent(out) ::alpha(np) !< Local coordinate using the central angles [rad]
252 real(rp), intent(out) ::beta (np) !< Local coordinate using the central angles [rad]
253
254 integer :: p
255 real(rp) :: tan_lat
256 real(rp) :: lon_(np)
257 !-----------------------------------------------------------------------------
258
259 select case(panelid)
260 case ( 1, 2, 3, 4 )
261 !$omp parallel
262 if ( panelid == 1 ) then
263 !$omp workshare
264 where (lon(:) >= 2.0_rp * pi - 0.25_rp * pi )
265 lon_(:) = lon(:) - 2.0_rp * pi
266 elsewhere
267 lon_(:) = lon(:)
268 end where
269 !$omp end workshare
270 else
271 !$omp workshare
272 where (lon(:) < 0.0_rp )
273 lon_(:) = lon(:) + 2.0_rp * pi
274 elsewhere
275 lon_(:) = lon(:)
276 end where
277 !$omp end workshare
278 end if
279 !$omp do
280 do p=1, np
281 alpha(p) = lon_(p) - 0.5_rp * pi * ( dble(panelid) - 1.0_rp )
282 beta(p) = atan( tan(lat(p)) / cos(alpha(p)) )
283 end do
284 !$omp end parallel
285 case ( 5 )
286 !$omp parallel private(tan_lat)
287 !$omp do
288 do p=1, np
289 tan_lat = tan(lat(p))
290 alpha(p) = + atan( sin(lon(p)) / tan_lat )
291 beta(p) = - atan( cos(lon(p)) / tan_lat )
292 end do
293 !$omp end parallel
294 case ( 6 )
295 !$omp parallel private(tan_lat)
296 !$omp do
297 do p=1, np
298 tan_lat = tan(lat(p))
299 alpha(p) = - atan( sin(lon(p)) / tan_lat )
300 beta(p) = - atan( cos(lon(p)) / tan_lat )
301 end do
302 !$omp end parallel
303 case default
304 log_error("CubedSphereCoordCnv_LonLat2CSPos",'(a,i2,a)') "panelID ", panelid, " is invalid. Check!"
305 call prc_abort
306 end select
307
308 return
310
311!> Covert the components of a vector in longitude and latitude coordinates to those in local coordinates with an equiangular gnomonic cubed-sphere projection
312!!
313!OCL SERIAL
315 panelID, alpha, beta, gam, Np, & ! (in)
316 veclon, veclat, & ! (in)
317 vecalpha, vecbeta, & ! (out)
318 lat ) ! (in, optional)
319
320 implicit none
321
322 integer, intent(in) :: panelid !< Panel ID of cubed-sphere coordinates
323 integer, intent(in) :: np !< Array size
324 real(rp), intent(in) :: alpha(np) !< Local coordinate using the central angles [rad]
325 real(rp), intent(in) :: beta (np) !< Local coordinate using the central angles [rad]
326 real(rp), intent(in) :: gam(np) !< A factor of RPlanet / r
327 real(dp), intent(in) :: veclon(np) !< A component of vector in the longitude-coordinate
328 real(dp), intent(in) :: veclat(np) !< A component of vector in the latitude-coordinate
329 real(dp), intent(out) :: vecalpha(np) !< A component of vector in the alpha-coordinate
330 real(dp), intent(out) :: vecbeta (np) !< A component of vector in the beta-coordinate
331 real(rp), intent(in), optional :: lat(np) !< latitude [rad]
332
333 integer :: p
334 real(rp) :: x ,y, del2
335 real(rp) :: s
336
337 real(rp) :: radius
338 real(rp) :: cos_lat(np)
339 real(rp) :: veclon_ov_coslat
340 !-----------------------------------------------------------------------------
341
342 if (present(lat)) then
343 !$omp parallel do
344 do p=1, np
345 cos_lat(p) = cos(lat(p))
346 end do
347 end if
348
349 select case( panelid )
350 case( 1, 2, 3, 4 )
351 if (.not. present(lat)) then
352 !$omp parallel do
353 do p=1, np
354 cos_lat(p) = cos( atan( tan( beta(p) ) * cos( alpha(p) ) ) )
355 end do
356 end if
357
358 !$omp parallel do private( X, Y, del2, radius, VecLon_ov_cosLat )
359 do p=1, np
360 x = tan( alpha(p) )
361 y = tan( beta(p) )
362 del2 = 1.0_rp + x**2 + y**2
363 radius = rplanet * gam(p)
364 veclon_ov_coslat = veclon(p) / cos_lat(p)
365
366 vecalpha(p) = veclon_ov_coslat / radius
367 vecbeta(p) = ( x * y * veclon_ov_coslat + del2 / sqrt( 1.0_rp + x**2 ) * veclat(p) ) &
368 / ( radius * (1.0_rp + y**2) )
369 end do
370 case ( 5 )
371 s = 1.0_rp
372 case ( 6 )
373 s = -1.0_rp
374 case default
375 log_error("CubedSphereCoordCnv_LonLat2CSVec",'(a,i2,a)') "panelID ", panelid, " is invalid. Check!"
376 call prc_abort
377 end select
378
379 select case( panelid )
380 case( 5, 6 )
381 !$omp parallel do private( X, Y, del2, radius, VecLon_ov_cosLat )
382 do p=1, np
383 x = tan( alpha(p) )
384 y = tan( beta(p) )
385 del2 = 1.0_rp + x**2 + y**2
386 radius = s * rplanet * gam(p)
387
388 if (.not. present(lat)) then
389 cos_lat(p) = cos( atan( s / sqrt(max(x**2 + y**2, eps)) ) )
390 end if
391 veclon_ov_coslat = veclon(p) / cos_lat(p)
392
393 vecalpha(p) = (- y * veclon_ov_coslat - del2 * x / sqrt( max(del2 - 1.0_rp,eps)) * veclat(p)) &
394 / ( radius * ( 1.0_rp + x**2 ) )
395 vecbeta(p) = ( x * veclon_ov_coslat - del2 * y / sqrt( max(del2 - 1.0_rp,eps)) * veclat(p)) &
396 / ( radius * ( 1.0_rp + y**2 ) )
397 end do
398 end select
399
400 return
402
403!> Calculate Cartesian coordinates from local coordinates using the central angles in an equiangular gnomonic cubed-sphere projection
404!!
405!OCL SERIAL
407 panelID, alpha, beta, gam, Np, & ! (in)
408 x, y, z ) ! (out)
409
410 implicit none
411 integer, intent(in) :: panelid !< Panel ID of cubed-sphere coordinates
412 integer, intent(in) :: np !< Array size
413 real(rp), intent(in) :: alpha(np) !< Local coordinate using the central angles [rad]
414 real(rp), intent(in) :: beta (np) !< Local coordinate using the central angles [rad]
415 real(rp), intent(in) :: gam(np) !< A factor of r/a
416 real(rp), intent(out) :: x(np) !< x-coordinate in the Cartesian coordinate
417 real(rp), intent(out) :: y(np) !< y-coordinate in the Cartesian coordinate
418 real(rp), intent(out) :: z(np) !< z-coordinate in the Cartesian coordinate
419
420 integer :: p
421 real(rp) :: x1, x2, fac
422
423 !-----------------------------------------------------------------------------
424
425 select case(panelid)
426 case(1)
427 !$omp parallel do private(x1, x2, fac)
428 do p=1, np
429 x1 = tan( alpha(p) )
430 x2 = tan( beta(p) )
431 fac = rplanet * gam(p) / sqrt( 1.0_rp + x1**2 + x2**2 )
432 x(p) = fac
433 y(p) = fac * x1
434 z(p) = fac * x2
435 end do
436 case(2)
437 !$omp parallel do private(x1, x2, fac)
438 do p=1, np
439 x1 = tan( alpha(p) )
440 x2 = tan( beta(p) )
441 fac = rplanet * gam(p) / sqrt( 1.0_rp + x1**2 + x2**2 )
442 x(p) = - fac * x1
443 y(p) = fac
444 z(p) = fac * x2
445 end do
446 case(3)
447 !$omp parallel do private(x1, x2, fac)
448 do p=1, np
449 x1 = tan( alpha(p) )
450 x2 = tan( beta(p) )
451 fac = rplanet * gam(p) / sqrt( 1.0_rp + x1**2 + x2**2 )
452 x(p) = - fac
453 y(p) = - fac * x1
454 z(p) = fac * x2
455 end do
456 case(4)
457 !$omp parallel do private(x1, x2, fac)
458 do p=1, np
459 x1 = tan( alpha(p) )
460 x2 = tan( beta(p) )
461 fac = rplanet * gam(p) / sqrt( 1.0_rp + x1**2 + x2**2 )
462 x(p) = fac * x1
463 y(p) = - fac
464 z(p) = fac * x2
465 end do
466 case(5)
467 !$omp parallel do private(x1, x2, fac)
468 do p=1, np
469 x1 = tan( alpha(p) )
470 x2 = tan( beta(p) )
471 fac = rplanet * gam(p) / sqrt( 1.0_rp + x1**2 + x2**2 )
472 x(p) = - fac * x2
473 y(p) = fac * x1
474 z(p) = fac
475 end do
476 case(6)
477 !$omp parallel do private(x1, x2, fac)
478 do p=1, np
479 x1 = tan( alpha(p) )
480 x2 = tan( beta(p) )
481 fac = rplanet * gam(p) / sqrt( 1.0_rp + x1**2 + x2**2 )
482 x(p) = fac * x2
483 y(p) = fac * x1
484 z(p) = - fac
485 end do
486 case default
487 log_error("CubedSphereCoordCnv_CS2CartPos",'(a,i2,a)') "panelID ", panelid, " is invalid. Check!"
488 call prc_abort
489 end select
490
491 return
492 end subroutine cubedspherecoordcnv_cs2cartpos
493
494!> Covert the components of a vector in local coordinates with an equiangular gnomonic cubed-sphere projection to those in the Cartesian coordinates
495!!
496!OCL SERIAL
498 panelID, alpha, beta, gam, Np, & ! (in)
499 vec_x, vec_y, vec_z, & ! (in)
500 vecalpha, vecbeta ) ! (out)
501
502 implicit none
503
504 integer, intent(in) :: panelid !< Panel ID of cubed-sphere coordinates
505 integer, intent(in) :: np !< Array size
506 real(rp), intent(in) :: alpha(np) !< Local coordinate using the central angles [rad]
507 real(rp), intent(in) :: beta (np) !< Local coordinate using the central angles [rad]
508 real(rp), intent(in) :: gam(np) !< A factor of RPlanet / r
509 real(dp), intent(in) :: vec_x(np) !< A component of vector in the x-coordinate with the Cartesian coordinate
510 real(dp), intent(in) :: vec_y(np) !< A component of vector in the y-coordinate with the Cartesian coordinate
511 real(dp), intent(in) :: vec_z(np) !< A component of vector in the z-coordinate with the Cartesian coordinate
512 real(rp), intent(out) :: vecalpha(np) !< A component of vector in the alpha-coordinate
513 real(rp), intent(out) :: vecbeta (np) !< A component of vector in the beta-coordinate
514
515 integer :: p
516 real(rp) :: x1, x2, fac
517 real(rp) :: r_sec2_alpha, r_sec2_beta
518 !-----------------------------------------------------------------------------
519
520 select case( panelid )
521 case(1)
522 !$omp parallel do private(x1, x2, r_sec2_alpha, r_sec2_beta, fac)
523 do p=1, np
524 x1 = tan( alpha(p) )
525 x2 = tan( beta(p) )
526 r_sec2_alpha = cos(alpha(p))**2
527 r_sec2_beta = cos(beta(p))**2
528 fac = sqrt( 1.0_rp + x1**2 + x2**2 ) / ( rplanet * gam(p) )
529
530 vecalpha(p) = fac * r_sec2_alpha * ( - x1 * vec_x(p) + vec_y(p) )
531 vecbeta(p) = fac * r_sec2_beta * ( - x2 * vec_x(p) + vec_z(p) )
532 end do
533 case(2)
534 !$omp parallel do private(x1, x2, r_sec2_alpha, r_sec2_beta, fac)
535 do p=1, np
536 x1 = tan( alpha(p) )
537 x2 = tan( beta(p) )
538 r_sec2_alpha = cos(alpha(p))**2
539 r_sec2_beta = cos(beta(p))**2
540 fac = sqrt( 1.0_rp + x1**2 + x2**2 ) / ( rplanet * gam(p) )
541
542 vecalpha(p) = - fac * r_sec2_alpha * ( vec_x(p) + x1 * vec_y(p) )
543 vecbeta(p) = fac * r_sec2_beta * ( - x2 * vec_y(p) + vec_z(p) )
544 end do
545 case(3)
546 !$omp parallel do private(x1, x2, r_sec2_alpha, r_sec2_beta, fac)
547 do p=1, np
548 x1 = tan( alpha(p) )
549 x2 = tan( beta(p) )
550 r_sec2_alpha = cos(alpha(p))**2
551 r_sec2_beta = cos(beta(p))**2
552 fac = sqrt( 1.0_rp + x1**2 + x2**2 ) / ( rplanet * gam(p) )
553
554 vecalpha(p) = fac * r_sec2_alpha * ( x1 * vec_x(p) - vec_y(p) )
555 vecbeta(p) = fac * r_sec2_beta * ( x2 * vec_x(p) + vec_z(p) )
556 end do
557 case(4)
558 !$omp parallel do private(x1, x2, r_sec2_alpha, r_sec2_beta, fac)
559 do p=1, np
560 x1 = tan( alpha(p) )
561 x2 = tan( beta(p) )
562 r_sec2_alpha = cos(alpha(p))**2
563 r_sec2_beta = cos(beta(p))**2
564 fac = sqrt( 1.0_rp + x1**2 + x2**2 ) / ( rplanet * gam(p) )
565
566 vecalpha(p) = fac * r_sec2_alpha * ( vec_x(p) + x1 * vec_y(p) )
567 vecbeta(p) = fac * r_sec2_beta * ( x2 * vec_y(p) + vec_z(p) )
568 end do
569 case ( 5 )
570 !$omp parallel do private(x1, x2, r_sec2_alpha, r_sec2_beta, fac)
571 do p=1, np
572 x1 = tan( alpha(p) )
573 x2 = tan( beta(p) )
574 r_sec2_alpha = cos(alpha(p))**2
575 r_sec2_beta = cos(beta(p))**2
576 fac = sqrt( 1.0_rp + x1**2 + x2**2 ) / ( rplanet * gam(p) )
577
578 vecalpha(p) = fac * r_sec2_alpha * ( vec_y(p) - x1 * vec_z(p) )
579 vecbeta(p) = - fac * r_sec2_beta * ( vec_x(p) + x2 * vec_z(p) )
580 end do
581 case ( 6 )
582 !$omp parallel do private(x1, x2, r_sec2_alpha, r_sec2_beta, fac)
583 do p=1, np
584 x1 = tan( alpha(p) )
585 x2 = tan( beta(p) )
586 r_sec2_alpha = cos(alpha(p))**2
587 r_sec2_beta = cos(beta(p))**2
588 fac = sqrt( 1.0_rp + x1**2 + x2**2 ) / ( rplanet * gam(p) )
589
590 vecalpha(p) = fac * r_sec2_alpha * ( vec_y(p) + x1 * vec_z(p) )
591 vecbeta(p) = fac * r_sec2_beta * ( vec_x(p) + x2 * vec_z(p) )
592 end do
593 case default
594 log_error("CubedSphereCoordCnv_Cart2CSVec",'(a,i2,a)') "panelID ", panelid, " is invalid. Check!"
595 call prc_abort
596 end select
597
598 return
599 end subroutine cubedspherecoordcnv_cart2csvec
600
601!OCL SERIAL
602 subroutine cubedspherecoordcnv_cs2localorthvec_alpha_0( &
603 alpha, beta, radius, Np, & ! (in)
604 vecorth1, vecorth2 ) ! (inout)
605
606 implicit none
607
608 integer, intent(in) :: Np
609 real(RP), intent(in) :: alpha(Np)
610 real(RP), intent(in) :: beta (Np)
611 real(RP), intent(in) :: radius(Np)
612 real(RP), intent(inout) :: VecOrth1(Np)
613 real(RP), intent(inout) :: VecOrth2(Np)
614
615 integer :: p
616 real(RP) :: x1, x2, del
617 real(RP) :: fac
618 real(RP) :: tmp
619 !-----------------------------------------------------------------------------
620
621 !$omp parallel do private(x1, x2, del, fac, tmp)
622 do p=1, np
623 x1 = tan( alpha(p) )
624 x2 = tan( beta(p) )
625 del = sqrt(1.0_rp + x1**2 + x2**2)
626 fac = radius(p) * sqrt(1.0_rp + x1**2) / del**2
627
628 tmp = vecorth1(p)
629 vecorth1(p) = fac * del * tmp
630 vecorth2(p) = fac * ( - x1 * x2 * tmp + (1.0_rp + x2**2) * vecorth2(p) )
631 end do
632
633 return
634 end subroutine cubedspherecoordcnv_cs2localorthvec_alpha_0
635
636!OCL SERIAL
637 subroutine cubedspherecoordcnv_cs2localorthvec_alpha_1( &
638 alpha, beta, radius, Np, & ! (in)
639 vecalpha, vecbeta, & ! (in)
640 vecorth1, vecorth2 ) ! (out)
641
642 implicit none
643
644 integer, intent(in) :: Np
645 real(RP), intent(in) :: alpha(Np)
646 real(RP), intent(in) :: beta (Np)
647 real(RP), intent(in) :: radius(Np)
648 real(RP), intent(in) :: VecAlpha(Np)
649 real(RP), intent(in) :: VecBeta (Np)
650 real(RP), intent(out) :: VecOrth1(Np)
651 real(RP), intent(out) :: VecOrth2(Np)
652
653 !-----------------------------------------------------------------------------
654
655 !$omp parallel workshare
656 vecorth1(:) = vecalpha(:)
657 vecorth2(:) = vecbeta(:)
658 !$omp end parallel workshare
659
660 call cubedspherecoordcnv_cs2localorthvec_alpha_0( &
661 alpha, beta, radius, np, & ! (in)
662 vecorth1, vecorth2 ) ! (inout)
663
664 return
665 end subroutine cubedspherecoordcnv_cs2localorthvec_alpha_1
666
667!OCL SERIAL
668 subroutine cubedspherecoordcnv_localorth2csvec_alpha_0( &
669 alpha, beta, radius, Np, & ! (in)
670 vecalpha, vecbeta ) ! (inout)
671
672 implicit none
673
674 integer, intent(in) :: Np
675 real(RP), intent(in) :: alpha(Np)
676 real(RP), intent(in) :: beta (Np)
677 real(RP), intent(in) :: radius(Np)
678 real(RP), intent(inout) :: VecAlpha(Np)
679 real(RP), intent(inout) :: VecBeta (Np)
680
681 integer :: p
682 real(RP) :: x1, x2, del
683 real(RP) :: fac
684 real(RP) :: tmp
685 !-----------------------------------------------------------------------------
686
687 !$omp parallel do private(x1, x2, del, fac, tmp)
688 do p=1, np
689 x1 = tan( alpha(p) )
690 x2 = tan( beta(p) )
691 del = sqrt(1.0_rp + x1**2 + x2**2)
692 fac = del / ( radius(p) * ( 1.0_rp + x2**2 ) * sqrt( 1.0_rp + x1**2 ) )
693
694 tmp = vecalpha(p)
695 vecalpha(p) = fac * ( 1.0_rp + x2**2 ) * tmp
696 vecbeta(p) = fac * ( x1 * x2 * tmp + del * vecbeta(p) )
697 end do
698
699 return
700 end subroutine cubedspherecoordcnv_localorth2csvec_alpha_0
701
702!OCL SERIAL
703 subroutine cubedspherecoordcnv_localorth2csvec_alpha_1( &
704 alpha, beta, radius, Np, & ! (in)
705 vecorth1, vecorth2, & ! (in)
706 vecalpha, vecbeta ) ! (out)
707
708 implicit none
709
710 integer, intent(in) :: Np
711 real(RP), intent(in) :: alpha(Np)
712 real(RP), intent(in) :: beta (Np)
713 real(RP), intent(in) :: radius(Np)
714 real(RP), intent(in) :: VecOrth1(Np)
715 real(RP), intent(in) :: VecOrth2(Np)
716 real(RP), intent(out) :: VecAlpha(Np)
717 real(RP), intent(out) :: VecBeta (Np)
718 !-----------------------------------------------------------------------------
719
720 !$omp parallel workshare
721 vecalpha(:) = vecorth1(:)
722 vecbeta(:) = vecorth2(:)
723 !$omp end parallel workshare
724
725 call cubedspherecoordcnv_localorth2csvec_alpha_0( &
726 alpha, beta, radius, np, & ! (in)
727 vecalpha, vecbeta ) ! (inout)
728
729 return
730 end subroutine cubedspherecoordcnv_localorth2csvec_alpha_1
731
732!> Calculate the metrics associated with an equiangular gnomonic cubed-sphere projection to those in longitude and latitude coordinates
733!!
734!OCL SERIAL
736 alpha, beta, Np, radius, & ! (in)
737 g_ij, gij, gsqrt ) ! (out)
738
739 implicit none
740
741 integer, intent(in) :: np !< Array size
742 real(rp), intent(in) :: alpha(np) !< Local coordinate using the central angles [rad]
743 real(rp), intent(in) :: beta (np) !< Local coordinate using the central angles [rad]
744 real(rp), intent(in) :: radius !< Planetary radius
745 real(rp), intent(out) :: g_ij(np,2,2) !< Horizontal covariant metric tensor
746 real(rp), intent(out) :: gij (np,2,2) !< Horizontal contravariant metric tensor
747 real(rp), intent(out) :: gsqrt(np) !< Horizontal Jacobian
748
749 real(rp) :: x, y
750 real(rp) :: r2
751 real(rp) :: fac
752 real(rp) :: g_ij_(2,2)
753 real(rp) :: gij_ (2,2)
754
755 integer :: p
756 real(rp) :: oneplusx2, oneplusy2
757 !-----------------------------------------------------------------------------
758
759 !$omp parallel do private( &
760 !$omp X, Y, r2, OnePlusX2, OnePlusY2, fac, &
761 !$omp G_ij_, GIJ_ )
762 do p=1, np
763 x = tan(alpha(p))
764 y = tan(beta(p))
765 r2 = 1.0_rp + x**2 + y**2
766 oneplusx2 = 1.0_rp + x**2
767 oneplusy2 = 1.0_rp + y**2
768
769 fac = oneplusx2 * oneplusy2 * ( radius / r2 )**2
770 g_ij_(1,1) = fac * oneplusx2
771 g_ij_(1,2) = - fac * (x * y)
772 g_ij_(2,1) = - fac * (x * y)
773 g_ij_(2,2) = fac * oneplusy2
774 g_ij(p,:,:) = g_ij_(:,:)
775
776 gsqrt(p) = radius**2 * oneplusx2 * oneplusy2 / ( r2 * sqrt(r2) )
777
778 fac = 1.0_rp / gsqrt(p)**2
779 gij_(1,1) = fac * g_ij_(2,2)
780 gij_(1,2) = - fac * g_ij_(1,2)
781 gij_(2,1) = - fac * g_ij_(2,1)
782 gij_(2,2) = fac * g_ij_(1,1)
783 gij(p,:,:) = gij_(:,:)
784 end do
785
786 return
787 end subroutine cubedspherecoordcnv_getmetric
788
Module common / Coordinate conversion with cubed-sphere projection.
subroutine, public cubedspherecoordcnv_cs2lonlatvec(panelid, alpha, beta, gam, np, vecalpha, vecbeta, veclon, veclat, lat)
Covert the components of a vector in local coordinates with an equiangular gnomonic cubed-sphere proj...
subroutine, public cubedspherecoordcnv_cs2lonlatpos(panelid, alpha, beta, gam, np, lon, lat)
Calculate longitude and latitude coordinates from local coordinates using the central angles in an eq...
subroutine, public cubedspherecoordcnv_lonlat2csvec(panelid, alpha, beta, gam, np, veclon, veclat, vecalpha, vecbeta, lat)
Covert the components of a vector in longitude and latitude coordinates to those in local coordinates...
subroutine, public cubedspherecoordcnv_getmetric(alpha, beta, np, radius, g_ij, gij, gsqrt)
Calculate the metrics associated with an equiangular gnomonic cubed-sphere projection to those in lon...
subroutine, public cubedspherecoordcnv_lonlat2cspos(panelid, lon, lat, np, alpha, beta)
Calculate local coordinates using the central angles in an equiangular gnomonic cubed-sphere projecti...
subroutine, public cubedspherecoordcnv_cart2csvec(panelid, alpha, beta, gam, np, vec_x, vec_y, vec_z, vecalpha, vecbeta)
Covert the components of a vector in local coordinates with an equiangular gnomonic cubed-sphere proj...
subroutine, public cubedspherecoordcnv_cs2cartpos(panelid, alpha, beta, gam, np, x, y, z)
Calculate Cartesian coordinates from local coordinates using the central angles in an equiangular gno...
Module common / Coordinate conversion with a geographic coordinate.
subroutine, public geographiccoordcnv_orth_to_geo_pos(orth_p, np, geo_p)