GRINS-0.7.0
parsed_source_term.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2016 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 // This class
27 
28 // GRINS
29 #include "grins/assembly_context.h"
30 
31 // libMesh
32 #include "libmesh/getpot.h"
33 #include "libmesh/fem_system.h"
34 #include "libmesh/quadrature.h"
35 
36 namespace GRINS
37 {
38  ParsedSourceTerm::ParsedSourceTerm( const std::string& physics_name, const GetPot& input )
39  : SourceTermBase(physics_name,input),
40  _value("")
41  {
42  this->set_parameter(_value, input,
43  "Physics/"+physics_name+"/Function/value",
44  "DIE!");
45  }
46 
48  {
49  return;
50  }
51 
53  {
54  for( std::vector<VariableIndex>::const_iterator v_it = _vars.begin();
55  v_it != _vars.end(); ++v_it )
56  {
57  VariableIndex var = *v_it;
58 
59  context.get_element_fe(var)->get_JxW();
60  context.get_element_fe(var)->get_phi();
61  context.get_element_fe(var)->get_xyz();
62  }
63  }
64 
65  void ParsedSourceTerm::element_time_derivative( bool /*compute_jacobian*/,
66  AssemblyContext& context,
67  CachedValues& /*cache*/ )
68  {
69  for( std::vector<VariableIndex>::const_iterator v_it = _vars.begin();
70  v_it != _vars.end(); ++v_it )
71  {
72  VariableIndex var = *v_it;
73 
74  // The number of local degrees of freedom in each variable.
75  const unsigned int n_dofs = context.get_dof_indices(var).size();
76 
77  // Element Jacobian * quadrature weights for interior integration.
78  const std::vector<libMesh::Real> &JxW =
79  context.get_element_fe(var)->get_JxW();
80 
81  // The temperature shape functions at interior quadrature points.
82  const std::vector<std::vector<libMesh::Real> >& phi =
83  context.get_element_fe(var)->get_phi();
84 
85  const std::vector<libMesh::Point>& x_qp = context.get_element_fe(var)->get_xyz();
86 
87  // Get residuals
88  libMesh::DenseSubVector<libMesh::Number> &F_var = context.get_elem_residual(var);
89 
90  libMesh::Real t = context.get_time();
91 
92  // Now we will build the element Jacobian and residual.
93  // Constructing the residual requires the solution and its
94  // gradient from the previous timestep. This must be
95  // calculated at each quadrature point by summing the
96  // solution degree-of-freedom values by the appropriate
97  // weight functions.
98  unsigned int n_qpoints = context.get_element_qrule().n_points();
99 
100  for (unsigned int qp=0; qp != n_qpoints; qp++)
101  {
102  libMesh::Real value = (this->_value)(x_qp[qp],t);
103 
104  for (unsigned int i=0; i != n_dofs; i++)
105  {
106  F_var(i) += value*phi[i][qp]*JxW[qp];
107  }
108  }
109 
110  } // Variable loop
111 
112  return;
113  }
114 
115 } // end namespace GRINS
std::vector< VariableIndex > _vars
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.
unsigned int VariableIndex
More descriptive name of the type used for variable indices.
Definition: var_typedefs.h:42
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
libMesh::ParsedFunction< libMesh::Real > _value
GRINS namespace.
Base class for generic source function term, f(x,y,z,t).
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.

Generated on Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10