Neko 1.99.2
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
207 abstract interface
209 import scalar_scheme_t
210 class(scalar_scheme_t), intent(inout) :: this
211 end subroutine scalar_scheme_free_intrf
212 end interface
213
215 abstract interface
216 subroutine scalar_scheme_step_intrf(this, time, ext_bdf, dt_controller, &
217 ksp_results)
218 import scalar_scheme_t
219 import time_state_t
222 import ksp_monitor_t
223 class(scalar_scheme_t), intent(inout) :: this
224 type(time_state_t), intent(in) :: time
225 type(time_scheme_controller_t), intent(in) :: ext_bdf
226 type(time_step_controller_t), intent(in) :: dt_controller
227 type(ksp_monitor_t), intent(inout) :: ksp_results
228 end subroutine scalar_scheme_step_intrf
229 end interface
230
231contains
232
241 subroutine scalar_scheme_init(this, msh, c_Xh, gs_Xh, params, scheme, user, &
242 rho)
243 class(scalar_scheme_t), target, intent(inout) :: this
244 type(mesh_t), target, intent(in) :: msh
245 type(coef_t), target, intent(in) :: c_Xh
246 type(gs_t), target, intent(inout) :: gs_Xh
247 type(json_file), target, intent(inout) :: params
248 character(len=*), intent(in) :: scheme
249 type(user_t), target, intent(in) :: user
250 type(field_t), target, intent(in) :: rho
251 ! IO buffer for log output
252 character(len=LOG_SIZE) :: log_buf
253 ! Variables for retrieving json parameters
254 logical :: logical_val
255 real(kind=rp) :: real_val, solver_abstol
256 integer :: integer_val, ierr
257 character(len=:), allocatable :: solver_type, solver_precon
258 type(json_file) :: precon_params
259 type(json_file) :: json_subdict
260 logical :: nut_dependency
261
262 this%u => neko_registry%get_field('u')
263 this%v => neko_registry%get_field('v')
264 this%w => neko_registry%get_field('w')
265 this%rho => rho
266
267 ! Assign a name
268 ! Note that the keyword is added by `scalars_t`, so there is always a
269 ! default.
270 call json_get(params, 'name', this%name)
271
272 ! Set the freeze flag
273 call json_get_or_default(params, 'freeze', this%freeze, .false.)
274
275 call neko_log%section('Scalar')
276 call json_get(params, 'solver.type', solver_type)
277 call json_get(params, 'solver.preconditioner.type', &
278 solver_precon)
279 call json_get(params, 'solver.preconditioner', precon_params)
280 call json_get(params, 'solver.absolute_tolerance', &
281 solver_abstol)
282
283 call json_get_or_default(params, &
284 'solver.projection_space_size', &
285 this%projection_dim, 0)
286 call json_get_or_default(params, &
287 'solver.projection_hold_steps', &
288 this%projection_activ_step, 5)
289
290
291 write(log_buf, '(A, A)') 'Type : ', trim(scheme)
292 call neko_log%message(log_buf)
293 write(log_buf, '(A, A)') 'Name : ', trim(this%name)
294 call neko_log%message(log_buf)
295 call neko_log%message('Ksp scalar : ('// trim(solver_type) // &
296 ', ' // trim(solver_precon) // ')')
297 write(log_buf, '(A,ES13.6)') ' `-abs tol :', solver_abstol
298 call neko_log%message(log_buf)
299
300 this%Xh => this%u%Xh
301 this%dm_Xh => this%u%dof
302 this%params => params
303 this%msh => msh
304
305 call neko_registry%add_field(this%dm_Xh, this%name, &
306 ignore_existing = .true.)
307
308 this%s => neko_registry%get_field(this%name)
309
310 call this%slag%init(this%s, 2)
311
312 this%gs_Xh => gs_xh
313 this%c_Xh => c_xh
314
315 !
316 ! Material properties
317 !
318 call this%set_material_properties(params, user)
319
320
321 !
322 ! Turbulence modelling
323 !
324 this%alphat_field_name = ""
325 this%nut_field_name = ""
326 if (params%valid_path('alphat')) then
327 call json_get(this%params, 'alphat', json_subdict)
328 call json_get(json_subdict, 'nut_dependency', nut_dependency)
329 if (nut_dependency) then
330 call json_get(json_subdict, 'Pr_t', this%pr_turb)
331 call json_get(json_subdict, 'nut_field', this%nut_field_name)
332 else
333 call json_get(json_subdict, 'alphat_field', this%alphat_field_name)
334 end if
335 end if
336
337 !
338 ! Setup right-hand side field.
339 !
340 allocate(this%f_Xh)
341 call this%f_Xh%init(this%dm_Xh, fld_name = "scalar_rhs")
342
343 ! Initialize the source term
344 call this%source_term%init(this%f_Xh, this%c_Xh, user, this%name)
345 call this%source_term%add(params, 'source_terms')
346
347 ! todo parameter file ksp tol should be added
348 call json_get_or_default(params, &
349 'solver.max_iterations', &
350 integer_val, ksp_max_iter)
351 call json_get_or_default(params, &
352 'solver.monitor', &
353 logical_val, .false.)
354 call scalar_scheme_solver_factory(this%ksp, this%dm_Xh%size(), &
355 solver_type, integer_val, solver_abstol, logical_val)
356 call scalar_scheme_precon_factory(this%pc, this%ksp, &
357 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs, &
358 solver_precon, precon_params)
359
360 call neko_log%end_section()
361
362 end subroutine scalar_scheme_init
363
364
366 subroutine scalar_scheme_free(this)
367 class(scalar_scheme_t), intent(inout) :: this
368
369 nullify(this%Xh)
370 nullify(this%dm_Xh)
371 nullify(this%gs_Xh)
372 nullify(this%c_Xh)
373 nullify(this%params)
374
375 if (allocated(this%ksp)) then
376 call this%ksp%free()
377 deallocate(this%ksp)
378 end if
379
380 if (allocated(this%pc)) then
381 call precon_destroy(this%pc)
382 deallocate(this%pc)
383 end if
384
385 if (allocated(this%name)) then
386 deallocate(this%name)
387 end if
388
389 call this%source_term%free()
390
391 call this%bcs%free()
392 call this%slag%free()
393
394 nullify(this%cp)
395 nullify(this%lambda)
396 nullify(this%lambda_tot)
397
398 end subroutine scalar_scheme_free
399
402 subroutine scalar_scheme_validate(this)
403 class(scalar_scheme_t), target, intent(inout) :: this
404
405 if ( (.not. allocated(this%u%x)) .or. &
406 (.not. allocated(this%v%x)) .or. &
407 (.not. allocated(this%w%x)) .or. &
408 (.not. allocated(this%s%x))) then
409 call neko_error('Fields are not allocated')
410 end if
411
412 if (.not. allocated(this%ksp)) then
413 call neko_error('No Krylov solver for velocity defined')
414 end if
415
416 if (.not. associated(this%Xh)) then
417 call neko_error('No function space defined')
418 end if
419
420 if (.not. associated(this%dm_Xh)) then
421 call neko_error('No dofmap defined')
422 end if
423
424 if (.not. associated(this%c_Xh)) then
425 call neko_error('No coefficients defined')
426 end if
427
428 if (.not. associated(this%f_Xh)) then
429 call neko_error('No rhs allocated')
430 end if
431
432 if (.not. associated(this%params)) then
433 call neko_error('No parameters defined')
434 end if
435
436 if (.not. associated(this%rho)) then
437 call neko_error('No density field defined')
438 end if
439
440 end subroutine scalar_scheme_validate
441
444 subroutine scalar_scheme_solver_factory(ksp, n, solver, max_iter, &
445 abstol, monitor)
446 class(ksp_t), allocatable, target, intent(inout) :: ksp
447 integer, intent(in), value :: n
448 integer, intent(in) :: max_iter
449 character(len=*), intent(in) :: solver
450 real(kind=rp) :: abstol
451 logical, intent(in) :: monitor
452
453 call krylov_solver_factory(ksp, n, solver, max_iter, &
454 abstol, monitor = monitor)
455
456 end subroutine scalar_scheme_solver_factory
457
459 subroutine scalar_scheme_precon_factory(pc, ksp, coef, dof, gs, bclst, &
460 pctype, pcparams)
461 class(pc_t), allocatable, target, intent(inout) :: pc
462 class(ksp_t), target, intent(inout) :: ksp
463 type(coef_t), target, intent(in) :: coef
464 type(dofmap_t), target, intent(in) :: dof
465 type(gs_t), target, intent(inout) :: gs
466 type(bc_list_t), target, intent(inout) :: bclst
467 character(len=*) :: pctype
468 type(json_file), intent(inout) :: pcparams
469
470 call precon_factory(pc, pctype)
471
472 select type (pcp => pc)
473 type is (jacobi_t)
474 call pcp%init(coef, dof, gs)
475 type is (sx_jacobi_t)
476 call pcp%init(coef, dof, gs)
477 type is (device_jacobi_t)
478 call pcp%init(coef, dof, gs)
479 type is (hsmg_t)
480 call pcp%init(coef, bclst, pcparams)
481 end select
482
483 call ksp%set_pc(pc)
484
485 end subroutine scalar_scheme_precon_factory
486
492 class(scalar_scheme_t), intent(inout) :: this
493 type(time_state_t), intent(in) :: time
494 type(field_t), pointer :: nut, alphat
495 integer :: index
496 ! Factor to transform nu_t to lambda_t
497 type(field_t), pointer :: lambda_factor
498
499 call this%user_material_properties(this%name, this%material_properties, &
500 time)
501
502 ! factor = rho * cp / pr_turb
503 if (len_trim(this%nut_field_name) .gt. 0 &
504 .and. len_trim(this%alphat_field_name) .eq. 0 ) then
505 nut => neko_registry%get_field(this%nut_field_name)
506
507 ! lambda_tot = lambda + rho * cp * nut / pr_turb
508 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
509 call field_cmult2(lambda_factor, nut, 1.0_rp / this%pr_turb)
510 call field_col2(lambda_factor, this%cp)
511 call field_col2(lambda_factor, this%rho)
512 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
513 call neko_scratch_registry%relinquish_field(index)
514
515 else if (len_trim(this%alphat_field_name) .gt. 0 &
516 .and. len_trim(this%nut_field_name) .eq. 0 ) then
517 alphat => neko_registry%get_field(this%alphat_field_name)
518
519 ! lambda_tot = lambda + rho * cp * alphat
520 call neko_scratch_registry%request_field(lambda_factor, index, .false.)
521 call field_col3(lambda_factor, this%cp, alphat)
522 call field_col2(lambda_factor, this%rho)
523 call field_add3(this%lambda_tot, this%lambda, lambda_factor)
524 call neko_scratch_registry%relinquish_field(index)
525
526 else if (len_trim(this%alphat_field_name) .gt. 0 &
527 .and. len_trim(this%nut_field_name) .gt. 0 ) then
528 call neko_error("Conflicting definition of eddy diffusivity " // &
529 "for the scalar equation")
530 end if
531
532 ! Since cp is a fields and we use the %x(1,1,1,1) of the
533 ! host array data to pass constant material properties
534 ! to some routines, we need to make sure that the host
535 ! values are also filled
536 if (neko_bcknd_device .eq. 1) then
537 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
538 device_to_host, sync=.false.)
539 end if
540
542
546 subroutine scalar_scheme_set_material_properties(this, params, user)
547 class(scalar_scheme_t), intent(inout) :: this
548 type(json_file), intent(inout) :: params
549 type(user_t), target, intent(in) :: user
550 character(len=LOG_SIZE) :: log_buf
551 ! A local pointer that is needed to make Intel happy
552 procedure(user_material_properties_intf), pointer :: dummy_mp_ptr
553 real(kind=rp) :: const_cp, const_lambda
554 ! Dummy time state set to 0
555 type(time_state_t) :: time
556
557 dummy_mp_ptr => dummy_user_material_properties
558
559 ! Fill lambda field with the physical value
560
561 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda")
562 call neko_registry%add_field(this%dm_Xh, this%name // "_lambda_tot")
563 call neko_registry%add_field(this%dm_Xh, this%name // "_cp")
564 this%lambda => neko_registry%get_field(this%name // "_lambda")
565 this%lambda_tot => neko_registry%get_field(this%name // "_lambda_tot")
566 this%cp => neko_registry%get_field(this%name // "_cp")
567
568 call this%material_properties%init(2)
569 call this%material_properties%assign(1, this%cp)
570 call this%material_properties%assign(2, this%lambda)
571
572 if (.not. associated(user%material_properties, dummy_mp_ptr)) then
573
574 write(log_buf, '(A)') "Material properties must be set in the user " // &
575 "file!"
576 call neko_log%message(log_buf)
577 this%user_material_properties => user%material_properties
578
579 call user%material_properties(this%name, this%material_properties, time)
580 else
581 this%user_material_properties => dummy_user_material_properties
582 if (params%valid_path('Pe') .and. &
583 (params%valid_path('lambda') .or. &
584 params%valid_path('cp'))) then
585 call neko_error("To set the material properties for the scalar, " // &
586 "either provide Pe OR lambda and cp in the case file.")
587 ! Non-dimensional case
588 else if (params%valid_path('Pe')) then
589 write(log_buf, '(A)') 'Non-dimensional scalar material properties' //&
590 ' input.'
591 call neko_log%message(log_buf, lvl = neko_log_verbose)
592 write(log_buf, '(A)') 'Specific heat capacity will be set to 1,'
593 call neko_log%message(log_buf, lvl = neko_log_verbose)
594 write(log_buf, '(A)') 'conductivity to 1/Pe. Assumes density is 1.'
595 call neko_log%message(log_buf, lvl = neko_log_verbose)
596
597 ! Read Pe into lambda for further manipulation.
598 call json_get(params, 'Pe', const_lambda)
599 write(log_buf, '(A,ES13.6)') 'Pe :', const_lambda
600 call neko_log%message(log_buf)
601
602 ! Set cp and rho to 1 since the setup is non-dimensional.
603 const_cp = 1.0_rp
604 ! Invert the Pe to get conductivity
605 const_lambda = 1.0_rp/const_lambda
606 ! Dimensional case
607 else
608 call json_get(params, 'lambda', const_lambda)
609 call json_get(params, 'cp', const_cp)
610 end if
611 end if
612 ! We need to fill the fields based on the parsed const values
613 ! if the user routine is not used.
614 if (associated(user%material_properties, dummy_mp_ptr)) then
615 ! Fill mu and rho field with the physical value
616 call field_cfill(this%lambda, const_lambda)
617 call field_cfill(this%cp, const_cp)
618
619 write(log_buf, '(A,ES13.6)') 'lambda :', const_lambda
620 call neko_log%message(log_buf)
621 write(log_buf, '(A,ES13.6)') 'cp :', const_cp
622 call neko_log%message(log_buf)
623 end if
624
625 ! Copy over material property to the total one
626 call field_copy(this%lambda_tot, this%lambda)
627
628 ! Since cp is a field and we use the %x(1,1,1,1) of the
629 ! host array data to pass constant material properties
630 ! to some routines, we need to make sure that the host
631 ! values are also filled
632 if (neko_bcknd_device .eq. 1) then
633 call device_memcpy(this%cp%x, this%cp%x_d, this%cp%size(), &
634 device_to_host, sync=.false.)
635 end if
637
638end 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...