Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
drag_torque.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2! Copyright (c) 2025, The Neko Authors
3!
4! The UChicago Argonne, LLC as Operator of Argonne National
5! Laboratory holds copyright in the Software. The copyright holder
6! reserves all rights except those expressly granted to licensees,
7! and U.S. Government license rights.
8!
9! Redistribution and use in source and binary forms, with or without
10! modification, are permitted provided that the following conditions
11! are met:
12!
13! 1. Redistributions of source code must retain the above copyright
14! notice, this list of conditions and the disclaimer below.
15!
16! 2. Redistributions in binary form must reproduce the above copyright
17! notice, this list of conditions and the disclaimer (as noted below)
18! in the documentation and/or other materials provided with the
19! distribution.
20!
21! 3. Neither the name of ANL nor the names of its contributors
22! may be used to endorse or promote products derived from this software
23! without specific prior written permission.
24!
25! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
26! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
27! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
28! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
29! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
30! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
31! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
32! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
33! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
34! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
35! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
36! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37!
38! Additional BSD Notice
39! ---------------------
40! 1. This notice is required to be provided under our contract with
41! the U.S. Department of Energy (DOE). This work was produced at
42! Argonne National Laboratory under Contract
43! No. DE-AC02-06CH11357 with the DOE.
44!
45! 2. Neither the United States Government nor UCHICAGO ARGONNE,
46! LLC nor any of their employees, makes any warranty,
47! express or implied, or assumes any liability or responsibility for the
48! accuracy, completeness, or usefulness of any information, apparatus,
49! product, or process disclosed, or represents that its use would not
50! infringe privately-owned rights.
51!
52! 3. Also, reference herein to any specific commercial products, process,
53! or services by trade name, trademark, manufacturer or otherwise does
54! not necessarily constitute or imply its endorsement, recommendation,
55! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
56! The views and opinions of authors expressed
57! herein do not necessarily state or reflect those of the United States
58! Government or UCHICAGO ARGONNE, LLC, and shall
59! not be used for advertising or product endorsement purposes.
60!
61
63 use field, only : field_t
64 use coefs, only : coef_t
65 use vector, only : vector_t
66 use facet_zone, only : facet_zone_t
67 use math, only : rzero, col3, vdot3, col2
68 use space, only : space_t
69 use num_types, only : rp
70 use utils, only : nonlinear_index
71 use iso_c_binding, only : c_ptr
73 use device, only : host_to_device
77 use mpi_f08, only : mpi_allreduce, mpi_in_place, mpi_sum
78 use utils, only : neko_error
79 implicit none
80 private
86
87contains
97 subroutine drag_torque_zone(dgtq, tstep, zone, center, s11, s22, s33, s12, &
98 s13, s23, p, coef, visc)
99 integer, intent(in) :: tstep
100 type(facet_zone_t) :: zone
101 type(coef_t), intent(inout) :: coef
102 real(kind=rp), intent(inout) :: s11(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
103 coef%msh%nelv)
104 real(kind=rp), intent(inout) :: s22(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
105 coef%msh%nelv)
106 real(kind=rp), intent(inout) :: s33(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
107 coef%msh%nelv)
108 real(kind=rp), intent(inout) :: s12(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
109 coef%msh%nelv)
110 real(kind=rp), intent(inout) :: s13(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
111 coef%msh%nelv)
112 real(kind=rp), intent(inout) :: s23(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
113 coef%msh%nelv)
114 type(field_t), intent(inout) :: p
115 real(kind=rp), intent(in) :: visc, center(3)
116 real(kind=rp) :: dgtq(3,4)
117 real(kind=rp) :: dragpx = 0.0_rp ! pressure
118 real(kind=rp) :: dragpy = 0.0_rp
119 real(kind=rp) :: dragpz = 0.0_rp
120 real(kind=rp) :: dragvx = 0.0_rp ! viscous
121 real(kind=rp) :: dragvy = 0.0_rp
122 real(kind=rp) :: dragvz = 0.0_rp
123 real(kind=rp) :: torqpx = 0.0_rp ! pressure
124 real(kind=rp) :: torqpy = 0.0_rp
125 real(kind=rp) :: torqpz = 0.0_rp
126 real(kind=rp) :: torqvx = 0.0_rp ! viscous
127 real(kind=rp) :: torqvy = 0.0_rp
128 real(kind=rp) :: torqvz = 0.0_rp
129 real(kind=rp) :: dragx, dragy, dragz
130 integer :: ie, ifc, mem, ierr
131 dragx = 0.0
132 dragy = 0.0
133 dragz = 0.0
134
135!
136! Fill up viscous array w/ default
137!
138 dragpx = 0.0
139 dragpy = 0.0
140 dragpz = 0.0
141 dragvx = 0.0
142 dragvy = 0.0
143 dragvz = 0.0
144 do mem = 1, zone%size
145 ie = zone%facet_el(mem)%x(2)
146 ifc = zone%facet_el(mem)%x(1)
147 call drag_torque_facet(dgtq, coef%dof%x, coef%dof%y, coef%dof%z, &
148 center, s11, s22, s33, s12, s13, s23, p%x, visc, ifc, ie, coef, &
149 coef%Xh)
150
151 dragpx = dragpx + dgtq(1,1) ! pressure
152 dragpy = dragpy + dgtq(2,1)
153 dragpz = dragpz + dgtq(3,1)
154
155 dragvx = dragvx + dgtq(1,2) ! viscous
156 dragvy = dragvy + dgtq(2,2)
157 dragvz = dragvz + dgtq(3,2)
158
159 torqpx = torqpx + dgtq(1,3) ! pressure
160 torqpy = torqpy + dgtq(2,3)
161 torqpz = torqpz + dgtq(3,3)
162
163 torqvx = torqvx + dgtq(1,4) ! viscous
164 torqvy = torqvy + dgtq(2,4)
165 torqvz = torqvz + dgtq(3,4)
166 end do
167!
168! Sum contributions from all processors
169!
170 call mpi_allreduce(mpi_in_place, dragpx, 1, &
171 mpi_real_precision, mpi_sum, neko_comm, ierr)
172 call mpi_allreduce(mpi_in_place, dragpy, 1, &
173 mpi_real_precision, mpi_sum, neko_comm, ierr)
174 call mpi_allreduce(mpi_in_place, dragpz, 1, &
175 mpi_real_precision, mpi_sum, neko_comm, ierr)
176 call mpi_allreduce(mpi_in_place, dragvx, 1, &
177 mpi_real_precision, mpi_sum, neko_comm, ierr)
178 call mpi_allreduce(mpi_in_place, dragvy, 1, &
179 mpi_real_precision, mpi_sum, neko_comm, ierr)
180 call mpi_allreduce(mpi_in_place, dragvz, 1, &
181 mpi_real_precision, mpi_sum, neko_comm, ierr)
182 !Torque
183 call mpi_allreduce(mpi_in_place, torqpx, 1, &
184 mpi_real_precision, mpi_sum, neko_comm, ierr)
185 call mpi_allreduce(mpi_in_place, torqpy, 1, &
186 mpi_real_precision, mpi_sum, neko_comm, ierr)
187 call mpi_allreduce(mpi_in_place, torqpz, 1, &
188 mpi_real_precision, mpi_sum, neko_comm, ierr)
189 call mpi_allreduce(mpi_in_place, torqvx, 1, &
190 mpi_real_precision, mpi_sum, neko_comm, ierr)
191 call mpi_allreduce(mpi_in_place, torqvy, 1, &
192 mpi_real_precision, mpi_sum, neko_comm, ierr)
193 call mpi_allreduce(mpi_in_place, torqvz, 1, &
194 mpi_real_precision, mpi_sum, neko_comm, ierr)
195
196 dgtq(1,1) = dragpx ! pressure
197 dgtq(2,1) = dragpy
198 dgtq(3,1) = dragpz
199
200 dgtq(1,2) = dragvx ! viscous
201 dgtq(2,2) = dragvy
202 dgtq(3,2) = dragvz
203
204 dgtq(1,3) = torqpx ! pressure
205 dgtq(2,3) = torqpy
206 dgtq(3,3) = torqpz
207
208 dgtq(1,4) = torqvx ! viscous
209 dgtq(2,4) = torqvy
210 dgtq(3,4) = torqvz
211
212 end subroutine drag_torque_zone
213
225 subroutine drag_torque_facet(dgtq, xm0, ym0, zm0, center, s11, s22, s33, &
226 s12, s13, s23, pm1, visc, f, e, coef, Xh)
227 type(coef_t), intent(in) :: coef
228 type(space_t), intent(in) :: xh
229 real(kind=rp), intent(out) :: dgtq(3,4)
230 real(kind=rp), intent(in) :: center(3)
231 real(kind=rp), intent(in) :: xm0(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
232 real(kind=rp), intent(in) :: ym0(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
233 real(kind=rp), intent(in) :: zm0(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
234 real(kind=rp), intent(in) :: s11(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
235 real(kind=rp), intent(in) :: s22(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
236 real(kind=rp), intent(in) :: s33(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
237 real(kind=rp), intent(in) :: s12(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
238 real(kind=rp), intent(in) :: s13(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
239 real(kind=rp), intent(in) :: s23(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
240 real(kind=rp), intent(in) :: pm1(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
241 real(kind=rp), intent(in) :: visc
242 integer, intent(in) :: f, e
243 integer :: pf, i, j1, j2
244 real(kind=rp) :: n1, n2, n3, a, v, dgtq_i(3,4)
245 integer :: skpdat(6,6), nx, ny, nz
246 integer :: js1
247 integer :: jf1
248 integer :: jskip1
249 integer :: js2
250 integer :: jf2
251 integer :: jskip2
252 real(kind=rp) :: s11_, s12_, s22_, s13_, s23_, s33_
253
254
255 nx = xh%lx
256 ny = xh%ly
257 nz = xh%lz
258 skpdat(1,1) = 1
259 skpdat(2,1) = nx*(ny-1) + 1
260 skpdat(3,1) = nx
261 skpdat(4,1) = 1
262 skpdat(5,1) = ny*(nz-1) + 1
263 skpdat(6,1) = ny
264
265 skpdat(1,2) = 1 + (nx-1)
266 skpdat(2,2) = nx*(ny-1)+1 + (nx-1)
267 skpdat(3,2) = nx
268 skpdat(4,2) = 1
269 skpdat(5,2) = ny*(nz-1)+1
270 skpdat(6,2) = ny
271
272 skpdat(1,3) = 1
273 skpdat(2,3) = nx
274 skpdat(3,3) = 1
275 skpdat(4,3) = 1
276 skpdat(5,3) = ny*(nz-1)+1
277 skpdat(6,3) = ny
278
279 skpdat(1,4) = 1 + nx*(ny-1)
280 skpdat(2,4) = nx + nx*(ny-1)
281 skpdat(3,4) = 1
282 skpdat(4,4) = 1
283 skpdat(5,4) = ny*(nz-1)+1
284 skpdat(6,4) = ny
285
286 skpdat(1,5) = 1
287 skpdat(2,5) = nx
288 skpdat(3,5) = 1
289 skpdat(4,5) = 1
290 skpdat(5,5) = ny
291 skpdat(6,5) = 1
292
293 skpdat(1,6) = 1 + nx*ny*(nz-1)
294 skpdat(2,6) = nx + nx*ny*(nz-1)
295 skpdat(3,6) = 1
296 skpdat(4,6) = 1
297 skpdat(5,6) = ny
298 skpdat(6,6) = 1
299 pf = f
300 js1 = skpdat(1, pf)
301 jf1 = skpdat(2, pf)
302 jskip1 = skpdat(3, pf)
303 js2 = skpdat(4, pf)
304 jf2 = skpdat(5, pf)
305 jskip2 = skpdat(6, pf)
306 call rzero(dgtq,12)
307 i = 0
308 a = 0
309 do j2 = js2, jf2, jskip2
310 do j1 = js1, jf1, jskip1
311 i = i + 1
312 n1 = coef%nx(i,1,f,e) * coef%area(i,1,f,e)
313 n2 = coef%ny(i,1,f,e) * coef%area(i,1,f,e)
314 n3 = coef%nz(i,1,f,e) * coef%area(i,1,f,e)
315 a = a + coef%area(i,1,f,e)
316 v = visc
317 s11_ = s11(j1, j2, 1, e)
318 s12_ = s12(j1, j2, 1, e)
319 s22_ = s22(j1, j2, 1, e)
320 s13_ = s13(j1, j2, 1, e)
321 s23_ = s23(j1, j2, 1, e)
322 s33_ = s33(j1, j2, 1, e)
323 call drag_torque_pt(dgtq_i, xm0(j1, j2, 1, e), ym0(j1, j2, 1, e), &
324 zm0(j1, j2, 1, e), center, s11_, s22_, s33_, s12_, s13_, s23_, &
325 pm1(j1, j2, 1, e), n1, n2, n3, v)
326 dgtq = dgtq + dgtq_i
327 end do
328 end do
329 end subroutine drag_torque_facet
330
343 subroutine drag_torque_pt(dgtq, x, y, z, center, s11, s22, s33, s12, s13, &
344 s23, p, n1, n2, n3, v)
345 real(kind=rp), intent(inout) :: dgtq(3, 4)
346 real(kind=rp), intent(in) :: x
347 real(kind=rp), intent(in) :: y
348 real(kind=rp), intent(in) :: z
349 real(kind=rp), intent(in) :: p
350 real(kind=rp), intent(in) :: v
351 real(kind=rp), intent(in) :: n1, n2, n3, center(3)
352 real(kind=rp), intent(in) :: s11, s12, s22, s13, s23, s33
353 real(kind=rp) :: s21, s31, s32, r1, r2, r3
354
355 call rzero(dgtq, 12)
356 s21 = s12
357 s32 = s23
358 s31 = s13
359 !pressure drag
360 dgtq(1,1) = p*n1
361 dgtq(2,1) = p*n2
362 dgtq(3,1) = p*n3
363 ! viscous drag
364 dgtq(1,2) = -2*v*(s11*n1 + s12*n2 + s13*n3)
365 dgtq(2,2) = -2*v*(s21*n1 + s22*n2 + s23*n3)
366 dgtq(3,2) = -2*v*(s31*n1 + s32*n2 + s33*n3)
367 r1 = x - center(1)
368 r2 = y - center(2)
369 r3 = z - center(3)
370 !pressure torque
371 dgtq(1,3) = r2*dgtq(3,1) - r3*dgtq(2,1)
372 dgtq(2,3) = r3*dgtq(1,1) - r1*dgtq(3,1)
373 dgtq(3,3) = r1*dgtq(2,1) - r2*dgtq(1,1)
374 !viscous torque
375 dgtq(1,4) = r2*dgtq(3,2) - r3*dgtq(2,2)
376 dgtq(2,4) = r3*dgtq(1,2) - r1*dgtq(3,2)
377 dgtq(3,4) = r1*dgtq(2,2) - r2*dgtq(1,2)
378 end subroutine drag_torque_pt
379
389 subroutine calc_force_array(force1, force2, force3, force4, force5, force6,&
390 s11, s22, s33, s12, s13, s23,&
391 p, n1, n2, n3, mu, n_pts)
392 integer :: n_pts
393 real(kind=rp), intent(inout), dimension(n_pts) :: force1, force2, force3
394 real(kind=rp), intent(inout), dimension(n_pts) :: force4, force5, force6
395 real(kind=rp), intent(in) :: p(n_pts)
396 real(kind=rp), intent(in) :: mu(n_pts)
397 real(kind=rp), intent(in) :: n1(n_pts), n2(n_pts), n3(n_pts)
398 real(kind=rp), intent(in), dimension(n_pts) :: s11, s12, s22, s13, s23, s33
399 real(kind=rp) :: v2(n_pts)
400 call rzero(force4, n_pts)
401 call rzero(force5, n_pts)
402 call rzero(force6, n_pts)
403 !pressure force
404 call col3(force1, p, n1, n_pts)
405 call col3(force2, p, n2, n_pts)
406 call col3(force3, p, n3, n_pts)
407 ! viscous force
408 v2 = -2.0_rp*mu
409 call vdot3(force4, s11, s12, s13, n1, n2, n3, n_pts)
410 call vdot3(force5, s12, s22, s23, n1, n2, n3, n_pts)
411 call vdot3(force6, s13, s23, s33, n1, n2, n3, n_pts)
412 call col2(force4, v2, n_pts)
413 call col2(force5, v2, n_pts)
414 call col2(force6, v2, n_pts)
415 end subroutine calc_force_array
416
426 subroutine device_calc_force_array(force1, force2, force3,&
427 force4, force5, force6,&
428 s11, s22, s33, s12, s13, s23,&
429 p, n1, n2, n3, mu, n_pts)
430 integer :: n_pts
431 type(c_ptr) :: force1, force2, force3
432 type(c_ptr) :: force4, force5, force6
433 type(c_ptr) :: p, n1, n2, n3
434 type(c_ptr) :: s11, s12, s22, s13, s23, s33
435 type(c_ptr) :: mu
436 if (neko_bcknd_device .eq. 1) then
437 call device_rzero(force4, n_pts)
438 call device_rzero(force5, n_pts)
439 call device_rzero(force6, n_pts)
440 !pressure force
441 call device_col3(force1, p, n1, n_pts)
442 call device_col3(force2, p, n2, n_pts)
443 call device_col3(force3, p, n3, n_pts)
444
445 call device_vdot3(force4, s11, s12, s13, n1, n2, n3, n_pts)
446 call device_vdot3(force5, s12, s22, s23, n1, n2, n3, n_pts)
447 call device_vdot3(force6, s13, s23, s33, n1, n2, n3, n_pts)
448 call device_cmult(force4, -2.0_rp, n_pts)
449 call device_cmult(force5, -2.0_rp, n_pts)
450 call device_cmult(force6, -2.0_rp, n_pts)
451 call device_col2(force4, mu, n_pts)
452 call device_col2(force5, mu, n_pts)
453 call device_col2(force6, mu, n_pts)
454 else
455 call neko_error('error in drag_torque, no device bcklnd configured')
456 end if
457 end subroutine device_calc_force_array
458
459
469 subroutine setup_normals(coef, mask, facets, n1, n2, n3, n_pts)
470 type(coef_t) :: coef
471 integer :: n_pts
472 type(vector_t), intent(inout) :: n1, n2, n3
473 integer :: mask(0:n_pts), facets(0:n_pts), fid, idx(4)
474 real(kind=rp) :: normal(3), area(3)
475 integer :: i
476
477 do i = 1, n_pts
478 fid = facets(i)
479 idx = nonlinear_index(mask(i), coef%Xh%lx, coef%Xh%lx, coef%Xh%lx)
480 normal = coef%get_normal(idx(1), idx(2), idx(3), idx(4), fid)
481 area = coef%get_area(idx(1), idx(2), idx(3), idx(4), fid)
482 n1%x(i) = normal(1)*area(1)
483 n2%x(i) = normal(2)*area(2)
484 n3%x(i) = normal(3)*area(3)
485 end do
486
487 if (neko_bcknd_device .eq. 1) then
488 call n1%copy_from(host_to_device, .false.)
489 call n2%copy_from(host_to_device, .false.)
490 call n3%copy_from(host_to_device, .true.)
491 end if
492 end subroutine setup_normals
493
494end module drag_torque
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:53
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:45
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
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_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:48
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 setup_normals(coef, mask, facets, n1, n2, n3, n_pts)
Computes the normals for a given set of boundary points accessed by the mask.
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 calc_force_array(force1, force2, force3, force4, force5, force6, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, mu, n_pts)
Calculate drag and torque from array of points.
subroutine, public device_calc_force_array(force1, force2, force3, force4, force5, force6, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, mu, n_pts)
Calculate drag and torque from array of points.
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 a zone as a subset of facets in a mesh.
Defines a field.
Definition field.f90:34
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:1044
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:1059
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:850
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:233
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
Utilities.
Definition utils.f90:35
Defines a vector.
Definition vector.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
The function space for the SEM solution fields.
Definition space.f90:63