Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
spectral_error.f90
Go to the documentation of this file.
1! Copyright (c) 2022, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only : rp
36 use field, only : field_t
37 use coefs, only : coef_t
38 use field_list, only : field_list_t
39 use math, only : rzero, copy
40 use file, only : file_t, file_free
41 use time_state, only : time_state_t
42 use tensor, only : tnsr3d
43 use device_math, only : device_copy
46 use logger, only : neko_log
48 use comm, only : pe_rank
52 use json_module, only : json_file
54 use case, only : case_t
55 use registry, only : neko_registry
56 use utils, only : neko_varname_len
57
58 use, intrinsic :: iso_c_binding
59 implicit none
60 private
61
67 type, public, extends(simulation_component_t) :: spectral_error_t
69 type(field_t), pointer :: u => null()
70 type(field_t), pointer :: v => null()
71 type(field_t), pointer :: w => null()
73 type(field_t), pointer :: u_hat => null()
74 type(field_t), pointer :: v_hat => null()
75 type(field_t), pointer :: w_hat => null()
77 type(field_t) :: wk
79 real(kind=rp) :: seri_small = 1.e-14
81 real(kind=rp) :: seri_smallr = 1.e-10
83 real(kind=rp) :: seri_smallg = 1.e-5
85 real(kind=rp) :: seri_smalls = 0.2
87 integer :: seri_np = 4
88 integer :: seri_np_max = 4
90 integer :: seri_elr = 0
92 real(kind=rp), allocatable :: eind_u(:), eind_v(:), eind_w(:)
94 real(kind=rp), allocatable :: sig_u(:), sig_v(:), sig_w(:)
96 type(field_list_t) :: speri_l
98 type(file_t) :: mf_speri
100 type(field_writer_t) :: writer
101 contains
103 procedure, pass(this) :: init => spectral_error_init
105 procedure, pass(this) :: free => spectral_error_free
107 procedure, pass(this) :: compute_ => spectral_error_compute
109 procedure, pass(this) :: get_indicators => spectral_error_get_indicators
110
111 end type spectral_error_t
112
113contains
114
116 subroutine spectral_error_init(this, json, case)
117 class(spectral_error_t), intent(inout), target :: this
118 type(json_file), intent(inout) :: json
119 class(case_t), intent(inout), target :: case
120
121 character(len=:), allocatable :: name
122 character(len=NEKO_VARNAME_LEN) :: fields(3)
123
124 call json_get_or_default(json, "name", name, "spectral_error")
125 this%name = name
127 ! picks it up. Will also add those fields to the registry.
128 fields(1) = "u_hat"
129 fields(2) = "v_hat"
130 fields(3) = "w_hat"
131 call json%add("fields", fields)
132
133 call this%init_base(json, case)
134 call this%writer%init(json, case)
135
136 call spectral_error_init_from_components(this, case%fluid%c_Xh)
137
138 end subroutine spectral_error_init
139
143 class(spectral_error_t), intent(inout) :: this
144 type(coef_t), intent(in) :: coef
145 integer :: il, jl, aa
146 character(len=NEKO_FNAME_LEN) :: fname_speri
147
148 this%u => neko_registry%get_field("u")
149 this%v => neko_registry%get_field("v")
150 this%w => neko_registry%get_field("w")
151 this%u_hat => neko_registry%get_field("u_hat")
152 this%v_hat => neko_registry%get_field("v_hat")
153 this%w_hat => neko_registry%get_field("w_hat")
154
156 this%wk = this%u
157
159 allocate(this%eind_u(coef%msh%nelv))
160 allocate(this%eind_v(coef%msh%nelv))
161 allocate(this%eind_w(coef%msh%nelv))
162 allocate(this%sig_u(coef%msh%nelv))
163 allocate(this%sig_v(coef%msh%nelv))
164 allocate(this%sig_w(coef%msh%nelv))
165
167 associate(lx1 => coef%Xh%lx, ly1 => coef%Xh%ly, &
168 lz1 => coef%Xh%lz, &
169 seri_small => this%SERI_SMALL, &
170 seri_smallr => this%SERI_SMALLR, &
171 seri_smallg => this%SERI_SMALLG, &
172 seri_smalls => this%SERI_SMALLS, &
173 seri_np => this%SERI_NP, &
174 seri_np_max => this%SERI_NP_MAX, &
175 seri_elr => this%SERI_ELR &
176 )
177 ! correctness check
178 if (seri_np .gt. seri_np_max) then
179 call neko_log%message('SETI_NP greater than SERI_NP_MAX')
180 end if
181 il = seri_np + seri_elr
182 jl = min(lx1, ly1)
183 jl = min(jl, lz1)
184 if (il .gt. jl) then
185 call neko_log%message('SERI_NP+SERI_ELR greater than L?1')
186 end if
187 end associate
188
190
192 subroutine spectral_error_free(this)
193 class(spectral_error_t), intent(inout) :: this
194
195 if (allocated(this%eind_u)) then
196 deallocate(this%eind_u)
197 end if
198
199 if (allocated(this%eind_v)) then
200 deallocate(this%eind_v)
201 end if
202
203 if (allocated(this%eind_w)) then
204 deallocate(this%eind_w)
205 end if
206
207 if (allocated(this%sig_u)) then
208 deallocate(this%sig_u)
209 end if
210
211 if (allocated(this%sig_v)) then
212 deallocate(this%sig_v)
213 end if
214
215 if (allocated(this%sig_w)) then
216 deallocate(this%sig_w)
217 end if
218
219 call this%wk%free()
220 call this%speri_l%free()
221
222 nullify(this%u)
223 nullify(this%v)
224 nullify(this%w)
225 nullify(this%u_hat)
226 nullify(this%v_hat)
227 nullify(this%w_hat)
228
229 call this%writer%free()
230 call this%free_base()
231
232 end subroutine spectral_error_free
233
235 subroutine spectral_error_compute(this, time)
236 class(spectral_error_t), intent(inout) :: this
237 type(time_state_t), intent(in) :: time
238
239 integer :: e, i, lx, ly, lz, nelv, n
240
241
242 call this%get_indicators(this%case%fluid%c_Xh)
243
244 lx = this%u_hat%Xh%lx
245 ly = this%u_hat%Xh%ly
246 lz = this%u_hat%Xh%lz
247 nelv = this%u_hat%msh%nelv
248
250 do e = 1, nelv
251 do i = 1, lx*ly*lx
252 this%u_hat%x(i, 1, 1, e) = this%eind_u(e)
253 this%v_hat%x(i, 1, 1, e) = this%eind_v(e)
254 this%w_hat%x(i, 1, 1, e) = this%eind_w(e)
255 end do
256 end do
257
258 ! We need this copy to GPU since the field writer does an internal copy
259 ! GPU -> CPU before writing the field files. This overwrites the *_hat
260 ! host arrays that contain the values.
261 if (neko_bcknd_device .eq. 1) then
262 call device_memcpy(this%u_hat%x, this%u_hat%x_d, lx*ly*lz*nelv, &
263 host_to_device, sync = .true.)
264 call device_memcpy(this%v_hat%x, this%v_hat%x_d, lx*ly*lz*nelv, &
265 host_to_device, sync = .true.)
266 call device_memcpy(this%w_hat%x, this%w_hat%x_d, lx*ly*lz*nelv, &
267 host_to_device, sync = .true.)
268 end if
269
270 end subroutine spectral_error_compute
271
280 subroutine transform_to_spec_or_phys(u_hat, u, wk, coef, space)
281 type(field_t), intent(inout) :: u_hat
282 type(field_t), intent(inout) :: u
283 type(field_t), intent(inout) :: wk
284 type(coef_t), intent(inout) :: coef
285 character(len=4), intent(in) :: space
286 integer :: i, j, k, e, nxyz, nelv, n
287
289 nxyz = coef%Xh%lx * coef%Xh%lx * coef%Xh%lx
290 nelv = coef%msh%nelv
291 n = nxyz*nelv
292
294 if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
295 (neko_bcknd_opencl .eq. 1)) then
296 call device_copy(wk%x_d, u%x_d, n)
297 else
298 call copy(wk%x, u%x, n)
299 end if
300
301 select case (space)
302 case ('spec')
303 call tnsr3d(u_hat%x, coef%Xh%lx, wk%x, &
304 coef%Xh%lx, coef%Xh%vinv, &
305 coef%Xh%vinvt, coef%Xh%vinvt, nelv)
306 case ('phys')
307 call tnsr3d(u_hat%x, coef%Xh%lx, wk%x, &
308 coef%Xh%lx, coef%Xh%v, &
309 coef%Xh%vt, coef%Xh%vt, nelv)
310 end select
311
312 ! Synchronize
313 if ((neko_bcknd_hip .eq. 1) .or. (neko_bcknd_cuda .eq. 1) .or. &
314 (neko_bcknd_opencl .eq. 1)) then
315
316 call device_memcpy(u_hat%x, u_hat%x_d, n, &
317 device_to_host, sync = .true.)
318 end if
319
320 end subroutine transform_to_spec_or_phys
321
324 subroutine spectral_error_get_indicators(this, coef)
325 class(spectral_error_t), intent(inout) :: this
326 type(coef_t), intent(inout) :: coef
327 integer :: i
328
329 ! Generate the uvwhat field (legendre coeff)
330 call transform_to_spec_or_phys(this%u_hat, this%u, this%wk, coef, 'spec')
331 call transform_to_spec_or_phys(this%v_hat, this%v, this%wk, coef, 'spec')
332 call transform_to_spec_or_phys(this%w_hat, this%w, this%wk, coef, 'spec')
333
334 ! Get the spectral error indicator
335 call calculate_indicators(this, coef, this%eind_u, this%sig_u, &
336 coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
337 this%u_hat%x)
338 call calculate_indicators(this, coef, this%eind_v, this%sig_v, &
339 coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
340 this%v_hat%x)
341 call calculate_indicators(this, coef, this%eind_w, this%sig_w, &
342 coef%msh%nelv, coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
343 this%w_hat%x)
344
345 end subroutine spectral_error_get_indicators
346
349 subroutine spectral_error_write(this, t)
350 class(spectral_error_t), intent(inout) :: this
351 real(kind=rp), intent(in) :: t
352
353 integer i, e
354 integer lx, ly, lz, nelv
355
359 call this%mf_speri%write(this%speri_l, t)
360
361 end subroutine spectral_error_write
362
372 subroutine calculate_indicators(this, coef, eind, sig, lnelt, LX1, LY1, LZ1, &
373 var)
374 type(spectral_error_t), intent(inout) :: this
375 type(coef_t), intent(in) :: coef
376 integer :: lnelt
377 integer :: LX1
378 integer :: LY1
379 integer :: LZ1
380 real(kind=rp) :: eind(lnelt)
381 real(kind=rp) :: sig(lnelt)
382 real(kind=rp) :: var(lx1, ly1, lz1, lnelt)
383
384 real(kind=rp) :: xa(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz)
385 real(kind=rp) :: xb(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz)
386 integer :: i, e
387
388 ! zero arrays
389 call rzero(eind, lnelt)
390 call rzero(sig, lnelt)
391
392 ! Get the indicator
393 call speri_var(this, eind, sig, var, lnelt, xa, xb, lx1, ly1, lz1)
394
395 end subroutine calculate_indicators
396
397
408 subroutine speri_var(this, est, sig, var, nell, xa, xb, LX1, LY1, LZ1)
409 type(spectral_error_t), intent(inout) :: this
410 integer :: nell
411 integer :: LX1
412 integer :: LY1
413 integer :: LZ1
414 real(kind=rp) :: est(nell)
415 real(kind=rp) :: sig(nell)
416 real(kind=rp) :: var(lx1, ly1, lz1, nell)
417 real(kind=rp) :: xa(lx1, ly1, lz1)
418 real(kind=rp) :: xb(lx1, ly1, lz1)
419
421 integer :: il, jl, kl, ll, j_st, j_en, ii
423 real(kind=rp) :: coeff(lx1, ly1, lz1)
425 real(kind=rp) :: coef11
427 real(kind=rp) :: coefx(this%SERI_NP_MAX, ly1, lz1), &
428 coefy(this%SERI_NP_MAX, lx1, lz1), &
429 coefz(this%SERI_NP_MAX, lx1, ly1)
431 real(kind=rp) :: estx, esty, estz
433 real(kind=rp) :: sigx, sigy, sigz
434 real(kind=rp) :: third
435 parameter(third = 1.0/3.0)
436
438 do il = 1, nell
441 do ii = 1, lx1*ly1*lz1
442 coeff(ii, 1, 1) = var(ii, 1, 1, il) * var(ii, 1, 1, il)
443 end do
444
446 coef11 = coeff(1, 1, 1)
447
449 if (coef11 .ge. this%SERI_SMALL) then
454 j_st = max(1, lx1 - this%SERI_NP + 1 - this%SERI_ELR)
455 j_en = max(1, lx1 - this%SERI_ELR)
456 do ll = 1, lz1
457 do kl = 1, ly1
458 do jl = j_st, j_en
459 coefx(j_en - jl + 1, kl, ll) = coeff(jl, kl, ll)
460 end do
461 end do
462 end do
464 call speri_extrap(this, estx, sigx, coef11, coefx, &
465 j_st, j_en, ly1, lz1)
466
470 j_st = max(1, ly1 - this%SERI_NP + 1 - this%SERI_ELR)
471 j_en = max(1, ly1 - this%SERI_ELR)
472 do ll = 1, lz1
473 do kl = j_st, j_en
474 do jl = 1, lx1
475 coefy(j_en - kl + 1, jl, ll) = coeff(jl, kl, ll)
476 end do
477 end do
478 end do
480 call speri_extrap(this, esty, sigy, coef11, coefy, &
481 j_st, j_en, lx1, lz1)
482
486 j_st = max(1, lz1 - this%SERI_NP + 1 - this%SERI_ELR)
487 j_en = max(1, lz1 - this%SERI_ELR)
488 do ll = j_st, j_en
489 do kl = 1, ly1
490 do jl = 1, lx1
491 coefz(j_en - ll + 1, jl, kl) = coeff(jl, kl, ll)
492 end do
493 end do
494 end do
496 call speri_extrap(this, estz, sigz, coef11, coefz, &
497 j_st, j_en, lx1, ly1)
498
500 est(il) = sqrt(estx + esty + estz)
501 sig(il) = third*(sigx + sigy + sigz)
502
503 else
505 estx = 0.0
506 esty = 0.0
507 estz = 0.0
508 sigx = -1.0
509 sigy = -1.0
510 sigz = -1.0
512
513 est(il) = 0.0
514 sig(il) = -1.0
515 end if
516
517 end do
518
519 end subroutine speri_var
520
530 subroutine speri_extrap(this, estx, sigx, coef11, coef, &
531 ix_st, ix_en, nyl, nzl)
532 implicit none
533 type(spectral_error_t), intent(inout) :: this
535 integer :: ix_st, ix_en, nyl, nzl
537 real(kind=rp) :: coef(this%SERI_NP_MAX, nyl, nzl)
539 real(kind=rp) :: coef11
541 real(kind=rp) :: estx, sigx
542
544 integer :: il, jl, kl, ll ! loop index
545 integer :: nsigt, pnr, nzlt
546 real(kind=rp) :: sigt, smallr, cmin, cmax, cnm, rtmp, rtmp2, rtmp3
547 real(kind=rp) :: sumtmp(4), cffl(this%SERI_NP_MAX)
548 real(kind=rp) :: stmp, estt, clog, ctmp, cave, erlog
549 logical :: cuse(this%seri_np_max)
550
551 associate(seri_small => this%SERI_SMALL, &
552 seri_smallr => this%SERI_SMALLR, &
553 seri_smallg => this%SERI_SMALLG, &
554 seri_smalls => this%SERI_SMALLS, &
555 seri_np => this%SERI_NP, &
556 seri_np_max => this%SERI_NP_MAX, &
557 seri_elr => this%SERI_ELR &
558 )
559 ! initial values
560 estx = 0.0
561 sigx = -1.0
562
563 ! relative cutoff
564 smallr = coef11*seri_smallr
565
566 ! number of points
567 pnr = ix_en - ix_st + 1
568
569 ! to few points to interpolate
570 !if ((ix_en - ix_st).le.1) return
571
572 ! for averaging, initial values
573 sigt = 0.0
574 nsigt = 0
575
576 ! loop over all face points
577 nzlt = max(1, nzl - seri_elr) ! for 2D runs
578 do il = 1, nzlt
579 ! weight
580 rtmp3 = 1.0/(2.0*(il - 1) + 1.0)
581 do jl = 1, nyl - seri_elr
582
583 ! find min and max coef along single row
584 cffl(1) = coef(1, jl, il)
585 cmin = cffl(1)
586 cmax = cmin
587 do kl = 2, pnr
588 cffl(kl) = coef(kl, jl, il)
589 cmin = min(cmin, cffl(kl))
590 cmax = max(cmax, cffl(kl))
591 end do
592
593 ! are coefficients sufficiently big
594 if ((cmin .gt. 0.0) .and. (cmax .gt. smallr)) then
595 ! mark array position we use in iterpolation
596 do kl = 1, pnr
597 cuse(kl) = .true.
598 end do
599 ! max n for polynomial order
600 cnm = real(ix_en)
601
602 ! check if all the points should be taken into account
603 ! in original code by Catherine Mavriplis this part is written
604 ! for 4 points, so I place if statement first
605 if (pnr .eq. 4) then
606 ! should we neglect last values
607 if ((cffl(1) .lt. smallr) .and. &
608 (cffl(2) .lt. smallr)) then
609 if (cffl(3) .lt. smallr) then
610 cuse(1) = .false.
611 cuse(2) = .false.
612 cnm = real(ix_en - 2)
613 else
614 cuse(1) = .false.
615 cnm = real(ix_en - 1)
616 end if
617 else
618 ! should we take stronger gradient
619 if ((cffl(1)/cffl(2) .lt. seri_smallg) .and. &
620 (cffl(3)/cffl(4) .lt. seri_smallg)) then
621 cuse(1) = .false.
622 cuse(3) = .false.
623 cnm = real(ix_en - 1)
624 elseif ((cffl(2)/cffl(1) .lt. seri_smallg) .and. &
625 (cffl(4)/cffl(3) .lt. seri_smallg)) then
626 cuse(2) = .false.
627 cuse(4) = .false.
628 end if
629 end if
630 end if
631
632 ! get sigma for given face point
633 do kl = 1, 4
634 sumtmp(kl) = 0.0
635 end do
636 ! find new min and count number of points
637 cmin = cmax
638 cmax = 0.0
639 do kl = 1, pnr
640 if (cuse(kl)) then
641 rtmp = real(ix_en - kl)
642 rtmp2 = log(cffl(kl))
643 sumtmp(1) = sumtmp(1) + rtmp2
644 sumtmp(2) = sumtmp(2) + rtmp
645 sumtmp(3) = sumtmp(3) + rtmp*rtmp
646 sumtmp(4) = sumtmp(4) + rtmp2*rtmp
647 ! find new min and count used points
648 cmin = min(cffl(kl), cmin)
649 cmax = cmax + 1.0
650 end if
651 end do
652 ! decay rate along single row
653 stmp = (sumtmp(1)*sumtmp(2) - sumtmp(4)*cmax)/ &
654 (sumtmp(3)*cmax - sumtmp(2)*sumtmp(2))
655 ! for averaging
656 sigt = sigt + stmp
657 nsigt = nsigt + 1
658
659 ! get error estimator depending on calculated decay rate
660 estt = 0.0
661 if (stmp .lt. seri_smalls) then
662 estt = cmin
663 else
664 ! get averaged constant in front of c*exp(-sig*n)
665 clog = (sumtmp(1)+stmp*sumtmp(2))/cmax
666 ctmp = exp(clog)
667 ! average exponent
668 cave = sumtmp(1)/cmax
669 ! check quality of approximation comparing is to the
670 ! constant cave
671 do kl = 1, 2
672 sumtmp(kl) = 0.0
673 end do
674 do kl = 1, pnr
675 if (cuse(kl)) then
676 erlog = clog - stmp*real(ix_en - kl)
677 sumtmp(1) = sumtmp(1) + &
678 (erlog - log(cffl(kl)))**2
679 sumtmp(2) = sumtmp(2) + &
680 (erlog - cave)**2
681 end if
682 end do
683 rtmp = 1.0 - sumtmp(1)/sumtmp(2)
684 if (rtmp .lt. seri_smalls) then
685 estt = cmin
686 else
687 ! last coefficient is not included in error estimator
688 estt = ctmp/stmp*exp(-stmp*cnm)
689 end if
690 end if
691 ! add contribution to error estimator; variable weight
692 estx = estx + estt/(2.0*(jl - 1) + 1.0)*rtmp3
693 end if ! if((cmin.gt.0.0).and.(cmax.gt.smallr))
694 end do
695 end do
696 ! constant weight
697 ! Multiplication by 4 in 2D / 8 in 3D
698 ! Normalization of the error by the volume of the reference element
699 ! which is equal to 4 in 2D / 8 in 3D
700 ! ==> Both operations cancel each other
701 estx = estx/(2.0*(ix_en - 1) + 1.0)
702
703 ! final everaging
704 ! sigt = 2*sigma so we divide by 2
705 if (nsigt .gt. 0) then
706 sigx = 0.5*sigt/nsigt
707 end if
708
709 end associate
710
711 end subroutine speri_extrap
712
713end module spectral_error
double real
Copy data between host and device (or device and device)
Definition device.F90:71
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Defines a simulation case.
Definition case.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
integer, public pe_rank
MPI rank.
Definition comm.F90:56
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
Implements the field_writer_t type.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
subroutine file_free(this)
File operation destructor.
Definition file.f90:160
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
Definition math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:252
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:208
Build configurations.
integer, parameter neko_bcknd_hip
integer, parameter neko_bcknd_device
integer, parameter neko_bcknd_opencl
integer, parameter neko_bcknd_cuda
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, time)
Dummy compute function.
Defines a function space.
Definition space.f90:34
Implements type spectral_error_t.
subroutine speri_extrap(this, estx, sigx, coef11, coef, ix_st, ix_en, nyl, nzl)
Extrapolate the Legendre spectrum from the last points.
subroutine spectral_error_write(this, t)
Write error indicators in a field file.
subroutine spectral_error_free(this)
Destructor.
subroutine spectral_error_get_indicators(this, coef)
Transform and get the spectral error indicators.
subroutine spectral_error_init_from_components(this, coef)
Actual constructor.
subroutine spectral_error_init(this, json, case)
Constructor.
subroutine spectral_error_compute(this, time)
Compute the spectral error indicator.
subroutine speri_var(this, est, sig, var, nell, xa, xb, lx1, ly1, lz1)
Calculate the indicator in a specified variable.
subroutine calculate_indicators(this, coef, eind, sig, lnelt, lx1, ly1, lz1, var)
Wrapper for old fortran 77 subroutines.
subroutine transform_to_spec_or_phys(u_hat, u, wk, coef, space)
Transform a field u to u_hat into physical or spectral space the result of the transformation is in u...
Tensor operations.
Definition tensor.f90:61
subroutine, public tnsr3d(v, nv, u, nu, a, bt, ct, nelv)
Tensor product performed on nelv elements.
Definition tensor.f90:234
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
integer, parameter, public neko_fname_len
Definition utils.f90:41
integer, parameter, public neko_varname_len
Definition utils.f90:42
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
field_list_t, To be able to group fields together
A simulation component that writes a 3d field to a file.
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:56
Base abstract class for simulation components.
Provides tools to calculate the spectral error indicator.
A struct that contains all info about the time, expand as needed.
#define max(a, b)
Definition tensor.cu:40