Skip to content

Commit

Permalink
1) used SSTKernelHex8Mesh class instead of a new calss for the transi…
Browse files Browse the repository at this point in the history
…tion model, 2) used gamma_intermittency instead of gamint in the unit test, 3) deleted some comments
  • Loading branch information
Bumseok Lee committed Sep 19, 2024
1 parent 450b8dc commit 2bcb163
Show file tree
Hide file tree
Showing 4 changed files with 30 additions and 93 deletions.
24 changes: 12 additions & 12 deletions unit_tests/kernels/UnitTestKernelUtils.C
Original file line number Diff line number Diff line change
Expand Up @@ -158,16 +158,16 @@ struct TrigFieldFunction
std::sin(a * pi * z));
}

void gamint(const double* coords, double* qField) const //added: transition
void gamma_intermittency(const double* coords, double* qField) const
{
double x = coords[0];
double y = coords[1];
double z = coords[2];

// Range should be from 0.02 to 1.0
qField[0] =
gamintnot + abs(std::cos(a * pi * x) * std::sin(a * pi * y) *
std::sin(a * pi * z))/(1.0 - gamintnot);
gamma_intermittencynot + abs(std::cos(a * pi * x) * std::sin(a * pi * y) *
std::sin(a * pi * z))/(1.0 - gamma_intermittencynot);
}

void dwalldistdx(const double* coords, double* qField) const
Expand Down Expand Up @@ -340,8 +340,8 @@ private:
/// Factor for sdr field
static constexpr double sdrnot{1.0};

/// Factor for gamint field
static constexpr double gamintnot{0.02};
/// Factor for gamma_intermittency field
static constexpr double gamma_intermittencynot{0.02};

