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 precon, only : pc_t, precon_factory, precon_destroy
51 use mesh, only : mesh_t
54 use registry, only : neko_registry
57 use json_module, only : json_file
60 use utils, only : neko_error
68 use time_state, only : time_state_t
70 implicit none
71
73 type, abstract :: scalar_scheme_t
75 character(len=:), allocatable :: name
77 type(field_t), pointer :: u
79 type(field_t), pointer :: v
81 type(field_t), pointer :: w
83 type(field_t), pointer :: s
85 type(field_series_t) :: slag
87 type(space_t), pointer :: xh
89 type(dofmap_t), pointer :: dm_xh
91 type(gs_t), pointer :: gs_xh
93 type(coef_t), pointer :: c_xh
95 type(field_t), pointer :: f_xh => null()
99 class(ksp_t), allocatable :: ksp
101 integer :: ksp_maxiter
103 integer :: projection_dim
104
105 integer :: projection_activ_step
107 class(pc_t), allocatable :: pc
109 type(bc_list_t) :: bcs
111 type(json_file), pointer :: params
113 type(mesh_t), pointer :: msh => null()
115 type(chkp_t), pointer :: chkp => null()
117 character(len=:), allocatable :: nut_field_name
119 character(len=:), allocatable :: alphat_field_name
121 type(field_t), pointer :: rho => null()
123 type(field_t), pointer :: lambda => null()
125 type(field_t), pointer :: cp => null()
127 type(field_t), pointer :: lambda_tot => null()
129 real(kind=rp) :: pr_turb
131 type(field_list_t) :: material_properties
132 procedure(user_material_properties_intf), nopass, pointer :: &
133 user_material_properties => null()
135 logical :: freeze = .false.
136 contains
138 procedure, pass(this) :: scheme_init => scalar_scheme_init
140 procedure, pass(this) :: scheme_free => scalar_scheme_free
142 procedure, pass(this) :: validate => scalar_scheme_validate
144 procedure, pass(this) :: set_material_properties => &
147 procedure, pass(this) :: update_material_properties => &
150 procedure(scalar_scheme_init_intrf), pass(this), deferred :: init
152 procedure(scalar_scheme_free_intrf), pass(this), deferred :: free
154 procedure(scalar_scheme_step_intrf), pass(this), deferred :: step
156 procedure(scalar_scheme_restart_intrf), pass(this), deferred :: restart
157 end type scalar_scheme_t
158
160 abstract interface
161 subroutine scalar_scheme_init_intrf(this, msh, coef, gs, params, &
162 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
163 import scalar_scheme_t
164 import json_file
165 import coef_t
166 import gs_t
167 import mesh_t
168 import user_t
169 import field_series_t, field_t
171 import rp
172 import chkp_t
173 class(scalar_scheme_t), target, intent(inout) :: this
174 type(mesh_t), target, intent(in) :: msh
175 type(coef_t), target, intent(in) :: coef
176 type(gs_t), target, intent(inout) :: gs
177 type(json_file), target, intent(inout) :: params
178 type(json_file), target, intent(inout) :: numerics_params
179 type(user_t), target, intent(in) :: user
180 type(chkp_t), target, intent(inout) :: chkp
181 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
182 type(time_scheme_controller_t), target, intent(in) :: time_scheme
183 type(field_t), target, intent(in) :: rho
184 end subroutine scalar_scheme_init_intrf
185 end interface
186
188 abstract interface
189 subroutine scalar_scheme_restart_intrf(this, chkp)
190 import scalar_scheme_t
191 import chkp_t
192 import rp
193 class(scalar_scheme_t), target, intent(inout) :: this
194 type(chkp_t), intent(inout) :: chkp
195 end subroutine scalar_scheme_restart_intrf
196 end interface
197
199 abstract interface
201 import scalar_scheme_t
202 class(scalar_scheme_t), intent(inout) :: this
203 end subroutine scalar_scheme_free_intrf
204 end interface
205
207 abstract interface
208 subroutine scalar_scheme_step_intrf(this, time, ext_bdf, dt_controller, &
209 ksp_results)
210 import scalar_scheme_t
211 import time_state_t
214 import ksp_monitor_t
215 class(scalar_scheme_t), intent(inout) :: this
216 type(time_state_t), intent(in) :: time
217 type(time_scheme_controller_t), intent(in) :: ext_bdf
218 type(time_step_controller_t), intent(in) :: dt_controller
219 type(ksp_monitor_t), intent(inout) :: ksp_results
220 end subroutine scalar_scheme_step_intrf
221 end interface
222
223 ! ========================================================================== !
224 ! Helper functions and types for scalar_scheme_t
225 ! ========================================================================== !
226
229 class(scalar_scheme_t), allocatable :: scalar
230 contains
232 procedure, pass(this) :: init => scalar_scheme_wrapper_init
234 procedure, pass(this) :: free => scalar_scheme_wrapper_free
237 procedure, pass(this) :: move_from => &
240 procedure, pass(this) :: is_allocated => &
243
244
245 interface
246
259 module subroutine scalar_scheme_factory(object, msh, coef, gs, params, &
260 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
261 class(scalar_scheme_t), allocatable, intent(inout) :: object
262 type(mesh_t), target, intent(in) :: msh
263 type(coef_t), target, intent(in) :: coef
264 type(gs_t), target, intent(inout) :: gs
265 type(json_file), target, intent(inout) :: params
266 type(json_file), target, intent(inout) :: numerics_params
267 type(user_t), target, intent(in) :: user
268 type(chkp_t), target, intent(inout) :: chkp
269 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
270 type(time_scheme_controller_t), target, intent(in) :: time_scheme
271 type(field_t), target, intent(in) :: rho
272 end subroutine scalar_scheme_factory
273 end interface
274
275 interface
276
279 module subroutine scalar_scheme_allocator(object, type_name)
280 class(scalar_scheme_t), allocatable, intent(inout) :: object
281 character(len=*), intent(in):: type_name
282 end subroutine scalar_scheme_allocator
283 end interface
284
285 !
286 ! Machinery for injecting user-defined types
287 !
288
292 abstract interface
293 subroutine scalar_scheme_allocate(obj)
294 import scalar_scheme_t
295 class(scalar_scheme_t), allocatable, intent(inout) :: obj
296 end subroutine scalar_scheme_allocate
297 end interface
298
299 interface
300
301 module subroutine register_scalar_scheme(type_name, allocator)
302 character(len=*), intent(in) :: type_name
303 procedure(scalar_scheme_allocate), pointer, intent(in) :: allocator
304 end subroutine register_scalar_scheme
305 end interface
306
307 ! A name-allocator pair for user-defined types. A helper type to define a
308 ! registry of custom allocators.
309 type scalar_scheme_allocator_entry
310 character(len=20) :: type_name
311 procedure(scalar_scheme_allocate), pointer, nopass :: allocator
312 end type scalar_scheme_allocator_entry
313
315 type(scalar_scheme_allocator_entry), allocatable, private :: &
316 scalar_scheme_registry(:)
317
319 integer, private :: scalar_scheme_registry_size = 0
320
321contains
322
331 subroutine scalar_scheme_init(this, msh, c_Xh, gs_Xh, params, scheme, user, &
332 rho)
333 class(scalar_scheme_t), target, intent(inout) :: this
334 type(mesh_t), target, intent(in) :: msh
335 type(coef_t), target, intent(in) :: c_Xh
336 type(gs_t), target, intent(inout) :: gs_Xh
337 type(json_file), target, intent(inout) :: params
338 character(len=*), intent(in) :: scheme
339 type(user_t), target, intent(in) :: user
340 type(field_t), target, intent(in) :: rho
341 ! IO buffer for log output
342 character(len=LOG_SIZE) :: log_buf
343 ! Variables for retrieving json parameters
344 logical :: logical_val
345 real(kind=rp) :: real_val, solver_abstol
346 integer :: integer_val, ierr
347 character(len=:), allocatable :: solver_type, solver_precon
348 type(json_file) :: precon_params
349 type(json_file) :: json_subdict
350 logical :: nut_dependency
351
352 this%u => neko_registry%get_field('u')
353 this%v => neko_registry%get_field('v')
354 this%w => neko_registry%get_field('w')
355 this%rho => rho
356
357 ! Assign a name
358 ! Note that the keyword is added by `scalars_t`, so there is always a
359 ! default.
360 call json_get(params, 'name', this%name)
361
362 ! Set the freeze flag
363 call json_get_or_default(params, 'freeze', this%freeze, .false.)
364
365 call neko_log%section('Scalar')
366 call json_get(params, 'solver.type', solver_type)
367 call json_get(params, 'solver.preconditioner.type', &
368 solver_precon)
369 call json_get(params, 'solver.preconditioner', precon_params)
370 call json_get_or_lookup(params, 'solver.absolute_tolerance', &
371 solver_abstol)
372
373 call json_get_or_lookup_or_default(params, &
374 'solver.projection_space_size', &
375 this%projection_dim, 0)
376 call json_get_or_lookup_or_default(params, &
377 'solver.projection_hold_steps', &
378 this%projection_activ_step, 5)
379
380
381 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
382 call neko_log%message(log_buf)
383 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
384 call neko_log%message(log_buf)
385 call neko_log%message('Ksp scalar : ('// trim(solver_type) // &
386 ', ' // trim(solver_precon) // ')')
387 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
388 call neko_log%message(log_buf)
389
390 this%Xh => this%u%Xh
391 this%dm_Xh => this%u%dof
392 this%params => params
393 this%msh => msh
394
395 call neko_registry%add_field(this%dm_Xh, this%name, &
396 ignore_existing = .true.)
397
398 this%s => neko_registry%get_field(this%name)
399
400 call this%slag%init(this%s, 2)
401
402 this%gs_Xh => gs_xh
403 this%c_Xh => c_xh
404
405 !
406 ! Material properties
407 !
408 call this%set_material_properties(params, user)
409
410
411 !
412 ! Turbulence modelling
413 !
414 this%alphat_field_name = ""
415 this%nut_field_name = ""
416 if (params%valid_path('alphat')) then
417 call json_get(this%params, 'alphat', json_subdict)
418 call json_get(json_subdict, 'nut_dependency', nut_dependency)
419 if (nut_dependency) then
420 call json_get(json_subdict, 'Pr_t', this%pr_turb)
421 call json_get(json_subdict, 'nut_field', this%nut_field_name)
422 else
423 call json_get(json_subdict, 'alphat_field', this%alphat_field_name)
424 end if
425 end if
426
427 !
428 ! Setup right-hand side field.
429 !
430 allocate(this%f_Xh)
431 call this%f_Xh%init(this%dm_Xh, fld_name = "scalar_rhs")
432
433 ! Initialize the source term
434 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
435 call this%source_term%add(params, 'source_terms')
436
437 ! todo parameter file ksp tol should be added
438 call json_get_or_lookup_or_default(params, &
439 'solver.max_iterations', &
440 integer_val, ksp_max_iter)
441 call json_get_or_default(params, &
442 'solver.monitor', &
443 logical_val, .false.)
444 call scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
445 solver_type, integer_val, solver_abstol, logical_val)
446 call scalar_scheme_precon_factory(this%pc, this%ksp, &
447 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
448 solver_precon, precon_params)
449
450 call neko_log%end_section()
451
452 end subroutine scalar_scheme_init
453
454
456 subroutine scalar_scheme_free(this)
457 class(scalar_scheme_t), intent(inout) :: this
458
459 nullify(this%Xh)
460 nullify(this%dm_Xh)
461 nullify(this%gs_Xh)
462 nullify(this%c_Xh)
463 nullify(this%params)
464
465 if (allocated(this%ksp)) then
466 call this%ksp%free()
467 deallocate(this%ksp)
468 end if
469
470 if (allocated(this%pc)) then
471 call precon_destroy(this%pc)
472 deallocate(this%pc)
473 end if
474
475 if (allocated(this%name)) then
476 deallocate(this%name)
477 end if
478
479 call this%source_term%free()
480
481 call this%bcs%free()
482 call this%slag%free()
483
484 nullify(this%cp)
485 nullify(this%lambda)
486 nullify(this%lambda_tot)
487
488 end subroutine scalar_scheme_free
489
492 subroutine scalar_scheme_validate(this)
493 class(scalar_scheme_t), target, intent(inout) :: this
494
495 if ( (.not. allocated(this%u%x)) .or. &
496 (.not. allocated(this%v%x)) .or. &
497 (.not. allocated(this%w%x)) .or. &
498 (.not. allocated(this%s%x))) then
499 call neko_error('Fields are not allocated')
500 end if
501
502 if (.not. allocated(this%ksp)) then
503 call neko_error('No Krylov solver for velocity defined')
504 end if
505
506 if (.not. associated(this%Xh)) then
507 call neko_error('No function space defined')
508 end if
509
510 if (.not. associated(this%dm_Xh)) then
511 call neko_error('No dofmap defined')
512 end if
513
514 if (.not. associated(this%c_Xh)) then
515 call neko_error('No coefficients defined')
516 end if
517
518 if (.not. associated(this%f_Xh)) then
519 call neko_error('No rhs allocated')
520 end if
521
522 if (.not. associated(this%params)) then
523 call neko_error('No parameters defined')
524 end if
525
526 if (.not. associated(this%rho)) then
527 call neko_error('No density field defined')
528 end if
529
530 end subroutine scalar_scheme_validate
531
534 subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
535 abstol, monitor)
536 class(ksp_t), allocatable, target, intent(inout) :: ksp
537 integer, intent(in), value :: n
538 integer, intent(in) :: max_iter
539 character(len=*), intent(in) :: solver
540 real(kind=rp) :: abstol
541 logical, intent(in) :: monitor
542
543 call krylov_solver_factory(ksp, n, solver, max_iter, &
544 abstol, monitor = monitor)
545
546 end subroutine scalar_scheme_solver_factory
547
549 subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
550 pctype, pcparams)
551 class(pc_t), allocatable, target, intent(inout) :: pc
552 class(ksp_t), target, intent(inout) :: ksp
553 type(coef_t), target, intent(in) :: coef
554 type(dofmap_t), target, intent(in) :: dof
555 type(gs_t), target, intent(inout) :: gs
556 type(bc_list_t), target, intent(inout) :: bclst
557 character(len=*) :: pctype
558 type(json_file), intent(inout) :: pcparams
559
560 call precon_factory(pc, pctype)
561
562 select type (pcp => pc)
563 type is (jacobi_t)
564 call pcp%init(coef, dof, gs)
565 type is (sx_jacobi_t)
566 call pcp%init(coef, dof, gs)
567 type is (device_jacobi_t)
568 call pcp%init(coef, dof, gs)
569 type is (hsmg_t)
570 call pcp%init(coef, bclst, pcparams)
571 end select
572
573 call ksp%set_pc(pc)
574
575 end subroutine scalar_scheme_precon_factory
576
581 subroutine scalar_scheme_update_material_properties(this, time)
582 class(scalar_scheme_t), intent(inout) :: this
583 type(time_state_t), intent(in) :: time
584 type(field_t), pointer :: nut, alphat
585 integer :: index
586 ! Factor to transform nu_t to lambda_t
587 type(field_t), pointer :: lambda_factor
588
589 call this%user_material_properties(this%name, this%material_properties, &
590 time)
591
592 ! factor = rho * cp / pr_turb
593 if (len_trim(this%nut_field_name) .gt. 0 &
594 .and. len_trim(this%alphat_field_name) .eq. 0 ) then
595 nut => neko_registry%get_field(this%nut_field_name)
596
597 ! lambda_tot = lambda + rho * cp * nut / pr_turb
598 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
599 call field_cmult2(lambda_factor, nut, 1.0_rp / this%pr_turb)
600 call field_col2(lambda_factor, this%cp)
601 call field_col2(lambda_factor, this%rho)
602 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
603 call neko_scratch_registry%relinquish_field(index)
604
605 else if (len_trim(this%alphat_field_name) .gt. 0 &
606 .and. len_trim(this%nut_field_name) .eq. 0 ) then
607 alphat => neko_registry%get_field(this%alphat_field_name)
608
609 ! lambda_tot = lambda + rho * cp * alphat
610 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
611 call field_col3(lambda_factor, this%cp, alphat)
612 call field_col2(lambda_factor, this%rho)
613 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
614 call neko_scratch_registry%relinquish_field(index)
615
616 else if (len_trim(this%alphat_field_name) .gt. 0 &
617 .and. len_trim(this%nut_field_name) .gt. 0 ) then
618 call neko_error("Conflicting definition of eddy diffusivity " // &
619 "for the scalar equation")
620 end if
621
622 ! Since cp is a fields and we use the %x(1,1,1,1) of the
623 ! host array data to pass constant material properties
624 ! to some routines, we need to make sure that the host
625 ! values are also filled
626 if (neko_bcknd_device .eq. 1) then
627 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
628 device_to_host, sync = .false.)
629 end if
630
631 end subroutine scalar_scheme_update_material_properties
632
636 subroutine scalar_scheme_set_material_properties(this, params, user)
637 class(scalar_scheme_t), intent(inout) :: this
638 type(json_file), intent(inout) :: params
639 type(user_t), target, intent(in) :: user
640 character(len=LOG_SIZE) :: log_buf
641 ! A local pointer that is needed to make Intel happy
642 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
643 real(kind=rp) :: const_cp, const_lambda
644 ! Dummy time state set to 0
645 type(time_state_t) :: time
646
647 dummy_mp_ptr => dummy_user_material_properties
648
649 ! Fill lambda field with the physical value
650
651 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda")
652 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda_tot")
653 call neko_registry%add_field(this%dm_Xh, this%name // "_cp")
654 this%lambda => neko_registry%get_field(this%name // "_lambda")
655 this%lambda_tot => neko_registry%get_field(this%name // "_lambda_tot")
656 this%cp => neko_registry%get_field(this%name // "_cp")
657
658 call this%material_properties%init(2)
659 call this%material_properties%assign(1, this%cp)
660 call this%material_properties%assign(2, this%lambda)
661
662 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
663
664 write(log_buf, '(A)') "Material properties must be set in the user " // &
665 "file!"
666 call neko_log%message(log_buf)
667 this%user_material_properties => user%material_properties
668
669 call user%material_properties(this%name, this%material_properties, time)
670 else
671 this%user_material_properties => dummy_user_material_properties
672 if (params%valid_path('Pe') .and. &
673 (params%valid_path('lambda') .or. &
674 params%valid_path('cp'))) then
675 call neko_error("To set the material properties for the scalar, " // &
676 "either provide Pe OR lambda and cp in the case file.")
677 ! Non-dimensional case
678 else if (params%valid_path('Pe')) then
679 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
680 ' input.'
681 call neko_log%message(log_buf, lvl = neko_log_verbose)
682 write(log_buf, '(A)') 'Specific heat capacity will be set to 1, '
683 call neko_log%message(log_buf, lvl = neko_log_verbose)
684 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
685 call neko_log%message(log_buf, lvl = neko_log_verbose)
686
687 ! Read Pe into lambda for further manipulation.
688 call json_get_or_lookup(params, 'Pe', const_lambda)
689 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
690 call neko_log%message(log_buf)
691
692 ! Set cp and rho to 1 since the setup is non-dimensional.
693 const_cp = 1.0_rp
694 ! Invert the Pe to get conductivity
695 const_lambda = 1.0_rp/const_lambda
696 ! Dimensional case
697 else
698 call json_get_or_lookup(params, 'lambda', const_lambda)
699 call json_get_or_lookup(params, 'cp', const_cp)
700 end if
701 end if
702 ! We need to fill the fields based on the parsed const values
703 ! if the user routine is not used.
704 if (associated(user%material_properties, dummy_mp_ptr)) then
705 ! Fill mu and rho field with the physical value
706 call field_cfill(this%lambda, const_lambda)
707 call field_cfill(this%cp, const_cp)
708
709 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
710 call neko_log%message(log_buf)
711 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
712 call neko_log%message(log_buf)
713 end if
714
715 ! Copy over material property to the total one
716 call field_copy(this%lambda_tot, this%lambda)
717
718 ! Since cp is a field and we use the %x(1,1,1,1) of the
719 ! host array data to pass constant material properties
720 ! to some routines, we need to make sure that the host
721 ! values are also filled
722 if (neko_bcknd_device .eq. 1) then
723 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
724 device_to_host, sync = .false.)
725 end if
726 end subroutine scalar_scheme_set_material_properties
727
728 ! ========================================================================== !
729 ! Scalar scheme wrapper type methods
730
732 subroutine scalar_scheme_wrapper_init(this, msh, coef, gs, params, &
733 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
734 class(scalar_scheme_wrapper_t), intent(inout) :: this
735 type(mesh_t), target, intent(in) :: msh
736 type(coef_t), target, intent(in) :: coef
737 type(gs_t), target, intent(inout) :: gs
738 type(json_file), target, intent(inout) :: params
739 type(json_file), target, intent(inout) :: numerics_params
740 type(user_t), target, intent(in) :: user
741 type(chkp_t), target, intent(inout) :: chkp
742 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
743 type(time_scheme_controller_t), target, intent(in) :: time_scheme
744 type(field_t), target, intent(in) :: rho
745
746 call this%free()
747 call scalar_scheme_factory(this%scalar, msh, coef, gs, params, &
748 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
749
750 end subroutine scalar_scheme_wrapper_init
751
753 subroutine scalar_scheme_wrapper_free(this)
754 class(scalar_scheme_wrapper_t), intent(inout) :: this
755
756 if (allocated(this%scalar)) then
757 call this%scalar%free()
758 deallocate(this%scalar)
759 end if
760
761 end subroutine scalar_scheme_wrapper_free
762
767 subroutine scalar_scheme_wrapper_move_from(this, other)
768 class(scalar_scheme_wrapper_t), intent(inout) :: this
769 class(scalar_scheme_wrapper_t), intent(inout) :: other
770
771 ! Move the pointer
772 call move_alloc(other%scalar, this%scalar)
773
774 end subroutine scalar_scheme_wrapper_move_from
775
778 function scalar_scheme_wrapper_is_allocated(this) result(is_alloc)
779 class(scalar_scheme_wrapper_t), intent(in) :: this
780 logical :: is_alloc
781 is_alloc = allocated(this%scalar)
782 end function scalar_scheme_wrapper_is_allocated
783
784end module scalar_scheme
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.
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 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:47
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:85
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
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
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:58
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...