-
Notifications
You must be signed in to change notification settings - Fork 19
Second Coupling Integral (times d tilde/dn) #1589
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
Open
int-ptr-ptr
wants to merge
20
commits into
PrincetonUniversity:devel
Choose a base branch
from
int-ptr-ptr:integral-dot-grad-shape
base: devel
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 all commits
Commits
Show all changes
20 commits
Select commit
Hold shift + click to select a range
0692e7d
coupling integral wip -- pre-devel-merge
int-ptr-ptr cd215a0
kernel up to validation
int-ptr-ptr 8b2e8e1
completed test -- passes
int-ptr-ptr 5d3153f
Merge remote-tracking branch 'public_repo/issue-1554' into integral-d…
int-ptr-ptr 7b84f12
1554 compatibility
int-ptr-ptr 7f6ea3d
compilation issues on dshape
int-ptr-ptr 52c5a71
purged snuck internal kokkos import
int-ptr-ptr 611e7cb
fixed unroll pragma location
int-ptr-ptr 1ccd299
fixed gpu compile and bug with accessor axis ordering
int-ptr-ptr 635c8ac
restored removed quadrature name/desc features
int-ptr-ptr c1c08ba
Merge branch 'devel' into integral-dot-grad-shape
int-ptr-ptr 8398c5c
Merge branch 'devel' into integral-dot-grad-shape
int-ptr-ptr 36242dc
Merge remote-tracking branch 'origin/devel' into integral-dot-grad-shape
int-ptr-ptr 31a8058
conflicts post-merge
int-ptr-ptr 27ec141
shape function deriv precompute kernel start
int-ptr-ptr 20f4259
completed precompute
int-ptr-ptr beda088
removed iterator patch, moistened edgeview code
int-ptr-ptr 72a7a38
missing pragma once
int-ptr-ptr 9945905
Merge remote-tracking branch 'origin/devel' into integral-dot-grad-shape
int-ptr-ptr 1be3df5
lagrange_derivative.hpp: specfem_setup -> specfem/setup
int-ptr-ptr 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
Some comments aren't visible on the classic Files Changed page.
There are no files selected for viewing
117 changes: 117 additions & 0 deletions
117
core/specfem/algorithms/coupling_integral1d_dnshape.hpp
This file contains hidden or 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,117 @@ | ||
| #pragma once | ||
|
|
||
| #include "specfem/assembly.hpp" | ||
| #include "specfem/element_coupling.hpp" | ||
| #include "specfem/execution.hpp" | ||
| #include "specfem/point.hpp" | ||
| #include <Kokkos_Core.hpp> | ||
|
|
||
| /** | ||
| * @brief Temporary file for integral in issue # 1550. this function should be | ||
| * merged with `coupling_integral1d.hpp:coupling_integral` once the intersection | ||
| * field typing is updated to store two fields. | ||
| * | ||
| */ | ||
|
|
||
| namespace specfem::algorithms { | ||
|
|
||
| /** | ||
| * @brief Takes a field `intersection_field` on the intersection and computes, | ||
| * for each self GLL point, the integral of `intersection_field` times the | ||
| * normal derivative of the shape function at that point. `intersection_field` | ||
| * should be call-accessible (e.g. Kokkos::View) with shape: | ||
| * | ||
| * (chunk_size, n_quad_intersection, self::components()) | ||
| * | ||
| * After handling any other intersection forces, boundary conditions, etc. the | ||
| * result can be `atomic_add`ed to the acceleration field. | ||
| * | ||
| * @tparam IndexType The chunk_edge iterator type | ||
| * @tparam IntersectionFieldViewType The type of `intersection_field` | ||
| * @tparam IntersectionFactor A nonconforming chunk_edge accessor | ||
| * holding `intersection_factor` | ||
| * @tparam ShapeFunctionNormalDerivativesType a view holding the normal | ||
| * derivatives of the shape function (implementation may change) | ||
| * @tparam CallableType The callback function, which will be given the point | ||
| * index and corresponding evaluated integral | ||
| * | ||
| * @param ngllz - number of quadrature points in z | ||
| * @param ngllx - number of quadrature points in x | ||
| * @param chunk_index - the outer index (chunk_edge) that gets iterated for | ||
| * points | ||
| * @param intersection_field - the field to integrate | ||
| * @param intersection_factor - nonconforming chunk_edge accessor holding | ||
| * `intersection_factor` | ||
| * @param shape_function_normal_derivatives - (TEMPORARY, TO BE LOADED LATER BY | ||
| * THIS FUNCTION) normal derivatives of the shape functions at each intersection | ||
| * point. | ||
| * @param callback - callback function to capture integral values | ||
| * @ingroup AlgorithmsIntegration | ||
| */ | ||
| template <typename IndexType, typename IntersectionFieldViewType, | ||
| typename IntersectionFactor, | ||
| typename ShapeFunctionNormalDerivativesType, typename CallableType> | ||
| KOKKOS_FUNCTION void coupling_integral_dnshape( | ||
| const int &ngllz, const int &ngllx, const IndexType &chunk_index, | ||
| const IntersectionFieldViewType &intersection_field, | ||
| const IntersectionFactor &intersection_factor, | ||
| const ShapeFunctionNormalDerivativesType &shape_function_normal_derivatives, | ||
| const CallableType &callback) { | ||
| constexpr specfem::element::dimension_tag dimension_tag = | ||
| specfem::element::dimension_tag::dim2; | ||
|
|
||
| constexpr auto self_medium_tag = specfem::element_coupling::attributes< | ||
| dimension_tag, IntersectionFactor::interface_tag>::self_medium(); | ||
|
|
||
| using PointIndexType = | ||
| typename IndexType::iterator_type::index_type::index_type; | ||
| using PointFieldType = | ||
| specfem::point::acceleration<dimension_tag, self_medium_tag, | ||
| IntersectionFieldViewType::using_simd>; | ||
|
|
||
| constexpr int ncomp = PointFieldType::components; | ||
| constexpr int nquad_intersection = IntersectionFactor::n_quad_intersection; | ||
|
|
||
| const auto &team = chunk_index.get_policy_index(); | ||
| const int &num_edges = chunk_index.nedges(); | ||
| specfem::execution::for_each_level( | ||
| specfem::execution::TeamThreadMDRangeIterator(team, num_edges, ngllz, | ||
| ngllx), | ||
| [&](const auto &index) { | ||
| const int iedge = index(0); | ||
| const int iz = index(1); | ||
| const int ix = index(2); | ||
|
|
||
| PointFieldType result; | ||
| #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_ENABLE_HIP) | ||
| #pragma unroll | ||
| #endif | ||
| for (int icomp = 0; icomp < ncomp; icomp++) { | ||
| result(icomp) = 0; | ||
| } | ||
|
|
||
| #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_ENABLE_HIP) | ||
| #pragma unroll | ||
| #endif | ||
| for (int iquad = 0; iquad < nquad_intersection; iquad++) { | ||
|
|
||
| #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_ENABLE_HIP) | ||
| #pragma unroll | ||
| #endif | ||
| for (int icomp = 0; icomp < ncomp; icomp++) { | ||
| result(icomp) += | ||
| intersection_field(iedge, iquad, icomp) * | ||
| intersection_factor(iedge, iquad) * | ||
| shape_function_normal_derivatives(iedge, iz, ix, iquad); | ||
| } | ||
| } | ||
| // this index will always have ipoint = 0, but we will not use it | ||
| const auto edge_index = chunk_index.get_iterator()(iedge).get_index(); | ||
|
|
||
| specfem::point::index<dimension_tag> self_index(edge_index.ispec, iz, | ||
| ix); | ||
| callback(self_index, result); | ||
| }); | ||
| } | ||
|
|
||
| } // namespace specfem::algorithms | ||
177 changes: 177 additions & 0 deletions
177
core/specfem/algorithms/shape_function_normal_derivative.hpp
This file contains hidden or 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,177 @@ | ||
| #pragma once | ||
|
|
||
| #include "specfem/assembly.hpp" | ||
| #include "specfem/assembly/edge_types.hpp" | ||
| #include "specfem/assembly/mesh.hpp" | ||
| #include "specfem/chunk_edge/nonconforming_interface.hpp" | ||
| #include "specfem/element/tags.hpp" | ||
| #include "specfem/element_coupling.hpp" | ||
| #include "specfem/execution.hpp" | ||
| #include "specfem/point.hpp" | ||
| #include <Kokkos_Core.hpp> | ||
|
|
||
| namespace specfem::algorithms { | ||
|
|
||
| /** | ||
| * @brief Takes a field `intersection_field` on the intersection and computes, | ||
| * for each self GLL point, the integral of `intersection_field` times the | ||
| * normal derivative of the shape function at that point. `intersection_field` | ||
| * should be call-accessible (e.g. Kokkos::View) with shape: | ||
| * | ||
| * (chunk_size, n_quad_intersection, self::components()) | ||
| * | ||
| * After handling any other intersection forces, boundary conditions, etc. the | ||
| * result can be `atomic_add`ed to the acceleration field. | ||
| * | ||
| * @tparam dimension_tag dimension of the simulation | ||
| * @param nonconforming_interfaces - assembly.nonconforming_interfaces struct | ||
| */ | ||
| template <specfem::element_coupling::interface_tag interface_tag, | ||
| specfem::element::boundary_tag boundary_tag, int nquad_element_, | ||
| int nquad_intersection_, typename SelfEdgeListView, | ||
| typename HPrimeView, typename TransformedIntersectionNormalView> | ||
| Kokkos::View<type_real ****, Kokkos::DefaultExecutionSpace> | ||
| shape_function_self_normal_derivatives( | ||
| const SelfEdgeListView &self_edges, | ||
| const specfem::assembly::nonconforming_interfaces< | ||
| specfem::element::dimension_tag::dim2> &nonconforming_interfaces, | ||
| const int &ngllz, const int &ngllx, const HPrimeView &hprime, | ||
| const TransformedIntersectionNormalView | ||
| &intersection_normal_contravariant_edgelocal) { | ||
| const auto dimension_tag = specfem::element::dimension_tag::dim2; | ||
| using ReturnViewType = | ||
| Kokkos::View<type_real ****, Kokkos::DefaultExecutionSpace>; | ||
|
|
||
| // TODO: get nquad_intersection from somewhere else | ||
| const int nquad_intersection = nquad_intersection_; | ||
| if (nquad_intersection != nquad_intersection_) { | ||
| throw std::runtime_error( | ||
| std::string("shape_function_self_normal_derivatives() kernel run with " | ||
| "nquad_intersection = ") + | ||
| std::to_string(nquad_intersection_) + | ||
| ", but assembly got nquad_intersection == " + | ||
| std::to_string(nquad_intersection)); | ||
| } | ||
| ReturnViewType normal_derivs("shape_function_self_normal_derivatives", | ||
| self_edges.n_edges, ngllz, ngllx, | ||
| nquad_intersection); | ||
|
|
||
| using parallel_config = | ||
| specfem::parallel_configuration::default_chunk_edge_config< | ||
| dimension_tag, Kokkos::DefaultExecutionSpace>; | ||
| specfem::execution::ChunkedEdgeIterator chunk(parallel_config(), self_edges); | ||
|
|
||
| using TransferFunctionSelf = specfem::chunk_edge::transfer_function_self< | ||
| dimension_tag, interface_tag, boundary_tag, parallel_config::chunk_size, | ||
| nquad_intersection_, nquad_element_>; | ||
|
|
||
| specfem::execution::for_each_level( | ||
| "specfem::compute::shape_function_normal_derivative", | ||
| chunk.set_scratch_size( | ||
| 0, Kokkos::PerTeam(TransferFunctionSelf::shmem_size())), | ||
| KOKKOS_LAMBDA( | ||
| const typename decltype(chunk)::index_type &chunk_iterator_index) { | ||
| const auto &chunk_index = chunk_iterator_index.get_index(); | ||
| const auto &team = chunk_index.get_policy_index(); | ||
| const int &num_edges = chunk_index.nedges(); | ||
|
|
||
| TransferFunctionSelf transfer_function_self(team); | ||
| specfem::assembly::load_on_device(chunk_index, nonconforming_interfaces, | ||
| transfer_function_self); | ||
|
|
||
| specfem::execution::for_each_level( | ||
| specfem::execution::TeamThreadMDRangeIterator(team, num_edges, | ||
| ngllz, ngllx), | ||
| [&](const auto &index) { | ||
| const int iedge = index(0); | ||
| const int iz = index(1); | ||
| const int ix = index(2); | ||
|
|
||
| // this index will always have ipoint = 0, but we will not use it | ||
| const auto edge_index = | ||
| chunk_index.get_iterator()(iedge).get_index(); | ||
|
|
||
| /* | ||
| * We convert local coordinates / indices ix, iz | ||
| * into an edge-local coordinate system: | ||
| * | ||
| * ================================== ipoint_s = ix for top | ||
| * | and bottom, iz for left | ||
| * | and right | ||
| * | | ||
| * ipoint_n = ngll-1 | ||
| * | ||
| * ipoint_n = 0 represents the edge. | ||
| */ | ||
|
|
||
| const bool ipoint_n_is_ix = | ||
| edge_index.edge_type == | ||
| specfem::mesh_entity::dim2::type::left || | ||
| edge_index.edge_type == | ||
| specfem::mesh_entity::dim2::type::right; | ||
| const int ngll_n = ipoint_n_is_ix ? ngllx : ngllz; | ||
| const int ngll_s = ipoint_n_is_ix ? ngllz : ngllx; | ||
|
|
||
| const int ipoint_s = ipoint_n_is_ix ? iz : ix; | ||
| int ipoint_n = ipoint_n_is_ix ? ix : iz; | ||
| if (edge_index.edge_type == | ||
| specfem::mesh_entity::dim2::type::right || | ||
| edge_index.edge_type == | ||
| specfem::mesh_entity::dim2::type::top) { | ||
| // flip ipoint_n coord so edge is 0 | ||
| ipoint_n = ngll_n - 1 - ipoint_n; | ||
| } | ||
|
|
||
| const type_real dlagrange_dn = -hprime(0, ipoint_n); | ||
|
|
||
| for (int iquad = 0; iquad < nquad_intersection; iquad++) { | ||
|
|
||
| // dshape_dn, local-normal derivative (derivative of | ||
| // normal-direction L, which used to be kronecker delta | ||
| // indicating on edge. Now has dlagrange_dn != 0 for all | ||
| // ipoint_n) | ||
|
|
||
| // first, in the local-normal component | ||
| type_real dshape_dn = | ||
| intersection_normal_contravariant_edgelocal( | ||
| edge_index.iedge, iquad, 0) * | ||
| (dlagrange_dn * | ||
| transfer_function_self(iedge, ipoint_s, iquad)); | ||
| // dshape_dn, local-tangential derivative (differentiate | ||
| // transfer_function_self instead of normal-direction L) | ||
| if (ipoint_n == 0) { | ||
| // the local-tangential is constant zero if we are not on the | ||
| // edge. | ||
|
|
||
| // recover L' (at point ipoint_s) on intersection quadrature | ||
| // points: use transfer function. | ||
| type_real transfer_function_derivative = 0; | ||
| for (int ipoint_s2 = 0; ipoint_s2 < ngll_s; ipoint_s2++) { | ||
| transfer_function_derivative += | ||
| hprime(ipoint_s2, ipoint_s) * | ||
| transfer_function_self(iedge, ipoint_s2, iquad); | ||
| } | ||
|
|
||
| dshape_dn += intersection_normal_contravariant_edgelocal( | ||
| edge_index.iedge, iquad, 1) * | ||
| transfer_function_derivative; | ||
| } | ||
|
|
||
| normal_derivs(edge_index.iedge, iz, ix, iquad) = dshape_dn; | ||
| } | ||
| }); | ||
| // specfem::execution::for_each_level( | ||
| // chunk_index.get_iterator(), | ||
| // [&](const typename ChunkIndexType::iterator_type::index_type | ||
| // &iterator_index) { | ||
| // const auto index = iterator_index.get_index(); | ||
| // const auto index_local = iterator_index.get_local_index(); | ||
|
|
||
| // index.ielem | ||
| // }); | ||
|
Comment on lines
+163
to
+171
Collaborator
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. ? |
||
| }); | ||
|
|
||
| return normal_derivs; | ||
| } | ||
|
|
||
| } // namespace specfem::algorithms | ||
This file contains hidden or 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 hidden or 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
Oops, something went wrong.
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.
Precomputing edge-related values and storing them in the assembly can clean up this code, improve readability, and facilitate integration into
coupling_integral.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.
Is it alright if we do this in a future PR? I was waiting for this one to go through before updating my work on #1342 in order to put the first
coupling_integralin the same test executable.Precomputation and storage would require changes to
nonconforming_interfacesthat will bloat this PR even further.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.
That sounds good. Can you create another issue for that refactor and link it here?
Also, see my other comment.