/// Factor for dwalldistdx field
static constexpr double dwalldistdxnot{1.0};
Expand Down Expand Up @@ -402,11 +402,11 @@ init_trigonometric_field(
funcPtr = &TrigFieldFunction::dkdx;
else if (fieldName == "specific_dissipation_rate")
funcPtr = &TrigFieldFunction::sdr;
else if (fieldName == "gamma_transition") // added: transition
funcPtr = &TrigFieldFunction::gamint;
else if (fieldName == "dwalldistdx") // added: transition
else if (fieldName == "gamma_transition")
funcPtr = &TrigFieldFunction::gamma_intermittency;
else if (fieldName == "dwalldistdx")
funcPtr = &TrigFieldFunction::dwalldistdx;
else if (fieldName == "dnDotVdx") // added: transition
else if (fieldName == "dnDotVdx")
funcPtr = &TrigFieldFunction::dnDotVdx;
else if (fieldName == "total_dissipation_rate")
funcPtr = &TrigFieldFunction::tdr;
Expand Down Expand Up @@ -564,12 +564,12 @@ sdr_test_function(
}

void
gamint_test_function( // added: transition
gamma_intermittency_test_function(
const stk::mesh::BulkData& bulk,
const sierra::nalu::VectorFieldType& coordinates,
sierra::nalu::ScalarFieldType& gamint)
sierra::nalu::ScalarFieldType& gamma_intermittency)
{
init_trigonometric_field(bulk, coordinates, gamint);
init_trigonometric_field(bulk, coordinates, gamma_intermittency);
}

void
Expand Down
96 changes: 17 additions & 79 deletions unit_tests/kernels/UnitTestKernelUtils.h
Original file line number Diff line number Diff line change
Expand Up @@ -83,10 +83,10 @@ void sdr_test_function(
const sierra::nalu::VectorFieldType& coordinates,
sierra::nalu::ScalarFieldType& sdr);

void gamint_test_function( //added
void gamma_intermittency_test_function(
const stk::mesh::BulkData& bulk,
const sierra::nalu::VectorFieldType& coordinates,
sierra::nalu::ScalarFieldType& gamint);
sierra::nalu::ScalarFieldType& gamma_intermittency);

void dwalldistdx_test_function(
const stk::mesh::BulkData& bulk,
Expand Down Expand Up @@ -746,8 +746,10 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
&meta_->declare_field<double>(stk::topology::NODE_RANK, "viscosity")),
tvisc_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "turbulent_viscosity")),
gamint_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "gamma_transition")), //added
gamma_intermittency_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "gamma_transition")),
dwalldistdx_(&meta_->declare_field<double>(stk::topology::NODE_RANK, "dwalldistdx")),
dnDotVdx_(&meta_->declare_field<double>(stk::topology::NODE_RANK, "dnDotVdx")),
maxLengthScale_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "sst_max_length_scale")),
minDistance_(&meta_->declare_field<double>(
Expand Down Expand Up @@ -776,7 +778,7 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
stk::mesh::put_field_on_mesh(*sdrbc_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*visc_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*tvisc_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*gamint_, meta_->universal_part(), nullptr); //added
stk::mesh::put_field_on_mesh(*gamma_intermittency_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(
*maxLengthScale_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(
Expand All @@ -794,6 +796,12 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
stk::mesh::put_field_on_mesh(
*dwdx_, meta_->universal_part(), spatialDim_, nullptr);
stk::io::set_field_output_type(*dwdx_, stk::io::FieldOutputType::VECTOR_3D);
stk::mesh::put_field_on_mesh(
*dwalldistdx_, meta_->universal_part(), spatialDim_, nullptr);
stk::io::set_field_output_type(*dwalldistdx_, stk::io::FieldOutputType::VECTOR_3D);
stk::mesh::put_field_on_mesh(
*dnDotVdx_, meta_->universal_part(), spatialDim_, nullptr);
stk::io::set_field_output_type(*dnDotVdx_, stk::io::FieldOutputType::VECTOR_3D);
double initOpenMassFlowRate[sierra::nalu::AlgTraitsQuad4::numScsIp_];
for (int i = 0; i < sierra::nalu::AlgTraitsQuad4::numScsIp_; ++i) {
initOpenMassFlowRate[i] = 10.0;
Expand Down Expand Up @@ -825,7 +833,7 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
*bulk_, *coordinates_, *density_);
unit_test_kernel_utils::tke_test_function(*bulk_, *coordinates_, *tke_);
unit_test_kernel_utils::sdr_test_function(*bulk_, *coordinates_, *sdr_);
unit_test_kernel_utils::gamint_test_function(*bulk_, *coordinates_, *gamint_); // added
unit_test_kernel_utils::gamma_intermittency_test_function(*bulk_, *coordinates_, *gamma_intermittency_);
unit_test_kernel_utils::minimum_distance_to_wall_test_function(
*bulk_, *coordinates_, *minDistance_);
unit_test_kernel_utils::sst_f_one_blending_test_function(
Expand All @@ -836,6 +844,8 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
stk::mesh::field_fill(0.0, *dkdx_);
stk::mesh::field_fill(0.0, *dwdx_);
stk::mesh::field_fill(0.0, *pecletFactor_);
stk::mesh::field_fill(0.0, *dwalldistdx_);
stk::mesh::field_fill(0.0, *dnDotVdx_);
}

sierra::nalu::ScalarFieldType* tke_{nullptr};
Expand All @@ -844,7 +854,7 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
sierra::nalu::ScalarFieldType* sdrbc_{nullptr};
sierra::nalu::ScalarFieldType* visc_{nullptr};
sierra::nalu::ScalarFieldType* tvisc_{nullptr};
sierra::nalu::ScalarFieldType* gamint_{nullptr}; // added
sierra::nalu::ScalarFieldType* gamma_intermittency_{nullptr};
sierra::nalu::ScalarFieldType* maxLengthScale_{nullptr};
sierra::nalu::ScalarFieldType* minDistance_{nullptr};
sierra::nalu::ScalarFieldType* fOneBlend_{nullptr};
Expand All @@ -858,78 +868,6 @@ class SSTKernelHex8Mesh : public LowMachKernelHex8Mesh
sierra::nalu::ScalarFieldType* sdrWallArea_{nullptr};
sierra::nalu::GenericFieldType* wallFricVel_{nullptr};
sierra::nalu::ScalarFieldType* pecletFactor_{nullptr};
};

/** Test Fixture for the BLT Gamma Kernels
*
*/
class BLTGammaM2015KernelHex8Mesh : public LowMachKernelHex8Mesh
{
public:
BLTGammaM2015KernelHex8Mesh()
: LowMachKernelHex8Mesh(),
tke_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "turbulent_ke")),
sdr_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "specific_dissipation_rate")),
visc_(
&meta_->declare_field<double>(stk::topology::NODE_RANK, "viscosity")),
tvisc_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "turbulent_viscosity")),
gamint_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "gamma_transition")), //added
minDistance_(&meta_->declare_field<double>(
stk::topology::NODE_RANK, "minimum_distance_to_wall")),
dudx_(&meta_->declare_field<double>(stk::topology::NODE_RANK, "dudx")),
dwalldistdx_(&meta_->declare_field<double>(stk::topology::NODE_RANK, "dwalldistdx")),
dnDotVdx_(&meta_->declare_field<double>(stk::topology::NODE_RANK, "dnDotVdx"))
{
stk::mesh::put_field_on_mesh(*tke_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*sdr_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*visc_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*tvisc_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(*gamint_, meta_->universal_part(), nullptr); //added
stk::mesh::put_field_on_mesh(
*minDistance_, meta_->universal_part(), nullptr);
stk::mesh::put_field_on_mesh(
*dudx_, meta_->universal_part(), spatialDim_ * spatialDim_, nullptr);
stk::io::set_field_output_type(
*dudx_, stk::io::FieldOutputType::FULL_TENSOR_36);
stk::mesh::put_field_on_mesh(
*dwalldistdx_, meta_->universal_part(), spatialDim_, nullptr);
stk::io::set_field_output_type(*dwalldistdx_, stk::io::FieldOutputType::VECTOR_3D);
stk::mesh::put_field_on_mesh(
*dnDotVdx_, meta_->universal_part(), spatialDim_, nullptr);
stk::io::set_field_output_type(*dnDotVdx_, stk::io::FieldOutputType::VECTOR_3D);
}
virtual ~BLTGammaM2015KernelHex8Mesh() {}

virtual void fill_mesh_and_init_fields(
const bool doPerturb = false, const bool generateSidesets = false) override
{
LowMachKernelHex8Mesh::fill_mesh_and_init_fields(
doPerturb, generateSidesets);
stk::mesh::field_fill(0.2, *visc_);
stk::mesh::field_fill(0.3, *tvisc_);
unit_test_kernel_utils::density_test_function(
*bulk_, *coordinates_, *density_);
unit_test_kernel_utils::tke_test_function(*bulk_, *coordinates_, *tke_);
unit_test_kernel_utils::sdr_test_function(*bulk_, *coordinates_, *sdr_);
unit_test_kernel_utils::gamint_test_function(*bulk_, *coordinates_, *gamint_); // added
unit_test_kernel_utils::minimum_distance_to_wall_test_function(
*bulk_, *coordinates_, *minDistance_);
unit_test_kernel_utils::dudx_test_function(*bulk_, *coordinates_, *dudx_);
stk::mesh::field_fill(0.0, *dwalldistdx_);
stk::mesh::field_fill(0.0, *dnDotVdx_);
}

sierra::nalu::ScalarFieldType* tke_{nullptr};
sierra::nalu::ScalarFieldType* sdr_{nullptr};
sierra::nalu::ScalarFieldType* visc_{nullptr};
sierra::nalu::ScalarFieldType* tvisc_{nullptr};
sierra::nalu::ScalarFieldType* gamint_{nullptr};
sierra::nalu::ScalarFieldType* minDistance_{nullptr};
sierra::nalu::TensorFieldType* dudx_{nullptr};
sierra::nalu::VectorFieldType* dwalldistdx_{nullptr};
sierra::nalu::VectorFieldType* dnDotVdx_{nullptr};
};
Expand Down
2 changes: 1 addition & 1 deletion unit_tests/node_kernels/UnitTestBLTGammaM2015NodeKernel.C
Original file line number Diff line number Diff line change
Expand Up @@ -114,7 +114,7 @@ static constexpr double lhs[8][8] = {
} // namespace hex8_golds
} // namespace

TEST_F(BLTGammaM2015KernelHex8Mesh, NGP_blt_gamma_node)
TEST_F(SSTKernelHex8Mesh, NGP_blt_gamma_node)
{
// Only execute for 1 processor runs
if (bulk_->parallel_size() > 1)
Expand Down
1 change: 0 additions & 1 deletion unit_tests/node_kernels/UnitTestSSTNodeKernel.C
Original file line number Diff line number Diff line change
Expand Up @@ -1368,7 +1368,6 @@ TEST_F(SSTKernelHex8Mesh, NGP_sdr_sst_trans_node)
helperObjs.nodeAlg->add_kernel<sierra::nalu::SDRSSTBLTM2015NodeKernel>(*meta_);

helperObjs.execute();
//helperObjs.print_lhs_and_rhs();

Kokkos::deep_copy(
helperObjs.linsys->hostNumSumIntoCalls_,
Expand Down

0 comments on commit 2bcb163

Please sign in to comment.