Skip to content
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

Unintuitive behavior of sum #943

Open
kmdeck opened this issue Sep 14, 2022 · 4 comments
Open

Unintuitive behavior of sum #943

kmdeck opened this issue Sep 14, 2022 · 4 comments
Labels
bug Something isn't working

Comments

@kmdeck
Copy link
Member

kmdeck commented Sep 14, 2022

Describe the bug

I am trying to use sum to compute integrals of different fields over different spaces. In particular, I am taking sum(ones(space)) - which should be volume or area or column or "point" integrals. But,

  • I get funny results when using level of the 3d space vs the horizontal space for the spherical shell. For example, I would expect that \int (ones(horizontal_space)) = 4pi R^2, and that \int ones(level) = 4pi(R+z)^2, where z is the height of the level. What \int ones(level) returns - for any level, at any height - is 4\pi R^2 dz/2, where dz is the discretization of the vertical.
  • I would expect sum(point_field) = 0, but I get a nonzero result.
  • I also found that using by column to take a column integral over each column of the spherical shell - returning something that lives on the horizontal space - also seemed to give unintuitive results.

To Reproduce

MWE below

using ClimaCore
using Test
@testset "Spherical Shell" begin
    FT = Float32
    radius = FT(10.0)
    height =FT(2.0)
    nelements = (10,10)
    npolynomial = 3
    vertdomain = ClimaCore.Domains.IntervalDomain(
        ClimaCore.Geometry.ZPoint(FT(0)),
        ClimaCore.Geometry.ZPoint(FT(height));
        boundary_tags = (:bottom, :top),
    )
    
    vertmesh = ClimaCore.Meshes.IntervalMesh(
        vertdomain,
        ClimaCore.Meshes.Uniform(),
        nelems = nelements[2],
    )
    vert_center_space = ClimaCore.Spaces.CenterFiniteDifferenceSpace(vertmesh)

    horzdomain = ClimaCore.Domains.SphereDomain(radius)
    horzmesh = ClimaCore.Meshes.EquiangularCubedSphere(horzdomain, nelements[1])
    horztopology = ClimaCore.Topologies.Topology2D(horzmesh)
    quad = ClimaCore.Spaces.Quadratures.GLL{npolynomial + 1}()
    horzspace = ClimaCore.Spaces.SpectralElementSpace2D(horztopology, quad)

    hv_center_space = ClimaCore.Spaces.ExtrudedFiniteDifferenceSpace(
        horzspace,
        vert_center_space,
    )
    hv_face_space = ClimaCore.Spaces.FaceExtrudedFiniteDifferenceSpace(hv_center_space)
    # Sum behavior
    ones_volume = ones(hv_center_space)
    ones_horizontal_space = ones(hv_center_space.horizontal_space)
    ones_bottom_space_level = ones(ClimaCore.Spaces.level(hv_face_space, ClimaCore.Utilities.PlusHalf(0)))
    ones_top_space_level = ones(ClimaCore.Spaces.level(hv_face_space, ClimaCore.Utilities.PlusHalf(10)))
    # Passing tests
    @test sum(ones_volume)  FT(4.0*π)*radius^2.0*height 
    @test sum(ones_horizontal_space)   FT(4.0*π)*radius^2.0

    # Not passing
    @test sum(ones_bottom_space_level)  FT(4.0*π)*radius^2.0 
    @test sum(ones_top_space_level)  FT(4.0*π)*(radius+height)^2.0
    column_integrals = zeros(horzspace)
    ClimaCore.Fields.bycolumn(horzspace) do colidx
        column_integrals[colidx] .= sum(ones_volume[colidx])
    end
    # Expect this should be ∫dz for each element? = height?
    @test sum(parent(column_integrals .- height)) FT(0.0)
end

@testset "Column" begin
    zlim = Float32.((0.0,10.0))
    nelements = 20
    boundary_tags = (:bottom, :top)
    column = ClimaCore.Domains.IntervalDomain(
        ClimaCore.Geometry.ZPoint{FT}(zlim[1]),
        ClimaCore.Geometry.ZPoint{FT}(zlim[2]);
        boundary_tags = boundary_tags,
    )
    mesh = ClimaCore.Meshes.IntervalMesh(column; nelems = nelements)
    center_space = ClimaCore.Spaces.CenterFiniteDifferenceSpace(mesh)
    face_space = ClimaCore.Spaces.FaceFiniteDifferenceSpace(mesh)
    bottom_space = ClimaCore.Spaces.level(face_space, ClimaCore.Utilities.PlusHalf(0))
    
    ones_column = ones(center_space)
    ones_point_bottom = ones(bottom_space)
    #Passes
    @test sum(ones_column)  FT(10)
    # Does not pass
    @test sum(ones_point_bottom)  FT(0) #∫dx over a point is zero
end

System details

Any relevant system information:

  • Julia version 1.7
  • operating system Mac
