-
Notifications
You must be signed in to change notification settings - Fork 8
/
Domains.jl
734 lines (661 loc) · 24.4 KB
/
Domains.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
module Domains
using ClimaCore
using ClimaComms
using IntervalSets
using DocStringExtensions
import ClimaCore: Meshes, Spaces, Topologies, Geometry
import ClimaCore.Meshes: Uniform
### General type and methods all model domains will need
"""
AbstractDomain{FT <:AbstractFloat}
An abstract type for domains.
The domain structs typically hold information
regarding the bounds of the domain, the boundary condition type
(periodic or not), and the spatial discretization.
Additionally, the domain struct holds the relevant spaces for that
domain. For example, a 3D domain holds the center space (in terms of
finite difference - the space corresponding to the centers of each
element), and the top face space where surface fluxes are computed.
"""
abstract type AbstractDomain{FT <: AbstractFloat} end
Base.eltype(::AbstractDomain{FT}) where {FT} = FT
"""
coordinates(domain::AbstractDomain)
Returns the coordinate fields for the domain as a NamedTuple.
The returned coordinates are stored with keys :surface, :subsurface, e.g.
as relevant for the domain.
"""
function coordinates(domain::AbstractDomain)
domain_names = propertynames(domain.space)
domain_spaces = getproperty.(Ref(domain.space), domain_names)
return NamedTuple{domain_names}(
ClimaCore.Fields.coordinate_field.(domain_spaces),
)
end
"""
Point{FT} <: AbstractDomain{FT}
A domain for single column surface variables.
For models such as ponds, snow, plant hydraulics, etc. Enables consistency
in variable initialization across all domains.
`space` is a NamedTuple holding the surface space (in this case,
the Point space).
# Fields
$(DocStringExtensions.FIELDS)
"""
struct Point{FT} <: AbstractDomain{FT}
"Surface elevation relative to a reference (m)"
z_sfc::FT
"A NamedTuple of associated ClimaCore spaces: in this case, the Point (surface) space"
space::NamedTuple
end
"""
Point(;z_sfc::FT,
comms = ClimaComms.SingletonCommsContext()
) where {FT}
Constructor for the `Point` domain using keyword arguments.
All other ClimaLand domains rely on default comms set internally
by ClimaCore. However, the Point space is unique in this context,
and does not have the same default defined in ClimaCore.
Because of this, we set the default here
in ClimaLand. In long term, we will repeat the same for all ClimaLand domains
and not rely on any internal defaults set in ClimaCore.
"""
function Point(;
z_sfc::FT,
comms = ClimaComms.SingletonCommsContext(),
) where {FT}
coord = ClimaCore.Geometry.ZPoint(z_sfc)
space = (; surface = ClimaCore.Spaces.PointSpace(comms, coord))
return Point{FT}(z_sfc, space)
end
"""
Column{FT} <: AbstractDomain{FT}
A struct holding the necessary information
to construct a domain, a mesh, a center and face
space, etc. for use when a finite difference in
1D is suitable, as for a soil column model.
`space` is a NamedTuple holding the surface space (in this case,
the top face space) and the center space for the subsurface.
These are stored using the keys :surface and :subsurface.
# Fields
$(DocStringExtensions.FIELDS)
"""
struct Column{FT} <: AbstractDomain{FT}
"Domain interval limits, (zmin, zmax), in meters"
zlim::Tuple{FT, FT}
"Number of elements used to discretize the interval"
nelements::Tuple{Int}
"Tuple for mesh stretching specifying *target* (dz_bottom, dz_top) (m). If nothing, no stretching is applied."
dz_tuple::Union{Tuple{FT, FT}, Nothing}
"Boundary face identifiers"
boundary_names::Tuple{Symbol, Symbol}
"A NamedTuple of associated ClimaCore spaces: in this case, the surface space and subsurface center space"
space::NamedTuple
end
"""
Column(;
zlim::Tuple{FT, FT},
nelements::Int,
dz_tuple::Union{Tuple{FT, FT}, Nothing} = nothing) where {FT}
Outer constructor for the `Column` type.
Using `ClimaCore` tools, the coordinate
mesh can be stretched such that the top of the domain has finer resolution
than the bottom of the domain. In order to activate this, a tuple with the
target dz_bottom and dz_top should be passed via keyword argument. The default is
uniform spacing. Please note that in order to use this feature, ClimaCore requires that
the elements of zlim be <=0. Additionally, the dz_tuple you supply may not be compatible
with the domain boundaries in some cases, in which case you may need to choose
different values.
The `boundary_names` field values are used to label the boundary faces
at the top and bottom of the domain.
"""
function Column(;
zlim::Tuple{FT, FT},
nelements::Int,
dz_tuple::Union{Tuple{FT, FT}, Nothing} = nothing,
) where {FT}
@assert zlim[1] < zlim[2]
boundary_names = (:bottom, :top)
column = ClimaCore.Domains.IntervalDomain(
ClimaCore.Geometry.ZPoint{FT}(zlim[1]),
ClimaCore.Geometry.ZPoint{FT}(zlim[2]);
boundary_names = boundary_names,
)
if dz_tuple isa Nothing
mesh = ClimaCore.Meshes.IntervalMesh(column; nelems = nelements)
else
@assert zlim[2] <= 0
mesh = ClimaCore.Meshes.IntervalMesh(
column,
ClimaCore.Meshes.GeneralizedExponentialStretching{FT}(
dz_tuple[1],
dz_tuple[2],
);
nelems = nelements,
reverse_mode = true,
)
end
subsurface_space = ClimaCore.Spaces.CenterFiniteDifferenceSpace(mesh)
surface_space = obtain_surface_space(subsurface_space)
space = (; surface = surface_space, subsurface = subsurface_space)
return Column{FT}(zlim, (nelements,), dz_tuple, boundary_names, space)
end
"""
Plane{FT} <: AbstractDomain{FT}
A struct holding the necessary information
to construct a domain, a mesh, a 2d spectral
element space, and the resulting coordinate field.
Note that only periodic domains are currently supported.
`space` is a NamedTuple holding the surface space (in this case,
the entire Plane space).
# Fields
$(DocStringExtensions.FIELDS)
"""
struct Plane{FT} <: AbstractDomain{FT}
"Domain interval limits along x axis, in meters"
xlim::Tuple{FT, FT}
"Domain interval limits along y axis, in meters"
ylim::Tuple{FT, FT}
"Number of elements to discretize interval, (nx, ny)"
nelements::Tuple{Int, Int}
"Flags for periodic boundaries; only true is supported"
periodic::Tuple{Bool, Bool}
"Polynomial order for both x and y"
npolynomial::Int
"A NamedTuple of associated ClimaCore spaces: in this case, the surface(Plane) space"
space::NamedTuple
end
"""
Plane(;
xlim::Tuple{FT,FT},
ylim::Tuple{FT,FT},
nelements::Tuple{Int,Int},
periodic::Tuple{Bool,Bool},
npolynomial::Int,
comms_ctx = ClimaComms.SingletonCommsContext(),
) where {FT}
Outer constructor for the `Plane` domain, using keyword arguments.
"""
function Plane(;
xlim::Tuple{FT, FT},
ylim::Tuple{FT, FT},
nelements::Tuple{Int, Int},
periodic::Tuple{Bool, Bool} = (true, true),
npolynomial::Int,
comms_ctx = ClimaComms.SingletonCommsContext(),
) where {FT}
@assert xlim[1] < xlim[2]
@assert ylim[1] < ylim[2]
@assert periodic == (true, true)
domain_x = ClimaCore.Domains.IntervalDomain(
ClimaCore.Geometry.XPoint(xlim[1]),
ClimaCore.Geometry.XPoint(xlim[2]);
periodic = periodic[1],
)
domain_y = ClimaCore.Domains.IntervalDomain(
ClimaCore.Geometry.YPoint(ylim[1]),
ClimaCore.Geometry.YPoint(ylim[2]);
periodic = periodic[2],
)
plane = ClimaCore.Domains.RectangleDomain(domain_x, domain_y)
mesh = ClimaCore.Meshes.RectilinearMesh(plane, nelements[1], nelements[2])
grid_topology = ClimaCore.Topologies.Topology2D(comms_ctx, mesh)
if npolynomial == 0
quad = ClimaCore.Spaces.Quadratures.GL{npolynomial + 1}()
else
quad = ClimaCore.Spaces.Quadratures.GLL{npolynomial + 1}()
end
space = ClimaCore.Spaces.SpectralElementSpace2D(grid_topology, quad)
space = (; surface = space)
return Plane{FT}(xlim, ylim, nelements, periodic, npolynomial, space)
end
"""
struct HybridBox{FT} <: AbstractDomain{FT}
xlim::Tuple{FT, FT}
ylim::Tuple{FT, FT}
zlim::Tuple{FT, FT}
dz_tuple::Union{Tuple{FT, FT}, Nothing}
nelements::Tuple{Int, Int, Int}
npolynomial::Int
periodic::Tuple{Bool, Bool}
end
A struct holding the necessary information to construct a domain, a mesh,
a 2d spectral element space (horizontal) x a 1d finite difference space
(vertical), and the resulting coordinate field.
This domain is not periodic along the z-axis. Note that
only periodic domains are supported
in the horizontal.
`space` is a NamedTuple holding the surface space (in this case,
the top face space) and the center space for the subsurface.
These are stored using the keys :surface and :subsurface.
# Fields
$(DocStringExtensions.FIELDS)
"""
struct HybridBox{FT} <: AbstractDomain{FT}
"Domain interval limits along x axis, in meters"
xlim::Tuple{FT, FT}
"Domain interval limits along y axis, in meters"
ylim::Tuple{FT, FT}
"Domain interval limits along z axis, in meters"
zlim::Tuple{FT, FT}
"Tuple for mesh stretching specifying *target* (dz_bottom, dz_top) (m). If nothing, no stretching is applied."
dz_tuple::Union{Tuple{FT, FT}, Nothing}
"Number of elements to discretize interval, (nx, ny,nz)"
nelements::Tuple{Int, Int, Int}
" Polynomial order for the horizontal directions"
npolynomial::Int
"Flag indicating periodic boundaries in horizontal. only true is supported"
periodic::Tuple{Bool, Bool}
"A NamedTuple of associated ClimaCore spaces: in this case, the surface space and subsurface center space"
space::NamedTuple
end
"""
HybridBox(;
xlim::Tuple{FT, FT},
ylim::Tuple{FT, FT},
zlim::Tuple{FT, FT},
nelements::Tuple{Int, Int, Int},
npolynomial::Int,
dz_tuple::Union{Tuple{FT, FT}, Nothing} = nothing,
periodic = (true, true),
) where {FT}
Constructs the `HybridBox` domain
with limits `xlim` `ylim` and `zlim
(where `xlim[1] < xlim[2]`,`ylim[1] < ylim[2]`, and `zlim[1] < zlim[2]`),
`nelements` must be a tuple with three values, with the first
value corresponding
to the x-axis, the second corresponding to the y-axis, and the third
corresponding to the z-axis. The domain is periodic at the (xy) boundaries,
and the function space is of polynomial order `npolynomial` in the
horizontal directions.
Using `ClimaCore` tools, the coordinate
mesh can be stretched such that the top of the domain has finer resolution
than the bottom of the domain. In order to activate this, a tuple with the
target dz_bottom and dz_top should be passed via keyword argument. The default is
uniform spacing. Please note that in order to use this feature, ClimaCore requires that
the elements of zlim be <=0. Additionally, the dz_tuple you supply may not be compatible
with the domain boundaries in some cases, in which case you may need to choose
different values.
"""
function HybridBox(;
xlim::Tuple{FT, FT},
ylim::Tuple{FT, FT},
zlim::Tuple{FT, FT},
nelements::Tuple{Int, Int, Int},
npolynomial::Int,
dz_tuple::Union{Tuple{FT, FT}, Nothing} = nothing,
periodic = (true, true),
) where {FT}
@assert xlim[1] < xlim[2]
@assert ylim[1] < ylim[2]
@assert zlim[1] < zlim[2]
@assert periodic == (true, true)
vertdomain = ClimaCore.Domains.IntervalDomain(
ClimaCore.Geometry.ZPoint(zlim[1]),
ClimaCore.Geometry.ZPoint(zlim[2]);
boundary_names = (:bottom, :top),
)
if dz_tuple isa Nothing
vertmesh =
ClimaCore.Meshes.IntervalMesh(vertdomain; nelems = nelements[3])
else
@assert zlim[2] <= 0
vertmesh = ClimaCore.Meshes.IntervalMesh(
vertdomain,
ClimaCore.Meshes.GeneralizedExponentialStretching{FT}(
dz_tuple[1],
dz_tuple[2],
);
nelems = nelements[3],
reverse_mode = true,
)
end
vert_center_space = ClimaCore.Spaces.CenterFiniteDifferenceSpace(vertmesh)
horzdomain = Plane(;
xlim = xlim,
ylim = ylim,
nelements = nelements[1:2],
periodic = periodic,
npolynomial = npolynomial,
)
horzspace = horzdomain.space.surface
subsurface_space = ClimaCore.Spaces.ExtrudedFiniteDifferenceSpace(
horzspace,
vert_center_space,
)
surface_space = obtain_surface_space(subsurface_space)
space = (; surface = surface_space, subsurface = subsurface_space)
return HybridBox{FT}(
xlim,
ylim,
zlim,
dz_tuple,
nelements,
npolynomial,
periodic,
space,
)
end
"""
struct SphericalShell{FT} <: AbstractDomain{FT}
radius::FT
depth::FT
dz_tuple::Union{Tuple{FT, FT}, Nothing}
nelements::Tuple{Int, Int}
npolynomial::Int
end
A struct holding the necessary information to construct a domain, a mesh,
a 2d spectral element space (non-radial directions)
x a 1d finite difference space (radial direction),
and the resulting coordinate field.
`space` is a NamedTuple holding the surface space (in this case,
the top face space) and the center space for the subsurface.
These are stored using the keys :surface and :subsurface.
# Fields
$(DocStringExtensions.FIELDS)
"""
struct SphericalShell{FT} <: AbstractDomain{FT}
"The radius of the shell"
radius::FT
"The radial extent of the shell"
depth::FT
"Tuple for mesh stretching specifying *target* (dz_bottom, dz_top) (m). If nothing, no stretching is applied."
dz_tuple::Union{Tuple{FT, FT}, Nothing}
"The number of elements to be used in the non-radial and radial directions"
nelements::Tuple{Int, Int}
"The polynomial order to be used in the non-radial directions"
npolynomial::Int
"A NamedTuple of associated ClimaCore spaces: in this case, the surface space and subsurface center space"
space::NamedTuple
end
"""
SphericalShell(;
radius::FT,
depth::FT,
nelements::Tuple{Int, Int},
npolynomial::Int,
dz_tuple::Union{Tuple{FT, FT}, Nothing} = nothing,
comms_ctx = ClimaComms.SingletonCommsContext(),
) where {FT}
Outer constructor for the `SphericalShell` domain, using keyword arguments.
Using `ClimaCore` tools, the coordinate
mesh can be stretched such that the top of the domain has finer resolution
than the bottom of the domain. In order to activate this, a tuple with the
target dz_bottom and dz_top should be passed via keyword argument. The default is
uniform spacing. Please note that the dz_tuple you supply may not be compatible
with the depth/nelements chosen, in which case you may need to choose
different values.
"""
function SphericalShell(;
radius::FT,
depth::FT,
nelements::Tuple{Int, Int},
npolynomial::Int,
dz_tuple::Union{Tuple{FT, FT}, Nothing} = nothing,
comms_ctx = ClimaComms.SingletonCommsContext(),
) where {FT}
@assert 0 < radius
@assert 0 < depth
vertdomain = ClimaCore.Domains.IntervalDomain(
ClimaCore.Geometry.ZPoint(FT(-depth)),
ClimaCore.Geometry.ZPoint(FT(0));
boundary_names = (:bottom, :top),
)
if dz_tuple isa Nothing
vertmesh =
ClimaCore.Meshes.IntervalMesh(vertdomain; nelems = nelements[2])
else
vertmesh = ClimaCore.Meshes.IntervalMesh(
vertdomain,
ClimaCore.Meshes.GeneralizedExponentialStretching{FT}(
dz_tuple[1],
dz_tuple[2],
);
nelems = nelements[2],
reverse_mode = true,
)
end
vert_center_space = ClimaCore.Spaces.CenterFiniteDifferenceSpace(vertmesh)
horzdomain = ClimaCore.Domains.SphereDomain(radius)
horzmesh = ClimaCore.Meshes.EquiangularCubedSphere(horzdomain, nelements[1])
horztopology = ClimaCore.Topologies.Topology2D(comms_ctx, horzmesh)
quad = ClimaCore.Spaces.Quadratures.GLL{npolynomial + 1}()
horzspace = ClimaCore.Spaces.SpectralElementSpace2D(horztopology, quad)
subsurface_space = ClimaCore.Spaces.ExtrudedFiniteDifferenceSpace(
horzspace,
vert_center_space,
)
surface_space = obtain_surface_space(subsurface_space)
space = (; surface = surface_space, subsurface = subsurface_space)
return SphericalShell{FT}(
radius,
depth,
dz_tuple,
nelements,
npolynomial,
space,
)
end
"""
struct SphericalSurface{FT} <: AbstractDomain{FT}
radius::FT
nelements::Tuple{Int, Int}
npolynomial::Int
end
A struct holding the necessary information to construct a domain, a mesh,
a 2d spectral element space (non-radial directions) and the resulting coordinate field.
`space` is a NamedTuple holding the surface space (in this case,
the entire SphericalSurface space).
# Fields
$(DocStringExtensions.FIELDS)
"""
struct SphericalSurface{FT} <: AbstractDomain{FT}
"The radius of the surface"
radius::FT
"The number of elements to be used in the non-radial directions"
nelements::Int
"The polynomial order to be used in the non-radial directions"
npolynomial::Int
"A NamedTuple of associated ClimaCore spaces: in this case, the surface (SphericalSurface) space"
space::NamedTuple
end
"""
SphericalSurface(;
radius::FT,
nelements::Int
npolynomial::Int,
comms_ctx = ClimaComms.SingletonCommsContext(),
) where {FT}
Outer constructor for the `SphericalSurface` domain, using keyword arguments.
"""
function SphericalSurface(;
radius::FT,
nelements::Int,
npolynomial::Int,
comms_ctx = ClimaComms.SingletonCommsContext(),
) where {FT}
@assert 0 < radius
horzdomain = ClimaCore.Domains.SphereDomain(radius)
horzmesh = Meshes.EquiangularCubedSphere(horzdomain, nelements)
horztopology = Topologies.Topology2D(comms_ctx, horzmesh)
quad = Spaces.Quadratures.GLL{npolynomial + 1}()
horzspace = Spaces.SpectralElementSpace2D(horztopology, quad)
space = (; surface = horzspace)
return SphericalSurface{FT}(radius, nelements, npolynomial, space)
end
"""
obtain_surface_domain(d::AbstractDomain) where {FT}
Default method throwing an error; any domain with a corresponding
domain should define a new method of this function.
"""
obtain_surface_domain(d::AbstractDomain) =
@error("No surface domain is defined for this domain.")
"""
obtain_surface_domain(c::Column{FT}) where {FT}
Returns the Point domain corresponding to the top face (surface) of the
Column domain `c`.
"""
function obtain_surface_domain(c::Column{FT}) where {FT}
surface_domain = Point{FT}(c.zlim[2], (; surface = c.space.surface))
return surface_domain
end
"""
obtain_surface_domain(b::HybridBox{FT}) where {FT}
Returns the Plane domain corresponding to the top face (surface) of the
HybridBox domain `b`.
"""
function obtain_surface_domain(b::HybridBox{FT}) where {FT}
surface_domain = Plane{FT}(
b.xlim,
b.ylim,
b.nelements[1:2],
b.periodic,
b.npolynomial,
(; surface = b.space.surface),
)
return surface_domain
end
"""
obtain_surface_domain(s::SphericalShell{FT}) where {FT}
Returns the SphericalSurface domain corresponding to the top face
(surface) of the SphericalShell domain `s`.
"""
function obtain_surface_domain(s::SphericalShell{FT}) where {FT}
surface_domain = SphericalSurface{FT}(
s.radius,
s.nelements[1],
s.npolynomial,
(; surface = s.space.surface),
)
return surface_domain
end
"""
obtain_face_space(cs::ClimaCore.Spaces.AbstractSpace)
Returns the face space, if applicable, for the center space `cs`.
"""
obtain_face_space(cs::ClimaCore.Spaces.AbstractSpace) =
@error("No face space is defined for this space.")
"""
obtain_surface_space(cs::ClimaCore.Spaces.AbstractSpace)
Returns the surface space, if applicable, for the center space `cs`.
"""
obtain_surface_space(cs::ClimaCore.Spaces.AbstractSpace) =
@error("No surface space is defined for this space.")
"""
obtain_face_space(cs::ClimaCore.Spaces.CenterExtrudedFiniteDifferenceSpace)
Returns the face space for the CenterExtrudedFiniteDifferenceSpace `cs`.
"""
function obtain_face_space(
cs::ClimaCore.Spaces.CenterExtrudedFiniteDifferenceSpace,
)
return ClimaCore.Spaces.FaceExtrudedFiniteDifferenceSpace(cs)
end
"""
obtain_face_space(cs::ClimaCore.Spaces.CenterFiniteDifferenceSpace)
Returns the face space corresponding to the CenterFiniteDifferenceSpace `cs`.
"""
function obtain_face_space(cs::ClimaCore.Spaces.CenterFiniteDifferenceSpace)
return ClimaCore.Spaces.FaceFiniteDifferenceSpace(cs)
end
"""
obtain_surface_space(cs::ClimaCore.Spaces.CenterExtrudedFiniteDifferenceSpace)
Returns the horizontal space for the CenterExtrudedFiniteDifferenceSpace `cs`.
"""
function obtain_surface_space(
cs::ClimaCore.Spaces.CenterExtrudedFiniteDifferenceSpace,
)
return Spaces.horizontal_space(cs)
end
"""
obtain_surface_space(cs::ClimaCore.Spaces.CenterFiniteDifferenceSpace)
Returns the top level of the face space corresponding to the CenterFiniteDifferenceSpace `cs`.
"""
function obtain_surface_space(cs::ClimaCore.Spaces.CenterFiniteDifferenceSpace)
fs = obtain_face_space(cs)
return ClimaCore.Spaces.level(
fs,
ClimaCore.Utilities.PlusHalf(ClimaCore.Spaces.nlevels(fs) - 1),
)
end
"""
top_center_to_surface(center_field::ClimaCore.Fields.Field)
Creates and returns a ClimaCore.Fields.Field defined on the space
corresponding to the surface of the space on which `center_field`
is defined, with values equal to the those at the level of the top
center.
For example, given a `center_field` defined on 1D center finite difference space,
this would return a field defined on the Point space of the surface of
the column. The value would be the value of the oroginal `center_field`
at the topmost location. Given a `center_field` defined on a 3D
extruded center finite difference space, this would return a 2D field
corresponding to the surface, with values equal to the topmost level.
"""
function top_center_to_surface(center_field::ClimaCore.Fields.Field)
cs = axes(center_field)
face_space = obtain_face_space(cs)
N = ClimaCore.Spaces.nlevels(face_space)
interp_c2f = ClimaCore.Operators.InterpolateC2F(
top = ClimaCore.Operators.Extrapolate(),
bottom = ClimaCore.Operators.Extrapolate(),
)
surface_space = obtain_surface_space(cs)
return ClimaCore.Fields.Field(
ClimaCore.Fields.field_values(
ClimaCore.Fields.level(
interp_c2f.(center_field),
ClimaCore.Utilities.PlusHalf(N - 1),
),
),
surface_space,
)
end
"""
linear_interpolation_to_surface!(sfc_field, center_field, z)
Linearly interpolate the center field `center_field` to the surface
defined by the top face coordinate of `z`; updates the `sfc_field`
on the surface (face) space in place.
"""
function linear_interpolation_to_surface!(sfc_field, center_field, z)
Δz_top, _ = get_Δz(z)
surface_space = obtain_surface_space(axes(z))
Δz_top = ClimaCore.Fields.field_values(Δz_top)
nz = Spaces.nlevels(axes(center_field))
f1 = ClimaCore.Fields.field_values(ClimaCore.Fields.level(center_field, nz))
f2 = ClimaCore.Fields.field_values(
ClimaCore.Fields.level(center_field, nz - 1),
)
z1 = ClimaCore.Fields.field_values(ClimaCore.Fields.level(z, nz))
z2 = ClimaCore.Fields.field_values(ClimaCore.Fields.level(z, nz - 1))
sfc_field .= ClimaCore.Fields.Field(
(@. (f1 - f2) / (z1 - z2) * (Δz_top + z1 - z2) + f2),
surface_space,
)
end
"""
get_Δz(z::ClimaCore.Fields.Field)
A function to return a tuple containing the distance between the top boundary
and its closest center, and the bottom boundary and its closest center,
both as Fields.
"""
function get_Δz(z::ClimaCore.Fields.Field)
# Extract the differences between levels of the face space
fs = obtain_face_space(axes(z))
z_face = ClimaCore.Fields.coordinate_field(fs).z
Δz = ClimaCore.Fields.Δz_field(z_face)
Δz_top = ClimaCore.Fields.level(
Δz,
ClimaCore.Utilities.PlusHalf(ClimaCore.Spaces.nlevels(fs) - 1),
)
Δz_bottom = ClimaCore.Fields.level(Δz, ClimaCore.Utilities.PlusHalf(0))
return Δz_top ./ 2, Δz_bottom ./ 2
end
export AbstractDomain
export Column, Plane, HybridBox, Point, SphericalShell, SphericalSurface
export coordinates,
obtain_face_space,
obtain_surface_space,
top_center_to_surface,
obtain_surface_domain,
linear_interpolation_to_surface!,
get_Δz
end