Neko 1.99.3
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 mask
40 use math, only : abscmp, rzero, izero, row_zero, rone, copy, cmult, cadd, &
47 use speclib
48 use dofmap, only : dofmap_t
49 use space, only : space_t, gl, gll, gj
50 use htable
51 use uset
52 use stack
53 use tuple
54 use mesh, only : mesh_t, neko_msh_max_zlbls
55 use point, only : point_t
56 use mesh_field, only : mesh_fld_t
57 use map
59 use mxm_wrapper, only : mxm
61 use file
62 use field, only : field_t, field_ptr_t
63 use field_math
66 use krylov
67 use coefs, only : coef_t
68 use bc, only : bc_t
70 use bc_list, only : bc_list_t
71 use dirichlet, only : dirichlet_t
72 use ax_product, only : ax_t, ax_helm_factory
74 use neko_config
75 use case, only : case_t
77 use output, only : output_t
79 use operators, only : dudxyz, opgrad, ortho, cdtp, conv1, curl, cfl,&
82 use projection
83 use user_intf
84 use signal
85 use time_state
88 use device
98 use map_1d, only : map_1d_t
99 use map_2d, only : map_2d_t
100 use cpr, only : cpr_t, cpr_init, cpr_free
101 use fluid_stats, only : fluid_stats_t
102 use field_list, only : field_list_t
104 use vector, only : vector_t, vector_ptr_t
105 use vector_list, only : vector_list_t
106 use matrix, only : matrix_t
107 use tensor
109 simulation_component_wrapper_t, simulation_component_factory, &
110 simulation_component_allocator, simulation_component_allocate, &
111 register_simulation_component
113 use probes, only : probes_t
125 use point_zone, only : point_zone_t
132 use runtime_stats, only : neko_rt_stats
133 use json_module, only : json_file
135 use bc_list, only : bc_list_t
136 use les_model, only : les_model_t, les_model_allocate, register_les_model, &
137 les_model_factory, les_model_allocator
138 use field_writer, only : field_writer_t
141 use curl_simcomp, only : curl_t
142 use gradient_simcomp, only : gradient_t
144 use force_torque, only : force_torque_t
145 use lambda2, only : lambda2_t
149 register_source_term, source_term_factory, source_term_allocator
151 use ale_manager, only : neko_ale
152 use, intrinsic :: iso_fortran_env
153 use mpi_f08
154 !$ use omp_lib
155 implicit none
156
157contains
158
160 subroutine neko_init(C, filename)
161 type(case_t), target, intent(inout), optional :: C
162 character(len=*), intent(in), optional :: filename
163 character(len=NEKO_FNAME_LEN) :: case_file, args
164 character(len=LOG_SIZE) :: log_buf
165 character(len=10) :: suffix
166 character(10) :: time
167 character(8) :: date
168 integer :: argc, i
169
170 call date_and_time(time = time, date = date)
171
172 call comm_init
174 call jobctrl_init
175 call device_init
176
177 call neko_log%init()
178 call neko_registry%init()
179 call neko_const_registry%init()
180 call neko_scratch_registry%init()
181
183
184 if (present(c)) then
185
186 if (present(filename)) then
187 case_file = filename
188 else
189 !
190 ! Command line arguments
191 !
192 argc = command_argument_count()
193
194 if (argc .lt. 1) then
195 if (pe_rank .eq. 0) write(*,*) 'Usage: ./neko <case file>'
196 stop
197 end if
198
199 call get_command_argument(1, case_file)
200 end if
201
202 call filename_suffix(case_file, suffix)
203
204 if (trim(suffix) .ne. 'case' .and. trim(suffix) .ne. 'json') then
205 call neko_error('Invalid case file')
206 end if
207
208 if (argc .gt. 1) then
209 write(log_buf, '(a)') 'Running with command line arguments: '
210 call neko_log%message(log_buf, neko_log_quiet)
211 do i = 2, argc
212 call get_command_argument(i, args)
213 call neko_log%message(args, neko_log_quiet)
214 end do
215 end if
216
217 !
218 ! Job information
219 !
220 call neko_job_info(date, time)
221
222 !
223 ! Create case
224 !
225 call c%init(case_file)
226
227 !
228 ! Setup runtime statistics
229 !
230 call neko_rt_stats%init(c%params)
231
232
233 !
234 ! Create simulation components
235 !
236 call neko_simcomps%init(c)
237
238 end if
239
240 end subroutine neko_init
241
243 subroutine neko_solve(C)
244 type(case_t), target, intent(inout) :: C
245 type(time_step_controller_t) :: dt_controller
246 type(json_file) :: dt_params
247 real(kind=dp) :: tstep_loop_start_time
248
249 call json_get(c%params, 'case.time', dt_params)
250 call dt_controller%init(dt_params)
251
252 call c%time%reset()
253 call simulation_init(c, dt_controller)
254
255 call profiler_start
256 tstep_loop_start_time = mpi_wtime()
257
258 do while (.not. (c%time%is_done() .or. jobctrl_time_limit()))
259 call simulation_step(c, dt_controller, tstep_loop_start_time)
260 end do
261 call profiler_stop
262
263 call simulation_finalize(c)
264
265 end subroutine neko_solve
266
268 subroutine neko_finalize(C)
269 type(case_t), intent(inout), optional :: c
270
271 call neko_rt_stats%report()
272 call neko_rt_stats%free()
273
274 call neko_scratch_registry%free()
275
276 if (present(c)) then
277 call c%free()
278 end if
279
280 call neko_simcomps%free()
281
282 call neko_registry%free()
283 call neko_user_access%free()
284 call neko_log%free()
285
286 call device_finalize
288 call comm_free
289 end subroutine neko_finalize
290
291
294 subroutine neko_job_info(date, time)
295 character(10), optional, intent(in) :: time
296 character(8), optional, intent(in) :: date
297 character(len=LOG_SIZE) :: log_buf
298 integer :: nthrds, rw, sw
299
300 call neko_log%section("Job Information")
301
302 if (present(time) .and. present(date)) then
303 write(log_buf, '(A,A,A,A,1x,A,1x,A,A,A,A,A)') 'Start time: ', &
304 time(1:2), ':', time(3:4), &
305 '/', date(1:4), '-', date(5:6), '-', date(7:8)
306 call neko_log%message(log_buf, neko_log_quiet)
307 end if
308 write(log_buf, '(a)') 'Running on: '
309 sw = 10
310 if (pe_size .lt. 1e1) then
311 write(log_buf(13:), '(i1,a)') pe_size, ' MPI '
312 if (pe_size .eq. 1) then
313 write(log_buf(19:), '(a)') 'rank'
314 sw = 9
315 else
316 write(log_buf(19:), '(a)') 'ranks'
317 end if
318 rw = 1
319 else if (pe_size .lt. 1e2) then
320 write(log_buf(13:), '(i2,a)') pe_size, ' MPI ranks'
321 rw = 2
322 else if (pe_size .lt. 1e3) then
323 write(log_buf(13:), '(i3,a)') pe_size, ' MPI ranks'
324 rw = 3
325 else if (pe_size .lt. 1e4) then
326 write(log_buf(13:), '(i4,a)') pe_size, ' MPI ranks'
327 rw = 4
328 else if (pe_size .lt. 1e5) then
329 write(log_buf(13:), '(i5,a)') pe_size, ' MPI ranks'
330 rw = 5
331 else
332 write(log_buf(13:), '(i6,a)') pe_size, ' MPI ranks'
333 rw = 6
334 end if
335 nthrds = 1
336 !$omp parallel
337 !$omp master
338 !$ nthrds = omp_get_num_threads()
339 !$omp end master
340 !$omp end parallel
341 if (nthrds .gt. 1) then
342 if (nthrds .lt. 1e1) then
343 write(log_buf(13 + rw + sw:), '(a,i1,a)') ', using ', &
344 nthrds, ' thrds each'
345 else if (nthrds .lt. 1e2) then
346 write(log_buf(13 + rw + sw:), '(a,i2,a)') ', using ', &
347 nthrds, ' thrds each'
348 else if (nthrds .lt. 1e3) then
349 write(log_buf(13 + rw + sw:), '(a,i3,a)') ', using ', &
350 nthrds, ' thrds each'
351 else if (nthrds .lt. 1e4) then
352 write(log_buf(13 + rw + sw:), '(a,i4,a)') ', using ', &
353 nthrds, ' thrds each'
354 end if
355 end if
356 call neko_log%message(log_buf, neko_log_quiet)
357
358 write(log_buf, '(a)') 'CPU type : '
359 call system_cpu_name(log_buf(13:))
360 call neko_log%message(log_buf, neko_log_quiet)
361
362 write(log_buf, '(a)') 'Bcknd type: '
363 if (neko_bcknd_sx .eq. 1) then
364 write(log_buf(13:), '(a)') 'SX-Aurora'
365 else if (neko_bcknd_xsmm .eq. 1) then
366 write(log_buf(13:), '(a)') 'CPU (libxsmm)'
367 else if (neko_bcknd_cuda .eq. 1) then
368 write(log_buf(13:), '(a)') 'Accelerator (CUDA)'
369 else if (neko_bcknd_hip .eq. 1) then
370 write(log_buf(13:), '(a)') 'Accelerator (HIP)'
371 else if (neko_bcknd_opencl .eq. 1) then
372 write(log_buf(13:), '(a)') 'Accelerator (OpenCL)'
373 else
374 write(log_buf(13:), '(a)') 'CPU'
375 end if
376 call neko_log%message(log_buf, neko_log_quiet)
377
378 if (neko_bcknd_hip .eq. 1 .or. neko_bcknd_cuda .eq. 1 .or. &
379 neko_bcknd_opencl .eq. 1) then
380 write(log_buf, '(a)') 'Dev. name : '
381 call device_name(log_buf(13:))
382 call neko_log%message(log_buf, neko_log_quiet)
383 end if
384 write(log_buf, '(a)') 'Real type : '
385 select case (rp)
386 case (real32)
387 write(log_buf(13:), '(a)') 'single precision'
388 case (real64)
389 write(log_buf(13:), '(a)') 'double precision'
390 case (real128)
391 write(log_buf(13:), '(a)') 'quad precision'
392 end select
393 call neko_log%message(log_buf, neko_log_quiet)
394 call neko_log%end()
395 end subroutine neko_job_info
396end 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
ALE Manager: Handles Mesh Motion.
type(ale_manager_t), pointer, public neko_ale
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
Implements boundary_operation_t.
Defines a simulation case.
Definition case.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
subroutine, public comm_free
Definition comm.F90:164
subroutine, public comm_init
Definition comm.F90:74
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
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)
real(kind=rp) function, public device_glmax(a_d, n, strm)
Max of a vector of length n.
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 .
real(kind=rp) function, public device_glmin(a_d, n, strm)
Min of a vector of length n.
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_finalize
Definition device.F90:156
subroutine, public device_name(name)
Definition device.F90:172
subroutine, public device_init
Definition device.F90:134
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 user neumann condition for a scalar field.
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:52
Importation of fields from fld files.
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:81
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
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
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:447
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:459
subroutine, public row_zero(a, m, n, e)
Sets row e to 0 in matrix a.
Definition math.f90:230
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:890
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
Definition math.f90:762
real(kind=rp), parameter, public pi
Definition math.f90:76
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1117
subroutine, public ascol5(a, b, c, d, e, n)
Returns .
Definition math.f90:1052
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition math.f90:689
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition math.f90:510
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:498
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:876
real(kind=rp) function, public glsc4(a, b, c, d, n)
Definition math.f90:1136
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:847
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1098
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:931
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:241
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:1083
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:789
integer function, public glimin(a, n)
Min of an integer vector of length n.
Definition math.f90:611
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:549
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:831
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:1022
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:776
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:523
subroutine, public absval(a, n)
Take the absolute value of an array.
Definition math.f90:1422
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:676
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:945
subroutine, public subcol4(a, b, c, d, n)
Returns .
Definition math.f90:993
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:1008
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:664
subroutine, public chsign(a, n)
Change sign of vector .
Definition math.f90:629
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:904
subroutine, public izero(a, n)
Zero an integer vector.
Definition math.f90:219
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:566
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:252
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:803
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:917
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:734
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:208
subroutine, public vdot2(dot, u1, u2, v1, v2, n)
Compute a dot product (2-d version) assuming vector components etc.
Definition math.f90:720
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:652
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:433
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:641
integer function, public glimax(a, n)
Max of an integer vector of length n.
Definition math.f90:581
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:818
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:862
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:596
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:271
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition math.f90:703
real(kind=rp) function, public vlsc3(u, v, w, n)
Compute multiplication sum .
Definition math.f90:748
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:1068
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:64
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:91
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:92
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.
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
type(registry_t), target, public neko_const_registry
This registry is used to store user-defined scalars and vectors, provided under the constants section...
Definition registry.f90:150
Runtime statistics.
type(runtime_stats_t), public neko_rt_stats
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
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:49
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:41
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:39
Utilities.
Definition utils.f90:35
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition utils.f90:118
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:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
A simulation component for boundary reductions on labelled zones.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
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:49
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:83
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
User defined neumann condition, for which the user can work with an entire field. The type stores a s...
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...
vector_list_t, To be able to group vectors together
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...