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

panzer: Fixed Dirichlet resid for edge basis and UVM #9670

Merged
merged 1 commit into from
Sep 7, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Original file line number Diff line number Diff line change
Expand Up @@ -88,18 +88,22 @@ template <typename EvalT,typename Traits>
void CurlSolution<EvalT,Traits>::evaluateFields(typename Traits::EvalData workset)
{
using panzer::index_t;
for (index_t cell = 0; cell < workset.num_cells; ++cell) {
for (int point = 0; point < solution.extent_int(1); ++point) {
auto ip_coordinates = this->wda(workset).int_rules[ir_index]->ip_coordinates.get_static_view();
auto solution_v = solution.get_static_view();
auto solution_curl_v = solution_curl.get_static_view();

const double & x = this->wda(workset).int_rules[ir_index]->ip_coordinates(cell,point,0);
const double & y = this->wda(workset).int_rules[ir_index]->ip_coordinates(cell,point,1);
Kokkos::parallel_for (workset.num_cells, KOKKOS_LAMBDA (const index_t cell) {
for (int point = 0; point < solution_v.extent_int(1); ++point) {

solution(cell,point,0) = -(y-1.0)*y + cos(2.0*M_PI*x)*sin(2.0*M_PI*y);
solution(cell,point,1) = -(x-1.0)*x + sin(2.0*M_PI*x)*cos(2.0*M_PI*y);
const double & x = ip_coordinates(cell,point,0);
const double & y = ip_coordinates(cell,point,1);

solution_curl(cell,point) = -2.0*x+2.0*y;
solution_v(cell,point,0) = -(y-1.0)*y + cos(2.0*M_PI*x)*sin(2.0*M_PI*y);
solution_v(cell,point,1) = -(x-1.0)*x + sin(2.0*M_PI*x)*cos(2.0*M_PI*y);

solution_curl_v(cell,point) = -2.0*x+2.0*y;
}
}
});
}

//**********************************************************************
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -88,18 +88,22 @@ template <typename EvalT,typename Traits>
void CurlSolution<EvalT,Traits>::evaluateFields(typename Traits::EvalData workset)
{
using panzer::index_t;
for (index_t cell = 0; cell < workset.num_cells; ++cell) {
for (int point = 0; point < solution.extent_int(1); ++point) {
auto ip_coordinates = this->wda(workset).int_rules[ir_index]->ip_coordinates.get_static_view();
auto solution_v = solution.get_static_view();
auto solution_curl_v = solution_curl.get_static_view();

const double & x = this->wda(workset).int_rules[ir_index]->ip_coordinates(cell,point,0);
const double & y = this->wda(workset).int_rules[ir_index]->ip_coordinates(cell,point,1);
Kokkos::parallel_for (workset.num_cells, KOKKOS_LAMBDA (const index_t cell) {
for (int point = 0; point < solution_v.extent_int(1); ++point) {

solution(cell,point,0) = -(y-1.0)*y + cos(2.0*M_PI*x)*sin(2.0*M_PI*y);
solution(cell,point,1) = -(x-1.0)*x + sin(2.0*M_PI*x)*cos(2.0*M_PI*y);
const double & x = ip_coordinates(cell,point,0);
const double & y = ip_coordinates(cell,point,1);

solution_curl(cell,point) = -2.0*x+2.0*y;
solution_v(cell,point,0) = -(y-1.0)*y + cos(2.0*M_PI*x)*sin(2.0*M_PI*y);
solution_v(cell,point,1) = -(x-1.0)*x + sin(2.0*M_PI*x)*cos(2.0*M_PI*y);

solution_curl_v(cell,point) = -2.0*x+2.0*y;
}
}
});
}

//**********************************************************************
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -40,6 +40,8 @@
// ***********************************************************************
// @HEADER

#include "Kokkos_View_Fad.hpp"

#include "PanzerDiscFE_config.hpp"

#include "Panzer_ExplicitTemplateInstantiation.hpp"
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -98,8 +98,8 @@ class DirichletResidual_EdgeBasis
PointValues2<double> pointValues;

Teuchos::RCP<const std::vector<Intrepid2::Orientation> > orientations;
Intrepid2::RefSubcellParametrization<PHX::Device>::ConstViewType edgeParam; //edge parametrization
Intrepid2::RefSubcellParametrization<PHX::Device>::ConstViewType faceParam; //face parametrization
Intrepid2::RefSubcellParametrization<Kokkos::HostSpace>::ConstViewType edgeParam; //edge parametrization
Intrepid2::RefSubcellParametrization<Kokkos::HostSpace>::ConstViewType faceParam; //face parametrization

}; // end of class DirichletResidual_EdgeBasis

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -122,10 +122,10 @@ postRegistrationSetup(
const int edgeDim = 1;
const int faceDim = 2;
if(cellTopo.getDimension() > edgeDim)
edgeParam = Intrepid2::RefSubcellParametrization<PHX::Device>::get(edgeDim, cellTopo.getKey());
edgeParam = Intrepid2::RefSubcellParametrization<Kokkos::HostSpace>::get(edgeDim, cellTopo.getKey());

if(cellTopo.getDimension() > faceDim)
faceParam = Intrepid2::RefSubcellParametrization<PHX::Device>::get(faceDim, cellTopo.getKey());
faceParam = Intrepid2::RefSubcellParametrization<Kokkos::HostSpace>::get(faceDim, cellTopo.getKey());
}

//**********************************************************************
Expand Down Expand Up @@ -157,9 +157,14 @@ evaluateFields(
const WorksetDetails & details = workset;

//const bool is_normalize = true;
auto work = Kokkos::createDynRankView(residual.get_static_view(),"work", 4, cellDim);
auto work = Kokkos::create_mirror_view(Kokkos::createDynRankView(residual.get_static_view(),"work", 4, cellDim));

// compute residual
auto residual_h = Kokkos::create_mirror_view(residual.get_static_view());
auto dof_h = Kokkos::create_mirror_view(dof.get_static_view());
auto value_h = Kokkos::create_mirror_view(value.get_static_view());
Kokkos::deep_copy(dof_h, dof.get_static_view());
Kokkos::deep_copy(value_h, value.get_static_view());
switch (subcellDim) {
case 1: { // 2D element Tri and Quad
if (intrepid_basis->getDofCount(1, subcellOrd)) {
Expand All @@ -169,6 +174,7 @@ evaluateFields(
const int ndofsEdge = intrepid_basis->getDofCount(1, subcellOrd);
const int numEdges = cellTopo.getEdgeCount();
/* */ int edgeOrts[4] = {};

for(index_t c=0;c<workset.num_cells;c++) {
orientations->at(details.cell_local_ids[c]).getEdgeOrientation(edgeOrts, numEdges);

Expand All @@ -180,11 +186,11 @@ evaluateFields(

for (int i=0;i<ndofsEdge;++i) {
const int b = intrepid_basis->getDofOrdinal(1, subcellOrd, i);
auto J = Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL());
auto J = Kokkos::create_mirror_view(Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL()));
Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);

for(int d=0;d<cellDim;d++) {
residual(c,b) += (dof(c,b,d)-value(c,b,d))*phyEdgeTan(d);
residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyEdgeTan(d);
}
}
}
Expand Down Expand Up @@ -217,11 +223,11 @@ evaluateFields(
edgeOrts[edgeOrd]);

{
auto J = Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL());
auto J = Kokkos::create_mirror_view(Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL()));
Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);

for(int d=0;d<dof.extent_int(2);d++) {
residual(c,b) += (dof(c,b,d)-value(c,b,d))*phyEdgeTan(d);
residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyEdgeTan(d);
}
}
}
Expand All @@ -245,13 +251,13 @@ evaluateFields(
faceOrts[subcellOrd]);

for(int b=0;b<dof.extent_int(1);b++) {
auto J = Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL());
auto J = Kokkos::create_mirror_view(Kokkos::subview(worksetJacobians, c, b, Kokkos::ALL(), Kokkos::ALL()));
Intrepid2::Kernels::Serial::matvec_product(phyFaceTanU, J, ortFaceTanU);
Intrepid2::Kernels::Serial::matvec_product(phyFaceTanV, J, ortFaceTanV);

for(int d=0;d<dof.extent_int(2);d++) {
residual(c,b) += (dof(c,b,d)-value(c,b,d))*phyFaceTanU(d);
residual(c,b) += (dof(c,b,d)-value(c,b,d))*phyFaceTanV(d);
residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyFaceTanU(d);
residual_h(c,b) += (dof_h(c,b,d)-value_h(c,b,d))*phyFaceTanV(d);
}
}
}
Expand All @@ -260,6 +266,7 @@ evaluateFields(
break;
}
}
Kokkos::deep_copy(residual.get_static_view(), residual_h);
}

}
Expand Down