Skip to content

Commit

Permalink
Added FlowComponentNS
Browse files Browse the repository at this point in the history
  • Loading branch information
joshuahansel committed Apr 13, 2023
1 parent c4fedb7 commit bf9fb0a
Show file tree
Hide file tree
Showing 8 changed files with 388 additions and 0 deletions.
@@ -0,0 +1,10 @@
# FlowComponentNS

This component inherits from [FileMeshComponent.md] to load a mesh and wraps
[NSFVAction.md] to add a Navier-Stokes flow model.

!syntax parameters /Components/FlowComponentNS

!syntax inputs /Components/FlowComponentNS

!syntax children /Components/FlowComponentNS
55 changes: 55 additions & 0 deletions modules/thermal_hydraulics/include/components/FlowComponentNS.h
@@ -0,0 +1,55 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#pragma once

#include "FileMeshComponent.h"
#include "NSFVBase.h"
#include "THMMesh.h"

/**
* Navier-Stokes flow component
*/
class FlowComponentNS : public NSFVBase<FileMeshComponent>
{
public:
static InputParameters validParams();

FlowComponentNS(const InputParameters & parameters);

virtual void addVariables() override;
virtual void addMooseObjects() override;

protected:
virtual void init() override;

virtual std::vector<SubdomainName> getBlocks() const override { return getSubdomainNames(); }
virtual Factory & getFactory() override { return getMooseApp().getFactory(); }
virtual FEProblemBase & getProblem() override { return getMooseApp().feProblem(); }
virtual const MooseMesh & getMesh() const override { return constMesh(); }
virtual void addNSNonlinearVariable(const std::string & var_type,
const std::string & var_name,
InputParameters & params) override
{
getTHMProblem().addSimVariable(true, var_type, var_name, params);
}
virtual void addNSAuxVariable(const std::string & var_type,
const std::string & var_name,
InputParameters & params) override
{
getTHMProblem().addSimVariable(false, var_type, var_name, params);
}
virtual void addNSInitialCondition(const std::string & type,
const std::string & name,
InputParameters & params) override
{
getTHMProblem().addSimInitialCondition(type, name, params);
}
virtual std::string prefix() const override { return name() + ":"; }
};
54 changes: 54 additions & 0 deletions modules/thermal_hydraulics/src/components/FlowComponentNS.C
@@ -0,0 +1,54 @@
//* This file is part of the MOOSE framework
//* https://www.mooseframework.org
//*
//* All rights reserved, see COPYRIGHT for full restrictions
//* https://github.com/idaholab/moose/blob/master/COPYRIGHT
//*
//* Licensed under LGPL 2.1, please see LICENSE for details
//* https://www.gnu.org/licenses/lgpl-2.1.html

#include "FlowComponentNS.h"

registerMooseObject("ThermalHydraulicsApp", FlowComponentNS);

InputParameters
FlowComponentNS::validParams()
{
InputParameters params = NSFVBase<FileMeshComponent>::validParams();

params.addClassDescription("Navier-Stokes flow component.");

return params;
}

FlowComponentNS::FlowComponentNS(const InputParameters & parameters)
: NSFVBase<FileMeshComponent>(parameters)
{
checkCopyNSNodalVariables();
}

void
FlowComponentNS::addVariables()
{
addNSVariables();
addNSInitialConditions();
}

void
FlowComponentNS::addMooseObjects()
{
addNSUserObjects();
addNSKernels();
addNSBoundaryConditions();
addNSMaterials();
addNSPostprocessors();
}

void
FlowComponentNS::init()
{
FileMeshComponent::init();

processMesh();
copyNSNodalVariables();
}
@@ -0,0 +1,225 @@
# This test was adapted from the following:
# modules/navier_stokes/test/tests/finite_volume/pwcns/channel-flow/2d-transient-action.i
# Mesh generation was moved outside of this file, and an additional channel was
# added (without heat transfer to solid region).

# Solid properties
cp_s = 2
rho_s = 4
k_s = 1e-2
h_fs = 10

# Operating conditions
u_inlet = 1
T_inlet = 200
p_outlet = 10
top_side_temperature = 150

[Variables]
[T_solid]
type = MooseVariableFVReal
initial_condition = 100
[]
[]

[AuxVariables]
[porosity]
type = MooseVariableFVReal
initial_condition = 0.5
[]
[velocity_norm]
type = MooseVariableFVReal
[]
[]

[FluidProperties]
[fp]
type = FlibeFluidProperties
[]
[]

[Components]
[comp1]
type = FlowComponentNS
file = rectangle.e
position = '0 0 0'

compressibility = 'weakly-compressible'
add_energy_equation = true
porous_medium_treatment = true

density = 'rho'
dynamic_viscosity = 'mu'
thermal_conductivity = 'k'
specific_heat = 'cp'

initial_velocity = '${u_inlet} 1e-6 0'
initial_pressure = '${p_outlet}'
initial_temperature = '${T_inlet}'

inlet_boundaries = 'comp1:left'
momentum_inlet_types = 'fixed-velocity'
momentum_inlet_function = '${u_inlet} 0'
energy_inlet_types = 'fixed-temperature'
energy_inlet_function = '${T_inlet}'

wall_boundaries = 'comp1:top comp1:bottom'
momentum_wall_types = 'noslip symmetry'
energy_wall_types = 'heatflux heatflux'
energy_wall_function = '0 0'