@kmdeck kmdeck added the bug Something isn't working label Sep 14, 2022
bors bot added a commit to CliMA/ClimaCoupler.jl that referenced this issue Sep 16, 2022
117: transition to using new bucket model version r=kmdeck a=kmdeck

# PULL REQUEST

## Purpose and Content
This PR updates the coupler driver file and other helper files so that the coupled runs work with the new version of the bucket model (with snow, and with a multi layer heat equation).

To do: followup on CliMA/ClimaCore.jl#943
## Benefits and Risks
The benefit of this is that our AMIP model will be slightly more realistic.

## Linked Issues
(Provide references to any link issues. Use closes #issuenum to automatically close an open issue)
- Fixes #118 
- Closes #118 

## PR Checklist
- [X] This PR has a corresponding issue OR is linked to an SDI.
- [X] I have followed CliMA's codebase [contribution](https://clima.github.io/ClimateMachine.jl/latest/Contributing/) and [style](https://clima.github.io/ClimateMachine.jl/latest/DevDocs/CodeStyle/) guidelines OR N/A.
- [X] I have followed CliMA's [documentation policy](https://github.com/CliMA/policies/wiki/Documentation-Policy).
- [X] I have checked all issues and PRs and I certify that this PR does not duplicate an open PR.
- [x] I linted my code on my local machine prior to submission OR N/A.
- [X] Unit tests are included OR N/A.
- [X] Code used in an integration test OR N/A.
- [x] All tests ran successfully on my local machine OR N/A.
- [X] All classes, modules, and function contain docstrings OR N/A.
- [X] Documentation has been added/updated OR N/A.


Co-authored-by: kmdeck <[email protected]>
@simonbyrne
Copy link
Member

This one has bothered me as well: I admit I don't have a good solution.

@kmdeck
Copy link
Member Author

kmdeck commented Sep 20, 2022

Is sum actually doing an integral over the domain? Maybe it would be clearer if there was a formula showing what sum was computing, and then a user could see how it reduces in some case to the integral and some cases not? Im not sure. I dont know what it is doing under the hood.

@valeriabarra
Copy link
Member

@kmdeck yes, my understanding is that sum has to use the Jacobian information of the underlying domain of the field you are trying to compute the integral of.
In essence, think about a change in coordinate transformation. Say, for instance, that you switch from Cartesian to polar coordinates, then any integral you compute has to account for the change of coordinates (and that is, in fact, the Jacobian that comes out of the chain rule). This is why we have to account for WJ and we tried to make this clear in the docstring (the W is due to the quadrature point weights needed to compute integrals, since we are using spectral spaces. These are just equal to 1 for Finite Difference spaces).

So to address your second bullet point, yes, we should expect sum of a point-field to return just zero.

For the first bullet point, is it due to the fact that we use the "shallow-atmosphere" approximation for our 3D sphere? See for instance this comment in this test.

@kmdeck
Copy link
Member Author

kmdeck commented Sep 22, 2022

Hi @valeriabarra - thank you! It does make sense that we need the Jacobian in the integral. That said, the results from sum still dont make sense in the spherical case (to me!).

I dont think that the issue in the first test/bullet point is because of the shallow-atmosphere approx. The sum I was showing in the test is over a 2D spectral element space, so the height is only involved in that it affects the radius of the sphere we are taking an area integral of. In the shallow case, we would expect 4pi(R+h)^2 \approx 4piR^2(1+2h/R), also, and not 4piR^2 dz/2...

as a comment - making the shallow approximation makes a lot of sense for the globe, but there is no restriction on what the user can pass in for h and r, so that assumption shouldnt be used in the formula for sum anywhere...

(Also note the funny behavior with bycolumn - seperate issue)

bors bot added a commit that referenced this issue Sep 28, 2022
962: Add support for definite column integrals r=simonbyrne a=charleskawczynski

This PR adds `column_integral_definite!`, a non-allocating definite integral for columns. The idea is that it will replace the existing (allocating) definite integral [here](https://github.com/CliMA/ClimaAtmos.jl/blob/44cee55def2a51433c5dcdcdae010b7777cd741b/examples/hybrid/sphere/baroclinic_wave_utilities.jl#L167-L186). However, this implementation has several advantages:
 - It's compatible with `bycolumn`, so the computation can occur in parallel
 - It's allocation-free
 - It's _tested_ on single column and sphere configurations

Looking at the flame graph for allocations (using the Julia 1.8 Allocs module with `sample_rate=1`), this is responsible for most of the remaining allocations in ClimaAtmos:
<img width="1993" alt="Screen Shot 2022-09-26 at 7 57 05 AM" src="https://user-images.githubusercontent.com/1880641/192353757-03101c41-2c0b-4ccb-a8b9-43faa78680f8.png">


The interface for the added function captures two cases:

```julia
function column_integral_definite!(col∫field::Field, field::Field)
    bycolumn(axes(field)) do colidx
        column_integral_definite!(col∫field[colidx], field[colidx])
        nothing
    end
    return nothing
end
```
and
```julia
function column_integral_definite!(
    col∫field::PointField,
    field::ColumnField,
)
    `@inbounds` col∫field[] = column_integral_definite(field)
    return nothing
end
```

A step towards closing #943, #748, [CA#686](CliMA/ClimaAtmos.jl#686).

## A note on an alternative approach
Ideally, we would write this function as `column_integral_definite!(fn, args...)` where we might be able to write a broadcast statement like:

```julia
`@.` f2D = column_integral_definite(f3D) do z
     f3D.x * z^2
end
```

however, this would require us to define broadcasting between planar and 3D domains, which is not trivial (or maybe not possible) because of ambiguities. The ambiguities arise because (2D, 3D) broadcasting may want different things in different cases, for example:

 - (f2D, f3D) -> f3D: mul full 3D field by planar surface value
 - (f2D, f3D) -> f2D: perform reduction over z-coordinate to yield 2D field

The situation is similar when thinking about what happens when we make views. For example,

```julia
Fields.bycolumn(axes(f3D)) do colidx
     `@.` f2D[colidx] = column_integral_definite(f3D[colidx]) do z
         f3D[colidx].x * z^2
      end
end
```

Now, we have to define how `DataF` data layouts broadcast with `VF`. Again, we have two cases:

 - (f0D, f1D) -> f1D: mul full 1D field by 0D field
 - (f0D, f1D) -> f0D: perform reduction over z-coordinate to yield 0D field

My vote/preference is to support the first cases (which is partially supported already) and write custom functions (e.g., reductions) that operate on single fields for the second case.

Co-authored-by: Charles Kawczynski <[email protected]>
Co-authored-by: Simon Byrne <[email protected]>
bors bot added a commit that referenced this issue Oct 1, 2022
962: Add support for definite column integrals r=charleskawczynski a=charleskawczynski

This PR adds `column_integral_definite!`, a non-allocating definite integral for columns. The idea is that it will replace the existing (allocating) definite integral [here](https://github.com/CliMA/ClimaAtmos.jl/blob/44cee55def2a51433c5dcdcdae010b7777cd741b/examples/hybrid/sphere/baroclinic_wave_utilities.jl#L167-L186). However, this implementation has several advantages:
 - It's compatible with `bycolumn`, so the computation can occur in parallel
 - It's allocation-free
 - It's _tested_ on single column and sphere configurations

Looking at the flame graph for allocations (using the Julia 1.8 Allocs module with `sample_rate=1`), this is responsible for most of the remaining allocations in ClimaAtmos:
<img width="1993" alt="Screen Shot 2022-09-26 at 7 57 05 AM" src="https://user-images.githubusercontent.com/1880641/192353757-03101c41-2c0b-4ccb-a8b9-43faa78680f8.png">


The interface for the added function captures two cases:

```julia
function column_integral_definite!(col∫field::Field, field::Field)
    bycolumn(axes(field)) do colidx
        column_integral_definite!(col∫field[colidx], field[colidx])
        nothing
    end
    return nothing
end
```
and
```julia
function column_integral_definite!(
    col∫field::PointField,
    field::ColumnField,
)
    `@inbounds` col∫field[] = column_integral_definite(field)
    return nothing
end
```

A step towards closing #943, #748, [CA#686](CliMA/ClimaAtmos.jl#686).

## A note on an alternative approach
Ideally, we would write this function as `column_integral_definite!(fn, args...)` where we might be able to write a broadcast statement like:

```julia
`@.` f2D = column_integral_definite(f3D) do z
     f3D.x * z^2
end
```

however, this would require us to define broadcasting between planar and 3D domains, which is not trivial (or maybe not possible) because of ambiguities. The ambiguities arise because (2D, 3D) broadcasting may want different things in different cases, for example:

 - (f2D, f3D) -> f3D: mul full 3D field by planar surface value
 - (f2D, f3D) -> f2D: perform reduction over z-coordinate to yield 2D field

The situation is similar when thinking about what happens when we make views. For example,

```julia
Fields.bycolumn(axes(f3D)) do colidx
     `@.` f2D[colidx] = column_integral_definite(f3D[colidx]) do z
         f3D[colidx].x * z^2
      end
end
```

Now, we have to define how `DataF` data layouts broadcast with `VF`. Again, we have two cases:

 - (f0D, f1D) -> f1D: mul full 1D field by 0D field
 - (f0D, f1D) -> f0D: perform reduction over z-coordinate to yield 0D field

My vote/preference is to support the first cases (which is partially supported already) and write custom functions (e.g., reductions) that operate on single fields for the second case.

A step towards closing [CA#686](CliMA/ClimaAtmos.jl#686)

Co-authored-by: Charles Kawczynski <[email protected]>
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
bug Something isn't working
Projects
None yet
Development

No branches or pull requests

3 participants