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