GRINS-0.8.0
List of all members | Public Member Functions | Private Member Functions
GRINS::HeatTransferAdjointStabilization< Conductivity > Class Template Reference

Adds VMS-based stabilization to LowMachNavierStokes physics class. More...

#include <heat_transfer_adjoint_stab.h>

Inheritance diagram for GRINS::HeatTransferAdjointStabilization< Conductivity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::HeatTransferAdjointStabilization< Conductivity >:
Collaboration graph
[legend]

Public Member Functions

 HeatTransferAdjointStabilization (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~HeatTransferAdjointStabilization ()
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
- Public Member Functions inherited from GRINS::HeatTransferStabilizationBase< Conductivity >
 HeatTransferStabilizationBase (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~HeatTransferStabilizationBase ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
- Public Member Functions inherited from GRINS::HeatTransferBase< Conductivity >
 HeatTransferBase (const std::string &physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~HeatTransferBase ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. More...
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Private Member Functions

 HeatTransferAdjointStabilization ()
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::Physics
static void set_is_axisymmetric (bool is_axisymmetric)
 Set whether we should treat the problem as axisymmetric. More...
 
static bool is_axisymmetric ()
 
- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::Physics
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void parse_enabled_subdomains (const GetPot &input, const std::string &physics_name)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. More...
 
- Protected Attributes inherited from GRINS::HeatTransferStabilizationBase< Conductivity >
HeatTransferStabilizationHelper _stab_helper
 
- Protected Attributes inherited from GRINS::HeatTransferBase< Conductivity >
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
PrimitiveTempFEVariables_temp_vars
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
libMesh::Number _Cp
 
Conductivity _k
 Conductivity. More...
 
- Protected Attributes inherited from GRINS::Physics
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
- Static Protected Attributes inherited from GRINS::Physics
static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 
static bool _is_axisymmetric = false
 Caches whether we are solving an axisymmetric problem or not. More...
 

Detailed Description

template<class Conductivity>
class GRINS::HeatTransferAdjointStabilization< Conductivity >

Adds VMS-based stabilization to LowMachNavierStokes physics class.

Definition at line 36 of file heat_transfer_adjoint_stab.h.

Constructor & Destructor Documentation

template<class K >
GRINS::HeatTransferAdjointStabilization< K >::HeatTransferAdjointStabilization ( const GRINS::PhysicsName physics_name,
const GetPot &  input 
)

Definition at line 38 of file heat_transfer_adjoint_stab.C.

40  : HeatTransferStabilizationBase<K>(physics_name,input)
41  {
42  return;
43  }

Definition at line 46 of file heat_transfer_adjoint_stab.C.

47  {
48  return;
49  }
template<class Conductivity >
GRINS::HeatTransferAdjointStabilization< Conductivity >::HeatTransferAdjointStabilization ( )
private

Member Function Documentation

template<class K >
void GRINS::HeatTransferAdjointStabilization< K >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

Time dependent part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 53 of file heat_transfer_adjoint_stab.C.

55  {
56  // The number of local degrees of freedom in each variable.
57  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
58  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
59 
60  // Element Jacobian * quadrature weights for interior integration.
61  const std::vector<libMesh::Real> &JxW =
62  context.get_element_fe(this->_temp_vars.T())->get_JxW();
63 
64  const std::vector<std::vector<libMesh::Real> >& u_phi =
65  context.get_element_fe(this->_flow_vars.u())->get_phi();
66 
67  const std::vector<std::vector<libMesh::Real> >& T_phi =
68  context.get_element_fe(this->_temp_vars.T())->get_phi();
69 
70  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
71  context.get_element_fe(this->_temp_vars.T())->get_dphi();
72 
73  const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
74  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
75 
76  /*
77  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
78 
79  const std::vector<std::vector<libMesh::Real> >& u_phi =
80  context.get_element_fe(this->_flow_vars.u())->get_phi();
81 
82  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{p}
83  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{p}
84  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
85  if(this->_flow_vars.dim() == 3)
86  {
87  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
88  }
89  */
90 
91  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
92  libMesh::DenseSubMatrix<libMesh::Number> &KTT =
93  context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T()); // J_{TT}
94  libMesh::DenseSubMatrix<libMesh::Number> &KTu =
95  context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.u()); // J_{Tu}
96  libMesh::DenseSubMatrix<libMesh::Number> &KTv =
97  context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.v()); // J_{Tv}
98  libMesh::DenseSubMatrix<libMesh::Number> *KTw = NULL;
99 
100  if(this->_flow_vars.dim() == 3)
101  {
102  KTw = &context.get_elem_jacobian
103  (this->_temp_vars.T(), this->_flow_vars.w()); // J_{Tw}
104  }
105 
106  unsigned int n_qpoints = context.get_element_qrule().n_points();
107 
108  for (unsigned int qp=0; qp != n_qpoints; qp++)
109  {
110  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
111 
112  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
113  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
114 
115  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
116  context.interior_value( this->_flow_vars.v(), qp ) );
117  if( this->_flow_vars.dim() == 3 )
118  {
119  U(2) = context.interior_value( this->_flow_vars.w(), qp );
120  }
121 
122  //libMesh::RealGradient grad_T = context.interior_gradient( this->_temp_vars.T(), qp );
123 
124  libMesh::Real tau_E, RE_s;
125  libMesh::Real d_tau_E_dT_rho, d_RE_s_dT;
126  libMesh::Gradient d_tau_E_dU, d_RE_s_dgradT, d_RE_s_dU;
127  libMesh::Tensor d_RE_s_dhessT;
128 
129  // Compute the conductivity at this qp
130  libMesh::Real _k_qp = this->_k(context, qp);
131 
132  if (compute_jacobian)
133  {
135  ( context, G, this->_rho, this->_Cp, _k_qp, U,
136  tau_E, d_tau_E_dT_rho, d_tau_E_dU, this->_is_steady );
138  ( context, qp, this->_rho, this->_Cp, _k_qp,
139  RE_s, d_RE_s_dT, d_RE_s_dgradT, d_RE_s_dhessT,
140  d_RE_s_dU );
141  }
142  else
143  {
144  tau_E = this->_stab_helper.compute_tau_energy
145  ( context, G, this->_rho, this->_Cp, _k_qp, U, this->_is_steady );
146 
148  ( context, qp, this->_rho, this->_Cp, _k_qp );
149  }
150 
151  /*
152  for (unsigned int i=0; i != n_u_dofs; i++)
153  {
154  Fu(i) += -tau_E*RE_s*this->_rho*this->_Cp*u_phi[i][qp]*grad_T(0)*JxW[qp];
155  Fv(i) += -tau_E*RE_s*this->_rho*this->_Cp*u_phi[i][qp]*grad_T(1)*JxW[qp];
156  if( this->_flow_vars.dim() == 3 )
157  {
158  (*Fw)(i) += -tau_E*RE_s*this->_rho*this->_Cp*u_phi[i][qp]*grad_T(2)*JxW[qp];
159  }
160  }
161  */
162 
163  for (unsigned int i=0; i != n_T_dofs; i++)
164  {
165  FT(i) += -tau_E*RE_s*( this->_rho*this->_Cp*U*T_gradphi[i][qp]
166  + _k_qp*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2))
167  )*JxW[qp];
168  if (compute_jacobian)
169  {
170  for (unsigned int j=0; j != n_T_dofs; ++j)
171  {
172  KTT(i,j) += -tau_E*
173  (d_RE_s_dT*T_phi[j][qp] +
174  d_RE_s_dgradT*T_gradphi[j][qp] +
175  d_RE_s_dhessT.contract(T_hessphi[j][qp])
176  ) *
177  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
178  + _k_qp*(T_hessphi[i][qp](0,0) +
179  T_hessphi[i][qp](1,1) +
180  T_hessphi[i][qp](2,2))
181  )*JxW[qp]
182  * context.get_fixed_solution_derivative();
183  }
184  for (unsigned int j=0; j != n_u_dofs; ++j)
185  {
186  KTu(i,j) += -tau_E*RE_s*
187  ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](0) )*JxW[qp]
188  * context.get_fixed_solution_derivative();
189  KTu(i,j) +=
190  -(tau_E*d_RE_s_dU(0)+d_tau_E_dU(0)*RE_s)*u_phi[j][qp]*
191  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
192  + _k_qp*(T_hessphi[i][qp](0,0) +
193  T_hessphi[i][qp](1,1) +
194  T_hessphi[i][qp](2,2))
195  )*JxW[qp]
196  * context.get_fixed_solution_derivative();
197  KTv(i,j) += -tau_E*RE_s*
198  ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](1) )*JxW[qp]
199  * context.get_fixed_solution_derivative();
200  KTv(i,j) +=
201  -(tau_E*d_RE_s_dU(1)+d_tau_E_dU(1)*RE_s)*u_phi[j][qp]*
202  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
203  + _k_qp*(T_hessphi[i][qp](0,0) +
204  T_hessphi[i][qp](1,1) +
205  T_hessphi[i][qp](2,2))
206  )*JxW[qp]
207  * context.get_fixed_solution_derivative();
208  }
209  if(this->_flow_vars.dim() == 3)
210  {
211  for (unsigned int j=0; j != n_u_dofs; ++j)
212  {
213  (*KTw)(i,j) += -tau_E*RE_s*
214  ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](2) )*JxW[qp]
215  * context.get_fixed_solution_derivative();
216  (*KTw)(i,j) +=
217  -(tau_E*d_RE_s_dU(2)+d_tau_E_dU(2)*RE_s)*u_phi[j][qp]*
218  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
219  + _k_qp*(T_hessphi[i][qp](0,0) +
220  T_hessphi[i][qp](1,1) +
221  T_hessphi[i][qp](2,2))
222  )*JxW[qp]
223  * context.get_fixed_solution_derivative();
224  }
225  }
226  }
227  }
228  }
229  }
VariableIndex T() const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
void compute_res_energy_steady_and_derivs(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp, const libMesh::Real k, libMesh::Real &res, libMesh::Real &d_res_dT, libMesh::Gradient &d_res_dgradT, libMesh::Tensor &d_res_dhessT, libMesh::Gradient &d_res_dU) const
VelocityVariable & _flow_vars
void compute_tau_energy_and_derivs(AssemblyContext &c, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Real cp, libMesh::Real k, libMesh::Gradient U, libMesh::Real &tau_E, libMesh::Real &d_tau_E_d_rho, libMesh::Gradient &d_tau_E_d_U, bool is_steady) const
libMesh::Number _rho
Material parameters, read from input.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
Conductivity _k
Conductivity.
HeatTransferStabilizationHelper _stab_helper
libMesh::Real compute_tau_energy(AssemblyContext &c, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Real cp, libMesh::Real k, libMesh::Gradient U, bool is_steady) const
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
libMesh::Real compute_res_energy_steady(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp, const libMesh::Real k) const
unsigned int dim() const
Number of components.
PrimitiveTempFEVariables & _temp_vars
template<class K >
void GRINS::HeatTransferAdjointStabilization< K >::mass_residual ( bool  ,
AssemblyContext  
)
virtual

Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part.

