-
Notifications
You must be signed in to change notification settings - Fork 196
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
(0.94.4) PolarBoundaryCondition
for tracer + views for AbstractOperations
#3953
Open
simone-silvestri
wants to merge
34
commits into
main
Choose a base branch
from
ss/top-boundary-condition
base: main
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
Open
Changes from 10 commits
Commits
Show all changes
34 commits
Select commit
Hold shift + click to select a range
9dddf39
add this check
simone-silvestri 95e6612
flux boundary conditions not working
simone-silvestri 3b03288
default polar BC for latitudelongitude grid
simone-silvestri ce3fef1
polar boundary conditions
simone-silvestri 83a23e1
update
simone-silvestri 09ce40b
update
simone-silvestri 755b3bc
will this work?
simone-silvestri e2dbac9
Merge branch 'ss/top-boundary-condition' of github.com:CliMA/Oceanani…
simone-silvestri 4c13313
update boundary conditions
simone-silvestri 8004cb5
add offsetarrays
simone-silvestri 575f8da
Merge branch 'main' into ss/top-boundary-condition
simone-silvestri f6f3dd1
indices in reductions
simone-silvestri 621f1e6
Merge branch 'ss/top-boundary-condition' of github.com:CliMA/Oceanani…
simone-silvestri 34dfe63
add a test
simone-silvestri 7f7b86b
hmmm
simone-silvestri 9add855
Merge branch 'main' into ss/top-boundary-condition
simone-silvestri 929bc24
bugfix
simone-silvestri 47b398f
another bugfix
simone-silvestri 091c6cb
extending `view`s for operations
simone-silvestri 6b0c234
should be ready to go
simone-silvestri 9706ad6
sprinkle more OneTo around
simone-silvestri 8c450e3
tests should pass
simone-silvestri 8b75606
this works unfortunate we cannot use fields
simone-silvestri 7c83f6b
revert what not needed
simone-silvestri 8b3db85
stricter conditions
simone-silvestri ae32aed
should be ready to go
simone-silvestri 391af30
another bugfix
simone-silvestri 4b50abe
allowscalar
simone-silvestri 44610fc
bump to 0.94.4
simone-silvestri 96f098d
support for Flat grids
simone-silvestri d4f1c10
Update multiary_operations.jl
simone-silvestri dd3817c
support for nothing
simone-silvestri bd38e5f
Merge branch 'ss/top-boundary-condition' of github.com:CliMA/Oceanani…
simone-silvestri 397bf28
bugfix
simone-silvestri File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,76 @@ | ||
const PolarBoundaryCondition{V} = BoundaryCondition{<:Value, <:PolarValue} | ||
|
||
condition_operand(data, grid, Loc, condition, mask) = data | ||
|
||
struct PolarValue{C} | ||
c :: C | ||
end | ||
|
||
PolarBoundaryCondition(field) = | ||
ValueBoundaryCondition(PolarValue(field)) | ||
|
||
@inline getbc(pv::PolarValue, args...) = getbc(pv.c, args...) | ||
|
||
# TODO: vectors should have a different treatment since vector components should account for the frame of reference | ||
# North - South flux boundary conditions are not valid on a Latitude-Longitude grid if the last / first rows represent the poles | ||
function regularize_north_boundary_condition(bc::DefaultBoundaryCondition, grid::LatitudeLongitudeGrid, loc, dim, args...) where C | ||
φmax = φnode(grid.Ny+1, grid, Face()) | ||
|
||
if φmax == 90 | ||
_, LY, LZ = loc | ||
field = Field{Nothing, LY, LZ}(grid; indices = (:, grid.Ny, :)) | ||
return PolarBoundaryCondition(field) | ||
else | ||
return regularize_boundary_condition(bc, grid, args...) | ||
end | ||
end | ||
|
||
# North - South flux boundary conditions are not valid on a Latitude-Longitude grid if the last / first rows represent the poles | ||
function regularize_south_boundary_condition(bc::DefaultBoundaryCondition, grid::LatitudeLongitudeGrid, loc, dim, args...) | ||
φmin = φnode(1, grid, Face()) | ||
|
||
if φmin == - 90 | ||
_, LY, LZ = loc | ||
field = Field{Nothing, LY, LZ}(grid; indices = (:, 1, :)) | ||
return PolarBoundaryCondition(field) | ||
else | ||
return regularize_boundary_condition(bc, grid, args...) | ||
end | ||
end | ||
|
||
function update_pole_value!(bc::PolarBoundaryCondition, c, grid, loc) | ||
operand = condition_operand(c, grid, loc, nothing, 0) | ||
mean!(bc.condition.c, operand) | ||
return nothing | ||
end | ||
|
||
function fill_south_halo!(c, bc::PolarBoundaryCondition, size, offset, loc, arch, grid, args...; only_local_halos = false, kwargs...) = | ||
update_pole_value!(bc, c, grid, loc) | ||
return launch!(arch, grid, KernelParameters(size, offset), | ||
_fill_only_south_halo!, c, bc, loc, grid, Tuple(args); kwargs...) | ||
end | ||
|
||
function fill_north_halo!(c, bc::PolarBoundaryCondition, size, offset, loc, arch, grid, args...; only_local_halos = false, kwargs...) = | ||
update_pole_value!(bc, c, grid, loc) | ||
return launch!(arch, grid, KernelParameters(size, offset), | ||
_fill_only_north_halo!, c, bc, loc, grid, Tuple(args); kwargs...) | ||
end | ||
|
||
function fill_south_and_north_halo!(c, south_bc::PolarBoundaryCondition, north_bc, size, offset, loc, arch, grid, args...; only_local_halos = false, kwargs...) | ||
update_pole_value!(south_bc, c, grid, loc) | ||
return launch!(arch, grid, KernelParameters(size, offset), | ||
_fill_south_and_north_halo!, c, south_bc, north_bc, loc, grid, Tuple(args); kwargs...) | ||
end | ||
|
||
function fill_south_and_north_halo!(c, south_bc, north_bc::PolarBoundaryCondition, size, offset, loc, arch, grid, args...; only_local_halos = false, kwargs...) | ||
update_pole_value!(north_bc, c, grid, loc) | ||
return launch!(arch, grid, KernelParameters(size, offset), | ||
_fill_south_and_north_halo!, c, south_bc, north_bc, loc, grid, Tuple(args); kwargs...) | ||
end | ||
|
||
function fill_south_and_north_halo!(c, south_bc::PolarBoundaryCondition, north_bc::PolarBoundaryCondition, size, offset, loc, arch, grid, args...; only_local_halos = false, kwargs...) | ||
update_pole_value!(south_bc, c, grid, loc) | ||
update_pole_value!(north_bc, c, grid, loc) | ||
return launch!(arch, grid, KernelParameters(size, offset), | ||
_fill_south_and_north_halo!, c, south_bc, north_bc, loc, grid, Tuple(args); kwargs...) | ||
end |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -25,3 +25,4 @@ function update_boundary_condition!(fields::Tuple, model) | |
|
||
return nothing | ||
end | ||
|
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
This does seem right
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
But should
LY
actually be flipped? Eg this is Face for a Center location. (It may not matter)There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Probably you are right. I have found out unfortunately that reductions on windowed fields do not work. I ll try to make it work so here we can do
mean!(windowed_field, full_field)