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