Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
neko.f90
Go to the documentation of this file.
1! Copyright (c) 2019-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 neko
35 use num_types, only : rp, sp, dp, qp, c_rp
36 use comm
37 use utils
38 use logger
39 use math, only : abscmp, rzero, izero, row_zero, rone, copy, cmult, cadd, &
46 use speclib
47 use dofmap, only : dofmap_t
48 use space, only : space_t, gl, gll, gj
49 use htable
50 use uset
51 use stack
52 use tuple
53 use mesh, only : mesh_t, neko_msh_max_zlbls
54 use point, only : point_t
55 use mesh_field, only : mesh_fld_t
56 use map
57 use mxm_wrapper, only : mxm
59 use file
60 use field, only : field_t, field_ptr_t
61 use field_math
64 use krylov
65 use coefs, only : coef_t
66 use bc, only : bc_t
68 use bc_list, only : bc_list_t
69 use dirichlet, only : dirichlet_t
70 use ax_product, only : ax_t, ax_helm_factory
72 use neko_config
73 use case, only : case_t, case_init, case_free
75 use output, only : output_t
77 use operators, only : dudxyz, opgrad, ortho, cdtp, conv1, curl, cfl,&
80 use projection
81 use user_intf
82 use signal
83 use time_state
86 use device
96 use map_1d, only : map_1d_t
97 use map_2d, only : map_2d_t
98 use cpr, only : cpr_t, cpr_init, cpr_free
99 use fluid_stats, only : fluid_stats_t
100 use field_list, only : field_list_t
102 use vector, only : vector_t, vector_ptr_t
103 use matrix, only : matrix_t
104 use tensor
106 simulation_component_wrapper_t, simulation_component_factory, &
107 simulation_component_allocator, simulation_component_allocate, &
108 register_simulation_component
109 use probes, only : probes_t
121 use point_zone, only: point_zone_t
127 use runtime_stats, only : neko_rt_stats
128 use json_module, only : json_file
130 use bc_list, only : bc_list_t
131 use les_model, only : les_model_t, les_model_allocate, register_les_model, &
132 les_model_factory, les_model_allocator
133 use field_writer, only : field_writer_t
136 use curl_simcomp, only : curl_t
137 use gradient_simcomp, only : gradient_t
139 use force_torque, only : force_torque_t
140 use lambda2, only : lambda2_t
144 register_source_term, source_term_factory, source_term_allocator
146 use, intrinsic :: iso_fortran_env
147 use mpi_f08
148 !$ use omp_lib
149 implicit none
150
151contains
152
154 subroutine neko_init(C)
155 type(case_t), target, intent(inout), optional :: C
156 character(len=NEKO_FNAME_LEN) :: case_file, args
157 character(len=LOG_SIZE) :: log_buf
158 character(len=10) :: suffix
159 character(10) :: time
160 character(8) :: date
161 integer :: argc, i
162
163 call date_and_time(time = time, date = date)
164
165 call comm_init
167 call jobctrl_init
168 call device_init
169
170 call neko_log%init()
171 call neko_field_registry%init()
172
174
175 if (present(c)) then
176
177 !
178 ! Command line arguments
179 !
180 argc = command_argument_count()
181
182 if (argc .lt. 1) then
183 if (pe_rank .eq. 0) write(*,*) 'Usage: ./neko <case file>'
184 stop
185 end if
186
187 call get_command_argument(1, case_file)
188
189 call filename_suffix(case_file, suffix)
190
191 if (trim(suffix) .ne. 'case' .and. trim(suffix) .ne. 'json') then
192 call neko_error('Invalid case file')
193 end if
194
195 if (argc .gt. 1) then
196 write(log_buf, '(a)') 'Running with command line arguments: '
197 call neko_log%message(log_buf, neko_log_quiet)
198 do i = 2, argc
199 call get_command_argument(i, args)
200 call neko_log%message(args, neko_log_quiet)
201 end do
202 end if
203
204 !
205 ! Job information
206 !
207 call neko_job_info(date, time)
208
209 !
210 ! Create case
211 !
212 call case_init(c, case_file)
213
214 !
215 ! Setup runtime statistics
216 !
217 call neko_rt_stats%init(c%params)
218
219
220 !
221 ! Create simulation components
222 !
223 call neko_simcomps%init(c)
224
225 end if
226
227 end subroutine neko_init
228
230 subroutine neko_solve(C)
231 type(case_t), target, intent(inout) :: C
232 type(time_step_controller_t) :: dt_controller
233 type(json_file) :: dt_params
234 real(kind=dp) :: tstep_loop_start_time
235
236 call json_get(c%params, 'case.time', dt_params)
237 call dt_controller%init(dt_params)
238
239 call c%time%reset()
240 call simulation_init(c, dt_controller)
241
242 call profiler_start
243 tstep_loop_start_time = mpi_wtime()
244
245 do while (.not. (c%time%is_done() .or. jobctrl_time_limit()))
246 call simulation_step(c, dt_controller, tstep_loop_start_time)
247 end do
248 call profiler_stop
249
250 call simulation_finalize(c)
251
252 end subroutine neko_solve
253
255 subroutine neko_finalize(C)
256 type(case_t), intent(inout), optional :: c
257
258 call neko_rt_stats%report()
259 call neko_rt_stats%free()
260
261 call neko_scratch_registry%free()
262
263 if (present(c)) then
264 call case_free(c)
265 end if
266
267 call neko_field_registry%free()
268 call neko_user_access%free()
269 call device_finalize
271 call comm_free
272 end subroutine neko_finalize
273
274
277 subroutine neko_job_info(date, time)
278 character(10), optional, intent(in) :: time
279 character(8), optional, intent(in) :: date
280 character(len=LOG_SIZE) :: log_buf
281 integer :: nthrds, rw, sw
282
283 call neko_log%section("Job Information")
284
285 if (present(time) .and. present(date)) then
286 write(log_buf, '(A,A,A,A,1x,A,1x,A,A,A,A,A)') 'Start time: ', &
287 time(1:2), ':', time(3:4), &
288 '/', date(1:4), '-', date(5:6), '-', date(7:8)
289 call neko_log%message(log_buf, neko_log_quiet)
290 end if
291 write(log_buf, '(a)') 'Running on: '
292 sw = 10
293 if (pe_size .lt. 1e1) then
294 write(log_buf(13:), '(i1,a)') pe_size, ' MPI '
295 if (pe_size .eq. 1) then
296 write(log_buf(19:), '(a)') 'rank'
297 sw = 9
298 else
299 write(log_buf(19:), '(a)') 'ranks'
300 end if
301 rw = 1
302 else if (pe_size .lt. 1e2) then
303 write(log_buf(13:), '(i2,a)') pe_size, ' MPI ranks'
304 rw = 2
305 else if (pe_size .lt. 1e3) then
306 write(log_buf(13:), '(i3,a)') pe_size, ' MPI ranks'
307 rw = 3
308 else if (pe_size .lt. 1e4) then
309 write(log_buf(13:), '(i4,a)') pe_size, ' MPI ranks'
310 rw = 4
311 else if (pe_size .lt. 1e5) then
312 write(log_buf(13:), '(i5,a)') pe_size, ' MPI ranks'
313 rw = 5
314 else
315 write(log_buf(13:), '(i6,a)') pe_size, ' MPI ranks'
316 rw = 6
317 end if
318 nthrds = 1
319 !$omp parallel
320 !$omp master
321 !$ nthrds = omp_get_num_threads()
322 !$omp end master
323 !$omp end parallel
324 if (nthrds .gt. 1) then
325 if (nthrds .lt. 1e1) then
326 write(log_buf(13 + rw + sw:), '(a,i1,a)') ', using ', &
327 nthrds, ' thrds each'
328 else if (nthrds .lt. 1e2) then
329 write(log_buf(13 + rw + sw:), '(a,i2,a)') ', using ', &
330 nthrds, ' thrds each'
331 else if (nthrds .lt. 1e3) then
332 write(log_buf(13 + rw + sw:), '(a,i3,a)') ', using ', &
333 nthrds, ' thrds each'
334 else if (nthrds .lt. 1e4) then
335 write(log_buf(13 + rw + sw:), '(a,i4,a)') ', using ', &
336 nthrds, ' thrds each'
337 end if
338 end if
339 call neko_log%message(log_buf, neko_log_quiet)
340
341 write(log_buf, '(a)') 'CPU type : '
342 call system_cpu_name(log_buf(13:))
343 call neko_log%message(log_buf, neko_log_quiet)
344
345 write(log_buf, '(a)') 'Bcknd type: '
346 if (neko_bcknd_sx .eq. 1) then
347 write(log_buf(13:), '(a)') 'SX-Aurora'
348 else if (neko_bcknd_xsmm .eq. 1) then
349 write(log_buf(13:), '(a)') 'CPU (libxsmm)'
350 else if (neko_bcknd_cuda .eq. 1) then
351 write(log_buf(13:), '(a)') 'Accelerator (CUDA)'
352 else if (neko_bcknd_hip .eq. 1) then
353 write(log_buf(13:), '(a)') 'Accelerator (HIP)'
354 else if (neko_bcknd_opencl .eq. 1) then
355 write(log_buf(13:), '(a)') 'Accelerator (OpenCL)'
356 else
357 write(log_buf(13:), '(a)') 'CPU'
358 end if
359 call neko_log%message(log_buf, neko_log_quiet)
360
361 if (neko_bcknd_hip .eq. 1 .or. neko_bcknd_cuda .eq. 1 .or. &
362 neko_bcknd_opencl .eq. 1) then
363 write(log_buf, '(a)') 'Dev. name : '
364 call device_name(log_buf(13:))
365 call neko_log%message(log_buf, neko_log_quiet)
366 end if
367 write(log_buf, '(a)') 'Real type : '
368 select case (rp)
369 case (real32)
370 write(log_buf(13:), '(a)') 'single precision'
371 case (real64)
372 write(log_buf(13:), '(a)') 'double precision'
373 case (real128)
374 write(log_buf(13:), '(a)') 'quad precision'
375 end select
376 call neko_log%message(log_buf, neko_log_quiet)
377 call neko_log%end()
378 end subroutine neko_job_info
379end module neko
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.
LES model allocator.
Simulation component factory. Both constructs and initializes the object.
Source term factory. Both constructs and initializes the object.
Interface to a C function to retrieve the CPU name (type).
Definition system.f90:44
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Defines a simulation case.
Definition case.f90:34
subroutine, public case_free(this)
Deallocate a case.
Definition case.f90:518
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
subroutine, public comm_free
Definition comm.F90:163
subroutine, public comm_init
Definition comm.F90:73
integer, public pe_size
MPI size of communicator.
Definition comm.F90:58
integer, public pe_rank
MPI rank.
Definition comm.F90:55
Compression.
Definition cpr.f90:34
subroutine, public cpr_free(cpr)
Deallocate coefficients.
Definition cpr.f90:107
Implements the curl_t type.
Implements type data_streamer_t.
Implements the derivative_t type.
subroutine, public device_add2s1(a_d, b_d, c1, n, strm)
subroutine, public device_add2s2_many(y_d, x_d_d, a_d, j, n, strm)
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_invcol1(a_d, n, strm)
Invert a vector .
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n, strm)
Returns .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_rone(a_d, n, strm)
Set all elements to one.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public device_glsc3_many(h, w_d, v_d_d, mult_d, j, n, strm)
subroutine, public device_cfill_mask(a_d, c, n, mask_d, n_mask, strm)
Fill a constant to a masked vector. .
subroutine, public device_cadd2(a_d, b_d, c, n, strm)
Add a scalar to vector .
subroutine, public device_sub2(a_d, b_d, n, strm)
Vector substraction .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
real(kind=rp) function, public device_vlsc3(u_d, v_d, w_d, n, strm)
Compute multiplication sum .
subroutine, public device_masked_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Copy a masked vector .
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_absval(a_d, n, strm)
subroutine, public device_invcol2(a_d, b_d, n, strm)
Vector division .
subroutine, public device_addsqr2s2(a_d, b_d, c1, n, strm)
Returns .
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n, strm)
Weighted inner product .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
subroutine, public device_cmult2(a_d, b_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
subroutine, public device_addcol4(a_d, b_d, c_d, d_d, n, strm)
Returns .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_add3(a_d, b_d, c_d, n, strm)
Vector addition .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_finalize
Definition device.F90:150
subroutine, public device_name(name)
Definition device.F90:166
subroutine, public device_init
Definition device.F90:128
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Implements the divergence_t type.
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
subroutine, public drag_torque_zone(dgtq, tstep, zone, center, s11, s22, s33, s12, s13, s23, p, coef, visc)
Some functions to calculate the lift/drag and torque Calculation can be done on a zone,...
subroutine, public drag_torque_facet(dgtq, xm0, ym0, zm0, center, s11, s22, s33, s12, s13, s23, pm1, visc, f, e, coef, xh)
Calculate drag and torque over a facet.
subroutine, public drag_torque_pt(dgtq, x, y, z, center, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, v)
Calculate drag and torque from one point.
Defines inflow dirichlet conditions.
Defines user dirichlet condition for a scalar field.
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Implements the field_writer_t type.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Computes various statistics for the fluid fields. We use the Reynolds decomposition for a field u = ...
Implements the force_torque_t type.
Gather-scatter.
Implements global_interpolation given a dofmap.
Implements the gradient_t type.
Implements a hash table ADT.
Definition htable.f90:36
Job control.
Definition jobctrl.f90:34
logical function, public jobctrl_time_limit()
Check if the job's time limit has been reached.
Definition jobctrl.f90:108
real(kind=rp) function, public jobctrl_jobtime()
Returns a job's time in seconds relative to the first call.
Definition jobctrl.f90:127
subroutine, public jobctrl_init()
Initialize jobctrl.
Definition jobctrl.f90:56
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
A simulation component that computes lambda2 The values are stored in the field registry under the na...
Definition lambda2.f90:37
Implements les_model_t.
Definition les_model.f90:35
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_quiet
Always logged.
Definition log.f90:72
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
Creates a 1d GLL point map along a specified direction based on the connectivity in the mesh.
Definition map_1d.f90:3
Maps a 3D dofmap to a 2D spectral element grid.
Definition map_2d.f90:3
NEKTON map.
Definition map.f90:3
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:411
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:423
subroutine, public row_zero(a, m, n, e)
Sets row e to 0 in matrix a.
Definition math.f90:227
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:840
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
Definition math.f90:712
real(kind=rp), parameter, public pi
Definition math.f90:75
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1067
subroutine, public ascol5(a, b, c, d, e, n)
Returns .
Definition math.f90:1002
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition math.f90:639
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition math.f90:474
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:462
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:826
real(kind=rp) function, public glsc4(a, b, c, d, n)
Definition math.f90:1086
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:797
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:881
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:238
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:1033
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:739
integer function, public glimin(a, n)
Min of an integer vector of length n.
Definition math.f90:561
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:499
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:781
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:972
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:487
subroutine, public absval(a, n)
Take the absolute value of an array.
Definition math.f90:1372
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:626
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:895
subroutine, public subcol4(a, b, c, d, n)
Returns .
Definition math.f90:943
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:958
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:614
subroutine, public chsign(a, n)
Change sign of vector .
Definition math.f90:579
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public izero(a, n)
Zero an integer vector.
Definition math.f90:216
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:753
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:867
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition math.f90:684
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:205
subroutine, public vdot2(dot, u1, u2, v1, v2, n)
Compute a dot product (2-d version) assuming vector components etc.
Definition math.f90:670
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:602
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:397
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:591
integer function, public glimax(a, n)
Max of an integer vector of length n.
Definition math.f90:531
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:768
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:546
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:268
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition math.f90:653
real(kind=rp) function, public vlsc3(u, v, w, n)
Compute multiplication sum .
Definition math.f90:698
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:1018
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition mathops.f90:65
subroutine, public opcolv(a1, a2, a3, c, gdim, n)
Definition mathops.f90:97
subroutine, public opadd2col(a1, a2, a3, b1, b2, b3, c, n, gdim)
Definition mathops.f90:165
subroutine, public opchsign(a1, a2, a3, gdim, n)
for and .
Definition mathops.f90:76
subroutine, public opadd2cm(a1, a2, a3, b1, b2, b3, c, n, gdim)
Definition mathops.f90:142
subroutine, public opcolv3c(a1, a2, a3, b1, b2, b3, c, d, n, gdim)
Definition mathops.f90:119
Defines a matrix.
Definition matrix.f90:34
Defines a mesh field.
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:62
Wrapper for all matrix-matrix product implementations.
subroutine, public mxm(a, n1, b, n2, c, n3)
Compute matrix-matrix product for contiguously packed matrices A,B, and C.
Build configurations.
integer, parameter neko_bcknd_sx
integer, parameter neko_bcknd_hip
character(len=80), parameter neko_build_info
character(len=10), parameter neko_version
integer, parameter neko_bcknd_opencl
integer, parameter neko_bcknd_cuda
integer, parameter neko_bcknd_xsmm
MPI derived types.
Definition mpi_types.f90:34
subroutine, public neko_mpi_types_init
Define all MPI derived types.
Definition mpi_types.f90:88
subroutine, public neko_mpi_types_free
Deallocate all derived MPI types.
Master module.
Definition neko.f90:34
integer, parameter, public qp
Definition num_types.f90:10
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public c_rp
Definition num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public ortho(x, glb_n_points, n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
real(kind=rp) function, public cfl(dt, u, v, w, xh, coef, nelv, gdim)
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
subroutine, public conv1(du, u, vx, vy, vz, xh, coef, es, ee)
Compute the advection term.
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef, event)
subroutine, public strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
subroutine, public lambda2op(lambda2, u, v, w, coef)
Compute the Lambda2 field for a given velocity field.
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
subroutine, public dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
Definition operators.f90:82
Implements output_controller_t
Defines an output.
Definition output.f90:34
Interface to ParMETIS.
Definition parmetis.F90:34
subroutine, public parmetis_partmeshkway(msh, parts, weights, nprts)
Compute a k-way partitioning of a mesh msh.
Definition parmetis.F90:112
subroutine, public parmetis_partgeom(msh, parts)
Compute a k-way partitioning of a mesh msh using a coordinated-based space-filing curves method.
Definition parmetis.F90:183
Routines to interpolate fields on a given element on a point in that element with given r,...
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Implements a point.
Definition point.f90:35
Implements probes.
Definition probes.F90:37
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start
Start profiling.
Definition profiler.F90:52
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:107
subroutine, public profiler_stop
Stop profiling.
Definition profiler.F90:65
Project x onto X, the space of old solutions and back again.
Runtime statistics.
type(runtime_stats_t), public neko_rt_stats
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Interface to signal handler.
Definition signal.f90:34
Contains the simcomp_executor_t type.
type(simcomp_executor_t), target, public neko_simcomps
Global variable for the simulation component driver.
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
Simulation driver.
subroutine, public simulation_step(c, dt_controller, tstep_loop_start_time)
Compute a single time-step of a case.
subroutine, public simulation_finalize(c)
Finalize a simulation of a case.
subroutine, public simulation_init(c, dt_controller)
Initialise a simulation of a case.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
integer, parameter, public gj
Definition space.f90:49
integer, parameter, public gl
Definition space.f90:49
LIBRARY ROUTINES FOR SPECTRAL METHODS.
Definition speclib.f90:148
Implements type spectral_error_t.
Implements a dynamic stack ADT.
Definition stack.f90:35
Interface to system information routines.
Definition system.f90:34
subroutine, public system_cpu_name(name)
Retrieve the system's CPU name (type)
Definition system.f90:59
Tensor operations.
Definition tensor.f90:61
Contains the time_based_controller_t type.
Implements type time_interpolator_t.
Module with things related to the simulation time.
Implements type time_step_controller.
Implements a n-tuple.
Definition tuple.f90:34
User access singleton Defines a singleton object available in the user file. Intended to allow unrest...
type(user_access_t), target, public neko_user_access
The singleton object.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Implements the user_source_term_t type.
Implements an unordered set ADT.
Definition uset.f90:35
Utilities.
Definition utils.f90:35
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition utils.f90:108
Defines a vector.
Definition vector.f90:34
Implements the weak_gradient_t type.
Defines a zero-valued Dirichlet boundary condition.
void neko_finalize()
void neko_init()
void neko_job_info()
void neko_solve(int **case_iptr)
Base type for a matrix-vector product providing .
Definition ax.f90:43
Base type for a boundary condition.
Definition bc.f90:61
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
include information needed for compressing fields
Definition cpr.f90:49
A simulation component that computes the curl of a vector field. Added to the field registry as curl_...
Provides access to data streaming by interfacing with c++ ADIOS2 subroutines.
A simulation component that computes a derivative of a field. Wraps the duxyz operator.
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:48
A simulation component that computes the divergence of a vector field. Added to the field registry as...
field_ptr_t, To easily obtain a pointer to a field
Definition field.f90:82
User defined dirichlet condition, for which the user can work with an entire field....
Extension of the user defined dirichlet condition field_dirichlet
field_list_t, To be able to group fields together
A simulation component that writes a 3d field to a file.
A simulation component that computes the force and torque on a given boundary zone.
A simulation component that computes the gradient of a field. Wraps the gradient operator.
Base abstract type for LES models based on the Boussinesq approximation.
Definition les_model.f90:64
Type that encapsulates a mapping from each gll point in the mesh to its corresponding (global) GLL po...
Definition map_1d.f90:27
Abstract type defining an output type.
Definition output.f90:41
Centralized controller for a list of outputs.
A point in with coordinates .
Definition point.f90:43
Field interpolator to arbitrary points within an element. Tailored for experimentation,...
Base abstract type for point zones.
Base abstract class for simulation components.
A helper type that is needed to have an array of polymorphic objects.
Base abstract type for source terms.
The function space for the SEM solution fields.
Definition space.f90:63
Provides tools to calculate the spectral error indicator.
A utility type for determining whether an action should be executed based on the current time value....
Provides a tool to perform interpolation in time.
A source term wrapping the user source term routine. Stores fields that are passed to the user routin...
A simulation component that computes the weak gradient of a field. Wraps the opgradient operator.
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...