Neko 1.0.0
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-2025, 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 dirichlet, only : dirichlet_t
46 use neumann, only : neumann_t
47 use jacobi, only : jacobi_t
49 use sx_jacobi, only : sx_jacobi_t
50 use hsmg, only : hsmg_t
51 use bc, only : bc_t
52 use bc_list, only : bc_list_t
53 use precon, only : pc_t, precon_factory, precon_destroy
56 use facet_zone, only : facet_zone_t
59 use registry, only : neko_registry
61 use json_module, only : json_file
64 use utils, only : neko_error, neko_warning
65 use comm, only: neko_comm
66 use mpi_f08, only : mpi_integer, mpi_sum
69 use math, only : cfill, add2s2
76 use time_state, only : time_state_t
78 implicit none
79
81 type, abstract :: scalar_scheme_t
83 character(len=:), allocatable :: name
85 type(field_t), pointer :: u
87 type(field_t), pointer :: v
89 type(field_t), pointer :: w
91 type(field_t), pointer :: s
93 type(field_series_t) :: slag
95 type(space_t), pointer :: xh
97 type(dofmap_t), pointer :: dm_xh
99 type(gs_t), pointer :: gs_xh
101 type(coef_t), pointer :: c_xh
103 type(field_t), pointer :: f_xh => null()
107 class(ksp_t), allocatable :: ksp
109 integer :: ksp_maxiter
111 integer :: projection_dim
112
113 integer :: projection_activ_step
115 class(pc_t), allocatable :: pc
117 type(bc_list_t) :: bcs
119 type(json_file), pointer :: params
121 type(mesh_t), pointer :: msh => null()
123 type(chkp_t), pointer :: chkp => null()
125 character(len=:), allocatable :: nut_field_name
127 character(len=:), allocatable :: alphat_field_name
129 type(field_t), pointer :: rho => null()
131 type(field_t), pointer :: lambda => null()
133 type(field_t), pointer :: cp => null()
135 type(field_t), pointer :: lambda_tot => null()
137 real(kind=rp) :: pr_turb
139 type(field_list_t) :: material_properties
140 procedure(user_material_properties_intf), nopass, pointer :: &
141 user_material_properties => null()
143 logical :: freeze = .false.
144 contains
146 procedure, pass(this) :: scheme_init => scalar_scheme_init
148 procedure, pass(this) :: scheme_free => scalar_scheme_free
150 procedure, pass(this) :: validate => scalar_scheme_validate
152 procedure, pass(this) :: set_material_properties => &
155 procedure, pass(this) :: update_material_properties => &
158 procedure(scalar_scheme_init_intrf), pass(this), deferred :: init
160 procedure(scalar_scheme_free_intrf), pass(this), deferred :: free
162 procedure(scalar_scheme_step_intrf), pass(this), deferred :: step
164 procedure(scalar_scheme_restart_intrf), pass(this), deferred :: restart
165 end type scalar_scheme_t
166
168 abstract interface
169 subroutine scalar_scheme_init_intrf(this, msh, coef, gs, params, &
170 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
171 import scalar_scheme_t
172 import json_file
173 import coef_t
174 import gs_t
175 import mesh_t
176 import user_t
177 import field_series_t, field_t
179 import rp
180 import chkp_t
181 class(scalar_scheme_t), target, intent(inout) :: this
182 type(mesh_t), target, intent(in) :: msh
183 type(coef_t), target, intent(in) :: coef
184 type(gs_t), target, intent(inout) :: gs
185 type(json_file), target, intent(inout) :: params
186 type(json_file), target, intent(inout) :: numerics_params
187 type(user_t), target, intent(in) :: user
188 type(chkp_t), target, intent(inout) :: chkp
189 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
190 type(time_scheme_controller_t), target, intent(in) :: time_scheme
191 type(field_t), target, intent(in) :: rho
192 end subroutine scalar_scheme_init_intrf
193 end interface
194
196 abstract interface
197 subroutine scalar_scheme_restart_intrf(this, chkp)
198 import scalar_scheme_t
199 import chkp_t
200 import rp
201 class(scalar_scheme_t), target, intent(inout) :: this
202 type(chkp_t), intent(inout) :: chkp
203 end subroutine scalar_scheme_restart_intrf
204 end interface
205
206
208 abstract interface
210 import scalar_scheme_t
211 class(scalar_scheme_t), intent(inout) :: this
212 end subroutine scalar_scheme_free_intrf
213 end interface
214
216 abstract interface
217 subroutine scalar_scheme_step_intrf(this, time, ext_bdf, dt_controller, &
218 ksp_results)
219 import scalar_scheme_t
220 import time_state_t
223 import ksp_monitor_t
224 class(scalar_scheme_t), intent(inout) :: this
225 type(time_state_t), intent(in) :: time
226 type(time_scheme_controller_t), intent(in) :: ext_bdf
227 type(time_step_controller_t), intent(in) :: dt_controller
228 type(ksp_monitor_t), intent(inout) :: ksp_results
229 end subroutine scalar_scheme_step_intrf
230 end interface
231
232contains
233
242 subroutine scalar_scheme_init(this, msh, c_Xh, gs_Xh, params, scheme, user, &
243 rho)
244 class(scalar_scheme_t), target, intent(inout) :: this
245 type(mesh_t), target, intent(in) :: msh
246 type(coef_t), target, intent(in) :: c_Xh
247 type(gs_t), target, intent(inout) :: gs_Xh
248 type(json_file), target, intent(inout) :: params
249 character(len=*), intent(in) :: scheme
250 type(user_t), target, intent(in) :: user
251 type(field_t), target, intent(in) :: rho
252 ! IO buffer for log output
253 character(len=LOG_SIZE) :: log_buf
254 ! Variables for retrieving json parameters
255 logical :: logical_val
256 real(kind=rp) :: real_val, solver_abstol
257 integer :: integer_val, ierr
258 character(len=:), allocatable :: solver_type, solver_precon
259 type(json_file) :: precon_params
260 type(json_file) :: json_subdict
261 logical :: nut_dependency
262
263 this%u => neko_registry%get_field('u')
264 this%v => neko_registry%get_field('v')
265 this%w => neko_registry%get_field('w')
266 this%rho => rho
267
268 ! Assign a name
269 ! Note that the keyword is added by `scalars_t`, so there is always a
270 ! default.
271 call json_get(params, 'name', this%name)
272
273 ! Set the freeze flag
274 call json_get_or_default(params, 'freeze', this%freeze, .false.)
275
276 call neko_log%section('Scalar')
277 call json_get(params, 'solver.type', solver_type)
278 call json_get(params, 'solver.preconditioner.type', &
279 solver_precon)
280 call json_get(params, 'solver.preconditioner', precon_params)
281 call json_get(params, 'solver.absolute_tolerance', &
282 solver_abstol)
283
284 call json_get_or_default(params, &
285 'solver.projection_space_size', &
286 this%projection_dim, 0)
287 call json_get_or_default(params, &
288 'solver.projection_hold_steps', &
289 this%projection_activ_step, 5)
290
291
292 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
293 call neko_log%message(log_buf)
294 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
295 call neko_log%message(log_buf)
296 call neko_log%message('Ksp scalar : ('// trim(solver_type) // &
297 ', ' // trim(solver_precon) // ')')
298 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
299 call neko_log%message(log_buf)
300
301 this%Xh => this%u%Xh
302 this%dm_Xh => this%u%dof
303 this%params => params
304 this%msh => msh
305
306 call neko_registry%add_field(this%dm_Xh, this%name, &
307 ignore_existing = .true.)
308
309 this%s => neko_registry%get_field(this%name)
310
311 call this%slag%init(this%s, 2)
312
313 this%gs_Xh => gs_xh
314 this%c_Xh => c_xh
315
316 !
317 ! Material properties
318 !
319 call this%set_material_properties(params, user)
320
321
322 !
323 ! Turbulence modelling
324 !
325 this%alphat_field_name = ""
326 this%nut_field_name = ""
327 if (params%valid_path('alphat')) then
328 call json_get(this%params, 'alphat', json_subdict)
329 call json_get(json_subdict, 'nut_dependency', nut_dependency)
330 if (nut_dependency) then
331 call json_get(json_subdict, 'Pr_t', this%pr_turb)
332 call json_get(json_subdict, 'nut_field', this%nut_field_name)
333 else
334 call json_get(json_subdict, 'alphat_field', this%alphat_field_name)
335 end if
336 end if
337
338 !
339 ! Setup right-hand side field.
340 !
341 allocate(this%f_Xh)
342 call this%f_Xh%init(this%dm_Xh, fld_name = "scalar_rhs")
343
344 ! Initialize the source term
345 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
346 call this%source_term%add(params, 'source_terms')
347
348 ! todo parameter file ksp tol should be added
349 call json_get_or_default(params, &
350 'solver.max_iterations', &
351 integer_val, ksp_max_iter)
352 call json_get_or_default(params, &
353 'solver.monitor', &
354 logical_val, .false.)
355 call scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
356 solver_type, integer_val, solver_abstol, logical_val)
357 call scalar_scheme_precon_factory(this%pc, this%ksp, &
358 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
359 solver_precon, precon_params)
360
361 call neko_log%end_section()
362
363 end subroutine scalar_scheme_init
364
365
367 subroutine scalar_scheme_free(this)
368 class(scalar_scheme_t), intent(inout) :: this
369
370 nullify(this%Xh)
371 nullify(this%dm_Xh)
372 nullify(this%gs_Xh)
373 nullify(this%c_Xh)
374 nullify(this%params)
375
376 if (allocated(this%ksp)) then
377 call this%ksp%free()
378 deallocate(this%ksp)
379 end if
380
381 if (allocated(this%pc)) then
382 call precon_destroy(this%pc)
383 deallocate(this%pc)
384 end if
385
386 if (allocated(this%name)) then
387 deallocate(this%name)
388 end if
389
390 call this%source_term%free()
391
392 call this%bcs%free()
393 call this%slag%free()
394
395 nullify(this%cp)
396 nullify(this%lambda)
397 nullify(this%lambda_tot)
398
399 end subroutine scalar_scheme_free
400
403 subroutine scalar_scheme_validate(this)
404 class(scalar_scheme_t), target, intent(inout) :: this
405
406 if ( (.not. allocated(this%u%x)) .or. &
407 (.not. allocated(this%v%x)) .or. &
408 (.not. allocated(this%w%x)) .or. &
409 (.not. allocated(this%s%x))) then
410 call neko_error('Fields are not allocated')
411 end if
412
413 if (.not. allocated(this%ksp)) then
414 call neko_error('No Krylov solver for velocity defined')
415 end if
416
417 if (.not. associated(this%Xh)) then
418 call neko_error('No function space defined')
419 end if
420
421 if (.not. associated(this%dm_Xh)) then
422 call neko_error('No dofmap defined')
423 end if
424
425 if (.not. associated(this%c_Xh)) then
426 call neko_error('No coefficients defined')
427 end if
428
429 if (.not. associated(this%f_Xh)) then
430 call neko_error('No rhs allocated')
431 end if
432
433 if (.not. associated(this%params)) then
434 call neko_error('No parameters defined')
435 end if
436
437 if (.not. associated(this%rho)) then
438 call neko_error('No density field defined')
439 end if
440
441 end subroutine scalar_scheme_validate
442
445 subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
446 abstol, monitor)
447 class(ksp_t), allocatable, target, intent(inout) :: ksp
448 integer, intent(in), value :: n
449 integer, intent(in) :: max_iter
450 character(len=*), intent(in) :: solver
451 real(kind=rp) :: abstol
452 logical, intent(in) :: monitor
453
454 call krylov_solver_factory(ksp, n, solver, max_iter, &
455 abstol, monitor = monitor)
456
457 end subroutine scalar_scheme_solver_factory
458
460 subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
461 pctype, pcparams)
462 class(pc_t), allocatable, target, intent(inout) :: pc
463 class(ksp_t), target, intent(inout) :: ksp
464 type(coef_t), target, intent(in) :: coef
465 type(dofmap_t), target, intent(in) :: dof
466 type(gs_t), target, intent(inout) :: gs
467 type(bc_list_t), target, intent(inout) :: bclst
468 character(len=*) :: pctype
469 type(json_file), intent(inout) :: pcparams
470
471 call precon_factory(pc, pctype)
472
473 select type (pcp => pc)
474 type is (jacobi_t)
475 call pcp%init(coef, dof, gs)
476 type is (sx_jacobi_t)
477 call pcp%init(coef, dof, gs)
478 type is (device_jacobi_t)
479 call pcp%init(coef, dof, gs)
480 type is (hsmg_t)
481 call pcp%init(coef, bclst, pcparams)
482 end select
483
484 call ksp%set_pc(pc)
485
486 end subroutine scalar_scheme_precon_factory
487
493 class(scalar_scheme_t), intent(inout) :: this
494 type(time_state_t), intent(in) :: time
495 type(field_t), pointer :: nut, alphat
496 integer :: index
497 ! Factor to transform nu_t to lambda_t
498 type(field_t), pointer :: lambda_factor
499
500 call this%user_material_properties(this%name, this%material_properties, &
501 time)
502
503 ! factor = rho * cp / pr_turb
504 if (len_trim(this%nut_field_name) .gt. 0 &
505 .and. len_trim(this%alphat_field_name) .eq. 0 ) then
506 nut => neko_registry%get_field(this%nut_field_name)
507
508 ! lambda_tot = lambda + rho * cp * nut / pr_turb
509 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
510 call field_cmult2(lambda_factor, nut, 1.0_rp / this%pr_turb)
511 call field_col2(lambda_factor, this%cp)
512 call field_col2(lambda_factor, this%rho)
513 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
514 call neko_scratch_registry%relinquish_field(index)
515
516 else if (len_trim(this%alphat_field_name) .gt. 0 &
517 .and. len_trim(this%nut_field_name) .eq. 0 ) then
518 alphat => neko_registry%get_field(this%alphat_field_name)
519
520 ! lambda_tot = lambda + rho * cp * alphat
521 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
522 call field_col3(lambda_factor, this%cp, alphat)
523 call field_col2(lambda_factor, this%rho)
524 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
525 call neko_scratch_registry%relinquish_field(index)
526
527 else if (len_trim(this%alphat_field_name) .gt. 0 &
528 .and. len_trim(this%nut_field_name) .gt. 0 ) then
529 call neko_error("Conflicting definition of eddy diffusivity " // &
530 "for the scalar equation")
531 end if
532
533 ! Since cp is a fields and we use the %x(1,1,1,1) of the
534 ! host array data to pass constant material properties
535 ! to some routines, we need to make sure that the host
536 ! values are also filled
537 if (neko_bcknd_device .eq. 1) then
538 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
539 device_to_host, sync=.false.)
540 end if
541
543
547 subroutine scalar_scheme_set_material_properties(this, params, user)
548 class(scalar_scheme_t), intent(inout) :: this
549 type(json_file), intent(inout) :: params
550 type(user_t), target, intent(in) :: user
551 character(len=LOG_SIZE) :: log_buf
552 ! A local pointer that is needed to make Intel happy
553 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
554 real(kind=rp) :: const_cp, const_lambda
555 ! Dummy time state set to 0
556 type(time_state_t) :: time
557
558 dummy_mp_ptr => dummy_user_material_properties
559
560 ! Fill lambda field with the physical value
561
562 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda")
563 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda_tot")
564 call neko_registry%add_field(this%dm_Xh, this%name // "_cp")
565 this%lambda => neko_registry%get_field(this%name // "_lambda")
566 this%lambda_tot => neko_registry%get_field(this%name // "_lambda_tot")
567 this%cp => neko_registry%get_field(this%name // "_cp")
568
569 call this%material_properties%init(2)
570 call this%material_properties%assign(1, this%cp)
571 call this%material_properties%assign(2, this%lambda)
572
573 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
574
575 write(log_buf, '(A)') "Material properties must be set in the user " // &
576 "file!"
577 call neko_log%message(log_buf)
578 this%user_material_properties => user%material_properties
579
580 call user%material_properties(this%name, this%material_properties, time)
581 else
582 this%user_material_properties => dummy_user_material_properties
583 if (params%valid_path('Pe') .and. &
584 (params%valid_path('lambda') .or. &
585 params%valid_path('cp'))) then
586 call neko_error("To set the material properties for the scalar, " // &
587 "either provide Pe OR lambda and cp in the case file.")
588 ! Non-dimensional case
589 else if (params%valid_path('Pe')) then
590 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
591 ' input.'
592 call neko_log%message(log_buf, lvl = neko_log_verbose)
593 write(log_buf, '(A)') 'Specific heat capacity will be set to 1,'
594 call neko_log%message(log_buf, lvl = neko_log_verbose)
595 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
596 call neko_log%message(log_buf, lvl = neko_log_verbose)
597
598 ! Read Pe into lambda for further manipulation.
599 call json_get(params, 'Pe', const_lambda)
600 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
601 call neko_log%message(log_buf)
602
603 ! Set cp and rho to 1 since the setup is non-dimensional.
604 const_cp = 1.0_rp
605 ! Invert the Pe to get conductivity
606 const_lambda = 1.0_rp/const_lambda
607 ! Dimensional case
608 else
609 call json_get(params, 'lambda', const_lambda)
610 call json_get(params, 'cp', const_cp)
611 end if
612 end if
613 ! We need to fill the fields based on the parsed const values
614 ! if the user routine is not used.
615 if (associated(user%material_properties, dummy_mp_ptr)) then
616 ! Fill mu and rho field with the physical value
617 call field_cfill(this%lambda, const_lambda)
618 call field_cfill(this%cp, const_cp)
619
620 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
621 call neko_log%message(log_buf)
622 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
623 call neko_log%message(log_buf)
624 end if
625
626 ! Copy over material property to the total one
627 call field_copy(this%lambda_tot, this%lambda)
628
629 ! Since cp is a field and we use the %x(1,1,1,1) of the
630 ! host array data to pass constant material properties
631 ! to some routines, we need to make sure that the host
632 ! values are also filled
633 if (neko_bcknd_device .eq. 1) then
634 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
635 device_to_host, sync=.false.)
636 end if
638
639end module scalar_scheme
Copy data between host and device (or device and device)
Definition device.F90:71
Abstract interface defining a dirichlet condition on a list of fields.
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
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Jacobi preconditioner accelerator backend.
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public device_to_host
Definition device.F90:47
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a zone as a subset of facets in a mesh.
Defines user dirichlet condition for a scalar field.
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:82
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:487
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:63
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:65
Build configurations.
integer, parameter neko_bcknd_device
Defines a Neumann boundary condition.
Definition neumann.f90:34
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:128
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.
subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype, pcparams)
Initialize a Krylov preconditioner.
subroutine scalar_scheme_free(this)
Deallocate a scalar formulation.
subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, abstol, monitor)
Initialize a linear solver.
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.
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.
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
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
Base type for a boundary condition.
Definition bc.f90:61
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:56
Defines a jacobi preconditioner.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:48
User defined dirichlet condition, for which the user can work with an entire field....
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
A Neumann boundary condition. Sets the flux of the field to the chosen values.
Definition neumann.f90:60
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Base type for a scalar advection-diffusion solver.
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...