Reimplemented from GRINS::Physics.

Definition at line 233 of file heat_transfer_adjoint_stab.C.

234  {
235  // The number of local degrees of freedom in each variable.
236  const unsigned int n_T_dofs = context.get_dof_indices(this->_temp_vars.T()).size();
237  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
238 
239  // Element Jacobian * quadrature weights for interior integration.
240  const std::vector<libMesh::Real> &JxW =
241  context.get_element_fe(this->_temp_vars.T())->get_JxW();
242 
243  const std::vector<std::vector<libMesh::Real> >& u_phi =
244  context.get_element_fe(this->_flow_vars.u())->get_phi();
245 
246  const std::vector<std::vector<libMesh::Real> >& T_phi =
247  context.get_element_fe(this->_temp_vars.T())->get_phi();
248 
249  const std::vector<std::vector<libMesh::RealGradient> >& T_gradphi =
250  context.get_element_fe(this->_temp_vars.T())->get_dphi();
251 
252  const std::vector<std::vector<libMesh::RealTensor> >& T_hessphi =
253  context.get_element_fe(this->_temp_vars.T())->get_d2phi();
254 
255  /*
256  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
257 
258  const std::vector<std::vector<libMesh::Real> >& u_phi =
259  context.get_element_fe(this->_flow_vars.u())->get_phi();
260 
261  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{p}
262  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{p}
263  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
264  if(this->_flow_vars.dim() == 3)
265  {
266  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
267  }
268  */
269 
270  libMesh::DenseSubVector<libMesh::Number> &FT = context.get_elem_residual(this->_temp_vars.T()); // R_{T}
271  libMesh::DenseSubMatrix<libMesh::Number> &KTT =
272  context.get_elem_jacobian(this->_temp_vars.T(), this->_temp_vars.T()); // J_{TT}
273  libMesh::DenseSubMatrix<libMesh::Number> &KTu =
274  context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.u()); // J_{Tu}
275  libMesh::DenseSubMatrix<libMesh::Number> &KTv =
276  context.get_elem_jacobian(this->_temp_vars.T(), this->_flow_vars.v()); // J_{Tv}
277  libMesh::DenseSubMatrix<libMesh::Number> *KTw = NULL;
278 
279  if(this->_flow_vars.dim() == 3)
280  {
281  KTw = &context.get_elem_jacobian
282  (this->_temp_vars.T(), this->_flow_vars.w()); // J_{Tw}
283  }
284 
285  unsigned int n_qpoints = context.get_element_qrule().n_points();
286 
287  for (unsigned int qp=0; qp != n_qpoints; qp++)
288  {
289  libMesh::FEBase* fe = context.get_element_fe(this->_temp_vars.T());
290 
291  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
292  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
293 
294  libMesh::RealGradient U( context.fixed_interior_value( this->_flow_vars.u(), qp ),
295  context.fixed_interior_value( this->_flow_vars.v(), qp ) );
296  if( this->_flow_vars.dim() == 3 )
297  U(2) = context.fixed_interior_value( this->_flow_vars.w(), qp );
298 
299  //libMesh::RealGradient grad_T = context.fixed_interior_gradient( this->_temp_vars.T(), qp );
300 
301  libMesh::Real tau_E, RE_t;
302  libMesh::Real d_tau_E_d_rho, d_RE_t_dT;
303  libMesh::Gradient d_tau_E_dU;
304 
305  // Compute the conductivity at this qp
306  libMesh::Real _k_qp = this->_k(context, qp);
307 
308  if (compute_jacobian)
309  {
311  ( context, G, this->_rho, this->_Cp, _k_qp, U,
312  tau_E, d_tau_E_d_rho, d_tau_E_dU, false );
314  ( context, qp, this->_rho, this->_Cp,
315  RE_t, d_RE_t_dT );
316  }
317  else
318  {
319  tau_E = this->_stab_helper.compute_tau_energy
320  ( context, G, this->_rho, this->_Cp, _k_qp, U, false );
321 
323  ( context, qp, this->_rho, this->_Cp );
324  }
325 
326 
327  /*
328  for (unsigned int i=0; i != n_u_dofs; i++)
329  {
330  Fu(i) += -tau_E*RE_t*this->_rho*this->_Cp*u_phi[i][qp]*grad_T(0)*JxW[qp];
331  Fv(i) += -tau_E*RE_t*this->_rho*this->_Cp*u_phi[i][qp]*grad_T(1)*JxW[qp];
332  if( this->_flow_vars.dim() == 3 )
333  {
334  (*Fw)(i) += -tau_E*RE_t*this->_rho*this->_Cp*u_phi[i][qp]*grad_T(2)*JxW[qp];
335  }
336  }
337  */
338 
339  for (unsigned int i=0; i != n_T_dofs; i++)
340  {
341  FT(i) -= tau_E*RE_t*( this->_rho*this->_Cp*U*T_gradphi[i][qp]
342  + _k_qp*(T_hessphi[i][qp](0,0) + T_hessphi[i][qp](1,1) + T_hessphi[i][qp](2,2))
343  )*JxW[qp];
344  if (compute_jacobian)
345  {
346  const libMesh::Real fixed_deriv =
347  context.get_fixed_solution_derivative();
348 
349  for (unsigned int j=0; j != n_T_dofs; ++j)
350  {
351  KTT(i,j) -=
352  (tau_E*d_RE_t_dT)*T_phi[j][qp]*
353  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
354  + _k_qp*(T_hessphi[i][qp](0,0) +
355  T_hessphi[i][qp](1,1) +
356  T_hessphi[i][qp](2,2))
357  )*JxW[qp];
358  }
359  for (unsigned int j=0; j != n_u_dofs; ++j)
360  {
361  KTu(i,j) -=
362  d_tau_E_dU(0)*u_phi[j][qp]*RE_t*
363  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
364  + _k_qp*(T_hessphi[i][qp](0,0) +
365  T_hessphi[i][qp](1,1) +
366  T_hessphi[i][qp](2,2))
367  )*fixed_deriv*JxW[qp];
368  KTu(i,j) -=
369  tau_E*RE_t*
370  ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](0)*fixed_deriv
371  )*JxW[qp];
372  KTv(i,j) -=
373  d_tau_E_dU(1)*u_phi[j][qp]*RE_t*
374  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
375  + _k_qp*(T_hessphi[i][qp](0,0) +
376  T_hessphi[i][qp](1,1) +
377  T_hessphi[i][qp](2,2))
378  )*fixed_deriv*JxW[qp];
379  KTv(i,j) -=
380  tau_E*RE_t*
381  ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](1)*fixed_deriv
382  )*JxW[qp];
383  }
384  if(this->_flow_vars.dim() == 3)
385  {
386  for (unsigned int j=0; j != n_u_dofs; ++j)
387  {
388  (*KTw)(i,j) -=
389  d_tau_E_dU(2)*u_phi[j][qp]*RE_t*
390  ( this->_rho*this->_Cp*U*T_gradphi[i][qp]
391  + _k_qp*(T_hessphi[i][qp](0,0) +
392  T_hessphi[i][qp](1,1) +
393  T_hessphi[i][qp](2,2))
394  )*fixed_deriv*JxW[qp];
395  (*KTw)(i,j) -=
396  tau_E*RE_t*
397  ( this->_rho*this->_Cp*u_phi[j][qp]*T_gradphi[i][qp](2)*fixed_deriv
398  )*JxW[qp];
399  }
400  }
401  }
402  }
403 
404  }
405  }
libMesh::Real compute_res_energy_transient(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp) const
VariableIndex T() const
void compute_res_energy_transient_and_derivs(AssemblyContext &context, unsigned int qp, const libMesh::Real rho, const libMesh::Real Cp, libMesh::Real &res, libMesh::Real &d_res_dTdot) const
libMesh::RealTensor compute_G(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:69
VelocityVariable & _flow_vars
void compute_tau_energy_and_derivs(AssemblyContext &c, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Real cp, libMesh::Real k, libMesh::Gradient U, libMesh::Real &tau_E, libMesh::Real &d_tau_E_d_rho, libMesh::Gradient &d_tau_E_d_U, bool is_steady) const
libMesh::Number _rho
Material parameters, read from input.
libMesh::RealGradient compute_g(libMesh::FEBase *fe, AssemblyContext &c, unsigned int qp) const
Definition: stab_helper.C:47
Conductivity _k
Conductivity.
HeatTransferStabilizationHelper _stab_helper
libMesh::Real compute_tau_energy(AssemblyContext &c, libMesh::RealTensor &G, libMesh::Real rho, libMesh::Real cp, libMesh::Real k, libMesh::Gradient U, bool is_steady) const
unsigned int dim() const
Number of components.
PrimitiveTempFEVariables & _temp_vars

The documentation for this class was generated from the following files:

Generated on Tue Dec 19 2017 12:47:30 for GRINS-0.8.0 by  doxygen 1.8.9.1