Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
flow_ic.f90
Go to the documentation of this file.
1! Copyright (c) 2021-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!
34module flow_ic
35 use num_types, only : rp
36 use logger, only : neko_log, log_size
37 use gather_scatter, only : gs_t, gs_op_add
42 use field, only : field_t
43 use utils, only : neko_error, filename_chsuffix, &
45 use coefs, only : coef_t
46 use math, only : col2, cfill, cfill_mask, abscmp
47 use device_math, only : device_col2
49 use json_module, only : json_file
52 use point_zone, only : point_zone_t
55 use fld_file, only : fld_file_t
56 use file, only : file_t
59 use space, only : space_t, gll
60 use field_list, only : field_list_t
61 use operators, only : rotate_cyc
62 implicit none
63 private
64
68 end interface set_flow_ic
69
71
72contains
73
75 subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
76 type(field_t), intent(inout) :: u
77 type(field_t), intent(inout) :: v
78 type(field_t), intent(inout) :: w
79 type(field_t), intent(inout) :: p
80 type(coef_t), intent(in) :: coef
81 type(gs_t), intent(inout) :: gs
82 character(len=*) :: type
83 type(json_file), intent(inout) :: params
84 real(kind=rp) :: delta, tol
85 real(kind=rp), allocatable :: uinf(:)
86 real(kind=rp), allocatable :: zone_value(:)
87 character(len=:), allocatable :: read_str
88 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
89 logical :: interpolate
90
91
92 !
93 ! Uniform (Uinf, Vinf, Winf)
94 !
95 if (trim(type) .eq. 'uniform') then
96
97 call json_get_or_lookup(params, 'value', uinf)
98 call set_flow_ic_uniform(u, v, w, uinf)
99
100 !
101 ! Blasius boundary layer
102 !
103 else if (trim(type) .eq. 'blasius') then
104
105 call json_get_or_lookup(params, 'delta', delta)
106 call json_get(params, 'approximation', read_str)
107 call json_get_or_lookup(params, 'freestream_velocity', uinf)
108
109 call set_flow_ic_blasius(u, v, w, delta, uinf, read_str)
110
111 !
112 ! Point zone initial condition
113 !
114 else if (trim(type) .eq. 'point_zone') then
115
116 call json_get_or_lookup(params, 'base_value', uinf)
117 call json_get(params, 'zone_name', read_str)
118 call json_get_or_lookup(params, 'zone_value', zone_value)
119
120 call set_flow_ic_point_zone(u, v, w, uinf, read_str, zone_value)
121
122 !
123 ! Field initial condition (from fld file)
124 !
125 else if (trim(type) .eq. 'field') then
126
127 call json_get(params, 'file_name', read_str)
128 fname = trim(read_str)
129 call json_get_or_default(params, 'interpolate', interpolate, &
130 .false.)
131 call json_get_or_lookup_or_default(params, 'tolerance', tol, &
132 0.000001_rp)
133 call json_get_or_default(params, 'mesh_file_name', read_str, "none")
134 mesh_fname = trim(read_str)
135
136 call set_flow_ic_fld(u, v, w, p, fname, interpolate, tol, mesh_fname)
137
138 else
139 call neko_error('Invalid initial condition')
140 end if
141
142 call set_flow_ic_common(u, v, w, p, coef, gs)
143
144 end subroutine set_flow_ic_int
145
147 subroutine set_flow_ic_usr(u, v, w, p, coef, gs, user_proc, scheme_name)
148 type(field_t), target, intent(inout) :: u
149 type(field_t), target, intent(inout) :: v
150 type(field_t), target, intent(inout) :: w
151 type(field_t), target, intent(inout) :: p
152 type(coef_t), intent(in) :: coef
153 type(gs_t), intent(inout) :: gs
154 procedure(user_initial_conditions_intf) :: user_proc
155 character(len=*), intent(in) :: scheme_name
156
157 type(field_list_t) :: fields
158
159
160 call neko_log%message("Type: user")
161
162 call fields%init(4)
163 call fields%assign_to_field(1, u)
164 call fields%assign_to_field(2, v)
165 call fields%assign_to_field(3, w)
166 call fields%assign_to_field(4, p)
167
168 call user_proc(scheme_name, fields)
169
170 call set_flow_ic_common(u, v, w, p, coef, gs)
171
172 end subroutine set_flow_ic_usr
173
176 subroutine set_compressible_flow_ic_usr(rho, u, v, w, p, coef, gs, &
177 user_proc, scheme_name)
178 type(field_t), target, intent(inout) :: rho
179 type(field_t), target, intent(inout) :: u
180 type(field_t), target, intent(inout) :: v
181 type(field_t), target, intent(inout) :: w
182 type(field_t), target, intent(inout) :: p
183 type(coef_t), intent(in) :: coef
184 type(gs_t), intent(inout) :: gs
185 procedure(user_initial_conditions_intf) :: user_proc
186 character(len=*), intent(in) :: scheme_name
187 integer :: n
188 type(field_list_t) :: fields
189
190
191 call neko_log%message("Type: user (compressible flows)")
192
193 call fields%init(5)
194 call fields%assign_to_field(1, rho)
195 call fields%assign_to_field(2, u)
196 call fields%assign_to_field(3, v)
197 call fields%assign_to_field(4, w)
198 call fields%assign_to_field(5, p)
199 call user_proc(scheme_name, fields)
200
201 call set_flow_ic_common(u, v, w, p, coef, gs)
202
203 n = u%dof%size()
204
205 if (neko_bcknd_device .eq. 1) then
206 call device_memcpy(p%x, p%x_d, n, host_to_device, sync = .false.)
207 call device_memcpy(rho%x, rho%x_d, n, host_to_device, sync = .false.)
208 end if
209
210 ! Ensure continuity across elements for initial conditions
211 ! These variables are not treated in the common constructor
212 call gs%op(p%x, p%dof%size(), gs_op_add)
213 call gs%op(rho%x, rho%dof%size(), gs_op_add)
214
215 if (neko_bcknd_device .eq. 1) then
216 call device_col2(rho%x_d, coef%mult_d, rho%dof%size())
217 call device_col2(p%x_d, coef%mult_d, p%dof%size())
218 else
219 call col2(rho%x, coef%mult, rho%dof%size())
220 call col2(p%x, coef%mult, p%dof%size())
221 end if
222
223 end subroutine set_compressible_flow_ic_usr
224
225 subroutine set_flow_ic_common(u, v, w, p, coef, gs)
226 type(field_t), intent(inout) :: u
227 type(field_t), intent(inout) :: v
228 type(field_t), intent(inout) :: w
229 type(field_t), intent(inout) :: p
230 type(coef_t), intent(in) :: coef
231 type(gs_t), intent(inout) :: gs
232 integer :: n
233
234 n = u%dof%size()
235
236 if (neko_bcknd_device .eq. 1) then
237 call device_memcpy(u%x, u%x_d, n, &
238 host_to_device, sync = .false.)
239 call device_memcpy(v%x, v%x_d, n, &
240 host_to_device, sync = .false.)
241 call device_memcpy(w%x, w%x_d, n, &
242 host_to_device, sync = .false.)
243 end if
244
245 ! Ensure continuity across elements for initial conditions
246 call rotate_cyc(u%x, v%x, w%x, 1, coef)
247 call gs%op(u%x, u%dof%size(), gs_op_add)
248 call gs%op(v%x, v%dof%size(), gs_op_add)
249 call gs%op(w%x, w%dof%size(), gs_op_add)
250 call rotate_cyc(u%x, v%x, w%x, 0, coef)
251
252 if (neko_bcknd_device .eq. 1) then
253 call device_col2(u%x_d, coef%mult_d, u%dof%size())
254 call device_col2(v%x_d, coef%mult_d, v%dof%size())
255 call device_col2(w%x_d, coef%mult_d, w%dof%size())
256 else
257 call col2(u%x, coef%mult, u%dof%size())
258 call col2(v%x, coef%mult, v%dof%size())
259 call col2(w%x, coef%mult, w%dof%size())
260 end if
261
262 end subroutine set_flow_ic_common
263
265 subroutine set_flow_ic_uniform(u, v, w, uinf)
266 type(field_t), intent(inout) :: u
267 type(field_t), intent(inout) :: v
268 type(field_t), intent(inout) :: w
269 real(kind=rp), intent(in) :: uinf(3)
270 integer :: n, i
271 character(len=LOG_SIZE) :: log_buf
272
273 call neko_log%message("Type : uniform")
274 write (log_buf, '(A, 3(ES12.6, A))') "Value: [", &
275 (uinf(i), ", ", i = 1, 2), uinf(3), "]"
276 call neko_log%message(log_buf)
277
278 u = uinf(1)
279 v = uinf(2)
280 w = uinf(3)
281 n = u%dof%size()
282 if (neko_bcknd_device .eq. 1) then
283 call cfill(u%x, uinf(1), n)
284 call cfill(v%x, uinf(2), n)
285 call cfill(w%x, uinf(3), n)
286 end if
287
288 end subroutine set_flow_ic_uniform
289
292 subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
293 type(field_t), intent(inout) :: u
294 type(field_t), intent(inout) :: v
295 type(field_t), intent(inout) :: w
296 real(kind=rp), intent(in) :: delta
297 real(kind=rp), intent(in) :: uinf(3)
298 character(len=*), intent(in) :: type
299 procedure(blasius_profile), pointer :: bla => null()
300 integer :: i
301 character(len=LOG_SIZE) :: log_buf
302
303 call neko_log%message("Type : blasius")
304 write (log_buf, '(A,ES12.6)') "delta : ", delta
305 call neko_log%message(log_buf)
306 call neko_log%message("Approximation : " // trim(type))
307 write (log_buf, '(A,"[",2(ES12.6,","),ES12.6,"]")') "Value : ", &
308 uinf(1), uinf(2), uinf(3)
309 call neko_log%message(log_buf)
310
311 select case (trim(type))
312 case ('linear')
313 bla => blasius_linear
314 case ('quadratic')
315 bla => blasius_quadratic
316 case ('cubic')
317 bla => blasius_cubic
318 case ('quartic')
319 bla => blasius_quartic
320 case ('sin')
321 bla => blasius_sin
322 case ('tanh')
323 bla => blasius_tanh
324 case default
325 call neko_error('Invalid Blasius approximation')
326 end select
327
328 if ((uinf(1) .gt. 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
329 .and. abscmp(uinf(3), 0.0_rp)) then
330 do i = 1, u%dof%size()
331 u%x(i,1,1,1) = bla(u%dof%z(i,1,1,1), delta, uinf(1))
332 v%x(i,1,1,1) = 0.0_rp
333 w%x(i,1,1,1) = 0.0_rp
334 end do
335 else if (abscmp(uinf(1), 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
336 .and. abscmp(uinf(3), 0.0_rp)) then
337 do i = 1, u%dof%size()
338 u%x(i,1,1,1) = 0.0_rp
339 v%x(i,1,1,1) = bla(u%dof%x(i,1,1,1), delta, uinf(2))
340 w%x(i,1,1,1) = 0.0_rp
341 end do
342 else if (abscmp(uinf(1), 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
343 .and. (uinf(3) .gt. 0.0_rp)) then
344 do i = 1, u%dof%size()
345 u%x(i,1,1,1) = 0.0_rp
346 v%x(i,1,1,1) = 0.0_rp
347 w%x(i,1,1,1) = bla(u%dof%y(i,1,1,1), delta, uinf(3))
348 end do
349 end if
350
351 end subroutine set_flow_ic_blasius
352
362 subroutine set_flow_ic_point_zone(u, v, w, base_value, zone_name, zone_value)
363 type(field_t), intent(inout) :: u
364 type(field_t), intent(inout) :: v
365 type(field_t), intent(inout) :: w
366 real(kind=rp), intent(in), dimension(3) :: base_value
367 character(len=*), intent(in) :: zone_name
368 real(kind=rp), intent(in) :: zone_value(:)
369 character(len=LOG_SIZE) :: log_buf
370
371 ! Internal variables
372 class(point_zone_t), pointer :: zone
373 integer :: size
374
375 call neko_log%message("Type : point_zone")
376 write (log_buf, '(A,ES12.6)') "Base value : ", base_value
377 call neko_log%message(log_buf)
378 call neko_log%message("Zone name : " // trim(zone_name))
379 write (log_buf, '(A,"[",2(ES12.6,","),ES12.6," ]")') "Value : ", &
380 zone_value(1), zone_value(2), zone_value(3)
381 call neko_log%message(log_buf)
382
383 call set_flow_ic_uniform(u, v, w, base_value)
384 size = u%dof%size()
385
386 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
387
388 call cfill_mask(u%x, zone_value(1), size, zone%mask%get(), zone%size)
389 call cfill_mask(v%x, zone_value(2), size, zone%mask%get(), zone%size)
390 call cfill_mask(w%x, zone_value(3), size, zone%mask%get(), zone%size)
391
392 end subroutine set_flow_ic_point_zone
393
411 subroutine set_flow_ic_fld(u, v, w, p, file_name, &
412 interpolate, tolerance, mesh_file_name)
413 type(field_t), intent(inout) :: u
414 type(field_t), intent(inout) :: v
415 type(field_t), intent(inout) :: w
416 type(field_t), intent(inout) :: p
417 character(len=*), intent(in) :: file_name
418 logical, intent(in) :: interpolate
419 real(kind=rp), intent(in) :: tolerance
420 character(len=*), intent(inout) :: mesh_file_name
421
422 character(len=LOG_SIZE) :: log_buf
423 integer :: sample_idx, sample_mesh_idx
424 integer :: last_index
425 type(fld_file_data_t) :: fld_data
426 type(file_t) :: f
427 logical :: mesh_mismatch
428
429 ! ---- For the mesh to mesh interpolation
430 type(global_interpolation_t) :: global_interp
431 ! -----
432
433 ! ---- For space to space interpolation
434 type(space_t) :: prev_xh
435 type(interpolator_t) :: space_interp
436 ! ----
437
438 call neko_log%message("Type : field")
439 call neko_log%message("File name : " // trim(file_name))
440 write (log_buf, '(A,L1)') "Interpolation : ", interpolate
441 call neko_log%message(log_buf)
442 if (interpolate) then
443 end if
444
445 ! Extract sample index from the file name
446 sample_idx = extract_fld_file_index(file_name, -1)
447
448 if (sample_idx .eq. -1) &
449 call neko_error("Invalid file name for the initial condition. The&
450 & file format must be e.g. 'mean0.f00001'")
451
452 ! Change from "field0.f000*" to "field0.fld" for the fld reader
453 call filename_chsuffix(file_name, file_name, 'fld')
454
455 call fld_data%init
456 call f%init(trim(file_name))
457
458 if (interpolate) then
459
460 ! If no mesh file is specified, use the default file name
461 if (mesh_file_name .eq. "none") then
462 mesh_file_name = trim(file_name)
463 sample_mesh_idx = sample_idx
464 else
465
466 ! Extract sample index from the mesh file name
467 sample_mesh_idx = extract_fld_file_index(mesh_file_name, -1)
468
469 if (sample_mesh_idx .eq. -1) then
470 call neko_error("Invalid file name for the initial condition. &
471 &The file format must be e.g. 'mean0.f00001'")
472 end if
473
474 write (log_buf, '(A,ES12.6)') "Tolerance : ", tolerance
475 call neko_log%message(log_buf)
476 write (log_buf, '(A,A)') "Mesh file : ", &
477 trim(mesh_file_name)
478 call neko_log%message(log_buf)
479
480 end if ! if mesh_file_name .eq. none
481
482 ! Read the mesh coordinates if they are not in our fld file
483 if (sample_mesh_idx .ne. sample_idx) then
484 call f%set_counter(sample_mesh_idx)
485 call f%read(fld_data)
486 end if
487
488 end if
489
490 ! Read the field file containing (u,v,w,p)
491 call f%set_counter(sample_idx)
492 call f%read(fld_data)
493
494 !
495 ! Check that the data in the fld file matches the current case.
496 ! Note that this is a safeguard and there are corner cases where
497 ! two different meshes have the same dimension and same # of elements
498 ! but this should be enough to cover obvious cases.
499 !
500 mesh_mismatch = (fld_data%glb_nelv .ne. u%msh%glb_nelv .or. &
501 fld_data%gdim .ne. u%msh%gdim)
502
503 if (mesh_mismatch .and. .not. interpolate) then
504 call neko_error("The fld file must match the current mesh! &
505 &Use 'interpolate': 'true' to enable interpolation.")
506 else if (.not. mesh_mismatch .and. interpolate) then
507 call neko_log%warning("You have activated interpolation but you might &
508 &still be using the same mesh.")
509 end if
510
511 ! Mesh interpolation if specified
512 if (interpolate) then
513
514 ! Issue a warning if the mesh is in single precision
515 select type (ft => f%file_type)
516 type is (fld_file_t)
517 if (.not. ft%dp_precision) then
518 call neko_warning("The coordinates read from the field file are &
519 &in single precision.")
520 call neko_log%message("It is recommended to use a mesh in double &
521 &precision for better interpolation results.")
522 call neko_log%message("If the interpolation does not work, you&
523 &can try to increase the tolerance.")
524 end if
525 class default
526 end select
527
528 ! Sync coordinates to device for the interpolation
529 if (neko_bcknd_device .eq. 1) then
530 call device_memcpy(fld_data%x%x, fld_data%x%x_d, fld_data%x%size(), &
531 host_to_device, sync = .false.)
532 call device_memcpy(fld_data%y%x, fld_data%y%x_d, fld_data%y%size(), &
533 host_to_device, sync = .false.)
534 call device_memcpy(fld_data%z%x, fld_data%z%x_d, fld_data%z%size(), &
535 host_to_device, sync = .true.)
536 end if
537
538 ! Generates an interpolator object and performs the point search
539 call fld_data%generate_interpolator(global_interp, u%dof, u%msh, &
540 tolerance)
541 call fld_data%u%copy_from(host_to_device, .true.)
542 call fld_data%v%copy_from(host_to_device, .true.)
543 call fld_data%w%copy_from(host_to_device, .true.)
544 call fld_data%p%copy_from(host_to_device, .true.)
545 ! Evaluate velocities and pressure
546
547 call global_interp%evaluate(u%x(:,1,1,1), fld_data%u%x, .false.)
548 call global_interp%evaluate(v%x(:,1,1,1), fld_data%v%x, .false.)
549 call global_interp%evaluate(w%x(:,1,1,1), fld_data%w%x, .false.)
550 call global_interp%evaluate(p%x(:,1,1,1), fld_data%p%x, .false.)
551 call u%copy_from(device_to_host, .true.)
552 call v%copy_from(device_to_host, .true.)
553 call w%copy_from(device_to_host, .true.)
554 call p%copy_from(device_to_host, .true.)
555
556 call global_interp%free
557
558 else ! No interpolation, but potentially just from different spaces
559
560 ! Build a space_t object from the data in the fld file
561 call prev_xh%init(gll, fld_data%lx, fld_data%ly, fld_data%lz)
562 call space_interp%init(u%Xh, prev_xh)
563
564 ! Do the space-to-space interpolation
565 call space_interp%map_host(u%x, fld_data%u%x, fld_data%nelv, u%Xh)
566 call space_interp%map_host(v%x, fld_data%v%x, fld_data%nelv, u%Xh)
567 call space_interp%map_host(w%x, fld_data%w%x, fld_data%nelv, u%Xh)
568 call space_interp%map_host(p%x, fld_data%p%x, fld_data%nelv, u%Xh)
569
570 call space_interp%free
571
572 end if
573
574 ! NOTE: we do not copy (u,v,w) to the GPU since `set_flow_ic_common` does
575 ! the copy for us
576 if (neko_bcknd_device .eq. 1) call device_memcpy(p%x, p%x_d, p%dof%size(), &
577 host_to_device, sync = .false.)
578
579 call fld_data%free
580 call prev_xh%free()
581
582 end subroutine set_flow_ic_fld
583
584end module flow_ic
Copy data between host and device (or device and device)
Definition device.F90:71
Synchronize a device or stream.
Definition device.F90:107
Abstract interface for computing a Blasius flow profile.
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 for user defined initial conditions.
Definition user_intf.f90:65
Coefficients.
Definition coef.f90:34
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
NEKTON fld file format.
Definition fld_file.f90:35
Initial flow condition.
Definition flow_ic.f90:34
subroutine set_flow_ic_usr(u, v, w, p, coef, gs, user_proc, scheme_name)
Set intial flow condition (user defined)
Definition flow_ic.f90:148
subroutine set_flow_ic_point_zone(u, v, w, base_value, zone_name, zone_value)
Set the initial condition of the flow based on a point zone.
Definition flow_ic.f90:363
subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
Set initial flow condition (builtin)
Definition flow_ic.f90:76
subroutine set_compressible_flow_ic_usr(rho, u, v, w, p, coef, gs, user_proc, scheme_name)
Set intial flow condition (user defined) for compressible flows.
Definition flow_ic.f90:178
subroutine set_flow_ic_uniform(u, v, w, uinf)
Uniform initial condition.
Definition flow_ic.f90:266
subroutine, public set_flow_ic_fld(u, v, w, p, file_name, interpolate, tolerance, mesh_file_name)
Set the initial condition of the flow based on a field. @detail The fields are read from an fld file....
Definition flow_ic.f90:413
subroutine set_flow_ic_common(u, v, w, p, coef, gs)
Definition flow_ic.f90:226
subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
Set a Blasius profile as initial condition.
Definition flow_ic.f90:293
Defines a flow profile.
real(kind=rp) function, public blasius_quadratic(y, delta, u)
Quadratic approximate Blasius Profile .
real(kind=rp) function, public blasius_quartic(y, delta, u)
Quartic approximate Blasius Profile .
real(kind=rp) function, public blasius_sin(y, delta, u)
Sinusoidal approximate Blasius Profile .
real(kind=rp) function, public blasius_cubic(y, delta, u)
Cubic approximate Blasius Profile .
real(kind=rp) function, public blasius_tanh(y, delta, u)
Hyperbolic tangent approximate Blasius Profile from O. Savas (2012) where is the 99 percent thickne...
real(kind=rp) function, public blasius_linear(y, delta, u)
Linear approximate Blasius profile .
Gather-scatter.
Implements global_interpolation given a dofmap.
Routines to interpolate between different spaces.
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:77
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 col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:397
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
integer function, public extract_fld_file_index(fld_filename, default_index)
Extracts the index of a field file. For example, "myfield.f00045" will return 45. If the suffix of th...
Definition utils.f90:157
integer, parameter, public neko_fname_len
Definition utils.f90:40
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:140
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
field_list_t, To be able to group fields together
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:55
Interface for NEKTON fld files.
Definition fld_file.f90:64
Implements global interpolation for arbitrary points in the domain.
Interpolation between two space::space_t.
Base abstract type for point zones.
The function space for the SEM solution fields.
Definition space.f90:63