outlet_boundaries = 'comp1:right'
momentum_outlet_types = 'fixed-pressure'
pressure_function = '${p_outlet}'

ambient_convection_alpha = 'h_cv'
ambient_temperature = 'T_solid'
[]

[comp2]
type = FlowComponentNS
file = rectangle.e
position = '0 2 0'

compressibility = 'weakly-compressible'
add_energy_equation = true
porous_medium_treatment = true

density = 'rho'
dynamic_viscosity = 'mu'
thermal_conductivity = 'k'
specific_heat = 'cp'

initial_velocity = '${u_inlet} 1e-6 0'
initial_pressure = '${p_outlet}'
initial_temperature = '${T_inlet}'

inlet_boundaries = 'comp2:left'
momentum_inlet_types = 'fixed-velocity'
momentum_inlet_function = '${u_inlet} 0'
energy_inlet_types = 'fixed-temperature'
energy_inlet_function = '${T_inlet}'

wall_boundaries = 'comp2:top comp2:bottom'
momentum_wall_types = 'noslip symmetry'
energy_wall_types = 'heatflux heatflux'
energy_wall_function = '0 0'

outlet_boundaries = 'comp2:right'
momentum_outlet_types = 'fixed-pressure'
pressure_function = '${p_outlet}'

ambient_convection_alpha = 'h_cv'
ambient_temperature = 'T_solid'
[]
[]

[FVKernels]
[solid_energy_time]
type = PINSFVEnergyTimeDerivative
variable = T_solid
cp = ${cp_s}
rho = ${rho_s}
is_solid = true
porosity = 'porosity'
[]
[solid_energy_diffusion]
type = FVDiffusion
variable = T_solid
coeff = ${k_s}
[]
[solid_energy_convection]
type = PINSFVEnergyAmbientConvection
variable = T_solid
is_solid = true
T_fluid = 'T_fluid'
T_solid = 'T_solid'
h_solid_fluid = 'h_cv'
[]
[]

[FVBCs]
[heated-side]
type = FVDirichletBC
boundary = 'comp1:top'
variable = 'T_solid'
value = ${top_side_temperature}
[]
[]

[Materials]
[const_functor]
type = ADGenericFunctorMaterial
prop_names = 'h_cv'
prop_values = '${h_fs}'
[]
[fluid_props_to_mat_props]
type = GeneralFunctorFluidProps
fp = fp
pressure = 'pressure'
T_fluid = 'T_fluid'
speed = 'velocity_norm'

# To initialize with a high viscosity
mu_rampdown = 'mu_rampdown'

# For porous flow
characteristic_length = 1
porosity = 'porosity'
[]
[]

[Functions]
[mu_rampdown]
type = PiecewiseLinear
x = '1 2 3 4'
y = '1e3 1e2 1e1 1'
[]
[]

[AuxKernels]
[speed]
type = ParsedAux
variable = 'velocity_norm'
coupled_variables = 'superficial_vel_x superficial_vel_y porosity'
expression = 'sqrt(superficial_vel_x*superficial_vel_x + superficial_vel_y*superficial_vel_y) / porosity'
[]
[]

[Executioner]
type = Transient
solve_type = 'NEWTON'
petsc_options_iname = '-pc_type -ksp_gmres_restart -sub_pc_type -sub_pc_factor_shift_type'
petsc_options_value = 'asm 100 lu NONZERO'
line_search = 'none'
nl_rel_tol = 1e-12

end_time = 3.0
[]

# Some basic Postprocessors to examine the solution
[Postprocessors]
[inlet-p]
type = SideAverageValue
variable = pressure
boundary = 'comp1:left'
[]
[outlet-u]
type = SideAverageValue
variable = superficial_vel_x
boundary = 'comp1:right'
[]
[outlet-temp]
type = SideAverageValue
variable = T_fluid
boundary = 'comp1:right'
[]
[solid-temp]
type = ElementAverageValue
variable = T_solid
[]
[]

[Outputs]
exodus = true
[]
Binary file not shown.
Binary file not shown.
@@ -0,0 +1,21 @@
# This input file is used to generate a rectangle mesh for other tests. It
# should be run with "--mesh-only rectangle.e".

[Mesh]
[gen]
type = GeneratedMeshGenerator
dim = 2
xmin = 0
xmax = 10
ymin = 0
ymax = 1
nx = 20
ny = 5
[]
[rename_block]
type = RenameBlockGenerator
input = gen
old_block = 0
new_block = 'body'
[]
[]
@@ -0,0 +1,23 @@
[Tests]
design = 'FlowComponentNS.md'
issues = '#23794'

[generate_mesh]
type = Exodiff
input = rectangle.i
exodiff = rectangle.e
cli_args = "--mesh-only rectangle.e"
requirement = 'The system shall generate a mesh for the FlowComponentNS tests.'
[]
[pwcns]
type = Exodiff
input = flow_component_ns_pwc.i
exodiff = flow_component_ns_pwc_out.e
abs_zero = 1e-9
prereq = 'generate_mesh'
method = "!dbg"
valgrind = HEAVY
recover = false # see #19126
requirement = 'The system shall model the porous, weakly compressible Navier-Stokes equations with a finite volume discretization, using a component.'
[]
[]

0 comments on commit bf9fb0a

Please sign in to comment.