GRINS-0.8.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-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 // 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  for( std::vector<VariableIndex>::const_iterator v_it = _vars.begin();
50  v_it != _vars.end(); ++v_it )
51  {
52  VariableIndex var = *v_it;
53 
54  context.get_element_fe(var)->get_JxW();
55  context.get_element_fe(var)->get_phi();
56  context.get_element_fe(var)->get_xyz();
57  }
58  }
59 
61  ( bool /*compute_jacobian*/, AssemblyContext & context )
62  {
63  for( std::vector<VariableIndex>::const_iterator v_it = _vars.begin();
64  v_it != _vars.end(); ++v_it )
65  {
66  VariableIndex var = *v_it;
67 
68  // The number of local degrees of freedom in each variable.
69  const unsigned int n_dofs = context.get_dof_indices(var).size();
70 
71  // Element Jacobian * quadrature weights for interior integration.
72  const std::vector<libMesh::Real> &JxW =
73  context.get_element_fe(var)->get_JxW();
74 
75  // The temperature shape functions at interior quadrature points.
76  const std::vector<std::vector<libMesh::Real> >& phi =
77  context.get_element_fe(var)->get_phi();
78 
79  const std::vector<libMesh::Point>& x_qp = context.get_element_fe(var)->get_xyz();
80 
81  // Get residuals
82  libMesh::DenseSubVector<libMesh::Number> &F_var = context.get_elem_residual(var);
83 
84  libMesh::Real t = context.get_time();
85 
86  // Now we will build the element Jacobian and residual.
87  // Constructing the residual requires the solution and its
88  // gradient from the previous timestep. This must be
89  // calculated at each quadrature point by summing the
90  // solution degree-of-freedom values by the appropriate
91  // weight functions.
92  unsigned int n_qpoints = context.get_element_qrule().n_points();
93 
94  for (unsigned int qp=0; qp != n_qpoints; qp++)
95  {
96  libMesh::Real value = (this->_value)(x_qp[qp],t);
97 
98  for (unsigned int i=0; i != n_dofs; i++)
99  {
100  F_var(i) += value*phi[i][qp]*JxW[qp];
101  }
102  }
103 
104  } // Variable loop
105 
106  return;
107  }
108 
109 } // 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
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.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.

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