Neko 1.99.1
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
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
77 use time_state, only : time_state_t
79 implicit none
80
82 type, abstract :: scalar_scheme_t
84 character(len=:), allocatable :: name
86 type(field_t), pointer :: u
88 type(field_t), pointer :: v
90 type(field_t), pointer :: w
92 type(field_t), pointer :: s
94 type(field_series_t) :: slag
96 type(space_t), pointer :: xh
98 type(dofmap_t), pointer :: dm_xh
100 type(gs_t), pointer :: gs_xh
102 type(coef_t), pointer :: c_xh
104 type(field_t), pointer :: f_xh => null()
108 class(ksp_t), allocatable :: ksp
110 integer :: ksp_maxiter
112 integer :: projection_dim
113
114 integer :: projection_activ_step
116 class(pc_t), allocatable :: pc
118 type(bc_list_t) :: bcs
120 type(json_file), pointer :: params
122 type(mesh_t), pointer :: msh => null()
124 type(chkp_t), pointer :: chkp => null()
126 character(len=:), allocatable :: nut_field_name
128 type(field_t), pointer :: rho => null()
130 type(field_t), pointer :: lambda => null()
132 type(field_t), pointer :: cp => null()
134 type(field_t), pointer :: lambda_tot => null()
136 real(kind=rp) :: pr_turb
138 type(field_list_t) :: material_properties
139 procedure(user_material_properties_intf), nopass, pointer :: &
140 user_material_properties => null()
141 contains
143 procedure, pass(this) :: scheme_init => scalar_scheme_init
145 procedure, pass(this) :: scheme_free => scalar_scheme_free
147 procedure, pass(this) :: validate => scalar_scheme_validate
149 procedure, pass(this) :: set_material_properties => &
152 procedure, pass(this) :: update_material_properties => &
155 procedure(scalar_scheme_init_intrf), pass(this), deferred :: init
157 procedure(scalar_scheme_free_intrf), pass(this), deferred :: free
159 procedure(scalar_scheme_step_intrf), pass(this), deferred :: step
161 procedure(scalar_scheme_restart_intrf), pass(this), deferred :: restart
162 end type scalar_scheme_t
163
165 abstract interface
166 subroutine scalar_scheme_init_intrf(this, msh, coef, gs, params, &
167 numerics_params, user, chkp, ulag, vlag, wlag, time_scheme, rho)
168 import scalar_scheme_t
169 import json_file
170 import coef_t
171 import gs_t
172 import mesh_t
173 import user_t
174 import field_series_t, field_t
176 import rp
177 import chkp_t
178 class(scalar_scheme_t), target, intent(inout) :: this
179 type(mesh_t), target, intent(in) :: msh
180 type(coef_t), target, intent(in) :: coef
181 type(gs_t), target, intent(inout) :: gs
182 type(json_file), target, intent(inout) :: params
183 type(json_file), target, intent(inout) :: numerics_params
184 type(user_t), target, intent(in) :: user
185 type(chkp_t), target, intent(inout) :: chkp
186 type(field_series_t), target, intent(in) :: ulag, vlag, wlag
187 type(time_scheme_controller_t), target, intent(in) :: time_scheme
188 type(field_t), target, intent(in) :: rho
189 end subroutine scalar_scheme_init_intrf
190 end interface
191
193 abstract interface
194 subroutine scalar_scheme_restart_intrf(this, chkp)
195 import scalar_scheme_t
196 import chkp_t
197 import rp
198 class(scalar_scheme_t), target, intent(inout) :: this
199 type(chkp_t), intent(inout) :: chkp
200 end subroutine scalar_scheme_restart_intrf
201 end interface
202
204 abstract interface
206 import scalar_scheme_t
207 class(scalar_scheme_t), intent(inout) :: this
208 end subroutine scalar_scheme_free_intrf
209 end interface
210
212 abstract interface
213 subroutine scalar_scheme_step_intrf(this, time, ext_bdf, dt_controller, &
214 ksp_results)
215 import scalar_scheme_t
216 import time_state_t
219 import ksp_monitor_t
220 class(scalar_scheme_t), intent(inout) :: this
221 type(time_state_t), intent(in) :: time
222 type(time_scheme_controller_t), intent(in) :: ext_bdf
223 type(time_step_controller_t), intent(in) :: dt_controller
224 type(ksp_monitor_t), intent(inout) :: ksp_results
225 end subroutine scalar_scheme_step_intrf
226 end interface
227
228contains
229
238 subroutine scalar_scheme_init(this, msh, c_Xh, gs_Xh, params, scheme, user, &
239 rho)
240 class(scalar_scheme_t), target, intent(inout) :: this
241 type(mesh_t), target, intent(in) :: msh
242 type(coef_t), target, intent(in) :: c_Xh
243 type(gs_t), target, intent(inout) :: gs_Xh
244 type(json_file), target, intent(inout) :: params
245 character(len=*), intent(in) :: scheme
246 type(user_t), target, intent(in) :: user
247 type(field_t), target, intent(in) :: rho
248 ! IO buffer for log output
249 character(len=LOG_SIZE) :: log_buf
250 ! Variables for retrieving json parameters
251 logical :: logical_val
252 real(kind=rp) :: real_val, solver_abstol
253 integer :: integer_val, ierr
254 character(len=:), allocatable :: solver_type, solver_precon
255 type(json_file) :: precon_params
256 real(kind=rp) :: gjp_param_a, gjp_param_b
257
258 this%u => neko_field_registry%get_field('u')
259 this%v => neko_field_registry%get_field('v')
260 this%w => neko_field_registry%get_field('w')
261 this%rho => rho
262
263 ! Assign a name
264 ! Note that the keyword is added by `scalars_t`, so there is always a
265 ! default.
266 call json_get(params, 'name', this%name)
267
268 call neko_log%section('Scalar')
269 call json_get(params, 'solver.type', solver_type)
270 call json_get(params, 'solver.preconditioner.type', &
271 solver_precon)
272 call json_get(params, 'solver.preconditioner', precon_params)
273 call json_get(params, 'solver.absolute_tolerance', &
274 solver_abstol)
275
276 call json_get_or_default(params, &
277 'solver.projection_space_size', &
278 this%projection_dim, 0)
279 call json_get_or_default(params, &
280 'solver.projection_hold_steps', &
281 this%projection_activ_step, 5)
282
283
284 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
285 call neko_log%message(log_buf)
286 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
287 call neko_log%message(log_buf)
288 call neko_log%message('Ksp scalar : ('// trim(solver_type) // &
289 ', ' // trim(solver_precon) // ')')
290 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
291 call neko_log%message(log_buf)
292
293 this%Xh => this%u%Xh
294 this%dm_Xh => this%u%dof
295 this%params => params
296 this%msh => msh
297
298 call neko_field_registry%add_field(this%dm_Xh, this%name, &
299 ignore_existing = .true.)
300
301 this%s => neko_field_registry%get_field(this%name)
302
303 call this%slag%init(this%s, 2)
304
305 this%gs_Xh => gs_xh
306 this%c_Xh => c_xh
307
308 !
309 ! Material properties
310 !
311 call this%set_material_properties(params, user)
312
313
314 !
315 ! Turbulence modelling
316 !
317 if (params%valid_path('nut_field')) then
318 call json_get(params, 'Pr_t', this%pr_turb)
319 call json_get(params, 'nut_field', this%nut_field_name)
320 else
321 this%nut_field_name = ""
322 end if
323
324 !
325 ! Setup right-hand side field.
326 !
327 allocate(this%f_Xh)
328 call this%f_Xh%init(this%dm_Xh, fld_name = "scalar_rhs")
329
330 ! Initialize the source term
331 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
332 call this%source_term%add(params, 'source_terms')
333
334 ! todo parameter file ksp tol should be added
335 call json_get_or_default(params, &
336 'solver.max_iterations', &
337 integer_val, ksp_max_iter)
338 call json_get_or_default(params, &
339 'solver.monitor', &
340 logical_val, .false.)
341 call scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
342 solver_type, integer_val, solver_abstol, logical_val)
343 call scalar_scheme_precon_factory(this%pc, this%ksp, &
344 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
345 solver_precon, precon_params)
346
347 call neko_log%end_section()
348
349 end subroutine scalar_scheme_init
350
351
353 subroutine scalar_scheme_free(this)
354 class(scalar_scheme_t), intent(inout) :: this
355
356 nullify(this%Xh)
357 nullify(this%dm_Xh)
358 nullify(this%gs_Xh)
359 nullify(this%c_Xh)
360 nullify(this%params)
361
362 if (allocated(this%ksp)) then
363 call this%ksp%free()
364 deallocate(this%ksp)
365 end if
366
367 if (allocated(this%pc)) then
368 call precon_destroy(this%pc)
369 deallocate(this%pc)
370 end if
371
372 if (allocated(this%name)) then
373 deallocate(this%name)
374 end if
375
376 call this%source_term%free()
377
378 call this%bcs%free()
379 call this%slag%free()
380
381 nullify(this%cp)
382 nullify(this%lambda)
383 nullify(this%lambda_tot)
384
385 end subroutine scalar_scheme_free
386
389 subroutine scalar_scheme_validate(this)
390 class(scalar_scheme_t), target, intent(inout) :: this
391
392 if ( (.not. allocated(this%u%x)) .or. &
393 (.not. allocated(this%v%x)) .or. &
394 (.not. allocated(this%w%x)) .or. &
395 (.not. allocated(this%s%x))) then
396 call neko_error('Fields are not allocated')
397 end if
398
399 if (.not. allocated(this%ksp)) then
400 call neko_error('No Krylov solver for velocity defined')
401 end if
402
403 if (.not. associated(this%Xh)) then
404 call neko_error('No function space defined')
405 end if
406
407 if (.not. associated(this%dm_Xh)) then
408 call neko_error('No dofmap defined')
409 end if
410
411 if (.not. associated(this%c_Xh)) then
412 call neko_error('No coefficients defined')
413 end if
414
415 if (.not. associated(this%f_Xh)) then
416 call neko_error('No rhs allocated')
417 end if
418
419 if (.not. associated(this%params)) then
420 call neko_error('No parameters defined')
421 end if
422
423 if (.not. associated(this%rho)) then
424 call neko_error('No density field defined')
425 end if
426
427 end subroutine scalar_scheme_validate
428
431 subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
432 abstol, monitor)
433 class(ksp_t), allocatable, target, intent(inout) :: ksp
434 integer, intent(in), value :: n
435 integer, intent(in) :: max_iter
436 character(len=*), intent(in) :: solver
437 real(kind=rp) :: abstol
438 logical, intent(in) :: monitor
439
440 call krylov_solver_factory(ksp, n, solver, max_iter, &
441 abstol, monitor = monitor)
442
443 end subroutine scalar_scheme_solver_factory
444
446 subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
447 pctype, pcparams)
448 class(pc_t), allocatable, target, intent(inout) :: pc
449 class(ksp_t), target, intent(inout) :: ksp
450 type(coef_t), target, intent(in) :: coef
451 type(dofmap_t), target, intent(in) :: dof
452 type(gs_t), target, intent(inout) :: gs
453 type(bc_list_t), target, intent(inout) :: bclst
454 character(len=*) :: pctype
455 type(json_file), intent(inout) :: pcparams
456
457 call precon_factory(pc, pctype)
458
459 select type (pcp => pc)
460 type is (jacobi_t)
461 call pcp%init(coef, dof, gs)
462 type is (sx_jacobi_t)
463 call pcp%init(coef, dof, gs)
464 type is (device_jacobi_t)
465 call pcp%init(coef, dof, gs)
466 type is (hsmg_t)
467 call pcp%init(coef, bclst, pcparams)
468 end select
469
470 call ksp%set_pc(pc)
471
472 end subroutine scalar_scheme_precon_factory
473
479 class(scalar_scheme_t), intent(inout) :: this
480 type(time_state_t), intent(in) :: time
481 type(field_t), pointer :: nut
482 integer :: index
483 ! Factor to transform nu_t to lambda_t
484 type(field_t), pointer :: lambda_factor
485
486 call this%user_material_properties(this%name, this%material_properties, &
487 time)
488
489 ! factor = rho * cp / pr_turb
490 if (len(trim(this%nut_field_name)) > 0) then
491 nut => neko_field_registry%get_field(this%nut_field_name)
492
493 ! lambda_tot = lambda + rho * cp * nut / pr_turb
494 call neko_scratch_registry%request_field(lambda_factor, index)
495 call field_col3(lambda_factor, this%cp, this%rho)
496 call field_cmult2(lambda_factor, nut, 1.0_rp / this%pr_turb)
497 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
498 call neko_scratch_registry%relinquish_field(index)
499 end if
500
501 ! Since cp is a fields and we use the %x(1,1,1,1) of the
502 ! host array data to pass constant material properties
503 ! to some routines, we need to make sure that the host
504 ! values are also filled
505 if (neko_bcknd_device .eq. 1) then
506 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
507 device_to_host, sync=.false.)
508 end if
509
511
515 subroutine scalar_scheme_set_material_properties(this, params, user)
516 class(scalar_scheme_t), intent(inout) :: this
517 type(json_file), intent(inout) :: params
518 type(user_t), target, intent(in) :: user
519 character(len=LOG_SIZE) :: log_buf
520 ! A local pointer that is needed to make Intel happy
521 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
522 real(kind=rp) :: const_cp, const_lambda
523 ! Dummy time state set to 0
524 type(time_state_t) :: time
525
526 dummy_mp_ptr => dummy_user_material_properties
527
528 ! Fill lambda field with the physical value
529
530 call neko_field_registry%add_field(this%dm_Xh, this%name // "_lambda")
531 call neko_field_registry%add_field(this%dm_Xh, this%name // "_lambda_tot")
532 call neko_field_registry%add_field(this%dm_Xh, this%name // "_cp")
533 this%lambda => neko_field_registry%get_field(this%name // "_lambda")
534 this%lambda_tot => neko_field_registry%get_field(this%name // "_lambda_tot")
535 this%cp => neko_field_registry%get_field(this%name // "_cp")
536
537 call this%material_properties%init(2)
538 call this%material_properties%assign(1, this%cp)
539 call this%material_properties%assign(2, this%lambda)
540
541 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
542
543 write(log_buf, '(A)') "Material properties must be set in the user " // &
544 "file!"
545 call neko_log%message(log_buf)
546 this%user_material_properties => user%material_properties
547
548 call user%material_properties(this%name, this%material_properties, time)
549 else
550 this%user_material_properties => dummy_user_material_properties
551 if (params%valid_path('Pe') .and. &
552 (params%valid_path('lambda') .or. &
553 params%valid_path('cp'))) then
554 call neko_error("To set the material properties for the scalar, " // &
555 "either provide Pe OR lambda and cp in the case file.")
556 ! Non-dimensional case
557 else if (params%valid_path('Pe')) then
558 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
559 ' input.'
560 call neko_log%message(log_buf, lvl = neko_log_verbose)
561 write(log_buf, '(A)') 'Specific heat capacity will be set to 1,'
562 call neko_log%message(log_buf, lvl = neko_log_verbose)
563 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
564 call neko_log%message(log_buf, lvl = neko_log_verbose)
565
566 ! Read Pe into lambda for further manipulation.
567 call json_get(params, 'Pe', const_lambda)
568 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
569 call neko_log%message(log_buf)
570
571 ! Set cp and rho to 1 since the setup is non-dimensional.
572 const_cp = 1.0_rp
573 ! Invert the Pe to get conductivity
574 const_lambda = 1.0_rp/const_lambda
575 ! Dimensional case
576 else
577 call json_get(params, 'lambda', const_lambda)
578 call json_get(params, 'cp', const_cp)
579 end if
580 end if
581 ! We need to fill the fields based on the parsed const values
582 ! if the user routine is not used.
583 if (associated(user%material_properties, dummy_mp_ptr)) then
584 ! Fill mu and rho field with the physical value
585 call field_cfill(this%lambda, const_lambda)
586 call field_cfill(this%cp, const_cp)
587
588 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
589 call neko_log%message(log_buf)
590 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
591 call neko_log%message(log_buf)
592 end if
593
594 ! Copy over material property to the total one
595 call field_copy(this%lambda_tot, this%lambda)
596
597 ! Since cp is a field and we use the %x(1,1,1,1) of the
598 ! host array data to pass constant material properties
599 ! to some routines, we need to make sure that the host
600 ! values are also filled
601 if (neko_bcknd_device .eq. 1) then
602 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
603 device_to_host, sync=.false.)
604 end if
606
607end 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:42
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_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 .
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
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:76
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
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:62
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:64
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
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 fields This can be used when you have a funct...
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:55
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 for scalar fields. Sets the flux of the field to the chosen value.
Definition neumann.f90:56
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...