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