FE-Project
Loading...
Searching...
No Matches
scale_meshfield_spectral_transform.F90
Go to the documentation of this file.
1!> module FElib / Data / Utility
2!!
3!! @par Description
4!! A module for providing utility routines associated with spectral transformation for DG fields
5!! @author Yuta Kawai, Team SCALE
6!!
7!<
8#include "scaleFElib.h"
10 !-----------------------------------------------------------------------------
11 !
12 !++ used modules
13 !
14 use scale_precision
15 use scale_io
16 use scale_const, only: &
17 pi => const_pi
18 use scale_prc, only: &
19 prc_abort
20
21 use scale_polynominal, only: &
23 use scale_element_base, only: &
25 use scale_element_line, only: &
32
33 use scale_meshfield_base, only: &
36
39 !-----------------------------------------------------------------------------
40 implicit none
41 private
42 !-----------------------------------------------------------------------------
43 !
44 !++ Public type & procedure
45 !
46
48 integer :: eval_type_id
49 integer :: ndim
50 integer :: var_num
51
52 real(rp), allocatable :: fftintrpmat(:,:)
53 real(rp), allocatable :: fft_xi(:)
54 integer :: nsampleptperelem
55
56 real(rp), allocatable :: intintrpmat(:,:)
57 real(rp), allocatable :: intweight(:)
58 real(rp), allocatable :: intxi(:)
59 integer :: nintglpt
61
62 !> A derived type for spectral transform of MeshField1D
64 integer :: kall
65 integer :: ks, ke
66
67 real(rp) :: xmin_gl, xmax_gl
68 real(rp) :: delx
69
70 real(rp), allocatable :: k(:)
71 real(rp), allocatable :: spectral_coef(:,:,:)
72
74 contains
75 procedure :: init => meshfield_spetraltransform1d_init
76 procedure :: final => meshfield_spetraltransform1d_final
77 procedure :: transform => meshfield_spetraltransform1d_transform
79
80 !> A derived type for spectral transform of MeshField2D
82 integer :: kall
83 integer :: ks, ke
84 integer :: lall
85 integer :: ls, le
86
87 real(rp) :: xmin_gl, xmax_gl
88 real(rp) :: delx
89 real(rp) :: ymin_gl, ymax_gl
90 real(rp) :: dely
91
92 integer :: nprcx
93 integer :: nprcy
94 integer :: negx, negy
95
96 integer :: nsampleptperelem1d
97
98 real(rp), allocatable :: k(:)
99 real(rp), allocatable :: l(:)
100 real(rp), allocatable :: spectral_coef(:,:,:,:)
101
104 contains
105 procedure :: init => meshfield_spetraltransform2d_init
106 procedure :: final => meshfield_spetraltransform2d_final
107 procedure :: transform => meshfield_spetraltransform2d_transform
109
110 !-----------------------------------------------------------------------------
111 !
112 !++ Public parameters & variables
113 !
114 integer, public, parameter :: st_evaltype_sample_uniform_pts = 1
115 integer, public, parameter :: st_evaltype_l2projection_1 = 2
116 integer, public, parameter :: st_evaltype_l2projection_2 = 3
117
118 !-----------------------------------------------------------------------------
119 !
120 !++ Private procedure
121 !
122
123 !-----------------------------------------------------------------------------
124 !
125 !++ Private parameters & variables
126 !
127 !-----------------------------------------------------------------------------
128
129contains
130!- Base type
131
132!OCL SERIAL
133 subroutine meshfield_spetraltransformbase_init( this, eval_type, ndim, var_num, NintGLpt, NsamplePtPerElem )
134 implicit none
135 class(meshfield_spetraltransformbase), intent(inout) :: this
136 integer, intent(in) :: eval_type
137 integer, intent(in) :: ndim
138 integer, intent(in) :: var_num
139 integer, intent(in), optional :: nintglpt
140 integer, intent(in), optional :: nsampleptperelem
141 !--------------------------------------------
142
143 this%eval_type_id = eval_type
144 call check_eval_type_id( eval_type )
145
146 this%ndim = ndim
147 this%var_num = var_num
148
149 if ( this%eval_type_id == st_evaltype_l2projection_1 &
150 .or. this%eval_type_id == st_evaltype_l2projection_2 ) then
151 if ( .not. present(nintglpt) ) then
152 log_info("MeshField_SpetralTransformBase_Init",*) "The order of Gaussian quadrature should be given. Check!"
153 call prc_abort
154 end if
155 this%NintGLPt = nintglpt
156 else
157 this%NintGLPt = -1
158 end if
159
160 if ( this%eval_type_id == st_evaltype_sample_uniform_pts ) then
161 if ( .not. present(nsampleptperelem) ) then
162 log_info("MeshField_SpetralTransformBase_Init",*) "NsamplePtPerElem shouled be given. Check!"
163 else
164 this%NsamplePtPerElem = nsampleptperelem
165 end if
166 else
167 this%NsamplePtPerElem = -1
168 end if
169 log_info("MeshField_SpetralTransformBase_Init",*) "NintGLPt=", this%NintGLpt, "NsamplePtPerElem=", this%NsamplePtPerElem
170
171 return
172 end subroutine meshfield_spetraltransformbase_init
173
174 subroutine check_eval_type_id( eval_type )
175 implicit none
176 integer, intent(in) :: eval_type
177 !-------------------------
178 select case(eval_type)
179! case(ST_EVALTYPE_L2PROJECTION_1,ST_EVALTYPE_L2PROJECTION_2,ST_EVALTYPE_SAMPLE_UNIFORM_PTS)
181 case default
182 log_info("MeshField_SpetralTransformBase_Init",*) "Unsupported evaluation type is specified. Check!"
183 call prc_abort
184 end select
185 return
186 end subroutine check_eval_type_id
187
188!OCL SERIAL
189 subroutine meshfield_spetraltransformbase_final( this )
190 implicit none
191 class(meshfield_spetraltransformbase), intent(inout) :: this
192 !--------------------------------------------
193
194 if ( this%NintGLpt > 0 ) &
195 deallocate( this%IntIntrpMat, this%IntXi, this%IntWeight )
196
197 if ( this%NsamplePtPerElem > 0 ) &
198 deallocate( this%FFTIntrpMat, this%FFT_xi )
199
200 return
201 end subroutine meshfield_spetraltransformbase_final
202
203!- 1D
204
205!OCL SERIAL
206 subroutine meshfield_spetraltransform1d_init( this, eval_type, ks, ke, mesh1D, var_num, &
207 GLQuadOrd, NsamplePtPerElem )
208 implicit none
209 class(meshfield_spetraltransform1d), intent(inout) :: this
210 integer, intent(in) :: eval_type
211 integer, intent(in) :: ks, ke
212 class(meshbase1d), intent(in), target :: mesh1d
213 integer, intent(in) :: var_num
214 integer, optional :: glquadord
215 integer, optional :: nsampleptperelem
216
217 integer :: i
218 integer :: lx
219 real(rp) :: dx
220
221 class(localmesh1d), pointer :: lmesh
222 class(elementbase1d), pointer :: elem
223
224 type(lineelement) :: elem_dummy
225 !--------------------------------------------
226
227 lmesh => mesh1d%lcmesh_list(1)
228 elem => lmesh%refElem1D
229
230 call meshfield_spetraltransformbase_init( this, eval_type, 1, var_num, &
231 glquadord, nsampleptperelem )
232
233 this%ks = ks; this%ke = ke
234 this%kall = ke - ks + 1
235
236 this%xmin_gl = mesh1d%xmin_gl
237 this%xmax_gl = mesh1d%xmax_gl
238 this%delx = 1.0_rp / real(mesh1d%NeG,kind=rp)
239
240 allocate( this%k(ks:ke) )
241 lx = this%xmax_gl - this%xmin_gl
242 do i=ks, ke
243 this%k(i) = i * 2.0_rp * pi / lx
244 end do
245
246 allocate( this%spectral_coef(ks:ke,2,var_num) )
247
248 call elem_dummy%Init(elem%PolyOrder, .false.)
249
250 if ( this%NintGLpt > 0 ) then
251 allocate( this%IntIntrpMat(this%NintGLpt,elem%Np) )
252 allocate( this%IntXi(this%NintGLpt), this%IntWeight(this%NintGLpt) )
253
254 this%IntIntrpMat(:,:) = elem_dummy%GenIntGaussLegendreIntrpMat( glquadord, this%IntWeight, this%IntXi )
255 end if
256
257 if ( this%NsamplePtPerElem > 0 ) then
258 allocate( this%FFTIntrpMat(this%NsamplePtPerElem,elem%Np) )
259 allocate( this%FFT_Xi(this%NsamplePtPerElem) )
260
261 dx = 2.0_rp / real(this%NsamplePtPerElem,kind=rp)
262 do i=1, this%NsamplePtPerElem
263 this%FFT_Xi(i) = -1.0_rp + real(i - 0.5_rp,kind=rp) * dx
264 end do
265 this%FFTIntrpMat(:,:) = polynominal_genlagrangepoly( elem_dummy%PolyOrder, elem_dummy%x1, this%FFT_xi )
266
267 call this%fft%Init( nsampleptperelem * mesh1d%NeG )
268 end if
269
270 call elem_dummy%Final()
271
272 return
273 end subroutine meshfield_spetraltransform1d_init
274
275!OCL SERIAL
276 subroutine meshfield_spetraltransform1d_final( this )
277 implicit none
278 class(meshfield_spetraltransform1d), intent(inout) :: this
279 !--------------------------------------------
280
281 call meshfield_spetraltransformbase_final( this )
282
283 if ( this%NsamplePtPerElem > 0 ) then
284 call this%fft%Final()
285 end if
286
287 deallocate( this%k, this%spectral_coef )
288
289 return
290 end subroutine meshfield_spetraltransform1d_final
291
292!OCL SERIAL
293 subroutine meshfield_spetraltransform1d_transform( this, q_list, mesh_num )
294 implicit none
295 class(meshfield_spetraltransform1d), intent(inout) :: this
296 integer, intent(in) :: mesh_num
297 type(meshfield1dlist), target :: q_list(this%var_num,mesh_num)
298
299 class(meshbase1d), pointer :: mesh1d
300 class(localmesh1d), pointer :: lmesh
301 !--------------------------------------------
302
303 mesh1d => q_list(1,1)%ptr%mesh
304 lmesh => mesh1d%lcmesh_list(1)
305
306 select case( this%eval_type_id )
308 call spectral_transform1d_dft( this%spectral_coef, & ! (out)
309 q_list, this%var_num, mesh_num, & ! (in)
310 this%ks, this%ke, lmesh%refElem1D%Np, lmesh%Ne, this%NsamplePtPerElem, & ! (in)
311 this%fft, this%FFTIntrpMat, this%FFT_xi, & ! (in)
312 this%xmin_gl, this%delx, this%xmax_gl - this%xmin_gl ) ! (in)
314 call spectral_transform1d_l2projection( this%spectral_coef, & ! (out)
315 q_list, this%var_num, mesh_num, & ! (in)
316 this%ks, this%ke, lmesh%refElem1D%Np, lmesh%Ne, & ! (in)
317 this%IntIntrpMat, this%IntXi, this%IntWeight, this%NintGLPt, & ! (in)
318 this%xmin_gl, this%delx, this%xmax_gl - this%xmin_gl ) ! (in)
319 end select
320
321 return
322 end subroutine meshfield_spetraltransform1d_transform
323
324!- 2D
325
326!OCL SERIAL
327 subroutine meshfield_spetraltransform2d_init( this, eval_type, ks, ke, ls, le, mesh2D, var_num, &
328 GLQuadOrd, NsamplePtPerElem1D )
329 use scale_element_line, only: &
331 use scale_mesh_rectdom2d, only: &
333 implicit none
334 class(meshfield_spetraltransform2d), intent(inout) :: this
335 integer, intent(in) :: eval_type
336 integer, intent(in) :: ks, ke
337 integer, intent(in) :: ls, le
338 class(meshbase2d), intent(in), target :: mesh2d
339 integer, intent(in) :: var_num
340 integer, optional :: glquadord
341 integer, optional :: nsampleptperelem1d
342
343 integer :: i, j
344 integer :: lx, ly
345 real(rp) :: dx
346
347 class(localmesh2d), pointer :: lmesh
348 class(elementbase2d), pointer :: elem
349
350 type(lineelement) :: elem_dummy
351 class(meshbase2d), pointer :: ptr_mesh2d
352 !--------------------------------------------
353
354 lmesh => mesh2d%lcmesh_list(1)
355 elem => lmesh%refElem2D
356
357 call meshfield_spetraltransformbase_init( this, eval_type, 1, var_num, glquadord, nsampleptperelem1d )
358
359 this%ks = ks; this%ke = ke
360 this%kall = ke - ks + 1
361
362 this%ls = ls; this%le = le
363 this%lall = le - ls + 1
364
365 select type(ptr_mesh2d => mesh2d)
366 class is (meshrectdom2d)
367 this%xmin_gl = ptr_mesh2d%xmin_gl
368 this%xmax_gl = ptr_mesh2d%xmax_gl
369 this%delx = ( this%xmax_gl - this%xmin_gl ) / real(ptr_mesh2d%NeGX, kind=rp)
370
371 this%ymin_gl = ptr_mesh2d%ymin_gl
372 this%ymax_gl = ptr_mesh2d%ymax_gl
373 this%dely = ( this%ymax_gl - this%ymin_gl ) / real(ptr_mesh2d%NeGY, kind=rp)
374
375 this%NeGX = ptr_mesh2d%NeGX
376 this%NeGY = ptr_mesh2d%NeGY
377 this%NprcX = ptr_mesh2d%NprcX
378 this%NprcY = ptr_mesh2d%NprcY
379 class default
380 log_info('MeshField_SpetralTransform2D_Init',*) 'Unexpected mesh type is given. Check!'
381 call prc_abort
382 end select
383
384 allocate( this%k(ks:ke) )
385 lx = this%xmax_gl - this%xmin_gl
386 do i=ks, ke
387 this%k(i) = i * 2.0_rp * pi / lx
388 end do
389
390 allocate( this%l(ls:le) )
391 ly = this%ymax_gl - this%ymin_gl
392 do j=ls, le
393 this%l(j) = j * 2.0_rp * pi / ly
394 end do
395
396 allocate( this%spectral_coef(ks:ke,ls:le,2,var_num) )
397
398 call elem_dummy%Init(elem%PolyOrder, .false.)
399
400 if ( this%NintGLpt > 0 ) then
401 allocate( this%IntIntrpMat(this%NintGLpt,elem%Nfp) )
402 allocate( this%IntXi(this%NintGLpt), this%IntWeight(this%NintGLpt) )
403
404 this%IntIntrpMat(:,:) = elem_dummy%GenIntGaussLegendreIntrpMat( glquadord, this%IntWeight, this%IntXi )
405 end if
406
407 if ( this%NsamplePtPerElem > 0 ) then
408 this%NsamplePtPerElem1D = nsampleptperelem1d
409 allocate( this%FFTIntrpMat(this%NsamplePtPerElem1D,elem%Nfp) )
410 allocate( this%FFT_Xi(this%NsamplePtPerElem1D) )
411
412 dx = 2.0_rp / real(this%NsamplePtPerElem1D,kind=rp)
413 do i=1, this%NsamplePtPerElem1D
414 this%FFT_Xi(i) = -1.0_rp + real(i - 0.5_rp,kind=rp) * dx
415 end do
416 this%FFTIntrpMat(:,:) = polynominal_genlagrangepoly( elem_dummy%PolyOrder, elem_dummy%x1, this%FFT_xi )
417
418 call this%fft_x%Init( this%NsamplePtPerElem1D * this%NeGX )
419 call this%fft_y%Init( this%NsamplePtPerElem1D * this%NeGY )
420 end if
421
422 call elem_dummy%Final()
423 return
424 end subroutine meshfield_spetraltransform2d_init
425
426!OCL SERIAL
427 subroutine meshfield_spetraltransform2d_final( this )
428 implicit none
429 class(meshfield_spetraltransform2d), intent(inout) :: this
430 !--------------------------------------------
431
432 call meshfield_spetraltransformbase_final( this )
433
434 if ( this%NsamplePtPerElem > 0 ) then
435 call this%fft_x%Final()
436 call this%fft_y%Final()
437 end if
438
439 deallocate( this%k, this%l, this%spectral_coef )
440 return
441 end subroutine meshfield_spetraltransform2d_final
442
443!OCL SERIAL
444 subroutine meshfield_spetraltransform2d_transform( this, q_list, mesh_num_x, mesh_num_y )
445 implicit none
446 class(meshfield_spetraltransform2d), intent(inout) :: this
447 integer, intent(in) :: mesh_num_x, mesh_num_y
448 type(meshfield2dlist), target :: q_list(this%var_num,mesh_num_x,mesh_num_y)
449
450 class(meshbase2d), pointer :: mesh2d
451 class(localmesh2d), pointer :: lmesh
452 !--------------------------------------------
453
454 mesh2d => q_list(1,1,1)%ptr%mesh
455 lmesh => mesh2d%lcmesh_list(1)
456
457 select case( this%eval_type_id )
459 call spectral_transform2d_dft( this%spectral_coef, &
460 q_list, this%var_num, mesh_num_x, mesh_num_y, & ! (in)
461 this%ks, this%ke, this%ls, this%le, lmesh%refElem2D%Nfp, lmesh%NeX, lmesh%NeY, & ! (in)
462 this%NeGY, this%NprcY, this%NsamplePtPerElem1D, & ! (in)
463 this%fft_x, this%fft_y, this%FFTIntrpMat )
465 call spectral_transform2d_l2projection( this%spectral_coef, & ! (out)
466 q_list, this%var_num, mesh_num_x*mesh_num_y, & ! (in)
467 this%ks, this%ke, this%ls, this%le, lmesh%refElem2D%Nfp, lmesh%NeX, lmesh%NeY, & ! (in)
468 this%IntIntrpMat, this%IntXi, this%IntWeight, this%NintGLPt, & ! (in)
469 this%xmin_gl, this%delx, this%xmax_gl - this%xmin_gl, & ! (in)
470 this%ymin_gl, this%dely, this%ymax_gl - this%ymin_gl ) ! (in)
471 end select
472
473 return
474 end subroutine meshfield_spetraltransform2d_transform
475
476!--- Private subroutines (1D)
477
478!OCL SERIAL
479 subroutine spectral_transform1d_dft( spectral_coef, &
480 q_list, var_num, mesh_num, ks, ke, Np, Ne, NsamplePerElem, &
481 fft, FFTIntrpMat, FFT_xi, &
482 xmin_gl, delx, Lx )
483
484 use mpi
485 use scale_prc, only: &
486 prc_local_comm_world
487 implicit none
488 integer, intent(in) :: ks, ke
489 integer, intent(in) :: np
490 integer, intent(in) :: ne
491 integer, intent(in) :: nsampleperelem
492 integer, intent(in) :: var_num
493 integer, intent(in) :: mesh_num
494 real(rp), intent(out) :: spectral_coef(ks:ke,2,var_num)
495 type(meshfield1dlist), target :: q_list(var_num,mesh_num)
496 class(fastfouriertransform1d), intent(in) :: fft
497 real(rp), intent(in) :: fftintrpmat(nsampleperelem,np)
498 real(rp), intent(in) :: fft_xi(nsampleperelem)
499 real(rp), intent(in) :: xmin_gl
500 real(rp), intent(in) :: delx
501 real(rp), intent(in) :: lx
502
503 real(rp) :: g_q(nsampleperelem,ne,mesh_num,var_num)
504 real(rp) :: q_tmp(np)
505 complex(RP) :: s_q(nsampleperelem*ne*mesh_num,var_num)
506
507 integer :: m
508 integer :: v
509 integer :: kel
510
511 integer :: nall
512 integer :: kk
513
514 real(rp) :: x_(nsampleperelem)
515 !-----------------------------------------------
516
517 !$omp parallel do collapse(3) private(v,m,kel,q_tmp, x_)
518 do v=1, var_num
519 do m=1, mesh_num
520 do kel=1, ne
521 q_tmp(:) = q_list(v,m)%ptr%local(1)%val(:,kel)
522 g_q(:,kel,m,v) = matmul( fftintrpmat, q_tmp )
523 end do
524 end do
525 end do
526
527 !$omp parallel do
528 do v=1, var_num
529 call fft%Forward_real( g_q(:,:,:,v), s_q(:,v) )
530 end do
531
532 nall = nsampleperelem * ne * mesh_num
533 !$omp parallel do private(kk)
534 do v=1, var_num
535 do kk=1, nall/2+1
536 spectral_coef(kk-1,1,v) = real(s_q(kk,v))
537 spectral_coef(kk-1,2,v) = aimag(s_q(kk,v))
538 end do
539 do kk=nall/2+1, nall
540 spectral_coef(kk-1-nall,1,v) = real(s_q(kk,v))
541 spectral_coef(kk-1-nall,2,v) = aimag(s_q(kk,v))
542 end do
543 end do
544 return
545 end subroutine spectral_transform1d_dft
546
547!OCL SERIAL
548 subroutine spectral_transform1d_l2projection( spectral_coef, &
549 q_list, var_num, mesh_num, ks, ke, Np, Ne, IntIntrpMat, IntXi, IntWeight, NintGLPt, &
550 xmin_gl, delx, Lx )
551
552 use mpi
553 use scale_prc, only: &
554 prc_local_comm_world
555 implicit none
556 integer, intent(in) :: ks, ke
557 integer, intent(in) :: np
558 integer, intent(in) :: ne
559 integer, intent(in) :: var_num
560 integer, intent(in) :: mesh_num
561 class(meshfield1dlist), target :: q_list(var_num,mesh_num)
562 real(rp), intent(out) :: spectral_coef(ks:ke,2,var_num)
563 integer, intent(in) :: nintglpt
564 real(rp), intent(in) :: intintrpmat(nintglpt,np)
565 real(rp), intent(in) :: intxi(nintglpt)
566 real(rp), intent(in) :: intweight(nintglpt)
567 real(rp), intent(in) :: xmin_gl
568 real(rp), intent(in) :: delx
569 real(rp), intent(in) :: lx
570
571 real(rp), allocatable :: s_coef_lc(:,:,:)
572 real(rp), allocatable :: s_coef(:,:,:)
573 real(rp), allocatable :: q_tmp(:,:,:)
574
575 integer :: vec_size
576 integer :: km
577
578 class(meshbase1d), pointer :: mesh1d
579 class(localmesh1d), pointer :: lmesh
580 integer :: ldom
581 integer :: meshid
582
583 integer :: kel, v
584 integer :: k, kk
585 integer :: ierr
586 !-----------------------------------------------
587
588 vec_size = var_num
589 km = ke
590
591 allocate( s_coef_lc(vec_size,2,0:km) )
592 allocate( s_coef(vec_size,2,0:km) )
593
594 !-
595 mesh1d => q_list(1,1)%ptr%mesh
596 lmesh => mesh1d%lcmesh_list(1)
597
598 !-
599
600 allocate( q_tmp(np,vec_size,ne) )
601
602 s_coef_lc(:,:,:) = 0.0_rp
603
604 do meshid=1, mesh_num
605 mesh1d => q_list(1,meshid)%ptr%mesh
606 do ldom=1, mesh1d%LOCAL_MESH_NUM
607 lmesh => mesh1d%lcmesh_list(ldom)
608 !$omp parallel do collapse(2)
609 do kel=lmesh%NeS, lmesh%NeE
610 do v=1, vec_size
611 q_tmp(:,v,kel) = q_list(v,meshid)%ptr%local(ldom)%val(:,kel)
612 end do
613 end do
614 call spectral_transform1d_l2projection_lc( s_coef_lc, &
615 q_tmp, 0, km, np, ne, vec_size, intintrpmat, intxi, intweight, nintglpt, &
616 (lmesh%xmin - xmin_gl)/lx-0.5_rp, delx/lx )
617 end do
618 end do
619
620 ! global sum
621 call mpi_allreduce( s_coef_lc, s_coef, vec_size * (km+1) * 2, &
622 mpi_double_precision, mpi_sum, prc_local_comm_world, ierr )
623
624 !-
625 !$omp parallel do private(kk)
626 do v=1, var_num
627 do k=ks, ke
628 kk = abs(k)
629 spectral_coef(k,1,v) = s_coef(v,1,kk)
630 spectral_coef(k,2,v) = s_coef(v,2,kk)
631 end do
632 end do
633
634 return
635 end subroutine spectral_transform1d_l2projection
636!OCL SERIAL
637 subroutine spectral_transform1d_l2projection_lc( spectral_coef, &
638 q, ks, ke, Np, Ne, vec_size, IntIntrpMat, IntXi, IntWeight, NintGLPt, &
639 xmin_lc, delx )
640 implicit none
641 integer, intent(in) :: ks, ke
642 integer, intent(in) :: np
643 integer, intent(in) :: ne
644 integer, intent(in) :: vec_size
645 real(rp), intent(in) :: q(np,vec_size,ne)
646 real(rp), intent(inout) :: spectral_coef(vec_size,2,ks:ke)
647 integer, intent(in) :: nintglpt
648 real(rp), intent(in) :: intintrpmat(nintglpt,np)
649 real(rp), intent(in) :: intxi(nintglpt)
650 real(rp), intent(in) :: intweight(nintglpt)
651 real(rp), intent(in) :: xmin_lc
652 real(rp), intent(in) :: delx
653
654 integer :: k
655 integer :: kel
656 integer :: v
657
658 real(rp) :: q_intrp(nintglpt,vec_size,ne)
659 real(rp) :: int_x(nintglpt,ne)
660
661 real(rp) :: phi(nintglpt)
662 real(rp) :: cos_kx(nintglpt)
663 real(rp) :: sin_kx(nintglpt)
664
665 real(rp) :: int_w(nintglpt)
666 !-------------------------------
667
668 int_w(:) = 0.5_rp * intweight(:) * delx
669
670 !$omp parallel private(kel,k,v,cos_kx,sin_kx,phi)
671 !$omp do
672 do kel=1, ne
673 int_x(:,kel) = 2.0_rp * pi * ( xmin_lc + delx * ( ( kel - 1.0_rp ) + 0.5_rp * ( 1.0_rp + intxi(:) ) ) )
674 q_intrp(:,:,kel) = matmul(intintrpmat, q(:,:,kel))
675 end do
676 !$omp do
677 do k=ks, ke
678 do kel=1, ne
679 phi(:) = mod( k * int_x(:,kel), 2.0_rp * pi )
680 cos_kx(:) = cos( phi )
681 sin_kx(:) = sin( phi )
682
683 do v=1, vec_size
684 spectral_coef(v,1,k) = spectral_coef(v,1,k) + sum(int_w(:) * cos_kx(:) * q_intrp(:,v,kel))
685 spectral_coef(v,2,k) = spectral_coef(v,2,k) - sum(int_w(:) * sin_kx(:) * q_intrp(:,v,kel))
686 end do
687 end do
688 end do
689 !$omp end parallel
690 return
691 end subroutine spectral_transform1d_l2projection_lc
692
693!--- Private subroutines (2D)
694
695!OCL SERIAL
696 subroutine spectral_transform2d_dft( spectral_coef, &
697 q_list, var_num, mesh_num_x, mesh_num_y, ks, ke, ls, le, Np1D, NeX, NeY, NeGY, NprcY, NsamplePerElem1D, &
698 fft_x, fft_y, FFTIntrpMat )
699
700 use mpi
701 use scale_prc, only: &
702 prc_local_comm_world
703 implicit none
704 integer, intent(in) :: ks, ke
705 integer, intent(in) :: ls, le
706 integer, intent(in) :: np1d
707 integer, intent(in) :: nex
708 integer, intent(in) :: ney, negy, nprcy
709 integer, intent(in) :: nsampleperelem1d
710 integer, intent(in) :: var_num
711 integer, intent(in) :: mesh_num_x
712 integer, intent(in) :: mesh_num_y
713 real(rp), intent(out) :: spectral_coef(ks:ke,ls:le,2,var_num)
714 class(meshfield2dlist), target :: q_list(var_num,mesh_num_x,mesh_num_y)
715 class(fastfouriertransform1d), intent(in) :: fft_x
716 class(fastfouriertransform1d), intent(in) :: fft_y
717 real(rp), intent(in) :: fftintrpmat(nsampleperelem1d,np1d)
718
719 real(rp) :: fftintrpmat_tr(np1d,nsampleperelem1d)
720
721 real(rp) :: g_q(nsampleperelem1d*nex*mesh_num_x,nsampleperelem1d*ney*mesh_num_y*var_num)
722 complex(RP) :: s_qx(nsampleperelem1d*nex*mesh_num_x,size(g_q,2))
723 complex(RP) :: g_qy(nsampleperelem1d*ney*mesh_num_y*nprcy,var_num,nsampleperelem1d*nex*mesh_num_x/nprcy)
724 complex(RP) :: s_q_lc(nsampleperelem1d*negy,var_num,nsampleperelem1d*nex*mesh_num_x/nprcy)
725 complex(RP) :: s_q(nsampleperelem1d*negy,var_num,nsampleperelem1d*nex*mesh_num_x)
726
727 integer :: v
728 integer :: xdim
729 integer :: prcy
730 integer :: j
731
732 integer :: nall_x, nall_y
733 integer :: kk, ll, kk_os
734
735 integer :: sendcount, recvcount
736 complex(RP) :: sendbuf(nsampleperelem1d*ney*mesh_num_y*var_num,nsampleperelem1d*nex*mesh_num_x)
737 complex(RP) :: recvbuf(nsampleperelem1d*ney*mesh_num_y,var_num,nsampleperelem1d*nex*mesh_num_x/nprcy,nprcy)
738 integer :: nx, local_nx
739 integer :: ny, local_ny
740 integer :: ierr
741 !-----------------------------------------------
742
743 !-
744 nall_x = nsampleperelem1d * nex * mesh_num_x
745 nall_y = nsampleperelem1d * negy
746
747 !-
748 fftintrpmat_tr(:,:) = transpose(fftintrpmat)
749 call spectral_transform2d_dft_sampling( g_q, &
750 q_list, var_num, mesh_num_x, mesh_num_y, np1d, nex, ney, nsampleperelem1d, &
751 fftintrpmat_tr )
752
753 !-
754 !$omp parallel do
755 do v=1, size(g_q,2)
756 call fft_x%Forward( g_q(:,v), s_qx(:,v) )
757 end do
758
759 !-
760 local_nx = nall_x
761 nx = local_nx
762
763 local_ny = nsampleperelem1d * ney * mesh_num_y * var_num
764 ny = nsampleperelem1d * negy * var_num
765
766 sendbuf(:,:) = transpose(s_qx)
767 sendcount = local_ny * (nx / nprcy)
768 recvcount = sendcount
769 call mpi_alltoall( sendbuf, sendcount, mpi_double_complex, &
770 recvbuf, recvcount, mpi_double_complex, prc_local_comm_world, ierr )
771
772 local_ny = nsampleperelem1d * ney * mesh_num_y
773 do v=1, var_num
774 do prcy=1, nprcy
775 do xdim=1, nsampleperelem1d*nex*mesh_num_x/nprcy
776 do j=1, local_ny
777 g_qy(j+(prcy-1)*local_ny,v,xdim) = recvbuf(j,v,xdim,prcy)
778 end do
779 end do
780 end do
781 end do
782
783 !-
784 !$omp parallel do collapse(2)
785 do v=1, var_num
786 do xdim=1, nx/nprcy
787 call fft_y%Forward( g_qy(:,v,xdim), s_q_lc(:,v,xdim) )
788 end do
789 end do
790
791 sendcount = size(s_q_lc)
792 recvcount = size(s_q_lc)
793 call mpi_allgather( s_q_lc, sendcount, mpi_double_complex, &
794 s_q, recvcount, mpi_double_complex, prc_local_comm_world, ierr )
795
796 !-
797 !$omp parallel do private(v,kk,ll)
798 do v=1, var_num
799 do ll=1, nall_y/2+1
800 do kk=1, nall_x/2+1
801 spectral_coef(kk-1,ll-1,1,v) = real(s_q(ll,v,kk))
802 spectral_coef(kk-1,ll-1,2,v) = aimag(s_q(ll,v,kk))
803 end do
804 do kk=nall_x/2+1, nall_x
805 spectral_coef(kk-1-nall_x,ll-1,1,v) = real(s_q(ll,v,kk))
806 spectral_coef(kk-1-nall_x,ll-1,2,v) = aimag(s_q(ll,v,kk))
807 end do
808 end do
809 do ll=nall_y/2+1, nall_y
810 do kk=1, nall_x/2+1
811 spectral_coef(kk-1,ll-1-nall_y,1,v) = real(s_q(ll,v,kk))
812 spectral_coef(kk-1,ll-1-nall_y,2,v) = aimag(s_q(ll,v,kk))
813 end do
814 do kk=nall_x/2+1, nall_x
815 spectral_coef(kk-1-nall_x,ll-1-nall_y,1,v) = real(s_q(ll,v,kk))
816 spectral_coef(kk-1-nall_x,ll-1-nall_y,2,v) = aimag(s_q(ll,v,kk))
817 end do
818 end do
819 end do
820 return
821 end subroutine spectral_transform2d_dft
822
823!OCL SERIAL
824 subroutine spectral_transform2d_dft_sampling( g_q, &
825 q_list, var_num, mesh_num_x, mesh_num_y, Np1D, NeX, NeY, NsamplePerElem1D, &
826 FFTIntrpMat_tr )
827
828 use mpi
829 use scale_prc, only: &
830 prc_local_comm_world
831 implicit none
832 integer, intent(in) :: np1d
833 integer, intent(in) :: nex, ney
834 integer, intent(in) :: nsampleperelem1d
835 integer, intent(in) :: var_num
836 integer, intent(in) :: mesh_num_x
837 integer, intent(in) :: mesh_num_y
838 real(rp), intent(out) :: g_q(nsampleperelem1d,nex,mesh_num_x,nsampleperelem1d,ney,mesh_num_y,var_num)
839 class(meshfield2dlist), intent(in) :: q_list(var_num,mesh_num_x,mesh_num_y)
840 real(rp), intent(in) :: fftintrpmat_tr(np1d,nsampleperelem1d)
841
842 real(rp) :: q_tmp(np1d,np1d)
843 real(rp) :: q_intrp1(nsampleperelem1d,np1d)
844 real(rp) :: q_intrp2(nsampleperelem1d,nsampleperelem1d)
845
846 integer :: mx, my, v
847 integer :: kel_x, kel_y, kel
848 integer :: p, p1, pp1, p2, pp2
849 !-----------------------------------------------
850
851 !$omp parallel do private(mx,my,v, kel_x,kel_y,kel, p,p1,pp1,p2,pp2, q_tmp,q_intrp1,q_intrp2) collapse(5)
852 do v=1, var_num
853 do my=1, mesh_num_y
854 do mx=1, mesh_num_x
855 do kel_y=1, ney
856 do kel_x=1, nex
857 kel = kel_x + (kel_y-1)*nex
858 do p2=1, np1d
859 do p1=1, np1d
860 p = p1 + (p2-1)*np1d
861 q_tmp(p1,p2) = q_list(v,mx,my)%ptr%local(1)%val(p,kel)
862 end do
863 end do
864
865 q_intrp1(:,:) = 0.0_rp
866 do p2=1, np1d
867 do p1=1, nsampleperelem1d
868 do pp1=1, np1d
869 q_intrp1(p1,p2) = q_intrp1(p1,p2) + fftintrpmat_tr(pp1,p1) * q_tmp(pp1,p2)
870 end do
871 end do
872 end do
873 q_intrp2(:,:) = 0.0_rp
874 do p2=1, nsampleperelem1d
875 do pp2=1, np1d
876 do p1=1, nsampleperelem1d
877 q_intrp2(p1,p2) = q_intrp2(p1,p2) + fftintrpmat_tr(pp2,p2) * q_intrp1(p1,pp2)
878 end do
879 end do
880 end do
881
882 do p2=1, nsampleperelem1d
883 g_q(:,kel_x,mx,p2,kel_y,my,v) = q_intrp2(:,p2)
884 end do
885 end do
886 end do
887 end do
888 end do
889 end do
890 return
891 end subroutine spectral_transform2d_dft_sampling
892
893!OCL SERIAL
894 subroutine spectral_transform2d_l2projection( spectral_coef, &
895 q_list, var_num, mesh_num, ks, ke, ls, le, Np1D, NeX, NeY, &
896 IntIntrpMat1D, IntXi1D, IntWeight1D, NintGLPt1D, &
897 xmin_gl, delx, Lx, ymin_gl, dely, Ly )
898
899 use mpi
900 use scale_prc, only: &
901 prc_local_comm_world
902 implicit none
903 integer, intent(in) :: ks, ke
904 integer, intent(in) :: ls, le
905 integer, intent(in) :: np1d
906 integer, intent(in) :: nex, ney
907 integer, intent(in) :: var_num
908 integer, intent(in) :: mesh_num
909 type(meshfield2dlist), target :: q_list(var_num,mesh_num)
910 real(rp), intent(out) :: spectral_coef(ks:ke,ls:le,2,var_num)
911 integer, intent(in) :: nintglpt1d
912 real(rp), intent(in) :: intintrpmat1d(nintglpt1d,np1d)
913 real(rp), intent(in) :: intxi1d(nintglpt1d)
914 real(rp), intent(in) :: intweight1d(nintglpt1d)
915 real(rp), intent(in) :: xmin_gl
916 real(rp), intent(in) :: delx
917 real(rp), intent(in) :: lx
918 real(rp), intent(in) :: ymin_gl
919 real(rp), intent(in) :: dely
920 real(rp), intent(in) :: ly
921
922 real(rp), allocatable :: s_coef_lc(:,:,:,:)
923 real(rp), allocatable :: s_coef(:,:,:,:)
924 real(rp), allocatable :: q_tmp(:,:,:)
925
926 real(rp) :: intintrpmat1d_tr(np1d,nintglpt1d)
927
928 integer :: vec_size
929 integer :: km
930
931 class(meshbase2d), pointer :: mesh2d
932 class(localmesh2d), pointer :: lmesh
933 integer :: ldom
934 integer :: meshid
935
936 integer :: kel, v
937 integer :: k, kk, l
938 integer :: ierr
939 !-----------------------------------------------
940
941 vec_size = var_num
942 km = ke
943
944 allocate( s_coef_lc(vec_size,2,0:km,ls:le) )
945 allocate( s_coef(vec_size,2,0:km,ls:le) )
946
947 !-
948 mesh2d => q_list(1,1)%ptr%mesh
949 lmesh => mesh2d%lcmesh_list(1)
950
951 intintrpmat1d_tr(:,:) = transpose(intintrpmat1d)
952 !-
953
954 allocate( q_tmp(np1d**2,vec_size,nex*ney) )
955
956 s_coef_lc(:,:,:,:) = 0.0_rp
957
958 do meshid=1, mesh_num
959 mesh2d => q_list(1,meshid)%ptr%mesh
960 do ldom=1, mesh2d%LOCAL_MESH_NUM
961 lmesh => mesh2d%lcmesh_list(ldom)
962
963 !$omp parallel do collapse(2)
964 do kel=lmesh%NeS, lmesh%NeE
965 do v=1, vec_size
966 q_tmp(:,v,kel) = q_list(v,meshid)%ptr%local(ldom)%val(:,kel)
967 end do
968 end do
969 call spectral_transform2d_l2projection_lc( s_coef_lc, &
970 q_tmp, 0, km, ls, le, np1d, nex, ney, vec_size, &
971 intintrpmat1d_tr, intxi1d, intweight1d, nintglpt1d, &
972 (lmesh%xmin - xmin_gl)/lx-0.5_rp, delx/lx, &
973 (lmesh%ymin - ymin_gl)/ly-0.5_rp, dely/ly )
974 end do
975 end do
976
977 ! global sum
978 call mpi_allreduce( s_coef_lc, s_coef, vec_size * (km+1)*(le-ls+1) * 2, &
979 mpi_double_precision, mpi_sum, prc_local_comm_world, ierr )
980
981 !-
982 !$omp parallel do collapse(2) private(kk)
983 do v=1, var_num
984 do l=ls, le
985 do k=ks, ke
986 kk = abs(k)
987 if ( k >= 0 ) then
988 spectral_coef(k,l,1,v) = s_coef(v,1,kk,l)
989 spectral_coef(k,l,2,v) = s_coef(v,2,kk,l)
990 else
991 spectral_coef(k,l,1,v) = s_coef(v,1,kk,-l)
992 spectral_coef(k,l,2,v) = - s_coef(v,2,kk,-l)
993 end if
994 end do
995 end do
996 end do
997
998 return
999 end subroutine spectral_transform2d_l2projection
1000!OCL SERIAL
1001 subroutine spectral_transform2d_l2projection_lc( spectral_coef, &
1002 q, ks, ke, ls, le, Np1D, NeX, NeY, vec_size, IntIntrpMat1D_tr, IntXi1D, IntWeight1D, NintGLpt1D, &
1003 xmin_lc, delx, ymin_lc, dely )
1004 implicit none
1005 integer, intent(in) :: ks, ke
1006 integer, intent(in) :: ls, le
1007 integer, intent(in) :: np1d
1008 integer, intent(in) :: nex
1009 integer, intent(in) :: ney
1010 integer, intent(in) :: vec_size
1011 real(rp), intent(in) :: q(np1d,np1d,vec_size,nex,ney)
1012 real(rp), intent(inout) :: spectral_coef(vec_size,2,ks:ke,ls:le)
1013 integer, intent(in) :: nintglpt1d
1014 real(rp), intent(in) :: intintrpmat1d_tr(np1d,nintglpt1d)
1015 real(rp), intent(in) :: intxi1d(nintglpt1d)
1016 real(rp), intent(in) :: intweight1d(nintglpt1d)
1017 real(rp), intent(in) :: xmin_lc
1018 real(rp), intent(in) :: delx
1019 real(rp), intent(in) :: ymin_lc
1020 real(rp), intent(in) :: dely
1021
1022 integer :: k, l
1023 integer :: kel_x, kel_y
1024 integer :: p, px, py
1025 integer :: v
1026
1027 real(rp) :: q_intrp_tmp(nintglpt1d,np1d,vec_size)
1028 real(rp) :: q_intrp_tmp2(nintglpt1d**2,vec_size)
1029 real(rp) :: q_intrp(nintglpt1d**2,vec_size,nex,ney)
1030 real(rp) :: int_x(nintglpt1d,nex)
1031 real(rp) :: int_y(nintglpt1d,ney)
1032
1033 real(rp) :: phi(nintglpt1d**2)
1034 real(rp) :: cos_phi(nintglpt1d**2)
1035 real(rp) :: sin_phi(nintglpt1d**2)
1036
1037 real(rp) :: int_w(nintglpt1d**2)
1038 !-------------------------------
1039
1040 !$omp parallel private(kel_x,kel_y,k,l,v,p,px,py,cos_phi,sin_phi,phi,q_intrp_tmp,q_intrp_tmp2)
1041 !$omp do
1042 do py=1, nintglpt1d
1043 do px=1, nintglpt1d
1044 p = px + (py-1)*nintglpt1d
1045 int_w(p) = 0.25_rp * delx * dely * intweight1d(px) * intweight1d(py)
1046 end do
1047 end do
1048 !$omp do
1049 do kel_x=1, nex
1050 int_x(:,kel_x) = 2.0_rp * pi * ( xmin_lc + delx * ( ( kel_x - 1.0_rp ) + 0.5_rp * ( 1.0_rp + intxi1d(:) ) ) )
1051 end do
1052 !$omp do
1053 do kel_y=1, ney
1054 int_y(:,kel_y) = 2.0_rp * pi * ( ymin_lc + dely * ( ( kel_y - 1.0_rp ) + 0.5_rp * ( 1.0_rp + intxi1d(:) ) ) )
1055 end do
1056 !$omp do collapse(2)
1057 do kel_y=1, ney
1058 do kel_x=1, nex
1059 q_intrp_tmp(:,:,:) = 0.0_rp
1060 do v=1, vec_size
1061 do py=1, np1d
1062 do px=1, nintglpt1d
1063 q_intrp_tmp(px,py,v) = q_intrp_tmp(px,py,v) &
1064 + sum(intintrpmat1d_tr(:,px) * q(:,py,v,kel_x,kel_y))
1065 end do
1066 end do
1067 end do
1068
1069 q_intrp_tmp2(:,:) = 0.0_rp
1070 do py=1, nintglpt1d
1071 do v=1, vec_size
1072 do px=1, nintglpt1d
1073 p = px + (py-1)*nintglpt1d
1074 q_intrp_tmp2(p,v) = q_intrp_tmp2(p,v) &
1075 + sum(intintrpmat1d_tr(:,py) * q_intrp_tmp(px,:,v))
1076 end do
1077 end do
1078 end do
1079 q_intrp(:,:,kel_x,kel_y) = q_intrp_tmp2(:,:)
1080 end do
1081 end do
1082 !$omp do collapse(2)
1083 do l=ls, le
1084 do k=ks, ke
1085 do kel_y=1, ney
1086 do kel_x=1, nex
1087 do py=1, nintglpt1d
1088 do px=1, nintglpt1d
1089 p = px + (py-1)*nintglpt1d
1090 phi(p) = mod( k * int_x(px,kel_x) + l * int_y(py,kel_y), 2.0_rp * pi )
1091 end do
1092 end do
1093 cos_phi(:) = cos( phi(:) )
1094 sin_phi(:) = sin( phi(:) )
1095
1096 do v=1, vec_size
1097 spectral_coef(v,1,k,l) = spectral_coef(v,1,k,l) + sum(int_w(:) * cos_phi(:) * q_intrp(:,v,kel_x,kel_y))
1098 spectral_coef(v,2,k,l) = spectral_coef(v,2,k,l) - sum(int_w(:) * sin_phi(:) * q_intrp(:,v,kel_x,kel_y))
1099 end do
1100 end do
1101 end do
1102 end do
1103 end do
1104 !$omp end parallel
1105 return
1106 end subroutine spectral_transform2d_l2projection_lc
1107
module FElib / Element / Base
module FElib / Element / line
module FElib / Mesh / Local 1D
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Base 1D
module FElib / Mesh / Base 2D
module FElib / Mesh / Rectangle 2D domain
module FElib / Data / base
Module common / Polynominal.
real(rp) function, dimension(size(x), nord+1), public polynominal_genlagrangepoly(nord, x_lgl, x)
A function to obtain the Lagrange basis functions related to the Gauss-Legendre-Lobatto (GLL) points.
Derived type representing a 1D reference element.
Derived type representing a 2D reference element.
Derived type representing a 3D reference element.
Derived type representing a line element.
Derived type to manage a local 3D computational domain.
Derived type representing a field with 1D mesh.
Derived type representing a field with 2D mesh.