GRINS-0.8.0
convection_diffusion.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2017 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
28 
29 // GRINS
30 #include "grins/common.h"
31 #include "grins/assembly_context.h"
35 #include "grins/single_variable.h"
36 
37 // libMesh
38 #include "libmesh/getpot.h"
39 #include "libmesh/string_to_enum.h"
40 #include "libmesh/quadrature.h"
41 #include "libmesh/fem_system.h"
42 
43 namespace GRINS
44 {
46  const GetPot& input )
47  : Physics(physics_name,input),
48  _v(3,libMesh::ParsedFunction<libMesh::Number>("0.0") ),
49  _kappa("0.0"),
50  _var(GRINSPrivate::VariableWarehouse::get_variable_subclass<SingleVariable>(VariablesParsing::single_variable_name(input,physics_name,VariablesParsing::PHYSICS)))
51  {
52  unsigned int n_v_comps = input.vector_variable_size("Physics/"+physics_name+"/velocity_field");
53 
54  for( unsigned int v = 0; v < n_v_comps; v++ )
55  _v[v].reparse( input("Physics/"+physics_name+"/velocity_field", "0.0", v) );
56 
57  std::string material_name = MaterialsParsing::material_name(input,physics_name);
58 
59  this->set_parameter(this->_kappa, input,
60  "Materials/"+material_name+"/Diffusivity/value",
61  "DIE!");
62 
63  _ic_handler = new GenericICHandler(physics_name,input);
64 
66  }
67 
68  void ConvectionDiffusion::set_time_evolving_vars( libMesh::FEMSystem* system )
69  {
70  system->time_evolving(_var.var(), 1);
71  }
72 
74  {
75  context.get_element_fe(_var.var())->get_JxW();
76  context.get_element_fe(_var.var())->get_phi();
77  context.get_element_fe(_var.var())->get_dphi();
78  context.get_element_fe(_var.var())->get_xyz();
79  }
80 
82  ( bool compute_jacobian,
83  AssemblyContext & context )
84  {
85  const unsigned int n_dofs = context.get_dof_indices(_var.var()).size();
86 
87  const std::vector<libMesh::Real> &JxW =
88  context.get_element_fe(_var.var())->get_JxW();
89 
90  const std::vector<std::vector<libMesh::Real> >& phi =
91  context.get_element_fe(_var.var())->get_phi();
92 
93  const std::vector<std::vector<libMesh::RealGradient> >& grad_phi =
94  context.get_element_fe(_var.var())->get_dphi();
95 
96  const std::vector<libMesh::Point>& x =
97  context.get_element_fe(_var.var())->get_xyz();
98 
99  libMesh::DenseSubVector<libMesh::Number> &F = context.get_elem_residual(_var.var());
100 
101  libMesh::DenseSubMatrix<libMesh::Number> &K = context.get_elem_jacobian(_var.var(), _var.var());
102 
103  unsigned int n_qpoints = context.get_element_qrule().n_points();
104 
105  for (unsigned int qp=0; qp != n_qpoints; qp++)
106  {
107  libMesh::Number c = context.interior_value(_var.var(), qp);
108 
109  libMesh::Gradient grad_c;
110  context.interior_gradient(_var.var(), qp, grad_c);
111 
112  libMesh::RealGradient v_qp;
113  v_qp(0) = this->_v[0](x[qp], context.get_time());
114  v_qp(1) = this->_v[1](x[qp], context.get_time());
115  v_qp(2) = this->_v[2](x[qp], context.get_time());
116 
117  libMesh::Real kappa_qp = this->_kappa(x[qp], context.time);
118 
119  for (unsigned int i=0; i != n_dofs; i++)
120  {
121  F(i) += JxW[qp]*( c*(v_qp*grad_phi[i][qp])
122  - kappa_qp*(grad_c*grad_phi[i][qp]) );
123 
124  if (compute_jacobian)
125  {
126  for (unsigned int j=0; j != n_dofs; j++)
127  {
128  K(i,j) += context.get_elem_solution_derivative()*
129  JxW[qp]*( phi[j][qp]*(v_qp*grad_phi[i][qp])
130  - kappa_qp*(grad_phi[j][qp]*grad_phi[i][qp]) );
131  }
132  }
133  }
134  }
135  }
136 
138  ( bool compute_jacobian, AssemblyContext & context )
139  {
140  const unsigned int n_dofs = context.get_dof_indices(_var.var()).size();
141 
142  const std::vector<libMesh::Real> &JxW =
143  context.get_element_fe(_var.var())->get_JxW();
144 
145  const std::vector<std::vector<libMesh::Real> >& phi =
146  context.get_element_fe(_var.var())->get_phi();
147 
148  libMesh::DenseSubVector<libMesh::Number> &F = context.get_elem_residual(_var.var());
149 
150  libMesh::DenseSubMatrix<libMesh::Number> &M = context.get_elem_jacobian(_var.var(), _var.var());
151 
152  unsigned int n_qpoints = context.get_element_qrule().n_points();
153 
154  for (unsigned int qp=0; qp != n_qpoints; qp++)
155  {
156  libMesh::Real c_dot;
157  context.interior_rate(_var.var(), qp, c_dot);
158 
159  for (unsigned int i=0; i != n_dofs; i++)
160  {
161  F(i) -= JxW[qp]*( c_dot*phi[i][qp] );
162 
163  if (compute_jacobian)
164  {
165  for (unsigned int j=0; j != n_dofs; j++)
166  {
167  M(i,j) -=
168  context.get_elem_solution_rate_derivative()*
169  JxW[qp]*( phi[j][qp]*phi[i][qp] );
170  }
171  }
172  }
173  }
174  }
175 
176 } // end namespace GRINS
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.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
Physics abstract base class. Defines API for physics to be added to MultiphysicsSystem.
Definition: physics.h:106
Variables with a single component.
void check_var_subdomain_consistency(const FEVariablesBase &var) const
Check that var is enabled on at least the subdomains this Physics is.
Definition: physics.C:174
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Set which variables are time evolving.
Base class for reading and handling initial conditions for physics classes.
VariableIndex var() const
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
libMesh::ParsedFunction< libMesh::Number > _kappa
Diffusivity, .
GRINS namespace.
std::string PhysicsName
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...
static std::string material_name(const GetPot &input, const std::string &physics)
Get the name of the material in the Physics/physics section.
std::vector< libMesh::ParsedFunction< libMesh::Number > > _v
Velocity field, ;.

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