Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_scheme.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2026, 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!
34
36 use gather_scatter, only : gs_t
37 use checkpoint, only : chkp_t
38 use num_types, only : rp
39 use field, only : field_t
40 use field_list, only : field_list_t
41 use space, only : space_t
42 use dofmap, only : dofmap_t
43 use krylov, only : ksp_t, krylov_solver_factory, ksp_max_iter, ksp_monitor_t
44 use coefs, only : coef_t
45 use jacobi, only : jacobi_t
47 use sx_jacobi, only : sx_jacobi_t
48 use hsmg, only : hsmg_t
49 use bc_list, only : bc_list_t
50 use bc, only : bc_t
51 use precon, only : pc_t, precon_factory, precon_destroy
52 use mesh, only : mesh_t
55 use registry, only : neko_registry
58 use json_module, only : json_file
61 use utils, only : neko_error
69 use time_state, only : time_state_t
71 implicit none
72
74 type, abstract :: scalar_scheme_t
76 character(len=:), allocatable :: name
78 type(field_t), pointer :: u
80 type(field_t), pointer :: v
82 type(field_t), pointer :: w
84 type(field_t), pointer :: s
86 type(field_series_t) :: slag
88 type(space_t), pointer :: xh
90 type(dofmap_t), pointer :: dm_xh
92 type(gs_t), pointer :: gs_xh
94 type(coef_t), pointer :: c_xh
96 type(field_t), pointer :: f_xh => null()
100 class(ksp_t), allocatable :: ksp
102 integer :: ksp_maxiter
104 integer :: projection_dim
105
106 integer :: projection_activ_step
108 class(pc_t), allocatable :: pc
110 type(bc_list_t) :: bcs
112 type(json_file), pointer :: params
114 type(mesh_t), pointer :: msh => null()
116 type(chkp_t), pointer :: chkp => null()
118 character(len=:), allocatable :: nut_field_name
120 character(len=:), allocatable :: alphat_field_name
122 type(field_t), pointer :: rho => null()
124 type(field_t), pointer :: lambda => null()
126 type(field_t), pointer :: cp => null()
128 type(field_t), pointer :: lambda_tot => null()
130 real(kind=rp) :: pr_turb
132 type(field_list_t) :: material_properties
133 procedure(user_material_properties_intf), nopass, pointer :: &
134 user_material_properties => null()
136 logical :: freeze = .false.
137 contains
139 procedure, pass(this) :: scheme_init => scalar_scheme_init
141 procedure, pass(this) :: scheme_free => scalar_scheme_free
143 procedure, pass(this) :: validate => scalar_scheme_validate
145 procedure, pass(this) :: set_material_properties => &
148 procedure, pass(this) :: update_material_properties => &
151 procedure(scalar_scheme_init_intrf), pass(this), deferred :: init
153 procedure(scalar_scheme_free_intrf), pass(this), deferred :: free
155 procedure(scalar_scheme_step_intrf), pass(this), deferred :: step
157 procedure(scalar_scheme_restart_intrf), pass(this), deferred :: restart
158 end type scalar_scheme_t
159
161 abstract interface
162 subroutine scalar_scheme_init_intrf(this, msh, coef, gs, params, &
163 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
164 import scalar_scheme_t
165 import json_file
166 import coef_t
167 import gs_t
168 import mesh_t
169 import user_t
170 import field_series_t, field_t
172 import rp
173 import chkp_t
174 class(scalar_scheme_t), target, intent(inout) :: this
175 type(mesh_t), target, intent(in) :: msh
176 type(coef_t), target, intent(in) :: coef
177 type(gs_t), target, intent(inout) :: gs
178 type(json_file), target, intent(inout) :: params
179 type(json_file), target, intent(inout) :: numerics_params
180 type(user_t), target, intent(in) :: user
181 type(chkp_t), target, intent(inout) :: chkp
182 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
183 type(time_scheme_controller_t), target, intent(in) :: time_scheme
184 type(field_t), target, intent(in) :: rho
185 end subroutine scalar_scheme_init_intrf
186 end interface
187
189 abstract interface
190 subroutine scalar_scheme_restart_intrf(this, chkp)
191 import scalar_scheme_t
192 import chkp_t
193 import rp
194 class(scalar_scheme_t), target, intent(inout) :: this
195 type(chkp_t), intent(inout) :: chkp
196 end subroutine scalar_scheme_restart_intrf
197 end interface
198
200 abstract interface
202 import scalar_scheme_t
203 class(scalar_scheme_t), intent(inout) :: this
204 end subroutine scalar_scheme_free_intrf
205 end interface
206
208 abstract interface
209 subroutine scalar_scheme_step_intrf(this, time, ext_bdf, dt_controller, &
210 ksp_results)
211 import scalar_scheme_t
212 import time_state_t
215 import ksp_monitor_t
216 class(scalar_scheme_t), intent(inout) :: this
217 type(time_state_t), intent(in) :: time
218 type(time_scheme_controller_t), intent(in) :: ext_bdf
219 type(time_step_controller_t), intent(in) :: dt_controller
220 type(ksp_monitor_t), intent(inout) :: ksp_results
221 end subroutine scalar_scheme_step_intrf
222 end interface
223
224 ! ========================================================================== !
225 ! Helper functions and types for scalar_scheme_t
226 ! ========================================================================== !
227
230 class(scalar_scheme_t), allocatable :: scalar
231 contains
233 procedure, pass(this) :: init => scalar_scheme_wrapper_init
235 procedure, pass(this) :: free => scalar_scheme_wrapper_free
238 procedure, pass(this) :: move_from => &
241 procedure, pass(this) :: is_allocated => &
244
245
246 interface
247
260 module subroutine scalar_scheme_factory(object, msh, coef, gs, params, &
261 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
262 class(scalar_scheme_t), allocatable, intent(inout) :: object
263 type(mesh_t), target, intent(in) :: msh
264 type(coef_t), target, intent(in) :: coef
265 type(gs_t), target, intent(inout) :: gs
266 type(json_file), target, intent(inout) :: params
267 type(json_file), target, intent(inout) :: numerics_params
268 type(user_t), target, intent(in) :: user
269 type(chkp_t), target, intent(inout) :: chkp
270 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
271 type(time_scheme_controller_t), target, intent(in) :: time_scheme
272 type(field_t), target, intent(in) :: rho
273 end subroutine scalar_scheme_factory
274 end interface
275
276 interface
277
280 module subroutine scalar_scheme_allocator(object, type_name)
281 class(scalar_scheme_t), allocatable, intent(inout) :: object
282 character(len=*), intent(in):: type_name
283 end subroutine scalar_scheme_allocator
284 end interface
285
286 !
287 ! Machinery for injecting user-defined types
288 !
289
293 abstract interface
294 subroutine scalar_scheme_allocate(obj)
295 import scalar_scheme_t
296 class(scalar_scheme_t), allocatable, intent(inout) :: obj
297 end subroutine scalar_scheme_allocate
298 end interface
299
300 interface
301
302 module subroutine register_scalar_scheme(type_name, allocator)
303 character(len=*), intent(in) :: type_name
304 procedure(scalar_scheme_allocate), pointer, intent(in) :: allocator
305 end subroutine register_scalar_scheme
306 end interface
307
308 ! A name-allocator pair for user-defined types. A helper type to define a
309 ! registry of custom allocators.
310 type scalar_scheme_allocator_entry
311 character(len=20) :: type_name
312 procedure(scalar_scheme_allocate), pointer, nopass :: allocator
313 end type scalar_scheme_allocator_entry
314
316 type(scalar_scheme_allocator_entry), allocatable, private :: &
317 scalar_scheme_registry(:)
318
320 integer, private :: scalar_scheme_registry_size = 0
321
322contains
323
332 subroutine scalar_scheme_init(this, msh, c_Xh, gs_Xh, params, scheme, user, &
333 rho)
334 class(scalar_scheme_t), target, intent(inout) :: this
335 type(mesh_t), target, intent(in) :: msh
336 type(coef_t), target, intent(in) :: c_Xh
337 type(gs_t), target, intent(inout) :: gs_Xh
338 type(json_file), target, intent(inout) :: params
339 character(len=*), intent(in) :: scheme
340 type(user_t), target, intent(in) :: user
341 type(field_t), target, intent(in) :: rho
342 ! IO buffer for log output
343 character(len=LOG_SIZE) :: log_buf
344 ! Variables for retrieving json parameters
345 logical :: logical_val
346 real(kind=rp) :: real_val, solver_abstol
347 integer :: integer_val, ierr
348 character(len=:), allocatable :: solver_type, solver_precon
349 type(json_file) :: precon_params
350 type(json_file) :: json_subdict
351 logical :: nut_dependency
352
353 this%u => neko_registry%get_field('u')
354 this%v => neko_registry%get_field('v')
355 this%w => neko_registry%get_field('w')
356 this%rho => rho
357
358 ! Assign a name
359 ! Note that the keyword is added by `scalars_t`, so there is always a
360 ! default.
361 call json_get(params, 'name', this%name)
362
363 ! Set the freeze flag
364 call json_get_or_default(params, 'freeze', this%freeze, .false.)
365
366 call neko_log%section('Scalar')
367 call json_get(params, 'solver.type', solver_type)
368 call json_get(params, 'solver.preconditioner.type', &
369 solver_precon)
370 call json_get(params, 'solver.preconditioner', precon_params)
371 call json_get_or_lookup(params, 'solver.absolute_tolerance', &
372 solver_abstol)
373
374 call json_get_or_lookup_or_default(params, &
375 'solver.projection_space_size', &
376 this%projection_dim, 0)
377 call json_get_or_lookup_or_default(params, &
378 'solver.projection_hold_steps', &
379 this%projection_activ_step, 5)
380
381
382 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
383 call neko_log%message(log_buf)
384 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
385 call neko_log%message(log_buf)
386 call neko_log%message('Ksp scalar : ('// trim(solver_type) // &
387 ', ' // trim(solver_precon) // ')')
388 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
389 call neko_log%message(log_buf)
390
391 this%Xh => this%u%Xh
392 this%dm_Xh => this%u%dof
393 this%params => params
394 this%msh => msh
395
396 call neko_registry%add_field(this%dm_Xh, this%name, &
397 ignore_existing = .true.)
398
399 this%s => neko_registry%get_field(this%name)
400
401 call this%slag%init(this%s, 2)
402
403 this%gs_Xh => gs_xh
404 this%c_Xh => c_xh
405
406 !
407 ! Material properties
408 !
409 call this%set_material_properties(params, user)
410
411
412 !
413 ! Turbulence modelling
414 !
415 this%alphat_field_name = ""
416 this%nut_field_name = ""
417 if (params%valid_path('alphat')) then
418 call json_get(this%params, 'alphat', json_subdict)
419 call json_get(json_subdict, 'nut_dependency', nut_dependency)
420 if (nut_dependency) then
421 call json_get_or_lookup(json_subdict, 'Pr_t', this%pr_turb)
422 call json_get(json_subdict, 'nut_field', this%nut_field_name)
423 else
424 call json_get(json_subdict, 'alphat_field', this%alphat_field_name)
425 end if
426 end if
427
428 !
429 ! Setup right-hand side field.
430 !
431 allocate(this%f_Xh)
432 call this%f_Xh%init(this%dm_Xh, fld_name = "scalar_rhs")
433
434 ! Initialize the source term
435 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
436 call this%source_term%add(params, 'source_terms')
437
438 ! todo parameter file ksp tol should be added
439 call json_get_or_lookup_or_default(params, &
440 'solver.max_iterations', &
441 integer_val, ksp_max_iter)
442 call json_get_or_default(params, &
443 'solver.monitor', &
444 logical_val, .false.)
445 call scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
446 solver_type, integer_val, solver_abstol, logical_val)
447 call scalar_scheme_precon_factory(this%pc, this%ksp, &
448 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
449 solver_precon, precon_params)
450
451 call neko_log%end_section()
452
453 end subroutine scalar_scheme_init
454
455
457 subroutine scalar_scheme_free(this)
458 class(scalar_scheme_t), intent(inout) :: this
459 class(bc_t), pointer :: bc
460 integer :: i
461
462 bc => null()
463
464 nullify(this%Xh)
465 nullify(this%dm_Xh)
466 nullify(this%gs_Xh)
467 nullify(this%c_Xh)
468 nullify(this%params)
469
470 if (allocated(this%ksp)) then
471 call this%ksp%free()
472 deallocate(this%ksp)
473 end if
474
475 if (allocated(this%pc)) then
476 call precon_destroy(this%pc)
477 deallocate(this%pc)
478 end if
479
480 if (allocated(this%name)) then
481 deallocate(this%name)
482 end if
483
484 call this%source_term%free()
485
486 if (associated(this%f_Xh)) then
487 call this%f_Xh%free()
488 deallocate(this%f_Xh)
489 end if
490
491 do i = 1, this%bcs%size()
492 bc => this%bcs%get(i)
493 if (associated(bc)) then
494 call bc%free()
495 deallocate(bc)
496 end if
497 end do
498
499 call this%bcs%free()
500 call this%slag%free()
501 call this%material_properties%free()
502
503 if (allocated(this%nut_field_name)) then
504 deallocate(this%nut_field_name)
505 end if
506
507 if (allocated(this%alphat_field_name)) then
508 deallocate(this%alphat_field_name)
509 end if
510
511 nullify(this%cp)
512 nullify(this%lambda)
513 nullify(this%lambda_tot)
514 nullify(bc)
515
516 end subroutine scalar_scheme_free
517
520 subroutine scalar_scheme_validate(this)
521 class(scalar_scheme_t), target, intent(inout) :: this
522
523 if ( (.not. allocated(this%u%x)) .or. &
524 (.not. allocated(this%v%x)) .or. &
525 (.not. allocated(this%w%x)) .or. &
526 (.not. allocated(this%s%x))) then
527 call neko_error('Fields are not allocated')
528 end if
529
530 if (.not. allocated(this%ksp)) then
531 call neko_error('No Krylov solver for velocity defined')
532 end if
533
534 if (.not. associated(this%Xh)) then
535 call neko_error('No function space defined')
536 end if
537
538 if (.not. associated(this%dm_Xh)) then
539 call neko_error('No dofmap defined')
540 end if
541
542 if (.not. associated(this%c_Xh)) then
543 call neko_error('No coefficients defined')
544 end if
545
546 if (.not. associated(this%f_Xh)) then
547 call neko_error('No rhs allocated')
548 end if
549
550 if (.not. associated(this%params)) then
551 call neko_error('No parameters defined')
552 end if
553
554 if (.not. associated(this%rho)) then
555 call neko_error('No density field defined')
556 end if
557
558 end subroutine scalar_scheme_validate
559
562 subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
563 abstol, monitor)
564 class(ksp_t), allocatable, target, intent(inout) :: ksp
565 integer, intent(in), value :: n
566 integer, intent(in) :: max_iter
567 character(len=*), intent(in) :: solver
568 real(kind=rp) :: abstol
569 logical, intent(in) :: monitor
570
571 call krylov_solver_factory(ksp, n, solver, max_iter, &
572 abstol, monitor = monitor)
573
574 end subroutine scalar_scheme_solver_factory
575
577 subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
578 pctype, pcparams)
579 class(pc_t), allocatable, target, intent(inout) :: pc
580 class(ksp_t), target, intent(inout) :: ksp
581 type(coef_t), target, intent(in) :: coef
582 type(dofmap_t), target, intent(in) :: dof
583 type(gs_t), target, intent(inout) :: gs
584 type(bc_list_t), target, intent(inout) :: bclst
585 character(len=*) :: pctype
586 type(json_file), intent(inout) :: pcparams
587
588 call precon_factory(pc, pctype)
589
590 select type (pcp => pc)
591 type is (jacobi_t)
592 call pcp%init(coef, dof, gs)
593 type is (sx_jacobi_t)
594 call pcp%init(coef, dof, gs)
595 type is (device_jacobi_t)
596 call pcp%init(coef, dof, gs)
597 type is (hsmg_t)
598 call pcp%init(coef, bclst, pcparams)
599 end select
600
601 call ksp%set_pc(pc)
602
603 end subroutine scalar_scheme_precon_factory
604
609 subroutine scalar_scheme_update_material_properties(this, time)
610 class(scalar_scheme_t), intent(inout) :: this
611 type(time_state_t), intent(in) :: time
612 type(field_t), pointer :: nut, alphat
613 integer :: index
614 ! Factor to transform nu_t to lambda_t
615 type(field_t), pointer :: lambda_factor
616
617 call this%user_material_properties(this%name, this%material_properties, &
618 time)
619
620 ! factor = rho * cp / pr_turb
621 if (len_trim(this%nut_field_name) .gt. 0 &
622 .and. len_trim(this%alphat_field_name) .eq. 0 ) then
623 nut => neko_registry%get_field(this%nut_field_name)
624
625 ! lambda_tot = lambda + rho * cp * nut / pr_turb
626 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
627 call field_cmult2(lambda_factor, nut, 1.0_rp / this%pr_turb)
628 call field_col2(lambda_factor, this%cp)
629 call field_col2(lambda_factor, this%rho)
630 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
631 call neko_scratch_registry%relinquish_field(index)
632
633 else if (len_trim(this%alphat_field_name) .gt. 0 &
634 .and. len_trim(this%nut_field_name) .eq. 0 ) then
635 alphat => neko_registry%get_field(this%alphat_field_name)
636
637 ! lambda_tot = lambda + rho * cp * alphat
638 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
639 call field_col3(lambda_factor, this%cp, alphat)
640 call field_col2(lambda_factor, this%rho)
641 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
642 call neko_scratch_registry%relinquish_field(index)
643
644 else if (len_trim(this%alphat_field_name) .gt. 0 &
645 .and. len_trim(this%nut_field_name) .gt. 0 ) then
646 call neko_error("Conflicting definition of eddy diffusivity " // &
647 "for the scalar equation")
648 end if
649
650 ! Since cp is a fields and we use the %x(1,1,1,1) of the
651 ! host array data to pass constant material properties
652 ! to some routines, we need to make sure that the host
653 ! values are also filled
654 if (neko_bcknd_device .eq. 1) then
655 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
656 device_to_host, sync = .false.)
657 end if
658
659 end subroutine scalar_scheme_update_material_properties
660
664 subroutine scalar_scheme_set_material_properties(this, params, user)
665 class(scalar_scheme_t), intent(inout) :: this
666 type(json_file), intent(inout) :: params
667 type(user_t), target, intent(in) :: user
668 character(len=LOG_SIZE) :: log_buf
669 ! A local pointer that is needed to make Intel happy
670 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
671 real(kind=rp) :: const_cp, const_lambda
672 ! Dummy time state set to 0
673 type(time_state_t) :: time
674
675 dummy_mp_ptr => dummy_user_material_properties
676
677 ! Fill lambda field with the physical value
678
679 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda")
680 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda_tot")
681 call neko_registry%add_field(this%dm_Xh, this%name // "_cp")
682 this%lambda => neko_registry%get_field(this%name // "_lambda")
683 this%lambda_tot => neko_registry%get_field(this%name // "_lambda_tot")
684 this%cp => neko_registry%get_field(this%name // "_cp")
685
686 call this%material_properties%init(2)
687 call this%material_properties%assign(1, this%cp)
688 call this%material_properties%assign(2, this%lambda)
689
690 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
691
692 write(log_buf, '(A)') "Material properties must be set in the user " // &
693 "file!"
694 call neko_log%message(log_buf)
695 this%user_material_properties => user%material_properties
696
697 call user%material_properties(this%name, this%material_properties, time)
698 else
699 this%user_material_properties => dummy_user_material_properties
700 if (params%valid_path('Pe') .and. &
701 (params%valid_path('lambda') .or. &
702 params%valid_path('cp'))) then
703 call neko_error("To set the material properties for the scalar, " // &
704 "either provide Pe OR lambda and cp in the case file.")
705 ! Non-dimensional case
706 else if (params%valid_path('Pe')) then
707 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
708 ' input.'
709 call neko_log%message(log_buf, lvl = neko_log_verbose)
710 write(log_buf, '(A)') 'Specific heat capacity will be set to 1, '
711 call neko_log%message(log_buf, lvl = neko_log_verbose)
712 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
713 call neko_log%message(log_buf, lvl = neko_log_verbose)
714
715 ! Read Pe into lambda for further manipulation.
716 call json_get_or_lookup(params, 'Pe', const_lambda)
717 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
718 call neko_log%message(log_buf)
719
720 ! Set cp and rho to 1 since the setup is non-dimensional.
721 const_cp = 1.0_rp
722 ! Invert the Pe to get conductivity
723 const_lambda = 1.0_rp/const_lambda
724 ! Dimensional case
725 else
726 call json_get_or_lookup(params, 'lambda', const_lambda)
727 call json_get_or_lookup(params, 'cp', const_cp)
728 end if
729 end if
730 ! We need to fill the fields based on the parsed const values
731 ! if the user routine is not used.
732 if (associated(user%material_properties, dummy_mp_ptr)) then
733 ! Fill mu and rho field with the physical value
734 call field_cfill(this%lambda, const_lambda)
735 call field_cfill(this%cp, const_cp)
736
737 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
738 call neko_log%message(log_buf)
739 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
740 call neko_log%message(log_buf)
741 end if
742
743 ! Copy over material property to the total one
744 call field_copy(this%lambda_tot, this%lambda)
745
746 ! Since cp is a field and we use the %x(1,1,1,1) of the
747 ! host array data to pass constant material properties
748 ! to some routines, we need to make sure that the host
749 ! values are also filled
750 if (neko_bcknd_device .eq. 1) then
751 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
752 device_to_host, sync = .false.)
753 end if
754 end subroutine scalar_scheme_set_material_properties
755
756 ! ========================================================================== !
757 ! Scalar scheme wrapper type methods
758
760 subroutine scalar_scheme_wrapper_init(this, msh, coef, gs, params, &
761 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
762 class(scalar_scheme_wrapper_t), intent(inout) :: this
763 type(mesh_t), target, intent(in) :: msh
764 type(coef_t), target, intent(in) :: coef
765 type(gs_t), target, intent(inout) :: gs
766 type(json_file), target, intent(inout) :: params
767 type(json_file), target, intent(inout) :: numerics_params
768 type(user_t), target, intent(in) :: user
769 type(chkp_t), target, intent(inout) :: chkp
770 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
771 type(time_scheme_controller_t), target, intent(in) :: time_scheme
772 type(field_t), target, intent(in) :: rho
773
774 call this%free()
775 call scalar_scheme_factory(this%scalar, msh, coef, gs, params, &
776 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
777
778 end subroutine scalar_scheme_wrapper_init
779
781 subroutine scalar_scheme_wrapper_free(this)
782 class(scalar_scheme_wrapper_t), intent(inout) :: this
783
784 if (allocated(this%scalar)) then
785 call this%scalar%free()
786 deallocate(this%scalar)
787 end if
788
789 end subroutine scalar_scheme_wrapper_free
790
795 subroutine scalar_scheme_wrapper_move_from(this, other)
796 class(scalar_scheme_wrapper_t), intent(inout) :: this
797 class(scalar_scheme_wrapper_t), intent(inout) :: other
798
799 ! Move the pointer
800 call move_alloc(other%scalar, this%scalar)
801
802 end subroutine scalar_scheme_wrapper_move_from
803
806 function scalar_scheme_wrapper_is_allocated(this) result(is_alloc)
807 class(scalar_scheme_wrapper_t), intent(in) :: this
808 logical :: is_alloc
809 is_alloc = allocated(this%scalar)
810 end function scalar_scheme_wrapper_is_allocated
811
812end module scalar_scheme
Copy data between host and device (or device and device)
Definition device.F90:72
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.
Abstract interface to dealocate a scalar formulation.
Abstract interface to initialize a scalar formulation.
Abstract interface to restart a scalar formulation.
Abstract interface to compute a time-step.
Abstract interface for setting material properties.
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Jacobi preconditioner accelerator backend.
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public device_to_host
Definition device.F90:48
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
subroutine, public field_col2(a, b, n)
Vector multiplication .
subroutine, public field_cmult2(a, b, c, n)
Multiplication by constant c .
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public field_copy(a, b, n)
Copy a vector .
subroutine, public field_add3(a, b, c, n)
Vector addition .
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Gather-scatter.
Krylov preconditioner.
Definition pc_hsmg.f90:61
Jacobi preconditioner.
Definition pc_jacobi.f90:34
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
integer, parameter, public ksp_max_iter
Maximum number of iters.
Definition krylov.f90:51
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_verbose
Verbose log level.
Definition log.f90:86
type(log_t), public neko_log
Global log stream.
Definition log.f90:80
integer, parameter, public log_size
Definition log.f90:46
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Krylov preconditioner.
Definition precon.f90:34
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
Contains the scalar_scheme_t type.
subroutine scalar_scheme_init(this, msh, c_xh, gs_xh, params, scheme, user, rho)
Initialize all related components of the current scheme.
logical function scalar_scheme_wrapper_is_allocated(this)
Return allocation status.
subroutine scalar_scheme_wrapper_free(this)
Destructor. Just deallocates the pointer.
subroutine scalar_scheme_free(this)
Deallocate a scalar formulation.
subroutine scalar_scheme_validate(this)
Validate that all fields, solvers etc necessary for performing time-stepping are defined.
subroutine scalar_scheme_update_material_properties(this, time)
Call user material properties routine and update the values of lambda if necessary.
subroutine scalar_scheme_set_material_properties(this, params, user)
Set lamdba and cp.
subroutine scalar_scheme_wrapper_init(this, msh, coef, gs, params, numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
Constructor. Initializes the object.
subroutine scalar_scheme_wrapper_move_from(this, other)
Move assignment operator for the wrapper, needed for storing schemes in lists and arrays.
Implements the scalar_source_term_t type.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Defines a function space.
Definition space.f90:34
Jacobi preconditioner SX-Aurora backend.
Compound scheme for the advection and diffusion operators in a transport equation.
Base class for time integration schemes.
Module with things related to the simulation time.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
subroutine, public dummy_user_material_properties(scheme_name, properties, time)
Utilities.
Definition utils.f90:35
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:49
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Defines a jacobi preconditioner.
field_list_t, To be able to group fields together
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
Defines a jacobi preconditioner.
Definition pc_jacobi.f90:45
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Base type for a scalar advection-diffusion solver.
A helper type that is needed to have an array of polymorphic objects.
Wrapper contaning and executing the scalar source terms.
The function space for the SEM solution fields.
Definition space.f90:63
Defines a jacobi preconditioner for SX-Aurora.
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
A struct that contains all info about the time, expand as needed.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...